23 integer(I4B),
public,
pointer :: zeromethod
41 method%subcell => subcell
42 method%name => method%subcell%type
43 method%delegates = .false.
49 deallocate (this%name)
56 real(DP),
intent(in) :: tmax
58 select type (subcell => this%subcell)
60 call this%track_subcell(subcell, particle, tmax)
70 real(DP),
intent(in) :: tmax
72 integer(I4B) :: exitFace
122 integer(I4B) :: izstatus
123 integer(I4B) :: itopbotexit
124 integer(I4B) :: ntmax
125 integer(I4B) :: isolv
126 integer(I4B) :: itrifaceenter
127 integer(I4B) :: itrifaceexit
132 integer(I4B) :: reason
136 if (particle%iexmeth == 0)
then
139 isolv = particle%iexmeth
165 vzbot = subcell%vzbot
166 vztop = subcell%vztop
179 v0x, v0y, v1x, v1y, v2x, v2y, &
181 rxx, rxy, ryx, ryy, &
183 alp0, bet0, alp1, bet1, alp2, bet2, alpi, beti)
191 zirel = (zi - zbot) / dz
193 az, dtexitz, izstatus, &
198 itrifaceenter = particle%iboundary(3) - 1
199 if (itrifaceenter == -1) itrifaceenter = 999
201 dtexitxy, alpexit, betexit, &
202 itrifaceenter, itrifaceexit, &
203 alp1, bet1, alp2, bet2, alpi, beti)
207 if (itopbotexit == 0 .and. itrifaceexit == 0)
then
209 particle%advancing = .false.
210 call this%save(particle, reason=3)
218 if (itrifaceexit /= 0)
then
220 exitface = itrifaceexit
222 else if (dtexitz < dtexitxy .and. dtexitz >= 0.0_dp)
then
224 if (itopbotexit == -1)
then
232 particle%advancing = .false.
233 call this%save(particle, reason=3)
236 texit = particle%ttrack + dtexit
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
251 izstatus, x0, y0, az, vzi, vzbot, &
252 ztop, zbot, zi, x, y, z)
258 call this%save(particle, reason=5)
265 if (texit .gt. tmax)
then
272 particle%advancing = .false.
282 izstatus, x0, y0, az, vzi, vzbot, &
283 ztop, zbot, zi, x, y, z, exitface)
288 particle%iboundary(3) = exitface
290 call this%save(particle, reason=reason)
301 dt, status, itopbotexit)
323 integer(I4B) :: status
324 integer(I4B) :: itopbotexit
325 logical(LGP) :: noOutflow
331 if (v2a .lt.
dzero) v2a = -v2a
333 if (v1a .lt.
dzero) v1a = -v1a
336 if (dva .lt.
dzero) dva = -dva
341 if ((v2a .lt. tol) .and. (v1a .lt. tol))
then
353 if (v2a .gt. vv) vv = v2a
355 if (vvv .lt. 1.0d-4)
then
360 if (v1 .gt. zro)
then
364 if (v1 .lt. zrom)
then
376 v = (
done - xl) * v1 + xl * v2
381 if (v1 .lt.
dzero) nooutflow = .false.
382 if (v2 .gt.
dzero) nooutflow = .false.
393 if ((v1 .le.
dzero) .and. (v2 .ge.
dzero))
then
394 if (abs(v) .le.
dzero)
then
396 if (v2 .le.
dzero) v = -v
406 if (vr .le.
dzero)
then
415 if (v1v2 .gt.
dzero)
then
416 if (v .gt.
dzero)
then
420 if (v .lt.
dzero)
then
427 dt = log(abs(vr)) / dvdx
433 izstatus, x0, y0, az, vzi, vzbot, &
434 ztop, zbot, zi, x, y, z, exitFace)
444 integer(I4B) :: izstatus
456 integer(I4B),
optional :: exitFace
458 integer(I4B) :: lexitface
459 real(DP) :: rot(2, 2), res(2), loc(2)
464 if (
present(exitface))
then
475 if (lexitface .eq. 1)
then
477 else if (lexitface .eq. 2)
then
479 else if (lexitface .eq. 3)
then
484 if (lexitface .eq. 4)
then
486 else if (lexitface .eq. 5)
then
490 if (izstatus .eq. 2)
then
493 else if (izstatus .eq. 1)
then
498 z = zbot + (vzi * dexp(az * dt) - vzbot) / az
504 loc =
skew(loc, (/sxx, sxy, syy/), invert=.true.)
505 rot = reshape((/rxx, rxy, ryx, ryy/), shape(rot))
506 res = matmul(rot, loc)
This module contains simulation constants.
real(dp), parameter dsame
real constant for values that are considered the same based on machine precision
real(dp), parameter dep3
real constant 1000
real(dp), parameter donethird
real constant 1/3
real(dp), parameter dhalf
real constant 1/2
real(dp), parameter dzero
real constant zero
real(dp), parameter dtwo
real constant 2
real(dp), parameter done
real constant 1
subroutine pstop(status, message)
Stop the program, optionally specifying an error status code.
subroutine, public clamp_bary(alpha, beta, gamma, pad)
Clamp barycentric coordinates to the interior of a triangle, with optional padding some minimum dista...
pure real(dp) function, dimension(2), public skew(v, s, invert)
Skew a 2D vector along the x-axis.
This module defines variable data types.
pure logical function, public is_close(a, b, rtol, atol, symmetric)
Check if a real value is approximately equal to another.
Particle tracking strategies.
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.
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...
Base type for particle tracking methods.
Ternary triangular subcell tracking method.
Particle tracked by the PRT model.