MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
timeseriesmanagermodule Module Reference

Data Types

type  timeseriesmanagertype
 

Functions/Subroutines

subroutine, public tsmanager_cr (this, iout, removeTsLinksOnCompletion, extendTsToEndOfSimulation)
 Create the tsmanager. More...
 
subroutine tsmanager_df (this)
 Define time series manager object. More...
 
subroutine add_tsfile (this, fname, inunit)
 Add a time series file to this manager. More...
 
subroutine tsmgr_ad (this)
 Time step (or subtime step) advance. Call this each time step or subtime step. More...
 
subroutine tsmgr_da (this)
 Deallocate memory. More...
 
subroutine reset (this, pkgName)
 Call this when a new BEGIN PERIOD block is read for a new stress period. More...
 
subroutine make_link (this, timeSeries, pkgName, auxOrBnd, bndElem, irow, jcol, iprpak, tsLink, text, bndName)
 Make link. More...
 
type(timeserieslinktype) function, pointer getlink (this, auxOrBnd, indx)
 Get link. More...
 
integer(i4b) function countlinks (this, auxOrBnd)
 Count links. More...
 
type(timeseriestype) function, pointer get_time_series (this, name)
 Get time series. More...
 
subroutine hashbndtimeseries (this)
 Store all boundary (stress) time series links in TsContainers and construct hash table BndTsHashTable. More...
 
subroutine, public read_value_or_time_series (textInput, ii, jj, bndElem, pkgName, auxOrBnd, tsManager, iprpak, tsLink)
 Call this subroutine if the time-series link is available or needed. More...
 
subroutine, public read_value_or_time_series_adv (textInput, ii, jj, bndElem, pkgName, auxOrBnd, tsManager, iprpak, varName)
 Call this subroutine from advanced packages to define timeseries link for a variable (varName). More...
 
logical function remove_existing_link (tsManager, ii, jj, pkgName, auxOrBnd, varName)
 Remove an existing timeseries link if it is defined. More...
 
logical function, public var_timeseries (tsManager, pkgName, varName, auxOrBnd)
 Determine if a timeseries link with varName is defined. More...
 

Function/Subroutine Documentation

◆ add_tsfile()

subroutine timeseriesmanagermodule::add_tsfile ( class(timeseriesmanagertype this,
character(len=*), intent(in)  fname,
integer(i4b), intent(in)  inunit 
)
private

Definition at line 96 of file TimeSeriesManager.f90.

97  ! -- modules
100  ! -- dummy
101  class(TimeSeriesManagerType) :: this
102  character(len=*), intent(in) :: fname
103  integer(I4B), intent(in) :: inunit
104  ! -- local
105  integer(I4B) :: isize
106  integer(I4B) :: i
107  class(TimeSeriesFileType), pointer :: tsfile => null()
108  !
109  ! -- Check for fname duplicates
110  if (this%numtsfiles > 0) then
111  do i = 1, this%numtsfiles
112  if (this%tsfiles(i) == fname) then
113  call store_error('Found duplicate time-series file name: '//trim(fname))
114  call store_error_unit(inunit)
115  end if
116  end do
117  end if
118  !
119  ! -- Save fname
120  this%numtsfiles = this%numtsfiles + 1
121  isize = size(this%tsfiles)
122  if (this%numtsfiles > isize) then
123  call expandarray(this%tsfiles, 1000)
124  end if
125  this%tsfiles(this%numtsfiles) = fname
126  !
127  ! --
128  call this%tsfileList%Add(fname, this%iout, tsfile)
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
subroutine, public store_error_unit(iunit, terminate)
Store the file unit number.
Definition: Sim.f90:168
Here is the call graph for this function:

◆ countlinks()

integer(i4b) function timeseriesmanagermodule::countlinks ( class(timeseriesmanagertype this,
character(len=3), intent(in)  auxOrBnd 
)
private

Definition at line 437 of file TimeSeriesManager.f90.

438  ! -- return
439  integer(I4B) :: CountLinks
440  ! -- dummy
441  class(TimeSeriesManagerType) :: this
442  character(len=3), intent(in) :: auxOrBnd
443  !
444  countlinks = 0
445  if (auxorbnd == 'BND') then
446  countlinks = this%boundTsLinks%Count()
447  elseif (auxorbnd == 'AUX') then
448  countlinks = this%auxvarTsLinks%count()
449  end if

◆ get_time_series()

type(timeseriestype) function, pointer timeseriesmanagermodule::get_time_series ( class(timeseriesmanagertype this,
character(len=*), intent(in)  name 
)
private

Definition at line 454 of file TimeSeriesManager.f90.

455  ! -- dummy
456  class(TimeSeriesManagerType) :: this
457  character(len=*), intent(in) :: name
458  ! -- result
459  type(TimeSeriesType), pointer :: res
460  ! -- local
461  integer(I4B) :: indx
462  !
463  ! Get index from hash table, get time series from TsContainers,
464  ! and assign result to time series contained in link.
465  res => null()
466  indx = this%BndTsHashTable%get(name)
467  if (indx > 0) then
468  res => this%TsContainers(indx)%timeSeries
469  end if

◆ getlink()

type(timeserieslinktype) function, pointer timeseriesmanagermodule::getlink ( class(timeseriesmanagertype this,
character(len=3), intent(in)  auxOrBnd,
integer(i4b), intent(in)  indx 
)
private

Definition at line 411 of file TimeSeriesManager.f90.

412  ! -- dummy
413  class(TimeSeriesManagerType) :: this
414  character(len=3), intent(in) :: auxOrBnd
415  integer(I4B), intent(in) :: indx
416  type(TimeSeriesLinkType), pointer :: tsLink
417  ! -- local
418  type(ListType), pointer :: list
419  !
420  list => null()
421  tslink => null()
422  !
423  select case (auxorbnd)
424  case ('AUX')
425  list => this%auxvarTsLinks
426  case ('BND')
427  list => this%boundTsLinks
428  end select
429  !
430  if (associated(list)) then
431  tslink => gettimeserieslinkfromlist(list, indx)
432  end if
Here is the call graph for this function:

◆ hashbndtimeseries()

subroutine timeseriesmanagermodule::hashbndtimeseries ( class(timeseriesmanagertype), intent(inout)  this)
private

Definition at line 475 of file TimeSeriesManager.f90.

476  ! -- dummy
477  class(TimeSeriesManagerType), intent(inout) :: this
478  ! -- local
479  integer(I4B) :: i, j, k, numtsfiles, numts
480  character(len=LENTIMESERIESNAME) :: name
481  type(TimeSeriesFileType), pointer :: tsfile => null()
482  !
483  ! Initialize the hash table
484  call hash_table_cr(this%BndTsHashTable)
485  !
486  ! Allocate the TsContainers array to accommodate all time-series links.
487  numts = this%tsfileList%CountTimeSeries()
488  allocate (this%TsContainers(numts))
489  !
490  ! Store a pointer to each time series in the TsContainers array
491  ! and put its key (time-series name) and index in the hash table.
492  numtsfiles = this%tsfileList%Counttsfiles()
493  k = 0
494  do i = 1, numtsfiles
495  tsfile => this%tsfileList%Gettsfile(i)
496  numts = tsfile%Count()
497  do j = 1, numts
498  k = k + 1
499  this%TsContainers(k)%timeSeries => tsfile%GetTimeSeries(j)
500  if (associated(this%TsContainers(k)%timeSeries)) then
501  name = this%TsContainers(k)%timeSeries%Name
502  call this%BndTsHashTable%add(name, k)
503  end if
504  end do
505  end do
Here is the call graph for this function:

◆ make_link()

subroutine timeseriesmanagermodule::make_link ( class(timeseriesmanagertype), intent(inout)  this,
type(timeseriestype), intent(inout), pointer  timeSeries,
character(len=*), intent(in)  pkgName,
character(len=3), intent(in)  auxOrBnd,
real(dp), intent(inout), pointer  bndElem,
integer(i4b), intent(in)  irow,
integer(i4b), intent(in)  jcol,
integer(i4b), intent(in)  iprpak,
type(timeserieslinktype), intent(inout), pointer  tsLink,
character(len=*), intent(in)  text,
character(len=*), intent(in)  bndName 
)
private

Definition at line 379 of file TimeSeriesManager.f90.

381  ! -- dummy
382  class(TimeSeriesManagerType), intent(inout) :: this
383  type(TimeSeriesType), pointer, intent(inout) :: timeSeries
384  character(len=*), intent(in) :: pkgName
385  character(len=3), intent(in) :: auxOrBnd
386  real(DP), pointer, intent(inout) :: bndElem
387  integer(I4B), intent(in) :: irow, jcol
388  integer(I4B), intent(in) :: iprpak
389  type(TimeSeriesLinkType), pointer, intent(inout) :: tsLink
390  character(len=*), intent(in) :: text
391  character(len=*), intent(in) :: bndName
392  !
393  tslink => null()
394  call constructtimeserieslink(tslink, timeseries, pkgname, &
395  auxorbnd, bndelem, irow, jcol, iprpak)
396  if (associated(tslink)) then
397  if (auxorbnd == 'BND') then
398  call addtimeserieslinktolist(this%boundTsLinks, tslink)
399  elseif (auxorbnd == 'AUX') then
400  call addtimeserieslinktolist(this%auxvarTsLinks, tslink)
401  else
402  call store_error('programmer error in make_link', terminate=.true.)
403  end if
404  tslink%Text = text
405  tslink%BndName = bndname
406  end if
Here is the call graph for this function:

◆ read_value_or_time_series()

subroutine, public timeseriesmanagermodule::read_value_or_time_series ( character(len=*), intent(in)  textInput,
integer(i4b), intent(in)  ii,
integer(i4b), intent(in)  jj,
real(dp), intent(inout), pointer  bndElem,
character(len=*), intent(in)  pkgName,
character(len=3), intent(in)  auxOrBnd,
type(timeseriesmanagertype), intent(inout)  tsManager,
integer(i4b), intent(in)  iprpak,
type(timeserieslinktype), intent(inout), pointer  tsLink 
)

This routine assumes that there is not an existing link for the specified package and array row and column. For the standard boundary package, all links are removed by calling the tsmanagerreset() method. For advanced packages, there is a separate routine called read_value_or_time_series_adv, which should be used instead of this one.

Definition at line 518 of file TimeSeriesManager.f90.

520  ! dummy
521  character(len=*), intent(in) :: textInput
522  integer(I4B), intent(in) :: ii
523  integer(I4B), intent(in) :: jj
524  real(DP), pointer, intent(inout) :: bndElem
525  character(len=*), intent(in) :: pkgName
526  character(len=3), intent(in) :: auxOrBnd
527  type(TimeSeriesManagerType), intent(inout) :: tsManager
528  integer(I4B), intent(in) :: iprpak
529  type(TimeSeriesLinkType), pointer, intent(inout) :: tsLink
530  ! local
531  type(TimeSeriesType), pointer :: timeseries => null()
532  integer(I4B) :: istat
533  real(DP) :: r
534  character(len=LINELENGTH) :: errmsg
535  character(len=LENTIMESERIESNAME) :: tsNameTemp
536 
537  read (textinput, *, iostat=istat) r
538  if (istat == 0) then
539  bndelem = r
540  else
541 
542  ! Check to see if this is a time series name
543  tsnametemp = textinput
544  call upcase(tsnametemp)
545  timeseries => tsmanager%get_time_series(tsnametemp)
546 
547  ! If this is a time series, then create a link between
548  ! the time series and the row and column position of
549  ! bndElem
550  if (associated(timeseries)) then
551  ! -- Assign value from time series to current
552  ! array element
553  r = timeseries%GetValue(totimsav, totim, &
554  tsmanager%extendTsToEndOfSimulation)
555  bndelem = r
556  ! Make a new link between the time series and the row and
557  ! column position in the array
558  call tsmanager%make_link(timeseries, pkgname, auxorbnd, bndelem, &
559  ii, jj, iprpak, tslink, '', '')
560  else
561  errmsg = 'Error in list input. Expected numeric value or '// &
562  "time-series name, but found '"//trim(textinput)//"'."
563  call store_error(errmsg)
564  end if
565  end if
566 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_value_or_time_series_adv()

subroutine, public timeseriesmanagermodule::read_value_or_time_series_adv ( character(len=*), intent(in)  textInput,
integer(i4b), intent(in)  ii,
integer(i4b), intent(in)  jj,
real(dp), intent(inout), pointer  bndElem,
character(len=*), intent(in)  pkgName,
character(len=3), intent(in)  auxOrBnd,
type(timeseriesmanagertype), intent(inout)  tsManager,
integer(i4b), intent(in)  iprpak,
character(len=*), intent(in)  varName 
)

Arguments are as follows: textInput : string that is either a float or a string name ii : column number jj : row number bndElem : pointer to a position in an array in package pkgName pkgName : package name auxOrBnd : 'AUX' or 'BND' keyword tsManager : timeseries manager object for package iprpak : integer flag indicating if interpolated timeseries values should be printed to package iout during TsManagerad() varName : variable name

Definition at line 584 of file TimeSeriesManager.f90.

586  ! -- dummy
587  character(len=*), intent(in) :: textInput
588  integer(I4B), intent(in) :: ii
589  integer(I4B), intent(in) :: jj
590  real(DP), pointer, intent(inout) :: bndElem
591  character(len=*), intent(in) :: pkgName
592  character(len=3), intent(in) :: auxOrBnd
593  type(TimeSeriesManagerType), intent(inout) :: tsManager
594  integer(I4B), intent(in) :: iprpak
595  character(len=*), intent(in) :: varName
596  ! -- local
597  integer(I4B) :: istat
598  real(DP) :: v
599  character(len=LINELENGTH) :: errmsg
600  character(len=LENTIMESERIESNAME) :: tsNameTemp
601  logical :: found
602  type(TimeSeriesType), pointer :: timeseries => null()
603  type(TimeSeriesLinkType), pointer :: tsLink => null()
604  !
605  ! -- attempt to read textInput as a real value
606  read (textinput, *, iostat=istat) v
607  !
608  ! -- numeric value
609  if (istat == 0) then
610  !
611  ! -- Numeric value was successfully read.
612  bndelem = v
613  !
614  ! -- remove existing link if it exists for this boundary element
615  found = remove_existing_link(tsmanager, ii, jj, pkgname, &
616  auxorbnd, varname)
617  !
618  ! -- timeseries
619  else
620  !
621  ! -- attempt to read numeric value from textInput failed.
622  ! Text should be a time-series name.
623  tsnametemp = textinput
624  call upcase(tsnametemp)
625  !
626  ! -- if textInput is a time-series name, get average value
627  ! from time series.
628  timeseries => tsmanager%get_time_series(tsnametemp)
629  !
630  ! -- create a time series link and add it to the package
631  ! list of time series links used by the array.
632  if (associated(timeseries)) then
633  !
634  ! -- Assign average value from time series to current array element
635  v = timeseries%GetValue(totimsav, totim, &
636  tsmanager%extendTsToEndOfSimulation)
637  bndelem = v
638  !
639  ! -- remove existing link if it exists for this boundary element
640  found = remove_existing_link(tsmanager, ii, jj, &
641  pkgname, auxorbnd, varname)
642  !
643  ! -- Add link to the list.
644  call tsmanager%make_link(timeseries, pkgname, auxorbnd, bndelem, &
645  ii, jj, iprpak, tslink, varname, '')
646  !
647  ! -- not a valid timeseries name
648  else
649  errmsg = 'Error in list input. Expected numeric value or '// &
650  "time-series name, but found '"//trim(textinput)//"'."
651  call store_error(errmsg)
652  end if
653  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ remove_existing_link()

logical function timeseriesmanagermodule::remove_existing_link ( type(timeseriesmanagertype), intent(inout)  tsManager,
integer(i4b), intent(in)  ii,
integer(i4b), intent(in)  jj,
character(len=*), intent(in)  pkgName,
character(len=3), intent(in)  auxOrBnd,
character(len=*), intent(in)  varName 
)
private

Arguments are as follows: tsManager : timeseries manager object for package ii : column number jj : row number pkgName : package name auxOrBnd : 'AUX' or 'BND' keyword varName : variable name

Definition at line 668 of file TimeSeriesManager.f90.

670  ! -- return
671  logical :: found
672  ! -- dummy
673  type(TimeSeriesManagerType), intent(inout) :: tsManager
674  integer(I4B), intent(in) :: ii
675  integer(I4B), intent(in) :: jj
676  character(len=*), intent(in) :: pkgName
677  character(len=3), intent(in) :: auxOrBnd
678  character(len=*), intent(in) :: varName
679  ! -- local
680  integer(I4B) :: i
681  integer(I4B) :: nlinks
682  integer(I4B) :: removeLink
683  type(TimeSeriesLinkType), pointer :: tslTemp => null()
684  !
685  ! -- determine if link exists
686  nlinks = tsmanager%CountLinks(auxorbnd)
687  found = .false.
688  removelink = -1
689  csearchlinks: do i = 1, nlinks
690  tsltemp => tsmanager%GetLink(auxorbnd, i)
691  !
692  ! -- Check ii against iRow, jj against jCol, and varName
693  ! against Text member of link
694  if (tsltemp%PackageName == pkgname) then
695  !
696  ! -- This array element is already linked to a time series.
697  if (tsltemp%IRow == ii .and. tsltemp%JCol == jj .and. &
698  same_word(tsltemp%Text, varname)) then
699  found = .true.
700  removelink = i
701  exit csearchlinks
702  end if
703  end if
704  end do csearchlinks
705  !
706  ! -- remove link if it was found
707  if (removelink > 0) then
708  if (auxorbnd == 'BND') then
709  call tsmanager%boundTsLinks%RemoveNode(removelink, .true.)
710  else if (auxorbnd == 'AUX') then
711  call tsmanager%auxvarTsLinks%RemoveNode(removelink, .true.)
712  end if
713  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ reset()

subroutine timeseriesmanagermodule::reset ( class(timeseriesmanagertype this,
character(len=*), intent(in)  pkgName 
)
private

Definition at line 331 of file TimeSeriesManager.f90.

332  ! -- dummy
333  class(TimeSeriesManagerType) :: this
334  character(len=*), intent(in) :: pkgName
335  ! -- local
336  integer(I4B) :: i, nlinks
337  type(TimeSeriesLinkType), pointer :: tslink
338  !
339  ! Zero out values for time-series controlled stresses.
340  ! Also deallocate all tslinks too.
341  ! Then when time series are
342  ! specified in this or another stress period,
343  ! a new tslink would be set up.
344  !
345  ! Reassign all linked elements to zero
346  nlinks = this%boundTsLinks%Count()
347  do i = 1, nlinks
348  tslink => gettimeserieslinkfromlist(this%boundTsLinks, i)
349  if (associated(tslink)) then
350  if (tslink%PackageName == pkgname) then
351  tslink%BndElement = dzero
352  end if
353  end if
354  end do
355  !
356  ! Remove links belonging to calling package
357  nlinks = this%boundTsLinks%Count()
358  do i = nlinks, 1, -1
359  tslink => gettimeserieslinkfromlist(this%boundTsLinks, i)
360  if (associated(tslink)) then
361  if (tslink%PackageName == pkgname) then
362  call this%boundTsLinks%RemoveNode(i, .true.)
363  end if
364  end if
365  end do
366  nlinks = this%auxvarTsLinks%Count()
367  do i = nlinks, 1, -1
368  tslink => gettimeserieslinkfromlist(this%auxvarTsLinks, i)
369  if (associated(tslink)) then
370  if (tslink%PackageName == pkgname) then
371  call this%auxvarTsLinks%RemoveNode(i, .true.)
372  end if
373  end if
374  end do
Here is the call graph for this function:

◆ tsmanager_cr()

subroutine, public timeseriesmanagermodule::tsmanager_cr ( type(timeseriesmanagertype this,
integer(i4b), intent(in)  iout,
logical, intent(in), optional  removeTsLinksOnCompletion,
logical, intent(in), optional  extendTsToEndOfSimulation 
)

Definition at line 62 of file TimeSeriesManager.f90.

64  ! -- dummy
65  type(TimeSeriesManagerType) :: this
66  integer(I4B), intent(in) :: iout
67  logical, intent(in), optional :: removeTsLinksOnCompletion
68  logical, intent(in), optional :: extendTsToEndOfSimulation
69  !
70  this%iout = iout
71  if (present(removetslinksoncompletion)) then
72  this%removeTsLinksOnCompletion = removetslinksoncompletion
73  end if
74  if (present(extendtstoendofsimulation)) then
75  this%extendTsToEndOfSimulation = extendtstoendofsimulation
76  end if
77  allocate (this%boundTsLinks)
78  allocate (this%auxvarTsLinks)
79  allocate (this%tsfileList)
80  allocate (this%tsfiles(1000))
Here is the caller graph for this function:

◆ tsmanager_df()

subroutine timeseriesmanagermodule::tsmanager_df ( class(timeseriesmanagertype this)
private

Definition at line 85 of file TimeSeriesManager.f90.

86  ! -- dummy
87  class(TimeSeriesManagerType) :: this
88  !
89  if (this%numtsfiles > 0) then
90  call this%HashBndTimeSeries()
91  end if

◆ tsmgr_ad()

subroutine timeseriesmanagermodule::tsmgr_ad ( class(timeseriesmanagertype this)

Definition at line 134 of file TimeSeriesManager.f90.

135  ! -- dummy
136  class(TimeSeriesManagerType) :: this
137  ! -- local
138  type(TimeSeriesLinkType), pointer :: tsLink => null()
139  type(TimeSeriesType), pointer :: timeseries => null()
140  integer(I4B) :: i, nlinks, nauxlinks
141  real(DP) :: begintime, endtime, tsendtime
142  character(len=LENPACKAGENAME + 2) :: pkgID
143  ! -- formats
144  character(len=*), parameter :: fmt5 = &
145  "(/,'Time-series controlled values in stress period: ', i0, &
146  &', time step ', i0, ':')"
147 10 format(a, ' package: Boundary ', i0, ', entry ', i0, &
148  ' value from time series "', a, '" = ', g12.5)
149 15 format(a, ' package: Boundary ', i0, ', entry ', i0, &
150  ' value from time series "', a, '" = ', g12.5, ' (', a, ')')
151 20 format(a, ' package: Boundary ', i0, ', ', a, &
152  ' value from time series "', a, '" = ', g12.5)
153 25 format(a, ' package: Boundary ', i0, ', ', a, &
154  ' value from time series "', a, '" = ', g12.5, ' (', a, ')')
155  !
156  ! -- Initialize time variables
157  begintime = totimc
158  endtime = begintime + delt
159  !
160  ! -- Determine number of ts links
161  nlinks = this%boundtslinks%Count()
162  nauxlinks = this%auxvartslinks%Count()
163  !
164  ! -- Iterate through auxvartslinks and replace specified
165  ! elements of auxvar with average value obtained from
166  ! appropriate time series. Need to do auxvartslinks
167  ! first because they may be a multiplier column
168  i = 1
169  do while (i <= nauxlinks)
170  tslink => gettimeserieslinkfromlist(this%auxvarTsLinks, i)
171  timeseries => tslink%timeSeries
172  !
173  ! -- Remove time series link once its end time has passed, if requested
174  if (this%removeTsLinksOnCompletion) then
175  tsendtime = timeseries%FindLatestTime(.true.)
176  if (tsendtime < begintime) then
177  call this%auxvarTsLinks%RemoveNode(i, .true.)
178  nauxlinks = this%auxvartslinks%Count()
179  cycle
180  end if
181  end if
182  !
183  if (i == 1) then
184  if (tslink%Iprpak == 1) then
185  write (this%iout, fmt5) kper, kstp
186  end if
187  end if
188  tslink%BndElement = timeseries%GetValue(begintime, endtime, &
189  this%extendTsToEndOfSimulation)
190  !
191  ! -- Write time series values to output file
192  if (tslink%Iprpak == 1) then
193  pkgid = '"'//trim(tslink%PackageName)//'"'
194  if (tslink%Text == '') then
195  if (tslink%BndName == '') then
196  write (this%iout, 10) trim(pkgid), tslink%IRow, tslink%JCol, &
197  trim(tslink%timeSeries%Name), &
198  tslink%BndElement
199  else
200  write (this%iout, 15) trim(pkgid), tslink%IRow, tslink%JCol, &
201  trim(tslink%timeSeries%Name), &
202  tslink%BndElement, trim(tslink%BndName)
203  end if
204  else
205  if (tslink%BndName == '') then
206  write (this%iout, 20) trim(pkgid), tslink%IRow, trim(tslink%Text), &
207  trim(tslink%timeSeries%Name), &
208  tslink%BndElement
209  else
210  write (this%iout, 25) trim(pkgid), tslink%IRow, trim(tslink%Text), &
211  trim(tslink%timeSeries%Name), &
212  tslink%BndElement, trim(tslink%BndName)
213  end if
214  end if
215  end if
216  !
217  i = i + 1
218  end do
219  !
220  ! -- Iterate through boundtslinks and replace specified
221  ! elements of bound with average value obtained from
222  ! appropriate time series. (For list-type packages)
223  i = 1
224  do while (i <= nlinks)
225  tslink => gettimeserieslinkfromlist(this%boundTsLinks, i)
226  timeseries => tslink%timeSeries
227  !
228  ! -- Remove time series link once its end time has passed, if requested
229  if (this%removeTsLinksOnCompletion) then
230  tsendtime = timeseries%FindLatestTime(.true.)
231  if (tsendtime < begintime) then
232  call this%boundTsLinks%RemoveNode(i, .true.)
233  nlinks = this%boundTsLinks%Count()
234  cycle
235  end if
236  end if
237  !
238  if (i == 1 .and. nauxlinks == 0) then
239  if (tslink%Iprpak == 1) then
240  write (this%iout, fmt5) kper, kstp
241  end if
242  end if
243  ! this part needs to be different for MAW because MAW does not use
244  ! bound array for well rate (although rate is stored in
245  ! this%bound(4,ibnd)), it uses this%mawwells(n)%rate%value
246  if (tslink%UseDefaultProc) then
247  timeseries => tslink%timeSeries
248  tslink%BndElement = timeseries%GetValue(begintime, endtime, &
249  this%extendTsToEndOfSimulation)
250  !
251  ! -- If multiplier is active and it applies to this element,
252  ! do the multiplication. This must be done after the auxlinks
253  ! have been calculated in case iauxmultcol is being used.
254  if (associated(tslink%RMultiplier)) then
255  tslink%BndElement = tslink%BndElement * tslink%RMultiplier
256  end if
257  !
258  ! -- Write time series values to output files
259  if (tslink%Iprpak == 1) then
260  pkgid = '"'//trim(tslink%PackageName)//'"'
261  if (tslink%Text == '') then
262  if (tslink%BndName == '') then
263  write (this%iout, 10) trim(pkgid), tslink%IRow, tslink%JCol, &
264  trim(tslink%timeSeries%Name), &
265  tslink%BndElement
266  else
267  write (this%iout, 15) trim(pkgid), tslink%IRow, tslink%JCol, &
268  trim(tslink%timeSeries%Name), &
269  tslink%BndElement, trim(tslink%BndName)
270  end if
271  else
272  if (tslink%BndName == '') then
273  write (this%iout, 20) trim(pkgid), tslink%IRow, trim(tslink%Text), &
274  trim(tslink%timeSeries%Name), &
275  tslink%BndElement
276  else
277  write (this%iout, 25) trim(pkgid), tslink%IRow, trim(tslink%Text), &
278  trim(tslink%timeSeries%Name), &
279  tslink%BndElement, trim(tslink%BndName)
280  end if
281  end if
282  end if
283  !
284  ! -- If conversion from flux to flow is required, multiply by cell area
285  if (tslink%ConvertFlux) then
286  tslink%BndElement = tslink%BndElement * tslink%CellArea
287  end if
288  end if
289  !
290  i = i + 1
291  end do
292  !
293  ! -- Finish with ending line
294  if (nlinks + nauxlinks > 0) then
295  if (tslink%Iprpak == 1) then
296  write (this%iout, '()')
297  end if
298  end if
Here is the call graph for this function:

◆ tsmgr_da()

subroutine timeseriesmanagermodule::tsmgr_da ( class(timeseriesmanagertype this)
private

Definition at line 303 of file TimeSeriesManager.f90.

304  ! -- dummy
305  class(TimeSeriesManagerType) :: this
306  ! -- local
307  !
308  ! -- Deallocate time-series links in boundTsLinks
309  call this%boundTsLinks%Clear(.true.)
310  deallocate (this%boundTsLinks)
311  !
312  ! -- Deallocate time-series links in auxvarTsLinks
313  call this%auxvarTsLinks%Clear(.true.)
314  deallocate (this%auxvarTsLinks)
315  !
316  ! -- Deallocate tsfileList
317  call this%tsfileList%da()
318  deallocate (this%tsfileList)
319  !
320  ! -- Deallocate the hash table
321  if (associated(this%BndTsHashTable)) then
322  call hash_table_da(this%BndTsHashTable)
323  end if
324  !
325  deallocate (this%tsfiles)
Here is the call graph for this function:

◆ var_timeseries()

logical function, public timeseriesmanagermodule::var_timeseries ( type(timeseriesmanagertype), intent(inout)  tsManager,
character(len=*), intent(in)  pkgName,
character(len=*), intent(in)  varName,
character(len=3), intent(in), optional  auxOrBnd 
)

Arguments are as follows: tsManager : timeseries manager object for package pkgName : package name varName : variable name auxOrBnd : optional 'AUX' or 'BND' keyword

Definition at line 724 of file TimeSeriesManager.f90.

725  ! -- return
726  logical :: tsexists
727  ! -- dummy
728  type(TimeSeriesManagerType), intent(inout) :: tsManager
729  character(len=*), intent(in) :: pkgName
730  character(len=*), intent(in) :: varName
731  character(len=3), intent(in), optional :: auxOrBnd
732  ! -- local
733  character(len=3) :: ctstype
734  integer(I4B) :: i
735  integer(I4B) :: nlinks
736  type(TimeSeriesLinkType), pointer :: tslTemp => null()
737  !
738  ! -- process optional variables
739  if (present(auxorbnd)) then
740  ctstype = auxorbnd
741  else
742  ctstype = 'BND'
743  end if
744  !
745  ! -- initialize the return variable and the number of timeseries links
746  tsexists = .false.
747  nlinks = tsmanager%CountLinks(ctstype)
748  !
749  ! -- determine if link exists
750  csearchlinks: do i = 1, nlinks
751  tsltemp => tsmanager%GetLink(ctstype, i)
752  if (tsltemp%PackageName == pkgname) then
753  !
754  ! -- Check varName against Text member of link
755  if (same_word(tsltemp%Text, varname)) then
756  tsexists = .true.
757  exit csearchlinks
758  end if
759  end if
760  end do csearchlinks
Here is the call graph for this function:
Here is the caller graph for this function: