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
45
47
49{
50
51///////////////////////////////////////////////////////////////////////////////
52// DIRK Order N
54{
55public:
56 DIRKTimeIntegrationScheme(std::string variant, size_t order,
57 std::vector<NekDouble> freeParams)
58 : TimeIntegrationSchemeGLM(variant, order, freeParams)
59 {
60 // Currently up to 4th order is implemented.
61 ASSERTL0(1 <= order && order <= 4,
62 "Diagonally Implicit Runge Kutta integration scheme bad order "
63 "(1-4): " +
64 std::to_string(order));
65
69
70 ASSERTL1((variant == "" || variant == "ES5" || variant == "ES6"),
71 "DIRK Time integration scheme bad variant: " + variant +
72 ". "
73 "Must blank, 'ES5', or 'ES6'");
74
75 if ("" == variant)
76 {
78 order);
79 }
80 else
81 {
83 m_integration_phases[0], variant, order, freeParams);
84 }
85 }
86
88 {
89 }
90
92 std::string variant, size_t order, std::vector<NekDouble> freeParams)
93 {
96 variant, order, freeParams);
97 return p;
98 }
99
100 static std::string className;
101
103 size_t order)
104 {
105 phase->m_schemeType = eDiagonallyImplicit;
106 phase->m_order = order;
107 phase->m_name =
108 std::string("DIRKOrder" + std::to_string(phase->m_order));
109
110 phase->m_numsteps = 1;
111 phase->m_numstages = phase->m_order;
112
113 phase->m_A = Array<OneD, Array<TwoD, NekDouble>>(1);
114 phase->m_B = Array<OneD, Array<TwoD, NekDouble>>(1);
115
116 phase->m_A[0] =
117 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numstages, 0.0);
118 phase->m_B[0] =
119 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numstages, 0.0);
120 phase->m_U =
121 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numsteps, 1.0);
122 phase->m_V =
123 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numsteps, 1.0);
124
125 switch (phase->m_order)
126 {
127 case 1:
128 {
129 // One-stage, 1st order, L-stable Diagonally Implicit
130 // Runge Kutta (backward Euler) method:
131
132 phase->m_A[0][0][0] = 1.0;
133
134 phase->m_B[0][0][0] = 1.0;
135 }
136 break;
137 case 2:
138 {
139 // Two-stage, 2nd order Diagonally Implicit Runge
140 // Kutta method. It is A-stable if and only if x ≥ 1/4
141 // and is L-stable if and only if x equals one of the
142 // roots of the polynomial x^2 - 2x + 1/2, i.e. if
143 // x = 1 ± sqrt(2)/2.
144 const NekDouble lambda = (2.0 - sqrt(2.0)) / 2.0;
145
146 phase->m_A[0][0][0] = lambda;
147 phase->m_A[0][1][0] = 1.0 - lambda;
148 phase->m_A[0][1][1] = lambda;
149
150 phase->m_B[0][0][0] = 1.0 - lambda;
151 phase->m_B[0][0][1] = lambda;
152 }
153 break;
154
155 case 3:
156 {
157 // Three-stage, 3rd order, L-stable Diagonally Implicit
158 // Runge Kutta method:
159 constexpr NekDouble lambda = 0.4358665215;
160
161 phase->m_A[0][0][0] = lambda;
162 phase->m_A[0][1][0] = 0.5 * (1.0 - lambda);
163 phase->m_A[0][2][0] =
164 0.25 * (-6.0 * lambda * lambda + 16.0 * lambda - 1.0);
165
166 phase->m_A[0][1][1] = lambda;
167
168 phase->m_A[0][2][1] =
169 0.25 * (6.0 * lambda * lambda - 20.0 * lambda + 5.0);
170 phase->m_A[0][2][2] = lambda;
171
172 phase->m_B[0][0][0] =
173 0.25 * (-6.0 * lambda * lambda + 16.0 * lambda - 1.0);
174 phase->m_B[0][0][1] =
175 0.25 * (6.0 * lambda * lambda - 20.0 * lambda + 5.0);
176 phase->m_B[0][0][2] = lambda;
177 }
178 break;
179 }
180
181 phase->m_numMultiStepValues = 1;
182 phase->m_numMultiStepImplicitDerivs = 0;
183 phase->m_numMultiStepExplicitDerivs = 0;
184 phase->m_timeLevelOffset = Array<OneD, size_t>(phase->m_numsteps);
185 phase->m_timeLevelOffset[0] = 0;
186
187 phase->CheckAndVerify();
188 }
189
191 TimeIntegrationAlgorithmGLMSharedPtr &phase, std::string variant,
192 size_t order, std::vector<NekDouble> freeParams)
193 {
194 phase->m_schemeType = eDiagonallyImplicit;
195 phase->m_order = order;
196 phase->m_name =
197 std::string("DIRKOrder" + std::to_string(phase->m_order) + variant);
198
199 phase->m_numsteps = 1;
200 char const &stage = variant.back();
201 phase->m_numstages = std::atoi(&stage);
202
203 phase->m_A = Array<OneD, Array<TwoD, NekDouble>>(1);
204 phase->m_B = Array<OneD, Array<TwoD, NekDouble>>(1);
205
206 phase->m_A[0] =
207 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numstages, 0.0);
208 phase->m_B[0] =
209 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numstages, 0.0);
210 phase->m_U =
211 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numsteps, 1.0);
212 phase->m_V =
213 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numsteps, 1.0);
214
215 const NekDouble ConstSqrt2 = sqrt(2.0);
216 switch (phase->m_order)
217 {
218 case 3:
219 {
220 // See: Kennedy, Christopher A., and Mark H. Carpenter.
221 // Diagonally implicit Runge-Kutta methods for ordinary
222 // differential equations. A review. No. NF1676L-19716. 2016.
223 ASSERTL0(5 == phase->m_numstages,
224 std::string("DIRKOrder3_ES" +
225 std::to_string(phase->m_numstages) +
226 " not defined"));
227
228 NekDouble lambda =
229 (freeParams.size()) ? freeParams[0] : 9.0 / 40.0;
230
231 phase->m_A[0][0][0] = 0.0;
232 phase->m_A[0][1][0] = lambda;
233 phase->m_A[0][2][0] = 9.0 * (1.0 + ConstSqrt2) / 80.0;
234 phase->m_A[0][3][0] =
235 (22.0 + 15.0 * ConstSqrt2) / (80.0 * (1.0 + ConstSqrt2));
236 phase->m_A[0][4][0] = (2398.0 + 1205.0 * ConstSqrt2) /
237 (2835.0 * (4.0 + 3.0 * ConstSqrt2));
238
239 phase->m_A[0][1][1] = phase->m_A[0][1][0];
240 phase->m_A[0][2][1] = phase->m_A[0][2][0];
241 phase->m_A[0][3][1] = phase->m_A[0][3][0];
242 phase->m_A[0][4][1] = phase->m_A[0][4][0];
243
244 phase->m_A[0][2][2] = lambda;
245 phase->m_A[0][3][2] = -7.0 / (40.0 * (1.0 + ConstSqrt2));
246 phase->m_A[0][4][2] = -2374.0 * (1.0 + 2.0 * ConstSqrt2) /
247 (2835.0 * (5.0 + 3.0 * ConstSqrt2));
248
249 phase->m_A[0][3][3] = lambda;
250 phase->m_A[0][4][3] = 5827.0 / 7560.0;
251
252 phase->m_A[0][4][4] = lambda;
253
254 phase->m_B[0][0][0] = phase->m_A[0][4][0];
255 phase->m_B[0][0][1] = phase->m_A[0][4][1];
256 phase->m_B[0][0][2] = phase->m_A[0][4][2];
257 phase->m_B[0][0][3] = phase->m_A[0][4][3];
258 phase->m_B[0][0][4] = phase->m_A[0][4][4];
259 }
260 break;
261
262 case 4:
263 {
264 // See: Kennedy, Christopher A., and Mark H. Carpenter.
265 // Diagonally implicit Runge-Kutta methods for ordinary
266 // differential equations. A review. No. NF1676L-19716. 2016.
267 ASSERTL0(6 == phase->m_numstages,
268 std::string("DIRKOrder4_ES" +
269 std::to_string(phase->m_numstages) +
270 " not defined"));
271
272 NekDouble lambda =
273 (freeParams.size()) ? freeParams[0] : 1.0 / 4.0;
274
275 phase->m_A[0][0][0] = 0.0;
276 phase->m_A[0][1][0] = lambda;
277 phase->m_A[0][2][0] = (1.0 - ConstSqrt2) / 8.0;
278 phase->m_A[0][3][0] = (5.0 - 7.0 * ConstSqrt2) / 64.0;
279 phase->m_A[0][4][0] =
280 (-13796.0 - 54539.0 * ConstSqrt2) / 125000.0;
281 phase->m_A[0][5][0] = (1181.0 - 987.0 * ConstSqrt2) / 13782.0;
282
283 phase->m_A[0][1][1] = phase->m_A[0][1][0];
284 phase->m_A[0][2][1] = phase->m_A[0][2][0];
285 phase->m_A[0][3][1] = phase->m_A[0][3][0];
286 phase->m_A[0][4][1] = phase->m_A[0][4][0];
287 phase->m_A[0][5][1] = phase->m_A[0][5][0];
288
289 phase->m_A[0][2][2] = lambda;
290 phase->m_A[0][3][2] = 7.0 * (1.0 + ConstSqrt2) / 32.0;
291 phase->m_A[0][4][2] =
292 (506605.0 + 132109.0 * ConstSqrt2) / 437500.0;
293 phase->m_A[0][5][2] =
294 47.0 * (-267.0 + 1783.0 * ConstSqrt2) / 273343.0;
295
296 phase->m_A[0][3][3] = lambda;
297 phase->m_A[0][4][3] =
298 166.0 * (-97.0 + 376.0 * ConstSqrt2) / 109375.0;
299 phase->m_A[0][5][3] =
300 -16.0 * (-22922.0 + 3525.0 * ConstSqrt2) / 571953.0;
301
302 phase->m_A[0][4][4] = lambda;
303 phase->m_A[0][5][4] =
304 -15625.0 * (97.0 + 376.0 * ConstSqrt2) / 90749876.0;
305
306 phase->m_A[0][5][5] = lambda;
307
308 phase->m_B[0][0][0] = phase->m_A[0][5][0];
309 phase->m_B[0][0][1] = phase->m_A[0][5][1];
310 phase->m_B[0][0][2] = phase->m_A[0][5][2];
311 phase->m_B[0][0][3] = phase->m_A[0][5][3];
312 phase->m_B[0][0][4] = phase->m_A[0][5][4];
313 phase->m_B[0][0][5] = phase->m_A[0][5][5];
314 }
315 break;
316 default:
317 {
318 ASSERTL0(false, std::string("ESDIRK of order" +
319 std::to_string(phase->m_order) +
320 " not defined"));
321 break;
322 }
323 }
324
325 phase->m_numMultiStepValues = 1;
326 phase->m_numMultiStepImplicitDerivs = 0;
327 phase->m_numMultiStepExplicitDerivs = 0;
328 phase->m_timeLevelOffset = Array<OneD, size_t>(phase->m_numsteps);
329 phase->m_timeLevelOffset[0] = 0;
330
331 phase->CheckAndVerify();
332 }
333
334protected:
335 LUE std::string v_GetName() const override
336 {
337 return std::string("DIRK");
338 }
339
341 {
342 return 1.0;
343 }
344
345}; // end class DIRKTimeIntegrationScheme
346
347////////////////////////////////////////////////////////////////////////////////
348// Backwards compatibility
350{
351public:
352 DIRKOrder1TimeIntegrationScheme(std::string variant, size_t order,
353 std::vector<NekDouble> freeParams)
354 : DIRKTimeIntegrationScheme("", 1, freeParams)
355 {
356 boost::ignore_unused(variant, order);
357 }
358
360 [[maybe_unused]] std::string variant, [[maybe_unused]] size_t order,
361 std::vector<NekDouble> freeParams)
362 {
365 "", 1, freeParams);
366
367 return p;
368 }
369
370 static std::string className;
371
372protected:
374
375}; // end class DIRKOrder1TimeIntegrationScheme
376
378{
379public:
380 DIRKOrder2TimeIntegrationScheme(std::string variant, size_t order,
381 std::vector<NekDouble> freeParams)
382 : DIRKTimeIntegrationScheme("", 2, freeParams)
383 {
384 boost::ignore_unused(variant, order);
385 }
386
388 [[maybe_unused]] std::string variant, [[maybe_unused]] size_t order,
389 std::vector<NekDouble> freeParams)
390 {
393 "", 2, freeParams);
394
395 return p;
396 }
397
398 static std::string className;
399
400protected:
402
403}; // end class DIRKOrder2TimeIntegrationScheme
404
406{
407public:
408 DIRKOrder3TimeIntegrationScheme(std::string variant, size_t order,
409 std::vector<NekDouble> freeParams)
410 : DIRKTimeIntegrationScheme("", 3, freeParams)
411 {
412 boost::ignore_unused(variant, order);
413 }
414
416 [[maybe_unused]] std::string variant, [[maybe_unused]] size_t order,
417 std::vector<NekDouble> freeParams)
418 {
421 "", 3, freeParams);
422
423 return p;
424 }
425
426 static std::string className;
427
428protected:
430
431}; // end class DIRKOrder3TimeIntegrationScheme
432
434{
435public:
436 DIRKOrder3_ES5TimeIntegrationScheme(std::string variant, size_t order,
437 std::vector<NekDouble> freeParams)
438 : DIRKTimeIntegrationScheme("ES5", 3, freeParams)
439 {
440 boost::ignore_unused(variant, order);
441 }
442
444 [[maybe_unused]] std::string variant, [[maybe_unused]] size_t order,
445 std::vector<NekDouble> freeParams)
446 {
449 "ES5", 3, freeParams);
450
451 return p;
452 }
453
454 static std::string className;
455
456protected:
458
459}; // end class DIRKOrder3_ES5TimeIntegrationScheme
460
462{
463public:
464 DIRKOrder4_ES6TimeIntegrationScheme(std::string variant, size_t order,
465 std::vector<NekDouble> freeParams)
466 : DIRKTimeIntegrationScheme("ES6", 4, freeParams)
467 {
468 boost::ignore_unused(variant, order);
469 }
470
472 [[maybe_unused]] std::string variant, [[maybe_unused]] size_t order,
473 std::vector<NekDouble> freeParams)
474 {
477 "ES6", 4, freeParams);
478
479 return p;
480 }
481
482 static std::string className;
483
484protected:
486}; // end class DIRKOrder4_ES6TimeIntegrationScheme
487
488} // namespace Nektar::LibUtilities
489
490#endif
#define LUE
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
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)
DIRKTimeIntegrationScheme(std::string variant, size_t order, std::vector< NekDouble > freeParams)
static TimeIntegrationSchemeSharedPtr create(std::string variant, size_t order, std::vector< NekDouble > freeParams)
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
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294