MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
methodsubcellpollockmodule Module Reference

Data Types

type  methodsubcellpollocktype
 Rectangular subcell tracking method. More...
 

Functions/Subroutines

subroutine, public create_method_subcell_pollock (method)
 Create a new Pollock's subcell-method object. More...
 
subroutine destroy_msp (this)
 Destructor for a Pollock's subcell-method object. More...
 
subroutine apply_msp (this, particle, tmax)
 Apply Pollock's method to a rectangular subcell. More...
 
subroutine track_subcell (this, subcell, particle, tmax)
 Track a particle across a rectangular subcell using Pollock's method. More...
 
integer(i4b) function, public calculate_dt (v1, v2, dx, xL, v, dvdx, dt)
 Calculate particle travel time to exit and exit status. More...
 
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. More...
 

Function/Subroutine Documentation

◆ apply_msp()

subroutine methodsubcellpollockmodule::apply_msp ( class(methodsubcellpollocktype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
real(dp), intent(in)  tmax 
)
private

Definition at line 50 of file MethodSubcellPollock.f90.

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

◆ calculate_dt()

integer(i4b) function, public methodsubcellpollockmodule::calculate_dt ( real(dp)  v1,
real(dp)  v2,
real(dp)  dx,
real(dp)  xL,
real(dp)  v,
real(dp)  dvdx,
real(dp)  dt 
)

This subroutine consists partly or entirely of code written by David W. Pollock of the USGS for MODPATH 7. The authors of the present code are responsible for its appropriate application in this context and for any modifications or errors.

Definition at line 279 of file MethodSubcellPollock.f90.

280  ! dummy
281  real(DP) :: v1
282  real(DP) :: v2
283  real(DP) :: dx
284  real(DP) :: xL
285  real(DP) :: v
286  real(DP) :: dvdx
287  real(DP) :: dt
288  ! result
289  integer(I4B) :: status
290  ! local
291  real(DP) :: v2a
292  real(DP) :: v1a
293  real(DP) :: dv
294  real(DP) :: dva
295  real(DP) :: vv
296  real(DP) :: vvv
297  real(DP) :: zro
298  real(DP) :: zrom
299  real(DP) :: x
300  real(DP) :: tol
301  real(DP) :: vr1
302  real(DP) :: vr2
303  real(DP) :: vr
304  real(DP) :: v1v2
305  logical(LGP) :: noOutflow
306 
307  ! -- Initialize variables.
308  status = -1
309  dt = 1.0d+20
310  v2a = v2
311  if (v2a .lt. 0d0) v2a = -v2a
312  v1a = v1
313  if (v1a .lt. 0d0) v1a = -v1a
314  dv = v2 - v1
315  dva = dv
316  if (dva .lt. 0d0) dva = -dva
317 
318  ! -- Check for a uniform zero velocity in this direction.
319  ! -- If so, set status = 2 and return (dt = 1.0d+20).
320  tol = 1.0d-15
321  if ((v2a .lt. tol) .and. (v1a .lt. tol)) then
322  v = 0d0
323  dvdx = 0d0
324  status = 2
325  return
326  end if
327 
328  ! -- Check for uniform non-zero velocity in this direction.
329  ! -- If so, set compute dt using the constant velocity,
330  ! -- set status = 1 and return.
331  vv = v1a
332  if (v2a .gt. vv) vv = v2a
333  vvv = dva / vv
334  if (vvv .lt. 1.0d-4) then
335  zro = tol
336  zrom = -zro
337  v = v1
338  x = xl * dx
339  if (v1 .gt. zro) dt = (dx - x) / v1
340  if (v1 .lt. zrom) dt = -x / v1
341  dvdx = 0d0
342  status = 1
343  return
344  end if
345 
346  ! -- Velocity has a linear variation.
347  ! -- Compute velocity corresponding to particle position.
348  dvdx = dv / dx
349  v = (1.0d0 - xl) * v1 + xl * v2
350 
351  ! -- If flow is into the cell from both sides there is no outflow.
352  ! -- In that case, set status = 3 and return.
353  nooutflow = .true.
354  if (v1 .lt. 0d0) nooutflow = .false.
355  if (v2 .gt. 0d0) nooutflow = .false.
356  if (nooutflow) then
357  status = 3
358  return
359  end if
360 
361  ! -- If there is a divide in the cell for this flow direction, check to
362  ! -- see if the particle is located exactly on the divide. If it is, move
363  ! -- it very slightly to get it off the divide. This avoids possible
364  ! -- numerical problems related to stagnation points.
365  if ((v1 .le. 0d0) .and. (v2 .ge. 0d0)) then
366  if (abs(v) .le. 0d0) then
367  v = 1.0d-20
368  if (v2 .le. 0d0) v = -v
369  end if
370  end if
371 
372  ! -- If there is a flow divide, this check finds out what side of the
373  ! -- divide the particle is on and sets the value of vr appropriately
374  ! -- to reflect that location.
375  vr1 = v1 / v
376  vr2 = v2 / v
377  vr = vr1
378  if (vr .le. 0d0) then
379  vr = vr2
380  end if
381 
382  ! -- If the product v1*v2 > 0, the velocity is in the same direction
383  ! -- throughout the cell (i.e. no flow divide). If so, set the value
384  ! -- of vr to reflect the appropriate direction.
385  v1v2 = v1 * v2
386  if (v1v2 .gt. 0d0) then
387  if (v .gt. 0d0) vr = vr2
388  if (v .lt. 0d0) vr = vr1
389  end if
390 
391  ! -- Check if vr is (very close to) zero.
392  ! -- If so, set status = 2 and return (dt = 1.0d+20).
393  if (dabs(vr) .lt. 1.0d-10) then
394  v = 0d0
395  dvdx = 0d0
396  status = 2
397  return
398  end if
399 
400  ! -- Compute travel time to exit face. Return with status = 0.
401  dt = log(vr) / dvdx
402  status = 0
403 
Here is the caller graph for this function:

◆ create_method_subcell_pollock()

subroutine, public methodsubcellpollockmodule::create_method_subcell_pollock ( type(methodsubcellpollocktype), pointer  method)

Definition at line 30 of file MethodSubcellPollock.f90.

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.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ destroy_msp()

subroutine methodsubcellpollockmodule::destroy_msp ( class(methodsubcellpollocktype), intent(inout)  this)
private

Definition at line 44 of file MethodSubcellPollock.f90.

45  class(MethodSubcellPollockType), intent(inout) :: this
46  deallocate (this%type)

◆ new_x()

pure real(dp) function methodsubcellpollockmodule::new_x ( real(dp), intent(in)  v,
real(dp), intent(in)  dvdx,
real(dp), intent(in)  v1,
real(dp), intent(in)  v2,
real(dp), intent(in)  dt,
real(dp), intent(in)  x,
real(dp), intent(in)  dx,
logical(lgp), intent(in), optional  velocity_profile 
)
private

This subroutine consists partly or entirely of code written by David W. Pollock of the USGS for MODPATH 7. The authors of the present code are responsible for its appropriate application in this context and for any modifications or errors.

Definition at line 413 of file MethodSubcellPollock.f90.

414  ! dummy
415  real(DP), intent(in) :: v
416  real(DP), intent(in) :: dvdx
417  real(DP), intent(in) :: v1
418  real(DP), intent(in) :: v2
419  real(DP), intent(in) :: dt
420  real(DP), intent(in) :: x
421  real(DP), intent(in) :: dx
422  logical(LGP), intent(in), optional :: velocity_profile
423  ! result
424  real(DP) :: newx
425  logical(LGP) :: lprofile
426 
427  ! -- process optional arguments
428  if (present(velocity_profile)) then
429  lprofile = velocity_profile
430  else
431  lprofile = .false.
432  end if
433 
434  ! -- recompute coordinate
435  newx = x
436  if (lprofile) then
437  newx = newx + (v1 * dt / dx)
438  else if (v .ne. 0d0) then
439  newx = newx + (v * (exp(dvdx * dt) - 1.0d0) / dvdx / dx)
440  end if
441 
442  ! -- clamp to [0, 1]
443  if (newx .lt. 0d0) newx = 0d0
444  if (newx .gt. 1.0d0) newx = 1.0d0
445 
Here is the caller graph for this function:

◆ track_subcell()

subroutine methodsubcellpollockmodule::track_subcell ( class(methodsubcellpollocktype), intent(inout)  this,
class(subcellrecttype), intent(in)  subcell,
type(particletype), intent(inout), pointer  particle,
real(dp), intent(in)  tmax 
)
private

This subroutine consists partly of code written by David W. Pollock of the USGS for MODPATH 7. PRT's authors take responsibility for its application in this context and for any modifications or errors.

Definition at line 86 of file MethodSubcellPollock.f90.

87  ! modules
89  ! dummy
90  class(MethodSubcellPollockType), intent(inout) :: this
91  class(SubcellRectType), intent(in) :: subcell
92  type(ParticleType), pointer, intent(inout) :: particle
93  real(DP), intent(in) :: tmax
94  ! local
95  real(DP) :: vx
96  real(DP) :: dvxdx
97  real(DP) :: vy
98  real(DP) :: dvydy
99  real(DP) :: vz
100  real(DP) :: dvzdz
101  real(DP) :: dtexitx
102  real(DP) :: dtexity
103  real(DP) :: dtexitz
104  real(DP) :: dtexit
105  real(DP) :: texit
106  real(DP) :: dt
107  real(DP) :: t
108  real(DP) :: t0
109  real(DP) :: x
110  real(DP) :: y
111  real(DP) :: z
112  integer(I4B) :: statusVX
113  integer(I4B) :: statusVY
114  integer(I4B) :: statusVZ
115  integer(I4B) :: i
116  real(DP) :: initialX
117  real(DP) :: initialY
118  real(DP) :: initialZ
119  integer(I4B) :: exitFace
120  integer(I4B) :: reason
121  integer(I4B) :: tslice(2) !< user-time slice for the current time step
122 
123  reason = -1
124 
125  ! -- Initial particle location in scaled subcell coordinates
126  initialx = particle%x / subcell%dx
127  initialy = particle%y / subcell%dy
128  initialz = particle%z / subcell%dz
129 
130  ! -- Compute time of travel to each possible exit face
131  statusvx = calculate_dt(subcell%vx1, subcell%vx2, subcell%dx, &
132  initialx, vx, dvxdx, dtexitx)
133  statusvy = calculate_dt(subcell%vy1, subcell%vy2, subcell%dy, &
134  initialy, vy, dvydy, dtexity)
135  statusvz = calculate_dt(subcell%vz1, subcell%vz2, subcell%dz, &
136  initialz, vz, dvzdz, dtexitz)
137 
138  ! -- Subcells should never be strong sinks, contact the developer situation
139  if ((statusvx .eq. 3) .and. (statusvy .eq. 3) .and. (statusvz .eq. 3)) then
140  print *, "Subcell with no exit face:", &
141  "particle", get_particle_id(particle), &
142  "cell", particle%idomain(2)
143  call pstop(1)
144  end if
145 
146  ! -- Determine (earliest) exit face and corresponding travel time to exit
147  exitface = 0
148  dtexit = 1.0d+30
149  if ((statusvx .lt. 2) .or. (statusvy .lt. 2) .or. (statusvz .lt. 2)) then
150  ! -- Consider x-oriented faces
151  dtexit = dtexitx
152  if (vx .lt. 0d0) then
153  exitface = 1
154  else if (vx .gt. 0) then
155  exitface = 2
156  end if
157  ! -- Consider y-oriented faces
158  if (dtexity .lt. dtexit) then
159  dtexit = dtexity
160  if (vy .lt. 0d0) then
161  exitface = 3
162  else if (vy .gt. 0d0) then
163  exitface = 4
164  end if
165  end if
166  ! -- Consider z-oriented faces
167  if (dtexitz .lt. dtexit) then
168  dtexit = dtexitz
169  if (vz .lt. 0d0) then
170  exitface = 5
171  else if (vz .gt. 0d0) then
172  exitface = 6
173  end if
174  end if
175  else
176  end if
177 
178  ! -- Compute exit time
179  texit = particle%ttrack + dtexit
180  t0 = particle%ttrack
181 
182  ! -- Select user tracking times to solve. If this is the first time step
183  ! of the simulation, include all times before it begins; if it is the
184  ! last time step include all times after it ends, otherwise the times
185  ! within the current period and time step only.
186  call this%tracktimes%try_advance()
187  tslice = this%tracktimes%selection
188 
189  if (all(tslice > 0)) then
190  do i = tslice(1), tslice(2)
191  t = this%tracktimes%times(i)
192  if (t < particle%ttrack .or. t >= texit .or. t >= tmax) cycle
193  dt = t - t0
194  x = new_x(vx, dvxdx, subcell%vx1, subcell%vx2, &
195  dt, initialx, subcell%dx, statusvx == 1)
196  y = new_x(vy, dvydy, subcell%vy1, subcell%vy2, &
197  dt, initialy, subcell%dy, statusvy == 1)
198  z = new_x(vz, dvzdz, subcell%vz1, subcell%vz2, &
199  dt, initialz, subcell%dz, statusvz == 1)
200  particle%x = x * subcell%dx
201  particle%y = y * subcell%dy
202  particle%z = z * subcell%dz
203  particle%ttrack = t
204  particle%istatus = 1
205  call this%save(particle, reason=5)
206  end do
207  end if
208 
209  if (texit .gt. tmax) then
210  ! -- The computed exit time is greater than the maximum time, so set
211  ! -- final time for particle trajectory equal to maximum time and
212  ! -- calculate particle location at that final time.
213  t = tmax
214  dt = t - t0
215  x = new_x(vx, dvxdx, subcell%vx1, subcell%vx2, &
216  dt, initialx, subcell%dx, statusvx == 1)
217  y = new_x(vy, dvydy, subcell%vy1, subcell%vy2, &
218  dt, initialy, subcell%dy, statusvy == 1)
219  z = new_x(vz, dvzdz, subcell%vz1, subcell%vz2, &
220  dt, initialz, subcell%dz, statusvz == 1)
221  exitface = 0
222  particle%istatus = 1
223  particle%advancing = .false.
224  reason = 2 ! timestep end
225  else
226  ! -- The computed exit time is less than or equal to the maximum time,
227  ! -- so set final time for particle trajectory equal to exit time and
228  ! -- calculate exit location.
229  t = texit
230  dt = dtexit
231  if ((exitface .eq. 1) .or. (exitface .eq. 2)) then
232  x = 0d0
233  y = new_x(vy, dvydy, subcell%vy1, subcell%vy2, &
234  dt, initialy, subcell%dy, statusvy == 1)
235  z = new_x(vz, dvzdz, subcell%vz1, subcell%vz2, &
236  dt, initialz, subcell%dz, statusvz == 1)
237  if (exitface .eq. 2) x = 1.0d0
238  else if ((exitface .eq. 3) .or. (exitface .eq. 4)) then
239  x = new_x(vx, dvxdx, subcell%vx1, subcell%vx2, dt, &
240  initialx, subcell%dx, statusvx == 1)
241  y = 0d0
242  z = new_x(vz, dvzdz, subcell%vz1, subcell%vz2, dt, &
243  initialz, subcell%dz, statusvz == 1)
244  if (exitface .eq. 4) y = 1.0d0
245  else if ((exitface .eq. 5) .or. (exitface .eq. 6)) then
246  x = new_x(vx, dvxdx, subcell%vx1, subcell%vx2, &
247  dt, initialx, subcell%dx, statusvx == 1)
248  y = new_x(vy, dvydy, subcell%vy1, subcell%vy2, &
249  dt, initialy, subcell%dy, statusvy == 1)
250  z = 0d0
251  if (exitface .eq. 6) z = 1.0d0
252  else
253  print *, "programmer error, invalid exit face", exitface
254  call pstop(1)
255  end if
256  reason = 1 ! cell transition
257  end if
258 
259  ! -- Set final particle location in local (unscaled) subcell coordinates,
260  ! -- final time for particle trajectory, and exit face
261  particle%x = x * subcell%dx
262  particle%y = y * subcell%dy
263  particle%z = z * subcell%dz
264  particle%ttrack = t
265  particle%iboundary(3) = exitface
266 
267  ! -- Save particle track record
268  if (reason >= 0) call this%save(particle, reason=reason)
269 
pure character(len=lenmempath) function, public get_particle_id(particle)
Return the particle's composite ID.
Definition: Particle.f90:354
Here is the call graph for this function: