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 
40 namespace Nektar
41 {
42 std::string RoeSolverSIMD::solverName =
44  "RoeOpt",
46  "Roe Riemann solver opt");
47 
49  :
50  CompressibleSolver(pSession)
51 {
52  m_requiresRotation = false;
53 }
54 
55 /// programmatic ctor
58 {
59  m_requiresRotation = false;
60 }
61 
63  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 auto nVars = fwd.size();
70  static auto 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.");
85  // ASSERTL1(CheckAuxVec("vecLocs"), "vecLocs not defined.");
87  m_vectors["N"]();
88  // const Array<OneD, const Array<OneD, NekDouble> > vecLocs =
89  // m_auxVec["vecLocs"]();
90 
91  // const unsigned int vx = (int)vecLocs[0][0];
92  // const unsigned int vy = (int)vecLocs[0][1];
93  // const unsigned int vz = (int)vecLocs[0][2];
94 
95  // Generate matrices if they don't already exist.
96  if (m_rotMat.size() == 0)
97  {
98  GenerateRotationMatrices(normals);
99  }
100 
101  // SIMD loop
102  size_t i = 0;
103  for (; i < sizeVec; i+=vec_t::width)
104  {
105  // load scalars
106  vec_t rhoL, rhoR, ER, EL;
107  rhoL.load(&(fwd[0][i]), is_not_aligned);
108  rhoR.load(&(bwd[0][i]), is_not_aligned);
109  ER.load(&(bwd[spaceDim+1][i]), is_not_aligned);
110  EL.load(&(fwd[spaceDim+1][i]), is_not_aligned);
111 
112  // load vectors left
113  vec_t tmpIn[3], tmpOut[3];
114  tmpIn[0].load(&(fwd[1][i]), is_not_aligned);
115  tmpIn[1].load(&(fwd[2][i]), is_not_aligned);
116  tmpIn[2].load(&(fwd[3][i]), is_not_aligned);
117 
118  // load rotation matrix
119  vec_t rotMat[9];
120  for (size_t j = 0; j < 9; ++j)
121  {
122  rotMat[j].load(&(m_rotMat[j][i]), is_not_aligned);
123  }
124 
125  // rotateTo kernel Fwd
126  rotateToNormalKernel(tmpIn, rotMat, tmpOut);
127 
128  vec_t rhouL = tmpOut[0];
129  vec_t rhovL = tmpOut[1];
130  vec_t rhowL = tmpOut[2];
131 
132  // load vectors right
133  tmpIn[0].load(&(bwd[1][i]), is_not_aligned);
134  tmpIn[1].load(&(bwd[2][i]), is_not_aligned);
135  tmpIn[2].load(&(bwd[3][i]), is_not_aligned);
136 
137  // rotateTo kernel Bwd
138  rotateToNormalKernel(tmpIn, rotMat, tmpOut);
139 
140  vec_t rhouR = tmpOut[0];
141  vec_t rhovR = tmpOut[1];
142  vec_t rhowR = tmpOut[2];
143 
144  // Roe kernel
145  vec_t rhof{}, Ef{};
146  RoeKernel(
147  rhoL, rhouL, rhovL, rhowL, EL,
148  rhoR, rhouR, rhovR, rhowR, ER,
149  rhof, tmpIn[0], tmpIn[1], tmpIn[2], Ef,
150  gamma);
151 
152  // rotateFrom kernel
153  rotateFromNormalKernel(tmpIn, rotMat, tmpOut);
154 
155  // store scalar
156  rhof.store(&(flux[0][i]), is_not_aligned);
157  Ef.store(&(flux[nVars-1][i]), is_not_aligned);
158 
159  // store vector 3D only
160  tmpOut[0].store(&(flux[1][i]), is_not_aligned);
161  tmpOut[1].store(&(flux[2][i]), is_not_aligned);
162  tmpOut[2].store(&(flux[3][i]), is_not_aligned);
163  }
164 
165  // spillover loop
166  for (; i < sizeScalar; ++i)
167  {
168  // load scalars
169  NekDouble rhoL = fwd[0][i];
170  NekDouble rhoR = bwd[0][i];
171  NekDouble EL = fwd[spaceDim+1][i];
172  NekDouble ER = bwd[spaceDim+1][i];
173 
174  // 3D case only
175  // load vectors left
176  NekDouble tmpIn[3], tmpOut[3];
177  tmpIn[0] = fwd[1][i];
178  tmpIn[1] = fwd[2][i];
179  tmpIn[2] = fwd[3][i];
180 
181  // load rotation matrix
182  NekDouble rotMat[9];
183  for (size_t j = 0; j < 9; ++j)
184  {
185  rotMat[j] = m_rotMat[j][i];
186  }
187 
188  // rotateTo kernel Fwd
189  rotateToNormalKernel(tmpIn, rotMat, tmpOut);
190 
191  NekDouble rhouL = tmpOut[0];
192  NekDouble rhovL = tmpOut[1];
193  NekDouble rhowL = tmpOut[2];
194 
195  // load vectors right
196  tmpIn[0] = bwd[1][i];
197  tmpIn[1] = bwd[2][i];
198  tmpIn[2] = bwd[3][i];
199 
200  // rotateTo kernel Bwd
201  rotateToNormalKernel(tmpIn, rotMat, tmpOut);
202 
203  NekDouble rhouR = tmpOut[0];
204  NekDouble rhovR = tmpOut[1];
205  NekDouble rhowR = tmpOut[2];
206 
207  // Roe kernel
208  NekDouble rhof{}, Ef{};
209  RoeKernel(
210  rhoL, rhouL, rhovL, rhowL, EL,
211  rhoR, rhouR, rhovR, rhowR, ER,
212  rhof, tmpIn[0], tmpIn[1], tmpIn[2], Ef,
213  gamma);
214 
215  // rotateFrom kernel
216  rotateFromNormalKernel(tmpIn, rotMat, tmpOut);
217 
218  // store scalar
219  flux[0][i] = rhof;
220  flux[nVars-1][i] = Ef;
221 
222  // store vector 3D only
223  flux[1][i] = tmpOut[0];
224  flux[2][i] = tmpOut[1];
225  flux[3][i] = tmpOut[2];
226 
227  }
228 }
229 
230 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:200
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:53
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:1
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:79
tinysimd::simd< NekDouble > vec_t
double NekDouble
static constexpr struct tinysimd::is_not_aligned_t is_not_aligned
typename abi< ScalarType >::type simd
Definition: tinysimd.hpp:83