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
57 std::vector<NekDouble> freeParams)
62 "Diagonally Implicit Runge Kutta integration scheme bad order "
64 std::to_string(order));
70 ASSERTL1((variant ==
"" || variant ==
"ES5" || variant ==
"ES6"),
71 "DIRK Time integration scheme bad variant: " + variant +
73 "Must blank, 'ES5', or 'ES6'");
92 std::string variant,
size_t order, std::vector<NekDouble> freeParams)
96 variant, order, freeParams);
106 phase->m_order = order;
108 std::string(
"DIRKOrder" + std::to_string(phase->m_order));
110 phase->m_numsteps = 1;
111 phase->m_numstages = phase->m_order;
125 switch (phase->m_order)
132 phase->m_A[0][0][0] = 1.0;
134 phase->m_B[0][0][0] = 1.0;
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;
150 phase->m_B[0][0][0] = 1.0 - lambda;
151 phase->m_B[0][0][1] = lambda;
159 constexpr NekDouble lambda = 0.4358665215;
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);
166 phase->m_A[0][1][1] = lambda;
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;
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;
181 phase->m_numMultiStepValues = 1;
182 phase->m_numMultiStepImplicitDerivs = 0;
183 phase->m_numMultiStepExplicitDerivs = 0;
185 phase->m_timeLevelOffset[0] = 0;
187 phase->CheckAndVerify();
192 size_t order, std::vector<NekDouble> freeParams)
195 phase->m_order = order;
197 std::string(
"DIRKOrder" + std::to_string(phase->m_order) + variant);
199 phase->m_numsteps = 1;
200 char const &stage = variant.back();
201 phase->m_numstages = std::atoi(&stage);
216 switch (phase->m_order)
224 std::string(
"DIRKOrder3_ES" +
225 std::to_string(phase->m_numstages) +
229 (freeParams.size()) ? freeParams[0] : 9.0 / 40.0;
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));
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];
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));
249 phase->m_A[0][3][3] = lambda;
250 phase->m_A[0][4][3] = 5827.0 / 7560.0;
252 phase->m_A[0][4][4] = lambda;
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];
268 std::string(
"DIRKOrder4_ES" +
269 std::to_string(phase->m_numstages) +
273 (freeParams.size()) ? freeParams[0] : 1.0 / 4.0;
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;
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];
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;
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;
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;
306 phase->m_A[0][5][5] = lambda;
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];
318 ASSERTL0(
false, std::string(
"ESDIRK of order" +
319 std::to_string(phase->m_order) +
325 phase->m_numMultiStepValues = 1;
326 phase->m_numMultiStepImplicitDerivs = 0;
327 phase->m_numMultiStepExplicitDerivs = 0;
329 phase->m_timeLevelOffset[0] = 0;
331 phase->CheckAndVerify();
337 return std::string(
"DIRK");
353 std::vector<NekDouble> freeParams)
356 boost::ignore_unused(variant, order);
360 [[maybe_unused]] std::string variant, [[maybe_unused]]
size_t order,
361 std::vector<NekDouble> freeParams)
381 std::vector<NekDouble> freeParams)
384 boost::ignore_unused(variant, order);
388 [[maybe_unused]] std::string variant, [[maybe_unused]]
size_t order,
389 std::vector<NekDouble> freeParams)
409 std::vector<NekDouble> freeParams)
412 boost::ignore_unused(variant, order);
416 [[maybe_unused]] std::string variant, [[maybe_unused]]
size_t order,
417 std::vector<NekDouble> freeParams)
437 std::vector<NekDouble> freeParams)
440 boost::ignore_unused(variant, order);
444 [[maybe_unused]] std::string variant, [[maybe_unused]]
size_t order,
445 std::vector<NekDouble> freeParams)
449 "ES5", 3, freeParams);
465 std::vector<NekDouble> freeParams)
468 boost::ignore_unused(variant, order);
472 [[maybe_unused]] std::string variant, [[maybe_unused]]
size_t order,
473 std::vector<NekDouble> freeParams)
477 "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)
LUE NekDouble v_GetTimeStability() const override
LUE std::string v_GetName() const override
DIRKTimeIntegrationScheme(std::string variant, size_t order, std::vector< NekDouble > freeParams)
~DIRKTimeIntegrationScheme() override
static TimeIntegrationSchemeSharedPtr create(std::string variant, size_t order, std::vector< NekDouble > freeParams)
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
scalarT< T > sqrt(scalarT< T > in)