Nektar++
UnsteadyInviscidBurgers.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: UnsteadyInviscidBurgers.cpp
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: Unsteady inviscid Burgers solve routines
32//
33///////////////////////////////////////////////////////////////////////////////
34
37
38namespace Nektar
39{
42 "UnsteadyInviscidBurgers", UnsteadyInviscidBurgers::create,
43 "Inviscid Burgers equation");
44
48 : UnsteadySystem(pSession, pGraph), AdvectionSystem(pSession, pGraph)
49{
50}
51
52/**
53 * @brief Initialisation object for the inviscid Burgers equation.
54 */
56{
57 AdvectionSystem::v_InitObject(DeclareFields);
58
59 // Type of advection class to be used
60 switch (m_projectionType)
61 {
62 // Discontinuous field
64 {
65 // Do not forwards transform initial condition
66 m_homoInitialFwd = false;
67
68 // Define the normal velocity fields
69 if (m_fields[0]->GetTrace())
70 {
72 }
73
74 // Advection term
75 std::string advName;
76 std::string riemName;
77 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
79 advName, advName);
81 this);
82 m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
85 riemName, m_session);
86 m_riemannSolver->SetScalar(
88 m_advObject->SetRiemannSolver(m_riemannSolver);
89 m_advObject->InitObject(m_session, m_fields);
90 break;
91 }
92 // Continuous field
94 {
95 std::string advName;
96 m_session->LoadSolverInfo("AdvectionType", advName,
97 "NonConservative");
99 advName, advName);
101 this);
102 break;
103 }
104 default:
105 {
106 ASSERTL0(false, "Unsupported projection type.");
107 break;
108 }
109 }
110
111 // Forcing terms
112 m_forcing = SolverUtils::Forcing::Load(m_session, shared_from_this(),
113 m_fields, m_fields.size());
114
116 {
119 }
120 else
121 {
122 ASSERTL0(false, "Implicit unsteady Advection not set up.");
123 }
124}
125
126/**
127 * @brief Get the normal velocity for the inviscid Burgers equation.
128 */
130{
131 // Number of trace (interface) points
132 int nTracePts = GetTraceNpoints();
133
134 // Number of solution points
135 int nSolutionPts = GetNpoints();
136
137 // Auxiliary variables to compute the normal velocity
138 Array<OneD, NekDouble> Fwd(nTracePts);
139 Array<OneD, NekDouble> Bwd(nTracePts);
140 Array<OneD, NekDouble> physfield(nSolutionPts);
141
142 // Reset the normal velocity
143 Vmath::Zero(nTracePts, m_traceVn, 1);
144
145 for (int i = 0; i < m_spacedim; ++i)
146 {
147 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), physfield);
148 m_fields[i]->GetFwdBwdTracePhys(physfield, Fwd, Bwd, true);
149 Vmath::Vadd(nTracePts, Fwd, 1, Bwd, 1, Fwd, 1);
150 Vmath::Smul(nTracePts, 0.5, Fwd, 1, Fwd, 1);
151 Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, Fwd, 1, m_traceVn, 1,
152 m_traceVn, 1);
153 }
154 Vmath::Smul(nTracePts, 0.5, m_traceVn, 1, m_traceVn, 1);
155
156 return m_traceVn;
157}
158
159/**
160 * @brief Compute the right-hand side for the inviscid Burgers equation.
161 *
162 * @param inarray Given fields.
163 * @param outarray Calculated solution.
164 * @param time Time.
165 */
167 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
168 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
169{
170 // Number of fields (variables of the problem)
171 int nVariables = inarray.size();
172
173 // Number of solution points
174 int nSolutionPts = GetNpoints();
175
177 timer.Start();
178 // RHS computation using the new advection base class
179 m_advObject->Advect(nVariables, m_fields, inarray, inarray, outarray, time);
180 timer.Stop();
181 // Elapsed time
182 timer.AccumulateRegion("Advect");
183
184 // Negate the RHS
185 for (int i = 0; i < nVariables; ++i)
186 {
187 Vmath::Neg(nSolutionPts, outarray[i], 1);
188 }
189
190 // Add forcing terms
191 for (auto &x : m_forcing)
192 {
193 // set up non-linear terms
194 x->Apply(m_fields, inarray, outarray, time);
195 }
196}
197
198/**
199 * @brief Compute the projection for the inviscid Burgers equation.
200 *
201 * @param inarray Given fields.
202 * @param outarray Calculated solution.
203 * @param time Time.
204 */
206 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
207 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
208{
209 // Number of variables of the problem
210 int nVariables = inarray.size();
211
212 // Set the boundary conditions
214
215 // Switch on the projection type (Discontinuous or Continuous)
216 switch (m_projectionType)
217 {
218 // Discontinuous projection
220 {
221 // Just copy over array
222 if (inarray != outarray)
223 {
224 int npoints = GetNpoints();
225
226 for (int i = 0; i < nVariables; ++i)
227 {
228 Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
229 }
230 }
231 break;
232 }
233 // Continuous projection
235 {
237
238 for (int i = 0; i < nVariables; ++i)
239 {
240 m_fields[i]->FwdTrans(inarray[i], coeffs);
241 m_fields[i]->BwdTrans(coeffs, outarray[i]);
242 }
243 break;
244 }
245 default:
246 {
247 ASSERTL0(false, "Unknown projection scheme");
248 break;
249 }
250 }
251}
252
253/**
254 * @brief Return the flux vector for the inviscid Burgers equation.
255 *
256 * @param physfield Fields.
257 * @param flux Resulting flux.
258 */
260 const Array<OneD, Array<OneD, NekDouble>> &physfield,
262{
263 const int nq = GetNpoints();
264
265 for (int i = 0; i < flux.size(); ++i)
266 {
267 for (int j = 0; j < flux[0].size(); ++j)
268 {
269 Vmath::Vmul(nq, physfield[i], 1, physfield[i], 1, flux[i][j], 1);
270 Vmath::Smul(nq, 0.5, flux[i][j], 1, flux[i][j], 1);
271 }
272 }
273}
274
276{
278}
279} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void AccumulateRegion(std::string, int iolevel=0)
Accumulate elapsed time for a region.
Definition: Timer.cpp:70
A base class for PDEs which include an advection component.
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Initialisation object for EquationSystem.
int m_spacedim
Spatial dimension (>= expansion dim).
SOLVER_UTILS_EXPORT int GetTraceNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
SOLVER_UTILS_EXPORT int GetNpoints()
SOLVER_UTILS_EXPORT int GetNcoeffs()
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
static SOLVER_UTILS_EXPORT std::vector< ForcingSharedPtr > Load(const LibUtilities::SessionReaderSharedPtr &pSession, const std::weak_ptr< EquationSystem > &pEquation, const Array< OneD, MultiRegions::ExpListSharedPtr > &pFields, const unsigned int &pNumForcingFields=0)
Definition: Forcing.cpp:118
Base class for unsteady solvers.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the projection.
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print Summary.
UnsteadyInviscidBurgers(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Session reader.
void v_InitObject(bool DeclareFields=true) override
Initialise the object.
Array< OneD, NekDouble > & GetNormalVelocity()
Get the normal velocity.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the RHS.
static std::string className
Name of class.
void GetFluxVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Evaluate the flux at each solution point.
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:43
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:46
EquationSystemFactory & GetEquationSystemFactory()
RiemannSolverFactory & GetRiemannSolverFactory()
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
double NekDouble
void Vmul(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Multiply vector z = x*y.
Definition: Vmath.hpp:72
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.hpp:292
void Vvtvp(int n, const T *w, const int incw, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
vvtvp (vector times vector plus vector): z = w*x + y
Definition: Vmath.hpp:366
void Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.hpp:180
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825