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  !-----------------------------------------------------------------------------
37  !
38  !++ Public parameters & variables
39  !
40  !-----------------------------------------------------------------------------
41  !
42  !++ Private procedure
43  !
44  private :: get_rapid
45 
46  !-----------------------------------------------------------------------------
47  !
48  !++ Private parameters & variables
49  !
50  integer, private, parameter :: prof_rapnlimit = 300
51  character(len=H_SHORT), private :: prof_prefix = ''
52  integer, private :: prof_rapnmax = 0
53  character(len=H_SHORT*2), private :: prof_rapname(prof_rapnlimit)
54  integer, private :: prof_grpnmax = 0
55  character(len=H_SHORT), private :: prof_grpname(prof_rapnlimit)
56  integer, private :: prof_grpid (prof_rapnlimit)
57  real(DP), private :: prof_raptstr(prof_rapnlimit)
58  real(DP), private :: prof_rapttot(prof_rapnlimit)
59  integer, private :: prof_rapnstr(prof_rapnlimit)
60  integer, private :: prof_rapnend(prof_rapnlimit)
61  integer, private :: prof_raplevel(prof_rapnlimit)
62 
63  integer, private, parameter :: prof_default_rap_level = 2
64  integer, private :: prof_rap_level = 2
65  logical, private :: prof_mpi_barrier = .false.
66 
67 #ifdef _PAPI_
68  integer(DP),private :: prof_papi_flops = 0
69  real(SP), private :: prof_papi_real_time = 0.0
70  real(SP), private :: prof_papi_proc_time = 0.0
71  real(SP), private :: prof_papi_mflops = 0.0
72  integer, private :: prof_papi_check
73 #endif
74 
75  !-----------------------------------------------------------------------------
76 contains
77  !-----------------------------------------------------------------------------
78  subroutine prof_setup
79  use scale_process, only: &
81  implicit none
82 
83  namelist / param_prof / &
84  prof_rap_level, &
85  prof_mpi_barrier
86 
87  integer :: ierr
88 
89  if( io_l ) write(io_fid_log,*)
90  if( io_l ) write(io_fid_log,*) '++++++ Module[PROF] / Categ[COMMON] / Origin[SCALElib]'
91 
92  !--- read namelist
93  rewind(io_fid_conf)
94  read(io_fid_conf,nml=param_prof,iostat=ierr)
95  if( ierr < 0 ) then !--- missing
96  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
97  elseif( ierr > 0 ) then !--- fatal error
98  write(*,*) 'xxx Not appropriate names in namelist PARAM_PROF. Check!'
99  call prc_mpistop
100  endif
101  if( io_lnml ) write(io_fid_log,nml=param_prof)
102 
103  if( io_l ) write(io_fid_log,*) '*** Rap output level = ', prof_rap_level
104  if( io_l ) write(io_fid_log,*) '*** Add MPI_barrier in every rap? = ', prof_mpi_barrier
105 
106  prof_prefix = ''
107 
108  return
109  end subroutine prof_setup
110 
111  !-----------------------------------------------------------------------------
112  subroutine prof_setprefx( &
113  prefxname )
114  implicit none
115 
116  character(len=*), intent(in) :: prefxname
117 
118  !---------------------------------------------------------------------------
119 
120  if ( prefxname == '' ) then !--- no prefix
121  prof_prefix = ''
122  else
123  prof_prefix = trim(prefxname)//'_'
124  endif
125 
126  return
127  end subroutine prof_setprefx
128 
129  !-----------------------------------------------------------------------------
131  subroutine prof_rapstart( rapname_base, level )
132  use scale_process, only: &
133  prc_mpibarrier, &
135  implicit none
136 
137  character(len=*), intent(in) :: rapname_base
138  integer, intent(in), optional :: level
139 
140  character(len=H_SHORT*2) :: rapname
141 
142  integer :: id
143  integer :: level_
144  !---------------------------------------------------------------------------
145 
146  if ( present(level) ) then
147  level_ = level
148  else
149  level_ = prof_default_rap_level
150  endif
151 
152  if( level_ > prof_rap_level ) return
153 
154  rapname = trim(prof_prefix)//trim(rapname_base)
155 
156  id = get_rapid( rapname, level_ )
157 
158  if(prof_mpi_barrier) call prc_mpibarrier
159 
160  prof_raptstr(id) = prc_mpitime()
161  prof_rapnstr(id) = prof_rapnstr(id) + 1
162 
163  !if( IO_L ) write(IO_FID_LOG,*) rapname, PROF_rapnstr(id)
164 
165 #ifdef _FAPP_
166  call fapp_start( trim(prof_grpname(get_grpid(rapname))), id, level_ )
167 #endif
168 #ifdef _FINEPA_
169  call start_collection( rapname )
170 #endif
171 
172  return
173  end subroutine prof_rapstart
174 
175  !-----------------------------------------------------------------------------
177  subroutine prof_rapend( rapname_base, level )
178  use scale_process, only: &
179  prc_mpibarrier, &
181  implicit none
182 
183  character(len=*), intent(in) :: rapname_base
184  integer, intent(in), optional :: level
185 
186  character(len=H_SHORT*2) :: rapname
187 
188  integer :: id
189  integer :: level_
190  !---------------------------------------------------------------------------
191 
192  if ( present(level) ) then
193  if( level > prof_rap_level ) return
194  endif
195 
196  rapname = trim(prof_prefix)//trim(rapname_base)
197 
198  id = get_rapid( rapname, level_ )
199 
200  if( level_ > prof_rap_level ) return
201 
202  if(prof_mpi_barrier) call prc_mpibarrier
203 
204  prof_rapttot(id) = prof_rapttot(id) + ( prc_mpitime()-prof_raptstr(id) )
205  prof_rapnend(id) = prof_rapnend(id) + 1
206 
207 #ifdef _FINEPA_
208  call stop_collection( rapname )
209 #endif
210 #ifdef _FAPP_
211  call fapp_stop( trim(prof_grpname(prof_grpid(id))), id, level_ )
212 #endif
213 
214  return
215  end subroutine prof_rapend
216 
217  !-----------------------------------------------------------------------------
219  subroutine prof_rapreport
220  use scale_process, only: &
221  prc_mpitimestat, &
223  implicit none
224 
225  real(DP) :: avgvar(prof_rapnlimit)
226  real(DP) :: maxvar(prof_rapnlimit)
227  real(DP) :: minvar(prof_rapnlimit)
228  integer :: maxidx(prof_rapnlimit)
229  integer :: minidx(prof_rapnlimit)
230 
231  integer :: id, gid
232  integer :: fid
233  !---------------------------------------------------------------------------
234 
235  do id = 1, prof_rapnmax
236  if ( prof_rapnstr(id) /= prof_rapnend(id) ) then
237  write(*,*) '*** Mismatch Report',id,prof_rapname(id),prof_rapnstr(id),prof_rapnend(id)
238  endif
239  enddo
240 
241  if( io_l ) write(io_fid_log,*)
242  if( io_l ) write(io_fid_log,*) '*** Computational Time Report'
243  if( io_l ) write(io_fid_log,*) '*** Rap level is ', prof_rap_level
244 
245  if ( io_log_allnode ) then ! report for each node
246 
247  do gid = 1, prof_rapnmax
248  do id = 1, prof_rapnmax
249  if ( prof_raplevel(id) <= prof_rap_level .AND. &
250  prof_grpid(id) == gid .AND. &
251  io_l ) then
252  write(io_fid_log,'(1x,A,I3.3,A,A,A,F10.3,A,I9)') &
253  '*** ID=',id,' : ',prof_rapname(id),' T=',prof_rapttot(id),' N=',prof_rapnstr(id)
254  endif
255  enddo
256  enddo
257 
258  else
259 
260  call prc_mpitimestat( avgvar(1:prof_rapnmax), &
261  maxvar(1:prof_rapnmax), &
262  minvar(1:prof_rapnmax), &
263  maxidx(1:prof_rapnmax), &
264  minidx(1:prof_rapnmax), &
265  prof_rapttot(1:prof_rapnmax) )
266 
267  fid = -1
268  if ( io_log_suppress ) then ! report to STDOUT
269  if ( prc_ismaster ) then
270  write(*,*) '*** Computational Time Report'
271  fid = 6 ! master node
272  endif
273  else
274  if ( io_l ) fid = io_fid_log
275  endif
276 
277  do gid = 1, prof_rapnmax
278  do id = 1, prof_rapnmax
279  if ( prof_raplevel(id) <= prof_rap_level .AND. &
280  prof_grpid(id) == gid .AND. &
281  fid > 0 ) then
282  write(io_fid_log,'(1x,A,I3.3,A,A,A,F10.3,A,F10.3,A,I5,A,A,F10.3,A,I5,A,A,I9)') &
283  '*** ID=',id,' : ',prof_rapname(id), &
284  ' T(avg)=',avgvar(id), &
285  ', T(max)=',maxvar(id),'[',maxidx(id),']', &
286  ', T(min)=',minvar(id),'[',minidx(id),']', &
287  ' N=',prof_rapnstr(id)
288  endif
289  enddo
290  enddo
291 
292  endif
293 
294  return
295  end subroutine prof_rapreport
296 
297 #ifdef _PAPI_
298  !-----------------------------------------------------------------------------
300  subroutine prof_papi_rapstart
301  implicit none
302  !---------------------------------------------------------------------------
303 
304  call papif_flops( prof_papi_real_time, prof_papi_proc_time, prof_papi_flops, prof_papi_mflops, prof_papi_check )
305 
306  return
307  end subroutine prof_papi_rapstart
308 
309  !-----------------------------------------------------------------------------
311  subroutine prof_papi_rapstop
312  implicit none
313  !---------------------------------------------------------------------------
314 
315  call papif_flops( prof_papi_real_time, prof_papi_proc_time, prof_papi_flops, prof_papi_mflops, prof_papi_check )
316 
317  return
318  end subroutine prof_papi_rapstop
319 
320  !-----------------------------------------------------------------------------
322  subroutine prof_papi_rapreport
323  use scale_process, only: &
324  prc_mpitimestat, &
325  prc_nprocs, &
327  implicit none
328 
329  real(DP) :: avgvar(3)
330  real(DP) :: maxvar(3)
331  real(DP) :: minvar(3)
332  integer :: maxidx(3)
333  integer :: minidx(3)
334 
335  real(DP) :: PROF_PAPI_gflop
336  real(DP) :: statistics(3)
337  !---------------------------------------------------------------------------
338 
339  prof_papi_gflop = real(PROF_PAPI_flops,kind=8) / 1024.0_DP**3
340 
341  if ( io_log_allnode ) then ! report for each node
342 
343  if( io_l ) write(io_fid_log,*)
344  if( io_l ) write(io_fid_log,*) '*** PAPI Report [Local PE information]'
345  if( io_l ) write(io_fid_log,'(1x,A,F15.3)') '*** Real time [sec] : ', prof_papi_real_time
346  if( io_l ) write(io_fid_log,'(1x,A,F15.3)') '*** CPU time [sec] : ', prof_papi_proc_time
347  if( io_l ) write(io_fid_log,'(1x,A,F15.3)') '*** FLOP [GFLOP] : ', prof_papi_gflop
348  if( io_l ) write(io_fid_log,'(1x,A,F15.3)') '*** FLOPS by PAPI [GFLOPS] : ', prof_papi_mflops/1024.0_dp
349  if( io_l ) write(io_fid_log,'(1x,A,F15.3)') '*** FLOP / CPU Time [GFLOPS] : ', prof_papi_gflop/prof_papi_proc_time
350 
351  else
352  statistics(1) = real(prof_papi_real_time,kind=8)
353  statistics(2) = real(prof_papi_proc_time,kind=8)
354  statistics(3) = prof_papi_gflop
355 
356  call prc_mpitimestat( avgvar(1:3), &
357  maxvar(1:3), &
358  minvar(1:3), &
359  maxidx(1:3), &
360  minidx(1:3), &
361  statistics(1:3) )
362 
363  if( io_l ) write(io_fid_log,*)
364  if( io_l ) write(io_fid_log,*) '*** PAPI Report'
365  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)') &
366  '*** Real time [sec]',' T(avg)=',avgvar(1), &
367  ', T(max)=',maxvar(1),'[',maxidx(1),']',', T(min)=',minvar(1),'[',minidx(1),']'
368  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)') &
369  '*** CPU time [sec]',' T(avg)=',avgvar(2), &
370  ', T(max)=',maxvar(2),'[',maxidx(2),']',', T(min)=',minvar(2),'[',minidx(2),']'
371  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)') &
372  '*** FLOP [GFLOP]',' N(avg)=',avgvar(3), &
373  ', N(max)=',maxvar(3),'[',maxidx(3),']',', N(min)=',minvar(3),'[',minidx(3),']'
374  if( io_l ) write(io_fid_log,*)
375  if( io_l ) write(io_fid_log,'(1x,A,F15.3,A,I6,A)') &
376  '*** TOTAL FLOP [GFLOP] : ', avgvar(3)*prc_nprocs, '(',prc_nprocs,' PEs)'
377  if( io_l ) write(io_fid_log,'(1x,A,F15.3)') &
378  '*** FLOPS [GFLOPS] : ', avgvar(3)*prc_nprocs/maxvar(2)
379  if( io_l ) write(io_fid_log,'(1x,A,F15.3)') &
380  '*** FLOPS per PE [GFLOPS] : ', avgvar(3)/maxvar(2)
381  if( io_l ) write(io_fid_log,*)
382 
383  if ( io_log_suppress ) then ! report to STDOUT
384  if ( prc_ismaster ) then ! master node
385  write(*,*)
386  write(*,*) '*** PAPI Report'
387  write(*,'(1x,A,A,F10.3,A,F10.3,A,I5,A,A,F10.3,A,I5,A,A,I7)') &
388  '*** Real time [sec]',' T(avg)=',avgvar(1), &
389  ', T(max)=',maxvar(1),'[',maxidx(1),']',', T(min)=',minvar(1),'[',minidx(1),']'
390  write(*,'(1x,A,A,F10.3,A,F10.3,A,I5,A,A,F10.3,A,I5,A,A,I7)') &
391  '*** CPU time [sec]',' T(avg)=',avgvar(2), &
392  ', T(max)=',maxvar(2),'[',maxidx(2),']',', T(min)=',minvar(2),'[',minidx(2),']'
393  write(*,'(1x,A,A,F10.3,A,F10.3,A,I5,A,A,F10.3,A,I5,A,A,I7)') &
394  '*** FLOP [GFLOP]',' N(avg)=',avgvar(3), &
395  ', N(max)=',maxvar(3),'[',maxidx(3),']',', N(min)=',minvar(3),'[',minidx(3),']'
396  write(*,*)
397  write(*,'(1x,A,F15.3,A,I6,A)') &
398  '*** TOTAL FLOP [GFLOP] : ', avgvar(3)*prc_nprocs, '(',prc_nprocs,' PEs)'
399  write(*,'(1x,A,F15.3)') &
400  '*** FLOPS [GFLOPS] : ', avgvar(3)*prc_nprocs/maxvar(2)
401  write(*,'(1x,A,F15.3)') &
402  '*** FLOPS per PE [GFLOPS] : ', avgvar(3)/maxvar(2)
403  endif
404  endif
405  endif
406 
407  return
408  end subroutine prof_papi_rapreport
409 #endif
410 
411  !-----------------------------------------------------------------------------
413  function get_rapid( rapname, level ) result(id)
414  implicit none
415 
416  character(len=*), intent(in) :: rapname
417  integer, intent(inout) :: level
418 
419  integer :: id
420  character (len=H_SHORT*2) :: trapname
421  character (len=H_SHORT) :: trapname2
422  !---------------------------------------------------------------------------
423 
424  trapname = trim(rapname)
425  trapname2 = trim(rapname)
426 
427  do id = 1, prof_rapnmax
428  if ( trapname == prof_rapname(id) ) then
429  level = prof_raplevel(id)
430  return
431  endif
432  enddo
433 
434  prof_rapnmax = prof_rapnmax + 1
435  id = prof_rapnmax
436  prof_rapname(id) = trapname
437 
438  prof_rapnstr(id) = 0
439  prof_rapnend(id) = 0
440  prof_rapttot(id) = 0.0_dp
441 
442  prof_grpid(id) = get_grpid(trapname2)
443  prof_raplevel(id) = level
444 
445  return
446  end function get_rapid
447 
448  !-----------------------------------------------------------------------------
450  function get_grpid( rapname ) result(gid)
451  implicit none
452 
453  character(len=*), intent(in) :: rapname
454 
455  integer :: gid
456  integer :: idx
457  character(len=H_SHORT) :: grpname
458  !---------------------------------------------------------------------------
459 
460  idx = index(rapname," ")
461  if ( idx > 1 ) then
462  grpname = rapname(1:idx-1)
463  else
464  grpname = rapname
465  endif
466 
467  do gid = 1, prof_grpnmax
468  if( grpname == prof_grpname(gid) ) return
469  enddo
470 
471  prof_grpnmax = prof_grpnmax + 1
472  gid = prof_grpnmax
473  prof_grpname(gid) = grpname
474 
475  return
476  end function get_grpid
477 
478 end module scale_prof
subroutine, public prof_setup
Definition: scale_prof.F90:79
logical, public prc_ismaster
master process in local communicator?
subroutine, public prc_mpistop
Abort MPI.
subroutine, public prof_setprefx(prefxname)
Definition: scale_prof.F90:114
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
module STDIO
Definition: scale_stdio.F90:12
logical, public io_log_suppress
suppress all of log output?
Definition: scale_stdio.F90:61
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:62
real(dp) function, public prc_mpitime()
Get MPI time.
module PROCESS
subroutine, public prof_rapstart(rapname_base, level)
Start raptime.
Definition: scale_prof.F90:132
module profiler
Definition: scale_prof.F90:10
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
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:178
integer, public prc_nprocs
myrank in local communicator
subroutine, public prof_rapreport
Report raptime.
Definition: scale_prof.F90:220