MODFLOW 6  version 6.8.0.dev0
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  use domainmodule, only: domaintype
13  use subcellmodule, only: subcelltype
14  use listmodule, only: listtype
19  implicit none
20  private
21  public :: methodsubcellpollocktype
23  public :: calculate_dt
24 
25  !> @brief Rectangular subcell tracking method
27  private
28  real(dp), allocatable, public :: qextl1(:), qextl2(:), qintl(:) !< external and internal subcell flows
29  type(linearexitsolutiontype), public :: exit_solutions(3) !< candidate exit solutions
30  contains
31  procedure, public :: find_exits
32  procedure, public :: pick_exit
33  procedure, public :: apply => apply_msp
34  procedure, public :: deallocate
35  procedure, private :: track_subcell
37 
38 contains
39 
40  !> @brief Create a new Pollock's subcell method
41  subroutine create_method_subcell_pollock(method)
42  ! dummy
43  type(methodsubcellpollocktype), pointer :: method
44  ! local
45  type(subcellrecttype), pointer :: subcell
46 
47  allocate (method)
48  call create_subcell_rect(subcell)
49  method%subcell => subcell
50  method%name => method%subcell%type
51  method%delegates = .false.
52  end subroutine create_method_subcell_pollock
53 
54  !> @brief Deallocate the Pollock's subcell method
55  subroutine deallocate (this)
56  class(methodsubcellpollocktype), intent(inout) :: this
57  deallocate (this%name)
58  end subroutine deallocate
59 
60  !> @brief Apply Pollock's method to a rectangular subcell
61  subroutine apply_msp(this, particle, tmax)
62  ! dummy
63  class(methodsubcellpollocktype), intent(inout) :: this
64  type(particletype), pointer, intent(inout) :: particle
65  real(DP), intent(in) :: tmax
66  ! local
67  real(DP) :: x_origin
68  real(DP) :: y_origin
69  real(DP) :: z_origin
70  real(DP) :: sinrot
71  real(DP) :: cosrot
72 
73  select type (subcell => this%subcell)
74  type is (subcellrecttype)
75  ! Transform particle position into local subcell coordinates,
76  ! track particle across subcell, convert back to model coords
77  ! (sinrot and cosrot should be 0 and 1, respectively, i.e. no
78  ! rotation, also no z translation; only x and y translations)
79  x_origin = subcell%xOrigin
80  y_origin = subcell%yOrigin
81  z_origin = subcell%zOrigin
82  sinrot = subcell%sinrot
83  cosrot = subcell%cosrot
84  call particle%transform(x_origin, y_origin)
85  call this%track_subcell(subcell, particle, tmax)
86  call particle%transform(x_origin, y_origin, invert=.true.)
87  end select
88  end subroutine apply_msp
89 
90  !> @brief Track a particle across a rectangular subcell using Pollock's method
91  !!
92  !! This subroutine consists partly of code written by
93  !! David W. Pollock of the USGS for MODPATH 7. PRT's
94  !! authors take responsibility for its application in
95  !! this context and for any modifications or errors.
96  !<
97  subroutine track_subcell(this, subcell, particle, tmax)
100  ! dummy
101  class(methodsubcellpollocktype), intent(inout) :: this
102  class(subcellrecttype), intent(in) :: subcell
103  type(particletype), pointer, intent(inout) :: particle
104  real(DP), intent(in) :: tmax
105  ! local
106  real(DP) :: dt, dtexit, texit
107  real(DP) :: t, x, y, z
108  real(DP) :: t0, x0, y0, z0
109  integer(I4B) :: i, exit_face, exit_soln
110  type(linearexitsolutiontype) :: exit_x, exit_y, exit_z
111 
112  t0 = particle%ttrack
113  x0 = particle%x / subcell%dx
114  y0 = particle%y / subcell%dy
115  z0 = particle%z / subcell%dz
116 
117  ! Find exit solution in each direction
118  call this%find_exits(particle, subcell)
119  exit_x = this%exit_solutions(1)
120  exit_y = this%exit_solutions(2)
121  exit_z = this%exit_solutions(3)
122 
123  ! Set solution, face, & travel time
124  exit_soln = this%pick_exit(particle)
125  if (exit_soln == 0) then
126  exit_face = 0
127  dtexit = 1.0d+30
128  else
129  exit_face = this%exit_solutions(exit_soln)%iboundary
130  dtexit = this%exit_solutions(exit_soln)%dt
131  end if
132  texit = particle%ttrack + dtexit
133 
134  ! Terminate if no valid exit solution was found.
135  ! MP7 is more nuanced in determining what to do
136  ! here. It considers whether this stress period
137  ! is steady state or transient, and whether the
138  ! flow is stationary, in determining whether to
139  ! terminate or allow the particle to survive to
140  ! the next time step. It does not compute paths
141  ! within the subcell even for particles it lets
142  ! remain active under this circumstance, though.
143  ! While we may consider that someday, we simply
144  ! terminate and sidestep the complexity for now.
145  if (all([this%exit_solutions%status] >= no_exit_stationary)) then
146  call this%terminate(particle, status=term_no_exits_sub)
147  return
148  end if
149 
150  ! Select user tracking times to solve. If this is the first time step
151  ! of the simulation, include all times before it begins; if it is the
152  ! last time step include all times after it ends only if the 'extend'
153  ! option is on, otherwise times in this period and time step only.
154  call this%tracktimes%advance()
155  if (this%tracktimes%any()) then
156  do i = this%tracktimes%selection(1), this%tracktimes%selection(2)
157  t = this%tracktimes%times(i)
158  if (t < particle%ttrack) cycle
159  if (t > texit .or. t > tmax) exit
160  dt = t - t0
161  x = new_x(exit_x%v, exit_x%dvdx, subcell%vx1, subcell%vx2, &
162  dt, x0, subcell%dx, exit_x%status == 1)
163  y = new_x(exit_y%v, exit_y%dvdx, subcell%vy1, subcell%vy2, &
164  dt, y0, subcell%dy, exit_y%status == 1)
165  z = new_x(exit_z%v, exit_z%dvdx, subcell%vz1, subcell%vz2, &
166  dt, z0, subcell%dz, exit_z%status == 1)
167  particle%x = x * subcell%dx
168  particle%y = y * subcell%dy
169  particle%z = z * subcell%dz
170  particle%ttrack = t
171  particle%istatus = active
172  call this%usertime(particle)
173  end do
174  end if
175 
176  ! Computed exit time greater than the maximum time? Set the
177  ! tracking time to tmax and calculate the particle location.
178  if (texit .gt. tmax) then
179  t = tmax
180  dt = t - t0
181  x = new_x(exit_x%v, exit_x%dvdx, subcell%vx1, subcell%vx2, &
182  dt, x0, subcell%dx, exit_x%status == 1)
183  y = new_x(exit_y%v, exit_y%dvdx, subcell%vy1, subcell%vy2, &
184  dt, y0, subcell%dy, exit_y%status == 1)
185  z = new_x(exit_z%v, exit_z%dvdx, subcell%vz1, subcell%vz2, &
186  dt, z0, subcell%dz, exit_z%status == 1)
187  exit_face = 0
188  particle%istatus = active
189  particle%advancing = .false.
190  particle%x = x * subcell%dx
191  particle%y = y * subcell%dy
192  particle%z = z * subcell%dz
193  particle%ttrack = t
194  particle%iboundary(level_subfeature) = exit_face
195  call this%timestep(particle)
196  return
197  end if
198 
199  ! If we get to here, the particle is exiting the subcell.
200  ! Its exit time is less than or equal to the maximum time.
201  ! Set tracking time to the exit time and set exit location.
202  t = texit
203  dt = dtexit
204  if ((exit_face .eq. 1) .or. (exit_face .eq. 2)) then
205  x = dzero
206  y = new_x(exit_y%v, exit_y%dvdx, subcell%vy1, subcell%vy2, &
207  dt, y0, subcell%dy, exit_y%status == 1)
208  z = new_x(exit_z%v, exit_z%dvdx, subcell%vz1, subcell%vz2, &
209  dt, z0, subcell%dz, exit_z%status == 1)
210  if (exit_face .eq. 2) x = done
211  else if ((exit_face .eq. 3) .or. (exit_face .eq. 4)) then
212  x = new_x(exit_x%v, exit_x%dvdx, subcell%vx1, subcell%vx2, dt, &
213  x0, subcell%dx, exit_x%status == 1)
214  y = dzero
215  z = new_x(exit_z%v, exit_z%dvdx, subcell%vz1, subcell%vz2, dt, &
216  z0, subcell%dz, exit_z%status == 1)
217  if (exit_face .eq. 4) y = done
218  else if ((exit_face .eq. 5) .or. (exit_face .eq. 6)) then
219  x = new_x(exit_x%v, exit_x%dvdx, subcell%vx1, subcell%vx2, &
220  dt, x0, subcell%dx, exit_x%status == 1)
221  y = new_x(exit_y%v, exit_y%dvdx, subcell%vy1, subcell%vy2, &
222  dt, y0, subcell%dy, exit_y%status == 1)
223  z = dzero
224  if (exit_face .eq. 6) z = done
225  else
226  print *, "programmer error, invalid exit face", exit_face
227  call pstop(1)
228  end if
229  particle%x = x * subcell%dx
230  particle%y = y * subcell%dy
231  particle%z = z * subcell%dz
232  particle%ttrack = t
233  particle%iboundary(level_subfeature) = exit_face
234  call this%subcellexit(particle)
235 
236  end subroutine track_subcell
237 
238  !> @brief Pick the exit solution with the shortest travel time
239  function pick_exit(this, particle) result(exit_soln)
240  class(methodsubcellpollocktype), intent(inout) :: this
241  type(particletype), pointer, intent(inout) :: particle
242  integer(I4B) :: exit_soln
243  ! local
244  real(dp) :: dtmin
245 
246  exit_soln = 0
247  dtmin = 1.0d+30
248 
249  if (this%exit_solutions(1)%status < 2) then
250  exit_soln = 1 ! x
251  dtmin = this%exit_solutions(1)%dt
252  end if
253  if (this%exit_solutions(2)%status < 2 .and. &
254  this%exit_solutions(2)%dt < dtmin) then
255  exit_soln = 2 ! y
256  dtmin = this%exit_solutions(2)%dt
257  end if
258  if (this%exit_solutions(3)%status < 2 .and. &
259  this%exit_solutions(3)%dt < dtmin) then
260  exit_soln = 3 ! z
261  end if
262 
263  end function pick_exit
264 
265  !> @brief Compute candidate exit solutions
266  subroutine find_exits(this, particle, domain)
267  class(methodsubcellpollocktype), intent(inout) :: this
268  type(particletype), pointer, intent(inout) :: particle
269  class(domaintype), intent(in) :: domain
270  ! local
271  real(DP) :: x0, y0, z0
272 
273  select type (domain)
274  type is (subcellrecttype)
275  ! Initial particle location in scaled subcell coordinates
276  x0 = particle%x / domain%dx
277  y0 = particle%y / domain%dy
278  z0 = particle%z / domain%dz
279 
280  ! Calculate exit solutions for each coordinate direction
281  this%exit_solutions = [ &
282  find_exit(domain%vx1, domain%vx2, domain%dx, x0), &
283  find_exit(domain%vy1, domain%vy2, domain%dy, y0), &
284  find_exit(domain%vz1, domain%vz2, domain%dz, z0) &
285  ]
286 
287  ! Set exit faces
288  if (this%exit_solutions(1)%v < dzero) then
289  this%exit_solutions(1)%iboundary = 1
290  else if (this%exit_solutions(1)%v > dzero) then
291  this%exit_solutions(1)%iboundary = 2
292  end if
293  if (this%exit_solutions(2)%v < dzero) then
294  this%exit_solutions(2)%iboundary = 3
295  else if (this%exit_solutions(2)%v > dzero) then
296  this%exit_solutions(2)%iboundary = 4
297  end if
298  if (this%exit_solutions(3)%v < dzero) then
299  this%exit_solutions(3)%iboundary = 5
300  else if (this%exit_solutions(3)%v > dzero) then
301  this%exit_solutions(3)%iboundary = 6
302  end if
303  end select
304  end subroutine find_exits
305 
306  !> @brief Find an exit solution for one dimension
307  function find_exit(v1, v2, dx, xL) result(solution)
308  ! dummy
309  real(dp), intent(in) :: v1
310  real(dp), intent(in) :: v2
311  real(dp), intent(in) :: dx
312  real(dp), intent(in) :: xl
313  type(linearexitsolutiontype) :: solution
314 
315  solution = linearexitsolutiontype()
316  solution%status = calculate_dt(v1, v2, dx, xl, &
317  solution%v, solution%dvdx, solution%dt)
318  end function find_exit
319 
320  !> @brief Calculate particle travel time to exit and exit status.
321  !!
322  !! This subroutine consists partly of code written by and/or adapted from
323  !! David W. Pollock of the USGS for MODPATH 7. The authors of the present
324  !! code are responsible for its appropriate application in this context
325  !! and for any modifications or errors.
326  !<
327  function calculate_dt(v1, v2, dx, xL, v, dvdx, dt) result(status)
328  ! dummy
329  real(dp) :: v1
330  real(dp) :: v2
331  real(dp) :: dx
332  real(dp) :: xl
333  real(dp) :: v
334  real(dp) :: dvdx
335  real(dp) :: dt
336  ! result
337  integer(I4B) :: status
338  ! local
339  real(dp) :: v2a
340  real(dp) :: v1a
341  real(dp) :: dv
342  real(dp) :: dva
343  real(dp) :: vv
344  real(dp) :: vvv
345  real(dp) :: zro
346  real(dp) :: zrom
347  real(dp) :: x
348  real(dp) :: tol
349  real(dp) :: vr1
350  real(dp) :: vr2
351  real(dp) :: vr
352  real(dp) :: v1v2
353  logical(LGP) :: nooutflow
354 
355  ! Initialize variables.
356  status = -1
357  dt = 1.0d+20
358  v2a = v2
359  if (v2a .lt. dzero) v2a = -v2a
360  v1a = v1
361  if (v1a .lt. dzero) v1a = -v1a
362  dv = v2 - v1
363  dva = dv
364  if (dva .lt. dzero) dva = -dva
365 
366  ! Check for a uniform zero velocity in this direction.
367  ! If so, set status = 2 and return (dt = 1.0d+20).
368  tol = 1.0d-15
369  if ((v2a .lt. tol) .and. (v1a .lt. tol)) then
370  v = dzero
371  dvdx = dzero
372  status = no_exit_stationary
373  return
374  end if
375 
376  ! Check for uniform non-zero velocity in this direction.
377  ! If so, set compute dt using the constant velocity,
378  ! set status = 1 and return.
379  vv = v1a
380  if (v2a .gt. vv) vv = v2a
381  vvv = dva / vv
382  if (vvv .lt. 1.0d-4) then
383  zro = tol
384  zrom = -zro
385  v = v1
386  x = xl * dx
387  if (v1 .gt. zro) dt = (dx - x) / v1
388  if (v1 .lt. zrom) dt = -x / v1
389  dvdx = dzero
390  status = ok_exit_constant
391  return
392  end if
393 
394  ! Velocity has a linear variation.
395  ! Compute velocity corresponding to particle position.
396  dvdx = dv / dx
397  v = (done - xl) * v1 + xl * v2
398 
399  ! If flow is into the cell from both sides there is no outflow.
400  ! In that case, set status = 3 and return.
401  nooutflow = .true.
402  if (v1 .lt. dzero) nooutflow = .false.
403  if (v2 .gt. dzero) nooutflow = .false.
404  if (nooutflow) then
405  status = no_exit_no_outflow
406  return
407  end if
408 
409  ! If there is a divide in the cell for this flow direction, check to
410  ! see if the particle is located exactly on the divide. If it is, move
411  ! it very slightly to get it off the divide. This avoids possible
412  ! numerical problems related to stagnation points.
413  if ((v1 .le. dzero) .and. (v2 .ge. dzero)) then
414  if (abs(v) .le. dzero) then
415  v = 1.0d-20
416  if (v2 .le. dzero) v = -v
417  end if
418  end if
419 
420  ! If there is a flow divide, this check finds out what side of the
421  ! divide the particle is on and sets the value of vr appropriately
422  ! to reflect that location.
423  vr1 = v1 / v
424  vr2 = v2 / v
425  vr = vr1
426  if (vr .le. dzero) then
427  vr = vr2
428  end if
429 
430  ! If the product v1*v2 > 0, the velocity is in the same direction
431  ! throughout the cell (i.e. no flow divide). If so, set the value
432  ! of vr to reflect the appropriate direction.
433  v1v2 = v1 * v2
434  if (v1v2 .gt. dzero) then
435  if (v .gt. dzero) vr = vr2
436  if (v .lt. dzero) vr = vr1
437  end if
438 
439  ! Check if vr is (very close to) zero.
440  ! If so, set status = 2 and return (dt = 1.0d+20).
441  if (dabs(vr) .lt. 1.0d-10) then
442  v = dzero
443  dvdx = dzero
444  status = no_exit_stationary
445  return
446  end if
447 
448  ! Compute travel time to exit face. Return with status = 0.
449  dt = log(vr) / dvdx
450  status = ok_exit
451 
452  end function calculate_dt
453 
454  !> @brief Update a cell-local coordinate based on a time increment.
455  !!
456  !! This subroutine consists partly or entirely of code written by
457  !! David W. Pollock of the USGS for MODPATH 7. The authors of the present
458  !! code are responsible for its appropriate application in this context
459  !! and for any modifications or errors.
460  !<
461  pure function new_x(v, dvdx, v1, v2, dt, x, dx, velocity_profile) result(newx)
462  ! dummy
463  real(dp), intent(in) :: v
464  real(dp), intent(in) :: dvdx
465  real(dp), intent(in) :: v1
466  real(dp), intent(in) :: v2
467  real(dp), intent(in) :: dt
468  real(dp), intent(in) :: x
469  real(dp), intent(in) :: dx
470  logical(LGP), intent(in), optional :: velocity_profile
471  ! result
472  real(dp) :: newx
473  logical(LGP) :: lprofile
474 
475  ! process optional arguments
476  if (present(velocity_profile)) then
477  lprofile = velocity_profile
478  else
479  lprofile = .false.
480  end if
481 
482  ! recompute coordinate
483  newx = x
484  if (lprofile) then
485  newx = newx + (v1 * dt / dx)
486  else if (v .ne. dzero) then
487  newx = newx + (v * (exp(dvdx * dt) - done) / dvdx / dx)
488  end if
489 
490  ! clamp to [0, 1]
491  if (newx .lt. dzero) newx = dzero
492  if (newx .gt. done) newx = done
493 
494  end function new_x
495 
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
@ ok_exit_constant
exit found, constant velocity
@ ok_exit
exit found using velocity interpolation
@ no_exit_stationary
no exit, zero velocity
@ no_exit_no_outflow
no exit, no outflow
This module defines variable data types.
Definition: kind.f90:8
Particle tracking strategies.
Definition: Method.f90:2
@, public level_subfeature
Definition: Method.f90:42
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 pick_exit(this, particle)
Pick the exit solution with the shortest travel time.
integer(i4b) function, public calculate_dt(v1, v2, dx, xL, v, dvdx, dt)
Calculate particle travel time to exit and exit status.
type(linearexitsolutiontype) function find_exit(v1, v2, dx, xL)
Find an exit solution for one dimension.
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 find_exits(this, particle, domain)
Compute candidate exit solutions.
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:38
@ term_no_exits_sub
terminated in a subcell with no exit face
Definition: Particle.f90:37
subroutine, public create_subcell_rect(subcell)
Create a new rectangular subcell.
Definition: SubcellRect.f90:28
Base type for grid cells of a concrete type. Contains a cell-definition which is information shared b...
Definition: Cell.f90:12
A tracking domain.
Definition: Domain.f90:8
Base type for exit solutions.
Linear velocity interpolation exit solution.
A generic heterogeneous doubly-linked list.
Definition: List.f90:14
Abstract base type for subcell tracking methods.
Particle tracked by the PRT model.
Definition: Particle.f90:62
A subcell of a cell.
Definition: Subcell.f90:10