Nektar++
RoeSolver.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: RoeSolver.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Roe Riemann solver.
32//
33///////////////////////////////////////////////////////////////////////////////
34
36
38
39namespace Nektar
40{
41std::string RoeSolver::solverName =
43 "Roe", RoeSolver::create, "Roe Riemann solver");
44
46 : CompressibleSolver(pSession)
47{
48 // m_pointSolve = false;
49}
50
51/// programmatic ctor
53{
54 // m_pointSolve = false;
55}
56
57/**
58 * @brief Roe Riemann solver.
59 *
60 * Stated equations numbers are from:
61 *
62 * "Riemann Solvers and Numerical Methods for Fluid Dynamics: A Practical
63 * Introduction", E. F. Toro (3rd edition, 2009).
64 *
65 * We follow the algorithm prescribed following equation 11.70.
66 *
67 * @param rhoL Density left state.
68 * @param rhoR Density right state.
69 * @param rhouL x-momentum component left state.
70 * @param rhouR x-momentum component right state.
71 * @param rhovL y-momentum component left state.
72 * @param rhovR y-momentum component right state.
73 * @param rhowL z-momentum component left state.
74 * @param rhowR z-momentum component right state.
75 * @param EL Energy left state.
76 * @param ER Energy right state.
77 * @param rhof Computed Riemann flux for density.
78 * @param rhouf Computed Riemann flux for x-momentum component
79 * @param rhovf Computed Riemann flux for y-momentum component
80 * @param rhowf Computed Riemann flux for z-momentum component
81 * @param Ef Computed Riemann flux for energy.
82 */
83void RoeSolver::v_PointSolve(double rhoL, double rhouL, double rhovL,
84 double rhowL, double EL, double rhoR, double rhouR,
85 double rhovR, double rhowR, double ER,
86 double &rhof, double &rhouf, double &rhovf,
87 double &rhowf, double &Ef)
88{
89 static NekDouble gamma = m_params["gamma"]();
90
91 RoeKernel(rhoL, rhouL, rhovL, rhowL, EL, rhoR, rhouR, rhovR, rhowR, ER,
92 rhof, rhouf, rhovf, rhowf, Ef, gamma);
93}
94
95/**
96 *
97 */
99 const Array<OneD, const Array<OneD, NekDouble>> &fwd,
100 const Array<OneD, const Array<OneD, NekDouble>> &bwd,
102{
103 static auto gamma = m_params["gamma"]();
104 static size_t nVars = fwd.size();
105 static size_t spaceDim = nVars - 2;
106
107 using namespace tinysimd;
108 using vec_t = simd<NekDouble>;
109
110 // get limit of vectorizable chunk
111 size_t sizeScalar = fwd[0].size();
112 size_t sizeVec = (sizeScalar / vec_t::width) * vec_t::width;
113
114 // SIMD loop
115 size_t i = 0;
116 for (; i < sizeVec; i += vec_t::width)
117 {
118 vec_t rhoL{}, rhouL{}, rhovL{}, rhowL{}, EL{};
119 vec_t rhoR{}, rhouR{}, rhovR{}, rhowR{}, ER{};
120
121 // load
122 rhoL.load(&(fwd[0][i]), is_not_aligned);
123 rhouL.load(&(fwd[1][i]), is_not_aligned);
124 EL.load(&(fwd[spaceDim + 1][i]), is_not_aligned);
125 rhoR.load(&(bwd[0][i]), is_not_aligned);
126 rhouR.load(&(bwd[1][i]), is_not_aligned);
127 ER.load(&(bwd[spaceDim + 1][i]), is_not_aligned);
128
129 if (spaceDim == 2)
130 {
131 rhovL.load(&(fwd[2][i]), is_not_aligned);
132 rhovR.load(&(bwd[2][i]), is_not_aligned);
133 }
134 else if (spaceDim == 3)
135 {
136 rhovL.load(&(fwd[2][i]), is_not_aligned);
137 rhowL.load(&(fwd[3][i]), is_not_aligned);
138 rhovR.load(&(bwd[2][i]), is_not_aligned);
139 rhowR.load(&(bwd[3][i]), is_not_aligned);
140 }
141
142 vec_t rhof{}, rhouf{}, rhovf{}, rhowf{}, Ef{};
143
144 RoeKernel(rhoL, rhouL, rhovL, rhowL, EL, rhoR, rhouR, rhovR, rhowR, ER,
145 rhof, rhouf, rhovf, rhowf, Ef, gamma);
146
147 // store
148 rhof.store(&(flux[0][i]), is_not_aligned);
149 rhouf.store(&(flux[1][i]), is_not_aligned);
150 Ef.store(&(flux[nVars - 1][i]), is_not_aligned);
151 if (spaceDim == 2)
152 {
153 rhovf.store(&(flux[2][i]), is_not_aligned);
154 }
155 else if (spaceDim == 3)
156 {
157 rhovf.store(&(flux[2][i]), is_not_aligned);
158 rhowf.store(&(flux[3][i]), is_not_aligned);
159 }
160
161 } // avx loop
162
163 // spillover loop
164 for (; i < sizeScalar; ++i)
165 {
166 NekDouble rhoL{}, rhouL{}, rhovL{}, rhowL{}, EL{};
167 NekDouble rhoR{}, rhouR{}, rhovR{}, rhowR{}, ER{};
168
169 // load
170 rhoL = fwd[0][i];
171 rhouL = fwd[1][i];
172 EL = fwd[spaceDim + 1][i];
173 rhoR = bwd[0][i];
174 rhouR = bwd[1][i];
175 ER = bwd[spaceDim + 1][i];
176
177 if (spaceDim == 2)
178 {
179 rhovL = fwd[2][i];
180 rhovR = bwd[2][i];
181 }
182 else if (spaceDim == 3)
183 {
184 rhovL = fwd[2][i];
185 rhowL = fwd[3][i];
186 rhovR = bwd[2][i];
187 rhowR = bwd[3][i];
188 }
189
190 NekDouble rhof{}, rhouf{}, rhovf{}, rhowf{}, Ef{};
191
192 RoeKernel(rhoL, rhouL, rhovL, rhowL, EL, rhoR, rhouR, rhovR, rhowR, ER,
193 rhof, rhouf, rhovf, rhowf, Ef, gamma);
194
195 // store
196 flux[0][i] = rhof;
197 flux[1][i] = rhouf;
198 flux[nVars - 1][i] = Ef;
199 if (spaceDim == 2)
200 {
201 flux[2][i] = rhovf;
202 }
203 else if (spaceDim == 3)
204 {
205 flux[2][i] = rhovf;
206 flux[3][i] = rhowf;
207 }
208
209 } // loop
210}
211
212} // namespace Nektar
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
void v_PointSolve(ND rhoL, ND rhouL, ND rhovL, ND rhowL, ND EL, ND rhoR, ND rhouR, ND rhovR, ND rhowR, ND ER, ND &rhof, ND &rhouf, ND &rhovf, ND &rhowf, ND &Ef) final
Roe Riemann solver.
Definition: RoeSolver.cpp:83
void v_ArraySolve(const Array< OneD, const Array< OneD, ND > > &Fwd, const Array< OneD, const Array< OneD, ND > > &Bwd, Array< OneD, Array< OneD, ND > > &flux) final
Definition: RoeSolver.cpp:98
RoeSolver()
programmatic ctor
Definition: RoeSolver.cpp:52
static RiemannSolverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
Definition: RoeSolver.h:45
static std::string solverName
Definition: RoeSolver.h:51
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
RiemannSolverFactory & GetRiemannSolverFactory()
void RoeKernel(T &rhoL, T &rhouL, T &rhovL, T &rhowL, T &EL, T &rhoR, T &rhouR, T &rhovR, T &rhowR, T &ER, T &rhof, T &rhouf, T &rhovf, T &rhowf, T &Ef, NekDouble gamma)
Definition: RoeSolver.h:73
tinysimd::simd< NekDouble > vec_t
double NekDouble
static constexpr struct tinysimd::is_not_aligned_t is_not_aligned
typename abi< ScalarType, width >::type simd
Definition: tinysimd.hpp:80