41 #ifndef NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_DIRK_TIME_INTEGRATION_SCHEME
42 #define NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_DIRK_TIME_INTEGRATION_SCHEME
44 #define LUE LIB_UTILITIES_EXPORT
52 namespace LibUtilities
61 std::vector<NekDouble> freeParams)
66 "Diagonally Implicit Runge Kutta integration scheme bad order "
68 std::to_string(order));
74 ASSERTL1((variant ==
"" || variant ==
"ES5" || variant ==
"ES6"),
75 "DIRK Time integration scheme bad variant: " + variant +
77 "Must blank, 'ES5', or 'ES6'");
81 DIRKTimeIntegrationScheme::SetupSchemeData(m_integration_phases[0],
86 DIRKTimeIntegrationScheme::SetupSchemeDataESDIRK(
87 m_integration_phases[0], variant, order, freeParams);
96 std::string variant,
unsigned int order,
97 std::vector<NekDouble> freeParams)
101 variant, order, freeParams);
111 phase->m_order = order;
113 std::string(
"DIRKOrder" + std::to_string(phase->m_order));
115 phase->m_numsteps = 1;
116 phase->m_numstages = phase->m_order;
130 switch (phase->m_order)
137 phase->m_A[0][0][0] = 1.0;
139 phase->m_B[0][0][0] = 1.0;
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;
155 phase->m_B[0][0][0] = 1.0 - lambda;
156 phase->m_B[0][0][1] = lambda;
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);
171 phase->m_A[0][1][1] = lambda;
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;
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;
186 phase->m_numMultiStepValues = 1;
187 phase->m_numMultiStepImplicitDerivs = 0;
188 phase->m_numMultiStepDerivs = 0;
190 phase->m_timeLevelOffset[0] = 0;
192 phase->CheckAndVerify();
197 unsigned int order, std::vector<NekDouble> freeParams)
200 phase->m_order = order;
202 std::string(
"DIRKOrder" + std::to_string(phase->m_order) + variant);
204 phase->m_numsteps = 1;
205 char const &stage = variant.back();
206 phase->m_numstages = std::atoi(&stage);
221 switch (phase->m_order)
229 std::string(
"DIRKOrder3_ES" +
230 std::to_string(phase->m_numstages) +
233 if (freeParams.size())
235 lambda = freeParams[0];
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));
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];
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));
260 phase->m_A[0][3][3] = lambda;
261 phase->m_A[0][4][3] = 5827.0 / 7560.0;
263 phase->m_A[0][4][4] = lambda;
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];
279 std::string(
"DIRKOrder4_ES" +
280 std::to_string(phase->m_numstages) +
283 if (freeParams.size())
285 lambda = freeParams[0];
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;
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];
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;
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;
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;
325 phase->m_A[0][5][5] = lambda;
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];
337 ASSERTL0(
false, std::string(
"ESDIRK of order" +
338 std::to_string(phase->m_order) +
344 phase->m_numMultiStepValues = 1;
345 phase->m_numMultiStepImplicitDerivs = 0;
346 phase->m_numMultiStepDerivs = 0;
348 phase->m_timeLevelOffset[0] = 0;
350 phase->CheckAndVerify();
356 return std::string(
"DIRK");
372 std::vector<NekDouble> freeParams)
375 boost::ignore_unused(variant);
376 boost::ignore_unused(order);
380 std::string variant,
unsigned int order,
381 std::vector<NekDouble> freeParams)
383 boost::ignore_unused(variant);
384 boost::ignore_unused(order);
401 std::vector<NekDouble> freeParams)
404 boost::ignore_unused(variant);
405 boost::ignore_unused(order);
409 std::string variant,
unsigned int order,
410 std::vector<NekDouble> freeParams)
412 boost::ignore_unused(variant);
413 boost::ignore_unused(order);
430 std::vector<NekDouble> freeParams)
433 boost::ignore_unused(variant);
434 boost::ignore_unused(order);
438 std::string variant,
unsigned int order,
439 std::vector<NekDouble> freeParams)
441 boost::ignore_unused(variant);
442 boost::ignore_unused(order);
459 std::vector<NekDouble> freeParams)
462 boost::ignore_unused(variant);
463 boost::ignore_unused(order);
467 std::string variant,
unsigned int order,
468 std::vector<NekDouble> freeParams)
470 boost::ignore_unused(variant);
471 boost::ignore_unused(order);
475 "ES5", 3, freeParams);
488 std::vector<NekDouble> freeParams)
491 boost::ignore_unused(variant);
492 boost::ignore_unused(order);
496 std::string variant,
unsigned int order,
497 std::vector<NekDouble> freeParams)
499 boost::ignore_unused(variant);
500 boost::ignore_unused(order);
504 "ES6", 4, freeParams);
#define ASSERTL0(condition, msg)
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
static std::string className
DIRKOrder1TimeIntegrationScheme(std::string variant, unsigned int order, std::vector< NekDouble > freeParams)
static TimeIntegrationSchemeSharedPtr create(std::string variant, unsigned int order, std::vector< NekDouble > freeParams)
static std::string className
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 std::string className
static std::string className
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)
static std::string className
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
virtual ~DIRKTimeIntegrationScheme()
static std::string className
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.
scalarT< T > sqrt(scalarT< T > in)