Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
CompressibleFlowSolver/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) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // License for the specific language governing rights and limitations under
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: LaxFriedrichs Riemann solver.
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
37 
38 namespace Nektar
39 {
42  "LaxFriedrichs",
44  "Lax-Friedrichs Riemann solver");
45 
46  LaxFriedrichsSolver::LaxFriedrichsSolver() : CompressibleSolver()
47  {
48 
49  }
50 
51  /**
52  * @brief Lax-Friedrichs Riemann solver
53  *
54  * @param rhoL Density left state.
55  * @param rhoR Density right state.
56  * @param rhouL x-momentum component left state.
57  * @param rhouR x-momentum component right state.
58  * @param rhovL y-momentum component left state.
59  * @param rhovR y-momentum component right state.
60  * @param rhowL z-momentum component left state.
61  * @param rhowR z-momentum component right state.
62  * @param EL Energy left state.
63  * @param ER Energy right state.
64  * @param rhof Computed Riemann flux for density.
65  * @param rhouf Computed Riemann flux for x-momentum component
66  * @param rhovf Computed Riemann flux for y-momentum component
67  * @param rhowf Computed Riemann flux for z-momentum component
68  * @param Ef Computed Riemann flux for energy.
69  */
71  double rhoL, double rhouL, double rhovL, double rhowL, double EL,
72  double rhoR, double rhouR, double rhovR, double rhowR, double ER,
73  double &rhof, double &rhouf, double &rhovf, double &rhowf, double &Ef)
74  {
75  static NekDouble gamma = m_params["gamma"]();
76 
77  // Left and right velocities
78  NekDouble uL = rhouL / rhoL;
79  NekDouble vL = rhovL / rhoL;
80  NekDouble wL = rhowL / rhoL;
81  NekDouble uR = rhouR / rhoR;
82  NekDouble vR = rhovR / rhoR;
83  NekDouble wR = rhowR / rhoR;
84 
85  // Left and right pressures
86  NekDouble pL = (gamma - 1.0) *
87  (EL - 0.5 * (rhouL * uL + rhovL * vL + rhowL * wL));
88  NekDouble pR = (gamma - 1.0) *
89  (ER - 0.5 * (rhouR * uR + rhovR * vR + rhowR * wR));
90 
91  // Left and right enthalpy
92  NekDouble hL = (EL + pL) / rhoL;
93  NekDouble hR = (ER + pR) / rhoR;
94 
95  // Square root of rhoL and rhoR
96  NekDouble srL = sqrt(rhoL);
97  NekDouble srR = sqrt(rhoR);
98  NekDouble srLR = srL + srR;
99 
100  // Roe average state
101  NekDouble uRoe = (srL * uL + srR * uR) / srLR;
102  NekDouble vRoe = (srL * vL + srR * vR) / srLR;
103  NekDouble wRoe = (srL * wL + srR * wR) / srLR;
104  NekDouble hRoe = (srL * hL + srR * hR) / srLR;
105  NekDouble URoe = (uRoe * uRoe + vRoe * vRoe + wRoe * wRoe);
106  NekDouble cRoe = sqrt((gamma - 1.0)*(hRoe - 0.5 * URoe));
107 
108  // Maximum eigenvalue
109  URoe = fabs(uRoe) + cRoe;
110 
111  // Lax-Friedrichs flux formula
112  rhof = 0.5*(rhouL + rhouR - URoe*(rhoR - rhoL));
113  rhouf = 0.5*(pL + rhouL*uL + pR + rhouR*uR - URoe*(rhouR - rhouL));
114  rhovf = 0.5*(rhouL*vL + rhouR*vR - URoe*(rhovR - rhovL));
115  rhowf = 0.5*(rhouL*wL + rhouR*wR - URoe*(rhowR - rhowL));
116  Ef = 0.5*(uL*(EL + pL) + uR*(ER + pR) - URoe*(ER - EL));
117  }
118 }
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.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, tDescription pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:215