SCALE-RM
scale_prof.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 #include "scalelib.h"
11 module scale_prof
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_io
18  !-----------------------------------------------------------------------------
19  implicit none
20  private
21  !-----------------------------------------------------------------------------
22  !
23  !++ Public procedure
24  !
25  public :: prof_setup
26  public :: prof_setprefx
27  public :: prof_rapstart
28  public :: prof_rapend
29  public :: prof_rapreport
30 
31 #ifdef PAPI
32  public :: prof_papi_rapstart
33  public :: prof_papi_rapstop
34  public :: prof_papi_rapreport
35 #endif
36 
37  public :: prof_valcheck
38 
39  interface prof_valcheck
40  module procedure prof_valcheck_sp_1d
41  module procedure prof_valcheck_sp_2d
42  module procedure prof_valcheck_sp_3d
43  module procedure prof_valcheck_sp_4d
44  module procedure prof_valcheck_sp_5d
45  module procedure prof_valcheck_sp_6d
46  module procedure prof_valcheck_dp_1d
47  module procedure prof_valcheck_dp_2d
48  module procedure prof_valcheck_dp_3d
49  module procedure prof_valcheck_dp_4d
50  module procedure prof_valcheck_dp_5d
51  module procedure prof_valcheck_dp_6d
52  end interface prof_valcheck
53 
54  !-----------------------------------------------------------------------------
55  !
56  !++ Public parameters & variables
57  !
58  !-----------------------------------------------------------------------------
59  !
60  !++ Private procedure
61  !
62  private :: get_rapid
63 
64  !-----------------------------------------------------------------------------
65  !
66  !++ Private parameters & variables
67  !
68  integer, private, parameter :: prof_rapnlimit = 300
69  character(len=H_SHORT), private :: prof_prefix = ''
70  integer, private :: prof_rapnmax = 0
71  character(len=H_SHORT), private :: prof_rapname(prof_rapnlimit)
72  integer, private :: prof_grpnmax = 0
73  character(len=H_SHORT), private :: prof_grpname(prof_rapnlimit)
74  integer, private :: prof_grpid (prof_rapnlimit)
75  real(DP), private :: prof_raptstr(prof_rapnlimit)
76  real(DP), private :: prof_rapttot(prof_rapnlimit)
77  integer, private :: prof_rapnstr(prof_rapnlimit)
78  integer, private :: prof_rapnend(prof_rapnlimit)
79  integer, private :: prof_raplevel(prof_rapnlimit)
80 
81  integer, private, parameter :: prof_default_rap_level = 2
82  integer, private :: prof_rap_level = 2
83  logical, private :: prof_mpi_barrier = .false.
84 
85 #ifdef PAPI
86  integer(DP), private :: prof_papi_flops = 0
87  real(SP), private :: prof_papi_real_time = 0.0
88  real(SP), private :: prof_papi_proc_time = 0.0
89  real(SP), private :: prof_papi_mflops = 0.0
90  integer, private :: prof_papi_check
91 #endif
92 
93  character(len=7), private :: prof_header
94  character(len=16), private :: prof_item
95  real(DP), private :: prof_max
96  real(DP), private :: prof_min
97  real(DP), private :: prof_sum
98 
99  !-----------------------------------------------------------------------------
100 contains
101  !-----------------------------------------------------------------------------
102  subroutine prof_setup
103  use scale_prc, only: &
104  prc_abort
105  implicit none
106 
107  namelist / param_prof / &
108  prof_rap_level, &
109  prof_mpi_barrier
110 
111  integer :: ierr
112 
113  log_newline
114  log_info("PRF_setup",*) 'Setup'
115 
116  !--- read namelist
117  rewind(io_fid_conf)
118  read(io_fid_conf,nml=param_prof,iostat=ierr)
119  if( ierr < 0 ) then !--- missing
120  log_info("PROF_setup",*) 'Not found namelist. Default used.'
121  elseif( ierr > 0 ) then !--- fatal error
122  log_error("PROF_setup",*) 'Not appropriate names in namelist PARAM_PROF. Check!'
123  call prc_abort
124  endif
125  log_nml(param_prof)
126 
127  log_newline
128  log_info("PROF_setup",*) 'Rap output level = ', prof_rap_level
129  log_info("PROF_setup",*) 'Add MPI_barrier in every rap? = ', prof_mpi_barrier
130 
131  prof_prefix = ''
132 
133  return
134  end subroutine prof_setup
135 
136  !-----------------------------------------------------------------------------
137  subroutine prof_setprefx( &
138  prefxname )
139  implicit none
140 
141  character(len=*), intent(in) :: prefxname
142 
143  !---------------------------------------------------------------------------
144 
145  if ( prefxname == '' ) then !--- no prefix
146  prof_prefix = ''
147  else
148  prof_prefix = trim(prefxname)//'_'
149  endif
150 
151  return
152  end subroutine prof_setprefx
153 
154  !-----------------------------------------------------------------------------
156  subroutine prof_rapstart( rapname_base, level )
157  use scale_prc, only: &
158  prc_mpibarrier, &
160  implicit none
161 
162  character(len=*), intent(in) :: rapname_base
163  integer, intent(in), optional :: level
164 
165  character(len=H_SHORT) :: rapname
166 
167  integer :: id
168  integer :: level_
169  integer :: tn
170  !$ integer :: omp_get_thread_num
171  !---------------------------------------------------------------------------
172 
173  tn = 0
174  !$ tn = omp_get_thread_num()
175  if ( tn > 0 ) return
176 
177  if ( present(level) ) then
178  level_ = level
179  else
180  level_ = prof_default_rap_level
181  endif
182 
183  if( level_ > prof_rap_level ) return
184 
185  rapname = trim(prof_prefix)//trim(rapname_base)
186 
187  id = get_rapid( rapname, level_ )
188 
189  if(prof_mpi_barrier) call prc_mpibarrier
190 
191  prof_raptstr(id) = prc_mpitime()
192  prof_rapnstr(id) = prof_rapnstr(id) + 1
193 
194  !LOG_INFO("PROF_rapstart",'(1x,A,I8)') rapname, PROF_rapnstr(id)
195  !call flush(IO_FID_LOG)
196 
197 #ifdef FAPP
198  call fapp_start( trim(prof_grpname(get_grpid(rapname))), id, level_ )
199 #endif
200 #ifdef FINEPA
201  call start_collection( rapname )
202 #endif
203 
204  return
205  end subroutine prof_rapstart
206 
207  !-----------------------------------------------------------------------------
209  subroutine prof_rapend( rapname_base, level )
210  use scale_prc, only: &
211  prc_mpibarrier, &
213  implicit none
214 
215  character(len=*), intent(in) :: rapname_base
216  integer, intent(in), optional :: level
217 
218  character(len=H_SHORT) :: rapname
219 
220  integer :: id
221  integer :: level_
222  integer :: tn
223  !$ integer :: omp_get_thread_num
224  !---------------------------------------------------------------------------
225 
226  tn = 0
227  !$ tn = omp_get_thread_num()
228  if ( tn > 0 ) return
229 
230  if ( present(level) ) then
231  if( level > prof_rap_level ) return
232  endif
233 
234  rapname = trim(prof_prefix)//trim(rapname_base)
235 
236  id = get_rapid( rapname, level_ )
237 
238  if( level_ > prof_rap_level ) return
239 
240  if(prof_mpi_barrier) call prc_mpibarrier
241 
242  prof_rapttot(id) = prof_rapttot(id) + ( prc_mpitime()-prof_raptstr(id) )
243  prof_rapnend(id) = prof_rapnend(id) + 1
244 
245 #ifdef FINEPA
246  call stop_collection( rapname )
247 #endif
248 #ifdef FAPP
249  call fapp_stop( trim(prof_grpname(prof_grpid(id))), id, level_ )
250 #endif
251 
252  return
253  end subroutine prof_rapend
254 
255  !-----------------------------------------------------------------------------
257  subroutine prof_rapreport
258  use scale_prc, only: &
259  prc_mpitimestat, &
261  implicit none
262 
263  real(DP) :: avgvar(prof_rapnlimit)
264  real(DP) :: maxvar(prof_rapnlimit)
265  real(DP) :: minvar(prof_rapnlimit)
266  integer :: maxidx(prof_rapnlimit)
267  integer :: minidx(prof_rapnlimit)
268 
269  integer :: id, gid
270  integer :: fid
271  !---------------------------------------------------------------------------
272 
273  do id = 1, prof_rapnmax
274  if ( prof_rapnstr(id) /= prof_rapnend(id) ) then
275  log_warn("PROF_rapreport",*) 'Mismatch Report',id,prof_rapname(id),prof_rapnstr(id),prof_rapnend(id)
276  endif
277  enddo
278 
279  log_newline
280  log_info("PROF_rapreport",'(1x,A,I2,A)') 'Computational Time Report (Rap level = ', prof_rap_level, ')'
281 
282  if ( io_log_allnode ) then ! report for each node
283 
284  do gid = 1, prof_rapnmax
285  do id = 1, prof_rapnmax
286  if ( prof_raplevel(id) <= prof_rap_level &
287  .AND. prof_grpid(id) == gid ) then
288  log_info_cont('(1x,A,I3.3,A,A,A,F10.3,A,I9)') &
289  'ID=',id,' : ',prof_rapname(id),' T=',prof_rapttot(id),' N=',prof_rapnstr(id)
290  endif
291  enddo
292  enddo
293 
294  else
295 
296  call prc_mpitimestat( avgvar(1:prof_rapnmax), &
297  maxvar(1:prof_rapnmax), &
298  minvar(1:prof_rapnmax), &
299  maxidx(1:prof_rapnmax), &
300  minidx(1:prof_rapnmax), &
301  prof_rapttot(1:prof_rapnmax) )
302 
303  fid = -1
304  if ( io_log_suppress ) then ! report to STDOUT
305  if ( prc_ismaster ) then
306  write(*,*) 'INFO [PROF_rapreport] Computational Time Report'
307  fid = io_fid_stdout ! master node
308  endif
309  else
310  if ( io_l ) fid = io_fid_log
311  endif
312 
313  do gid = 1, prof_rapnmax
314  do id = 1, prof_rapnmax
315  if ( prof_raplevel(id) <= prof_rap_level &
316  .AND. prof_grpid(id) == gid &
317  .AND. fid > 0 ) then
318  write(fid,'(6x,A,I3.3,3A,F10.3,A,F10.3,A,I5,2A,F10.3,A,I5,2A,I9)') &
319  'ID=',id,' : ',prof_rapname(id), &
320  ' T(avg)=',avgvar(id), &
321  ', T(max)=',maxvar(id),'[',maxidx(id),']', &
322  ', T(min)=',minvar(id),'[',minidx(id),']', &
323  ', N=',prof_rapnstr(id)
324  endif
325  enddo
326  enddo
327 
328  endif
329 
330  return
331  end subroutine prof_rapreport
332 
333 #ifdef PAPI
334  !-----------------------------------------------------------------------------
336  subroutine prof_papi_rapstart
337  implicit none
338  !---------------------------------------------------------------------------
339 
340  call papif_flops( prof_papi_real_time, prof_papi_proc_time, prof_papi_flops, prof_papi_mflops, prof_papi_check )
341 
342  return
343  end subroutine prof_papi_rapstart
344 
345  !-----------------------------------------------------------------------------
347  subroutine prof_papi_rapstop
348  implicit none
349  !---------------------------------------------------------------------------
350 
351  call papif_flops( prof_papi_real_time, prof_papi_proc_time, prof_papi_flops, prof_papi_mflops, prof_papi_check )
352 
353  return
354  end subroutine prof_papi_rapstop
355 
356  !-----------------------------------------------------------------------------
358  subroutine prof_papi_rapreport
359  use scale_prc, only: &
360  prc_mpitimestat, &
361  prc_nprocs, &
363  implicit none
364 
365  real(DP) :: avgvar(3)
366  real(DP) :: maxvar(3)
367  real(DP) :: minvar(3)
368  integer :: maxidx(3)
369  integer :: minidx(3)
370 
371  real(DP) :: zerosw
372  real(DP) :: PROF_PAPI_gflop
373  real(DP) :: statistics(3)
374  !---------------------------------------------------------------------------
375 
376  prof_papi_gflop = real(PROF_PAPI_flops,kind=8) / 1024.0_DP**3
377 
378  if ( io_log_allnode ) then ! report for each node
379 
380  log_newline
381  log_info("PROF_PAPI_rapreport",*) 'PAPI Report [Local PE information]'
382  log_info("PROF_PAPI_rapreport",'(1x,A,F15.3)') 'Real time [sec] : ', prof_papi_real_time
383  log_info("PROF_PAPI_rapreport",'(1x,A,F15.3)') 'CPU time [sec] : ', prof_papi_proc_time
384  log_info("PROF_PAPI_rapreport",'(1x,A,F15.3)') 'FLOP [GFLOP] : ', prof_papi_gflop
385  log_info("PROF_PAPI_rapreport",'(1x,A,F15.3)') 'FLOPS by PAPI [GFLOPS] : ', prof_papi_mflops/1024.0_dp
386  log_info("PROF_PAPI_rapreport",'(1x,A,F15.3)') 'FLOP / CPU Time [GFLOPS] : ', prof_papi_gflop/prof_papi_proc_time
387 
388  else
389  statistics(1) = real(prof_papi_real_time,kind=8)
390  statistics(2) = real(prof_papi_proc_time,kind=8)
391  statistics(3) = prof_papi_gflop
392 
393  call prc_mpitimestat( avgvar(1:3), &
394  maxvar(1:3), &
395  minvar(1:3), &
396  maxidx(1:3), &
397  minidx(1:3), &
398  statistics(1:3) )
399 
400  zerosw = 0.5_dp - sign(0.5_dp,maxvar(2)-1.d-12) ! if maxvar(2) = 0 then zerosw = 1
401 
402  log_newline
403  log_info("PROF_PAPI_rapreport",*) 'PAPI Report'
404  log_info("PROF_PAPI_rapreport",'(1x,A,A,F10.3,A,F10.3,A,I5,A,A,F10.3,A,I5,A,A,I7)') &
405  'Real time [sec]',' T(avg)=',avgvar(1), &
406  ', T(max)=',maxvar(1),'[',maxidx(1),']',', T(min)=',minvar(1),'[',minidx(1),']'
407  log_info("PROF_PAPI_rapreport",'(1x,A,A,F10.3,A,F10.3,A,I5,A,A,F10.3,A,I5,A,A,I7)') &
408  'CPU time [sec]',' T(avg)=',avgvar(2), &
409  ', T(max)=',maxvar(2),'[',maxidx(2),']',', T(min)=',minvar(2),'[',minidx(2),']'
410  log_info("PROF_PAPI_rapreport",'(1x,A,A,F10.3,A,F10.3,A,I5,A,A,F10.3,A,I5,A,A,I7)') &
411  'FLOP [GFLOP]',' N(avg)=',avgvar(3), &
412  ', N(max)=',maxvar(3),'[',maxidx(3),']',', N(min)=',minvar(3),'[',minidx(3),']'
413  log_newline
414  log_info("PROF_PAPI_rapreport",'(1x,A,F15.3,A,I6,A)') &
415  'TOTAL FLOP [GFLOP] : ', avgvar(3)*prc_nprocs, '(',prc_nprocs,' PEs)'
416  log_info("PROF_PAPI_rapreport",'(1x,A,F15.3)') &
417  'FLOPS [GFLOPS] : ', avgvar(3)*prc_nprocs * ( 1.0_dp-zerosw ) / ( maxvar(2)+zerosw )
418  log_info("PROF_PAPI_rapreport",'(1x,A,F15.3)') &
419  'FLOPS per PE [GFLOPS] : ', avgvar(3) * ( 1.0_dp-zerosw ) / ( maxvar(2)+zerosw )
420  log_newline
421 
422  if ( io_log_suppress ) then ! report to STDOUT
423  if ( prc_ismaster ) then ! master node
424  write(*,*)
425  write(*,*) '*** PAPI Report'
426  write(*,'(1x,A,A,F10.3,A,F10.3,A,I5,A,A,F10.3,A,I5,A,A,I7)') &
427  '*** Real time [sec]',' T(avg)=',avgvar(1), &
428  ', T(max)=',maxvar(1),'[',maxidx(1),']',', T(min)=',minvar(1),'[',minidx(1),']'
429  write(*,'(1x,A,A,F10.3,A,F10.3,A,I5,A,A,F10.3,A,I5,A,A,I7)') &
430  '*** CPU time [sec]',' T(avg)=',avgvar(2), &
431  ', T(max)=',maxvar(2),'[',maxidx(2),']',', T(min)=',minvar(2),'[',minidx(2),']'
432  write(*,'(1x,A,A,F10.3,A,F10.3,A,I5,A,A,F10.3,A,I5,A,A,I7)') &
433  '*** FLOP [GFLOP]',' N(avg)=',avgvar(3), &
434  ', N(max)=',maxvar(3),'[',maxidx(3),']',', N(min)=',minvar(3),'[',minidx(3),']'
435  write(*,*)
436  write(*,'(1x,A,F15.3,A,I6,A)') &
437  '*** TOTAL FLOP [GFLOP] : ', avgvar(3)*prc_nprocs, '(',prc_nprocs,' PEs)'
438  write(*,'(1x,A,F15.3)') &
439  '*** FLOPS [GFLOPS] : ', avgvar(3)*prc_nprocs * ( 1.0_dp-zerosw ) / ( maxvar(2)+zerosw )
440  write(*,'(1x,A,F15.3)') &
441  '*** FLOPS per PE [GFLOPS] : ', avgvar(3) * ( 1.0_dp-zerosw ) / ( maxvar(2)+zerosw )
442  endif
443  endif
444  endif
445 
446  return
447  end subroutine prof_papi_rapreport
448 #endif
449 
450  !-----------------------------------------------------------------------------
452  function get_rapid( rapname, level ) result(id)
453  implicit none
454 
455  character(len=*), intent(in) :: rapname
456  integer, intent(inout) :: level
457 
458  character (len=H_SHORT) :: trapname
459 
460  integer :: id
461  !---------------------------------------------------------------------------
462 
463  trapname = trim(rapname)
464 
465  do id = 1, prof_rapnmax
466  if ( trapname == prof_rapname(id) ) then
467  level = prof_raplevel(id)
468  return
469  endif
470  enddo
471 
472  prof_rapnmax = prof_rapnmax + 1
473  id = prof_rapnmax
474  prof_rapname(id) = trapname
475 
476  prof_rapnstr(id) = 0
477  prof_rapnend(id) = 0
478  prof_rapttot(id) = 0.0_dp
479 
480  prof_grpid(id) = get_grpid(trapname)
481  prof_raplevel(id) = level
482 
483  return
484  end function get_rapid
485 
486  !-----------------------------------------------------------------------------
488  function get_grpid( rapname ) result(gid)
489  implicit none
490 
491  character(len=*), intent(in) :: rapname
492 
493  character(len=H_SHORT) :: grpname
494 
495  integer :: gid
496  integer :: idx
497  !---------------------------------------------------------------------------
498 
499  idx = index(rapname," ")
500  if ( idx > 1 ) then
501  grpname = rapname(1:idx-1)
502  else
503  grpname = rapname
504  endif
505 
506  do gid = 1, prof_grpnmax
507  if( grpname == prof_grpname(gid) ) return
508  enddo
509 
510  prof_grpnmax = prof_grpnmax + 1
511  gid = prof_grpnmax
512  prof_grpname(gid) = grpname
513 
514  return
515  end function get_grpid
516 
517  !-----------------------------------------------------------------------------
518  subroutine prof_valcheck_sp_1d( &
519  header, &
520  varname, &
521  var )
522  implicit none
523 
524  character(len=*), intent(in) :: header
525  character(len=*), intent(in) :: varname
526  real(SP), intent(in) :: var(:)
527  !---------------------------------------------------------------------------
528 
529  prof_header = trim(header)
530  prof_item = trim(varname)
531  prof_max = real(maxval(var),kind=dp)
532  prof_min = real(minval(var),kind=dp)
533  prof_sum = real(sum (var),kind=dp)
534  log_info("PROF_valcheck_SP_1D",'(1x,A,A7,A,A16,3(A,ES24.16))') &
535  '+',prof_header,'[',prof_item,'] max=',prof_max,',min=',prof_min,',sum=',prof_sum
536 
537  return
538  end subroutine prof_valcheck_sp_1d
539 
540  !-----------------------------------------------------------------------------
541  subroutine prof_valcheck_sp_2d( &
542  header, &
543  varname, &
544  var )
545  implicit none
546 
547  character(len=*), intent(in) :: header
548  character(len=*), intent(in) :: varname
549  real(SP), intent(in) :: var(:,:)
550  !---------------------------------------------------------------------------
551 
552  prof_header = trim(header)
553  prof_item = trim(varname)
554  prof_max = real(maxval(var),kind=dp)
555  prof_min = real(minval(var),kind=dp)
556  prof_sum = real(sum (var),kind=dp)
557  log_info("PROF_valcheck_SP_2D",'(1x,A,A7,A,A16,3(A,ES24.16))') &
558  '+',prof_header,'[',prof_item,'] max=',prof_max,',min=',prof_min,',sum=',prof_sum
559 
560  return
561  end subroutine prof_valcheck_sp_2d
562 
563  !-----------------------------------------------------------------------------
564  subroutine prof_valcheck_sp_3d( &
565  header, &
566  varname, &
567  var )
568  implicit none
569 
570  character(len=*), intent(in) :: header
571  character(len=*), intent(in) :: varname
572  real(SP), intent(in) :: var(:,:,:)
573  !---------------------------------------------------------------------------
574 
575  prof_header = trim(header)
576  prof_item = trim(varname)
577  prof_max = real(maxval(var),kind=dp)
578  prof_min = real(minval(var),kind=dp)
579  prof_sum = real(sum (var),kind=dp)
580  log_info("PROF_valcheck_SP_3D",'(1x,A,A7,A,A16,3(A,ES24.16))') &
581  '+',prof_header,'[',prof_item,'] max=',prof_max,',min=',prof_min,',sum=',prof_sum
582 
583  return
584  end subroutine prof_valcheck_sp_3d
585 
586  !-----------------------------------------------------------------------------
587  subroutine prof_valcheck_sp_4d( &
588  header, &
589  varname, &
590  var )
591  implicit none
592 
593  character(len=*), intent(in) :: header
594  character(len=*), intent(in) :: varname
595  real(SP), intent(in) :: var(:,:,:,:)
596  !---------------------------------------------------------------------------
597 
598  prof_header = trim(header)
599  prof_item = trim(varname)
600  prof_max = real(maxval(var),kind=dp)
601  prof_min = real(minval(var),kind=dp)
602  prof_sum = real(sum (var),kind=dp)
603  log_info("PROF_valcheck_SP_4D",'(1x,A,A7,A,A16,3(A,ES24.16))') &
604  '+',prof_header,'[',prof_item,'] max=',prof_max,',min=',prof_min,',sum=',prof_sum
605 
606  return
607  end subroutine prof_valcheck_sp_4d
608 
609  !-----------------------------------------------------------------------------
610  subroutine prof_valcheck_sp_5d( &
611  header, &
612  varname, &
613  var )
614  implicit none
615 
616  character(len=*), intent(in) :: header
617  character(len=*), intent(in) :: varname
618  real(SP), intent(in) :: var(:,:,:,:,:)
619  !---------------------------------------------------------------------------
620 
621  prof_header = trim(header)
622  prof_item = trim(varname)
623  prof_max = real(maxval(var),kind=dp)
624  prof_min = real(minval(var),kind=dp)
625  prof_sum = real(sum (var),kind=dp)
626  log_info("PROF_valcheck_SP_5D",'(1x,A,A7,A,A16,3(A,ES24.16))') &
627  '+',prof_header,'[',prof_item,'] max=',prof_max,',min=',prof_min,',sum=',prof_sum
628 
629  return
630  end subroutine prof_valcheck_sp_5d
631 
632  !-----------------------------------------------------------------------------
633  subroutine prof_valcheck_sp_6d( &
634  header, &
635  varname, &
636  var )
637  implicit none
638 
639  character(len=*), intent(in) :: header
640  character(len=*), intent(in) :: varname
641  real(SP), intent(in) :: var(:,:,:,:,:,:)
642  !---------------------------------------------------------------------------
643 
644  prof_header = trim(header)
645  prof_item = trim(varname)
646  prof_max = real(maxval(var),kind=dp)
647  prof_min = real(minval(var),kind=dp)
648  prof_sum = real(sum (var),kind=dp)
649  log_info("PROF_valcheck_SP_6D",'(1x,A,A7,A,A16,3(A,ES24.16))') &
650  '+',prof_header,'[',prof_item,'] max=',prof_max,',min=',prof_min,',sum=',prof_sum
651 
652  return
653  end subroutine prof_valcheck_sp_6d
654 
655  !-----------------------------------------------------------------------------
656  subroutine prof_valcheck_dp_1d( &
657  header, &
658  varname, &
659  var )
660  implicit none
661 
662  character(len=*), intent(in) :: header
663  character(len=*), intent(in) :: varname
664  real(DP), intent(in) :: var(:)
665  !---------------------------------------------------------------------------
666 
667  prof_header = trim(header)
668  prof_item = trim(varname)
669  prof_max = real(maxval(var),kind=dp)
670  prof_min = real(minval(var),kind=dp)
671  prof_sum = real(sum (var),kind=dp)
672  log_info("PROF_valcheck_DP_1D",'(1x,A,A7,A,A16,3(A,ES24.16))') &
673  '+',prof_header,'[',prof_item,'] max=',prof_max,',min=',prof_min,',sum=',prof_sum
674 
675  return
676  end subroutine prof_valcheck_dp_1d
677 
678  !-----------------------------------------------------------------------------
679  subroutine prof_valcheck_dp_2d( &
680  header, &
681  varname, &
682  var )
683  implicit none
684 
685  character(len=*), intent(in) :: header
686  character(len=*), intent(in) :: varname
687  real(DP), intent(in) :: var(:,:)
688  !---------------------------------------------------------------------------
689 
690  prof_header = trim(header)
691  prof_item = trim(varname)
692  prof_max = real(maxval(var),kind=dp)
693  prof_min = real(minval(var),kind=dp)
694  prof_sum = real(sum (var),kind=dp)
695  log_info("PROF_valcheck_DP_2D",'(1x,A,A7,A,A16,3(A,ES24.16))') &
696  '+',prof_header,'[',prof_item,'] max=',prof_max,',min=',prof_min,',sum=',prof_sum
697 
698  return
699  end subroutine prof_valcheck_dp_2d
700 
701  !-----------------------------------------------------------------------------
702  subroutine prof_valcheck_dp_3d( &
703  header, &
704  varname, &
705  var )
706  implicit none
707 
708  character(len=*), intent(in) :: header
709  character(len=*), intent(in) :: varname
710  real(DP), intent(in) :: var(:,:,:)
711  !---------------------------------------------------------------------------
712 
713  prof_header = trim(header)
714  prof_item = trim(varname)
715  prof_max = real(maxval(var),kind=dp)
716  prof_min = real(minval(var),kind=dp)
717  prof_sum = real(sum (var),kind=dp)
718  log_info("PROF_valcheck_DP_3D",'(1x,A,A7,A,A16,3(A,ES24.16))') &
719  '+',prof_header,'[',prof_item,'] max=',prof_max,',min=',prof_min,',sum=',prof_sum
720 
721  return
722  end subroutine prof_valcheck_dp_3d
723 
724  !-----------------------------------------------------------------------------
725  subroutine prof_valcheck_dp_4d( &
726  header, &
727  varname, &
728  var )
729  implicit none
730 
731  character(len=*), intent(in) :: header
732  character(len=*), intent(in) :: varname
733  real(DP), intent(in) :: var(:,:,:,:)
734  !---------------------------------------------------------------------------
735 
736  prof_header = trim(header)
737  prof_item = trim(varname)
738  prof_max = real(maxval(var),kind=dp)
739  prof_min = real(minval(var),kind=dp)
740  prof_sum = real(sum (var),kind=dp)
741  log_info("PROF_valcheck_DP_4D",'(1x,A,A7,A,A16,3(A,ES24.16))') &
742  '+',prof_header,'[',prof_item,'] max=',prof_max,',min=',prof_min,',sum=',prof_sum
743 
744  return
745  end subroutine prof_valcheck_dp_4d
746 
747  !-----------------------------------------------------------------------------
748  subroutine prof_valcheck_dp_5d( &
749  header, &
750  varname, &
751  var )
752  implicit none
753 
754  character(len=*), intent(in) :: header
755  character(len=*), intent(in) :: varname
756  real(DP), intent(in) :: var(:,:,:,:,:)
757  !---------------------------------------------------------------------------
758 
759  prof_header = trim(header)
760  prof_item = trim(varname)
761  prof_max = real(maxval(var),kind=dp)
762  prof_min = real(minval(var),kind=dp)
763  prof_sum = real(sum (var),kind=dp)
764  log_info("PROF_valcheck_DP_5D",'(1x,A,A7,A,A16,3(A,ES24.16))') &
765  '+',prof_header,'[',prof_item,'] max=',prof_max,',min=',prof_min,',sum=',prof_sum
766 
767  return
768  end subroutine prof_valcheck_dp_5d
769 
770  !-----------------------------------------------------------------------------
771  subroutine prof_valcheck_dp_6d( &
772  header, &
773  varname, &
774  var )
775  implicit none
776 
777  character(len=*), intent(in) :: header
778  character(len=*), intent(in) :: varname
779  real(DP), intent(in) :: var(:,:,:,:,:,:)
780  !---------------------------------------------------------------------------
781 
782  prof_header = trim(header)
783  prof_item = trim(varname)
784  prof_max = real(maxval(var),kind=dp)
785  prof_min = real(minval(var),kind=dp)
786  prof_sum = real(sum (var),kind=dp)
787  log_info("PROF_valcheck_DP_6D",'(1x,A,A7,A,A16,3(A,ES24.16))') &
788  '+',prof_header,'[',prof_item,'] max=',prof_max,',min=',prof_min,',sum=',prof_sum
789 
790  return
791  end subroutine prof_valcheck_dp_6d
792 
793 end module scale_prof
subroutine, public prof_setup
Definition: scale_prof.F90:103
subroutine, public prof_setprefx(prefxname)
Definition: scale_prof.F90:139
logical, public io_log_suppress
suppress all of log output?
Definition: scale_io.F90:63
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
integer, public prc_nprocs
myrank in local communicator
Definition: scale_prc.F90:88
logical, public io_l
output log or not? (this process)
Definition: scale_io.F90:61
integer, parameter, public io_fid_stdout
Definition: scale_io.F90:54
logical, public io_log_allnode
output log for each node?
Definition: scale_io.F90:65
module PROCESS
Definition: scale_prc.F90:11
real(dp) function, public prc_mpitime()
Get MPI time.
Definition: scale_prc.F90:813
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:157
subroutine, public prc_mpibarrier
Barrier MPI.
Definition: scale_prc.F90:798
module profiler
Definition: scale_prof.F90:11
logical, public prc_ismaster
master process in local communicator?
Definition: scale_prc.F90:90
module PRECISION
module STDIO
Definition: scale_io.F90:10
subroutine, public prof_rapend(rapname_base, level)
Save raptime.
Definition: scale_prof.F90:210
integer, public io_fid_log
Log file ID.
Definition: scale_io.F90:56
subroutine, public prc_mpitimestat(avgvar, maxvar, minvar, maxidx, minidx, var)
Calc global statistics for timer.
Definition: scale_prc.F90:835
integer, parameter, public dp
subroutine, public prof_rapreport
Report raptime.
Definition: scale_prof.F90:258