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
37
38namespace Nektar
39{
40
43 "APELaxFriedrichs", APELaxFriedrichsSolver::create,
44 "Lax-Friedrichs Solver");
45
46/**
47 *
48 */
51 : AcousticSolver(pSession)
52{
53}
54
55/**
56 * @brief Lax-Friedrichs Riemann solver
57 *
58 * @param pL Perturbation pressure left state
59 * @param rhoL Perturbation density left state
60 * @param pR Perturbation pressure right state
61 * @param rhoR Perturbation density right state
62 * @param uL x perturbation velocity component left state
63 * @param uR x perturbation velocity component right state
64 * @param vL y perturbation velocity component left state
65 * @param vR y perturbation velocity component right state
66 * @param wL z perturbation velocity component left state
67 * @param wR z perturbation velocity component right state
68 * @param c0sqL Base pressure left state
69 * @param c0sqR Base pressure right state
70 * @param rho0L Base density left state
71 * @param rho0R Base density right state
72 * @param u0L Base x velocity component left state
73 * @param u0R Base x velocity component right state
74 * @param v0L Base y velocity component left state
75 * @param v0R Base y velocity component right state
76 * @param w0L Base z velocity component left state
77 * @param w0R Base z velocity component right state
78 * @param pF Computed Riemann flux for perturbation pressure
79 * @param rhoF Computed Riemann flux for perturbation density
80 * @param uF Computed Riemann flux for x perturbation velocity component
81 * @param vF Computed Riemann flux for y perturbation velocity component
82 * @param wF Computed Riemann flux for z perturbation velocity component
83 */
85 NekDouble pL, [[maybe_unused]] NekDouble rhoL, NekDouble uL, NekDouble vL,
86 NekDouble wL, NekDouble pR, [[maybe_unused]] NekDouble rhoR, NekDouble uR,
87 NekDouble vR, NekDouble wR, NekDouble c0sqL, NekDouble rho0L, NekDouble u0L,
88 NekDouble v0L, NekDouble w0L, NekDouble c0sqR, NekDouble rho0R,
89 NekDouble u0R, NekDouble v0R, NekDouble w0R, NekDouble &pF,
90 [[maybe_unused]] NekDouble &rhoF, NekDouble &uF, NekDouble &vF,
91 NekDouble &wF)
92{
93 // Speed of sound
94 NekDouble cL = sqrt(c0sqL);
95 NekDouble cR = sqrt(c0sqR);
96
97 // max absolute eigenvalue of the jacobian of F_n1
98 NekDouble a_1_max = 0;
99 a_1_max = std::max(a_1_max, std::abs(u0L - cL));
100 a_1_max = std::max(a_1_max, std::abs(u0R - cR));
101 a_1_max = std::max(a_1_max, std::abs(u0L + cL));
102 a_1_max = std::max(a_1_max, std::abs(u0R + cR));
103
104 NekDouble pFL = rho0L * c0sqL * uL + pL * u0L;
105 NekDouble uFL = pL / rho0L + u0L * uL + v0L * vL + w0L * wL;
106 NekDouble vFL = 0;
107 NekDouble wFL = 0;
108
109 NekDouble pFR = rho0R * c0sqR * uR + pR * u0R;
110 NekDouble uFR = pR / rho0R + u0R * uR + v0R * vR + w0R * wR;
111 NekDouble vFR = 0;
112 NekDouble wFR = 0;
113
114 // assemble the face-normal fluxes
115 pF = 0.5 * (pFL + pFR - a_1_max * (pR - pL));
116 uF = 0.5 * (uFL + uFR - a_1_max * (uR - uL));
117 vF = 0.5 * (vFL + vFR - a_1_max * (vR - vL));
118 wF = 0.5 * (wFL + wFR - a_1_max * (wR - wL));
119}
120
121} // namespace Nektar
static RiemannSolverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
APELaxFriedrichsSolver(const LibUtilities::SessionReaderSharedPtr &pSession)
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) override
Lax-Friedrichs Riemann solver.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
RiemannSolverFactory & GetRiemannSolverFactory()
double NekDouble
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:289
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:285