MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
TimeSeriesLink.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b
6  use inputoutputmodule, only: upcase
7  use listmodule, only: listtype
9 
10  implicit none
11 
12  private
15  private :: castastimeserieslinktype
16 
18  ! -- Public members
19  integer(I4B), public :: irow = 0 ! row index (2nd dim) in bound or auxval array
20  integer(I4B), public :: jcol = 0 ! column index (1st dim) in bound or auxval array
21  integer(I4B), public :: iprpak = 1
22  ! BndElement can point to an element in either the bound or auxval
23  ! array of BndType, or any other double precision variable or array
24  ! element that contains a value that could be controlled by a time series.
25  real(dp), pointer, public :: bndelement => null()
26  real(dp), pointer, public :: rmultiplier => null()
27  real(dp), public :: cellarea = dzero
28  character(len=LENPACKAGENAME), public :: packagename = ''
29  character(len=3), public :: auxorbnd = ''
30  character(len=LENTIMESERIESTEXT), public :: text = ''
31  character(len=LENBOUNDNAME), public :: bndname = ''
32  logical, public :: active = .true.
33  logical, public :: usedefaultproc = .true.
34  logical, public :: convertflux = .false.
35  type(timeseriestype), pointer, public :: timeseries => null()
36  end type timeserieslinktype
37 
38 contains
39 
40  !> @brief Construct time series link
41  !<
42  subroutine constructtimeserieslink(newTsLink, timeSeries, pkgName, &
43  auxOrBnd, bndElem, iRow, jCol, iprpak, &
44  text)
45  implicit none
46  ! -- dummy
47  type(timeserieslinktype), pointer, intent(out) :: newtslink
48  type(timeseriestype), pointer, intent(in) :: timeseries
49  integer(I4B), intent(in) :: irow, jcol
50  character(len=*), intent(in) :: pkgname
51  character(len=3), intent(in) :: auxorbnd
52  real(dp), pointer, intent(in) :: bndelem
53  integer(I4B), intent(in) :: iprpak
54  character(len=*), optional :: text
55  ! -- local
56  character(len=LENPACKAGENAME) :: pkgnametemp
57  !
58  allocate (newtslink)
59  !
60  ! Store package name as all caps
61  pkgnametemp = pkgname
62  call upcase(pkgnametemp)
63  newtslink%PackageName = pkgnametemp
64  newtslink%AuxOrBnd = auxorbnd
65  newtslink%timeSeries => timeseries
66  newtslink%iRow = irow
67  newtslink%jCol = jcol
68  newtslink%BndElement => bndelem
69  newtslink%Iprpak = iprpak
70  !
71  if (present(text)) then
72  newtslink%Text = text
73  end if
74  !
75  ! -- Return
76  return
77  end subroutine constructtimeserieslink
78 
79  !> @brief Cast an unlimited polymorphic object as TimeSeriesLinkType
80  !<
81  function castastimeserieslinktype(obj) result(res)
82  implicit none
83  ! -- dummy
84  class(*), pointer, intent(inout) :: obj
85  ! -- return
86  type(timeserieslinktype), pointer :: res
87  !
88  res => null()
89  if (.not. associated(obj)) return
90  !
91  select type (obj)
92  type is (timeserieslinktype)
93  res => obj
94  class default
95  continue
96  end select
97  !
98  ! -- Return
99  return
100  end function castastimeserieslinktype
101 
102  !> @brief Get time series link from a list
103  !<
104  function gettimeserieslinkfromlist(list, indx) result(tsLink)
105  implicit none
106  ! -- dummy
107  type(listtype), intent(inout) :: list
108  integer(I4B), intent(in) :: indx
109  ! -- return
110  type(timeserieslinktype), pointer :: tslink
111  ! -- local
112  class(*), pointer :: obj
113  !
114  tslink => null()
115  obj => list%GetItem(indx)
116  tslink => castastimeserieslinktype(obj)
117  !
118  ! -- Return
119  return
120  end function gettimeserieslinkfromlist
121 
122  !> @brief Add time series link to a list
123  !<
124  subroutine addtimeserieslinktolist(list, tslink)
125  implicit none
126  ! -- dummy
127  type(listtype), intent(inout) :: list
128  type(timeserieslinktype), pointer, intent(inout) :: tslink
129  ! -- local
130  class(*), pointer :: obj
131  !
132  obj => tslink
133  call list%Add(obj)
134  !
135  ! -- Return
136  return
137  end subroutine addtimeserieslinktolist
138 
139 end module timeserieslinkmodule
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
integer(i4b), parameter lenboundname
maximum length of a bound name
Definition: Constants.f90:35
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
subroutine, public upcase(word)
Convert to upper case.
This module defines variable data types.
Definition: kind.f90:8
subroutine, public constructtimeserieslink(newTsLink, timeSeries, pkgName, auxOrBnd, bndElem, iRow, jCol, iprpak, text)
Construct time series link.
subroutine, public addtimeserieslinktolist(list, tslink)
Add time series link to a list.
type(timeserieslinktype) function, pointer, private castastimeserieslinktype(obj)
Cast an unlimited polymorphic object as TimeSeriesLinkType.
type(timeserieslinktype) function, pointer, public gettimeserieslinkfromlist(list, indx)
Get time series link from a list.
A generic heterogeneous doubly-linked list.
Definition: List.f90:10