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

Data Types

type  celldefntype
 Base grid cell definition. More...
 

Functions/Subroutines

subroutine, public create_defn (cellDefn)
 Create a new cell definition object. More...
 
integer(i4b) function, public get_iatop (ncpl, icu)
 Get the index corresponding to top elevation of a cell in the grid. This is -1 if the cell is in the top layer and 1 otherwise. More...
 
logical(lgp) function get_ispv180 (this, m)
 Return 180-degree indicator for a vertex. More...
 
real(dp) function get_botflow (this)
 Return the bottom flow. More...
 
real(dp) function get_topflow (this)
 Return the top flow. More...
 
real(dp) function get_distflow (this)
 Return the distributed flow. More...
 
real(dp) function get_faceflow (this, m)
 Return a face flow. More...
 

Function/Subroutine Documentation

◆ create_defn()

subroutine, public celldefnmodule::create_defn ( type(celldefntype), pointer  cellDefn)

Definition at line 40 of file CellDefn.f90.

41  type(CellDefnType), pointer :: cellDefn
42  allocate (celldefn)
43  ! Initially, allocate arrays to size for structured grid tracking method.
44  ! They can be (lazily) expanded as necessary for the unstructured method.
45  allocate (celldefn%ispv180(5))
46  allocate (celldefn%facenbr(7))
47  allocate (celldefn%faceflow(7))
Here is the caller graph for this function:

◆ get_botflow()

real(dp) function celldefnmodule::get_botflow ( class(celldefntype), intent(inout)  this)
private

Definition at line 72 of file CellDefn.f90.

73  class(CellDefnType), intent(inout) :: this
74  real(DP) :: botflow
75  botflow = this%faceflow(this%npolyverts + 2)

◆ get_distflow()

real(dp) function celldefnmodule::get_distflow ( class(celldefntype), intent(inout)  this)
private

Definition at line 86 of file CellDefn.f90.

87  class(CellDefnType), intent(inout) :: this
88  real(DP) :: distflow
89  distflow = this%distflow

◆ get_faceflow()

real(dp) function celldefnmodule::get_faceflow ( class(celldefntype), intent(inout)  this,
integer(i4b)  m 
)
private

Definition at line 93 of file CellDefn.f90.

94  class(CellDefnType), intent(inout) :: this
95  integer(I4B) :: m
96  real(DP) :: faceflow
97  faceflow = this%faceflow(m)

◆ get_iatop()

integer(i4b) function, public celldefnmodule::get_iatop ( integer(i4b), intent(in)  ncpl,
integer(i4b), intent(in)  icu 
)

Definition at line 52 of file CellDefn.f90.

53  integer(I4B), intent(in) :: ncpl, icu
54  integer(I4B) :: iatop
55 
56  if (icu .le. ncpl) then
57  iatop = -1
58  else
59  iatop = 1
60  end if
Here is the caller graph for this function:

◆ get_ispv180()

logical(lgp) function celldefnmodule::get_ispv180 ( class(celldefntype), intent(inout)  this,
integer(i4b)  m 
)
private

Definition at line 64 of file CellDefn.f90.

65  class(CellDefnType), intent(inout) :: this
66  integer(I4B) :: m
67  logical(LGP) :: ispv180
68  ispv180 = this%ispv180(m)

◆ get_topflow()

real(dp) function celldefnmodule::get_topflow ( class(celldefntype), intent(inout)  this)
private

Definition at line 79 of file CellDefn.f90.

80  class(CellDefnType), intent(inout) :: this
81  real(DP) :: topflow
82  topflow = this%faceflow(this%npolyverts + 3)