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
41using namespace std;
42
43namespace Nektar
44{
47 "UnsteadyViscousBurgers", UnsteadyViscousBurgers::create);
48
52 : UnsteadySystem(pSession, pGraph), AdvectionSystem(pSession, pGraph),
53 m_varCoeffLap(StdRegions::NullVarCoeffMap)
54{
55 m_planeNumber = 0;
56}
57
58/**
59 * @brief Initialisation object for the unsteady linear advection
60 * diffusion equation.
61 */
63{
64 AdvectionSystem::v_InitObject(DeclareFields);
65
66 m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
67 m_session->LoadParameter("epsilon", m_epsilon, 0.0);
68
69 m_session->MatchSolverInfo("SpectralVanishingViscosity", "True",
70 m_useSpecVanVisc, false);
71 m_session->MatchSolverInfo("SpectralVanishingViscosity", "VarDiff",
74 {
75 m_useSpecVanVisc = true;
76 }
77
79 {
80 m_session->LoadParameter("SVVCutoffRatio", m_sVVCutoffRatio, 0.75);
81 m_session->LoadParameter("SVVDiffCoeff", m_sVVDiffCoeff, 0.1);
82 }
83
84 // Type of advection and diffusion classes to be used
85 switch (m_projectionType)
86 {
87 // Discontinuous field
89 {
90 ASSERTL0(false, "Need to implement for DG");
91 // Do not forwards transform initial condition
92 m_homoInitialFwd = false;
93
94 // Advection term
95 string advName;
96 string riemName;
97 m_session->LoadSolverInfo("AdvectionType", advName, "WeakDG");
99 advName, advName);
100 m_advObject->SetFluxVector(
102 m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
105 riemName, m_session);
106 m_advObject->SetRiemannSolver(m_riemannSolver);
107 m_advObject->InitObject(m_session, m_fields);
108
109 // Diffusion term
110 std::string diffName;
111 m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
113 diffName, diffName);
114 m_diffusion->SetFluxVector(
116 m_diffusion->InitObject(m_session, m_fields);
117 break;
118 }
119 // Continuous field
122 {
123 // Advection term
124 std::string advName;
125 m_session->LoadSolverInfo("AdvectionType", advName,
126 "NonConservative");
128 advName, advName);
129 m_advObject->SetFluxVector(
131
133 {
135 for (int i = 0; i < m_fields.size(); ++i)
136 {
137 vel[i] = m_fields[i]->UpdatePhys();
138 }
140 }
141
142 // In case of Galerkin explicit diffusion gives an error
144 {
145 ASSERTL0(false, "Explicit Galerkin diffusion not set up.");
146 }
147 // In case of Galerkin implicit diffusion: do nothing
148 break;
149 }
150 default:
151 {
152 ASSERTL0(false, "Unsupported projection type.");
153 break;
154 }
155 }
156
157 // Forcing terms
158 m_forcing = SolverUtils::Forcing::Load(m_session, shared_from_this(),
159 m_fields, m_fields.size());
160
163
166 {
168 }
169}
170
171/**
172 * @brief Unsteady linear advection diffusion equation destructor.
173 */
175{
176}
177
178/**
179 * @brief Get the normal velocity for the unsteady linear advection
180 * diffusion equation.
181 */
184{
185 // Number of trace (interface) points
186 int i;
187 int nTracePts = GetTraceNpoints();
188
189 // Auxiliary variable to compute the normal velocity
190 Array<OneD, NekDouble> tmp(nTracePts);
191 m_traceVn = Array<OneD, NekDouble>(nTracePts, 0.0);
192
193 // Reset the normal velocity
194 Vmath::Zero(nTracePts, m_traceVn, 1);
195
196 for (i = 0; i < inarray.size(); ++i)
197 {
198 m_fields[0]->ExtractTracePhys(inarray[i], tmp);
199
200 Vmath::Vvtvp(nTracePts, m_traceNormals[i], 1, tmp, 1, m_traceVn, 1,
201 m_traceVn, 1);
202 }
203
204 return m_traceVn;
205}
206
207/**
208 * @brief Compute the right-hand side for the unsteady linear advection
209 * diffusion problem.
210 *
211 * @param inarray Given fields.
212 * @param outarray Calculated solution.
213 * @param time Time.
214 */
216 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
217 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
218{
219 // Number of fields (variables of the problem)
220 int nVariables = inarray.size();
221
222 // Number of solution points
223 int nSolutionPts = GetNpoints();
224
225 Array<OneD, Array<OneD, NekDouble>> outarrayDiff(nVariables);
226
227 for (int i = 0; i < nVariables; ++i)
228 {
229 outarrayDiff[i] = Array<OneD, NekDouble>(nSolutionPts, 0.0);
230 }
231
232 // RHS computation using the new advection base class
233 m_advObject->Advect(nVariables, m_fields, inarray, inarray, outarray, time);
234
235 // Negate the RHS
236 for (int i = 0; i < nVariables; ++i)
237 {
238 Vmath::Neg(nSolutionPts, outarray[i], 1);
239 }
240
241 // No explicit diffusion for CG
243 {
244 m_diffusion->Diffuse(nVariables, m_fields, inarray, outarrayDiff);
245
246 for (int i = 0; i < nVariables; ++i)
247 {
248 Vmath::Vadd(nSolutionPts, &outarray[i][0], 1, &outarrayDiff[i][0],
249 1, &outarray[i][0], 1);
250 }
251 }
252
253 // Add forcing terms
254 for (auto &x : m_forcing)
255 {
256 // set up non-linear terms
257 x->Apply(m_fields, inarray, outarray, time);
258 }
259}
260
261/**
262 * @brief Compute the projection for the unsteady advection
263 * diffusion problem.
264 *
265 * @param inarray Given fields.
266 * @param outarray Calculated solution.
267 * @param time Time.
268 */
270 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
271 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
272{
273 int i;
274 int nvariables = inarray.size();
276 switch (m_projectionType)
277 {
279 {
280 // Just copy over array
281 if (inarray != outarray)
282 {
283 int npoints = GetNpoints();
284
285 for (i = 0; i < nvariables; ++i)
286 {
287 Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
288 }
289 }
290 break;
291 }
294 {
295 // Do nothing for the moment.
296 }
297 default:
298 {
299 ASSERTL0(false, "Unknown projection scheme");
300 break;
301 }
302 }
303}
304
305/* @brief Compute the diffusion term implicitly.
306 *
307 * @param inarray Given fields.
308 * @param outarray Calculated solution.
309 * @param time Time.
310 * @param lambda Diffusion coefficient.
311 */
313 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
314 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time,
315 const NekDouble lambda)
316{
317 int nvariables = inarray.size();
318 int nq = m_fields[0]->GetNpoints();
319
322
324 {
327 }
328
330 F[0] = Array<OneD, NekDouble>(nq * nvariables);
331
332 for (int n = 1; n < nvariables; ++n)
333 {
334 F[n] = F[n - 1] + nq;
335 }
336
337 // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
338 // inarray = input: \hat{rhs} -> output: \hat{Y}
339 // outarray = output: nabla^2 \hat{Y}
340 // where \hat = modal coeffs
341 for (int i = 0; i < nvariables; ++i)
342 {
343 // Multiply 1.0/timestep/lambda
344 Vmath::Smul(nq, -factors[StdRegions::eFactorLambda], inarray[i], 1,
345 F[i], 1);
346 }
347
348 // Setting boundary conditions
350
352 {
353 static int cnt = 0;
354
355 if (cnt % 10 == 0)
356 {
358 for (int i = 0; i < m_fields.size(); ++i)
359 {
360 m_fields[i]->ClearGlobalLinSysManager();
361 vel[i] = m_fields[i]->UpdatePhys();
362 }
364 }
365 ++cnt;
366 }
367 for (int i = 0; i < nvariables; ++i)
368 {
369 // Solve a system of equations with Helmholtz solver
370 m_fields[i]->HelmSolve(F[i], m_fields[i]->UpdateCoeffs(), factors,
372
373 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
374 }
375}
376
377/**
378 * @brief Return the flux vector for the advection part.
379 *
380 * @param physfield Fields.
381 * @param flux Resulting flux.
382 */
384 const Array<OneD, Array<OneD, NekDouble>> &physfield,
386{
387
388 const int nq = m_fields[0]->GetNpoints();
389
390 for (int i = 0; i < flux.size(); ++i)
391 {
392 for (int j = 0; j < flux[0].size(); ++j)
393 {
394 Vmath::Vmul(nq, physfield[i], 1, physfield[j], 1, flux[i][j], 1);
395 }
396 }
397}
398
399/**
400 * @brief Return the flux vector for the diffusion part.
401 *
402 * @param i Equation number.
403 * @param j Spatial direction.
404 * @param physfield Fields.
405 * @param derivatives First order derivatives.
406 * @param flux Resulting flux.
407 */
409 const Array<OneD, Array<OneD, NekDouble>> &inarray,
410 const Array<OneD, Array<OneD, Array<OneD, NekDouble>>> &qfield,
411 Array<OneD, Array<OneD, Array<OneD, NekDouble>>> &viscousTensor)
412{
413 boost::ignore_unused(inarray);
414
415 unsigned int nDim = qfield.size();
416 unsigned int nConvectiveFields = qfield[0].size();
417 unsigned int nPts = qfield[0][0].size();
418
419 for (unsigned int j = 0; j < nDim; ++j)
420 {
421 for (unsigned int i = 0; i < nConvectiveFields; ++i)
422 {
423 Vmath::Smul(nPts, m_epsilon, qfield[j][i], 1, viscousTensor[j][i],
424 1);
425 }
426 }
427}
428
430{
433 {
434 stringstream ss;
435 ss << "SVV (cut off = " << m_sVVCutoffRatio
436 << ", coeff = " << m_sVVDiffCoeff << ")";
437 AddSummaryItem(s, "Smoothing", ss.str());
438 }
439}
440} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
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 DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
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.
SOLVER_UTILS_EXPORT int GetTraceNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
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.
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_explicitDiffusion
Indicates if explicit or implicit treatment of diffusion is used.
virtual SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.
SOLVER_UTILS_EXPORT void SVVVarDiffCoeff(const Array< OneD, Array< OneD, NekDouble > > vel, StdRegions::VarCoeffMap &varCoeffMap)
Evaluate the SVV diffusion coefficient according to Moura's paper where it should proportional to h t...
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
Array< OneD, NekDouble > m_traceVn
static std::string className
Name of class.
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.
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
virtual void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print Summary.
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.
Array< OneD, NekDouble > & GetNormalVelocity(Array< OneD, Array< OneD, NekDouble > > &inarray)
Get the normal velocity.
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.
SolverUtils::DiffusionSharedPtr m_diffusion
std::vector< SolverUtils::ForcingSharedPtr > m_forcing
Forcing terms.
UnsteadyViscousBurgers(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Session reader.
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Perform the projection.
virtual void v_InitObject(bool DeclareFields=true) override
Initialise the object.
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the RHS.
StdRegions::VarCoeffMap m_varCoeffLap
Variable Coefficient map for the Laplacian which can be activated as part of SVV or otherwise.
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:48
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:41
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:176
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:408
static VarCoeffMap NullVarCoeffMap
Definition: StdRegions.hpp:353
StdRegions::ConstFactorMap factors
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
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:207
void Neg(int n, T *x, const int incx)
Negate x = -x.
Definition: Vmath.cpp:513
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:569
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:354
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.cpp:245
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.cpp:487
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191