MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
MethodSubcellTernary.f90
Go to the documentation of this file.
2  use kindmodule, only: dp, i4b, lgp
3  use errorutilmodule, only: pstop
4  use geomutilmodule, only: skew
5  use methodmodule, only: methodtype
6  use cellmodule, only: celltype
7  use subcellmodule, only: subcelltype
12  use prtfmimodule, only: prtfmitype
13  use basedismodule, only: disbasetype
14  implicit none
15 
16  private
17  public :: methodsubcellternarytype
19 
20  !> @brief Ternary triangular subcell tracking method
22  integer(I4B), public, pointer :: zeromethod
23  contains
24  procedure, public :: apply => apply_mst
25  procedure, public :: destroy
26  procedure, private :: track_subcell
28 
29 contains
30 
31  !> @brief Create a new ternary subcell-method object
32  subroutine create_method_subcell_ternary(method)
33  ! -- dummy
34  type(methodsubcellternarytype), pointer :: method
35  ! -- local
36  type(subcelltritype), pointer :: subcell
37 
38  allocate (method)
39  allocate (method%zeromethod)
40  call create_subcell_tri(subcell)
41  method%subcell => subcell
42  method%type => method%subcell%type
43  method%delegates = .false.
44  method%zeromethod = 0
45  end subroutine create_method_subcell_ternary
46 
47  !> @brief Destructor for a ternary subcell-method object
48  subroutine destroy(this)
49  class(methodsubcellternarytype), intent(inout) :: this
50  deallocate (this%type)
51  end subroutine destroy
52 
53  !> @brief Apply the ternary subcell method
54  subroutine apply_mst(this, particle, tmax)
55  class(methodsubcellternarytype), intent(inout) :: this
56  type(particletype), pointer, intent(inout) :: particle
57  real(DP), intent(in) :: tmax
58 
59  select type (subcell => this%subcell)
60  type is (subcelltritype)
61  call this%track_subcell(subcell, particle, tmax)
62  end select
63  end subroutine apply_mst
64 
65  !> @brief Track a particle across a triangular subcell using the ternary method
66  subroutine track_subcell(this, subcell, particle, tmax)
67  ! dummy
68  class(methodsubcellternarytype), intent(inout) :: this
69  class(subcelltritype), intent(in) :: subcell
70  type(particletype), pointer, intent(inout) :: particle
71  real(DP), intent(in) :: tmax
72  ! local
73  integer(I4B) :: exitFace
74  logical(LGP) :: lbary ! kluge
75  real(DP) :: x0
76  real(DP) :: y0
77  real(DP) :: x1
78  real(DP) :: y1
79  real(DP) :: x2
80  real(DP) :: y2
81  real(DP) :: v0x
82  real(DP) :: v0y
83  real(DP) :: v1x
84  real(DP) :: v1y
85  real(DP) :: v2x
86  real(DP) :: v2y
87  real(DP) :: xi
88  real(DP) :: yi
89  real(DP) :: zi
90  real(DP) :: zirel
91  real(DP) :: ztop
92  real(DP) :: zbot
93  real(DP) :: dz
94  real(DP) :: rxx
95  real(DP) :: rxy
96  real(DP) :: ryx
97  real(DP) :: ryy
98  real(DP) :: sxx
99  real(DP) :: sxy
100  real(DP) :: syy
101  real(DP) :: rot(2, 2), res(2), loc(2)
102  real(DP) :: alp
103  real(DP) :: bet
104  real(DP) :: alp0
105  real(DP) :: bet0
106  real(DP) :: alp1
107  real(DP) :: bet1
108  real(DP) :: alp2
109  real(DP) :: bet2
110  real(DP) :: alpi
111  real(DP) :: beti
112  real(DP) :: vzbot
113  real(DP) :: vztop
114  real(DP) :: vzi
115  real(DP) :: vziodz
116  real(DP) :: az
117  real(DP) :: dtexitz
118  real(DP) :: dt
119  real(DP) :: t
120  real(DP) :: t0
121  real(DP) :: dtexitxy
122  real(DP) :: texit
123  real(DP) :: x
124  real(DP) :: y
125  real(DP) :: z
126  integer(I4B) :: izstatus
127  integer(I4B) :: itopbotexit
128  integer(I4B) :: ntmax
129  integer(I4B) :: nsave
130  integer(I4B) :: isolv
131  integer(I4B) :: itrifaceenter
132  integer(I4B) :: itrifaceexit
133  real(DP) :: tol
134  real(DP) :: step
135  real(DP) :: dtexit
136  real(DP) :: alpexit
137  real(DP) :: betexit
138  integer(I4B) :: ntdebug ! kluge
139  integer(I4B) :: reason
140  integer(I4B) :: i
141  integer(I4B) :: tslice(2)
142 
143  lbary = .true. ! kluge
144  ntmax = 10000
145  nsave = 1 ! needed???
146  isolv = this%zeromethod
147  tol = 1d-7
148  step = 1e-3 ! needed only for euler
149  reason = -1
150 
151  ! -- Set some local variables for convenience
152  xi = particle%x
153  yi = particle%y
154  zi = particle%z
155  x0 = subcell%x0
156  y0 = subcell%y0
157  x1 = subcell%x1
158  y1 = subcell%y1
159  x2 = subcell%x2
160  y2 = subcell%y2
161  v0x = subcell%v0x
162  v0y = subcell%v0y
163  v1x = subcell%v1x
164  v1y = subcell%v1y
165  v2x = subcell%v2x
166  v2y = subcell%v2y
167  zbot = subcell%zbot
168  ztop = subcell%ztop
169  dz = subcell%dz
170  vzbot = subcell%vzbot
171  vztop = subcell%vztop
172 
173  ! -- Translate and rotate coordinates to "canonical" configuration
174  call canonical(x0, y0, x1, y1, x2, y2, &
175  v0x, v0y, v1x, v1y, v2x, v2y, &
176  xi, yi, &
177  rxx, rxy, ryx, ryy, &
178  sxx, sxy, syy, &
179  alp0, bet0, alp1, bet1, alp2, bet2, alpi, beti, &
180  lbary)
181 
182  ! -- Do calculations related to analytical z solution, which can be done
183  ! -- after traverse_triangle call if results not needed for adaptive time
184  ! -- stepping during triangle (subcell) traversal
185  ! kluge note: actually, can probably do z calculation just once for each cell
186  zirel = (zi - zbot) / dz
187  call calculate_dt(vzbot, vztop, dz, zirel, vzi, &
188  az, dtexitz, izstatus, &
189  itopbotexit)
190  vziodz = vzi / dz
191 
192  ! -- Traverse triangular subcell
193  ntdebug = -999 ! kluge debug bludebug
194  itrifaceenter = particle%iboundary(3) - 1
195  if (itrifaceenter .eq. -1) itrifaceenter = 999
196 
197  ! kluge note: can probably avoid calculating alpexit
198  ! here in many cases and wait to calculate it later,
199  ! once the final trajectory time is known
200  call traverse_triangle(isolv, tol, step, &
201  dtexitxy, alpexit, betexit, &
202  itrifaceenter, itrifaceexit, &
203  rxx, rxy, ryx, ryy, &
204  alp0, bet0, alp1, bet1, alp2, bet2, alpi, beti, &
205  vziodz, az, lbary)
206 
207  ! -- Subcell with no exit face, terminate the particle
208  if ((itopbotexit .eq. 0) .and. (itrifaceexit .eq. 0)) then
209  call this%save(particle, reason=9)
210  return
211  end if
212 
213  ! -- Determine (earliest) exit face and corresponding travel time to exit
214  if (itopbotexit .eq. 0) then
215  ! -- Exits through triangle face first
216  exitface = itrifaceexit
217  dtexit = dtexitxy
218  else if (itrifaceexit .eq. 0) then
219  ! -- Exits through top/bottom first
220  exitface = 45
221  dtexit = dtexitz
222  else if (dtexitz .lt. dtexitxy) then
223  ! -- Exits through top/bottom first
224  exitface = 45
225  dtexit = dtexitz
226  else
227  ! -- Exits through triangle face first
228  exitface = itrifaceexit
229  dtexit = dtexitxy
230  end if
231  if (exitface .eq. 45) then
232  if (itopbotexit .eq. -1) then
233  exitface = 4
234  else
235  exitface = 5
236  end if
237  end if
238 
239  ! -- Compute exit time
240  texit = particle%ttrack + dtexit
241  t0 = particle%ttrack
242 
243  ! -- Select user tracking times to solve. If this is the first time step
244  ! of the simulation, include all times before it begins; if it is the
245  ! last time step, include all times after it ends. Otherwise take the
246  ! times within the current period and time step only.
247  call this%tracktimes%try_advance()
248  tslice = this%tracktimes%selection
249  if (all(tslice > 0)) then
250  do i = tslice(1), tslice(2)
251  t = this%tracktimes%times(i)
252  if (t < particle%ttrack .or. t >= texit .or. t >= tmax) cycle
253  dt = t - t0
254  call step_analytical(dt, alp, bet)
255  loc = (/alp, bet/)
256  if (lbary) loc = skew(loc, (/sxx, sxy, syy/), invert=.true.)
257  rot = reshape((/rxx, rxy, ryx, ryy/), shape(rot))
258  res = matmul(rot, loc) ! rotate vector
259  x = res(1) + x0
260  y = res(2) + y0
261  ! kluge note: make this into a function
262  if (izstatus .eq. 2) then
263  ! -- vz uniformly zero
264  z = zi
265  else if (izstatus .eq. 1) then
266  ! -- vz uniform, nonzero
267  z = zi + vzi * dt
268  else
269  ! -- vz nonuniform
270  z = zbot + (vzi * dexp(az * dt) - vzbot) / az
271  end if
272  particle%x = x
273  particle%y = y
274  particle%z = z
275  particle%ttrack = t
276  particle%istatus = 1
277  call this%save(particle, reason=5)
278  end do
279  end if
280 
281  if (texit .gt. tmax) then
282  ! -- The computed exit time is greater than the maximum time, so set
283  ! -- final time for particle trajectory equal to maximum time.
284  t = tmax
285  dt = t - t0
286  exitface = 0
287  particle%istatus = 1
288  particle%advancing = .false.
289  reason = 2 ! timestep end
290  else
291  ! -- The computed exit time is less than or equal to the maximum time,
292  ! -- so set final time for particle trajectory equal to exit time.
293  t = texit
294  dt = dtexit
295  reason = 1 ! cell transition
296  end if
297 
298  ! -- Calculate final particle location
299  ! -- kluge note: need to evaluate both alpha and beta here only
300  ! -- for exitFace=0, otherwise just one or the other
301  call step_analytical(dt, alp, bet)
302  if (exitface .eq. 1) then
303  bet = 0d0
304  else if (exitface .eq. 2) then
305  alp = 1d0 - bet
306  else if (exitface .eq. 3) then
307  alp = 0d0
308  end if
309  loc = (/alp, bet/)
310  if (lbary) loc = skew(loc, (/sxx, sxy, syy/), invert=.true.)
311  rot = reshape((/rxx, rxy, ryx, ryy/), shape(rot))
312  res = matmul(rot, loc) ! rotate vector
313  x = res(1) + x0
314  y = res(2) + y0
315  if (exitface .eq. 4) then
316  z = zbot
317  else if (exitface .eq. 5) then
318  z = ztop
319  else
320  if (izstatus .eq. 2) then ! kluge note: make this into a function
321  ! -- vz uniformly zero
322  z = zi
323  else if (izstatus .eq. 1) then
324  ! -- vz uniform, nonzero
325  z = zi + vzi * dt
326  else
327  ! -- vz nonuniform
328  z = zbot + (vzi * dexp(az * dt) - vzbot) / az
329  end if
330  end if
331 
332  ! -- Set final particle location in local (unscaled) subcell coordinates,
333  ! -- final time for particle trajectory, and exit face
334  particle%x = x
335  particle%y = y
336  particle%z = z
337  particle%ttrack = t
338  particle%iboundary(3) = exitface
339 
340  ! -- Save particle track record
341  if (reason > -1) &
342  call this%save(particle, reason=reason) ! reason=2: timestep
343  end subroutine track_subcell
344 
345  !> @brief Do calculations related to analytical z solution
346  !!
347  !! This subroutine consists partly or entirely of code written by
348  !! David W. Pollock of the USGS for MODPATH 7. The authors of the present
349  !! code are responsible for its appropriate application in this context
350  !! and for any modifications or errors.
351  !<
352  subroutine calculate_dt(v1, v2, dx, xL, v, dvdx, &
353  dt, status, itopbotexit)
354  real(DP) :: v1
355  real(DP) :: v2
356  real(DP) :: dx
357  real(DP) :: xL
358  real(DP) :: v
359  real(DP) :: dvdx
360  real(DP) :: dt
361  real(DP) :: v2a
362  real(DP) :: v1a
363  real(DP) :: dv
364  real(DP) :: dva
365  real(DP) :: vv
366  real(DP) :: vvv
367  real(DP) :: zro
368  real(DP) :: zrom
369  real(DP) :: x
370  real(DP) :: tol
371  real(DP) :: vr1
372  real(DP) :: vr2
373  real(DP) :: vr
374  real(DP) :: v1v2
375  integer(I4B) :: status
376  integer(I4B) :: itopbotexit
377  logical(LGP) :: noOutflow
378 
379  ! Initialize variables
380  status = -1
381  dt = 1.0d+20
382  v2a = v2
383  if (v2a .lt. 0d0) v2a = -v2a
384  v1a = v1
385  if (v1a .lt. 0d0) v1a = -v1a
386  dv = v2 - v1
387  dva = dv
388  if (dva .lt. 0d0) dva = -dva
389 
390  ! Check for a uniform zero velocity in this direction.
391  ! If so, set status = 2 and return (dt = 1.0d+20).
392  tol = 1.0d-15
393  if ((v2a .lt. tol) .and. (v1a .lt. tol)) then
394  v = 0d0
395  dvdx = 0d0
396  status = 2
397  itopbotexit = 0
398  return
399  end if
400 
401  ! Check for uniform non-zero velocity in this direction.
402  ! If so, set compute dt using the constant velocity,
403  ! set status = 1 and return.
404  vv = v1a
405  if (v2a .gt. vv) vv = v2a
406  vvv = dva / vv
407  if (vvv .lt. 1.0d-4) then
408  zro = tol
409  zrom = -zro
410  v = v1
411  x = xl * dx
412  if (v1 .gt. zro) then
413  dt = (dx - x) / v1
414  itopbotexit = -2
415  end if
416  if (v1 .lt. zrom) then
417  dt = -x / v1
418  itopbotexit = -1
419  end if
420  dvdx = 0d0
421  status = 1
422  return
423  end if
424 
425  ! Velocity has a linear variation.
426  ! Compute velocity corresponding to particle position
427  dvdx = dv / dx
428  v = (1.0d0 - xl) * v1 + xl * v2
429 
430  ! If flow is into the cell from both sides there is no outflow.
431  ! In that case, set status = 3 and return
432  nooutflow = .true.
433  if (v1 .lt. 0d0) nooutflow = .false.
434  if (v2 .gt. 0d0) nooutflow = .false.
435  if (nooutflow) then
436  status = 3
437  itopbotexit = 0
438  return
439  end if
440 
441  ! If there is a divide in the cell for this flow direction, check to see if the
442  ! particle is located exactly on the divide. If it is, move it very slightly to
443  ! get it off the divide. This avoids possible numerical problems related to
444  ! stagnation points.
445  if ((v1 .le. 0d0) .and. (v2 .ge. 0d0)) then
446  if (abs(v) .le. 0d0) then
447  v = 1.0d-20
448  if (v2 .le. 0d0) v = -v
449  end if
450  end if
451 
452  ! If there is a flow divide, find out what side of the divide the particle
453  ! is on and set the value of vr appropriately to reflect that location.
454  vr1 = v1 / v
455  vr2 = v2 / v
456  vr = vr1
457  itopbotexit = -1
458  if (vr .le. 0d0) then
459  vr = vr2
460  itopbotexit = -2
461  end if
462 
463  ! Check if velocity is in the same direction throughout cell (i.e. no flow divide).
464  ! Check if product v1*v2 > 0 then the velocity is in the same direction throughout
465  ! the cell (i.e. no flow divide). If so, set vr to reflect appropriate direction.
466  v1v2 = v1 * v2
467  if (v1v2 .gt. 0d0) then
468  if (v .gt. 0d0) then
469  vr = vr2
470  itopbotexit = -2
471  end if
472  if (v .lt. 0d0) then
473  vr = vr1
474  itopbotexit = -1
475  end if
476  end if
477 
478  ! Compute travel time to exit face. Return with status = 0
479  dt = log(vr) / dvdx
480  status = 0
481  end subroutine calculate_dt
482 
subroutine pstop(status, message)
Stop the program, optionally specifying an error status code.
Definition: ErrorUtil.f90:24
pure real(dp) function, dimension(2), public skew(v, s, invert)
Skew a 2D vector along the x-axis.
Definition: GeomUtil.f90:147
This module defines variable data types.
Definition: kind.f90:8
Particle tracking strategies.
Definition: Method.f90:2
subroutine track_subcell(this, subcell, particle, tmax)
Track a particle across a triangular subcell using the ternary method.
subroutine apply_mst(this, particle, tmax)
Apply the ternary subcell method.
subroutine, public create_method_subcell_ternary(method)
Create a new ternary subcell-method object.
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 traverse_triangle(isolv, tol, step, texit, alpexit, betexit, itrifaceenter, itrifaceexit, rxx, rxy, ryx, ryy, alp0, bet0, alp1, bet1, alp2, bet2, alpi, beti, vziodz, az, bary)
Traverse triangular cell.
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, bary)
Set coordinates to "canonical" configuration.
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:10
Base type for particle tracking methods.
Definition: Method.f90:29
Ternary triangular subcell tracking method.
A particle tracked by the PRT model.
Definition: Particle.f90:31
A subcell of a cell.
Definition: Subcell.f90:9
Manages particle track (i.e. pathline) files.
Definition: TrackData.f90:38