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 
54 
55 namespace Nektar
56 {
57 /**
58  *
59  */
62 {
63 public:
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  virtual void v_GetPressure(
77  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
79 
80  virtual void v_GetDensity(
81  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
82  Array<OneD, NekDouble> &density) override;
83 
84  virtual bool v_HasConstantDensity() override
85  {
86  return false;
87  }
88 
89  virtual void v_GetVelocity(
90  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
91  Array<OneD, Array<OneD, NekDouble>> &velocity) override;
92 
93 protected:
98  std::string m_shockCaptureType;
99 
100  // Parameters for exponential filtering
105 
106  // Parameters for local time-stepping
108 
109  // Store physical artificial viscosity
111 
112  // Store physical artificial viscosity
114 
115  // Auxiliary object to convert variables
117 
118  // User defined boundary conditions
119  std::vector<CFSBndCondSharedPtr> m_bndConds;
120 
122 
123  // Forcing term
124  std::vector<SolverUtils::ForcingSharedPtr> m_forcing;
125 
127  const SpatialDomains::MeshGraphSharedPtr &pGraph);
128 
129  virtual void v_InitObject(bool DeclareFields = true) override;
130 
131  void InitialiseParameters();
132 
133  void InitAdvection();
134 
135  void DoOdeRhs(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
136  Array<OneD, Array<OneD, NekDouble>> &outarray,
137  const NekDouble time);
138  void DoOdeProjection(
139  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
140  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time);
141 
142  void DoAdvection(const Array<OneD, Array<OneD, NekDouble>> &inarray,
143  Array<OneD, Array<OneD, NekDouble>> &outarray,
144  const NekDouble time,
145  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
146  const Array<OneD, Array<OneD, NekDouble>> &pBwd);
147 
148  void DoDiffusion(const Array<OneD, 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, const Array<OneD, NekDouble>> &physfield,
157  const Array<OneD, const Array<OneD, NekDouble>> &physfield,
159 
161  NekDouble time);
162 
164 
165  void GetElmtTimeStep(
166  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
167  Array<OneD, NekDouble> &tstep);
168 
169  virtual NekDouble v_GetTimeStep(
170  const Array<OneD, const Array<OneD, NekDouble>> &inarray) override;
171 
172  virtual void v_SetInitialConditions(NekDouble initialtime = 0.0,
173  bool dumpInitialConditions = true,
174  const int domain = 0) override;
175 
177  {
178  return m_gamma;
179  }
180 
182  {
183  return m_vecLocs;
184  }
185 
187  {
188  return m_traceNormals;
189  }
190 
192  {
193  ASSERTL0(false, "This function is not valid for this class");
195  return null;
196  }
197 
198  virtual void v_ExtraFldOutput(
199  std::vector<Array<OneD, NekDouble>> &fieldcoeffs,
200  std::vector<std::string> &variables) override;
201 
202  virtual void v_DoDiffusion(
203  const Array<OneD, Array<OneD, NekDouble>> &inarray,
204  Array<OneD, Array<OneD, NekDouble>> &outarray,
205  const Array<OneD, Array<OneD, NekDouble>> &pFwd,
206  const Array<OneD, Array<OneD, NekDouble>> &pBwd) = 0;
207 
209  const NekDouble SpeedSoundFactor) override;
210 
211  virtual void v_SteadyStateResidual(int step,
212  Array<OneD, NekDouble> &L2) override;
213 
214  // Virtual function that returns true if derived class supports a given
215  // shock capturing method
216  virtual bool v_SupportsShockCaptType(const std::string type) const = 0;
217 };
218 } // namespace Nektar
219 #endif
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
virtual void v_DoDiffusion(const Array< OneD, Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const Array< OneD, Array< OneD, NekDouble >> &pFwd, const Array< OneD, Array< OneD, NekDouble >> &pBwd)=0
CompressibleFlowSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
virtual MultiRegions::ExpListSharedPtr v_GetPressure() override
virtual NekDouble v_GetTimeStep(const Array< OneD, const Array< OneD, NekDouble >> &inarray) override
Calculate the maximum timestep subject to CFL restrictions.
void DoAdvection(const Array< OneD, 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.
virtual void v_GetVelocity(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, Array< OneD, NekDouble >> &velocity) override
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
Compute the right-hand side.
virtual void v_GetDensity(const Array< OneD, const Array< OneD, NekDouble >> &physfield, Array< OneD, NekDouble > &density) override
void DoDiffusion(const Array< OneD, 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.
Array< OneD, NekDouble > m_muavTrace
virtual bool v_SupportsShockCaptType(const std::string type) const =0
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.
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) override
Set up logic for residual calculation.
virtual void v_InitObject(bool DeclareFields=true) override
Initialization object for CompressibleFlowSystem class.
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.
virtual bool v_HasConstantDensity() override
Array< OneD, Array< OneD, NekDouble > > m_vecLocs
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.
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
ArtificialDiffusionSharedPtr m_artificialDiffusion
virtual Array< OneD, NekDouble > v_GetMaxStdVelocity(const NekDouble SpeedSoundFactor) override
Compute the advection velocity in the standard space for each element of the expansion.
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals()
void InitialiseParameters()
Load CFS parameters from the session file.
virtual void v_ExtraFldOutput(std::vector< Array< OneD, NekDouble >> &fieldcoeffs, std::vector< std::string > &variables) override
virtual void v_SteadyStateResidual(int step, Array< OneD, NekDouble > &L2) override
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs()
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.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< ExpList > ExpListSharedPtr
Shared pointer to an ExpList object.
std::shared_ptr< SolverUtils::Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:550
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< ArtificialDiffusion > ArtificialDiffusionSharedPtr
A shared pointer to a artificial diffusion object.
double NekDouble