36 #include <boost/core/ignore_unused.hpp>
45 std::string LEEUpwindSolver::solverName =
47 "LEEUpwind", LEEUpwindSolver::create,
"Upwind Solver for LEE");
52 LEEUpwindSolver::LEEUpwindSolver(
95 boost::ignore_unused(rho0L, v0L, w0L, rho0R, v0R, w0R);
113 rhoF = rhoF + u0L * (c0sqL * rhoL - pL) / c0sqL;
114 rhovF = rhovF + rhovL * u0L;
115 rhowF = rhowF + rhowL * u0L;
119 rhoF = rhoF + u0R * (c0sqR * rhoR - pR) / c0sqR;
120 rhovF = rhovF + rhovR * u0R;
121 rhowF = rhowF + rhowR * u0R;
127 pF = pF + 0.5 * (c0L - u0L) * (rhouL * c0L - pL);
128 rhoF = rhoF + 0.5 * (c0L - u0L) * (rhouL * c0sqL - c0L * pL) /
129 pow(c0sqL, 3.0 / 2.0);
130 rhouF = rhouF + 0.5 * (c0L - u0L) * (-rhouL * c0L + pL) / c0L;
134 pF = pF + 0.5 * (c0R - u0R) * (rhouR * c0R - pR);
135 rhoF = rhoF + 0.5 * (c0R - u0R) * (rhouR * c0sqR - c0R * pR) /
136 pow(c0sqR, 3.0 / 2.0);
137 rhouF = rhouF + 0.5 * (c0R - u0R) * (-rhouR * c0R + pR) / c0R;
143 pF = pF + 0.5 * (c0L + u0L) * (rhouL * c0L + pL);
144 rhoF = rhoF + 0.5 * (c0L + u0L) * (rhouL * c0sqL + c0L * pL) /
145 pow(c0sqL, 3.0 / 2.0);
146 rhouF = rhouF + 0.5 * (c0L + u0L) * (rhouL * c0L + pL) / c0L;
150 pF = pF + 0.5 * (c0R + u0R) * (rhouR * c0R + pR);
151 rhoF = rhoF + 0.5 * (c0R + u0R) * (rhouR * c0sqR + c0R * pR) /
152 pow(c0sqR, 3.0 / 2.0);
153 rhouF = rhouF + 0.5 * (c0R + u0R) * (rhouR * c0R + pR) / c0R;
virtual void v_PointSolve(NekDouble pL, NekDouble rhoL, NekDouble rhouL, NekDouble rhovL, NekDouble rhowL, NekDouble pR, NekDouble rhoR, NekDouble rhouR, NekDouble rhovR, NekDouble rhowR, 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 &rhouF, NekDouble &rhovF, NekDouble &rhowF) override
Lax-Friedrichs Riemann solver.
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)