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 52 of file MethodSubcellTernary.f90.

53  class(MethodSubcellTernaryType), intent(inout) :: this
54  type(ParticleType), pointer, intent(inout) :: particle
55  real(DP), intent(in) :: tmax
56 
57  select type (subcell => this%subcell)
58  type is (subcelltritype)
59  call this%track_subcell(subcell, particle, tmax)
60  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 293 of file MethodSubcellTernary.f90.

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

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

33  ! -- dummy
34  type(MethodSubcellTernaryType), pointer :: method
35  ! -- local
36  type(SubcellTriType), pointer :: subcell
37 
38  allocate (method)
39  call create_subcell_tri(subcell)
40  method%subcell => subcell
41  method%name => method%subcell%type
42  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 46 of file MethodSubcellTernary.f90.

47  class(MethodSubcellTernaryType), intent(inout) :: this
48  deallocate (this%name)

◆ 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 64 of file MethodSubcellTernary.f90.

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