Nektar++
UnsteadyInviscidBurger.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: UnsteadyInviscidBurger.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 Burger solve routines
32//
33///////////////////////////////////////////////////////////////////////////////
34
36
37using namespace std;
38
39namespace Nektar
40{
43 "UnsteadyInviscidBurger", UnsteadyInviscidBurger::create,
44 "Inviscid Burger equation");
45
49 : UnsteadySystem(pSession, pGraph), AdvectionSystem(pSession, pGraph)
50{
51}
52
53/**
54 * @brief Initialisation object for the inviscid Burger equation.
55 */
57{
58 // Call to the initialisation object of UnsteadySystem
59 AdvectionSystem::v_InitObject(DeclareFields);
60
61 // Define the normal velocity fields
62 if (m_fields[0]->GetTrace())
63 {
65 }
66
67 // Type of advection class to be used
68 switch (m_projectionType)
69 {
70 // Continuous field
72 {
73 string advName;
74 m_session->LoadSolverInfo("AdvectionType", advName,
75 "NonConservative");
77 advName, advName);
79 this);
80 break;
81 }
82 // Discontinuous field
84 {
85 string advName;
86 string riemName;
87
88 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
90 advName, advName);
92 this);
93
94 m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
97 riemName, m_session);
98 m_riemannSolver->SetScalar(
100
101 m_advObject->SetRiemannSolver(m_riemannSolver);
102 m_advObject->InitObject(m_session, m_fields);
103 break;
104 }
105 default:
106 {
107 ASSERTL0(false, "Unsupported projection type.");
108 break;
109 }
110 }
111
112 // Forcing terms
113 m_forcing = SolverUtils::Forcing::Load(m_session, shared_from_this(),
114 m_fields, m_fields.size());
115
116 // If explicit it computes RHS and PROJECTION for the time integration
118 {
121 }
122 // Otherwise it gives an error because (no implicit integration)
123 else
124 {
125 ASSERTL0(false, "Implicit unsteady Advection not set up.");
126 }
127}
128
129/**
130 * @brief Inviscid Burger equation destructor.
131 */
133{
134}
135
136/**
137 * @brief Get the normal velocity for the inviscid Burger equation.
138 */
140{
141 // Counter variable
142 int i;
143
144 // Number of trace (interface) points
145 int nTracePts = GetTraceNpoints();
146
147 // Number of solution points
148 int nSolutionPts = GetNpoints();
149
150 // Number of fields (variables of the problem)
151 int nVariables = m_fields.size();
152
153 // Auxiliary variables to compute the normal velocity
154 Array<OneD, NekDouble> Fwd(nTracePts);
155 Array<OneD, NekDouble> Bwd(nTracePts);
156 Array<OneD, Array<OneD, NekDouble>> physfield(nVariables);
157
158 // Reset the normal velocity
159 Vmath::Zero(nTracePts, m_traceVn, 1);
160
161 // The TimeIntegration Class does not update the physical values of the
162 // solution. It is thus necessary to transform back the coefficient into
163 // the physical space and save them in physfield to compute the normal
164 // advection velocity properly. However it remains a critical point.
165 for (i = 0; i < nVariables; ++i)
166 {
167 physfield[i] = Array<OneD, NekDouble>(nSolutionPts);
168 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), physfield[i]);
169 }
170
171 /// Extract the physical values at the trace space
172 m_fields[0]->GetFwdBwdTracePhys(physfield[0], Fwd, Bwd, true);
173 Vmath::Vadd(nTracePts, Fwd, 1, Bwd, 1, Fwd, 1);
174 Vmath::Smul(nTracePts, 0.5, Fwd, 1, Fwd, 1);
175
176 /// Compute the normal velocity
177 for (i = 0; i < m_spacedim; ++i)
178 {
179 Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, Fwd, 1, m_traceVn, 1,
180 m_traceVn, 1);
181
182 Vmath::Smul(nTracePts, 0.5, m_traceVn, 1, m_traceVn, 1);
183 }
184 return m_traceVn;
185}
186
187/**
188 * @brief Compute the right-hand side for the inviscid Burger equation.
189 *
190 * @param inarray Given fields.
191 * @param outarray Calculated solution.
192 * @param time Time.
193 */
195 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
196 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
197{
198 // Counter variable
199 int i;
200
201 // Number of fields (variables of the problem)
202 int nVariables = inarray.size();
203
204 // Number of solution points
205 int nSolutionPts = GetNpoints();
206
207 // Unused variable for WeakDG and FR
208 Array<OneD, Array<OneD, NekDouble>> advVel(nVariables);
209 for (int i = 0; i < nVariables; ++i)
210 {
211 advVel[i] = inarray[i];
212 }
213
214 // RHS computation using the new advection base class
215 m_advObject->Advect(nVariables, m_fields, advVel, inarray, outarray, time);
216
217 // Negate the RHS
218 for (i = 0; i < nVariables; ++i)
219 {
220 Vmath::Neg(nSolutionPts, outarray[i], 1);
221 }
222
223 // Add forcing terms
224 for (auto &x : m_forcing)
225 {
226 // set up non-linear terms
227 x->Apply(m_fields, inarray, outarray, time);
228 }
229}
230
231/**
232 * @brief Compute the projection for the inviscid Burger equation.
233 *
234 * @param inarray Given fields.
235 * @param outarray Calculated solution.
236 * @param time Time.
237 */
239 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
240 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
241{
242 // Counter variable
243 int i;
244
245 // Number of variables of the problem
246 int nVariables = inarray.size();
247
248 // Set the boundary conditions
250
251 // Switch on the projection type (Discontinuous or Continuous)
252 switch (m_projectionType)
253 {
254 // Discontinuous projection
256 {
257 // Number of quadrature points
258 int nQuadraturePts = GetNpoints();
259
260 // Just copy over array
261 if (inarray != outarray)
262 {
263 for (i = 0; i < nVariables; ++i)
264 {
265 Vmath::Vcopy(nQuadraturePts, inarray[i], 1, outarray[i], 1);
266 }
267 }
268 break;
269 }
270
271 // Continuous projection
274 {
276
277 for (i = 0; i < nVariables; ++i)
278 {
279 m_fields[i]->FwdTrans(inarray[i], coeffs);
280 m_fields[i]->BwdTrans(coeffs, outarray[i]);
281 }
282 break;
283 }
284 default:
285 ASSERTL0(false, "Unknown projection scheme");
286 break;
287 }
288}
289
290/**
291 * @brief Return the flux vector for the inviscid Burger equation.
292 *
293 * @param i Component of the flux vector to calculate.
294 * @param physfield Fields.
295 * @param flux Resulting flux.
296 */
298 const Array<OneD, Array<OneD, NekDouble>> &physfield,
300{
301 const int nq = GetNpoints();
302
303 for (int i = 0; i < flux.size(); ++i)
304 {
305 for (int j = 0; j < flux[0].size(); ++j)
306 {
307 Vmath::Vmul(nq, physfield[i], 1, physfield[i], 1, flux[i][j], 1);
308 Vmath::Smul(nq, 0.5, flux[i][j], 1, flux[i][j], 1);
309 }
310 }
311}
312} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:198
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
A base class for PDEs which include an advection component.
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
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:120
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.
static std::string className
Name of class.
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the projection.
Array< OneD, NekDouble > m_traceVn
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this 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.
virtual 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.
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
UnsteadyInviscidBurger(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Session reader.
std::shared_ptr< SessionReader > SessionReaderSharedPtr
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
EquationSystemFactory & GetEquationSystemFactory()
RiemannSolverFactory & GetRiemannSolverFactory()
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:176
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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.cpp:207
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:513
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.cpp:569
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.cpp:354
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.cpp:245
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191