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