Nektar++
TimeIntegrationSchemeSDC.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: TimeIntegrationSchemeSDC.h
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: Header file of time integration scheme SDC base class
32//
33///////////////////////////////////////////////////////////////////////////////
34
35#ifndef NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_TIME_INTEGRATION_SCHEME_SDC
36#define NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_TIME_INTEGRATION_SCHEME_SDC
37
38#define LUE LIB_UTILITIES_EXPORT
39
40#include <string>
41
48
49namespace Nektar
50{
51namespace LibUtilities
52{
53
54///////////////////////////////////////////////////////////////////////////////
55/// Class for spectral deferred correction integration.
57{
58public:
59 TimeIntegrationSchemeSDC(std::string variant, size_t order,
60 std::vector<NekDouble> freeParams)
61 : TimeIntegrationScheme(variant, order, freeParams),
62 m_name("SpectralDeferredCorrection")
63 {
64 ASSERTL0(freeParams.size() == 2,
65 "SDC Time integration scheme invalid number "
66 "of free parameters, expected two "
67 "<theta, number of quadrature>, received " +
68 std::to_string(freeParams.size()));
69
70 ASSERTL0(0.0 <= freeParams[0] && freeParams[0] <= 1.0,
71 "Spectral Deferred Correction Time integration "
72 "scheme bad parameter numbers (0.0 - 1.0): " +
73 std::to_string(freeParams[0]));
74
75 // Set scheme properties
76 m_variant = variant;
77 m_order = order;
78 m_freeParams = freeParams;
81
82 if (m_variant == "Equidistant")
83 {
85 m_variant +
86 " quadrature require quadrature "
87 "numbers (>=1" +
88 "): " + std::to_string(freeParams[1]));
89
90 m_first_quadrature = (m_nQuadPts == 1) ? false : true;
91 m_last_quadrature = (m_nQuadPts == 1) ? false : true;
92 m_ordermin = 1;
96 }
97 else if (m_variant == "GaussLobattoLegendre")
98 {
100 m_variant +
101 " quadrature require quadrature "
102 "numbers (>=2" +
103 "): " + std::to_string(freeParams[1]));
104
105 m_first_quadrature = true;
106 m_last_quadrature = true;
107 m_ordermin = 1;
108 m_ordermax = 2 * m_nQuadPts - 2;
111 }
112 else if (m_variant == "GaussRadauLegendre")
113 {
114 ASSERTL0(2 <= m_nQuadPts,
115 m_variant +
116 " quadrature require quadrature "
117 "numbers (>=1" +
118 "): " + std::to_string(freeParams[1]));
119
120 m_first_quadrature = false;
121 m_last_quadrature = true;
122 m_ordermin = 1;
123 m_ordermax = 2 * m_nQuadPts - 1;
126 }
127 else if (m_variant == "GaussGaussLegendre")
128 {
129 ASSERTL0(1 <= m_nQuadPts,
130 m_variant +
131 " quadrature require quadrature "
132 "numbers (>=1" +
133 "): " + std::to_string(freeParams[1]));
134
135 m_first_quadrature = false;
136 m_last_quadrature = false;
137 m_ordermin = 1;
141 }
142 else
143 {
144 ASSERTL0(false, "unknow variant (quadrature) type");
145 }
146
148 "Spectral Deferred Correction Time integration "
149 "scheme bad order numbers (>=" +
150 std::to_string(m_ordermin) +
151 "): " + std::to_string(m_order));
152
154 "Spectral Deferred Correction Time integration "
155 "scheme bad order numbers (<=" +
156 std::to_string(m_ordermax) +
157 "): " + std::to_string(m_order));
158
159 // Add one extra quadrature points for i.c., if necessary
161 {
163 }
164
165 // Get quadrature points and rescale to [0, 1]
167 size_t offset = m_first_quadrature ? 0 : 1;
168 for (size_t i = offset; i < m_nQuadPts; i++)
169 {
170 NekDouble tau = PointsManager()[m_pointsKey]->GetZ()[i - offset];
171 m_tau[i] = (tau + 1.0) / 2.0;
172 }
173 }
174
175 /// Destructor
177 {
178 }
179
181 std::string variant, size_t order, std::vector<NekDouble> freeParams)
182 {
185 variant, order, freeParams);
186
187 return p;
188 }
189
190 static std::string className;
191
192 LUE void SetPFASST(bool pfasst)
193 {
194 m_PFASST = pfasst;
195 }
196
197 LUE void SetTime(double time)
198 {
199 m_time = time;
200 }
201
202 LUE size_t GetMaxOrder() const
203 {
204 return m_ordermax;
205 }
206
208 {
209 return m_first_quadrature;
210 }
211
213 {
214 return m_last_quadrature;
215 }
216
217 LUE size_t GetQuadPtsNumber() const
218 {
219 return m_nQuadPts;
220 }
221
222 LUE size_t GetNpoints() const
223 {
224 return m_npoints;
225 }
226
227 LUE size_t GetNvars() const
228 {
229 return m_nvars;
230 }
231
233 {
234 return m_pointsKey;
235 }
236
238 {
239 return m_Y[0];
240 }
241
243 {
244 return m_Y[0];
245 }
246
248 {
249 return m_last_quadrature ? m_Y[m_nQuadPts - 1] : m_Y_f;
250 }
251
253 {
254 return m_last_quadrature ? m_Y[m_nQuadPts - 1] : m_Y_f;
255 }
256
258 {
259 return m_F;
260 }
261
263 {
264 return m_F;
265 }
266
268 {
269 return m_QFint;
270 }
271
273 {
274 return m_QFint;
275 }
276
278 {
279 return m_FAScorr;
280 }
281
283 {
284 return m_FAScorr;
285 }
286
287 LUE void ResidualEval(const NekDouble &delta_t, const size_t n)
288 {
289 v_ResidualEval(delta_t, n);
290 }
291
292 LUE void ResidualEval(const NekDouble &delta_t)
293 {
294 v_ResidualEval(delta_t);
295 }
296
297 LUE void ComputeInitialGuess(const NekDouble &delta_t)
298 {
299 v_ComputeInitialGuess(delta_t);
300 }
301
302 LUE void SDCIterationLoop(const NekDouble &delta_t)
303 {
304 v_SDCIterationLoop(delta_t);
305 }
306
307 LUE void UpdateFirstQuadrature(void);
308 LUE void UpdateLastQuadrature(void);
309 LUE void AddFASCorrectionToSFint(void);
310 LUE void UpdateIntegratedResidualSFint(const NekDouble &delta_t);
311 LUE void UpdateIntegratedResidualQFint(const NekDouble &delta_t);
312
313protected:
314 LUE virtual std::string v_GetName() const override;
315 LUE virtual std::string v_GetVariant() const override;
316 LUE virtual size_t v_GetOrder() const override;
317 LUE virtual std::vector<NekDouble> v_GetFreeParams() const override;
319 const override;
320 LUE virtual NekDouble v_GetTimeStability() const override;
321 LUE virtual size_t v_GetNumIntegrationPhases() const override;
322
323 /**
324 * \brief Gets the solution vector of the ODE
325 */
326 LUE virtual const TripleArray &v_GetSolutionVector() const override;
327 LUE virtual TripleArray &v_UpdateSolutionVector() override;
328
329 /**
330 * \brief Sets the solution vector of the ODE
331 */
332 LUE virtual void v_SetSolutionVector(const size_t Offset,
333 const DoubleArray &y) override;
334
335 // The worker methods from the base class that are virtual
336 LUE virtual void v_InitializeScheme(
337 const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time,
338 const TimeIntegrationSchemeOperators &op) override;
339
341 const size_t timestep, const NekDouble delta_t) override;
342
343 LUE virtual void v_ResidualEval(const NekDouble &delta_t, const size_t n)
344 {
345 ASSERTL0(false, "Specific version of spectral deferred correction "
346 "not implemented");
347 boost::ignore_unused(delta_t, n);
348 }
349
350 LUE virtual void v_ResidualEval(const NekDouble &delta_t)
351 {
352 ASSERTL0(false, "Specific version of spectral deferred correction "
353 "not implemented");
354 boost::ignore_unused(delta_t);
355 }
356
357 LUE virtual void v_ComputeInitialGuess(const NekDouble &delta_t)
358 {
359 ASSERTL0(false, "Specific version of spectral deferred correction "
360 "not implemented");
361 boost::ignore_unused(delta_t);
362 }
363
364 LUE virtual void v_SDCIterationLoop(const NekDouble &delta_t)
365 {
366 ASSERTL0(false, "Specific version of spectral deferred correction "
367 "not implemented");
368 boost::ignore_unused(delta_t);
369 }
370
371 LUE virtual void v_print(std::ostream &os) const override;
372 LUE virtual void v_printFull(std::ostream &os) const override;
373
374 // Variables common to all schemes
376 std::string m_name;
377 std::string m_variant;
378 std::vector<NekDouble> m_freeParams;
381
382 // Storage of states and associated timesteps
383 PointsKey m_pointsKey; /// Object containing quadrature data
384 SingleArray m_tau; /// Array containing the quadrature points
385 DoubleArray m_Y_f; /// Array containing the last stage values
386 TripleArray m_Y; /// Array containing the stage values
387 TripleArray m_F; /// Array containing the stage derivatives
388 TripleArray m_FAScorr; /// Array containing the FAS correction term
389 TripleArray m_SFint; /// Array containing the integrated residual term
390 TripleArray m_QFint; /// Array containing the integrated residual term
391 SingleArray m_QMat; /// Array containing the integration matrix
392 SingleArray m_interp; /// Array containing the interpolation coefficients
393
394 // SDC parameter
395 NekDouble m_theta{1.0}; /// SDC parameter
396 size_t m_ordermin{0}; /// Minimum order of the integration scheme
397 size_t m_ordermax{0}; /// Maximum order of the integration scheme
398 size_t m_order{0}; /// Order of the integration scheme
399 size_t m_nQuadPts{0}; /// Number of quadrature points
400 size_t m_nvars{0}; /// Number of variables in the integration scheme
401 size_t m_npoints{0}; /// Number of points in the integration scheme
404 bool m_initialized{false};
405 bool m_PFASST{false};
406
407}; // end class TimeIntegrationSchemeSDC
408
409} // end namespace LibUtilities
410} // end namespace Nektar
411
412#endif
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define LUE
Defines a specification for a set of points.
Definition: Points.h:55
Base class for time integration schemes.
Binds a set of functions for use by time integration schemes.
Class for spectral deferred correction integration.
size_t m_ordermax
Minimum order of the integration scheme.
static TimeIntegrationSchemeSharedPtr create(std::string variant, size_t order, std::vector< NekDouble > freeParams)
LUE void UpdateIntegratedResidualSFint(const NekDouble &delta_t)
Worker method that compute residual integral.
virtual LUE size_t v_GetNumIntegrationPhases() const override
LUE const DoubleArray & GetLastQuadratureSolutionVector() const
SingleArray m_interp
Array containing the integration matrix.
virtual LUE void v_ComputeInitialGuess(const NekDouble &delta_t)
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.
LUE const TripleArray & GetIntegratedResidualQFintVector() const
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.
virtual LUE void v_SDCIterationLoop(const NekDouble &delta_t)
DoubleArray m_Y_f
Array containing the quadrature points.
NekDouble m_theta
Array containing the interpolation coefficients.
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.
LUE const DoubleArray & GetFirstQuadratureSolutionVector() const
virtual LUE void v_ResidualEval(const NekDouble &delta_t, const size_t n)
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 void v_ResidualEval(const NekDouble &delta_t)
virtual LUE NekDouble v_GetTimeStability() const override
LUE void ResidualEval(const NekDouble &delta_t, const size_t n)
LUE void SDCIterationLoop(const NekDouble &delta_t)
TimeIntegrationSchemeSDC(std::string variant, size_t order, std::vector< NekDouble > freeParams)
virtual LUE const TripleArray & v_GetSolutionVector() const override
Gets the solution vector of the ODE.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
PointsManagerT & PointsManager(void)
std::shared_ptr< TimeIntegrationScheme > TimeIntegrationSchemeSharedPtr
@ eGaussLobattoLegendre
1D Gauss-Lobatto-Legendre quadrature points
Definition: PointsType.h:53
@ ePolyEvenlySpaced
1D Evenly-spaced points using Lagrange polynomial
Definition: PointsType.h:75
@ eGaussGaussLegendre
1D Gauss-Gauss-Legendre quadrature points
Definition: PointsType.h:48
@ eGaussRadauPLegendre
1D Gauss-Radau-Legendre quadrature points, pinned at x=1
Definition: PointsType.h:51
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble