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
37namespace Nektar
38{
39namespace LibUtilities
40{
41
43{
44 return m_name;
45}
46
48{
49 return m_variant;
50}
51
53{
54 return m_order;
55}
56
57std::vector<NekDouble> TimeIntegrationSchemeSDC::v_GetFreeParams() const
58{
59 return m_freeParams;
60}
61
63 const
64{
65 return m_schemeType;
66}
67
69{
70 return 1.0;
71}
72
74{
75 return 1;
76}
77
79{
80 return m_Y;
81}
82
84{
85 return m_Y;
86}
87
89 const DoubleArray &y)
90{
91 m_Y[Offset] = y;
92}
93
94/**
95 * @brief Worker method to initialize the integration scheme.
96 */
98 const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time,
100{
101 boost::ignore_unused(deltaT);
102
103 if (m_initialized)
104 {
105 m_time = time;
106 for (size_t i = 0; i < m_nvars; ++i)
107 {
108 // Store the initial values as the first previous state.
109 Vmath::Vcopy(m_npoints, y_0[i], 1, m_Y[0][i], 1);
110 }
111 }
112 else
113 {
114 m_op = op;
115 m_time = time;
116 m_nvars = y_0.size();
117 m_npoints = y_0[0].size();
118
119 // Compute integration matrix.
120 size_t colOffset = m_first_quadrature ? 0 : 1;
121 size_t rowOffset = m_first_quadrature ? 0 : m_nQuadPts - 1;
122 size_t nCols = m_nQuadPts - colOffset;
123 size_t nRows = m_nQuadPts;
124 m_QMat = SingleArray(nRows * nCols, 0.0);
125 Polylib::Qg(&m_QMat[rowOffset], &m_tau[colOffset], nCols);
126
127 // Compute intepolation coefficient.
129 for (size_t i = 0; i < m_nQuadPts; i++)
130 {
131 m_interp[i] = Polylib::hgj(i, 1.0, &m_tau[0], m_nQuadPts, 0.0, 0.0);
132 }
133
134 // Storage of previous states and associated timesteps.
139 for (size_t m = 0; m < m_nQuadPts; ++m)
140 {
141 m_Y[m] = DoubleArray(m_nvars);
142 m_F[m] = DoubleArray(m_nvars);
145 for (size_t i = 0; i < m_nvars; ++i)
146 {
147 m_Y[m][i] = SingleArray(m_npoints, 0.0);
148 m_F[m][i] = SingleArray(m_npoints, 0.0);
149 m_SFint[m][i] = SingleArray(m_npoints, 0.0);
150 m_QFint[m][i] = SingleArray(m_npoints, 0.0);
151 // Store the initial values as the first previous state.
152 if (m == 0)
153 {
154 Vmath::Vcopy(m_npoints, y_0[i], 1, m_Y[m][i], 1);
155 }
156 }
157 }
158
160 {
162 for (size_t i = 0; i < m_nvars; ++i)
163 {
164 m_Y_f[i] = SingleArray(m_npoints, 0.0);
165 }
166 }
167
168 if (m_PFASST)
169 {
171 for (size_t m = 0; m < m_nQuadPts; ++m)
172 {
174 for (size_t i = 0; i < m_nvars; ++i)
175 {
176 m_FAScorr[m][i] = SingleArray(m_npoints, 0.0);
177 }
178 }
179 }
180
181 m_initialized = true;
182 }
183}
184
185/**
186 * @brief Worker method that performs the time integration.
187 */
189 const size_t timestep, const NekDouble delta_t)
190{
191 boost::ignore_unused(timestep);
192
193 for (size_t k = 0; k < m_order; ++k)
194 {
195 // Compute initial guess
196 if (k == 0)
197 {
198 ComputeInitialGuess(delta_t);
199 }
200 // Apply SDC correction loop
201 else
202 {
203 SDCIterationLoop(delta_t);
204 }
205 }
206
207 // Update last quadrature
209
210 // Update first quadrature
212
213 // Update time step
214 m_time += delta_t;
215
216 // Return solution
217 return m_Y[0];
218}
219
220/**
221 * @brief Worker method that update the first quadrature.
222 */
224{
226 for (size_t i = 0; i < m_nvars; ++i)
227 {
228 Vmath::Vcopy(m_npoints, Y_f[i], 1, m_Y[0][i], 1);
229 }
230}
231
232/**
233 * @brief Worker method that update the last quadrature.
234 */
236{
238 {
239 for (size_t i = 0; i < m_nvars; ++i)
240 {
242 for (size_t n = 0; n < m_nQuadPts; ++n)
243 {
244 Vmath::Svtvp(m_npoints, m_interp[n], m_Y[n][i], 1, m_Y_f[i], 1,
245 m_Y_f[i], 1);
246 }
247 }
248 }
249}
250
251/**
252 * @brief Worker method that add the FASCorrection.
253 */
255{
256 // Add PFASST correction term
257 if (m_PFASST)
258 {
259 for (size_t n = 1; n < m_nQuadPts; ++n)
260 {
261 for (size_t i = 0; i < m_nvars; ++i)
262 {
263 Vmath::Vadd(m_npoints, m_SFint[n][i], 1, m_FAScorr[n][i], 1,
264 m_SFint[n][i], 1);
265 Vmath::Vsub(m_npoints, m_SFint[n][i], 1, m_FAScorr[n - 1][i], 1,
266 m_SFint[n][i], 1);
267 }
268 }
269 }
270}
271
272/**
273 * @brief Worker method that compute residual integral.
274 */
276 const NekDouble &delta_t)
277{
278 // Zeroing integrated flux
279 for (size_t n = 0; n < m_nQuadPts; ++n)
280 {
281 for (size_t i = 0; i < m_nvars; ++i)
282 {
283 Vmath::Zero(m_npoints, m_SFint[n][i], 1);
284 }
285 }
286
287 // Update integrated flux
288 size_t offset = m_first_quadrature ? 0 : 1;
289 for (size_t n = 1; n < m_nQuadPts; ++n)
290 {
291 for (size_t p = 0; p < m_nQuadPts - offset; ++p)
292 {
293 for (size_t i = 0; i < m_nvars; ++i)
294 {
296 m_npoints,
297 delta_t * (m_QMat[n * (m_nQuadPts - offset) + p] -
298 m_QMat[(n - 1) * (m_nQuadPts - offset) + p]),
299 m_F[p + offset][i], 1, m_SFint[n][i], 1, m_SFint[n][i], 1);
300 }
301 }
302 }
303}
304
306 const NekDouble &delta_t)
307{
308 // Zeroing integrated flux
309 for (size_t n = 0; n < m_nQuadPts; ++n)
310 {
311 for (size_t i = 0; i < m_nvars; ++i)
312 {
313 Vmath::Zero(m_npoints, m_QFint[n][i], 1);
314 }
315 }
316
317 // Update integrated flux
318 size_t offset = m_first_quadrature ? 0 : 1;
319 for (size_t n = 0; n < m_nQuadPts; ++n)
320 {
321 for (size_t p = 0; p < m_nQuadPts - offset; ++p)
322 {
323 for (size_t i = 0; i < m_nvars; ++i)
324 {
326 m_npoints, delta_t * m_QMat[n * (m_nQuadPts - offset) + p],
327 m_F[p + offset][i], 1, m_QFint[n][i], 1, m_QFint[n][i], 1);
328 }
329 }
330 }
331}
332
333/**
334 * @brief Worker method to print details on the integration scheme
335 */
336void TimeIntegrationSchemeSDC::v_print(std::ostream &os) const
337{
338 os << "Time Integration Scheme: " << GetFullName() << std::endl;
339}
340
341void TimeIntegrationSchemeSDC::v_printFull(std::ostream &os) const
342{
343 os << "Time Integration Scheme: " << GetFullName() << std::endl;
344}
345
346// Friend Operators
347std::ostream &operator<<(std::ostream &os, const TimeIntegrationSchemeSDC &rhs)
348{
349 rhs.print(os);
350
351 return os;
352}
353
354std::ostream &operator<<(std::ostream &os,
356{
357 os << *rhs.get();
358
359 return os;
360}
361
362} // end namespace LibUtilities
363} // namespace Nektar
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.
virtual LUE size_t v_GetNumIntegrationPhases() const override
SingleArray m_interp
Array containing the integration matrix.
virtual LUE TimeIntegrationSchemeType v_GetIntegrationSchemeType() const override
TripleArray m_Y
Array containing the last stage values.
virtual LUE void v_print(std::ostream &os) const override
Worker method to print details on the integration scheme.
virtual LUE std::string v_GetName() const override
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.
virtual LUE TripleArray & v_UpdateSolutionVector() override
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)
virtual 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.
virtual 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.
virtual LUE void v_SetSolutionVector(const size_t Offset, const DoubleArray &y) override
Sets the solution vector of the ODE.
virtual 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.
virtual LUE ConstDoubleArray & v_TimeIntegrate(const size_t timestep, const NekDouble delta_t) override
Worker method that performs the time integration.
virtual LUE std::string v_GetVariant() const override
virtual LUE NekDouble v_GetTimeStability() const override
LUE void SDCIterationLoop(const NekDouble &delta_t)
virtual 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)
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble
void Qg(double *Q, const double *z, const int np)
Compute the Integration Matrix.
Definition: Polylib.cpp:605
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:878
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.cpp:617
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 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
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.cpp:414