Nektar++
Static Public Member Functions | Static Public Attributes | Protected Member Functions | List of all members
Nektar::HLLCSolver Class Reference

#include <HLLCSolver.h>

Inheritance diagram for Nektar::HLLCSolver:
Inheritance graph
[legend]
Collaboration diagram for Nektar::HLLCSolver:
Collaboration graph
[legend]

Static Public Member Functions

static RiemannSolverSharedPtr create ()
 
static RiemannSolverSharedPtr create ()
 

Static Public Attributes

static std::string solverName
 

Protected Member Functions

 HLLCSolver ()
 
virtual void v_PointSolve (NekDouble rhoL, NekDouble rhouL, NekDouble rhovL, NekDouble rhowL, NekDouble EL, NekDouble rhoR, NekDouble rhouR, NekDouble rhovR, NekDouble rhowR, NekDouble ER, NekDouble &rhof, NekDouble &rhouf, NekDouble &rhovf, NekDouble &rhowf, NekDouble &Ef)
 HLLC Riemann solver. More...
 
virtual void v_PointSolveVisc (NekDouble rhoL, NekDouble rhouL, NekDouble rhovL, NekDouble rhowL, NekDouble EL, NekDouble EpsL, NekDouble rhoR, NekDouble rhouR, NekDouble rhovR, NekDouble rhowR, NekDouble ER, NekDouble EpsR, NekDouble &rhof, NekDouble &rhouf, NekDouble &rhovf, NekDouble &rhowf, NekDouble &Ef, NekDouble &Epsf)
 
 HLLCSolver ()
 
virtual void v_PointSolve (double hL, double huL, double hvL, double hR, double huR, double hvR, double &hf, double &huf, double &hvf)
 HLLC Riemann solver for the Nonlinear Shallow Water Equations. More...
 
- Protected Member Functions inherited from Nektar::NonlinearSWESolver
 NonlinearSWESolver ()
 
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)
 
virtual void v_ArraySolve (const Array< OneD, const Array< OneD, NekDouble > > &Fwd, const Array< OneD, const Array< OneD, NekDouble > > &Bwd, Array< OneD, Array< OneD, NekDouble > > &flux)
 
- Protected Member Functions inherited from Nektar::SolverUtils::RiemannSolver
SOLVER_UTILS_EXPORT RiemannSolver ()
 
void GenerateRotationMatrices (const Array< OneD, const Array< OneD, NekDouble > > &normals)
 Generate rotation matrices for 3D expansions. More...
 
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. More...
 
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. More...
 
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. More...
 
SOLVER_UTILS_EXPORT bool CheckScalars (std::string name)
 Determine whether a scalar has been defined in m_scalars. More...
 
SOLVER_UTILS_EXPORT bool CheckVectors (std::string name)
 Determine whether a vector has been defined in m_vectors. More...
 
SOLVER_UTILS_EXPORT bool CheckParams (std::string name)
 Determine whether a parameter has been defined in m_params. More...
 
SOLVER_UTILS_EXPORT bool CheckAuxScal (std::string name)
 Determine whether a scalar has been defined in m_auxScal. More...
 
SOLVER_UTILS_EXPORT bool CheckAuxVec (std::string name)
 Determine whether a vector has been defined in m_auxVec. More...
 
- Protected Member Functions inherited from Nektar::CompressibleSolver
 CompressibleSolver ()
 
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)
 
virtual void v_ArraySolve (const Array< OneD, const Array< OneD, NekDouble > > &Fwd, const Array< OneD, const Array< OneD, NekDouble > > &Bwd, Array< OneD, Array< OneD, NekDouble > > &flux)
 

Additional Inherited Members

- Public Member Functions inherited from Nektar::SolverUtils::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. More...
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetScalar (std::string name, FuncPointerT func, ObjectPointerT obj)
 
void SetScalar (std::string name, RSScalarFuncType fp)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetVector (std::string name, FuncPointerT func, ObjectPointerT obj)
 
void SetVector (std::string name, RSVecFuncType fp)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetParam (std::string name, FuncPointerT func, ObjectPointerT obj)
 
void SetParam (std::string name, RSParamFuncType fp)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetAuxScal (std::string name, FuncPointerT func, ObjectPointerT obj)
 
template<typename FuncPointerT , typename ObjectPointerT >
void SetAuxVec (std::string name, FuncPointerT func, ObjectPointerT obj)
 
std::map< std::string, RSScalarFuncType > & GetScalars ()
 
std::map< std::string, RSVecFuncType > & GetVectors ()
 
std::map< std::string, RSParamFuncType > & GetParams ()
 
- Public Attributes inherited from Nektar::SolverUtils::RiemannSolver
int m_spacedim
 
- Protected Attributes inherited from Nektar::NonlinearSWESolver
bool m_pointSolve
 
- Protected Attributes inherited from Nektar::SolverUtils::RiemannSolver
bool m_requiresRotation
 Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields. More...
 
std::map< std::string, RSScalarFuncTypem_scalars
 Map of scalar function types. More...
 
std::map< std::string, RSVecFuncTypem_vectors
 Map of vector function types. More...
 
std::map< std::string, RSParamFuncTypem_params
 Map of parameter function types. More...
 
std::map< std::string, RSScalarFuncTypem_auxScal
 Map of auxiliary scalar function types. More...
 
std::map< std::string, RSVecFuncTypem_auxVec
 Map of auxiliary vector function types. More...
 
Array< OneD, Array< OneD, NekDouble > > m_rotMat
 Rotation matrices for each trace quadrature point. More...
 
Array< OneD, Array< OneD, Array< OneD, NekDouble > > > m_rotStorage
 Rotation storage. More...
 
- Protected Attributes inherited from Nektar::CompressibleSolver
bool m_pointSolve
 

Detailed Description

Definition at line 43 of file CompressibleFlowSolver/RiemannSolvers/HLLCSolver.h.

Constructor & Destructor Documentation

Nektar::HLLCSolver::HLLCSolver ( )
protected

Definition at line 44 of file CompressibleFlowSolver/RiemannSolvers/HLLCSolver.cpp.

Referenced by create().

45  {
46 
47  }
Nektar::HLLCSolver::HLLCSolver ( )
protected

Member Function Documentation

static RiemannSolverSharedPtr Nektar::HLLCSolver::create ( )
inlinestatic

Definition at line 46 of file CompressibleFlowSolver/RiemannSolvers/HLLCSolver.h.

References HLLCSolver().

47  {
48  return RiemannSolverSharedPtr(new HLLCSolver());
49  }
boost::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.
static RiemannSolverSharedPtr Nektar::HLLCSolver::create ( )
inlinestatic

Definition at line 46 of file ShallowWaterSolver/RiemannSolvers/HLLCSolver.h.

References HLLCSolver().

47  {
49  new HLLCSolver());
50  }
boost::shared_ptr< RiemannSolver > RiemannSolverSharedPtr
A shared pointer to an EquationSystem object.
void Nektar::HLLCSolver::v_PointSolve ( NekDouble  rhoL,
NekDouble  rhouL,
NekDouble  rhovL,
NekDouble  rhowL,
NekDouble  EL,
NekDouble  rhoR,
NekDouble  rhouR,
NekDouble  rhovR,
NekDouble  rhowR,
NekDouble  ER,
NekDouble rhof,
NekDouble rhouf,
NekDouble rhovf,
NekDouble rhowf,
NekDouble Ef 
)
protectedvirtual

HLLC Riemann solver.

Parameters
rhoLDensity left state.
rhoRDensity right state.
rhouLx-momentum component left state.
rhouRx-momentum component right state.
rhovLy-momentum component left state.
rhovRy-momentum component right state.
rhowLz-momentum component left state.
rhowRz-momentum component right state.
ELEnergy left state.
EREnergy right state.
rhofComputed Riemann flux for density.
rhoufComputed Riemann flux for x-momentum component
rhovfComputed Riemann flux for y-momentum component
rhowfComputed Riemann flux for z-momentum component
EfComputed Riemann flux for energy.

Reimplemented from Nektar::CompressibleSolver.

Definition at line 68 of file CompressibleFlowSolver/RiemannSolvers/HLLCSolver.cpp.

References Nektar::SolverUtils::RiemannSolver::m_params.

72  {
73  static NekDouble gamma = m_params["gamma"]();
74 
75  // Left and Right velocities
76  NekDouble uL = rhouL / rhoL;
77  NekDouble vL = rhovL / rhoL;
78  NekDouble wL = rhowL / rhoL;
79  NekDouble uR = rhouR / rhoR;
80  NekDouble vR = rhovR / rhoR;
81  NekDouble wR = rhowR / rhoR;
82 
83  // Left and right pressure, sound speed and enthalpy.
84  NekDouble pL = (gamma - 1.0) *
85  (EL - 0.5 * (rhouL * uL + rhovL * vL + rhowL * wL));
86  NekDouble pR = (gamma - 1.0) *
87  (ER - 0.5 * (rhouR * uR + rhovR * vR + rhowR * wR));
88  NekDouble cL = sqrt(gamma * pL / rhoL);
89  NekDouble cR = sqrt(gamma * pR / rhoR);
90  NekDouble hL = (EL + pL) / rhoL;
91  NekDouble hR = (ER + pR) / rhoR;
92 
93  // Square root of rhoL and rhoR.
94  NekDouble srL = sqrt(rhoL);
95  NekDouble srR = sqrt(rhoR);
96  NekDouble srLR = srL + srR;
97 
98  // Velocity Roe averages
99  NekDouble uRoe = (srL * uL + srR * uR) / srLR;
100  NekDouble vRoe = (srL * vL + srR * vR) / srLR;
101  NekDouble wRoe = (srL * wL + srR * wR) / srLR;
102  NekDouble hRoe = (srL * hL + srR * hR) / srLR;
103  NekDouble cRoe = sqrt((gamma - 1.0)*(hRoe - 0.5 *
104  (uRoe * uRoe + vRoe * vRoe + wRoe * wRoe)));
105 
106  // Maximum wave speeds
107  NekDouble SL = std::min(uL-cL, uRoe-cRoe);
108  NekDouble SR = std::max(uR+cR, uRoe+cRoe);
109 
110  // HLLC Riemann fluxes (positive case)
111  if (SL >= 0)
112  {
113  rhof = rhouL;
114  rhouf = rhouL * uL + pL;
115  rhovf = rhouL * vL;
116  rhowf = rhouL * wL;
117  Ef = uL * (EL + pL);
118  }
119  // HLLC Riemann fluxes (negative case)
120  else if (SR <= 0)
121  {
122  rhof = rhouR;
123  rhouf = rhouR * uR + pR;
124  rhovf = rhouR * vR;
125  rhowf = rhouR * wR;
126  Ef = uR * (ER + pR);
127  }
128  // HLLC Riemann fluxes (general case (SL < 0 | SR > 0)
129  else
130  {
131  NekDouble SM = (pR - pL + rhouL * (SL - uL) - rhouR * (SR - uR)) /
132  (rhoL * (SL - uL) - rhoR * (SR - uR));
133  NekDouble rhoML = rhoL * (SL - uL) / (SL - SM);
134  NekDouble rhouML = rhoML * SM;
135  NekDouble rhovML = rhoML * vL;
136  NekDouble rhowML = rhoML * wL;
137  NekDouble EML = rhoML * (EL / rhoL +
138  (SM - uL) * (SM + pL / (rhoL * (SL - uL))));
139 
140  NekDouble rhoMR = rhoR * (SR - uR) / (SR - SM);
141  NekDouble rhouMR = rhoMR * SM;
142  NekDouble rhovMR = rhoMR * vR;
143  NekDouble rhowMR = rhoMR * wR;
144  NekDouble EMR = rhoMR * (ER / rhoR +
145  (SM - uR) * (SM + pR / (rhoR * (SR - uR))));
146 
147  if (SL < 0.0 && SM >= 0.0)
148  {
149  rhof = rhouL + SL * (rhoML - rhoL);
150  rhouf = rhouL * uL + pL + SL * (rhouML - rhouL);
151  rhovf = rhouL * vL + SL * (rhovML - rhovL);
152  rhowf = rhouL * wL + SL * (rhowML - rhowL);
153  Ef = uL * (EL + pL) + SL * (EML - EL);
154  }
155  else if(SM < 0.0 && SR > 0.0)
156  {
157  rhof = rhouR + SR * (rhoMR - rhoR);
158  rhouf = rhouR * uR + pR + SR * (rhouMR - rhouR);
159  rhovf = rhouR * vR + SR * (rhovMR - rhovR);
160  rhowf = rhouR * wR + SR * (rhowMR - rhowR);
161  Ef = uR * (ER + pR) + SR * (EMR - ER);
162  }
163  }
164  }
double NekDouble
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
void Nektar::HLLCSolver::v_PointSolve ( double  hL,
double  huL,
double  hvL,
double  hR,
double  huR,
double  hvR,
double &  hf,
double &  huf,
double &  hvf 
)
protectedvirtual

HLLC Riemann solver for the Nonlinear Shallow Water Equations.

Parameters
hLWater depth left state.
hRWater depth right state.
huLx-momentum component left state.
huRx-momentum component right state.
hvLy-momentum component left state.
hvRy-momentum component right state.
hfComputed Riemann flux for density.
hufComputed Riemann flux for x-momentum component
hvfComputed Riemann flux for y-momentum component

Reimplemented from Nektar::NonlinearSWESolver.

Definition at line 64 of file ShallowWaterSolver/RiemannSolvers/HLLCSolver.cpp.

References ASSERTL0, and Nektar::SolverUtils::RiemannSolver::m_params.

68  {
69  static NekDouble g = m_params["gravity"]();
70 
71  // Left and Right velocities
72  NekDouble uL = huL / hL;
73  NekDouble vL = hvL / hL;
74  NekDouble uR = huR / hR;
75  NekDouble vR = hvR / hR;
76 
77 
78  // Left and right wave speeds
79  NekDouble cL = sqrt(g * hL);
80  NekDouble cR = sqrt(g * hR);
81 
82  // the two-rarefaction wave assumption
83  NekDouble hC,huC,hvC,SL,SR,hstar,ustar,Sstar;
84  hstar = 0.5*(cL + cR) + 0.25*(uL - uR);
85  hstar *= hstar;
86  hstar *= (1.0/g);
87  ustar = 0.5*(uL + uR) + cL - cR;
88 
89 
90  // Compute SL
91  if (hstar > hL)
92  SL = uL - cL * sqrt(0.5*((hstar*hstar + hstar*hL)/(hL*hL)));
93  else
94  SL = uL - cL;
95 
96  // Compute SR
97  if (hstar > hR)
98  SR = uR + cR * sqrt(0.5*((hstar*hstar + hstar*hR)/(hR*hR)));
99  else
100  SR = uR + cR;
101 
102  if (fabs(hR*(uR-SR)-hL*(uL-SL)) <= 1.0e-10)
103  Sstar = ustar;
104  else
105  Sstar = (SL*hR*(uR-SR)-SR*hL*(uL-SL))/(hR*(uR-SR)-hL*(uL-SL));
106 
107  if (SL >= 0)
108  {
109  hf = hL * uL;
110  huf = uL * uL * hL + 0.5 * g * hL * hL;
111  hvf = hL * uL * vL;
112  }
113  else if (SR <= 0)
114  {
115  hf = hR * uR;
116  huf = uR * uR * hR + 0.5 * g * hR * hR;
117  hvf = hR * uR *vR;
118  }
119  else if ((SL < 0) && (Sstar >= 0))
120  {
121  hC = hL * ((SL - uL) / (SL - Sstar));
122  huC = hC * Sstar;
123  hvC = hC * vL;
124 
125  hf = hL*uL + SL * (hC - hL);
126  huf = (uL*uL*hL+0.5*g*hL*hL) + SL * (huC - hL*uL);
127  hvf = (uL*vL*hL) + SL * (hvC - hL*vL);
128  }
129  else if ((SR > 0) && (Sstar <= 0))
130  {
131  hC = hR * ((SR - uR) / (SR - Sstar));
132  huC = hC * Sstar;
133  hvC = hC * vR;
134 
135  hf = hR*uR + SR * (hC - hR);
136  huf = (uR*uR*hR+0.5*g*hR*hR) + SR * (huC - hR*uR);
137  hvf = (uR*vR*hR) + SR * (hvC - hR*vR);
138  }
139  else
140  {
141  ASSERTL0(false,"Error in HLLC solver -- non physical combination of SR, SL and Sstar");
142  }
143  }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:161
double NekDouble
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
void Nektar::HLLCSolver::v_PointSolveVisc ( NekDouble  rhoL,
NekDouble  rhouL,
NekDouble  rhovL,
NekDouble  rhowL,
NekDouble  EL,
NekDouble  EpsL,
NekDouble  rhoR,
NekDouble  rhouR,
NekDouble  rhovR,
NekDouble  rhowR,
NekDouble  ER,
NekDouble  EpsR,
NekDouble rhof,
NekDouble rhouf,
NekDouble rhovf,
NekDouble rhowf,
NekDouble Ef,
NekDouble Epsf 
)
protectedvirtual

Reimplemented from Nektar::CompressibleSolver.

Definition at line 166 of file CompressibleFlowSolver/RiemannSolvers/HLLCSolver.cpp.

References Nektar::SolverUtils::RiemannSolver::m_params.

170  {
171  static NekDouble gamma = m_params["gamma"]();
172 
173  // Left and Right velocities
174  NekDouble uL = rhouL / rhoL;
175  NekDouble vL = rhovL / rhoL;
176  NekDouble wL = rhowL / rhoL;
177  NekDouble uR = rhouR / rhoR;
178  NekDouble vR = rhovR / rhoR;
179  NekDouble wR = rhowR / rhoR;
180 
181  // Left and right pressure, sound speed and enthalpy.
182  NekDouble pL = (gamma - 1.0) *
183  (EL - 0.5 * (rhouL * uL + rhovL * vL + rhowL * wL));
184  NekDouble pR = (gamma - 1.0) *
185  (ER - 0.5 * (rhouR * uR + rhovR * vR + rhowR * wR));
186  NekDouble cL = sqrt(gamma * pL / rhoL);
187  NekDouble cR = sqrt(gamma * pR / rhoR);
188  NekDouble hL = (EL + pL) / rhoL;
189  NekDouble hR = (ER + pR) / rhoR;
190 
191  // Square root of rhoL and rhoR.
192  NekDouble srL = sqrt(rhoL);
193  NekDouble srR = sqrt(rhoR);
194  NekDouble srLR = srL + srR;
195 
196  // Velocity Roe averages
197  NekDouble uRoe = (srL * uL + srR * uR) / srLR;
198  NekDouble vRoe = (srL * vL + srR * vR) / srLR;
199  NekDouble wRoe = (srL * wL + srR * wR) / srLR;
200  NekDouble hRoe = (srL * hL + srR * hR) / srLR;
201  NekDouble cRoe = sqrt((gamma - 1.0)*(hRoe - 0.5 *
202  (uRoe * uRoe + vRoe * vRoe + wRoe * wRoe)));
203 
204  // Maximum wave speeds
205  NekDouble SL = std::min(uL-cL, uRoe-cRoe);
206  NekDouble SR = std::max(uR+cR, uRoe+cRoe);
207 
208  // HLLC Riemann fluxes (positive case)
209  if (SL >= 0)
210  {
211  rhof = rhouL;
212  rhouf = rhouL * uL + pL;
213  rhovf = rhouL * vL;
214  rhowf = rhouL * wL;
215  Ef = uL * (EL + pL);
216  Epsf = 0.0;
217  }
218  // HLLC Riemann fluxes (negative case)
219  else if (SR <= 0)
220  {
221  rhof = rhouR;
222  rhouf = rhouR * uR + pR;
223  rhovf = rhouR * vR;
224  rhowf = rhouR * wR;
225  Ef = uR * (ER + pR);
226  Epsf = 0.0;
227  }
228  // HLLC Riemann fluxes (general case (SL < 0 | SR > 0)
229  else
230  {
231  NekDouble SM = (pR - pL + rhouL * (SL - uL) - rhouR * (SR - uR)) /
232  (rhoL * (SL - uL) - rhoR * (SR - uR));
233  NekDouble rhoML = rhoL * (SL - uL) / (SL - SM);
234  NekDouble rhouML = rhoML * SM;
235  NekDouble rhovML = rhoML * vL;
236  NekDouble rhowML = rhoML * wL;
237  NekDouble EML = rhoML * (EL / rhoL +
238  (SM - uL) * (SM + pL / (rhoL * (SL - uL))));
239  NekDouble EpsML = EpsL * (SL - uL) / (SL - SM);
240 
241  NekDouble rhoMR = rhoR * (SR - uR) / (SR - SM);
242  NekDouble rhouMR = rhoMR * SM;
243  NekDouble rhovMR = rhoMR * vR;
244  NekDouble rhowMR = rhoMR * wR;
245  NekDouble EMR = rhoMR * (ER / rhoR +
246  (SM - uR) * (SM + pR / (rhoR * (SR - uR))));
247  NekDouble EpsMR = EpsR * (SL - uR) / (SL - SM);
248 
249  if (SL < 0.0 && SM >= 0.0)
250  {
251  rhof = rhouL + SL * (rhoML - rhoL);
252  rhouf = rhouL * uL + pL + SL * (rhouML - rhouL);
253  rhovf = rhouL * vL + SL * (rhovML - rhovL);
254  rhowf = rhouL * wL + SL * (rhowML - rhowL);
255  Ef = uL * (EL + pL) + SL * (EML - EL);
256  Epsf = 0.0 + SL * (EpsML - EpsL);
257  }
258  else if(SM < 0.0 && SR > 0.0)
259  {
260  rhof = rhouR + SR * (rhoMR - rhoR);
261  rhouf = rhouR * uR + pR + SR * (rhouMR - rhouR);
262  rhovf = rhouR * vR + SR * (rhovMR - rhovR);
263  rhowf = rhouR * wR + SR * (rhowMR - rhowR);
264  Ef = uR * (ER + pR) + SR * (EMR - ER);
265  Epsf = 0.0 + SR * (EpsMR - EpsR);
266  }
267  }
268  }
double NekDouble
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.

Member Data Documentation

static std::string Nektar::HLLCSolver::solverName
static