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

#include <DIRKTimeIntegrationSchemes.h>

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

Public Member Functions

 DIRKTimeIntegrationScheme (std::string variant, size_t order, std::vector< NekDouble > freeParams)
 
 ~DIRKTimeIntegrationScheme () 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)
 
static LUE void SetupSchemeDataESDIRK (TimeIntegrationAlgorithmGLMSharedPtr &phase, std::string variant, size_t order, std::vector< NekDouble > freeParams)
 

Static Public Attributes

static std::string className
 

Protected Member Functions

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 54 of file DIRKTimeIntegrationSchemes.h.

Constructor & Destructor Documentation

◆ DIRKTimeIntegrationScheme()

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

Definition at line 57 of file DIRKTimeIntegrationSchemes.h.

59 : TimeIntegrationSchemeGLM(variant, order, freeParams)
60 {
61 // Currently up to 4th order is implemented.
62 ASSERTL0(1 <= order && order <= 4,
63 "Diagonally Implicit Runge Kutta integration scheme bad order "
64 "(1-4): " +
65 std::to_string(order));
66
69 new TimeIntegrationAlgorithmGLM(this));
70
71 ASSERTL1((variant == "" || variant == "ES5" || variant == "ES6"),
72 "DIRK Time integration scheme bad variant: " + variant +
73 ". "
74 "Must blank, 'ES5', or 'ES6'");
75
76 if ("" == variant)
77 {
79 order);
80 }
81 else
82 {
84 m_integration_phases[0], variant, order, freeParams);
85 }
86 }
#define ASSERTL0(condition, msg)
Definition: ErrorUtil.hpp:208
#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 SetupSchemeDataESDIRK(TimeIntegrationAlgorithmGLMSharedPtr &phase, std::string variant, size_t order, std::vector< NekDouble > freeParams)
static LUE void SetupSchemeData(TimeIntegrationAlgorithmGLMSharedPtr &phase, size_t order)
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 ASSERTL0, ASSERTL1, Nektar::LibUtilities::TimeIntegrationSchemeGLM::m_integration_phases, SetupSchemeData(), and SetupSchemeDataESDIRK().

◆ ~DIRKTimeIntegrationScheme()

Nektar::LibUtilities::DIRKTimeIntegrationScheme::~DIRKTimeIntegrationScheme ( )
inlineoverride

Definition at line 88 of file DIRKTimeIntegrationSchemes.h.

89 {
90 }

Member Function Documentation

◆ create()

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

Definition at line 92 of file DIRKTimeIntegrationSchemes.h.

94 {
97 variant, order, freeParams);
98 return p;
99 }
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::DIRKTimeIntegrationScheme::SetupSchemeData ( TimeIntegrationAlgorithmGLMSharedPtr phase,
size_t  order 
)
inlinestatic

Definition at line 103 of file DIRKTimeIntegrationSchemes.h.

105 {
106 phase->m_schemeType = eDiagonallyImplicit;
107 phase->m_order = order;
108 phase->m_name =
109 std::string("DIRKOrder" + std::to_string(phase->m_order));
110
111 phase->m_numsteps = 1;
112 phase->m_numstages = phase->m_order;
113
114 phase->m_A = Array<OneD, Array<TwoD, NekDouble>>(1);
115 phase->m_B = Array<OneD, Array<TwoD, NekDouble>>(1);
116
117 phase->m_A[0] =
118 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numstages, 0.0);
119 phase->m_B[0] =
120 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numstages, 0.0);
121 phase->m_U =
122 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numsteps, 1.0);
123 phase->m_V =
124 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numsteps, 1.0);
125
126 switch (phase->m_order)
127 {
128 case 1:
129 {
130 // One-stage, 1st order, L-stable Diagonally Implicit
131 // Runge Kutta (backward Euler) method:
132
133 phase->m_A[0][0][0] = 1.0;
134
135 phase->m_B[0][0][0] = 1.0;
136 }
137 break;
138 case 2:
139 {
140 // Two-stage, 2nd order Diagonally Implicit Runge
141 // Kutta method. It is A-stable if and only if x ≥ 1/4
142 // and is L-stable if and only if x equals one of the
143 // roots of the polynomial x^2 - 2x + 1/2, i.e. if
144 // x = 1 ± sqrt(2)/2.
145 const NekDouble lambda = (2.0 - sqrt(2.0)) / 2.0;
146
147 phase->m_A[0][0][0] = lambda;
148 phase->m_A[0][1][0] = 1.0 - lambda;
149 phase->m_A[0][1][1] = lambda;
150
151 phase->m_B[0][0][0] = 1.0 - lambda;
152 phase->m_B[0][0][1] = lambda;
153 }
154 break;
155
156 case 3:
157 {
158 // Three-stage, 3rd order, L-stable Diagonally Implicit
159 // Runge Kutta method:
160 constexpr NekDouble lambda = 0.4358665215;
161
162 phase->m_A[0][0][0] = lambda;
163 phase->m_A[0][1][0] = 0.5 * (1.0 - lambda);
164 phase->m_A[0][2][0] =
165 0.25 * (-6.0 * lambda * lambda + 16.0 * lambda - 1.0);
166
167 phase->m_A[0][1][1] = lambda;
168
169 phase->m_A[0][2][1] =
170 0.25 * (6.0 * lambda * lambda - 20.0 * lambda + 5.0);
171 phase->m_A[0][2][2] = lambda;
172
173 phase->m_B[0][0][0] =
174 0.25 * (-6.0 * lambda * lambda + 16.0 * lambda - 1.0);
175 phase->m_B[0][0][1] =
176 0.25 * (6.0 * lambda * lambda - 20.0 * lambda + 5.0);
177 phase->m_B[0][0][2] = lambda;
178 }
179 break;
180 }
181
182 phase->m_numMultiStepValues = 1;
183 phase->m_numMultiStepImplicitDerivs = 0;
184 phase->m_numMultiStepExplicitDerivs = 0;
185 phase->m_timeLevelOffset = Array<OneD, size_t>(phase->m_numsteps);
186 phase->m_timeLevelOffset[0] = 0;
187
188 phase->CheckAndVerify();
189 }
@ eDiagonallyImplicit
Diagonally implicit scheme (e.g. the DIRK schemes)
double NekDouble
scalarT< T > sqrt(scalarT< T > in)
Definition: scalar.hpp:294

References Nektar::LibUtilities::eDiagonallyImplicit, and tinysimd::sqrt().

Referenced by Nektar::LibUtilities::AdamsMoultonTimeIntegrationScheme::AdamsMoultonTimeIntegrationScheme(), Nektar::LibUtilities::BDFImplicitTimeIntegrationScheme::BDFImplicitTimeIntegrationScheme(), and DIRKTimeIntegrationScheme().

◆ SetupSchemeDataESDIRK()

static LUE void Nektar::LibUtilities::DIRKTimeIntegrationScheme::SetupSchemeDataESDIRK ( TimeIntegrationAlgorithmGLMSharedPtr phase,
std::string  variant,
size_t  order,
std::vector< NekDouble freeParams 
)
inlinestatic

Definition at line 191 of file DIRKTimeIntegrationSchemes.h.

194 {
195 phase->m_schemeType = eDiagonallyImplicit;
196 phase->m_order = order;
197 phase->m_name =
198 std::string("DIRKOrder" + std::to_string(phase->m_order) + variant);
199
200 phase->m_numsteps = 1;
201 char const &stage = variant.back();
202 phase->m_numstages = std::atoi(&stage);
203
204 phase->m_A = Array<OneD, Array<TwoD, NekDouble>>(1);
205 phase->m_B = Array<OneD, Array<TwoD, NekDouble>>(1);
206
207 phase->m_A[0] =
208 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numstages, 0.0);
209 phase->m_B[0] =
210 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numstages, 0.0);
211 phase->m_U =
212 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numsteps, 1.0);
213 phase->m_V =
214 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numsteps, 1.0);
215
216 const NekDouble ConstSqrt2 = sqrt(2.0);
217 switch (phase->m_order)
218 {
219 case 3:
220 {
221 // See: Kennedy, Christopher A., and Mark H. Carpenter.
222 // Diagonally implicit Runge-Kutta methods for ordinary
223 // differential equations. A review. No. NF1676L-19716. 2016.
224 ASSERTL0(5 == phase->m_numstages,
225 std::string("DIRKOrder3_ES" +
226 std::to_string(phase->m_numstages) +
227 " not defined"));
228
229 NekDouble lambda =
230 (freeParams.size()) ? freeParams[0] : 9.0 / 40.0;
231
232 phase->m_A[0][0][0] = 0.0;
233 phase->m_A[0][1][0] = lambda;
234 phase->m_A[0][2][0] = 9.0 * (1.0 + ConstSqrt2) / 80.0;
235 phase->m_A[0][3][0] =
236 (22.0 + 15.0 * ConstSqrt2) / (80.0 * (1.0 + ConstSqrt2));
237 phase->m_A[0][4][0] = (2398.0 + 1205.0 * ConstSqrt2) /
238 (2835.0 * (4.0 + 3.0 * ConstSqrt2));
239
240 phase->m_A[0][1][1] = phase->m_A[0][1][0];
241 phase->m_A[0][2][1] = phase->m_A[0][2][0];
242 phase->m_A[0][3][1] = phase->m_A[0][3][0];
243 phase->m_A[0][4][1] = phase->m_A[0][4][0];
244
245 phase->m_A[0][2][2] = lambda;
246 phase->m_A[0][3][2] = -7.0 / (40.0 * (1.0 + ConstSqrt2));
247 phase->m_A[0][4][2] = -2374.0 * (1.0 + 2.0 * ConstSqrt2) /
248 (2835.0 * (5.0 + 3.0 * ConstSqrt2));
249
250 phase->m_A[0][3][3] = lambda;
251 phase->m_A[0][4][3] = 5827.0 / 7560.0;
252
253 phase->m_A[0][4][4] = lambda;
254
255 phase->m_B[0][0][0] = phase->m_A[0][4][0];
256 phase->m_B[0][0][1] = phase->m_A[0][4][1];
257 phase->m_B[0][0][2] = phase->m_A[0][4][2];
258 phase->m_B[0][0][3] = phase->m_A[0][4][3];
259 phase->m_B[0][0][4] = phase->m_A[0][4][4];
260 }
261 break;
262
263 case 4:
264 {
265 // See: Kennedy, Christopher A., and Mark H. Carpenter.
266 // Diagonally implicit Runge-Kutta methods for ordinary
267 // differential equations. A review. No. NF1676L-19716. 2016.
268 ASSERTL0(6 == phase->m_numstages,
269 std::string("DIRKOrder4_ES" +
270 std::to_string(phase->m_numstages) +
271 " not defined"));
272
273 NekDouble lambda =
274 (freeParams.size()) ? freeParams[0] : 1.0 / 4.0;
275
276 phase->m_A[0][0][0] = 0.0;
277 phase->m_A[0][1][0] = lambda;
278 phase->m_A[0][2][0] = (1.0 - ConstSqrt2) / 8.0;
279 phase->m_A[0][3][0] = (5.0 - 7.0 * ConstSqrt2) / 64.0;
280 phase->m_A[0][4][0] =
281 (-13796.0 - 54539.0 * ConstSqrt2) / 125000.0;
282 phase->m_A[0][5][0] = (1181.0 - 987.0 * ConstSqrt2) / 13782.0;
283
284 phase->m_A[0][1][1] = phase->m_A[0][1][0];
285 phase->m_A[0][2][1] = phase->m_A[0][2][0];
286 phase->m_A[0][3][1] = phase->m_A[0][3][0];
287 phase->m_A[0][4][1] = phase->m_A[0][4][0];
288 phase->m_A[0][5][1] = phase->m_A[0][5][0];
289
290 phase->m_A[0][2][2] = lambda;
291 phase->m_A[0][3][2] = 7.0 * (1.0 + ConstSqrt2) / 32.0;
292 phase->m_A[0][4][2] =
293 (506605.0 + 132109.0 * ConstSqrt2) / 437500.0;
294 phase->m_A[0][5][2] =
295 47.0 * (-267.0 + 1783.0 * ConstSqrt2) / 273343.0;
296
297 phase->m_A[0][3][3] = lambda;
298 phase->m_A[0][4][3] =
299 166.0 * (-97.0 + 376.0 * ConstSqrt2) / 109375.0;
300 phase->m_A[0][5][3] =
301 -16.0 * (-22922.0 + 3525.0 * ConstSqrt2) / 571953.0;
302
303 phase->m_A[0][4][4] = lambda;
304 phase->m_A[0][5][4] =
305 -15625.0 * (97.0 + 376.0 * ConstSqrt2) / 90749876.0;
306
307 phase->m_A[0][5][5] = lambda;
308
309 phase->m_B[0][0][0] = phase->m_A[0][5][0];
310 phase->m_B[0][0][1] = phase->m_A[0][5][1];
311 phase->m_B[0][0][2] = phase->m_A[0][5][2];
312 phase->m_B[0][0][3] = phase->m_A[0][5][3];
313 phase->m_B[0][0][4] = phase->m_A[0][5][4];
314 phase->m_B[0][0][5] = phase->m_A[0][5][5];
315 }
316 break;
317 default:
318 {
319 ASSERTL0(false, std::string("ESDIRK of order" +
320 std::to_string(phase->m_order) +
321 " not defined"));
322 break;
323 }
324 }
325
326 phase->m_numMultiStepValues = 1;
327 phase->m_numMultiStepImplicitDerivs = 0;
328 phase->m_numMultiStepExplicitDerivs = 0;
329 phase->m_timeLevelOffset = Array<OneD, size_t>(phase->m_numsteps);
330 phase->m_timeLevelOffset[0] = 0;
331
332 phase->CheckAndVerify();
333 }

References ASSERTL0, Nektar::LibUtilities::eDiagonallyImplicit, and tinysimd::sqrt().

Referenced by DIRKTimeIntegrationScheme().

◆ v_GetName()

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

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 336 of file DIRKTimeIntegrationSchemes.h.

337 {
338 return std::string("DIRK");
339 }

◆ v_GetTimeStability()

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

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 341 of file DIRKTimeIntegrationSchemes.h.

342 {
343 return 1.0;
344 }

Member Data Documentation

◆ className

std::string Nektar::LibUtilities::DIRKTimeIntegrationScheme::className
static

Definition at line 101 of file DIRKTimeIntegrationSchemes.h.