MODFLOW 6  version 6.8.0.dev0
USGS Modular Hydrologic Model
TimeArray.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b
4  use listmodule, only: listtype
5  use simvariablesmodule, only: errmsg
6  use simmodule, only: store_error
7 
8  implicit none
9  private
10  public :: timearraytype, constructtimearray, &
13 
14  type :: timearraytype
15  ! -- Public members
16  real(dp), public :: tatime
17  real(dp), dimension(:), pointer, contiguous, public :: taarray => null()
18 
19  contains
20 
21  ! -- Public procedures
22  ! -- When gfortran adds support for finalization, the
23  ! following declaration could be: final :: finalize
24  procedure, public :: da => ta_da
25  end type timearraytype
26 
27 contains
28 
29  !> @brief Construct time array
30  !!
31  !! Allocate and assign members of a new TimeArrayType object. Allocate space
32  !! for the array so that this subroutine can be called repeatedly with the
33  !! same array (but with different contents).
34  !<
35  subroutine constructtimearray(newTa, modelname)
36  ! -- modules
37  use constantsmodule, only: lenmempath
41  ! -- dummy
42  type(timearraytype), pointer, intent(out) :: newta
43  character(len=*), intent(in) :: modelname
44  ! -- local
45  integer(I4B), dimension(:), contiguous, &
46  pointer :: mshape
47  character(len=LENMEMPATH) :: mempath
48  integer(I4B) :: isize
49  !
50  ! -- initialize
51  nullify (mshape)
52  !
53  ! -- create mempath to model input context where MODEL_SHAPE is
54  ! -- stored during static DIS loading
55  mempath = create_mem_path(component=modelname, context=idm_context)
56  !
57  ! -- set mshape pointer
58  call mem_setptr(mshape, 'MODEL_SHAPE', mempath)
59  !
60  ! Get dimensions for supported discretization type
61  if (size(mshape) == 2) then
62  isize = mshape(2)
63  else if (size(mshape) == 3) then
64  isize = mshape(2) * mshape(3)
65  else
66  errmsg = 'Time array series is not supported for discretization type'
67  call store_error(errmsg, terminate=.true.)
68  end if
69  !
70  allocate (newta)
71  allocate (newta%taArray(isize))
72  end subroutine constructtimearray
73 
74  !> @brief Cast an unlimited polymorphic object as TimeArrayType
75  !<
76  function castastimearraytype(obj) result(res)
77  ! -- dummy
78  class(*), pointer, intent(inout) :: obj
79  ! -- return
80  type(timearraytype), pointer :: res
81  !
82  res => null()
83  if (.not. associated(obj)) return
84  !
85  select type (obj)
86  type is (timearraytype)
87  res => obj
88  end select
89  end function castastimearraytype
90 
91  !> @brief Add a time array to a to list
92  !<
93  subroutine addtimearraytolist(list, timearray)
94  ! -- dummy
95  type(listtype), intent(inout) :: list
96  type(timearraytype), pointer, intent(inout) :: timearray
97  ! -- local
98  class(*), pointer :: obj
99  !
100  obj => timearray
101  call list%Add(obj)
102  end subroutine addtimearraytolist
103 
104  !> @brief Retrieve a time array from a list
105  !<
106  function gettimearrayfromlist(list, indx) result(res)
107  ! -- dummy
108  type(listtype), intent(inout) :: list
109  integer(I4B), intent(in) :: indx
110  ! -- return
111  type(timearraytype), pointer :: res
112  ! -- local
113  class(*), pointer :: obj
114  !
115  obj => list%GetItem(indx)
116  res => castastimearraytype(obj)
117  end function gettimearrayfromlist
118 
119  !> @brief Deallocate memory
120  !<
121  subroutine ta_da(this)
122  ! -- dummy
123  class(timearraytype) :: this
124  !
125  deallocate (this%taArray)
126  this%taArray => null()
127  end subroutine ta_da
128 
129 end module timearraymodule
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27
This module defines variable data types.
Definition: kind.f90:8
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains simulation methods.
Definition: Sim.f90:10
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=maxcharlen) errmsg
error message string
character(len=linelength) idm_context
type(timearraytype) function, pointer, public castastimearraytype(obj)
Cast an unlimited polymorphic object as TimeArrayType.
Definition: TimeArray.f90:77
type(timearraytype) function, pointer, public gettimearrayfromlist(list, indx)
Retrieve a time array from a list.
Definition: TimeArray.f90:107
subroutine, public constructtimearray(newTa, modelname)
Construct time array.
Definition: TimeArray.f90:36
subroutine, public addtimearraytolist(list, timearray)
Add a time array to a to list.
Definition: TimeArray.f90:94
subroutine ta_da(this)
Deallocate memory.
Definition: TimeArray.f90:122
A generic heterogeneous doubly-linked list.
Definition: List.f90:14