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  ASSERTL1(2 <= order && order <= 4,
66  "Diagonally Implicit Runge Kutta integration scheme bad order "
67  "(2-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 
107  LUE virtual std::string GetName() const
108  {
109  return std::string("DIRK");
110  }
111 
113  {
114  return 1.0;
115  }
116 
118  unsigned int order)
119  {
120  phase->m_schemeType = eDiagonallyImplicit;
121  phase->m_order = order;
122  phase->m_name =
123  std::string("DIRKOrder" + std::to_string(phase->m_order));
124 
125  phase->m_numsteps = 1;
126  phase->m_numstages = phase->m_order;
127 
128  phase->m_A = Array<OneD, Array<TwoD, NekDouble>>(1);
129  phase->m_B = Array<OneD, Array<TwoD, NekDouble>>(1);
130 
131  phase->m_A[0] =
132  Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numstages, 0.0);
133  phase->m_B[0] =
134  Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numstages, 0.0);
135  phase->m_U =
136  Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numsteps, 1.0);
137  phase->m_V =
138  Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numsteps, 1.0);
139 
140  switch (phase->m_order)
141  {
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_numMultiStepDerivs = 0;
188  phase->m_timeLevelOffset = Array<OneD, unsigned int>(phase->m_numsteps);
189  phase->m_timeLevelOffset[0] = 0;
190 
191  phase->CheckAndVerify();
192  }
193 
195  TimeIntegrationAlgorithmGLMSharedPtr &phase, std::string variant,
196  unsigned int order, std::vector<NekDouble> freeParams)
197  {
198  phase->m_schemeType = eDiagonallyImplicit;
199  phase->m_order = order;
200  phase->m_name =
201  std::string("DIRKOrder" + std::to_string(phase->m_order) + variant);
202 
203  phase->m_numsteps = 1;
204  char const &stage = variant.back();
205  phase->m_numstages = std::atoi(&stage);
206 
207  phase->m_A = Array<OneD, Array<TwoD, NekDouble>>(1);
208  phase->m_B = Array<OneD, Array<TwoD, NekDouble>>(1);
209 
210  phase->m_A[0] =
211  Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numstages, 0.0);
212  phase->m_B[0] =
213  Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numstages, 0.0);
214  phase->m_U =
215  Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numsteps, 1.0);
216  phase->m_V =
217  Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numsteps, 1.0);
218 
219  const NekDouble ConstSqrt2 = sqrt(2.0);
220  switch (phase->m_order)
221  {
222  case 3:
223  {
224  ASSERTL0(5 == phase->m_numstages,
225  std::string("DIRKOrder3_ES" +
226  std::to_string(phase->m_numstages) +
227  " not defined"));
228  NekDouble lambda;
229  if (freeParams.size())
230  {
231  lambda = freeParams[0];
232  }
233  else
234  {
235  lambda = 9.0 / 40.0;
236  }
237 
238  phase->m_A[0][0][0] = 0.0;
239  phase->m_A[0][1][0] = lambda;
240  phase->m_A[0][2][0] = 9.0 * (1.0 + ConstSqrt2) / 80.0;
241  phase->m_A[0][3][0] =
242  (22.0 + 15.0 * ConstSqrt2) / (80.0 * (1 + ConstSqrt2));
243  phase->m_A[0][4][0] = (2398.0 + 1205.0 * ConstSqrt2) /
244  (2835.0 * (4 + 3.0 * ConstSqrt2));
245 
246  phase->m_A[0][1][1] = phase->m_A[0][1][0];
247  phase->m_A[0][2][1] = phase->m_A[0][2][0];
248  phase->m_A[0][3][1] = phase->m_A[0][3][0];
249  phase->m_A[0][4][1] = phase->m_A[0][4][0];
250 
251  phase->m_A[0][2][2] = lambda;
252  phase->m_A[0][3][2] = -7.0 / (40.0 * (1.0 + ConstSqrt2));
253  phase->m_A[0][4][2] = -2374 * (2.0 + ConstSqrt2) /
254  (2835.0 * (4.0 + 3.0 * ConstSqrt2));
255 
256  phase->m_A[0][3][3] = lambda;
257  phase->m_A[0][4][3] = 5827.0 / 7560.0;
258 
259  phase->m_A[0][4][4] = lambda;
260 
261  phase->m_B[0][0][0] = phase->m_A[0][4][0];
262  phase->m_B[0][0][1] = phase->m_A[0][4][1];
263  phase->m_B[0][0][2] = phase->m_A[0][4][2];
264  phase->m_B[0][0][3] = phase->m_A[0][4][3];
265  phase->m_B[0][0][4] = phase->m_A[0][4][4];
266  }
267  break;
268 
269  case 4:
270  {
271  ASSERTL0(6 == phase->m_numstages,
272  std::string("DIRKOrder4_ES" +
273  std::to_string(phase->m_numstages) +
274  " not defined"));
275  NekDouble lambda;
276  if (freeParams.size())
277  {
278  lambda = freeParams[0];
279  }
280  else
281  {
282  lambda = 1.0 / 4.0;
283  }
284 
285  const NekDouble ConstSqrt2 = sqrt(2.0);
286 
287  phase->m_A[0][0][0] = 0.0;
288  phase->m_A[0][1][0] = lambda;
289  phase->m_A[0][2][0] = (1.0 - ConstSqrt2) / 8.0;
290  phase->m_A[0][3][0] = (5.0 - 7.0 * ConstSqrt2) / 64.0;
291  phase->m_A[0][4][0] =
292  (-13796.0 - 54539 * ConstSqrt2) / 125000.0;
293  phase->m_A[0][5][0] = (1181.0 - 987.0 * ConstSqrt2) / 13782.0;
294 
295  phase->m_A[0][1][1] = phase->m_A[0][1][0];
296  phase->m_A[0][2][1] = phase->m_A[0][2][0];
297  phase->m_A[0][3][1] = phase->m_A[0][3][0];
298  phase->m_A[0][4][1] = phase->m_A[0][4][0];
299  phase->m_A[0][5][1] = phase->m_A[0][5][0];
300 
301  phase->m_A[0][2][2] = lambda;
302  phase->m_A[0][3][2] = 7.0 * (1.0 + ConstSqrt2) / 32.0;
303  phase->m_A[0][4][2] =
304  (506605.0 + 132109.0 * ConstSqrt2) / 437500.0;
305  phase->m_A[0][5][2] =
306  47.0 * (-267.0 + 1783.0 * ConstSqrt2) / 273343.0;
307 
308  phase->m_A[0][3][3] = lambda;
309  phase->m_A[0][4][3] =
310  166.0 * (-97.0 + 376.0 * ConstSqrt2) / 109375.0;
311  phase->m_A[0][5][3] =
312  -16.0 * (-22922.0 + 3525.0 * ConstSqrt2) / 571953.0;
313 
314  phase->m_A[0][4][4] = lambda;
315  phase->m_A[0][5][4] =
316  -15625.0 * (97.0 + 376.0 * ConstSqrt2) / 90749876.0;
317 
318  phase->m_A[0][5][5] = lambda;
319  }
320  break;
321  default:
322  {
323  ASSERTL0(false, std::string("ESDIRK of order" +
324  std::to_string(phase->m_order) +
325  " not defined"));
326  break;
327  }
328  }
329 
330  phase->m_numMultiStepValues = 1;
331  phase->m_numMultiStepDerivs = 0;
332  phase->m_timeLevelOffset = Array<OneD, unsigned int>(phase->m_numsteps);
333  phase->m_timeLevelOffset[0] = 0;
334 
335  phase->CheckAndVerify();
336  }
337 
338 }; // end class DIRKTimeIntegrationScheme
339 
340 ////////////////////////////////////////////////////////////////////////////////
341 // Backwards compatibility
343 {
344 public:
345  DIRKOrder2TimeIntegrationScheme(std::string variant, unsigned int order,
346  std::vector<NekDouble> freeParams)
347  : DIRKTimeIntegrationScheme("", 2, freeParams)
348  {
349  boost::ignore_unused(variant);
350  boost::ignore_unused(order);
351  }
352 
354  std::string variant, unsigned int order,
355  std::vector<NekDouble> freeParams)
356  {
357  boost::ignore_unused(variant);
358  boost::ignore_unused(order);
359 
362  "", 2, freeParams);
363 
364  return p;
365  }
366 
367  static std::string className;
368 
369 }; // end class DIRKOrder2TimeIntegrationScheme
370 
372 {
373 public:
374  DIRKOrder3TimeIntegrationScheme(std::string variant, unsigned int order,
375  std::vector<NekDouble> freeParams)
376  : DIRKTimeIntegrationScheme("", 3, freeParams)
377  {
378  boost::ignore_unused(variant);
379  boost::ignore_unused(order);
380  }
381 
383  std::string variant, unsigned int order,
384  std::vector<NekDouble> freeParams)
385  {
386  boost::ignore_unused(variant);
387  boost::ignore_unused(order);
388 
391  "", 3, freeParams);
392 
393  return p;
394  }
395 
396  static std::string className;
397 
398 }; // end class DIRKOrder3TimeIntegrationScheme
399 
401 {
402 public:
403  DIRKOrder3_ES5TimeIntegrationScheme(std::string variant, unsigned int order,
404  std::vector<NekDouble> freeParams)
405  : DIRKTimeIntegrationScheme("ES5", 3, freeParams)
406  {
407  boost::ignore_unused(variant);
408  boost::ignore_unused(order);
409  }
410 
412  std::string variant, unsigned int order,
413  std::vector<NekDouble> freeParams)
414  {
415  boost::ignore_unused(variant);
416  boost::ignore_unused(order);
417 
420  "ES5", 3, freeParams);
421 
422  return p;
423  }
424 
425  static std::string className;
426 
427 }; // end class DIRKOrder3_ES5TimeIntegrationScheme
428 
430 {
431 public:
432  DIRKOrder4_ES6TimeIntegrationScheme(std::string variant, unsigned int order,
433  std::vector<NekDouble> freeParams)
434  : DIRKTimeIntegrationScheme("ES6", 4, freeParams)
435  {
436  boost::ignore_unused(variant);
437  boost::ignore_unused(order);
438  }
439 
441  std::string variant, unsigned int order,
442  std::vector<NekDouble> freeParams)
443  {
444  boost::ignore_unused(variant);
445  boost::ignore_unused(order);
446 
449  "ES6", 4, freeParams);
450 
451  return p;
452  }
453 
454  static std::string className;
455 };
456 
457 } // end namespace LibUtilities
458 } // end namespace Nektar
459 
460 #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
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)
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)
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:1
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:291