Nektar++
TimeIntegrationAlgorithmGLM.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: TimeIntegrationAlgorithmGLM.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 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: Implementation of time integration scheme data class.
33//
34///////////////////////////////////////////////////////////////////////////////
35
37
40
41#include <iostream>
42
43#include <cmath>
44
46{
47
48////////////////////////////////////////////////////////////////////////////////
51{
52 m_op = op;
53}
54
56 const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time)
57{
58 // Create a TimeIntegrationSolutionGLM object based upon the initial value
59 // y_0. Initialise all other multi-step values and derivatives to zero.
62 this, y_0, time, deltaT);
63
65 {
66 // Ensure initial solution is in correct space.
67 m_op.DoProjection(y_out->GetSolution(), y_out->UpdateSolution(), time);
68 }
70 fabs(A(0, 0)) < NekConstants::kNekZeroTol)
71 {
72 // Explicit first-stage when first diagonal coefficient is equal to
73 // zero (EDIRK/ESDIRK schemes).
74 m_op.DoProjection(y_out->GetSolution(), y_out->UpdateSolution(), time);
75 }
76
77 // Calculate the initial derivative, if is part of the solution vector of
78 // the current scheme.
80 {
82
84 {
85 size_t nvar = y_0.size();
86 size_t npoints = y_0[0].size();
87 DoubleArray f_y_0(nvar);
88 for (size_t i = 0; i < nvar; i++)
89 {
90 f_y_0[i] = Array<OneD, NekDouble>(npoints);
91 }
92
93 // Calculate the derivative of the initial value.
94 m_op.DoOdeRhs(y_out->GetSolution(), f_y_0, time);
95
96 // Multiply by the step size.
97 for (size_t i = 0; i < nvar; i++)
98 {
99 Vmath::Smul(npoints, deltaT, f_y_0[i], 1, f_y_0[i], 1);
100 }
101 y_out->SetExplicitDerivative(0, f_y_0, deltaT);
102 }
103 }
104
105 return y_out;
106}
107
109 const NekDouble deltaT, TimeIntegrationSolutionGLMSharedPtr &solvector)
110{
111 size_t nvar = solvector->GetFirstDim();
112 size_t npoints = solvector->GetSecondDim();
113
114 if (solvector->GetIntegrationSchemeData() != this)
115 {
116 // This branch will be taken when the solution vector (solvector) is set
117 // up for a different scheme than the object this method is called from.
118 // (typically needed to calculate the first time-levels of a multi-step
119 // scheme).
120
121 // To do this kind of 'non-matching' integration, we perform the
122 // following three steps:
123 //
124 // 1: Copy the required input information from the solution vector of
125 // the master scheme to the input solution vector of the current
126 // scheme.
127 //
128 // 2: Time-integrate for one step using the current scheme.
129 //
130 // 3: Copy the information contained in the output vector of the current
131 // scheme to the solution vector of the master scheme.
132
133 // STEP 1: Copy the required input information from the solution
134 // vector of the master scheme to the input solution vector
135 // of the current scheme.
136
137 // 1.1 Determine which information is required for the current scheme.
138
139 // Number of required values of the current scheme.
140 size_t nCurSchemeVals = GetNmultiStepValues();
141
142 // Number of required implicit derivs of the current scheme.
143 size_t nCurSchemeImpDers = GetNmultiStepImplicitDerivs();
144
145 // Number of required explicit derivs of the current scheme.
146 size_t nCurSchemeExpDers = GetNmultiStepExplicitDerivs();
147
148 // Number of steps in the current scheme.
149 size_t nCurSchemeSteps = GetNsteps();
150
151 // Number of values of the master scheme.
152 size_t nMasterSchemeVals = solvector->GetNvalues();
153
154 // Number of implicit derivs of the master scheme.
155 size_t nMasterSchemeImpDers = solvector->GetNimplicitderivs();
156
157 // Number of explicit derivs of the master scheme.
158 size_t nMasterSchemeExpDers = solvector->GetNexplicitderivs();
159
160 // Array indicating which time-level the values and derivatives of the
161 // current schemes belong.
162 const Array<OneD, const size_t> &curTimeLevels = GetTimeLevelOffset();
163
164 // Array indicating which time-level the values and derivatives of the
165 // master schemes belong.
166 const Array<OneD, const size_t> &masterTimeLevels =
167 solvector->GetTimeLevelOffset();
168
169 // 1.2 Copy the required information from the master solution vector to
170 // the input solution vector of the current scheme.
171
172 NekDouble t_n = 0.0;
173 DoubleArray y_n;
174 DoubleArray dtFy_n;
175
176 // Input solution vector of the current scheme.
179
180 // Set solution.
181 for (size_t n = 0; n < nCurSchemeVals; n++)
182 {
183 // Get the required value out of the master solution vector.
184 y_n = solvector->GetValue(curTimeLevels[n]);
185 t_n = solvector->GetValueTime(curTimeLevels[n]);
186
187 // Set the required value in the input solution vector of the
188 // current scheme.
189 solvector_in->SetValue(curTimeLevels[n], y_n, t_n);
190 }
191
192 // Set implicit derivatives.
193 for (size_t n = nCurSchemeVals; n < nCurSchemeVals + nCurSchemeImpDers;
194 n++)
195 {
196 // Get the required derivative out of the master solution vector.
197 dtFy_n = solvector->GetImplicitDerivative(curTimeLevels[n]);
198
199 // Set the required derivative in the input solution vector of the
200 // current scheme.
201 solvector_in->SetImplicitDerivative(curTimeLevels[n], dtFy_n,
202 deltaT);
203 }
204
205 // Set explicit derivatives.
206 for (size_t n = nCurSchemeVals + nCurSchemeImpDers; n < nCurSchemeSteps;
207 n++)
208 {
209 // Get the required derivative out of the master solution vector.
210 dtFy_n = solvector->GetExplicitDerivative(curTimeLevels[n]);
211
212 // Set the required derivative in the input solution vector of the
213 // current scheme.
214 solvector_in->SetExplicitDerivative(curTimeLevels[n], dtFy_n,
215 deltaT);
216 }
217
218 // STEP 2: Time-integrate for one step using the current scheme.
219
220 // Output solution vector of the current scheme.
223 this, nvar, npoints);
224
225 // Integrate one step.
226 TimeIntegrate(deltaT, solvector_in->GetSolutionVector(),
227 solvector_in->GetTimeVector(),
228 solvector_out->UpdateSolutionVector(),
229 solvector_out->UpdateTimeVector());
230
231 // STEP 3: Copy the information contained in the output vector of the
232 // current scheme to the solution vector of the master scheme.
233
234 // 3.1 Check whether the current time scheme updates the most recent
235 // derivative that should be updated in the master scheme. If not,
236 // calculate the derivative. This can be done based upon the
237 // corresponding value and the DoOdeRhs operator.
238
239 // This flag indicates whether the new implicit derivative is available
240 // in the output of the current scheme or whether it should be
241 // calculated.
242 bool CalcNewImpDeriv = false;
243
244 if (nMasterSchemeImpDers > 0)
245 {
246 if (nCurSchemeImpDers == 0 || (masterTimeLevels[nMasterSchemeVals] <
247 curTimeLevels[nCurSchemeVals]))
248 {
249 CalcNewImpDeriv = true;
250 }
251 }
252
253 // This flag indicates whether the new explicit derivative is available
254 // in the output of the current scheme or whether it should be
255 // calculated.
256 bool CalcNewExpDeriv = false;
257
258 if (nMasterSchemeExpDers > 0)
259 {
260 if (nCurSchemeExpDers == 0 ||
261 (masterTimeLevels[nMasterSchemeVals + nMasterSchemeImpDers] <
262 curTimeLevels[nCurSchemeVals + nCurSchemeImpDers]))
263 {
264 CalcNewExpDeriv = true;
265 }
266 }
267
268 // Indicates the time level at which the implicit derivative of the
269 // master scheme has to be updated.
270 size_t newImpDerivTimeLevel =
271 (masterTimeLevels.size() > nMasterSchemeVals)
272 ? masterTimeLevels[nMasterSchemeVals]
273 : 0;
274
275 DoubleArray f_impn(nvar);
276
277 // Calculate implicit derivatives.
278 if (CalcNewImpDeriv)
279 {
280 if (newImpDerivTimeLevel == 0 || newImpDerivTimeLevel == 1)
281 {
282 y_n = solvector->GetValue(newImpDerivTimeLevel);
283 t_n = solvector->GetValueTime(newImpDerivTimeLevel);
284 }
285 else
286 {
287 ASSERTL1(false, "Problems with initialising scheme");
288 }
289
290 for (size_t j = 0; j < nvar; j++)
291 {
292 f_impn[j] = Array<OneD, NekDouble>(npoints);
293 }
294
295 // Calculate the derivative of the initial value.
296 m_op.DoImplicitSolve(y_n, f_impn, t_n + deltaT, deltaT);
297 for (size_t j = 0; j < nvar; j++)
298 {
299 Vmath::Vsub(m_npoints, f_impn[j], 1, y_n[j], 1, f_impn[j], 1);
300 }
301 }
302
303 // Indiciates the time level at which the explicit derivative of the
304 // master scheme has to be updated.
305 size_t newExpDerivTimeLevel =
306 (masterTimeLevels.size() > nMasterSchemeVals + nMasterSchemeImpDers)
307 ? masterTimeLevels[nMasterSchemeVals + nMasterSchemeImpDers]
308 : 0;
309
310 DoubleArray f_expn(nvar);
311
312 // Calculate explicit derivatives.
313 if (CalcNewExpDeriv)
314 {
315 // If the time level corresponds to 0, calculate the derivative
316 // based upon the solution value at the new time-level.
317 if (newExpDerivTimeLevel == 0)
318 {
319 y_n = solvector_out->GetValue(0);
320 t_n = solvector_out->GetValueTime(0);
321 }
322 // If the time level corresponds to 1, calculate the derivative
323 // based upon the solution value at the old time-level.
324 else if (newExpDerivTimeLevel == 1)
325 {
326 y_n = solvector->GetValue(0);
327 t_n = solvector->GetValueTime(0);
328 }
329 else
330 {
331 ASSERTL1(false, "Problems with initialising scheme");
332 }
333
334 for (size_t j = 0; j < nvar; j++)
335 {
336 f_expn[j] = Array<OneD, NekDouble>(npoints);
337 }
338
339 // Ensure solution is in correct space.
340 if (newExpDerivTimeLevel == 1)
341 {
342 m_op.DoProjection(y_n, y_n, t_n);
343 }
344
345 // Calculate the derivative.
346 m_op.DoOdeRhs(y_n, f_expn, t_n);
347
348 // Multiply by dt (as required by the General Linear Method
349 // framework).
350 for (size_t j = 0; j < nvar; j++)
351 {
352 Vmath::Smul(npoints, deltaT, f_expn[j], 1, f_expn[j], 1);
353 }
354 }
355
356 // Rotate the solution vector (i.e. updating without
357 // calculating/inserting new values).
358 solvector->RotateSolutionVector();
359
360 // 3.2 Copy the information calculated using the current scheme from
361 // the output solution vector to the master solution vector.
362
363 // Set solution.
364 y_n = solvector_out->GetValue(0);
365 t_n = solvector_out->GetValueTime(0);
366
367 // Ensure solution is in correct space.
368 if (newExpDerivTimeLevel == 1)
369 {
370 m_op.DoProjection(y_n, y_n, t_n);
371 }
372
373 solvector->SetValue(0, y_n, t_n);
374
375 // Set implicit derivative.
376 if (CalcNewImpDeriv)
377 {
378 // Set the calculated derivative in the master solution vector.
379 solvector->SetImplicitDerivative(newImpDerivTimeLevel, f_impn,
380 deltaT);
381 }
382 else if (nCurSchemeImpDers > 0 && nMasterSchemeImpDers > 0)
383 {
384 // Get the calculated derivative out of the output solution vector
385 // of the current scheme.
386 dtFy_n = solvector_out->GetImplicitDerivative(newImpDerivTimeLevel);
387
388 // Set the calculated derivative in the master solution vector.
389 solvector->SetImplicitDerivative(newImpDerivTimeLevel, dtFy_n,
390 deltaT);
391 }
392
393 // Set explicit derivative.
394 if (CalcNewExpDeriv)
395 {
396 // Set the calculated derivative in the master solution vector.
397 solvector->SetExplicitDerivative(newExpDerivTimeLevel, f_expn,
398 deltaT);
399 }
400 else if (nCurSchemeExpDers > 0 && nMasterSchemeExpDers > 0)
401 {
402 // Get the calculated derivative out of the output solution vector
403 // of the current scheme.
404 dtFy_n = solvector_out->GetExplicitDerivative(newExpDerivTimeLevel);
405
406 // Set the calculated derivative in the master solution vector.
407 solvector->SetExplicitDerivative(newExpDerivTimeLevel, dtFy_n,
408 deltaT);
409 }
410 }
411 else
412 {
415 this, nvar, npoints);
416
417 TimeIntegrate(deltaT, solvector->GetSolutionVector(),
418 solvector->GetTimeVector(),
419 solvector_new->UpdateSolutionVector(),
420 solvector_new->UpdateTimeVector());
421
422 solvector = solvector_new;
423 }
424
425 return solvector->GetSolution();
426
427} // end TimeIntegrate()
428
429// Does the actual multi-stage multi-step integration.
431 ConstTripleArray &y_old,
432 ConstSingleArray &t_old,
433 TripleArray &y_new,
434 SingleArray &t_new)
435{
436 ASSERTL1(CheckTimeIntegrateArguments(y_old, t_old, y_new, t_new),
437 "Arguments not well defined");
438
440
441 // Check if storage has already been initialised. If so, we just zero the
442 // temporary storage.
443 if (m_initialised && m_nvars == GetFirstDim(y_old) &&
444 m_npoints == GetSecondDim(y_old))
445 {
446 for (size_t j = 0; j < m_nvars; j++)
447 {
449 }
450 }
451 else
452 {
453 m_nvars = GetFirstDim(y_old);
454 m_npoints = GetSecondDim(y_old);
455
456 // First, calculate the various stage values and stage derivatives
457 // (this is the multi-stage part of the method)
458 // - m_Y corresponds to the stage values
459 // - m_F corresponds to the stage derivatives
460 // - m_T corresponds to the time at the different stages
461 // - m_tmp corresponds to the explicit right hand side of each stage
462 // equation
463 // (for explicit schemes, this correspond to m_Y)
464
465 // Allocate memory for the arrays m_Y and m_F and m_tmp. The same
466 // storage will be used for every stage -> m_Y is a DoubleArray.
468 for (size_t j = 0; j < m_nvars; j++)
469 {
471 }
472
473 // The same storage will be used for every stage -> m_tmp is a
474 // DoubleArray.
475 if (type == eExplicit || type == eExponential)
476 {
477 m_Y = m_tmp;
478 }
479 else
480 {
482 for (size_t j = 0; j < m_nvars; j++)
483 {
485 }
486 }
487
488 // Different storage for every stage derivative as the data
489 // will be re-used to update the solution -> m_F is a TripleArray.
491 for (size_t i = 0; i < m_numstages; ++i)
492 {
493 m_F[i] = DoubleArray(m_nvars);
494 for (size_t j = 0; j < m_nvars; j++)
495 {
497 }
498 }
499
500 if (type == eIMEX)
501 {
503 for (size_t i = 0; i < m_numstages; ++i)
504 {
506 for (size_t j = 0; j < m_nvars; j++)
507 {
509 }
510 }
511 }
512
513 // Finally, flag indicating that the memory has been initialised.
514 m_initialised = true;
515
516 } // end else
517
518 // For an exponential integrator if the time increment or the number of
519 // variables has changed then the exponenial matrices must be recomputed.
520
521 // NOTE: The eigenvalues are not yet seeded!
522 if (type == eExponential)
523 {
524 if (m_lastDeltaT != deltaT || m_lastNVars != GetFirstDim(y_old))
525 {
526 m_parent->InitializeSecondaryData(this, deltaT);
527
528 m_lastDeltaT = deltaT;
529 m_lastNVars = GetFirstDim(y_old);
530 }
531 }
532
533 // The loop below calculates the stage values and derivatives.
534 for (size_t stage = 0; stage < m_numstages; stage++)
535 {
536 if ((stage == 0) && m_firstStageEqualsOldSolution)
537 {
538 for (size_t k = 0; k < m_nvars; k++)
539 {
540 Vmath::Vcopy(m_npoints, y_old[0][k], 1, m_Y[k], 1);
541 }
542 }
543 else
544 {
545 // The stage values m_Y are a linear combination of:
546 // 1: The stage derivatives:
547 if (stage != 0)
548 {
549 for (size_t k = 0; k < m_nvars; k++)
550 {
551 Vmath::Smul(m_npoints, deltaT * A(k, stage, 0), m_F[0][k],
552 1, m_tmp[k], 1);
553
554 if (type == eIMEX)
555 {
556 Vmath::Svtvp(m_npoints, deltaT * A_IMEX(stage, 0),
557 m_F_IMEX[0][k], 1, m_tmp[k], 1, m_tmp[k],
558 1);
559 }
560 }
561 }
562
563 for (size_t j = 1; j < stage; j++)
564 {
565 for (size_t k = 0; k < m_nvars; k++)
566 {
567 Vmath::Svtvp(m_npoints, deltaT * A(k, stage, j), m_F[j][k],
568 1, m_tmp[k], 1, m_tmp[k], 1);
569
570 if (type == eIMEX)
571 {
572 Vmath::Svtvp(m_npoints, deltaT * A_IMEX(stage, j),
573 m_F_IMEX[j][k], 1, m_tmp[k], 1, m_tmp[k],
574 1);
575 }
576 }
577 }
578
579 // 2: The imported multi-step solution of the previous time level:
580 for (size_t j = 0; j < m_numsteps; j++)
581 {
582 for (size_t k = 0; k < m_nvars; k++)
583 {
584 Vmath::Svtvp(m_npoints, U(k, stage, j), y_old[j][k], 1,
585 m_tmp[k], 1, m_tmp[k], 1);
586 }
587 }
588 } // end else
589
590 // Calculate the stage derivative based upon the stage value.
591 if (type == eDiagonallyImplicit)
592 {
593 if (m_numstages == 1)
594 {
595 m_T = t_old[0] + deltaT;
596 }
597 else
598 {
599 m_T = t_old[0];
600 for (size_t j = 0; j <= stage; ++j)
601 {
602 m_T += A(stage, j) * deltaT;
603 }
604 }
605
606 // Explicit first-stage when first diagonal coefficient is equal to
607 // zero (EDIRK/ESDIRK schemes).
608 if ((stage == 0) && (fabs(A(0, 0)) < NekConstants::kNekZeroTol))
609 {
610 m_op.DoOdeRhs(m_Y, m_F[stage], t_old[0]);
611 }
612 // Otherwise, the stage is updated implicitly.
613 else
614 {
615 m_op.DoImplicitSolve(m_tmp, m_Y, m_T, A(stage, stage) * deltaT);
616
617 for (size_t k = 0; k < m_nvars; ++k)
618 {
619 Vmath::Vsub(m_npoints, m_Y[k], 1, m_tmp[k], 1,
620 m_F[stage][k], 1);
621 Vmath::Smul(m_npoints, 1.0 / (A(stage, stage) * deltaT),
622 m_F[stage][k], 1, m_F[stage][k], 1);
623 }
624 }
625 }
626 else if (type == eIMEX)
627 {
628 if (m_numstages == 1)
629 {
630 m_T = t_old[0] + deltaT;
631 }
632 else
633 {
634 m_T = t_old[0];
635 for (size_t j = 0; j <= stage; ++j)
636 {
637 m_T += A(stage, j) * deltaT;
638 }
639 }
640
641 if (fabs(A(stage, stage)) > NekConstants::kNekZeroTol)
642 {
643 m_op.DoImplicitSolve(m_tmp, m_Y, m_T, A(stage, stage) * deltaT);
644
645 for (size_t k = 0; k < m_nvars; k++)
646 {
647 Vmath::Vsub(m_npoints, m_Y[k], 1, m_tmp[k], 1,
648 m_F[stage][k], 1);
649 Vmath::Smul(m_npoints, 1.0 / (A(stage, stage) * deltaT),
650 m_F[stage][k], 1, m_F[stage][k], 1);
651 }
652 }
653
654 m_op.DoOdeRhs(m_Y, m_F_IMEX[stage], m_T);
655 }
656 else if (type == eExplicit || type == eExponential)
657 {
658 if (m_numstages == 1)
659 {
660 m_T = t_old[0];
661 }
662 else
663 {
664 m_T = t_old[0];
665 for (size_t j = 0; j < stage; ++j)
666 {
667 m_T += A(stage, j) * deltaT;
668 }
669 }
670
671 // Avoid projecting the same solution twice.
672 if (!((stage == 0) && m_firstStageEqualsOldSolution))
673 {
674 // Ensure solution is in correct space.
676 }
677
678 m_op.DoOdeRhs(m_Y, m_F[stage], m_T);
679 }
680 }
681
682 // Next, the solution vector y at the new time level will be calculated.
683 //
684 // For multi-step methods, this includes updating the values of the
685 // auxiliary parameters.
686 //
687 // The loop below calculates the solution at the new time level.
688 //
689 // If last stage equals the new solution, the new solution needs not be
690 // calculated explicitly but can simply be copied. This saves a solve.
691 size_t i_start = 0;
693 {
694 for (size_t k = 0; k < m_nvars; k++)
695 {
696 Vmath::Vcopy(m_npoints, m_Y[k], 1, y_new[0][k], 1);
697 }
698
699 t_new[0] = t_old[0] + deltaT;
700
701 i_start = 1;
702 }
703
704 for (size_t i = i_start; i < m_numsteps; i++)
705 {
706 // The solution at the new time level is a linear combination of:
707 // 1: The stage derivatives:
708 for (size_t k = 0; k < m_nvars; k++)
709 {
710 Vmath::Smul(m_npoints, deltaT * B(k, i, 0), m_F[0][k], 1,
711 y_new[i][k], 1);
712
713 if (type == eIMEX)
714 {
715 Vmath::Svtvp(m_npoints, deltaT * B_IMEX(i, 0), m_F_IMEX[0][k],
716 1, y_new[i][k], 1, y_new[i][k], 1);
717 }
718 }
719
720 if (m_numstages != 1 || type != eIMEX)
721 {
722 t_new[i] = B(i, 0) * deltaT;
723 }
724
725 for (size_t j = 1; j < m_numstages; j++)
726 {
727 for (size_t k = 0; k < m_nvars; k++)
728 {
729 Vmath::Svtvp(m_npoints, deltaT * B(k, i, j), m_F[j][k], 1,
730 y_new[i][k], 1, y_new[i][k], 1);
731
732 if (type == eIMEX)
733 {
734 Vmath::Svtvp(m_npoints, deltaT * B_IMEX(i, j),
735 m_F_IMEX[j][k], 1, y_new[i][k], 1, y_new[i][k],
736 1);
737 }
738 }
739
740 if (m_numstages != 1 || type != eIMEX)
741 {
742 t_new[i] += B(i, j) * deltaT;
743 }
744 }
745
746 // 2: The imported multi-step solution of the previous time level:
747 for (size_t j = 0; j < m_numsteps; j++)
748 {
749 for (size_t k = 0; k < m_nvars; k++)
750 {
751 Vmath::Svtvp(m_npoints, V(k, i, j), y_old[j][k], 1, y_new[i][k],
752 1, y_new[i][k], 1);
753 }
754
755 if (m_numstages != 1 || type != eIMEX)
756 {
757 t_new[i] += V(i, j) * t_old[j];
758 }
759 }
760 }
761
762 // Ensure that the new solution is projected if necessary.
763 if (type == eExplicit || type == eExponential ||
764 fabs(m_T - t_new[0]) > NekConstants::kNekZeroTol)
765 {
766 m_op.DoProjection(y_new[0], y_new[0], t_new[0]);
767 }
768
769} // end TimeIntegrate()
770
772{
773 // First stage equals old solution if:
774 // 1. the first row of the coefficient matrix A consists of zeros
775 // 2. U[0][0] is equal to one and all other first row entries of U are zero
776
777 // 1. Check first condition
778 for (size_t m = 0; m < m_A.size(); m++)
779 {
780 for (size_t i = 0; i < m_numstages; i++)
781 {
782 if (fabs(m_A[m][0][i]) > NekConstants::kNekZeroTol)
783 {
785 return;
786 }
787 }
788 }
789
790 // 2. Check second condition
791 if (fabs(m_U[0][0] - 1.0) > NekConstants::kNekZeroTol)
792 {
794 return;
795 }
796
797 for (size_t i = 1; i < m_numsteps; i++)
798 {
799 if (fabs(m_U[0][i]) > NekConstants::kNekZeroTol)
800 {
802 return;
803 }
804 }
805
807}
808
810{
811 // Last stage equals new solution if:
812 // 1. the last row of the coefficient matrix A is equal to the first row of
813 // matrix B
814 // 2. the last row of the coefficient matrix U is equal to the first row of
815 // matrix V
816
817 // 1. Check first condition
818 for (size_t m = 0; m < m_A.size(); m++)
819 {
820 for (size_t i = 0; i < m_numstages; i++)
821 {
822 if (fabs(m_A[m][m_numstages - 1][i] - m_B[m][0][i]) >
824 {
826 return;
827 }
828 }
829 }
830
831 // 2. Check second condition
832 for (size_t i = 0; i < m_numsteps; i++)
833 {
834 if (fabs(m_U[m_numstages - 1][i] - m_V[0][i]) >
836 {
838 return;
839 }
840 }
841
843}
844
846{
847#ifdef NEKTAR_DEBUG
848 size_t IMEXdim = m_A.size();
849 size_t dim = m_A[0].GetRows();
850
852
854 {
855 vertype[0] = eExponential;
856 }
857
858 for (size_t m = 0; m < IMEXdim; m++)
859 {
860 for (size_t i = 0; i < dim; i++)
861 {
862 if (fabs(m_A[m][i][i]) > NekConstants::kNekZeroTol)
863 {
864 vertype[m] = eDiagonallyImplicit;
865 }
866 }
867
868 for (size_t i = 0; i < dim; i++)
869 {
870 for (size_t j = i + 1; j < dim; j++)
871 {
872 if (fabs(m_A[m][i][j]) > NekConstants::kNekZeroTol)
873 {
874 vertype[m] = eImplicit;
875 ASSERTL1(false, "Fully Implicit schemes cannnot be handled "
876 "by the TimeIntegrationScheme class");
877 }
878 }
879 }
880 }
881
882 if (IMEXdim == 2)
883 {
884 ASSERTL1(m_B.size() == 2,
885 "Coefficient Matrix B should have an "
886 "implicit and explicit part for IMEX schemes");
887
888 if ((vertype[0] == eDiagonallyImplicit) && (vertype[1] == eExplicit))
889 {
890 vertype[0] = eIMEX;
891 }
892 else
893 {
894 ASSERTL1(false, "This is not a proper IMEX scheme");
895 }
896 }
897
898 ASSERTL1(vertype[0] == m_schemeType,
899 "Time integration scheme coefficients do not match its type");
900#endif
901}
902
904 [[maybe_unused]] ConstTripleArray &y_old,
905 [[maybe_unused]] ConstSingleArray &t_old,
906 [[maybe_unused]] TripleArray &y_new,
907 [[maybe_unused]] SingleArray &t_new) const
908{
909 // Check if arrays are all of consistent size
910
911 ASSERTL1(y_old.size() == m_numsteps, "Non-matching number of steps.");
912 ASSERTL1(y_new.size() == m_numsteps, "Non-matching number of steps.");
913
914 ASSERTL1(y_old[0].size() == y_new[0].size(),
915 "Non-matching number of variables.");
916 ASSERTL1(y_old[0][0].size() == y_new[0][0].size(),
917 "Non-matching number of coefficients.");
918
919 ASSERTL1(t_old.size() == m_numsteps, "Non-matching number of steps.");
920 ASSERTL1(t_new.size() == m_numsteps, "Non-matching number of steps.");
921
922 return true;
923}
924
925// Friend Operators
926std::ostream &operator<<(std::ostream &os,
928{
929 return operator<<(os, *rhs);
930}
931
932std::ostream &operator<<(std::ostream &os,
934{
935 size_t r = rhs.m_numsteps;
936 size_t s = rhs.m_numstages;
938
939 size_t oswidth = 9;
940 size_t osprecision = 6;
941
942 os << "Time Integration Scheme (Master): " << rhs.m_parent->GetFullName()
943 << "\n"
944 << "Time Integration Phase : " << rhs.m_name << "\n"
945 << "- number of steps: " << r << "\n"
946 << "- number of stages: " << s << "\n"
947 << "General linear method tableau:\n";
948
949 for (size_t i = 0; i < s; i++)
950 {
951 for (size_t j = 0; j < s; j++)
952 {
953 os.width(oswidth);
954 os.precision(osprecision);
955 os << std::right << rhs.A(i, j) << " ";
956 }
957 if (type == eIMEX)
958 {
959 os << " '";
960 for (size_t j = 0; j < s; j++)
961 {
962 os.width(oswidth);
963 os.precision(osprecision);
964 os << std::right << rhs.A_IMEX(i, j) << " ";
965 }
966 }
967
968 os << " |";
969
970 for (size_t j = 0; j < r; j++)
971 {
972 os.width(oswidth);
973 os.precision(osprecision);
974 os << std::right << rhs.U(i, j);
975 }
976 os << std::endl;
977 }
978
979 size_t imexflag = (type == eIMEX) ? 2 : 1;
980 for (size_t i = 0;
981 i < (r + imexflag * s) * (oswidth + 1) + imexflag * 2 - 1; i++)
982 {
983 os << "-";
984 }
985 os << std::endl;
986
987 for (size_t i = 0; i < r; i++)
988 {
989 for (size_t j = 0; j < s; j++)
990 {
991 os.width(oswidth);
992 os.precision(osprecision);
993 os << std::right << rhs.B(i, j) << " ";
994 }
995 if (type == eIMEX)
996 {
997 os << " '";
998 for (size_t j = 0; j < s; j++)
999 {
1000 os.width(oswidth);
1001 os.precision(osprecision);
1002 os << std::right << rhs.B_IMEX(i, j) << " ";
1003 }
1004 }
1005
1006 os << " |";
1007
1008 for (size_t j = 0; j < r; j++)
1009 {
1010 os.width(oswidth);
1011 os.precision(osprecision);
1012 os << std::right << rhs.V(i, j);
1013 }
1014
1015 os << " |";
1016
1017 os.width(oswidth);
1018 os.precision(osprecision);
1019 os << std::right << rhs.m_timeLevelOffset[i];
1020
1021 if (i < rhs.m_numMultiStepValues)
1022 {
1023 os << std::right << " value";
1024 }
1025 else
1026 {
1027 os << std::right << " derivative";
1028 }
1029
1030 os << std::endl;
1031 }
1032
1033 if (type == eExponential)
1034 {
1035 for (size_t k = 0; k < rhs.m_nvars; k++)
1036 {
1037 os << std::endl
1038 << "General linear method exponential tableau for variable " << k
1039 << ":\n";
1040
1041 for (size_t i = 0; i < s; i++)
1042 {
1043 for (size_t j = 0; j < s; j++)
1044 {
1045 os.width(oswidth);
1046 os.precision(osprecision);
1047 os << std::right << rhs.A(k, i, j) << " ";
1048 }
1049
1050 os << " |";
1051
1052 for (size_t j = 0; j < r; j++)
1053 {
1054 os.width(oswidth);
1055 os.precision(osprecision);
1056 os << std::right << rhs.B(k, i, j);
1057 }
1058 os << std::endl;
1059 }
1060
1061 for (size_t i = 0;
1062 i < (r + imexflag * s) * (oswidth + 1) + imexflag * 2 - 1; i++)
1063 {
1064 os << "-";
1065 }
1066 os << std::endl;
1067
1068 for (size_t i = 0; i < r; i++)
1069 {
1070 for (size_t j = 0; j < s; j++)
1071 {
1072 os.width(oswidth);
1073 os.precision(osprecision);
1074 os << std::right << rhs.B(k, i, j) << " ";
1075 }
1076
1077 os << " |";
1078
1079 for (size_t j = 0; j < r; j++)
1080 {
1081 os.width(oswidth);
1082 os.precision(osprecision);
1083 os << std::right << rhs.V(k, i, j);
1084 }
1085 os << std::endl;
1086 }
1087 }
1088 }
1089
1090 return os;
1091} // end function operator<<
1092
1093} // namespace Nektar::LibUtilities
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
NekDouble V(const size_t i, const size_t j) const
const TimeIntegrationSchemeGLM * m_parent
Parent scheme object.
TripleArray m_F
Explicit RHS of each stage equation.
DoubleArray m_tmp
Array containing the stage values.
NekDouble A(const size_t i, const size_t j) const
NekDouble B_IMEX(const size_t i, const size_t j) const
bool CheckTimeIntegrateArguments(ConstTripleArray &y_old, ConstSingleArray &t_old, TripleArray &y_new, SingleArray &t_new) const
NekDouble A_IMEX(const size_t i, const size_t j) const
NekDouble B(const size_t i, const size_t j) const
NekDouble U(const size_t i, const size_t j) const
NekDouble m_T
Used to store the Explicit stage derivative of IMEX schemes.
const Array< OneD, const size_t > & GetTimeLevelOffset() const
size_t m_npoints
The number of variables in integration scheme.
void InitializeScheme(const TimeIntegrationSchemeOperators &op)
LUE ConstDoubleArray & TimeIntegrate(const NekDouble deltaT, TimeIntegrationSolutionGLMSharedPtr &y)
Explicit integration of an ODE.
LUE TimeIntegrationSolutionGLMSharedPtr InitializeData(const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time)
This function initialises the time integration scheme.
TripleArray m_F_IMEX
Array corresponding to the stage Derivatives.
LUE void InitializeSecondaryData(TimeIntegrationAlgorithmGLM *phase, NekDouble deltaT) const
Binds a set of functions for use by time integration schemes.
void DoImplicitSolve(InArrayType &inarray, OutArrayType &outarray, const NekDouble time, const NekDouble lambda) const
void DoProjection(InArrayType &inarray, OutArrayType &outarray, const NekDouble time) const
void DoOdeRhs(InArrayType &inarray, OutArrayType &outarray, const NekDouble time) const
static std::shared_ptr< DataType > AllocateSharedPtr(const Args &...args)
Allocate a shared pointer from the memory pool.
AT< AT< NekDouble > > DoubleArray
std::shared_ptr< TimeIntegrationAlgorithmGLM > TimeIntegrationAlgorithmGLMSharedPtr
@ eImplicit
Fully implicit scheme.
@ eExplicit
Formally explicit scheme.
@ eExponential
Exponential scheme.
@ eDiagonallyImplicit
Diagonally implicit scheme (e.g. the DIRK schemes)
@ eIMEX
Implicit Explicit General Linear Method.
std::shared_ptr< TimeIntegrationSolutionGLM > TimeIntegrationSolutionGLMSharedPtr
AT< AT< AT< NekDouble > > > TripleArray
std::ostream & operator<<(std::ostream &os, const BasisKey &rhs)
static const NekDouble kNekZeroTol
double NekDouble
void Svtvp(int n, const T alpha, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Svtvp (scalar times vector plus vector): z = alpha*x + y.
Definition: Vmath.hpp:396
void Smul(int n, const T alpha, const T *x, const int incx, T *y, const int incy)
Scalar multiply y = alpha*x.
Definition: Vmath.hpp:100
void Zero(int n, T *x, const int incx)
Zero vector.
Definition: Vmath.hpp:273
void Vcopy(int n, const T *x, const int incx, T *y, const int incy)
Definition: Vmath.hpp:825
void Vsub(int n, const T *x, const int incx, const T *y, const int incy, T *z, const int incz)
Subtract vector z = x-y.
Definition: Vmath.hpp:220