MODFLOW 6  version 6.6.0.dev0
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 204 of file CellRectQuad.f90.

205  ! -- dummy
206  class(CellRectQuadType), intent(inout) :: this
207  integer(I4B) :: i !< face index
208  logical(LGP) :: is_refined
209 
210  if (this%ipv4irv(2, i) .ne. 0) then
211  is_refined = .true.
212  else
213  is_refined = .false.
214  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 154 of file CellRectQuad.f90.

155  ! -- dummy
156  class(CellRectQuadType), intent(inout) :: this
157  ! -- local
158  integer(I4B) :: irv2, irv4, ipv1, ipv2, ipv4
159  integer(I4B), dimension(4) :: irvnxt = (/2, 3, 4, 1/)
160  real(DP) :: x1, y1, x2, y2, x4, y4, dx2, dy2, dx4, dy4
161 
162  ! -- Get rectangle vertex neighbors irv2 and irv4
163  irv2 = irvnxt(this%irvOrigin)
164  irv4 = irvnxt(irvnxt(irv2))
165 
166  ! -- Get model coordinates at irv1, irv2, and irv4
167  ipv1 = this%irectvert(this%irvOrigin)
168  x1 = this%defn%polyvert(1, ipv1)
169  y1 = this%defn%polyvert(2, ipv1)
170  ipv2 = this%irectvert(irv2)
171  x2 = this%defn%polyvert(1, ipv2)
172  y2 = this%defn%polyvert(2, ipv2)
173  ipv4 = this%irectvert(irv4)
174  x4 = this%defn%polyvert(1, ipv4)
175  y4 = this%defn%polyvert(2, ipv4)
176 
177  ! -- Compute rectangle dimensions
178  this%xOrigin = x1
179  this%yOrigin = y1
180  this%zOrigin = this%defn%bot
181  dx2 = x2 - this%xOrigin
182  dy2 = y2 - this%yOrigin
183  dx4 = x4 - this%xOrigin
184  dy4 = y4 - this%yOrigin
185  this%dx = dsqrt(dx4 * dx4 + dy4 * dy4)
186  this%dy = dsqrt(dx2 * dx2 + dy2 * dy2)
187  this%dz = this%defn%top - this%zOrigin
188 
189  ! -- Compute sine and cosine of rotation angle (angle between "southern"
190  ! -- rectangle side irv1-irv4 and the model x axis)
191  this%sinrot = dy4 / this%dx
192  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 196 of file CellRectQuad.f90.

197  class(CellRectQuadType), intent(inout) :: this
198  integer(I4B) :: iq, irv
199  real(DP) :: rectflow
200  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 109 of file CellRectQuad.f90.

110  ! -- dummy
111  class(CellRectQuadType), intent(inout) :: this
112  integer(I4B) :: irv1
113  ! -- local
114  integer(I4B) :: irv, irv2, irv4, ipv1, ipv2, ipv4
115  integer(I4B), dimension(4) :: irvnxt = (/2, 3, 4, 1/)
116  real(DP) :: x1, y1, x2, y2, x4, y4
117 
118  ! -- Find the "southwest" rectangle vertex by finding the vertex formed
119  ! -- either by (1) a rectangle edge over which x decreases (going
120  ! -- clockwise) followed by an edge over which x does not increase, or by
121  ! -- (2) a rectangle edge over which y does not decrease (again going
122  ! -- clockwise) followed by a rectangle edge over which y increases. In
123  ! -- the end, ipv1 is the index (1, 2, 3, or 4) of the southwest
124  ! -- rectangle vertex.
125  do irv = 1, 4
126  irv4 = irv
127  irv1 = irvnxt(irv4)
128  ipv4 = this%irectvert(irv4)
129  ipv1 = this%irectvert(irv1)
130  x4 = this%defn%polyvert(1, ipv4)
131  y4 = this%defn%polyvert(2, ipv4)
132  x1 = this%defn%polyvert(1, ipv1)
133  y1 = this%defn%polyvert(2, ipv1)
134  if (x1 .lt. x4) then
135  irv2 = irvnxt(irv1)
136  ipv2 = this%irectvert(irv2)
137  x2 = this%defn%polyvert(1, ipv2)
138  if (x2 .le. x1) return
139  else if (y1 .ge. y4) then
140  irv2 = irvnxt(irv1)
141  ipv2 = this%irectvert(irv2)
142  y2 = this%defn%polyvert(2, ipv2)
143  if (y2 .gt. y1) return
144  end if
145  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 
83  n = 0
84  do m = 1, this%defn%npolyverts
85  if (.not. this%defn%get_ispv180(m)) then
86  n = n + 1
87  this%irectvert(n) = m
88  this%ipv4irv(1, n) = m
89  this%rectflow(1, n) = this%defn%get_faceflow(m)
90  this%ipv4irv(2, n) = 0
91  this%rectflow(2, n) = dzero
92  else
93  if (n .ne. 0) then
94  this%ipv4irv(2, n) = m
95  this%rectflow(2, n) = this%defn%get_faceflow(m)
96  end if
97  end if
98  end do
99 
100  ! Wrap around for convenience
101  this%irectvert(5) = this%irectvert(1)