Nektar++
Protected Member Functions | Protected Attributes | List of all members
Nektar::CompressibleSolver Class Reference

#include <CompressibleSolver.h>

Inheritance diagram for Nektar::CompressibleSolver:
[legend]

Protected Member Functions

 CompressibleSolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
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)
 
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)
 
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)
 
NekDouble GetRoeSoundSpeed (NekDouble rhoL, NekDouble pL, NekDouble eL, NekDouble HL, NekDouble srL, NekDouble rhoR, NekDouble pR, NekDouble eR, NekDouble HR, NekDouble srR, NekDouble HRoe, NekDouble URoe2, NekDouble srLR)
 
- Protected Member Functions inherited from Nektar::SolverUtils::RiemannSolver
SOLVER_UTILS_EXPORT RiemannSolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
virtual 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 Attributes

bool m_pointSolve
 
EquationOfStateSharedPtr m_eos
 
bool m_idealGas
 
- 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...
 

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
 

Detailed Description

Definition at line 47 of file CompressibleSolver.h.

Constructor & Destructor Documentation

◆ CompressibleSolver()

Nektar::CompressibleSolver::CompressibleSolver ( const LibUtilities::SessionReaderSharedPtr pSession)
protected

Definition at line 41 of file CompressibleSolver.cpp.

References Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::GetEquationOfStateFactory(), m_eos, m_idealGas, and Nektar::SolverUtils::RiemannSolver::m_requiresRotation.

43  : RiemannSolver(pSession), m_pointSolve(true)
44  {
45  m_requiresRotation = true;
46 
47  // Create equation of state object
48  std::string eosType;
49  pSession->LoadSolverInfo("EquationOfState",
50  eosType, "IdealGas");
52  .CreateInstance(eosType, pSession);
53  // Check if using ideal gas
54  m_idealGas = boost::iequals(eosType,"IdealGas");
55  }
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
EquationOfStateFactory & GetEquationOfStateFactory()
Declaration of the equation of state factory singleton.
SOLVER_UTILS_EXPORT RiemannSolver(const LibUtilities::SessionReaderSharedPtr &pSession)
bool m_requiresRotation
Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields...
EquationOfStateSharedPtr m_eos

Member Function Documentation

◆ GetRoeSoundSpeed()

NekDouble Nektar::CompressibleSolver::GetRoeSoundSpeed ( NekDouble  rhoL,
NekDouble  pL,
NekDouble  eL,
NekDouble  HL,
NekDouble  srL,
NekDouble  rhoR,
NekDouble  pR,
NekDouble  eR,
NekDouble  HR,
NekDouble  srR,
NekDouble  HRoe,
NekDouble  URoe2,
NekDouble  srLR 
)
protected

Definition at line 133 of file CompressibleSolver.cpp.

References Nektar::NekConstants::kNekZeroTol, m_eos, m_idealGas, and Nektar::SolverUtils::RiemannSolver::m_params.

Referenced by Nektar::HLLCSolver::v_PointSolve(), Nektar::HLLSolver::v_PointSolve(), Nektar::LaxFriedrichsSolver::v_PointSolve(), and Nektar::HLLCSolver::v_PointSolveVisc().

137  {
138  boost::ignore_unused(HL, srL, HR, srR, srLR);
139 
140  static NekDouble gamma = m_params["gamma"]();
141  NekDouble cRoe;
142  if(m_idealGas)
143  {
144  cRoe = sqrt((gamma - 1.0)*(HRoe - 0.5 * URoe2));
145  }
146  else
147  {
148  // Calculate static enthalpy of left and right states
149  NekDouble hL = eL + pL/rhoL;
150  NekDouble hR = eR + pR/rhoR;
151 
152  // Get partial derivatives of P(rho,e)
153  NekDouble dpdeL = m_eos->GetDPDe_rho(rhoL,eL);
154  NekDouble dpdeR = m_eos->GetDPDe_rho(rhoR,eR);
155  NekDouble dpdrhoL = m_eos->GetDPDrho_e(rhoL,eL);
156  NekDouble dpdrhoR = m_eos->GetDPDrho_e(rhoR,eR);
157 
158  // Evaluate chi and kappa parameters
159  NekDouble chiL = dpdrhoL - eL / rhoL * dpdeL;
160  NekDouble kappaL = dpdeL / rhoL;
161  NekDouble chiR = dpdrhoR - eR / rhoR * dpdeR;
162  NekDouble kappaR = dpdeR / rhoR;
163 
164  //
165  // Calculate interface speed of sound using procedure from
166  // Vinokur, M.; Montagné, J.-L. "Generalized Flux-Vector
167  // Splitting and Roe Average for an Equilibrium Real Gas",
168  // JCP (1990).
169  //
170 
171  // Calculate averages
172  NekDouble avgChi = 0.5 * (chiL + chiR);
173  NekDouble avgKappa = 0.5 * (kappaL + kappaR);
174  NekDouble avgKappaH = 0.5 * (kappaL*hL + kappaR*hR);
175 
176  // Calculate jumps
177  NekDouble deltaP = pR - pL;
178  NekDouble deltaRho = rhoR - rhoL;
179  NekDouble deltaRhoe = rhoR*eR - rhoL*eL;
180 
181  // Evaluate dP: equation (64) from Vinokur-Montagné
182  NekDouble dP = deltaP - avgChi * deltaRho - avgKappa * deltaRhoe;
183  // s (eq 66)
184  NekDouble s = avgChi + avgKappaH;
185  // D (eq 65)
186  NekDouble D = (s*deltaRho)*(s*deltaRho) + deltaP*deltaP;
187  // chiRoe and kappaRoe (eq 66)
188  NekDouble chiRoe, kappaRoe;
189  NekDouble fac = D - deltaP*deltaRho;
190  if( std::abs(fac) > NekConstants::kNekZeroTol)
191  {
192  chiRoe = (D*avgChi + s*s*deltaRho*dP) / fac;
193  kappaRoe = D*avgKappa / fac;
194  }
195  else
196  {
197  chiRoe = avgChi;
198  kappaRoe = avgKappa;
199  }
200  // Speed of sound (eq 53)
201  cRoe = sqrt( chiRoe + kappaRoe*(HRoe - 0.5 * URoe2));
202  }
203  return cRoe;
204  }
static const NekDouble kNekZeroTol
double NekDouble
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
EquationOfStateSharedPtr m_eos

◆ v_ArraySolve()

virtual void Nektar::CompressibleSolver::v_ArraySolve ( const Array< OneD, const Array< OneD, NekDouble > > &  Fwd,
const Array< OneD, const Array< OneD, NekDouble > > &  Bwd,
Array< OneD, Array< OneD, NekDouble > > &  flux 
)
inlineprotectedvirtual

Reimplemented in Nektar::AverageSolver.

Definition at line 63 of file CompressibleSolver.h.

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by v_Solve().

67  {
68  boost::ignore_unused(Fwd, Bwd, flux);
70  "This function should be defined by subclasses.");
71  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:209

◆ v_PointSolve()

virtual void Nektar::CompressibleSolver::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 
)
inlineprotectedvirtual

Reimplemented in Nektar::AUSM2Solver, Nektar::ExactSolverToro, Nektar::LaxFriedrichsSolver, Nektar::AUSM0Solver, Nektar::AUSM1Solver, Nektar::AUSM3Solver, Nektar::HLLSolver, Nektar::RoeSolver, and Nektar::HLLCSolver.

Definition at line 73 of file CompressibleSolver.h.

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by v_Solve().

77  {
78  boost::ignore_unused(rhoL, rhouL, rhovL, rhowL, EL,
79  rhoR, rhouR, rhovR, rhowR, ER,
80  rhof, rhouf, rhovf, rhowf, Ef);
82  "This function should be defined by subclasses.");
83  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:209

◆ v_PointSolveVisc()

virtual void Nektar::CompressibleSolver::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 
)
inlineprotectedvirtual

Reimplemented in Nektar::HLLCSolver.

Definition at line 85 of file CompressibleSolver.h.

References Nektar::ErrorUtil::efatal, and NEKERROR.

Referenced by v_Solve().

89  {
90  boost::ignore_unused(rhoL, rhouL, rhovL, rhowL, EL, EpsL,
91  rhoR, rhouR, rhovR, rhowR, ER, EpsR,
92  rhof, rhouf, rhovf, rhowf, Ef, Epsf);
94  "This function should be defined by subclasses.");
95  }
#define NEKERROR(type, msg)
Assert Level 0 – Fundamental assert which is used whether in FULLDEBUG, DEBUG or OPT compilation mod...
Definition: ErrorUtil.hpp:209

◆ v_Solve()

void Nektar::CompressibleSolver::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 
)
protectedvirtual

Implements Nektar::SolverUtils::RiemannSolver.

Definition at line 57 of file CompressibleSolver.cpp.

References m_pointSolve, v_ArraySolve(), v_PointSolve(), and v_PointSolveVisc().

62  {
63  if (m_pointSolve)
64  {
65  int expDim = nDim;
66  int nvariables = Fwd.num_elements();
67 
68  NekDouble rhouf, rhovf;
69 
70  // Check if PDE-based SC is used
71  if (expDim == 1)
72  {
73  for (int i = 0; i < Fwd[0].num_elements(); ++i)
74  {
76  Fwd [0][i], Fwd [1][i], 0.0, 0.0, Fwd [2][i],
77  Bwd [0][i], Bwd [1][i], 0.0, 0.0, Bwd [2][i],
78  flux[0][i], flux[1][i], rhouf, rhovf, flux[2][i]);
79  }
80  }
81  else if (expDim == 2)
82  {
83  if (nvariables == expDim+2)
84  {
85  for (int i = 0; i < Fwd[0].num_elements(); ++i)
86  {
88  Fwd [0][i], Fwd [1][i], Fwd [2][i], 0.0, Fwd [3][i],
89  Bwd [0][i], Bwd [1][i], Bwd [2][i], 0.0, Bwd [3][i],
90  flux[0][i], flux[1][i], flux[2][i], rhovf, flux[3][i]);
91  }
92  }
93 
94  if (nvariables > expDim+2)
95  {
96  for (int i = 0; i < Fwd[0].num_elements(); ++i)
97  {
99  Fwd [0][i], Fwd [1][i], Fwd [2][i], 0.0, Fwd [3][i], Fwd [4][i],
100  Bwd [0][i], Bwd [1][i], Bwd [2][i], 0.0, Bwd [3][i], Bwd [4][i],
101  flux[0][i], flux[1][i], flux[2][i], rhovf, flux[3][i], flux[4][i]);
102  }
103  }
104 
105  }
106  else if (expDim == 3)
107  {
108  for (int i = 0; i < Fwd[0].num_elements(); ++i)
109  {
110  v_PointSolve(
111  Fwd [0][i], Fwd [1][i], Fwd [2][i], Fwd [3][i], Fwd [4][i],
112  Bwd [0][i], Bwd [1][i], Bwd [2][i], Bwd [3][i], Bwd [4][i],
113  flux[0][i], flux[1][i], flux[2][i], flux[3][i], flux[4][i]);
114  }
115  if (nvariables > expDim+2)
116  {
117  for (int i = 0; i < Fwd[0].num_elements(); ++i)
118  {
120  Fwd [0][i], Fwd [1][i], Fwd [2][i], Fwd [3][i], Fwd [4][i], Fwd [5][i],
121  Bwd [0][i], Bwd [1][i], Bwd [2][i], Bwd [3][i], Bwd [4][i], Bwd [5][i],
122  flux[0][i], flux[1][i], flux[2][i], flux[3][i], flux[4][i], flux[5][i]);
123  }
124  }
125  }
126  }
127  else
128  {
129  v_ArraySolve(Fwd, Bwd, flux);
130  }
131  }
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)
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)
double NekDouble
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)

Member Data Documentation

◆ m_eos

EquationOfStateSharedPtr Nektar::CompressibleSolver::m_eos
protected

◆ m_idealGas

bool Nektar::CompressibleSolver::m_idealGas
protected

Definition at line 52 of file CompressibleSolver.h.

Referenced by CompressibleSolver(), and GetRoeSoundSpeed().

◆ m_pointSolve

bool Nektar::CompressibleSolver::m_pointSolve
protected

Definition at line 50 of file CompressibleSolver.h.

Referenced by Nektar::AverageSolver::AverageSolver(), and v_Solve().