MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
tdis.f90
Go to the documentation of this file.
1 !stress periods and time stepping is handled by these routines
2 !convert this to a derived type? May not be necessary since only
3 !one of them is needed.
4 
5 module tdismodule
6 
7  use kindmodule, only: dp, i4b, lgp
11  use atsfactorymodule, only: create_ats
12  !
13  implicit none
14  !
15  private
16  public :: tdis_cr
17  public :: tdis_set_counters
18  public :: tdis_set_timestep
19  public :: tdis_delt_reset
20  public :: tdis_ot
21  public :: tdis_da
22  public :: ats
23  !
24  integer(I4B), public, pointer :: nper => null() !< number of stress period
25  integer(I4B), public, pointer :: itmuni => null() !< flag indicating time units
26  integer(I4B), public, pointer :: kper => null() !< current stress period number
27  integer(I4B), public, pointer :: kstp => null() !< current time step number
28  integer(I4B), public, pointer :: inats => null() !< flag indicating ats active for simulation
29  logical(LGP), public, pointer :: readnewdata => null() !< flag indicating time to read new data
30  logical(LGP), public, pointer :: endofperiod => null() !< flag indicating end of stress period
31  logical(LGP), public, pointer :: endofsimulation => null() !< flag indicating end of simulation
32  real(dp), public, pointer :: delt => null() !< length of the current time step
33  real(dp), public, pointer :: pertim => null() !< time relative to start of stress period
34  real(dp), public, pointer :: topertim => null() !< simulation time at start of stress period
35  real(dp), public, pointer :: totim => null() !< time relative to start of simulation
36  real(dp), public, pointer :: totimc => null() !< simulation time at start of time step
37  real(dp), public, pointer :: deltsav => null() !< saved value for delt, used for subtiming
38  real(dp), public, pointer :: totimsav => null() !< saved value for totim, used for subtiming
39  real(dp), public, pointer :: pertimsav => null() !< saved value for pertim, used for subtiming
40  real(dp), public, pointer :: totalsimtime => null() !< time at end of simulation
41  real(dp), public, dimension(:), pointer, contiguous :: perlen => null() !< length of each stress period
42  integer(I4B), public, dimension(:), pointer, contiguous :: nstp => null() !< number of time steps in each stress period
43  real(dp), public, dimension(:), pointer, contiguous :: tsmult => null() !< time step multiplier for each stress period
44  character(len=LENDATETIME), public, pointer :: datetime0 => null() !< starting date and time for the simulation
45  character(len=LENMEMPATH), pointer :: input_mempath => null() !< input context mempath for tdis
46  character(len=LINELENGTH), pointer :: input_fname => null() !< input filename for tdis
47 
48  class(atstype), pointer :: ats => null()
49  !
50 contains
51 
52  !> @brief Create temporal discretization
53  !<
54  subroutine tdis_cr(fname, inmempath)
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
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
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
91  end subroutine tdis_cr
92 
93  !> @brief Set kstp and kper
94  !<
95  subroutine tdis_set_counters()
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
151  end subroutine tdis_set_counters
152 
153  !> @brief Set time step length
154  !<
155  subroutine tdis_set_timestep()
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
205  end subroutine tdis_set_timestep
206 
207  !> @brief Reset delt and update timing variables and indicators
208  !!
209  !! This routine is called when a timestep fails to converge, and so it is
210  !! retried using a smaller time step (deltnew).
211  !<
212  subroutine tdis_delt_reset(deltnew)
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.
240  end if
241  end subroutine tdis_delt_reset
242 
243  !> @brief Set time step length
244  !<
245  subroutine tdis_set_delt()
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
266  end subroutine tdis_set_delt
267 
268  !> @brief Print simulation time
269  !<
270  subroutine tdis_ot(iout)
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
337  end subroutine tdis_ot
338 
339  !> @brief Deallocate memory
340  !<
341  subroutine tdis_da()
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)
360  call mem_deallocate(delt)
361  call mem_deallocate(pertim)
363  call mem_deallocate(totim)
364  call mem_deallocate(totimc)
365  call mem_deallocate(deltsav)
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)
379  end subroutine tdis_da
380 
381  !> @brief Source the timing discretization options
382  !<
383  subroutine tdis_source_options()
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.')
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'
458  end subroutine tdis_source_options
459 
460  !> @brief Allocate tdis scalars
461  !<
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
509  datetime0 = ''
510  end subroutine tdis_allocate_scalars
511 
512  !> @brief Allocate tdis arrays
513  !<
514  subroutine tdis_allocate_arrays()
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')
521  end subroutine tdis_allocate_arrays
522 
523  !> @brief Source dimension NPER
524  !<
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'
548  end subroutine tdis_source_dimensions
549 
550  !> @brief Source timing information
551  !<
552  subroutine tdis_source_timing()
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
576  !
577  ! -- Check for errors
578  if (count_errors() > 0) then
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)
589  end do
590  !
591  write (iout, '(1x,a)') 'END OF TDIS PERIODDATA'
592  end subroutine tdis_source_timing
593 
594  !> @brief Check the tdis timing information
595  !!
596  !! Return back to tdis_read_timing if an error condition is found and let the
597  !! ustop routine be called there instead so the StoreErrorUnit routine can be
598  !! called to assign the correct file name.
599  !<
600  subroutine check_tdis_timing(nper, perlen, nstp, tsmult)
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
680  end subroutine check_tdis_timing
681 
682 end module tdismodule
683 
class(atstype) function, pointer, public create_ats()
Factory method to create ATS object, to hide.
Definition: AtsFactory.F90:17
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:45
@ mvalidate
validation mode - do not run time steps
Definition: Constants.f90:205
@ mnormal
normal output mode
Definition: Constants.f90:206
real(dp), parameter dsixty
real constant 60
Definition: Constants.f90:85
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:95
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
integer(i4b), parameter lendatetime
maximum length of a date time string
Definition: Constants.f90:44
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 dzero
real constant zero
Definition: Constants.f90:65
@ vall
write all simulation notes and warnings
Definition: Constants.f90:189
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27
real(dp), parameter dhrperday
real constant representing number of hours per day (used in tdis)
Definition: Constants.f90:98
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
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
This module defines variable data types.
Definition: kind.f90:8
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 methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
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
This module contains simulation variables.
Definition: SimVariables.f90:9
integer(i4b) isim_level
simulation output level
integer(i4b) iout
file unit number for simulation output
integer(i4b) isim_mode
simulation mode
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
real(dp), dimension(:), pointer, public, contiguous tsmult
time step multiplier for each stress period
Definition: tdis.f90:43
real(dp), pointer, public pertim
time relative to start of stress period
Definition: tdis.f90:33
character(len=lenmempath), pointer input_mempath
input context mempath for tdis
Definition: tdis.f90:45
logical(lgp), pointer, public endofperiod
flag indicating end of stress period
Definition: tdis.f90:30
logical(lgp), pointer, public endofsimulation
flag indicating end of simulation
Definition: tdis.f90:31
subroutine, public tdis_ot(iout)
Print simulation time.
Definition: tdis.f90:271
integer(i4b), pointer, public itmuni
flag indicating time units
Definition: tdis.f90:25
integer(i4b), dimension(:), pointer, public, contiguous nstp
number of time steps in each stress period
Definition: tdis.f90:42
real(dp), dimension(:), pointer, public, contiguous perlen
length of each stress period
Definition: tdis.f90:41
subroutine tdis_source_dimensions()
Source dimension NPER.
Definition: tdis.f90:526
subroutine tdis_source_timing()
Source timing information.
Definition: tdis.f90:553
real(dp), pointer, public totim
time relative to start of simulation
Definition: tdis.f90:35
class(atstype), pointer, public ats
Definition: tdis.f90:48
character(len=linelength), pointer input_fname
input filename for tdis
Definition: tdis.f90:46
character(len=lendatetime), pointer, public datetime0
starting date and time for the simulation
Definition: tdis.f90:44
logical(lgp), pointer, public readnewdata
flag indicating time to read new data
Definition: tdis.f90:29
real(dp), pointer, public pertimsav
saved value for pertim, used for subtiming
Definition: tdis.f90:39
subroutine, public tdis_cr(fname, inmempath)
Create temporal discretization.
Definition: tdis.f90:55
real(dp), pointer, public topertim
simulation time at start of stress period
Definition: tdis.f90:34
subroutine tdis_allocate_arrays()
Allocate tdis arrays.
Definition: tdis.f90:515
integer(i4b), pointer, public inats
flag indicating ats active for simulation
Definition: tdis.f90:28
subroutine, public tdis_set_timestep()
Set time step length.
Definition: tdis.f90:156
subroutine, public tdis_set_counters()
Set kstp and kper.
Definition: tdis.f90:96
real(dp), pointer, public totimc
simulation time at start of time step
Definition: tdis.f90:36
real(dp), pointer, public totalsimtime
time at end of simulation
Definition: tdis.f90:40
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:27
real(dp), pointer, public totimsav
saved value for totim, used for subtiming
Definition: tdis.f90:38
subroutine, public tdis_da()
Deallocate memory.
Definition: tdis.f90:342
subroutine tdis_set_delt()
Set time step length.
Definition: tdis.f90:246
subroutine, public tdis_delt_reset(deltnew)
Reset delt and update timing variables and indicators.
Definition: tdis.f90:213
subroutine tdis_source_options()
Source the timing discretization options.
Definition: tdis.f90:384
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:26
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:32
subroutine check_tdis_timing(nper, perlen, nstp, tsmult)
Check the tdis timing information.
Definition: tdis.f90:601
subroutine tdis_allocate_scalars()
Allocate tdis scalars.
Definition: tdis.f90:463
real(dp), pointer, public deltsav
saved value for delt, used for subtiming
Definition: tdis.f90:37
integer(i4b), pointer, public nper
number of stress period
Definition: tdis.f90:24