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

#include <RungeKuttaTimeIntegrationSchemes.h>

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

Public Member Functions

 RungeKuttaTimeIntegrationScheme (std::string variant, size_t order, std::vector< NekDouble > freeParams)
 
 ~RungeKuttaTimeIntegrationScheme () 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, 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 56 of file RungeKuttaTimeIntegrationSchemes.h.

Constructor & Destructor Documentation

◆ RungeKuttaTimeIntegrationScheme()

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

Definition at line 59 of file RungeKuttaTimeIntegrationSchemes.h.

61 : TimeIntegrationSchemeGLM(variant, order, freeParams)
62 {
63 ASSERTL1(variant == "" || variant == "SSP",
64 "Runge Kutta Time integration scheme unknown variant: " +
65 variant + ". Must be blank or 'SSP'");
66
67 // Std - Currently up to 5th order is implemented.
68 // SSP - Currently 1st through 3rd order is implemented.
69 ASSERTL1((variant == "" && 1 <= order && order <= 5) ||
70 (variant == "SSP" && 1 <= order && order <= 3),
71 "Runge Kutta Time integration scheme bad order, "
72 "Std (1-5) or SSP (1-3): " +
73 std::to_string(order));
74
77 new TimeIntegrationAlgorithmGLM(this));
78
80 m_integration_phases[0], variant, order, freeParams);
81 }
#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, std::string variant, size_t order, std::vector< NekDouble > freeParams)
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, and SetupSchemeData().

◆ ~RungeKuttaTimeIntegrationScheme()

Nektar::LibUtilities::RungeKuttaTimeIntegrationScheme::~RungeKuttaTimeIntegrationScheme ( )
inlineoverride

Definition at line 83 of file RungeKuttaTimeIntegrationSchemes.h.

84 {
85 }

Member Function Documentation

◆ create()

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

Definition at line 87 of file RungeKuttaTimeIntegrationSchemes.h.

89 {
92 variant, order, freeParams);
93
94 return p;
95 }
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::RungeKuttaTimeIntegrationScheme::SetupSchemeData ( TimeIntegrationAlgorithmGLMSharedPtr phase,
std::string  variant,
size_t  order,
std::vector< NekDouble freeParams 
)
inlinestatic

Definition at line 99 of file RungeKuttaTimeIntegrationSchemes.h.

102 {
103 constexpr size_t nStages[6] = {0, 1, 2, 3, 4, 6};
104
105 // A Coefficients for the lower diagonal quadrant stored in a
106 // contiguous fashion. For the fourth order, six coefficients
107 // from the Butcher tableau would be stored as the following.
108 //
109 // 0 0 0 0
110 // Butcher a 0 0 0 Stored as a
111 // Tableau b c 0 0 b c
112 // d e f 0 d e f 0 ... 0
113
114 // clang-format off
115 constexpr NekDouble Acoefficients[2][6][15] =
116 { { { 0., 0., 0., 0., 0.,
117 0., 0., 0., 0., 0.,
118 0., 0., 0., 0., 0. },
119 // 1st Order
120 { 0., 0., 0., 0., 0.,
121 0., 0., 0., 0., 0.,
122 0., 0., 0., 0., 0. },
123 // 2nd Order - midpoint
124 { 1./2, // Last entry
125 0., 0., 0., 0.,
126 0., 0., 0., 0., 0.,
127 0., 0., 0., 0., 0. },
128 // 3rd Order - Ralston's
129 { 1./2.,
130 0., 3./4., // Last entry
131 0., 0.,
132 0., 0., 0., 0., 0.,
133 0., 0., 0., 0., 0. },
134 // 4th Order - Classic
135 { 1./2.,
136 0., 1./2.,
137 0., 0., 1., // Last entry
138 0., 0., 0., 0.,
139 0., 0., 0., 0., 0. },
140 // 5th Order - 6 stages
141 // Rabiei, Faranak, and Fudziah Ismail. "Fifth order improved
142 // Runge-Kutta method for solving ordinary differential
143 // equations." In Proceedings of the 11th WSEAS international
144 // conference on Applied computer science, pp. 129-133. 2011.
145 { 1./4.,
146 1./8., 1./8.,
147 0., -1./2., 1.,
148 3./16., 0., 0., 9./16.,
149 -3./7., 2./7., 12./7., -12./7., 8./7. } },
150 // Strong Stability Preserving
151 { { 0., 0., 0., 0., 0.,
152 0., 0., 0., 0., 0.,
153 0., 0., 0., 0., 0. },
154 // 1st Order
155 { 0., 0., 0., 0., 0.,
156 0., 0., 0., 0., 0.,
157 0., 0., 0., 0., 0. },
158 // 2nd Order - strong scaling - improved
159 { 1., // Last entry
160 0., 0., 0., 0.,
161 0., 0., 0., 0., 0.,
162 0., 0., 0., 0., 0. },
163 // 3rd Order - strong scaling
164 { 1.,
165 1./4., 1./4., // Last entry
166 0, 0.,
167 0., 0., 0., 0., 0.,
168 0., 0., 0., 0., 0. },
169 // 4th Order - Classic - not used
170 { 1./2.,
171 0., 1./2.,
172 0., 0., 1., // Last entry
173 0., 0., 0., 0.,
174 0., 0., 0., 0., 0. },
175 // 5th Order - 6 stages - not used
176 // Rabiei, Faranak, and Fudziah Ismail. "Fifth order improved
177 // Runge-Kutta method for solving ordinary differential
178 // equations." In Proceedings of the 11th WSEAS international
179 // conference on Applied computer science, pp. 129-133. 2011.
180 { 1./4.,
181 1./8., 1./8.,
182 0., -1./2., 1.,
183 3./16., 0., 0., 9./16.,
184 -3./7., 2./7., 12./7., -12./7., 8./7. } } };
185 // clang-format on
186
187 // B Coefficients for the final summing.
188
189 // clang-format off
190 constexpr NekDouble Bcoefficients[2][6][6] =
191 { { { 0., 0., 0., 0., 0., 0. },
192 // 1st Order
193 { 1., 0., 0., 0., 0., 0. },
194 // 2nd Order - midpoint
195 { 0., 1., 0., 0., 0., 0. },
196 // 3rd Order - Ralston's
197 { 2./9., 3./9., 4./9., 0., 0., 0. },
198 // 4th Order - Classic
199 { 1./6., 2./6., 2./6., 1./6., 0., 0. },
200 // 5th Order - 6 stages
201 // Rabiei, Faranak, and Fudziah Ismail. "Fifth order improved
202 // Runge-Kutta method for solving ordinary differential
203 // equations." In Proceedings of the 11th WSEAS international
204 // conference on Applied computer science, pp. 129-133. 2011.
205 { 7./90., 0., 32./90., 12./90., 32./90., 7./90.} },
206 // Strong Stability Preserving
207 { { 0., 0., 0., 0., 0., 0. },
208 // 1st Order
209 { 1., 0., 0., 0., 0., 0. },
210 // 2nd Order - improved
211 { 1./2., 1./2., 0., 0., 0., 0. },
212 // 3rd Order - strong scaling
213 { 1./6., 1./6., 4./6., 0., 0., 0. },
214 // 4th Order - Classic - not used
215 { 1./6., 2./6., 2./6., 1./6., 0., 0. },
216 // 5th Order - 6 stages - not used
217 // Rabiei, Faranak, and Fudziah Ismail. "Fifth order improved
218 // Runge-Kutta method for solving ordinary differential
219 // equations." In Proceedings of the 11th WSEAS international
220 // conference on Applied computer science, pp. 129-133. 2011.
221 { 7./90., 0., 32./90., 12./90., 32./90., 7./90. } } };
222 // clang-format on
223
224 size_t index = (variant == "SSP" || variant == "ImprovedEuler");
225
226 phase->m_schemeType = eExplicit;
227 phase->m_variant = variant;
228 phase->m_order = order;
229 phase->m_name = std::string("RungeKutta") + phase->m_variant +
230 std::string("Order") + std::to_string(phase->m_order);
231
232 phase->m_numsteps = 1;
233 phase->m_numstages = nStages[phase->m_order];
234
235 phase->m_A = Array<OneD, Array<TwoD, NekDouble>>(1);
236 phase->m_B = Array<OneD, Array<TwoD, NekDouble>>(1);
237
238 phase->m_A[0] =
239 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numstages, 0.0);
240 phase->m_B[0] =
241 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numstages, 0.0);
242 phase->m_U =
243 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numsteps, 1.0);
244 phase->m_V =
245 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numsteps, 1.0);
246
247 // Coefficients
248
249 // A Coefficients for each stages along the lower diagonal quadrant.
250 size_t cc = 0;
251
252 for (size_t s = 1; s < phase->m_numstages; ++s)
253 {
254 for (size_t i = 0; i < s; ++i)
255 {
256 phase->m_A[0][s][i] =
257 Acoefficients[index][phase->m_order][cc++];
258 }
259 }
260
261 // B Coefficients for the finial summing.
262 for (size_t n = 0; n < phase->m_numstages; ++n)
263 {
264 phase->m_B[0][0][n] = Bcoefficients[index][phase->m_order][n];
265 }
266
267 phase->m_numMultiStepValues = 1;
268 phase->m_numMultiStepImplicitDerivs = 0;
269 phase->m_numMultiStepExplicitDerivs = 0;
270 phase->m_timeLevelOffset = Array<OneD, size_t>(phase->m_numsteps);
271 phase->m_timeLevelOffset[0] = 0;
272
273 phase->CheckAndVerify();
274 }
@ eExplicit
Formally explicit scheme.
double NekDouble

References Nektar::LibUtilities::eExplicit.

Referenced by Nektar::LibUtilities::AdamsBashforthTimeIntegrationScheme::AdamsBashforthTimeIntegrationScheme(), and RungeKuttaTimeIntegrationScheme().

◆ v_GetName()

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

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 277 of file RungeKuttaTimeIntegrationSchemes.h.

278 {
279 return std::string("RungeKutta");
280 }

◆ v_GetTimeStability()

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

Implements Nektar::LibUtilities::TimeIntegrationScheme.

Definition at line 282 of file RungeKuttaTimeIntegrationSchemes.h.

283 {
284 if (GetOrder() == 1 || GetOrder() == 2)
285 {
286 return 2.0;
287 }
288 else if (GetOrder() == 3)
289 {
290 return 2.51274532661833;
291 }
292 else if (GetOrder() == 4)
293 {
294 // return 2.78529356340528;
295 return 2.784;
296 }
297 else if (GetOrder() == 5)
298 {
299 return 3.21704786664011;
300 }
301 else
302 {
303 return 2.0;
304 }
305 }

References Nektar::LibUtilities::TimeIntegrationScheme::GetOrder().

Member Data Documentation

◆ className

std::string Nektar::LibUtilities::RungeKuttaTimeIntegrationScheme::className
static

Definition at line 97 of file RungeKuttaTimeIntegrationSchemes.h.