Nektar++
APELaxFriedrichsSolver.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File: APELaxFriedrichsSolver.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 APE equations.
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #include <boost/core/ignore_unused.hpp>
37 
39 
40 using namespace std;
41 
42 namespace Nektar
43 {
44 
45 std::string APELaxFriedrichsSolver::solverName =
47  "APELaxFriedrichs", APELaxFriedrichsSolver::create,
48  "Lax-Friedrichs Solver");
49 
50 /**
51  *
52  */
53 APELaxFriedrichsSolver::APELaxFriedrichsSolver(
55  : AcousticSolver(pSession)
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 c0sqL Base pressure left state
73  * @param c0sqR Base pressure right state
74  * @param rho0L Base density left state
75  * @param rho0R Base density right state
76  * @param u0L Base x velocity component left state
77  * @param u0R Base x velocity component right state
78  * @param v0L Base y velocity component left state
79  * @param v0R Base y velocity component right state
80  * @param w0L Base z velocity component left state
81  * @param w0R Base z velocity component right state
82  * @param pF Computed Riemann flux for perturbation pressure
83  * @param rhoF Computed Riemann flux for perturbation density
84  * @param uF Computed Riemann flux for x perturbation velocity component
85  * @param vF Computed Riemann flux for y perturbation velocity component
86  * @param wF Computed Riemann flux for z perturbation velocity component
87  */
89  NekDouble pL, NekDouble rhoL, NekDouble uL, NekDouble vL, NekDouble wL,
90  NekDouble pR, NekDouble rhoR, NekDouble uR, NekDouble vR, NekDouble wR,
91  NekDouble c0sqL, NekDouble rho0L, NekDouble u0L, NekDouble v0L, NekDouble w0L,
92  NekDouble c0sqR, NekDouble rho0R, NekDouble u0R, NekDouble v0R, NekDouble w0R,
93  NekDouble &pF, NekDouble &rhoF, NekDouble &uF, NekDouble &vF, NekDouble &wF)
94 {
95  boost::ignore_unused(rhoL, rhoR, rhoF);
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 = rho0L * c0sqL * uL + pL * u0L;
109  NekDouble uFL = pL / rho0L + u0L * uL + v0L * vL + w0L * wL;
110  NekDouble vFL = 0;
111  NekDouble wFL = 0;
112 
113  NekDouble pFR = rho0R * c0sqR * uR + pR * u0R;
114  NekDouble uFR = pR / rho0R + u0R * uR + v0R * vR + w0R * wR;
115  NekDouble vFR = 0;
116  NekDouble wFR = 0;
117 
118  // assemble the face-normal fluxes
119  pF = 0.5 * (pFL + pFR - a_1_max * (pR - pL));
120  uF = 0.5 * (uFL + uFR - a_1_max * (uR - uL));
121  vF = 0.5 * (vFL + vFR - a_1_max * (vR - vL));
122  wF = 0.5 * (wFL + wFR - a_1_max * (wR - wL));
123 }
124 
125 } // namespace Nektar
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 c0sqL, NekDouble rho0L, NekDouble u0L, NekDouble v0L, NekDouble w0L, NekDouble c0sqR, NekDouble rho0R, NekDouble u0R, NekDouble v0R, NekDouble w0R, NekDouble &pF, NekDouble &rhoF, NekDouble &uF, NekDouble &vF, NekDouble &wF)
Lax-Friedrichs Riemann solver.
RiemannSolverFactory & GetRiemannSolverFactory()
double NekDouble
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:199
std::shared_ptr< SessionReader > SessionReaderSharedPtr