36#include <boost/core/ignore_unused.hpp>
48 "Upwind solver for the APE equation");
96 boost::ignore_unused(rhoL, rhoR, rhoF);
113 pF = pF + 0.5 * (c0L - u0L) *
114 (-rho0L * v0L * vL * (c0L * u0L + c0sqL) -
115 rho0L * w0L * wL * (c0L * u0L + c0sqL) +
116 (c0sqL - pow(u0L, 2)) * (rho0L * c0L * uL - pL)) /
117 (c0sqL - pow(u0L, 2));
121 (rho0L * c0L * (c0L * u0L + c0sqL) * (v0L * vL + w0L * wL) -
122 rho0L * c0sqL * uL * (c0sqL - pow(u0L, 2)) +
123 c0L * pL * (c0sqL - pow(u0L, 2))) /
124 (rho0L * c0sqL * (c0sqL - pow(u0L, 2)));
128 pF = pF + 0.5 * (c0R - u0R) *
129 (-rho0R * v0R * vR * (c0R * u0R + c0sqR) -
130 rho0R * w0R * wR * (c0R * u0R + c0sqR) +
131 (c0sqR - pow(u0R, 2)) * (rho0R * c0R * uR - pR)) /
132 (c0sqR - pow(u0R, 2));
136 (rho0R * c0R * (c0R * u0R + c0sqR) * (v0R * vR + w0R * wR) -
137 rho0R * c0sqR * uR * (c0sqR - pow(u0R, 2)) +
138 c0R * pR * (c0sqR - pow(u0R, 2))) /
139 (rho0R * c0sqR * (c0sqR - pow(u0R, 2)));
145 pF = pF + 0.5 * (c0L + u0L) *
146 (-rho0L * v0L * vL * (c0L * u0L - c0sqL) -
147 rho0L * w0L * wL * (c0L * u0L - c0sqL) +
148 (c0sqL - pow(u0L, 2)) * (rho0L * c0L * uL + pL)) /
149 (c0sqL - pow(u0L, 2));
153 (-rho0L * c0L * (c0L * u0L - c0sqL) * (v0L * vL + w0L * wL) +
154 rho0L * c0sqL * uL * (c0sqL - pow(u0L, 2)) +
155 c0L * pL * (c0sqL - pow(u0L, 2))) /
156 (rho0L * c0sqL * (c0sqL - pow(u0L, 2)));
160 pF = pF + 0.5 * (c0R + u0R) *
161 (-rho0R * v0R * vR * (c0R * u0R - c0sqR) -
162 rho0R * w0R * wR * (c0R * u0R - c0sqR) +
163 (c0sqR - pow(u0R, 2)) * (rho0R * c0R * uR + pR)) /
164 (c0sqR - pow(u0R, 2));
168 (-rho0R * c0R * (c0R * u0R - c0sqR) * (v0R * vR + w0R * wR) +
169 rho0R * c0sqR * uR * (c0sqR - pow(u0R, 2)) +
170 c0R * pR * (c0sqR - pow(u0R, 2))) /
171 (rho0R * c0sqR * (c0sqR - pow(u0R, 2)));
APEUpwindSolver(const LibUtilities::SessionReaderSharedPtr &pSession)
static std::string solverName
virtual void v_PointSolve(NekDouble pL, NekDouble rhoL, NekDouble uL, NekDouble vL, NekDouble wL, NekDouble pR, NekDouble rhoR, NekDouble uR, NekDouble vR, NekDouble wR, NekDouble c0sqL, NekDouble rho0L, NekDouble u0L, NekDouble v0L, NekDouble w0L, NekDouble c0sqR, NekDouble rho0R, NekDouble u0R, NekDouble v0R, NekDouble w0R, NekDouble &pF, NekDouble &rhoF, NekDouble &uF, NekDouble &vF, NekDouble &wF) override
Upwind Riemann solver.
static RiemannSolverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
RiemannSolverFactory & GetRiemannSolverFactory()
The above copyright notice and this permission notice shall be included.
scalarT< T > sqrt(scalarT< T > in)