35 #include <boost/core/ignore_unused.hpp>
40 #define EPSILON 0.000001
42 #define CROSS(dest, v1, v2) \
44 dest[0] = v1[1] * v2[2] - v1[2] * v2[1]; \
45 dest[1] = v1[2] * v2[0] - v1[0] * v2[2]; \
46 dest[2] = v1[0] * v2[1] - v1[1] * v2[0]; \
49 #define DOT(v1, v2) (v1[0] * v2[0] + v1[1] * v2[1] + v1[2] * v2[2])
51 #define SUB(dest, v1, v2) \
53 dest[0] = v1[0] - v2[0]; \
54 dest[1] = v1[1] - v2[1]; \
55 dest[2] = v1[2] - v2[2]; \
80 RiemannSolver::RiemannSolver() : m_requiresRotation(false), m_rotStorage(3)
86 : m_requiresRotation(false), m_rotStorage(3)
88 boost::ignore_unused(pSession);
121 int nFields = Fwd.size();
122 int nPts = Fwd[0].size();
127 for (
int i = 0; i < 3; ++i)
130 for (
int j = 0; j < nFields; ++j)
173 for (
int i = 0; i < inarray.size(); ++i)
175 Vmath::Vcopy(inarray[i].size(), inarray[i], 1, outarray[i], 1);
178 for (
int i = 0; i < vecLocs.size(); i++)
180 ASSERTL1(vecLocs[i].size() == normals.size(),
181 "vecLocs[i] element count mismatch");
183 switch (normals.size())
187 const int nq = inarray[0].size();
188 const int vx = (int)vecLocs[i][0];
189 Vmath::Vmul(nq, inarray[vx], 1, normals[0], 1, outarray[vx], 1);
194 const int nq = inarray[0].size();
195 const int vx = (int)vecLocs[i][0];
196 const int vy = (int)vecLocs[i][1];
198 Vmath::Vmul(nq, inarray[vx], 1, normals[0], 1, outarray[vx], 1);
199 Vmath::Vvtvp(nq, inarray[vy], 1, normals[1], 1, outarray[vx], 1,
201 Vmath::Vmul(nq, inarray[vx], 1, normals[1], 1, outarray[vy], 1);
202 Vmath::Vvtvm(nq, inarray[vy], 1, normals[0], 1, outarray[vy], 1,
209 const int nq = inarray[0].size();
210 const int vx = (int)vecLocs[i][0];
211 const int vy = (int)vecLocs[i][1];
212 const int vz = (int)vecLocs[i][2];
222 1,
m_rotMat[1], 1, outarray[vx], 1);
226 1,
m_rotMat[4], 1, outarray[vy], 1);
230 1,
m_rotMat[7], 1, outarray[vz], 1);
237 ASSERTL1(
false,
"Invalid space dimension.");
257 for (
int i = 0; i < inarray.size(); ++i)
259 Vmath::Vcopy(inarray[i].size(), inarray[i], 1, outarray[i], 1);
262 for (
int i = 0; i < vecLocs.size(); i++)
264 ASSERTL1(vecLocs[i].size() == normals.size(),
265 "vecLocs[i] element count mismatch");
267 switch (normals.size())
271 const int nq = normals[0].size();
272 const int vx = (int)vecLocs[i][0];
273 Vmath::Vmul(nq, inarray[vx], 1, normals[0], 1, outarray[vx], 1);
278 const int nq = normals[0].size();
279 const int vx = (int)vecLocs[i][0];
280 const int vy = (int)vecLocs[i][1];
282 Vmath::Vmul(nq, inarray[vy], 1, normals[1], 1, outarray[vx], 1);
283 Vmath::Vvtvm(nq, inarray[vx], 1, normals[0], 1, outarray[vx], 1,
285 Vmath::Vmul(nq, inarray[vx], 1, normals[1], 1, outarray[vy], 1);
286 Vmath::Vvtvp(nq, inarray[vy], 1, normals[0], 1, outarray[vy], 1,
293 const int nq = normals[0].size();
294 const int vx = (int)vecLocs[i][0];
295 const int vy = (int)vecLocs[i][1];
296 const int vz = (int)vecLocs[i][2];
299 1,
m_rotMat[3], 1, outarray[vx], 1);
303 1,
m_rotMat[4], 1, outarray[vy], 1);
307 1,
m_rotMat[5], 1, outarray[vz], 1);
314 ASSERTL1(
false,
"Invalid space dimension.");
379 const int nq = normals[0].size();
386 for (i = 0; i < 9; ++i)
390 for (i = 0; i < normals[0].size(); ++i)
394 tn[0] = normals[0][i];
395 tn[1] = normals[1][i];
396 tn[2] = normals[2][i];
399 for (j = 0; j < 9; ++j)
427 f = (e < 0) ? -e : e;
435 x[0] = (from[0] > 0.0) ? from[0] : -from[0];
436 x[1] = (from[1] > 0.0) ? from[1] : -from[1];
437 x[2] = (from[2] > 0.0) ? from[2] : -from[2];
466 u[0] = x[0] - from[0];
467 u[1] = x[1] - from[1];
468 u[2] = x[2] - from[2];
473 c1 = 2.0 /
DOT(u, u);
474 c2 = 2.0 /
DOT(v, v);
475 c3 = c1 * c2 *
DOT(u, v);
477 for (i = 0; i < 3; i++)
479 for (j = 0; j < 3; j++)
482 -c1 * u[i] * u[j] - c2 * v[i] * v[j] + c3 * v[i] * u[j];
484 mat[i + 3 * i] += 1.0;
496 mat[0] = e + hvx * v[0];
497 mat[1] = hvxy - v[2];
498 mat[2] = hvxz + v[1];
499 mat[3] = hvxy + v[2];
500 mat[4] = e + h * v[1] * v[1];
501 mat[5] = hvyz - v[0];
502 mat[6] = hvxz - v[1];
503 mat[7] = hvyz + v[0];
504 mat[8] = e + hvz * v[2];
520 int nPts = Fwd[0].size();
536 for (
int i = 0; i < nDim; i++)
552 boost::ignore_unused(nDim, Fwd, Bwd, normals, FJac, BJac);
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mode...
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
#define CROSS(dest, v1, v2)
Provides a generic Factory class.
SOLVER_UTILS_EXPORT bool CheckAuxVec(std::string name)
Determine whether a vector has been defined in m_auxVec.
SOLVER_UTILS_EXPORT void rotateToNormal(const Array< OneD, const Array< OneD, NekDouble >> &inarray, const Array< OneD, const Array< OneD, NekDouble >> &normals, const Array< OneD, const Array< OneD, NekDouble >> &vecLocs, Array< OneD, Array< OneD, NekDouble >> &outarray)
Rotate a vector field to trace normal.
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_rotStorage
Rotation storage.
bool m_requiresRotation
Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields.
SOLVER_UTILS_EXPORT void rotateFromNormal(const Array< OneD, const Array< OneD, NekDouble >> &inarray, const Array< OneD, const Array< OneD, NekDouble >> &normals, const Array< OneD, const Array< OneD, NekDouble >> &vecLocs, Array< OneD, Array< OneD, NekDouble >> &outarray)
Rotate a vector field from trace normal.
SOLVER_UTILS_EXPORT void GenerateRotationMatrices(const Array< OneD, const Array< OneD, NekDouble >> &normals)
Generate rotation matrices for 3D expansions.
SOLVER_UTILS_EXPORT bool CheckParams(std::string name)
Determine whether a parameter has been defined in m_params.
std::map< std::string, RSScalarFuncType > m_auxScal
Map of auxiliary scalar function types.
SOLVER_UTILS_EXPORT void CalcFluxJacobian(const int nDim, const Array< OneD, const Array< OneD, NekDouble >> &Fwd, const Array< OneD, const Array< OneD, NekDouble >> &Bwd, DNekBlkMatSharedPtr &FJac, DNekBlkMatSharedPtr &BJac)
Calculate the flux jacobian of Fwd and Bwd.
SOLVER_UTILS_EXPORT void Solve(const int nDim, const Array< OneD, const Array< OneD, NekDouble >> &Fwd, const Array< OneD, const Array< OneD, NekDouble >> &Bwd, Array< OneD, Array< OneD, NekDouble >> &flux)
Perform the Riemann solve given the forwards and backwards spaces.
virtual SOLVER_UTILS_EXPORT void v_CalcFluxJacobian(const int nDim, const Array< OneD, const Array< OneD, NekDouble >> &Fwd, const Array< OneD, const Array< OneD, NekDouble >> &Bwd, const Array< OneD, const Array< OneD, NekDouble >> &normals, DNekBlkMatSharedPtr &FJac, DNekBlkMatSharedPtr &BJac)
Array< OneD, Array< OneD, NekDouble > > m_rotMat
Rotation matrices for each trace quadrature point.
SOLVER_UTILS_EXPORT bool CheckScalars(std::string name)
Determine whether a scalar has been defined in m_scalars.
std::map< std::string, RSVecFuncType > m_vectors
Map of vector function types.
std::map< std::string, RSScalarFuncType > m_scalars
Map of scalar function types.
SOLVER_UTILS_EXPORT bool CheckVectors(std::string name)
Determine whether a vector has been defined in m_vectors.
void FromToRotation(Array< OneD, const NekDouble > &from, Array< OneD, const NekDouble > &to, NekDouble *mat)
A function for creating a rotation matrix that rotates a vector from into another vector to.
virtual void v_Solve(const int nDim, const Array< OneD, const Array< OneD, NekDouble >> &Fwd, const Array< OneD, const Array< OneD, NekDouble >> &Bwd, Array< OneD, Array< OneD, NekDouble >> &flux)=0
SOLVER_UTILS_EXPORT bool CheckAuxScal(std::string name)
Determine whether a scalar has been defined in m_auxScal.
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
SOLVER_UTILS_EXPORT RiemannSolver()
std::map< std::string, RSVecFuncType > m_auxVec
Map of auxiliary vector function types.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
RiemannSolverFactory & GetRiemannSolverFactory()
The above copyright notice and this permission notice shall be included.
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
void Vvtvm(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvm (vector times vector minus vector): z = w*x - y
void Fill(int n, const T alpha, T *x, const int incx)
Fill a vector with a constant value.
void Vvtvvtp(int n, const T *v, int incv, const T *w, int incw, const T *x, int incx, const T *y, int incy, T *z, int incz)
vvtvvtp (vector times vector plus vector times vector):
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)