Nektar++
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
UnsteadyDiffusion.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: UnsteadyDiffusion.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 diffusion solve routines
32//
33///////////////////////////////////////////////////////////////////////////////
34
37#include <iomanip>
38#include <iostream>
39
40#include <boost/core/ignore_unused.hpp>
41
42using namespace std;
43
44namespace Nektar
45{
48 "UnsteadyDiffusion", UnsteadyDiffusion::create);
49
53 : UnsteadySystem(pSession, pGraph)
54{
55}
56
57/**
58 * @brief Initialisation object for the unsteady diffusion problem.
59 */
60void UnsteadyDiffusion::v_InitObject(bool DeclareFields)
61{
62 UnsteadySystem::v_InitObject(DeclareFields);
63
64 m_session->LoadParameter("wavefreq", m_waveFreq, 0.0);
65 m_session->LoadParameter("epsilon", m_epsilon, 1.0);
66
67 m_session->MatchSolverInfo("SpectralVanishingViscosity", "True",
68 m_useSpecVanVisc, false);
69
71 {
72 m_session->LoadParameter("SVVCutoffRatio", m_sVVCutoffRatio, 0.75);
73 m_session->LoadParameter("SVVDiffCoeff", m_sVVDiffCoeff, 0.1);
74 }
75
76 int npoints = m_fields[0]->GetNpoints();
77
78 if (m_session->DefinesParameter("d00"))
79 {
81 Array<OneD, NekDouble>(npoints, m_session->GetParameter("d00"));
82 }
83 if (m_session->DefinesParameter("d11"))
84 {
86 Array<OneD, NekDouble>(npoints, m_session->GetParameter("d11"));
87 }
88 if (m_session->DefinesParameter("d22"))
89 {
91 Array<OneD, NekDouble>(npoints, m_session->GetParameter("d22"));
92 }
93
94 switch (m_projectionType)
95 {
97 {
98 std::string diffName;
99
100 // Do not forwards transform initial condition
101 m_homoInitialFwd = false;
102
103 m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
105 diffName, diffName);
106 m_diffusion->SetFluxVector(&UnsteadyDiffusion::GetFluxVector, this);
107 m_diffusion->InitObject(m_session, m_fields);
108 break;
109 }
110
113 {
114 // In case of Galerkin explicit diffusion gives an error
116 {
117 ASSERTL0(false, "Explicit Galerkin diffusion not set up.");
118 }
119 // In case of Galerkin implicit diffusion: do nothing
120 }
121 }
122
124 {
127 }
128 else
129 {
133 }
134}
135
136/**
137 * @brief Unsteady diffusion problem destructor.
138 */
140{
141}
142
144{
147 {
148 stringstream ss;
149 ss << "SVV (cut off = " << m_sVVCutoffRatio
150 << ", coeff = " << m_sVVDiffCoeff << ")";
151 AddSummaryItem(s, "Smoothing", ss.str());
152 }
153}
154
155/* @brief Compute the right-hand side for the unsteady diffusion problem.
156 *
157 * @param inarray Given fields.
158 * @param outarray Calculated solution.
159 * @param time Time.
160 */
162 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
163 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
164{
165 boost::ignore_unused(time);
166
167 // Number of fields (variables of the problem)
168 int nVariables = inarray.size();
169
170 // RHS computation using the new diffusion base class
171 m_diffusion->Diffuse(nVariables, m_fields, inarray, outarray);
172}
173
174/**
175 * @brief Compute the projection for the unsteady diffusion problem.
176 *
177 * @param inarray Given fields.
178 * @param outarray Calculated solution.
179 * @param time Time.
180 */
182 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
183 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
184{
185 int i;
186 int nvariables = inarray.size();
188
189 switch (m_projectionType)
190 {
192 {
193 // Just copy over array
194 if (inarray != outarray)
195 {
196 int npoints = GetNpoints();
197
198 for (i = 0; i < nvariables; ++i)
199 {
200 Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
201 }
202 }
203 break;
204 }
207 {
209
210 for (i = 0; i < nvariables; ++i)
211 {
212 m_fields[i]->FwdTrans(inarray[i], coeffs);
213 m_fields[i]->BwdTrans(coeffs, outarray[i]);
214 }
215 break;
216 }
217 default:
218 {
219 ASSERTL0(false, "Unknown projection scheme");
220 break;
221 }
222 }
223}
224
225/**
226 * @brief Implicit solution of the unsteady diffusion problem.
227 */
229 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
230 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time,
231 const NekDouble lambda)
232{
233 boost::ignore_unused(time);
234
236
237 int nvariables = inarray.size();
238 int npoints = m_fields[0]->GetNpoints();
241
243 {
246 }
247
248 // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
249 // inarray = input: \hat{rhs} -> output: \hat{Y}
250 // outarray = output: nabla^2 \hat{Y}
251 // where \hat = modal coeffs
252 for (int i = 0; i < nvariables; ++i)
253 {
254 // Multiply 1.0/timestep/lambda
255 Vmath::Smul(npoints, -factors[StdRegions::eFactorLambda], inarray[i], 1,
256 outarray[i], 1);
257
258 // Solve a system of equations with Helmholtz solver
259 m_fields[i]->HelmSolve(outarray[i], m_fields[i]->UpdateCoeffs(),
261
262 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
263
264 m_fields[i]->SetPhysState(false);
265 }
266}
267
268/**
269 * @brief Return the flux vector for the unsteady diffusion problem.
270 */
272 const Array<OneD, Array<OneD, NekDouble>> &inarray,
273 const Array<OneD, Array<OneD, Array<OneD, NekDouble>>> &qfield,
274 Array<OneD, Array<OneD, Array<OneD, NekDouble>>> &viscousTensor)
275{
276 boost::ignore_unused(inarray);
277
278 unsigned int nDim = qfield.size();
279 unsigned int nConvectiveFields = qfield[0].size();
280 unsigned int nPts = qfield[0][0].size();
281
282 for (unsigned int j = 0; j < nDim; ++j)
283 {
284 for (unsigned int i = 0; i < nConvectiveFields; ++i)
285 {
286 Vmath::Smul(nPts, m_epsilon, qfield[j][i], 1, viscousTensor[j][i],
287 1);
288 }
289 }
290}
291} // 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)
Array< OneD, MultiRegions::ExpListSharedPtr > m_fields
Array holding all dependent variables.
LibUtilities::SessionReaderSharedPtr m_session
The session reader.
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.
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.
virtual SOLVER_UTILS_EXPORT void v_InitObject(bool DeclareField=true) override
Init object for UnsteadySystem class.
bool m_homoInitialFwd
Flag to determine if simulation should start in homogeneous forward transformed state.
void DoImplicitSolve(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, NekDouble time, NekDouble lambda)
Implicit solution of the unsteady diffusion problem.
void GetFluxVector(const Array< OneD, Array< OneD, NekDouble > > &inarray, const Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &qfield, Array< OneD, Array< OneD, Array< OneD, NekDouble > > > &viscousTensor)
Return the flux vector for the unsteady diffusion problem.
void DoOdeRhs(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
virtual void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.
StdRegions::VarCoeffMap m_varcoeff
static EquationSystemSharedPtr create(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
Creates an instance of this class.
virtual void v_InitObject(bool DeclareFields=true) override
Initialisation object for the unsteady diffusion problem.
UnsteadyDiffusion(const LibUtilities::SessionReaderSharedPtr &pSession, const SpatialDomains::MeshGraphSharedPtr &pGraph)
void DoOdeProjection(const Array< OneD, const Array< OneD, NekDouble > > &inarray, Array< OneD, Array< OneD, NekDouble > > &outarray, const NekDouble time)
Compute the projection for the unsteady diffusion problem.
static std::string className
Name of class.
virtual ~UnsteadyDiffusion()
Destructor.
SolverUtils::DiffusionSharedPtr m_diffusion
std::shared_ptr< SessionReader > SessionReaderSharedPtr
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
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:176
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:408
StdRegions::ConstFactorMap factors
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble
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 Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.cpp:1191