MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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  !
49  ! -- Return
50  return
51  end subroutine tasl_da
52 
53  !> @brief Construct a time series of arrays that are linked
54  !<
55  subroutine constructtimearrayserieslink(newTasLink, timeArraySeries, &
56  pkgName, bndArray, iprpak, text)
57  ! -- dummy
58  type(timearrayserieslinktype), pointer, intent(out) :: newtaslink
59  type(timearrayseriestype), pointer, intent(in) :: timearrayseries
60  character(len=*), intent(in) :: pkgname
61  real(dp), dimension(:), pointer, intent(in) :: bndarray
62  integer(I4B), intent(in) :: iprpak
63  character(len=*), intent(in) :: text
64  ! -- local
65  character(len=LENPACKAGENAME) :: pkgnametemp
66  !
67  allocate (newtaslink)
68  ! Store package name as all caps
69  pkgnametemp = pkgname
70  call upcase(pkgnametemp)
71  newtaslink%PackageName = pkgnametemp
72  newtaslink%timeArraySeries => timearrayseries
73  newtaslink%BndArray => bndarray
74  newtaslink%Iprpak = iprpak
75  newtaslink%Text = text
76  !
77  ! -- Return
78  return
79  end subroutine constructtimearrayserieslink
80 
81  !> @brief Cast an unlimited polymorphic object as TimeArraySeriesLinkType
82  !<
83  function castastimearrayserieslinktype(obj) result(res)
84  ! -- dummy
85  class(*), pointer, intent(inout) :: obj
86  ! -- return
87  type(timearrayserieslinktype), pointer :: res
88  !
89  res => null()
90  if (.not. associated(obj)) return
91  !
92  select type (obj)
93  type is (timearrayserieslinktype)
94  res => obj
95  end select
96  !
97  ! -- Return
98  return
100 
101  !> @brief Add time-array series to list
102  !<
103  subroutine addtimearrayserieslinktolist(list, tasLink)
104  ! -- dummy
105  type(listtype), intent(inout) :: list
106  ! -- return
107  type(timearrayserieslinktype), pointer, intent(inout) :: taslink
108  ! -- local
109  class(*), pointer :: obj
110  !
111  obj => taslink
112  call list%Add(obj)
113  !
114  ! -- Return
115  return
116  end subroutine addtimearrayserieslinktolist
117 
118  !> @brief Get time-array series from a list and return
119  !<
120  function gettimearrayserieslinkfromlist(list, idx) result(res)
121  ! -- dummy
122  type(listtype), intent(inout) :: list
123  integer(I4B), intent(in) :: idx
124  ! -- return
125  type(timearrayserieslinktype), pointer :: res
126  ! -- local
127  class(*), pointer :: obj
128  !
129  obj => list%GetItem(idx)
131  !
132  ! -- Return
133  return
134  end function gettimearrayserieslinkfromlist
135 
136 end module timearrayserieslinkmodule
137 
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:22
integer(i4b), parameter lentimeseriestext
maximum length of a time series text
Definition: Constants.f90:42
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:10