34 itrifaceenter, itrifaceexit, &
35 alp1, bet1, alp2, bet2, alpi, beti)
37 integer(I4B),
intent(in) :: isolv
38 real(dp),
intent(in) :: tol
39 real(dp),
intent(out) :: texit
42 integer(I4B) :: itrifaceenter
43 integer(I4B) :: itrifaceexit
70 texit0, alpexit0, betexit0)
73 texit1, alpexit1, betexit1)
76 texit2, alpexit2, betexit2)
77 texit = min(texit0, texit1, texit2)
81 if (texit .eq. texit0)
then
85 else if (texit .eq. texit1)
then
89 else if (texit .eq. texit2)
then
94 if (texit .eq. huge(1d0)) itrifaceexit = 0
100 v0x, v0y, v1x, v1y, v2x, v2y, &
102 rxx, rxy, ryx, ryy, &
104 alp0, bet0, alp1, bet1, alp2, bet2, alpi, beti)
124 real(dp),
intent(inout) :: sxx, sxy, syy
135 real(dp) :: oobaselen
144 real(dp) :: rot(2, 2), res(2)
151 baselen = dsqrt(x1diff * x1diff + y1diff * y1diff)
152 oobaselen =
done / baselen
153 cosomega = x1diff * oobaselen
154 sinomega = y1diff * oobaselen
164 rot = reshape((/rxx, ryx, rxy, ryy/), shape(rot))
165 res = matmul(rot, (/x2diff, y2diff/))
169 cv0 = matmul(rot, (/v0x, v0y/))
170 cv1 = matmul(rot, (/v1x, v1y/))
171 cv2 = matmul(rot, (/v2x, v2y/))
175 res = matmul(rot, (/xidiff, yidiff/))
181 sxy = -alp2 * sxx * syy
189 res =
skew(res, (/sxx, sxy, syy/))
197 alp1, bet1, alp2, bet2, &
209 real(dp) :: v1alpdiff
210 real(dp) :: v2alpdiff
211 real(dp) :: v2betdiff
215 v1alpdiff =
cv1(1) -
cv0(1)
216 v2alpdiff =
cv2(1) -
cv0(1)
217 v2betdiff =
cv2(2) -
cv0(2)
241 if (dabs(
wbb) .gt. zerotol)
then
243 acoef =
cv0(1) - wratv
244 bcoef = wratv +
wab * beti
251 if (dabs(
waa) .gt. zerotol)
then
253 if (dabs(
wbb -
waa) .gt. zerotol)
then
255 bfact = bcoef / (
wbb -
waa)
257 ca2 = alpi + afact - bfact
279 if (dabs(
waa) .gt. zerotol)
then
282 vfact = (
wab * oowaa) *
cv0(2)
283 ca1 = -oowaa * (
cv0(1) +
wab * beti + vfact)
301 real(dp),
intent(in) :: t
305 if (
icase .eq. 1)
then
308 else if (
icase .eq. -1)
then
311 else if (
icase .eq. 2)
then
314 else if (
icase .eq. 3)
then
317 else if (
icase .eq. 4)
then
327 texit, alpexit, betexit)
329 integer(I4B) :: isolv
330 integer(I4B) :: itriface
331 integer(I4B) :: itrifaceenter
348 real(dp) :: v0alpstar
359 integer(I4B) :: ibettrend
365 if (itriface .eq. 0)
then
367 if (itrifaceenter .eq. 0)
then
380 vbeti =
cv0(2) +
wbb * beti
381 if (vbeti .ge.
dzero)
then
387 if ((alpt .ge.
dzero) .and. (alpt .le.
done))
then
402 if (itriface .eq. 1)
then
411 if ((v1n .ge.
dzero) .and. (v2n .ge.
dzero))
then
420 if (ibettrend .eq. 0)
then
423 if ((beti .gt. betouthi) .or. (beti .lt. betoutlo))
then
430 v0alpstar =
cv0(1) +
wab * beti
431 valpi = v0alpstar +
waa * alpi
432 if ((itriface .eq. 1) .and. (valpi .le.
dzero))
then
435 else if ((itriface .eq. 2) .and. (valpi .ge.
dzero))
then
440 if (itriface .eq. 1)
then
441 alpexit =
done - beti
446 if (dabs(
waa) > zerotol)
then
447 alplim = -v0alpstar /
waa
448 texit = dlog(alpexit - alplim / (alpi - alplim)) /
waa
450 texit = (alpexit - alpi) / v0alpstar
457 bethi = min(betouthi, betsolhi)
458 betlo = max(betoutlo, betsollo)
459 if (betlo .gt. bethi)
then
466 if (itriface .eq. 1)
then
467 fax =
done - betlo - alplo
468 fbx =
done - bethi - alphi
473 if (fax * fbx .gt.
dzero)
then
477 if (isolv .eq. 1)
then
480 call soln_brent(itriface, betlo, bethi, tol, texit, &
482 else if (isolv .eq. 2)
then
485 call soln_chand(itriface, betlo, bethi, tol, texit, &
488 call pstop(1,
"Invalid isolv, expected one of: 1, 2")
498 if (texit .ne. huge(
done) .and. texit .lt.
dzero)
then
499 call pstop(1,
"texit is negative (unexpected)")
507 real(dp),
intent(in) :: bet
515 fb =
done - alpt - bet
521 real(dp),
intent(in) :: bet
535 real(dp),
intent(in) :: bet
552 if (
icase .eq. 1)
then
553 term = max(term, zerotol)
556 if (waat > waatlim)
then
557 alp = sign(ewaatlim,
ca2)
559 alp =
ca1 +
ca2 * dexp(waat) +
ca3 * term
561 else if (
icase .eq. -1)
then
562 term = max(term, zerotol)
566 if (waat > waatlim)
then
567 alp = sign(ewaatlim, coef)
569 alp =
ca1 + coef * dexp(waat)
571 else if (
icase .eq. 2)
then
572 term = max(term, zerotol)
575 else if (
icase .eq. 3)
then
578 if (waat > waatlim)
then
579 alp = sign(ewaatlim,
ca3)
583 else if (
icase .eq. 4)
then
601 if (vn1 .lt.
dzero)
then
604 if (vn2 .le.
dzero)
then
609 betouthi = -vn1 / vndiff
614 if (vn1 .le.
dzero)
then
619 betoutlo = -vn1 / vndiff
628 real(dp),
intent(in) :: beti
631 integer(I4B),
intent(inout) :: ibettrend
637 betsolhi = huge(
done)
642 betsollo = -huge(
done)
651 if (beti .gt. betlim)
then
652 betsolhi = huge(
done)
655 else if (beti .lt. betlim)
then
657 betsollo = -huge(
done)
666 if (beti .gt. betlim)
then
670 else if (beti .lt. betlim)
then
685 texit, alpexit, betexit)
687 integer(I4B),
intent(in) :: itriface
690 real(dp),
intent(in) :: tol
697 procedure(
f1d),
pointer :: f
702 if (itriface .eq. 1)
then
704 betexit =
zero_br(blo, bhi, f, tol)
707 betexit =
zero_br(blo, bhi, f, tol)
715 texit, alpexit, betexit)
717 integer(I4B),
intent(in) :: itriface
720 real(dp),
intent(in) :: tol
727 procedure(
f1d),
pointer :: f
732 if (itriface .eq. 1)
then
734 betexit =
zero_ch(blo, bhi, f, tol)
737 betexit =
zero_ch(blo, bhi, f, tol)
This module contains simulation constants.
real(dp), parameter dzero
real constant zero
real(dp), parameter dprec
real constant machine precision
real(dp), parameter done
real constant 1
subroutine pstop(status, message)
Stop the program, optionally specifying an error status code.
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.
real(dp) function, public zero_ch(x0, x1, f, epsa)
Compute zeros on an interval using Chadrupatla's method.
real(dp) function, public zero_br(ax, bx, f, tol)
Compute a zero of the function f(x) in the interval (x0, x1).
real(dp) function fbary1(bet)
Brent's method applied to canonical face 1 (gamma = 0)
real(dp), dimension(2) cv2
"Canonical" velocity components at corners of triangular subcell
integer(i4b) icase
Case index for analytical solution.
real(dp) wbb
Elements of the "velocity matrix," W.
subroutine, public get_w(alp1, bet1, alp2, bet2, waa, wab, wba, wbb)
Compute elements of W matrix.
real(dp), dimension(2) cv1
subroutine, public get_bet_outflow_bary(vn1, vn2, betoutlo, betouthi)
Find outflow interval.
subroutine, public get_t_alpt(bet, t, alp)
Given beta evaluate t and alpha depending on case.
subroutine, public get_bet_soln_limits(beti, betsollo, betsolhi, ibettrend)
Find trend of and limits on beta from beta{t} solution.
subroutine, public solve_coefs(alpi, beti)
Compute analytical solution coefficients depending on case.
subroutine, public find_exit_bary(isolv, itriface, itrifaceenter, alpi, beti, tol, texit, alpexit, betexit)
Find the exit time and location in barycentric coordinates.
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 soln_brent(itriface, betlo, bethi, tol, texit, alpexit, betexit)
Use Brent's method with initial bounds on beta of betlo and bethi.
real(dp) cb2
Analytical solution coefficients.
real(dp), dimension(2) cv0
real(dp) function fbary2(bet)
Brent's method applied to canonical face 2 (alpha = 0)
subroutine, public soln_chand(itriface, betlo, bethi, tol, texit, alpexit, betexit)
Use Chandrupatla's method with initial bounds on beta of betlo and bethi.
subroutine, public step_analytical(t, alp, bet)
Step (evaluate) analytically depending on case.