MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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 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 691 of file tdis.f90.

692  ! -- modules
693  use constantsmodule, only: linelength, dzero, done
694  use simmodule, only: store_error
695  ! -- dummy
696  integer(I4B), intent(in) :: nper
697  real(DP), dimension(:), contiguous, intent(in) :: perlen
698  integer(I4B), dimension(:), contiguous, intent(in) :: nstp
699  real(DP), dimension(:), contiguous, intent(in) :: tsmult
700  ! -- local
701  integer(I4B) :: kper, kstp
702  real(DP) :: tstart, tend, dt
703  character(len=LINELENGTH) :: errmsg
704  ! -- formats
705  character(len=*), parameter :: fmtpwarn = &
706  "(1X,/1X,'PERLEN is zero for stress period ', I0, &
707  &'. PERLEN must not be zero for transient periods.')"
708  character(len=*), parameter :: fmtsperror = &
709  &"(A,' for stress period ', I0)"
710  character(len=*), parameter :: fmtdterror = &
711  "('Time step length of ', G0, ' is too small in period ', I0, &
712  &' and time step ', I0)"
713  !
714  ! -- Initialize
715  tstart = dzero
716  !
717  ! -- Go through and check each stress period
718  do kper = 1, nper
719  !
720  ! -- Error if nstp less than or equal to zero
721  if (nstp(kper) <= 0) then
722  write (errmsg, fmtsperror) 'Number of time steps less than one ', kper
723  call store_error(errmsg)
724  return
725  end if
726  !
727  ! -- Warn if perlen is zero
728  if (perlen(kper) == dzero) then
729  write (iout, fmtpwarn) kper
730  return
731  end if
732  !
733  ! -- Error if tsmult is less than zero
734  if (tsmult(kper) <= dzero) then
735  write (errmsg, fmtsperror) 'TSMULT must be greater than 0.0 ', kper
736  call store_error(errmsg)
737  return
738  end if
739  !
740  ! -- Error if negative period length
741  if (perlen(kper) < dzero) then
742  write (errmsg, fmtsperror) 'PERLEN cannot be less than 0.0 ', kper
743  call store_error(errmsg)
744  return
745  end if
746  !
747  ! -- Go through all time step lengths and make sure they are valid
748  do kstp = 1, nstp(kper)
749  if (kstp == 1) then
750  dt = perlen(kper) / float(nstp(kper))
751  if (tsmult(kper) /= done) &
752  dt = perlen(kper) * (done - tsmult(kper)) / &
753  (done - tsmult(kper)**nstp(kper))
754  else
755  dt = dt * tsmult(kper)
756  end if
757  tend = tstart + dt
758  !
759  ! -- Error condition if tstart == tend
760  if (tstart == tend) then
761  write (errmsg, fmtdterror) dt, kper, kstp
762  call store_error(errmsg)
763  return
764  end if
765  end do
766  !
767  ! -- reset tstart = tend
768  tstart = tend
769  !
770  end do
771  ! -- Return
772  return
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
real(dp), parameter done
real constant 1
Definition: Constants.f90:75
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 596 of file tdis.f90.

597  ! -- modules
599  !
600  call mem_allocate(perlen, nper, 'PERLEN', 'TDIS')
601  call mem_allocate(nstp, nper, 'NSTP', 'TDIS')
602  call mem_allocate(tsmult, nper, 'TSMULT', 'TDIS')
603  !
604  ! -- Return
605  return
Here is the caller graph for this function:

◆ tdis_allocate_scalars()

subroutine tdismodule::tdis_allocate_scalars

Definition at line 541 of file tdis.f90.

542  ! -- modules
544  use constantsmodule, only: dzero
545  !
546  ! -- memory manager variables
547  call mem_allocate(nper, 'NPER', 'TDIS')
548  call mem_allocate(itmuni, 'ITMUNI', 'TDIS')
549  call mem_allocate(kper, 'KPER', 'TDIS')
550  call mem_allocate(kstp, 'KSTP', 'TDIS')
551  call mem_allocate(inats, 'INATS', 'TDIS')
552  call mem_allocate(readnewdata, 'READNEWDATA', 'TDIS')
553  call mem_allocate(endofperiod, 'ENDOFPERIOD', 'TDIS')
554  call mem_allocate(endofsimulation, 'ENDOFSIMULATION', 'TDIS')
555  call mem_allocate(delt, 'DELT', 'TDIS')
556  call mem_allocate(pertim, 'PERTIM', 'TDIS')
557  call mem_allocate(topertim, 'TOPERTIM', 'TDIS')
558  call mem_allocate(totim, 'TOTIM', 'TDIS')
559  call mem_allocate(totimc, 'TOTIMC', 'TDIS')
560  call mem_allocate(deltsav, 'DELTSAV', 'TDIS')
561  call mem_allocate(totimsav, 'TOTIMSAV', 'TDIS')
562  call mem_allocate(pertimsav, 'PERTIMSAV', 'TDIS')
563  call mem_allocate(totalsimtime, 'TOTALSIMTIME', 'TDIS')
564  !
565  ! -- strings
566  allocate (datetime0)
567  allocate (input_mempath)
568  allocate (input_fname)
569  !
570  ! -- Initialize variables
571  nper = 0
572  itmuni = 0
573  kper = 0
574  kstp = 0
575  inats = 0
576  readnewdata = .true.
577  endofperiod = .true.
578  endofsimulation = .false.
579  delt = dzero
580  pertim = dzero
581  topertim = dzero
582  totim = dzero
583  totimc = dzero
584  deltsav = dzero
585  totimsav = dzero
586  pertimsav = dzero
587  totalsimtime = dzero
588  datetime0 = ''
589  !
590  ! -- Return
591  return
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
86  !
87  ! -- Return
88  return
subroutine, public ats_cr(inunit, nper_tdis)
@ brief Create ATS object
Definition: ats.f90:62
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 422 of file tdis.f90.

423  ! -- modules
425  use adaptivetimestepmodule, only: ats_da
426  !
427  ! -- ats
428  if (inats > 0) call ats_da()
429  !
430  ! -- Scalars
431  call mem_deallocate(nper)
432  call mem_deallocate(itmuni)
433  call mem_deallocate(kper)
434  call mem_deallocate(kstp)
435  call mem_deallocate(inats)
436  call mem_deallocate(readnewdata)
437  call mem_deallocate(endofperiod)
438  call mem_deallocate(endofsimulation)
439  call mem_deallocate(delt)
440  call mem_deallocate(pertim)
441  call mem_deallocate(topertim)
442  call mem_deallocate(totim)
443  call mem_deallocate(totimc)
444  call mem_deallocate(deltsav)
445  call mem_deallocate(totimsav)
446  call mem_deallocate(pertimsav)
447  call mem_deallocate(totalsimtime)
448  !
449  ! -- strings
450  deallocate (datetime0)
451  deallocate (input_mempath)
452  deallocate (input_fname)
453  !
454  ! -- Arrays
455  call mem_deallocate(perlen)
456  call mem_deallocate(nstp)
457  call mem_deallocate(tsmult)
458  !
459  ! -- Return
460  return
subroutine, public ats_da()
@ brief Deallocate variables
Definition: ats.f90:177
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 221 of file tdis.f90.

222  ! -- modules
223  use constantsmodule, only: done, dzero
225  ats_set_delt, &
227  ! -- dummy
228  real(DP), intent(in) :: deltnew
229  ! -- local
230  logical(LGP) :: adaptivePeriod
231  !
232  ! -- Set values
233  adaptiveperiod = isadaptiveperiod(kper)
234  delt = deltnew
235  totim = totimsav + delt
236  pertim = pertimsav + delt
237  !
238  ! -- Set end of period indicator
239  endofperiod = .false.
240  if (adaptiveperiod) then
241  call ats_set_endofperiod(kper, pertim, perlen(kper), endofperiod)
242  else
243  if (kstp == nstp(kper)) then
244  endofperiod = .true.
245  end if
246  end if
247  !
248  ! -- Set end of simulation indicator
249  if (endofperiod .and. kper == nper) then
250  endofsimulation = .true.
251  totim = totalsimtime
252  end if
253  !
254  ! -- Return
255  return
subroutine, public ats_set_delt(kstp, kper, pertim, perlencurrent, delt)
@ brief Set time step
Definition: ats.f90:566
subroutine, public ats_set_endofperiod(kper, pertim, perlencurrent, endofperiod)
@ brief Set end of period indicator
Definition: ats.f90:674
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 348 of file tdis.f90.

349  ! -- modules
352  ! -- dummy
353  integer(I4B), intent(in) :: iout
354  ! -- local
355  real(DP) :: cnv, delsec, totsec, persec, delmn, delhr, totmn, tothr, &
356  totdy, totyr, permn, perhr, perdy, peryr, deldy, delyr
357  ! -- format
358  character(len=*), parameter :: fmttmsmry = "(1X, ///9X, &
359  &'TIME SUMMARY AT END OF TIME STEP', I5,' IN STRESS PERIOD ', I4)"
360  character(len=*), parameter :: fmttmstpmsg = &
361  &"(21X, ' TIME STEP LENGTH =', G15.6 / &
362  & 21X, ' STRESS PERIOD TIME =', G15.6 / &
363  & 21X, 'TOTAL SIMULATION TIME =', G15.6)"
364  character(len=*), parameter :: fmttottmmsg = &
365  &"(19X, ' SECONDS MINUTES HOURS', 7X, &
366  &'DAYS YEARS'/20X, 59('-'))"
367  character(len=*), parameter :: fmtdelttm = &
368  &"(1X, ' TIME STEP LENGTH', 1P, 5G12.5)"
369  character(len=*), parameter :: fmtpertm = &
370  &"(1X, 'STRESS PERIOD TIME', 1P, 5G12.5)"
371  character(len=*), parameter :: fmttottm = &
372  &"(1X, ' TOTAL TIME', 1P, 5G12.5,/)"
373  !
374  ! -- Write header message for the information that follows
375  write (iout, fmttmsmry) kstp, kper
376  !
377  ! -- Use time unit indicator to get factor to convert to seconds
378  cnv = dzero
379  if (itmuni == 1) cnv = done
380  if (itmuni == 2) cnv = dsixty
381  if (itmuni == 3) cnv = dsecperhr
382  if (itmuni == 4) cnv = dsecperdy
383  if (itmuni == 5) cnv = dsecperyr
384  !
385  ! -- If FACTOR=0 then time units are non-standard
386  if (cnv == dzero) then
387  ! -- Print times in non-standard time units
388  write (iout, fmttmstpmsg) delt, pertim, totim
389  else
390  ! -- Calculate length of time step & elapsed time in seconds
391  delsec = cnv * delt
392  totsec = cnv * totim
393  persec = cnv * pertim
394  !
395  ! -- Calculate times in minutes, hours, days, and years
396  delmn = delsec / dsixty
397  delhr = delmn / dsixty
398  deldy = delhr / dhrperday
399  delyr = deldy / ddyperyr
400  totmn = totsec / dsixty
401  tothr = totmn / dsixty
402  totdy = tothr / dhrperday
403  totyr = totdy / ddyperyr
404  permn = persec / dsixty
405  perhr = permn / dsixty
406  perdy = perhr / dhrperday
407  peryr = perdy / ddyperyr
408  !
409  ! -- Print time step length and elapsed times in all time units
410  write (iout, fmttottmmsg)
411  write (iout, fmtdelttm) delsec, delmn, delhr, deldy, delyr
412  write (iout, fmtpertm) persec, permn, perhr, perdy, peryr
413  write (iout, fmttottm) totsec, totmn, tothr, totdy, totyr
414  end if
415  !
416  ! -- Return
417  return
real(dp), parameter dsixty
real constant 60
Definition: Constants.f90:84
real(dp), parameter dsecperdy
real constant representing the number of seconds per day (used in tdis)
Definition: Constants.f90:99
real(dp), parameter dsecperyr
real constant representing the average number of seconds per year (used in tdis)
Definition: Constants.f90:100
real(dp), parameter dsecperhr
real constant representing number of seconds per hour (used in tdis)
Definition: Constants.f90:96
real(dp), parameter ddyperyr
real constant representing the average number of days per year (used in tdis)
Definition: Constants.f90:98
real(dp), parameter dhrperday
real constant representing number of hours per day (used in tdis)
Definition: Constants.f90:97
Here is the caller graph for this function:

◆ tdis_set_counters()

subroutine, public tdismodule::tdis_set_counters

Definition at line 93 of file tdis.f90.

94  ! -- modules
96  use simvariablesmodule, only: isim_mode
97  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) 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 (isadaptiveperiod(kper)) then
146  call ats_period_message(kper)
147  else
148  write (iout, fmtspits) nstp(kper), tsmult(kper)
149  end if
150  end if
151  !
152  ! -- Return
153  return
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:492
@ mvalidate
validation mode - do not run time steps
Definition: Constants.f90:204
@ mnormal
normal output mode
Definition: Constants.f90:205
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:94
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 260 of file tdis.f90.

261  ! -- modules
262  use constantsmodule, only: done
263  !
264  if (kstp == 1) then
265  ! -- Calculate the first value of delt for this stress period
266  topertim = totim
267  if (tsmult(kper) /= done) then
268  ! -- Timestep length has a geometric progression
269  delt = perlen(kper) * (done - tsmult(kper)) / &
270  (done - tsmult(kper)**nstp(kper))
271  else
272  ! -- Timestep length is constant
273  delt = perlen(kper) / float(nstp(kper))
274  end if
275  elseif (kstp == nstp(kper)) then
276  ! -- Calculate exact last delt to avoid accumulation errors
277  delt = topertim + perlen(kper) - totim
278  else
279  delt = tsmult(kper) * delt
280  end if
281  !
282  ! -- Return
283  return
Here is the caller graph for this function:

◆ tdis_set_timestep()

subroutine, public tdismodule::tdis_set_timestep

Definition at line 158 of file tdis.f90.

159  ! -- modules
160  use constantsmodule, only: done, dzero
162  ats_set_delt, &
164  ! -- local
165  logical(LGP) :: adaptivePeriod
166  ! -- format
167  character(len=*), parameter :: fmttsi = &
168  "(1X,'INITIAL TIME STEP SIZE =',G15.7)"
169  !
170  ! -- Initialize
171  adaptiveperiod = isadaptiveperiod(kper)
172  if (kstp == 1) then
173  pertim = dzero
174  topertim = dzero
175  end if
176  !
177  ! -- Set delt
178  if (adaptiveperiod) then
179  call ats_set_delt(kstp, kper, pertim, perlen(kper), delt)
180  else
181  call tdis_set_delt()
182  if (kstp == 1) then
183  write (iout, fmttsi) delt
184  end if
185  end if
186  !
187  ! -- Advance timers and update totim and pertim based on delt
188  totimsav = totim
189  pertimsav = pertim
190  totimc = totimsav
191  totim = totimsav + delt
192  pertim = pertimsav + delt
193  !
194  ! -- Set end of period indicator
195  endofperiod = .false.
196  if (adaptiveperiod) then
197  call ats_set_endofperiod(kper, pertim, perlen(kper), endofperiod)
198  else
199  if (kstp == nstp(kper)) then
200  endofperiod = .true.
201  end if
202  end if
203  if (endofperiod) then
204  pertim = perlen(kper)
205  end if
206  !
207  ! -- Set end of simulation indicator
208  if (endofperiod .and. kper == nper) then
209  endofsimulation = .true.
210  end if
211  !
212  ! -- Return
213  return
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 610 of file tdis.f90.

611  ! -- modules
612  use constantsmodule, only: linelength
616  ! -- local
617  type(SimTdisParamFoundType) :: found
618  ! -- formats
619  character(len=*), parameter :: fmtnper = &
620  "(1X,I4,' STRESS PERIOD(S) IN SIMULATION')"
621  !
622  ! -- source dimensions from input context
623  call mem_set_value(nper, 'NPER', input_mempath, found%nper)
624  !
625  ! -- log values to list file
626  write (iout, '(1x,a)') 'PROCESSING TDIS DIMENSIONS'
627  !
628  if (found%nper) then
629  write (iout, fmtnper) nper
630  end if
631  !
632  write (iout, '(1x,a)') 'END OF TDIS DIMENSIONS'
633  !
634  ! -- Return
635  return
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 465 of file tdis.f90.

466  ! -- modules
467  use constantsmodule, only: linelength
472  ! -- local
473  type(SimTdisParamFoundType) :: found
474  character(len=LINELENGTH), dimension(6) :: time_units = &
475  &[character(len=LINELENGTH) :: 'UNDEFINED', 'SECONDS', 'MINUTES', 'HOURS', &
476  'DAYS', 'YEARS']
477  character(len=LINELENGTH) :: fname
478  ! -- formats
479  character(len=*), parameter :: fmtitmuni = &
480  &"(4x,'SIMULATION TIME UNIT IS ',A)"
481  character(len=*), parameter :: fmtdatetime0 = &
482  &"(4x,'SIMULATION STARTING DATE AND TIME IS ',A)"
483  !
484  ! -- initialize time unit to undefined
485  itmuni = 0
486  !
487  ! -- source options from input context
488  call mem_set_value(itmuni, 'TIME_UNITS', input_mempath, time_units, &
489  found%time_units)
490  call mem_set_value(datetime0, 'START_DATE_TIME', input_mempath, &
491  found%start_date_time)
492  !
493  if (found%time_units) then
494  !
495  ! -- adjust to 0-based indexing for itmuni
496  itmuni = itmuni - 1
497  end if
498  !
499  ! -- enforce 0 or 1 ATS6_FILENAME entries in option block
500  if (filein_fname(fname, 'ATS6_FILENAME', input_mempath, &
501  input_fname)) then
502  inats = getunit()
503  call openfile(inats, iout, fname, 'ATS')
504  end if
505  !
506  ! -- log values to list file
507  write (iout, '(1x,a)') 'PROCESSING TDIS OPTIONS'
508  !
509  if (found%time_units) then
510  select case (itmuni)
511  case (0)
512  write (iout, fmtitmuni) 'UNDEFINED'
513  case (1)
514  write (iout, fmtitmuni) 'SECONDS'
515  case (2)
516  write (iout, fmtitmuni) 'MINUTES'
517  case (3)
518  write (iout, fmtitmuni) 'HOURS'
519  case (4)
520  write (iout, fmtitmuni) 'DAYS'
521  case (5)
522  write (iout, fmtitmuni) 'YEARS'
523  case default
524  end select
525  else
526  write (iout, fmtitmuni) 'UNDEFINED'
527  end if
528  !
529  if (found%start_date_time) then
530  write (iout, fmtdatetime0) datetime0
531  end if
532  !
533  write (iout, '(1x,a)') 'END OF TDIS OPTIONS'
534  !
535  ! -- Return
536  return
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 640 of file tdis.f90.

641  ! -- modules
642  use constantsmodule, only: linelength, dzero
647  ! -- local
648  type(SimTdisParamFoundType) :: found
649  integer(I4B) :: n
650  ! -- formats
651  character(len=*), parameter :: fmtheader = &
652  "(1X,//1X,'STRESS PERIOD LENGTH TIME STEPS', &
653  &' MULTIPLIER FOR DELT',/1X,76('-'))"
654  character(len=*), parameter :: fmtrow = &
655  "(1X,I8,1PG21.7,I7,0PF25.3)"
656  !
657  ! -- source perioddata from input context
658  call mem_set_value(perlen, 'PERLEN', input_mempath, found%perlen)
659  call mem_set_value(nstp, 'NSTP', input_mempath, found%nstp)
660  call mem_set_value(tsmult, 'TSMULT', input_mempath, found%tsmult)
661  !
662  ! -- Check timing information
663  call check_tdis_timing(nper, perlen, nstp, tsmult)
664  !
665  ! -- Check for errors
666  if (count_errors() > 0) then
667  call store_error_filename(input_fname)
668  end if
669  !
670  ! -- log timing
671  write (iout, '(1x,a)') 'PROCESSING TDIS PERIODDATA'
672  write (iout, fmtheader)
673  !
674  do n = 1, size(perlen)
675  write (iout, fmtrow) n, perlen(n), nstp(n), tsmult(n)
676  totalsimtime = totalsimtime + perlen(n)
677  end do
678  !
679  write (iout, '(1x,a)') 'END OF TDIS PERIODDATA'
680  !
681  ! -- Return
682  return
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 tdismodule::datetime0 => null()
private

Definition at line 41 of file tdis.f90.

41  character(len=LENDATETIME), 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