35 #include <boost/core/ignore_unused.hpp> 40 #define EPSILON 0.000001 42 #define CROSS(dest, v1, v2){ \ 43 dest[0] = v1[1] * v2[2] - v1[2] * v2[1]; \ 44 dest[1] = v1[2] * v2[0] - v1[0] * v2[2]; \ 45 dest[2] = v1[0] * v2[1] - v1[1] * v2[0];} 47 #define DOT(v1, v2) (v1[0] * v2[0] + v1[1] * v2[1] + v1[2] * v2[2]) 49 #define SUB(dest, v1, v2){ \ 50 dest[0] = v1[0] - v2[0]; \ 51 dest[1] = v1[1] - v2[1]; \ 52 dest[2] = v1[2] - v2[2];} 76 RiemannSolver::RiemannSolver(
78 : m_requiresRotation(false), m_rotStorage (3)
80 boost::ignore_unused(pSession);
114 int nFields = Fwd .num_elements();
115 int nPts = Fwd[0].num_elements();
120 for (
int i = 0; i < 3; ++i)
124 for (
int j = 0; j < nFields; ++j)
168 for (
int i = 0; i < inarray.num_elements(); ++i)
174 for (
int i = 0; i < vecLocs.num_elements(); i++)
176 ASSERTL1(vecLocs[i].num_elements() == normals.num_elements(),
177 "vecLocs[i] element count mismatch");
179 switch (normals.num_elements())
183 const int nq = inarray[0].num_elements();
184 const int vx = (int)vecLocs[i][0];
191 const int nq = inarray[0].num_elements();
192 const int vx = (int)vecLocs[i][0];
193 const int vy = (int)vecLocs[i][1];
198 outarray[vx], 1, outarray[vx], 1);
202 outarray[vy], 1, outarray[vy], 1);
208 const int nq = inarray[0].num_elements();
209 const int vx = (int)vecLocs[i][0];
210 const int vy = (int)vecLocs[i][1];
211 const int vz = (int)vecLocs[i][2];
224 outarray[vx], 1, outarray[vx], 1);
229 outarray[vy], 1, outarray[vy], 1);
234 outarray[vz], 1, outarray[vz], 1);
239 ASSERTL1(
false,
"Invalid space dimension.");
259 for (
int i = 0; i < inarray.num_elements(); ++i)
265 for (
int i = 0; i < vecLocs.num_elements(); i++)
267 ASSERTL1(vecLocs[i].num_elements() == normals.num_elements(),
268 "vecLocs[i] element count mismatch");
270 switch (normals.num_elements())
274 const int nq = normals[0].num_elements();
275 const int vx = (int)vecLocs[i][0];
282 const int nq = normals[0].num_elements();
283 const int vx = (int)vecLocs[i][0];
284 const int vy = (int)vecLocs[i][1];
289 outarray[vx], 1, outarray[vx], 1);
293 outarray[vy], 1, outarray[vy], 1);
299 const int nq = normals[0].num_elements();
300 const int vx = (int)vecLocs[i][0];
301 const int vy = (int)vecLocs[i][1];
302 const int vz = (int)vecLocs[i][2];
308 outarray[vx], 1, outarray[vx], 1);
313 outarray[vy], 1, outarray[vy], 1);
318 outarray[vz], 1, outarray[vz], 1);
323 ASSERTL1(
false,
"Invalid space dimension.");
388 const int nq = normals[0].num_elements();
395 for (i = 0; i < 9; ++i)
399 for (i = 0; i < normals[0].num_elements(); ++i)
403 tn[0] = normals[0][i];
404 tn[1] = normals[1][i];
405 tn[2] = normals[2][i];
408 for (j = 0; j < 9; ++j)
445 x[0] = (from[0] > 0.0)? from[0] : -from[0];
446 x[1] = (from[1] > 0.0)? from[1] : -from[1];
447 x[2] = (from[2] > 0.0)? from[2] : -from[2];
453 x[0] = 1.0; x[1] = x[2] = 0.0;
457 x[2] = 1.0; x[0] = x[1] = 0.0;
464 x[1] = 1.0; x[0] = x[2] = 0.0;
468 x[2] = 1.0; x[0] = x[1] = 0.0;
472 u[0] = x[0] - from[0];
473 u[1] = x[1] - from[1];
474 u[2] = x[2] - from[2];
475 v[0] = x[0] - to [0];
476 v[1] = x[1] - to [1];
477 v[2] = x[2] - to [2];
479 c1 = 2.0 /
DOT(u, u);
480 c2 = 2.0 /
DOT(v, v);
481 c3 = c1 * c2 *
DOT(u, v);
483 for (i = 0; i < 3; i++) {
484 for (j = 0; j < 3; j++) {
485 mat[3*i+j] = - c1 * u[i] * u[j]
501 mat[0] = e + hvx * v[0];
502 mat[1] = hvxy - v[2];
503 mat[2] = hvxz + v[1];
504 mat[3] = hvxy + v[2];
505 mat[4] = e + h * v[1] * v[1];
506 mat[5] = hvyz - v[0];
507 mat[6] = hvxz - v[1];
508 mat[7] = hvyz + v[0];
509 mat[8] = e + hvz * v[2];
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...
SOLVER_UTILS_EXPORT bool CheckScalars(std::string name)
Determine whether a scalar has been defined in m_scalars.
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
SOLVER_UTILS_EXPORT bool CheckVectors(std::string name)
Determine whether a vector has been defined in m_vectors.
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.
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
#define CROSS(dest, v1, v2)
std::map< std::string, RSScalarFuncType > m_scalars
Map of scalar function types.
void GenerateRotationMatrices(const Array< OneD, const Array< OneD, NekDouble > > &normals)
Generate rotation matrices for 3D expansions.
RiemannSolverFactory & GetRiemannSolverFactory()
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.
std::map< std::string, RSScalarFuncType > m_auxScal
Map of auxiliary scalar function types.
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_rotStorage
Rotation storage.
SOLVER_UTILS_EXPORT bool CheckAuxScal(std::string name)
Determine whether a scalar has been defined in m_auxScal.
Array< OneD, Array< OneD, NekDouble > > m_rotMat
Rotation matrices for each trace quadrature point.
SOLVER_UTILS_EXPORT bool CheckAuxVec(std::string name)
Determine whether a vector has been defined in m_auxVec.
std::map< std::string, RSVecFuncType > m_auxVec
Map of auxiliary vector function types.
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):
bool m_requiresRotation
Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields...
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 plus vector): z = w*x - y
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
SOLVER_UTILS_EXPORT bool CheckParams(std::string name)
Determine whether a parameter has been defined in m_params.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
std::map< std::string, RSVecFuncType > m_vectors
Map of vector function types.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
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.
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.
Provides a generic Factory class.