Nektar++
APEUpwindSolver.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: APEUpwindSolver.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: Upwind Riemann solver for the APE equations.
33//
34///////////////////////////////////////////////////////////////////////////////
35
37
38namespace Nektar
39{
40
43 "APEUpwind", APEUpwindSolver::create,
44 "Upwind solver for the APE equation");
45
48 : AcousticSolver(pSession)
49{
50}
51
52/**
53 * @brief Upwind Riemann solver
54 *
55 * The fluxes are straight out of sympy, so lets just hope the compiler
56 * optimizes them for us.
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 c0L = sqrt(c0sqL);
95 NekDouble c0R = sqrt(c0sqR);
96 NekDouble c0M = (c0L + c0R) / 2.0;
97
98 NekDouble u0M = (u0L + u0R) / 2.0;
99
100 pF = 0.0;
101 uF = 0.0;
102 vF = 0.0;
103 wF = 0.0;
104
105 // lambda_3
106 if (u0M - c0M > 0)
107 {
108 pF = pF + 0.5 * (c0L - u0L) *
109 (-rho0L * v0L * vL * (c0L * u0L + c0sqL) -
110 rho0L * w0L * wL * (c0L * u0L + c0sqL) +
111 (c0sqL - pow(u0L, 2)) * (rho0L * c0L * uL - pL)) /
112 (c0sqL - pow(u0L, 2));
113
114 uF = uF +
115 0.5 * (c0L - u0L) *
116 (rho0L * c0L * (c0L * u0L + c0sqL) * (v0L * vL + w0L * wL) -
117 rho0L * c0sqL * uL * (c0sqL - pow(u0L, 2)) +
118 c0L * pL * (c0sqL - pow(u0L, 2))) /
119 (rho0L * c0sqL * (c0sqL - pow(u0L, 2)));
120 }
121 else
122 {
123 pF = pF + 0.5 * (c0R - u0R) *
124 (-rho0R * v0R * vR * (c0R * u0R + c0sqR) -
125 rho0R * w0R * wR * (c0R * u0R + c0sqR) +
126 (c0sqR - pow(u0R, 2)) * (rho0R * c0R * uR - pR)) /
127 (c0sqR - pow(u0R, 2));
128
129 uF = uF +
130 0.5 * (c0R - u0R) *
131 (rho0R * c0R * (c0R * u0R + c0sqR) * (v0R * vR + w0R * wR) -
132 rho0R * c0sqR * uR * (c0sqR - pow(u0R, 2)) +
133 c0R * pR * (c0sqR - pow(u0R, 2))) /
134 (rho0R * c0sqR * (c0sqR - pow(u0R, 2)));
135 }
136
137 // lambda_4
138 if (u0M + c0M > 0)
139 {
140 pF = pF + 0.5 * (c0L + u0L) *
141 (-rho0L * v0L * vL * (c0L * u0L - c0sqL) -
142 rho0L * w0L * wL * (c0L * u0L - c0sqL) +
143 (c0sqL - pow(u0L, 2)) * (rho0L * c0L * uL + pL)) /
144 (c0sqL - pow(u0L, 2));
145
146 uF = uF +
147 0.5 * (c0L + u0L) *
148 (-rho0L * c0L * (c0L * u0L - c0sqL) * (v0L * vL + w0L * wL) +
149 rho0L * c0sqL * uL * (c0sqL - pow(u0L, 2)) +
150 c0L * pL * (c0sqL - pow(u0L, 2))) /
151 (rho0L * c0sqL * (c0sqL - pow(u0L, 2)));
152 }
153 else
154 {
155 pF = pF + 0.5 * (c0R + u0R) *
156 (-rho0R * v0R * vR * (c0R * u0R - c0sqR) -
157 rho0R * w0R * wR * (c0R * u0R - c0sqR) +
158 (c0sqR - pow(u0R, 2)) * (rho0R * c0R * uR + pR)) /
159 (c0sqR - pow(u0R, 2));
160
161 uF = uF +
162 0.5 * (c0R + u0R) *
163 (-rho0R * c0R * (c0R * u0R - c0sqR) * (v0R * vR + w0R * wR) +
164 rho0R * c0sqR * uR * (c0sqR - pow(u0R, 2)) +
165 c0R * pR * (c0sqR - pow(u0R, 2))) /
166 (rho0R * c0sqR * (c0sqR - pow(u0R, 2)));
167 }
168}
169
170} // namespace Nektar
APEUpwindSolver(const LibUtilities::SessionReaderSharedPtr &pSession)
static std::string solverName
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
Upwind Riemann solver.
static RiemannSolverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
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 > sqrt(scalarT< T > in)
Definition: scalar.hpp:285