MODFLOW 6  version 6.7.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...
 

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 595 of file tdis.f90.

596  ! -- modules
597  use constantsmodule, only: linelength, dzero, done
598  use simmodule, only: store_error
599  ! -- dummy
600  integer(I4B), intent(in) :: nper
601  real(DP), dimension(:), contiguous, intent(in) :: perlen
602  integer(I4B), dimension(:), contiguous, intent(in) :: nstp
603  real(DP), dimension(:), contiguous, intent(in) :: tsmult
604  ! -- local
605  integer(I4B) :: kper, kstp
606  real(DP) :: tstart, tend, dt
607  character(len=LINELENGTH) :: errmsg
608  ! -- formats
609  character(len=*), parameter :: fmtpwarn = &
610  "(1X,/1X,'PERLEN is zero for stress period ', I0, &
611  &'. PERLEN must not be zero for transient periods.')"
612  character(len=*), parameter :: fmtsperror = &
613  &"(A,' for stress period ', I0)"
614  character(len=*), parameter :: fmtdterror = &
615  "('Time step length of ', G0, ' is too small in period ', I0, &
616  &' and time step ', I0)"
617  !
618  ! -- Initialize
619  tstart = dzero
620  !
621  ! -- Go through and check each stress period
622  do kper = 1, nper
623  !
624  ! -- Error if nstp less than or equal to zero
625  if (nstp(kper) <= 0) then
626  write (errmsg, fmtsperror) 'Number of time steps less than one ', kper
627  call store_error(errmsg)
628  return
629  end if
630  !
631  ! -- Warn if perlen is zero
632  if (perlen(kper) == dzero) then
633  write (iout, fmtpwarn) kper
634  return
635  end if
636  !
637  ! -- Error if tsmult is less than zero
638  if (tsmult(kper) <= dzero) then
639  write (errmsg, fmtsperror) 'TSMULT must be greater than 0.0 ', kper
640  call store_error(errmsg)
641  return
642  end if
643  !
644  ! -- Error if negative period length
645  if (perlen(kper) < dzero) then
646  write (errmsg, fmtsperror) 'PERLEN cannot be less than 0.0 ', kper
647  call store_error(errmsg)
648  return
649  end if
650  !
651  ! -- Go through all time step lengths and make sure they are valid
652  do kstp = 1, nstp(kper)
653  if (kstp == 1) then
654  dt = perlen(kper) / float(nstp(kper))
655  if (tsmult(kper) /= done) &
656  dt = perlen(kper) * (done - tsmult(kper)) / &
657  (done - tsmult(kper)**nstp(kper))
658  else
659  dt = dt * tsmult(kper)
660  end if
661  tend = tstart + dt
662  !
663  ! -- Error condition if tstart == tend
664  if (tstart == tend) then
665  write (errmsg, fmtdterror) dt, kper, kstp
666  call store_error(errmsg)
667  return
668  end if
669  end do
670  !
671  ! -- reset tstart = tend
672  tstart = tend
673  !
674  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 509 of file tdis.f90.

510  ! -- modules
512  !
513  call mem_allocate(perlen, nper, 'PERLEN', 'TDIS')
514  call mem_allocate(nstp, nper, 'NSTP', 'TDIS')
515  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 457 of file tdis.f90.

458  ! -- modules
460  use constantsmodule, only: dzero
461  !
462  ! -- memory manager variables
463  call mem_allocate(nper, 'NPER', 'TDIS')
464  call mem_allocate(itmuni, 'ITMUNI', 'TDIS')
465  call mem_allocate(kper, 'KPER', 'TDIS')
466  call mem_allocate(kstp, 'KSTP', 'TDIS')
467  call mem_allocate(inats, 'INATS', 'TDIS')
468  call mem_allocate(readnewdata, 'READNEWDATA', 'TDIS')
469  call mem_allocate(endofperiod, 'ENDOFPERIOD', 'TDIS')
470  call mem_allocate(endofsimulation, 'ENDOFSIMULATION', 'TDIS')
471  call mem_allocate(delt, 'DELT', 'TDIS')
472  call mem_allocate(pertim, 'PERTIM', 'TDIS')
473  call mem_allocate(topertim, 'TOPERTIM', 'TDIS')
474  call mem_allocate(totim, 'TOTIM', 'TDIS')
475  call mem_allocate(totimc, 'TOTIMC', 'TDIS')
476  call mem_allocate(deltsav, 'DELTSAV', 'TDIS')
477  call mem_allocate(totimsav, 'TOTIMSAV', 'TDIS')
478  call mem_allocate(pertimsav, 'PERTIMSAV', 'TDIS')
479  call mem_allocate(totalsimtime, 'TOTALSIMTIME', 'TDIS')
480  !
481  ! -- strings
482  allocate (datetime0)
483  allocate (input_mempath)
484  allocate (input_fname)
485  !
486  ! -- Initialize variables
487  nper = 0
488  itmuni = 0
489  kper = 0
490  kstp = 0
491  inats = 0
492  readnewdata = .true.
493  endofperiod = .true.
494  endofsimulation = .false.
495  delt = dzero
496  pertim = dzero
497  topertim = dzero
498  totim = dzero
499  totimc = dzero
500  deltsav = dzero
501  totimsav = dzero
502  pertimsav = dzero
503  totalsimtime = dzero
504  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 49 of file tdis.f90.

50  ! -- modules
52  use constantsmodule, only: linelength, dzero
54  ! -- dummy
55  character(len=*), intent(in) :: fname
56  character(len=*), intent(in) :: inmempath
57  ! -- local
58  ! -- formats
59  character(len=*), parameter :: fmtheader = &
60  "(1X,/1X,'TDIS -- TEMPORAL DISCRETIZATION PACKAGE,', / &
61  &' VERSION 1 : 11/13/2014 - INPUT READ FROM MEMPATH: ', A)"
62  !
63  ! -- Allocate the scalar variables
64  call tdis_allocate_scalars()
65  !
66  ! -- set input context and fname
67  input_fname = fname
68  input_mempath = inmempath
69  !
70  ! -- Identify package
71  write (iout, fmtheader) input_mempath
72  !
73  ! -- Source options
74  call tdis_source_options()
75  !
76  ! -- Source dimensions and then allocate arrays
77  call tdis_source_dimensions()
78  call tdis_allocate_arrays()
79  !
80  ! -- Source timing
81  call tdis_source_timing()
82  !
83  if (inats > 0) then
84  call ats_cr(inats, nper)
85  end if
subroutine, public ats_cr(inunit, nper_tdis)
@ brief Create ATS object
Definition: ats.f90:61
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 344 of file tdis.f90.

345  ! -- modules
347  use adaptivetimestepmodule, only: ats_da
348  !
349  ! -- ats
350  if (inats > 0) call ats_da()
351  !
352  ! -- Scalars
353  call mem_deallocate(nper)
354  call mem_deallocate(itmuni)
355  call mem_deallocate(kper)
356  call mem_deallocate(kstp)
357  call mem_deallocate(inats)
358  call mem_deallocate(readnewdata)
359  call mem_deallocate(endofperiod)
360  call mem_deallocate(endofsimulation)
361  call mem_deallocate(delt)
362  call mem_deallocate(pertim)
363  call mem_deallocate(topertim)
364  call mem_deallocate(totim)
365  call mem_deallocate(totimc)
366  call mem_deallocate(deltsav)
367  call mem_deallocate(totimsav)
368  call mem_deallocate(pertimsav)
369  call mem_deallocate(totalsimtime)
370  !
371  ! -- strings
372  deallocate (datetime0)
373  deallocate (input_mempath)
374  deallocate (input_fname)
375  !
376  ! -- Arrays
377  call mem_deallocate(perlen)
378  call mem_deallocate(nstp)
379  call mem_deallocate(tsmult)
subroutine, public ats_da()
@ brief Deallocate variables
Definition: ats.f90:167
Here is the call graph for this function:
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
216  ats_set_delt, &
218  ! -- dummy
219  real(DP), intent(in) :: deltnew
220  ! -- local
221  logical(LGP) :: adaptivePeriod
222  !
223  ! -- Set values
224  adaptiveperiod = isadaptiveperiod(kper)
225  delt = deltnew
226  totim = totimsav + delt
227  pertim = pertimsav + delt
228  !
229  ! -- Set end of period indicator
230  endofperiod = .false.
231  if (adaptiveperiod) then
232  call ats_set_endofperiod(kper, pertim, perlen(kper), endofperiod)
233  else
234  if (kstp == nstp(kper)) then
235  endofperiod = .true.
236  end if
237  end if
238  !
239  ! -- Set end of simulation indicator
240  if (endofperiod .and. kper == nper) then
241  endofsimulation = .true.
242  totim = totalsimtime
243  end if
subroutine, public ats_set_delt(kstp, kper, pertim, perlencurrent, delt)
@ brief Set time step
Definition: ats.f90:541
subroutine, public ats_set_endofperiod(kper, pertim, perlencurrent, endofperiod)
@ brief Set end of period indicator
Definition: ats.f90:646
logical(lgp) function, public isadaptiveperiod(kper)
@ brief Determine if period is adaptive
Definition: ats.f90:45
Here is the call graph for this function:
Here is the caller graph for this function:

◆ tdis_ot()

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

Definition at line 273 of file tdis.f90.

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

91  ! -- modules
93  use simvariablesmodule, only: isim_mode
94  use messagemodule, only: write_message
97  ! -- local
98  character(len=LINELENGTH) :: line
99  character(len=4) :: cpref
100  character(len=10) :: cend
101  ! -- formats
102  character(len=*), parameter :: fmtspts = &
103  &"(a, 'Solving: Stress period: ',i5,4x, 'Time step: ',i5,4x, a)"
104  character(len=*), parameter :: fmtvspts = &
105  &"(' Validating: Stress period: ',i5,4x,'Time step: ',i5,4x)"
106  character(len=*), parameter :: fmtspi = &
107  "('1',/1X,'STRESS PERIOD NO. ',I0,', LENGTH =',G15.7,/ &
108  &1X,42('-'))"
109  character(len=*), parameter :: fmtspits = &
110  "(1X,'NUMBER OF TIME STEPS = ',I0,/ &
111  &1X,'MULTIPLIER FOR DELT =',F10.3)"
112  !
113  ! -- Initialize variables for this step
114  if (inats > 0) dtstable = dnodata
115  readnewdata = .false.
116  cpref = ' '
117  cend = ''
118  !
119  ! -- Increment kstp and kper
120  if (endofperiod) then
121  kstp = 1
122  kper = kper + 1
123  readnewdata = .true.
124  else
125  kstp = kstp + 1
126  end if
127  !
128  ! -- Print stress period and time step to console
129  select case (isim_mode)
130  case (mvalidate)
131  write (line, fmtvspts) kper, kstp
132  case (mnormal)
133  write (line, fmtspts) cpref, kper, kstp, trim(cend)
134  end select
135  if (isim_level >= vall) &
136  call write_message(line)
137  call write_message(line, iunit=iout, skipbefore=1, skipafter=1)
138  !
139  ! -- Write message if first time step
140  if (kstp == 1) then
141  write (iout, fmtspi) kper, perlen(kper)
142  if (isadaptiveperiod(kper)) then
143  call ats_period_message(kper)
144  else
145  write (iout, fmtspits) nstp(kper), tsmult(kper)
146  end if
147  end if
real(dp), pointer, public dtstable
delt value required for stability
Definition: ats.f90:26
subroutine, public ats_period_message(kper)
@ brief Write period message
Definition: ats.f90:469
@ 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 248 of file tdis.f90.

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

◆ tdis_set_timestep()

subroutine, public tdismodule::tdis_set_timestep

Definition at line 152 of file tdis.f90.

153  ! -- modules
154  use constantsmodule, only: done, dzero
156  ats_set_delt, &
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 = 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_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_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 520 of file tdis.f90.

521  ! -- modules
522  use constantsmodule, only: linelength
526  ! -- local
527  type(SimTdisParamFoundType) :: found
528  ! -- formats
529  character(len=*), parameter :: fmtnper = &
530  "(1X,I4,' STRESS PERIOD(S) IN SIMULATION')"
531  !
532  ! -- source dimensions from input context
533  call mem_set_value(nper, 'NPER', input_mempath, found%nper)
534  !
535  ! -- log values to list file
536  write (iout, '(1x,a)') 'PROCESSING TDIS DIMENSIONS'
537  !
538  if (found%nper) then
539  write (iout, fmtnper) nper
540  end if
541  !
542  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 384 of file tdis.f90.

385  ! -- modules
386  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  !
414  ! -- adjust to 0-based indexing for itmuni
415  itmuni = itmuni - 1
416  end if
417  !
418  ! -- enforce 0 or 1 ATS6_FILENAME entries in option block
419  if (filein_fname(fname, 'ATS6_FILENAME', input_mempath, &
420  input_fname)) then
421  inats = getunit()
422  call openfile(inats, iout, fname, 'ATS')
423  end if
424  !
425  ! -- log values to list file
426  write (iout, '(1x,a)') 'PROCESSING TDIS OPTIONS'
427  !
428  if (found%time_units) then
429  select case (itmuni)
430  case (0)
431  write (iout, fmtitmuni) 'UNDEFINED'
432  case (1)
433  write (iout, fmtitmuni) 'SECONDS'
434  case (2)
435  write (iout, fmtitmuni) 'MINUTES'
436  case (3)
437  write (iout, fmtitmuni) 'HOURS'
438  case (4)
439  write (iout, fmtitmuni) 'DAYS'
440  case (5)
441  write (iout, fmtitmuni) 'YEARS'
442  case default
443  end select
444  else
445  write (iout, fmtitmuni) 'UNDEFINED'
446  end if
447  !
448  if (found%start_date_time) then
449  write (iout, fmtdatetime0) datetime0
450  end if
451  !
452  write (iout, '(1x,a)') 'END OF TDIS OPTIONS'
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 547 of file tdis.f90.

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

Variable Documentation

◆ datetime0

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

Definition at line 41 of file tdis.f90.

41  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 29 of file tdis.f90.

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

◆ deltsav

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

Definition at line 34 of file tdis.f90.

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

◆ endofperiod

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

Definition at line 27 of file tdis.f90.

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

◆ endofsimulation

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

Definition at line 28 of file tdis.f90.

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

◆ inats

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

Definition at line 25 of file tdis.f90.

25  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 43 of file tdis.f90.

43  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 42 of file tdis.f90.

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

◆ itmuni

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

Definition at line 22 of file tdis.f90.

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

◆ kper

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

Definition at line 23 of file tdis.f90.

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

◆ kstp

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

Definition at line 24 of file tdis.f90.

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

◆ nper

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

Definition at line 21 of file tdis.f90.

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

◆ nstp

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

Definition at line 39 of file tdis.f90.

39  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 38 of file tdis.f90.

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

◆ pertim

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

Definition at line 30 of file tdis.f90.

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

◆ pertimsav

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

Definition at line 36 of file tdis.f90.

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

◆ readnewdata

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

Definition at line 26 of file tdis.f90.

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

◆ topertim

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

Definition at line 31 of file tdis.f90.

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

◆ totalsimtime

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

Definition at line 37 of file tdis.f90.

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

◆ totim

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

Definition at line 32 of file tdis.f90.

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

◆ totimc

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

Definition at line 33 of file tdis.f90.

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

◆ totimsav

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

Definition at line 35 of file tdis.f90.

35  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 40 of file tdis.f90.

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