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
60 std::vector<NekDouble> freeParams)
65 "Diagonally Implicit Runge Kutta integration scheme bad order "
67 std::to_string(order));
73 ASSERTL1((variant ==
"" || variant ==
"ES5" || variant ==
"ES6"),
74 "DIRK Time integration scheme bad variant: " + variant +
76 "Must blank, 'ES5', or 'ES6'");
95 std::string variant,
size_t order, std::vector<NekDouble> freeParams)
99 variant, order, freeParams);
109 phase->m_order = order;
111 std::string(
"DIRKOrder" + std::to_string(phase->m_order));
113 phase->m_numsteps = 1;
114 phase->m_numstages = phase->m_order;
128 switch (phase->m_order)
135 phase->m_A[0][0][0] = 1.0;
137 phase->m_B[0][0][0] = 1.0;
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;
153 phase->m_B[0][0][0] = 1.0 - lambda;
154 phase->m_B[0][0][1] = lambda;
162 constexpr NekDouble lambda = 0.4358665215;
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);
169 phase->m_A[0][1][1] = lambda;
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;
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;
184 phase->m_numMultiStepValues = 1;
185 phase->m_numMultiStepImplicitDerivs = 0;
186 phase->m_numMultiStepExplicitDerivs = 0;
188 phase->m_timeLevelOffset[0] = 0;
190 phase->CheckAndVerify();
195 size_t order, std::vector<NekDouble> freeParams)
198 phase->m_order = order;
200 std::string(
"DIRKOrder" + std::to_string(phase->m_order) + variant);
202 phase->m_numsteps = 1;
203 char const &stage = variant.back();
204 phase->m_numstages = std::atoi(&stage);
219 switch (phase->m_order)
227 std::string(
"DIRKOrder3_ES" +
228 std::to_string(phase->m_numstages) +
232 (freeParams.size()) ? freeParams[0] : 9.0 / 40.0;
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));
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];
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));
252 phase->m_A[0][3][3] = lambda;
253 phase->m_A[0][4][3] = 5827.0 / 7560.0;
255 phase->m_A[0][4][4] = lambda;
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];
271 std::string(
"DIRKOrder4_ES" +
272 std::to_string(phase->m_numstages) +
276 (freeParams.size()) ? freeParams[0] : 1.0 / 4.0;
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;
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];
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;
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;
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;
309 phase->m_A[0][5][5] = lambda;
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];
321 ASSERTL0(
false, std::string(
"ESDIRK of order" +
322 std::to_string(phase->m_order) +
328 phase->m_numMultiStepValues = 1;
329 phase->m_numMultiStepImplicitDerivs = 0;
330 phase->m_numMultiStepExplicitDerivs = 0;
332 phase->m_timeLevelOffset[0] = 0;
334 phase->CheckAndVerify();
340 return std::string(
"DIRK");
356 std::vector<NekDouble> freeParams)
359 boost::ignore_unused(variant);
360 boost::ignore_unused(order);
364 std::string variant,
size_t order, std::vector<NekDouble> freeParams)
366 boost::ignore_unused(variant);
367 boost::ignore_unused(order);
387 std::vector<NekDouble> freeParams)
390 boost::ignore_unused(variant);
391 boost::ignore_unused(order);
395 std::string variant,
size_t order, std::vector<NekDouble> freeParams)
397 boost::ignore_unused(variant);
398 boost::ignore_unused(order);
418 std::vector<NekDouble> freeParams)
421 boost::ignore_unused(variant);
422 boost::ignore_unused(order);
426 std::string variant,
size_t order, std::vector<NekDouble> freeParams)
428 boost::ignore_unused(variant);
429 boost::ignore_unused(order);
449 std::vector<NekDouble> freeParams)
452 boost::ignore_unused(variant);
453 boost::ignore_unused(order);
457 std::string variant,
size_t order, std::vector<NekDouble> freeParams)
459 boost::ignore_unused(variant);
460 boost::ignore_unused(order);
464 "ES5", 3, freeParams);
480 std::vector<NekDouble> freeParams)
483 boost::ignore_unused(variant);
484 boost::ignore_unused(order);
488 std::string variant,
size_t order, std::vector<NekDouble> freeParams)
490 boost::ignore_unused(variant);
491 boost::ignore_unused(order);
495 "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
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 std::string TimeIntegrationMethodLookupId
static TimeIntegrationSchemeSharedPtr create(std::string variant, size_t order, std::vector< NekDouble > freeParams)
static std::string className
DIRKOrder2TimeIntegrationScheme(std::string variant, size_t order, std::vector< NekDouble > freeParams)
static std::string TimeIntegrationMethodLookupId
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)
static std::string TimeIntegrationMethodLookupId
static std::string className
static std::string TimeIntegrationMethodLookupId
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 std::string className
static TimeIntegrationSchemeSharedPtr create(std::string variant, size_t order, std::vector< NekDouble > freeParams)
static std::string className
static std::string TimeIntegrationMethodLookupId
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
virtual ~DIRKTimeIntegrationScheme()
static std::string className
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.
scalarT< T > sqrt(scalarT< T > in)