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"
40 #include <MultiRegions/ContField.h>
42 
43 namespace Nektar
44 {
45 // Forward declarations
46 class VariableConverter;
47 typedef std::shared_ptr<VariableConverter> VariableConverterSharedPtr;
48 /**
49  *
50  */
52 {
53 public:
56  const int spaceDim,
57  const SpatialDomains::MeshGraphSharedPtr &pGraph = nullptr);
58 
60 
61  // Variable manipulations valid for all fluids
62  void GetDynamicEnergy(
63  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
64  Array<OneD, NekDouble> &energy);
65  void GetInternalEnergy(
66  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
67  Array<OneD, NekDouble> &energy);
68  template <class T, typename = typename std::enable_if<
69  std::is_floating_point<T>::value ||
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  std::is_floating_point<T>::value ||
100  inline T GetDynamicViscosity(T &temperature)
101  {
102  const NekDouble onePlusC = 1.0 + m_TRatioSutherland;
103 
104  NekDouble mu_star = m_mu;
105 
106  T ratio = temperature * m_oneOverT_star;
107  return mu_star * ratio * sqrt(ratio) * onePlusC /
108  (ratio + m_TRatioSutherland);
109  }
110 
111  void GetAbsoluteVelocity(
112  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
113  Array<OneD, NekDouble> &Vtot);
114 
115  // Transformations depending on the equation of state
116  void GetTemperature(
117  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
118  Array<OneD, NekDouble> &temperature);
119  template <class T, typename = typename std::enable_if<
120  std::is_floating_point<T>::value ||
122  inline T GetTemperature(T *physfield)
123  {
124  T energy = GetInternalEnergy(physfield);
125  return m_eos->GetTemperature(physfield[0], energy);
126  }
127  //
128  void GetPressure(const Array<OneD, const Array<OneD, NekDouble>> &physfield,
130  template <class T, typename = typename std::enable_if<
131  std::is_floating_point<T>::value ||
133  inline T GetPressure(T *physfield)
134  {
135  T energy = GetInternalEnergy(physfield);
136  return m_eos->GetPressure(physfield[0], energy);
137  }
138 
139  void GetSoundSpeed(
140  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
141  Array<OneD, NekDouble> &soundspeed);
142  void GetEntropy(const Array<OneD, const Array<OneD, NekDouble>> &physfield,
143  Array<OneD, NekDouble> &entropy);
144  void GetEFromRhoP(const Array<OneD, NekDouble> &rho,
146  Array<OneD, NekDouble> &energy);
148  const Array<OneD, NekDouble> &temperature,
150  void GetDmuDT(const Array<OneD, const NekDouble> &temperature,
152  Array<OneD, NekDouble> &DmuDT);
153 
155  {
156  return m_eos;
157  }
158 
159  // Shock sensor methods
160  void SetAv(
162  const Array<OneD, const Array<OneD, NekDouble>> &consVar,
164  const Array<OneD, NekDouble> &curlSquared = NullNekDouble1DArray);
165 
167 
169 
170  bool GetFlagCalcDivCurl(void) const
171  {
172  return m_flagCalcDivCurl;
173  }
174 
175  void SetElmtMinHP(
177 
179 
180  void GetSensor(const MultiRegions::ExpListSharedPtr &field,
181  const Array<OneD, const Array<OneD, NekDouble>> &physarray,
182  Array<OneD, NekDouble> &Sensor,
183  Array<OneD, NekDouble> &SensorKappa, int offset = 1);
184 
186  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
187  Array<OneD, NekDouble> &muAv);
188 
190  const Array<OneD, const Array<OneD, NekDouble>> &consVar,
191  const Array<OneD, NekDouble> &div,
192  Array<OneD, NekDouble> &muAv);
193 
194  void ApplyDucros(const Array<OneD, NekDouble> &div,
195  const Array<OneD, NekDouble> &curlSquare,
196  Array<OneD, NekDouble> &muAv);
197 
199 
200 protected:
203  size_t m_spacedim;
213 
214  /// Shock sensor
216  std::string m_shockCaptureType;
217  std::string m_shockSensorType;
218  std::string m_ducrosSensor;
219  std::string m_smoothing;
221 
222  /// h/p scaling
224  /// storage
227  bool m_flagCalcDivCurl = false;
228 };
229 
230 } // namespace Nektar
231 #endif
void ApplyC0Smooth(Array< OneD, NekDouble > &field)
Make field C0.
Array< OneD, NekDouble > & GetAvTrace()
EquationOfStateSharedPtr m_eos
void GetDynamicEnergy(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &energy)
Compute the dynamic energy .
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 GetTemperature(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &temperature)
Compute the temperature using the equation of state.
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 ,.
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.
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
LibUtilities::SessionReaderSharedPtr m_session
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)
NekDouble m_mu0
Shock sensor.
void GetPressure(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &pressure)
Calculate the pressure using the equation of state.
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 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)
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 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 .
bool GetFlagCalcDivCurl(void) const
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()
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.
~VariableConverter()
Destructor for VariableConverter class.
Array< OneD, NekDouble > & GetElmtMinHP()
Array< OneD, NekDouble > m_muAv
storage
Array< OneD, NekDouble > & GetAv()
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 , 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
Definition: ContField.h:289
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
std::shared_ptr< VariableConverter > VariableConverterSharedPtr
std::shared_ptr< EquationOfState > EquationOfStateSharedPtr
A shared pointer to an equation of state object.
static Array< OneD, NekDouble > NullNekDouble1DArray
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294