Nektar++
UnsteadyAdvection.cpp
Go to the documentation of this file.
1/////////////////////////////////////////////////////////////////////////////
2//
3// File: UnsteadyAdvection.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 linear advection solve routines
32//
33///////////////////////////////////////////////////////////////////////////////
34
38
39using namespace std;
40
41namespace Nektar
42{
45 "UnsteadyAdvection", UnsteadyAdvection::create,
46 "Unsteady Advection equation.");
47
51 : UnsteadySystem(pSession, pGraph), AdvectionSystem(pSession, pGraph)
52{
53}
54
55/**
56 * @brief Initialisation object for the unsteady linear advection equation.
57 */
58void UnsteadyAdvection::v_InitObject(bool DeclareFields)
59{
60 // Call to the initialisation object of UnsteadySystem
61 AdvectionSystem::v_InitObject(DeclareFields);
62
63 // Read the advection velocities from session file
64 m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
65
66 // check to see if it is explicity turned off
67 m_session->MatchSolverInfo("GJPStabilisation", "False",
69
70 // if GJPStabilisation set to False bool will be true and
71 // if not false so negate/revese bool
73
74 m_session->LoadParameter("GJPJumpScale", m_GJPJumpScale, 1.0);
75
76 // Define Velocity fields
77 std::vector<std::string> vel;
78 vel.push_back("Vx");
79 vel.push_back("Vy");
80 vel.push_back("Vz");
81
82 // Resize the advection velocities vector to dimension of the problem
83 vel.resize(m_spacedim);
84
85 // Store in the global variable m_velocity the advection velocities
87 GetFunction("AdvectionVelocity")->Evaluate(vel, m_velocity);
88
89 // Type of advection class to be used
90 switch (m_projectionType)
91 {
92 // Discontinuous field
94 {
95 // Do not forwards transform initial condition
96 m_homoInitialFwd = false;
97
98 // Define the normal velocity fields
99 if (m_fields[0]->GetTrace())
100 {
102 }
103
104 string advName;
105 string riemName;
106 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
108 advName, advName);
110 {
111 m_advObject->SetFluxVector(
113 }
114 else
115 {
117 this);
118 }
119 m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
122 riemName, m_session);
123 m_riemannSolver->SetScalar(
125 m_advObject->SetRiemannSolver(m_riemannSolver);
126 m_advObject->InitObject(m_session, m_fields);
127 break;
128 }
129 // Continuous field
132 {
133 string advName;
134 m_session->LoadSolverInfo("AdvectionType", advName,
135 "NonConservative");
137 advName, advName);
139 {
140 m_advObject->SetFluxVector(
142 }
143 else
144 {
146 this);
147 }
148 break;
149 }
150 default:
151 {
152 ASSERTL0(false, "Unsupported projection type.");
153 break;
154 }
155 }
156
157 // Forcing terms
158 m_forcing = SolverUtils::Forcing::Load(m_session, shared_from_this(),
159 m_fields, m_fields.size());
160
161 // If explicit it computes RHS and PROJECTION for the time integration
163 {
166 }
167 // Otherwise it gives an error (no implicit integration)
168 else
169 {
170 ASSERTL0(false, "Implicit unsteady Advection not set up.");
171 }
172}
173
174/**
175 * @brief Get the normal velocity for the linear advection equation.
176 */
178{
180 return m_traceVn;
181}
182
184 const Array<OneD, const Array<OneD, NekDouble>> &velfield)
185{
186 // Number of trace (interface) points
187 int nTracePts = GetTraceNpoints();
188
189 // Auxiliary variable to compute the normal velocity
190 Array<OneD, NekDouble> tmp(nTracePts);
191
192 // Reset the normal velocity
193 Vmath::Zero(nTracePts, m_traceVn, 1);
194
195 for (int i = 0; i < velfield.size(); ++i)
196 {
197 m_fields[0]->ExtractTracePhys(velfield[i], tmp);
198
199 Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, tmp, 1, m_traceVn, 1,
200 m_traceVn, 1);
201 }
202
203 return m_traceVn;
204}
205
206/**
207 * @brief Compute the right-hand side for the linear advection equation.
208 *
209 * @param inarray Given fields.
210 * @param outarray Calculated solution.
211 * @param time Time.
212 */
214 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
215 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
216{
217 // Number of fields (variables of the problem)
218 int nVariables = inarray.size();
219
220 // Number of solution points
221 int nSolutionPts = GetNpoints();
222
224 timer.Start();
225 // RHS computation using the new advection base class
226 m_advObject->Advect(nVariables, m_fields, m_velocity, inarray, outarray,
227 time);
228 timer.Stop();
229 // Elapsed time
230 timer.AccumulateRegion("Advect");
231
232 // Negate the RHS
233 for (int i = 0; i < nVariables; ++i)
234 {
235 Vmath::Neg(nSolutionPts, outarray[i], 1);
236 }
237
238 // Add forcing terms
239 for (auto &x : m_forcing)
240 {
241 // set up non-linear terms
242 x->Apply(m_fields, inarray, outarray, time);
243 }
244}
245
246/**
247 * @brief Compute the projection for the linear advection equation.
248 *
249 * @param inarray Given fields.
250 * @param outarray Calculated solution.
251 * @param time Time.
252 */
254 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
255 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
256{
257 // Number of fields (variables of the problem)
258 int nVariables = inarray.size();
259
260 // Set the boundary conditions
262
263 // Switch on the projection type (Discontinuous or Continuous)
264 switch (m_projectionType)
265 {
266 // Discontinuous projection
268 {
269 // Just copy over array
270 if (inarray != outarray)
271 {
272 int npoints = GetNpoints();
273
274 for (int i = 0; i < nVariables; ++i)
275 {
276 Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
277 }
278 }
279 break;
280 }
281 // Continuous projection
284 {
285 int ncoeffs = m_fields[0]->GetNcoeffs();
286 Array<OneD, NekDouble> coeffs(ncoeffs, 0.0);
288 {
291
292 Array<OneD, NekDouble> wsp(ncoeffs);
293
294 for (int i = 0; i < nVariables; ++i)
295 {
297 std::dynamic_pointer_cast<MultiRegions::ContField>(
298 m_fields[i]);
299
300 m_fields[i]->IProductWRTBase(inarray[i], wsp);
301
303 cfield->GetGJPForcing();
304
307
308 if (GJPData->IsSemiImplicit())
309 {
310 mtype = StdRegions::eMassGJP;
311 }
312
313 // to set up forcing need initial guess in
314 // physical space
316
317 GJPData->Apply(inarray[i], wsp, NullNekDouble1DArray,
318 scale);
319
320 // Solve the system
322 mtype, cfield->GetLocalToGlobalMap(), factors);
323
324 cfield->GlobalSolve(key, wsp, coeffs, NullNekDouble1DArray);
325
326 m_fields[i]->BwdTrans(coeffs, outarray[i]);
327 }
328 }
329 else
330 {
331 for (int i = 0; i < nVariables; ++i)
332 {
333 m_fields[i]->FwdTrans(inarray[i], coeffs);
334 m_fields[i]->BwdTrans(coeffs, outarray[i]);
335 }
336 }
337 break;
338 }
339 default:
340 ASSERTL0(false, "Unknown projection scheme");
341 break;
342 }
343}
344
345/**
346 * @brief Return the flux vector for the linear advection equation.
347 *
348 * @param physfield Fields.
349 * @param flux Resulting flux.
350 */
352 const Array<OneD, Array<OneD, NekDouble>> &physfield,
354{
355 ASSERTL1(flux[0].size() == m_velocity.size(),
356 "Dimension of flux array and velocity array do not match");
357
358 const int nq = m_fields[0]->GetNpoints();
359
360 for (int i = 0; i < flux.size(); ++i)
361 {
362 for (int j = 0; j < flux[0].size(); ++j)
363 {
364 Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1, flux[i][j], 1);
365 }
366 }
367}
368
369/**
370 * @brief Return the flux vector for the linear advection equation using
371 * the dealiasing technique.
372 *
373 * @param physfield Fields.
374 * @param flux Resulting flux.
375 */
377 const Array<OneD, Array<OneD, NekDouble>> &physfield,
379{
380 ASSERTL1(flux[0].size() == m_velocity.size(),
381 "Dimension of flux array and velocity array do not match");
382
383 int nq = physfield[0].size();
384 int nVariables = physfield.size();
385
386 // Factor to rescale 1d points in dealiasing
387 NekDouble OneDptscale = 2;
388
390
391 // Get number of points to dealias a cubic non-linearity
392 nq = m_fields[0]->Get1DScaledTotPoints(OneDptscale);
393
394 // Initialisation of higher-space variables
395 Array<OneD, Array<OneD, NekDouble>> physfieldInterp(nVariables);
397 Array<OneD, Array<OneD, Array<OneD, NekDouble>>> fluxInterp(nVariables);
398
399 // Interpolation to higher space of physfield
400 for (int i = 0; i < nVariables; ++i)
401 {
402 physfieldInterp[i] = Array<OneD, NekDouble>(nq);
404 for (int j = 0; j < m_expdim; ++j)
405 {
406 fluxInterp[i][j] = Array<OneD, NekDouble>(nq);
407 }
408
409 m_fields[0]->PhysInterp1DScaled(OneDptscale, physfield[i],
410 physfieldInterp[i]);
411 }
412
413 // Interpolation to higher space of velocity
414 for (int j = 0; j < m_expdim; ++j)
415 {
416 velocityInterp[j] = Array<OneD, NekDouble>(nq);
417
418 m_fields[0]->PhysInterp1DScaled(OneDptscale, m_velocity[j],
419 velocityInterp[j]);
420 }
421
422 // Evaluation of flux vector in the higher space
423 for (int i = 0; i < flux.size(); ++i)
424 {
425 for (int j = 0; j < flux[0].size(); ++j)
426 {
427 Vmath::Vmul(nq, physfieldInterp[i], 1, velocityInterp[j], 1,
428 fluxInterp[i][j], 1);
429 }
430 }
431
432 // Galerkin project solution back to original space
433 for (int i = 0; i < nVariables; ++i)
434 {
435 for (int j = 0; j < m_spacedim; ++j)
436 {
437 m_fields[0]->PhysGalerkinProjection1DScaled(
438 OneDptscale, fluxInterp[i][j], flux[i][j]);
439 }
440 }
441}
442
444{
447 {
449 s, "GJP Stab. Impl. ",
450 m_session->GetSolverInfo("GJPStabilisation"));
451 SolverUtils::AddSummaryItem(s, "GJP Stab. JumpScale", m_GJPJumpScale);
452 }
453}
454} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:197
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
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()
int m_expdim
Expansion dimension.
NekDouble m_timestep
Time step size.
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
bool m_specHP_dealiasing
Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
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()
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.
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
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.
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
Array< OneD, NekDouble > m_traceVn
Array< OneD, NekDouble > & GetNormalVel(const Array< OneD, const Array< OneD, NekDouble > > &velfield)
Get the normal velocity based on input velfield.
Array< OneD, NekDouble > & GetNormalVelocity()
Get the normal velocity.
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print Summary.
Array< OneD, Array< OneD, NekDouble > > m_velocity
Advection velocity.
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the RHS.
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the projection.
void GetFluxVectorDeAlias(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Evaluate the flux at each solution point using dealiasing.
UnsteadyAdvection(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Session reader.
void v_InitObject(bool DeclareFields=true) override
Initialise the object.
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.
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< GJPStabilisation > GJPStabilisationSharedPtr
std::shared_ptr< ContField > ContFieldSharedPtr
Definition: ContField.h:268
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()
void AddSummaryItem(SummaryList &l, const std::string &name, const std::string &value)
Adds a summary item to the summary info list.
Definition: Misc.cpp:47
RiemannSolverFactory & GetRiemannSolverFactory()
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:402
StdRegions::ConstFactorMap factors
static Array< OneD, NekDouble > NullNekDouble1DArray
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 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