SCALE-RM
Data Types | Functions/Subroutines
scale_prof Module Reference

module profiler More...

Functions/Subroutines

subroutine, public prof_setup
 
subroutine, public prof_setprefx (prefxname)
 
subroutine, public prof_rapstart (rapname_base, level)
 Start raptime. More...
 
subroutine, public prof_rapend (rapname_base, level)
 Save raptime. More...
 
subroutine, public prof_rapreport
 Report raptime. More...
 

Detailed Description

module profiler

Description
Time counter & FLOP counter(PAPI) toolbox
Author
Team SCALE
NAMELIST
  • PARAM_PROF
    nametypedefault valuecomment
    PROF_RAP_LEVEL integer 2
    PROF_MPI_BARRIER logical .false.

History Output
No history output

Function/Subroutine Documentation

◆ prof_setup()

subroutine, public scale_prof::prof_setup ( )

Definition at line 102 of file scale_prof.F90.

References scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_fid_nml, scale_stdio::io_l, scale_stdio::io_nml, and scale_process::prc_mpistop().

Referenced by mod_rm_driver::scalerm(), and mod_rm_prep::scalerm_prep().

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
subroutine, public prc_mpistop
Abort MPI.
module PROCESS
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
Here is the call graph for this function:
Here is the caller graph for this function:

◆ prof_setprefx()

subroutine, public scale_prof::prof_setprefx ( character(len=*), intent(in)  prefxname)
Parameters
[in]prefxnameprefix

Definition at line 138 of file scale_prof.F90.

Referenced by mod_rm_driver::scalerm(), and mod_rm_prep::scalerm_prep().

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
Here is the caller graph for this function:

◆ prof_rapstart()

subroutine, public scale_prof::prof_rapstart ( character(len=*), intent(in)  rapname_base,
integer, intent(in), optional  level 
)

Start raptime.

Parameters
[in]rapname_basename of item
[in]levellevel of item

Definition at line 156 of file scale_prof.F90.

References scale_process::prc_mpibarrier(), and scale_process::prc_mpitime().

Referenced by scale_atmos_phy_ae_kajino13::aerosol_zerochem(), mod_atmos_driver::atmos_driver(), mod_atmos_driver::atmos_driver_resume1(), mod_atmos_driver::atmos_driver_resume2(), mod_atmos_driver::atmos_driver_setup(), scale_atmos_dyn::atmos_dyn(), scale_atmos_dyn_common::atmos_dyn_divergence(), scale_atmos_dyn_common::atmos_dyn_numfilter_coef(), scale_atmos_dyn_common::atmos_dyn_numfilter_coef_q(), scale_atmos_dyn_tinteg_large_rk3::atmos_dyn_tinteg_large_rk3(), scale_atmos_dyn_tinteg_short_rk3::atmos_dyn_tinteg_short_rk3(), scale_atmos_dyn_tinteg_short_rk4::atmos_dyn_tinteg_short_rk4(), scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve(), mod_atmos_phy_ae_driver::atmos_phy_ae_driver_resume(), mod_atmos_phy_ch_driver::atmos_phy_ch_driver_resume(), mod_atmos_phy_cp_driver::atmos_phy_cp_driver_resume(), scale_atmos_phy_cp_kf::atmos_phy_cp_kf(), mod_atmos_phy_mp_driver::atmos_phy_mp_driver_resume(), scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler(), scale_atmos_phy_mp_common::atmos_phy_mp_negative_fixer(), scale_atmos_phy_mp_common::atmos_phy_mp_precipitation(), scale_atmos_phy_mp_common::atmos_phy_mp_saturation_adjustment(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08(), mod_atmos_phy_rd_driver::atmos_phy_rd_driver_resume(), scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx(), mod_atmos_phy_sf_driver::atmos_phy_sf_driver_resume(), mod_atmos_phy_tb_driver::atmos_phy_tb_driver_resume(), mod_atmos_vars::atmos_vars_restart_check(), scale_atmos_dyn_common::calc_numdiff(), scale_debug::check(), scale_comm::comm_bcast_1d(), scale_comm::comm_bcast_2d(), scale_comm::comm_bcast_3d(), scale_comm::comm_bcast_4d(), scale_comm::comm_bcast_int_1d(), scale_comm::comm_bcast_int_2d(), scale_comm::comm_bcast_int_scr(), scale_comm::comm_bcast_logical_scr(), scale_comm::comm_bcast_scr(), scale_comm::comm_horizontal_max_2d(), scale_comm::comm_horizontal_mean(), scale_comm::comm_horizontal_min_2d(), scale_comm::comm_vars8_init(), scale_comm::comm_vars_init(), scale_fileio::fileio_close(), scale_fileio::fileio_create(), scale_fileio::fileio_def_var(), scale_fileio::fileio_enddef(), scale_fileio::fileio_flush(), scale_fileio::fileio_open(), scale_fileio::fileio_read_var_1d(), scale_fileio::fileio_read_var_2d(), scale_fileio::fileio_read_var_3d(), scale_fileio::fileio_read_var_4d(), scale_fileio::fileio_write_var_1d(), scale_fileio::fileio_write_var_2d(), scale_fileio::fileio_write_var_3d(), scale_fileio::fileio_write_var_3d_t(), scale_fileio::fileio_write_var_4d(), scale_history::hist_get_2d(), scale_history::hist_get_3d(), scale_history::hist_in_0d(), scale_history::hist_put_0d(), scale_history::hist_put_1d(), scale_history::hist_put_2d(), scale_history::hist_put_3d(), scale_history::hist_query(), scale_history::hist_reg(), scale_history::hist_setup(), scale_history::hist_write(), mod_land_driver::land_driver(), mod_land_driver::land_driver_resume(), mod_land_phy_driver::land_phy_driver_resume(), scale_monitor::monit_write(), scale_atmos_phy_mp_sn14::mp_negativefilter(), scale_atmos_phy_mp_suzuki10::mp_suzuki10(), scale_grid_nest::nest_comm_nestdown(), scale_grid_nest::nest_comm_recvwait_issue(), scale_grid_nest::nest_comm_test(), mod_ocean_driver::ocean_driver(), mod_ocean_driver::ocean_driver_resume(), mod_ocean_phy_driver::ocean_phy_driver_resume(), scale_atmos_phy_mp_suzuki10::r_collcoag(), mod_rm_driver::scalerm(), mod_rm_prep::scalerm_prep(), scale_rm_statistics::stat_detail(), scale_rm_statistics::stat_total_2d(), scale_rm_statistics::stat_total_3d(), mod_urban_driver::urban_driver(), mod_urban_driver::urban_driver_resume(), mod_urban_phy_driver::urban_phy_driver_resume(), scale_debug::valcheck_1d(), scale_debug::valcheck_2d(), scale_debug::valcheck_3d(), and scale_comm::wait_3d_mpi().

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
real(dp) function, public prc_mpitime()
Get MPI time.
module PROCESS
subroutine, public prc_mpibarrier
Barrier MPI.
Here is the call graph for this function:

◆ prof_rapend()

subroutine, public scale_prof::prof_rapend ( character(len=*), intent(in)  rapname_base,
integer, intent(in), optional  level 
)

Save raptime.

Parameters
[in]rapname_basename of item
[in]levellevel of item

Definition at line 204 of file scale_prof.F90.

References scale_process::prc_mpibarrier(), and scale_process::prc_mpitime().

Referenced by scale_atmos_phy_ae_kajino13::aerosol_zerochem(), mod_atmos_driver::atmos_driver(), mod_atmos_driver::atmos_driver_resume1(), mod_atmos_driver::atmos_driver_resume2(), mod_atmos_driver::atmos_driver_setup(), scale_atmos_dyn::atmos_dyn(), scale_atmos_dyn_common::atmos_dyn_divergence(), scale_atmos_dyn_common::atmos_dyn_numfilter_coef(), scale_atmos_dyn_common::atmos_dyn_numfilter_coef_q(), scale_atmos_dyn_tinteg_large_rk3::atmos_dyn_tinteg_large_rk3(), scale_atmos_dyn_tinteg_short_rk3::atmos_dyn_tinteg_short_rk3(), scale_atmos_dyn_tinteg_short_rk4::atmos_dyn_tinteg_short_rk4(), scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve(), mod_atmos_phy_ae_driver::atmos_phy_ae_driver_resume(), mod_atmos_phy_ch_driver::atmos_phy_ch_driver_resume(), mod_atmos_phy_cp_driver::atmos_phy_cp_driver_resume(), scale_atmos_phy_cp_kf::atmos_phy_cp_kf(), mod_atmos_phy_mp_driver::atmos_phy_mp_driver_resume(), scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler(), scale_atmos_phy_mp_common::atmos_phy_mp_negative_fixer(), scale_atmos_phy_mp_common::atmos_phy_mp_precipitation(), scale_atmos_phy_mp_common::atmos_phy_mp_saturation_adjustment(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08(), mod_atmos_phy_rd_driver::atmos_phy_rd_driver_resume(), scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx(), mod_atmos_phy_sf_driver::atmos_phy_sf_driver_resume(), mod_atmos_phy_tb_driver::atmos_phy_tb_driver_resume(), mod_atmos_vars::atmos_vars_restart_check(), scale_atmos_dyn_common::calc_numdiff(), scale_debug::check(), scale_comm::comm_bcast_1d(), scale_comm::comm_bcast_2d(), scale_comm::comm_bcast_3d(), scale_comm::comm_bcast_4d(), scale_comm::comm_bcast_int_1d(), scale_comm::comm_bcast_int_2d(), scale_comm::comm_bcast_int_scr(), scale_comm::comm_bcast_logical_scr(), scale_comm::comm_bcast_scr(), scale_comm::comm_horizontal_max_2d(), scale_comm::comm_horizontal_mean(), scale_comm::comm_horizontal_min_2d(), scale_comm::comm_vars8_init(), scale_comm::comm_vars_init(), scale_fileio::fileio_close(), scale_fileio::fileio_create(), scale_fileio::fileio_def_var(), scale_fileio::fileio_enddef(), scale_fileio::fileio_flush(), scale_fileio::fileio_open(), scale_fileio::fileio_read_var_1d(), scale_fileio::fileio_read_var_2d(), scale_fileio::fileio_read_var_3d(), scale_fileio::fileio_read_var_4d(), scale_fileio::fileio_write_var_1d(), scale_fileio::fileio_write_var_2d(), scale_fileio::fileio_write_var_3d(), scale_fileio::fileio_write_var_3d_t(), scale_fileio::fileio_write_var_4d(), scale_history::hist_get_2d(), scale_history::hist_get_3d(), scale_history::hist_in_0d(), scale_history::hist_put_0d(), scale_history::hist_put_1d(), scale_history::hist_put_2d(), scale_history::hist_put_3d(), scale_history::hist_query(), scale_history::hist_reg(), scale_history::hist_setup(), scale_history::hist_write(), mod_land_driver::land_driver(), mod_land_driver::land_driver_resume(), mod_land_phy_driver::land_phy_driver_resume(), scale_monitor::monit_write(), scale_atmos_phy_mp_sn14::mp_negativefilter(), scale_atmos_phy_mp_suzuki10::mp_suzuki10(), scale_grid_nest::nest_comm_nestdown(), scale_grid_nest::nest_comm_recvwait_issue(), scale_grid_nest::nest_comm_test(), mod_ocean_driver::ocean_driver(), mod_ocean_driver::ocean_driver_resume(), mod_ocean_phy_driver::ocean_phy_driver_resume(), scale_atmos_phy_mp_suzuki10::r_collcoag(), mod_rm_driver::scalerm(), mod_rm_prep::scalerm_prep(), scale_rm_statistics::stat_detail(), scale_rm_statistics::stat_total_2d(), scale_rm_statistics::stat_total_3d(), mod_urban_driver::urban_driver(), mod_urban_driver::urban_driver_resume(), mod_urban_phy_driver::urban_phy_driver_resume(), scale_debug::valcheck_1d(), scale_debug::valcheck_2d(), scale_debug::valcheck_3d(), and scale_comm::wait_3d_mpi().

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
real(dp) function, public prc_mpitime()
Get MPI time.
module PROCESS
subroutine, public prc_mpibarrier
Barrier MPI.
Here is the call graph for this function:

◆ prof_rapreport()

subroutine, public scale_prof::prof_rapreport ( )

Report raptime.

Definition at line 246 of file scale_prof.F90.

References scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_log_allnode, scale_stdio::io_log_suppress, scale_process::prc_ismaster, scale_process::prc_mpitimestat(), and scale_process::prc_nprocs.

Referenced by mod_rm_driver::scalerm(), and mod_rm_prep::scalerm_prep().

246  use scale_process, only: &
247  prc_mpitimestat, &
248  prc_ismaster
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
subroutine, public prc_mpitimestat(avgvar, maxvar, minvar, maxidx, minidx, var)
Calc global statistics for timer.
module PROCESS
Here is the call graph for this function:
Here is the caller graph for this function: