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

103  ! -- modules
106  ! -- dummy
107  class(TimeSeriesManagerType) :: this
108  character(len=*), intent(in) :: fname
109  integer(I4B), intent(in) :: inunit
110  ! -- local
111  integer(I4B) :: isize
112  integer(I4B) :: i
113  class(TimeSeriesFileType), pointer :: tsfile => null()
114  !
115  ! -- Check for fname duplicates
116  if (this%numtsfiles > 0) then
117  do i = 1, this%numtsfiles
118  if (this%tsfiles(i) == fname) then
119  call store_error('Found duplicate time-series file name: '//trim(fname))
120  call store_error_unit(inunit)
121  end if
122  end do
123  end if
124  !
125  ! -- Save fname
126  this%numtsfiles = this%numtsfiles + 1
127  isize = size(this%tsfiles)
128  if (this%numtsfiles > isize) then
129  call expandarray(this%tsfiles, 1000)
130  end if
131  this%tsfiles(this%numtsfiles) = fname
132  !
133  ! --
134  call this%tsfileList%Add(fname, this%iout, tsfile)
135  !
136  ! -- Return
137  return
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 461 of file TimeSeriesManager.f90.

462  ! -- return
463  integer(I4B) :: CountLinks
464  ! -- dummy
465  class(TimeSeriesManagerType) :: this
466  character(len=3), intent(in) :: auxOrBnd
467  !
468  countlinks = 0
469  if (auxorbnd == 'BND') then
470  countlinks = this%boundTsLinks%Count()
471  elseif (auxorbnd == 'AUX') then
472  countlinks = this%auxvarTsLinks%count()
473  end if
474  !
475  ! -- Return
476  return

◆ get_time_series()

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

Definition at line 481 of file TimeSeriesManager.f90.

482  ! -- dummy
483  class(TimeSeriesManagerType) :: this
484  character(len=*), intent(in) :: name
485  ! -- result
486  type(TimeSeriesType), pointer :: res
487  ! -- local
488  integer(I4B) :: indx
489  !
490  ! Get index from hash table, get time series from TsContainers,
491  ! and assign result to time series contained in link.
492  res => null()
493  indx = this%BndTsHashTable%get(name)
494  if (indx > 0) then
495  res => this%TsContainers(indx)%timeSeries
496  end if
497  !
498  ! -- Return
499  return

◆ 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 432 of file TimeSeriesManager.f90.

433  ! -- dummy
434  class(TimeSeriesManagerType) :: this
435  character(len=3), intent(in) :: auxOrBnd
436  integer(I4B), intent(in) :: indx
437  type(TimeSeriesLinkType), pointer :: tsLink
438  ! -- local
439  type(ListType), pointer :: list
440  !
441  list => null()
442  tslink => null()
443  !
444  select case (auxorbnd)
445  case ('AUX')
446  list => this%auxvarTsLinks
447  case ('BND')
448  list => this%boundTsLinks
449  end select
450  !
451  if (associated(list)) then
452  tslink => gettimeserieslinkfromlist(list, indx)
453  end if
454  !
455  ! -- Return
456  return
Here is the call graph for this function:

◆ hashbndtimeseries()

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

Definition at line 505 of file TimeSeriesManager.f90.

506  ! -- dummy
507  class(TimeSeriesManagerType), intent(inout) :: this
508  ! -- local
509  integer(I4B) :: i, j, k, numtsfiles, numts
510  character(len=LENTIMESERIESNAME) :: name
511  type(TimeSeriesFileType), pointer :: tsfile => null()
512  !
513  ! Initialize the hash table
514  call hash_table_cr(this%BndTsHashTable)
515  !
516  ! Allocate the TsContainers array to accommodate all time-series links.
517  numts = this%tsfileList%CountTimeSeries()
518  allocate (this%TsContainers(numts))
519  !
520  ! Store a pointer to each time series in the TsContainers array
521  ! and put its key (time-series name) and index in the hash table.
522  numtsfiles = this%tsfileList%Counttsfiles()
523  k = 0
524  do i = 1, numtsfiles
525  tsfile => this%tsfileList%Gettsfile(i)
526  numts = tsfile%Count()
527  do j = 1, numts
528  k = k + 1
529  this%TsContainers(k)%timeSeries => tsfile%GetTimeSeries(j)
530  if (associated(this%TsContainers(k)%timeSeries)) then
531  name = this%TsContainers(k)%timeSeries%Name
532  call this%BndTsHashTable%add(name, k)
533  end if
534  end do
535  end do
536  !
537  ! -- Return
538  return
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 397 of file TimeSeriesManager.f90.

399  ! -- dummy
400  class(TimeSeriesManagerType), intent(inout) :: this
401  type(TimeSeriesType), pointer, intent(inout) :: timeSeries
402  character(len=*), intent(in) :: pkgName
403  character(len=3), intent(in) :: auxOrBnd
404  real(DP), pointer, intent(inout) :: bndElem
405  integer(I4B), intent(in) :: irow, jcol
406  integer(I4B), intent(in) :: iprpak
407  type(TimeSeriesLinkType), pointer, intent(inout) :: tsLink
408  character(len=*), intent(in) :: text
409  character(len=*), intent(in) :: bndName
410  !
411  tslink => null()
412  call constructtimeserieslink(tslink, timeseries, pkgname, &
413  auxorbnd, bndelem, irow, jcol, iprpak)
414  if (associated(tslink)) then
415  if (auxorbnd == 'BND') then
416  call addtimeserieslinktolist(this%boundTsLinks, tslink)
417  elseif (auxorbnd == 'AUX') then
418  call addtimeserieslinktolist(this%auxvarTsLinks, tslink)
419  else
420  call store_error('programmer error in make_link', terminate=.true.)
421  end if
422  tslink%Text = text
423  tslink%BndName = bndname
424  end if
425  !
426  ! -- Return
427  return
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 551 of file TimeSeriesManager.f90.

553  ! dummy
554  character(len=*), intent(in) :: textInput
555  integer(I4B), intent(in) :: ii
556  integer(I4B), intent(in) :: jj
557  real(DP), pointer, intent(inout) :: bndElem
558  character(len=*), intent(in) :: pkgName
559  character(len=3), intent(in) :: auxOrBnd
560  type(TimeSeriesManagerType), intent(inout) :: tsManager
561  integer(I4B), intent(in) :: iprpak
562  type(TimeSeriesLinkType), pointer, intent(inout) :: tsLink
563  ! local
564  type(TimeSeriesType), pointer :: timeseries => null()
565  integer(I4B) :: istat
566  real(DP) :: r
567  character(len=LINELENGTH) :: errmsg
568  character(len=LENTIMESERIESNAME) :: tsNameTemp
569 
570  read (textinput, *, iostat=istat) r
571  if (istat == 0) then
572  bndelem = r
573  else
574 
575  ! Check to see if this is a time series name
576  tsnametemp = textinput
577  call upcase(tsnametemp)
578  timeseries => tsmanager%get_time_series(tsnametemp)
579 
580  ! If this is a time series, then create a link between
581  ! the time series and the row and column position of
582  ! bndElem
583  if (associated(timeseries)) then
584  ! -- Assign value from time series to current
585  ! array element
586  r = timeseries%GetValue(totimsav, totim, &
587  tsmanager%extendTsToEndOfSimulation)
588  bndelem = r
589  ! Make a new link between the time series and the row and
590  ! column position in the array
591  call tsmanager%make_link(timeseries, pkgname, auxorbnd, bndelem, &
592  ii, jj, iprpak, tslink, '', '')
593  else
594  errmsg = 'Error in list input. Expected numeric value or '// &
595  "time-series name, but found '"//trim(textinput)//"'."
596  call store_error(errmsg)
597  end if
598  end if
599 
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 617 of file TimeSeriesManager.f90.

619  ! -- dummy
620  character(len=*), intent(in) :: textInput
621  integer(I4B), intent(in) :: ii
622  integer(I4B), intent(in) :: jj
623  real(DP), pointer, intent(inout) :: bndElem
624  character(len=*), intent(in) :: pkgName
625  character(len=3), intent(in) :: auxOrBnd
626  type(TimeSeriesManagerType), intent(inout) :: tsManager
627  integer(I4B), intent(in) :: iprpak
628  character(len=*), intent(in) :: varName
629  ! -- local
630  integer(I4B) :: istat
631  real(DP) :: v
632  character(len=LINELENGTH) :: errmsg
633  character(len=LENTIMESERIESNAME) :: tsNameTemp
634  logical :: found
635  type(TimeSeriesType), pointer :: timeseries => null()
636  type(TimeSeriesLinkType), pointer :: tsLink => null()
637  !
638  ! -- attempt to read textInput as a real value
639  read (textinput, *, iostat=istat) v
640  !
641  ! -- numeric value
642  if (istat == 0) then
643  !
644  ! -- Numeric value was successfully read.
645  bndelem = v
646  !
647  ! -- remove existing link if it exists for this boundary element
648  found = remove_existing_link(tsmanager, ii, jj, pkgname, &
649  auxorbnd, varname)
650  !
651  ! -- timeseries
652  else
653  !
654  ! -- attempt to read numeric value from textInput failed.
655  ! Text should be a time-series name.
656  tsnametemp = textinput
657  call upcase(tsnametemp)
658  !
659  ! -- if textInput is a time-series name, get average value
660  ! from time series.
661  timeseries => tsmanager%get_time_series(tsnametemp)
662  !
663  ! -- create a time series link and add it to the package
664  ! list of time series links used by the array.
665  if (associated(timeseries)) then
666  !
667  ! -- Assign average value from time series to current array element
668  v = timeseries%GetValue(totimsav, totim, &
669  tsmanager%extendTsToEndOfSimulation)
670  bndelem = v
671  !
672  ! -- remove existing link if it exists for this boundary element
673  found = remove_existing_link(tsmanager, ii, jj, &
674  pkgname, auxorbnd, varname)
675  !
676  ! -- Add link to the list.
677  call tsmanager%make_link(timeseries, pkgname, auxorbnd, bndelem, &
678  ii, jj, iprpak, tslink, varname, '')
679  !
680  ! -- not a valid timeseries name
681  else
682  errmsg = 'Error in list input. Expected numeric value or '// &
683  "time-series name, but found '"//trim(textinput)//"'."
684  call store_error(errmsg)
685  end if
686  end if
687  !
688  ! -- Return
689  return
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 704 of file TimeSeriesManager.f90.

706  ! -- return
707  logical :: found
708  ! -- dummy
709  type(TimeSeriesManagerType), intent(inout) :: tsManager
710  integer(I4B), intent(in) :: ii
711  integer(I4B), intent(in) :: jj
712  character(len=*), intent(in) :: pkgName
713  character(len=3), intent(in) :: auxOrBnd
714  character(len=*), intent(in) :: varName
715  ! -- local
716  integer(I4B) :: i
717  integer(I4B) :: nlinks
718  integer(I4B) :: removeLink
719  type(TimeSeriesLinkType), pointer :: tslTemp => null()
720  !
721  ! -- determine if link exists
722  nlinks = tsmanager%CountLinks(auxorbnd)
723  found = .false.
724  removelink = -1
725  csearchlinks: do i = 1, nlinks
726  tsltemp => tsmanager%GetLink(auxorbnd, i)
727  !
728  ! -- Check ii against iRow, jj against jCol, and varName
729  ! against Text member of link
730  if (tsltemp%PackageName == pkgname) then
731  !
732  ! -- This array element is already linked to a time series.
733  if (tsltemp%IRow == ii .and. tsltemp%JCol == jj .and. &
734  same_word(tsltemp%Text, varname)) then
735  found = .true.
736  removelink = i
737  exit csearchlinks
738  end if
739  end if
740  end do csearchlinks
741  !
742  ! -- remove link if it was found
743  if (removelink > 0) then
744  if (auxorbnd == 'BND') then
745  call tsmanager%boundTsLinks%RemoveNode(removelink, .true.)
746  else if (auxorbnd == 'AUX') then
747  call tsmanager%auxvarTsLinks%RemoveNode(removelink, .true.)
748  end if
749  end if
750  !
751  ! -- Return
752  return
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 346 of file TimeSeriesManager.f90.

347  ! -- dummy
348  class(TimeSeriesManagerType) :: this
349  character(len=*), intent(in) :: pkgName
350  ! -- local
351  integer(I4B) :: i, nlinks
352  type(TimeSeriesLinkType), pointer :: tslink
353  !
354  ! Zero out values for time-series controlled stresses.
355  ! Also deallocate all tslinks too.
356  ! Then when time series are
357  ! specified in this or another stress period,
358  ! a new tslink would be set up.
359  !
360  ! Reassign all linked elements to zero
361  nlinks = this%boundTsLinks%Count()
362  do i = 1, nlinks
363  tslink => gettimeserieslinkfromlist(this%boundTsLinks, i)
364  if (associated(tslink)) then
365  if (tslink%PackageName == pkgname) then
366  tslink%BndElement = dzero
367  end if
368  end if
369  end do
370  !
371  ! Remove links belonging to calling package
372  nlinks = this%boundTsLinks%Count()
373  do i = nlinks, 1, -1
374  tslink => gettimeserieslinkfromlist(this%boundTsLinks, i)
375  if (associated(tslink)) then
376  if (tslink%PackageName == pkgname) then
377  call this%boundTsLinks%RemoveNode(i, .true.)
378  end if
379  end if
380  end do
381  nlinks = this%auxvarTsLinks%Count()
382  do i = nlinks, 1, -1
383  tslink => gettimeserieslinkfromlist(this%auxvarTsLinks, i)
384  if (associated(tslink)) then
385  if (tslink%PackageName == pkgname) then
386  call this%auxvarTsLinks%RemoveNode(i, .true.)
387  end if
388  end if
389  end do
390  !
391  ! -- Return
392  return
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))
81  !
82  ! -- Return
83  return
Here is the caller graph for this function:

◆ tsmanager_df()

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

Definition at line 88 of file TimeSeriesManager.f90.

89  ! -- dummy
90  class(TimeSeriesManagerType) :: this
91  !
92  if (this%numtsfiles > 0) then
93  call this%HashBndTimeSeries()
94  end if
95  !
96  ! -- Return
97  return

◆ tsmgr_ad()

subroutine timeseriesmanagermodule::tsmgr_ad ( class(timeseriesmanagertype this)

Definition at line 143 of file TimeSeriesManager.f90.

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

◆ tsmgr_da()

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

Definition at line 315 of file TimeSeriesManager.f90.

316  ! -- dummy
317  class(TimeSeriesManagerType) :: this
318  ! -- local
319  !
320  ! -- Deallocate time-series links in boundTsLinks
321  call this%boundTsLinks%Clear(.true.)
322  deallocate (this%boundTsLinks)
323  !
324  ! -- Deallocate time-series links in auxvarTsLinks
325  call this%auxvarTsLinks%Clear(.true.)
326  deallocate (this%auxvarTsLinks)
327  !
328  ! -- Deallocate tsfileList
329  call this%tsfileList%da()
330  deallocate (this%tsfileList)
331  !
332  ! -- Deallocate the hash table
333  if (associated(this%BndTsHashTable)) then
334  call hash_table_da(this%BndTsHashTable)
335  end if
336  !
337  deallocate (this%tsfiles)
338  !
339  ! -- Return
340  return
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 763 of file TimeSeriesManager.f90.

764  ! -- return
765  logical :: tsexists
766  ! -- dummy
767  type(TimeSeriesManagerType), intent(inout) :: tsManager
768  character(len=*), intent(in) :: pkgName
769  character(len=*), intent(in) :: varName
770  character(len=3), intent(in), optional :: auxOrBnd
771  ! -- local
772  character(len=3) :: ctstype
773  integer(I4B) :: i
774  integer(I4B) :: nlinks
775  type(TimeSeriesLinkType), pointer :: tslTemp => null()
776  !
777  ! -- process optional variables
778  if (present(auxorbnd)) then
779  ctstype = auxorbnd
780  else
781  ctstype = 'BND'
782  end if
783  !
784  ! -- initialize the return variable and the number of timeseries links
785  tsexists = .false.
786  nlinks = tsmanager%CountLinks(ctstype)
787  !
788  ! -- determine if link exists
789  csearchlinks: do i = 1, nlinks
790  tsltemp => tsmanager%GetLink(ctstype, i)
791  if (tsltemp%PackageName == pkgname) then
792  !
793  ! -- Check varName against Text member of link
794  if (same_word(tsltemp%Text, varname)) then
795  tsexists = .true.
796  exit csearchlinks
797  end if
798  end if
799  end do csearchlinks
800  !
801  ! -- Return
802  return
Here is the call graph for this function:
Here is the caller graph for this function: