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