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