Nektar++
UnsteadyViscousBurgers.cpp
Go to the documentation of this file.
1 ///////////////////////////////////////////////////////////////////////////////
2 //
3 // File UnsteadyViscousBurgers.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 advection-diffusion solve routines
32 //
33 ///////////////////////////////////////////////////////////////////////////////
34 
35 #include <iostream>
36 
37 #include <boost/core/ignore_unused.hpp>
38 
40 #include <StdRegions/StdQuadExp.h>
41 
42 using namespace std;
43 
44 namespace Nektar
45 {
46  string UnsteadyViscousBurgers::className
48  "UnsteadyViscousBurgers",
49  UnsteadyViscousBurgers::create);
50 
51  UnsteadyViscousBurgers::UnsteadyViscousBurgers(
54  : UnsteadySystem(pSession, pGraph),
55  AdvectionSystem(pSession, pGraph),
56  m_varCoeffLap(StdRegions::NullVarCoeffMap)
57  {
58  m_planeNumber = 0;
59  }
60 
61  /**
62  * @brief Initialisation object for the unsteady linear advection
63  * diffusion equation.
64  */
66  {
68 
69  m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
70  m_session->LoadParameter("epsilon", m_epsilon, 0.0);
71 
72  m_session->MatchSolverInfo(
73  "SpectralVanishingViscosity", "True", m_useSpecVanVisc, false);
74  m_session->MatchSolverInfo(
75  "SpectralVanishingViscosity", "VarDiff", m_useSpecVanViscVarDiff, false);
77  {
78  m_useSpecVanVisc = true;
79  }
80 
82  {
83  m_session->LoadParameter("SVVCutoffRatio",m_sVVCutoffRatio,0.75);
84  m_session->LoadParameter("SVVDiffCoeff",m_sVVDiffCoeff,0.1);
85  }
86 
87  // Type of advection and diffusion classes to be used
88  switch(m_projectionType)
89  {
90  // Discontinuous field
92  {
93  ASSERTL0(false,"Need to implement for DG");
94  // Do not forwards transform initial condition
95  m_homoInitialFwd = false;
96 
97  // Advection term
98  string advName;
99  string riemName;
100  m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
102  CreateInstance(advName, advName);
103  m_advObject->SetFluxVector(&UnsteadyViscousBurgers::
104  GetFluxVectorAdv, this);
105  m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
107  CreateInstance(riemName, m_session);
108  m_advObject->SetRiemannSolver(m_riemannSolver);
109  m_advObject->InitObject (m_session, m_fields);
110 
111  // Diffusion term
112  std::string diffName;
113  m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
115  CreateInstance(diffName, diffName);
116  m_diffusion->SetFluxVector(&UnsteadyViscousBurgers::
117  GetFluxVectorDiff, this);
118  m_diffusion->InitObject(m_session, m_fields);
119  break;
120  }
121  // Continuous field
124  {
125  // Advection term
126  std::string advName;
127  m_session->LoadSolverInfo("AdvectionType", advName,
128  "NonConservative");
130  CreateInstance(advName, advName);
131  m_advObject->SetFluxVector(&UnsteadyViscousBurgers::
132  GetFluxVectorAdv, this);
133 
135  {
136  Array<OneD, Array<OneD, NekDouble> > vel(m_fields.num_elements());
137  for(int i = 0; i < m_fields.num_elements(); ++i)
138  {
139  vel[i] = m_fields[i]->UpdatePhys();
140  }
142  }
143 
144  // In case of Galerkin explicit diffusion gives an error
146  {
147  ASSERTL0(false, "Explicit Galerkin diffusion not set up.");
148  }
149  // In case of Galerkin implicit diffusion: do nothing
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.num_elements());
162 
165 
167  m_explicitDiffusion == 1)
168  {
170  }
171  }
172 
173  /**
174  * @brief Unsteady linear advection diffusion equation destructor.
175  */
177  {
178  }
179 
180  /**
181  * @brief Get the normal velocity for the unsteady linear advection
182  * diffusion equation.
183  */
185  Array<OneD, Array<OneD, NekDouble> >&inarray)
186  {
187  // Number of trace (interface) points
188  int i;
189  int nTracePts = GetTraceNpoints();
190 
191  // Auxiliary variable to compute the normal velocity
192  Array<OneD, NekDouble> tmp(nTracePts);
193  m_traceVn = Array<OneD, NekDouble>(nTracePts, 0.0);
194 
195  // Reset the normal velocity
196  Vmath::Zero(nTracePts, m_traceVn, 1);
197 
198  for (i = 0; i < inarray.num_elements(); ++i)
199  {
200  m_fields[0]->ExtractTracePhys(inarray[i], tmp);
201 
202  Vmath::Vvtvp(nTracePts,
203  m_traceNormals[i], 1,
204  tmp, 1,
205  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 unsteady linear advection
214  * diffusion problem.
215  *
216  * @param inarray Given fields.
217  * @param outarray Calculated solution.
218  * @param time Time.
219  */
221  const Array<OneD, const Array<OneD, NekDouble> >&inarray,
222  Array<OneD, Array<OneD, NekDouble> >&outarray,
223  const NekDouble time)
224  {
225  // Number of fields (variables of the problem)
226  int nVariables = inarray.num_elements();
227 
228  // Number of solution points
229  int nSolutionPts = GetNpoints();
230 
231  Array<OneD, Array<OneD, NekDouble> > outarrayDiff(nVariables);
232 
233  for (int i = 0; i < nVariables; ++i)
234  {
235  outarrayDiff[i] = Array<OneD, NekDouble>(nSolutionPts, 0.0);
236  }
237 
238  // RHS computation using the new advection base class
239  m_advObject->Advect(nVariables, m_fields, inarray,
240  inarray, outarray, time);
241 
242  // Negate the RHS
243  for (int i = 0; i < nVariables; ++i)
244  {
245  Vmath::Neg(nSolutionPts, outarray[i], 1);
246  }
247 
248  // No explicit diffusion for CG
250  {
251  m_diffusion->Diffuse(nVariables, m_fields, inarray, outarrayDiff);
252 
253  for (int i = 0; i < nVariables; ++i)
254  {
255  Vmath::Vadd(nSolutionPts, &outarray[i][0], 1,
256  &outarrayDiff[i][0], 1, &outarray[i][0], 1);
257  }
258  }
259 
260  // Add forcing terms
261  for (auto &x : m_forcing)
262  {
263  // set up non-linear terms
264  x->Apply(m_fields, inarray, outarray, time);
265  }
266  }
267 
268  /**
269  * @brief Compute the projection for the unsteady advection
270  * diffusion problem.
271  *
272  * @param inarray Given fields.
273  * @param outarray Calculated solution.
274  * @param time Time.
275  */
277  const Array<OneD, const Array<OneD, NekDouble> > &inarray,
278  Array<OneD, Array<OneD, NekDouble> > &outarray,
279  const NekDouble time)
280  {
281  int i;
282  int nvariables = inarray.num_elements();
283  SetBoundaryConditions(time);
284  switch(m_projectionType)
285  {
287  {
288  // Just copy over array
289  int npoints = GetNpoints();
290 
291  for(i = 0; i < nvariables; ++i)
292  {
293  Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
294  }
295  break;
296  }
299  {
300  // Do nothing for the moment.
301  }
302  default:
303  {
304  ASSERTL0(false, "Unknown projection scheme");
305  break;
306  }
307  }
308  }
309 
310 
311  /* @brief Compute the diffusion term implicitly.
312  *
313  * @param inarray Given fields.
314  * @param outarray Calculated solution.
315  * @param time Time.
316  * @param lambda Diffusion coefficient.
317  */
319  const Array<OneD, const Array<OneD, NekDouble> >&inarray,
320  Array<OneD, Array<OneD, NekDouble> >&outarray,
321  const NekDouble time,
322  const NekDouble lambda)
323  {
324  int nvariables = inarray.num_elements();
325  int nq = m_fields[0]->GetNpoints();
326 
328  factors[StdRegions::eFactorLambda] = 1.0/lambda/m_epsilon;
329 
330  if(m_useSpecVanVisc)
331  {
334  }
335 
336  Array<OneD, Array< OneD, NekDouble> > F(nvariables);
337  F[0] = Array<OneD, NekDouble> (nq*nvariables);
338 
339  for (int n = 1; n < nvariables; ++n)
340  {
341  F[n] = F[n-1] + nq;
342  }
343 
344  // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
345  // inarray = input: \hat{rhs} -> output: \hat{Y}
346  // outarray = output: nabla^2 \hat{Y}
347  // where \hat = modal coeffs
348  for (int i = 0; i < nvariables; ++i)
349  {
350  // Multiply 1.0/timestep/lambda
352  inarray[i], 1, F[i], 1);
353  }
354 
355  //Setting boundary conditions
356  SetBoundaryConditions(time);
357 
359  {
360  static int cnt = 0;
361 
362  if(cnt %10 == 0)
363  {
364  Array<OneD, Array<OneD, NekDouble> > vel(m_fields.num_elements());
365  for(int i = 0; i < m_fields.num_elements(); ++i)
366  {
367  m_fields[i]->ClearGlobalLinSysManager();
368  vel[i] = m_fields[i]->UpdatePhys();
369  }
371  }
372  ++cnt;
373  }
374  for (int i = 0; i < nvariables; ++i)
375  {
376  // Solve a system of equations with Helmholtz solver
377  m_fields[i]->HelmSolve(F[i], m_fields[i]->UpdateCoeffs(),
378  NullFlagList, factors, m_varCoeffLap);
379 
380  m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
381  }
382  }
383 
384  /**
385  * @brief Return the flux vector for the advection part.
386  *
387  * @param physfield Fields.
388  * @param flux Resulting flux.
389  */
391  const Array<OneD, Array<OneD, NekDouble> > &physfield,
393  {
394 
395  const int nq = m_fields[0]->GetNpoints();
396 
397  for (int i = 0; i < flux.num_elements(); ++i)
398  {
399  for (int j = 0; j < flux[0].num_elements(); ++j)
400  {
401  Vmath::Vmul(nq, physfield[i], 1, physfield[j], 1,
402  flux[i][j], 1);
403  }
404  }
405  }
406 
407  /**
408  * @brief Return the flux vector for the diffusion part.
409  *
410  * @param i Equation number.
411  * @param j Spatial direction.
412  * @param physfield Fields.
413  * @param derivatives First order derivatives.
414  * @param flux Resulting flux.
415  */
417  const Array<OneD, Array<OneD, NekDouble> > &inarray,
418  const Array<OneD, Array<OneD, Array<OneD, NekDouble> > >&qfield,
419  Array<OneD, Array<OneD, Array<OneD, NekDouble> > >&viscousTensor)
420  {
421  boost::ignore_unused(inarray);
422 
423  unsigned int nDim = qfield.num_elements();
424  unsigned int nConvectiveFields = qfield[0].num_elements();
425  unsigned int nPts = qfield[0][0].num_elements();
426 
427  for (unsigned int j = 0; j < nDim; ++j)
428  {
429  for (unsigned int i = 0; i < nConvectiveFields; ++i)
430  {
431  Vmath::Smul(nPts, m_epsilon, qfield[j][i], 1,
432  viscousTensor[j][i], 1 );
433  }
434  }
435  }
436 
439  {
441  if(m_useSpecVanVisc)
442  {
443  stringstream ss;
444  ss << "SVV (cut off = " << m_sVVCutoffRatio
445  << ", coeff = " << m_sVVDiffCoeff << ")";
446  AddSummaryItem(s, "Smoothing", ss.str());
447  }
448  }
449 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:216
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:163
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
bool m_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
LibUtilities::TimeIntegrationSchemeOperators m_ode
The time integration scheme operators to use.
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:46
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:41
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.
Array< OneD, NekDouble > m_traceVn
STL namespace.
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:294
StdRegions::VarCoeffMap m_varCoeffLap
Variable Coefficient map for the Laplacian which can be activated as part of SVV or otherwise...
Array< OneD, Array< OneD, NekDouble > > m_traceNormals
Array holding trace normals for DG simulations in the forwards direction.
virtual ~UnsteadyViscousBurgers()
Destructor.
void DefineProjection(FuncPointerT func, ObjectPointerT obj)
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s)
Print a summary of time stepping parameters.
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*y.
Definition: Vmath.cpp:216
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
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:85
Base class for unsteady solvers.
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()
void GetFluxVectorAdv(const Array< OneD, Array< OneD, NekDouble > > &physfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &flux)
Evaluate the flux at each solution point for the advection part.
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Perform the projection.
AdvectionFactory & GetAdvectionFactory()
Gets the factory for initialising advection objects.
Definition: Advection.cpp:47
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:399
double NekDouble
virtual void v_GenerateSummary(SolverUtils::SummaryList &s)
Print Summary.
EquationSystemFactory & GetEquationSystemFactory()
SOLVER_UTILS_EXPORT void SetBoundaryConditions(NekDouble time)
Evaluates the boundary conditions at the given time.
void GetFluxVectorDiff(const Array< OneD, Array< OneD, NekDouble > > &inarray, const Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &viscousTensor)
Evaluate the flux at each solution point for the diffusion part.
virtual 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 > & GetNormalVelocity(Array< OneD, Array< OneD, NekDouble > > &inarray)
Get the normal velocity.
SolverUtils::DiffusionSharedPtr m_diffusion
SOLVER_UTILS_EXPORT int GetNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
SOLVER_UTILS_EXPORT int GetTraceNpoints()
virtual void v_InitObject()
Initialise the object.
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
Definition: NekFactory.hpp:199
virtual void DoImplicitSolve(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble lambda)
Solve implicitly the diffusion term.
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff(const Array< OneD, Array< OneD, NekDouble >> vel, StdRegions::VarCoeffMap &varCoeffMap)
Evaluate the SVV diffusion coefficient according to Moura&#39;s paper where it should proportional to h t...
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.
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 Vadd(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Add vector z = x+y.
Definition: Vmath.cpp:302
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
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
static FlagList NullFlagList
An empty flag list.
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:265