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 const int nDim, const Array<OneD, const Array<OneD, NekDouble>> &fwd,
64 const Array<OneD, const Array<OneD, NekDouble>> &bwd,
66{
67 boost::ignore_unused(nDim);
68
69 static auto gamma = m_params["gamma"]();
70 static size_t nVars = fwd.size();
71 static size_t spaceDim = nVars - 2;
72
73 // 3D case only so far
74 ASSERTL0(spaceDim == 3, "SIMD Roe implemented only for 3D case...");
75
76 using namespace tinysimd;
77 using vec_t = simd<NekDouble>;
78 // using vec_t = typename tinysimd::abi::scalar<NekDouble>::type;
79
80 // get limit of vectorizable chunk
81 size_t sizeScalar = fwd[0].size();
82 size_t sizeVec = (sizeScalar / vec_t::width) * vec_t::width;
83
84 // get normal, vellocs
85 ASSERTL1(CheckVectors("N"), "N not defined.");
87
88 // Generate matrices if they don't already exist.
89 if (m_rotMat.size() == 0)
90 {
92 }
93
94 // SIMD loop
95 size_t i = 0;
96 for (; i < sizeVec; i += vec_t::width)
97 {
98 // load scalars
99 vec_t rhoL, rhoR, ER, EL;
100 rhoL.load(&(fwd[0][i]), is_not_aligned);
101 rhoR.load(&(bwd[0][i]), is_not_aligned);
102 ER.load(&(bwd[spaceDim + 1][i]), is_not_aligned);
103 EL.load(&(fwd[spaceDim + 1][i]), is_not_aligned);
104
105 // load vectors left
106 vec_t tmpIn[3], tmpOut[3];
107 tmpIn[0].load(&(fwd[1][i]), is_not_aligned);
108 tmpIn[1].load(&(fwd[2][i]), is_not_aligned);
109 tmpIn[2].load(&(fwd[3][i]), is_not_aligned);
110
111 // load rotation matrix
112 vec_t rotMat[9];
113 for (size_t j = 0; j < 9; ++j)
114 {
115 rotMat[j].load(&(m_rotMat[j][i]), is_not_aligned);
116 }
117
118 // rotateTo kernel Fwd
119 rotateToNormalKernel(tmpIn, rotMat, tmpOut);
120
121 vec_t rhouL = tmpOut[0];
122 vec_t rhovL = tmpOut[1];
123 vec_t rhowL = tmpOut[2];
124
125 // load vectors right
126 tmpIn[0].load(&(bwd[1][i]), is_not_aligned);
127 tmpIn[1].load(&(bwd[2][i]), is_not_aligned);
128 tmpIn[2].load(&(bwd[3][i]), is_not_aligned);
129
130 // rotateTo kernel Bwd
131 rotateToNormalKernel(tmpIn, rotMat, tmpOut);
132
133 vec_t rhouR = tmpOut[0];
134 vec_t rhovR = tmpOut[1];
135 vec_t rhowR = tmpOut[2];
136
137 // Roe kernel
138 vec_t rhof{}, Ef{};
139 RoeKernel(rhoL, rhouL, rhovL, rhowL, EL, rhoR, rhouR, rhovR, rhowR, ER,
140 rhof, tmpIn[0], tmpIn[1], tmpIn[2], Ef, gamma);
141
142 // rotateFrom kernel
143 rotateFromNormalKernel(tmpIn, rotMat, tmpOut);
144
145 // store scalar
146 rhof.store(&(flux[0][i]), is_not_aligned);
147 Ef.store(&(flux[nVars - 1][i]), is_not_aligned);
148
149 // store vector 3D only
150 tmpOut[0].store(&(flux[1][i]), is_not_aligned);
151 tmpOut[1].store(&(flux[2][i]), is_not_aligned);
152 tmpOut[2].store(&(flux[3][i]), is_not_aligned);
153 }
154
155 // spillover loop
156 for (; i < sizeScalar; ++i)
157 {
158 // load scalars
159 NekDouble rhoL = fwd[0][i];
160 NekDouble rhoR = bwd[0][i];
161 NekDouble EL = fwd[spaceDim + 1][i];
162 NekDouble ER = bwd[spaceDim + 1][i];
163
164 // 3D case only
165 // load vectors left
166 NekDouble tmpIn[3], tmpOut[3];
167 tmpIn[0] = fwd[1][i];
168 tmpIn[1] = fwd[2][i];
169 tmpIn[2] = fwd[3][i];
170
171 // load rotation matrix
172 NekDouble rotMat[9];
173 for (size_t j = 0; j < 9; ++j)
174 {
175 rotMat[j] = m_rotMat[j][i];
176 }
177
178 // rotateTo kernel Fwd
179 rotateToNormalKernel(tmpIn, rotMat, tmpOut);
180
181 NekDouble rhouL = tmpOut[0];
182 NekDouble rhovL = tmpOut[1];
183 NekDouble rhowL = tmpOut[2];
184
185 // load vectors right
186 tmpIn[0] = bwd[1][i];
187 tmpIn[1] = bwd[2][i];
188 tmpIn[2] = bwd[3][i];
189
190 // rotateTo kernel Bwd
191 rotateToNormalKernel(tmpIn, rotMat, tmpOut);
192
193 NekDouble rhouR = tmpOut[0];
194 NekDouble rhovR = tmpOut[1];
195 NekDouble rhowR = tmpOut[2];
196
197 // Roe kernel
198 NekDouble rhof{}, Ef{};
199 RoeKernel(rhoL, rhouL, rhovL, rhowL, EL, rhoR, rhouR, rhovR, rhowR, ER,
200 rhof, tmpIn[0], tmpIn[1], tmpIn[2], Ef, gamma);
201
202 // rotateFrom kernel
203 rotateFromNormalKernel(tmpIn, rotMat, tmpOut);
204
205 // store scalar
206 flux[0][i] = rhof;
207 flux[nVars - 1][i] = Ef;
208
209 // store vector 3D only
210 flux[1][i] = tmpOut[0];
211 flux[2][i] = tmpOut[1];
212 flux[3][i] = tmpOut[2];
213 }
214}
215
216} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
static RiemannSolverSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession)
Definition: RoeSolverSIMD.h:46
virtual 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) override final
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()
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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:75
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