Nektar++
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
Nektar::EquationOfState Class Referenceabstract

Encapsulates equations of state allowing us to obtain thermodynamic properties: most relations are in the form X(rho,e) More...

#include <EquationOfState.h>

Inheritance diagram for Nektar::EquationOfState:
[legend]

Public Member Functions

virtual ~EquationOfState ()
 
NekDouble GetTemperature (const NekDouble &rho, const NekDouble &e)
 Calculate the temperature. More...
 
vec_t GetTemperature (const vec_t &rho, const vec_t &e)
 
NekDouble GetPressure (const NekDouble &rho, const NekDouble &e)
 Calculate the pressure. More...
 
vec_t GetPressure (const vec_t &rho, const vec_t &e)
 
NekDouble GetSoundSpeed (const NekDouble &rho, const NekDouble &e)
 Calculate the sound speed. More...
 
NekDouble GetEntropy (const NekDouble &rho, const NekDouble &e)
 Calculate the entropy. More...
 
NekDouble GetDPDrho_e (const NekDouble &rho, const NekDouble &e)
 Calculate the partial derivative of P(rho,e) with respect to rho. More...
 
NekDouble GetDPDe_rho (const NekDouble &rho, const NekDouble &e)
 Calculate the partial derivative of P(rho,e) with respect to e. More...
 
NekDouble GetEFromRhoP (const NekDouble &rho, const NekDouble &p)
 Obtain the internal energy from rho and P. More...
 
NekDouble GetRhoFromPT (const NekDouble &p, const NekDouble &T)
 Obtain the density from P and T. More...
 

Protected Member Functions

 EquationOfState (const LibUtilities::SessionReaderSharedPtr &pSession)
 Constructor. More...
 
 EquationOfState (const NekDouble &gamma, const NekDouble &gasConstant)
 Programmatic Constructor. More...
 
virtual NekDouble v_GetTemperature (const NekDouble &rho, const NekDouble &e)=0
 
virtual vec_t v_GetTemperature (const vec_t &rho, const vec_t &e)=0
 
virtual NekDouble v_GetPressure (const NekDouble &rho, const NekDouble &e)=0
 
virtual vec_t v_GetPressure (const vec_t &rho, const vec_t &e)=0
 
virtual NekDouble v_GetSoundSpeed (const NekDouble &rho, const NekDouble &e)
 
virtual NekDouble v_GetEntropy (const NekDouble &rho, const NekDouble &e)=0
 
virtual NekDouble v_GetDPDrho_e (const NekDouble &rho, const NekDouble &e)=0
 
virtual NekDouble v_GetDPDe_rho (const NekDouble &rho, const NekDouble &e)=0
 
virtual NekDouble v_GetEFromRhoP (const NekDouble &rho, const NekDouble &p)=0
 
virtual NekDouble v_GetRhoFromPT (const NekDouble &rho, const NekDouble &p)=0
 

Protected Attributes

NekDouble m_gamma
 
NekDouble m_gasConstant
 
NekDouble m_gammaMone
 
NekDouble m_gammaMoneOgasConst
 

Detailed Description

Encapsulates equations of state allowing us to obtain thermodynamic properties: most relations are in the form X(rho,e)

Definition at line 66 of file EquationOfState.h.

Constructor & Destructor Documentation

◆ ~EquationOfState()

virtual Nektar::EquationOfState::~EquationOfState ( )
inlinevirtual

Definition at line 69 of file EquationOfState.h.

70 {
71 }

◆ EquationOfState() [1/2]

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

Constructor.

Definition at line 47 of file EquationOfState.cpp.

49{
50 pSession->LoadParameter("Gamma", m_gamma, 1.4);
51 pSession->LoadParameter("GasConstant", m_gasConstant, 287.058);
52
53 m_gammaMone = m_gamma - 1.0;
55}

References m_gamma, m_gammaMone, m_gammaMoneOgasConst, and m_gasConstant.

◆ EquationOfState() [2/2]

Nektar::EquationOfState::EquationOfState ( const NekDouble gamma,
const NekDouble gasConstant 
)
protected

Programmatic Constructor.

Definition at line 57 of file EquationOfState.cpp.

59 : m_gamma{gamma}, m_gasConstant{gasConstant}
60{
61}

Member Function Documentation

◆ GetDPDe_rho()

NekDouble Nektar::EquationOfState::GetDPDe_rho ( const NekDouble rho,
const NekDouble e 
)

Calculate the partial derivative of P(rho,e) with respect to e.

Definition at line 100 of file EquationOfState.cpp.

101{
102 return v_GetDPDe_rho(rho, e);
103}
virtual NekDouble v_GetDPDe_rho(const NekDouble &rho, const NekDouble &e)=0

References v_GetDPDe_rho().

Referenced by Nektar::PengRobinsonEoS::v_GetDPDrho_e(), Nektar::RedlichKwongEoS::v_GetDPDrho_e(), and v_GetSoundSpeed().

◆ GetDPDrho_e()

NekDouble Nektar::EquationOfState::GetDPDrho_e ( const NekDouble rho,
const NekDouble e 
)

Calculate the partial derivative of P(rho,e) with respect to rho.

Definition at line 95 of file EquationOfState.cpp.

96{
97 return v_GetDPDrho_e(rho, e);
98}
virtual NekDouble v_GetDPDrho_e(const NekDouble &rho, const NekDouble &e)=0

References v_GetDPDrho_e().

Referenced by v_GetSoundSpeed().

◆ GetEFromRhoP()

NekDouble Nektar::EquationOfState::GetEFromRhoP ( const NekDouble rho,
const NekDouble p 
)

Obtain the internal energy from rho and P.

Definition at line 105 of file EquationOfState.cpp.

107{
108 return v_GetEFromRhoP(rho, p);
109}
virtual NekDouble v_GetEFromRhoP(const NekDouble &rho, const NekDouble &p)=0

References CellMLToNektar.cellml_metadata::p, and v_GetEFromRhoP().

◆ GetEntropy()

NekDouble Nektar::EquationOfState::GetEntropy ( const NekDouble rho,
const NekDouble e 
)

Calculate the entropy.

Definition at line 90 of file EquationOfState.cpp.

91{
92 return v_GetEntropy(rho, e);
93}
virtual NekDouble v_GetEntropy(const NekDouble &rho, const NekDouble &e)=0

References v_GetEntropy().

◆ GetPressure() [1/2]

NekDouble Nektar::EquationOfState::GetPressure ( const NekDouble rho,
const NekDouble e 
)

Calculate the pressure.

Definition at line 74 of file EquationOfState.cpp.

75{
76 return v_GetPressure(rho, e);
77}
virtual NekDouble v_GetPressure(const NekDouble &rho, const NekDouble &e)=0

References v_GetPressure().

Referenced by v_GetSoundSpeed().

◆ GetPressure() [2/2]

vec_t Nektar::EquationOfState::GetPressure ( const vec_t rho,
const vec_t e 
)

Definition at line 79 of file EquationOfState.cpp.

80{
81 return v_GetPressure(rho, e);
82}

References v_GetPressure().

◆ GetRhoFromPT()

NekDouble Nektar::EquationOfState::GetRhoFromPT ( const NekDouble p,
const NekDouble T 
)

Obtain the density from P and T.

Definition at line 111 of file EquationOfState.cpp.

112{
113 return v_GetRhoFromPT(p, T);
114}
virtual NekDouble v_GetRhoFromPT(const NekDouble &rho, const NekDouble &p)=0

References CellMLToNektar.cellml_metadata::p, and v_GetRhoFromPT().

◆ GetSoundSpeed()

NekDouble Nektar::EquationOfState::GetSoundSpeed ( const NekDouble rho,
const NekDouble e 
)

Calculate the sound speed.

Definition at line 84 of file EquationOfState.cpp.

86{
87 return v_GetSoundSpeed(rho, e);
88}
virtual NekDouble v_GetSoundSpeed(const NekDouble &rho, const NekDouble &e)

References v_GetSoundSpeed().

◆ GetTemperature() [1/2]

NekDouble Nektar::EquationOfState::GetTemperature ( const NekDouble rho,
const NekDouble e 
)

◆ GetTemperature() [2/2]

vec_t Nektar::EquationOfState::GetTemperature ( const vec_t rho,
const vec_t e 
)

Definition at line 69 of file EquationOfState.cpp.

70{
71 return v_GetTemperature(rho, e);
72}

References v_GetTemperature().

◆ v_GetDPDe_rho()

virtual NekDouble Nektar::EquationOfState::v_GetDPDe_rho ( const NekDouble rho,
const NekDouble e 
)
protectedpure virtual

◆ v_GetDPDrho_e()

virtual NekDouble Nektar::EquationOfState::v_GetDPDrho_e ( const NekDouble rho,
const NekDouble e 
)
protectedpure virtual

◆ v_GetEFromRhoP()

virtual NekDouble Nektar::EquationOfState::v_GetEFromRhoP ( const NekDouble rho,
const NekDouble p 
)
protectedpure virtual

◆ v_GetEntropy()

virtual NekDouble Nektar::EquationOfState::v_GetEntropy ( const NekDouble rho,
const NekDouble e 
)
protectedpure virtual

◆ v_GetPressure() [1/2]

virtual NekDouble Nektar::EquationOfState::v_GetPressure ( const NekDouble rho,
const NekDouble e 
)
protectedpure virtual

◆ v_GetPressure() [2/2]

virtual vec_t Nektar::EquationOfState::v_GetPressure ( const vec_t rho,
const vec_t e 
)
protectedpure virtual

◆ v_GetRhoFromPT()

virtual NekDouble Nektar::EquationOfState::v_GetRhoFromPT ( const NekDouble rho,
const NekDouble p 
)
protectedpure virtual

◆ v_GetSoundSpeed()

NekDouble Nektar::EquationOfState::v_GetSoundSpeed ( const NekDouble rho,
const NekDouble e 
)
protectedvirtual

Reimplemented in Nektar::IdealGasEoS.

Definition at line 118 of file EquationOfState.cpp.

120{
121 NekDouble p = GetPressure(rho, e);
122 NekDouble dpde = GetDPDe_rho(rho, e);
123 NekDouble dpdrho = GetDPDrho_e(rho, e);
124
125 NekDouble enthalpy = e + p / rho;
126
127 NekDouble chi = dpdrho - e / rho * dpde;
128 NekDouble kappa = dpde / rho;
129
130 return sqrt(chi + kappa * enthalpy);
131}
NekDouble GetPressure(const NekDouble &rho, const NekDouble &e)
Calculate the pressure.
NekDouble GetDPDrho_e(const NekDouble &rho, const NekDouble &e)
Calculate the partial derivative of P(rho,e) with respect to rho.
NekDouble GetDPDe_rho(const NekDouble &rho, const NekDouble &e)
Calculate the partial derivative of P(rho,e) with respect to e.
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References GetDPDe_rho(), GetDPDrho_e(), GetPressure(), CellMLToNektar.cellml_metadata::p, and tinysimd::sqrt().

Referenced by GetSoundSpeed().

◆ v_GetTemperature() [1/2]

virtual NekDouble Nektar::EquationOfState::v_GetTemperature ( const NekDouble rho,
const NekDouble e 
)
protectedpure virtual

◆ v_GetTemperature() [2/2]

virtual vec_t Nektar::EquationOfState::v_GetTemperature ( const vec_t rho,
const vec_t e 
)
protectedpure virtual

Member Data Documentation

◆ m_gamma

NekDouble Nektar::EquationOfState::m_gamma
protected

◆ m_gammaMone

NekDouble Nektar::EquationOfState::m_gammaMone
protected

Definition at line 104 of file EquationOfState.h.

Referenced by EquationOfState(), and Nektar::IdealGasEoS::GetPressureKernel().

◆ m_gammaMoneOgasConst

NekDouble Nektar::EquationOfState::m_gammaMoneOgasConst
protected

◆ m_gasConstant

NekDouble Nektar::EquationOfState::m_gasConstant
protected