Nektar++
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 
46 
47 namespace Nektar
48 {
49  template <typename Dim, typename DataType>
50  class Array;
51 
52  namespace SolverUtils
53  {
54  typedef std::function<
56  typedef std::function<
58  typedef std::function<
60 
62  {
63  public:
65  const int nDim,
66  const Array<OneD, const Array<OneD, NekDouble> > &Fwd,
67  const Array<OneD, const Array<OneD, NekDouble> > &Bwd,
69 
70  template<typename FuncPointerT, typename ObjectPointerT>
71  void SetScalar(std::string name,
72  FuncPointerT func,
73  ObjectPointerT obj)
74  {
75  m_scalars[name] = std::bind(func, obj);
76  }
77 
78  void SetScalar(std::string name, RSScalarFuncType fp)
79  {
80  m_scalars[name] = fp;
81  }
82 
83  template<typename FuncPointerT, typename ObjectPointerT>
84  void SetVector(std::string name,
85  FuncPointerT func,
86  ObjectPointerT obj)
87  {
88  m_vectors[name] = std::bind(func, obj);
89  }
90 
91  void SetVector(std::string name, RSVecFuncType fp)
92  {
93  m_vectors[name] = fp;
94  }
95 
96  template<typename FuncPointerT, typename ObjectPointerT>
97  void SetParam(std::string name,
98  FuncPointerT func,
99  ObjectPointerT obj)
100  {
101  m_params[name] = std::bind(func, obj);
102  }
103 
104  void SetParam(std::string name, RSParamFuncType fp)
105  {
106  m_params[name] = fp;
107  }
108 
109  template<typename FuncPointerT, typename ObjectPointerT>
110  void SetAuxScal(std::string name,
111  FuncPointerT func,
112  ObjectPointerT obj)
113  {
114  m_auxScal[name] = std::bind(func, obj);
115  }
116 
117  template<typename FuncPointerT, typename ObjectPointerT>
118  void SetAuxVec(std::string name,
119  FuncPointerT func,
120  ObjectPointerT obj)
121  {
122  m_auxVec[name] = std::bind(func, obj);
123  }
124 
125  void SetAuxVec(std::string name, RSVecFuncType fp)
126  {
127  m_auxVec[name] = fp;
128  }
129 
130  std::map<std::string, RSScalarFuncType> &GetScalars()
131  {
132  return m_scalars;
133  }
134 
135  std::map<std::string, RSVecFuncType> &GetVectors()
136  {
137  return m_vectors;
138  }
139 
140  std::map<std::string, RSParamFuncType> &GetParams()
141  {
142  return m_params;
143  }
144 
146 
148  const int nDim,
149  const Array<OneD, const Array<OneD, NekDouble> > &Fwd,
150  const Array<OneD, const Array<OneD, NekDouble> > &Bwd,
151  DNekBlkMatSharedPtr &FJac,
152  DNekBlkMatSharedPtr &BJac);
153 
154  protected:
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  /// Rotation matrices for each trace quadrature point.
170  /// Rotation storage
172 
175  const LibUtilities::SessionReaderSharedPtr& pSession);
176 
178  {};
179 
180  virtual void v_Solve(
181  const int nDim,
182  const Array<OneD, const Array<OneD, NekDouble> > &Fwd,
183  const Array<OneD, const Array<OneD, NekDouble> > &Bwd,
184  Array<OneD, Array<OneD, NekDouble> > &flux) = 0;
185 
187  const Array<OneD, const Array<OneD, NekDouble> > &normals);
188  void FromToRotation(
191  NekDouble *mat);
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,
210  const Array<OneD, const Array<OneD, NekDouble> > &Fwd,
211  const Array<OneD, const Array<OneD, NekDouble> > &Bwd,
212  const Array<OneD, const Array<OneD, NekDouble> > &normals,
213  DNekBlkMatSharedPtr &FJac,
214  DNekBlkMatSharedPtr &BJac);
215 
216  };
217 
218  /// A shared pointer to an EquationSystem object
219  typedef std::shared_ptr<RiemannSolver> RiemannSolverSharedPtr;
220  /// Datatype of the NekFactory used to instantiate classes derived
221  /// from the RiemannSolver class.
222  typedef LibUtilities::NekFactory<std::string, RiemannSolver,
226 
227 
228  template <class T, typename = typename std::enable_if
229  <
230  std::is_floating_point<T>::value ||
232  >::type
233  >
234  inline void rotateToNormalKernel(T* in, T* rotMat, T* out)
235  {
236 
237  // Apply rotation matrices.
238  out[0] = in[0] * rotMat[0]
239  + in[1] * rotMat[1]
240  + in[2] * rotMat[2];
241 
242  out[1] = in[0] * rotMat[3]
243  + in[1] * rotMat[4]
244  + in[2] * rotMat[5];
245 
246  out[2] = in[0] * rotMat[6]
247  + in[1] * rotMat[7]
248  + in[2] * rotMat[8];
249  }
250 
251  template <class T, typename = typename std::enable_if
252  <
253  std::is_floating_point<T>::value ||
255  >::type
256  >
257  inline void rotateFromNormalKernel(T* in, T* rotMat, T* out)
258  {
259 
260  // Apply rotation matrices.
261  out[0] = in[0] * rotMat[0]
262  + in[1] * rotMat[3]
263  + in[2] * rotMat[6];
264 
265  out[1] = in[0] * rotMat[1]
266  + in[1] * rotMat[4]
267  + in[2] * rotMat[7];
268 
269  out[2] = in[0] * rotMat[2]
270  + in[1] * rotMat[5]
271  + in[2] * rotMat[8];
272  }
273 
274 
275 
276 
277 
278  } // namespace SolverUtils
279 }
280 
281 #endif
#define SOLVER_UTILS_EXPORT
Provides a generic Factory class.
Definition: NekFactory.hpp:105
The RiemannSolver class provides an abstract interface under which solvers for various Riemann proble...
Definition: RiemannSolver.h:62
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)
Definition: RiemannSolver.h:71
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)
Definition: RiemannSolver.h:91
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)
Definition: RiemannSolver.h:97
void SetScalar(std::string name, RSScalarFuncType fp)
Definition: RiemannSolver.h:78
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)
Definition: RiemannSolver.h:84
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
Definition: RiemannSolver.h:55
std::function< const Array< OneD, const Array< OneD, NekDouble > > &()> RSVecFuncType
Definition: RiemannSolver.h:57
std::function< NekDouble()> RSParamFuncType
Definition: RiemannSolver.h:59
RiemannSolverFactory & GetRiemannSolverFactory()
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:71
double NekDouble