MODFLOW 6  version 6.7.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
40  ! -- dummy
41  type(timearraytype), pointer, intent(out) :: newta
42  character(len=*), intent(in) :: modelname
43  ! -- local
44  integer(I4B), dimension(:), contiguous, &
45  pointer :: mshape
46  character(len=LENMEMPATH) :: mempath
47  integer(I4B) :: isize
48  !
49  ! -- initialize
50  nullify (mshape)
51  !
52  ! -- create mempath
53  mempath = create_mem_path(component=modelname, subcomponent='DIS')
54  !
55  ! -- set mshape pointer
56  call mem_setptr(mshape, 'MSHAPE', mempath)
57  !
58  ! Get dimensions for supported discretization type
59  if (size(mshape) == 2) then
60  isize = mshape(2)
61  else if (size(mshape) == 3) then
62  isize = mshape(2) * mshape(3)
63  else
64  errmsg = 'Time array series is not supported for discretization type'
65  call store_error(errmsg, terminate=.true.)
66  end if
67  !
68  allocate (newta)
69  allocate (newta%taArray(isize))
70  end subroutine constructtimearray
71 
72  !> @brief Cast an unlimited polymorphic object as TimeArrayType
73  !<
74  function castastimearraytype(obj) result(res)
75  ! -- dummy
76  class(*), pointer, intent(inout) :: obj
77  ! -- return
78  type(timearraytype), pointer :: res
79  !
80  res => null()
81  if (.not. associated(obj)) return
82  !
83  select type (obj)
84  type is (timearraytype)
85  res => obj
86  end select
87  end function castastimearraytype
88 
89  !> @brief Add a time array to a to list
90  !<
91  subroutine addtimearraytolist(list, timearray)
92  ! -- dummy
93  type(listtype), intent(inout) :: list
94  type(timearraytype), pointer, intent(inout) :: timearray
95  ! -- local
96  class(*), pointer :: obj
97  !
98  obj => timearray
99  call list%Add(obj)
100  end subroutine addtimearraytolist
101 
102  !> @brief Retrieve a time array from a list
103  !<
104  function gettimearrayfromlist(list, indx) result(res)
105  ! -- dummy
106  type(listtype), intent(inout) :: list
107  integer(I4B), intent(in) :: indx
108  ! -- return
109  type(timearraytype), pointer :: res
110  ! -- local
111  class(*), pointer :: obj
112  !
113  obj => list%GetItem(indx)
114  res => castastimearraytype(obj)
115  end function gettimearrayfromlist
116 
117  !> @brief Deallocate memory
118  !<
119  subroutine ta_da(this)
120  ! -- dummy
121  class(timearraytype) :: this
122  !
123  deallocate (this%taArray)
124  this%taArray => null()
125  end subroutine ta_da
126 
127 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
type(timearraytype) function, pointer, public castastimearraytype(obj)
Cast an unlimited polymorphic object as TimeArrayType.
Definition: TimeArray.f90:75
type(timearraytype) function, pointer, public gettimearrayfromlist(list, indx)
Retrieve a time array from a list.
Definition: TimeArray.f90:105
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:92
subroutine ta_da(this)
Deallocate memory.
Definition: TimeArray.f90:120
A generic heterogeneous doubly-linked list.
Definition: List.f90:14