Nektar++
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
36
37using namespace std;
38
39namespace Nektar
40{
43 "UnsteadyDiffusion", UnsteadyDiffusion::create);
44
48 : UnsteadySystem(pSession, pGraph)
49{
50}
51
52/**
53 * @brief Initialisation object for the unsteady diffusion problem.
54 */
55void UnsteadyDiffusion::v_InitObject(bool DeclareFields)
56{
57 UnsteadySystem::v_InitObject(DeclareFields);
58
59 // Load diffusion parameter
60 m_session->LoadParameter("epsilon", m_epsilon, 1.0);
61
62 m_session->MatchSolverInfo("SpectralVanishingViscosity", "True",
63 m_useSpecVanVisc, false);
64
66 {
67 m_session->LoadParameter("SVVCutoffRatio", m_sVVCutoffRatio, 0.75);
68 m_session->LoadParameter("SVVDiffCoeff", m_sVVDiffCoeff, 0.1);
69 }
70
71 int npoints = m_fields[0]->GetNpoints();
72
73 if (m_session->DefinesParameter("d00"))
74 {
75 m_d00 = m_session->GetParameter("d00");
77 Array<OneD, NekDouble>(npoints, m_session->GetParameter("d00"));
78 }
79 if (m_session->DefinesParameter("d11"))
80 {
81 m_d11 = m_session->GetParameter("d11");
83 Array<OneD, NekDouble>(npoints, m_session->GetParameter("d11"));
84 }
85 if (m_session->DefinesParameter("d22"))
86 {
87 m_d22 = m_session->GetParameter("d22");
89 Array<OneD, NekDouble>(npoints, m_session->GetParameter("d22"));
90 }
91
92 switch (m_projectionType)
93 {
95 {
96 std::string diffName;
97
98 // Do not forwards transform initial condition
99 m_homoInitialFwd = false;
100
101 m_session->LoadSolverInfo("DiffusionType", diffName, "LDG");
103 diffName, diffName);
104 m_diffusion->SetFluxVector(&UnsteadyDiffusion::GetFluxVector, this);
105 m_diffusion->InitObject(m_session, m_fields);
106 break;
107 }
109 {
110 // In case of Galerkin explicit diffusion gives an error
112 {
113 ASSERTL0(false, "Explicit Galerkin diffusion not set up.");
114 }
115 // In case of Galerkin implicit diffusion: do nothing
116 break;
117 }
118 default:
119 {
120 ASSERTL0(false, "Unknown projection scheme");
121 break;
122 }
123 }
124
128}
129
131{
134 {
135 stringstream ss;
136 ss << "SVV (cut off = " << m_sVVCutoffRatio
137 << ", coeff = " << m_sVVDiffCoeff << ")";
138 AddSummaryItem(s, "Smoothing", ss.str());
139 }
140}
141
142/* @brief Compute the right-hand side for the unsteady diffusion problem.
143 *
144 * @param inarray Given fields.
145 * @param outarray Calculated solution.
146 * @param time Time.
147 */
149 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
151 [[maybe_unused]] const NekDouble time)
152{
153 // Number of fields (variables of the problem)
154 int nVariables = inarray.size();
155
156 // RHS computation using the new diffusion base class
157 m_diffusion->Diffuse(nVariables, m_fields, inarray, outarray);
158}
159
160/**
161 * @brief Compute the projection for the unsteady diffusion problem.
162 *
163 * @param inarray Given fields.
164 * @param outarray Calculated solution.
165 * @param time Time.
166 */
168 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
169 Array<OneD, Array<OneD, NekDouble>> &outarray, const NekDouble time)
170{
171 int i;
172 int nvariables = inarray.size();
174
175 switch (m_projectionType)
176 {
178 {
179 // Just copy over array
180 if (inarray != outarray)
181 {
182 int npoints = GetNpoints();
183
184 for (i = 0; i < nvariables; ++i)
185 {
186 Vmath::Vcopy(npoints, inarray[i], 1, outarray[i], 1);
187 }
188 }
189 break;
190 }
192 {
194
195 for (i = 0; i < nvariables; ++i)
196 {
197 m_fields[i]->FwdTrans(inarray[i], coeffs);
198 m_fields[i]->BwdTrans(coeffs, outarray[i]);
199 }
200 break;
201 }
202 default:
203 {
204 ASSERTL0(false, "Unknown projection scheme");
205 break;
206 }
207 }
208}
209
210/**
211 * @brief Implicit solution of the unsteady diffusion problem.
212 */
214 const Array<OneD, const Array<OneD, NekDouble>> &inarray,
216 [[maybe_unused]] const NekDouble time, const NekDouble lambda)
217{
219
220 int nvariables = inarray.size();
221 int npoints = m_fields[0]->GetNpoints();
224
226 {
229 }
230
231 // We solve ( \nabla^2 - HHlambda ) Y[i] = rhs [i]
232 // inarray = input: \hat{rhs} -> output: \hat{Y}
233 // outarray = output: nabla^2 \hat{Y}
234 // where \hat = modal coeffs
235 for (int i = 0; i < nvariables; ++i)
236 {
237 // Multiply 1.0/timestep/lambda
238 Vmath::Smul(npoints, -factors[StdRegions::eFactorLambda], inarray[i], 1,
239 outarray[i], 1);
240
241 // Solve a system of equations with Helmholtz solver
242 m_fields[i]->HelmSolve(outarray[i], m_fields[i]->UpdateCoeffs(),
244
245 m_fields[i]->BwdTrans(m_fields[i]->GetCoeffs(), outarray[i]);
246
247 m_fields[i]->SetPhysState(false);
248 }
249}
250
251/**
252 * @brief Return the flux vector for the unsteady diffusion problem.
253 */
255 [[maybe_unused]] const Array<OneD, Array<OneD, NekDouble>> &inarray,
256 const Array<OneD, Array<OneD, Array<OneD, NekDouble>>> &qfield,
257 Array<OneD, Array<OneD, Array<OneD, NekDouble>>> &viscousTensor)
258{
259 unsigned int nDim = qfield.size();
260 unsigned int nConvectiveFields = qfield[0].size();
261 unsigned int nPts = qfield[0][0].size();
262
263 NekDouble d[3] = {m_d00, m_d11, m_d22};
264
265 for (unsigned int j = 0; j < nDim; ++j)
266 {
267 for (unsigned int i = 0; i < nConvectiveFields; ++i)
268 {
269 Vmath::Smul(nPts, m_epsilon * d[j], qfield[j][i], 1,
270 viscousTensor[j][i], 1);
271 }
272 }
273}
274} // 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.
Definition: NekFactory.hpp:197
tBaseSharedPtr CreateInstance(tKey idKey, tParam... args)
Create an instance of the class referred to by idKey.
Definition: NekFactory.hpp:143
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.
SOLVER_UTILS_EXPORT void v_GenerateSummary(SummaryList &s) override
Print a summary of time stepping parameters.
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)
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.
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.
SolverUtils::DiffusionSharedPtr m_diffusion
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
std::shared_ptr< MeshGraph > MeshGraphSharedPtr
Definition: MeshGraph.h:174
std::map< ConstFactorType, NekDouble > ConstFactorMap
Definition: StdRegions.hpp:402
std::vector< double > d(NPUPPER *NPUPPER)
StdRegions::ConstFactorMap factors
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.hpp:100
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825