Nektar++
|
#include <VariableConverter.h>
Public Member Functions | |
VariableConverter (const LibUtilities::SessionReaderSharedPtr &pSession, const int spaceDim, const SpatialDomains::MeshGraphSharedPtr &pGraph=nullptr) | |
~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) \), C : 110. /Tref Tref : the reference temperature, Tref, should always given in Kelvin, if non-dimensional should be the reference for non-dimensionalizing muref : the dynamic viscosity or the 1/Re corresponding to Tref T_star : m_pInf / (m_rhoInf * m_gasConstant),non-dimensional or dimensional. 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 | 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) \),. More... | |
const EquationOfStateSharedPtr | Geteos () |
void | SetAv (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const Array< OneD, NekDouble > > &consVar, const Array< OneD, NekDouble > &div=NullNekDouble1DArray, const Array< OneD, NekDouble > &curlSquared=NullNekDouble1DArray) |
Array< OneD, NekDouble > & | GetAv () |
Array< OneD, NekDouble > & | GetAvTrace () |
bool | GetFlagCalcDivCurl (void) const |
void | SetElmtMinHP (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields) |
Compute an estimate of minimum h/p for each element of the expansion. More... | |
Array< OneD, NekDouble > & | GetElmtMinHP () |
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 | GetMuAv (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &muAv) |
Calculate the physical artificial viscosity based on modal sensor. More... | |
void | GetMuAv (const Array< OneD, MultiRegions::ExpListSharedPtr > &fields, const Array< OneD, const Array< OneD, NekDouble > > &consVar, const Array< OneD, NekDouble > &div, Array< OneD, NekDouble > &muAv) |
Calculate the physical artificial viscosity based on dilatation of velocity vector. More... | |
void | ApplyDucros (const Array< OneD, NekDouble > &div, const Array< OneD, NekDouble > &curlSquare, Array< OneD, NekDouble > &muAv) |
Apply Ducros (anti-vorticity) sensor averaged over the element. More... | |
void | ApplyC0Smooth (Array< OneD, NekDouble > &field) |
Make field C0. More... | |
Protected Attributes | |
LibUtilities::SessionReaderSharedPtr | m_session |
EquationOfStateSharedPtr | m_eos |
size_t | m_spacedim |
NekDouble | m_pInf |
NekDouble | m_rhoInf |
NekDouble | m_gasConstant |
NekDouble | m_mu |
NekDouble | m_Skappa |
NekDouble | m_Kappa |
NekDouble | m_oneOverT_star |
NekDouble | m_Tref |
NekDouble | m_TRatioSutherland |
NekDouble | m_mu0 |
Shock sensor. More... | |
std::string | m_shockCaptureType |
std::string | m_shockSensorType |
std::string | m_ducrosSensor |
std::string | m_smoothing |
MultiRegions::ContFieldSharedPtr | m_C0ProjectExp = nullptr |
Array< OneD, NekDouble > | m_hOverP |
h/p scaling More... | |
Array< OneD, NekDouble > | m_muAv |
storage More... | |
Array< OneD, NekDouble > | m_muAvTrace |
bool | m_flagCalcDivCurl = false |
Definition at line 51 of file VariableConverter.h.
Nektar::VariableConverter::VariableConverter | ( | const LibUtilities::SessionReaderSharedPtr & | pSession, |
const int | spaceDim, | ||
const SpatialDomains::MeshGraphSharedPtr & | pGraph = nullptr |
||
) |
Definition at line 47 of file VariableConverter.cpp.
References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), Nektar::LibUtilities::NekFactory< tKey, tBase, tParam >::CreateInstance(), Nektar::GetEquationOfStateFactory(), m_C0ProjectExp, m_ducrosSensor, m_eos, m_flagCalcDivCurl, m_gasConstant, m_hOverP, m_Kappa, m_mu, m_mu0, m_muAv, m_muAvTrace, m_oneOverT_star, m_pInf, m_rhoInf, m_session, m_shockCaptureType, m_shockSensorType, m_Skappa, m_smoothing, m_TRatioSutherland, m_Tref, Nektar::NullNekDouble1DArray, and WARNINGL0.
Nektar::VariableConverter::~VariableConverter | ( | ) |
Destructor for VariableConverter class.
Definition at line 124 of file VariableConverter.cpp.
Make field C0.
field | Input Field |
Definition at line 782 of file VariableConverter.cpp.
References ASSERTL0, FilterPython_Function::field, m_C0ProjectExp, Vmath::Vcopy(), and Vmath::Vdiv().
Referenced by SetAv().
void Nektar::VariableConverter::ApplyDucros | ( | const Array< OneD, NekDouble > & | div, |
const Array< OneD, NekDouble > & | curlSquare, | ||
Array< OneD, NekDouble > & | muAv | ||
) |
Apply Ducros (anti-vorticity) sensor averaged over the element.
field | Input Field |
Definition at line 756 of file VariableConverter.cpp.
References CellMLToNektar.cellml_metadata::p.
Referenced by SetAv().
void Nektar::VariableConverter::GetAbsoluteVelocity | ( | const Array< OneD, const Array< OneD, NekDouble > > & | physfield, |
Array< OneD, NekDouble > & | Vtot | ||
) |
Definition at line 283 of file VariableConverter.cpp.
References GetVelocityVector(), m_spacedim, Vmath::Vsqrt(), Vmath::Vvtvp(), and Vmath::Zero().
Referenced by GetMuAv().
Definition at line 478 of file VariableConverter.cpp.
References ASSERTL1, m_muAv, and Nektar::NullNekDouble1DArray.
Definition at line 484 of file VariableConverter.cpp.
References ASSERTL1, m_muAvTrace, and Nektar::NullNekDouble1DArray.
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) \),.
Definition at line 266 of file VariableConverter.cpp.
References m_oneOverT_star, and m_TRatioSutherland.
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 132 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) \), C : 110. /Tref Tref : the reference temperature, Tref, should always given in Kelvin, if non-dimensional should be the reference for non-dimensionalizing muref : the dynamic viscosity or the 1/Re corresponding to Tref T_star : m_pInf / (m_rhoInf * m_gasConstant),non-dimensional or dimensional.
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. |
mu | The resulting dynamic viscosity. |
Definition at line 251 of file VariableConverter.cpp.
References GetDynamicViscosity().
Referenced by GetDynamicViscosity().
|
inline |
Definition at line 100 of file VariableConverter.h.
References m_mu, m_oneOverT_star, m_TRatioSutherland, 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 400 of file VariableConverter.cpp.
References m_eos, and CG_Iterations::pressure.
Definition at line 557 of file VariableConverter.cpp.
References ASSERTL1, m_hOverP, and Nektar::NullNekDouble1DArray.
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 172 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 378 of file VariableConverter.cpp.
References GetInternalEnergy(), and m_eos.
|
inline |
|
inline |
Definition at line 170 of file VariableConverter.h.
References m_flagCalcDivCurl.
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 154 of file VariableConverter.cpp.
References GetDynamicEnergy(), m_spacedim, Vmath::Vdiv(), and Vmath::Vsub().
Referenced by GetEnthalpy(), GetEntropy(), GetPressure(), GetSoundSpeed(), and GetTemperature().
|
inline |
Definition at line 71 of file VariableConverter.h.
References Nektar::UnitTests::d(), and m_spacedim.
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 215 of file VariableConverter.cpp.
References m_spacedim, Vmath::Vdiv(), Vmath::Vmul(), Vmath::Vsqrt(), and Vmath::Vvtvp().
void Nektar::VariableConverter::GetMuAv | ( | const Array< OneD, MultiRegions::ExpListSharedPtr > & | fields, |
const Array< OneD, const Array< OneD, NekDouble > > & | consVar, | ||
const Array< OneD, NekDouble > & | div, | ||
Array< OneD, NekDouble > & | muAv | ||
) |
Calculate the physical artificial viscosity based on dilatation of velocity vector.
Definition at line 707 of file VariableConverter.cpp.
References GetAbsoluteVelocity(), GetSoundSpeed(), m_hOverP, m_mu0, and CellMLToNektar.cellml_metadata::p.
void Nektar::VariableConverter::GetMuAv | ( | const Array< OneD, MultiRegions::ExpListSharedPtr > & | fields, |
const Array< OneD, const Array< OneD, NekDouble > > & | consVar, | ||
Array< OneD, NekDouble > & | muAv | ||
) |
Calculate the physical artificial viscosity based on modal sensor.
consVar | Input field. |
Definition at line 659 of file VariableConverter.cpp.
References GetAbsoluteVelocity(), GetSensor(), GetSoundSpeed(), m_hOverP, m_mu0, Smath::Smax(), Vmath::Smul(), Vmath::Vadd(), Vmath::Vmax(), and Vmath::Vsum().
Referenced by SetAv().
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 315 of file VariableConverter.cpp.
References GetInternalEnergy(), m_eos, and CG_Iterations::pressure.
Referenced by GetEnthalpy().
|
inline |
Definition at line 133 of file VariableConverter.h.
References GetInternalEnergy(), and m_eos.
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 419 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 563 of file VariableConverter.cpp.
References Vmath::Dot(), FilterPython_Function::field, Vmath::Fill(), Nektar::NekConstants::kNekMachineEpsilon, m_Kappa, m_Skappa, tinysimd::sqrt(), and Vmath::Vsub().
Referenced by GetMuAv().
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 357 of file VariableConverter.cpp.
References GetInternalEnergy(), and m_eos.
Referenced by GetMuAv().
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 336 of file VariableConverter.cpp.
References GetInternalEnergy(), and m_eos.
|
inline |
Definition at line 122 of file VariableConverter.h.
References GetInternalEnergy(), and m_eos.
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 196 of file VariableConverter.cpp.
References m_spacedim, and Vmath::Vdiv().
Referenced by GetAbsoluteVelocity().
void Nektar::VariableConverter::SetAv | ( | const Array< OneD, MultiRegions::ExpListSharedPtr > & | fields, |
const Array< OneD, const Array< OneD, NekDouble > > & | consVar, | ||
const Array< OneD, NekDouble > & | div = NullNekDouble1DArray , |
||
const Array< OneD, NekDouble > & | curlSquared = NullNekDouble1DArray |
||
) |
Definition at line 431 of file VariableConverter.cpp.
References ApplyC0Smooth(), ApplyDucros(), GetMuAv(), m_ducrosSensor, m_muAv, m_muAvTrace, m_shockSensorType, m_smoothing, Nektar::NullNekDouble1DArray, CellMLToNektar.cellml_metadata::p, and SetElmtMinHP().
void Nektar::VariableConverter::SetElmtMinHP | ( | const Array< OneD, MultiRegions::ExpListSharedPtr > & | fields | ) |
Compute an estimate of minimum h/p for each element of the expansion.
Definition at line 494 of file VariableConverter.cpp.
References ASSERTL0, Nektar::LocalRegions::Expansion1D::GetGeom1D(), m_hOverP, and Nektar::NullNekDouble1DArray.
Referenced by SetAv().
|
protected |
Definition at line 220 of file VariableConverter.h.
Referenced by ApplyC0Smooth(), and VariableConverter().
|
protected |
Definition at line 218 of file VariableConverter.h.
Referenced by SetAv(), and VariableConverter().
|
protected |
Definition at line 202 of file VariableConverter.h.
Referenced by GetEFromRhoP(), GetEntropy(), Geteos(), GetPressure(), GetRhoFromPT(), GetSoundSpeed(), GetTemperature(), and VariableConverter().
|
protected |
Definition at line 227 of file VariableConverter.h.
Referenced by GetFlagCalcDivCurl(), and VariableConverter().
|
protected |
Definition at line 206 of file VariableConverter.h.
Referenced by VariableConverter().
h/p scaling
Definition at line 223 of file VariableConverter.h.
Referenced by GetElmtMinHP(), GetMuAv(), SetElmtMinHP(), and VariableConverter().
|
protected |
Definition at line 209 of file VariableConverter.h.
Referenced by GetSensor(), and VariableConverter().
|
protected |
Definition at line 207 of file VariableConverter.h.
Referenced by GetDynamicViscosity(), and VariableConverter().
|
protected |
Shock sensor.
Definition at line 215 of file VariableConverter.h.
Referenced by GetMuAv(), and VariableConverter().
storage
Definition at line 225 of file VariableConverter.h.
Referenced by GetAv(), SetAv(), and VariableConverter().
Definition at line 226 of file VariableConverter.h.
Referenced by GetAvTrace(), SetAv(), and VariableConverter().
|
protected |
Definition at line 210 of file VariableConverter.h.
Referenced by GetDmuDT(), GetDynamicViscosity(), and VariableConverter().
|
protected |
Definition at line 204 of file VariableConverter.h.
Referenced by VariableConverter().
|
protected |
Definition at line 205 of file VariableConverter.h.
Referenced by VariableConverter().
|
protected |
Definition at line 201 of file VariableConverter.h.
Referenced by VariableConverter().
|
protected |
Definition at line 216 of file VariableConverter.h.
Referenced by VariableConverter().
|
protected |
Definition at line 217 of file VariableConverter.h.
Referenced by SetAv(), and VariableConverter().
|
protected |
Definition at line 208 of file VariableConverter.h.
Referenced by GetSensor(), and VariableConverter().
|
protected |
Definition at line 219 of file VariableConverter.h.
Referenced by SetAv(), and VariableConverter().
|
protected |
Definition at line 203 of file VariableConverter.h.
Referenced by GetAbsoluteVelocity(), GetDynamicEnergy(), GetInternalEnergy(), GetMach(), and GetVelocityVector().
|
protected |
Definition at line 212 of file VariableConverter.h.
Referenced by GetDmuDT(), GetDynamicViscosity(), and VariableConverter().
|
protected |
Definition at line 211 of file VariableConverter.h.
Referenced by VariableConverter().