MODFLOW 6  version 6.8.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 93 of file TimeArray.f90.

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

◆ castastimearraytype()

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

Definition at line 76 of file TimeArray.f90.

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
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
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))
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
This module contains simulation variables.
Definition: SimVariables.f90:9
character(len=linelength) idm_context
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 106 of file TimeArray.f90.

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)
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 121 of file TimeArray.f90.

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