Nektar++
Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | List of all members
Nektar::LibUtilities::IMEXdirkTimeIntegrationScheme Class Reference

#include <IMEXdirkTimeIntegrationSchemes.h>

Inheritance diagram for Nektar::LibUtilities::IMEXdirkTimeIntegrationScheme:
[legend]

Public Member Functions

 IMEXdirkTimeIntegrationScheme (std::string variant, size_t order, std::vector< NekDouble > freeParams)
 
 ~IMEXdirkTimeIntegrationScheme () override
 
- Public Member Functions inherited from Nektar::LibUtilities::TimeIntegrationSchemeGLM
LUE void InitializeSecondaryData (TimeIntegrationAlgorithmGLM *phase, NekDouble deltaT) const
 
- Public Member Functions inherited from Nektar::LibUtilities::TimeIntegrationScheme
LUE std::string GetFullName () const
 
LUE std::string GetName () const
 
LUE std::string GetVariant () const
 
LUE size_t GetOrder () const
 
LUE std::vector< NekDoubleGetFreeParams ()
 
LUE TimeIntegrationSchemeType GetIntegrationSchemeType ()
 
LUE NekDouble GetTimeStability () const
 
LUE size_t GetNumIntegrationPhases ()
 
LUE const TripleArrayGetSolutionVector () const
 Gets the solution vector of the ODE. More...
 
LUE TripleArrayUpdateSolutionVector ()
 
LUE void SetSolutionVector (const size_t Offset, const DoubleArray &y)
 Sets the solution vector of the ODE. More...
 
LUE void InitializeScheme (const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time, const TimeIntegrationSchemeOperators &op)
 Explicit integration of an ODE. More...
 
LUE ConstDoubleArrayTimeIntegrate (const size_t timestep, const NekDouble delta_t)
 
LUE void print (std::ostream &os) const
 
LUE void printFull (std::ostream &os) const
 

Static Public Member Functions

static TimeIntegrationSchemeSharedPtr create (std::string variant, size_t order, std::vector< NekDouble > freeParams)
 
static LUE void SetupSchemeData (TimeIntegrationAlgorithmGLMSharedPtr &phase, size_t order, std::vector< NekDouble > freeParams)
 
static LUE void SetupSchemeData_1_1_1 (TimeIntegrationAlgorithmGLMSharedPtr &phase)
 
static LUE void SetupSchemeData_1_2_1 (TimeIntegrationAlgorithmGLMSharedPtr &phase)
 
static LUE void SetupSchemeData_1_2_2 (TimeIntegrationAlgorithmGLMSharedPtr &phase)
 
static LUE void SetupSchemeData_2_2_2 (TimeIntegrationAlgorithmGLMSharedPtr &phase)
 
static LUE void SetupSchemeData_2_3_2 (TimeIntegrationAlgorithmGLMSharedPtr &phase)
 
static LUE void SetupSchemeData_2_3_3 (TimeIntegrationAlgorithmGLMSharedPtr &phase)
 
static LUE void SetupSchemeData_3_4_3 (TimeIntegrationAlgorithmGLMSharedPtr &phase)
 
static LUE void SetupSchemeData_4_4_3 (TimeIntegrationAlgorithmGLMSharedPtr &phase)
 

Static Public Attributes

static std::string className
 

Protected Member Functions

LUE std::string v_GetFullName () const override
 
LUE std::string v_GetName () const override
 
LUE NekDouble v_GetTimeStability () const override
 
- Protected Member Functions inherited from Nektar::LibUtilities::TimeIntegrationSchemeGLM
LUE std::string v_GetVariant () const override
 
LUE size_t v_GetOrder () const override
 
LUE std::vector< NekDoublev_GetFreeParams () const override
 
LUE TimeIntegrationSchemeType v_GetIntegrationSchemeType () const override
 
LUE size_t v_GetNumIntegrationPhases () const override
 
LUE const TripleArrayv_GetSolutionVector () const override
 
LUE TripleArrayv_UpdateSolutionVector () override
 
LUE void v_SetSolutionVector (const size_t Offset, const DoubleArray &y) override
 
LUE void v_InitializeScheme (const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time, const TimeIntegrationSchemeOperators &op) override
 Worker method to initialize the integration scheme. More...
 
LUE ConstDoubleArrayv_TimeIntegrate (const size_t timestep, const NekDouble delta_t) override
 Worker method that actually does the time integration. More...
 
virtual LUE void v_InitializeSecondaryData (TimeIntegrationAlgorithmGLM *phase, NekDouble deltaT) const
 
LUE void v_print (std::ostream &os) const override
 Worker method to print details on the integration scheme. More...
 
LUE void v_printFull (std::ostream &os) const override
 
LUE TimeIntegrationSchemeGLM (std::string variant, size_t order, std::vector< NekDouble > freeParams)
 
 ~TimeIntegrationSchemeGLM () override
 
- Protected Member Functions inherited from Nektar::LibUtilities::TimeIntegrationScheme
virtual LUE std::string v_GetFullName () const
 
virtual LUE std::string v_GetName () const =0
 
virtual LUE std::string v_GetVariant () const =0
 
virtual LUE size_t v_GetOrder () const =0
 
virtual LUE std::vector< NekDoublev_GetFreeParams () const =0
 
virtual LUE TimeIntegrationSchemeType v_GetIntegrationSchemeType () const =0
 
virtual LUE NekDouble v_GetTimeStability () const =0
 
virtual LUE size_t v_GetNumIntegrationPhases () const =0
 
virtual LUE const TripleArrayv_GetSolutionVector () const =0
 
virtual LUE TripleArrayv_UpdateSolutionVector ()=0
 
virtual LUE void v_SetSolutionVector (const size_t Offset, const DoubleArray &y)=0
 
virtual LUE void v_InitializeScheme (const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time, const TimeIntegrationSchemeOperators &op)=0
 
virtual LUE ConstDoubleArrayv_TimeIntegrate (const size_t timestep, const NekDouble delta_t)=0
 
virtual LUE void v_print (std::ostream &os) const =0
 
virtual LUE void v_printFull (std::ostream &os) const =0
 
LUE TimeIntegrationScheme (std::string variant, size_t order, std::vector< NekDouble > freeParams)
 
LUE TimeIntegrationScheme (const TimeIntegrationScheme &in)=delete
 
virtual ~TimeIntegrationScheme ()=default
 

Additional Inherited Members

- Protected Attributes inherited from Nektar::LibUtilities::TimeIntegrationSchemeGLM
TimeIntegrationAlgorithmGLMVector m_integration_phases
 
TimeIntegrationSolutionGLMSharedPtr m_solVector
 

Detailed Description

Definition at line 58 of file IMEXdirkTimeIntegrationSchemes.h.

Constructor & Destructor Documentation

◆ IMEXdirkTimeIntegrationScheme()

Nektar::LibUtilities::IMEXdirkTimeIntegrationScheme::IMEXdirkTimeIntegrationScheme ( std::string  variant,
size_t  order,
std::vector< NekDouble freeParams 
)
inline

Definition at line 61 of file IMEXdirkTimeIntegrationSchemes.h.

63 : TimeIntegrationSchemeGLM(variant, order, freeParams)
64 {
65 ASSERTL1(freeParams.size() == 2,
66 "IMEX DIRK Time integration scheme invalid number "
67 "of free parameters, expected two "
68 "<implicit stages, explicit stages>, received " +
69 std::to_string(freeParams.size()));
70
71 size_t s = freeParams[0];
72 size_t sigma = freeParams[1];
73
76 new TimeIntegrationAlgorithmGLM(this));
77
78 if (order == 1 && s == 1 && sigma == 1)
79 {
80 // This phase is Forward Backward Euler which has two steps.
83 }
84 else
85 {
87 m_integration_phases[0], order, freeParams);
88 }
89 }
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
static LUE void SetupSchemeData(TimeIntegrationAlgorithmGLMSharedPtr &phase, size_t order, std::vector< NekDouble > freeParams)
static LUE void SetupSchemeData_1_1_1(TimeIntegrationAlgorithmGLMSharedPtr &phase)
TimeIntegrationAlgorithmGLMVector m_integration_phases
LUE TimeIntegrationSchemeGLM(std::string variant, size_t order, std::vector< NekDouble > freeParams)
std::shared_ptr< TimeIntegrationAlgorithmGLM > TimeIntegrationAlgorithmGLMSharedPtr
std::vector< TimeIntegrationAlgorithmGLMSharedPtr > TimeIntegrationAlgorithmGLMVector

References ASSERTL1, Nektar::LibUtilities::TimeIntegrationSchemeGLM::m_integration_phases, SetupSchemeData(), and SetupSchemeData_1_1_1().

◆ ~IMEXdirkTimeIntegrationScheme()

Nektar::LibUtilities::IMEXdirkTimeIntegrationScheme::~IMEXdirkTimeIntegrationScheme ( )
inlineoverride

Definition at line 91 of file IMEXdirkTimeIntegrationSchemes.h.

92 {
93 }

Member Function Documentation

◆ create()

static TimeIntegrationSchemeSharedPtr Nektar::LibUtilities::IMEXdirkTimeIntegrationScheme::create ( std::string  variant,
size_t  order,
std::vector< NekDouble freeParams 
)
inlinestatic

Definition at line 95 of file IMEXdirkTimeIntegrationSchemes.h.

97 {
100 variant, order, freeParams);
101
102 return p;
103 }
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< TimeIntegrationScheme > TimeIntegrationSchemeSharedPtr

References Nektar::MemoryManager< DataType >::AllocateSharedPtr(), and CellMLToNektar.cellml_metadata::p.

◆ SetupSchemeData()

static LUE void Nektar::LibUtilities::IMEXdirkTimeIntegrationScheme::SetupSchemeData ( TimeIntegrationAlgorithmGLMSharedPtr phase,
size_t  order,
std::vector< NekDouble freeParams 
)
inlinestatic

Definition at line 107 of file IMEXdirkTimeIntegrationSchemes.h.

110 {
111 ASSERTL1(freeParams.size() == 2,
112 "IMEX DIRK Time integration scheme invalid number "
113 "of free parameters, expected two "
114 "<implicit stages, explicit stages>, received " +
115 std::to_string(freeParams.size()) + ".");
116
117 size_t s = freeParams[0];
118 size_t sigma = freeParams[1];
119
120 phase->m_schemeType = eIMEX;
121 phase->m_variant = "DIRK";
122 phase->m_order = order;
123 phase->m_name = "IMEX_DIRK"
124 "_" +
125 std::to_string(s) + "_" + std::to_string(sigma) + "_" +
126 std::to_string(phase->m_order);
127
128 phase->m_numsteps = 1;
129 phase->m_numstages = s + 1;
130
131 phase->m_A = Array<OneD, Array<TwoD, NekDouble>>(2);
132 phase->m_B = Array<OneD, Array<TwoD, NekDouble>>(2);
133
134 phase->m_A[0] =
135 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numstages, 0.0);
136 phase->m_B[0] =
137 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numstages, 0.0);
138 phase->m_A[1] =
139 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numstages, 0.0);
140 phase->m_B[1] =
141 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numstages, 0.0);
142 phase->m_U =
143 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numsteps, 1.0);
144 phase->m_V =
145 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numsteps, 1.0);
146
147 if (s == 1 && sigma == 2 && order == 1)
148 {
150 }
151 else if (s == 1 && sigma == 2 && order == 2)
152 {
154 }
155 else if (s == 2 && sigma == 2 && order == 2)
156 {
158 }
159 else if (s == 2 && sigma == 3 && order == 2)
160 {
162 }
163 else if (s == 2 && sigma == 3 && order == 3)
164 {
166 }
167 else if (s == 3 && sigma == 4 && order == 3)
168 {
170 }
171 else if (s == 4 && sigma == 4 && order == 3)
172 {
174 }
175 else
176 {
177 ASSERTL1(false,
178 "IMEX DIRK Time integration scheme bad type "
179 "(s, sigma, order) : (" +
180 std::to_string(s) + "," + std::to_string(sigma) + "," +
181 std::to_string(order) + "). " +
182 "Allowed combinations: (1,1,1), (1,2,1), (1,2,2), "
183 "(2,2,2), (2,3,2), (2,3,3), (3,4,3), (4,4,3)");
184 }
185
186 phase->m_numMultiStepValues = 1;
187 phase->m_numMultiStepImplicitDerivs = 0;
188 phase->m_numMultiStepExplicitDerivs = 0;
189 phase->m_timeLevelOffset = Array<OneD, size_t>(phase->m_numsteps);
190 phase->m_timeLevelOffset[0] = 0;
191
192 phase->CheckAndVerify();
193 }
static LUE void SetupSchemeData_2_2_2(TimeIntegrationAlgorithmGLMSharedPtr &phase)
static LUE void SetupSchemeData_1_2_2(TimeIntegrationAlgorithmGLMSharedPtr &phase)
static LUE void SetupSchemeData_2_3_3(TimeIntegrationAlgorithmGLMSharedPtr &phase)
static LUE void SetupSchemeData_2_3_2(TimeIntegrationAlgorithmGLMSharedPtr &phase)
static LUE void SetupSchemeData_3_4_3(TimeIntegrationAlgorithmGLMSharedPtr &phase)
static LUE void SetupSchemeData_4_4_3(TimeIntegrationAlgorithmGLMSharedPtr &phase)
static LUE void SetupSchemeData_1_2_1(TimeIntegrationAlgorithmGLMSharedPtr &phase)
@ eIMEX
Implicit Explicit General Linear Method.

References ASSERTL1, Nektar::LibUtilities::eIMEX, SetupSchemeData_1_2_1(), SetupSchemeData_1_2_2(), SetupSchemeData_2_2_2(), SetupSchemeData_2_3_2(), SetupSchemeData_2_3_3(), SetupSchemeData_3_4_3(), and SetupSchemeData_4_4_3().

Referenced by Nektar::LibUtilities::CNABTimeIntegrationScheme::CNABTimeIntegrationScheme(), IMEXdirkTimeIntegrationScheme(), Nektar::LibUtilities::IMEXGearTimeIntegrationScheme::IMEXGearTimeIntegrationScheme(), Nektar::LibUtilities::IMEXTimeIntegrationScheme::IMEXTimeIntegrationScheme(), and Nektar::LibUtilities::MCNABTimeIntegrationScheme::MCNABTimeIntegrationScheme().

◆ SetupSchemeData_1_1_1()

static LUE void Nektar::LibUtilities::IMEXdirkTimeIntegrationScheme::SetupSchemeData_1_1_1 ( TimeIntegrationAlgorithmGLMSharedPtr phase)
inlinestatic

Definition at line 197 of file IMEXdirkTimeIntegrationSchemes.h.

199 {
200 phase->m_schemeType = eIMEX;
201 phase->m_order = 1;
202 phase->m_name =
203 std::string("IMEXdirk_1_1_" + std::to_string(phase->m_order));
204
205 phase->m_numsteps = 2;
206 phase->m_numstages = 1;
207
208 phase->m_A = Array<OneD, Array<TwoD, NekDouble>>(2);
209 phase->m_B = Array<OneD, Array<TwoD, NekDouble>>(2);
210
211 phase->m_A[0] =
212 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numstages, 0.0);
213 phase->m_B[0] =
214 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numstages, 0.0);
215 phase->m_A[1] =
216 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numstages, 0.0);
217 phase->m_B[1] =
218 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numstages, 0.0);
219 phase->m_U =
220 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numsteps, 0.0);
221 phase->m_V =
222 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numsteps, 0.0);
223
224 phase->m_A[0][0][0] = 1.0;
225
226 phase->m_B[0][0][0] = 1.0;
227 phase->m_B[1][1][0] = 1.0;
228
229 phase->m_U[0][0] = 1.0;
230 phase->m_U[0][1] = 1.0;
231
232 phase->m_V[0][0] = 1.0;
233 phase->m_V[0][1] = 1.0;
234
235 phase->m_numMultiStepValues = 1;
236 phase->m_numMultiStepImplicitDerivs = 0;
237 phase->m_numMultiStepExplicitDerivs = 1;
238
239 phase->m_timeLevelOffset = Array<OneD, size_t>(phase->m_numsteps);
240 phase->m_timeLevelOffset[0] = 0;
241 phase->m_timeLevelOffset[1] = 0;
242
243 phase->CheckAndVerify();
244 }

References Nektar::LibUtilities::eIMEX.

Referenced by IMEXdirkTimeIntegrationScheme(), and Nektar::LibUtilities::IMEXTimeIntegrationScheme::IMEXTimeIntegrationScheme().

◆ SetupSchemeData_1_2_1()

static LUE void Nektar::LibUtilities::IMEXdirkTimeIntegrationScheme::SetupSchemeData_1_2_1 ( TimeIntegrationAlgorithmGLMSharedPtr phase)
inlinestatic

Definition at line 249 of file IMEXdirkTimeIntegrationSchemes.h.

251 {
252 phase->m_A[0][1][1] = 1.0;
253 phase->m_B[0][0][1] = 1.0;
254
255 phase->m_A[1][1][0] = 1.0;
256 phase->m_B[1][0][1] = 1.0;
257
258 // U and V set to 1 when allocated.
259 }

Referenced by SetupSchemeData().

◆ SetupSchemeData_1_2_2()

static LUE void Nektar::LibUtilities::IMEXdirkTimeIntegrationScheme::SetupSchemeData_1_2_2 ( TimeIntegrationAlgorithmGLMSharedPtr phase)
inlinestatic

Definition at line 263 of file IMEXdirkTimeIntegrationSchemes.h.

265 {
266 phase->m_A[0][1][1] = 1.0 / 2.0;
267 phase->m_B[0][0][1] = 1.0;
268
269 phase->m_A[1][1][0] = 1.0 / 2.0;
270 phase->m_B[1][0][1] = 1.0;
271
272 // U and V set to 1 when allocated.
273 }

Referenced by SetupSchemeData().

◆ SetupSchemeData_2_2_2()

static LUE void Nektar::LibUtilities::IMEXdirkTimeIntegrationScheme::SetupSchemeData_2_2_2 ( TimeIntegrationAlgorithmGLMSharedPtr phase)
inlinestatic

Definition at line 277 of file IMEXdirkTimeIntegrationSchemes.h.

279 {
280 const NekDouble glambda = 1.0 - sqrt(2.0) / 2.0;
281 const NekDouble gdelta = -sqrt(2.0) / 2.0;
282
283 phase->m_A[0][1][1] = glambda;
284 phase->m_A[0][2][1] = 1.0 - glambda;
285 phase->m_A[0][2][2] = glambda;
286
287 phase->m_B[0][0][1] = 1.0 - glambda;
288 phase->m_B[0][0][2] = glambda;
289
290 phase->m_A[1][1][0] = glambda;
291 phase->m_A[1][2][0] = gdelta;
292 phase->m_A[1][2][1] = 1.0 - gdelta;
293
294 phase->m_B[1][0][0] = gdelta;
295 phase->m_B[1][0][1] = 1.0 - gdelta;
296
297 // U and V set to 1 when allocated.
298 }
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References tinysimd::sqrt().

Referenced by SetupSchemeData().

◆ SetupSchemeData_2_3_2()

static LUE void Nektar::LibUtilities::IMEXdirkTimeIntegrationScheme::SetupSchemeData_2_3_2 ( TimeIntegrationAlgorithmGLMSharedPtr phase)
inlinestatic

Definition at line 302 of file IMEXdirkTimeIntegrationSchemes.h.

304 {
305 const NekDouble lambda = (2.0 - sqrt(2.0)) / 2.0;
306 const NekDouble delta = -2.0 * sqrt(2.0) / 3.0;
307
308 phase->m_A[0][1][1] = lambda;
309 phase->m_A[0][2][1] = 1.0 - lambda;
310 phase->m_A[0][2][2] = lambda;
311
312 phase->m_B[0][0][1] = 1.0 - lambda;
313 phase->m_B[0][0][2] = lambda;
314
315 phase->m_A[1][1][0] = lambda;
316 phase->m_A[1][2][0] = delta;
317 phase->m_A[1][2][1] = 1.0 - delta;
318
319 phase->m_B[1][0][1] = 1.0 - lambda;
320 phase->m_B[1][0][2] = lambda;
321
322 // U and V set to 1 when allocated.
323 }

References tinysimd::sqrt().

Referenced by SetupSchemeData().

◆ SetupSchemeData_2_3_3()

static LUE void Nektar::LibUtilities::IMEXdirkTimeIntegrationScheme::SetupSchemeData_2_3_3 ( TimeIntegrationAlgorithmGLMSharedPtr phase)
inlinestatic

Definition at line 327 of file IMEXdirkTimeIntegrationSchemes.h.

329 {
330 const NekDouble glambda = (3.0 + sqrt(3.0)) / 6.0;
331
332 phase->m_A[0][1][1] = glambda;
333 phase->m_A[0][2][1] = 1.0 - 2.0 * glambda;
334 phase->m_A[0][2][2] = glambda;
335
336 phase->m_B[0][0][1] = 0.5;
337 phase->m_B[0][0][2] = 0.5;
338
339 phase->m_A[1][1][0] = glambda;
340 phase->m_A[1][2][0] = glambda - 1.0;
341 phase->m_A[1][2][1] = 2.0 * (1.0 - glambda);
342
343 phase->m_B[1][0][1] = 0.5;
344 phase->m_B[1][0][2] = 0.5;
345
346 // U and V set to 1 when allocated.
347 }

References tinysimd::sqrt().

Referenced by SetupSchemeData().

◆ SetupSchemeData_3_4_3()

static LUE void Nektar::LibUtilities::IMEXdirkTimeIntegrationScheme::SetupSchemeData_3_4_3 ( TimeIntegrationAlgorithmGLMSharedPtr phase)
inlinestatic

Definition at line 351 of file IMEXdirkTimeIntegrationSchemes.h.

353 {
354 constexpr NekDouble lambda = 0.4358665215;
355
356 phase->m_A[0][1][1] = lambda;
357 phase->m_A[0][2][1] = 0.5 * (1.0 - lambda);
358 phase->m_A[0][2][2] = lambda;
359 phase->m_A[0][3][1] =
360 0.25 * (-6.0 * lambda * lambda + 16.0 * lambda - 1.0);
361 phase->m_A[0][3][2] =
362 0.25 * (6.0 * lambda * lambda - 20.0 * lambda + 5.0);
363 phase->m_A[0][3][3] = lambda;
364
365 phase->m_B[0][0][1] =
366 0.25 * (-6.0 * lambda * lambda + 16.0 * lambda - 1.0);
367 phase->m_B[0][0][2] =
368 0.25 * (6.0 * lambda * lambda - 20.0 * lambda + 5.0);
369 phase->m_B[0][0][3] = lambda;
370
371 phase->m_A[1][1][0] = 0.4358665215;
372 phase->m_A[1][2][0] = 0.3212788860;
373 phase->m_A[1][2][1] = 0.3966543747;
374 phase->m_A[1][3][0] = -0.105858296;
375 phase->m_A[1][3][1] = 0.5529291479;
376 phase->m_A[1][3][2] = 0.5529291479;
377
378 phase->m_B[1][0][1] =
379 0.25 * (-6.0 * lambda * lambda + 16.0 * lambda - 1.0);
380 phase->m_B[1][0][2] =
381 0.25 * (6.0 * lambda * lambda - 20.0 * lambda + 5.0);
382 phase->m_B[1][0][3] = lambda;
383
384 // U and V set to 1 when allocated.
385 }

Referenced by SetupSchemeData().

◆ SetupSchemeData_4_4_3()

static LUE void Nektar::LibUtilities::IMEXdirkTimeIntegrationScheme::SetupSchemeData_4_4_3 ( TimeIntegrationAlgorithmGLMSharedPtr phase)
inlinestatic

Definition at line 389 of file IMEXdirkTimeIntegrationSchemes.h.

391 {
392 phase->m_A[0][1][1] = 1.0 / 2.0;
393 phase->m_A[0][2][1] = 1.0 / 6.0;
394 phase->m_A[0][2][2] = 1.0 / 2.0;
395 phase->m_A[0][3][1] = -1.0 / 2.0;
396 phase->m_A[0][3][2] = 1.0 / 2.0;
397 phase->m_A[0][3][3] = 1.0 / 2.0;
398 phase->m_A[0][4][1] = 3.0 / 2.0;
399 phase->m_A[0][4][2] = -3.0 / 2.0;
400 phase->m_A[0][4][3] = 1.0 / 2.0;
401 phase->m_A[0][4][4] = 1.0 / 2.0;
402
403 phase->m_B[0][0][1] = 3.0 / 2.0;
404 phase->m_B[0][0][2] = -3.0 / 2.0;
405 phase->m_B[0][0][3] = 1.0 / 2.0;
406 phase->m_B[0][0][4] = 1.0 / 2.0;
407
408 phase->m_A[1][1][0] = 1.0 / 2.0;
409 phase->m_A[1][2][0] = 11.0 / 18.0;
410 phase->m_A[1][2][1] = 1.0 / 18.0;
411 phase->m_A[1][3][0] = 5.0 / 6.0;
412 phase->m_A[1][3][1] = -5.0 / 6.0;
413 phase->m_A[1][3][2] = 1.0 / 2.0;
414 phase->m_A[1][4][0] = 1.0 / 4.0;
415 phase->m_A[1][4][1] = 7.0 / 4.0;
416 phase->m_A[1][4][2] = 3.0 / 4.0;
417 phase->m_A[1][4][3] = -7.0 / 4.0;
418
419 phase->m_B[1][0][0] = 1.0 / 4.0;
420 phase->m_B[1][0][1] = 7.0 / 4.0;
421 phase->m_B[1][0][2] = 3.0 / 4.0;
422 phase->m_B[1][0][3] = -7.0 / 4.0;
423
424 // U and V set to 1 when allocated.
425 }

Referenced by SetupSchemeData().

◆ v_GetFullName()

LUE std::string Nektar::LibUtilities::IMEXdirkTimeIntegrationScheme::v_GetFullName ( ) const
inlineoverrideprotectedvirtual

◆ v_GetName()

LUE std::string Nektar::LibUtilities::IMEXdirkTimeIntegrationScheme::v_GetName ( ) const
inlineoverrideprotectedvirtual

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 433 of file IMEXdirkTimeIntegrationSchemes.h.

434 {
435 return std::string("IMEX");
436 }

◆ v_GetTimeStability()

LUE NekDouble Nektar::LibUtilities::IMEXdirkTimeIntegrationScheme::v_GetTimeStability ( ) const
inlineoverrideprotectedvirtual

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 438 of file IMEXdirkTimeIntegrationSchemes.h.

439 {
440 return 1.0;
441 }

Member Data Documentation

◆ className

std::string Nektar::LibUtilities::IMEXdirkTimeIntegrationScheme::className
static

Definition at line 105 of file IMEXdirkTimeIntegrationSchemes.h.