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, dimension(6) (/ -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 112 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, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_fid_nml, scale_stdio::io_l, scale_stdio::io_nml, 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().

112  use gtool_file, only: &
113  filegetdatainfo
114  use scale_process, only: &
115  prc_myrank, &
116  prc_mpistop, &
118  use scale_const, only: &
119  undef8 => const_undef8
120  use scale_calendar, only: &
128  use scale_time, only: &
129  time_dtsec, &
130  time_nowdate, &
131  time_nowms, &
132  time_nowday, &
133  time_nowsec, &
134  time_nowdaysec, &
135  time_nowstep, &
136  time_nstep, &
138  nstep_dyn => time_nstep_atmos_dyn, &
147  time_dtsec_land, &
159  time_dstep_land, &
164  use mod_atmos_vars, only: &
165  restart_in_basename => atmos_restart_in_basename
166  implicit none
167 
168  logical, intent(in) :: setup_TimeIntegration
169 
170  real(DP) :: TIME_DURATION = undef8
171  character(len=H_SHORT) :: TIME_DURATION_UNIT = "SEC"
172  real(DP) :: TIME_DT = undef8
173  character(len=H_SHORT) :: TIME_DT_UNIT = "SEC"
174 
175  real(DP) :: TIME_DT_ATMOS_DYN = undef8
176  character(len=H_SHORT) :: TIME_DT_ATMOS_DYN_UNIT = "SEC"
177  integer :: TIME_NSTEP_ATMOS_DYN = -1
178  real(DP) :: TIME_DT_ATMOS_PHY_CP = undef8
179  character(len=H_SHORT) :: TIME_DT_ATMOS_PHY_CP_UNIT = ""
180  real(DP) :: TIME_DT_ATMOS_PHY_MP = undef8
181  character(len=H_SHORT) :: TIME_DT_ATMOS_PHY_MP_UNIT = ""
182  real(DP) :: TIME_DT_ATMOS_PHY_RD = undef8
183  character(len=H_SHORT) :: TIME_DT_ATMOS_PHY_RD_UNIT = ""
184  real(DP) :: TIME_DT_ATMOS_PHY_SF = undef8
185  character(len=H_SHORT) :: TIME_DT_ATMOS_PHY_SF_UNIT = ""
186  real(DP) :: TIME_DT_ATMOS_PHY_TB = undef8
187  character(len=H_SHORT) :: TIME_DT_ATMOS_PHY_TB_UNIT = ""
188  real(DP) :: TIME_DT_ATMOS_PHY_CH = undef8
189  character(len=H_SHORT) :: TIME_DT_ATMOS_PHY_CH_UNIT = ""
190  real(DP) :: TIME_DT_ATMOS_PHY_AE = undef8
191  character(len=H_SHORT) :: TIME_DT_ATMOS_PHY_AE_UNIT = ""
192  real(DP) :: TIME_DT_ATMOS_RESTART = undef8
193  character(len=H_SHORT) :: TIME_DT_ATMOS_RESTART_UNIT = ""
194 
195  real(DP) :: TIME_DT_OCEAN = undef8
196  character(len=H_SHORT) :: TIME_DT_OCEAN_UNIT = ""
197  real(DP) :: TIME_DT_OCEAN_RESTART = undef8
198  character(len=H_SHORT) :: TIME_DT_OCEAN_RESTART_UNIT = ""
199 
200  real(DP) :: TIME_DT_LAND = undef8
201  character(len=H_SHORT) :: TIME_DT_LAND_UNIT = ""
202  real(DP) :: TIME_DT_LAND_RESTART = undef8
203  character(len=H_SHORT) :: TIME_DT_LAND_RESTART_UNIT = ""
204 
205  real(DP) :: TIME_DT_URBAN = undef8
206  character(len=H_SHORT) :: TIME_DT_URBAN_UNIT = ""
207  real(DP) :: TIME_DT_URBAN_RESTART = undef8
208  character(len=H_SHORT) :: TIME_DT_URBAN_RESTART_UNIT = ""
209 
210  real(DP) :: TIME_DT_RESUME = undef8
211  character(len=H_SHORT) :: TIME_DT_RESUME_UNIT = ""
212 
213  real(DP) :: TIME_DT_WALLCLOCK_CHECK = undef8
214  character(len=H_SHORT) :: TIME_DT_WALLCLOCK_CHECK_UNIT = ""
215 
216  namelist / param_time / &
217  time_startdate, &
218  time_startms, &
219  time_duration, &
220  time_duration_unit, &
221  time_dt, &
222  time_dt_unit, &
223  time_dt_atmos_dyn, &
224  time_dt_atmos_dyn_unit, &
226  time_dt_atmos_phy_cp, &
227  time_dt_atmos_phy_cp_unit, &
228  time_dt_atmos_phy_mp, &
229  time_dt_atmos_phy_mp_unit, &
230  time_dt_atmos_phy_rd, &
231  time_dt_atmos_phy_rd_unit, &
232  time_dt_atmos_phy_sf, &
233  time_dt_atmos_phy_sf_unit, &
234  time_dt_atmos_phy_tb, &
235  time_dt_atmos_phy_tb_unit, &
236  time_dt_atmos_phy_ch, &
237  time_dt_atmos_phy_ch_unit, &
238  time_dt_atmos_phy_ae, &
239  time_dt_atmos_phy_ae_unit, &
240  time_dt_atmos_restart, &
241  time_dt_atmos_restart_unit, &
242  time_dt_ocean, &
243  time_dt_ocean_unit, &
244  time_dt_ocean_restart, &
245  time_dt_ocean_restart_unit, &
246  time_dt_land, &
247  time_dt_land_unit, &
248  time_dt_land_restart, &
249  time_dt_land_restart_unit, &
250  time_dt_urban, &
251  time_dt_urban_unit, &
252  time_dt_urban_restart, &
253  time_dt_urban_restart_unit, &
254  time_dt_wallclock_check, &
255  time_dt_wallclock_check_unit, &
256  time_dt_resume, &
257  time_dt_resume_unit, &
258  time_wallclock_limit, &
259  time_wallclock_safe
260 
261  integer :: dateday
262  real(DP) :: datesec
263  real(DP) :: cftime
264  character(len=H_MID) :: cfunits
265 
266  real(DP) :: TIME_DURATIONSEC
267  character(len=27) :: startchardate
268  character(len=27) :: endchardate
269 
270  integer :: ierr
271  !---------------------------------------------------------------------------
272 
273  if( io_l ) write(io_fid_log,*)
274  if( io_l ) write(io_fid_log,*) '++++++ Module[TIME] / Categ[ADMIN] / Origin[SCALE-RM]'
275 
277 
278  !--- read namelist
279  rewind(io_fid_conf)
280  read(io_fid_conf,nml=param_time,iostat=ierr)
281  if( ierr < 0 ) then !--- missing
282  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
283  elseif( ierr > 0 ) then !--- fatal error
284  write(*,*) 'xxx Not appropriate names in namelist PARAM_TIME. Check!'
285  call prc_mpistop
286  endif
287  if( io_nml ) write(io_fid_nml,nml=param_time)
288 
289  ! check time setting
290  if ( setup_timeintegration ) then
291 
292  if( io_l ) write(io_fid_log,*)
293  if( io_l ) write(io_fid_log,*) '*** Check time interval and unit for each component ***'
294 
295  if ( time_dt == undef8 ) then
296  write(*,*) 'xxx Not found TIME_DT. STOP.'
297  call prc_mpistop
298  endif
299  if ( time_duration == undef8 ) then
300  write(*,*) 'xxx Not found TIME_DURATION. STOP.'
301  call prc_mpistop
302  endif
303 
304  ! DYN
305  if ( time_dt_atmos_dyn == undef8 ) then
306  if ( time_nstep_atmos_dyn < 0 ) then
307  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_ATMOS_DYN. ', &
308  'TIME_DT is used.'
309  time_dt_atmos_dyn = time_dt
310  endif
311  endif
312  if ( time_dt_atmos_dyn_unit == '' ) then
313  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_ATMOS_DYN_UNIT. ', &
314  'TIME_DT_UNIT is used.'
315  time_dt_atmos_dyn_unit = time_dt_unit
316  endif
317  ! PHY_CP
318  if ( time_dt_atmos_phy_cp == undef8 ) then
319  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_ATMOS_PHY_CP. ', &
320  'TIME_DT is used.'
321  time_dt_atmos_phy_cp = time_dt
322  endif
323  if ( time_dt_atmos_phy_cp_unit == '' ) then
324  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_ATMOS_PHY_CP_UNIT. ', &
325  'TIME_DT_UNIT is used.'
326  time_dt_atmos_phy_cp_unit = time_dt_unit
327  endif
328  ! PHY_MP
329  if ( time_dt_atmos_phy_mp == undef8 ) then
330  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_ATMOS_PHY_MP. ', &
331  'TIME_DT is used.'
332  time_dt_atmos_phy_mp = time_dt
333  endif
334  if ( time_dt_atmos_phy_mp_unit == '' ) then
335  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_ATMOS_PHY_MP_UNIT. ', &
336  'TIME_DT_UNIT is used.'
337  time_dt_atmos_phy_mp_unit = time_dt_unit
338  endif
339  ! PHY_RD
340  if ( time_dt_atmos_phy_rd == undef8 ) then
341  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_ATMOS_PHY_RD. ', &
342  'TIME_DT is used.'
343  time_dt_atmos_phy_rd = time_dt
344  endif
345  if ( time_dt_atmos_phy_rd_unit == '' ) then
346  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_ATMOS_PHY_RD_UNIT. ', &
347  'TIME_DT_UNIT is used.'
348  time_dt_atmos_phy_rd_unit = time_dt_unit
349  endif
350  ! PHY_SF
351  if ( time_dt_atmos_phy_sf == undef8 ) then
352  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_ATMOS_PHY_SF. ', &
353  'TIME_DT is used.'
354  time_dt_atmos_phy_sf = time_dt
355  endif
356  if ( time_dt_atmos_phy_sf_unit == '' ) then
357  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_ATMOS_PHY_SF_UNIT. ', &
358  'TIME_DT_UNIT is used.'
359  time_dt_atmos_phy_sf_unit = time_dt_unit
360  endif
361  ! PHY_TB
362  if ( time_dt_atmos_phy_tb == undef8 ) then
363  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_ATMOS_PHY_TB. ', &
364  'TIME_DT is used.'
365  time_dt_atmos_phy_tb = time_dt
366  endif
367  if ( time_dt_atmos_phy_tb_unit == '' ) then
368  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_ATMOS_PHY_TB_UNIT. ', &
369  'TIME_DT_UNIT is used.'
370  time_dt_atmos_phy_tb_unit = time_dt_unit
371  endif
372  ! PHY_CH
373  if ( time_dt_atmos_phy_ch == undef8 ) then
374  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_ATMOS_PHY_CH. ', &
375  'TIME_DT is used.'
376  time_dt_atmos_phy_ch = time_dt
377  endif
378  if ( time_dt_atmos_phy_ch_unit == '' ) then
379  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_ATMOS_PHY_CH_UNIT. ', &
380  'TIME_DT_UNIT is used.'
381  time_dt_atmos_phy_ch_unit = time_dt_unit
382  endif
383  ! PHY_AE
384  if ( time_dt_atmos_phy_ae == undef8 ) then
385  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_ATMOS_PHY_AE. ', &
386  'TIME_DT is used.'
387  time_dt_atmos_phy_ae = time_dt
388  endif
389  if ( time_dt_atmos_phy_ae_unit == '' ) then
390  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_ATMOS_PHY_AE_UNIT. ', &
391  'TIME_DT_UNIT is used.'
392  time_dt_atmos_phy_ae_unit = time_dt_unit
393  endif
394  ! ATMOS RESTART
395  if ( time_dt_atmos_restart == undef8 ) then
396  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_ATMOS_RESTART. ', &
397  'TIME_DURATION is used.'
398  time_dt_atmos_restart = time_duration
399  endif
400  if ( time_dt_atmos_restart_unit == '' ) then
401  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_ATMOS_RESTART_UNIT. ', &
402  'TIME_DURATION_UNIT is used.'
403  time_dt_atmos_restart_unit = time_duration_unit
404  endif
405  ! OCEAN
406  if ( time_dt_ocean == undef8 ) then
407  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_OCEAN. ', &
408  'TIME_DT is used.'
409  time_dt_ocean = time_dt
410  endif
411  if ( time_dt_ocean_unit == '' ) then
412  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_OCEAN_UNIT. ', &
413  'TIME_DT_UNIT is used.'
414  time_dt_ocean_unit = time_dt_unit
415  endif
416  ! OCEAN RESTART
417  if ( time_dt_ocean_restart == undef8 ) then
418  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_OCEAN_RESTART. ', &
419  'TIME_DURATION is used.'
420  time_dt_ocean_restart = time_duration
421  endif
422  if ( time_dt_ocean_restart_unit == '' ) then
423  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_OCEAN_RESTART_UNIT. ', &
424  'TIME_DURATION_UNIT is used.'
425  time_dt_ocean_restart_unit = time_duration_unit
426  endif
427  ! LAND
428  if ( time_dt_land == undef8 ) then
429  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_LAND. ', &
430  'TIME_DT is used.'
431  time_dt_land = time_dt
432  endif
433  if ( time_dt_land_unit == '' ) then
434  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_LAND_UNIT. ', &
435  'TIME_DT_UNIT is used.'
436  time_dt_land_unit = time_dt_unit
437  endif
438  ! LAND RESTART
439  if ( time_dt_land_restart == undef8 ) then
440  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_LAND_RESTART. ', &
441  'TIME_DURATION is used.'
442  time_dt_land_restart = time_duration
443  endif
444  if ( time_dt_land_restart_unit == '' ) then
445  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_LAND_RESTART_UNIT. ', &
446  'TIME_DURATION_UNIT is used.'
447  time_dt_land_restart_unit = time_duration_unit
448  endif
449  ! URBAN
450  if ( time_dt_urban == undef8 ) then
451  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_URBAN. ', &
452  'TIME_DT is used.'
453  time_dt_urban = time_dt
454  endif
455  if ( time_dt_urban_unit == '' ) then
456  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_URBAN_UNIT. ', &
457  'TIME_DT_UNIT is used.'
458  time_dt_urban_unit = time_dt_unit
459  endif
460  ! URBAN RESTART
461  if ( time_dt_urban_restart == undef8 ) then
462  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_URBAN_RESTART. ', &
463  'TIME_DURATION is used.'
464  time_dt_urban_restart = time_duration
465  endif
466  if ( time_dt_urban_restart_unit == '' ) then
467  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_URBAN_RESTART_UNIT. ', &
468  'TIME_DURATION_UNIT is used.'
469  time_dt_urban_restart_unit = time_duration_unit
470  endif
471  ! Resume
472  if ( time_dt_resume == undef8 ) then
473  time_dt_resume = time_duration
474  endif
475  if ( time_dt_resume_unit == '' ) then
476  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_RESUME_UNIT. ', &
477  'TIME_DURATION_UNIT is used.'
478  time_dt_resume_unit = time_duration_unit
479  endif
480  endif
481 
482  !--- calculate time
483  if ( time_startdate(1) == -999 ) then
484  if ( restart_in_basename /= '' ) then ! read start time from the restart data
485  call filegetdatainfo( restart_in_basename, & ! [IN]
486  'DENS', & ! [IN]
487  prc_myrank, & ! [IN]
488  0, & ! [IN] step
489  time_start = cftime, & ! [OUT]
490  time_units = cfunits ) ! [OUT]
491 
492  dateday = 0
493  datesec = calendar_cfunits2sec( cftime, cfunits, 0 )
494 
495  call calendar_adjust_daysec( dateday, datesec )
496 
497  call calendar_daysec2date( time_startdate, & ! [OUT]
498  time_startms, & ! [OUT]
499  dateday, & ! [IN]
500  datesec, & ! [IN]
501  0 ) ! [IN]
502  else
503  time_startdate = (/ 0, 1, 1, 0, 0, 0 /)
504  time_startms = 0.0_dp
505  endif
506  else
507  time_startms = time_startms * 1.e-3_dp
508  endif
509 
510  time_offset_year = time_startdate(1)
511 
512  call calendar_date2daysec( time_startday, & ! [OUT]
513  time_startsec, & ! [OUT]
514  time_startdate(:), & ! [IN]
515  time_startms, & ! [IN]
516  time_offset_year ) ! [IN]
517 
518  call calendar_date2char( startchardate, & ! [OUT]
519  time_startdate(:), & ! [IN]
520  time_startms ) ! [IN]
521 
522  time_startdaysec = calendar_combine_daysec( time_startday, time_startsec )
523 
524  time_nowdate(:) = time_startdate(:)
525  time_nowms = time_startms
526  time_nowday = time_startday
527  time_nowsec = time_startsec
529 
530  time_endday = time_startday
531 
532  if ( setup_timeintegration ) then
533  call calendar_unit2sec( time_durationsec, time_duration, time_duration_unit )
534  time_endsec = time_startsec + time_durationsec
535  else
536  time_endsec = time_startsec
537  endif
538 
539  call calendar_adjust_daysec( time_endday, time_endsec ) ! [INOUT]
540 
541  call calendar_daysec2date( time_enddate(:), & ! [OUT]
542  time_endms, & ! [OUT]
543  time_endday, & ! [IN]
544  time_endsec, & ! [IN]
545  time_offset_year ) ! [IN]
546 
547  call calendar_date2char( endchardate, & ! [OUT]
548  time_enddate(:), & ! [IN]
549  time_endms ) ! [IN]
550 
551  if( io_l ) write(io_fid_log,*)
552  if( io_l ) write(io_fid_log,*) '*** Global date / time setting ***'
553  if( io_l ) write(io_fid_log,'(1x,A,A)') '*** START Date : ', startchardate
554  if( io_l ) write(io_fid_log,'(1x,A,A)') '*** END Date : ', endchardate
555 
556  if ( setup_timeintegration ) then
557 
558  call calendar_unit2sec( time_dtsec, time_dt, time_dt_unit )
559 
560  time_nstep = int( time_durationsec / time_dtsec )
561  time_nowstep = 1
562 
563  if( io_l ) write(io_fid_log,'(1x,A,F10.3)') '*** delta t (sec.) :', time_dtsec
564  if( io_l ) write(io_fid_log,'(1x,A,I10)' ) '*** No. of steps :', time_nstep
565 
566  !--- calculate intervals for atmosphere
567  if ( time_dt_atmos_dyn /= undef8 ) then
568  call calendar_unit2sec( time_dtsec_atmos_dyn, time_dt_atmos_dyn, time_dt_atmos_dyn_unit )
569  else
571  endif
572  call calendar_unit2sec( time_dtsec_atmos_phy_cp, time_dt_atmos_phy_cp, time_dt_atmos_phy_cp_unit )
573  call calendar_unit2sec( time_dtsec_atmos_phy_mp, time_dt_atmos_phy_mp, time_dt_atmos_phy_mp_unit )
574  call calendar_unit2sec( time_dtsec_atmos_phy_rd, time_dt_atmos_phy_rd, time_dt_atmos_phy_rd_unit )
575  call calendar_unit2sec( time_dtsec_atmos_phy_sf, time_dt_atmos_phy_sf, time_dt_atmos_phy_sf_unit )
576  call calendar_unit2sec( time_dtsec_atmos_phy_tb, time_dt_atmos_phy_tb, time_dt_atmos_phy_tb_unit )
577  call calendar_unit2sec( time_dtsec_atmos_phy_ch, time_dt_atmos_phy_ch, time_dt_atmos_phy_ch_unit )
578  call calendar_unit2sec( time_dtsec_atmos_phy_ae, time_dt_atmos_phy_ae, time_dt_atmos_phy_ae_unit )
579  call calendar_unit2sec( time_dtsec_atmos_restart, time_dt_atmos_restart, time_dt_atmos_restart_unit )
580  call calendar_unit2sec( time_dtsec_ocean, time_dt_ocean, time_dt_ocean_unit )
581  call calendar_unit2sec( time_dtsec_ocean_restart, time_dt_ocean_restart, time_dt_ocean_restart_unit )
582  call calendar_unit2sec( time_dtsec_land, time_dt_land, time_dt_land_unit )
583  call calendar_unit2sec( time_dtsec_land_restart, time_dt_land_restart, time_dt_land_restart_unit )
584  call calendar_unit2sec( time_dtsec_urban, time_dt_urban, time_dt_urban_unit )
585  call calendar_unit2sec( time_dtsec_urban_restart, time_dt_urban_restart, time_dt_urban_restart_unit )
586  call calendar_unit2sec( time_dtsec_resume, time_dt_resume, time_dt_resume_unit )
587 
589  nstep_dyn = time_nstep_atmos_dyn
590 
599  time_dtsec_atmos_restart = max( time_dtsec_atmos_restart, time_dtsec_atmos_dyn*time_nstep_atmos_dyn )
601  time_dtsec_ocean_restart = max( time_dtsec_ocean_restart, time_dtsec_atmos_dyn*time_nstep_atmos_dyn )
603  time_dtsec_land_restart = max( time_dtsec_land_restart, time_dtsec_atmos_dyn*time_nstep_atmos_dyn )
605  time_dtsec_urban_restart = max( time_dtsec_urban_restart, time_dtsec_atmos_dyn*time_nstep_atmos_dyn )
606  time_dtsec_resume = max( time_dtsec_resume, time_dtsec_atmos_dyn*time_nstep_atmos_dyn )
607 
619  time_dstep_atmos_restart = nint( time_dtsec_atmos_restart / time_dtsec )
620  time_dstep_ocean_restart = nint( time_dtsec_ocean_restart / time_dtsec )
621  time_dstep_land_restart = nint( time_dtsec_land_restart / time_dtsec )
622  time_dstep_urban_restart = nint( time_dtsec_urban_restart / time_dtsec )
623  time_dstep_resume = nint( time_dtsec_resume / time_dtsec )
624 
625  time_res_resume = time_dstep_resume - 1
626 
627  if ( abs( real(TIME_NSTEP_ATMOS_DYN,kind=dp)*time_dtsec_atmos_dyn &
628  - real(TIME_DSTEP_ATMOS_DYN,kind=dp)*time_dtsec ) > eps ) then
629  write(*,*) 'xxx delta t(ATMOS_DYN) must be a multiple of delta t ', &
630  time_dtsec_atmos_dyn, real(TIME_DSTEP_ATMOS_DYN,kind=dp)*time_dtsec
631  call prc_mpistop
632  endif
633  if ( abs(time_dtsec_atmos_phy_cp-real(TIME_DSTEP_ATMOS_PHY_CP,kind=dp)*time_dtsec) > eps ) then
634  write(*,*) 'xxx delta t(ATMOS_PHY_CP) must be a multiple of delta t ', &
635  time_dtsec_atmos_phy_cp, real(TIME_DSTEP_ATMOS_PHY_CP,kind=dp)*time_dtsec
636  call prc_mpistop
637  endif
638  if ( abs(time_dtsec_atmos_phy_mp-real(TIME_DSTEP_ATMOS_PHY_MP,kind=dp)*time_dtsec) > eps ) then
639  write(*,*) 'xxx delta t(ATMOS_PHY_MP) must be a multiple of delta t ', &
640  time_dtsec_atmos_phy_mp, real(TIME_DSTEP_ATMOS_PHY_MP,kind=dp)*time_dtsec
641  call prc_mpistop
642  endif
643  if ( abs(time_dtsec_atmos_phy_rd-real(TIME_DSTEP_ATMOS_PHY_RD,kind=dp)*time_dtsec) > eps ) then
644  write(*,*) 'xxx delta t(ATMOS_PHY_RD) must be a multiple of delta t ', &
645  time_dtsec_atmos_phy_rd, real(TIME_DSTEP_ATMOS_PHY_RD,kind=dp)*time_dtsec
646  call prc_mpistop
647  endif
648  if ( abs(time_dtsec_atmos_phy_sf-real(TIME_DSTEP_ATMOS_PHY_SF,kind=dp)*time_dtsec) > eps ) then
649  write(*,*) 'xxx delta t(ATMOS_PHY_SF) must be a multiple of delta t ', &
650  time_dtsec_atmos_phy_sf, real(TIME_DSTEP_ATMOS_PHY_SF,kind=dp)*time_dtsec
651  call prc_mpistop
652  endif
653  if ( abs(time_dtsec_atmos_phy_tb-real(TIME_DSTEP_ATMOS_PHY_TB,kind=dp)*time_dtsec) > eps ) then
654  write(*,*) 'xxx delta t(ATMOS_PHY_TB) must be a multiple of delta t ', &
655  time_dtsec_atmos_phy_tb, real(TIME_DSTEP_ATMOS_PHY_TB,kind=dp)*time_dtsec
656  call prc_mpistop
657  endif
658  if ( abs(time_dtsec_atmos_phy_ch-real(TIME_DSTEP_ATMOS_PHY_CH,kind=dp)*time_dtsec) > eps ) then
659  write(*,*) 'xxx delta t(ATMOS_PHY_CH) must be a multiple of delta t ', &
660  time_dtsec_atmos_phy_ch, real(TIME_DSTEP_ATMOS_PHY_CH,kind=dp)*time_dtsec
661  call prc_mpistop
662  endif
663  if ( abs(time_dtsec_atmos_phy_ae-real(TIME_DSTEP_ATMOS_PHY_AE,kind=dp)*time_dtsec) > eps ) then
664  write(*,*) 'xxx delta t(ATMOS_PHY_AE) must be a multiple of delta t ', &
665  time_dtsec_atmos_phy_ae, real(TIME_DSTEP_ATMOS_PHY_AE,kind=dp)*time_dtsec
666  call prc_mpistop
667  endif
668  if ( abs(time_dtsec_ocean-real(TIME_DSTEP_OCEAN,kind=dp)*time_dtsec) > eps ) then
669  write(*,*) 'xxx delta t(OCEAN) must be a multiple of delta t ', &
670  time_dtsec_ocean, real(TIME_DSTEP_OCEAN,kind=dp)*time_dtsec
671  call prc_mpistop
672  endif
673  if ( abs(time_dtsec_land-real(TIME_DSTEP_LAND,kind=dp)*time_dtsec) > eps ) then
674  write(*,*) 'xxx delta t(LAND) must be a multiple of delta t ', &
675  time_dtsec_land, real(TIME_DSTEP_LAND,kind=dp)*time_dtsec
676  call prc_mpistop
677  endif
678  if ( abs(time_dtsec_urban-real(TIME_DSTEP_URBAN,kind=dp)*time_dtsec) > eps ) then
679  write(*,*) 'xxx delta t(URBAN) must be a multiple of delta t ', &
680  time_dtsec_urban, real(TIME_DSTEP_URBAN,kind=dp)*time_dtsec
681  call prc_mpistop
682  endif
683  if ( abs(time_dtsec_atmos_restart-real(TIME_DSTEP_ATMOS_RESTART,kind=dp)*time_dtsec) > eps ) then
684  write(*,*) 'xxx delta t(ATMOS_RESTART) must be a multiple of delta t ', &
685  time_dtsec_atmos_restart, real(TIME_DSTEP_ATMOS_RESTART,kind=dp)*time_dtsec
686  call prc_mpistop
687  endif
688  if ( abs(time_dtsec_ocean_restart-real(TIME_DSTEP_OCEAN_RESTART,kind=dp)*time_dtsec) > eps ) then
689  write(*,*) 'xxx delta t(OCEAN_RESTART) must be a multiple of delta t ', &
690  time_dtsec_ocean_restart, real(TIME_DSTEP_OCEAN_RESTART,kind=dp)*time_dtsec
691  call prc_mpistop
692  endif
693  if ( abs(time_dtsec_land_restart-real(TIME_DSTEP_LAND_RESTART,kind=dp)*time_dtsec) > eps ) then
694  write(*,*) 'xxx delta t(LAND_RESTART) must be a multiple of delta t ', &
695  time_dtsec_land_restart, real(TIME_DSTEP_LAND_RESTART,kind=dp)*time_dtsec
696  call prc_mpistop
697  endif
698  if ( abs(time_dtsec_urban_restart-real(TIME_DSTEP_URBAN_RESTART,kind=dp)*time_dtsec) > eps ) then
699  write(*,*) 'xxx delta t(URBAN_RESTART) must be a multiple of delta t ', &
700  time_dtsec_urban_restart, real(TIME_DSTEP_URBAN_RESTART,kind=dp)*time_dtsec
701  call prc_mpistop
702  endif
703  if ( abs(time_dtsec_urban_restart-real(TIME_DSTEP_URBAN_RESTART,kind=dp)*time_dtsec) > eps ) then
704  write(*,*) 'xxx delta t(URBAN_RESTART) must be a multiple of delta t ', &
705  time_dtsec_urban_restart, real(TIME_DSTEP_URBAN_RESTART,kind=dp)*time_dtsec
706  call prc_mpistop
707  endif
708  if ( abs(time_dtsec_resume-real(TIME_DSTEP_RESUME,kind=dp)*time_dtsec) > eps ) then
709  write(*,*) 'xxx delta t(RESUME) must be a multiple of delta t ', &
710  time_dtsec_resume, real(TIME_DSTEP_RESUME,kind=dp)*time_dtsec
711  call prc_mpistop
712  endif
713 
714  if( io_l ) write(io_fid_log,*)
715  if( io_l ) write(io_fid_log,*) '*** Time interval for each component (sec.)'
716  if( io_l ) write(io_fid_log,*) '*** Atmosphere'
717  if( io_l ) write(io_fid_log,'(1x,A,F10.3)') '*** Dynamics (time) : ', time_dtsec_atmos_dyn
718  if( io_l ) write(io_fid_log,'(1x,A,I10,A,I8,A)') '*** (step) : ', time_nstep_atmos_dyn, &
719  ' (step interval=', time_dstep_atmos_dyn, ')'
720  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Physics, Cumulus : ', time_dtsec_atmos_phy_cp, &
721  ' (step interval=', time_dstep_atmos_phy_cp, ')'
722  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Physics, Cloud Microphysics : ', time_dtsec_atmos_phy_mp, &
723  ' (step interval=', time_dstep_atmos_phy_mp, ')'
724  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Physics, Radiation : ', time_dtsec_atmos_phy_rd, &
725  ' (step interval=', time_dstep_atmos_phy_rd, ')'
726  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Physics, Surface Flux : ', time_dtsec_atmos_phy_sf, &
727  ' (step interval=', time_dstep_atmos_phy_sf, ')'
728  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Physics, Turbulence : ', time_dtsec_atmos_phy_tb, &
729  ' (step interval=', time_dstep_atmos_phy_tb, ')'
730  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Physics, Chemistry : ', time_dtsec_atmos_phy_ch, &
731  ' (step interval=', time_dstep_atmos_phy_ch, ')'
732  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Physics, Aerosol : ', time_dtsec_atmos_phy_ae, &
733  ' (step interval=', time_dstep_atmos_phy_ae, ')'
734  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Ocean : ', time_dtsec_ocean, &
735  ' (step interval=', time_dstep_ocean, ')'
736  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Land : ', time_dtsec_land, &
737  ' (step interval=', time_dstep_land, ')'
738  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Urban : ', time_dtsec_urban, &
739  ' (step interval=', time_dstep_urban, ')'
740  if( io_l ) write(io_fid_log,*)
741  if( io_l ) write(io_fid_log,*) '*** Time interval for restart (sec.)'
742  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Atmospheric Variables : ', time_dtsec_atmos_restart, &
743  ' (step interval=', time_dstep_atmos_restart, ')'
744  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Ocean Variables : ', time_dtsec_ocean_restart, &
745  ' (step interval=', time_dstep_ocean_restart, ')'
746  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Land Variables : ', time_dtsec_land_restart, &
747  ' (step interval=', time_dstep_land_restart, ')'
748  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Urban Variables : ', time_dtsec_urban_restart, &
749  ' (step interval=', time_dstep_urban_restart, ')'
750  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Resume : ', time_dtsec_resume, &
751  ' (step interval=', time_dstep_resume, ')'
752  else
753  time_dtsec = 1.0_rp
754  endif
755 
756  time_wallclock_start = prc_mpitime()
757 
758  ! WALLCLOCK TERMINATOR SETUP
759  if ( time_wallclock_limit > 0.0_dp ) then
760  if( io_l ) write(io_fid_log,*)
761  if( io_l ) write(io_fid_log,*) '*** Wall clock time limit of execution is specified.'
762 
763  if ( time_dt_wallclock_check == undef8 ) then
764  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_WALLCLOCK_CHECK. ', &
765  ' largest time step interval is used.'
775  time_dtsec_land, &
777  else
778  if ( time_dt_wallclock_check_unit == '' ) then
779  if( io_l ) write(io_fid_log,'(1x,2A)') '*** Not found TIME_DT_WALLCLOCK_CHECK_UNIT. ', &
780  'TIME_DURATION_UNIT is used.'
781  time_dt_wallclock_check_unit = time_duration_unit
782  endif
783  call calendar_unit2sec( time_dtsec_wallclock_check, time_dt_wallclock_check, time_dt_wallclock_check_unit )
785  endif
786 
788 
789  time_wallclock_safe = max( min( time_wallclock_safe, 1.0_dp ), 0.0_dp )
790  time_wallclock_safelim = time_wallclock_limit * time_wallclock_safe
791 
792  if( io_l ) write(io_fid_log,'(1x,A,F10.1,A)') '*** This job stops after ', time_wallclock_safelim, ' seconds.'
793  if( io_l ) write(io_fid_log,'(1x,A,F10.3,A,I8,A)') '*** Time interval for check : ', time_dtsec_wallclock_check, &
794  ' (step interval=', time_dstep_wallclock_check, ')'
795  endif
796 
797  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
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 input 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 803 of file mod_admin_time.f90.

References scale_calendar::calendar_date2char(), scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_step_to_stdout, scale_process::prc_mpitime(), scale_process::prc_universal_ismaster, 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().

803  use scale_process, only: &
806  use scale_calendar, only: &
808  use scale_time, only: &
809  time_nowdate, &
810  time_nowms, &
811  time_nowstep, &
812  time_nstep, &
822  time_dstep_land, &
824  implicit none
825 
826  real(DP) :: WALLCLOCK_elapse
827  character(len=27) :: nowchardate
828  logical :: TO_STDOUT
829  !---------------------------------------------------------------------------
830 
831  time_doatmos_step = .false.
832  time_doatmos_dyn = .false.
833  time_doatmos_phy_cp = .false.
834  time_doatmos_phy_mp = .false.
835  time_doatmos_phy_rd = .false.
836  time_doatmos_phy_sf = .false.
837  time_doatmos_phy_tb = .false.
838  time_doatmos_phy_ch = .false.
839  time_doatmos_phy_ae = .false.
840  time_doocean_step = .false.
841  time_doland_step = .false.
842  time_dourban_step = .false.
843  time_doresume = .false.
844 
845  time_res_atmos_dyn = time_res_atmos_dyn + 1
846  time_res_atmos_phy_cp = time_res_atmos_phy_cp + 1
847  time_res_atmos_phy_mp = time_res_atmos_phy_mp + 1
848  time_res_atmos_phy_rd = time_res_atmos_phy_rd + 1
849  time_res_atmos_phy_sf = time_res_atmos_phy_sf + 1
850  time_res_atmos_phy_tb = time_res_atmos_phy_tb + 1
851  time_res_atmos_phy_ch = time_res_atmos_phy_ch + 1
852  time_res_atmos_phy_ae = time_res_atmos_phy_ae + 1
853  time_res_ocean = time_res_ocean + 1
854  time_res_land = time_res_land + 1
855  time_res_urban = time_res_urban + 1
856  time_res_resume = time_res_resume + 1
857 
858  if ( time_res_atmos_dyn == time_dstep_atmos_dyn ) then
859  time_doatmos_step = .true.
860  time_doatmos_dyn = .true.
861  time_res_atmos_dyn = 0
862  endif
863  if ( time_res_atmos_phy_cp == time_dstep_atmos_phy_cp ) then
864  time_doatmos_step = .true.
865  time_doatmos_phy_cp = .true.
866  time_res_atmos_phy_cp = 0
867  endif
868  if ( time_res_atmos_phy_mp == time_dstep_atmos_phy_mp ) then
869  time_doatmos_step = .true.
870  time_doatmos_phy_mp = .true.
871  time_res_atmos_phy_mp = 0
872  endif
873  if ( time_res_atmos_phy_rd == time_dstep_atmos_phy_rd ) then
874  time_doatmos_step = .true.
875  time_doatmos_phy_rd = .true.
876  time_res_atmos_phy_rd = 0
877  endif
878  if ( time_res_atmos_phy_sf == time_dstep_atmos_phy_sf ) then
879  time_doatmos_step = .true.
880  time_doatmos_phy_sf = .true.
881  time_res_atmos_phy_sf = 0
882  endif
883  if ( time_res_atmos_phy_tb == time_dstep_atmos_phy_tb ) then
884  time_doatmos_step = .true.
885  time_doatmos_phy_tb = .true.
886  time_res_atmos_phy_tb = 0
887  endif
888  if ( time_res_atmos_phy_ch == time_dstep_atmos_phy_ch ) then
889  time_doatmos_step = .true.
890  time_doatmos_phy_ch = .true.
891  time_res_atmos_phy_ch = 0
892  endif
893  if ( time_res_atmos_phy_ae == time_dstep_atmos_phy_ae ) then
894  time_doatmos_step = .true.
895  time_doatmos_phy_ae = .true.
896  time_res_atmos_phy_ae = 0
897  endif
898 
899  if ( time_res_ocean == time_dstep_ocean ) then
900  time_doocean_step = .true.
901  time_res_ocean = 0
902  endif
903  if ( time_res_land == time_dstep_land ) then
904  time_doland_step = .true.
905  time_res_land = 0
906  endif
907  if ( time_res_urban == time_dstep_urban ) then
908  time_dourban_step = .true.
909  time_res_urban = 0
910  endif
911  if ( time_res_resume == time_dstep_resume ) then
912  time_doresume = .true.
913  time_res_resume = 0
914  endif
915 
916  to_stdout = .false.
917  if ( io_step_to_stdout > 0 ) then
918  if( mod(time_nowstep-1,io_step_to_stdout) == 0 ) to_stdout = .true.
919  endif
920 
921  call calendar_date2char( nowchardate, & ! [OUT]
922  time_nowdate(:), & ! [IN]
923  time_nowms ) ! [IN]
924 
925  wallclock_elapse = prc_mpitime() - time_wallclock_start
926 
927  if ( time_wallclock_limit > 0.0_dp ) then
928  if( io_l ) write(io_fid_log,'(1x,2A,2(A,I7),2(A,F10.1))') '*** TIME: ', nowchardate,' STEP:',time_nowstep, '/', time_nstep, &
929  ' WCLOCK:', wallclock_elapse, '/', time_wallclock_safelim
930  if ( prc_universal_ismaster .AND. to_stdout ) then ! universal master node
931  if( io_l ) write(*,'(1x,2A,2(A,I7),2(A,F10.1))') '*** TIME: ', nowchardate,' STEP:',time_nowstep, '/', time_nstep, &
932  ' WCLOCK:', wallclock_elapse, '/', time_wallclock_safelim
933  endif
934  else
935  if( io_l ) write(io_fid_log,'(1x,2A,2(A,I7),A,F10.1)') '*** TIME: ', nowchardate,' STEP:',time_nowstep, '/', time_nstep, &
936  ' WCLOCK:', wallclock_elapse
937  if ( prc_universal_ismaster .AND. to_stdout ) then ! universal master node
938  if( io_l ) write(*,'(1x,2A,2(A,I7),A,F10.1)') '*** TIME: ', nowchardate,' STEP:',time_nowstep, '/', time_nstep, &
939  ' WCLOCK:', wallclock_elapse
940  endif
941  endif
942 
943  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
logical, public prc_universal_ismaster
master process in universal communicator?
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 949 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().

949  use scale_process, only: &
950  prc_ismaster, &
952  use scale_calendar, only: &
956  use scale_time, only: &
957  time_dtsec, &
958  time_nowdate, &
959  time_nowms, &
960  time_nowday, &
961  time_nowsec, &
962  time_nowdaysec, &
963  time_nowstep, &
964  time_nstep, &
967  implicit none
968 
969  real(DP) :: WALLCLOCK_elapse
970  logical :: exists
971  !---------------------------------------------------------------------------
972 
973  time_doend = .false.
974 
976  time_nowday = time_startday
977  time_nowsec = time_startsec + real(TIME_NOWSTEP-1,kind=DP) * TIME_DTSEC
978 
979  ! reallocate day & sub-day
981 
982  call calendar_daysec2date( time_nowdate(:), & ! [OUT]
983  time_nowms, & ! [OUT]
984  time_nowday, & ! [IN]
985  time_nowsec, & ! [IN]
986  time_offset_year ) ! [IN]
987 
989 
990  if ( time_nowstep > time_nstep ) then
991  time_doend = .true.
992  endif
993 
994  if ( time_wallclock_limit > 0.0_dp &
995  .AND. mod(time_nowstep-1,time_dstep_wallclock_check) == 0 ) then
996  wallclock_elapse = prc_mpitime() - time_wallclock_start
997 
998  if ( wallclock_elapse > time_wallclock_safelim ) then
999  if( io_l ) write(io_fid_log,*)
1000  if( io_l ) write(io_fid_log,*) '********************************************************************'
1001  if( io_l ) write(io_fid_log,*) '*** Elapse time limit is detected. Termination operation starts. ***'
1002  if( io_l ) write(io_fid_log,*) '********************************************************************'
1003  if( io_l ) write(io_fid_log,*)
1004  time_doend = .true.
1005  endif
1006  endif
1007 
1008  ! QUIT file control
1009  if ( prc_ismaster ) then ! master node
1010  inquire(file='QUIT', exist=exists)
1011 
1012  if( exists ) then
1013  if( io_l ) write(io_fid_log,*)
1014  if( io_l ) write(io_fid_log,*) '*********************************************************'
1015  if( io_l ) write(io_fid_log,*) '*** QUIT file is found. Termination operation starts. ***'
1016  if( io_l ) write(io_fid_log,*) '*********************************************************'
1017  if( io_l ) write(io_fid_log,*)
1018  time_doend = .true.
1019  endif
1020  endif
1021 
1022  time_doatmos_restart = .false.
1023  time_doocean_restart = .false.
1024  time_doland_restart = .false.
1025  time_dourban_restart = .false.
1026 
1027  time_res_atmos_restart = time_res_atmos_restart + 1
1028  time_res_ocean_restart = time_res_ocean_restart + 1
1029  time_res_land_restart = time_res_land_restart + 1
1030  time_res_urban_restart = time_res_urban_restart + 1
1031 
1032  if ( time_res_atmos_restart == time_dstep_atmos_restart ) then
1033  time_doatmos_restart = .true.
1034  time_res_atmos_restart = 0
1035  elseif( time_doend ) then
1036  time_doatmos_restart = .true.
1037  endif
1038 
1039  if ( time_res_ocean_restart == time_dstep_ocean_restart ) then
1040  time_doocean_restart = .true.
1041  time_res_ocean_restart = 0
1042  elseif( time_doend ) then
1043  time_doocean_restart = .true.
1044  endif
1045 
1046  if ( time_res_land_restart == time_dstep_land_restart ) then
1047  time_doland_restart = .true.
1048  time_res_land_restart = 0
1049  elseif( time_doend ) then
1050  time_doland_restart = .true.
1051  endif
1052 
1053  if ( time_res_urban_restart == time_dstep_urban_restart ) then
1054  time_dourban_restart = .true.
1055  time_res_urban_restart = 0
1056  elseif( time_doend ) then
1057  time_dourban_restart = .true.
1058  endif
1059 
1060  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_write(), admin_time_advance(), and mod_mkinit::mkinit().

54  logical, public :: TIME_DOATMOS_restart

◆ 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_write(), admin_time_advance(), and mod_mkinit::mkinit().

56  logical, public :: TIME_DOOCEAN_restart

◆ 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_write(), admin_time_advance(), and mod_mkinit::mkinit().

58  logical, public :: TIME_DOLAND_restart

◆ 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_write(), admin_time_advance(), and mod_mkinit::mkinit().

60  logical, public :: TIME_DOURBAN_restart

◆ 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