MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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 function get_npolyverts (this)
 Return the number of polygon vertices. More...
 
logical function get_ispv180 (this, m)
 Return 180-degree indicator for a vertex. More...
 
double precision function get_botflow (this)
 Return the bottom flow. More...
 
double precision function get_topflow (this)
 Return the top flow. More...
 
double precision function get_distflow (this)
 Return the distributed flow. More...
 
double precision 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 41 of file CellDefn.f90.

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

◆ get_botflow()

double precision function celldefnmodule::get_botflow ( class(celldefntype), intent(inout)  this)
private

Definition at line 67 of file CellDefn.f90.

68  class(CellDefnType), intent(inout) :: this
69  double precision :: botflow
70  botflow = this%faceflow(this%npolyverts + 2)

◆ get_distflow()

double precision function celldefnmodule::get_distflow ( class(celldefntype), intent(inout)  this)
private

Definition at line 81 of file CellDefn.f90.

82  class(CellDefnType), intent(inout) :: this
83  double precision :: distflow
84  distflow = this%distflow

◆ get_faceflow()

double precision function celldefnmodule::get_faceflow ( class(celldefntype), intent(inout)  this,
integer  m 
)
private

Definition at line 88 of file CellDefn.f90.

89  class(CellDefnType), intent(inout) :: this
90  integer :: m
91  double precision :: faceflow
92  faceflow = this%faceflow(m)

◆ get_ispv180()

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

Definition at line 59 of file CellDefn.f90.

60  class(CellDefnType), intent(inout) :: this
61  integer :: m
62  logical :: ispv180
63  ispv180 = this%ispv180(m)

◆ get_npolyverts()

integer function celldefnmodule::get_npolyverts ( class(celldefntype), intent(inout)  this)
private

Definition at line 52 of file CellDefn.f90.

53  class(CellDefnType), intent(inout) :: this
54  integer :: npolyverts
55  npolyverts = this%npolyverts

◆ get_topflow()

double precision function celldefnmodule::get_topflow ( class(celldefntype), intent(inout)  this)
private

Definition at line 74 of file CellDefn.f90.

75  class(CellDefnType), intent(inout) :: this
76  double precision :: topflow
77  topflow = this%faceflow(this%npolyverts + 3)