MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
timearraymodule Module Reference

Data Types

type  timearraytype
 

Functions/Subroutines

subroutine, public constructtimearray (newTa, modelname)
 Construct time array. More...
 
type(timearraytype) function, pointer, public castastimearraytype (obj)
 Cast an unlimited polymorphic object as TimeArrayType. More...
 
subroutine, public addtimearraytolist (list, timearray)
 Add a time array to a to list. More...
 
type(timearraytype) function, pointer, public gettimearrayfromlist (list, indx)
 Retrieve a time array from a list. More...
 
subroutine ta_da (this)
 Deallocate memory. More...
 

Function/Subroutine Documentation

◆ addtimearraytolist()

subroutine, public timearraymodule::addtimearraytolist ( type(listtype), intent(inout)  list,
type(timearraytype), intent(inout), pointer  timearray 
)

Definition at line 91 of file TimeArray.f90.

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)
Here is the caller graph for this function:

◆ castastimearraytype()

type(timearraytype) function, pointer, public timearraymodule::castastimearraytype ( class(*), intent(inout), pointer  obj)

Definition at line 74 of file TimeArray.f90.

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
Here is the caller graph for this function:

◆ constructtimearray()

subroutine, public timearraymodule::constructtimearray ( type(timearraytype), intent(out), pointer  newTa,
character(len=*), intent(in)  modelname 
)

Allocate and assign members of a new TimeArrayType object. Allocate space for the array so that this subroutine can be called repeatedly with the same array (but with different contents).

Definition at line 35 of file TimeArray.f90.

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))
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:27
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
Here is the call graph for this function:
Here is the caller graph for this function:

◆ gettimearrayfromlist()

type(timearraytype) function, pointer, public timearraymodule::gettimearrayfromlist ( type(listtype), intent(inout)  list,
integer(i4b), intent(in)  indx 
)

Definition at line 104 of file TimeArray.f90.

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)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ta_da()

subroutine timearraymodule::ta_da ( class(timearraytype this)
private

Definition at line 119 of file TimeArray.f90.

120  ! -- dummy
121  class(TimeArrayType) :: this
122  !
123  deallocate (this%taArray)
124  this%taArray => null()