Nektar++
LEELaxFriedrichsSolver.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: LEELaxFriedrichsSolver.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2015 Kilian Lackhove
10// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
11// Department of Aeronautics, Imperial College London (UK), and Scientific
12// Computing and Imaging Institute, University of Utah (USA).
13//
14// Permission is hereby granted, free of charge, to any person obtaining a
15// copy of this software and associated documentation files (the "Software"),
16// to deal in the Software without restriction, including without limitation
17// the rights to use, copy, modify, merge, publish, distribute, sublicense,
18// and/or sell copies of the Software, and to permit persons to whom the
19// Software is furnished to do so, subject to the following conditions:
20//
21// The above copyright notice and this permission notice shall be included
22// in all copies or substantial portions of the Software.
23//
24// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30// DEALINGS IN THE SOFTWARE.
31//
32// Description: Lax-Friedrichs solver for the LEE equations.
33//
34///////////////////////////////////////////////////////////////////////////////
35
37
38using namespace std;
39
40namespace Nektar
41{
42
45 "LEELaxFriedrichs", LEELaxFriedrichsSolver::create,
46 "Lax-Friedrichs Solver for LEE");
47
48/**
49 *
50 */
53 : LEESolver(pSession)
54{
55}
56
57/**
58 * @brief Lax-Friedrichs Riemann solver
59 *
60 * @param pL Perturbation pressure left state
61 * @param rhoL Perturbation density left state
62 * @param pR Perturbation pressure right state
63 * @param rhoR Perturbation density right state
64 * @param rhouL x perturbation velocity component left state
65 * @param rhouR x perturbation velocity component right state
66 * @param rhovL y perturbation velocity component left state
67 * @param rhovR y perturbation velocity component right state
68 * @param rhowL z perturbation velocity component left state
69 * @param rhowR z perturbation velocity component right state
70 * @param c0sqL Base pressure left state
71 * @param c0sqR Base pressure right state
72 * @param rho0L Base density left state
73 * @param rho0R Base density right state
74 * @param u0L Base x velocity component left state
75 * @param u0R Base x velocity component right state
76 * @param v0L Base y velocity component left state
77 * @param v0R Base y velocity component right state
78 * @param w0L Base z velocity component left state
79 * @param w0R Base z velocity component right state
80 * @param pF Computed Riemann flux for perturbation pressure
81 * @param rhoF Computed Riemann flux for perturbation density
82 * @param rhouF Computed Riemann flux for x perturbation velocity component
83 * @param rhovF Computed Riemann flux for y perturbation velocity component
84 * @param rhowF Computed Riemann flux for z perturbation velocity component
85 */
87 NekDouble pL, NekDouble rhoL, NekDouble rhouL, NekDouble rhovL,
88 NekDouble rhowL, NekDouble pR, NekDouble rhoR, NekDouble rhouR,
89 NekDouble rhovR, NekDouble rhowR, NekDouble c0sqL,
90 [[maybe_unused]] NekDouble rho0L, NekDouble u0L,
91 [[maybe_unused]] NekDouble v0L, [[maybe_unused]] NekDouble w0L,
92 NekDouble c0sqR, [[maybe_unused]] NekDouble rho0R, NekDouble u0R,
93 [[maybe_unused]] NekDouble v0R, [[maybe_unused]] NekDouble w0R,
94 NekDouble &pF, NekDouble &rhoF, NekDouble &rhouF, NekDouble &rhovF,
95 NekDouble &rhowF)
96{
97 // Speed of sound
98 NekDouble cL = sqrt(c0sqL);
99 NekDouble cR = sqrt(c0sqR);
100
101 // max absolute eigenvalue of the jacobian of F_n1
102 NekDouble a_1_max = 0;
103 a_1_max = std::max(a_1_max, std::abs(u0L - cL));
104 a_1_max = std::max(a_1_max, std::abs(u0R - cR));
105 a_1_max = std::max(a_1_max, std::abs(u0L + cL));
106 a_1_max = std::max(a_1_max, std::abs(u0R + cR));
107
108 NekDouble pFL = rhouL * c0sqL + u0L * pL;
109 NekDouble rhoFL = rhouL + rhoL * u0L;
110 NekDouble rhouFL = pL + rhouL * u0L;
111 NekDouble rhovFL = 0;
112 NekDouble rhowFL = 0;
113
114 NekDouble pFR = rhouR * c0sqR + u0R * pR;
115 NekDouble rhoFR = rhouR + rhoR * u0R;
116 NekDouble rhouFR = pR + rhouR * u0R;
117 NekDouble rhovFR = 0;
118 NekDouble rhowFR = 0;
119
120 // assemble the face-normal fluxes
121 pF = 0.5 * (pFL + pFR - a_1_max * (pR - pL));
122 rhoF = 0.5 * (rhoFL + rhoFR - a_1_max * (rhoR - rhoL));
123 rhouF = 0.5 * (rhouFL + rhouFR - a_1_max * (rhouR - rhouL));
124 rhovF = 0.5 * (rhovFL + rhovFR - a_1_max * (rhovR - rhovL));
125 rhowF = 0.5 * (rhowFL + rhowFR - a_1_max * (rhowR - rhowL));
126}
127
128} // namespace Nektar
static RiemannSolverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
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.
LEELaxFriedrichsSolver(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()
double NekDouble
STL namespace.
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:298
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294