Nektar++
|
#include <VariableConverter.h>
Public Member Functions | |
VariableConverter (const LibUtilities::SessionReaderSharedPtr &pSession, const int spaceDim) | |
~VariableConverter () | |
Destructor for VariableConverter class. More... | |
void | GetDynamicEnergy (const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &energy) |
Compute the dynamic energy \( e = rho*V^2/2 \). More... | |
void | GetInternalEnergy (const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &energy) |
Compute the specific internal energy \( e = (E - rho*V^2/2)/rho \). More... | |
template<class T , typename = typename std::enable_if < std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value >::type> | |
T | GetInternalEnergy (T *physfield) |
void | GetEnthalpy (const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &enthalpy) |
Compute the specific enthalpy \( h = e + p/rho \). More... | |
void | GetVelocityVector (const Array< OneD, Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, NekDouble >> &velocity) |
Compute the velocity field \( \mathbf{v} \) given the momentum \( \rho\mathbf{v} \). More... | |
void | GetMach (Array< OneD, Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &soundspeed, Array< OneD, NekDouble > &mach) |
Compute the mach number \( M = \| \mathbf{v} \|^2 / c \). More... | |
void | GetDynamicViscosity (const Array< OneD, const NekDouble > &temperature, Array< OneD, NekDouble > &mu) |
Compute the dynamic viscosity using the Sutherland's law \( \mu = \mu_star * (T / T_star)^3/2 * (1 + C) / (T / T_star + C) \), where: \mu_star = 1.7894 * 10^-5 Kg / (m * s) T_star = 288.15 K C = 110. / 288.15. More... | |
template<class T , typename = typename std::enable_if < std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value >::type> | |
T | GetDynamicViscosity (T &temperature) |
void | GetAbsoluteVelocity (const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &Vtot) |
void | GetSensor (const MultiRegions::ExpListSharedPtr &field, const Array< OneD, const Array< OneD, NekDouble >> &physarray, Array< OneD, NekDouble > &Sensor, Array< OneD, NekDouble > &SensorKappa, int offset=1) |
void | GetTemperature (const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &temperature) |
Compute the temperature using the equation of state. More... | |
template<class T , typename = typename std::enable_if < std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value >::type> | |
T | GetTemperature (T *physfield) |
void | GetPressure (const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &pressure) |
Calculate the pressure using the equation of state. More... | |
template<class T , typename = typename std::enable_if < std::is_floating_point<T>::value || tinysimd::is_vector_floating_point<T>::value >::type> | |
T | GetPressure (T *physfield) |
void | GetSoundSpeed (const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &soundspeed) |
Compute the sound speed using the equation of state. More... | |
void | GetEntropy (const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &entropy) |
Compute the entropy using the equation of state. More... | |
void | GetEFromRhoP (const Array< OneD, NekDouble > &rho, const Array< OneD, NekDouble > &pressure, Array< OneD, NekDouble > &energy) |
Compute \( e(rho,p) \) using the equation of state. More... | |
void | GetRhoFromPT (const Array< OneD, NekDouble > &pressure, const Array< OneD, NekDouble > &temperature, Array< OneD, NekDouble > &rho) |
Compute \( rho(p,T) \) using the equation of state. More... | |
void | GetDmuDT (const Array< OneD, const NekDouble > &temperature, const Array< OneD, const NekDouble > &mu, Array< OneD, NekDouble > &DmuDT) |
Compute the dynamic viscosity using the Sutherland's law \( \mu = \mu_star * (T / T_star)^3/2 * (T_star + 110) / (T + 110) \), where: \mu_star = 1.7894 * 10^-5 Kg / (m * s) T_star = 288.15 K. More... | |
const EquationOfStateSharedPtr | Geteos () |
Definition at line 50 of file VariableConverter.h.
Nektar::VariableConverter::VariableConverter | ( | const LibUtilities::SessionReaderSharedPtr & | pSession, |
const int | spaceDim | ||
) |
Definition at line 44 of file VariableConverter.cpp.
References Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::GetEquationOfStateFactory(), m_eos, m_gasConstant, m_Kappa, m_mu, m_oneOverT_star, m_pInf, m_rhoInf, m_session, and m_Skappa.
Nektar::VariableConverter::~VariableConverter | ( | ) |
Destructor for VariableConverter class.
Definition at line 70 of file VariableConverter.cpp.
void Nektar::VariableConverter::GetAbsoluteVelocity | ( | const Array< OneD, const Array< OneD, NekDouble >> & | physfield, |
Array< OneD, NekDouble > & | Vtot | ||
) |
Definition at line 231 of file VariableConverter.cpp.
References GetVelocityVector(), m_spacedim, Vmath::Vsqrt(), Vmath::Vvtvp(), and Vmath::Zero().
void Nektar::VariableConverter::GetDmuDT | ( | const Array< OneD, const NekDouble > & | temperature, |
const Array< OneD, const NekDouble > & | mu, | ||
Array< OneD, NekDouble > & | DmuDT | ||
) |
Compute the dynamic viscosity using the Sutherland's law \( \mu = \mu_star * (T / T_star)^3/2 * (T_star + 110) / (T + 110) \), where: \mu_star = 1.7894 * 10^-5 Kg / (m * s) T_star = 288.15 K.
physfield | Input physical field. |
mu | The resulting dynamic viscosity. |
Definition at line 215 of file VariableConverter.cpp.
void Nektar::VariableConverter::GetDynamicEnergy | ( | const Array< OneD, const Array< OneD, NekDouble >> & | physfield, |
Array< OneD, NekDouble > & | energy | ||
) |
Compute the dynamic energy \( e = rho*V^2/2 \).
Definition at line 78 of file VariableConverter.cpp.
References m_spacedim, Vmath::Smul(), Vmath::Vdiv(), Vmath::Vvtvp(), and Vmath::Zero().
Referenced by GetInternalEnergy().
void Nektar::VariableConverter::GetDynamicViscosity | ( | const Array< OneD, const NekDouble > & | temperature, |
Array< OneD, NekDouble > & | mu | ||
) |
Compute the dynamic viscosity using the Sutherland's law \( \mu = \mu_star * (T / T_star)^3/2 * (1 + C) / (T / T_star + C) \), where: \mu_star = 1.7894 * 10^-5 Kg / (m * s) T_star = 288.15 K C = 110. / 288.15.
WARNING, if this routine is modified the same must be done in the FieldConvert utility ProcessWSS.cpp (this class should be restructured).
temperature | Input temperature field. |
mu | The resulting dynamic viscosity. |
Definition at line 195 of file VariableConverter.cpp.
|
inline |
Definition at line 103 of file VariableConverter.h.
References m_mu, m_oneOverT_star, and tinysimd::sqrt().
void Nektar::VariableConverter::GetEFromRhoP | ( | const Array< OneD, NekDouble > & | rho, |
const Array< OneD, NekDouble > & | pressure, | ||
Array< OneD, NekDouble > & | energy | ||
) |
Compute \( e(rho,p) \) using the equation of state.
rho | Input density |
pressure | Input pressure |
energy | The resulting internal energy. |
Definition at line 441 of file VariableConverter.cpp.
References m_eos, and CG_Iterations::pressure.
void Nektar::VariableConverter::GetEnthalpy | ( | const Array< OneD, const Array< OneD, NekDouble >> & | physfield, |
Array< OneD, NekDouble > & | enthalpy | ||
) |
Compute the specific enthalpy \( h = e + p/rho \).
Definition at line 118 of file VariableConverter.cpp.
References GetInternalEnergy(), GetPressure(), CG_Iterations::pressure, Vmath::Vadd(), and Vmath::Vdiv().
void Nektar::VariableConverter::GetEntropy | ( | const Array< OneD, const Array< OneD, NekDouble >> & | physfield, |
Array< OneD, NekDouble > & | entropy | ||
) |
Compute the entropy using the equation of state.
physfield | Input physical field |
soundspeed | The resulting sound speed \( c \). |
Definition at line 419 of file VariableConverter.cpp.
References GetInternalEnergy(), and m_eos.
|
inline |
void Nektar::VariableConverter::GetInternalEnergy | ( | const Array< OneD, const Array< OneD, NekDouble >> & | physfield, |
Array< OneD, NekDouble > & | energy | ||
) |
Compute the specific internal energy \( e = (E - rho*V^2/2)/rho \).
Definition at line 100 of file VariableConverter.cpp.
References GetDynamicEnergy(), m_spacedim, Vmath::Vdiv(), and Vmath::Vsub().
Referenced by GetEnthalpy(), GetEntropy(), GetPressure(), GetSoundSpeed(), and GetTemperature().
|
inline |
void Nektar::VariableConverter::GetMach | ( | Array< OneD, Array< OneD, NekDouble >> & | physfield, |
Array< OneD, NekDouble > & | soundspeed, | ||
Array< OneD, NekDouble > & | mach | ||
) |
Compute the mach number \( M = \| \mathbf{v} \|^2 / c \).
physfield | Input physical field. |
soundfield | The speed of sound corresponding to physfield. |
mach | The resulting mach number \( M \). |
Definition at line 161 of file VariableConverter.cpp.
References m_spacedim, Vmath::Vdiv(), Vmath::Vmul(), Vmath::Vsqrt(), and Vmath::Vvtvp().
void Nektar::VariableConverter::GetPressure | ( | const Array< OneD, const Array< OneD, NekDouble >> & | physfield, |
Array< OneD, NekDouble > & | pressure | ||
) |
Calculate the pressure using the equation of state.
physfield | Input momentum. |
pressure | Computed pressure field. |
Definition at line 356 of file VariableConverter.cpp.
References GetInternalEnergy(), m_eos, and CG_Iterations::pressure.
Referenced by GetEnthalpy().
|
inline |
void Nektar::VariableConverter::GetRhoFromPT | ( | const Array< OneD, NekDouble > & | pressure, |
const Array< OneD, NekDouble > & | temperature, | ||
Array< OneD, NekDouble > & | rho | ||
) |
Compute \( rho(p,T) \) using the equation of state.
pressure | Input pressure |
temperature | Input temperature |
rho | The resulting density |
Definition at line 460 of file VariableConverter.cpp.
References m_eos, and CG_Iterations::pressure.
void Nektar::VariableConverter::GetSensor | ( | const MultiRegions::ExpListSharedPtr & | field, |
const Array< OneD, const Array< OneD, NekDouble >> & | physarray, | ||
Array< OneD, NekDouble > & | Sensor, | ||
Array< OneD, NekDouble > & | SensorKappa, | ||
int | offset = 1 |
||
) |
Definition at line 257 of file VariableConverter.cpp.
References Vmath::Dot(), Vmath::Fill(), Nektar::NekConstants::kNekSqrtTol, m_Kappa, m_Skappa, tinysimd::sqrt(), and Vmath::Vsub().
void Nektar::VariableConverter::GetSoundSpeed | ( | const Array< OneD, const Array< OneD, NekDouble >> & | physfield, |
Array< OneD, NekDouble > & | soundspeed | ||
) |
Compute the sound speed using the equation of state.
physfield | Input physical field |
soundspeed | The resulting sound speed \( c \). |
Definition at line 398 of file VariableConverter.cpp.
References GetInternalEnergy(), and m_eos.
void Nektar::VariableConverter::GetTemperature | ( | const Array< OneD, const Array< OneD, NekDouble >> & | physfield, |
Array< OneD, NekDouble > & | temperature | ||
) |
Compute the temperature using the equation of state.
physfield | Input physical field. |
temperature | The resulting temperature \( T \). |
Definition at line 377 of file VariableConverter.cpp.
References GetInternalEnergy(), and m_eos.
|
inline |
void Nektar::VariableConverter::GetVelocityVector | ( | const Array< OneD, Array< OneD, NekDouble >> & | physfield, |
Array< OneD, Array< OneD, NekDouble >> & | velocity | ||
) |
Compute the velocity field \( \mathbf{v} \) given the momentum \( \rho\mathbf{v} \).
physfield | Momentum field. |
velocity | Velocity field. |
Definition at line 142 of file VariableConverter.cpp.
References m_spacedim, and Vmath::Vdiv().
Referenced by GetAbsoluteVelocity().
|
protected |
Definition at line 175 of file VariableConverter.h.
Referenced by GetEFromRhoP(), GetEntropy(), Geteos(), GetPressure(), GetRhoFromPT(), GetSoundSpeed(), GetTemperature(), and VariableConverter().
|
protected |
Definition at line 179 of file VariableConverter.h.
Referenced by VariableConverter().
|
protected |
Definition at line 182 of file VariableConverter.h.
Referenced by GetSensor(), and VariableConverter().
|
protected |
Definition at line 180 of file VariableConverter.h.
Referenced by GetDynamicViscosity(), and VariableConverter().
|
protected |
Definition at line 183 of file VariableConverter.h.
Referenced by GetDynamicViscosity(), and VariableConverter().
|
protected |
Definition at line 177 of file VariableConverter.h.
Referenced by VariableConverter().
|
protected |
Definition at line 178 of file VariableConverter.h.
Referenced by VariableConverter().
|
protected |
Definition at line 174 of file VariableConverter.h.
Referenced by VariableConverter().
|
protected |
Definition at line 181 of file VariableConverter.h.
Referenced by GetSensor(), and VariableConverter().
|
protected |
Definition at line 176 of file VariableConverter.h.
Referenced by GetAbsoluteVelocity(), GetDynamicEnergy(), GetInternalEnergy(), GetMach(), and GetVelocityVector().