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
36#include <boost/core/ignore_unused.hpp>
37
39
40using namespace std;
41
42namespace Nektar
43{
44
47 "APEUpwind", APEUpwindSolver::create,
48 "Upwind solver for the APE equation");
49
52 : AcousticSolver(pSession)
53{
54}
55
56/**
57 * @brief Upwind Riemann solver
58 *
59 * The fluxes are straight out of sympy, so lets just hope the compiler
60 * optimizes them for us.
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,
92 NekDouble w0L, NekDouble c0sqR, NekDouble rho0R, NekDouble u0R,
93 NekDouble v0R, NekDouble w0R, NekDouble &pF, NekDouble &rhoF, NekDouble &uF,
94 NekDouble &vF, NekDouble &wF)
95{
96 boost::ignore_unused(rhoL, rhoR, rhoF);
97
98 // Speed of sound
99 NekDouble c0L = sqrt(c0sqL);
100 NekDouble c0R = sqrt(c0sqR);
101 NekDouble c0M = (c0L + c0R) / 2.0;
102
103 NekDouble u0M = (u0L + u0R) / 2.0;
104
105 pF = 0.0;
106 uF = 0.0;
107 vF = 0.0;
108 wF = 0.0;
109
110 // lambda_3
111 if (u0M - c0M > 0)
112 {
113 pF = pF + 0.5 * (c0L - u0L) *
114 (-rho0L * v0L * vL * (c0L * u0L + c0sqL) -
115 rho0L * w0L * wL * (c0L * u0L + c0sqL) +
116 (c0sqL - pow(u0L, 2)) * (rho0L * c0L * uL - pL)) /
117 (c0sqL - pow(u0L, 2));
118
119 uF = uF +
120 0.5 * (c0L - u0L) *
121 (rho0L * c0L * (c0L * u0L + c0sqL) * (v0L * vL + w0L * wL) -
122 rho0L * c0sqL * uL * (c0sqL - pow(u0L, 2)) +
123 c0L * pL * (c0sqL - pow(u0L, 2))) /
124 (rho0L * c0sqL * (c0sqL - pow(u0L, 2)));
125 }
126 else
127 {
128 pF = pF + 0.5 * (c0R - u0R) *
129 (-rho0R * v0R * vR * (c0R * u0R + c0sqR) -
130 rho0R * w0R * wR * (c0R * u0R + c0sqR) +
131 (c0sqR - pow(u0R, 2)) * (rho0R * c0R * uR - pR)) /
132 (c0sqR - pow(u0R, 2));
133
134 uF = uF +
135 0.5 * (c0R - u0R) *
136 (rho0R * c0R * (c0R * u0R + c0sqR) * (v0R * vR + w0R * wR) -
137 rho0R * c0sqR * uR * (c0sqR - pow(u0R, 2)) +
138 c0R * pR * (c0sqR - pow(u0R, 2))) /
139 (rho0R * c0sqR * (c0sqR - pow(u0R, 2)));
140 }
141
142 // lambda_4
143 if (u0M + c0M > 0)
144 {
145 pF = pF + 0.5 * (c0L + u0L) *
146 (-rho0L * v0L * vL * (c0L * u0L - c0sqL) -
147 rho0L * w0L * wL * (c0L * u0L - c0sqL) +
148 (c0sqL - pow(u0L, 2)) * (rho0L * c0L * uL + pL)) /
149 (c0sqL - pow(u0L, 2));
150
151 uF = uF +
152 0.5 * (c0L + u0L) *
153 (-rho0L * c0L * (c0L * u0L - c0sqL) * (v0L * vL + w0L * wL) +
154 rho0L * c0sqL * uL * (c0sqL - pow(u0L, 2)) +
155 c0L * pL * (c0sqL - pow(u0L, 2))) /
156 (rho0L * c0sqL * (c0sqL - pow(u0L, 2)));
157 }
158 else
159 {
160 pF = pF + 0.5 * (c0R + u0R) *
161 (-rho0R * v0R * vR * (c0R * u0R - c0sqR) -
162 rho0R * w0R * wR * (c0R * u0R - c0sqR) +
163 (c0sqR - pow(u0R, 2)) * (rho0R * c0R * uR + pR)) /
164 (c0sqR - pow(u0R, 2));
165
166 uF = uF +
167 0.5 * (c0R + u0R) *
168 (-rho0R * c0R * (c0R * u0R - c0sqR) * (v0R * vR + w0R * wR) +
169 rho0R * c0sqR * uR * (c0sqR - pow(u0R, 2)) +
170 c0R * pR * (c0sqR - pow(u0R, 2))) /
171 (rho0R * c0sqR * (c0sqR - pow(u0R, 2)));
172 }
173}
174
175} // namespace Nektar
APEUpwindSolver(const LibUtilities::SessionReaderSharedPtr &pSession)
static std::string solverName
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) 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.
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