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 
35 #include <iostream>
37 
38 using namespace std;
39 
40 namespace Nektar
41 {
42  string UnsteadyAdvection::className = SolverUtils::GetEquationSystemFactory().
43  RegisterCreatorFunction("UnsteadyAdvection",
44  UnsteadyAdvection::create,
45  "Unsteady Advection equation.");
46 
47  UnsteadyAdvection::UnsteadyAdvection(
50  : UnsteadySystem(pSession, pGraph),
51  AdvectionSystem(pSession, pGraph)
52  {
53  m_planeNumber = 0;
54  }
55 
56  /**
57  * @brief Initialisation object for the unsteady linear advection equation.
58  */
60  {
61  // Call to the initialisation object of UnsteadySystem
63 
64  m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
65  // Read the advection velocities from session file
66 
67  std::vector<std::string> vel;
68  vel.push_back("Vx");
69  vel.push_back("Vy");
70  vel.push_back("Vz");
71 
72  // Resize the advection velocities vector to dimension of the problem
73  vel.resize(m_spacedim);
74 
75  // Store in the global variable m_velocity the advection velocities
77  GetFunction( "AdvectionVelocity")->Evaluate(vel, m_velocity);
78 
79  // Type of advection class to be used
80  switch(m_projectionType)
81  {
82  // Continuous field
84  {
85  string advName;
86  m_session->LoadSolverInfo(
87  "AdvectionType", advName, "NonConservative");
89  GetAdvectionFactory().CreateInstance(advName, advName);
91  {
92  m_advObject->SetFluxVector(
94  }
95  else
96  {
97  m_advObject->SetFluxVector(
99  }
100  break;
101  }
102  // Discontinuous field
104  {
105  // Do not forwards transform initial condition
106  m_homoInitialFwd = false;
107 
108  // Define the normal velocity fields
109  if (m_fields[0]->GetTrace())
110  {
112  }
113 
114  string advName;
115  string riemName;
116  m_session->LoadSolverInfo(
117  "AdvectionType", advName, "WeakDG");
119  GetAdvectionFactory().CreateInstance(advName, advName);
121  {
122  m_advObject->SetFluxVector(
124  }
125  else
126  {
127  m_advObject->SetFluxVector(
129  }
130  m_session->LoadSolverInfo(
131  "UpwindType", riemName, "Upwind");
134  riemName, m_session);
135  m_riemannSolver->SetScalar(
137 
138  m_advObject->SetRiemannSolver(m_riemannSolver);
139  m_advObject->InitObject(m_session, m_fields);
140  break;
141  }
142  default:
143  {
144  ASSERTL0(false, "Unsupported projection type.");
145  break;
146  }
147  }
148 
149  // If explicit it computes RHS and PROJECTION for the time integration
151  {
154  }
155  // Otherwise it gives an error (no implicit integration)
156  else
157  {
158  ASSERTL0(false, "Implicit unsteady Advection not set up.");
159  }
160  }
161 
162  /**
163  * @brief Unsteady linear advection equation destructor.
164  */
166  {
167  }
168 
169  /**
170  * @brief Get the normal velocity for the linear advection equation.
171  */
173  {
174  // Number of trace (interface) points
175  int i;
176  int nTracePts = GetTraceNpoints();
177 
178  // Auxiliary variable to compute the normal velocity
179  Array<OneD, NekDouble> tmp(nTracePts);
180 
181  // Reset the normal velocity
182  Vmath::Zero(nTracePts, m_traceVn, 1);
183 
184  for (i = 0; i < m_velocity.num_elements(); ++i)
185  {
186  m_fields[0]->ExtractTracePhys(m_velocity[i], tmp);
187 
188  Vmath::Vvtvp(nTracePts,
189  m_traceNormals[i], 1,
190  tmp, 1,
191  m_traceVn, 1,
192  m_traceVn, 1);
193  }
194 
195  return m_traceVn;
196  }
197 
198  /**
199  * @brief Compute the right-hand side for the linear advection 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,
208  const NekDouble time)
209  {
210  // Counter variable
211  int i;
212 
213  // Number of fields (variables of the problem)
214  int nVariables = inarray.num_elements();
215 
216  // Number of solution points
217  int nSolutionPts = GetNpoints();
218 
219  // RHS computation using the new advection base class
220  m_advObject->Advect(nVariables, m_fields, m_velocity, inarray,
221  outarray, time);
222 
223  // Negate the RHS
224  for (i = 0; i < nVariables; ++i)
225  {
226  Vmath::Neg(nSolutionPts, outarray[i], 1);
227  }
228  }
229 
230  /**
231  * @brief Compute the projection for the linear advection equation.
232  *
233  * @param inarray Given fields.
234  * @param outarray Calculated solution.
235  * @param time Time.
236  */
238  const Array<OneD, const Array<OneD, NekDouble> >&inarray,
239  Array<OneD, Array<OneD, NekDouble> >&outarray,
240  const NekDouble time)
241  {
242  // Counter variable
243  int i;
244 
245  // Number of fields (variables of the problem)
246  int nVariables = inarray.num_elements();
247 
248  // Set the boundary conditions
249  SetBoundaryConditions(time);
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  for(i = 0; i < nVariables; ++i)
262  {
263  Vmath::Vcopy(nQuadraturePts, inarray[i], 1, outarray[i], 1);
264  }
265  break;
266  }
267 
268  // Continuous projection
271  {
272  Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs(),0.0);
273  for(i = 0; i < nVariables; ++i)
274  {
275  m_fields[i]->FwdTrans(inarray[i], coeffs);
276  m_fields[i]->BwdTrans_IterPerExp(coeffs, outarray[i]);
277  }
278  break;
279  }
280 
281  default:
282  ASSERTL0(false,"Unknown projection scheme");
283  break;
284  }
285  }
286 
287  /**
288  * @brief Return the flux vector for the linear advection equation.
289  *
290  * @param i Component of the flux vector to calculate.
291  * @param physfield Fields.
292  * @param flux Resulting flux.
293  */
295  const Array<OneD, Array<OneD, NekDouble> > &physfield,
297  {
298  ASSERTL1(flux[0].num_elements() == m_velocity.num_elements(),
299  "Dimension of flux array and velocity array do not match");
300 
301  int i , j;
302  int nq = physfield[0].num_elements();
303 
304  for (i = 0; i < flux.num_elements(); ++i)
305  {
306  for (j = 0; j < flux[0].num_elements(); ++j)
307  {
308  Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1,
309  flux[i][j], 1);
310  }
311  }
312  }
313 
314  /**
315  * @brief Return the flux vector for the linear advection equation using
316  * the dealiasing technique.
317  *
318  * @param i Component of the flux vector to calculate.
319  * @param physfield Fields.
320  * @param flux Resulting flux.
321  */
323  const Array<OneD, Array<OneD, NekDouble> > &physfield,
325  {
326  ASSERTL1(flux[0].num_elements() == m_velocity.num_elements(),
327  "Dimension of flux array and velocity array do not match");
328 
329  int i, j;
330  int nq = physfield[0].num_elements();
331  int nVariables = physfield.num_elements();
332 
333  // Factor to rescale 1d points in dealiasing
334  NekDouble OneDptscale = 2;
335 
337  advVel_plane(m_velocity.num_elements());
338 
339  // Get number of points to dealias a cubic non-linearity
340  nq = m_fields[0]->Get1DScaledTotPoints(OneDptscale);
341 
342  // Initialisation of higher-space variables
343  Array<OneD, Array<OneD, NekDouble> >physfieldInterp(nVariables);
345  Array<OneD, Array<OneD, Array<OneD, NekDouble> > >fluxInterp(nVariables);
346 
347  // Interpolation to higher space of physfield
348  for (i = 0; i < nVariables; ++i)
349  {
350  physfieldInterp[i] = Array<OneD, NekDouble>(nq);
351  fluxInterp[i] = Array<OneD, Array<OneD, NekDouble> >(m_expdim);
352  for (j = 0; j < m_expdim; ++j)
353  {
354  fluxInterp[i][j] = Array<OneD, NekDouble>(nq);
355  }
356 
357  m_fields[0]->PhysInterp1DScaled(
358  OneDptscale, physfield[i], physfieldInterp[i]);
359  }
360 
361  // Interpolation to higher space of velocity
362  for (j = 0; j < m_expdim; ++j)
363  {
364  velocityInterp[j] = Array<OneD, NekDouble>(nq);
365 
366  m_fields[0]->PhysInterp1DScaled(
367  OneDptscale, m_velocity[j], velocityInterp[j]);
368  }
369 
370  // Evaluation of flux vector in the higher space
371  for (i = 0; i < flux.num_elements(); ++i)
372  {
373  for (j = 0; j < flux[0].num_elements(); ++j)
374  {
375  Vmath::Vmul(nq, physfieldInterp[i], 1, velocityInterp[j], 1,
376  fluxInterp[i][j], 1);
377  }
378  }
379 
380  // Galerkin project solution back to original space
381  for (i = 0; i < nVariables; ++i)
382  {
383  for (j = 0; j < m_spacedim; ++j)
384  {
385  m_fields[0]->PhysGalerkinProjection1DScaled(
386  OneDptscale, fluxInterp[i][j], flux[i][j]);
387  }
388  }
389  }
390 
392  {
394  }
395 }
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the projection.
virtual void v_InitObject()
Initialise the object.
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:163
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:46
void GetFluxVector(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Evaluate the flux at each solution point.
int m_expdim
Expansion dimension.
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:445
enum MultiRegions::ProjectionType m_projectionType
Type of projection; e.g continuous or discontinuous.
STL namespace.
bool m_specHP_dealiasing
Flag to determine if dealisising is usde for the Spectral/hp element discretisation.
virtual void v_GenerateSummary(SolverUtils::SummaryList &s)
Print Summary.
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:144
bool m_explicitAdvection
Indicates if explicit or implicit treatment of advection is used.
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
Array< OneD, Array< OneD, NekDouble > > m_velocity
Advection velocity.
Base class for unsteady solvers.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the RHS.
RiemannSolverFactory & GetRiemannSolverFactory()
int m_spacedim
Spatial dimension (>= expansion dim).
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:399
double NekDouble
EquationSystemFactory & GetEquationSystemFactory()
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
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.
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 GetNpoints()
Array< OneD, NekDouble > m_traceVn
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT int GetTraceNpoints()
SOLVER_UTILS_EXPORT int GetNcoeffs()
Array< OneD, NekDouble > & GetNormalVelocity()
Get the normal velocity.
virtual ~UnsteadyAdvection()
Destructor.
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Init object for UnsteadySystem class.
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:376
A base class for PDEs which include an advection component.
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode...
Definition: ErrorUtil.hpp:250
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1064
std::shared_ptr< SessionReader > SessionReaderSharedPtr
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:186