MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
TimeArraySeriesLink.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b
5  use inputoutputmodule, only: upcase
6  use listmodule, only: listtype
8 
9  implicit none
10 
11  private
15 
17  ! -- Public members
18  character(len=LENPACKAGENAME), public :: packagename = ''
19  character(len=LENTIMESERIESTEXT), public :: text = ''
20  integer(I4B), public :: iprpak = 1
21  logical, public :: usedefaultproc = .true.
22  logical, public :: convertflux = .false.
23  integer(I4B), dimension(:), pointer, contiguous, public :: nodelist => null()
24  ! BndArray can point to an array in either the bound or auxval
25  ! array of BndType, or any other double precision variable or array
26  ! element that contains a value that could be controlled by a time series.
27  real(dp), dimension(:), pointer, public :: bndarray => null()
28  real(dp), dimension(:), pointer, public :: rmultarray => null()
29  type(timearrayseriestype), pointer, public :: timearrayseries => null()
30 
31  contains
32 
33  procedure, public :: da => tasl_da
35 
36 contains
37 
38  !> @brief Deallocate
39  !<
40  subroutine tasl_da(this)
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()
48  end subroutine tasl_da
49 
50  !> @brief Construct a time series of arrays that are linked
51  !<
52  subroutine constructtimearrayserieslink(newTasLink, timeArraySeries, &
53  pkgName, bndArray, iprpak, text)
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
73  end subroutine constructtimearrayserieslink
74 
75  !> @brief Cast an unlimited polymorphic object as TimeArraySeriesLinkType
76  !<
77  function castastimearrayserieslinktype(obj) result(res)
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
91 
92  !> @brief Add time-array series to list
93  !<
94  subroutine addtimearrayserieslinktolist(list, tasLink)
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)
104  end subroutine addtimearrayserieslinktolist
105 
106  !> @brief Get time-array series from a list and return
107  !<
108  function gettimearrayserieslinkfromlist(list, idx) result(res)
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)
119  end function gettimearrayserieslinkfromlist
120 
121 end module timearrayserieslinkmodule
122 
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:23
integer(i4b), parameter lentimeseriestext
maximum length of a time series text
Definition: Constants.f90:43
subroutine, public upcase(word)
Convert to upper case.
This module defines variable data types.
Definition: kind.f90:8
subroutine tasl_da(this)
Deallocate.
type(timearrayserieslinktype) function, pointer, public gettimearrayserieslinkfromlist(list, idx)
Get time-array series from a list and return.
type(timearrayserieslinktype) function, pointer, private castastimearrayserieslinktype(obj)
Cast an unlimited polymorphic object as TimeArraySeriesLinkType.
subroutine, public addtimearrayserieslinktolist(list, tasLink)
Add time-array series to list.
subroutine, public constructtimearrayserieslink(newTasLink, timeArraySeries, pkgName, bndArray, iprpak, text)
Construct a time series of arrays that are linked.
A generic heterogeneous doubly-linked list.
Definition: List.f90:14