MODFLOW 6  version 6.7.0.dev3
USGS Modular Hydrologic Model
cellrectquadmodule Module Reference

Data Types

type  cellrectquadtype
 

Functions/Subroutines

subroutine, public create_cell_rect_quad (cell)
 Create a new rectangular-quad cell. More...
 
subroutine destroy_cell_rect_quad (this)
 Destroy the rectangular-quad cell. More...
 
subroutine init_from (this, defn)
 Initialize a rectangular-quad cell from cell definition. More...
 
subroutine load_rect_verts_flows (this)
 Load local polygon vertex indices and rectangular face flows. More...
 
integer(i4b) function get_rect_ivert_sw (this)
 Get index of SW rectangle vertex. More...
 
subroutine get_rect_dim_rot (this)
 Get rectangular cell dimensions and rotation. More...
 
real(dp) function get_rect_flow (this, iq, irv)
 Return a rectangle face flow. More...
 
logical(lgp) function face_is_refined (this, i)
 Return whether a rectangle face is refined. More...
 

Function/Subroutine Documentation

◆ create_cell_rect_quad()

subroutine, public cellrectquadmodule::create_cell_rect_quad ( type(cellrectquadtype), pointer  cell)

Definition at line 44 of file CellRectQuad.f90.

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'
Here is the call graph for this function:
Here is the caller graph for this function:

◆ destroy_cell_rect_quad()

subroutine cellrectquadmodule::destroy_cell_rect_quad ( class(cellrectquadtype), intent(inout)  this)
private

Definition at line 56 of file CellRectQuad.f90.

57  class(CellRectQuadType), intent(inout) :: this
58  deallocate (this%defn)
59  deallocate (this%irectvert)
60  deallocate (this%type)

◆ face_is_refined()

logical(lgp) function cellrectquadmodule::face_is_refined ( class(cellrectquadtype), intent(inout)  this,
integer(i4b)  i 
)
private
Parameters
iface index

Definition at line 212 of file CellRectQuad.f90.

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

◆ get_rect_dim_rot()

subroutine cellrectquadmodule::get_rect_dim_rot ( class(cellrectquadtype), intent(inout)  this)
private

Compute rectangular dimensions and rotation of the cell using the specified rectangle vertex as the origin

Definition at line 162 of file CellRectQuad.f90.

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

◆ get_rect_flow()

real(dp) function cellrectquadmodule::get_rect_flow ( class(cellrectquadtype), intent(inout)  this,
integer(i4b)  iq,
integer(i4b)  irv 
)
private

Definition at line 204 of file CellRectQuad.f90.

205  class(CellRectQuadType), intent(inout) :: this
206  integer(I4B) :: iq, irv
207  real(DP) :: rectflow
208  rectflow = this%rectflow(iq, irv)

◆ get_rect_ivert_sw()

integer(i4b) function cellrectquadmodule::get_rect_ivert_sw ( class(cellrectquadtype), intent(inout)  this)
private

Return the index (1, 2, 3, or 4) of the southwest rectangle vertex of a rectangular-quad cell

Definition at line 117 of file CellRectQuad.f90.

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

◆ init_from()

subroutine cellrectquadmodule::init_from ( class(cellrectquadtype), intent(inout)  this,
type(celldefntype), pointer  defn 
)
private

Definition at line 64 of file CellRectQuad.f90.

65  class(CellRectQuadType), intent(inout) :: this
66  type(CellDefnType), pointer :: defn
67  this%defn => defn
68  call this%load_rect_verts_flows()

◆ load_rect_verts_flows()

subroutine cellrectquadmodule::load_rect_verts_flows ( class(cellrectquadtype), intent(inout)  this)
private

Loads local polygon vertex indices of the four rectangle vertices and face flows of a rectangular-quad cell.

Definition at line 77 of file CellRectQuad.f90.

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)