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