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