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
45 class VariableConverter;
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 GetDynamicEnergy(
60  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
61  Array<OneD, NekDouble> &energy);
62  void GetInternalEnergy(
63  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
64  Array<OneD, NekDouble> &energy);
65  template <class T, typename = typename std::enable_if
66  <
67  std::is_floating_point<T>::value ||
69  >::type
70  >
71  inline T GetInternalEnergy(T* physfield)
72  {
73  // get dynamic energy
74  T oneOrho = 1.0 / physfield[0];
75  T dynEne{};
76  for (size_t d = 1; d < m_spacedim + 1; ++d)
77  {
78  T tmp = physfield[d]; //load 1x
79  dynEne += tmp * tmp;
80  }
81  dynEne = 0.5 * dynEne * oneOrho;
82 
83  // Calculate rhoe = E - rho*V^2/2
84  T energy = physfield[m_spacedim + 1] - dynEne;
85  return energy * oneOrho;
86  }
87  void GetEnthalpy(const Array<OneD, const Array<OneD, NekDouble>> &physfield,
88  Array<OneD, NekDouble> &enthalpy);
89  void GetVelocityVector(const Array<OneD, Array<OneD, NekDouble>> &physfield,
90  Array<OneD, Array<OneD, NekDouble>> &velocity);
91  void GetMach(Array<OneD, Array<OneD, NekDouble>> &physfield,
92  Array<OneD, NekDouble> &soundspeed,
94  void GetDynamicViscosity(const Array<OneD, const NekDouble> &temperature,
96 
97  template <class T, typename = typename std::enable_if
98  <
99  std::is_floating_point<T>::value ||
101  >::type
102  >
103  inline T GetDynamicViscosity(T &temperature)
104  {
105  constexpr NekDouble C = .38175;
106  constexpr NekDouble onePlusC = 1.0 + C;
107 
108  NekDouble mu_star = m_mu;
109 
110  T ratio = temperature * m_oneOverT_star;
111  return mu_star * ratio * sqrt(ratio) * onePlusC / (ratio + C);
112  }
113 
114  void GetAbsoluteVelocity(
115  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
116  Array<OneD, NekDouble> &Vtot);
117  void GetSensor(const MultiRegions::ExpListSharedPtr &field,
118  const Array<OneD, const Array<OneD, NekDouble>> &physarray,
119  Array<OneD, NekDouble> &Sensor,
120  Array<OneD, NekDouble> &SensorKappa, int offset = 1);
121 
122  // Transformations depending on the equation of state
123  void GetTemperature(
124  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
125  Array<OneD, NekDouble> &temperature);
126  template <class T, typename = typename std::enable_if
127  <
128  std::is_floating_point<T>::value ||
130  >::type
131  >
132  inline T GetTemperature(T* physfield)
133  {
134  T energy = GetInternalEnergy(physfield);
135  return m_eos->GetTemperature(physfield[0], energy);
136  }
137  //
138  void GetPressure(const Array<OneD, const Array<OneD, NekDouble>> &physfield,
140  template <class T, typename = typename std::enable_if
141  <
142  std::is_floating_point<T>::value ||
144  >::type
145  >
146  inline T GetPressure(T* physfield)
147  {
148  T energy = GetInternalEnergy(physfield);
149  return m_eos->GetPressure(physfield[0], energy);
150  }
151 
152  void GetSoundSpeed(
153  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
154  Array<OneD, NekDouble> &soundspeed);
155  void GetEntropy(const Array<OneD, const Array<OneD, NekDouble>> &physfield,
156  Array<OneD, NekDouble> &entropy);
157  void GetEFromRhoP(const Array<OneD, NekDouble> &rho,
159  Array<OneD, NekDouble> &energy);
161  const Array<OneD, NekDouble> &temperature,
163  void GetDmuDT(
164  const Array<OneD, const NekDouble> &temperature,
165  const Array<OneD, const NekDouble> &mu,
166  Array<OneD, NekDouble> &DmuDT);
167 
169  {
170  return m_eos;
171  }
172 
173 protected:
184 };
185 
186 
187 
188 
189 
190 
191 
192 }
193 #endif
EquationOfStateSharedPtr m_eos
void GetDynamicEnergy(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &energy)
Compute the dynamic energy .
void GetRhoFromPT(const Array< OneD, NekDouble > &pressure, const Array< OneD, NekDouble > &temperature, Array< OneD, NekDouble > &rho)
Compute using the equation of state.
void GetTemperature(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &temperature)
Compute the temperature using the equation of state.
VariableConverter(const LibUtilities::SessionReaderSharedPtr &pSession, const int spaceDim)
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 , where: \mu_star = 1.7894 * 10^-5 Kg / (m *...
void GetEnthalpy(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &enthalpy)
Compute the specific enthalpy .
void GetEFromRhoP(const Array< OneD, NekDouble > &rho, const Array< OneD, NekDouble > &pressure, Array< OneD, NekDouble > &energy)
Compute using the equation of state.
LibUtilities::SessionReaderSharedPtr m_session
void GetPressure(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &pressure)
Calculate the pressure using the equation of state.
void GetSoundSpeed(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &soundspeed)
Compute the sound speed using the equation of state.
void GetAbsoluteVelocity(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &Vtot)
void GetEntropy(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &entropy)
Compute the entropy 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 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)
const EquationOfStateSharedPtr Geteos()
~VariableConverter()
Destructor for VariableConverter class.
void GetInternalEnergy(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &energy)
Compute the specific internal energy .
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 , where: \mu_star = 1.7894 * 10^-5 Kg / (m *...
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
std::shared_ptr< VariableConverter > VariableConverterSharedPtr
std::shared_ptr< EquationOfState > EquationOfStateSharedPtr
A shared pointer to an equation of state object.
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:267