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
53
54namespace Nektar
55{
56/**
57 *
58 */
61{
62public:
64
65 ~CompressibleFlowSystem() override = default;
66
67 /// Function to calculate the stability limit for DG/CG.
69
70 /// Function to calculate the stability limit for DG/CG
71 /// (a vector of them).
73 const Array<OneD, int> &ExpOrder);
74
75protected:
80 std::string m_shockCaptureType;
81
82 // Parameters for exponential filtering
87
88 // Parameters for local time-stepping
90
91 // Store physical artificial viscosity
93
94 // Store physical artificial viscosity
96
97 // Auxiliary object to convert variables
99
100 // User defined boundary conditions
101 std::vector<CFSBndCondSharedPtr> m_bndConds;
102
104
105 // Forcing term
106 std::vector<SolverUtils::ForcingSharedPtr> m_forcing;
107
110
111 void v_InitObject(bool DeclareFields = true) override;
112
113 void v_GetPressure(
114 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
116
117 void v_GetDensity(
118 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
119 Array<OneD, NekDouble> &density) override;
120
121 bool v_HasConstantDensity() override
122 {
123 return false;
124 }
125
126 void v_GetVelocity(
127 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
129
131
132 void InitAdvection();
133
134 void DoOdeRhs(const Array<OneD, const Array<OneD, NekDouble>> &inarray,
136 const NekDouble time);
137 void DoOdeProjection(
138 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
139 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time);
140
141 void DoAdvection(const Array<OneD, Array<OneD, NekDouble>> &inarray,
143 const NekDouble time,
144 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
145 const Array<OneD, Array<OneD, NekDouble>> &pBwd);
146
147 void DoDiffusion(const Array<OneD, Array<OneD, NekDouble>> &inarray,
149 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
150 const Array<OneD, Array<OneD, NekDouble>> &pBwd);
151
152 void GetFluxVector(
153 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
156 const Array<OneD, const Array<OneD, NekDouble>> &physfield,
158
160 NekDouble time);
161
163
164 void GetElmtTimeStep(
165 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
167
169 const Array<OneD, const Array<OneD, NekDouble>> &inarray) override;
170
171 /// Print a summary of time stepping parameters.
173
174 void v_SetInitialConditions(NekDouble initialtime = 0.0,
175 bool dumpInitialConditions = true,
176 const int domain = 0) override;
177
178 void v_EvaluateExactSolution(unsigned int field,
179 Array<OneD, NekDouble> &outfield,
180 const NekDouble time = 0.0) override;
182 {
183 return m_gamma;
184 }
185
187 {
188 return m_vecLocs;
189 }
190
192 {
193 return m_traceNormals;
194 }
195
197 {
198 ASSERTL0(false, "This function is not valid for this class");
200 return null;
201 }
202
203 void v_ExtraFldOutput(std::vector<Array<OneD, NekDouble>> &fieldcoeffs,
204 std::vector<std::string> &variables) override;
205
206 virtual void v_DoDiffusion(
207 const Array<OneD, Array<OneD, NekDouble>> &inarray,
209 const Array<OneD, Array<OneD, NekDouble>> &pFwd,
210 const Array<OneD, Array<OneD, NekDouble>> &pBwd) = 0;
211
213 const NekDouble SpeedSoundFactor) override;
214
215 void v_SteadyStateResidual(int step, Array<OneD, NekDouble> &L2) override;
216
217 // Virtual function that returns true if derived class supports a given
218 // shock capturing method
219 virtual bool v_SupportsShockCaptType(const std::string type) const = 0;
220
221private:
222 /// Isentropic Vortex Test Case.
223 void EvaluateIsentropicVortex(unsigned int field,
224 Array<OneD, NekDouble> &outfield,
225 NekDouble time, const int o = 0);
226
227 /// Ringleb Flow Test Case.
228 void GetExactRinglebFlow(int field, Array<OneD, NekDouble> &outarray);
229};
230} // namespace Nektar
231#endif
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
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)
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
MultiRegions::ExpListSharedPtr v_GetPressure() override
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.
NekDouble v_GetTimeStep(const Array< OneD, const Array< OneD, NekDouble > > &inarray) override
Calculate the maximum timestep subject to CFL restrictions.
const Array< OneD, const Array< OneD, NekDouble > > & GetNormals()
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
void v_SetInitialConditions(NekDouble initialtime=0.0, bool dumpInitialConditions=true, const int domain=0) override
Set up logic for residual calculation.
~CompressibleFlowSystem() override=default
void GetExactRinglebFlow(int field, Array< OneD, NekDouble > &outarray)
Ringleb Flow Test Case.
void v_InitObject(bool DeclareFields=true) override
Initialization object for CompressibleFlowSystem class.
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
Array< OneD, Array< OneD, NekDouble > > m_vecLocs
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
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.
void v_SteadyStateResidual(int step, Array< OneD, NekDouble > &L2) override
void v_EvaluateExactSolution(unsigned int field, Array< OneD, NekDouble > &outfield, const NekDouble time=0.0) override
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
const std::vector< NekDouble > velocity
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:54
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:46
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::shared_ptr< VariableConverter > VariableConverterSharedPtr
std::shared_ptr< ArtificialDiffusion > ArtificialDiffusionSharedPtr
A shared pointer to a artificial diffusion object.
double NekDouble