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