MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
tdismodule Module Reference

Functions/Subroutines

subroutine, public tdis_cr (fname, inmempath)
 Create temporal discretization. More...
 
subroutine, public tdis_set_counters ()
 Set kstp and kper. More...
 
subroutine, public tdis_set_timestep ()
 Set time step length. More...
 
subroutine, public tdis_delt_reset (deltnew)
 Reset delt and update timing variables and indicators. More...
 
subroutine tdis_set_delt ()
 Set time step length. More...
 
subroutine, public tdis_ot (iout)
 Print simulation time. More...
 
subroutine, public tdis_da ()
 Deallocate memory. More...
 
subroutine tdis_source_options ()
 Source the timing discretization options. More...
 
subroutine tdis_allocate_scalars ()
 Allocate tdis scalars. More...
 
subroutine tdis_allocate_arrays ()
 Allocate tdis arrays. More...
 
subroutine tdis_source_dimensions ()
 Source dimension NPER. More...
 
subroutine tdis_source_timing ()
 Source timing information. More...
 
subroutine check_tdis_timing (nper, perlen, nstp, tsmult)
 Check the tdis timing information. More...
 

Variables

integer(i4b), pointer, public nper => null()
 number of stress period More...
 
integer(i4b), pointer, public itmuni => null()
 flag indicating time units More...
 
integer(i4b), pointer, public kper => null()
 current stress period number More...
 
integer(i4b), pointer, public kstp => null()
 current time step number More...
 
integer(i4b), pointer, public inats => null()
 flag indicating ats active for simulation More...
 
logical(lgp), pointer, public readnewdata => null()
 flag indicating time to read new data More...
 
logical(lgp), pointer, public endofperiod => null()
 flag indicating end of stress period More...
 
logical(lgp), pointer, public endofsimulation => null()
 flag indicating end of simulation More...
 
real(dp), pointer, public delt => null()
 length of the current time step More...
 
real(dp), pointer, public pertim => null()
 time relative to start of stress period More...
 
real(dp), pointer, public topertim => null()
 simulation time at start of stress period More...
 
real(dp), pointer, public totim => null()
 time relative to start of simulation More...
 
real(dp), pointer, public totimc => null()
 simulation time at start of time step More...
 
real(dp), pointer, public deltsav => null()
 saved value for delt, used for subtiming More...
 
real(dp), pointer, public totimsav => null()
 saved value for totim, used for subtiming More...
 
real(dp), pointer, public pertimsav => null()
 saved value for pertim, used for subtiming More...
 
real(dp), pointer, public totalsimtime => null()
 time at end of simulation More...
 
real(dp), dimension(:), pointer, public, contiguous perlen => null()
 length of each stress period More...
 
integer(i4b), dimension(:), pointer, public, contiguous nstp => null()
 number of time steps in each stress period More...
 
real(dp), dimension(:), pointer, public, contiguous tsmult => null()
 time step multiplier for each stress period More...
 
character(len=lendatetime), pointer, public datetime0 => null()
 starting date and time for the simulation More...
 
character(len=lenmempath), pointer input_mempath => null()
 input context mempath for tdis More...
 
character(len=linelength), pointer input_fname => null()
 input filename for tdis More...
 
class(atstype), pointer, public ats => null()
 

Function/Subroutine Documentation

◆ check_tdis_timing()

subroutine tdismodule::check_tdis_timing ( integer(i4b), intent(in)  nper,
real(dp), dimension(:), intent(in), contiguous  perlen,
integer(i4b), dimension(:), intent(in), contiguous  nstp,
real(dp), dimension(:), intent(in), contiguous  tsmult 
)

Return back to tdis_read_timing if an error condition is found and let the ustop routine be called there instead so the StoreErrorUnit routine can be called to assign the correct file name.

Definition at line 600 of file tdis.f90.

601  ! -- modules
602  use constantsmodule, only: linelength, dzero, done
603  use simmodule, only: store_error
604  ! -- dummy
605  integer(I4B), intent(in) :: nper
606  real(DP), dimension(:), contiguous, intent(in) :: perlen
607  integer(I4B), dimension(:), contiguous, intent(in) :: nstp
608  real(DP), dimension(:), contiguous, intent(in) :: tsmult
609  ! -- local
610  integer(I4B) :: kper, kstp
611  real(DP) :: tstart, tend, dt
612  character(len=LINELENGTH) :: errmsg
613  ! -- formats
614  character(len=*), parameter :: fmtpwarn = &
615  "(1X,/1X,'PERLEN is zero for stress period ', I0, &
616  &'. PERLEN must not be zero for transient periods.')"
617  character(len=*), parameter :: fmtsperror = &
618  &"(A,' for stress period ', I0)"
619  character(len=*), parameter :: fmtdterror = &
620  "('Time step length of ', G0, ' is too small in period ', I0, &
621  &' and time step ', I0)"
622  !
623  ! -- Initialize
624  tstart = dzero
625  !
626  ! -- Go through and check each stress period
627  do kper = 1, nper
628  !
629  ! -- Error if nstp less than or equal to zero
630  if (nstp(kper) <= 0) then
631  write (errmsg, fmtsperror) 'Number of time steps less than one ', kper
632  call store_error(errmsg)
633  return
634  end if
635  !
636  ! -- Warn if perlen is zero
637  if (perlen(kper) == dzero) then
638  write (iout, fmtpwarn) kper
639  return
640  end if
641  !
642  ! -- Error if tsmult is less than zero
643  if (tsmult(kper) <= dzero) then
644  write (errmsg, fmtsperror) 'TSMULT must be greater than 0.0 ', kper
645  call store_error(errmsg)
646  return
647  end if
648  !
649  ! -- Error if negative period length
650  if (perlen(kper) < dzero) then
651  write (errmsg, fmtsperror) 'PERLEN cannot be less than 0.0 ', kper
652  call store_error(errmsg)
653  return
654  end if
655  !
656  ! -- Go through all time step lengths and make sure they are valid
657  do kstp = 1, nstp(kper)
658  if (kstp == 1) then
659  dt = perlen(kper) / float(nstp(kper))
660  if (tsmult(kper) /= done) &
661  dt = perlen(kper) * (done - tsmult(kper)) / &
662  (done - tsmult(kper)**nstp(kper))
663  else
664  dt = dt * tsmult(kper)
665  end if
666  tend = tstart + dt
667  !
668  ! -- Error condition if tstart == tend
669  if (tstart == tend) then
670  write (errmsg, fmtdterror) dt, kper, kstp
671  call store_error(errmsg)
672  return
673  end if
674  end do
675  !
676  ! -- reset tstart = tend
677  tstart = tend
678  !
679  end do
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
Here is the call graph for this function:
Here is the caller graph for this function:

◆ tdis_allocate_arrays()

subroutine tdismodule::tdis_allocate_arrays

Definition at line 514 of file tdis.f90.

515  ! -- modules
517  !
518  call mem_allocate(perlen, nper, 'PERLEN', 'TDIS')
519  call mem_allocate(nstp, nper, 'NSTP', 'TDIS')
520  call mem_allocate(tsmult, nper, 'TSMULT', 'TDIS')
Here is the caller graph for this function:

◆ tdis_allocate_scalars()

subroutine tdismodule::tdis_allocate_scalars

Definition at line 462 of file tdis.f90.

463  ! -- modules
465  use constantsmodule, only: dzero
466  !
467  ! -- memory manager variables
468  call mem_allocate(nper, 'NPER', 'TDIS')
469  call mem_allocate(itmuni, 'ITMUNI', 'TDIS')
470  call mem_allocate(kper, 'KPER', 'TDIS')
471  call mem_allocate(kstp, 'KSTP', 'TDIS')
472  call mem_allocate(inats, 'INATS', 'TDIS')
473  call mem_allocate(readnewdata, 'READNEWDATA', 'TDIS')
474  call mem_allocate(endofperiod, 'ENDOFPERIOD', 'TDIS')
475  call mem_allocate(endofsimulation, 'ENDOFSIMULATION', 'TDIS')
476  call mem_allocate(delt, 'DELT', 'TDIS')
477  call mem_allocate(pertim, 'PERTIM', 'TDIS')
478  call mem_allocate(topertim, 'TOPERTIM', 'TDIS')
479  call mem_allocate(totim, 'TOTIM', 'TDIS')
480  call mem_allocate(totimc, 'TOTIMC', 'TDIS')
481  call mem_allocate(deltsav, 'DELTSAV', 'TDIS')
482  call mem_allocate(totimsav, 'TOTIMSAV', 'TDIS')
483  call mem_allocate(pertimsav, 'PERTIMSAV', 'TDIS')
484  call mem_allocate(totalsimtime, 'TOTALSIMTIME', 'TDIS')
485  !
486  ! -- strings
487  allocate (datetime0)
488  allocate (input_mempath)
489  allocate (input_fname)
490  !
491  ! -- Initialize variables
492  nper = 0
493  itmuni = 0
494  kper = 0
495  kstp = 0
496  inats = 0
497  readnewdata = .true.
498  endofperiod = .true.
499  endofsimulation = .false.
500  delt = dzero
501  pertim = dzero
502  topertim = dzero
503  totim = dzero
504  totimc = dzero
505  deltsav = dzero
506  totimsav = dzero
507  pertimsav = dzero
508  totalsimtime = dzero
509  datetime0 = ''
Here is the caller graph for this function:

◆ tdis_cr()

subroutine, public tdismodule::tdis_cr ( character(len=*), intent(in)  fname,
character(len=*), intent(in)  inmempath 
)

Definition at line 54 of file tdis.f90.

55  ! -- modules
57  use constantsmodule, only: linelength, dzero
58  ! -- dummy
59  character(len=*), intent(in) :: fname
60  character(len=*), intent(in) :: inmempath
61  ! -- local
62  ! -- formats
63  character(len=*), parameter :: fmtheader = &
64  "(1X,/1X,'TDIS -- TEMPORAL DISCRETIZATION PACKAGE,', / &
65  &' VERSION 1 : 11/13/2014 - INPUT READ FROM MEMPATH: ', A)"
66  !
67  ! -- Allocate the scalar variables
68  call tdis_allocate_scalars()
69  !
70  ! -- set input context and fname
71  input_fname = fname
72  input_mempath = inmempath
73  !
74  ! -- Identify package
75  write (iout, fmtheader) input_mempath
76  !
77  ! -- Source options
78  call tdis_source_options()
79  !
80  ! -- Source dimensions and then allocate arrays
81  call tdis_source_dimensions()
82  call tdis_allocate_arrays()
83  !
84  ! -- Source timing
85  call tdis_source_timing()
86  !
87  ats => create_ats()
88  if (inats > 0) then
89  call ats%ats_init(inats, nper)
90  end if
integer(i4b) function, public getunit()
Get a free unit number.
subroutine, public openfile(iu, iout, fname, ftype, fmtarg_opt, accarg_opt, filstat_opt, mode_opt)
Open a file.
Definition: InputOutput.f90:30
Here is the call graph for this function:
Here is the caller graph for this function:

◆ tdis_da()

subroutine, public tdismodule::tdis_da

Definition at line 341 of file tdis.f90.

342  ! -- modules
344  !
345  ! -- ats
346  if (inats > 0) then
347  call ats%ats_da()
348  end if
349  deallocate (ats)
350  !
351  ! -- Scalars
352  call mem_deallocate(nper)
353  call mem_deallocate(itmuni)
354  call mem_deallocate(kper)
355  call mem_deallocate(kstp)
356  call mem_deallocate(inats)
357  call mem_deallocate(readnewdata)
358  call mem_deallocate(endofperiod)
359  call mem_deallocate(endofsimulation)
360  call mem_deallocate(delt)
361  call mem_deallocate(pertim)
362  call mem_deallocate(topertim)
363  call mem_deallocate(totim)
364  call mem_deallocate(totimc)
365  call mem_deallocate(deltsav)
366  call mem_deallocate(totimsav)
367  call mem_deallocate(pertimsav)
368  call mem_deallocate(totalsimtime)
369  !
370  ! -- strings
371  deallocate (datetime0)
372  deallocate (input_mempath)
373  deallocate (input_fname)
374  !
375  ! -- Arrays
376  call mem_deallocate(perlen)
377  call mem_deallocate(nstp)
378  call mem_deallocate(tsmult)
Here is the caller graph for this function:

◆ tdis_delt_reset()

subroutine, public tdismodule::tdis_delt_reset ( real(dp), intent(in)  deltnew)

This routine is called when a timestep fails to converge, and so it is retried using a smaller time step (deltnew).

Definition at line 212 of file tdis.f90.

213  ! -- modules
214  use constantsmodule, only: done, dzero
215  ! -- dummy
216  real(DP), intent(in) :: deltnew
217  ! -- local
218  logical(LGP) :: adaptivePeriod
219  !
220  ! -- Set values
221  adaptiveperiod = ats%isAdaptivePeriod(kper)
222  delt = deltnew
223  totim = totimsav + delt
224  pertim = pertimsav + delt
225  !
226  ! -- Set end of period indicator
227  endofperiod = .false.
228  if (adaptiveperiod) then
229  call ats%ats_set_endofperiod(kper, pertim, perlen(kper), endofperiod)
230  else
231  if (kstp == nstp(kper)) then
232  endofperiod = .true.
233  end if
234  end if
235  !
236  ! -- Set end of simulation indicator
237  if (endofperiod .and. kper == nper) then
238  endofsimulation = .true.
239  totim = totalsimtime
240  end if
Here is the caller graph for this function:

◆ tdis_ot()

subroutine, public tdismodule::tdis_ot ( integer(i4b), intent(in)  iout)

Definition at line 270 of file tdis.f90.

271  ! -- modules
274  ! -- dummy
275  integer(I4B), intent(in) :: iout
276  ! -- local
277  real(DP) :: cnv, delsec, totsec, persec, delmn, delhr, totmn, tothr, &
278  totdy, totyr, permn, perhr, perdy, peryr, deldy, delyr
279  ! -- format
280  character(len=*), parameter :: fmttmsmry = "(1X, ///9X, &
281  &'TIME SUMMARY AT END OF TIME STEP', I5,' IN STRESS PERIOD ', I4)"
282  character(len=*), parameter :: fmttmstpmsg = &
283  &"(21X, ' TIME STEP LENGTH =', G15.6 / &
284  & 21X, ' STRESS PERIOD TIME =', G15.6 / &
285  & 21X, 'TOTAL SIMULATION TIME =', G15.6)"
286  character(len=*), parameter :: fmttottmmsg = &
287  &"(19X, ' SECONDS MINUTES HOURS', 7X, &
288  &'DAYS YEARS'/20X, 59('-'))"
289  character(len=*), parameter :: fmtdelttm = &
290  &"(1X, ' TIME STEP LENGTH', 1P, 5G12.5)"
291  character(len=*), parameter :: fmtpertm = &
292  &"(1X, 'STRESS PERIOD TIME', 1P, 5G12.5)"
293  character(len=*), parameter :: fmttottm = &
294  &"(1X, ' TOTAL TIME', 1P, 5G12.5,/)"
295  !
296  ! -- Write header message for the information that follows
297  write (iout, fmttmsmry) kstp, kper
298  !
299  ! -- Use time unit indicator to get factor to convert to seconds
300  cnv = dzero
301  if (itmuni == 1) cnv = done
302  if (itmuni == 2) cnv = dsixty
303  if (itmuni == 3) cnv = dsecperhr
304  if (itmuni == 4) cnv = dsecperdy
305  if (itmuni == 5) cnv = dsecperyr
306  !
307  ! -- If FACTOR=0 then time units are non-standard
308  if (cnv == dzero) then
309  ! -- Print times in non-standard time units
310  write (iout, fmttmstpmsg) delt, pertim, totim
311  else
312  ! -- Calculate length of time step & elapsed time in seconds
313  delsec = cnv * delt
314  totsec = cnv * totim
315  persec = cnv * pertim
316  !
317  ! -- Calculate times in minutes, hours, days, and years
318  delmn = delsec / dsixty
319  delhr = delmn / dsixty
320  deldy = delhr / dhrperday
321  delyr = deldy / ddyperyr
322  totmn = totsec / dsixty
323  tothr = totmn / dsixty
324  totdy = tothr / dhrperday
325  totyr = totdy / ddyperyr
326  permn = persec / dsixty
327  perhr = permn / dsixty
328  perdy = perhr / dhrperday
329  peryr = perdy / ddyperyr
330  !
331  ! -- Print time step length and elapsed times in all time units
332  write (iout, fmttottmmsg)
333  write (iout, fmtdelttm) delsec, delmn, delhr, deldy, delyr
334  write (iout, fmtpertm) persec, permn, perhr, perdy, peryr
335  write (iout, fmttottm) totsec, totmn, tothr, totdy, totyr
336  end if
real(dp), parameter dsixty
real constant 60
Definition: Constants.f90:85
real(dp), parameter dsecperdy
real constant representing the number of seconds per day (used in tdis)
Definition: Constants.f90:100
real(dp), parameter dsecperyr
real constant representing the average number of seconds per year (used in tdis)
Definition: Constants.f90:101
real(dp), parameter dsecperhr
real constant representing number of seconds per hour (used in tdis)
Definition: Constants.f90:97
real(dp), parameter ddyperyr
real constant representing the average number of days per year (used in tdis)
Definition: Constants.f90:99
real(dp), parameter dhrperday
real constant representing number of hours per day (used in tdis)
Definition: Constants.f90:98
Here is the caller graph for this function:

◆ tdis_set_counters()

subroutine, public tdismodule::tdis_set_counters

Definition at line 95 of file tdis.f90.

96  ! -- modules
98  use simvariablesmodule, only: isim_mode
99  use messagemodule, only: write_message
100  ! -- local
101  character(len=LINELENGTH) :: line
102  character(len=4) :: cpref
103  character(len=10) :: cend
104  ! -- formats
105  character(len=*), parameter :: fmtspts = &
106  &"(a, 'Solving: Stress period: ',i5,4x, 'Time step: ',i5,4x, a)"
107  character(len=*), parameter :: fmtvspts = &
108  &"(' Validating: Stress period: ',i5,4x,'Time step: ',i5,4x)"
109  character(len=*), parameter :: fmtspi = &
110  "('1',/1X,'STRESS PERIOD NO. ',I0,', LENGTH =',G15.7,/ &
111  &1X,42('-'))"
112  character(len=*), parameter :: fmtspits = &
113  "(1X,'NUMBER OF TIME STEPS = ',I0,/ &
114  &1X,'MULTIPLIER FOR DELT =',F10.3)"
115  !
116  ! -- Initialize variables for this step
117  if (inats > 0) ats%dtstable = dnodata
118  readnewdata = .false.
119  cpref = ' '
120  cend = ''
121  !
122  ! -- Increment kstp and kper
123  if (endofperiod) then
124  kstp = 1
125  kper = kper + 1
126  readnewdata = .true.
127  else
128  kstp = kstp + 1
129  end if
130  !
131  ! -- Print stress period and time step to console
132  select case (isim_mode)
133  case (mvalidate)
134  write (line, fmtvspts) kper, kstp
135  case (mnormal)
136  write (line, fmtspts) cpref, kper, kstp, trim(cend)
137  end select
138  if (isim_level >= vall) &
139  call write_message(line)
140  call write_message(line, iunit=iout, skipbefore=1, skipafter=1)
141  !
142  ! -- Write message if first time step
143  if (kstp == 1) then
144  write (iout, fmtspi) kper, perlen(kper)
145  if (ats%isAdaptivePeriod(kper)) then
146  call ats%ats_period_message(kper)
147  else
148  write (iout, fmtspits) nstp(kper), tsmult(kper)
149  end if
150  end if
@ mvalidate
validation mode - do not run time steps
Definition: Constants.f90:205
@ mnormal
normal output mode
Definition: Constants.f90:206
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:95
Store and issue logging messages to output units.
Definition: Message.f90:2
subroutine, public write_message(text, iunit, fmt, skipbefore, skipafter, advance)
Write a message to an output unit.
Definition: Message.f90:210
This module contains simulation variables.
Definition: SimVariables.f90:9
integer(i4b) isim_mode
simulation mode
Here is the call graph for this function:
Here is the caller graph for this function:

◆ tdis_set_delt()

subroutine tdismodule::tdis_set_delt

Definition at line 245 of file tdis.f90.

246  ! -- modules
247  use constantsmodule, only: done
248  !
249  if (kstp == 1) then
250  ! -- Calculate the first value of delt for this stress period
251  topertim = totim
252  if (tsmult(kper) /= done) then
253  ! -- Timestep length has a geometric progression
254  delt = perlen(kper) * (done - tsmult(kper)) / &
255  (done - tsmult(kper)**nstp(kper))
256  else
257  ! -- Timestep length is constant
258  delt = perlen(kper) / float(nstp(kper))
259  end if
260  elseif (kstp == nstp(kper)) then
261  ! -- Calculate exact last delt to avoid accumulation errors
262  delt = topertim + perlen(kper) - totim
263  else
264  delt = tsmult(kper) * delt
265  end if
Here is the caller graph for this function:

◆ tdis_set_timestep()

subroutine, public tdismodule::tdis_set_timestep

Definition at line 155 of file tdis.f90.

156  ! -- modules
157  use constantsmodule, only: done, dzero
158  ! -- local
159  logical(LGP) :: adaptivePeriod
160  ! -- format
161  character(len=*), parameter :: fmttsi = &
162  "(1X,'INITIAL TIME STEP SIZE =',G15.7)"
163  !
164  ! -- Initialize
165  adaptiveperiod = ats%isAdaptivePeriod(kper)
166  if (kstp == 1) then
167  pertim = dzero
168  topertim = dzero
169  end if
170  !
171  ! -- Set delt
172  if (adaptiveperiod) then
173  call ats%ats_set_delt(kstp, kper, pertim, perlen(kper), delt)
174  else
175  call tdis_set_delt()
176  if (kstp == 1) then
177  write (iout, fmttsi) delt
178  end if
179  end if
180  !
181  ! -- Advance timers and update totim and pertim based on delt
182  totimsav = totim
183  pertimsav = pertim
184  totimc = totimsav
185  totim = totimsav + delt
186  pertim = pertimsav + delt
187  !
188  ! -- Set end of period indicator
189  endofperiod = .false.
190  if (adaptiveperiod) then
191  call ats%ats_set_endofperiod(kper, pertim, perlen(kper), endofperiod)
192  else
193  if (kstp == nstp(kper)) then
194  endofperiod = .true.
195  end if
196  end if
197  if (endofperiod) then
198  pertim = perlen(kper)
199  end if
200  !
201  ! -- Set end of simulation indicator
202  if (endofperiod .and. kper == nper) then
203  endofsimulation = .true.
204  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ tdis_source_dimensions()

subroutine tdismodule::tdis_source_dimensions

Definition at line 525 of file tdis.f90.

526  ! -- modules
527  use constantsmodule, only: linelength
531  ! -- local
532  type(SimTdisParamFoundType) :: found
533  ! -- formats
534  character(len=*), parameter :: fmtnper = &
535  "(1X,I4,' STRESS PERIOD(S) IN SIMULATION')"
536  !
537  ! -- source dimensions from input context
538  call mem_set_value(nper, 'NPER', input_mempath, found%nper)
539  !
540  ! -- log values to list file
541  write (iout, '(1x,a)') 'PROCESSING TDIS DIMENSIONS'
542  !
543  if (found%nper) then
544  write (iout, fmtnper) nper
545  end if
546  !
547  write (iout, '(1x,a)') 'END OF TDIS DIMENSIONS'
This module contains the SourceCommonModule.
Definition: SourceCommon.f90:7
logical(lgp) function, public filein_fname(filename, tagname, input_mempath, input_fname)
enforce and set a single input filename provided via FILEIN keyword
Here is the call graph for this function:
Here is the caller graph for this function:

◆ tdis_source_options()

subroutine tdismodule::tdis_source_options

Definition at line 383 of file tdis.f90.

384  ! -- modules
385  use constantsmodule, only: linelength
391  ! -- local
392  type(SimTdisParamFoundType) :: found
393  character(len=LINELENGTH), dimension(6) :: time_units = &
394  &[character(len=LINELENGTH) :: 'UNDEFINED', 'SECONDS', 'MINUTES', 'HOURS', &
395  'DAYS', 'YEARS']
396  character(len=LINELENGTH) :: fname
397  ! -- formats
398  character(len=*), parameter :: fmtitmuni = &
399  &"(4x,'SIMULATION TIME UNIT IS ',A)"
400  character(len=*), parameter :: fmtdatetime0 = &
401  &"(4x,'SIMULATION STARTING DATE AND TIME IS ',A)"
402  !
403  ! -- initialize time unit to undefined
404  itmuni = 0
405  !
406  ! -- source options from input context
407  call mem_set_value(itmuni, 'TIME_UNITS', input_mempath, time_units, &
408  found%time_units)
409  call mem_set_value(datetime0, 'START_DATE_TIME', input_mempath, &
410  found%start_date_time)
411  !
412  if (found%time_units) then
413  if (itmuni == 0) then
414  call store_error('Unrecognized input value for TIME_UNITS option.')
415  call store_error_filename(input_fname)
416  else
417  !
418  ! -- adjust to 0-based indexing for itmuni
419  itmuni = itmuni - 1
420  end if
421  end if
422  !
423  ! -- enforce 0 or 1 ATS6_FILENAME entries in option block
424  if (filein_fname(fname, 'ATS6_FILENAME', input_mempath, &
425  input_fname)) then
426  inats = getunit()
427  call openfile(inats, iout, fname, 'ATS')
428  end if
429  !
430  ! -- log values to list file
431  write (iout, '(1x,a)') 'PROCESSING TDIS OPTIONS'
432  !
433  if (found%time_units) then
434  select case (itmuni)
435  case (0)
436  write (iout, fmtitmuni) 'UNDEFINED'
437  case (1)
438  write (iout, fmtitmuni) 'SECONDS'
439  case (2)
440  write (iout, fmtitmuni) 'MINUTES'
441  case (3)
442  write (iout, fmtitmuni) 'HOURS'
443  case (4)
444  write (iout, fmtitmuni) 'DAYS'
445  case (5)
446  write (iout, fmtitmuni) 'YEARS'
447  case default
448  end select
449  else
450  write (iout, fmtitmuni) 'UNDEFINED'
451  end if
452  !
453  if (found%start_date_time) then
454  write (iout, fmtdatetime0) datetime0
455  end if
456  !
457  write (iout, '(1x,a)') 'END OF TDIS OPTIONS'
subroutine, public store_error_filename(filename, terminate)
Store the erroring file name.
Definition: Sim.f90:203
Here is the call graph for this function:
Here is the caller graph for this function:

◆ tdis_source_timing()

subroutine tdismodule::tdis_source_timing

Definition at line 552 of file tdis.f90.

553  ! -- modules
554  use constantsmodule, only: linelength, dzero
559  ! -- local
560  type(SimTdisParamFoundType) :: found
561  integer(I4B) :: n
562  ! -- formats
563  character(len=*), parameter :: fmtheader = &
564  "(1X,//1X,'STRESS PERIOD LENGTH TIME STEPS', &
565  &' MULTIPLIER FOR DELT',/1X,76('-'))"
566  character(len=*), parameter :: fmtrow = &
567  "(1X,I8,1PG21.7,I7,0PF25.3)"
568  !
569  ! -- source perioddata from input context
570  call mem_set_value(perlen, 'PERLEN', input_mempath, found%perlen)
571  call mem_set_value(nstp, 'NSTP', input_mempath, found%nstp)
572  call mem_set_value(tsmult, 'TSMULT', input_mempath, found%tsmult)
573  !
574  ! -- Check timing information
575  call check_tdis_timing(nper, perlen, nstp, tsmult)
576  !
577  ! -- Check for errors
578  if (count_errors() > 0) then
579  call store_error_filename(input_fname)
580  end if
581  !
582  ! -- log timing
583  write (iout, '(1x,a)') 'PROCESSING TDIS PERIODDATA'
584  write (iout, fmtheader)
585  !
586  do n = 1, size(perlen)
587  write (iout, fmtrow) n, perlen(n), nstp(n), tsmult(n)
588  totalsimtime = totalsimtime + perlen(n)
589  end do
590  !
591  write (iout, '(1x,a)') 'END OF TDIS PERIODDATA'
integer(i4b) function, public count_errors()
Return number of errors.
Definition: Sim.f90:59
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ ats

class(atstype), pointer, public tdismodule::ats => null()

Definition at line 48 of file tdis.f90.

48  class(AtsType), pointer :: ats => null()

◆ datetime0

character(len=lendatetime), pointer, public tdismodule::datetime0 => null()

Definition at line 44 of file tdis.f90.

44  character(len=LENDATETIME), public, pointer :: datetime0 => null() !< starting date and time for the simulation

◆ delt

real(dp), pointer, public tdismodule::delt => null()

Definition at line 32 of file tdis.f90.

32  real(DP), public, pointer :: delt => null() !< length of the current time step

◆ deltsav

real(dp), pointer, public tdismodule::deltsav => null()

Definition at line 37 of file tdis.f90.

37  real(DP), public, pointer :: deltsav => null() !< saved value for delt, used for subtiming

◆ endofperiod

logical(lgp), pointer, public tdismodule::endofperiod => null()

Definition at line 30 of file tdis.f90.

30  logical(LGP), public, pointer :: endofperiod => null() !< flag indicating end of stress period

◆ endofsimulation

logical(lgp), pointer, public tdismodule::endofsimulation => null()

Definition at line 31 of file tdis.f90.

31  logical(LGP), public, pointer :: endofsimulation => null() !< flag indicating end of simulation

◆ inats

integer(i4b), pointer, public tdismodule::inats => null()

Definition at line 28 of file tdis.f90.

28  integer(I4B), public, pointer :: inats => null() !< flag indicating ats active for simulation

◆ input_fname

character(len=linelength), pointer tdismodule::input_fname => null()
private

Definition at line 46 of file tdis.f90.

46  character(len=LINELENGTH), pointer :: input_fname => null() !< input filename for tdis

◆ input_mempath

character(len=lenmempath), pointer tdismodule::input_mempath => null()
private

Definition at line 45 of file tdis.f90.

45  character(len=LENMEMPATH), pointer :: input_mempath => null() !< input context mempath for tdis

◆ itmuni

integer(i4b), pointer, public tdismodule::itmuni => null()

Definition at line 25 of file tdis.f90.

25  integer(I4B), public, pointer :: itmuni => null() !< flag indicating time units

◆ kper

integer(i4b), pointer, public tdismodule::kper => null()

Definition at line 26 of file tdis.f90.

26  integer(I4B), public, pointer :: kper => null() !< current stress period number

◆ kstp

integer(i4b), pointer, public tdismodule::kstp => null()

Definition at line 27 of file tdis.f90.

27  integer(I4B), public, pointer :: kstp => null() !< current time step number

◆ nper

integer(i4b), pointer, public tdismodule::nper => null()

Definition at line 24 of file tdis.f90.

24  integer(I4B), public, pointer :: nper => null() !< number of stress period

◆ nstp

integer(i4b), dimension(:), pointer, public, contiguous tdismodule::nstp => null()

Definition at line 42 of file tdis.f90.

42  integer(I4B), public, dimension(:), pointer, contiguous :: nstp => null() !< number of time steps in each stress period

◆ perlen

real(dp), dimension(:), pointer, public, contiguous tdismodule::perlen => null()

Definition at line 41 of file tdis.f90.

41  real(DP), public, dimension(:), pointer, contiguous :: perlen => null() !< length of each stress period

◆ pertim

real(dp), pointer, public tdismodule::pertim => null()

Definition at line 33 of file tdis.f90.

33  real(DP), public, pointer :: pertim => null() !< time relative to start of stress period

◆ pertimsav

real(dp), pointer, public tdismodule::pertimsav => null()

Definition at line 39 of file tdis.f90.

39  real(DP), public, pointer :: pertimsav => null() !< saved value for pertim, used for subtiming

◆ readnewdata

logical(lgp), pointer, public tdismodule::readnewdata => null()

Definition at line 29 of file tdis.f90.

29  logical(LGP), public, pointer :: readnewdata => null() !< flag indicating time to read new data

◆ topertim

real(dp), pointer, public tdismodule::topertim => null()

Definition at line 34 of file tdis.f90.

34  real(DP), public, pointer :: topertim => null() !< simulation time at start of stress period

◆ totalsimtime

real(dp), pointer, public tdismodule::totalsimtime => null()

Definition at line 40 of file tdis.f90.

40  real(DP), public, pointer :: totalsimtime => null() !< time at end of simulation

◆ totim

real(dp), pointer, public tdismodule::totim => null()

Definition at line 35 of file tdis.f90.

35  real(DP), public, pointer :: totim => null() !< time relative to start of simulation

◆ totimc

real(dp), pointer, public tdismodule::totimc => null()

Definition at line 36 of file tdis.f90.

36  real(DP), public, pointer :: totimc => null() !< simulation time at start of time step

◆ totimsav

real(dp), pointer, public tdismodule::totimsav => null()

Definition at line 38 of file tdis.f90.

38  real(DP), public, pointer :: totimsav => null() !< saved value for totim, used for subtiming

◆ tsmult

real(dp), dimension(:), pointer, public, contiguous tdismodule::tsmult => null()

Definition at line 43 of file tdis.f90.

43  real(DP), public, dimension(:), pointer, contiguous :: tsmult => null() !< time step multiplier for each stress period