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

#include <CompressibleSolver.h>

Inheritance diagram for Nektar::CompressibleSolver:
[legend]

Protected Types

using ND = NekDouble
 

Protected Member Functions

 CompressibleSolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 Session ctor. More...
 
 CompressibleSolver ()
 Programmatic ctor. More...
 
virtual void v_Solve (const int nDim, const Array< OneD, const Array< OneD, ND >> &Fwd, const Array< OneD, const Array< OneD, ND >> &Bwd, Array< OneD, Array< OneD, ND >> &flux) override
 
virtual void v_ArraySolve (const Array< OneD, const Array< OneD, ND >> &Fwd, const Array< OneD, const Array< OneD, ND >> &Bwd, Array< OneD, Array< OneD, ND >> &flux)
 
virtual void v_PointSolve (ND rhoL, ND rhouL, ND rhovL, ND rhowL, ND EL, ND rhoR, ND rhouR, ND rhovR, ND rhowR, ND ER, ND &rhof, ND &rhouf, ND &rhovf, ND &rhowf, ND &Ef)
 
ND GetRoeSoundSpeed (ND rhoL, ND pL, ND eL, ND HL, ND srL, ND rhoR, ND pR, ND eR, ND HR, ND srR, ND HRoe, ND URoe2, ND srLR)
 
- Protected Member Functions inherited from Nektar::SolverUtils::RiemannSolver
SOLVER_UTILS_EXPORT RiemannSolver ()
 
SOLVER_UTILS_EXPORT RiemannSolver (const LibUtilities::SessionReaderSharedPtr &pSession)
 
virtual SOLVER_UTILS_EXPORT ~RiemannSolver ()
 
SOLVER_UTILS_EXPORT 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...
 
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)
 

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)
 
void SetAuxVec (std::string name, RSVecFuncType fp)
 
std::map< std::string, RSScalarFuncType > & GetScalars ()
 
std::map< std::string, RSVecFuncType > & GetVectors ()
 
std::map< std::string, RSParamFuncType > & GetParams ()
 
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. More...
 
- Public Attributes inherited from Nektar::SolverUtils::RiemannSolver
int m_spacedim
 

Detailed Description

Definition at line 47 of file CompressibleSolver.h.

Member Typedef Documentation

◆ ND

Definition at line 60 of file CompressibleSolver.h.

Constructor & Destructor Documentation

◆ CompressibleSolver() [1/2]

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

Session ctor.

Definition at line 41 of file CompressibleSolver.cpp.

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", eosType, "IdealGas");
50  m_eos = GetEquationOfStateFactory().CreateInstance(eosType, pSession);
51  // Check if using ideal gas
52  m_idealGas = boost::iequals(eosType, "IdealGas");
53 }
EquationOfStateSharedPtr m_eos
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
bool m_requiresRotation
Indicates whether the Riemann solver requires a rotation to be applied to the velocity fields.
SOLVER_UTILS_EXPORT RiemannSolver()
EquationOfStateFactory & GetEquationOfStateFactory()
Declaration of the equation of state factory singleton.

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

◆ CompressibleSolver() [2/2]

Nektar::CompressibleSolver::CompressibleSolver ( )
protected

Programmatic ctor.

Definition at line 55 of file CompressibleSolver.cpp.

56  : RiemannSolver(), m_pointSolve(true), m_idealGas(true)
57 {
58  m_requiresRotation = true;
59 }

References Nektar::SolverUtils::RiemannSolver::m_requiresRotation.

Member Function Documentation

◆ GetRoeSoundSpeed()

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

Definition at line 107 of file CompressibleSolver.cpp.

111 {
112  boost::ignore_unused(HL, srL, HR, srR, srLR);
113 
114  static NekDouble gamma = m_params["gamma"]();
115  NekDouble cRoe;
116  if (m_idealGas)
117  {
118  cRoe = sqrt((gamma - 1.0) * (HRoe - 0.5 * URoe2));
119  }
120  else
121  {
122  // Calculate static enthalpy of left and right states
123  NekDouble hL = eL + pL / rhoL;
124  NekDouble hR = eR + pR / rhoR;
125 
126  // Get partial derivatives of P(rho,e)
127  NekDouble dpdeL = m_eos->GetDPDe_rho(rhoL, eL);
128  NekDouble dpdeR = m_eos->GetDPDe_rho(rhoR, eR);
129  NekDouble dpdrhoL = m_eos->GetDPDrho_e(rhoL, eL);
130  NekDouble dpdrhoR = m_eos->GetDPDrho_e(rhoR, eR);
131 
132  // Evaluate chi and kappa parameters
133  NekDouble chiL = dpdrhoL - eL / rhoL * dpdeL;
134  NekDouble kappaL = dpdeL / rhoL;
135  NekDouble chiR = dpdrhoR - eR / rhoR * dpdeR;
136  NekDouble kappaR = dpdeR / rhoR;
137 
138  //
139  // Calculate interface speed of sound using procedure from
140  // Vinokur, M.; Montagné, J.-L. "Generalized Flux-Vector
141  // Splitting and Roe Average for an Equilibrium Real Gas",
142  // JCP (1990).
143  //
144 
145  // Calculate averages
146  NekDouble avgChi = 0.5 * (chiL + chiR);
147  NekDouble avgKappa = 0.5 * (kappaL + kappaR);
148  NekDouble avgKappaH = 0.5 * (kappaL * hL + kappaR * hR);
149 
150  // Calculate jumps
151  NekDouble deltaP = pR - pL;
152  NekDouble deltaRho = rhoR - rhoL;
153  NekDouble deltaRhoe = rhoR * eR - rhoL * eL;
154 
155  // Evaluate dP: equation (64) from Vinokur-Montagné
156  NekDouble dP = deltaP - avgChi * deltaRho - avgKappa * deltaRhoe;
157  // s (eq 66)
158  NekDouble s = avgChi + avgKappaH;
159  // D (eq 65)
160  NekDouble D = (s * deltaRho) * (s * deltaRho) + deltaP * deltaP;
161  // chiRoe and kappaRoe (eq 66)
162  NekDouble chiRoe, kappaRoe;
163  NekDouble fac = D - deltaP * deltaRho;
165  {
166  chiRoe = (D * avgChi + s * s * deltaRho * dP) / fac;
167  kappaRoe = D * avgKappa / fac;
168  }
169  else
170  {
171  chiRoe = avgChi;
172  kappaRoe = avgKappa;
173  }
174  // Speed of sound (eq 53)
175  cRoe = sqrt(chiRoe + kappaRoe * (HRoe - 0.5 * URoe2));
176  }
177  return cRoe;
178 }
std::map< std::string, RSParamFuncType > m_params
Map of parameter function types.
static const NekDouble kNekZeroTol
double NekDouble
scalarT< T > abs(scalarT< T > in)
Definition: scalar.hpp:298
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References tinysimd::abs(), Nektar::NekConstants::kNekZeroTol, m_eos, m_idealGas, Nektar::SolverUtils::RiemannSolver::m_params, and tinysimd::sqrt().

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

◆ v_ArraySolve()

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

Reimplemented in Nektar::AverageSolver, and Nektar::RoeSolver.

Definition at line 67 of file CompressibleSolver.h.

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

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

Referenced by v_Solve().

◆ v_PointSolve()

virtual void Nektar::CompressibleSolver::v_PointSolve ( ND  rhoL,
ND  rhouL,
ND  rhovL,
ND  rhowL,
ND  EL,
ND  rhoR,
ND  rhouR,
ND  rhovR,
ND  rhowR,
ND  ER,
ND rhof,
ND rhouf,
ND rhovf,
ND rhowf,
ND Ef 
)
inlineprotectedvirtual

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

Definition at line 76 of file CompressibleSolver.h.

79  {
80  boost::ignore_unused(rhoL, rhouL, rhovL, rhowL, EL, rhoR, rhouR, rhovR,
81  rhowR, ER, rhof, rhouf, rhovf, rhowf, Ef);
83  "This function should be defined by subclasses.");
84  }

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

Referenced by v_Solve().

◆ v_Solve()

void Nektar::CompressibleSolver::v_Solve ( const int  nDim,
const Array< OneD, const Array< OneD, ND >> &  Fwd,
const Array< OneD, const Array< OneD, ND >> &  Bwd,
Array< OneD, Array< OneD, ND >> &  flux 
)
overrideprotectedvirtual

Implements Nektar::SolverUtils::RiemannSolver.

Reimplemented in Nektar::RoeSolverSIMD.

Definition at line 61 of file CompressibleSolver.cpp.

65 {
66  if (m_pointSolve)
67  {
68  size_t expDim = nDim;
69  NekDouble rhouf{}, rhovf{};
70 
71  if (expDim == 1)
72  {
73  for (size_t i = 0; i < Fwd[0].size(); ++i)
74  {
75  v_PointSolve(Fwd[0][i], Fwd[1][i], 0.0, 0.0, Fwd[2][i],
76  Bwd[0][i], Bwd[1][i], 0.0, 0.0, Bwd[2][i],
77  flux[0][i], flux[1][i], rhouf, rhovf, flux[2][i]);
78  }
79  }
80  else if (expDim == 2)
81  {
82  for (size_t i = 0; i < Fwd[0].size(); ++i)
83  {
84  v_PointSolve(Fwd[0][i], Fwd[1][i], Fwd[2][i], 0.0, Fwd[3][i],
85  Bwd[0][i], Bwd[1][i], Bwd[2][i], 0.0, Bwd[3][i],
86  flux[0][i], flux[1][i], flux[2][i], rhovf,
87  flux[3][i]);
88  }
89  }
90  else if (expDim == 3)
91  {
92  for (size_t i = 0; i < Fwd[0].size(); ++i)
93  {
94  v_PointSolve(Fwd[0][i], Fwd[1][i], Fwd[2][i], Fwd[3][i],
95  Fwd[4][i], Bwd[0][i], Bwd[1][i], Bwd[2][i],
96  Bwd[3][i], Bwd[4][i], flux[0][i], flux[1][i],
97  flux[2][i], flux[3][i], flux[4][i]);
98  }
99  }
100  }
101  else
102  {
103  v_ArraySolve(Fwd, Bwd, flux);
104  }
105 }
virtual void v_PointSolve(ND rhoL, ND rhouL, ND rhovL, ND rhowL, ND EL, ND rhoR, ND rhouR, ND rhovR, ND rhowR, ND ER, ND &rhof, ND &rhouf, ND &rhovf, ND &rhowf, ND &Ef)
virtual void v_ArraySolve(const Array< OneD, const Array< OneD, ND >> &Fwd, const Array< OneD, const Array< OneD, ND >> &Bwd, Array< OneD, Array< OneD, ND >> &flux)

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

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().