MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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), 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  !
87  ! -- Return
88  return
89  end subroutine tdis_cr
90 
91  !> @brief Set kstp and kper
92  !<
93  subroutine tdis_set_counters()
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
147  else
148  write (iout, fmtspits) nstp(kper), tsmult(kper)
149  end if
150  end if
151  !
152  ! -- Return
153  return
154  end subroutine tdis_set_counters
155 
156  !> @brief Set time step length
157  !<
158  subroutine tdis_set_timestep()
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
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
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
214  end subroutine tdis_set_timestep
215 
216  !> @brief Reset delt and update timing variables and indicators
217  !!
218  !! This routine is called when a timestep fails to converge, and so it is
219  !! retried using a smaller time step (deltnew).
220  !<
221  subroutine tdis_delt_reset(deltnew)
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
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.
252  end if
253  !
254  ! -- Return
255  return
256  end subroutine tdis_delt_reset
257 
258  !> @brief Set time step length
259  !<
260  subroutine tdis_set_delt()
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
284  end subroutine tdis_set_delt
285 
286 ! subroutine tdis_set_delt_std()
287 !! ******************************************************************************
288 !! tdis_tu_std -- Standard non-adaptive time update
289 !! ******************************************************************************
290 !!
291 !! SPECIFICATIONS:
292 !! ------------------------------------------------------------------------------
293 ! ! -- modules
294 ! use ConstantsModule, only: DONE, DZERO
295 ! ! -- formats
296 ! character(len=*),parameter :: fmttsi = &
297 ! "(28X,'INITIAL TIME STEP SIZE =',G15.7)"
298 !! ------------------------------------------------------------------------------
299 ! !
300 ! ! -- Setup new stress period if kstp is 1
301 ! if(kstp == 1) then
302 ! !
303 ! ! -- Calculate the first value of delt for this stress period
304 ! delt = perlen(kper) / float(nstp(kper))
305 ! if(tsmult(kper) /= DONE) &
306 ! delt = perlen(kper) * (DONE-tsmult(kper)) / &
307 ! (DONE - tsmult(kper) ** nstp(kper))
308 ! !
309 ! ! -- Print length of first time step
310 ! write(iout, fmttsi) delt
311 ! !
312 ! ! -- Initialize pertim (Elapsed time within stress period)
313 ! pertim = DZERO
314 ! !
315 ! ! -- Clear flag that indicates last time step of a stress period
316 ! endofperiod = .false.
317 ! endif
318 ! !
319 ! ! -- Calculate delt for kstp > 1
320 ! if (kstp /= 1) then
321 ! delt = tsmult(kper) * delt
322 ! end if
323 ! !
324 ! ! -- Store totim and pertim, which are times at end of previous time step
325 ! totimsav = totim
326 ! pertimsav = pertim
327 ! totimc = totim
328 ! !
329 ! ! -- Update totim and pertim
330 ! totim = totimsav + delt
331 ! pertim = pertimsav + delt
332 ! !
333 ! ! -- End of stress period and/or simulation?
334 ! if (kstp == nstp(kper)) then
335 ! endofperiod = .true.
336 ! end if
337 ! if (endofperiod .and. kper==nper) then
338 ! endofsimulation = .true.
339 ! totim = totalsimtime
340 ! end if
341 ! !
342 ! ! -- Return
343 ! return
344 ! end subroutine tdis_set_delt_std
345 
346  !> @brief Print simulation time
347  !<
348  subroutine tdis_ot(iout)
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
418  end subroutine tdis_ot
419 
420  !> @brief Deallocate memory
421  !<
422  subroutine tdis_da()
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)
439  call mem_deallocate(delt)
440  call mem_deallocate(pertim)
442  call mem_deallocate(totim)
443  call mem_deallocate(totimc)
444  call mem_deallocate(deltsav)
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
461  end subroutine tdis_da
462 
463  !> @brief Source the timing discretization options
464  !<
465  subroutine tdis_source_options()
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
537  end subroutine tdis_source_options
538 
539  !> @brief Allocate tdis scalars
540  !<
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
588  datetime0 = ''
589  !
590  ! -- Return
591  return
592  end subroutine tdis_allocate_scalars
593 
594  !> @brief Allocate tdis arrays
595  !<
596  subroutine tdis_allocate_arrays()
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
606  end subroutine tdis_allocate_arrays
607 
608  !> @brief Source dimension NPER
609  !<
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
636  end subroutine tdis_source_dimensions
637 
638  !> @brief Source timing information
639  !<
640  subroutine tdis_source_timing()
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
664  !
665  ! -- Check for errors
666  if (count_errors() > 0) then
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)
677  end do
678  !
679  write (iout, '(1x,a)') 'END OF TDIS PERIODDATA'
680  !
681  ! -- Return
682  return
683  end subroutine tdis_source_timing
684 
685  !> @brief Check the tdis timing information
686  !!
687  !! Return back to tdis_read_timing if an error condition is found and let the
688  !! ustop routine be called there instead so the StoreErrorUnit routine can be
689  !! called to assign the correct file name.
690  !<
691  subroutine check_tdis_timing(nper, perlen, nstp, tsmult)
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
773  end subroutine check_tdis_timing
774 
775 end module tdismodule
776 
subroutine, public ats_set_delt(kstp, kper, pertim, perlencurrent, delt)
@ brief Set time step
Definition: ats.f90:566
subroutine, public ats_cr(inunit, nper_tdis)
@ brief Create ATS object
Definition: ats.f90:62
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
real(dp), pointer, public dtstable
delt value required for stability
Definition: ats.f90:26
subroutine, public ats_da()
@ brief Deallocate variables
Definition: ats.f90:177
subroutine, public ats_period_message(kper)
@ brief Write period message
Definition: ats.f90:492
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
@ mvalidate
validation mode - do not run time steps
Definition: Constants.f90:204
@ mnormal
normal output mode
Definition: Constants.f90:205
real(dp), parameter dsixty
real constant 60
Definition: Constants.f90:84
real(dp), parameter dnodata
real no data constant
Definition: Constants.f90:94
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
integer(i4b), parameter lendatetime
maximum length of a date time string
Definition: Constants.f90:43
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 dzero
real constant zero
Definition: Constants.f90:64
@ vall
write all simulation notes and warnings
Definition: Constants.f90:188
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:26
real(dp), parameter dhrperday
real constant representing number of hours per day (used in tdis)
Definition: Constants.f90:97
real(dp), parameter done
real constant 1
Definition: Constants.f90:75
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:349
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:611
subroutine tdis_source_timing()
Source timing information.
Definition: tdis.f90:641
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
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:597
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:159
subroutine, public tdis_set_counters()
Set kstp and kper.
Definition: tdis.f90:94
character(len=lendatetime), pointer datetime0
starting date and time for the simulation
Definition: tdis.f90:41
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:423
subroutine tdis_set_delt()
Set time step length.
Definition: tdis.f90:261
subroutine, public tdis_delt_reset(deltnew)
Reset delt and update timing variables and indicators.
Definition: tdis.f90:222
subroutine tdis_source_options()
Source the timing discretization options.
Definition: tdis.f90:466
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:692
subroutine tdis_allocate_scalars()
Allocate tdis scalars.
Definition: tdis.f90:542
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