Nektar++
CompressibleFlowSystem.h
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File CompressibleFlowSystem.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 for the compressible flow system
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #ifndef NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_EQUATIONSYSTEMS_COMPRESSIBLEFLOWSYSTEM_H
36 #define NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_EQUATIONSYSTEMS_COMPRESSIBLEFLOWSYSTEM_H
37 
38 #include <boost/core/ignore_unused.hpp>
39 
55 
56 namespace Nektar
57 {
58  /**
59  *
60  */
63  {
64  public:
65 
67 
68  virtual ~CompressibleFlowSystem();
69 
70  /// Function to calculate the stability limit for DG/CG.
72 
73  /// Function to calculate the stability limit for DG/CG
74  /// (a vector of them).
76  const Array<OneD,int> &ExpOrder);
77 
78  /// Function to get estimate of min h/p factor per element
80 
81  virtual void GetPressure(
82  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
84 
85  virtual void GetDensity(
86  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
87  Array<OneD, NekDouble> &density);
88 
89  virtual bool HasConstantDensity()
90  {
91  return false;
92  }
93 
94  virtual void GetVelocity(
95  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
96  Array<OneD, Array<OneD, NekDouble>> &velocity);
97 
98  protected:
103  std::string m_shockCaptureType;
104 
105  // Parameters for exponential filtering
110 
111  // Parameters for local time-stepping
113 
114  // Store physical artificial viscosity
116 
117  // Store physical artificial viscosity
119 
120  // Auxiliary object to convert variables
122 
123  // User defined boundary conditions
124  std::vector<CFSBndCondSharedPtr> m_bndConds;
125 
127 
128  // Forcing term
129  std::vector<SolverUtils::ForcingSharedPtr> m_forcing;
130 
132  const LibUtilities::SessionReaderSharedPtr& pSession,
133  const SpatialDomains::MeshGraphSharedPtr& pGraph);
134 
135  virtual void v_InitObject();
136 
137  void InitialiseParameters();
138 
139  void InitAdvection();
140 
141  void DoOdeRhs(
142  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
143  Array<OneD, Array<OneD, NekDouble>> &outarray,
144  const NekDouble time);
145  void DoOdeProjection(
146  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
147  Array<OneD, Array<OneD, NekDouble>> &outarray,
148  const NekDouble time);
149 
150  void DoAdvection(
151  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
152  Array<OneD, Array<OneD, NekDouble>> &outarray,
153  const NekDouble time,
154  const Array<OneD, const Array<OneD, NekDouble>> &pFwd,
155  const Array<OneD, const Array<OneD, NekDouble>> &pBwd);
156 
157  void DoDiffusion(
158  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
159  Array<OneD, Array<OneD, NekDouble>> &outarray,
160  const Array<OneD, const Array<OneD, NekDouble>> &pFwd,
161  const Array<OneD, const Array<OneD, NekDouble>> &pBwd);
162 
163  void GetFluxVector(
164  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
167  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
169 
171  Array<OneD, Array<OneD, NekDouble>> &physarray,
172  NekDouble time);
173 
175 
176  void GetElmtTimeStep(
177  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
178  Array<OneD, NekDouble> &tstep);
179 
180  virtual NekDouble v_GetTimeStep(
181  const Array<OneD, const Array<OneD, NekDouble> > &inarray);
182  virtual void v_SetInitialConditions(
183  NekDouble initialtime = 0.0,
184  bool dumpInitialConditions = true,
185  const int domain = 0);
186 
188  {
189  return m_gamma;
190  }
191 
193  {
194  return m_vecLocs;
195  }
196 
198  {
199  return m_traceNormals;
200  }
201 
202  virtual void v_ExtraFldOutput(
203  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
204  std::vector<std::string> &variables);
205 
206  virtual void v_DoDiffusion(
207  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
208  Array<OneD, Array<OneD, NekDouble>> &outarray,
209  const Array<OneD, const Array<OneD, NekDouble>> &pFwd,
210  const Array<OneD, const Array<OneD, NekDouble>> &pBwd)
211  {
212  boost::ignore_unused(inarray, outarray, pFwd, pBwd);
213  // Do nothing by default
214  }
215 
217  const NekDouble SpeedSoundFactor);
218 
219  virtual void v_SteadyStateResidual(
220  int step,
222  };
223 }
224 #endif
CompressibleFlowSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
virtual void v_DoDiffusion(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, const Array< OneD, NekDouble >> &pFwd, const Array< OneD, const Array< OneD, NekDouble >> &pBwd)
virtual void GetVelocity(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, NekDouble >> &velocity)
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
Compute the right-hand side.
Array< OneD, NekDouble > m_muavTrace
void DoDiffusion(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, const Array< OneD, NekDouble >> &pFwd, const Array< OneD, const Array< OneD, NekDouble >> &pBwd)
Add the diffusions terms to the right-hand side.
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble >> &physfield, TensorOfArray3D< NekDouble > &flux)
Return the flux vector for the compressible Euler equations.
Array< OneD, NekDouble > GetStabilityLimitVector(const Array< OneD, int > &ExpOrder)
Function to calculate the stability limit for DG/CG (a vector of them).
void InitAdvection()
Create advection and diffusion objects for CFS.
virtual ~CompressibleFlowSystem()
Destructor for CompressibleFlowSystem class.
virtual Array< OneD, NekDouble > v_GetMaxStdVelocity(const NekDouble SpeedSoundFactor)
Compute the advection velocity in the standard space for each element of the expansion.
void GetFluxVectorDeAlias(const Array< OneD, const Array< OneD, NekDouble >> &physfield, TensorOfArray3D< NekDouble > &flux)
Return the flux vector for the compressible Euler equations by using the de-aliasing technique.
void SetBoundaryConditions(Array< OneD, Array< OneD, NekDouble >> &physarray, NekDouble time)
SolverUtils::DiffusionSharedPtr m_diffusion
virtual void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
Set up logic for residual calculation.
NekDouble GetStabilityLimit(int n)
Function to calculate the stability limit for DG/CG.
std::vector< CFSBndCondSharedPtr > m_bndConds
void GetElmtTimeStep(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, NekDouble > &tstep)
Calculate the maximum timestep on each element subject to CFL restrictions.
Array< OneD, Array< OneD, NekDouble > > m_vecLocs
Array< OneD, NekDouble > GetElmtMinHP(void)
Function to get estimate of min h/p factor per element.
virtual void v_SteadyStateResidual(int step, Array< OneD, NekDouble > &L2)
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
Compute the projection and call the method for imposing the boundary conditions in case of discontinu...
VariableConverterSharedPtr m_varConv
void SetBoundaryConditionsBwdWeight()
Set up a weight on physical boundaries for boundary condition applications.
virtual void v_InitObject()
Initialization object for CompressibleFlowSystem class.
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
ArtificialDiffusionSharedPtr m_artificialDiffusion
virtual void GetDensity(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &density)
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals()
void InitialiseParameters()
Load CFS parameters from the session file.
void DoAdvection(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time, const Array< OneD, const Array< OneD, NekDouble >> &pFwd, const Array< OneD, const Array< OneD, NekDouble >> &pBwd)
Compute the advection terms for the right-hand side.
virtual void v_ExtraFldOutput(std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs()
virtual NekDouble v_GetTimeStep(const Array< OneD, const Array< OneD, NekDouble > > &inarray)
Calculate the maximum timestep subject to CFL restrictions.
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
A base class for PDEs which include an advection component.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr GetPressure()
Get pressure field if available.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< SolverUtils::Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:627
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
std::shared_ptr< VariableConverter > VariableConverterSharedPtr
std::shared_ptr< ArtificialDiffusion > ArtificialDiffusionSharedPtr
A shared pointer to a artificial diffusion object.
double NekDouble