MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
MethodCellPollockQuad.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b
4  use errorutilmodule, only: pstop
5  use constantsmodule, only: done, dzero
6  use methodmodule, only: methodtype
11  use particlemodule, only: particletype
12  implicit none
13 
14  private
16  public :: create_method_cell_quad
17 
19  contains
20  procedure, public :: apply => apply_mcpq
21  procedure, public :: deallocate
22  procedure, public :: load => load_mcpq
23  procedure, public :: load_subcell
24  procedure, public :: pass => pass_mcpq
26 
27 contains
28 
29  !> @brief Create a new Pollock quad-refined cell method
30  subroutine create_method_cell_quad(method)
31  ! dummy
32  type(methodcellpollockquadtype), pointer :: method
33  ! local
34  type(cellrectquadtype), pointer :: cell
35  type(subcellrecttype), pointer :: subcell
36 
37  allocate (method)
38  call create_cell_rect_quad(cell)
39  method%cell => cell
40  method%type => method%cell%type
41  method%delegates = .true.
42  call create_subcell_rect(subcell)
43  method%subcell => subcell
44  end subroutine create_method_cell_quad
45 
46  !> @brief Deallocate the Pollock quad-refined cell method
47  subroutine deallocate (this)
48  class(methodcellpollockquadtype), intent(inout) :: this
49  deallocate (this%type)
50  end subroutine deallocate
51 
52  !> @brief Load subcell into tracking method
53  subroutine load_mcpq(this, particle, next_level, submethod)
54  class(methodcellpollockquadtype), intent(inout) :: this
55  type(particletype), pointer, intent(inout) :: particle
56  integer, intent(in) :: next_level
57  class(methodtype), pointer, intent(inout) :: submethod
58 
59  select type (subcell => this%subcell)
60  type is (subcellrecttype)
61  call this%load_subcell(particle, subcell)
62  end select
63  call method_subcell_plck%init( &
64  fmi=this%fmi, &
65  cell=this%cell, &
66  subcell=this%subcell, &
67  trackctl=this%trackctl, &
68  tracktimes=this%tracktimes)
69  submethod => method_subcell_plck
70  end subroutine load_mcpq
71 
72  !> @brief Pass particle to next subcell if there is one, or to the cell face
73  subroutine pass_mcpq(this, particle)
74  ! dummy
75  class(methodcellpollockquadtype), intent(inout) :: this
76  type(particletype), pointer, intent(inout) :: particle
77  ! local
78  integer(I4B) :: isc, exitFace, npolyverts, inface, infaceoff
79 
80  select type (cell => this%cell)
81  type is (cellrectquadtype)
82  exitface = particle%iboundary(3)
83  isc = particle%idomain(3)
84  npolyverts = cell%defn%npolyverts
85 
86  ! exitFace uses MODPATH 7 iface convention here
87  select case (exitface)
88  case (0)
89  ! Subcell interior (cell interior)
90  inface = -1
91  case (1)
92  select case (isc)
93  case (1)
94  ! W face, subcell 1 --> E face, subcell 4 (cell interior)
95  particle%idomain(3) = 4
96  particle%iboundary(3) = 2
97  inface = 0 ! want Domain(2) unchanged; Boundary(2) = 0
98  case (2)
99  ! W face, subcell 2 --> E face, subcell 3 (cell interior)
100  particle%idomain(3) = 3
101  particle%iboundary(3) = 2
102  inface = 0 ! want Domain(2) unchanged; Boundary(2) = 0
103  case (3)
104  ! W face, subcell 3 (cell face)
105  inface = 1 ! want Domain(2) = -Domain(2); Boundary(2) = inface
106  infaceoff = 0
107  case (4)
108  ! W face, subcell 4 (cell face)
109  inface = 2 ! want Domain(2) = -Domain(2); Boundary(2) = inface
110  infaceoff = -1
111  end select
112  case (2)
113  select case (isc)
114  case (1)
115  ! E face, subcell 1 (cell face)
116  inface = 3 ! want Domain(2) = -Domain(2); Boundary(2) = inface
117  infaceoff = 0
118  case (2)
119  ! E face, subcell 2 (cell face)
120  inface = 4 ! want Domain(2) = -Domain(2); Boundary(2) = inface
121  infaceoff = -1
122  case (3)
123  ! E face, subcell 3 --> W face, subcell 2 (cell interior)
124  particle%idomain(3) = 2
125  particle%iboundary(3) = 1
126  inface = 0 ! want Domain(2) unchanged; Boundary(2) = 0
127  case (4)
128  ! E face, subcell 4 --> W face subcell 1 (cell interior)
129  particle%idomain(3) = 1
130  particle%iboundary(3) = 1
131  inface = 0 ! want Domain(2) unchanged; Boundary(2) = 0
132  end select
133  case (3)
134  select case (isc)
135  case (1)
136  ! S face, subcell 1 --> N face, subcell 2 (cell interior)
137  particle%idomain(3) = 2
138  particle%iboundary(3) = 4
139  inface = 0 ! want Domain(2) unchanged; Boundary(2) = 0
140  case (2)
141  ! S face, subcell 2 (cell face)
142  inface = 4 ! want Domain(2) = -Domain(2); Boundary(2) = inface
143  infaceoff = 0
144  case (3)
145  ! S face, subcell 3 (cell face)
146  inface = 1 ! want Domain(2) = -Domain(2); Boundary(2) = inface
147  infaceoff = -1
148  case (4)
149  ! S face, subcell 4 --> N face, subcell 3 (cell interior)
150  particle%idomain(3) = 3
151  particle%iboundary(3) = 4
152  inface = 0 ! want Domain(2) unchanged; Boundary(2) = 0
153  end select
154  case (4)
155  select case (isc)
156  case (1)
157  ! N face, subcell 1 (cell face)
158  inface = 3 ! want Domain(2) = -Domain(2); Boundary(2) = inface
159  infaceoff = -1
160  case (2)
161  ! N face, subcell 2 --> S face, subcell 1 (cell interior)
162  particle%idomain(3) = 1
163  particle%iboundary(3) = 3
164  inface = 0 ! want Domain(2) unchanged; Boundary(2) = 0
165  case (3)
166  ! N face, subcell 3 --> S face, subcell 4 (cell interior)
167  particle%idomain(3) = 4
168  particle%iboundary(3) = 3
169  inface = 0 ! want Domain(2) unchanged; Boundary(2) = 0
170  case (4)
171  ! N face, subcell 4 (cell face)
172  inface = 2 ! want Domain(2) = -Domain(2); Boundary(2) = inface
173  infaceoff = 0
174  end select
175  case (5)
176  ! Subcell bottom (cell bottom)
177  inface = npolyverts + 2 ! want Domain(2) = -Domain(2); Boundary(2) = inface
178  case (6)
179  ! Subcell top (cell top)
180  inface = npolyverts + 3 ! want Domain(2) = -Domain(2); Boundary(2) = inface
181  end select
182 
183  if (inface .eq. -1) then
184  particle%iboundary(2) = 0
185  else if (inface .eq. 0) then
186  particle%iboundary(2) = 0
187  else
188  if ((inface .ge. 1) .and. (inface .le. 4)) then
189  ! Account for local cell rotation
190  inface = inface + cell%irvOrigin - 1
191  if (inface .gt. 4) inface = inface - 4
192  inface = cell%irectvert(inface) + infaceoff
193  if (inface .lt. 1) inface = inface + npolyverts
194  end if
195  particle%iboundary(2) = inface
196  end if
197  end select
198  end subroutine pass_mcpq
199 
200  !> @brief Solve the quad-rectangular cell via Pollock's method
201  subroutine apply_mcpq(this, particle, tmax)
202  ! dummy
203  class(methodcellpollockquadtype), intent(inout) :: this
204  type(particletype), pointer, intent(inout) :: particle
205  real(DP), intent(in) :: tmax
206  ! local
207  double precision :: xOrigin, yOrigin, zOrigin, sinrot, cosrot
208 
209  select type (cell => this%cell)
210  type is (cellrectquadtype)
211  ! Update particle state, return early if done advancing
212  call this%update(particle, cell%defn)
213  if (.not. particle%advancing) return
214 
215  ! If the particle is above the top of the cell (presumed water table)
216  ! pass it vertically and instantaneously to the top
217  if (particle%z > cell%defn%top) then
218  particle%z = cell%defn%top
219  call this%save(particle, reason=1) ! reason=1: cell transition
220  end if
221 
222  ! Transform particle location into local cell coordinates
223  ! (translated and rotated but not scaled relative to model).
224  xorigin = cell%xOrigin
225  yorigin = cell%yOrigin
226  zorigin = cell%zOrigin
227  sinrot = cell%sinrot
228  cosrot = cell%cosrot
229  call particle%transform(xorigin, yorigin, zorigin, &
230  sinrot, cosrot)
231 
232  ! Track the particle across the cell.
233  call this%track(particle, 2, tmax)
234 
235  ! Transform particle location back to model coordinates, then
236  ! reset transformation and eliminate accumulated roundoff error.
237  call particle%transform(xorigin, yorigin, zorigin, &
238  sinrot, cosrot, invert=.true.)
239  call particle%transform(reset=.true.)
240  end select
241  end subroutine apply_mcpq
242 
243  !> @brief Load the rectangular subcell from the rectangular cell
244  subroutine load_subcell(this, particle, subcell)
245  ! dummy
246  class(methodcellpollockquadtype), intent(inout) :: this
247  type(particletype), pointer, intent(inout) :: particle
248  class(subcellrecttype), intent(inout) :: subcell
249  ! local
250  real(DP) :: dx, dy, dz, areax, areay, areaz
251  real(DP) :: dxprel, dyprel
252  integer(I4B) :: isc, npolyverts, m1, m2
253  real(DP) :: qextl1, qextl2, qintl1, qintl2
254  real(DP) :: factor, term
255 
256  select type (cell => this%cell)
257  type is (cellrectquadtype)
258  factor = done / cell%defn%retfactor
259  factor = factor / cell%defn%porosity
260  npolyverts = cell%defn%npolyverts
261 
262  isc = particle%idomain(3)
263  ! Subcells 1, 2, 3, and 4 are Pollock's subcells A, B, C, and D,
264  ! respectively
265 
266  dx = cell%dx
267  dy = cell%dy
268  ! If not already known, determine subcell number
269  if (isc .le. 0) then
270  dxprel = particle%x / dx
271  dyprel = particle%y / dy
272 
273  if (dyprel .ge. 5d-1) then
274  if (dxprel .le. 5d-1) then
275  isc = 4
276  else
277  isc = 1
278  end if
279  else
280  if (dxprel .le. 5d-1) then
281  isc = 3
282  else
283  isc = 2
284  end if
285  end if
286 
287  subcell%isubcell = isc
288  particle%idomain(3) = isc
289  end if
290  dx = 5d-1 * dx
291  dy = 5d-1 * dy
292  dz = cell%defn%top - &
293  cell%defn%bot
294  areax = dy * dz
295  areay = dx * dz
296  areaz = dx * dy
297  qintl1 = cell%qintl(isc)
298  ! qintl list wraps around, so isc+1=5 is ok
299  qintl2 = cell%qintl(isc + 1)
300  qextl1 = cell%qextl1(isc)
301  qextl2 = cell%qextl2(isc)
302  !
303  subcell%dx = dx
304  subcell%dy = dy
305  subcell%dz = dz
306  subcell%sinrot = dzero
307  subcell%cosrot = done
308  subcell%zOrigin = dzero
309  select case (isc)
310  case (1)
311  subcell%xOrigin = dx
312  subcell%yOrigin = dy
313  term = factor / areax
314  subcell%vx1 = qintl1 * term
315  subcell%vx2 = -qextl2 * term
316  term = factor / areay
317  subcell%vy1 = -qintl2 * term
318  subcell%vy2 = -qextl1 * term
319  case (2)
320  subcell%xOrigin = dx
321  subcell%yOrigin = dzero
322  term = factor / areax
323  subcell%vx1 = -qintl2 * term
324  subcell%vx2 = -qextl1 * term
325  term = factor / areay
326  subcell%vy1 = qextl2 * term
327  subcell%vy2 = -qintl1 * term
328  case (3)
329  subcell%xOrigin = dzero
330  subcell%yOrigin = dzero
331  term = factor / areax
332  subcell%vx1 = qextl2 * term
333  subcell%vx2 = -qintl1 * term
334  term = factor / areay
335  subcell%vy1 = qextl1 * term
336  subcell%vy2 = qintl2 * term
337  case (4)
338  subcell%xOrigin = dzero
339  subcell%yOrigin = dy
340  term = factor / areax
341  subcell%vx1 = qextl1 * term
342  subcell%vx2 = qintl2 * term
343  term = factor / areay
344  subcell%vy1 = qintl1 * term
345  subcell%vy2 = -qextl2 * term
346  end select
347  m1 = npolyverts + 2
348  m2 = m1 + 1
349  term = factor / areaz
350  subcell%vz1 = 2.5d-1 * cell%defn%faceflow(m1) * term
351  subcell%vz2 = -2.5d-1 * cell%defn%faceflow(m2) * term
352  end select
353  end subroutine load_subcell
354 
subroutine, public create_cell_rect_quad(cell)
Create a new rectangular-quad cell.
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
subroutine pstop(status, message)
Stop the program, optionally specifying an error status code.
Definition: ErrorUtil.f90:24
This module defines variable data types.
Definition: kind.f90:8
procedure subroutine, public create_method_cell_quad(method)
Create a new Pollock quad-refined cell method.
subroutine load_mcpq(this, particle, next_level, submethod)
Load subcell into tracking method.
subroutine pass_mcpq(this, particle)
Pass particle to next subcell if there is one, or to the cell face.
subroutine load_subcell(this, particle, subcell)
Load the rectangular subcell from the rectangular cell.
subroutine apply_mcpq(this, particle, tmax)
Solve the quad-rectangular cell via Pollock's method.
Particle tracking strategies.
Definition: Method.f90:2
Subcell-level tracking methods.
type(methodsubcellpollocktype), pointer, public method_subcell_plck
subroutine, public create_subcell_rect(subcell)
Create a new rectangular subcell.
Definition: SubcellRect.f90:28
Base grid cell definition.
Definition: CellDefn.f90:10
Base type for particle tracking methods.
Definition: Method.f90:29
Particle tracked by the PRT model.
Definition: Particle.f90:32