Nektar++
VariableConverter.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File VariableConverter.h
4 //
5 // For more information, please see: http://www.nektar.info
6 //
7 // The MIT License
8 //
9 // Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10 // Department of Aeronautics, Imperial College London (UK), and Scientific
11 // Computing and Imaging Institute, University of Utah (USA).
12 //
13 // Permission is hereby granted, free of charge, to any person obtaining a
14 // copy of this software and associated documentation files (the "Software"),
15 // to deal in the Software without restriction, including without limitation
16 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
17 // and/or sell copies of the Software, and to permit persons to whom the
18 // Software is furnished to do so, subject to the following conditions:
19 //
20 // The above copyright notice and this permission notice shall be included
21 // in all copies or substantial portions of the Software.
22 //
23 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29 // DEALINGS IN THE SOFTWARE.
30 //
31 // Description: Auxiliary functions to convert variables in
32 // the compressible flow system
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #ifndef NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_MISC_VARCONVERT_H
37 #define NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_MISC_VARCONVERT_H
38 
39 #include "EquationOfState.h"
41 
42 namespace Nektar
43 {
44 // Forward declarations
46 typedef std::shared_ptr<VariableConverter> VariableConverterSharedPtr;
47 /**
48  *
49  */
51 {
52 public:
54  const int spaceDim);
55 
57 
58  // Variable manipulations valid for all fluids
59  void GetInternalEnergy(
60  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
61  Array<OneD, NekDouble> &energy);
62  void GetEnthalpy(const Array<OneD, const Array<OneD, NekDouble>> &physfield,
63  Array<OneD, NekDouble> &enthalpy);
64  void GetVelocityVector(const Array<OneD, Array<OneD, NekDouble>> &physfield,
65  Array<OneD, Array<OneD, NekDouble>> &velocity);
66  void GetMach(Array<OneD, Array<OneD, NekDouble>> &physfield,
67  Array<OneD, NekDouble> &soundspeed,
69  void GetDynamicViscosity(const Array<OneD, const NekDouble> &temperature,
72  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
74  void GetSensor(const MultiRegions::ExpListSharedPtr &field,
75  const Array<OneD, const Array<OneD, NekDouble>> &physarray,
76  Array<OneD, NekDouble> &Sensor,
77  Array<OneD, NekDouble> &SensorKappa, int offset = 1);
78 
79  // Transformations depending on the equation of state
80  void GetTemperature(
81  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
82  Array<OneD, NekDouble> &temperature);
83  void GetPressure(const Array<OneD, const Array<OneD, NekDouble>> &physfield,
85  void GetSoundSpeed(
86  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
87  Array<OneD, NekDouble> &soundspeed);
88  void GetEntropy(const Array<OneD, const Array<OneD, NekDouble>> &physfield,
89  Array<OneD, NekDouble> &entropy);
90  void GetEFromRhoP(const Array<OneD, NekDouble> &rho,
91  const Array<OneD, NekDouble> &pressure,
92  Array<OneD, NekDouble> &energy);
93  void GetRhoFromPT(const Array<OneD, NekDouble> &pressure,
94  const Array<OneD, NekDouble> &temperature,
96 
97 protected:
107 };
108 }
109 #endif
std::shared_ptr< VariableConverter > VariableConverterSharedPtr
EquationOfStateSharedPtr m_eos
VariableConverter(const LibUtilities::SessionReaderSharedPtr &pSession, const int spaceDim)
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
void GetEFromRhoP(const Array< OneD, NekDouble > &rho, const Array< OneD, NekDouble > &pressure, Array< OneD, NekDouble > &energy)
Compute using the equation of state.
void GetPressure(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &pressure)
Calculate the pressure using the equation of state.
void GetEntropy(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &entropy)
Compute the entropy using the equation of state.
void GetAbsoluteVelocity(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &Vtot)
LibUtilities::SessionReaderSharedPtr m_session
void GetEnthalpy(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &enthalpy)
Compute the specific enthalpy .
void GetDynamicViscosity(const Array< OneD, const NekDouble > &temperature, Array< OneD, NekDouble > &mu)
Compute the dynamic viscosity using the Sutherland&#39;s law , where: = 1.7894 * 10^-5 Kg / (m * s) T_st...
void GetMach(Array< OneD, Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &soundspeed, Array< OneD, NekDouble > &mach)
Compute the mach number .
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)
~VariableConverter()
Destructor for VariableConverter class.
void GetTemperature(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &temperature)
Compute the temperature using the equation of state.
double NekDouble
std::shared_ptr< EquationOfState > EquationOfStateSharedPtr
A shared pointer to an equation of state object.
void GetSoundSpeed(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &soundspeed)
Compute the sound speed using the equation of state.
void GetVelocityVector(const Array< OneD, Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, NekDouble >> &velocity)
Compute the velocity field given the momentum .
void GetInternalEnergy(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &energy)
Compute the specific internal energy .
void GetRhoFromPT(const Array< OneD, NekDouble > &pressure, const Array< OneD, NekDouble > &temperature, Array< OneD, NekDouble > &rho)
Compute using the equation of state.
std::shared_ptr< SessionReader > SessionReaderSharedPtr