MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
MethodSubcellPollock.f90
Go to the documentation of this file.
2  use kindmodule, only: dp, i4b, lgp
3  use errorutilmodule, only: pstop
8  use prtfmimodule, only: prtfmitype
9  use basedismodule, only: disbasetype
10  use cellmodule, only: celltype
11  use constantsmodule, only: dzero, done
12  implicit none
13  private
14  public :: methodsubcellpollocktype
16  public :: calculate_dt
17 
18  !> @brief Rectangular subcell tracking method
20  private
21  real(dp), allocatable, public :: qextl1(:), qextl2(:), qintl(:) !< external and internal subcell flows
22  contains
23  procedure, public :: apply => apply_msp
24  procedure, public :: deallocate
25  procedure, private :: track_subcell
27 
28 contains
29 
30  !> @brief Create a new Pollock's subcell method
31  subroutine create_method_subcell_pollock(method)
32  ! dummy
33  type(methodsubcellpollocktype), pointer :: method
34  ! local
35  type(subcellrecttype), pointer :: subcell
36 
37  allocate (method)
38  call create_subcell_rect(subcell)
39  method%subcell => subcell
40  method%name => method%subcell%type
41  method%delegates = .false.
42  end subroutine create_method_subcell_pollock
43 
44  !> @brief Deallocate the Pollock's subcell method
45  subroutine deallocate (this)
46  class(methodsubcellpollocktype), intent(inout) :: this
47  deallocate (this%name)
48  end subroutine deallocate
49 
50  !> @brief Apply Pollock's method to a rectangular subcell
51  subroutine apply_msp(this, particle, tmax)
52  ! dummy
53  class(methodsubcellpollocktype), intent(inout) :: this
54  type(particletype), pointer, intent(inout) :: particle
55  real(DP), intent(in) :: tmax
56  ! local
57  real(DP) :: xOrigin
58  real(DP) :: yOrigin
59  real(DP) :: zOrigin
60  real(DP) :: sinrot
61  real(DP) :: cosrot
62 
63  select type (subcell => this%subcell)
64  type is (subcellrecttype)
65  ! Transform particle position into local subcell coordinates,
66  ! track particle across subcell, convert back to model coords
67  ! (sinrot and cosrot should be 0 and 1, respectively, i.e. no
68  ! rotation, also no z translation; only x and y translations)
69  xorigin = subcell%xOrigin
70  yorigin = subcell%yOrigin
71  zorigin = subcell%zOrigin
72  sinrot = subcell%sinrot
73  cosrot = subcell%cosrot
74  call particle%transform(xorigin, yorigin)
75  call this%track_subcell(subcell, particle, tmax)
76  call particle%transform(xorigin, yorigin, invert=.true.)
77  end select
78  end subroutine apply_msp
79 
80  !> @brief Track a particle across a rectangular subcell using Pollock's method
81  !!
82  !! This subroutine consists partly of code written by
83  !! David W. Pollock of the USGS for MODPATH 7. PRT's
84  !! authors take responsibility for its application in
85  !! this context and for any modifications or errors.
86  !<
87  subroutine track_subcell(this, subcell, particle, tmax)
88  use tdismodule, only: endofsimulation
91  ! dummy
92  class(methodsubcellpollocktype), intent(inout) :: this
93  class(subcellrecttype), intent(in) :: subcell
94  type(particletype), pointer, intent(inout) :: particle
95  real(DP), intent(in) :: tmax
96  ! local
97  real(DP) :: vx
98  real(DP) :: dvxdx
99  real(DP) :: vy
100  real(DP) :: dvydy
101  real(DP) :: vz
102  real(DP) :: dvzdz
103  real(DP) :: dtexitx
104  real(DP) :: dtexity
105  real(DP) :: dtexitz
106  real(DP) :: dtexit
107  real(DP) :: texit
108  real(DP) :: dt
109  real(DP) :: t
110  real(DP) :: t0
111  real(DP) :: x
112  real(DP) :: y
113  real(DP) :: z
114  integer(I4B) :: statusVX
115  integer(I4B) :: statusVY
116  integer(I4B) :: statusVZ
117  integer(I4B) :: i
118  real(DP) :: initialX
119  real(DP) :: initialY
120  real(DP) :: initialZ
121  integer(I4B) :: exitFace
122  integer(I4B) :: event_code
123 
124  event_code = -1
125 
126  ! Initial particle location in scaled subcell coordinates
127  initialx = particle%x / subcell%dx
128  initialy = particle%y / subcell%dy
129  initialz = particle%z / subcell%dz
130 
131  ! Compute time of travel to each possible exit face
132  statusvx = calculate_dt(subcell%vx1, subcell%vx2, subcell%dx, &
133  initialx, vx, dvxdx, dtexitx)
134  statusvy = calculate_dt(subcell%vy1, subcell%vy2, subcell%dy, &
135  initialy, vy, dvydy, dtexity)
136  statusvz = calculate_dt(subcell%vz1, subcell%vz2, subcell%dz, &
137  initialz, vz, dvzdz, dtexitz)
138 
139  ! Subcell has no exit face, terminate the particle
140  ! todo: after initial release, consider ramifications
141  if ((statusvx .eq. 3) .and. (statusvy .eq. 3) .and. (statusvz .eq. 3)) then
142  call this%terminate(particle, status=term_no_exits_sub)
143  return
144  end if
145 
146  ! If particle stationary and extended tracking is on, terminate here if it's
147  ! the last timestep. TODO: temporary solution, consider where to catch this?
148  ! Should we really have to special case this here? We do diverge from MP7 in
149  ! guaranteeing that every particle terminates at the end of the simulation..
150  ! ideally that would be handled at a higher scope but with extended tracking
151  ! tmax is not the end of the simulation, it's just a wildly high upper bound.
152  if ((statusvx .eq. 2) .and. (statusvy .eq. 2) .and. (statusvz .eq. 2) .and. &
153  particle%iextend > 0 .and. endofsimulation) then
154  call this%terminate(particle, status=term_timeout)
155  return
156  end if
157 
158  ! Determine (earliest) exit face and corresponding travel time to exit
159  exitface = 0
160  dtexit = 1.0d+30
161  if ((statusvx .lt. 2) .or. (statusvy .lt. 2) .or. (statusvz .lt. 2)) then
162  ! Consider x-oriented faces
163  dtexit = dtexitx
164  if (vx .lt. dzero) then
165  exitface = 1
166  else if (vx .gt. 0) then
167  exitface = 2
168  end if
169  ! Consider y-oriented faces
170  if (dtexity .lt. dtexit) then
171  dtexit = dtexity
172  if (vy .lt. dzero) then
173  exitface = 3
174  else if (vy .gt. dzero) then
175  exitface = 4
176  end if
177  end if
178  ! Consider z-oriented faces
179  if (dtexitz .lt. dtexit) then
180  dtexit = dtexitz
181  if (vz .lt. dzero) then
182  exitface = 5
183  else if (vz .gt. dzero) then
184  exitface = 6
185  end if
186  end if
187  else
188  end if
189 
190  ! Compute exit time
191  texit = particle%ttrack + dtexit
192  t0 = particle%ttrack
193 
194  ! Select user tracking times to solve. If this is the first time step
195  ! of the simulation, include all times before it begins; if it is the
196  ! last time step include all times after it ends only if the 'extend'
197  ! option is on, otherwise times in this period and time step only.
198  call this%tracktimes%advance()
199  if (this%tracktimes%any()) then
200  do i = this%tracktimes%selection(1), this%tracktimes%selection(2)
201  t = this%tracktimes%times(i)
202  if (t < particle%ttrack) cycle
203  if (t >= texit .or. t >= tmax) exit
204  dt = t - t0
205  x = new_x(vx, dvxdx, subcell%vx1, subcell%vx2, &
206  dt, initialx, subcell%dx, statusvx == 1)
207  y = new_x(vy, dvydy, subcell%vy1, subcell%vy2, &
208  dt, initialy, subcell%dy, statusvy == 1)
209  z = new_x(vz, dvzdz, subcell%vz1, subcell%vz2, &
210  dt, initialz, subcell%dz, statusvz == 1)
211  particle%x = x * subcell%dx
212  particle%y = y * subcell%dy
213  particle%z = z * subcell%dz
214  particle%ttrack = t
215  particle%istatus = active
216  call this%usertime(particle)
217  end do
218  end if
219 
220  if (texit .gt. tmax) then
221  ! The computed exit time is greater than the maximum time, so set
222  ! final time for particle trajectory equal to maximum time and
223  ! calculate particle location at that final time.
224  t = tmax
225  dt = t - t0
226  x = new_x(vx, dvxdx, subcell%vx1, subcell%vx2, &
227  dt, initialx, subcell%dx, statusvx == 1)
228  y = new_x(vy, dvydy, subcell%vy1, subcell%vy2, &
229  dt, initialy, subcell%dy, statusvy == 1)
230  z = new_x(vz, dvzdz, subcell%vz1, subcell%vz2, &
231  dt, initialz, subcell%dz, statusvz == 1)
232  exitface = 0
233  particle%istatus = active
234  particle%advancing = .false.
235  event_code = timestep ! timestep end
236  else
237  ! The computed exit time is less than or equal to the maximum time,
238  ! so set final time for particle trajectory equal to exit time and
239  ! calculate exit location.
240  t = texit
241  dt = dtexit
242  if ((exitface .eq. 1) .or. (exitface .eq. 2)) then
243  x = dzero
244  y = new_x(vy, dvydy, subcell%vy1, subcell%vy2, &
245  dt, initialy, subcell%dy, statusvy == 1)
246  z = new_x(vz, dvzdz, subcell%vz1, subcell%vz2, &
247  dt, initialz, subcell%dz, statusvz == 1)
248  if (exitface .eq. 2) x = done
249  else if ((exitface .eq. 3) .or. (exitface .eq. 4)) then
250  x = new_x(vx, dvxdx, subcell%vx1, subcell%vx2, dt, &
251  initialx, subcell%dx, statusvx == 1)
252  y = dzero
253  z = new_x(vz, dvzdz, subcell%vz1, subcell%vz2, dt, &
254  initialz, subcell%dz, statusvz == 1)
255  if (exitface .eq. 4) y = done
256  else if ((exitface .eq. 5) .or. (exitface .eq. 6)) then
257  x = new_x(vx, dvxdx, subcell%vx1, subcell%vx2, &
258  dt, initialx, subcell%dx, statusvx == 1)
259  y = new_x(vy, dvydy, subcell%vy1, subcell%vy2, &
260  dt, initialy, subcell%dy, statusvy == 1)
261  z = dzero
262  if (exitface .eq. 6) z = done
263  else
264  print *, "programmer error, invalid exit face", exitface
265  call pstop(1)
266  end if
267  event_code = featexit
268  end if
269 
270  ! Set final particle location in local (unscaled) subcell coordinates,
271  ! final time for particle trajectory, and exit face
272  particle%x = x * subcell%dx
273  particle%y = y * subcell%dy
274  particle%z = z * subcell%dz
275  particle%ttrack = t
276  particle%iboundary(level_subfeature) = exitface
277 
278  ! Save particle track record
279  if (event_code == timestep) then
280  call this%timestep(particle)
281  else if (event_code == featexit) then
282  call this%subcellexit(particle)
283  end if
284 
285  end subroutine track_subcell
286 
287  !> @brief Calculate particle travel time to exit and exit status.
288  !!
289  !! This subroutine consists partly or entirely of code written by
290  !! David W. Pollock of the USGS for MODPATH 7. The authors of the present
291  !! code are responsible for its appropriate application in this context
292  !! and for any modifications or errors.
293  !<
294  function calculate_dt(v1, v2, dx, xL, v, dvdx, dt) result(status)
295  ! dummy
296  real(dp) :: v1
297  real(dp) :: v2
298  real(dp) :: dx
299  real(dp) :: xl
300  real(dp) :: v
301  real(dp) :: dvdx
302  real(dp) :: dt
303  ! result
304  integer(I4B) :: status
305  ! local
306  real(dp) :: v2a
307  real(dp) :: v1a
308  real(dp) :: dv
309  real(dp) :: dva
310  real(dp) :: vv
311  real(dp) :: vvv
312  real(dp) :: zro
313  real(dp) :: zrom
314  real(dp) :: x
315  real(dp) :: tol
316  real(dp) :: vr1
317  real(dp) :: vr2
318  real(dp) :: vr
319  real(dp) :: v1v2
320  logical(LGP) :: nooutflow
321 
322  ! Initialize variables.
323  status = -1
324  dt = 1.0d+20
325  v2a = v2
326  if (v2a .lt. dzero) v2a = -v2a
327  v1a = v1
328  if (v1a .lt. dzero) v1a = -v1a
329  dv = v2 - v1
330  dva = dv
331  if (dva .lt. dzero) dva = -dva
332 
333  ! Check for a uniform zero velocity in this direction.
334  ! If so, set status = 2 and return (dt = 1.0d+20).
335  tol = 1.0d-15
336  if ((v2a .lt. tol) .and. (v1a .lt. tol)) then
337  v = dzero
338  dvdx = dzero
339  status = 2
340  return
341  end if
342 
343  ! Check for uniform non-zero velocity in this direction.
344  ! If so, set compute dt using the constant velocity,
345  ! set status = 1 and return.
346  vv = v1a
347  if (v2a .gt. vv) vv = v2a
348  vvv = dva / vv
349  if (vvv .lt. 1.0d-4) then
350  zro = tol
351  zrom = -zro
352  v = v1
353  x = xl * dx
354  if (v1 .gt. zro) dt = (dx - x) / v1
355  if (v1 .lt. zrom) dt = -x / v1
356  dvdx = dzero
357  status = 1
358  return
359  end if
360 
361  ! Velocity has a linear variation.
362  ! Compute velocity corresponding to particle position.
363  dvdx = dv / dx
364  v = (done - xl) * v1 + xl * v2
365 
366  ! If flow is into the cell from both sides there is no outflow.
367  ! In that case, set status = 3 and return.
368  nooutflow = .true.
369  if (v1 .lt. dzero) nooutflow = .false.
370  if (v2 .gt. dzero) nooutflow = .false.
371  if (nooutflow) then
372  status = 3
373  return
374  end if
375 
376  ! If there is a divide in the cell for this flow direction, check to
377  ! see if the particle is located exactly on the divide. If it is, move
378  ! it very slightly to get it off the divide. This avoids possible
379  ! numerical problems related to stagnation points.
380  if ((v1 .le. dzero) .and. (v2 .ge. dzero)) then
381  if (abs(v) .le. dzero) then
382  v = 1.0d-20
383  if (v2 .le. dzero) v = -v
384  end if
385  end if
386 
387  ! If there is a flow divide, this check finds out what side of the
388  ! divide the particle is on and sets the value of vr appropriately
389  ! to reflect that location.
390  vr1 = v1 / v
391  vr2 = v2 / v
392  vr = vr1
393  if (vr .le. dzero) then
394  vr = vr2
395  end if
396 
397  ! If the product v1*v2 > 0, the velocity is in the same direction
398  ! throughout the cell (i.e. no flow divide). If so, set the value
399  ! of vr to reflect the appropriate direction.
400  v1v2 = v1 * v2
401  if (v1v2 .gt. dzero) then
402  if (v .gt. dzero) vr = vr2
403  if (v .lt. dzero) vr = vr1
404  end if
405 
406  ! Check if vr is (very close to) zero.
407  ! If so, set status = 2 and return (dt = 1.0d+20).
408  if (dabs(vr) .lt. 1.0d-10) then
409  v = dzero
410  dvdx = dzero
411  status = 2
412  return
413  end if
414 
415  ! Compute travel time to exit face. Return with status = 0.
416  dt = log(vr) / dvdx
417  status = 0
418 
419  end function calculate_dt
420 
421  !> @brief Update a cell-local coordinate based on a time increment.
422  !!
423  !! This subroutine consists partly or entirely of code written by
424  !! David W. Pollock of the USGS for MODPATH 7. The authors of the present
425  !! code are responsible for its appropriate application in this context
426  !! and for any modifications or errors.
427  !<
428  pure function new_x(v, dvdx, v1, v2, dt, x, dx, velocity_profile) result(newx)
429  ! dummy
430  real(dp), intent(in) :: v
431  real(dp), intent(in) :: dvdx
432  real(dp), intent(in) :: v1
433  real(dp), intent(in) :: v2
434  real(dp), intent(in) :: dt
435  real(dp), intent(in) :: x
436  real(dp), intent(in) :: dx
437  logical(LGP), intent(in), optional :: velocity_profile
438  ! result
439  real(dp) :: newx
440  logical(LGP) :: lprofile
441 
442  ! process optional arguments
443  if (present(velocity_profile)) then
444  lprofile = velocity_profile
445  else
446  lprofile = .false.
447  end if
448 
449  ! recompute coordinate
450  newx = x
451  if (lprofile) then
452  newx = newx + (v1 * dt / dx)
453  else if (v .ne. dzero) then
454  newx = newx + (v * (exp(dvdx * dt) - done) / dvdx / dx)
455  end if
456 
457  ! clamp to [0, 1]
458  if (newx .lt. dzero) newx = dzero
459  if (newx .gt. done) newx = done
460 
461  end function new_x
462 
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
Particle tracking strategies.
Definition: Method.f90:2
@, public level_subfeature
Definition: Method.f90:38
pure real(dp) function new_x(v, dvdx, v1, v2, dt, x, dx, velocity_profile)
Update a cell-local coordinate based on a time increment.
integer(i4b) function, public calculate_dt(v1, v2, dx, xL, v, dvdx, dt)
Calculate particle travel time to exit and exit status.
subroutine track_subcell(this, subcell, particle, tmax)
Track a particle across a rectangular subcell using Pollock's method.
subroutine, public create_method_subcell_pollock(method)
Create a new Pollock's subcell method.
subroutine apply_msp(this, particle, tmax)
Apply Pollock's method to a rectangular subcell.
@, public featexit
particle exited a grid feature
@, public timestep
time step ended
@ term_timeout
terminated at stop time or end of simulation
Definition: Particle.f90:37
@ term_no_exits_sub
terminated in a subcell with no exit face
Definition: Particle.f90:36
subroutine, public create_subcell_rect(subcell)
Create a new rectangular subcell.
Definition: SubcellRect.f90:28
logical(lgp), pointer, public endofsimulation
flag indicating end of simulation
Definition: tdis.f90:28
Base type for grid cells of a concrete type. Contains a cell-definition which is information shared b...
Definition: Cell.f90:11
Particle tracked by the PRT model.
Definition: Particle.f90:56