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
36#include <boost/core/ignore_unused.hpp>
37
39
40using namespace std;
41
42namespace Nektar
43{
44
47 "LEEUpwind", LEEUpwindSolver::create, "Upwind Solver for LEE");
48
49/**
50 *
51 */
54 : LEESolver(pSession)
55{
56}
57
58/**
59 * @brief Lax-Friedrichs Riemann solver
60 *
61 * @param pL Perturbation pressure left state
62 * @param rhoL Perturbation density left state
63 * @param pR Perturbation pressure right state
64 * @param rhoR Perturbation density right state
65 * @param rhouL x perturbation velocity component left state
66 * @param rhouR x perturbation velocity component right state
67 * @param rhovL y perturbation velocity component left state
68 * @param rhovR y perturbation velocity component right state
69 * @param rhowL z perturbation velocity component left state
70 * @param rhowR z perturbation velocity component right state
71 * @param c0sqL Base pressure left state
72 * @param c0sqR Base pressure right state
73 * @param rho0L Base density left state
74 * @param rho0R Base density right state
75 * @param u0L Base x velocity component left state
76 * @param u0R Base x velocity component right state
77 * @param v0L Base y velocity component left state
78 * @param v0R Base y velocity component right state
79 * @param w0L Base z velocity component left state
80 * @param w0R Base z velocity component right state
81 * @param pF Computed Riemann flux for perturbation pressure
82 * @param rhoF Computed Riemann flux for perturbation density
83 * @param rhouF Computed Riemann flux for x perturbation velocity component
84 * @param rhovF Computed Riemann flux for y perturbation velocity component
85 * @param rhowF Computed Riemann flux for z perturbation velocity component
86 */
88 NekDouble pL, NekDouble rhoL, NekDouble rhouL, NekDouble rhovL,
89 NekDouble rhowL, NekDouble pR, NekDouble rhoR, NekDouble rhouR,
90 NekDouble rhovR, NekDouble rhowR, NekDouble c0sqL, NekDouble rho0L,
91 NekDouble u0L, NekDouble v0L, NekDouble w0L, NekDouble c0sqR,
92 NekDouble rho0R, NekDouble u0R, NekDouble v0R, NekDouble w0R, NekDouble &pF,
93 NekDouble &rhoF, NekDouble &rhouF, NekDouble &rhovF, NekDouble &rhowF)
94{
95 boost::ignore_unused(rho0L, v0L, w0L, rho0R, v0R, w0R);
96
97 // Speed of sound
98 NekDouble c0L = sqrt(c0sqL);
99 NekDouble c0R = sqrt(c0sqR);
100 NekDouble c0M = (c0L + c0R) / 2.0;
101
102 NekDouble u0M = (u0L + u0R) / 2.0;
103
104 pF = 0.0;
105 rhoF = 0.0;
106 rhouF = 0.0;
107 rhovF = 0.0;
108 rhowF = 0.0;
109
110 // lambda_1,2,3
111 if (u0M > 0)
112 {
113 rhoF = rhoF + u0L * (c0sqL * rhoL - pL) / c0sqL;
114 rhovF = rhovF + rhovL * u0L;
115 rhowF = rhowF + rhowL * u0L;
116 }
117 else
118 {
119 rhoF = rhoF + u0R * (c0sqR * rhoR - pR) / c0sqR;
120 rhovF = rhovF + rhovR * u0R;
121 rhowF = rhowF + rhowR * u0R;
122 }
123
124 // lambda_4
125 if (u0M - c0M > 0)
126 {
127 pF = pF + 0.5 * (c0L - u0L) * (rhouL * c0L - pL);
128 rhoF = rhoF + 0.5 * (c0L - u0L) * (rhouL * c0sqL - c0L * pL) /
129 pow(c0sqL, 3.0 / 2.0);
130 rhouF = rhouF + 0.5 * (c0L - u0L) * (-rhouL * c0L + pL) / c0L;
131 }
132 else
133 {
134 pF = pF + 0.5 * (c0R - u0R) * (rhouR * c0R - pR);
135 rhoF = rhoF + 0.5 * (c0R - u0R) * (rhouR * c0sqR - c0R * pR) /
136 pow(c0sqR, 3.0 / 2.0);
137 rhouF = rhouF + 0.5 * (c0R - u0R) * (-rhouR * c0R + pR) / c0R;
138 }
139
140 // lambda_5
141 if (u0M + c0M > 0)
142 {
143 pF = pF + 0.5 * (c0L + u0L) * (rhouL * c0L + pL);
144 rhoF = rhoF + 0.5 * (c0L + u0L) * (rhouL * c0sqL + c0L * pL) /
145 pow(c0sqL, 3.0 / 2.0);
146 rhouF = rhouF + 0.5 * (c0L + u0L) * (rhouL * c0L + pL) / c0L;
147 }
148 else
149 {
150 pF = pF + 0.5 * (c0R + u0R) * (rhouR * c0R + pR);
151 rhoF = rhoF + 0.5 * (c0R + u0R) * (rhouR * c0sqR + c0R * pR) /
152 pow(c0sqR, 3.0 / 2.0);
153 rhouF = rhouF + 0.5 * (c0R + u0R) * (rhouR * c0R + pR) / c0R;
154 }
155}
156} // namespace Nektar
static RiemannSolverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
virtual 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:198
std::shared_ptr< SessionReader > SessionReaderSharedPtr
RiemannSolverFactory & GetRiemannSolverFactory()
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294