Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
APESolver/RiemannSolvers/LaxFriedrichsSolver.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: LaxFriedrichsSolver.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 // License for the specific language governing rights and limitations under
15 // Permission is hereby granted, free of charge, to any person obtaining a
16 // copy of this software and associated documentation files (the "Software"),
17 // to deal in the Software without restriction, including without limitation
18 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
19 // and/or sell copies of the Software, and to permit persons to whom the
20 // Software is furnished to do so, subject to the following conditions:
21 //
22 // The above copyright notice and this permission notice shall be included
23 // in all copies or substantial portions of the Software.
24 //
25 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
26 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
27 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
28 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
29 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
30 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
31 // DEALINGS IN THE SOFTWARE.
32 //
33 // Description: Lax-Friedrichs solver for the APE equations.
34 //
35 ///////////////////////////////////////////////////////////////////////////////
36 
38 
39 using namespace std;
40 
41 namespace Nektar
42 {
43 
44 std::string LaxFriedrichsSolver::solverName =
46  RegisterCreatorFunction("LaxFriedrichs", LaxFriedrichsSolver::create,
47  "Lax-Friedrichs Solver");
48 
49 
50 /**
51 *
52 */
53 LaxFriedrichsSolver::LaxFriedrichsSolver() :
54  APESolver()
55 {
56 
57 }
58 
59 /**
60  * @brief Lax-Friedrichs Riemann solver
61  *
62  * @param pL Perturbation pressure left state
63  * @param rhoL Perturbation density left state
64  * @param pR Perturbation pressure right state
65  * @param rhoR Perturbation density right state
66  * @param uL x perturbation velocity component left state
67  * @param uR x perturbation velocity component right state
68  * @param vL y perturbation velocity component left state
69  * @param vR y perturbation velocity component right state
70  * @param wL z perturbation velocity component left state
71  * @param wR z perturbation velocity component right state
72  * @param p0 Base pressure
73  * @param rho0 Base density
74  * @param u0 Base x velocity component
75  * @param v0 Base y velocity component
76  * @param w0 Base z velocity component
77  * @param pF Computed Riemann flux for perturbation pressure
78  * @param rhoF Computed Riemann flux for perturbation density
79  * @param uF Computed Riemann flux for x perturbation velocity component
80  * @param vF Computed Riemann flux for y perturbation velocity component
81  * @param wF Computed Riemann flux for z perturbation velocity component
82  */
84  NekDouble pL, NekDouble rhoL, NekDouble uL, NekDouble vL, NekDouble wL,
85  NekDouble pR, NekDouble rhoR, NekDouble uR, NekDouble vR, NekDouble wR,
86  NekDouble p0, NekDouble rho0, NekDouble u0, NekDouble v0, NekDouble w0,
87  NekDouble &pF, NekDouble &rhoF, NekDouble &uF, NekDouble &vF, NekDouble &wF)
88 {
89  ASSERTL1(CheckParams("Gamma"), "Gamma not defined.");
90  const NekDouble &gamma = m_params["Gamma"]();
91 
92  // Speed of sound
93  NekDouble c = sqrt(gamma * p0 / rho0);
94 
95  // max absolute eigenvalue of the jacobian of F_n1
96  NekDouble a_1_max = std::max(std::abs(u0 - c), std::abs(u0 + c));
97 
98  NekDouble pFL = rho0*uL + u0*pL/(c*c);
99  NekDouble uFL = pL/rho0 + u0*uL + v0*vL + w0*wL;
100  NekDouble vFL = 0;
101  NekDouble wFL = 0;
102 
103  NekDouble pFR = rho0*uR + u0*pR/(c*c);
104  NekDouble uFR = pR/rho0 + u0*uR + v0*vR + w0*wR;
105  NekDouble vFR = 0;
106  NekDouble wFR = 0;
107 
108  // assemble the face-normal fluxes
109  pF = 0.5*(pFL + pFR - a_1_max*(pR/(c*c) - pL/(c*c)));
110  uF = 0.5*(uFL + uFR - a_1_max*(uR - uL));
111  vF = 0.5*(vFL + vFR - a_1_max*(vR - vL));
112  wF = 0.5*(wFL + wFR - a_1_max*(wR - wL));
113 }
114 
115 }
STL namespace.
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 p0, NekDouble rho0, NekDouble u0, NekDouble v0, NekDouble w0, NekDouble &pF, NekDouble &rhoF, NekDouble &uF, NekDouble &vF, NekDouble &wF)
Lax-Friedrichs Riemann solver.
RiemannSolverFactory & GetRiemannSolverFactory()
double NekDouble
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
SOLVER_UTILS_EXPORT bool CheckParams(std::string name)
Determine whether a parameter has been defined in m_params.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:191