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);
109 return std::string(
"DIRK");
121 phase->m_order = order;
123 std::string(
"DIRKOrder" + std::to_string(phase->m_order));
125 phase->m_numsteps = 1;
126 phase->m_numstages = phase->m_order;
140 switch (phase->m_order)
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_numMultiStepDerivs = 0;
189 phase->m_timeLevelOffset[0] = 0;
191 phase->CheckAndVerify();
196 unsigned int order, std::vector<NekDouble> freeParams)
199 phase->m_order = order;
201 std::string(
"DIRKOrder" + std::to_string(phase->m_order) + variant);
203 phase->m_numsteps = 1;
204 char const &stage = variant.back();
205 phase->m_numstages = std::atoi(&stage);
220 switch (phase->m_order)
225 std::string(
"DIRKOrder3_ES" +
226 std::to_string(phase->m_numstages) +
229 if (freeParams.size())
231 lambda = freeParams[0];
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));
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];
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));
256 phase->m_A[0][3][3] = lambda;
257 phase->m_A[0][4][3] = 5827.0 / 7560.0;
259 phase->m_A[0][4][4] = lambda;
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];
272 std::string(
"DIRKOrder4_ES" +
273 std::to_string(phase->m_numstages) +
276 if (freeParams.size())
278 lambda = freeParams[0];
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;
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];
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;
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;
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;
318 phase->m_A[0][5][5] = lambda;
323 ASSERTL0(
false, std::string(
"ESDIRK of order" +
324 std::to_string(phase->m_order) +
330 phase->m_numMultiStepValues = 1;
331 phase->m_numMultiStepDerivs = 0;
333 phase->m_timeLevelOffset[0] = 0;
335 phase->CheckAndVerify();
346 std::vector<NekDouble> freeParams)
349 boost::ignore_unused(variant);
350 boost::ignore_unused(order);
354 std::string variant,
unsigned int order,
355 std::vector<NekDouble> freeParams)
357 boost::ignore_unused(variant);
358 boost::ignore_unused(order);
375 std::vector<NekDouble> freeParams)
378 boost::ignore_unused(variant);
379 boost::ignore_unused(order);
383 std::string variant,
unsigned int order,
384 std::vector<NekDouble> freeParams)
386 boost::ignore_unused(variant);
387 boost::ignore_unused(order);
404 std::vector<NekDouble> freeParams)
407 boost::ignore_unused(variant);
408 boost::ignore_unused(order);
412 std::string variant,
unsigned int order,
413 std::vector<NekDouble> freeParams)
415 boost::ignore_unused(variant);
416 boost::ignore_unused(order);
420 "ES5", 3, freeParams);
433 std::vector<NekDouble> freeParams)
436 boost::ignore_unused(variant);
437 boost::ignore_unused(order);
441 std::string variant,
unsigned int order,
442 std::vector<NekDouble> freeParams)
444 boost::ignore_unused(variant);
445 boost::ignore_unused(order);
449 "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
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 std::string GetName() const
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 ~DIRKTimeIntegrationScheme()
virtual LUE NekDouble GetTimeStability() const
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)