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 
52 
53 
54 namespace Nektar
55 {
56  /**
57  *
58  */
61  {
62  public:
63 
65 
66  virtual ~CompressibleFlowSystem();
67 
68  /// Function to calculate the stability limit for DG/CG.
70 
71  /// Function to calculate the stability limit for DG/CG
72  /// (a vector of them).
74  const Array<OneD,int> &ExpOrder);
75 
76  /// Function to get estimate of min h/p factor per element
78 
79  virtual void GetPressure(
80  const Array<OneD, const Array<OneD, NekDouble> > &physfield,
82 
83  virtual void GetDensity(
84  const Array<OneD, const Array<OneD, NekDouble> > &physfield,
85  Array<OneD, NekDouble> &density);
86 
87  virtual bool HasConstantDensity()
88  {
89  return false;
90  }
91 
92  virtual void GetVelocity(
93  const Array<OneD, const Array<OneD, NekDouble> > &physfield,
94  Array<OneD, Array<OneD, NekDouble> > &velocity);
95 
96  protected:
101  std::string m_shockCaptureType;
102 
103  // Parameters for exponential filtering
108 
109  // Parameters for local time-stepping
111 
112  // Auxiliary object to convert variables
114 
115  // User defined boundary conditions
116  std::vector<CFSBndCondSharedPtr> m_bndConds;
117 
118  // Forcing term
119  std::vector<SolverUtils::ForcingSharedPtr> m_forcing;
120 
122  const LibUtilities::SessionReaderSharedPtr& pSession,
123  const SpatialDomains::MeshGraphSharedPtr& pGraph);
124 
125  virtual void v_InitObject();
126 
127  void InitialiseParameters();
128 
129  void InitAdvection();
130 
131  void DoOdeRhs(
132  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
133  Array<OneD, Array<OneD, NekDouble> > &outarray,
134  const NekDouble time);
135  void DoOdeProjection(
136  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
137  Array<OneD, Array<OneD, NekDouble> > &outarray,
138  const NekDouble time);
139 
140  void DoAdvection(
141  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
142  Array<OneD, Array<OneD, NekDouble> > &outarray,
143  const NekDouble time,
144  const Array<OneD, Array<OneD, NekDouble> > &pFwd,
145  const Array<OneD, Array<OneD, NekDouble> > &pBwd);
146 
147  void DoDiffusion(
148  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
149  Array<OneD, Array<OneD, NekDouble> > &outarray,
150  const Array<OneD, Array<OneD, NekDouble> > &pFwd,
151  const Array<OneD, Array<OneD, NekDouble> > &pBwd);
152 
153  void GetFluxVector(
154  const Array<OneD, Array<OneD, NekDouble> > &physfield,
157  const Array<OneD, Array<OneD, NekDouble> > &physfield,
159 
161  Array<OneD, Array<OneD, NekDouble> > &physarray,
162  NekDouble time);
163 
164  void GetElmtTimeStep(
165  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
166  Array<OneD, NekDouble> &tstep);
167 
168  virtual NekDouble v_GetTimeStep(
169  const Array<OneD, const Array<OneD, NekDouble> > &inarray);
170  virtual void v_SetInitialConditions(
171  NekDouble initialtime = 0.0,
172  bool dumpInitialConditions = true,
173  const int domain = 0);
174 
176  {
177  return m_gamma;
178  }
179 
181  {
182  return m_vecLocs;
183  }
184 
186  {
187  return m_traceNormals;
188  }
189 
190  virtual void v_ExtraFldOutput(
191  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
192  std::vector<std::string> &variables);
193 
194  virtual void v_DoDiffusion(
195  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
196  Array<OneD, Array<OneD, NekDouble> > &outarray,
197  const Array<OneD, Array<OneD, NekDouble> > &pFwd,
198  const Array<OneD, Array<OneD, NekDouble> > &pBwd)
199  {
200  boost::ignore_unused(inarray, outarray, pFwd, pBwd);
201  // Do nothing by default
202  }
203 
205 
206  };
207 }
208 #endif
std::shared_ptr< Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:183
void DoDiffusion(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
Add the diffusions terms to the right-hand side.
std::shared_ptr< VariableConverter > VariableConverterSharedPtr
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:163
virtual void v_InitObject()
Initialization object for CompressibleFlowSystem class.
virtual ~CompressibleFlowSystem()
Destructor for CompressibleFlowSystem class.
General purpose memory allocation routines with the ability to allocate from thread specific memory p...
VariableConverterSharedPtr m_varConv
virtual void v_DoDiffusion(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
NekDouble GetStabilityLimit(int n)
Function to calculate the stability limit for DG/CG.
virtual void v_ExtraFldOutput(std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
std::shared_ptr< ArtificialDiffusion > ArtificialDiffusionSharedPtr
A shared pointer to a artificial diffusion object.
Array< OneD, Array< OneD, NekDouble > > m_vecLocs
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.
CompressibleFlowSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
virtual NekDouble v_GetTimeStep(const Array< OneD, const Array< OneD, NekDouble > > &inarray)
Calculate the maximum timestep subject to CFL restrictions.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the right-hand side.
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...
virtual void GetDensity(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &density)
Extract array with density from physfield.
virtual void GetVelocity(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity)
Extract array with velocity from physfield.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs()
void GetFluxVectorDeAlias(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Return the flux vector for the compressible Euler equations by using the de-aliasing technique...
virtual Array< OneD, NekDouble > v_GetMaxStdVelocity()
Compute the advection velocity in the standard space for each element of the expansion.
Array< OneD, NekDouble > GetStabilityLimitVector(const Array< OneD, int > &ExpOrder)
Function to calculate the stability limit for DG/CG (a vector of them).
void SetBoundaryConditions(Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
void InitialiseParameters()
Load CFS parameters from the session file.
double NekDouble
virtual void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0)
Set up logic for residual calculation.
void DoAdvection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time, const Array< OneD, Array< OneD, NekDouble > > &pFwd, const Array< OneD, Array< OneD, NekDouble > > &pBwd)
Compute the advection terms for the right-hand side.
void InitAdvection()
Create advection and diffusion objects for CFS.
SolverUtils::DiffusionSharedPtr m_diffusion
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals()
ArtificialDiffusionSharedPtr m_artificialDiffusion
SOLVER_UTILS_EXPORT MultiRegions::ExpListSharedPtr GetPressure()
Get pressure field if available.
Array< OneD, NekDouble > GetElmtMinHP(void)
Function to get estimate of min h/p factor per element.
void GetFluxVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Return the flux vector for the compressible Euler equations.
A base class for PDEs which include an advection component.
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
std::vector< CFSBndCondSharedPtr > m_bndConds
std::shared_ptr< SessionReader > SessionReaderSharedPtr