Nektar++
Loading...
Searching...
No Matches
RiemannSolver.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: RiemannSolver.h
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: Abstract base class for Riemann solvers with factory.
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#ifndef NEKTAR_SOLVERUTILS_RIEMANNSOLVER
36#define NEKTAR_SOLVERUTILS_RIEMANNSOLVER
37
43
44#include <string>
45
46namespace Nektar
47{
48template <typename Dim, typename DataType> class Array;
49
50namespace SolverUtils
51{
52typedef std::function<const Array<OneD, const NekDouble> &()> RSScalarFuncType;
53typedef std::function<const Array<OneD, const Array<OneD, NekDouble>> &()>
55typedef std::function<NekDouble()> RSParamFuncType;
56typedef std::function<bool()> RSFlagype;
57
59{
60public:
62 const int nDim, const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
63 const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
65
66 template <typename FuncPointerT, typename ObjectPointerT>
67 void SetScalar(std::string name, FuncPointerT func, ObjectPointerT obj)
68 {
69 m_scalars[name] = std::bind(func, obj);
70 }
71
72 void SetScalar(std::string name, RSScalarFuncType fp)
73 {
74 m_scalars[name] = fp;
75 }
76
77 template <typename FuncPointerT, typename ObjectPointerT>
78 void SetVector(std::string name, FuncPointerT func, ObjectPointerT obj)
79 {
80 m_vectors[name] = std::bind(func, obj);
81 }
82
83 void SetVector(std::string name, RSVecFuncType fp)
84 {
85 m_vectors[name] = fp;
86 }
87
88 template <typename FuncPointerT, typename ObjectPointerT>
89 void SetParam(std::string name, FuncPointerT func, ObjectPointerT obj)
90 {
91 m_params[name] = std::bind(func, obj);
92 }
93
94 template <typename FuncPointerT, typename ObjectPointerT>
95 void SetUpdateNormalsFlag(FuncPointerT func, ObjectPointerT obj)
96 {
97 m_updateNormals = std::bind(func, obj);
98 }
99
100 void SetALEFlag(bool &ALE)
101 {
102 m_ALESolver = ALE;
103 }
104
105 void SetParam(std::string name, RSParamFuncType fp)
106 {
107 m_params[name] = fp;
108 }
109
110 template <typename FuncPointerT, typename ObjectPointerT>
111 void SetAuxScal(std::string name, FuncPointerT func, ObjectPointerT obj)
112 {
113 m_auxScal[name] = std::bind(func, obj);
114 }
115
116 template <typename FuncPointerT, typename ObjectPointerT>
117 void SetAuxVec(std::string name, FuncPointerT func, ObjectPointerT obj)
118 {
119 m_auxVec[name] = std::bind(func, obj);
120 }
121
122 void SetAuxVec(std::string name, RSVecFuncType fp)
123 {
124 m_auxVec[name] = fp;
125 }
126
127 std::map<std::string, RSScalarFuncType> &GetScalars()
128 {
129 return m_scalars;
130 }
131
132 std::map<std::string, RSVecFuncType> &GetVectors()
133 {
134 return m_vectors;
135 }
136
137 std::map<std::string, RSParamFuncType> &GetParams()
138 {
139 return m_params;
140 }
141
143
145 const int nDim, const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
146 const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
148
153
154protected:
155 /// Indicates whether the Riemann solver requires a rotation to be
156 /// applied to the velocity fields.
158 /// Map of scalar function types.
159 std::map<std::string, RSScalarFuncType> m_scalars;
160 /// Map of vector function types.
161 std::map<std::string, RSVecFuncType> m_vectors;
162 /// Map of parameter function types.
163 std::map<std::string, RSParamFuncType> m_params;
164 /// Map of auxiliary scalar function types.
165 std::map<std::string, RSScalarFuncType> m_auxScal;
166 /// Map of auxiliary vector function types.
167 std::map<std::string, RSVecFuncType> m_auxVec;
168 /// Function of normals updated flag
169 std::function<bool()> m_updateNormals;
170 /// Rotation matrices for each trace quadrature point.
172 /// Rotation storage
174 /// Flag if using the ALE formulation
175 bool m_ALESolver = false;
176
180
182
183 virtual void v_Solve(const int nDim,
184 const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
185 const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
186 Array<OneD, Array<OneD, NekDouble>> &flux) = 0;
187
189 const Array<OneD, const Array<OneD, NekDouble>> &normals);
193 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
194 const Array<OneD, const Array<OneD, NekDouble>> &normals,
195 const Array<OneD, const Array<OneD, NekDouble>> &vecLocs,
196 Array<OneD, Array<OneD, NekDouble>> &outarray);
198 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
199 const Array<OneD, const Array<OneD, NekDouble>> &normals,
200 const Array<OneD, const Array<OneD, NekDouble>> &vecLocs,
201 Array<OneD, Array<OneD, NekDouble>> &outarray);
202 SOLVER_UTILS_EXPORT bool CheckScalars(std::string name);
203 SOLVER_UTILS_EXPORT bool CheckVectors(std::string name);
204 SOLVER_UTILS_EXPORT bool CheckParams(std::string name);
205 SOLVER_UTILS_EXPORT bool CheckAuxScal(std::string name);
206 SOLVER_UTILS_EXPORT bool CheckAuxVec(std::string name);
207
209 const int nDim, const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
210 const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
211 const Array<OneD, const Array<OneD, NekDouble>> &normals,
213};
214
215/// A shared pointer to an EquationSystem object
216typedef std::shared_ptr<RiemannSolver> RiemannSolverSharedPtr;
217/// Datatype of the NekFactory used to instantiate classes derived
218/// from the RiemannSolver class.
219typedef LibUtilities::NekFactory<std::string, RiemannSolver,
223
224template <class T, typename = typename std::enable_if<
225 std::is_floating_point_v<T> ||
226 tinysimd::is_vector_floating_point_v<T>>::type>
227inline void rotateToNormalKernel(T *in, T *rotMat, T *out)
228{
229
230 // Apply rotation matrices.
231 out[0] = in[0] * rotMat[0] + in[1] * rotMat[1] + in[2] * rotMat[2];
232
233 out[1] = in[0] * rotMat[3] + in[1] * rotMat[4] + in[2] * rotMat[5];
234
235 out[2] = in[0] * rotMat[6] + in[1] * rotMat[7] + in[2] * rotMat[8];
236}
237
238template <class T, typename = typename std::enable_if<
239 std::is_floating_point_v<T> ||
240 tinysimd::is_vector_floating_point_v<T>>::type>
241inline void rotateFromNormalKernel(T *in, T *rotMat, T *out)
242{
243
244 // Apply rotation matrices.
245 out[0] = in[0] * rotMat[0] + in[1] * rotMat[3] + in[2] * rotMat[6];
246
247 out[1] = in[0] * rotMat[1] + in[1] * rotMat[4] + in[2] * rotMat[7];
248
249 out[2] = in[0] * rotMat[2] + in[1] * rotMat[5] + in[2] * rotMat[8];
250}
251
252} // namespace SolverUtils
253} // namespace Nektar
254
255#endif
#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()
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::function< bool()> m_updateNormals
Function of normals updated flag.
std::map< std::string, RSScalarFuncType > & GetScalars()
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.
SOLVER_UTILS_EXPORT void ResetRotMatrix(void)
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)
std::map< std::string, RSParamFuncType > & GetParams()
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, 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.
std::map< std::string, RSVecFuncType > & GetVectors()
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.
void SetUpdateNormalsFlag(FuncPointerT func, ObjectPointerT obj)
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
SOLVER_UTILS_EXPORT RiemannSolver()
bool m_ALESolver
Flag if using the ALE formulation.
void SetAuxVec(std::string name, RSVecFuncType fp)
std::map< std::string, RSVecFuncType > m_auxVec
Map of auxiliary vector function types.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.
std::function< bool()> RSFlagype
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()
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
static Array< OneD, Array< OneD, NekDouble > > NullNekDoubleArrayOfArray