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