File riemann_2shock_solvers.H
-
namespace TwoShock
Functions
- AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void wsqge (const amrex::Real p, const amrex::Real v, const amrex::Real gam, const amrex::Real gdot, amrex::Real &gstar, const amrex::Real gmin, const amrex::Real gmax, const amrex::Real csq, const amrex::Real pstar, amrex::Real &wsq)
- AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void pstar_bisection (amrex::Real &pstar_lo, amrex::Real &pstar_hi, const amrex::Real ul, const amrex::Real pl, const amrex::Real taul, const amrex::Real gamel, const amrex::Real clsql, const amrex::Real ur, const amrex::Real pr, const amrex::Real taur, const amrex::Real gamer, const amrex::Real clsqr, const amrex::Real gdot, const amrex::Real gmin, const amrex::Real gmax, const int lriemann_shock_maxiter, const amrex::Real lriemann_pstar_tol, amrex::Real &pstar, amrex::Real &gamstar, bool &converged, amrex::GpuArray< amrex::Real, riemann_constants::PSTAR_BISECT_FACTOR *riemann_constants::HISTORY_SIZE > &pstar_hist_extra)
- AMREX_GPU_HOST_DEVICE AMREX_INLINE void riemanncg (const RiemannState &ql, const RiemannState &qr, const RiemannAux &raux, RiemannState &qint)
The Colella-Glaz Riemann solver for pure hydrodynamics. This is a two shock approximate state Riemann solver.
- Parameters:
bx – the box to operate over
ql – the left interface state
qr – the right interface state
qaux_arr – the auxiliary state
qint – the full Godunov state on the interface
idir – coordinate direction for the solve (0 = x, 1 = y, 2 = z)
- AMREX_GPU_HOST_DEVICE AMREX_INLINE void riemannus (const RiemannState &ql, const RiemannState &qr, const RiemannAux &raux, RiemannState &qint)
The Colella-Glaz-Ferguson Riemann solver for hydrodynamics and radiation hydrodynamics. This is a two shock approximate state Riemann solver.
- Parameters:
bx – the box to operate over
ql – the left interface state
qr – the right interface state
qaux_arr – the auxiliary state
qint – the full Godunov state on the interface
lambda_int – the radiation flux limiter on the interface