SCALE-RM
Functions/Subroutines | Variables
mod_admin_time Module Reference

module ADMIN TIME More...

Functions/Subroutines

subroutine, public admin_time_setup (setup_TimeIntegration)
 Setup. More...
 
subroutine, public admin_time_checkstate
 Evaluate component execution. More...
 
subroutine, public admin_time_advance
 Advance the time & evaluate restart & stop. More...
 

Variables

real(dp), public time_dtsec_atmos_restart
 time interval of atmosphere restart [sec] More...
 
real(dp), public time_dtsec_ocean_restart
 time interval of ocean restart [sec] More...
 
real(dp), public time_dtsec_land_restart
 time interval of land restart [sec] More...
 
real(dp), public time_dtsec_urban_restart
 time interval of urban restart [sec] More...
 
real(dp), public time_dtsec_resume
 time interval for resume [sec] More...
 
integer, public time_dstep_atmos_restart
 interval of atmosphere restart [step] More...
 
integer, public time_dstep_ocean_restart
 interval of ocean restart [step] More...
 
integer, public time_dstep_land_restart
 interval of land restart [step] More...
 
integer, public time_dstep_urban_restart
 interval of urban restart [step] More...
 
integer, public time_dstep_resume
 interval for resume [step] More...
 
logical, public time_doatmos_step
 execute atmosphere component in this step? More...
 
logical, public time_doatmos_dyn
 execute dynamics in this step? More...
 
logical, public time_doatmos_phy_cp
 execute physics in this step? (cumulus ) More...
 
logical, public time_doatmos_phy_mp
 execute physics in this step? (microphysics) More...
 
logical, public time_doatmos_phy_rd
 execute physics in this step? (radiation ) More...
 
logical, public time_doatmos_phy_sf
 execute physics in this step? (surface flux) More...
 
logical, public time_doatmos_phy_tb
 execute physics in this step? (turbulence ) More...
 
logical, public time_doatmos_phy_ch
 execute physics in this step? (chemistry ) More...
 
logical, public time_doatmos_phy_ae
 execute physics in this step? (aerosol ) More...
 
logical, public time_doatmos_restart
 execute atmosphere restart output in this step? More...
 
logical, public time_doocean_step
 execute ocean component in this step? More...
 
logical, public time_doocean_restart
 execute ocean restart output in this step? More...
 
logical, public time_doland_step
 execute land component in this step? More...
 
logical, public time_doland_restart
 execute land restart output in this step? More...
 
logical, public time_dourban_step
 execute urban component in this step? More...
 
logical, public time_dourban_restart
 execute urban restart output in this step? More...
 
logical, public time_doresume
 resume in this step? More...
 
logical, public time_doend
 finish program in this step? More...
 

Detailed Description

module ADMIN TIME

Description
Time management for SCALE-RM
Author
Team SCALE
NAMELIST
  • PARAM_TIME
    nametypedefault valuecomment
    TIME_STARTDATE integer (/ -999, 1, 1, 0, 0, 0 /)
    TIME_STARTMS real(DP) 0.0_DP [millisec]
    TIME_DURATION real(DP) UNDEF8
    TIME_DURATION_UNIT character(len=H_SHORT) "SEC"
    TIME_DT real(DP) UNDEF8
    TIME_DT_UNIT character(len=H_SHORT) "SEC"
    TIME_DT_ATMOS_DYN real(DP) UNDEF8
    TIME_DT_ATMOS_DYN_UNIT character(len=H_SHORT) "SEC"
    TIME_NSTEP_ATMOS_DYN integer -1
    TIME_DT_ATMOS_PHY_CP real(DP) UNDEF8
    TIME_DT_ATMOS_PHY_CP_UNIT character(len=H_SHORT) ""
    TIME_DT_ATMOS_PHY_MP real(DP) UNDEF8
    TIME_DT_ATMOS_PHY_MP_UNIT character(len=H_SHORT) ""
    TIME_DT_ATMOS_PHY_RD real(DP) UNDEF8
    TIME_DT_ATMOS_PHY_RD_UNIT character(len=H_SHORT) ""
    TIME_DT_ATMOS_PHY_SF real(DP) UNDEF8
    TIME_DT_ATMOS_PHY_SF_UNIT character(len=H_SHORT) ""
    TIME_DT_ATMOS_PHY_TB real(DP) UNDEF8
    TIME_DT_ATMOS_PHY_TB_UNIT character(len=H_SHORT) ""
    TIME_DT_ATMOS_PHY_CH real(DP) UNDEF8
    TIME_DT_ATMOS_PHY_CH_UNIT character(len=H_SHORT) ""
    TIME_DT_ATMOS_PHY_AE real(DP) UNDEF8
    TIME_DT_ATMOS_PHY_AE_UNIT character(len=H_SHORT) ""
    TIME_DT_ATMOS_RESTART real(DP) UNDEF8
    TIME_DT_ATMOS_RESTART_UNIT character(len=H_SHORT) ""
    TIME_DT_OCEAN real(DP) UNDEF8
    TIME_DT_OCEAN_UNIT character(len=H_SHORT) ""
    TIME_DT_OCEAN_RESTART real(DP) UNDEF8
    TIME_DT_OCEAN_RESTART_UNIT character(len=H_SHORT) ""
    TIME_DT_LAND real(DP) UNDEF8
    TIME_DT_LAND_UNIT character(len=H_SHORT) ""
    TIME_DT_LAND_RESTART real(DP) UNDEF8
    TIME_DT_LAND_RESTART_UNIT character(len=H_SHORT) ""
    TIME_DT_URBAN real(DP) UNDEF8
    TIME_DT_URBAN_UNIT character(len=H_SHORT) ""
    TIME_DT_URBAN_RESTART real(DP) UNDEF8
    TIME_DT_URBAN_RESTART_UNIT character(len=H_SHORT) ""
    TIME_DT_WALLCLOCK_CHECK real(DP) UNDEF8
    TIME_DT_WALLCLOCK_CHECK_UNIT character(len=H_SHORT) ""
    TIME_DT_RESUME real(DP) UNDEF8
    TIME_DT_RESUME_UNIT character(len=H_SHORT) ""
    TIME_WALLCLOCK_LIMIT real(DP) -1.0_DP Elapse time limit of wall clock time [sec]
    TIME_WALLCLOCK_SAFE real(DP) 0.9_DP Safety coefficient for elapse time limit

History Output
No history output

Function/Subroutine Documentation

◆ admin_time_setup()

subroutine, public mod_admin_time::admin_time_setup ( logical, intent(in)  setup_TimeIntegration)

Setup.

Definition at line 111 of file mod_admin_time.f90.

References mod_atmos_vars::atmos_restart_in_basename, scale_calendar::calendar_adjust_daysec(), scale_calendar::calendar_cfunits2sec(), scale_calendar::calendar_combine_daysec(), scale_calendar::calendar_date2char(), scale_calendar::calendar_date2daysec(), scale_calendar::calendar_daysec2date(), scale_calendar::calendar_unit2sec(), scale_const::const_undef8, gtool_file::filegetdatainfo(), scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, scale_process::prc_mpistop(), scale_process::prc_mpitime(), scale_process::prc_myrank, scale_time::time_dstep_atmos_dyn, scale_time::time_dstep_atmos_phy_ae, scale_time::time_dstep_atmos_phy_ch, scale_time::time_dstep_atmos_phy_cp, scale_time::time_dstep_atmos_phy_mp, scale_time::time_dstep_atmos_phy_rd, scale_time::time_dstep_atmos_phy_sf, scale_time::time_dstep_atmos_phy_tb, time_dstep_atmos_restart, scale_time::time_dstep_land, time_dstep_land_restart, scale_time::time_dstep_ocean, time_dstep_ocean_restart, time_dstep_resume, scale_time::time_dstep_urban, time_dstep_urban_restart, scale_time::time_dstep_wallclock_check, scale_time::time_dtsec, scale_time::time_dtsec_atmos_dyn, scale_time::time_dtsec_atmos_phy_ae, scale_time::time_dtsec_atmos_phy_ch, scale_time::time_dtsec_atmos_phy_cp, scale_time::time_dtsec_atmos_phy_mp, scale_time::time_dtsec_atmos_phy_rd, scale_time::time_dtsec_atmos_phy_sf, scale_time::time_dtsec_atmos_phy_tb, time_dtsec_atmos_restart, scale_time::time_dtsec_land, time_dtsec_land_restart, scale_time::time_dtsec_ocean, time_dtsec_ocean_restart, time_dtsec_resume, scale_time::time_dtsec_urban, time_dtsec_urban_restart, scale_time::time_dtsec_wallclock_check, scale_time::time_nowdate, scale_time::time_nowday, scale_time::time_nowdaysec, scale_time::time_nowms, scale_time::time_nowsec, scale_time::time_nowstep, scale_time::time_nstep, scale_time::time_nstep_atmos_dyn, scale_time::time_offset_year, and scale_time::time_startdaysec.

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

111  use gtool_file, only: &
113  use scale_process, only: &
114  prc_myrank, &
115  prc_mpistop, &
117  use scale_const, only: &
118  undef8 => const_undef8
119  use scale_calendar, only: &
127  use scale_time, only: &
128  time_dtsec, &
129  time_nowdate, &
130  time_nowms, &
131  time_nowday, &
132  time_nowsec, &
133  time_nowdaysec, &
134  time_nowstep, &
135  time_nstep, &
137  nstep_dyn => time_nstep_atmos_dyn, &
146  time_dtsec_land, &
158  time_dstep_land, &
163  use mod_atmos_vars, only: &
164  restart_in_basename => atmos_restart_in_basename
165  implicit none
166 
167  logical, intent(in) :: setup_timeintegration
168 
169  real(DP) :: time_duration = undef8
170  character(len=H_SHORT) :: time_duration_unit = "SEC"
171  real(DP) :: time_dt = undef8
172  character(len=H_SHORT) :: time_dt_unit = "SEC"
173 
174  real(DP) :: time_dt_atmos_dyn = undef8
175  character(len=H_SHORT) :: time_dt_atmos_dyn_unit = "SEC"
176  integer :: time_nstep_atmos_dyn = -1
177  real(DP) :: time_dt_atmos_phy_cp = undef8
178  character(len=H_SHORT) :: time_dt_atmos_phy_cp_unit = ""
179  real(DP) :: time_dt_atmos_phy_mp = undef8
180  character(len=H_SHORT) :: time_dt_atmos_phy_mp_unit = ""
181  real(DP) :: time_dt_atmos_phy_rd = undef8
182  character(len=H_SHORT) :: time_dt_atmos_phy_rd_unit = ""
183  real(DP) :: time_dt_atmos_phy_sf = undef8
184  character(len=H_SHORT) :: time_dt_atmos_phy_sf_unit = ""
185  real(DP) :: time_dt_atmos_phy_tb = undef8
186  character(len=H_SHORT) :: time_dt_atmos_phy_tb_unit = ""
187  real(DP) :: time_dt_atmos_phy_ch = undef8
188  character(len=H_SHORT) :: time_dt_atmos_phy_ch_unit = ""
189  real(DP) :: time_dt_atmos_phy_ae = undef8
190  character(len=H_SHORT) :: time_dt_atmos_phy_ae_unit = ""
191  real(DP) :: time_dt_atmos_restart = undef8
192  character(len=H_SHORT) :: time_dt_atmos_restart_unit = ""
193 
194  real(DP) :: time_dt_ocean = undef8
195  character(len=H_SHORT) :: time_dt_ocean_unit = ""
196  real(DP) :: time_dt_ocean_restart = undef8
197  character(len=H_SHORT) :: time_dt_ocean_restart_unit = ""
198 
199  real(DP) :: time_dt_land = undef8
200  character(len=H_SHORT) :: time_dt_land_unit = ""
201  real(DP) :: time_dt_land_restart = undef8
202  character(len=H_SHORT) :: time_dt_land_restart_unit = ""
203 
204  real(DP) :: time_dt_urban = undef8
205  character(len=H_SHORT) :: time_dt_urban_unit = ""
206  real(DP) :: time_dt_urban_restart = undef8
207  character(len=H_SHORT) :: time_dt_urban_restart_unit = ""
208 
209  real(DP) :: time_dt_resume = undef8
210  character(len=H_SHORT) :: time_dt_resume_unit = ""
211 
212  real(DP) :: time_dt_wallclock_check = undef8
213  character(len=H_SHORT) :: time_dt_wallclock_check_unit = ""
214 
215  namelist / param_time / &
216  time_startdate, &
217  time_startms, &
218  time_duration, &
219  time_duration_unit, &
220  time_dt, &
221  time_dt_unit, &
222  time_dt_atmos_dyn, &
223  time_dt_atmos_dyn_unit, &
225  time_dt_atmos_phy_cp, &
226  time_dt_atmos_phy_cp_unit, &
227  time_dt_atmos_phy_mp, &
228  time_dt_atmos_phy_mp_unit, &
229  time_dt_atmos_phy_rd, &
230  time_dt_atmos_phy_rd_unit, &
231  time_dt_atmos_phy_sf, &
232  time_dt_atmos_phy_sf_unit, &
233  time_dt_atmos_phy_tb, &
234  time_dt_atmos_phy_tb_unit, &
235  time_dt_atmos_phy_ch, &
236  time_dt_atmos_phy_ch_unit, &
237  time_dt_atmos_phy_ae, &
238  time_dt_atmos_phy_ae_unit, &
239  time_dt_atmos_restart, &
240  time_dt_atmos_restart_unit, &
241  time_dt_ocean, &
242  time_dt_ocean_unit, &
243  time_dt_ocean_restart, &
244  time_dt_ocean_restart_unit, &
245  time_dt_land, &
246  time_dt_land_unit, &
247  time_dt_land_restart, &
248  time_dt_land_restart_unit, &
249  time_dt_urban, &
250  time_dt_urban_unit, &
251  time_dt_urban_restart, &
252  time_dt_urban_restart_unit, &
253  time_dt_wallclock_check, &
254  time_dt_wallclock_check_unit, &
255  time_dt_resume, &
256  time_dt_resume_unit, &
257  time_wallclock_limit, &
258  time_wallclock_safe
259 
260  integer :: dateday
261  real(DP) :: datesec
262  real(DP) :: cftime
263  character(len=H_MID) :: cfunits
264 
265  real(DP) :: time_durationsec
266  character(len=27) :: startchardate
267  character(len=27) :: endchardate
268 
269  integer :: ierr
270  !---------------------------------------------------------------------------
271 
272  if( io_l ) write(io_fid_log,*)
273  if( io_l ) write(io_fid_log,*) '++++++ Module[TIME] / Categ[COMMON] / Origin[SCALElib]'
274 
276 
277  !--- read namelist
278  rewind(io_fid_conf)
279  read(io_fid_conf,nml=param_time,iostat=ierr)
280  if( ierr < 0 ) then !--- missing
281  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
282  elseif( ierr > 0 ) then !--- fatal error
283  write(*,*) 'xxx Not appropriate names in namelist PARAM_TIME. Check!'
284  call prc_mpistop
285  endif
286  if( io_lnml ) write(io_fid_log,nml=param_time)
287 
288  ! check time setting
289  if ( setup_timeintegration ) then
290  if ( time_dt == undef8 ) then
291  write(*,*) 'xxx Not found TIME_DT.'
292  call prc_mpistop
293  endif
294  if ( time_duration == undef8 ) then
295  write(*,*) 'xxx Not found TIME_DURATION.'
296  call prc_mpistop
297  endif
298 
299  ! DYN
300  if ( time_dt_atmos_dyn == undef8 ) then
301  if ( time_nstep_atmos_dyn < 0 ) then
302  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_ATMOS_DYN. TIME_DT is used.'
303  time_dt_atmos_dyn = time_dt
304  endif
305  endif
306  if ( time_dt_atmos_dyn_unit == '' ) then
307  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_ATMOS_DYN_UNIT. TIME_DT_UNIT is used.'
308  time_dt_atmos_dyn_unit = time_dt_unit
309  endif
310  ! PHY_CP
311  if ( time_dt_atmos_phy_cp == undef8 ) then
312  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_ATMOS_PHY_CP. TIME_DT is used.'
313  time_dt_atmos_phy_cp = time_dt
314  endif
315  if ( time_dt_atmos_phy_cp_unit == '' ) then
316  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_ATMOS_PHY_CP_UNIT. TIME_DT_UNIT is used.'
317  time_dt_atmos_phy_cp_unit = time_dt_unit
318  endif
319  ! PHY_MP
320  if ( time_dt_atmos_phy_mp == undef8 ) then
321  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_ATMOS_PHY_MP. TIME_DT is used.'
322  time_dt_atmos_phy_mp = time_dt
323  endif
324  if ( time_dt_atmos_phy_mp_unit == '' ) then
325  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_ATMOS_PHY_MP_UNIT. TIME_DT_UNIT is used.'
326  time_dt_atmos_phy_mp_unit = time_dt_unit
327  endif
328  ! PHY_RD
329  if ( time_dt_atmos_phy_rd == undef8 ) then
330  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_ATMOS_PHY_RD. TIME_DT is used.'
331  time_dt_atmos_phy_rd = time_dt
332  endif
333  if ( time_dt_atmos_phy_rd_unit == '' ) then
334  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_ATMOS_PHY_RD_UNIT. TIME_DT_UNIT is used.'
335  time_dt_atmos_phy_rd_unit = time_dt_unit
336  endif
337  ! PHY_SF
338  if ( time_dt_atmos_phy_sf == undef8 ) then
339  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_ATMOS_PHY_SF. TIME_DT is used.'
340  time_dt_atmos_phy_sf = time_dt
341  endif
342  if ( time_dt_atmos_phy_sf_unit == '' ) then
343  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_ATMOS_PHY_SF_UNIT. TIME_DT_UNIT is used.'
344  time_dt_atmos_phy_sf_unit = time_dt_unit
345  endif
346  ! PHY_TB
347  if ( time_dt_atmos_phy_tb == undef8 ) then
348  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_ATMOS_PHY_TB. TIME_DT is used.'
349  time_dt_atmos_phy_tb = time_dt
350  endif
351  if ( time_dt_atmos_phy_tb_unit == '' ) then
352  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_ATMOS_PHY_TB_UNIT. TIME_DT_UNIT is used.'
353  time_dt_atmos_phy_tb_unit = time_dt_unit
354  endif
355  ! PHY_CH
356  if ( time_dt_atmos_phy_ch == undef8 ) then
357  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_ATMOS_PHY_CH. TIME_DT is used.'
358  time_dt_atmos_phy_ch = time_dt
359  endif
360  if ( time_dt_atmos_phy_ch_unit == '' ) then
361  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_ATMOS_PHY_CH_UNIT. TIME_DT_UNIT is used.'
362  time_dt_atmos_phy_ch_unit = time_dt_unit
363  endif
364  ! PHY_AE
365  if ( time_dt_atmos_phy_ae == undef8 ) then
366  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_ATMOS_PHY_AE. TIME_DT is used.'
367  time_dt_atmos_phy_ae = time_dt
368  endif
369  if ( time_dt_atmos_phy_ae_unit == '' ) then
370  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_ATMOS_PHY_AE_UNIT. TIME_DT_UNIT is used.'
371  time_dt_atmos_phy_ae_unit = time_dt_unit
372  endif
373  ! ATMOS RESTART
374  if ( time_dt_atmos_restart == undef8 ) then
375  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_ATMOS_RESTART. TIME_DURATION is used.'
376  time_dt_atmos_restart = time_duration
377  endif
378  if ( time_dt_atmos_restart_unit == '' ) then
379  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_ATMOS_RESTART_UNIT. TIME_DURATION_UNIT is used.'
380  time_dt_atmos_restart_unit = time_duration_unit
381  endif
382  ! OCEAN
383  if ( time_dt_ocean == undef8 ) then
384  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_OCEAN. TIME_DT is used.'
385  time_dt_ocean = time_dt
386  endif
387  if ( time_dt_ocean_unit == '' ) then
388  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_OCEAN_UNIT. TIME_DT_UNIT is used.'
389  time_dt_ocean_unit = time_dt_unit
390  endif
391  ! OCEAN RESTART
392  if ( time_dt_ocean_restart == undef8 ) then
393  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_OCEAN_RESTART. TIME_DURATION is used.'
394  time_dt_ocean_restart = time_duration
395  endif
396  if ( time_dt_ocean_restart_unit == '' ) then
397  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_OCEAN_RESTART_UNIT. TIME_DURATION_UNIT is used.'
398  time_dt_ocean_restart_unit = time_duration_unit
399  endif
400  ! LAND
401  if ( time_dt_land == undef8 ) then
402  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_LAND. TIME_DT is used.'
403  time_dt_land = time_dt
404  endif
405  if ( time_dt_land_unit == '' ) then
406  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_LAND_UNIT. TIME_DT_UNIT is used.'
407  time_dt_land_unit = time_dt_unit
408  endif
409  ! LAND RESTART
410  if ( time_dt_land_restart == undef8 ) then
411  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_LAND_RESTART. TIME_DURATION is used.'
412  time_dt_land_restart = time_duration
413  endif
414  if ( time_dt_land_restart_unit == '' ) then
415  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_LAND_RESTART_UNIT. TIME_DURATION_UNIT is used.'
416  time_dt_land_restart_unit = time_duration_unit
417  endif
418  ! URBAN
419  if ( time_dt_urban == undef8 ) then
420  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_URBAN. TIME_DT is used.'
421  time_dt_urban = time_dt
422  endif
423  if ( time_dt_urban_unit == '' ) then
424  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_URBAN_UNIT. TIME_DT_UNIT is used.'
425  time_dt_urban_unit = time_dt_unit
426  endif
427  ! URBAN RESTART
428  if ( time_dt_urban_restart == undef8 ) then
429  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_URBAN_RESTART. TIME_DURATION is used.'
430  time_dt_urban_restart = time_duration
431  endif
432  if ( time_dt_urban_restart_unit == '' ) then
433  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_URBAN_RESTART_UNIT. TIME_DURATION_UNIT is used.'
434  time_dt_urban_restart_unit = time_duration_unit
435  endif
436  ! Resume
437  if ( time_dt_resume == undef8 ) then
438  time_dt_resume = time_duration
439  endif
440  if ( time_dt_resume_unit == '' ) then
441  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_RESUME_UNIT. TIME_DURATION_UNIT is used.'
442  time_dt_resume_unit = time_duration_unit
443  endif
444  endif
445 
446  !--- calculate time
447  if ( time_startdate(1) == -999 ) then
448  if ( restart_in_basename /= '' ) then ! read start time from the restart data
449  call filegetdatainfo( restart_in_basename, & ! [IN]
450  'DENS', & ! [IN]
451  prc_myrank, & ! [IN]
452  0, & ! [IN] step
453  time_start = cftime, & ! [OUT]
454  time_units = cfunits ) ! [OUT]
455 
456  dateday = 0
457  datesec = calendar_cfunits2sec( cftime, cfunits, 0 )
458 
459  call calendar_adjust_daysec( dateday, datesec )
460 
461  call calendar_daysec2date( time_startdate, & ! [OUT]
462  time_startms, & ! [OUT]
463  dateday, & ! [IN]
464  datesec, & ! [IN]
465  0 ) ! [IN]
466  else
467  time_startdate = (/ 0, 1, 1, 0, 0, 0 /)
468  time_startms = 0.0_dp
469  endif
470  else
471  time_startms = time_startms * 1.e-3_dp
472  endif
473 
474  time_offset_year = time_startdate(1)
475 
476  call calendar_date2daysec( time_startday, & ! [OUT]
477  time_startsec, & ! [OUT]
478  time_startdate(:), & ! [IN]
479  time_startms, & ! [IN]
480  time_offset_year ) ! [IN]
481 
482  call calendar_date2char( startchardate, & ! [OUT]
483  time_startdate(:), & ! [IN]
484  time_startms ) ! [IN]
485 
486  time_startdaysec = calendar_combine_daysec( time_startday, time_startsec )
487 
488  time_nowdate(:) = time_startdate(:)
489  time_nowms = time_startms
490  time_nowday = time_startday
491  time_nowsec = time_startsec
493 
494  time_endday = time_startday
495 
496  if ( setup_timeintegration ) then
497  call calendar_unit2sec( time_durationsec, time_duration, time_duration_unit )
498  time_endsec = time_startsec + time_durationsec
499  else
500  time_endsec = time_startsec
501  endif
502 
503  call calendar_adjust_daysec( time_endday, time_endsec ) ! [INOUT]
504 
505  call calendar_daysec2date( time_enddate(:), & ! [OUT]
506  time_endms, & ! [OUT]
507  time_endday, & ! [IN]
508  time_endsec, & ! [IN]
509  time_offset_year ) ! [IN]
510 
511  call calendar_date2char( endchardate, & ! [OUT]
512  time_enddate(:), & ! [IN]
513  time_endms ) ! [IN]
514 
515  if( io_l ) write(io_fid_log,*)
516  if( io_l ) write(io_fid_log,*) '*** Date/time setting ***'
517  if( io_l ) write(io_fid_log,'(1x,A,A)') '*** START Date : ', startchardate
518  if( io_l ) write(io_fid_log,'(1x,A,A)') '*** END Date : ', endchardate
519 
520  if ( setup_timeintegration ) then
521 
522  call calendar_unit2sec( time_dtsec, time_dt, time_dt_unit )
523 
524  time_nstep = int( time_durationsec / time_dtsec )
525  time_nowstep = 1
526 
527  if( io_l ) write(io_fid_log,'(1x,A,F10.3)') '*** delta t (sec.) :', time_dtsec
528  if( io_l ) write(io_fid_log,'(1x,A,I10)' ) '*** No. of steps :', time_nstep
529 
530  !--- calculate intervals for atmosphere
531  if ( time_dt_atmos_dyn /= undef8 ) then
532  call calendar_unit2sec( time_dtsec_atmos_dyn, time_dt_atmos_dyn, time_dt_atmos_dyn_unit )
533  else
535  endif
536  call calendar_unit2sec( time_dtsec_atmos_phy_cp, time_dt_atmos_phy_cp, time_dt_atmos_phy_cp_unit )
537  call calendar_unit2sec( time_dtsec_atmos_phy_mp, time_dt_atmos_phy_mp, time_dt_atmos_phy_mp_unit )
538  call calendar_unit2sec( time_dtsec_atmos_phy_rd, time_dt_atmos_phy_rd, time_dt_atmos_phy_rd_unit )
539  call calendar_unit2sec( time_dtsec_atmos_phy_sf, time_dt_atmos_phy_sf, time_dt_atmos_phy_sf_unit )
540  call calendar_unit2sec( time_dtsec_atmos_phy_tb, time_dt_atmos_phy_tb, time_dt_atmos_phy_tb_unit )
541  call calendar_unit2sec( time_dtsec_atmos_phy_ch, time_dt_atmos_phy_ch, time_dt_atmos_phy_ch_unit )
542  call calendar_unit2sec( time_dtsec_atmos_phy_ae, time_dt_atmos_phy_ae, time_dt_atmos_phy_ae_unit )
543  call calendar_unit2sec( time_dtsec_atmos_restart, time_dt_atmos_restart, time_dt_atmos_restart_unit )
544  call calendar_unit2sec( time_dtsec_ocean, time_dt_ocean, time_dt_ocean_unit )
545  call calendar_unit2sec( time_dtsec_ocean_restart, time_dt_ocean_restart, time_dt_ocean_restart_unit )
546  call calendar_unit2sec( time_dtsec_land, time_dt_land, time_dt_land_unit )
547  call calendar_unit2sec( time_dtsec_land_restart, time_dt_land_restart, time_dt_land_restart_unit )
548  call calendar_unit2sec( time_dtsec_urban, time_dt_urban, time_dt_urban_unit )
549  call calendar_unit2sec( time_dtsec_urban_restart, time_dt_urban_restart, time_dt_urban_restart_unit )
550  call calendar_unit2sec( time_dtsec_resume, time_dt_resume, time_dt_resume_unit )
551 
553  nstep_dyn = time_nstep_atmos_dyn
554 
563  time_dtsec_atmos_restart = max( time_dtsec_atmos_restart, time_dtsec_atmos_dyn*time_nstep_atmos_dyn )
565  time_dtsec_ocean_restart = max( time_dtsec_ocean_restart, time_dtsec_atmos_dyn*time_nstep_atmos_dyn )
567  time_dtsec_land_restart = max( time_dtsec_land_restart, time_dtsec_atmos_dyn*time_nstep_atmos_dyn )
569  time_dtsec_urban_restart = max( time_dtsec_urban_restart, time_dtsec_atmos_dyn*time_nstep_atmos_dyn )
570  time_dtsec_resume = max( time_dtsec_resume, time_dtsec_atmos_dyn*time_nstep_atmos_dyn )
571 
583  time_dstep_atmos_restart = nint( time_dtsec_atmos_restart / time_dtsec )
584  time_dstep_ocean_restart = nint( time_dtsec_ocean_restart / time_dtsec )
585  time_dstep_land_restart = nint( time_dtsec_land_restart / time_dtsec )
586  time_dstep_urban_restart = nint( time_dtsec_urban_restart / time_dtsec )
587  time_dstep_resume = nint( time_dtsec_resume / time_dtsec )
588 
589  time_res_resume = time_dstep_resume - 1
590 
591  if ( abs( real(time_nstep_atmos_dyn,kind=dp)*time_dtsec_atmos_dyn &
592  - real(time_dstep_atmos_dyn,kind=dp)*time_dtsec ) > eps ) then
593  write(*,*) 'xxx delta t(ATMOS_DYN) must be a multiple of delta t ', &
595  call prc_mpistop
596  endif
597  if ( abs(time_dtsec_atmos_phy_cp-real(time_dstep_atmos_phy_cp,kind=dp)*time_dtsec) > eps ) then
598  write(*,*) 'xxx delta t(ATMOS_PHY_CP) must be a multiple of delta t ', &
600  call prc_mpistop
601  endif
602  if ( abs(time_dtsec_atmos_phy_mp-real(time_dstep_atmos_phy_mp,kind=dp)*time_dtsec) > eps ) then
603  write(*,*) 'xxx delta t(ATMOS_PHY_MP) must be a multiple of delta t ', &
605  call prc_mpistop
606  endif
607  if ( abs(time_dtsec_atmos_phy_rd-real(time_dstep_atmos_phy_rd,kind=dp)*time_dtsec) > eps ) then
608  write(*,*) 'xxx delta t(ATMOS_PHY_RD) must be a multiple of delta t ', &
610  call prc_mpistop
611  endif
612  if ( abs(time_dtsec_atmos_phy_sf-real(time_dstep_atmos_phy_sf,kind=dp)*time_dtsec) > eps ) then
613  write(*,*) 'xxx delta t(ATMOS_PHY_SF) must be a multiple of delta t ', &
615  call prc_mpistop
616  endif
617  if ( abs(time_dtsec_atmos_phy_tb-real(time_dstep_atmos_phy_tb,kind=dp)*time_dtsec) > eps ) then
618  write(*,*) 'xxx delta t(ATMOS_PHY_TB) must be a multiple of delta t ', &
620  call prc_mpistop
621  endif
622  if ( abs(time_dtsec_atmos_phy_ch-real(time_dstep_atmos_phy_ch,kind=dp)*time_dtsec) > eps ) then
623  write(*,*) 'xxx delta t(ATMOS_PHY_CH) must be a multiple of delta t ', &
625  call prc_mpistop
626  endif
627  if ( abs(time_dtsec_atmos_phy_ae-real(time_dstep_atmos_phy_ae,kind=dp)*time_dtsec) > eps ) then
628  write(*,*) 'xxx delta t(ATMOS_PHY_AE) must be a multiple of delta t ', &
630  call prc_mpistop
631  endif
632  if ( abs(time_dtsec_ocean-real(time_dstep_ocean,kind=dp)*time_dtsec) > eps ) then
633  write(*,*) 'xxx delta t(OCEAN) must be a multiple of delta t ', &
635  call prc_mpistop
636  endif
637  if ( abs(time_dtsec_land-real(time_dstep_land,kind=dp)*time_dtsec) > eps ) then
638  write(*,*) 'xxx delta t(LAND) must be a multiple of delta t ', &
640  call prc_mpistop
641  endif
642  if ( abs(time_dtsec_urban-real(time_dstep_urban,kind=dp)*time_dtsec) > eps ) then
643  write(*,*) 'xxx delta t(URBAN) must be a multiple of delta t ', &
645  call prc_mpistop
646  endif
647  if ( abs(time_dtsec_atmos_restart-real(time_dstep_atmos_restart,kind=dp)*time_dtsec) > eps ) then
648  write(*,*) 'xxx delta t(ATMOS_RESTART) must be a multiple of delta t ', &
649  time_dtsec_atmos_restart, real(time_dstep_atmos_restart,kind=dp)*time_dtsec
650  call prc_mpistop
651  endif
652  if ( abs(time_dtsec_ocean_restart-real(time_dstep_ocean_restart,kind=dp)*time_dtsec) > eps ) then
653  write(*,*) 'xxx delta t(OCEAN_RESTART) must be a multiple of delta t ', &
654  time_dtsec_ocean_restart, real(time_dstep_ocean_restart,kind=dp)*time_dtsec
655  call prc_mpistop
656  endif
657  if ( abs(time_dtsec_land_restart-real(time_dstep_land_restart,kind=dp)*time_dtsec) > eps ) then
658  write(*,*) 'xxx delta t(LAND_RESTART) must be a multiple of delta t ', &
659  time_dtsec_land_restart, real(time_dstep_land_restart,kind=dp)*time_dtsec
660  call prc_mpistop
661  endif
662  if ( abs(time_dtsec_urban_restart-real(time_dstep_urban_restart,kind=dp)*time_dtsec) > eps ) then
663  write(*,*) 'xxx delta t(URBAN_RESTART) must be a multiple of delta t ', &
664  time_dtsec_urban_restart, real(time_dstep_urban_restart,kind=dp)*time_dtsec
665  call prc_mpistop
666  endif
667  if ( abs(time_dtsec_urban_restart-real(time_dstep_urban_restart,kind=dp)*time_dtsec) > eps ) then
668  write(*,*) 'xxx delta t(URBAN_RESTART) must be a multiple of delta t ', &
669  time_dtsec_urban_restart, real(time_dstep_urban_restart,kind=dp)*time_dtsec
670  call prc_mpistop
671  endif
672  if ( abs(time_dtsec_resume-real(time_dstep_resume,kind=dp)*time_dtsec) > eps ) then
673  write(*,*) 'xxx delta t(RESUME) must be a multiple of delta t ', &
674  time_dtsec_resume, real(time_dstep_resume,kind=dp)*time_dtsec
675  call prc_mpistop
676  endif
677 
678  if( io_l ) write(io_fid_log,*)
679  if( io_l ) write(io_fid_log,*) '*** Time interval for each processes (sec.)'
680  if( io_l ) write(io_fid_log,*) '*** Atmosphere'
681  if( io_l ) write(io_fid_log,'(1x,A,F10.3)') '*** Dynamics (time) : ', time_dtsec_atmos_dyn
682  if( io_l ) write(io_fid_log,'(1x,A,I10,A,I8,A)') '*** (step) : ', time_nstep_atmos_dyn, &
683  ' (step interval=', time_dstep_atmos_dyn, ')'
684  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Physics, Cumulus : ', time_dtsec_atmos_phy_cp, &
685  ' (step interval=', time_dstep_atmos_phy_cp, ')'
686  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Physics, Cloud Microphysics : ', time_dtsec_atmos_phy_mp, &
687  ' (step interval=', time_dstep_atmos_phy_mp, ')'
688  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Physics, Radiation : ', time_dtsec_atmos_phy_rd, &
689  ' (step interval=', time_dstep_atmos_phy_rd, ')'
690  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Physics, Surface Flux : ', time_dtsec_atmos_phy_sf, &
691  ' (step interval=', time_dstep_atmos_phy_sf, ')'
692  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Physics, Turbulence : ', time_dtsec_atmos_phy_tb, &
693  ' (step interval=', time_dstep_atmos_phy_tb, ')'
694  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Physics, Chemistry : ', time_dtsec_atmos_phy_ch, &
695  ' (step interval=', time_dstep_atmos_phy_ch, ')'
696  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Physics, Aerosol : ', time_dtsec_atmos_phy_ae, &
697  ' (step interval=', time_dstep_atmos_phy_ae, ')'
698  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Ocean : ', time_dtsec_ocean, &
699  ' (step interval=', time_dstep_ocean, ')'
700  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Land : ', time_dtsec_land, &
701  ' (step interval=', time_dstep_land, ')'
702  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Urban : ', time_dtsec_urban, &
703  ' (step interval=', time_dstep_urban, ')'
704  if( io_l ) write(io_fid_log,*)
705  if( io_l ) write(io_fid_log,*) '*** Time interval for restart (sec.)'
706  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Atmospheric Variables : ', time_dtsec_atmos_restart, &
707  ' (step interval=', time_dstep_atmos_restart, ')'
708  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Ocean Variables : ', time_dtsec_ocean_restart, &
709  ' (step interval=', time_dstep_ocean_restart, ')'
710  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Land Variables : ', time_dtsec_land_restart, &
711  ' (step interval=', time_dstep_land_restart, ')'
712  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Urban Variables : ', time_dtsec_urban_restart, &
713  ' (step interval=', time_dstep_urban_restart, ')'
714  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Resume : ', time_dtsec_resume, &
715  ' (step interval=', time_dstep_resume, ')'
716  else
717  time_dtsec = 1.0_rp
718  endif
719 
720  ! WALLCLOCK TERMINATOR SETUP
721  if ( time_wallclock_limit > 0.0_dp ) then
722  if( io_l ) write(io_fid_log,*)
723  if( io_l ) write(io_fid_log,*) '*** Wall clock time limit of execution is specified.'
724 
725  if ( time_dt_wallclock_check == undef8 ) then
726  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_WALLCLOCK_CHECK. largest time step interval is used.'
736  time_dtsec_land, &
738  else
739  if ( time_dt_wallclock_check_unit == '' ) then
740  if( io_l ) write(io_fid_log,*) '*** Not found TIME_DT_WALLCLOCK_CHECK_UNIT. TIME_DURATION_UNIT is used.'
741  time_dt_wallclock_check_unit = time_duration_unit
742  endif
743  call calendar_unit2sec( time_dtsec_wallclock_check, time_dt_wallclock_check, time_dt_wallclock_check_unit )
745  endif
746 
748 
749  time_wallclock_safe = max( min( time_wallclock_safe, 1.0_dp ), 0.0_dp )
750  time_wallclock_start = prc_mpitime()
751 
752  if( io_l ) write(io_fid_log,'(1x,A,F10.1,A)') '*** This job stops after ', &
753  time_wallclock_limit * time_wallclock_safe, ' seconds.'
754  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Step interval for check : ', time_dtsec_wallclock_check, &
755  ' (step interval=', time_dstep_wallclock_check, ')'
756  endif
757 
758  return
integer, public time_nowstep
current step [number]
Definition: scale_time.F90:70
real(dp), public time_dtsec_ocean
time interval of ocean step [sec]
Definition: scale_time.F90:47
integer, public time_dstep_atmos_phy_mp
step interval of physics(microphysics)
Definition: scale_time.F90:54
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public calendar_unit2sec(second, value, unit)
Convert several units to second.
subroutine, public prc_mpistop
Abort MPI.
real(dp) function, public calendar_combine_daysec(absday, abssec)
Combine day and second.
real(dp), public time_dtsec_atmos_phy_rd
time interval of physics(radiation ) [sec]
Definition: scale_time.F90:42
integer, public time_nstep_atmos_dyn
small step of dynamics
Definition: scale_time.F90:39
real(dp), public time_nowms
subsecond part of current time [millisec]
Definition: scale_time.F90:66
real(dp), public time_nowsec
subday part of current time [sec]
Definition: scale_time.F90:68
subroutine, public calendar_date2char(chardate, ymdhms, subsec)
Convert from gregorian date to absolute day/second.
integer, public time_nstep
total steps [number]
Definition: scale_time.F90:71
real(dp), public time_dtsec_atmos_phy_mp
time interval of physics(microphysics) [sec]
Definition: scale_time.F90:41
real(dp), public time_dtsec_atmos_phy_ch
time interval of physics(chemistry ) [sec]
Definition: scale_time.F90:45
module ATMOSPHERIC Variables
integer, public time_dstep_atmos_dyn
step interval of dynamics
Definition: scale_time.F90:52
real(dp), public time_nowdaysec
second of current time [sec]
Definition: scale_time.F90:69
real(dp), public time_startdaysec
second of start time [sec]
Definition: scale_time.F90:74
integer, public time_nowday
absolute day of current time [day]
Definition: scale_time.F90:67
integer, public time_dstep_urban
step interval of urban step
Definition: scale_time.F90:62
real(dp), public time_dtsec_land
time interval of land step [sec]
Definition: scale_time.F90:48
integer, public time_dstep_wallclock_check
step interval of wallclock terminator
Definition: scale_time.F90:63
subroutine, public calendar_adjust_daysec(absday, abssec)
Adjust day and second.
integer, public time_dstep_ocean
step interval of ocean step
Definition: scale_time.F90:60
real(dp), public time_dtsec
time interval of model [sec]
Definition: scale_time.F90:36
real(dp) function, public calendar_cfunits2sec(cftime, cfunits, offset_year, startdaysec)
Convert time in units of the CF convention to second.
integer, public time_dstep_atmos_phy_sf
step interval of physics(surface flux)
Definition: scale_time.F90:56
integer, public time_offset_year
time offset [year]
Definition: scale_time.F90:73
integer, public time_dstep_atmos_phy_tb
step interval of physics(turbulence )
Definition: scale_time.F90:57
real(dp), public time_dtsec_wallclock_check
time interval of wallclock terminator [sec]
Definition: scale_time.F90:50
integer, public time_dstep_atmos_phy_rd
step interval of physics(radiation )
Definition: scale_time.F90:55
module TIME
Definition: scale_time.F90:15
real(dp), public time_dtsec_atmos_dyn
time interval of dynamics [sec]
Definition: scale_time.F90:38
real(dp) function, public prc_mpitime()
Get MPI time.
module PROCESS
real(dp), public time_dtsec_atmos_phy_sf
time interval of physics(surface flux) [sec]
Definition: scale_time.F90:43
real(dp), public time_dtsec_atmos_phy_tb
time interval of physics(turbulence ) [sec]
Definition: scale_time.F90:44
integer, public time_dstep_land
step interval of land step
Definition: scale_time.F90:61
real(dp), public time_dtsec_atmos_phy_ae
time interval of physics(aerosol ) [sec]
Definition: scale_time.F90:46
subroutine, public filegetdatainfo(basename, varname, myrank, istep, single, description, units, datatype, dim_rank, dim_name, dim_size, time_start, time_end, time_units)
module CONSTANT
Definition: scale_const.F90:14
integer, public time_dstep_atmos_phy_cp
step interval of physics(cumulus )
Definition: scale_time.F90:53
real(dp), public time_dtsec_atmos_phy_cp
time interval of physics(cumulus ) [sec]
Definition: scale_time.F90:40
integer, public prc_myrank
process num in local communicator
real(dp), parameter, public const_undef8
undefined value (REAL8)
Definition: scale_const.F90:42
integer, public time_dstep_atmos_phy_ae
step interval of physics(aerosol )
Definition: scale_time.F90:59
subroutine, public calendar_daysec2date(ymdhms, subsec, absday, abssec, offset_year)
Convert from gregorian date to absolute day/second.
integer, public time_dstep_atmos_phy_ch
step interval of physics(chemistry )
Definition: scale_time.F90:58
character(len=h_long), public atmos_restart_in_basename
basename of the restart file
module CALENDAR
integer, dimension(6), public time_nowdate
current time [YYYY MM DD HH MM SS]
Definition: scale_time.F90:65
real(dp), public time_dtsec_urban
time interval of urban step [sec]
Definition: scale_time.F90:49
subroutine, public calendar_date2daysec(absday, abssec, ymdhms, subsec, offset_year)
Convert from gregorian date to absolute day/second.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ admin_time_checkstate()

subroutine, public mod_admin_time::admin_time_checkstate ( )

Evaluate component execution.

Definition at line 764 of file mod_admin_time.f90.

References scale_calendar::calendar_date2char(), scale_stdio::io_fid_log, scale_stdio::io_l, scale_process::prc_mpitime(), time_doatmos_dyn, time_doatmos_phy_ae, time_doatmos_phy_ch, time_doatmos_phy_cp, time_doatmos_phy_mp, time_doatmos_phy_rd, time_doatmos_phy_sf, time_doatmos_phy_tb, time_doatmos_step, time_doland_step, time_doocean_step, time_doresume, time_dourban_step, scale_time::time_dstep_atmos_dyn, scale_time::time_dstep_atmos_phy_ae, scale_time::time_dstep_atmos_phy_ch, scale_time::time_dstep_atmos_phy_cp, scale_time::time_dstep_atmos_phy_mp, scale_time::time_dstep_atmos_phy_rd, scale_time::time_dstep_atmos_phy_sf, scale_time::time_dstep_atmos_phy_tb, scale_time::time_dstep_land, scale_time::time_dstep_ocean, time_dstep_resume, scale_time::time_dstep_urban, scale_time::time_nowdate, scale_time::time_nowms, scale_time::time_nowstep, and scale_time::time_nstep.

Referenced by mod_rm_driver::scalerm().

764  use scale_process, only: &
766  use scale_calendar, only: &
768  use scale_time, only: &
769  time_nowdate, &
770  time_nowms, &
771  time_nowstep, &
772  time_nstep, &
782  time_dstep_land, &
784  implicit none
785 
786  real(DP) :: wallclock_elapse
787  character(len=27) :: nowchardate
788  !---------------------------------------------------------------------------
789 
790  time_doatmos_step = .false.
791  time_doatmos_dyn = .false.
792  time_doatmos_phy_cp = .false.
793  time_doatmos_phy_mp = .false.
794  time_doatmos_phy_rd = .false.
795  time_doatmos_phy_sf = .false.
796  time_doatmos_phy_tb = .false.
797  time_doatmos_phy_ch = .false.
798  time_doatmos_phy_ae = .false.
799  time_doocean_step = .false.
800  time_doland_step = .false.
801  time_dourban_step = .false.
802  time_doresume = .false.
803 
804  time_res_atmos_dyn = time_res_atmos_dyn + 1
805  time_res_atmos_phy_cp = time_res_atmos_phy_cp + 1
806  time_res_atmos_phy_mp = time_res_atmos_phy_mp + 1
807  time_res_atmos_phy_rd = time_res_atmos_phy_rd + 1
808  time_res_atmos_phy_sf = time_res_atmos_phy_sf + 1
809  time_res_atmos_phy_tb = time_res_atmos_phy_tb + 1
810  time_res_atmos_phy_ch = time_res_atmos_phy_ch + 1
811  time_res_atmos_phy_ae = time_res_atmos_phy_ae + 1
812  time_res_ocean = time_res_ocean + 1
813  time_res_land = time_res_land + 1
814  time_res_urban = time_res_urban + 1
815  time_res_resume = time_res_resume + 1
816 
817  if ( time_res_atmos_dyn == time_dstep_atmos_dyn ) then
818  time_doatmos_step = .true.
819  time_doatmos_dyn = .true.
820  time_res_atmos_dyn = 0
821  endif
822  if ( time_res_atmos_phy_cp == time_dstep_atmos_phy_cp ) then
823  time_doatmos_step = .true.
824  time_doatmos_phy_cp = .true.
825  time_res_atmos_phy_cp = 0
826  endif
827  if ( time_res_atmos_phy_mp == time_dstep_atmos_phy_mp ) then
828  time_doatmos_step = .true.
829  time_doatmos_phy_mp = .true.
830  time_res_atmos_phy_mp = 0
831  endif
832  if ( time_res_atmos_phy_rd == time_dstep_atmos_phy_rd ) then
833  time_doatmos_step = .true.
834  time_doatmos_phy_rd = .true.
835  time_res_atmos_phy_rd = 0
836  endif
837  if ( time_res_atmos_phy_sf == time_dstep_atmos_phy_sf ) then
838  time_doatmos_step = .true.
839  time_doatmos_phy_sf = .true.
840  time_res_atmos_phy_sf = 0
841  endif
842  if ( time_res_atmos_phy_tb == time_dstep_atmos_phy_tb ) then
843  time_doatmos_step = .true.
844  time_doatmos_phy_tb = .true.
845  time_res_atmos_phy_tb = 0
846  endif
847  if ( time_res_atmos_phy_ch == time_dstep_atmos_phy_ch ) then
848  time_doatmos_step = .true.
849  time_doatmos_phy_ch = .true.
850  time_res_atmos_phy_ch = 0
851  endif
852  if ( time_res_atmos_phy_ae == time_dstep_atmos_phy_ae ) then
853  time_doatmos_step = .true.
854  time_doatmos_phy_ae = .true.
855  time_res_atmos_phy_ae = 0
856  endif
857 
858  if ( time_res_ocean == time_dstep_ocean ) then
859  time_doocean_step = .true.
860  time_res_ocean = 0
861  endif
862  if ( time_res_land == time_dstep_land ) then
863  time_doland_step = .true.
864  time_res_land = 0
865  endif
866  if ( time_res_urban == time_dstep_urban ) then
867  time_dourban_step = .true.
868  time_res_urban = 0
869  endif
870  if ( time_res_resume == time_dstep_resume ) then
871  time_doresume = .true.
872  time_res_resume = 0
873  endif
874 
875  call calendar_date2char( nowchardate, & ! [OUT]
876  time_nowdate(:), & ! [IN]
877  time_nowms ) ! [IN]
878 
879  if ( time_wallclock_limit > 0.0_dp ) then
880  wallclock_elapse = prc_mpitime() - time_wallclock_start
881  if( io_l ) write(io_fid_log,'(1x,3A,I7,A,I7,A,F10.1)') '*** TIME: ', nowchardate,' STEP:',time_nowstep, '/', time_nstep, &
882  ' WCLOCK:', wallclock_elapse
883  else
884  if( io_l ) write(io_fid_log,'(1x,3A,I7,A,I7)') '*** TIME: ', nowchardate,' STEP:',time_nowstep, '/', time_nstep
885  endif
886 
887  return
integer, public time_nowstep
current step [number]
Definition: scale_time.F90:70
integer, public time_dstep_atmos_phy_mp
step interval of physics(microphysics)
Definition: scale_time.F90:54
real(dp), public time_nowms
subsecond part of current time [millisec]
Definition: scale_time.F90:66
subroutine, public calendar_date2char(chardate, ymdhms, subsec)
Convert from gregorian date to absolute day/second.
integer, public time_nstep
total steps [number]
Definition: scale_time.F90:71
integer, public time_dstep_atmos_dyn
step interval of dynamics
Definition: scale_time.F90:52
integer, public time_dstep_urban
step interval of urban step
Definition: scale_time.F90:62
integer, public time_dstep_ocean
step interval of ocean step
Definition: scale_time.F90:60
integer, public time_dstep_atmos_phy_sf
step interval of physics(surface flux)
Definition: scale_time.F90:56
integer, public time_dstep_atmos_phy_tb
step interval of physics(turbulence )
Definition: scale_time.F90:57
integer, public time_dstep_atmos_phy_rd
step interval of physics(radiation )
Definition: scale_time.F90:55
module TIME
Definition: scale_time.F90:15
real(dp) function, public prc_mpitime()
Get MPI time.
module PROCESS
integer, public time_dstep_land
step interval of land step
Definition: scale_time.F90:61
integer, public time_dstep_atmos_phy_cp
step interval of physics(cumulus )
Definition: scale_time.F90:53
integer, public time_dstep_atmos_phy_ae
step interval of physics(aerosol )
Definition: scale_time.F90:59
integer, public time_dstep_atmos_phy_ch
step interval of physics(chemistry )
Definition: scale_time.F90:58
module CALENDAR
integer, dimension(6), public time_nowdate
current time [YYYY MM DD HH MM SS]
Definition: scale_time.F90:65
Here is the call graph for this function:
Here is the caller graph for this function:

◆ admin_time_advance()

subroutine, public mod_admin_time::admin_time_advance ( )

Advance the time & evaluate restart & stop.

Definition at line 893 of file mod_admin_time.f90.

References scale_calendar::calendar_adjust_daysec(), scale_calendar::calendar_combine_daysec(), scale_calendar::calendar_daysec2date(), scale_stdio::io_fid_log, scale_stdio::io_l, scale_process::prc_ismaster, scale_process::prc_mpitime(), time_doatmos_restart, time_doend, time_doland_restart, time_doocean_restart, time_dourban_restart, time_dstep_atmos_restart, time_dstep_land_restart, time_dstep_ocean_restart, time_dstep_urban_restart, scale_time::time_dstep_wallclock_check, scale_time::time_dtsec, scale_time::time_nowdate, scale_time::time_nowday, scale_time::time_nowdaysec, scale_time::time_nowms, scale_time::time_nowsec, scale_time::time_nowstep, scale_time::time_nstep, and scale_time::time_offset_year.

Referenced by mod_rm_driver::scalerm().

893  use scale_process, only: &
894  prc_ismaster, &
896  use scale_calendar, only: &
900  use scale_time, only: &
901  time_dtsec, &
902  time_nowdate, &
903  time_nowms, &
904  time_nowday, &
905  time_nowsec, &
906  time_nowdaysec, &
907  time_nowstep, &
908  time_nstep, &
911  implicit none
912 
913  real(DP) :: wallclock_elapse
914  logical :: exists
915  !---------------------------------------------------------------------------
916 
917  time_doend = .false.
918 
920  time_nowday = time_startday
921  time_nowsec = time_startsec + real(TIME_NOWSTEP-1,kind=DP) * time_dtsec
922 
923  ! reallocate day & sub-day
925 
926  call calendar_daysec2date( time_nowdate(:), & ! [OUT]
927  time_nowms, & ! [OUT]
928  time_nowday, & ! [IN]
929  time_nowsec, & ! [IN]
930  time_offset_year ) ! [IN]
931 
933 
934  if ( time_nowstep > time_nstep ) then
935  time_doend = .true.
936  endif
937 
938  if ( time_wallclock_limit > 0.0_dp &
939  .AND. mod(time_nowstep-1,time_dstep_wallclock_check) == 0 ) then
940  wallclock_elapse = prc_mpitime() - time_wallclock_start
941 
942  if ( wallclock_elapse > time_wallclock_limit * time_wallclock_safe ) then
943  if( io_l ) write(io_fid_log,*) '*** Elapse time limit is detected. Termination operation starts.'
944  time_doend = .true.
945  endif
946  endif
947 
948  ! QUIT file control
949  if ( prc_ismaster ) then ! master node
950  inquire(file='QUIT', exist=exists)
951 
952  if( exists ) then
953  if( io_l ) write(io_fid_log,*) '*** QUIT file is found. Termination operation starts.'
954  time_doend = .true.
955  endif
956  endif
957 
958  time_doatmos_restart = .false.
959  time_doocean_restart = .false.
960  time_doland_restart = .false.
961  time_dourban_restart = .false.
962 
963  time_res_atmos_restart = time_res_atmos_restart + 1
964  time_res_ocean_restart = time_res_ocean_restart + 1
965  time_res_land_restart = time_res_land_restart + 1
966  time_res_urban_restart = time_res_urban_restart + 1
967 
968  if ( time_res_atmos_restart == time_dstep_atmos_restart ) then
969  time_doatmos_restart = .true.
970  time_res_atmos_restart = 0
971  elseif( time_doend ) then
972  time_doatmos_restart = .true.
973  endif
974 
975  if ( time_res_ocean_restart == time_dstep_ocean_restart ) then
976  time_doocean_restart = .true.
977  time_res_ocean_restart = 0
978  elseif( time_doend ) then
979  time_doocean_restart = .true.
980  endif
981 
982  if ( time_res_land_restart == time_dstep_land_restart ) then
983  time_doland_restart = .true.
984  time_res_land_restart = 0
985  elseif( time_doend ) then
986  time_doland_restart = .true.
987  endif
988 
989  if ( time_res_urban_restart == time_dstep_urban_restart ) then
990  time_dourban_restart = .true.
991  time_res_urban_restart = 0
992  elseif( time_doend ) then
993  time_dourban_restart = .true.
994  endif
995 
996  return
integer, public time_nowstep
current step [number]
Definition: scale_time.F90:70
logical, public prc_ismaster
master process in local communicator?
real(dp) function, public calendar_combine_daysec(absday, abssec)
Combine day and second.
real(dp), public time_nowms
subsecond part of current time [millisec]
Definition: scale_time.F90:66
real(dp), public time_nowsec
subday part of current time [sec]
Definition: scale_time.F90:68
integer, public time_nstep
total steps [number]
Definition: scale_time.F90:71
real(dp), public time_nowdaysec
second of current time [sec]
Definition: scale_time.F90:69
logical, public time_doatmos_restart
execute atmosphere restart output in this step?
logical, public time_doocean_restart
execute ocean restart output in this step?
integer, public time_nowday
absolute day of current time [day]
Definition: scale_time.F90:67
integer, public time_dstep_wallclock_check
step interval of wallclock terminator
Definition: scale_time.F90:63
subroutine, public calendar_adjust_daysec(absday, abssec)
Adjust day and second.
real(dp), public time_dtsec
time interval of model [sec]
Definition: scale_time.F90:36
integer, public time_offset_year
time offset [year]
Definition: scale_time.F90:73
module TIME
Definition: scale_time.F90:15
real(dp) function, public prc_mpitime()
Get MPI time.
module PROCESS
logical, public time_dourban_restart
execute urban restart output in this step?
subroutine, public calendar_daysec2date(ymdhms, subsec, absday, abssec, offset_year)
Convert from gregorian date to absolute day/second.
logical, public time_doland_restart
execute land restart output in this step?
module CALENDAR
integer, dimension(6), public time_nowdate
current time [YYYY MM DD HH MM SS]
Definition: scale_time.F90:65
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ time_dtsec_atmos_restart

real(dp), public mod_admin_time::time_dtsec_atmos_restart

time interval of atmosphere restart [sec]

Definition at line 33 of file mod_admin_time.f90.

Referenced by admin_time_setup().

33  real(DP), public :: time_dtsec_atmos_restart

◆ time_dtsec_ocean_restart

real(dp), public mod_admin_time::time_dtsec_ocean_restart

time interval of ocean restart [sec]

Definition at line 34 of file mod_admin_time.f90.

Referenced by admin_time_setup().

34  real(DP), public :: time_dtsec_ocean_restart

◆ time_dtsec_land_restart

real(dp), public mod_admin_time::time_dtsec_land_restart

time interval of land restart [sec]

Definition at line 35 of file mod_admin_time.f90.

Referenced by admin_time_setup().

35  real(DP), public :: time_dtsec_land_restart

◆ time_dtsec_urban_restart

real(dp), public mod_admin_time::time_dtsec_urban_restart

time interval of urban restart [sec]

Definition at line 36 of file mod_admin_time.f90.

Referenced by admin_time_setup().

36  real(DP), public :: time_dtsec_urban_restart

◆ time_dtsec_resume

real(dp), public mod_admin_time::time_dtsec_resume

time interval for resume [sec]

Definition at line 37 of file mod_admin_time.f90.

Referenced by admin_time_setup().

37  real(DP), public :: time_dtsec_resume

◆ time_dstep_atmos_restart

integer, public mod_admin_time::time_dstep_atmos_restart

interval of atmosphere restart [step]

Definition at line 39 of file mod_admin_time.f90.

Referenced by admin_time_advance(), and admin_time_setup().

39  integer, public :: time_dstep_atmos_restart

◆ time_dstep_ocean_restart

integer, public mod_admin_time::time_dstep_ocean_restart

interval of ocean restart [step]

Definition at line 40 of file mod_admin_time.f90.

Referenced by admin_time_advance(), and admin_time_setup().

40  integer, public :: time_dstep_ocean_restart

◆ time_dstep_land_restart

integer, public mod_admin_time::time_dstep_land_restart

interval of land restart [step]

Definition at line 41 of file mod_admin_time.f90.

Referenced by admin_time_advance(), and admin_time_setup().

41  integer, public :: time_dstep_land_restart

◆ time_dstep_urban_restart

integer, public mod_admin_time::time_dstep_urban_restart

interval of urban restart [step]

Definition at line 42 of file mod_admin_time.f90.

Referenced by admin_time_advance(), and admin_time_setup().

42  integer, public :: time_dstep_urban_restart

◆ time_dstep_resume

integer, public mod_admin_time::time_dstep_resume

interval for resume [step]

Definition at line 43 of file mod_admin_time.f90.

Referenced by admin_time_checkstate(), and admin_time_setup().

43  integer, public :: time_dstep_resume

◆ time_doatmos_step

logical, public mod_admin_time::time_doatmos_step

execute atmosphere component in this step?

Definition at line 45 of file mod_admin_time.f90.

Referenced by admin_time_checkstate(), and mod_rm_driver::scalerm().

45  logical, public :: time_doatmos_step

◆ time_doatmos_dyn

logical, public mod_admin_time::time_doatmos_dyn

execute dynamics in this step?

Definition at line 46 of file mod_admin_time.f90.

Referenced by admin_time_checkstate(), and mod_atmos_driver::atmos_driver().

46  logical, public :: time_doatmos_dyn

◆ time_doatmos_phy_cp

logical, public mod_admin_time::time_doatmos_phy_cp

execute physics in this step? (cumulus )

Definition at line 47 of file mod_admin_time.f90.

Referenced by admin_time_checkstate(), and mod_atmos_driver::atmos_driver().

47  logical, public :: time_doatmos_phy_cp

◆ time_doatmos_phy_mp

logical, public mod_admin_time::time_doatmos_phy_mp

execute physics in this step? (microphysics)

Definition at line 48 of file mod_admin_time.f90.

Referenced by admin_time_checkstate(), and mod_atmos_driver::atmos_driver().

48  logical, public :: time_doatmos_phy_mp

◆ time_doatmos_phy_rd

logical, public mod_admin_time::time_doatmos_phy_rd

execute physics in this step? (radiation )

Definition at line 49 of file mod_admin_time.f90.

Referenced by admin_time_checkstate(), and mod_atmos_driver::atmos_driver().

49  logical, public :: time_doatmos_phy_rd

◆ time_doatmos_phy_sf

logical, public mod_admin_time::time_doatmos_phy_sf

execute physics in this step? (surface flux)

Definition at line 50 of file mod_admin_time.f90.

Referenced by admin_time_checkstate(), and mod_atmos_driver::atmos_driver().

50  logical, public :: time_doatmos_phy_sf

◆ time_doatmos_phy_tb

logical, public mod_admin_time::time_doatmos_phy_tb

execute physics in this step? (turbulence )

Definition at line 51 of file mod_admin_time.f90.

Referenced by admin_time_checkstate(), and mod_atmos_driver::atmos_driver().

51  logical, public :: time_doatmos_phy_tb

◆ time_doatmos_phy_ch

logical, public mod_admin_time::time_doatmos_phy_ch

execute physics in this step? (chemistry )

Definition at line 52 of file mod_admin_time.f90.

Referenced by admin_time_checkstate(), and mod_atmos_driver::atmos_driver().

52  logical, public :: time_doatmos_phy_ch

◆ time_doatmos_phy_ae

logical, public mod_admin_time::time_doatmos_phy_ae

execute physics in this step? (aerosol )

Definition at line 53 of file mod_admin_time.f90.

Referenced by admin_time_checkstate(), and mod_atmos_driver::atmos_driver().

53  logical, public :: time_doatmos_phy_ae

◆ time_doatmos_restart

logical, public mod_admin_time::time_doatmos_restart

execute atmosphere restart output in this step?

Definition at line 54 of file mod_admin_time.f90.

Referenced by mod_admin_restart::admin_restart(), admin_time_advance(), mod_mkinit::mkinit(), and mod_rm_driver::scalerm().

54  logical, public :: time_doatmos_restart
logical, public time_doatmos_restart
execute atmosphere restart output in this step?

◆ time_doocean_step

logical, public mod_admin_time::time_doocean_step

execute ocean component in this step?

Definition at line 55 of file mod_admin_time.f90.

Referenced by admin_time_checkstate(), and mod_rm_driver::scalerm().

55  logical, public :: time_doocean_step

◆ time_doocean_restart

logical, public mod_admin_time::time_doocean_restart

execute ocean restart output in this step?

Definition at line 56 of file mod_admin_time.f90.

Referenced by mod_admin_restart::admin_restart(), admin_time_advance(), mod_mkinit::mkinit(), and mod_rm_driver::scalerm().

56  logical, public :: time_doocean_restart
logical, public time_doocean_restart
execute ocean restart output in this step?

◆ time_doland_step

logical, public mod_admin_time::time_doland_step

execute land component in this step?

Definition at line 57 of file mod_admin_time.f90.

Referenced by admin_time_checkstate(), and mod_rm_driver::scalerm().

57  logical, public :: time_doland_step

◆ time_doland_restart

logical, public mod_admin_time::time_doland_restart

execute land restart output in this step?

Definition at line 58 of file mod_admin_time.f90.

Referenced by mod_admin_restart::admin_restart(), admin_time_advance(), mod_mkinit::mkinit(), and mod_rm_driver::scalerm().

58  logical, public :: time_doland_restart
logical, public time_doland_restart
execute land restart output in this step?

◆ time_dourban_step

logical, public mod_admin_time::time_dourban_step

execute urban component in this step?

Definition at line 59 of file mod_admin_time.f90.

Referenced by admin_time_checkstate(), and mod_rm_driver::scalerm().

59  logical, public :: time_dourban_step

◆ time_dourban_restart

logical, public mod_admin_time::time_dourban_restart

execute urban restart output in this step?

Definition at line 60 of file mod_admin_time.f90.

Referenced by mod_admin_restart::admin_restart(), admin_time_advance(), mod_mkinit::mkinit(), and mod_rm_driver::scalerm().

60  logical, public :: time_dourban_restart
logical, public time_dourban_restart
execute urban restart output in this step?

◆ time_doresume

logical, public mod_admin_time::time_doresume

resume in this step?

Definition at line 61 of file mod_admin_time.f90.

Referenced by admin_time_checkstate(), and mod_rm_driver::scalerm().

61  logical, public :: time_doresume

◆ time_doend

logical, public mod_admin_time::time_doend

finish program in this step?

Definition at line 62 of file mod_admin_time.f90.

Referenced by admin_time_advance(), and mod_rm_driver::scalerm().

62  logical, public :: time_doend