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
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;
56
58{
59public:
61 const int nDim, const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
62 const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
64
65 template <typename FuncPointerT, typename ObjectPointerT>
66 void SetScalar(std::string name, FuncPointerT func, ObjectPointerT obj)
67 {
68 m_scalars[name] = std::bind(func, obj);
69 }
70
71 void SetScalar(std::string name, RSScalarFuncType fp)
72 {
73 m_scalars[name] = fp;
74 }
75
76 template <typename FuncPointerT, typename ObjectPointerT>
77 void SetVector(std::string name, FuncPointerT func, ObjectPointerT obj)
78 {
79 m_vectors[name] = std::bind(func, obj);
80 }
81
82 void SetVector(std::string name, RSVecFuncType fp)
83 {
84 m_vectors[name] = fp;
85 }
86
87 template <typename FuncPointerT, typename ObjectPointerT>
88 void SetParam(std::string name, FuncPointerT func, ObjectPointerT obj)
89 {
90 m_params[name] = std::bind(func, obj);
91 }
92
93 void SetParam(std::string name, RSParamFuncType fp)
94 {
95 m_params[name] = fp;
96 }
97
98 template <typename FuncPointerT, typename ObjectPointerT>
99 void SetAuxScal(std::string name, FuncPointerT func, ObjectPointerT obj)
100 {
101 m_auxScal[name] = std::bind(func, obj);
102 }
103
104 template <typename FuncPointerT, typename ObjectPointerT>
105 void SetAuxVec(std::string name, FuncPointerT func, ObjectPointerT obj)
106 {
107 m_auxVec[name] = std::bind(func, obj);
108 }
109
110 void SetAuxVec(std::string name, RSVecFuncType fp)
111 {
112 m_auxVec[name] = fp;
113 }
114
115 std::map<std::string, RSScalarFuncType> &GetScalars()
116 {
117 return m_scalars;
118 }
119
120 std::map<std::string, RSVecFuncType> &GetVectors()
121 {
122 return m_vectors;
123 }
124
125 std::map<std::string, RSParamFuncType> &GetParams()
126 {
127 return m_params;
128 }
129
131
133 const int nDim, const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
134 const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
136
137protected:
138 /// Indicates whether the Riemann solver requires a rotation to be
139 /// applied to the velocity fields.
141 /// Map of scalar function types.
142 std::map<std::string, RSScalarFuncType> m_scalars;
143 /// Map of vector function types.
144 std::map<std::string, RSVecFuncType> m_vectors;
145 /// Map of parameter function types.
146 std::map<std::string, RSParamFuncType> m_params;
147 /// Map of auxiliary scalar function types.
148 std::map<std::string, RSScalarFuncType> m_auxScal;
149 /// Map of auxiliary vector function types.
150 std::map<std::string, RSVecFuncType> m_auxVec;
151 /// Rotation matrices for each trace quadrature point.
153 /// Rotation storage
155
159
161
162 virtual void v_Solve(const int nDim,
163 const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
164 const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
165 Array<OneD, Array<OneD, NekDouble>> &flux) = 0;
166
168 const Array<OneD, const Array<OneD, NekDouble>> &normals);
172 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
173 const Array<OneD, const Array<OneD, NekDouble>> &normals,
174 const Array<OneD, const Array<OneD, NekDouble>> &vecLocs,
175 Array<OneD, Array<OneD, NekDouble>> &outarray);
177 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
178 const Array<OneD, const Array<OneD, NekDouble>> &normals,
179 const Array<OneD, const Array<OneD, NekDouble>> &vecLocs,
180 Array<OneD, Array<OneD, NekDouble>> &outarray);
181 SOLVER_UTILS_EXPORT bool CheckScalars(std::string name);
182 SOLVER_UTILS_EXPORT bool CheckVectors(std::string name);
183 SOLVER_UTILS_EXPORT bool CheckParams(std::string name);
184 SOLVER_UTILS_EXPORT bool CheckAuxScal(std::string name);
185 SOLVER_UTILS_EXPORT bool CheckAuxVec(std::string name);
186
188 const int nDim, const Array<OneD, const Array<OneD, NekDouble>> &Fwd,
189 const Array<OneD, const Array<OneD, NekDouble>> &Bwd,
190 const Array<OneD, const Array<OneD, NekDouble>> &normals,
192};
193
194/// A shared pointer to an EquationSystem object
195typedef std::shared_ptr<RiemannSolver> RiemannSolverSharedPtr;
196/// Datatype of the NekFactory used to instantiate classes derived
197/// from the RiemannSolver class.
198typedef LibUtilities::NekFactory<std::string, RiemannSolver,
202
203template <class T, typename = typename std::enable_if<
204 std::is_floating_point<T>::value ||
206inline void rotateToNormalKernel(T *in, T *rotMat, T *out)
207{
208
209 // Apply rotation matrices.
210 out[0] = in[0] * rotMat[0] + in[1] * rotMat[1] + in[2] * rotMat[2];
211
212 out[1] = in[0] * rotMat[3] + in[1] * rotMat[4] + in[2] * rotMat[5];
213
214 out[2] = in[0] * rotMat[6] + in[1] * rotMat[7] + in[2] * rotMat[8];
215}
216
217template <class T, typename = typename std::enable_if<
218 std::is_floating_point<T>::value ||
220inline void rotateFromNormalKernel(T *in, T *rotMat, T *out)
221{
222
223 // Apply rotation matrices.
224 out[0] = in[0] * rotMat[0] + in[1] * rotMat[3] + in[2] * rotMat[6];
225
226 out[1] = in[0] * rotMat[1] + in[1] * rotMat[4] + in[2] * rotMat[7];
227
228 out[2] = in[0] * rotMat[2] + in[1] * rotMat[5] + in[2] * rotMat[8];
229}
230
231} // namespace SolverUtils
232} // namespace Nektar
233
234#endif
#define SOLVER_UTILS_EXPORT
Provides a generic Factory class.
Definition: NekFactory.hpp:104
The RiemannSolver class provides an abstract interface under which solvers for various Riemann proble...
Definition: RiemannSolver.h:58
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:66
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::map< std::string, RSScalarFuncType > & GetScalars()
void SetVector(std::string name, RSVecFuncType fp)
Definition: RiemannSolver.h:82
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)
Definition: RiemannSolver.h:93
void SetParam(std::string name, FuncPointerT func, ObjectPointerT obj)
Definition: RiemannSolver.h:88
void SetScalar(std::string name, RSScalarFuncType fp)
Definition: RiemannSolver.h:71
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:77
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)
Definition: RiemannSolver.h:99
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.
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
SOLVER_UTILS_EXPORT RiemannSolver()
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.
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:52
std::function< const Array< OneD, const Array< OneD, NekDouble > > &()> RSVecFuncType
Definition: RiemannSolver.h:54
std::function< NekDouble()> RSParamFuncType
Definition: RiemannSolver.h:55
RiemannSolverFactory & GetRiemannSolverFactory()
std::shared_ptr< DNekBlkMat > DNekBlkMatSharedPtr
Definition: NekTypeDefs.hpp:77
double NekDouble