MODFLOW 6  version 6.7.0.dev3
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
10  !
11  implicit none
12  !
13  private
14  public :: tdis_cr
15  public :: tdis_set_counters
16  public :: tdis_set_timestep
17  public :: tdis_delt_reset
18  public :: tdis_ot
19  public :: tdis_da
20  !
21  integer(I4B), public, pointer :: nper => null() !< number of stress period
22  integer(I4B), public, pointer :: itmuni => null() !< flag indicating time units
23  integer(I4B), public, pointer :: kper => null() !< current stress period number
24  integer(I4B), public, pointer :: kstp => null() !< current time step number
25  integer(I4B), public, pointer :: inats => null() !< flag indicating ats active for simulation
26  logical(LGP), public, pointer :: readnewdata => null() !< flag indicating time to read new data
27  logical(LGP), public, pointer :: endofperiod => null() !< flag indicating end of stress period
28  logical(LGP), public, pointer :: endofsimulation => null() !< flag indicating end of simulation
29  real(dp), public, pointer :: delt => null() !< length of the current time step
30  real(dp), public, pointer :: pertim => null() !< time relative to start of stress period
31  real(dp), public, pointer :: topertim => null() !< simulation time at start of stress period
32  real(dp), public, pointer :: totim => null() !< time relative to start of simulation
33  real(dp), public, pointer :: totimc => null() !< simulation time at start of time step
34  real(dp), public, pointer :: deltsav => null() !< saved value for delt, used for subtiming
35  real(dp), public, pointer :: totimsav => null() !< saved value for totim, used for subtiming
36  real(dp), public, pointer :: pertimsav => null() !< saved value for pertim, used for subtiming
37  real(dp), public, pointer :: totalsimtime => null() !< time at end of simulation
38  real(dp), public, dimension(:), pointer, contiguous :: perlen => null() !< length of each stress period
39  integer(I4B), public, dimension(:), pointer, contiguous :: nstp => null() !< number of time steps in each stress period
40  real(dp), public, dimension(:), pointer, contiguous :: tsmult => null() !< time step multiplier for each stress period
41  character(len=LENDATETIME), public, pointer :: datetime0 => null() !< starting date and time for the simulation
42  character(len=LENMEMPATH), pointer :: input_mempath => null() !< input context mempath for tdis
43  character(len=LINELENGTH), pointer :: input_fname => null() !< input filename for tdis
44  !
45 contains
46 
47  !> @brief Create temporal discretization
48  !<
49  subroutine tdis_cr(fname, inmempath)
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
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
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  end subroutine tdis_cr
87 
88  !> @brief Set kstp and kper
89  !<
90  subroutine tdis_set_counters()
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
144  else
145  write (iout, fmtspits) nstp(kper), tsmult(kper)
146  end if
147  end if
148  end subroutine tdis_set_counters
149 
150  !> @brief Set time step length
151  !<
152  subroutine tdis_set_timestep()
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
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
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
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
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.
243  end if
244  end subroutine tdis_delt_reset
245 
246  !> @brief Set time step length
247  !<
248  subroutine tdis_set_delt()
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
269  end subroutine tdis_set_delt
270 
271  !> @brief Print simulation time
272  !<
273  subroutine tdis_ot(iout)
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
340  end subroutine tdis_ot
341 
342  !> @brief Deallocate memory
343  !<
344  subroutine tdis_da()
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)
361  call mem_deallocate(delt)
362  call mem_deallocate(pertim)
364  call mem_deallocate(totim)
365  call mem_deallocate(totimc)
366  call mem_deallocate(deltsav)
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)
380  end subroutine tdis_da
381 
382  !> @brief Source the timing discretization options
383  !<
384  subroutine tdis_source_options()
385  ! -- modules
386  use constantsmodule, only: linelength
392  ! -- local
393  type(simtdisparamfoundtype) :: found
394  character(len=LINELENGTH), dimension(6) :: time_units = &
395  &[character(len=LINELENGTH) :: 'UNDEFINED', 'SECONDS', 'MINUTES', 'HOURS', &
396  'DAYS', 'YEARS']
397  character(len=LINELENGTH) :: fname
398  ! -- formats
399  character(len=*), parameter :: fmtitmuni = &
400  &"(4x,'SIMULATION TIME UNIT IS ',A)"
401  character(len=*), parameter :: fmtdatetime0 = &
402  &"(4x,'SIMULATION STARTING DATE AND TIME IS ',A)"
403  !
404  ! -- initialize time unit to undefined
405  itmuni = 0
406  !
407  ! -- source options from input context
408  call mem_set_value(itmuni, 'TIME_UNITS', input_mempath, time_units, &
409  found%time_units)
410  call mem_set_value(datetime0, 'START_DATE_TIME', input_mempath, &
411  found%start_date_time)
412  !
413  if (found%time_units) then
414  if (itmuni == 0) then
415  call store_error('Unrecognized input value for TIME_UNITS option.')
417  else
418  !
419  ! -- adjust to 0-based indexing for itmuni
420  itmuni = itmuni - 1
421  end if
422  end if
423  !
424  ! -- enforce 0 or 1 ATS6_FILENAME entries in option block
425  if (filein_fname(fname, 'ATS6_FILENAME', input_mempath, &
426  input_fname)) then
427  inats = getunit()
428  call openfile(inats, iout, fname, 'ATS')
429  end if
430  !
431  ! -- log values to list file
432  write (iout, '(1x,a)') 'PROCESSING TDIS OPTIONS'
433  !
434  if (found%time_units) then
435  select case (itmuni)
436  case (0)
437  write (iout, fmtitmuni) 'UNDEFINED'
438  case (1)
439  write (iout, fmtitmuni) 'SECONDS'
440  case (2)
441  write (iout, fmtitmuni) 'MINUTES'
442  case (3)
443  write (iout, fmtitmuni) 'HOURS'
444  case (4)
445  write (iout, fmtitmuni) 'DAYS'
446  case (5)
447  write (iout, fmtitmuni) 'YEARS'
448  case default
449  end select
450  else
451  write (iout, fmtitmuni) 'UNDEFINED'
452  end if
453  !
454  if (found%start_date_time) then
455  write (iout, fmtdatetime0) datetime0
456  end if
457  !
458  write (iout, '(1x,a)') 'END OF TDIS OPTIONS'
459  end subroutine tdis_source_options
460 
461  !> @brief Allocate tdis scalars
462  !<
464  ! -- modules
466  use constantsmodule, only: dzero
467  !
468  ! -- memory manager variables
469  call mem_allocate(nper, 'NPER', 'TDIS')
470  call mem_allocate(itmuni, 'ITMUNI', 'TDIS')
471  call mem_allocate(kper, 'KPER', 'TDIS')
472  call mem_allocate(kstp, 'KSTP', 'TDIS')
473  call mem_allocate(inats, 'INATS', 'TDIS')
474  call mem_allocate(readnewdata, 'READNEWDATA', 'TDIS')
475  call mem_allocate(endofperiod, 'ENDOFPERIOD', 'TDIS')
476  call mem_allocate(endofsimulation, 'ENDOFSIMULATION', 'TDIS')
477  call mem_allocate(delt, 'DELT', 'TDIS')
478  call mem_allocate(pertim, 'PERTIM', 'TDIS')
479  call mem_allocate(topertim, 'TOPERTIM', 'TDIS')
480  call mem_allocate(totim, 'TOTIM', 'TDIS')
481  call mem_allocate(totimc, 'TOTIMC', 'TDIS')
482  call mem_allocate(deltsav, 'DELTSAV', 'TDIS')
483  call mem_allocate(totimsav, 'TOTIMSAV', 'TDIS')
484  call mem_allocate(pertimsav, 'PERTIMSAV', 'TDIS')
485  call mem_allocate(totalsimtime, 'TOTALSIMTIME', 'TDIS')
486  !
487  ! -- strings
488  allocate (datetime0)
489  allocate (input_mempath)
490  allocate (input_fname)
491  !
492  ! -- Initialize variables
493  nper = 0
494  itmuni = 0
495  kper = 0
496  kstp = 0
497  inats = 0
498  readnewdata = .true.
499  endofperiod = .true.
500  endofsimulation = .false.
501  delt = dzero
502  pertim = dzero
503  topertim = dzero
504  totim = dzero
505  totimc = dzero
506  deltsav = dzero
507  totimsav = dzero
508  pertimsav = dzero
510  datetime0 = ''
511  end subroutine tdis_allocate_scalars
512 
513  !> @brief Allocate tdis arrays
514  !<
515  subroutine tdis_allocate_arrays()
516  ! -- modules
518  !
519  call mem_allocate(perlen, nper, 'PERLEN', 'TDIS')
520  call mem_allocate(nstp, nper, 'NSTP', 'TDIS')
521  call mem_allocate(tsmult, nper, 'TSMULT', 'TDIS')
522  end subroutine tdis_allocate_arrays
523 
524  !> @brief Source dimension NPER
525  !<
527  ! -- modules
528  use constantsmodule, only: linelength
532  ! -- local
533  type(simtdisparamfoundtype) :: found
534  ! -- formats
535  character(len=*), parameter :: fmtnper = &
536  "(1X,I4,' STRESS PERIOD(S) IN SIMULATION')"
537  !
538  ! -- source dimensions from input context
539  call mem_set_value(nper, 'NPER', input_mempath, found%nper)
540  !
541  ! -- log values to list file
542  write (iout, '(1x,a)') 'PROCESSING TDIS DIMENSIONS'
543  !
544  if (found%nper) then
545  write (iout, fmtnper) nper
546  end if
547  !
548  write (iout, '(1x,a)') 'END OF TDIS DIMENSIONS'
549  end subroutine tdis_source_dimensions
550 
551  !> @brief Source timing information
552  !<
553  subroutine tdis_source_timing()
554  ! -- modules
555  use constantsmodule, only: linelength, dzero
560  ! -- local
561  type(simtdisparamfoundtype) :: found
562  integer(I4B) :: n
563  ! -- formats
564  character(len=*), parameter :: fmtheader = &
565  "(1X,//1X,'STRESS PERIOD LENGTH TIME STEPS', &
566  &' MULTIPLIER FOR DELT',/1X,76('-'))"
567  character(len=*), parameter :: fmtrow = &
568  "(1X,I8,1PG21.7,I7,0PF25.3)"
569  !
570  ! -- source perioddata from input context
571  call mem_set_value(perlen, 'PERLEN', input_mempath, found%perlen)
572  call mem_set_value(nstp, 'NSTP', input_mempath, found%nstp)
573  call mem_set_value(tsmult, 'TSMULT', input_mempath, found%tsmult)
574  !
575  ! -- Check timing information
577  !
578  ! -- Check for errors
579  if (count_errors() > 0) then
581  end if
582  !
583  ! -- log timing
584  write (iout, '(1x,a)') 'PROCESSING TDIS PERIODDATA'
585  write (iout, fmtheader)
586  !
587  do n = 1, size(perlen)
588  write (iout, fmtrow) n, perlen(n), nstp(n), tsmult(n)
590  end do
591  !
592  write (iout, '(1x,a)') 'END OF TDIS PERIODDATA'
593  end subroutine tdis_source_timing
594 
595  !> @brief Check the tdis timing information
596  !!
597  !! Return back to tdis_read_timing if an error condition is found and let the
598  !! ustop routine be called there instead so the StoreErrorUnit routine can be
599  !! called to assign the correct file name.
600  !<
601  subroutine check_tdis_timing(nper, perlen, nstp, tsmult)
602  ! -- modules
603  use constantsmodule, only: linelength, dzero, done
604  use simmodule, only: store_error
605  ! -- dummy
606  integer(I4B), intent(in) :: nper
607  real(DP), dimension(:), contiguous, intent(in) :: perlen
608  integer(I4B), dimension(:), contiguous, intent(in) :: nstp
609  real(DP), dimension(:), contiguous, intent(in) :: tsmult
610  ! -- local
611  integer(I4B) :: kper, kstp
612  real(DP) :: tstart, tend, dt
613  character(len=LINELENGTH) :: errmsg
614  ! -- formats
615  character(len=*), parameter :: fmtpwarn = &
616  "(1X,/1X,'PERLEN is zero for stress period ', I0, &
617  &'. PERLEN must not be zero for transient periods.')"
618  character(len=*), parameter :: fmtsperror = &
619  &"(A,' for stress period ', I0)"
620  character(len=*), parameter :: fmtdterror = &
621  "('Time step length of ', G0, ' is too small in period ', I0, &
622  &' and time step ', I0)"
623  !
624  ! -- Initialize
625  tstart = dzero
626  !
627  ! -- Go through and check each stress period
628  do kper = 1, nper
629  !
630  ! -- Error if nstp less than or equal to zero
631  if (nstp(kper) <= 0) then
632  write (errmsg, fmtsperror) 'Number of time steps less than one ', kper
633  call store_error(errmsg)
634  return
635  end if
636  !
637  ! -- Warn if perlen is zero
638  if (perlen(kper) == dzero) then
639  write (iout, fmtpwarn) kper
640  return
641  end if
642  !
643  ! -- Error if tsmult is less than zero
644  if (tsmult(kper) <= dzero) then
645  write (errmsg, fmtsperror) 'TSMULT must be greater than 0.0 ', kper
646  call store_error(errmsg)
647  return
648  end if
649  !
650  ! -- Error if negative period length
651  if (perlen(kper) < dzero) then
652  write (errmsg, fmtsperror) 'PERLEN cannot be less than 0.0 ', kper
653  call store_error(errmsg)
654  return
655  end if
656  !
657  ! -- Go through all time step lengths and make sure they are valid
658  do kstp = 1, nstp(kper)
659  if (kstp == 1) then
660  dt = perlen(kper) / float(nstp(kper))
661  if (tsmult(kper) /= done) &
662  dt = perlen(kper) * (done - tsmult(kper)) / &
663  (done - tsmult(kper)**nstp(kper))
664  else
665  dt = dt * tsmult(kper)
666  end if
667  tend = tstart + dt
668  !
669  ! -- Error condition if tstart == tend
670  if (tstart == tend) then
671  write (errmsg, fmtdterror) dt, kper, kstp
672  call store_error(errmsg)
673  return
674  end if
675  end do
676  !
677  ! -- reset tstart = tend
678  tstart = tend
679  !
680  end do
681  end subroutine check_tdis_timing
682 
683 end module tdismodule
684 
subroutine, public ats_set_delt(kstp, kper, pertim, perlencurrent, delt)
@ brief Set time step
Definition: ats.f90:541
subroutine, public ats_cr(inunit, nper_tdis)
@ brief Create ATS object
Definition: ats.f90:61
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
real(dp), pointer, public dtstable
delt value required for stability
Definition: ats.f90:26
subroutine, public ats_da()
@ brief Deallocate variables
Definition: ats.f90:167
subroutine, public ats_period_message(kper)
@ brief Write period message
Definition: ats.f90:469
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:40
real(dp), pointer, public pertim
time relative to start of stress period
Definition: tdis.f90:30
character(len=lenmempath), pointer input_mempath
input context mempath for tdis
Definition: tdis.f90:42
logical(lgp), pointer, public endofperiod
flag indicating end of stress period
Definition: tdis.f90:27
logical(lgp), pointer, public endofsimulation
flag indicating end of simulation
Definition: tdis.f90:28
subroutine, public tdis_ot(iout)
Print simulation time.
Definition: tdis.f90:274
integer(i4b), pointer, public itmuni
flag indicating time units
Definition: tdis.f90:22
integer(i4b), dimension(:), pointer, public, contiguous nstp
number of time steps in each stress period
Definition: tdis.f90:39
real(dp), dimension(:), pointer, public, contiguous perlen
length of each stress period
Definition: tdis.f90:38
subroutine tdis_source_dimensions()
Source dimension NPER.
Definition: tdis.f90:527
subroutine tdis_source_timing()
Source timing information.
Definition: tdis.f90:554
real(dp), pointer, public totim
time relative to start of simulation
Definition: tdis.f90:32
character(len=linelength), pointer input_fname
input filename for tdis
Definition: tdis.f90:43
character(len=lendatetime), pointer, public datetime0
starting date and time for the simulation
Definition: tdis.f90:41
logical(lgp), pointer, public readnewdata
flag indicating time to read new data
Definition: tdis.f90:26
real(dp), pointer, public pertimsav
saved value for pertim, used for subtiming
Definition: tdis.f90:36
subroutine, public tdis_cr(fname, inmempath)
Create temporal discretization.
Definition: tdis.f90:50
real(dp), pointer, public topertim
simulation time at start of stress period
Definition: tdis.f90:31
subroutine tdis_allocate_arrays()
Allocate tdis arrays.
Definition: tdis.f90:516
integer(i4b), pointer, public inats
flag indicating ats active for simulation
Definition: tdis.f90:25
subroutine, public tdis_set_timestep()
Set time step length.
Definition: tdis.f90:153
subroutine, public tdis_set_counters()
Set kstp and kper.
Definition: tdis.f90:91
real(dp), pointer, public totimc
simulation time at start of time step
Definition: tdis.f90:33
real(dp), pointer, public totalsimtime
time at end of simulation
Definition: tdis.f90:37
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24
real(dp), pointer, public totimsav
saved value for totim, used for subtiming
Definition: tdis.f90:35
subroutine, public tdis_da()
Deallocate memory.
Definition: tdis.f90:345
subroutine tdis_set_delt()
Set time step length.
Definition: tdis.f90:249
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:385
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23
real(dp), pointer, public delt
length of the current time step
Definition: tdis.f90:29
subroutine check_tdis_timing(nper, perlen, nstp, tsmult)
Check the tdis timing information.
Definition: tdis.f90:602
subroutine tdis_allocate_scalars()
Allocate tdis scalars.
Definition: tdis.f90:464
real(dp), pointer, public deltsav
saved value for delt, used for subtiming
Definition: tdis.f90:34
integer(i4b), pointer, public nper
number of stress period
Definition: tdis.f90:21