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 43 of file CellDefn.f90.

44  type(CellDefnType), pointer :: cellDefn
45  allocate (celldefn)
46  ! Initially, allocate arrays to size for structured grid tracking method.
47  ! They can be (lazily) expanded as necessary for the unstructured method.
48  allocate (celldefn%ispv180(5))
49  allocate (celldefn%facenbr(7))
50  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 75 of file CellDefn.f90.

76  class(CellDefnType), intent(inout) :: this
77  real(DP) :: botflow
78  botflow = this%faceflow(this%npolyverts + 2)

◆ get_distflow()

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

Definition at line 89 of file CellDefn.f90.

90  class(CellDefnType), intent(inout) :: this
91  real(DP) :: distflow
92  distflow = this%distflow

◆ get_faceflow()

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

Definition at line 96 of file CellDefn.f90.

97  class(CellDefnType), intent(inout) :: this
98  integer(I4B) :: m
99  real(DP) :: faceflow
100  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 55 of file CellDefn.f90.

56  integer(I4B), intent(in) :: ncpl, icu
57  integer(I4B) :: iatop
58 
59  if (icu .le. ncpl) then
60  iatop = -1
61  else
62  iatop = 1
63  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 67 of file CellDefn.f90.

68  class(CellDefnType), intent(inout) :: this
69  integer(I4B) :: m
70  logical(LGP) :: ispv180
71  ispv180 = this%ispv180(m)

◆ get_topflow()

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

Definition at line 82 of file CellDefn.f90.

83  class(CellDefnType), intent(inout) :: this
84  real(DP) :: topflow
85  topflow = this%faceflow(this%npolyverts + 3)