MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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
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_npolyverts !< returns the number of polygon vertices
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 Return the number of polygon vertices
52  function get_npolyverts(this) result(npolyverts)
53  class(celldefntype), intent(inout) :: this
54  integer :: npolyverts
55  npolyverts = this%npolyverts
56  end function get_npolyverts
57 
58  !> @brief Return 180-degree indicator for a vertex
59  function get_ispv180(this, m) result(ispv180)
60  class(celldefntype), intent(inout) :: this
61  integer :: m
62  logical :: ispv180
63  ispv180 = this%ispv180(m)
64  end function get_ispv180
65 
66  !> @brief Return the bottom flow
67  function get_botflow(this) result(botflow)
68  class(celldefntype), intent(inout) :: this
69  double precision :: botflow
70  botflow = this%faceflow(this%npolyverts + 2)
71  end function get_botflow
72 
73  !> @brief Return the top flow
74  function get_topflow(this) result(topflow)
75  class(celldefntype), intent(inout) :: this
76  double precision :: topflow
77  topflow = this%faceflow(this%npolyverts + 3)
78  end function get_topflow
79 
80  !> @brief Return the distributed flow
81  function get_distflow(this) result(distflow)
82  class(celldefntype), intent(inout) :: this
83  double precision :: distflow
84  distflow = this%distflow
85  end function get_distflow
86 
87  !> @brief Return a face flow
88  function get_faceflow(this, m) result(faceflow)
89  class(celldefntype), intent(inout) :: this
90  integer :: m
91  double precision :: faceflow
92  faceflow = this%faceflow(m)
93  end function get_faceflow
94 
95 end module celldefnmodule
logical function get_ispv180(this, m)
Return 180-degree indicator for a vertex.
Definition: CellDefn.f90:60
subroutine, public create_defn(cellDefn)
Create a new cell definition object.
Definition: CellDefn.f90:42
double precision function get_topflow(this)
Return the top flow.
Definition: CellDefn.f90:75
double precision function get_distflow(this)
Return the distributed flow.
Definition: CellDefn.f90:82
double precision function get_botflow(this)
Return the bottom flow.
Definition: CellDefn.f90:68
integer function get_npolyverts(this)
Return the number of polygon vertices.
Definition: CellDefn.f90:53
double precision function get_faceflow(this, m)
Return a face flow.
Definition: CellDefn.f90:89
This module defines variable data types.
Definition: kind.f90:8
Base grid cell definition.
Definition: CellDefn.f90:10