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) 2014 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 = SolverUtils::GetRiemannSolverFactory().
45  RegisterCreatorFunction("LaxFriedrichs", LaxFriedrichsSolver::create, "Lax-Friedrichs Solver");
46 
47 
48 /**
49 *
50 */
51 LaxFriedrichsSolver::LaxFriedrichsSolver() :
52  APESolver()
53 {
54 
55 }
56 
57 /**
58  * @brief Lax-Friedrichs Riemann solver
59  *
60  * @param pL Perturbation pressure left state.
61  * @param pR Perturbation pressure right state.
62  * @param uL x perturbation verlocity component left state.
63  * @param uR x perturbation verlocity component right state.
64  * @param vL y perturbation verlocity component left state.
65  * @param vR y perturbation verlocity component right state.
66  * @param wL z perturbation verlocity component left state.
67  * @param wR z perturbation verlocity component right state.
68  * @param p0 Base pressure.
69  * @param u0 Base x verlocity component
70  * @param v0 Base y verlocity component
71  * @param w0 Base z verlocity component
72  * @param pF Computed Riemann flux for perturbation pressure.
73  * @param uF Computed Riemann flux for x perturbation verlocity component
74  * @param vF Computed Riemann flux for y perturbation verlocity component
75  * @param wF Computed Riemann flux for z perturbation verlocity component
76  */
78  NekDouble pL, NekDouble uL, NekDouble vL, NekDouble wL,
79  NekDouble pR, NekDouble uR, NekDouble vR, NekDouble wR,
80  NekDouble p0, NekDouble u0, NekDouble v0, NekDouble w0,
81  NekDouble &pF, NekDouble &uF, NekDouble &vF, NekDouble &wF)
82 {
83  ASSERTL1(CheckParams("Gamma"), "Gamma not defined.");
84  ASSERTL1(CheckParams("Rho"), "Rho not defined.");
85  const NekDouble &gamma= m_params["Gamma"]();
86  const NekDouble &rho = m_params["Rho"]();
87 
88  // Speed of sound
89  NekDouble c = sqrt(gamma*p0/rho);
90 
91  // max absolute eigenvalue of the jacobian of F_n1
92  NekDouble a_1_max = std::max(std::abs(u0 - c), std::abs(u0 + c));
93 
94  NekDouble pFL = gamma*p0*uL + pL*u0;
95  NekDouble uFL = pL/rho + u0*uL + v0*vL + w0*wL;
96  NekDouble vFL = 0;
97  NekDouble wFL = 0;
98 
99  NekDouble pFR = gamma*p0*uR + pR*u0;
100  NekDouble uFR = pR/rho + u0*uR + v0*vR + w0*wR;
101  NekDouble vFR = 0;
102  NekDouble wFR = 0;
103 
104  // assemble the face-normal fluxes
105  pF = 0.5*(pFL + pFR - a_1_max*(pR - pL));
106  uF = 0.5*(uFL + uFR - a_1_max*(uR - uL));
107  vF = 0.5*(vFL + vFR - a_1_max*(vR - vL));
108  wF = 0.5*(wFL + wFR - a_1_max*(wR - wL));
109 }
110 
111 }