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
56namespace Nektar
57{
58/**
59 *
60 */
63{
64public:
66
68
69 /// Function to calculate the stability limit for DG/CG.
71
72 /// Function to calculate the stability limit for DG/CG
73 /// (a vector of them).
75 const Array<OneD, int> &ExpOrder);
76
77protected:
82 std::string m_shockCaptureType;
83
84 // Parameters for exponential filtering
89
90 // Parameters for local time-stepping
92
93 // Store physical artificial viscosity
95
96 // Store physical artificial viscosity
98
99 // Auxiliary object to convert variables
101
102 // User defined boundary conditions
103 std::vector<CFSBndCondSharedPtr> m_bndConds;
104
106
107 // Forcing term
108 std::vector<SolverUtils::ForcingSharedPtr> m_forcing;
109
112
113 virtual void v_InitObject(bool DeclareFields = true) override;
114
115 virtual void v_GetPressure(
116 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
118
119 virtual void v_GetDensity(
120 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
121 Array<OneD, NekDouble> &density) override;
122
123 virtual bool v_HasConstantDensity() override
124 {
125 return false;
126 }
127
128 virtual void v_GetVelocity(
129 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
130 Array<OneD, Array<OneD, NekDouble>> &velocity) override;
131
133
134 void InitAdvection();
135
136 void DoOdeRhs(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
138 const NekDouble time);
139 void DoOdeProjection(
140 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
141 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time);
142
143 void DoAdvection(const Array<OneD, Array<OneD, NekDouble>> &inarray,
145 const NekDouble time,
146 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
147 const Array<OneD, Array<OneD, NekDouble>> &pBwd);
148
149 void DoDiffusion(const Array<OneD, Array<OneD, NekDouble>> &inarray,
151 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
152 const Array<OneD, Array<OneD, NekDouble>> &pBwd);
153
154 void GetFluxVector(
155 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
158 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
160
162 NekDouble time);
163
165
166 void GetElmtTimeStep(
167 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
169
170 virtual NekDouble v_GetTimeStep(
171 const Array<OneD, const Array<OneD, NekDouble>> &inarray) override;
172
173 /// Print a summary of time stepping parameters.
174 virtual void v_GenerateSummary(SolverUtils::SummaryList &s) override;
175
176 virtual void v_SetInitialConditions(NekDouble initialtime = 0.0,
177 bool dumpInitialConditions = true,
178 const int domain = 0) override;
179
180 virtual void v_EvaluateExactSolution(unsigned int field,
181 Array<OneD, NekDouble> &outfield,
182 const NekDouble time = 0.0) override;
184 {
185 return m_gamma;
186 }
187
189 {
190 return m_vecLocs;
191 }
192
194 {
195 return m_traceNormals;
196 }
197
199 {
200 ASSERTL0(false, "This function is not valid for this class");
202 return null;
203 }
204
205 virtual void v_ExtraFldOutput(
206 std::vector<Array<OneD, NekDouble>> &fieldcoeffs,
207 std::vector<std::string> &variables) override;
208
209 virtual void v_DoDiffusion(
210 const Array<OneD, Array<OneD, NekDouble>> &inarray,
212 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
213 const Array<OneD, Array<OneD, NekDouble>> &pBwd) = 0;
214
216 const NekDouble SpeedSoundFactor) override;
217
218 virtual void v_SteadyStateResidual(int step,
219 Array<OneD, NekDouble> &L2) override;
220
221 // Virtual function that returns true if derived class supports a given
222 // shock capturing method
223 virtual bool v_SupportsShockCaptType(const std::string type) const = 0;
224
225private:
226 /// Isentropic Vortex Test Case.
227 void EvaluateIsentropicVortex(unsigned int field,
228 Array<OneD, NekDouble> &outfield,
229 NekDouble time, const int o = 0);
230
231 /// Ringleb Flow Test Case.
232 void GetExactRinglebFlow(int field, Array<OneD, NekDouble> &outarray);
233};
234} // namespace Nektar
235#endif
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
void GetFluxVector(const Array< OneD, const Array< OneD, NekDouble > > &physfield, TensorOfArray3D< NekDouble > &flux)
Return the flux vector for the compressible Euler equations.
CompressibleFlowSystem(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
virtual MultiRegions::ExpListSharedPtr v_GetPressure() override
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.
void SetBoundaryConditions(Array< OneD, Array< OneD, NekDouble > > &physarray, NekDouble time)
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, NekDouble > m_muavTrace
virtual bool v_SupportsShockCaptType(const std::string type) const =0
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 NekDouble v_GetTimeStep(const Array< OneD, const Array< OneD, NekDouble > > &inarray) override
Calculate the maximum timestep subject to CFL restrictions.
virtual ~CompressibleFlowSystem()
Destructor for CompressibleFlowSystem class.
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals()
virtual void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print a summary of time stepping parameters.
SolverUtils::DiffusionSharedPtr m_diffusion
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
virtual void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0) override
Set up logic for residual calculation.
void GetExactRinglebFlow(int field, Array< OneD, NekDouble > &outarray)
Ringleb Flow Test Case.
virtual void v_InitObject(bool DeclareFields=true) override
Initialization object for CompressibleFlowSystem class.
virtual void v_GetDensity(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, NekDouble > &density) override
void EvaluateIsentropicVortex(unsigned int field, Array< OneD, NekDouble > &outfield, NekDouble time, const int o=0)
Isentropic Vortex Test Case.
NekDouble GetStabilityLimit(int n)
Function to calculate the stability limit for DG/CG.
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.
std::vector< CFSBndCondSharedPtr > m_bndConds
virtual bool v_HasConstantDensity() override
Array< OneD, Array< OneD, NekDouble > > m_vecLocs
virtual void v_ExtraFldOutput(std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables) override
VariableConverterSharedPtr m_varConv
const Array< OneD, const Array< OneD, NekDouble > > & GetVecLocs()
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.
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.
void InitialiseParameters()
Load CFS parameters from the session file.
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.
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_SteadyStateResidual(int step, Array< OneD, NekDouble > &L2) override
virtual void v_EvaluateExactSolution(unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time=0.0) override
virtual void v_GetVelocity(const Array< OneD, const Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, NekDouble > > &velocity) override
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< Diffusion > DiffusionSharedPtr
A shared pointer to an EquationSystem object.
Definition: Diffusion.h:58
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:48
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:176
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