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 
37 #include <MultiRegions/ContField.h>
38 #include <iostream>
39 
40 using namespace std;
41 
42 namespace Nektar
43 {
44 string UnsteadyAdvection::className =
46  "UnsteadyAdvection", UnsteadyAdvection::create,
47  "Unsteady Advection equation.");
48 
49 UnsteadyAdvection::UnsteadyAdvection(
52  : UnsteadySystem(pSession, pGraph), AdvectionSystem(pSession, pGraph)
53 {
54  m_planeNumber = 0;
55 }
56 
57 /**
58  * @brief Initialisation object for the unsteady linear advection equation.
59  */
60 void UnsteadyAdvection::v_InitObject(bool DeclareFields)
61 {
62  // Call to the initialisation object of UnsteadySystem
63  AdvectionSystem::v_InitObject(DeclareFields);
64 
65  m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
66  // Read the advection velocities from session file
67 
68  // check to see if it is explicity turned off
69  m_session->MatchSolverInfo("GJPStabilisation", "False",
70  m_useGJPStabilisation, true);
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  std::vector<std::string> vel;
79  vel.push_back("Vx");
80  vel.push_back("Vy");
81  vel.push_back("Vz");
82 
83  // Resize the advection velocities vector to dimension of the problem
84  vel.resize(m_spacedim);
85 
86  // Store in the global variable m_velocity the advection velocities
88  GetFunction("AdvectionVelocity")->Evaluate(vel, m_velocity);
89 
90  // Type of advection class to be used
91  switch (m_projectionType)
92  {
93  // Continuous field
96  {
97  string advName;
98  m_session->LoadSolverInfo("AdvectionType", advName,
99  "NonConservative");
101  advName, advName);
103  {
104  m_advObject->SetFluxVector(
106  }
107  else
108  {
110  this);
111  }
112  break;
113  }
114  // Discontinuous field
116  {
117  // Do not forwards transform initial condition
118  m_homoInitialFwd = false;
119 
120  // Define the normal velocity fields
121  if (m_fields[0]->GetTrace())
122  {
124  }
125 
126  string advName;
127  string riemName;
128  m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
130  advName, advName);
132  {
133  m_advObject->SetFluxVector(
135  }
136  else
137  {
139  this);
140  }
141  m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
144  riemName, m_session);
145  m_riemannSolver->SetScalar(
147 
148  m_advObject->SetRiemannSolver(m_riemannSolver);
149  m_advObject->InitObject(m_session, m_fields);
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 Unsteady linear advection equation destructor.
178  */
180 {
181 }
182 
183 /**
184  * @brief Get the normal velocity for the linear advection equation.
185  */
187 {
188  // Number of trace (interface) points
189  int i;
190  int nTracePts = GetTraceNpoints();
191 
192  // Auxiliary variable to compute the normal velocity
193  Array<OneD, NekDouble> tmp(nTracePts);
194 
195  // Reset the normal velocity
196  Vmath::Zero(nTracePts, m_traceVn, 1);
197 
198  for (i = 0; i < m_velocity.size(); ++i)
199  {
200  m_fields[0]->ExtractTracePhys(m_velocity[i], tmp);
201 
202  Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, tmp, 1, m_traceVn, 1,
203  m_traceVn, 1);
204  }
205 
206  return m_traceVn;
207 }
208 
209 /**
210  * @brief Compute the right-hand side for the linear advection equation.
211  *
212  * @param inarray Given fields.
213  * @param outarray Calculated solution.
214  * @param time Time.
215  */
217  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
218  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
219 {
220  // Counter variable
221  int i;
222 
223  // Number of fields (variables of the problem)
224  int nVariables = inarray.size();
225 
226  // Number of solution points
227  int nSolutionPts = GetNpoints();
228 
229  LibUtilities::Timer timer;
230  timer.Start();
231  // RHS computation using the new advection base class
232  m_advObject->Advect(nVariables, m_fields, m_velocity, inarray, outarray,
233  time);
234  timer.Stop();
235  // Elapsed time
236  timer.AccumulateRegion("Advect");
237 
238  // Negate the RHS
239  for (i = 0; i < nVariables; ++i)
240  {
241  Vmath::Neg(nSolutionPts, outarray[i], 1);
242  }
243 
244  for (auto &x : m_forcing)
245  {
246  // set up non-linear terms
247  x->Apply(m_fields, inarray, outarray, time);
248  }
249 }
250 
251 /**
252  * @brief Compute the projection for the linear advection equation.
253  *
254  * @param inarray Given fields.
255  * @param outarray Calculated solution.
256  * @param time Time.
257  */
259  const Array<OneD, const Array<OneD, NekDouble>> &inarray,
260  Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
261 {
262  // Counter variable
263  int i;
264 
265  // Number of fields (variables of the problem)
266  int nVariables = inarray.size();
267 
268  // Set the boundary conditions
269  SetBoundaryConditions(time);
270 
271  // Switch on the projection type (Discontinuous or Continuous)
272  switch (m_projectionType)
273  {
274  // Discontinuous projection
276  {
277  // Number of quadrature points
278  int nQuadraturePts = GetNpoints();
279 
280  // Just copy over array
281  for (i = 0; i < nVariables; ++i)
282  {
283  Vmath::Vcopy(nQuadraturePts, inarray[i], 1, outarray[i], 1);
284  }
285  break;
286  }
287 
288  // Continuous projection
291  {
292  int ncoeffs = m_fields[0]->GetNcoeffs();
293  Array<OneD, NekDouble> coeffs(ncoeffs, 0.0);
294 
295 #if 0
296  for(i = 0; i < nVariables; ++i)
297  {
298  m_fields[i]->FwdTrans(inarray[i], coeffs);
299  m_fields[i]->BwdTrans_IterPerExp(coeffs, outarray[i]);
300  }
301 #else
304 
305  Array<OneD, NekDouble> wsp(ncoeffs);
306 
307  for (i = 0; i < nVariables; ++i)
308  {
310  std::dynamic_pointer_cast<MultiRegions::ContField>(
311  m_fields[i]);
312 
313  // copy inarray
314  Array<OneD, NekDouble> in = inarray[i];
315 
316  m_fields[i]->IProductWRTBase(in, wsp);
317 
319  {
321  cfield->GetGJPForcing();
322 
323  factors[StdRegions::eFactorGJP] =
325 
326  if (GJPData->IsSemiImplicit())
327  {
328  mtype = StdRegions::eMassGJP;
329  }
330 
331  // to set up forcing need initial guess in
332  // physical space
333  NekDouble scale = -1.0 * factors[StdRegions::eFactorGJP];
334 
335  GJPData->Apply(inarray[i], wsp, NullNekDouble1DArray,
336  scale);
337  }
338 
339  // Solve the system
341  mtype, cfield->GetLocalToGlobalMap(), factors);
342 
343  cfield->GlobalSolve(key, wsp, coeffs, NullNekDouble1DArray);
344 
345  m_fields[i]->BwdTrans(coeffs, outarray[i]);
346  }
347 #endif
348  break;
349  }
350  default:
351  ASSERTL0(false, "Unknown projection scheme");
352  break;
353  }
354 }
355 
356 /**
357  * @brief Return the flux vector for the linear advection equation.
358  *
359  * @param i Component of the flux vector to calculate.
360  * @param physfield Fields.
361  * @param flux Resulting flux.
362  */
364  const Array<OneD, Array<OneD, NekDouble>> &physfield,
366 {
367  ASSERTL1(flux[0].size() == m_velocity.size(),
368  "Dimension of flux array and velocity array do not match");
369 
370  int i, j;
371  int nq = physfield[0].size();
372 
373  for (i = 0; i < flux.size(); ++i)
374  {
375  for (j = 0; j < flux[0].size(); ++j)
376  {
377  Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1, flux[i][j], 1);
378  }
379  }
380 }
381 
382 /**
383  * @brief Return the flux vector for the linear advection equation using
384  * the dealiasing technique.
385  *
386  * @param i Component of the flux vector to calculate.
387  * @param physfield Fields.
388  * @param flux Resulting flux.
389  */
391  const Array<OneD, Array<OneD, NekDouble>> &physfield,
393 {
394  ASSERTL1(flux[0].size() == m_velocity.size(),
395  "Dimension of flux array and velocity array do not match");
396 
397  int i, j;
398  int nq = physfield[0].size();
399  int nVariables = physfield.size();
400 
401  // Factor to rescale 1d points in dealiasing
402  NekDouble OneDptscale = 2;
403 
404  Array<OneD, Array<OneD, NekDouble>> advVel_plane(m_velocity.size());
405 
406  // Get number of points to dealias a cubic non-linearity
407  nq = m_fields[0]->Get1DScaledTotPoints(OneDptscale);
408 
409  // Initialisation of higher-space variables
410  Array<OneD, Array<OneD, NekDouble>> physfieldInterp(nVariables);
412  Array<OneD, Array<OneD, Array<OneD, NekDouble>>> fluxInterp(nVariables);
413 
414  // Interpolation to higher space of physfield
415  for (i = 0; i < nVariables; ++i)
416  {
417  physfieldInterp[i] = Array<OneD, NekDouble>(nq);
419  for (j = 0; j < m_expdim; ++j)
420  {
421  fluxInterp[i][j] = Array<OneD, NekDouble>(nq);
422  }
423 
424  m_fields[0]->PhysInterp1DScaled(OneDptscale, physfield[i],
425  physfieldInterp[i]);
426  }
427 
428  // Interpolation to higher space of velocity
429  for (j = 0; j < m_expdim; ++j)
430  {
431  velocityInterp[j] = Array<OneD, NekDouble>(nq);
432 
433  m_fields[0]->PhysInterp1DScaled(OneDptscale, m_velocity[j],
434  velocityInterp[j]);
435  }
436 
437  // Evaluation of flux vector in the higher space
438  for (i = 0; i < flux.size(); ++i)
439  {
440  for (j = 0; j < flux[0].size(); ++j)
441  {
442  Vmath::Vmul(nq, physfieldInterp[i], 1, velocityInterp[j], 1,
443  fluxInterp[i][j], 1);
444  }
445  }
446 
447  // Galerkin project solution back to original space
448  for (i = 0; i < nVariables; ++i)
449  {
450  for (j = 0; j < m_spacedim; ++j)
451  {
452  m_fields[0]->PhysGalerkinProjection1DScaled(
453  OneDptscale, fluxInterp[i][j], flux[i][j]);
454  }
455  }
456 }
457 
459 {
462  {
464  s, "GJP Stab. Impl. ",
465  m_session->GetSolverInfo("GJPStabilisation"));
466  SolverUtils::AddSummaryItem(s, "GJP Stab. JumpScale", m_GJPJumpScale);
467  }
468 }
469 } // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
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)
void AccumulateRegion(std::string, int iolevel=0)
Accumulate elapsed time for a region.
Definition: Timer.cpp:73
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()
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: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.
virtual 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.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
Compute the RHS.
Array< OneD, NekDouble > m_traceVn
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble >> &inarray, Array< OneD, Array< OneD, NekDouble >> &outarray, const NekDouble time)
Compute the projection.
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 ~UnsteadyAdvection()
Destructor.
Array< OneD, NekDouble > & GetNormalVelocity()
Get the normal velocity.
virtual void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print Summary.
Array< OneD, Array< OneD, NekDouble > > m_velocity
Advection velocity.
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.
virtual void v_InitObject(bool DeclareFields=true) override
Initialise the object.
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
std::shared_ptr< SessionReader > SessionReaderSharedPtr
std::shared_ptr< GJPStabilisation > GJPStabilisationSharedPtr
std::shared_ptr< ContField > ContFieldSharedPtr
Definition: ContField.h:289
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:48
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:49
RiemannSolverFactory & GetRiemannSolverFactory()
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:172
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:399
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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.cpp:209
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:518
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:574
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:492
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1255