Nektar++
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 // License for the specific language governing rights and limitations under
14 // Permission is hereby granted, free of charge, to any person obtaining a
15 // copy of this software and associated documentation files (the "Software"),
16 // to deal in the Software without restriction, including without limitation
17 // the rights to use, copy, modify, merge, publish, distribute, sublicense,
18 // and/or sell copies of the Software, and to permit persons to whom the
19 // Software is furnished to do so, subject to the following conditions:
20 //
21 // The above copyright notice and this permission notice shall be included
22 // in all copies or substantial portions of the Software.
23 //
24 // THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25 // OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26 // FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27 // THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28 // LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29 // FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30 // DEALINGS IN THE SOFTWARE.
31 //
32 // Description: Auxiliary functions for the compressible flow system
33 //
34 ///////////////////////////////////////////////////////////////////////////////
35 
36 #ifndef NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_EQUATIONSYSTEMS_COMPRESSIBLEFLOWSYSTEM_H
37 #define NEKTAR_SOLVERS_COMPRESSIBLEFLOWSOLVER_EQUATIONSYSTEMS_COMPRESSIBLEFLOWSYSTEM_H
38 
47 
48 namespace Nektar
49 {
50  /**
51  *
52  */
54  {
55  public:
56 
58 
59  virtual ~CompressibleFlowSystem();
60 
61  /// Function to calculate the stability limit for DG/CG.
63 
64  /// Function to calculate the stability limit for DG/CG
65  /// (a vector of them).
67  const Array<OneD,int> &ExpOrder);
68 
69  protected:
78  std::string m_ViscosityType;
79  std::string m_shockCaptureType;
84 
85  // Auxiliary object to convert variables
87 
88  // User defined boundary conditions
89  std::vector<CFSBndCondSharedPtr> m_bndConds;
90 
91  // L2 error file
92  std::ofstream m_errFile;
93 
94  // Tolerance to which steady state should be evaluated at
96 
97  // Forcing term
98  std::vector<SolverUtils::ForcingSharedPtr> m_forcing;
99 
100  // Storage for L2 norm error
102 
104  const LibUtilities::SessionReaderSharedPtr& pSession);
105 
106  virtual void v_InitObject();
107 
108  void InitialiseParameters();
109 
110  void InitAdvection();
111 
112  void DoOdeRhs(
113  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
114  Array<OneD, Array<OneD, NekDouble> > &outarray,
115  const NekDouble time);
116  void DoOdeProjection(
117  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
118  Array<OneD, Array<OneD, NekDouble> > &outarray,
119  const NekDouble time);
120 
121  void DoAdvection(
122  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
123  Array<OneD, Array<OneD, NekDouble> > &outarray,
124  const NekDouble time,
125  const Array<OneD, Array<OneD, NekDouble> > &pFwd,
126  const Array<OneD, Array<OneD, NekDouble> > &pBwd);
127 
128  void DoDiffusion(
129  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
130  Array<OneD, Array<OneD, NekDouble> > &outarray,
131  const Array<OneD, Array<OneD, NekDouble> > &pFwd,
132  const Array<OneD, Array<OneD, NekDouble> > &pBwd);
133 
134  void GetFluxVector(
135  const Array<OneD, Array<OneD, NekDouble> > &physfield,
138  const Array<OneD, Array<OneD, NekDouble> > &physfield,
140 
141  void InitializeSteadyState();
142 
144  Array<OneD, Array<OneD, NekDouble> > &physarray,
145  NekDouble time);
146 
147  void GetStdVelocity(
148  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
149  Array<OneD, NekDouble> &stdV);
150 
151  virtual bool v_PostIntegrate(int step);
152  bool CalcSteadyState(bool output);
153 
154  virtual NekDouble v_GetTimeStep(
155  const Array<OneD, const Array<OneD, NekDouble> > &inarray);
156  virtual void v_SetInitialConditions(
157  NekDouble initialtime = 0.0,
158  bool dumpInitialConditions = true,
159  const int domain = 0);
160 
162  {
163  return m_gamma;
164  }
165 
167  {
168  return m_vecLocs;
169  }
170 
172  {
173  return m_traceNormals;
174  }
175 
176  virtual void v_ExtraFldOutput(
177  std::vector<Array<OneD, NekDouble> > &fieldcoeffs,
178  std::vector<std::string> &variables);
179 
180  virtual void v_DoDiffusion(
181  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
182  Array<OneD, Array<OneD, NekDouble> > &outarray,
183  const Array<OneD, Array<OneD, NekDouble> > &pFwd,
184  const Array<OneD, Array<OneD, NekDouble> > &pBwd)
185  {
186  // Do nothing by default
187  }
188  };
189 }
190 #endif
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.
virtual void v_InitObject()
Initialization object for CompressibleFlowSystem class.
virtual ~CompressibleFlowSystem()
Destructor for CompressibleFlowSystem class.
void GetStdVelocity(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, NekDouble > &stdV)
Compute the advection velocity in the standard space for each element of the expansion.
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)
boost::shared_ptr< ArtificialDiffusion > ArtificialDiffusionSharedPtr
A shared pointer to a artificial diffusion object.
boost::shared_ptr< Advection > AdvectionSharedPtr
A shared pointer to an Advection object.
Definition: Advection.h:165
Array< OneD, Array< OneD, NekDouble > > m_vecLocs
Array< OneD, Array< OneD, NekDouble > > m_un
boost::shared_ptr< SessionReader > SessionReaderSharedPtr
Definition: MeshPartition.h:51
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...
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...
Base class for unsteady solvers.
Array< OneD, NekDouble > GetStabilityLimitVector(const Array< OneD, int > &ExpOrder)
Function to calculate the stability limit for DG/CG (a vector of them).
boost::shared_ptr< Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:162
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.
virtual bool v_PostIntegrate(int step)
Perform post-integration checks, presently just to check steady state behaviour.
bool CalcSteadyState(bool output)
Calculate whether the system has reached a steady state by observing residuals to a user-defined tole...
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.
SolverUtils::AdvectionSharedPtr m_advection
void InitAdvection()
Create advection and diffusion objects for CFS.
SolverUtils::DiffusionSharedPtr m_diffusion
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals()
ArtificialDiffusionSharedPtr m_artificialDiffusion
boost::shared_ptr< VariableConverter > VariableConverterSharedPtr
CompressibleFlowSystem(const LibUtilities::SessionReaderSharedPtr &pSession)
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.
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
std::vector< CFSBndCondSharedPtr > m_bndConds