63 [[maybe_unused]]
const int nDim,
68 static auto gamma =
m_params[
"gamma"]();
69 static size_t nVars = fwd.size();
70 static size_t spaceDim = nVars - 2;
73 ASSERTL0(spaceDim == 3,
"SIMD Roe implemented only for 3D case...");
80 size_t sizeScalar = fwd[0].size();
81 size_t sizeVec = (sizeScalar / vec_t::width) * vec_t::width;
95 for (; i < sizeVec; i += vec_t::width)
98 vec_t rhoL, rhoR, ER, EL;
105 vec_t tmpIn[3], tmpOut[3];
112 for (
size_t j = 0; j < 9; ++j)
120 vec_t rhouL = tmpOut[0];
121 vec_t rhovL = tmpOut[1];
122 vec_t rhowL = tmpOut[2];
132 vec_t rhouR = tmpOut[0];
133 vec_t rhovR = tmpOut[1];
134 vec_t rhowR = tmpOut[2];
138 RoeKernel(rhoL, rhouL, rhovL, rhowL, EL, rhoR, rhouR, rhovR, rhowR, ER,
139 rhof, tmpIn[0], tmpIn[1], tmpIn[2], Ef, gamma);
155 for (; i < sizeScalar; ++i)
166 tmpIn[0] = fwd[1][i];
167 tmpIn[1] = fwd[2][i];
168 tmpIn[2] = fwd[3][i];
172 for (
size_t j = 0; j < 9; ++j)
185 tmpIn[0] = bwd[1][i];
186 tmpIn[1] = bwd[2][i];
187 tmpIn[2] = bwd[3][i];
198 RoeKernel(rhoL, rhouL, rhovL, rhowL, EL, rhoR, rhouR, rhovR, rhowR, ER,
199 rhof, tmpIn[0], tmpIn[1], tmpIn[2], Ef, gamma);
206 flux[nVars - 1][i] = Ef;
209 flux[1][i] = tmpOut[0];
210 flux[2][i] = tmpOut[1];
211 flux[3][i] = tmpOut[2];
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
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)
static std::string solverName
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)
tinysimd::simd< NekDouble > vec_t
static constexpr struct tinysimd::is_not_aligned_t is_not_aligned
typename abi< ScalarType, width >::type simd