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
39#include <iostream>
40
41using namespace std;
42
43namespace Nektar
44{
47 "UnsteadyAdvection", UnsteadyAdvection::create,
48 "Unsteady Advection equation.");
49
53 : UnsteadySystem(pSession, pGraph), AdvectionSystem(pSession, pGraph)
54{
55}
56
57/**
58 * @brief Initialisation object for the unsteady linear advection equation.
59 */
60void UnsteadyAdvection::v_InitObject(bool DeclareFields)
61{
62 // Call to the initialisation object of UnsteadySystem
63 AdvectionSystem::v_InitObject(DeclareFields);
64
65 // Read the advection velocities from session file
66 m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
67
68 // check to see if it is explicity turned off
69 m_session->MatchSolverInfo("GJPStabilisation", "False",
71
72 // if GJPStabilisation set to False bool will be true and
73 // if not false so negate/revese bool
75
76 m_session->LoadParameter("GJPJumpScale", m_GJPJumpScale, 1.0);
77
78 // Define Velocity fields
79 std::vector<std::string> vel;
80 vel.push_back("Vx");
81 vel.push_back("Vy");
82 vel.push_back("Vz");
83
84 // Resize the advection velocities vector to dimension of the problem
85 vel.resize(m_spacedim);
86
87 // Store in the global variable m_velocity the advection velocities
89 GetFunction("AdvectionVelocity")->Evaluate(vel, m_velocity);
90
91 // Type of advection class to be used
92 switch (m_projectionType)
93 {
94 // Discontinuous field
96 {
97 // Do not forwards transform initial condition
98 m_homoInitialFwd = false;
99
100 // Define the normal velocity fields
101 if (m_fields[0]->GetTrace())
102 {
104 }
105
106 string advName;
107 string riemName;
108 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
110 advName, advName);
112 {
113 m_advObject->SetFluxVector(
115 }
116 else
117 {
119 this);
120 }
121 m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
124 riemName, m_session);
125 m_riemannSolver->SetScalar(
127 m_advObject->SetRiemannSolver(m_riemannSolver);
128 m_advObject->InitObject(m_session, m_fields);
129 break;
130 }
131 // Continuous field
134 {
135 string advName;
136 m_session->LoadSolverInfo("AdvectionType", advName,
137 "NonConservative");
139 advName, advName);
141 {
142 m_advObject->SetFluxVector(
144 }
145 else
146 {
148 this);
149 }
150 break;
151 }
152 default:
153 {
154 ASSERTL0(false, "Unsupported projection type.");
155 break;
156 }
157 }
158
159 // Forcing terms
160 m_forcing = SolverUtils::Forcing::Load(m_session, shared_from_this(),
161 m_fields, m_fields.size());
162
163 // If explicit it computes RHS and PROJECTION for the time integration
165 {
168 }
169 // Otherwise it gives an error (no implicit integration)
170 else
171 {
172 ASSERTL0(false, "Implicit unsteady Advection not set up.");
173 }
174}
175
176/**
177 * @brief Get the normal velocity for the linear advection equation.
178 */
180{
182 return m_traceVn;
183}
184
186 const Array<OneD, const Array<OneD, NekDouble>> &velfield)
187{
188 // Number of trace (interface) points
189 int nTracePts = GetTraceNpoints();
190 int nPts = m_velocity[0].size();
191
192 // Auxiliary variable to compute the normal velocity
193 Array<OneD, NekDouble> tmp(nPts), tmp2(nTracePts);
194
195 // Reset the normal velocity
196 Vmath::Zero(nTracePts, m_traceVn, 1);
197
198 for (int i = 0; i < velfield.size(); ++i)
199 {
200 // velocity - grid velocity for ALE before getting trace velocity
201 Vmath::Vsub(nPts, velfield[i], 1, m_gridVelocity[i], 1, tmp, 1);
202
203 m_fields[0]->ExtractTracePhys(tmp, tmp2);
204
205 Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, tmp2, 1, m_traceVn, 1,
206 m_traceVn, 1);
207 }
208
209 return m_traceVn;
210}
211
212/**
213 * @brief Compute the right-hand side for the linear advection equation.
214 *
215 * @param inarray Given fields.
216 * @param outarray Calculated solution.
217 * @param time Time.
218 */
220 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
221 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
222{
223 // Number of fields (variables of the problem)
224 int nVariables = inarray.size();
225
227 if (m_ALESolver)
228 {
229 timer.Start();
230 Array<OneD, Array<OneD, NekDouble>> tmpIn(nVariables);
231 // If ALE we must take Mu coefficient space to u physical space
233 auto advWeakDGObject =
234 std::dynamic_pointer_cast<SolverUtils::AdvectionWeakDG>(
236 advWeakDGObject->AdvectCoeffs(nVariables, m_fields, m_velocity, tmpIn,
237 outarray, time);
238 timer.Stop();
239 }
240 else
241 {
242 timer.Start();
243 m_advObject->Advect(nVariables, m_fields, m_velocity, inarray, outarray,
244 time);
245 timer.Stop();
246 }
247
248 // Elapsed time
249 timer.AccumulateRegion("Advect");
250
251 // Negate the RHS
252 for (int i = 0; i < nVariables; ++i)
253 {
254 Vmath::Neg(outarray[i].size(), outarray[i], 1);
255 }
256
257 // Add forcing terms
258 for (auto &x : m_forcing)
259 {
260 // set up non-linear terms
261 x->Apply(m_fields, inarray, outarray, time);
262 }
263}
264
265/**
266 * @brief Compute the projection for the linear advection equation.
267 *
268 * @param inarray Given fields.
269 * @param outarray Calculated solution.
270 * @param time Time.
271 */
273 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
274 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
275{
276 // Number of fields (variables of the problem)
277 int nVariables = inarray.size();
278
279 // Perform ALE movement
280 if (m_ALESolver)
281 {
283 }
284
285 // Set the boundary conditions
287
288 // Switch on the projection type (Discontinuous or Continuous)
289 switch (m_projectionType)
290 {
291 // Discontinuous projection
293 {
294 // Just copy over array
295 if (inarray != outarray)
296 {
297 int npoints = GetNpoints();
298
299 for (int i = 0; i < nVariables; ++i)
300 {
301 Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
302 }
303 }
304 break;
305 }
306 // Continuous projection
309 {
310 int ncoeffs = m_fields[0]->GetNcoeffs();
311 Array<OneD, NekDouble> coeffs(ncoeffs, 0.0);
313 {
316
317 Array<OneD, NekDouble> wsp(ncoeffs);
318
319 for (int i = 0; i < nVariables; ++i)
320 {
322 std::dynamic_pointer_cast<MultiRegions::ContField>(
323 m_fields[i]);
324
325 m_fields[i]->IProductWRTBase(inarray[i], wsp);
326
328 cfield->GetGJPForcing();
329
332
333 if (GJPData->IsSemiImplicit())
334 {
335 mtype = StdRegions::eMassGJP;
336 }
337
338 // to set up forcing need initial guess in
339 // physical space
341
342 GJPData->Apply(inarray[i], wsp, NullNekDouble1DArray,
343 scale);
344
345 // Solve the system
347 mtype, cfield->GetLocalToGlobalMap(), factors);
348
349 cfield->GlobalSolve(key, wsp, coeffs, NullNekDouble1DArray);
350
351 m_fields[i]->BwdTrans(coeffs, outarray[i]);
352 }
353 }
354 else
355 {
356 for (int i = 0; i < nVariables; ++i)
357 {
358 m_fields[i]->FwdTrans(inarray[i], coeffs);
359 m_fields[i]->BwdTrans(coeffs, outarray[i]);
360 }
361 }
362 break;
363 }
364 default:
365 ASSERTL0(false, "Unknown projection scheme");
366 break;
367 }
368}
369
370/**
371 * @brief Return the flux vector for the linear advection equation.
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 const int nq = m_fields[0]->GetNpoints();
384
385 for (int i = 0; i < flux.size(); ++i)
386 {
387 for (int j = 0; j < flux[0].size(); ++j)
388 {
389 for (int k = 0; k < nq; ++k)
390 {
391 // If ALE we need to take off the grid velocity
392 flux[i][j][k] =
393 physfield[i][k] * (m_velocity[j][k] - m_gridVelocity[j][k]);
394 }
395 }
396 }
397}
398
399/**
400 * @brief Return the flux vector for the linear advection equation using
401 * the dealiasing technique.
402 *
403 * @param physfield Fields.
404 * @param flux Resulting flux.
405 */
407 const Array<OneD, Array<OneD, NekDouble>> &physfield,
409{
410 ASSERTL1(flux[0].size() == m_velocity.size(),
411 "Dimension of flux array and velocity array do not match");
412
413 int nq = physfield[0].size();
414 int nVariables = physfield.size();
415
416 // Factor to rescale 1d points in dealiasing
417 NekDouble OneDptscale = 2;
418
420
421 // Get number of points to dealias a cubic non-linearity
422 nq = m_fields[0]->Get1DScaledTotPoints(OneDptscale);
423
424 // Initialisation of higher-space variables
425 Array<OneD, Array<OneD, NekDouble>> physfieldInterp(nVariables);
427 Array<OneD, Array<OneD, Array<OneD, NekDouble>>> fluxInterp(nVariables);
428
429 // Interpolation to higher space of physfield
430 for (int i = 0; i < nVariables; ++i)
431 {
432 physfieldInterp[i] = Array<OneD, NekDouble>(nq);
434 for (int j = 0; j < m_expdim; ++j)
435 {
436 fluxInterp[i][j] = Array<OneD, NekDouble>(nq);
437 }
438
439 m_fields[0]->PhysInterp1DScaled(OneDptscale, physfield[i],
440 physfieldInterp[i]);
441 }
442
443 // Interpolation to higher space of velocity
444 for (int j = 0; j < m_expdim; ++j)
445 {
446 velocityInterp[j] = Array<OneD, NekDouble>(nq);
447
448 m_fields[0]->PhysInterp1DScaled(OneDptscale, m_velocity[j],
449 velocityInterp[j]);
450 }
451
452 // Evaluation of flux vector in the higher space
453 for (int i = 0; i < flux.size(); ++i)
454 {
455 for (int j = 0; j < flux[0].size(); ++j)
456 {
457 Vmath::Vmul(nq, physfieldInterp[i], 1, velocityInterp[j], 1,
458 fluxInterp[i][j], 1);
459 }
460 }
461
462 // Galerkin project solution back to original space
463 for (int i = 0; i < nVariables; ++i)
464 {
465 for (int j = 0; j < m_spacedim; ++j)
466 {
467 m_fields[0]->PhysGalerkinProjection1DScaled(
468 OneDptscale, fluxInterp[i][j], flux[i][j]);
469 }
470 }
471}
472
474{
477 {
479 s, "GJP Stab. Impl. ",
480 m_session->GetSolverInfo("GJPStabilisation"));
481 SolverUtils::AddSummaryItem(s, "GJP Stab. JumpScale", m_GJPJumpScale);
482 }
483}
484
486{
487 boost::ignore_unused(step);
488 return false;
489}
490
492 std::vector<Array<OneD, NekDouble>> &fieldcoeffs,
493 std::vector<std::string> &variables)
494{
495 bool extraFields;
496 m_session->MatchSolverInfo("OutputExtraFields", "True", extraFields, true);
497
498 if (extraFields && m_ALESolver)
499 {
500 ExtraFldOutputGridVelocity(fieldcoeffs, variables);
501 }
502}
503
506{
508 {
509 m_spaceDim = spaceDim;
510 m_fieldsALE = fields;
511
512 // Initialise grid velocities as 0s
515 for (int i = 0; i < spaceDim; ++i)
516 {
517 m_gridVelocity[i] =
518 Array<OneD, NekDouble>(fields[0]->GetTotPoints(), 0.0);
520 fields[0]->GetTrace()->GetTotPoints(), 0.0);
521 }
522 }
523 ALEHelper::InitObject(spaceDim, fields);
524}
525
526} // 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.
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
Array< OneD, MultiRegions::ExpListSharedPtr > m_fieldsALE
Definition: ALEHelper.h:89
SOLVER_UTILS_EXPORT void InitObject(int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields)
Definition: ALEHelper.cpp:48
Array< OneD, Array< OneD, NekDouble > > m_gridVelocityTrace
Definition: ALEHelper.h:91
SOLVER_UTILS_EXPORT void ALEDoElmtInvMassBwdTrans(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray)
Definition: ALEHelper.cpp:149
SOLVER_UTILS_EXPORT void ExtraFldOutputGridVelocity(std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables)
Definition: ALEHelper.cpp:392
SOLVER_UTILS_EXPORT void MoveMesh(const NekDouble &time, Array< OneD, Array< OneD, NekDouble > > &traceNormals)
Definition: ALEHelper.cpp:169
Array< OneD, Array< OneD, NekDouble > > m_gridVelocity
Definition: ALEHelper.h:90
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.
SOLVER_UTILS_EXPORT int GetTotPoints()
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.
void v_ALEInitObject(int spaceDim, Array< OneD, MultiRegions::ExpListSharedPtr > &fields) override
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
void v_ExtraFldOutput(std::vector< Array< OneD, NekDouble > > &fieldcoeffs, std::vector< std::string > &variables) override
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
bool v_PreIntegrate(int step) override
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:430
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
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.hpp:220
STL namespace.