MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
cellutilmodule Module Reference

Functions/Subroutines

subroutine, public cell_poly_to_rect (poly, rect)
 Convert CellPoly representation to CellRect. Assumes the conversion is possible. More...
 
subroutine, public cell_poly_to_quad (poly, quad)
 Convert CellPoly representation to CellRectQuad. Assumes the conversion is possible. More...
 

Function/Subroutine Documentation

◆ cell_poly_to_quad()

subroutine, public cellutilmodule::cell_poly_to_quad ( type(cellpolytype), intent(in), pointer  poly,
type(cellrectquadtype), intent(inout), pointer  quad 
)

Definition at line 114 of file CellUtil.f90.

116  use cellpolymodule, only: cellpolytype
117  use mathutilmodule, only: mod_offset
118  ! -- dummy
119  type(CellPolyType), intent(in), pointer :: poly
120  type(CellRectQuadType), intent(inout), pointer :: quad
121  ! -- local
122  integer(I4B) :: i, irv, isc
123  real(DP) :: qhalf, qdisttopbot, q1, q2, q4
124 
125  call create_cell_rect_quad(quad)
126  call quad%init_from(poly%defn)
127  ! -- Translate and rotate the rect-quad cell into local coordinates with
128  ! -- x varying from 0 to dx and y varying from 0 to dy. Choose the "south-
129  ! -- west" rectangle vertex to be the local origin so that the rotation
130  ! -- angle is zero if the cell already aligns with the model x and y
131  ! -- coordinates.
132  quad%irvOrigin = quad%get_rect_ivert_sw()
133  call quad%get_rect_dim_rot()
134 
135  ! -- Set the external and internal face flows used for subcells
136  do i = 0, 3
137  irv = mod_offset(i + quad%irvOrigin, 4, 1)
138  isc = mod_offset(i + 3, 4, 1)
139  if (.not. quad%face_is_refined(irv)) then
140  qhalf = 5d-1 * quad%get_rect_flow(1, irv)
141  quad%qextl2(isc) = qhalf
142  isc = mod_offset(isc + 1, 4, 1)
143  quad%qextl1(isc) = qhalf
144  else
145  quad%qextl2(isc) = quad%get_rect_flow(1, irv)
146  isc = mod_offset(isc + 1, 4, 1)
147  quad%qextl1(isc) = quad%get_rect_flow(2, irv)
148  end if
149  end do
150  qdisttopbot = 2.5d-1 * (quad%defn%get_distflow() &
151  + quad%defn%get_botflow() &
152  + quad%defn%get_topflow())
153  q1 = qdisttopbot + quad%qextl1(1) + quad%qextl2(1)
154  q2 = qdisttopbot + quad%qextl1(2) + quad%qextl2(2)
155  q4 = qdisttopbot + quad%qextl1(4) + quad%qextl2(4)
156  quad%qintl(1) = -5d-1 * (q1 + 5d-1 * (q2 - q4))
157  quad%qintl(2) = quad%qintl(1) + q1
158  quad%qintl(3) = quad%qintl(2) + q2
159  quad%qintl(4) = quad%qintl(1) - q4
160  quad%qintl(5) = quad%qintl(1)
subroutine, public create_cell_rect_quad(cell)
Create a new rectangular-quad cell.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ cell_poly_to_rect()

subroutine, public cellutilmodule::cell_poly_to_rect ( type(cellpolytype), intent(in), pointer  poly,
type(cellrecttype), intent(inout), pointer  rect 
)

Definition at line 13 of file CellUtil.f90.

14  use constantsmodule, only: done
16  use cellpolymodule, only: cellpolytype
17  use celldefnmodule, only: celldefntype
18  ! -- dummy
19  type(CellPolyType), intent(in), pointer :: poly
20  type(CellRectType), intent(inout), pointer :: rect
21  ! -- local
22  type(CellDefnType), pointer :: defn
23  integer(I4B) :: ipv, ipv1, ipv2, ipv3, ipv4
24  integer(I4B), dimension(4) :: ipvnxt = (/2, 3, 4, 1/)
25  real(DP) :: x1, y1, x2, y2, x4, y4
26  real(DP) :: dx2, dy2, dx4, dy4, areax, areay, areaz
27  real(DP) :: xOrigin, yOrigin, zOrigin, dx, dy, dz, sinrot, cosrot
28  real(DP) :: factor, term
29 
30  call create_cell_rect(rect)
31  defn => poly%defn
32  ! -- Translate and rotate the rectangular cell into local coordinates
33  ! -- with x varying from 0 to dx and y varying from 0 to dy. Choose the
34  ! -- "south-west" vertex to be the local origin so that the rotation
35  ! -- angle is zero if the cell already aligns with the model x and y
36  ! -- coordinates. The "southwest" vertex is found by finding the vertex
37  ! -- formed either by (1) an edge over which x decreases (going
38  ! -- clockwise) followed by an edge over which x does not increase, or
39  ! -- by (2) an edge over which y does not decrease (again going
40  ! -- clockwise) followed by an edge over which y increases. In the end,
41  ! -- ipv1 is the local vertex number (within the cell, taking a value
42  ! -- of 1, 2, 3, or 4) of the southwest vertex, and ipv2, ipv3, and
43  ! -- ipv4 are the local vertex numbers of the remaining three vertices
44  ! -- going clockwise.
45  do ipv = 1, 4
46  ipv4 = ipv
47  ipv1 = ipvnxt(ipv4)
48  x4 = defn%polyvert(1, ipv4)
49  y4 = defn%polyvert(2, ipv4)
50  x1 = defn%polyvert(1, ipv1)
51  y1 = defn%polyvert(2, ipv1)
52  if (x1 .lt. x4) then
53  ipv2 = ipvnxt(ipv1)
54  x2 = defn%polyvert(1, ipv2)
55  if (x2 .le. x1) then
56  y2 = defn%polyvert(2, ipv2)
57  exit
58  end if
59  else if (y1 .ge. y4) then
60  ipv2 = ipvnxt(ipv1)
61  y2 = defn%polyvert(2, ipv2)
62  if (y2 .gt. y1) then
63  x2 = defn%polyvert(1, ipv2)
64  exit
65  end if
66  end if
67  end do
68  ipv3 = ipvnxt(ipv2)
69 
70  ! -- Compute upper bounds on the local coordinates (the rectangular
71  ! -- dimensions of the cell) and the sine and cosine of the rotation
72  ! -- angle, and store local origin information
73  xorigin = x1
74  yorigin = y1
75  zorigin = defn%bot
76  dx2 = x2 - xorigin
77  dy2 = y2 - yorigin
78  dx4 = x4 - xorigin
79  dy4 = y4 - yorigin
80  dx = dsqrt(dx4 * dx4 + dy4 * dy4)
81  dy = dsqrt(dx2 * dx2 + dy2 * dy2)
82  dz = defn%top - zorigin
83  sinrot = dy4 / dx
84  cosrot = dx4 / dx
85  rect%defn = poly%defn
86  rect%dx = dx
87  rect%dy = dy
88  rect%dz = dz
89  rect%sinrot = sinrot
90  rect%cosrot = cosrot
91  rect%xOrigin = xorigin
92  rect%yOrigin = yorigin
93  rect%zOrigin = zorigin
94  rect%ipvOrigin = ipv1
95 
96  ! -- Compute (unscaled) cell edge velocities from face flows
97  areax = dx * dz
98  areay = dy * dz
99  areaz = dx * dy
100  factor = done / (defn%retfactor * defn%porosity)
101  term = factor / areax
102  rect%vx1 = defn%faceflow(ipv1) * term
103  rect%vx2 = -defn%faceflow(ipv3) * term
104  term = factor / areay
105  rect%vy1 = defn%faceflow(ipv4) * term
106  rect%vy2 = -defn%faceflow(ipv2) * term
107  term = factor / areaz
108  rect%vz1 = defn%faceflow(6) * term
109  rect%vz2 = -defn%faceflow(7) * term
subroutine, public create_cell_rect(cell)
Create a new rectangular cell.
Definition: CellRect.f90:40
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
Base grid cell definition.
Definition: CellDefn.f90:12
Here is the call graph for this function:
Here is the caller graph for this function: