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 viscous Burgers solve routines
32//
33///////////////////////////////////////////////////////////////////////////////
34
36
37namespace Nektar
38{
41 "UnsteadyViscousBurgers", UnsteadyViscousBurgers::create,
42 "Viscous Burgers equation");
43
47 : UnsteadySystem(pSession, pGraph),
48 UnsteadyInviscidBurgers(pSession, pGraph)
49{
50}
51
52/**
53 * @brief Initialisation object for the viscous Burgers equation.
54 */
56{
58
59 m_session->LoadParameter("epsilon", m_epsilon, 1.0);
60 m_session->MatchSolverInfo("SpectralVanishingViscosity", "True",
61 m_useSpecVanVisc, false);
62 m_session->MatchSolverInfo("SpectralVanishingViscosity", "VarDiff",
64
66 {
67 m_useSpecVanVisc = true;
68 }
69
71 {
72 m_session->LoadParameter("SVVCutoffRatio", m_sVVCutoffRatio, 0.75);
73 m_session->LoadParameter("SVVDiffCoeff", m_sVVDiffCoeff, 0.1);
74 }
75
76 // Type of diffusion classes to be used
77 switch (m_projectionType)
78 {
79 // Discontinuous field
81 {
83 {
84 std::string diffName;
85 m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
87 diffName, diffName);
88 m_diffusion->SetFluxVector(
90 m_diffusion->InitObject(m_session, m_fields);
91 }
92 break;
93 }
94 // Continuous field
96 {
97 std::string advName;
98 m_session->LoadSolverInfo("AdvectionType", advName,
99 "NonConservative");
100 if (advName.compare("WeakDG") == 0)
101 {
102 // Define the normal velocity fields
103 if (m_fields[0]->GetTrace())
104 {
106 }
107
108 std::string riemName;
109 m_session->LoadSolverInfo("UpwindType", riemName, "Upwind");
112 riemName, m_session);
113 m_riemannSolver->SetScalar(
115 m_advObject->SetRiemannSolver(m_riemannSolver);
116 m_advObject->InitObject(m_session, m_fields);
117 }
118
119 // In case of Galerkin explicit diffusion gives an error
121 {
122 ASSERTL0(false, "Explicit Galerkin diffusion not set up.");
123 }
124 break;
125 }
126 default:
127 {
128 ASSERTL0(false, "Unsupported projection type.");
129 break;
130 }
131 }
132
135}
136
137/**
138 * @brief Compute the right-hand side for the viscous Burgers equation.
139 *
140 * @param inarray Given fields.
141 * @param outarray Calculated solution.
142 * @param time Time.
143 */
145 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
146 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
147{
148 // Number of fields (variables of the problem)
149 int nVariables = inarray.size();
150
151 // Number of solution points
152 int nSolutionPts = GetNpoints();
153
154 UnsteadyInviscidBurgers::DoOdeRhs(inarray, outarray, time);
155
157 {
158 Array<OneD, Array<OneD, NekDouble>> outarrayDiff(nVariables);
159
160 for (int i = 0; i < nVariables; ++i)
161 {
162 outarrayDiff[i] = Array<OneD, NekDouble>(nSolutionPts, 0.0);
163 }
164
165 m_diffusion->Diffuse(nVariables, m_fields, inarray, outarrayDiff);
166
167 for (int i = 0; i < nVariables; ++i)
168 {
169 Vmath::Vadd(nSolutionPts, &outarrayDiff[i][0], 1, &outarray[i][0],
170 1, &outarray[i][0], 1);
171 }
172 }
173}
174
175/* @brief Compute the diffusion term implicitly.
176 *
177 * @param inarray Given fields.
178 * @param outarray Calculated solution.
179 * @param time Time.
180 * @param lambda Diffusion coefficient.
181 */
183 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
184 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time,
185 const NekDouble lambda)
186{
187 int nvariables = inarray.size();
188 int nq = m_fields[0]->GetNpoints();
189
192
194 {
197 }
198
200 {
202 }
203
205 for (int n = 0; n < nvariables; ++n)
206 {
207 F[n] = Array<OneD, NekDouble>(nq);
208 }
209
210 // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
211 // inarray = input: \hat{rhs} -> output: \hat{Y}
212 // outarray = output: nabla^2 \hat{Y}
213 // where \hat = modal coeffs
214 for (int i = 0; i < nvariables; ++i)
215 {
216 // Multiply 1.0/timestep/lambda
217 Vmath::Smul(nq, -factors[StdRegions::eFactorLambda], inarray[i], 1,
218 F[i], 1);
219 }
220
221 // Setting boundary conditions
223
225 {
226 static int cnt = 0;
227
228 if (cnt % 10 == 0)
229 {
231 for (int i = 0; i < m_fields.size(); ++i)
232 {
233 m_fields[i]->ClearGlobalLinSysManager();
234 vel[i] = m_fields[i]->UpdatePhys();
235 }
237 }
238 ++cnt;
239 }
240
241 for (int i = 0; i < nvariables; ++i)
242 {
243 // Solve a system of equations with Helmholtz solver
244 m_fields[i]->HelmSolve(F[i], m_fields[i]->UpdateCoeffs(), factors,
246
247 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
248 }
249}
250
251/**
252 * @brief Return the flux vector for the diffusion part.
253 *
254 */
256 [[maybe_unused]] const Array<OneD, Array<OneD, NekDouble>> &inarray,
257 const Array<OneD, Array<OneD, Array<OneD, NekDouble>>> &qfield,
258 Array<OneD, Array<OneD, Array<OneD, NekDouble>>> &viscousTensor)
259{
260 unsigned int nDim = qfield.size();
261 unsigned int nConvectiveFields = qfield[0].size();
262 unsigned int nPts = qfield[0][0].size();
263 for (unsigned int j = 0; j < nDim; ++j)
264 {
265 for (unsigned int i = 0; i < nConvectiveFields; ++i)
266 {
267 Vmath::Smul(nPts, m_epsilon, qfield[j][i], 1, viscousTensor[j][i],
268 1);
269 }
270 }
271}
272
274{
277 {
278 std::stringstream ss;
279 ss << "SVV (cut off = " << m_sVVCutoffRatio
280 << ", coeff = " << m_sVVDiffCoeff << ")";
281 SolverUtils::AddSummaryItem(s, "Smoothing", ss.str());
282 }
283}
284} // namespace Nektar
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
tKey RegisterCreatorFunction(tKey idKey, CreatorFunction classCreator, std::string pDesc="")
Register a class with the factory.
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
void DefineOdeRhs(FuncPointerT func, ObjectPointerT obj)
void DefineImplicitSolve(FuncPointerT func, ObjectPointerT obj)
SolverUtils::AdvectionSharedPtr m_advObject
Advection term.
SOLVER_UTILS_EXPORT int GetTraceNpoints()
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
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.
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.
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...
SolverUtils::RiemannSolverSharedPtr m_riemannSolver
void v_GenerateSummary(SolverUtils::SummaryList &s) override
Print Summary.
void v_InitObject(bool DeclareFields=true) override
Initialise the object.
Array< OneD, NekDouble > & GetNormalVelocity()
Get the normal velocity.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the RHS.
static SolverUtils::EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
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.
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.
static std::string className
Name of class.
SolverUtils::DiffusionSharedPtr m_diffusion
UnsteadyViscousBurgers(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
void v_InitObject(bool DeclareFields=true) override
Initialise the object.
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
std::vector< std::pair< std::string, std::string > > SummaryList
Definition: Misc.h:46
DiffusionFactory & GetDiffusionFactory()
Definition: Diffusion.cpp:39
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:47
RiemannSolverFactory & GetRiemannSolverFactory()
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:430
StdRegions::ConstFactorMap factors
double NekDouble
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.hpp:180
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.hpp:100