36#ifndef NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_MISC_VARCONVERT_H
37#define NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_MISC_VARCONVERT_H
46class VariableConverter;
68 template <
class T,
typename =
typename std::enable_if<
69 std::is_floating_point<T>::value ||
74 T oneOrho = 1.0 / physfield[0];
81 dynEne = 0.5 * dynEne * oneOrho;
85 return energy * oneOrho;
97 template <
class T,
typename =
typename std::enable_if<
98 std::is_floating_point<T>::value ||
107 return mu_star * ratio *
sqrt(ratio) * onePlusC /
119 template <
class T,
typename =
typename std::enable_if<
120 std::is_floating_point<T>::value ||
125 return m_eos->GetTemperature(physfield[0], energy);
130 template <
class T,
typename =
typename std::enable_if<
131 std::is_floating_point<T>::value ||
136 return m_eos->GetPressure(physfield[0], energy);
void ApplyC0Smooth(Array< OneD, NekDouble > &field)
Make field C0.
Array< OneD, NekDouble > & GetAvTrace()
void GetInternalEnergy(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &energy)
Compute the specific internal energy .
EquationOfStateSharedPtr m_eos
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)
Array< OneD, NekDouble > m_hOverP
h/p scaling
void GetRhoFromPT(const Array< OneD, NekDouble > &pressure, const Array< OneD, NekDouble > &temperature, Array< OneD, NekDouble > &rho)
Compute using the equation of state.
void GetAbsoluteVelocity(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &Vtot)
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.
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)
T GetTemperature(T *physfield)
NekDouble m_TRatioSutherland
Array< OneD, NekDouble > m_muAvTrace
T GetInternalEnergy(T *physfield)
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 ,.
std::string m_shockSensorType
std::string m_shockCaptureType
void GetEFromRhoP(const Array< OneD, NekDouble > &rho, const Array< OneD, NekDouble > &pressure, Array< OneD, NekDouble > &energy)
Compute using the equation of state.
void SetElmtMinHP(const Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
Compute an estimate of minimum h/p for each element of the expansion.
MultiRegions::ContFieldSharedPtr m_C0ProjectExp
T GetPressure(T *physfield)
void GetPressure(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &pressure)
Calculate the pressure using the equation of state.
LibUtilities::SessionReaderSharedPtr m_session
NekDouble m_mu0
Shock sensor.
void GetTemperature(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &temperature)
Compute the temperature using the equation of state.
NekDouble m_oneOverT_star
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.
void GetDynamicEnergy(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &energy)
Compute the dynamic energy .
void GetSoundSpeed(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &soundspeed)
Compute the sound speed using the equation of state.
VariableConverter(const LibUtilities::SessionReaderSharedPtr &pSession, const int spaceDim, const SpatialDomains::MeshGraphSharedPtr &pGraph=nullptr)
void GetEntropy(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &entropy)
Compute the entropy using the equation of state.
void GetMach(Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &soundspeed, Array< OneD, NekDouble > &mach)
Compute the mach number .
bool GetFlagCalcDivCurl(void) const
const EquationOfStateSharedPtr Geteos()
void GetVelocityVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)
Compute the velocity field given the momentum .
std::string m_ducrosSensor
~VariableConverter()
Destructor for VariableConverter class.
void GetEnthalpy(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &enthalpy)
Compute the specific enthalpy .
Array< OneD, NekDouble > & GetElmtMinHP()
Array< OneD, NekDouble > m_muAv
storage
Array< OneD, NekDouble > & GetAv()
T GetDynamicViscosity(T &temperature)
void GetDynamicViscosity(const Array< OneD, const NekDouble > &temperature, Array< OneD, NekDouble > &mu)
Compute the dynamic viscosity using the Sutherland's law , C : 110. /Tref Tref : the reference temper...
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< ContField > ContFieldSharedPtr
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
std::vector< double > d(NPUPPER *NPUPPER)
std::shared_ptr< VariableConverter > VariableConverterSharedPtr
std::shared_ptr< EquationOfState > EquationOfStateSharedPtr
A shared pointer to an equation of state object.
static Array< OneD, NekDouble > NullNekDouble1DArray
scalarT< T > sqrt(scalarT< T > in)