MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
CellRectQuad.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b, lgp
4  use cellmodule, only: celltype
6  use constantsmodule, only: dzero
7  implicit none
8 
9  private
10  public :: cellrectquadtype
11  public :: create_cell_rect_quad
12 
13  type, extends(celltype) :: cellrectquadtype
14  real(dp) :: dx !< dimension of cell in local x direction
15  real(dp) :: dy !< dimension of cell in local y direction
16  real(dp) :: dz !< dimension of cell in z direction
17 
18  real(dp) :: sinrot !< sine of rotation angle for local (x, y)
19  real(dp) :: cosrot !< cosine of rotation angle for local (x, y)
20 
21  integer(I4B) :: irvorigin !< origin rectangle vertex
22  real(dp) :: xorigin !< model x origin for local (x, y)
23  real(dp) :: yorigin !< model y origin for local (x, y)
24  real(dp) :: zorigin !< model z origin for local z
25 
26  real(dp) :: qextl1(4), qextl2(4), qintl(5) !< external and internal subcell flows for the cell
27  integer(I4B), allocatable :: irectvert(:) !< list of indices of the rectangle vertices
28  integer(I4B), allocatable :: ipv4irv(:, :) !< list of the polygon vertex indices that correspond to the rectangle vertex indices
29  real(dp), allocatable :: rectflow(:, :) !< flow(s) for each rectangle face
30  contains
31  procedure :: destroy => destroy_cell_rect_quad !< destructor for the cell
32  procedure :: init_from !< initializes the cell from an existing cell
33 
34  procedure :: load_rect_verts_flows ! loads list of indices of the rectangle vertices and face flows
35  procedure :: get_rect_ivert_sw ! gets index of southwest rectangle vertex
36  procedure :: get_rect_dim_rot ! gets rectangular dimensions and rotation
37  procedure :: get_rect_flow !< returns a rectangle face flow
38  procedure :: face_is_refined !< returns whether a rectangle face is refined
39  end type cellrectquadtype
40 
41 contains
42 
43  !> @brief Create a new rectangular-quad cell
44  subroutine create_cell_rect_quad(cell)
45  type(cellrectquadtype), pointer :: cell
46  allocate (cell)
47  call create_defn(cell%defn)
48  allocate (cell%irectvert(5))
49  allocate (cell%ipv4irv(2, 4))
50  allocate (cell%rectflow(2, 4))
51  allocate (cell%type)
52  cell%type = 'rectquad'
53  end subroutine create_cell_rect_quad
54 
55  !> @brief Destroy the rectangular-quad cell
56  subroutine destroy_cell_rect_quad(this)
57  class(cellrectquadtype), intent(inout) :: this
58  deallocate (this%defn)
59  deallocate (this%irectvert)
60  deallocate (this%type)
61  end subroutine destroy_cell_rect_quad
62 
63  !> @brief Initialize a rectangular-quad cell from cell definition
64  subroutine init_from(this, defn)
65  class(cellrectquadtype), intent(inout) :: this
66  type(celldefntype), pointer :: defn
67  this%defn => defn
68  call this%load_rect_verts_flows()
69  end subroutine init_from
70 
71  !> @brief Load local polygon vertex indices and rectangular
72  !> face flows
73  !!
74  !! Loads local polygon vertex indices of the four rectangle
75  !! vertices and face flows of a rectangular-quad cell.
76  !<
77  subroutine load_rect_verts_flows(this)
78  ! dummy
79  class(cellrectquadtype), intent(inout) :: this
80  ! local
81  integer(I4B) :: n, m
82  logical(LGP) :: istart180
83 
84  n = 0
85  istart180 = .false.
86  do m = 1, this%defn%npolyverts
87  if (.not. this%defn%get_ispv180(m)) then
88  n = n + 1
89  this%irectvert(n) = m
90  this%ipv4irv(1, n) = m
91  this%rectflow(1, n) = this%defn%get_faceflow(m)
92  this%ipv4irv(2, n) = 0
93  this%rectflow(2, n) = dzero
94  else
95  if (n .ne. 0) then
96  this%ipv4irv(2, n) = m
97  this%rectflow(2, n) = this%defn%get_faceflow(m)
98  else
99  istart180 = .true.
100  end if
101  end if
102  end do
103  if (istart180) then
104  this%ipv4irv(2, 4) = 1
105  this%rectflow(2, 4) = this%defn%get_faceflow(1)
106  end if
107 
108  ! Wrap around for convenience
109  this%irectvert(5) = this%irectvert(1)
110  end subroutine load_rect_verts_flows
111 
112  !> @brief Get index of SW rectangle vertex
113  !!
114  !! Return the index (1, 2, 3, or 4) of the southwest
115  !! rectangle vertex of a rectangular-quad cell
116  !<
117  function get_rect_ivert_sw(this) result(irv1)
118  ! dummy
119  class(cellrectquadtype), intent(inout) :: this
120  integer(I4B) :: irv1
121  ! local
122  integer(I4B) :: irv, irv2, irv4, ipv1, ipv2, ipv4
123  integer(I4B), dimension(4) :: irvnxt = (/2, 3, 4, 1/)
124  real(dp) :: x1, y1, x2, y2, x4, y4
125 
126  ! Find the "southwest" rectangle vertex by finding the vertex formed
127  ! either by (1) a rectangle edge over which x decreases (going
128  ! clockwise) followed by an edge over which x does not increase, or by
129  ! (2) a rectangle edge over which y does not decrease (again going
130  ! clockwise) followed by a rectangle edge over which y increases. In
131  ! the end, ipv1 is the index (1, 2, 3, or 4) of the southwest
132  ! rectangle vertex.
133  do irv = 1, 4
134  irv4 = irv
135  irv1 = irvnxt(irv4)
136  ipv4 = this%irectvert(irv4)
137  ipv1 = this%irectvert(irv1)
138  x4 = this%defn%polyvert(1, ipv4)
139  y4 = this%defn%polyvert(2, ipv4)
140  x1 = this%defn%polyvert(1, ipv1)
141  y1 = this%defn%polyvert(2, ipv1)
142  if (x1 .lt. x4) then
143  irv2 = irvnxt(irv1)
144  ipv2 = this%irectvert(irv2)
145  x2 = this%defn%polyvert(1, ipv2)
146  if (x2 .le. x1) return
147  else if (y1 .ge. y4) then
148  irv2 = irvnxt(irv1)
149  ipv2 = this%irectvert(irv2)
150  y2 = this%defn%polyvert(2, ipv2)
151  if (y2 .gt. y1) return
152  end if
153  end do
154  end function get_rect_ivert_sw
155 
156  !> @brief Get rectangular cell dimensions and rotation
157  !!
158  !! Compute rectangular dimensions and rotation of
159  !! the cell using the specified rectangle vertex
160  !! as the origin
161  !<
162  subroutine get_rect_dim_rot(this)
163  ! dummy
164  class(cellrectquadtype), intent(inout) :: this
165  ! local
166  integer(I4B) :: irv2, irv4, ipv1, ipv2, ipv4
167  integer(I4B), dimension(4) :: irvnxt = (/2, 3, 4, 1/)
168  real(DP) :: x1, y1, x2, y2, x4, y4, dx2, dy2, dx4, dy4
169 
170  ! Get rectangle vertex neighbors irv2 and irv4
171  irv2 = irvnxt(this%irvOrigin)
172  irv4 = irvnxt(irvnxt(irv2))
173 
174  ! Get model coordinates at irv1, irv2, and irv4
175  ipv1 = this%irectvert(this%irvOrigin)
176  x1 = this%defn%polyvert(1, ipv1)
177  y1 = this%defn%polyvert(2, ipv1)
178  ipv2 = this%irectvert(irv2)
179  x2 = this%defn%polyvert(1, ipv2)
180  y2 = this%defn%polyvert(2, ipv2)
181  ipv4 = this%irectvert(irv4)
182  x4 = this%defn%polyvert(1, ipv4)
183  y4 = this%defn%polyvert(2, ipv4)
184 
185  ! Compute rectangle dimensions
186  this%xOrigin = x1
187  this%yOrigin = y1
188  this%zOrigin = this%defn%bot
189  dx2 = x2 - this%xOrigin
190  dy2 = y2 - this%yOrigin
191  dx4 = x4 - this%xOrigin
192  dy4 = y4 - this%yOrigin
193  this%dx = dsqrt(dx4 * dx4 + dy4 * dy4)
194  this%dy = dsqrt(dx2 * dx2 + dy2 * dy2)
195  this%dz = this%defn%top - this%zOrigin
196 
197  ! Compute sine and cosine of rotation angle (angle between "southern"
198  ! rectangle side irv1-irv4 and the model x axis)
199  this%sinrot = dy4 / this%dx
200  this%cosrot = dx4 / this%dx
201  end subroutine get_rect_dim_rot
202 
203  !> @brief Return a rectangle face flow
204  function get_rect_flow(this, iq, irv) result(rectflow)
205  class(cellrectquadtype), intent(inout) :: this
206  integer(I4B) :: iq, irv
207  real(dp) :: rectflow
208  rectflow = this%rectflow(iq, irv)
209  end function get_rect_flow
210 
211  !> @brief Return whether a rectangle face is refined
212  function face_is_refined(this, i) result(is_refined)
213  ! dummy
214  class(cellrectquadtype), intent(inout) :: this
215  integer(I4B) :: i !< face index
216  logical(LGP) :: is_refined
217 
218  if (this%ipv4irv(2, i) .ne. 0) then
219  is_refined = .true.
220  else
221  is_refined = .false.
222  end if
223  end function face_is_refined
224 
225 end module cellrectquadmodule
subroutine, public create_defn(cellDefn)
Create a new cell definition object.
Definition: CellDefn.f90:58
subroutine, public create_cell_rect_quad(cell)
Create a new rectangular-quad cell.
subroutine init_from(this, defn)
Initialize a rectangular-quad cell from cell definition.
subroutine get_rect_dim_rot(this)
Get rectangular cell dimensions and rotation.
subroutine destroy_cell_rect_quad(this)
Destroy the rectangular-quad cell.
real(dp) function get_rect_flow(this, iq, irv)
Return a rectangle face flow.
subroutine load_rect_verts_flows(this)
Load local polygon vertex indices and rectangular face flows.
logical(lgp) function face_is_refined(this, i)
Return whether a rectangle face is refined.
integer(i4b) function get_rect_ivert_sw(this)
Get index of SW rectangle vertex.
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
Base grid cell definition.
Definition: CellDefn.f90:25
Base type for grid cells of a concrete type. Contains a cell-definition which is information shared b...
Definition: Cell.f90:12