35 #ifndef NEKTAR_SOLVERUTILS_RIEMANNSOLVER
36 #define NEKTAR_SOLVERUTILS_RIEMANNSOLVER
49 template <
typename Dim,
typename DataType>
54 typedef std::function<
56 typedef std::function<
58 typedef std::function<
70 template<
typename FuncPo
interT,
typename ObjectPo
interT>
83 template<
typename FuncPo
interT,
typename ObjectPo
interT>
96 template<
typename FuncPo
interT,
typename ObjectPo
interT>
109 template<
typename FuncPo
interT,
typename ObjectPo
interT>
117 template<
typename FuncPo
interT,
typename ObjectPo
interT>
228 template <
class T,
typename =
typename std::enable_if
230 std::is_floating_point<T>::value ||
238 out[0] = in[0] * rotMat[0]
242 out[1] = in[0] * rotMat[3]
246 out[2] = in[0] * rotMat[6]
251 template <
class T,
typename =
typename std::enable_if
253 std::is_floating_point<T>::value ||
261 out[0] = in[0] * rotMat[0]
265 out[1] = in[0] * rotMat[1]
269 out[2] = in[0] * rotMat[2]
#define SOLVER_UTILS_EXPORT
Provides a generic Factory class.
The RiemannSolver class provides an abstract interface under which solvers for various Riemann proble...
SOLVER_UTILS_EXPORT bool CheckAuxVec(std::string name)
Determine whether a vector has been defined in m_auxVec.
void SetScalar(std::string name, FuncPointerT func, ObjectPointerT obj)
virtual SOLVER_UTILS_EXPORT ~RiemannSolver()
std::map< std::string, RSScalarFuncType > & GetScalars()
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.
void SetVector(std::string name, RSVecFuncType fp)
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.
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)
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.
void SetParam(std::string name, RSParamFuncType fp)
void SetParam(std::string name, FuncPointerT func, ObjectPointerT obj)
void SetScalar(std::string name, RSScalarFuncType fp)
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 SetAuxVec(std::string name, FuncPointerT func, ObjectPointerT obj)
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 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.
void SetVector(std::string name, FuncPointerT func, ObjectPointerT obj)
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.
void SetAuxScal(std::string name, FuncPointerT func, ObjectPointerT obj)
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.
std::map< std::string, RSVecFuncType > m_vectors
Map of vector function types.
std::map< std::string, RSParamFuncType > & GetParams()
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.
SOLVER_UTILS_EXPORT void GenerateRotationMatrices(const Array< OneD, const Array< OneD, NekDouble > > &normals)
Generate rotation matrices for 3D expansions.
SOLVER_UTILS_EXPORT bool CheckAuxScal(std::string name)
Determine whether a scalar has been defined in m_auxScal.
SOLVER_UTILS_EXPORT RiemannSolver()
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
void SetAuxVec(std::string name, RSVecFuncType fp)
std::map< std::string, RSVecFuncType > m_auxVec
Map of auxiliary vector function types.
std::map< std::string, RSVecFuncType > & GetVectors()
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.
void rotateFromNormalKernel(T *in, T *rotMat, T *out)
void rotateToNormalKernel(T *in, T *rotMat, T *out)
LibUtilities::NekFactory< std::string, RiemannSolver, const LibUtilities::SessionReaderSharedPtr & > RiemannSolverFactory
Datatype of the NekFactory used to instantiate classes derived from the RiemannSolver class.
std::function< const Array< OneD, const NekDouble > &()> RSScalarFuncType
std::function< const Array< OneD, const Array< OneD, NekDouble > > &()> RSVecFuncType
std::function< NekDouble()> RSParamFuncType
RiemannSolverFactory & GetRiemannSolverFactory()
The above copyright notice and this permission notice shall be included.
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr