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

Data Types

type  timearrayserieslinktype
 

Functions/Subroutines

subroutine tasl_da (this)
 Deallocate. More...
 
subroutine, public constructtimearrayserieslink (newTasLink, timeArraySeries, pkgName, bndArray, iprpak, text)
 Construct a time series of arrays that are linked. More...
 
type(timearrayserieslinktype) function, pointer, private castastimearrayserieslinktype (obj)
 Cast an unlimited polymorphic object as TimeArraySeriesLinkType. More...
 
subroutine, public addtimearrayserieslinktolist (list, tasLink)
 Add time-array series to list. More...
 
type(timearrayserieslinktype) function, pointer, public gettimearrayserieslinkfromlist (list, idx)
 Get time-array series from a list and return. More...
 

Function/Subroutine Documentation

◆ addtimearrayserieslinktolist()

subroutine, public timearrayserieslinkmodule::addtimearrayserieslinktolist ( type(listtype), intent(inout)  list,
type(timearrayserieslinktype), intent(inout), pointer  tasLink 
)

Definition at line 94 of file TimeArraySeriesLink.f90.

95  ! -- dummy
96  type(ListType), intent(inout) :: list
97  ! -- return
98  type(TimeArraySeriesLinkType), pointer, intent(inout) :: tasLink
99  ! -- local
100  class(*), pointer :: obj
101  !
102  obj => taslink
103  call list%Add(obj)
Here is the caller graph for this function:

◆ castastimearrayserieslinktype()

type(timearrayserieslinktype) function, pointer, private timearrayserieslinkmodule::castastimearrayserieslinktype ( class(*), intent(inout), pointer  obj)
private

Definition at line 77 of file TimeArraySeriesLink.f90.

78  ! -- dummy
79  class(*), pointer, intent(inout) :: obj
80  ! -- return
81  type(TimeArraySeriesLinkType), pointer :: res
82  !
83  res => null()
84  if (.not. associated(obj)) return
85  !
86  select type (obj)
87  type is (timearrayserieslinktype)
88  res => obj
89  end select
Here is the caller graph for this function:

◆ constructtimearrayserieslink()

subroutine, public timearrayserieslinkmodule::constructtimearrayserieslink ( type(timearrayserieslinktype), intent(out), pointer  newTasLink,
type(timearrayseriestype), intent(in), pointer  timeArraySeries,
character(len=*), intent(in)  pkgName,
real(dp), dimension(:), intent(in), pointer  bndArray,
integer(i4b), intent(in)  iprpak,
character(len=*), intent(in)  text 
)

Definition at line 52 of file TimeArraySeriesLink.f90.

54  ! -- dummy
55  type(TimeArraySeriesLinkType), pointer, intent(out) :: newTasLink
56  type(TimeArraySeriesType), pointer, intent(in) :: timeArraySeries
57  character(len=*), intent(in) :: pkgName
58  real(DP), dimension(:), pointer, intent(in) :: bndArray
59  integer(I4B), intent(in) :: iprpak
60  character(len=*), intent(in) :: text
61  ! -- local
62  character(len=LENPACKAGENAME) :: pkgNameTemp
63  !
64  allocate (newtaslink)
65  ! Store package name as all caps
66  pkgnametemp = pkgname
67  call upcase(pkgnametemp)
68  newtaslink%PackageName = pkgnametemp
69  newtaslink%timeArraySeries => timearrayseries
70  newtaslink%BndArray => bndarray
71  newtaslink%Iprpak = iprpak
72  newtaslink%Text = text
Here is the call graph for this function:
Here is the caller graph for this function:

◆ gettimearrayserieslinkfromlist()

type(timearrayserieslinktype) function, pointer, public timearrayserieslinkmodule::gettimearrayserieslinkfromlist ( type(listtype), intent(inout)  list,
integer(i4b), intent(in)  idx 
)

Definition at line 108 of file TimeArraySeriesLink.f90.

109  ! -- dummy
110  type(ListType), intent(inout) :: list
111  integer(I4B), intent(in) :: idx
112  ! -- return
113  type(TimeArraySeriesLinkType), pointer :: res
114  ! -- local
115  class(*), pointer :: obj
116  !
117  obj => list%GetItem(idx)
118  res => castastimearrayserieslinktype(obj)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ tasl_da()

subroutine timearrayserieslinkmodule::tasl_da ( class(timearrayserieslinktype), intent(inout)  this)
private

Definition at line 40 of file TimeArraySeriesLink.f90.

41  ! -- dummy
42  class(TimeArraySeriesLinkType), intent(inout) :: this
43  !
44  this%nodelist => null()
45  this%bndarray => null()
46  this%rmultarray => null()
47  this%TimeArraySeries => null()