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 
42 
43 #include <string>
44 
45 namespace Nektar
46 {
47  template <typename Dim, typename DataType>
48  class Array;
49 
50  namespace SolverUtils
51  {
52  typedef std::function<
54  typedef std::function<
56  typedef std::function<
58 
60  {
61  public:
63  const int nDim,
64  const Array<OneD, const Array<OneD, NekDouble> > &Fwd,
65  const Array<OneD, const Array<OneD, NekDouble> > &Bwd,
67 
68  template<typename FuncPointerT, typename ObjectPointerT>
69  void SetScalar(std::string name,
70  FuncPointerT func,
71  ObjectPointerT obj)
72  {
73  m_scalars[name] = std::bind(func, obj);
74  }
75 
76  void SetScalar(std::string name, RSScalarFuncType fp)
77  {
78  m_scalars[name] = fp;
79  }
80 
81  template<typename FuncPointerT, typename ObjectPointerT>
82  void SetVector(std::string name,
83  FuncPointerT func,
84  ObjectPointerT obj)
85  {
86  m_vectors[name] = std::bind(func, obj);
87  }
88 
89  void SetVector(std::string name, RSVecFuncType fp)
90  {
91  m_vectors[name] = fp;
92  }
93 
94  template<typename FuncPointerT, typename ObjectPointerT>
95  void SetParam(std::string name,
96  FuncPointerT func,
97  ObjectPointerT obj)
98  {
99  m_params[name] = std::bind(func, obj);
100  }
101 
102  void SetParam(std::string name, RSParamFuncType fp)
103  {
104  m_params[name] = fp;
105  }
106 
107  template<typename FuncPointerT, typename ObjectPointerT>
108  void SetAuxScal(std::string name,
109  FuncPointerT func,
110  ObjectPointerT obj)
111  {
112  m_auxScal[name] = std::bind(func, obj);
113  }
114 
115  template<typename FuncPointerT, typename ObjectPointerT>
116  void SetAuxVec(std::string name,
117  FuncPointerT func,
118  ObjectPointerT obj)
119  {
120  m_auxVec[name] = std::bind(func, obj);
121  }
122 
123  std::map<std::string, RSScalarFuncType> &GetScalars()
124  {
125  return m_scalars;
126  }
127 
128  std::map<std::string, RSVecFuncType> &GetVectors()
129  {
130  return m_vectors;
131  }
132 
133  std::map<std::string, RSParamFuncType> &GetParams()
134  {
135  return m_params;
136  }
137 
139 
140  protected:
141  /// Indicates whether the Riemann solver requires a rotation to be
142  /// applied to the velocity fields.
144  /// Map of scalar function types.
145  std::map<std::string, RSScalarFuncType> m_scalars;
146  /// Map of vector function types.
147  std::map<std::string, RSVecFuncType> m_vectors;
148  /// Map of parameter function types.
149  std::map<std::string, RSParamFuncType > m_params;
150  /// Map of auxiliary scalar function types.
151  std::map<std::string, RSScalarFuncType> m_auxScal;
152  /// Map of auxiliary vector function types.
153  std::map<std::string, RSVecFuncType> m_auxVec;
154  /// Rotation matrices for each trace quadrature point.
156  /// Rotation storage
158 
160  const LibUtilities::SessionReaderSharedPtr& pSession);
161 
163  {};
164 
165  virtual void v_Solve(
166  const int nDim,
167  const Array<OneD, const Array<OneD, NekDouble> > &Fwd,
168  const Array<OneD, const Array<OneD, NekDouble> > &Bwd,
169  Array<OneD, Array<OneD, NekDouble> > &flux) = 0;
170 
172  const Array<OneD, const Array<OneD, NekDouble> > &normals);
173  void FromToRotation(
176  NekDouble *mat);
178  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
179  const Array<OneD, const Array<OneD, NekDouble> > &normals,
180  const Array<OneD, const Array<OneD, NekDouble> > &vecLocs,
181  Array<OneD, Array<OneD, NekDouble> > &outarray);
183  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
184  const Array<OneD, const Array<OneD, NekDouble> > &normals,
185  const Array<OneD, const Array<OneD, NekDouble> > &vecLocs,
186  Array<OneD, Array<OneD, NekDouble> > &outarray);
187  SOLVER_UTILS_EXPORT bool CheckScalars (std::string name);
188  SOLVER_UTILS_EXPORT bool CheckVectors (std::string name);
189  SOLVER_UTILS_EXPORT bool CheckParams (std::string name);
190  SOLVER_UTILS_EXPORT bool CheckAuxScal (std::string name);
191  SOLVER_UTILS_EXPORT bool CheckAuxVec (std::string name);
192  };
193 
194  /// A shared pointer to an EquationSystem object
195  typedef std::shared_ptr<RiemannSolver> RiemannSolverSharedPtr;
196  /// Datatype of the NekFactory used to instantiate classes derived
197  /// from the RiemannSolver class.
198  typedef LibUtilities::NekFactory<std::string, RiemannSolver,
202  }
203 }
204 
205 #endif
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 SOLVER_UTILS_EXPORT ~RiemannSolver()
SOLVER_UTILS_EXPORT bool CheckScalars(std::string name)
Determine whether a scalar has been defined in m_scalars.
void SetParam(std::string name, FuncPointerT func, ObjectPointerT obj)
Definition: RiemannSolver.h:95
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.
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
std::map< std::string, RSParamFuncType > & GetParams()
void SetVector(std::string name, FuncPointerT func, ObjectPointerT obj)
Definition: RiemannSolver.h:82
std::map< std::string, RSScalarFuncType > m_scalars
Map of scalar function types.
std::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.
std::function< const Array< OneD, const Array< OneD, NekDouble > > &()> RSVecFuncType
Definition: RiemannSolver.h:55
void GenerateRotationMatrices(const Array< OneD, const Array< OneD, NekDouble > > &normals)
Generate rotation matrices for 3D expansions.
void SetParam(std::string name, RSParamFuncType fp)
void SetScalar(std::string name, RSScalarFuncType fp)
Definition: RiemannSolver.h:76
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.
SOLVER_UTILS_EXPORT RiemannSolver(const LibUtilities::SessionReaderSharedPtr &pSession)
void SetAuxScal(std::string name, FuncPointerT func, ObjectPointerT obj)
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_rotStorage
Rotation storage.
double NekDouble
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::function< const Array< OneD, const NekDouble > &()> RSScalarFuncType
Definition: RiemannSolver.h:53
std::map< std::string, RSVecFuncType > m_auxVec
Map of auxiliary vector function types.
LibUtilities::NekFactory< std::string, RiemannSolver, const LibUtilities::SessionReaderSharedPtr & > RiemannSolverFactory
Datatype of the NekFactory used to instantiate classes derived from the RiemannSolver class...
bool m_requiresRotation
Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields...
The RiemannSolver class provides an abstract interface under which solvers for various Riemann proble...
Definition: RiemannSolver.h:59
#define SOLVER_UTILS_EXPORT
std::function< NekDouble()> RSParamFuncType
Definition: RiemannSolver.h:57
void SetScalar(std::string name, FuncPointerT func, ObjectPointerT obj)
Definition: RiemannSolver.h:69
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.
void SetVector(std::string name, RSVecFuncType fp)
Definition: RiemannSolver.h:89
std::map< std::string, RSVecFuncType > & GetVectors()
std::map< std::string, RSScalarFuncType > & GetScalars()
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.
Provides a generic Factory class.
Definition: NekFactory.hpp:103