Nektar++
TimeIntegrationSchemeSDC.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: TimeIntegrationSchemeSDC.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: implementation of time integration scheme SDC class
32//
33///////////////////////////////////////////////////////////////////////////////
34
36
38{
39
41{
42 return m_name;
43}
44
46{
47 return m_variant;
48}
49
51{
52 return m_order;
53}
54
55std::vector<NekDouble> TimeIntegrationSchemeSDC::v_GetFreeParams() const
56{
57 return m_freeParams;
58}
59
61 const
62{
63 return m_schemeType;
64}
65
67{
68 return 1.0;
69}
70
72{
73 return 1;
74}
75
77{
78 return m_Y;
79}
80
82{
83 return m_Y;
84}
85
87 const DoubleArray &y)
88{
89 m_Y[Offset] = y;
90}
91
92/**
93 * @brief Worker method to initialize the integration scheme.
94 */
96 [[maybe_unused]] const NekDouble deltaT, ConstDoubleArray &y_0,
97 const NekDouble time, const TimeIntegrationSchemeOperators &op)
98{
99 if (m_initialized)
100 {
101 m_time = time;
102 for (size_t i = 0; i < m_nvars; ++i)
103 {
104 // Store the initial values as the first previous state.
105 Vmath::Vcopy(m_npoints, y_0[i], 1, m_Y[0][i], 1);
106 }
107 }
108 else
109 {
110 m_op = op;
111 m_time = time;
112 m_nvars = y_0.size();
113 m_npoints = y_0[0].size();
114
115 // Compute integration matrix.
116 size_t colOffset = m_first_quadrature ? 0 : 1;
117 size_t rowOffset = m_first_quadrature ? 0 : m_nQuadPts - 1;
118 size_t nCols = m_nQuadPts - colOffset;
119 size_t nRows = m_nQuadPts;
120 m_QMat = SingleArray(nRows * nCols, 0.0);
121 Polylib::Qg(&m_QMat[rowOffset], &m_tau[colOffset], nCols);
122
123 // Compute intepolation coefficient.
125 for (size_t i = 0; i < m_nQuadPts; i++)
126 {
127 m_interp[i] = Polylib::hgj(i, 1.0, &m_tau[0], m_nQuadPts, 0.0, 0.0);
128 }
129
130 // Storage of previous states and associated timesteps.
135 for (size_t m = 0; m < m_nQuadPts; ++m)
136 {
137 m_Y[m] = DoubleArray(m_nvars);
138 m_F[m] = DoubleArray(m_nvars);
141 for (size_t i = 0; i < m_nvars; ++i)
142 {
143 m_Y[m][i] = SingleArray(m_npoints, 0.0);
144 m_F[m][i] = SingleArray(m_npoints, 0.0);
145 m_SFint[m][i] = SingleArray(m_npoints, 0.0);
146 m_QFint[m][i] = SingleArray(m_npoints, 0.0);
147 // Store the initial values as the first previous state.
148 if (m == 0)
149 {
150 Vmath::Vcopy(m_npoints, y_0[i], 1, m_Y[m][i], 1);
151 }
152 }
153 }
154
156 {
158 for (size_t i = 0; i < m_nvars; ++i)
159 {
160 m_Y_f[i] = SingleArray(m_npoints, 0.0);
161 }
162 }
163
164 if (m_PFASST)
165 {
167 for (size_t m = 0; m < m_nQuadPts; ++m)
168 {
170 for (size_t i = 0; i < m_nvars; ++i)
171 {
172 m_FAScorr[m][i] = SingleArray(m_npoints, 0.0);
173 }
174 }
175 }
176
177 m_initialized = true;
178 }
179}
180
181/**
182 * @brief Worker method that performs the time integration.
183 */
185 [[maybe_unused]] const size_t timestep, const NekDouble delta_t)
186{
187 for (size_t k = 0; k < m_order; ++k)
188 {
189 // Compute initial guess
190 if (k == 0)
191 {
192 ComputeInitialGuess(delta_t);
193 }
194 // Apply SDC correction loop
195 else
196 {
197 SDCIterationLoop(delta_t);
198 }
199 }
200
201 // Update last quadrature
203
204 // Update first quadrature
206
207 // Update time step
208 m_time += delta_t;
209
210 // Return solution
211 return m_Y[0];
212}
213
214/**
215 * @brief Worker method that update the first quadrature.
216 */
218{
220 for (size_t i = 0; i < m_nvars; ++i)
221 {
222 Vmath::Vcopy(m_npoints, Y_f[i], 1, m_Y[0][i], 1);
223 }
224}
225
226/**
227 * @brief Worker method that update the last quadrature.
228 */
230{
232 {
233 for (size_t i = 0; i < m_nvars; ++i)
234 {
236 for (size_t n = 0; n < m_nQuadPts; ++n)
237 {
238 Vmath::Svtvp(m_npoints, m_interp[n], m_Y[n][i], 1, m_Y_f[i], 1,
239 m_Y_f[i], 1);
240 }
241 }
242 }
243}
244
245/**
246 * @brief Worker method that add the FASCorrection.
247 */
249{
250 // Add PFASST correction term
251 if (m_PFASST)
252 {
253 for (size_t n = 1; n < m_nQuadPts; ++n)
254 {
255 for (size_t i = 0; i < m_nvars; ++i)
256 {
257 Vmath::Vadd(m_npoints, m_SFint[n][i], 1, m_FAScorr[n][i], 1,
258 m_SFint[n][i], 1);
259 Vmath::Vsub(m_npoints, m_SFint[n][i], 1, m_FAScorr[n - 1][i], 1,
260 m_SFint[n][i], 1);
261 }
262 }
263 }
264}
265
266/**
267 * @brief Worker method that compute residual integral.
268 */
270 const NekDouble &delta_t)
271{
272 // Zeroing integrated flux
273 for (size_t n = 0; n < m_nQuadPts; ++n)
274 {
275 for (size_t i = 0; i < m_nvars; ++i)
276 {
277 Vmath::Zero(m_npoints, m_SFint[n][i], 1);
278 }
279 }
280
281 // Update integrated flux
282 size_t offset = m_first_quadrature ? 0 : 1;
283 for (size_t n = 1; n < m_nQuadPts; ++n)
284 {
285 for (size_t p = 0; p < m_nQuadPts - offset; ++p)
286 {
287 for (size_t i = 0; i < m_nvars; ++i)
288 {
290 m_npoints,
291 delta_t * (m_QMat[n * (m_nQuadPts - offset) + p] -
292 m_QMat[(n - 1) * (m_nQuadPts - offset) + p]),
293 m_F[p + offset][i], 1, m_SFint[n][i], 1, m_SFint[n][i], 1);
294 }
295 }
296 }
297}
298
300 const NekDouble &delta_t)
301{
302 // Zeroing integrated flux
303 for (size_t n = 0; n < m_nQuadPts; ++n)
304 {
305 for (size_t i = 0; i < m_nvars; ++i)
306 {
307 Vmath::Zero(m_npoints, m_QFint[n][i], 1);
308 }
309 }
310
311 // Update integrated flux
312 size_t offset = m_first_quadrature ? 0 : 1;
313 for (size_t n = 0; n < m_nQuadPts; ++n)
314 {
315 for (size_t p = 0; p < m_nQuadPts - offset; ++p)
316 {
317 for (size_t i = 0; i < m_nvars; ++i)
318 {
320 m_npoints, delta_t * m_QMat[n * (m_nQuadPts - offset) + p],
321 m_F[p + offset][i], 1, m_QFint[n][i], 1, m_QFint[n][i], 1);
322 }
323 }
324 }
325}
326
327/**
328 * @brief Worker method to print details on the integration scheme
329 */
330void TimeIntegrationSchemeSDC::v_print(std::ostream &os) const
331{
332 os << "Time Integration Scheme: " << GetFullName() << std::endl;
333}
334
335void TimeIntegrationSchemeSDC::v_printFull(std::ostream &os) const
336{
337 os << "Time Integration Scheme: " << GetFullName() << std::endl;
338}
339
340// Friend Operators
341std::ostream &operator<<(std::ostream &os, const TimeIntegrationSchemeSDC &rhs)
342{
343 rhs.print(os);
344
345 return os;
346}
347
348std::ostream &operator<<(std::ostream &os,
350{
351 os << *rhs.get();
352
353 return os;
354}
355
356} // namespace Nektar::LibUtilities
Binds a set of functions for use by time integration schemes.
Class for spectral deferred correction integration.
LUE void UpdateIntegratedResidualSFint(const NekDouble &delta_t)
Worker method that compute residual integral.
SingleArray m_interp
Array containing the integration matrix.
LUE TimeIntegrationSchemeType v_GetIntegrationSchemeType() const override
TripleArray m_Y
Array containing the last stage values.
LUE void v_print(std::ostream &os) const override
Worker method to print details on the integration scheme.
TripleArray m_QFint
Array containing the integrated residual term.
LUE void UpdateIntegratedResidualQFint(const NekDouble &delta_t)
size_t m_nQuadPts
Order of the integration scheme.
LUE void AddFASCorrectionToSFint(void)
Worker method that add the FASCorrection.
bool m_first_quadrature
Number of points in the integration scheme.
size_t m_npoints
Number of variables in the integration scheme.
SingleArray m_QMat
Array containing the integrated residual term.
LUE void ComputeInitialGuess(const NekDouble &delta_t)
LUE void v_printFull(std::ostream &os) const override
TripleArray m_FAScorr
Array containing the stage derivatives.
TripleArray m_SFint
Array containing the FAS correction term.
DoubleArray m_Y_f
Array containing the quadrature points.
LUE std::vector< NekDouble > v_GetFreeParams() const override
size_t m_order
Maximum order of the integration scheme.
SingleArray m_tau
Object containing quadrature data.
LUE void UpdateFirstQuadrature(void)
Worker method that update the first quadrature.
TripleArray m_F
Array containing the stage values.
LUE void v_SetSolutionVector(const size_t Offset, const DoubleArray &y) override
Sets the solution vector of the ODE.
LUE void v_InitializeScheme(const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time, const TimeIntegrationSchemeOperators &op) override
Worker method to initialize the integration scheme.
LUE void UpdateLastQuadrature(void)
Worker method that update the last quadrature.
LUE ConstDoubleArray & v_TimeIntegrate(const size_t timestep, const NekDouble delta_t) override
Worker method that performs the time integration.
LUE void SDCIterationLoop(const NekDouble &delta_t)
LUE const TripleArray & v_GetSolutionVector() const override
Gets the solution vector of the ODE.
AT< AT< NekDouble > > DoubleArray
std::shared_ptr< TimeIntegrationSchemeSDC > TimeIntegrationSchemeSDCSharedPtr
AT< AT< AT< NekDouble > > > TripleArray
std::ostream & operator<<(std::ostream &os, const BasisKey &rhs)
double NekDouble
void Qg(double *Q, const double *z, const int np)
Compute the Integration Matrix.
Definition: Polylib.cpp:611
double hgj(const int i, const double z, const double *zgj, const int np, const double alpha, const double beta)
Compute the value of the i th Lagrangian interpolant through the np Gauss-Jacobi points zgj at the ar...
Definition: Polylib.cpp:914
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvp (scalar times vector plus vector): z = alpha*x + y.
Definition: Vmath.hpp:396
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 Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.hpp:220