Nektar++
LEEUpwindSolver.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: LEEUpwindSolver.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2017 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 LEE equations.
33//
34///////////////////////////////////////////////////////////////////////////////
35
37
38using namespace std;
39
40namespace Nektar
41{
42
45 "LEEUpwind", LEEUpwindSolver::create, "Upwind Solver for LEE");
46
47/**
48 *
49 */
52 : LEESolver(pSession)
53{
54}
55
56/**
57 * @brief Lax-Friedrichs Riemann solver
58 *
59 * @param pL Perturbation pressure left state
60 * @param rhoL Perturbation density left state
61 * @param pR Perturbation pressure right state
62 * @param rhoR Perturbation density right state
63 * @param rhouL x perturbation velocity component left state
64 * @param rhouR x perturbation velocity component right state
65 * @param rhovL y perturbation velocity component left state
66 * @param rhovR y perturbation velocity component right state
67 * @param rhowL z perturbation velocity component left state
68 * @param rhowR z perturbation velocity component right state
69 * @param c0sqL Base pressure left state
70 * @param c0sqR Base pressure right state
71 * @param rho0L Base density left state
72 * @param rho0R Base density right state
73 * @param u0L Base x velocity component left state
74 * @param u0R Base x velocity component right state
75 * @param v0L Base y velocity component left state
76 * @param v0R Base y velocity component right state
77 * @param w0L Base z velocity component left state
78 * @param w0R Base z velocity component right state
79 * @param pF Computed Riemann flux for perturbation pressure
80 * @param rhoF Computed Riemann flux for perturbation density
81 * @param rhouF Computed Riemann flux for x perturbation velocity component
82 * @param rhovF Computed Riemann flux for y perturbation velocity component
83 * @param rhowF Computed Riemann flux for z perturbation velocity component
84 */
86 NekDouble pL, NekDouble rhoL, NekDouble rhouL, NekDouble rhovL,
87 NekDouble rhowL, NekDouble pR, NekDouble rhoR, NekDouble rhouR,
88 NekDouble rhovR, NekDouble rhowR, NekDouble c0sqL,
89 [[maybe_unused]] NekDouble rho0L, NekDouble u0L,
90 [[maybe_unused]] NekDouble v0L, [[maybe_unused]] NekDouble w0L,
91 NekDouble c0sqR, [[maybe_unused]] NekDouble rho0R, NekDouble u0R,
92 [[maybe_unused]] NekDouble v0R, [[maybe_unused]] NekDouble w0R,
93 NekDouble &pF, NekDouble &rhoF, NekDouble &rhouF, NekDouble &rhovF,
94 NekDouble &rhowF)
95{
96 // Speed of sound
97 NekDouble c0L = sqrt(c0sqL);
98 NekDouble c0R = sqrt(c0sqR);
99 NekDouble c0M = (c0L + c0R) / 2.0;
100
101 NekDouble u0M = (u0L + u0R) / 2.0;
102
103 pF = 0.0;
104 rhoF = 0.0;
105 rhouF = 0.0;
106 rhovF = 0.0;
107 rhowF = 0.0;
108
109 // lambda_1,2,3
110 if (u0M > 0)
111 {
112 rhoF = rhoF + u0L * (c0sqL * rhoL - pL) / c0sqL;
113 rhovF = rhovF + rhovL * u0L;
114 rhowF = rhowF + rhowL * u0L;
115 }
116 else
117 {
118 rhoF = rhoF + u0R * (c0sqR * rhoR - pR) / c0sqR;
119 rhovF = rhovF + rhovR * u0R;
120 rhowF = rhowF + rhowR * u0R;
121 }
122
123 // lambda_4
124 if (u0M - c0M > 0)
125 {
126 pF = pF + 0.5 * (c0L - u0L) * (rhouL * c0L - pL);
127 rhoF = rhoF + 0.5 * (c0L - u0L) * (rhouL * c0sqL - c0L * pL) /
128 pow(c0sqL, 3.0 / 2.0);
129 rhouF = rhouF + 0.5 * (c0L - u0L) * (-rhouL * c0L + pL) / c0L;
130 }
131 else
132 {
133 pF = pF + 0.5 * (c0R - u0R) * (rhouR * c0R - pR);
134 rhoF = rhoF + 0.5 * (c0R - u0R) * (rhouR * c0sqR - c0R * pR) /
135 pow(c0sqR, 3.0 / 2.0);
136 rhouF = rhouF + 0.5 * (c0R - u0R) * (-rhouR * c0R + pR) / c0R;
137 }
138
139 // lambda_5
140 if (u0M + c0M > 0)
141 {
142 pF = pF + 0.5 * (c0L + u0L) * (rhouL * c0L + pL);
143 rhoF = rhoF + 0.5 * (c0L + u0L) * (rhouL * c0sqL + c0L * pL) /
144 pow(c0sqL, 3.0 / 2.0);
145 rhouF = rhouF + 0.5 * (c0L + u0L) * (rhouL * c0L + pL) / c0L;
146 }
147 else
148 {
149 pF = pF + 0.5 * (c0R + u0R) * (rhouR * c0R + pR);
150 rhoF = rhoF + 0.5 * (c0R + u0R) * (rhouR * c0sqR + c0R * pR) /
151 pow(c0sqR, 3.0 / 2.0);
152 rhouF = rhouF + 0.5 * (c0R + u0R) * (rhouR * c0R + pR) / c0R;
153 }
154}
155} // namespace Nektar
static RiemannSolverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
void v_PointSolve(NekDouble pL, NekDouble rhoL, NekDouble rhouL, NekDouble rhovL, NekDouble rhowL, NekDouble pR, NekDouble rhoR, NekDouble rhouR, NekDouble rhovR, NekDouble rhowR, 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 &rhouF, NekDouble &rhovF, NekDouble &rhowF) override
Lax-Friedrichs Riemann solver.
static std::string solverName
LEEUpwindSolver(const LibUtilities::SessionReaderSharedPtr &pSession)
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
std::shared_ptr< SessionReader > SessionReaderSharedPtr
RiemannSolverFactory & GetRiemannSolverFactory()
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294