Nektar++
RoeSolverSIMD.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: RoeSolverSIMD.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 using simd types.
32//
33///////////////////////////////////////////////////////////////////////////////
34
37
39
40namespace Nektar
41{
42std::string RoeSolverSIMD::solverName =
44 "RoeOpt", RoeSolverSIMD::create, "Roe Riemann solver opt");
45
48 : CompressibleSolver(pSession)
49{
50 m_requiresRotation = false;
51}
52
53/// programmatic ctor
55{
56 m_requiresRotation = false;
57}
58
59/**
60 *
61 */
63 [[maybe_unused]] const int nDim,
64 const Array<OneD, const Array<OneD, NekDouble>> &fwd,
65 const Array<OneD, const Array<OneD, NekDouble>> &bwd,
67{
68 static auto gamma = m_params["gamma"]();
69 static size_t nVars = fwd.size();
70 static size_t spaceDim = nVars - 2;
71
72 // 3D case only so far
73 ASSERTL0(spaceDim == 3, "SIMD Roe implemented only for 3D case...");
74
75 using namespace tinysimd;
76 using vec_t = simd<NekDouble>;
77 // using vec_t = typename tinysimd::abi::scalar<NekDouble>::type;
78
79 // get limit of vectorizable chunk
80 size_t sizeScalar = fwd[0].size();
81 size_t sizeVec = (sizeScalar / vec_t::width) * vec_t::width;
82
83 // get normal, vellocs
84 ASSERTL1(CheckVectors("N"), "N not defined.");
86
87 // Generate matrices if they don't already exist.
88 if (m_rotMat.size() == 0)
89 {
91 }
92
93 // SIMD loop
94 size_t i = 0;
95 for (; i < sizeVec; i += vec_t::width)
96 {
97 // load scalars
98 vec_t rhoL, rhoR, ER, EL;
99 rhoL.load(&(fwd[0][i]), is_not_aligned);
100 rhoR.load(&(bwd[0][i]), is_not_aligned);
101 ER.load(&(bwd[spaceDim + 1][i]), is_not_aligned);
102 EL.load(&(fwd[spaceDim + 1][i]), is_not_aligned);
103
104 // load vectors left
105 vec_t tmpIn[3], tmpOut[3];
106 tmpIn[0].load(&(fwd[1][i]), is_not_aligned);
107 tmpIn[1].load(&(fwd[2][i]), is_not_aligned);
108 tmpIn[2].load(&(fwd[3][i]), is_not_aligned);
109
110 // load rotation matrix
111 vec_t rotMat[9];
112 for (size_t j = 0; j < 9; ++j)
113 {
114 rotMat[j].load(&(m_rotMat[j][i]), is_not_aligned);
115 }
116
117 // rotateTo kernel Fwd
118 rotateToNormalKernel(tmpIn, rotMat, tmpOut);
119
120 vec_t rhouL = tmpOut[0];
121 vec_t rhovL = tmpOut[1];
122 vec_t rhowL = tmpOut[2];
123
124 // load vectors right
125 tmpIn[0].load(&(bwd[1][i]), is_not_aligned);
126 tmpIn[1].load(&(bwd[2][i]), is_not_aligned);
127 tmpIn[2].load(&(bwd[3][i]), is_not_aligned);
128
129 // rotateTo kernel Bwd
130 rotateToNormalKernel(tmpIn, rotMat, tmpOut);
131
132 vec_t rhouR = tmpOut[0];
133 vec_t rhovR = tmpOut[1];
134 vec_t rhowR = tmpOut[2];
135
136 // Roe kernel
137 vec_t rhof{}, Ef{};
138 RoeKernel(rhoL, rhouL, rhovL, rhowL, EL, rhoR, rhouR, rhovR, rhowR, ER,
139 rhof, tmpIn[0], tmpIn[1], tmpIn[2], Ef, gamma);
140
141 // rotateFrom kernel
142 rotateFromNormalKernel(tmpIn, rotMat, tmpOut);
143
144 // store scalar
145 rhof.store(&(flux[0][i]), is_not_aligned);
146 Ef.store(&(flux[nVars - 1][i]), is_not_aligned);
147
148 // store vector 3D only
149 tmpOut[0].store(&(flux[1][i]), is_not_aligned);
150 tmpOut[1].store(&(flux[2][i]), is_not_aligned);
151 tmpOut[2].store(&(flux[3][i]), is_not_aligned);
152 }
153
154 // spillover loop
155 for (; i < sizeScalar; ++i)
156 {
157 // load scalars
158 NekDouble rhoL = fwd[0][i];
159 NekDouble rhoR = bwd[0][i];
160 NekDouble EL = fwd[spaceDim + 1][i];
161 NekDouble ER = bwd[spaceDim + 1][i];
162
163 // 3D case only
164 // load vectors left
165 NekDouble tmpIn[3], tmpOut[3];
166 tmpIn[0] = fwd[1][i];
167 tmpIn[1] = fwd[2][i];
168 tmpIn[2] = fwd[3][i];
169
170 // load rotation matrix
171 NekDouble rotMat[9];
172 for (size_t j = 0; j < 9; ++j)
173 {
174 rotMat[j] = m_rotMat[j][i];
175 }
176
177 // rotateTo kernel Fwd
178 rotateToNormalKernel(tmpIn, rotMat, tmpOut);
179
180 NekDouble rhouL = tmpOut[0];
181 NekDouble rhovL = tmpOut[1];
182 NekDouble rhowL = tmpOut[2];
183
184 // load vectors right
185 tmpIn[0] = bwd[1][i];
186 tmpIn[1] = bwd[2][i];
187 tmpIn[2] = bwd[3][i];
188
189 // rotateTo kernel Bwd
190 rotateToNormalKernel(tmpIn, rotMat, tmpOut);
191
192 NekDouble rhouR = tmpOut[0];
193 NekDouble rhovR = tmpOut[1];
194 NekDouble rhowR = tmpOut[2];
195
196 // Roe kernel
197 NekDouble rhof{}, Ef{};
198 RoeKernel(rhoL, rhouL, rhovL, rhowL, EL, rhoR, rhouR, rhovR, rhowR, ER,
199 rhof, tmpIn[0], tmpIn[1], tmpIn[2], Ef, gamma);
200
201 // rotateFrom kernel
202 rotateFromNormalKernel(tmpIn, rotMat, tmpOut);
203
204 // store scalar
205 flux[0][i] = rhof;
206 flux[nVars - 1][i] = Ef;
207
208 // store vector 3D only
209 flux[1][i] = tmpOut[0];
210 flux[2][i] = tmpOut[1];
211 flux[3][i] = tmpOut[2];
212 }
213}
214
215} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
void v_Solve(const int nDim, const Array< OneD, const Array< OneD, ND > > &Fwd, const Array< OneD, const Array< OneD, ND > > &Bwd, Array< OneD, Array< OneD, ND > > &flux) final
static RiemannSolverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
Definition: RoeSolverSIMD.h:46
static std::string solverName
Definition: RoeSolverSIMD.h:52
RoeSolverSIMD()
programmatic ctor
bool m_requiresRotation
Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields.
Array< OneD, Array< OneD, NekDouble > > m_rotMat
Rotation matrices for each trace quadrature point.
std::map< std::string, RSVecFuncType > m_vectors
Map of vector function types.
SOLVER_UTILS_EXPORT bool CheckVectors(std::string name)
Determine whether a vector has been defined in m_vectors.
SOLVER_UTILS_EXPORT void GenerateRotationMatrices(const Array< OneD, const Array< OneD, NekDouble > > &normals)
Generate rotation matrices for 3D expansions.
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
void rotateFromNormalKernel(T *in, T *rotMat, T *out)
void rotateToNormalKernel(T *in, T *rotMat, T *out)
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