MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
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 another cell. More...
 
subroutine load_irectvert (this)
 Load local polygon vertex indices. More...
 
integer function get_irectvertsw (this)
 Get index of SW rectangle vertex. More...
 
subroutine get_rectdimensionsrotation (this, irv1, xOrigin, yOrigin, zOrigin, dx, dy, dz, sinrot, cosrot)
 Get rectangular cell dimensions and rotation. More...
 
double precision function get_rectflow (this, iq, irv)
 Return a rectangle face flow. More...
 
logical 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 43 of file CellRectQuad.f90.

44  type(CellRectQuadType), pointer :: cell
45  allocate (cell)
46  call create_defn(cell%defn)
47  allocate (cell%irectvert(5))
48  allocate (cell%ipv4irv(2, 4))
49  allocate (cell%rectflow(2, 4))
50  allocate (cell%type)
51  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 55 of file CellRectQuad.f90.

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

◆ face_is_refined()

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

Definition at line 207 of file CellRectQuad.f90.

208  ! -- dummy
209  class(CellRectQuadType), intent(inout) :: this
210  integer :: i !< face index
211  logical :: is_refined
212 
213  if (this%ipv4irv(2, i) .ne. 0) then
214  is_refined = .true.
215  else
216  is_refined = .false.
217  end if

◆ get_irectvertsw()

integer function cellrectquadmodule::get_irectvertsw ( 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 :: irv1
113  ! -- local
114  integer :: irv, irv2, irv4, ipv1, ipv2, ipv4
115  integer, dimension(4) :: irvnxt = (/2, 3, 4, 1/) ! kluge???
116  double precision :: 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

◆ get_rectdimensionsrotation()

subroutine cellrectquadmodule::get_rectdimensionsrotation ( class(cellrectquadtype), intent(inout)  this,
integer  irv1,
double precision  xOrigin,
double precision  yOrigin,
double precision  zOrigin,
double precision  dx,
double precision  dy,
double precision  dz,
double precision  sinrot,
double precision  cosrot 
)
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.

156  ! -- dummy
157  class(CellRectQuadType), intent(inout) :: this
158  integer :: irv1
159  double precision :: xOrigin, yOrigin, zOrigin, dx, dy, dz, sinrot, cosrot
160  ! -- local
161  integer :: irv2, irv4, ipv1, ipv2, ipv4
162  integer, dimension(4) :: irvnxt = (/2, 3, 4, 1/) ! kluge???
163  double precision :: x1, y1, x2, y2, x4, y4, dx2, dy2, dx4, dy4
164 
165  ! -- Get rectangle vertex neighbors irv2 and irv4
166  irv2 = irvnxt(irv1)
167  irv4 = irvnxt(irvnxt(irv2)) ! kluge
168 
169  ! -- Get model coordinates at irv1, irv2, and irv4
170  ipv1 = this%irectvert(irv1)
171  x1 = this%defn%polyvert(1, ipv1)
172  y1 = this%defn%polyvert(2, ipv1)
173  ipv2 = this%irectvert(irv2)
174  x2 = this%defn%polyvert(1, ipv2)
175  y2 = this%defn%polyvert(2, ipv2)
176  ipv4 = this%irectvert(irv4)
177  x4 = this%defn%polyvert(1, ipv4)
178  y4 = this%defn%polyvert(2, ipv4)
179 
180  ! -- Compute rectangle dimensions
181  xorigin = x1
182  yorigin = y1
183  zorigin = this%defn%bot
184  dx2 = x2 - xorigin
185  dy2 = y2 - yorigin
186  dx4 = x4 - xorigin
187  dy4 = y4 - yorigin
188  dx = dsqrt(dx4 * dx4 + dy4 * dy4)
189  dy = dsqrt(dx2 * dx2 + dy2 * dy2)
190  dz = this%defn%top - zorigin ! kluge note: need to account for partial saturation
191 
192  ! -- Compute sine and cosine of rotation angle (angle between "southern"
193  ! -- rectangle side irv1-irv4 and the model x axis)
194  sinrot = dy4 / dx
195  cosrot = dx4 / dx

◆ get_rectflow()

double precision function cellrectquadmodule::get_rectflow ( class(cellrectquadtype), intent(inout)  this,
integer  iq,
integer  irv 
)
private

Definition at line 199 of file CellRectQuad.f90.

200  class(CellRectQuadType), intent(inout) :: this
201  integer :: iq, irv
202  double precision :: rectflow
203  rectflow = this%rectflow(iq, irv)

◆ init_from()

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

Definition at line 63 of file CellRectQuad.f90.

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

◆ load_irectvert()

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

Loads local polygon vertex indices of the four rectangle vertices of a rectangular-quad cell. Todo: rename?

Definition at line 75 of file CellRectQuad.f90.

76  ! -- dummy
77  class(CellRectQuadType), intent(inout) :: this
78  ! -- local
79  integer :: npolyverts, n, m
80 
81  npolyverts = this%defn%get_npolyverts()
82 
83  n = 0
84  do m = 1, 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) = 0d0
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)