Nektar++
Loading...
Searching...
No Matches
RungeKuttaTimeIntegrationSchemes.h
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: RungeKuttaTimeIntegrationSchemes.h
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2018 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// License for the specific language governing rights and limitations under
14// Permission is hereby granted, free of charge, to any person obtaining a
15// copy of this software and associated documentation files (the "Software"),
16// to deal in the Software without restriction, including without limitation
17// the rights to use, copy, modify, merge, publish, distribute, sublicense,
18// and/or sell copies of the Software, and to permit persons to whom the
19// Software is furnished to do so, subject to the following conditions:
20//
21// The above copyright notice and this permission notice shall be included
22// in all copies or substantial portions of the Software.
23//
24// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
25// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
26// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
27// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
28// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
29// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
30// DEALINGS IN THE SOFTWARE.
31//
32// Description: Combined header file for all Runge Kutta based time integration
33// schemes.
34//
35///////////////////////////////////////////////////////////////////////////////
36
37// Note : If adding a new integrator be sure to register the
38// integrator with the Time Integration Scheme Facatory in
39// SchemeInitializor.cpp.
40
41#ifndef NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_RK_TIME_INTEGRATION_SCHEME
42#define NEKTAR_LIB_UTILITIES_TIME_INTEGRATION_RK_TIME_INTEGRATION_SCHEME
43
44#define LUE LIB_UTILITIES_EXPORT
45
47
49{
50
51////////////////////////////////////////////////////////////////////////////////
52// Runge Kutta Order N where the number of stages == order
53
55{
56public:
57 RungeKuttaTimeIntegrationScheme(std::string variant, size_t order,
58 std::vector<NekDouble> freeParams)
59 : TimeIntegrationSchemeGLM(variant, order, freeParams)
60 {
61 ASSERTL1(variant == "" || variant == "SSP",
62 "Runge Kutta Time integration scheme unknown variant: " +
63 variant + ". Must be blank or 'SSP'");
64
65 // Std - Currently up to 5th order is implemented.
66 // SSP - Currently 1st through 3rd order is implemented.
67 ASSERTL1((variant == "" && 1 <= order && order <= 5) ||
68 (variant == "SSP" && 1 <= order && order <= 3),
69 "Runge Kutta Time integration scheme bad order, "
70 "Std (1-5) or SSP (1-3): " +
71 std::to_string(order));
72
76
78 m_integration_phases[0], variant, order, freeParams);
79 }
80
82
84 std::string variant, size_t order, std::vector<NekDouble> freeParams)
85 {
88 variant, order, freeParams);
89
90 return p;
91 }
92
93 static std::string className;
94
95 LUE static void SetupSchemeData(
96 TimeIntegrationAlgorithmGLMSharedPtr &phase, std::string variant,
97 size_t order, [[maybe_unused]] std::vector<NekDouble> freeParams)
98 {
99 constexpr size_t nStages[6] = {0, 1, 2, 3, 4, 6};
100
101 // A Coefficients for the lower diagonal quadrant stored in a
102 // contiguous fashion. For the fourth order, six coefficients
103 // from the Butcher tableau would be stored as the following.
104 //
105 // 0 0 0 0
106 // Butcher a 0 0 0 Stored as a
107 // Tableau b c 0 0 b c
108 // d e f 0 d e f 0 ... 0
109
110 // clang-format off
111 constexpr NekDouble Acoefficients[2][6][15] =
112 { { { 0., 0., 0., 0., 0.,
113 0., 0., 0., 0., 0.,
114 0., 0., 0., 0., 0. },
115 // 1st Order
116 { 0., 0., 0., 0., 0.,
117 0., 0., 0., 0., 0.,
118 0., 0., 0., 0., 0. },
119 // 2nd Order - midpoint
120 { 1./2, // Last entry
121 0., 0., 0., 0.,
122 0., 0., 0., 0., 0.,
123 0., 0., 0., 0., 0. },
124 // 3rd Order - Ralston's
125 { 1./2.,
126 0., 3./4., // Last entry
127 0., 0.,
128 0., 0., 0., 0., 0.,
129 0., 0., 0., 0., 0. },
130 // 4th Order - Classic
131 { 1./2.,
132 0., 1./2.,
133 0., 0., 1., // Last entry
134 0., 0., 0., 0.,
135 0., 0., 0., 0., 0. },
136 // 5th Order - 6 stages
137 // Rabiei, Faranak, and Fudziah Ismail. "Fifth order improved
138 // Runge-Kutta method for solving ordinary differential
139 // equations." In Proceedings of the 11th WSEAS international
140 // conference on Applied computer science, pp. 129-133. 2011.
141 { 1./4.,
142 1./8., 1./8.,
143 0., -1./2., 1.,
144 3./16., 0., 0., 9./16.,
145 -3./7., 2./7., 12./7., -12./7., 8./7. } },
146 // Strong Stability Preserving
147 { { 0., 0., 0., 0., 0.,
148 0., 0., 0., 0., 0.,
149 0., 0., 0., 0., 0. },
150 // 1st Order
151 { 0., 0., 0., 0., 0.,
152 0., 0., 0., 0., 0.,
153 0., 0., 0., 0., 0. },
154 // 2nd Order - strong scaling - improved
155 { 1., // Last entry
156 0., 0., 0., 0.,
157 0., 0., 0., 0., 0.,
158 0., 0., 0., 0., 0. },
159 // 3rd Order - strong scaling
160 { 1.,
161 1./4., 1./4., // Last entry
162 0, 0.,
163 0., 0., 0., 0., 0.,
164 0., 0., 0., 0., 0. },
165 // 4th Order - Classic - not used
166 { 1./2.,
167 0., 1./2.,
168 0., 0., 1., // Last entry
169 0., 0., 0., 0.,
170 0., 0., 0., 0., 0. },
171 // 5th Order - 6 stages - not used
172 // Rabiei, Faranak, and Fudziah Ismail. "Fifth order improved
173 // Runge-Kutta method for solving ordinary differential
174 // equations." In Proceedings of the 11th WSEAS international
175 // conference on Applied computer science, pp. 129-133. 2011.
176 { 1./4.,
177 1./8., 1./8.,
178 0., -1./2., 1.,
179 3./16., 0., 0., 9./16.,
180 -3./7., 2./7., 12./7., -12./7., 8./7. } } };
181 // clang-format on
182
183 // B Coefficients for the final summing.
184
185 // clang-format off
186 constexpr NekDouble Bcoefficients[2][6][6] =
187 { { { 0., 0., 0., 0., 0., 0. },
188 // 1st Order
189 { 1., 0., 0., 0., 0., 0. },
190 // 2nd Order - midpoint
191 { 0., 1., 0., 0., 0., 0. },
192 // 3rd Order - Ralston's
193 { 2./9., 3./9., 4./9., 0., 0., 0. },
194 // 4th Order - Classic
195 { 1./6., 2./6., 2./6., 1./6., 0., 0. },
196 // 5th Order - 6 stages
197 // Rabiei, Faranak, and Fudziah Ismail. "Fifth order improved
198 // Runge-Kutta method for solving ordinary differential
199 // equations." In Proceedings of the 11th WSEAS international
200 // conference on Applied computer science, pp. 129-133. 2011.
201 { 7./90., 0., 32./90., 12./90., 32./90., 7./90.} },
202 // Strong Stability Preserving
203 { { 0., 0., 0., 0., 0., 0. },
204 // 1st Order
205 { 1., 0., 0., 0., 0., 0. },
206 // 2nd Order - improved
207 { 1./2., 1./2., 0., 0., 0., 0. },
208 // 3rd Order - strong scaling
209 { 1./6., 1./6., 4./6., 0., 0., 0. },
210 // 4th Order - Classic - not used
211 { 1./6., 2./6., 2./6., 1./6., 0., 0. },
212 // 5th Order - 6 stages - not used
213 // Rabiei, Faranak, and Fudziah Ismail. "Fifth order improved
214 // Runge-Kutta method for solving ordinary differential
215 // equations." In Proceedings of the 11th WSEAS international
216 // conference on Applied computer science, pp. 129-133. 2011.
217 { 7./90., 0., 32./90., 12./90., 32./90., 7./90. } } };
218 // clang-format on
219
220 size_t index = (variant == "SSP" || variant == "ImprovedEuler");
221
222 phase->m_schemeType = eExplicit;
223 phase->m_variant = variant;
224 phase->m_order = order;
225 phase->m_name = std::string("RungeKutta") + phase->m_variant +
226 std::string("Order") + std::to_string(phase->m_order);
227
228 phase->m_numsteps = 1;
229 phase->m_numstages = nStages[phase->m_order];
230
231 phase->m_A = Array<OneD, Array<TwoD, NekDouble>>(1);
232 phase->m_B = Array<OneD, Array<TwoD, NekDouble>>(1);
233
234 phase->m_A[0] =
235 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numstages, 0.0);
236 phase->m_B[0] =
237 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numstages, 0.0);
238 phase->m_U =
239 Array<TwoD, NekDouble>(phase->m_numstages, phase->m_numsteps, 1.0);
240 phase->m_V =
241 Array<TwoD, NekDouble>(phase->m_numsteps, phase->m_numsteps, 1.0);
242
243 // Coefficients
244
245 // A Coefficients for each stages along the lower diagonal quadrant.
246 size_t cc = 0;
247
248 for (size_t s = 1; s < phase->m_numstages; ++s)
249 {
250 for (size_t i = 0; i < s; ++i)
251 {
252 phase->m_A[0][s][i] =
253 Acoefficients[index][phase->m_order][cc++];
254 }
255 }
256
257 // B Coefficients for the finial summing.
258 for (size_t n = 0; n < phase->m_numstages; ++n)
259 {
260 phase->m_B[0][0][n] = Bcoefficients[index][phase->m_order][n];
261 }
262
263 phase->m_numMultiStepValues = 1;
264 phase->m_numMultiStepImplicitDerivs = 0;
265 phase->m_numMultiStepExplicitDerivs = 0;
266 phase->m_timeLevelOffset = Array<OneD, size_t>(phase->m_numsteps);
267 phase->m_timeLevelOffset[0] = 0;
268
269 phase->CheckAndVerify();
270 }
271
272protected:
273 LUE std::string v_GetName() const override
274 {
275 return std::string("RungeKutta");
276 }
277
279 {
280 if (GetOrder() == 1 || GetOrder() == 2)
281 {
282 return 2.0;
283 }
284 else if (GetOrder() == 3)
285 {
286 return 2.51274532661833;
287 }
288 else if (GetOrder() == 4)
289 {
290 // return 2.78529356340528;
291 return 2.784;
292 }
293 else if (GetOrder() == 5)
294 {
295 return 3.21704786664011;
296 }
297 else
298 {
299 return 2.0;
300 }
301 }
302
303}; // end class RungeKuttaTimeIntegrator
304
305} // namespace Nektar::LibUtilities
306
307#endif
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
static LUE void SetupSchemeData(TimeIntegrationAlgorithmGLMSharedPtr &phase, std::string variant, size_t order, std::vector< NekDouble > freeParams)
static TimeIntegrationSchemeSharedPtr create(std::string variant, size_t order, std::vector< NekDouble > freeParams)
RungeKuttaTimeIntegrationScheme(std::string variant, size_t order, std::vector< NekDouble > freeParams)
Base class for GLM time integration schemes.
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
std::shared_ptr< TimeIntegrationAlgorithmGLM > TimeIntegrationAlgorithmGLMSharedPtr
@ eExplicit
Formally explicit scheme.
std::vector< TimeIntegrationAlgorithmGLMSharedPtr > TimeIntegrationAlgorithmGLMVector
std::shared_ptr< TimeIntegrationScheme > TimeIntegrationSchemeSharedPtr