Nektar++
DIRKTimeIntegrationSchemes.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: DIRKTimeIntegrationSchemes.h
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2018 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// License for the specific language governing rights and limitations under
14// Permission is hereby granted, free of charge, to any person obtaining a
15// copy of this software and associated documentation files (the "Software"),
16// to deal in the Software without restriction, including without limitation
17// the rights to use, copy, modify, merge, publish, distribute, sublicense,
18// and/or sell copies of the Software, and to permit persons to whom the
19// Software is furnished to do so, subject to the following conditions:
20//
21// The above copyright notice and this permission notice shall be included
22// in all copies or substantial portions of the Software.
23//
24// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30// DEALINGS IN THE SOFTWARE.
31//
32// Description: Combined header file for all basic DIRK based time integration
33// schemes.
34//
35///////////////////////////////////////////////////////////////////////////////
36
37// Note : If adding a new integrator be sure to register the
38// integrator with the Time Integration Scheme Facatory in
39// SchemeInitializor.cpp.
40
41#ifndef NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_DIRK_TIME_INTEGRATION_SCHEME
42#define NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_DIRK_TIME_INTEGRATION_SCHEME
43
44#define LUE LIB_UTILITIES_EXPORT
45using namespace std;
46
48
49namespace Nektar
50{
51namespace LibUtilities
52{
53
54///////////////////////////////////////////////////////////////////////////////
55// DIRK Order N
57{
58public:
59 DIRKTimeIntegrationScheme(std::string variant, size_t order,
60 std::vector<NekDouble> freeParams)
61 : TimeIntegrationSchemeGLM(variant, order, freeParams)
62 {
63 // Currently up to 4th order is implemented.
64 ASSERTL0(1 <= order && order <= 4,
65 "Diagonally Implicit Runge Kutta integration scheme bad order "
66 "(1-4): " +
67 std::to_string(order));
68
72
73 ASSERTL1((variant == "" || variant == "ES5" || variant == "ES6"),
74 "DIRK Time integration scheme bad variant: " + variant +
75 ". "
76 "Must blank, 'ES5', or 'ES6'");
77
78 if ("" == variant)
79 {
81 order);
82 }
83 else
84 {
86 m_integration_phases[0], variant, order, freeParams);
87 }
88 }
89
91 {
92 }
93
95 std::string variant, size_t order, std::vector<NekDouble> freeParams)
96 {
99 variant, order, freeParams);
100 return p;
101 }
102
103 static std::string className;
104
106 size_t order)
107 {
108 phase->m_schemeType = eDiagonallyImplicit;
109 phase->m_order = order;
110 phase->m_name =
111 std::string("DIRKOrder" + std::to_string(phase->m_order));
112
113 phase->m_numsteps = 1;
114 phase->m_numstages = phase->m_order;
115
116 phase->m_A = Array<OneD, Array<TwoD, NekDouble>>(1);
117 phase->m_B = Array<OneD, Array<TwoD, NekDouble>>(1);
118
119 phase->m_A[0] =
120 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numstages, 0.0);
121 phase->m_B[0] =
122 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numstages, 0.0);
123 phase->m_U =
124 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numsteps, 1.0);
125 phase->m_V =
126 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numsteps, 1.0);
127
128 switch (phase->m_order)
129 {
130 case 1:
131 {
132 // One-stage, 1st order, L-stable Diagonally Implicit
133 // Runge Kutta (backward Euler) method:
134
135 phase->m_A[0][0][0] = 1.0;
136
137 phase->m_B[0][0][0] = 1.0;
138 }
139 break;
140 case 2:
141 {
142 // Two-stage, 2nd order Diagonally Implicit Runge
143 // Kutta method. It is A-stable if and only if x ≥ 1/4
144 // and is L-stable if and only if x equals one of the
145 // roots of the polynomial x^2 - 2x + 1/2, i.e. if
146 // x = 1 ± sqrt(2)/2.
147 const NekDouble lambda = (2.0 - sqrt(2.0)) / 2.0;
148
149 phase->m_A[0][0][0] = lambda;
150 phase->m_A[0][1][0] = 1.0 - lambda;
151 phase->m_A[0][1][1] = lambda;
152
153 phase->m_B[0][0][0] = 1.0 - lambda;
154 phase->m_B[0][0][1] = lambda;
155 }
156 break;
157
158 case 3:
159 {
160 // Three-stage, 3rd order, L-stable Diagonally Implicit
161 // Runge Kutta method:
162 constexpr NekDouble lambda = 0.4358665215;
163
164 phase->m_A[0][0][0] = lambda;
165 phase->m_A[0][1][0] = 0.5 * (1.0 - lambda);
166 phase->m_A[0][2][0] =
167 0.25 * (-6.0 * lambda * lambda + 16.0 * lambda - 1.0);
168
169 phase->m_A[0][1][1] = lambda;
170
171 phase->m_A[0][2][1] =
172 0.25 * (6.0 * lambda * lambda - 20.0 * lambda + 5.0);
173 phase->m_A[0][2][2] = lambda;
174
175 phase->m_B[0][0][0] =
176 0.25 * (-6.0 * lambda * lambda + 16.0 * lambda - 1.0);
177 phase->m_B[0][0][1] =
178 0.25 * (6.0 * lambda * lambda - 20.0 * lambda + 5.0);
179 phase->m_B[0][0][2] = lambda;
180 }
181 break;
182 }
183
184 phase->m_numMultiStepValues = 1;
185 phase->m_numMultiStepImplicitDerivs = 0;
186 phase->m_numMultiStepExplicitDerivs = 0;
187 phase->m_timeLevelOffset = Array<OneD, size_t>(phase->m_numsteps);
188 phase->m_timeLevelOffset[0] = 0;
189
190 phase->CheckAndVerify();
191 }
192
194 TimeIntegrationAlgorithmGLMSharedPtr &phase, std::string variant,
195 size_t order, std::vector<NekDouble> freeParams)
196 {
197 phase->m_schemeType = eDiagonallyImplicit;
198 phase->m_order = order;
199 phase->m_name =
200 std::string("DIRKOrder" + std::to_string(phase->m_order) + variant);
201
202 phase->m_numsteps = 1;
203 char const &stage = variant.back();
204 phase->m_numstages = std::atoi(&stage);
205
206 phase->m_A = Array<OneD, Array<TwoD, NekDouble>>(1);
207 phase->m_B = Array<OneD, Array<TwoD, NekDouble>>(1);
208
209 phase->m_A[0] =
210 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numstages, 0.0);
211 phase->m_B[0] =
212 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numstages, 0.0);
213 phase->m_U =
214 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numsteps, 1.0);
215 phase->m_V =
216 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numsteps, 1.0);
217
218 const NekDouble ConstSqrt2 = sqrt(2.0);
219 switch (phase->m_order)
220 {
221 case 3:
222 {
223 // See: Kennedy, Christopher A., and Mark H. Carpenter.
224 // Diagonally implicit Runge-Kutta methods for ordinary
225 // differential equations. A review. No. NF1676L-19716. 2016.
226 ASSERTL0(5 == phase->m_numstages,
227 std::string("DIRKOrder3_ES" +
228 std::to_string(phase->m_numstages) +
229 " not defined"));
230
231 NekDouble lambda =
232 (freeParams.size()) ? freeParams[0] : 9.0 / 40.0;
233
234 phase->m_A[0][0][0] = 0.0;
235 phase->m_A[0][1][0] = lambda;
236 phase->m_A[0][2][0] = 9.0 * (1.0 + ConstSqrt2) / 80.0;
237 phase->m_A[0][3][0] =
238 (22.0 + 15.0 * ConstSqrt2) / (80.0 * (1.0 + ConstSqrt2));
239 phase->m_A[0][4][0] = (2398.0 + 1205.0 * ConstSqrt2) /
240 (2835.0 * (4.0 + 3.0 * ConstSqrt2));
241
242 phase->m_A[0][1][1] = phase->m_A[0][1][0];
243 phase->m_A[0][2][1] = phase->m_A[0][2][0];
244 phase->m_A[0][3][1] = phase->m_A[0][3][0];
245 phase->m_A[0][4][1] = phase->m_A[0][4][0];
246
247 phase->m_A[0][2][2] = lambda;
248 phase->m_A[0][3][2] = -7.0 / (40.0 * (1.0 + ConstSqrt2));
249 phase->m_A[0][4][2] = -2374.0 * (1.0 + 2.0 * ConstSqrt2) /
250 (2835.0 * (5.0 + 3.0 * ConstSqrt2));
251
252 phase->m_A[0][3][3] = lambda;
253 phase->m_A[0][4][3] = 5827.0 / 7560.0;
254
255 phase->m_A[0][4][4] = lambda;
256
257 phase->m_B[0][0][0] = phase->m_A[0][4][0];
258 phase->m_B[0][0][1] = phase->m_A[0][4][1];
259 phase->m_B[0][0][2] = phase->m_A[0][4][2];
260 phase->m_B[0][0][3] = phase->m_A[0][4][3];
261 phase->m_B[0][0][4] = phase->m_A[0][4][4];
262 }
263 break;
264
265 case 4:
266 {
267 // See: Kennedy, Christopher A., and Mark H. Carpenter.
268 // Diagonally implicit Runge-Kutta methods for ordinary
269 // differential equations. A review. No. NF1676L-19716. 2016.
270 ASSERTL0(6 == phase->m_numstages,
271 std::string("DIRKOrder4_ES" +
272 std::to_string(phase->m_numstages) +
273 " not defined"));
274
275 NekDouble lambda =
276 (freeParams.size()) ? freeParams[0] : 1.0 / 4.0;
277
278 phase->m_A[0][0][0] = 0.0;
279 phase->m_A[0][1][0] = lambda;
280 phase->m_A[0][2][0] = (1.0 - ConstSqrt2) / 8.0;
281 phase->m_A[0][3][0] = (5.0 - 7.0 * ConstSqrt2) / 64.0;
282 phase->m_A[0][4][0] =
283 (-13796.0 - 54539.0 * ConstSqrt2) / 125000.0;
284 phase->m_A[0][5][0] = (1181.0 - 987.0 * ConstSqrt2) / 13782.0;
285
286 phase->m_A[0][1][1] = phase->m_A[0][1][0];
287 phase->m_A[0][2][1] = phase->m_A[0][2][0];
288 phase->m_A[0][3][1] = phase->m_A[0][3][0];
289 phase->m_A[0][4][1] = phase->m_A[0][4][0];
290 phase->m_A[0][5][1] = phase->m_A[0][5][0];
291
292 phase->m_A[0][2][2] = lambda;
293 phase->m_A[0][3][2] = 7.0 * (1.0 + ConstSqrt2) / 32.0;
294 phase->m_A[0][4][2] =
295 (506605.0 + 132109.0 * ConstSqrt2) / 437500.0;
296 phase->m_A[0][5][2] =
297 47.0 * (-267.0 + 1783.0 * ConstSqrt2) / 273343.0;
298
299 phase->m_A[0][3][3] = lambda;
300 phase->m_A[0][4][3] =
301 166.0 * (-97.0 + 376.0 * ConstSqrt2) / 109375.0;
302 phase->m_A[0][5][3] =
303 -16.0 * (-22922.0 + 3525.0 * ConstSqrt2) / 571953.0;
304
305 phase->m_A[0][4][4] = lambda;
306 phase->m_A[0][5][4] =
307 -15625.0 * (97.0 + 376.0 * ConstSqrt2) / 90749876.0;
308
309 phase->m_A[0][5][5] = lambda;
310
311 phase->m_B[0][0][0] = phase->m_A[0][5][0];
312 phase->m_B[0][0][1] = phase->m_A[0][5][1];
313 phase->m_B[0][0][2] = phase->m_A[0][5][2];
314 phase->m_B[0][0][3] = phase->m_A[0][5][3];
315 phase->m_B[0][0][4] = phase->m_A[0][5][4];
316 phase->m_B[0][0][5] = phase->m_A[0][5][5];
317 }
318 break;
319 default:
320 {
321 ASSERTL0(false, std::string("ESDIRK of order" +
322 std::to_string(phase->m_order) +
323 " not defined"));
324 break;
325 }
326 }
327
328 phase->m_numMultiStepValues = 1;
329 phase->m_numMultiStepImplicitDerivs = 0;
330 phase->m_numMultiStepExplicitDerivs = 0;
331 phase->m_timeLevelOffset = Array<OneD, size_t>(phase->m_numsteps);
332 phase->m_timeLevelOffset[0] = 0;
333
334 phase->CheckAndVerify();
335 }
336
337protected:
338 LUE virtual std::string v_GetName() const override
339 {
340 return std::string("DIRK");
341 }
342
343 LUE virtual NekDouble v_GetTimeStability() const override
344 {
345 return 1.0;
346 }
347
348}; // end class DIRKTimeIntegrationScheme
349
350////////////////////////////////////////////////////////////////////////////////
351// Backwards compatibility
353{
354public:
355 DIRKOrder1TimeIntegrationScheme(std::string variant, size_t order,
356 std::vector<NekDouble> freeParams)
357 : DIRKTimeIntegrationScheme("", 1, freeParams)
358 {
359 boost::ignore_unused(variant);
360 boost::ignore_unused(order);
361 }
362
364 std::string variant, size_t order, std::vector<NekDouble> freeParams)
365 {
366 boost::ignore_unused(variant);
367 boost::ignore_unused(order);
368
371 "", 1, freeParams);
372
373 return p;
374 }
375
376 static std::string className;
377
378protected:
380
381}; // end class DIRKOrder1TimeIntegrationScheme
382
384{
385public:
386 DIRKOrder2TimeIntegrationScheme(std::string variant, size_t order,
387 std::vector<NekDouble> freeParams)
388 : DIRKTimeIntegrationScheme("", 2, freeParams)
389 {
390 boost::ignore_unused(variant);
391 boost::ignore_unused(order);
392 }
393
395 std::string variant, size_t order, std::vector<NekDouble> freeParams)
396 {
397 boost::ignore_unused(variant);
398 boost::ignore_unused(order);
399
402 "", 2, freeParams);
403
404 return p;
405 }
406
407 static std::string className;
408
409protected:
411
412}; // end class DIRKOrder2TimeIntegrationScheme
413
415{
416public:
417 DIRKOrder3TimeIntegrationScheme(std::string variant, size_t order,
418 std::vector<NekDouble> freeParams)
419 : DIRKTimeIntegrationScheme("", 3, freeParams)
420 {
421 boost::ignore_unused(variant);
422 boost::ignore_unused(order);
423 }
424
426 std::string variant, size_t order, std::vector<NekDouble> freeParams)
427 {
428 boost::ignore_unused(variant);
429 boost::ignore_unused(order);
430
433 "", 3, freeParams);
434
435 return p;
436 }
437
438 static std::string className;
439
440protected:
442
443}; // end class DIRKOrder3TimeIntegrationScheme
444
446{
447public:
448 DIRKOrder3_ES5TimeIntegrationScheme(std::string variant, size_t order,
449 std::vector<NekDouble> freeParams)
450 : DIRKTimeIntegrationScheme("ES5", 3, freeParams)
451 {
452 boost::ignore_unused(variant);
453 boost::ignore_unused(order);
454 }
455
457 std::string variant, size_t order, std::vector<NekDouble> freeParams)
458 {
459 boost::ignore_unused(variant);
460 boost::ignore_unused(order);
461
464 "ES5", 3, freeParams);
465
466 return p;
467 }
468
469 static std::string className;
470
471protected:
473
474}; // end class DIRKOrder3_ES5TimeIntegrationScheme
475
477{
478public:
479 DIRKOrder4_ES6TimeIntegrationScheme(std::string variant, size_t order,
480 std::vector<NekDouble> freeParams)
481 : DIRKTimeIntegrationScheme("ES6", 4, freeParams)
482 {
483 boost::ignore_unused(variant);
484 boost::ignore_unused(order);
485 }
486
488 std::string variant, size_t order, std::vector<NekDouble> freeParams)
489 {
490 boost::ignore_unused(variant);
491 boost::ignore_unused(order);
492
495 "ES6", 4, freeParams);
496
497 return p;
498 }
499
500 static std::string className;
501
502protected:
504}; // end class DIRKOrder4_ES6TimeIntegrationScheme
505
506} // end namespace LibUtilities
507} // end namespace Nektar
508
509#endif
#define LUE
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:215
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:249
static TimeIntegrationSchemeSharedPtr create(std::string variant, size_t order, std::vector< NekDouble > freeParams)
DIRKOrder1TimeIntegrationScheme(std::string variant, size_t order, std::vector< NekDouble > freeParams)
static TimeIntegrationSchemeSharedPtr create(std::string variant, size_t order, std::vector< NekDouble > freeParams)
DIRKOrder2TimeIntegrationScheme(std::string variant, size_t order, std::vector< NekDouble > freeParams)
DIRKOrder3_ES5TimeIntegrationScheme(std::string variant, size_t order, std::vector< NekDouble > freeParams)
static TimeIntegrationSchemeSharedPtr create(std::string variant, size_t order, std::vector< NekDouble > freeParams)
DIRKOrder3TimeIntegrationScheme(std::string variant, size_t order, std::vector< NekDouble > freeParams)
static TimeIntegrationSchemeSharedPtr create(std::string variant, size_t order, std::vector< NekDouble > freeParams)
static TimeIntegrationSchemeSharedPtr create(std::string variant, size_t order, std::vector< NekDouble > freeParams)
DIRKOrder4_ES6TimeIntegrationScheme(std::string variant, size_t order, std::vector< NekDouble > freeParams)
static LUE void SetupSchemeDataESDIRK(TimeIntegrationAlgorithmGLMSharedPtr &phase, std::string variant, size_t order, std::vector< NekDouble > freeParams)
static LUE void SetupSchemeData(TimeIntegrationAlgorithmGLMSharedPtr &phase, size_t order)
virtual LUE NekDouble v_GetTimeStability() const override
DIRKTimeIntegrationScheme(std::string variant, size_t order, std::vector< NekDouble > freeParams)
static TimeIntegrationSchemeSharedPtr create(std::string variant, size_t order, std::vector< NekDouble > freeParams)
virtual LUE std::string v_GetName() const override
Base class for GLM time integration schemes.
TimeIntegrationAlgorithmGLMVector m_integration_phases
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< TimeIntegrationAlgorithmGLM > TimeIntegrationAlgorithmGLMSharedPtr
@ eDiagonallyImplicit
Diagonally implicit scheme (e.g. the DIRK schemes)
std::vector< TimeIntegrationAlgorithmGLMSharedPtr > TimeIntegrationAlgorithmGLMVector
std::shared_ptr< TimeIntegrationScheme > TimeIntegrationSchemeSharedPtr
The above copyright notice and this permission notice shall be included.
Definition: CoupledSolver.h:2
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294