MODFLOW 6  version 6.7.0.dev0
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
8 
9  !> @brief Base grid cell definition.
11  private
12  integer(I4B), public :: icell !< index of cell in source grid
13  logical(LGP), public :: can_be_rect !< whether cell is representable as a rectangular cell
14  logical(LGP), public :: can_be_quad !< whether cell is representable as a rectangular quad cell
15  integer(I4B), public :: npolyverts !< number of vertices for cell polygon
16  real(dp), public :: porosity !< cell porosity
17  real(dp), public :: retfactor !< cell retardation factor
18  integer(I4B), public :: ilay !< layer number
19  integer(I4B), public :: izone !< cell zone number
20  integer(I4B), public :: iweaksink !< weak sink indicator
21  integer(I4B), public :: inoexitface !< no exit face indicator
22  integer(I4B), public :: iatop !< index of cell top in grid's top/bot arrays (<0 => top array)
23  real(dp), public :: top, bot !< top and bottom elevations of cell
24  real(dp), public :: sat !< cell saturation
25  real(dp), allocatable, public :: polyvert(:, :) !< vertices for cell polygon
26  logical(LGP), allocatable, public :: ispv180(:) !< indicator of 180-degree vertices (.true. = 180-degree angle at vertex)
27  integer(I4B), allocatable, public :: facenbr(:) !< neighbors that correspond to faces(/vertices)
28  real(dp), allocatable, public :: faceflow(:) !< flows that correspond to faces(/vertices)
29  real(dp), public :: distflow !< net distributed flow into cell
30  contains
31  procedure, public :: get_ispv180 !< returns 180-degree indicator for a vertex
32  procedure, public :: get_botflow !< returns bottom flow
33  procedure, public :: get_topflow !< returns top flow
34  procedure, public :: get_distflow !< returns distributed flow
35  procedure, public :: get_faceflow !< returns a face flow
36  end type celldefntype
37 
38 contains
39 
40  !> @brief Create a new cell definition object
41  subroutine create_defn(cellDefn)
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))
49  end subroutine create_defn
50 
51  !> @brief Get the index corresponding to top elevation of a cell in the grid.
52  !! This is -1 if the cell is in the top layer and 1 otherwise.
53  function get_iatop(ncpl, icu) result(iatop)
54  integer(I4B), intent(in) :: ncpl, icu
55  integer(I4B) :: iatop
56 
57  if (icu .le. ncpl) then
58  iatop = -1
59  else
60  iatop = 1
61  end if
62  end function get_iatop
63 
64  !> @brief Return 180-degree indicator for a vertex
65  function get_ispv180(this, m) result(ispv180)
66  class(celldefntype), intent(inout) :: this
67  integer(I4B) :: m
68  logical(LGP) :: ispv180
69  ispv180 = this%ispv180(m)
70  end function get_ispv180
71 
72  !> @brief Return the bottom flow
73  function get_botflow(this) result(botflow)
74  class(celldefntype), intent(inout) :: this
75  real(dp) :: botflow
76  botflow = this%faceflow(this%npolyverts + 2)
77  end function get_botflow
78 
79  !> @brief Return the top flow
80  function get_topflow(this) result(topflow)
81  class(celldefntype), intent(inout) :: this
82  real(dp) :: topflow
83  topflow = this%faceflow(this%npolyverts + 3)
84  end function get_topflow
85 
86  !> @brief Return the distributed flow
87  function get_distflow(this) result(distflow)
88  class(celldefntype), intent(inout) :: this
89  real(dp) :: distflow
90  distflow = this%distflow
91  end function get_distflow
92 
93  !> @brief Return a face flow
94  function get_faceflow(this, m) result(faceflow)
95  class(celldefntype), intent(inout) :: this
96  integer(I4B) :: m
97  real(dp) :: faceflow
98  faceflow = this%faceflow(m)
99  end function get_faceflow
100 
101 end module celldefnmodule
real(dp) function get_faceflow(this, m)
Return a face flow.
Definition: CellDefn.f90:95
real(dp) function get_topflow(this)
Return the top flow.
Definition: CellDefn.f90:81
subroutine, public create_defn(cellDefn)
Create a new cell definition object.
Definition: CellDefn.f90:42
real(dp) function get_botflow(this)
Return the bottom flow.
Definition: CellDefn.f90:74
real(dp) function get_distflow(this)
Return the distributed flow.
Definition: CellDefn.f90:88
logical(lgp) function get_ispv180(this, m)
Return 180-degree indicator for a vertex.
Definition: CellDefn.f90:66
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:54
This module defines variable data types.
Definition: kind.f90:8
Base grid cell definition.
Definition: CellDefn.f90:10