MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
CellDefn.f90
Go to the documentation of this file.
2  use kindmodule, only: dp, i4b, lgp
3  implicit none
4 
5  private
6  public :: celldefntype
7  public :: create_defn, get_iatop
9 
10  !> @brief Cell saturation status.
11  !!
12  !! Status 0 indicates that the cell is dry.
13  !! Status 1 indicates that the cell contains a water table.
14  !! Since a water table might exist exactly at the top face,
15  !! status 1 includes that case. Status 2 is reserved for a
16  !! fully saturated cell with no water table.
17  !<
18  enum, bind(C)
19  enumerator :: saturation_dry = 0
20  enumerator :: saturation_watertable = 1
21  enumerator :: saturation_saturated = 2
22  end enum
23 
24  !> @brief Base grid cell definition.
26  private
27  integer(I4B), public :: icell !< index of cell in source grid
28  logical(LGP), public :: can_be_rect !< whether cell is representable as a rectangular cell
29  logical(LGP), public :: can_be_quad !< whether cell is representable as a rectangular quad cell
30  integer(I4B), public :: npolyverts !< number of vertices for cell polygon
31  real(dp), public :: porosity !< cell porosity
32  real(dp), public :: retfactor !< cell retardation factor
33  integer(I4B), public :: ilay !< layer number
34  integer(I4B), public :: izone !< cell zone number
35  integer(I4B), public :: iweaksink !< weak sink indicator
36  integer(I4B), public :: inoexitface !< no exit face indicator
37  integer(I4B), public :: iatop !< index of cell top in grid's top/bot arrays (<0 => top array)
38  integer(I4B), public :: isatstat !< saturation status. 0 = dry, 1 = contains a water table, 2 = fully saturated w/ no water table
39  real(dp), public :: top, bot !< top and bottom elevations of cell
40  real(dp), public :: sat !< cell saturation
41  real(dp), allocatable, public :: polyvert(:, :) !< vertices for cell polygon
42  logical(LGP), allocatable, public :: ispv180(:) !< indicator of 180-degree vertices (.true. = 180-degree angle at vertex)
43  integer(I4B), allocatable, public :: facenbr(:) !< neighbors that correspond to faces(/vertices)
44  real(dp), allocatable, public :: faceflow(:) !< flows that correspond to faces(/vertices)
45  real(dp), public :: distflow !< net distributed flow into cell
46  contains
47  procedure, public :: get_ispv180 !< returns 180-degree indicator for a vertex
48  procedure, public :: get_botflow !< returns bottom flow
49  procedure, public :: get_topflow !< returns top flow
50  procedure, public :: get_distflow !< returns distributed flow
51  procedure, public :: get_faceflow !< returns a face flow
52  end type celldefntype
53 
54 contains
55 
56  !> @brief Create a new cell definition object
57  subroutine create_defn(cellDefn)
58  type(celldefntype), pointer :: celldefn
59  allocate (celldefn)
60  ! Initially, allocate arrays to size for structured grid tracking method.
61  ! They can be (lazily) expanded as necessary for the unstructured method.
62  allocate (celldefn%ispv180(5))
63  allocate (celldefn%facenbr(7))
64  allocate (celldefn%faceflow(7))
65  end subroutine create_defn
66 
67  !> @brief Get the index corresponding to top elevation of a cell in the grid.
68  !! This is -1 if the cell is in the top layer and 1 otherwise.
69  function get_iatop(ncpl, icu) result(iatop)
70  integer(I4B), intent(in) :: ncpl, icu
71  integer(I4B) :: iatop
72 
73  if (icu .le. ncpl) then
74  iatop = -1
75  else
76  iatop = 1
77  end if
78  end function get_iatop
79 
80  !> @brief Return 180-degree indicator for a vertex
81  function get_ispv180(this, m) result(ispv180)
82  class(celldefntype), intent(inout) :: this
83  integer(I4B) :: m
84  logical(LGP) :: ispv180
85  ispv180 = this%ispv180(m)
86  end function get_ispv180
87 
88  !> @brief Return the bottom flow
89  function get_botflow(this) result(botflow)
90  class(celldefntype), intent(inout) :: this
91  real(dp) :: botflow
92  botflow = this%faceflow(this%npolyverts + 2)
93  end function get_botflow
94 
95  !> @brief Return the top flow
96  function get_topflow(this) result(topflow)
97  class(celldefntype), intent(inout) :: this
98  real(dp) :: topflow
99  topflow = this%faceflow(this%npolyverts + 3)
100  end function get_topflow
101 
102  !> @brief Return the distributed flow
103  function get_distflow(this) result(distflow)
104  class(celldefntype), intent(inout) :: this
105  real(dp) :: distflow
106  distflow = this%distflow
107  end function get_distflow
108 
109  !> @brief Return a face flow
110  function get_faceflow(this, m) result(faceflow)
111  class(celldefntype), intent(inout) :: this
112  integer(I4B) :: m
113  real(dp) :: faceflow
114  faceflow = this%faceflow(m)
115  end function get_faceflow
116 
117 end module celldefnmodule
real(dp) function get_faceflow(this, m)
Return a face flow.
Definition: CellDefn.f90:111
@, public saturation_saturated
Definition: CellDefn.f90:21
real(dp) function get_topflow(this)
Return the top flow.
Definition: CellDefn.f90:97
subroutine, public create_defn(cellDefn)
Create a new cell definition object.
Definition: CellDefn.f90:58
@, public saturation_dry
Definition: CellDefn.f90:19
real(dp) function get_botflow(this)
Return the bottom flow.
Definition: CellDefn.f90:90
real(dp) function get_distflow(this)
Return the distributed flow.
Definition: CellDefn.f90:104
logical(lgp) function get_ispv180(this, m)
Return 180-degree indicator for a vertex.
Definition: CellDefn.f90:82
@, public saturation_watertable
Definition: CellDefn.f90:20
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 ...
Definition: CellDefn.f90:70
This module defines variable data types.
Definition: kind.f90:8
Base grid cell definition.
Definition: CellDefn.f90:25