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