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>
38 
39 using namespace std;
40 
41 namespace Nektar
42 {
43  string UnsteadyAdvection::className = SolverUtils::GetEquationSystemFactory().
44  RegisterCreatorFunction("UnsteadyAdvection",
45  UnsteadyAdvection::create,
46  "Unsteady Advection equation.");
47 
48  UnsteadyAdvection::UnsteadyAdvection(
51  : UnsteadySystem(pSession, pGraph),
52  AdvectionSystem(pSession, pGraph)
53  {
54  m_planeNumber = 0;
55  }
56 
57  /**
58  * @brief Initialisation object for the unsteady linear advection equation.
59  */
61  {
62  // Call to the initialisation object of UnsteadySystem
64 
65  m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
66  // Read the advection velocities from session file
67 
68  std::vector<std::string> vel;
69  vel.push_back("Vx");
70  vel.push_back("Vy");
71  vel.push_back("Vz");
72 
73  // Resize the advection velocities vector to dimension of the problem
74  vel.resize(m_spacedim);
75 
76  // Store in the global variable m_velocity the advection velocities
78  GetFunction( "AdvectionVelocity")->Evaluate(vel, m_velocity);
79 
80  // Type of advection class to be used
81  switch(m_projectionType)
82  {
83  // Continuous field
85  {
86  string advName;
87  m_session->LoadSolverInfo(
88  "AdvectionType", advName, "NonConservative");
90  GetAdvectionFactory().CreateInstance(advName, advName);
92  {
93  m_advObject->SetFluxVector(
95  }
96  else
97  {
98  m_advObject->SetFluxVector(
100  }
101  break;
102  }
103  // Discontinuous field
105  {
106  // Do not forwards transform initial condition
107  m_homoInitialFwd = false;
108 
109  // Define the normal velocity fields
110  if (m_fields[0]->GetTrace())
111  {
113  }
114 
115  string advName;
116  string riemName;
117  m_session->LoadSolverInfo(
118  "AdvectionType", advName, "WeakDG");
120  GetAdvectionFactory().CreateInstance(advName, advName);
122  {
123  m_advObject->SetFluxVector(
125  }
126  else
127  {
128  m_advObject->SetFluxVector(
130  }
131  m_session->LoadSolverInfo(
132  "UpwindType", riemName, "Upwind");
135  riemName, m_session);
136  m_riemannSolver->SetScalar(
138 
139  m_advObject->SetRiemannSolver(m_riemannSolver);
140  m_advObject->InitObject(m_session, m_fields);
141  break;
142  }
143  default:
144  {
145  ASSERTL0(false, "Unsupported projection type.");
146  break;
147  }
148  }
149 
150  // If explicit it computes RHS and PROJECTION for the time integration
152  {
155  }
156  // Otherwise it gives an error (no implicit integration)
157  else
158  {
159  ASSERTL0(false, "Implicit unsteady Advection not set up.");
160  }
161  }
162 
163  /**
164  * @brief Unsteady linear advection equation destructor.
165  */
167  {
168  }
169 
170  /**
171  * @brief Get the normal velocity for the linear advection equation.
172  */
174  {
175  // Number of trace (interface) points
176  int i;
177  int nTracePts = GetTraceNpoints();
178 
179  // Auxiliary variable to compute the normal velocity
180  Array<OneD, NekDouble> tmp(nTracePts);
181 
182  // Reset the normal velocity
183  Vmath::Zero(nTracePts, m_traceVn, 1);
184 
185  for (i = 0; i < m_velocity.size(); ++i)
186  {
187  m_fields[0]->ExtractTracePhys(m_velocity[i], tmp);
188 
189  Vmath::Vvtvp(nTracePts,
190  m_traceNormals[i], 1,
191  tmp, 1,
192  m_traceVn, 1,
193  m_traceVn, 1);
194  }
195 
196  return m_traceVn;
197  }
198 
199  /**
200  * @brief Compute the right-hand side for the linear advection equation.
201  *
202  * @param inarray Given fields.
203  * @param outarray Calculated solution.
204  * @param time Time.
205  */
207  const Array<OneD, const Array<OneD, NekDouble> >&inarray,
208  Array<OneD, Array<OneD, NekDouble> >&outarray,
209  const NekDouble time)
210  {
211  // Counter variable
212  int i;
213 
214  // Number of fields (variables of the problem)
215  int nVariables = inarray.size();
216 
217  // Number of solution points
218  int nSolutionPts = GetNpoints();
219 
220 LibUtilities::Timer timer;
221 timer.Start();
222  // RHS computation using the new advection base class
223  m_advObject->Advect(nVariables, m_fields, m_velocity, inarray,
224  outarray, time);
225 timer.Stop();
226 // Elapsed time
227 timer.AccumulateRegion("Advect");
228 
229  // Negate the RHS
230  for (i = 0; i < nVariables; ++i)
231  {
232  Vmath::Neg(nSolutionPts, outarray[i], 1);
233  }
234  }
235 
236  /**
237  * @brief Compute the projection for the linear advection equation.
238  *
239  * @param inarray Given fields.
240  * @param outarray Calculated solution.
241  * @param time Time.
242  */
244  const Array<OneD, const Array<OneD, NekDouble> >&inarray,
245  Array<OneD, Array<OneD, NekDouble> >&outarray,
246  const NekDouble time)
247  {
248  // Counter variable
249  int i;
250 
251  // Number of fields (variables of the problem)
252  int nVariables = inarray.size();
253 
254  // Set the boundary conditions
255  SetBoundaryConditions(time);
256 
257  // Switch on the projection type (Discontinuous or Continuous)
258  switch(m_projectionType)
259  {
260  // Discontinuous projection
262  {
263  // Number of quadrature points
264  int nQuadraturePts = GetNpoints();
265 
266  // Just copy over array
267  for(i = 0; i < nVariables; ++i)
268  {
269  Vmath::Vcopy(nQuadraturePts, inarray[i], 1, outarray[i], 1);
270  }
271  break;
272  }
273 
274  // Continuous projection
277  {
278  Array<OneD, NekDouble> coeffs(m_fields[0]->GetNcoeffs(),0.0);
279  for(i = 0; i < nVariables; ++i)
280  {
281  m_fields[i]->FwdTrans(inarray[i], coeffs);
282  m_fields[i]->BwdTrans_IterPerExp(coeffs, outarray[i]);
283  }
284  break;
285  }
286 
287  default:
288  ASSERTL0(false,"Unknown projection scheme");
289  break;
290  }
291  }
292 
293  /**
294  * @brief Return the flux vector for the linear advection equation.
295  *
296  * @param i Component of the flux vector to calculate.
297  * @param physfield Fields.
298  * @param flux Resulting flux.
299  */
301  const Array<OneD, Array<OneD, NekDouble> > &physfield,
303  {
304  ASSERTL1(flux[0].size() == m_velocity.size(),
305  "Dimension of flux array and velocity array do not match");
306 
307  int i , j;
308  int nq = physfield[0].size();
309 
310  for (i = 0; i < flux.size(); ++i)
311  {
312  for (j = 0; j < flux[0].size(); ++j)
313  {
314  Vmath::Vmul(nq, physfield[i], 1, m_velocity[j], 1,
315  flux[i][j], 1);
316  }
317  }
318  }
319 
320  /**
321  * @brief Return the flux vector for the linear advection equation using
322  * the dealiasing technique.
323  *
324  * @param i Component of the flux vector to calculate.
325  * @param physfield Fields.
326  * @param flux Resulting flux.
327  */
329  const Array<OneD, Array<OneD, NekDouble> > &physfield,
331  {
332  ASSERTL1(flux[0].size() == m_velocity.size(),
333  "Dimension of flux array and velocity array do not match");
334 
335  int i, j;
336  int nq = physfield[0].size();
337  int nVariables = physfield.size();
338 
339  // Factor to rescale 1d points in dealiasing
340  NekDouble OneDptscale = 2;
341 
343  advVel_plane(m_velocity.size());
344 
345  // Get number of points to dealias a cubic non-linearity
346  nq = m_fields[0]->Get1DScaledTotPoints(OneDptscale);
347 
348  // Initialisation of higher-space variables
349  Array<OneD, Array<OneD, NekDouble> >physfieldInterp(nVariables);
351  Array<OneD, Array<OneD, Array<OneD, NekDouble> > >fluxInterp(nVariables);
352 
353  // Interpolation to higher space of physfield
354  for (i = 0; i < nVariables; ++i)
355  {
356  physfieldInterp[i] = Array<OneD, NekDouble>(nq);
357  fluxInterp[i] = Array<OneD, Array<OneD, NekDouble> >(m_expdim);
358  for (j = 0; j < m_expdim; ++j)
359  {
360  fluxInterp[i][j] = Array<OneD, NekDouble>(nq);
361  }
362 
363  m_fields[0]->PhysInterp1DScaled(
364  OneDptscale, physfield[i], physfieldInterp[i]);
365  }
366 
367  // Interpolation to higher space of velocity
368  for (j = 0; j < m_expdim; ++j)
369  {
370  velocityInterp[j] = Array<OneD, NekDouble>(nq);
371 
372  m_fields[0]->PhysInterp1DScaled(
373  OneDptscale, m_velocity[j], velocityInterp[j]);
374  }
375 
376  // Evaluation of flux vector in the higher space
377  for (i = 0; i < flux.size(); ++i)
378  {
379  for (j = 0; j < flux[0].size(); ++j)
380  {
381  Vmath::Vmul(nq, physfieldInterp[i], 1, velocityInterp[j], 1,
382  fluxInterp[i][j], 1);
383  }
384  }
385 
386  // Galerkin project solution back to original space
387  for (i = 0; i < nVariables; ++i)
388  {
389  for (j = 0; j < m_spacedim; ++j)
390  {
391  m_fields[0]->PhysGalerkinProjection1DScaled(
392  OneDptscale, fluxInterp[i][j], flux[i][j]);
393  }
394  }
395  }
396 
398  {
400  }
401 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:250
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:145
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:67
A base class for PDEs which include an advection component.
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
virtual SOLVER_UTILS_EXPORT void v_InitObject()
Init object for UnsteadySystem class.
int m_spacedim
Spatial dimension (>= expansion dim).
SOLVER_UTILS_EXPORT int GetTraceNpoints()
int m_expdim
Expansion dimension.
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()
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.
SOLVER_UTILS_EXPORT SessionFunctionSharedPtr GetFunction(std::string name, const MultiRegions::ExpListSharedPtr &field=MultiRegions::NullExpListSharedPtr, bool cache=false)
Get a SessionFunction by name.
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)
Print a summary of time stepping parameters.
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
virtual void v_InitObject()
Initialise the object.
Array< OneD, NekDouble > m_traceVn
virtual void v_GenerateSummary(SolverUtils::SummaryList &s)
Print Summary.
virtual ~UnsteadyAdvection()
Destructor.
Array< OneD, NekDouble > & GetNormalVelocity()
Get the normal velocity.
Array< OneD, Array< OneD, NekDouble > > m_velocity
Advection velocity.
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.
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
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:46
EquationSystemFactory & GetEquationSystemFactory()
RiemannSolverFactory & GetRiemannSolverFactory()
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:1
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:192
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:461
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:513
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:436
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1199