MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
methodsubcellternarymodule Module Reference

Data Types

type  methodsubcellternarytype
 Ternary triangular subcell tracking method. More...
 

Functions/Subroutines

subroutine, public create_method_subcell_ternary (method)
 Create a new ternary subcell tracking method. More...
 
subroutine deallocate (this)
 Deallocate the ternary subcell tracking method. More...
 
subroutine apply_mst (this, particle, tmax)
 Apply the ternary subcell tracking method. More...
 
subroutine track_subcell (this, subcell, particle, tmax)
 Track a particle across a triangular subcell. More...
 
subroutine calculate_dt (v1, v2, dx, xL, v, dvdx, dt, status, itopbotexit)
 Do calculations related to analytical z solution. More...
 
subroutine calculate_xyz_position (dt, rxx, rxy, ryx, ryy, sxx, sxy, syy, izstatus, x0, y0, az, vzi, vzbot, ztop, zbot, zi, x, y, z, exitFace)
 Calculate the particle's local unscaled xyz coordinates after dt. More...
 

Function/Subroutine Documentation

◆ apply_mst()

subroutine methodsubcellternarymodule::apply_mst ( class(methodsubcellternarytype), intent(inout)  this,
type(particletype), intent(inout), pointer  particle,
real(dp), intent(in)  tmax 
)
private

Definition at line 53 of file MethodSubcellTernary.f90.

54  class(MethodSubcellTernaryType), intent(inout) :: this
55  type(ParticleType), pointer, intent(inout) :: particle
56  real(DP), intent(in) :: tmax
57 
58  select type (subcell => this%subcell)
59  type is (subcelltritype)
60  call this%track_subcell(subcell, particle, tmax)
61  end select

◆ calculate_dt()

subroutine methodsubcellternarymodule::calculate_dt ( real(dp)  v1,
real(dp)  v2,
real(dp)  dx,
real(dp)  xL,
real(dp)  v,
real(dp)  dvdx,
real(dp)  dt,
integer(i4b)  status,
integer(i4b)  itopbotexit 
)
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 294 of file MethodSubcellTernary.f90.

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  real(DP) :: v2a
304  real(DP) :: v1a
305  real(DP) :: dv
306  real(DP) :: dva
307  real(DP) :: vv
308  real(DP) :: vvv
309  real(DP) :: zro
310  real(DP) :: zrom
311  real(DP) :: x
312  real(DP) :: tol
313  real(DP) :: vr1
314  real(DP) :: vr2
315  real(DP) :: vr
316  real(DP) :: v1v2
317  integer(I4B) :: status
318  integer(I4B) :: itopbotexit
319  logical(LGP) :: noOutflow
320 
321  ! Initialize variables
322  status = -1
323  dt = 1.0d+20
324  v2a = v2
325  if (v2a .lt. dzero) v2a = -v2a
326  v1a = v1
327  if (v1a .lt. dzero) v1a = -v1a
328  dv = v2 - v1
329  dva = dv
330  if (dva .lt. dzero) dva = -dva
331 
332  ! Check for a uniform zero velocity in this direction.
333  ! If so, set status = 2 and return (dt = 1.0d+20).
334  tol = 1.0d-15
335  if ((v2a .lt. tol) .and. (v1a .lt. tol)) then
336  v = dzero
337  dvdx = dzero
338  status = 2
339  itopbotexit = 0
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) then
355  dt = (dx - x) / v1
356  itopbotexit = -2
357  end if
358  if (v1 .lt. zrom) then
359  dt = -x / v1
360  itopbotexit = -1
361  end if
362  dvdx = dzero
363  status = 1
364  return
365  end if
366 
367  ! Velocity has a linear variation.
368  ! Compute velocity corresponding to particle position
369  dvdx = dv / dx
370  v = (done - xl) * v1 + xl * v2
371 
372  ! If flow is into the cell from both sides there is no outflow.
373  ! In that case, set status = 3 and return
374  nooutflow = .true.
375  if (v1 .lt. dzero) nooutflow = .false.
376  if (v2 .gt. dzero) nooutflow = .false.
377  if (nooutflow) then
378  status = 3
379  itopbotexit = 0
380  return
381  end if
382 
383  ! If there is a divide in the cell for this flow direction, check to see if the
384  ! particle is located exactly on the divide. If it is, move it very slightly to
385  ! get it off the divide. This avoids possible numerical problems related to
386  ! stagnation points.
387  if ((v1 .le. dzero) .and. (v2 .ge. dzero)) then
388  if (abs(v) .le. dzero) then
389  v = 1.0d-20
390  if (v2 .le. dzero) v = -v
391  end if
392  end if
393 
394  ! If there is a flow divide, find out what side of the divide the particle
395  ! is on and set the value of vr appropriately to reflect that location.
396  vr1 = v1 / v
397  vr2 = v2 / v
398  vr = vr1
399  itopbotexit = -1
400  if (vr .le. dzero) then
401  vr = vr2
402  itopbotexit = -2
403  end if
404 
405  ! Check if velocity is in the same direction throughout cell (i.e. no flow divide).
406  ! Check if product v1*v2 > 0 then the velocity is in the same direction throughout
407  ! the cell (i.e. no flow divide). If so, set vr to reflect appropriate direction.
408  v1v2 = v1 * v2
409  if (v1v2 .gt. dzero) then
410  if (v .gt. dzero) then
411  vr = vr2
412  itopbotexit = -2
413  end if
414  if (v .lt. dzero) then
415  vr = vr1
416  itopbotexit = -1
417  end if
418  end if
419 
420  ! Compute travel time to exit face. Return with status = 0
421  dt = log(vr) / dvdx
422  status = 0
Here is the caller graph for this function:

◆ calculate_xyz_position()

subroutine methodsubcellternarymodule::calculate_xyz_position ( real(dp)  dt,
real(dp)  rxx,
real(dp)  rxy,
real(dp)  ryx,
real(dp)  ryy,
real(dp)  sxx,
real(dp)  sxy,
real(dp)  syy,
integer(i4b)  izstatus,
real(dp)  x0,
real(dp)  y0,
real(dp)  az,
real(dp)  vzi,
real(dp)  vzbot,
real(dp)  ztop,
real(dp)  zbot,
real(dp)  zi,
real(dp)  x,
real(dp)  y,
real(dp)  z,
integer(i4b), optional  exitFace 
)
private

Definition at line 426 of file MethodSubcellTernary.f90.

429  ! dummy
430  real(DP) :: dt
431  real(DP) :: rxx
432  real(DP) :: rxy
433  real(DP) :: ryx
434  real(DP) :: ryy
435  real(DP) :: sxx
436  real(DP) :: sxy
437  real(DP) :: syy
438  integer(I4B) :: izstatus
439  real(DP) :: x0
440  real(DP) :: y0
441  real(DP) :: az
442  real(DP) :: vzi
443  real(DP) :: vzbot
444  real(DP) :: ztop
445  real(DP) :: zbot
446  real(DP) :: zi
447  real(DP) :: x
448  real(DP) :: y
449  real(DP) :: z
450  integer(I4B), optional :: exitFace
451  ! local
452  integer(I4B) :: lexitface
453  real(DP) :: rot(2, 2), res(2), loc(2)
454  real(DP) :: alp
455  real(DP) :: bet
456 
457  ! process optional exit face argument
458  if (present(exitface)) then
459  lexitface = exitface
460  else
461  lexitface = 0
462  end if
463 
464  ! calculate alpha and beta
465  call step_analytical(dt, alp, bet)
466 
467  ! if exit face is known, set alpha or beta coordinate
468  ! corresponding to the exit face exactly.
469  if (lexitface .eq. 1) then
470  bet = dzero
471  else if (lexitface .eq. 2) then
472  alp = done - bet
473  else if (lexitface .eq. 3) then
474  alp = dzero
475  end if
476 
477  ! if exit face is top or bottom, set z coordinate exactly.
478  if (lexitface .eq. 4) then
479  z = zbot
480  else if (lexitface .eq. 5) then
481  z = ztop
482  else
483  ! otherwise calculate z.
484  if (izstatus .eq. 2) then
485  ! -- vz uniformly zero
486  z = zi
487  else if (izstatus .eq. 1) then
488  ! -- vz uniform, nonzero
489  z = zi + vzi * dt
490  else
491  ! -- vz nonuniform
492  z = zbot + (vzi * dexp(az * dt) - vzbot) / az
493  end if
494  end if
495 
496  ! transform (alp, beta) to (x, y)
497  loc = (/alp, bet/)
498  loc = skew(loc, (/sxx, sxy, syy/), invert=.true.)
499  rot = reshape((/rxx, rxy, ryx, ryy/), shape(rot))
500  res = matmul(rot, loc) ! rotate vector
501  x = res(1) + x0
502  y = res(2) + y0
503 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ create_method_subcell_ternary()

subroutine, public methodsubcellternarymodule::create_method_subcell_ternary ( type(methodsubcellternarytype), pointer  method)

Definition at line 33 of file MethodSubcellTernary.f90.

34  ! -- dummy
35  type(MethodSubcellTernaryType), pointer :: method
36  ! -- local
37  type(SubcellTriType), pointer :: subcell
38 
39  allocate (method)
40  call create_subcell_tri(subcell)
41  method%subcell => subcell
42  method%type => method%subcell%type
43  method%delegates = .false.
Here is the call graph for this function:
Here is the caller graph for this function:

◆ deallocate()

subroutine methodsubcellternarymodule::deallocate ( class(methodsubcellternarytype), intent(inout)  this)
private

Definition at line 47 of file MethodSubcellTernary.f90.

48  class(MethodSubcellTernaryType), intent(inout) :: this
49  deallocate (this%type)

◆ track_subcell()

subroutine methodsubcellternarymodule::track_subcell ( class(methodsubcellternarytype), intent(inout)  this,
class(subcelltritype), intent(in)  subcell,
type(particletype), intent(inout), pointer  particle,
real(dp), intent(in)  tmax 
)
private

Definition at line 65 of file MethodSubcellTernary.f90.

66  ! dummy
67  class(MethodSubcellTernaryType), intent(inout) :: this
68  class(SubcellTriType), intent(in) :: subcell
69  type(ParticleType), pointer, intent(inout) :: particle
70  real(DP), intent(in) :: tmax
71  ! local
72  integer(I4B) :: exitFace
73  real(DP) :: x0
74  real(DP) :: y0
75  real(DP) :: x1
76  real(DP) :: y1
77  real(DP) :: x2
78  real(DP) :: y2
79  real(DP) :: v0x
80  real(DP) :: v0y
81  real(DP) :: v1x
82  real(DP) :: v1y
83  real(DP) :: v2x
84  real(DP) :: v2y
85  real(DP) :: xi
86  real(DP) :: yi
87  real(DP) :: zi
88  real(DP) :: zirel
89  real(DP) :: ztop
90  real(DP) :: zbot
91  real(DP) :: dz
92  real(DP) :: rxx
93  real(DP) :: rxy
94  real(DP) :: ryx
95  real(DP) :: ryy
96  real(DP) :: sxx
97  real(DP) :: sxy
98  real(DP) :: syy
99  real(DP) :: alp0
100  real(DP) :: bet0
101  real(DP) :: alp1
102  real(DP) :: bet1
103  real(DP) :: alp2
104  real(DP) :: bet2
105  real(DP) :: alpi
106  real(DP) :: beti
107  real(DP) :: gami
108  real(DP) :: vzbot
109  real(DP) :: vztop
110  real(DP) :: vzi
111  real(DP) :: vziodz
112  real(DP) :: az
113  real(DP) :: dtexitz
114  real(DP) :: dt
115  real(DP) :: t
116  real(DP) :: t0
117  real(DP) :: dtexitxy
118  real(DP) :: texit
119  real(DP) :: x
120  real(DP) :: y
121  real(DP) :: z
122  integer(I4B) :: izstatus
123  integer(I4B) :: itopbotexit
124  integer(I4B) :: ntmax
125  integer(I4B) :: isolv
126  integer(I4B) :: itrifaceenter
127  integer(I4B) :: itrifaceexit
128  real(DP) :: tol
129  real(DP) :: dtexit
130  real(DP) :: alpexit
131  real(DP) :: betexit
132  integer(I4B) :: reason
133  integer(I4B) :: i
134 
135  ! Set solution method
136  if (particle%iexmeth == 0) then
137  isolv = 1 ! default to Brent's
138  else
139  isolv = particle%iexmeth
140  end if
141 
142  ntmax = 10000
143  tol = particle%extol
144  reason = -1
145 
146  ! Set some local variables for convenience.
147  xi = particle%x
148  yi = particle%y
149  zi = particle%z
150  x0 = subcell%x0
151  y0 = subcell%y0
152  x1 = subcell%x1
153  y1 = subcell%y1
154  x2 = subcell%x2
155  y2 = subcell%y2
156  v0x = subcell%v0x
157  v0y = subcell%v0y
158  v1x = subcell%v1x
159  v1y = subcell%v1y
160  v2x = subcell%v2x
161  v2y = subcell%v2y
162  zbot = subcell%zbot
163  ztop = subcell%ztop
164  dz = subcell%dz
165  vzbot = subcell%vzbot
166  vztop = subcell%vztop
167 
168  ! Transform coordinates to the "canonical" configuration:
169  ! barycentric in two dimensions with alpha, beta & gamma
170  ! such that at f2 alpha = 0, f0 beta = 0, f1 gamma = 0.
171  !
172  ! v2
173  ! |\
174  ! f2| \f1
175  ! |__\
176  ! v0 f0 v1
177  !
178  call canonical(x0, y0, x1, y1, x2, y2, &
179  v0x, v0y, v1x, v1y, v2x, v2y, &
180  xi, yi, &
181  rxx, rxy, ryx, ryy, &
182  sxx, sxy, syy, &
183  alp0, bet0, alp1, bet1, alp2, bet2, alpi, beti)
184 
185  ! Clamp particle coordinates to the canonical triangular
186  ! subcell and nudge it ever so slightly inside if needed.
187  call clamp_bary(alpi, beti, gami, pad=dsame * dep3)
188 
189  ! Do calculations related to the analytical z solution.
190  ! todo: just once for each cell? store at cell-level?
191  zirel = (zi - zbot) / dz
192  call calculate_dt(vzbot, vztop, dz, zirel, vzi, &
193  az, dtexitz, izstatus, &
194  itopbotexit)
195  vziodz = vzi / dz
196 
197  ! If possible, track the particle across the subcell.
198  itrifaceenter = particle%iboundary(3) - 1
199  if (itrifaceenter == -1) itrifaceenter = 999
200  call traverse_triangle(isolv, tol, &
201  dtexitxy, alpexit, betexit, &
202  itrifaceenter, itrifaceexit, &
203  alp1, bet1, alp2, bet2, alpi, beti)
204 
205  ! If the subcell has no exit face, terminate the particle.
206  ! todo: after initial release, consider ramifications
207  if (itopbotexit == 0 .and. itrifaceexit == 0) then
208  particle%istatus = 9
209  particle%advancing = .false.
210  call this%save(particle, reason=3)
211  return
212  end if
213 
214  ! Determine the particle's exit face and travel time to exit.
215  ! The exit face is the face through which it would exit first,
216  ! considering only the velocity component in the direction of
217  ! the face. Then compute the particle's exit time.
218  if (itrifaceexit /= 0) then
219  ! Exit through lateral subcell face
220  exitface = itrifaceexit
221  dtexit = dtexitxy
222  else if (dtexitz < dtexitxy) then
223  ! Exit through top or bottom
224  if (itopbotexit == -1) then
225  exitface = 4
226  else
227  exitface = 5
228  end if
229  dtexit = dtexitz
230  end if
231  texit = particle%ttrack + dtexit
232  t0 = particle%ttrack
233 
234  ! Select user tracking times to solve. If this is the first time step
235  ! of the simulation, include all times before it begins; if it is the
236  ! last time step include all times after it ends only if the 'extend'
237  ! option is on, otherwise times in this period and time step only.
238  call this%tracktimes%advance()
239  if (this%tracktimes%any()) then
240  do i = this%tracktimes%selection(1), this%tracktimes%selection(2)
241  t = this%tracktimes%times(i)
242  if (t < particle%ttrack) cycle
243  if (t >= texit .or. t >= tmax) exit
244  dt = t - t0
245  call calculate_xyz_position(dt, rxx, rxy, ryx, ryy, sxx, sxy, syy, &
246  izstatus, x0, y0, az, vzi, vzbot, &
247  ztop, zbot, zi, x, y, z)
248  particle%x = x
249  particle%y = y
250  particle%z = z
251  particle%ttrack = t
252  particle%istatus = 1
253  call this%save(particle, reason=5)
254  end do
255  end if
256 
257  ! Compute exit time and face and update the particle's coordinates
258  ! (local, unscaled) and other properties. The particle may at this
259  ! point lie on a boundary of the subcell or may still be within it.
260  if (texit .gt. tmax) then
261  ! -- The computed exit time is greater than the maximum time, so set
262  ! -- final time for particle trajectory equal to maximum time.
263  t = tmax
264  dt = t - t0
265  exitface = 0
266  particle%istatus = 1
267  particle%advancing = .false.
268  reason = 2 ! timestep end
269  else
270  ! -- The computed exit time is less than or equal to the maximum time,
271  ! -- so set final time for particle trajectory equal to exit time.
272  t = texit
273  dt = dtexit
274  reason = 1 ! (sub)cell transition
275  end if
276  call calculate_xyz_position(dt, rxx, rxy, ryx, ryy, sxx, sxy, syy, &
277  izstatus, x0, y0, az, vzi, vzbot, &
278  ztop, zbot, zi, x, y, z, exitface)
279  particle%x = x
280  particle%y = y
281  particle%z = z
282  particle%ttrack = t
283  particle%iboundary(3) = exitface
284  call this%save(particle, reason=reason)
Here is the call graph for this function: