MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
MethodSubcellTernary.f90
Go to the documentation of this file.
2  use kindmodule, only: dp, i4b, lgp
4  use errorutilmodule, only: pstop
5  use geomutilmodule, only: clamp_bary, skew
6  use methodmodule, only: methodtype
7  use cellmodule, only: celltype
8  use subcellmodule, only: subcelltype
10  use particlemodule, only: particletype
12  use prtfmimodule, only: prtfmitype
13  use basedismodule, only: disbasetype
14  use mathutilmodule, only: is_close
15  implicit none
16 
17  private
18  public :: methodsubcellternarytype
20 
21  !> @brief Ternary triangular subcell tracking method.
23  integer(I4B), public, pointer :: zeromethod
24  contains
25  procedure, public :: apply => apply_mst
26  procedure, public :: deallocate
27  procedure, private :: track_subcell
29 
30 contains
31 
32  !> @brief Create a new ternary subcell tracking method.
33  subroutine create_method_subcell_ternary(method)
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%name => method%subcell%type
43  method%delegates = .false.
44  end subroutine create_method_subcell_ternary
45 
46  !> @brief Deallocate the ternary subcell tracking method.
47  subroutine deallocate (this)
48  class(methodsubcellternarytype), intent(inout) :: this
49  deallocate (this%name)
50  end subroutine deallocate
51 
52  !> @brief Apply the ternary subcell tracking method.
53  subroutine apply_mst(this, particle, tmax)
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
62  end subroutine apply_mst
63 
64  !> @brief Track a particle across a triangular subcell.
65  subroutine track_subcell(this, subcell, particle, tmax)
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 .and. dtexitz >= 0.0_dp) 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  else
231  particle%istatus = 9
232  particle%advancing = .false.
233  call this%save(particle, reason=3)
234  return
235  end if
236  texit = particle%ttrack + dtexit
237  t0 = particle%ttrack
238 
239  ! Select user tracking times to solve. If this is the first time step
240  ! of the simulation, include all times before it begins; if it is the
241  ! last time step include all times after it ends only if the 'extend'
242  ! option is on, otherwise times in this period and time step only.
243  call this%tracktimes%advance()
244  if (this%tracktimes%any()) then
245  do i = this%tracktimes%selection(1), this%tracktimes%selection(2)
246  t = this%tracktimes%times(i)
247  if (t < particle%ttrack) cycle
248  if (t >= texit .or. t >= tmax) exit
249  dt = t - t0
250  call calculate_xyz_position(dt, rxx, rxy, ryx, ryy, sxx, sxy, syy, &
251  izstatus, x0, y0, az, vzi, vzbot, &
252  ztop, zbot, zi, x, y, z)
253  particle%x = x
254  particle%y = y
255  particle%z = z
256  particle%ttrack = t
257  particle%istatus = 1
258  call this%save(particle, reason=5)
259  end do
260  end if
261 
262  ! Compute exit time and face and update the particle's coordinates
263  ! (local, unscaled) and other properties. The particle may at this
264  ! point lie on a boundary of the subcell or may still be within it.
265  if (texit .gt. tmax) then
266  ! The computed exit time is greater than the maximum time, so set
267  ! final time for particle trajectory equal to maximum time.
268  t = tmax
269  dt = t - t0
270  exitface = 0
271  particle%istatus = 1
272  particle%advancing = .false.
273  reason = 2 ! timestep end
274  else
275  ! The computed exit time is less than or equal to the maximum time,
276  ! so set final time for particle trajectory equal to exit time.
277  t = texit
278  dt = dtexit
279  reason = 1 ! (sub)cell transition
280  end if
281  call calculate_xyz_position(dt, rxx, rxy, ryx, ryy, sxx, sxy, syy, &
282  izstatus, x0, y0, az, vzi, vzbot, &
283  ztop, zbot, zi, x, y, z, exitface)
284  particle%x = x
285  particle%y = y
286  particle%z = z
287  particle%ttrack = t
288  particle%iboundary(3) = exitface
289 
290  call this%save(particle, reason=reason)
291  end subroutine track_subcell
292 
293  !> @brief Do calculations related to analytical z solution
294  !!
295  !! This subroutine consists partly or entirely of code written by
296  !! David W. Pollock of the USGS for MODPATH 7. The authors of the present
297  !! code are responsible for its appropriate application in this context
298  !! and for any modifications or errors.
299  !<
300  subroutine calculate_dt(v1, v2, dx, xL, v, dvdx, &
301  dt, status, itopbotexit)
302  real(DP) :: v1
303  real(DP) :: v2
304  real(DP) :: dx
305  real(DP) :: xL
306  real(DP) :: v
307  real(DP) :: dvdx
308  real(DP) :: dt
309  real(DP) :: v2a
310  real(DP) :: v1a
311  real(DP) :: dv
312  real(DP) :: dva
313  real(DP) :: vv
314  real(DP) :: vvv
315  real(DP) :: zro
316  real(DP) :: zrom
317  real(DP) :: x
318  real(DP) :: tol
319  real(DP) :: vr1
320  real(DP) :: vr2
321  real(DP) :: vr
322  real(DP) :: v1v2
323  integer(I4B) :: status
324  integer(I4B) :: itopbotexit
325  logical(LGP) :: noOutflow
326 
327  ! Initialize variables
328  status = -1
329  dt = 1.0d+20
330  v2a = v2
331  if (v2a .lt. dzero) v2a = -v2a
332  v1a = v1
333  if (v1a .lt. dzero) v1a = -v1a
334  dv = v2 - v1
335  dva = dv
336  if (dva .lt. dzero) dva = -dva
337 
338  ! Check for a uniform zero velocity in this direction.
339  ! If so, set status = 2 and return (dt = 1.0d+20).
340  tol = 1.0d-15
341  if ((v2a .lt. tol) .and. (v1a .lt. tol)) then
342  v = dzero
343  dvdx = dzero
344  status = 2
345  itopbotexit = 0
346  return
347  end if
348 
349  ! Check for uniform non-zero velocity in this direction.
350  ! If so, set compute dt using the constant velocity,
351  ! set status = 1 and return.
352  vv = v1a
353  if (v2a .gt. vv) vv = v2a
354  vvv = dva / vv
355  if (vvv .lt. 1.0d-4) then
356  zro = tol
357  zrom = -zro
358  v = v1
359  x = xl * dx
360  if (v1 .gt. zro) then
361  dt = (dx - x) / v1
362  itopbotexit = -2
363  end if
364  if (v1 .lt. zrom) then
365  dt = -x / v1
366  itopbotexit = -1
367  end if
368  dvdx = dzero
369  status = 1
370  return
371  end if
372 
373  ! Velocity has a linear variation.
374  ! Compute velocity corresponding to particle position
375  dvdx = dv / dx
376  v = (done - xl) * v1 + xl * v2
377 
378  ! If flow is into the cell from both sides there is no outflow.
379  ! In that case, set status = 3 and return
380  nooutflow = .true.
381  if (v1 .lt. dzero) nooutflow = .false.
382  if (v2 .gt. dzero) nooutflow = .false.
383  if (nooutflow) then
384  status = 3
385  itopbotexit = 0
386  return
387  end if
388 
389  ! If there is a divide in the cell for this flow direction, check to see if the
390  ! particle is located exactly on the divide. If it is, move it very slightly to
391  ! get it off the divide. This avoids possible numerical problems related to
392  ! stagnation points.
393  if ((v1 .le. dzero) .and. (v2 .ge. dzero)) then
394  if (abs(v) .le. dzero) then
395  v = 1.0d-20
396  if (v2 .le. dzero) v = -v
397  end if
398  end if
399 
400  ! If there is a flow divide, find out what side of the divide the particle
401  ! is on and set the value of vr appropriately to reflect that location.
402  vr1 = v1 / v
403  vr2 = v2 / v
404  vr = vr1
405  itopbotexit = -1
406  if (vr .le. dzero) then
407  vr = vr2
408  itopbotexit = -2
409  end if
410 
411  ! Check if velocity is in the same direction throughout cell (i.e. no flow divide).
412  ! Check if product v1*v2 > 0 then the velocity is in the same direction throughout
413  ! the cell (i.e. no flow divide). If so, set vr to reflect appropriate direction.
414  v1v2 = v1 * v2
415  if (v1v2 .gt. dzero) then
416  if (v .gt. dzero) then
417  vr = vr2
418  itopbotexit = -2
419  end if
420  if (v .lt. dzero) then
421  vr = vr1
422  itopbotexit = -1
423  end if
424  end if
425 
426  ! Compute travel time to exit face. Return with status = 0
427  dt = log(abs(vr)) / dvdx
428  status = 0
429  end subroutine calculate_dt
430 
431  !> @brief Calculate the particle's local unscaled xyz coordinates after dt.
432  subroutine calculate_xyz_position(dt, rxx, rxy, ryx, ryy, sxx, sxy, syy, &
433  izstatus, x0, y0, az, vzi, vzbot, &
434  ztop, zbot, zi, x, y, z, exitFace)
435  ! dummy
436  real(DP) :: dt
437  real(DP) :: rxx
438  real(DP) :: rxy
439  real(DP) :: ryx
440  real(DP) :: ryy
441  real(DP) :: sxx
442  real(DP) :: sxy
443  real(DP) :: syy
444  integer(I4B) :: izstatus
445  real(DP) :: x0
446  real(DP) :: y0
447  real(DP) :: az
448  real(DP) :: vzi
449  real(DP) :: vzbot
450  real(DP) :: ztop
451  real(DP) :: zbot
452  real(DP) :: zi
453  real(DP) :: x
454  real(DP) :: y
455  real(DP) :: z
456  integer(I4B), optional :: exitFace
457  ! local
458  integer(I4B) :: lexitface
459  real(DP) :: rot(2, 2), res(2), loc(2)
460  real(DP) :: alp
461  real(DP) :: bet
462 
463  ! process optional exit face argument
464  if (present(exitface)) then
465  lexitface = exitface
466  else
467  lexitface = 0
468  end if
469 
470  ! calculate alpha and beta
471  call step_analytical(dt, alp, bet)
472 
473  ! if exit face is known, set alpha or beta coordinate
474  ! corresponding to the exit face exactly.
475  if (lexitface .eq. 1) then
476  bet = dzero
477  else if (lexitface .eq. 2) then
478  alp = done - bet
479  else if (lexitface .eq. 3) then
480  alp = dzero
481  end if
482 
483  ! if exit face is top or bottom, set z coordinate exactly.
484  if (lexitface .eq. 4) then
485  z = zbot
486  else if (lexitface .eq. 5) then
487  z = ztop
488  else
489  ! otherwise calculate z.
490  if (izstatus .eq. 2) then
491  ! vz uniformly zero
492  z = zi
493  else if (izstatus .eq. 1) then
494  ! vz uniform, nonzero
495  z = zi + vzi * dt
496  else
497  ! vz nonuniform
498  z = zbot + (vzi * dexp(az * dt) - vzbot) / az
499  end if
500  end if
501 
502  ! transform (alp, beta) to (x, y)
503  loc = (/alp, bet/)
504  loc = skew(loc, (/sxx, sxy, syy/), invert=.true.)
505  rot = reshape((/rxx, rxy, ryx, ryy/), shape(rot))
506  res = matmul(rot, loc) ! rotate vector
507  x = res(1) + x0
508  y = res(2) + y0
509 
510  end subroutine calculate_xyz_position
511 
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dsame
real constant for values that are considered the same based on machine precision
Definition: Constants.f90:122
real(dp), parameter dep3
real constant 1000
Definition: Constants.f90:88
real(dp), parameter donethird
real constant 1/3
Definition: Constants.f90:67
real(dp), parameter dhalf
real constant 1/2
Definition: Constants.f90:68
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter dtwo
real constant 2
Definition: Constants.f90:79
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
subroutine, public clamp_bary(alpha, beta, gamma, pad)
Clamp barycentric coordinates to the interior of a triangle, with optional padding some minimum dista...
Definition: GeomUtil.f90:466
pure real(dp) function, dimension(2), public skew(v, s, invert)
Skew a 2D vector along the x-axis.
Definition: GeomUtil.f90:149
This module defines variable data types.
Definition: kind.f90:8
pure logical function, public is_close(a, b, rtol, atol, symmetric)
Check if a real value is approximately equal to another.
Definition: MathUtil.f90:46
Particle tracking strategies.
Definition: Method.f90:2
subroutine track_subcell(this, subcell, particle, tmax)
Track a particle across a triangular subcell.
subroutine apply_mst(this, particle, tmax)
Apply the ternary subcell tracking method.
subroutine, public create_method_subcell_ternary(method)
Create a new ternary subcell tracking method.
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.
subroutine calculate_dt(v1, v2, dx, xL, v, dvdx, dt, status, itopbotexit)
Do calculations related to analytical z solution.
subroutine, public create_subcell_tri(subcell)
Create a new triangular subcell.
Definition: SubcellTri.f90:26
subroutine, public canonical(x0, y0, x1, y1, x2, y2, v0x, v0y, v1x, v1y, v2x, v2y, xi, yi, rxx, rxy, ryx, ryy, sxx, sxy, syy, alp0, bet0, alp1, bet1, alp2, bet2, alpi, beti)
Set coordinates to "canonical" configuration.
subroutine, public traverse_triangle(isolv, tol, texit, alpexit, betexit, itrifaceenter, itrifaceexit, alp1, bet1, alp2, bet2, alpi, beti)
Traverse triangular cell.
subroutine, public step_analytical(t, alp, bet)
Step (evaluate) analytically depending on case.
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:31
Ternary triangular subcell tracking method.
Particle tracked by the PRT model.
Definition: Particle.f90:32
A subcell of a cell.
Definition: Subcell.f90:9