Nektar++
TimeIntegrationSchemeFIT.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: TimeIntegrationSchemeFIT.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// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: implementation of time integration scheme FIT class
32//
33///////////////////////////////////////////////////////////////////////////////
34
35// Note: The file is named TimeIntegrationSchemeFIT to parallel the
36// TimeIntegrationSchemeGLM file but the class is named
37// FractionalInTimeIntegrationScheme so keep with the factory naming
38// convention.
39
41
43{
44/**
45 * @class FractionalInTimeIntegrationScheme
46 *
47 * A fast convolution algorithm for computing solutions to (Caputo)
48 * time-fractional differential equations. This is an explicit solver
49 * that expresses the solution as an integral over a Talbot curve,
50 * which is discretized with quadrature. First-order quadrature is
51 * currently implemented (Soon be expanded to forth order).
52 */
54 std::string variant, size_t order, std::vector<NekDouble> freeParams)
55 : TimeIntegrationScheme(variant, order, freeParams),
56 m_name("FractionalInTime")
57{
58 m_variant = variant;
59 m_order = order;
60 m_freeParams = freeParams;
61
62 // Currently up to 4th order is implemented.
63 ASSERTL1(1 <= order && order <= 4,
64 "FractionalInTime Time integration scheme bad order: " +
65 std::to_string(order));
66
67 ASSERTL1(freeParams.size() == 0 || // Use defaults
68 freeParams.size() == 1 || // Alpha
69 freeParams.size() == 2 || // Base
70 freeParams.size() == 6, // Talbot quadrature rule
71 "FractionalInTime Time integration scheme invalid number "
72 "of free parameters, expected zero, one <alpha>, "
73 "two <alpha, base>, or "
74 "six <alpha, base, nQuadPts, sigma, mu0, nu> received " +
75 std::to_string(freeParams.size()));
76
77 if (freeParams.size() >= 1)
78 {
79 m_alpha = freeParams[0]; // Value for exp integration.
80 }
81
82 if (freeParams.size() >= 2)
83 {
84 m_base = freeParams[1]; // "Base" of the algorithm.
85 }
86
87 if (freeParams.size() == 6)
88 {
89 m_nQuadPts = freeParams[2]; // Number of Talbot quadrature rule points
90 m_sigma = freeParams[3];
91 m_mu0 = freeParams[4];
92 m_nu = freeParams[5];
93 }
94}
95
96/**
97 * @brief Worker method to initialize the integration scheme.
98 */
100 const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time,
102
103{
104 m_op = op;
105 m_nvars = y_0.size();
106 m_npoints = y_0[0].size();
107
108 m_deltaT = deltaT;
109
110 m_T = time; // Finial time;
112
113 // The +2 below is a buffer, and keeps +2 extra rectangle groups
114 // in case T needs to be increased later.
116
117 // Demarcation integers - one array that is re-used
118 m_qml = Array<OneD, size_t>(m_Lmax - 1, size_t(0));
119
120 // Demarcation interval markers - one array that is re-used
121 m_taus = Array<OneD, size_t>(m_Lmax + 1, size_t(0));
122
123 // Storage of the initial values.
124 m_u0 = y_0;
125
126 // Storage for the exponential factor in the integral
127 // contribution. One array that is re-used
129
130 // Storage of previous states and associated timesteps.
131 m_u = TripleArray(m_order + 1);
132
133 for (size_t m = 0; m <= m_order; ++m)
134 {
135 m_u[m] = DoubleArray(m_nvars);
136
137 for (size_t i = 0; i < m_nvars; ++i)
138 {
139 m_u[m][i] = SingleArray(m_npoints, 0.0);
140
141 for (size_t j = 0; j < m_npoints; ++j)
142 {
143 // Store the initial values as the first previous state.
144 if (m == 0)
145 {
146 m_u[m][i][j] = m_u0[i][j];
147 }
148 else
149 {
150 m_u[m][i][j] = 0;
151 }
152 }
153 }
154 }
155
156 // Storage for the stage derivative as the data will be re-used to
157 // update the solution.
159 // Storage of the next solution from the final increment.
161 // Storage for the integral contribution.
163
164 for (size_t i = 0; i < m_nvars; ++i)
165 {
166 m_F[i] = SingleArray(m_npoints, 0.0);
167 m_uNext[i] = SingleArray(m_npoints, 0.0);
169 }
170
171 // J
172 m_J = SingleArray(m_order, 0.0);
173
174 m_J[0] = pow(m_deltaT, m_alpha) / std::tgamma(m_alpha + 1.);
175
176 for (size_t m = 1, m_1 = 0; m < m_order; ++m, ++m_1)
177 {
178 m_J[m] = m_J[m_1] * NekDouble(m) / (NekDouble(m) + m_alpha);
179 }
180
181 // Ahat array, one for each order.
182 // These are elements in a multi-step exponential integrator tableau
184
185 for (size_t m = 1; m <= m_order; ++m)
186 {
187 m_Ahats[m] = DoubleArray(m);
188
189 for (size_t n = 0; n < m; ++n)
190 {
191 m_Ahats[m][n] = SingleArray(m, 0.0);
192 }
193
194 switch (m)
195 {
196 case 1:
197 m_Ahats[m][0][0] = 1.;
198 break;
199
200 case 2:
201 m_Ahats[m][0][0] = 1.;
202 m_Ahats[m][0][1] = 0.;
203 m_Ahats[m][1][0] = 1.;
204 m_Ahats[m][1][1] = -1.;
205 break;
206
207 case 3:
208 m_Ahats[m][0][0] = 1.;
209 m_Ahats[m][0][1] = 0.;
210 m_Ahats[m][0][2] = 0;
211 m_Ahats[m][1][0] = 3. / 2.;
212 m_Ahats[m][1][1] = -2.;
213 m_Ahats[m][1][2] = 1. / 2.;
214 m_Ahats[m][2][0] = 1. / 2.;
215 m_Ahats[m][2][1] = -1.;
216 m_Ahats[m][2][2] = 1. / 2.;
217 break;
218
219 case 4:
220 m_Ahats[m][0][0] = 1.;
221 m_Ahats[m][0][1] = 0.;
222 m_Ahats[m][0][2] = 0.;
223 m_Ahats[m][0][3] = 0.;
224
225 m_Ahats[m][1][0] = 11. / 6.;
226 m_Ahats[m][1][1] = -3;
227 m_Ahats[m][1][2] = 3. / 2.;
228 m_Ahats[m][1][3] = -1. / 3.;
229
230 m_Ahats[m][2][0] = 1.;
231 m_Ahats[m][2][1] = -5. / 2.;
232 m_Ahats[m][2][2] = 2.;
233 m_Ahats[m][2][3] = -1. / 2.;
234
235 m_Ahats[m][3][0] = 1. / 6.;
236 m_Ahats[m][3][1] = -1. / 2.;
237 m_Ahats[m][3][2] = 1. / 2.;
238 m_Ahats[m][3][3] = -1. / 6.;
239 break;
240
241 default:
242
243 m_Ahats[m][0][0] = 1;
244
245 for (size_t j = 2; j <= m; ++j)
246 {
247 for (size_t i = 0; i < m; ++i)
248 {
249 m_Ahats[m][j - 1][i] = pow((1 - j), i);
250 }
251 }
252
253 ASSERTL1(false, "No matrix inverse.");
254
255 // Future code: m_Ahats[m] = inv(m_Ahats[m]);
256
257 break;
258 }
259 }
260
261 // Mulitply the last Ahat array, transposed, by J
263
264 for (size_t i = 0; i < m_order; ++i)
265 {
266 for (size_t j = 0; j < m_order; ++j)
267 {
268 m_AhattJ[i] += m_Ahats[m_order][j][i] * m_J[j];
269 }
270 }
271
273
274 for (size_t l = 0; l < m_Lmax; ++l)
275 {
277 }
278}
279
280/**
281 * @brief Worker method that performs the time integration.
282 */
284 const size_t timestep, [[maybe_unused]] const NekDouble delta_t)
285{
286 ASSERTL1(delta_t == m_deltaT,
287 "Delta T has changed which is not permitted.");
288
289 // The Fractional in Time works via the logical? time step value.
290 size_t timeStep = timestep + 1;
291
292 // Update the storage and counters for integral classes. Performs
293 // staging for updating u.
294 for (size_t l = 0; l < m_Lmax; ++l)
295 {
296 updateStage(timeStep, m_integral_classes[l]);
297 }
298
299 // Compute u update to time timeStep * m_deltaT. Stored in
300 // m_uNext.
301 finalIncrement(timeStep);
302
303 // Contributions to the current integral
304 size_t L = computeTaus(m_base, timeStep);
305
306 for (size_t l = 0; l < L; ++l)
307 {
308 // Integral contribution over [taus(i+1) taus(i)]. Stored in
309 // m_uInt.
311
312 for (size_t i = 0; i < m_nvars; ++i)
313 {
314 for (size_t j = 0; j < m_npoints; ++j)
315 {
316 m_uNext[i][j] += m_uInt[i][j].real();
317 }
318 }
319 }
320
321 // Shuffle the previous solutions back one in the history.
322 for (size_t m = m_order; m > 0; --m)
323 {
324 for (size_t i = 0; i < m_nvars; ++i)
325 {
326 for (size_t j = 0; j < m_npoints; ++j)
327 {
328 m_u[m][i][j] = m_u[m - 1][i][j];
329 }
330 }
331 }
332
333 // Get the current solution.
334 for (size_t i = 0; i < m_nvars; ++i)
335 {
336 for (size_t j = 0; j < m_npoints; ++j)
337 {
338 m_u[0][i][j] = m_uNext[i][j] + m_u0[i][j];
339
340 m_uNext[i][j] = 0; // Zero out for the next itereation.
341 }
342 }
343
344 // Update the storage and counters for integral classes to
345 // time timeStep * m_deltaT. Also time-steps the sandboxes and stashes.
346 for (size_t i = 0; i < m_Lmax; ++i)
347 {
349 }
350
351 return m_u[0];
352}
353
354/**
355 * @brief Method that increments the counter then performs mod
356 * calculation.
357 */
359 const size_t base) const
360{
361 return (counter + 1) % base;
362}
363
364/**
365 * @brief Method to compute the smallest integer L such that base < 2
366 * * base^l.
367 */
369 const size_t l) const
370{
371 size_t L = ceil(log(l / 2.0) / log(base));
372
373 if (l % (size_t)(2 * pow(base, L)) == 0)
374 {
375 ++L;
376 }
377
378 return L;
379}
380
381/**
382 * @brief Method to compute the demarcation integers q_{m, ell}.
383 *
384 * Returns a length-(L-1) vector qml such that h*taus are interval
385 * boundaries for a partition of [0, m h]. The value of h is not
386 * needed to compute this vector.
387 */
389 const size_t m)
390{
391 size_t L = computeL(base, m);
392
393 // m_qml is set in InitializeScheme to be the largest length expected.
394 // qml = Array<OneD, size_t>( L-1, 0 );
395
396 for (size_t i = 0; i < L - 1; ++i)
397 {
398 m_qml[i] = floor(m / pow(base, i + 1)) - 1;
399 }
400
401 return L;
402}
403
404/**
405 * @brief Method to compute the demarcation interval marker tau_{m, ell}.
406 *
407 * Returns a length-(L+1) vector tauml such that h*taus are interval
408 * boundaries for a partition of [0, m h]. The value of h is not
409 * needed to compute this vector.
410 */
412 const size_t m)
413{
414 if (m == 1)
415 {
416 m_taus[0] = 0;
417
418 return 0;
419 }
420 else
421 {
422 size_t L = computeQML(base, m);
423
424 // m_taus is set in InitializeScheme to be the largest length
425 // expected.
426
427 m_taus[0] = m - 1;
428
429 for (size_t i = 1; i < L; ++i)
430 {
431 m_taus[i] = m_qml[i - 1] * pow(base, i);
432 }
433
434 m_taus[L] = 0;
435
436 return L;
437 }
438}
439
440/**
441 * @brief Method to compute the quadrature rule over Tablot contour
442 *
443 * Returns a quadrature rule over the Tablot contour defined by the
444 * parameterization.
445 *
446 * gamma(th) = sigma + mu * ( th*cot(th) + i*nu*th ), -pi < th < pi
447 *
448 * An N-point rule is returned, equidistant in the parameter theta. The
449 * returned quadrature rule approximes an integral over the contour.
450 */
452 const size_t nQuadPts, const NekDouble mu, const NekDouble nu,
453 const NekDouble sigma, ComplexSingleArray &lamb,
454 ComplexSingleArray &w) const
455{
456 lamb = ComplexSingleArray(nQuadPts, 0.0);
457 w = ComplexSingleArray(nQuadPts, 0.0);
458
459 for (size_t q = 0; q < nQuadPts; ++q)
460 {
461 NekDouble th =
462 (NekDouble(q) + 0.5) / NekDouble(nQuadPts) * 2.0 * M_PI - M_PI;
463
464 lamb[q] = sigma + mu * th * std::complex<NekDouble>(1. / tan(th), nu);
465
466 w[q] = std::complex<NekDouble>(0, -1. / NekDouble(nQuadPts)) * mu *
467 std::complex<NekDouble>(1. / tan(th) - th / (sin(th) * sin(th)),
468 nu);
469 }
470
471 // Special case for th = 0 which happens when there is an odd
472 // number of quadrature points.
473 if (nQuadPts % 2 == 1)
474 {
475 size_t q = (nQuadPts + 1) / 2;
476
477 lamb[q] = std::complex<NekDouble>(sigma + mu, 0);
478
479 w[q] = std::complex<NekDouble>(nu * mu / nQuadPts, 0);
480 }
481}
482
483/**
484 * @brief Method to initialize the integral class
485 */
487 const size_t index, Instance &instance) const
488{
489 /**
490 * /brief
491 *
492 * This object stores information for performing integration over
493 * an interval [a, b]. (Defined by taus in the parent calling
494 * function.)
495 *
496 * The "main" object stores information about [a,b]. In
497 * particular, main.ind identifies [a,b] via multiples of h.
498 *
499 * Periodically the values of [a,b] need to be incremented. The
500 * necessary background storage to accomplish this increment
501 * depends whether a or b is being incremented.
502 *
503 * The objects with "f" ("Floor") modifiers are associated with
504 * increments of the interval floor a.
505 *
506 * The objects with "c" ("Ceiling") modifiers are associated with
507 * increments of the interval ceiling b.
508 *
509 * Items on the "stage" are stored for use in computing u at the
510 * current time. Items in the "stash" are stored for use for
511 * future staging. Items in the "sandbox" are being actively
512 * updated at the current time for future stashing. Only items in
513 * the sandbox are time-stepped. the stage and stash locations are
514 * for storage only.
515 *
516 * This is the same for all integral classes, so there's probably
517 * a better way to engineer this. And technically, all that's
518 * needed is the array K(instance.z) anyway.
519 */
520
521 instance.base = m_base;
522 instance.index = index; // Index of this instance
523 instance.active = false; // Used to determine if active
524 instance.activecounter = 0; // Counter used to flip active bit
525 instance.activebase = 2. * pow(m_base, (index - 1));
526
527 // Storage for values of y currently used to update u
533
534 for (size_t q = 0; q < m_nvars; ++q)
535 {
541
542 for (size_t i = 0; i < m_npoints; ++i)
543 {
544 instance.stage_y[q][i] = ComplexSingleArray(m_nQuadPts, 0.0);
545 instance.cstash_y[q][i] = ComplexSingleArray(m_nQuadPts, 0.0);
546 instance.csandbox_y[q][i] = ComplexSingleArray(m_nQuadPts, 0.0);
547 instance.fstash_y[q][i] = ComplexSingleArray(m_nQuadPts, 0.0);
548 instance.fsandbox_y[q][i] = ComplexSingleArray(m_nQuadPts, 0.0);
549 }
550 }
551
552 // Major storage for auxilliary ODE solutions.
553 instance.stage_ind =
554 std::pair<size_t, size_t>(0, 0); // Time-step counters
555 // indicating the interval
556 // ymain is associated with
557
558 // Staging allocation
559 instance.stage_active = false;
560 instance.stage_ccounter = 0;
561 instance.stage_cbase = pow(m_base, index - 1); // This base is halved
562 // after the first cycle
563 instance.stage_fcounter = 0;
564 instance.stage_fbase = pow(m_base, index); // This base is halved
565 // after the first cycle
566
567 // Ceiling stash allocation
568 instance.cstash_counter = 0; // Counter used to determine
569 // when to stash
570
571 instance.cstash_base = pow(m_base, index - 1); // base for counter ind(1)
572 instance.cstash_ind =
573 std::pair<size_t, size_t>(0, 0); // is never used: it always
574 // matches main.ind(1)
575
576 // Ceiling sandbox allocation
577 instance.csandbox_active = false; // Flag to determine when stash 2
578 // is utilized
579 instance.csandbox_counter = 0;
580 instance.csandbox_ind = std::pair<size_t, size_t>(0, 0);
581
582 // Floor stash
583 instance.fstash_base = 2 * pow(m_base, index);
584 instance.fstash_ind = std::pair<size_t, size_t>(0, 0);
585
586 // Floor sandbox
587 instance.fsandbox_active = false;
588 instance.fsandbox_activebase = pow(m_base, index);
589 instance.fsandbox_stashincrement = (m_base - 1) * pow(m_base, index - 1);
590 instance.fsandbox_ind = std::pair<size_t, size_t>(0, 0);
591
592 // Defining parameters of the Talbot contour quadrature rule
593 NekDouble Tl =
594 m_deltaT * (2. * pow(m_base, index) - 1. - pow(m_base, index - 1));
595 NekDouble mu = m_mu0 / Tl;
596
597 // Talbot quadrature rule
598 talbotQuadrature(m_nQuadPts, mu, m_nu, m_sigma, instance.z, instance.w);
599
600 /**
601 * /brief
602 *
603 * With sigma == 0, the dependence of z and w on index is just a
604 * multiplicative scaling factor (mu). So technically we'll only
605 * need one instance of this N-point rule and can scale it
606 * accordingly inside each integral_class instance. Not sure if
607 * this optimization is worth it. Cumulative memory savings would
608 * only be about 4*N*Lmax floats.
609
610 * Below: precomputation for time integration of auxiliary
611 * variables. Everything below here is independent of the
612 * instance index index. Therefore, we could actually just
613 * generate and store one copy of this stuff and use it
614 * everywhere.
615 */
616
617 // 'As' array - one for each order.
618 TripleArray &As = instance.As;
619
620 As = TripleArray(m_order + 2);
621
622 for (size_t m = 1; m <= m_order + 1; ++m)
623 {
624 As[m] = DoubleArray(m);
625
626 for (size_t n = 0; n < m; ++n)
627 {
628 As[m][n] = SingleArray(m, 0.0);
629 }
630
631 switch (m)
632 {
633 case 1:
634 As[m][0][0] = 1.;
635 break;
636
637 case 2:
638 As[m][0][0] = 0.;
639 As[m][0][1] = 1.;
640 As[m][1][0] = 1.;
641 As[m][1][1] = -1.;
642 break;
643
644 case 3:
645 As[m][0][0] = 0.;
646 As[m][0][1] = 1.;
647 As[m][0][2] = 0;
648 As[m][1][0] = 1. / 2.;
649 As[m][1][1] = 0.;
650 As[m][1][2] = -1. / 2.;
651 As[m][2][0] = 1. / 2.;
652 As[m][2][1] = -1.;
653 As[m][2][2] = 1. / 2.;
654 break;
655
656 case 4:
657 As[m][0][0] = 0.;
658 As[m][0][1] = 1.;
659 As[m][0][2] = 0.;
660 As[m][0][3] = 0.;
661
662 As[m][1][0] = 1. / 3.;
663 As[m][1][1] = 1. / 2.;
664 As[m][1][2] = -1.;
665 As[m][1][3] = 1. / 6.;
666
667 As[m][2][0] = 1. / 2.;
668 As[m][2][1] = -1.;
669 As[m][2][2] = 1. / 2.;
670 As[m][2][3] = 0.;
671
672 As[m][3][0] = 1. / 6.;
673 As[m][3][1] = -1. / 2.;
674 As[m][3][2] = 1. / 2.;
675 As[m][3][3] = -1. / 6.;
676 break;
677
678 case 5:
679 As[m][0][0] = 0.;
680 As[m][0][1] = 1.;
681 As[m][0][2] = 0.;
682 As[m][0][3] = 0.;
683 As[m][0][4] = 0.;
684
685 As[m][1][0] = 1. / 4.;
686 As[m][1][1] = 5. / 6.;
687 As[m][1][2] = -3. / 2.;
688 As[m][1][3] = 1. / 2.;
689 As[m][1][4] = -1. / 12.;
690
691 As[m][2][0] = 11. / 24.;
692 As[m][2][1] = -5. / 6.;
693 As[m][2][2] = 1. / 4.;
694 As[m][2][3] = 1. / 6.;
695 As[m][2][4] = -1. / 24.;
696
697 As[m][3][0] = 1. / 4.;
698 As[m][3][1] = -5. / 6.;
699 As[m][3][2] = 1.;
700 As[m][3][3] = -1. / 2.;
701 As[m][3][4] = 1. / 12.;
702
703 As[m][4][0] = 1. / 24.;
704 As[m][4][1] = -1. / 6.;
705 As[m][4][2] = 1. / 4.;
706 As[m][4][3] = -1. / 6.;
707 As[m][4][4] = 1. / 24.;
708 break;
709
710 // The default is a general formula, but the matrix inversion
711 // involved is ill-conditioned, so the special cases above are
712 // epxlicitly given to combat roundoff error in the most-used
713 // scenarios.
714 default:
715 ASSERTL1(false, "No matrix inverse.");
716 break;
717 }
718 }
719
720 // Initialize the exponenetial integrators.
721 instance.E = ComplexSingleArray(m_nQuadPts, 0.0);
722
723 for (size_t q = 0; q < m_nQuadPts; ++q)
724 {
725 instance.E[q] = exp(instance.z[q] * m_deltaT);
726 }
727
728 instance.Eh = ComplexDoubleArray(m_order + 1);
729
730 for (size_t m = 0; m < m_order + 1; ++m)
731 {
732 instance.Eh[m] = ComplexSingleArray(m_nQuadPts, 0.0);
733
734 for (size_t q = 0; q < m_nQuadPts; ++q)
735 {
736 if (m == 0)
737 {
738 instance.Eh[0][q] =
739 1. / instance.z[q] * (exp(instance.z[q] * m_deltaT) - 1.0);
740 }
741 else
742 {
743 instance.Eh[m][q] = -1. / instance.z[q] +
744 NekDouble(m) / (instance.z[q] * m_deltaT) *
745 instance.Eh[m - 1][q];
746 }
747 }
748 }
749
750 // 'AtEh' is set for the primary order. If a lower order method is
751 // needed for initializing it will be changed in time_advance then
752 // restored.
753 instance.AtEh = ComplexDoubleArray(m_order + 1);
754
755 for (size_t m = 0; m <= m_order; ++m)
756 {
757 instance.AtEh[m] = ComplexSingleArray(m_nQuadPts, 0.0);
758
759 for (size_t q = 0; q < m_nQuadPts; ++q)
760 {
761 for (size_t i = 0; i <= m_order; ++i)
762 {
763 instance.AtEh[m][q] +=
764 instance.As[m_order + 1][m][i] * instance.Eh[i][q];
765 }
766 }
767 }
768}
769
770/**
771 * @brief Method to rearrange of staging/stashing for current time
772 *
773 * (1) activates ceiling staging
774 * (2) moves ceiling stash ---> stage
775 * (3) moves floor stash --> stage (+ updates all ceiling data)
776 */
778 Instance &instance)
779{
780 // Counter to flip active bit
781 if (!instance.active)
782 {
783 instance.active = (timeStep % instance.activebase == 0);
784 }
785
786 // Determine if staging is necessary
787 if (instance.active)
788 {
789 // Floor staging superscedes ceiling staging
790 if (timeStep % instance.fstash_base == 0)
791 {
792 // Here a swap of the contents can be done because values
793 // will copied into the stash and the f sandbox values will
794 // cleared next.
795 std::swap(instance.stage_y, instance.fstash_y);
796 instance.stage_ind = instance.fstash_ind;
797
798 std::swap(instance.csandbox_y, instance.fsandbox_y);
799 instance.csandbox_ind = instance.fsandbox_ind;
800
801 // After floor staging happens once, new base is base^index
802 instance.fstash_base = pow(instance.base, instance.index);
803
804 // Restart floor sandbox
805 instance.fsandbox_ind = std::pair<size_t, size_t>(0, 0);
806 instance.fsandbox_active = false;
807
808 // Clear the floor sandbox values.
809 for (size_t i = 0; i < m_nvars; ++i)
810 {
811 for (size_t j = 0; j < m_npoints; ++j)
812 {
813 for (size_t q = 0; q < m_nQuadPts; ++q)
814 {
815 instance.fsandbox_y[i][j][q] = 0;
816 }
817 }
818 }
819 }
820
821 // Check for ceiling staging
822 else if (timeStep % instance.stage_cbase == 0)
823 {
824 instance.stage_ind = instance.cstash_ind;
825
826 // A swap of the contents can be done because values will
827 // copied into the stash.
828 std::swap(instance.stage_y, instance.cstash_y);
829 }
830 }
831}
832
833/**
834 * @brief Method to approximate the integral
835 *
836 * \int_{(m-1) h}^{m h} k(m*h -s) f(u, s) dx{s}
837 *
838 * Using a time-stepping scheme of a particular order. Here, k depends
839 * on alpha, the derivative order.
840 */
842{
843 // Note: m_uNext is initialized to zero and then reset to zero
844 // after it is used to update the current solution in TimeIntegrate.
845 for (size_t m = 0; m < m_order; ++m)
846 {
847 m_op.DoOdeRhs(m_u[m], m_F, m_deltaT * (timeStep - m));
848
849 for (size_t i = 0; i < m_nvars; ++i)
850 {
851 for (size_t j = 0; j < m_npoints; ++j)
852 {
853 m_uNext[i][j] += m_F[i][j] * m_AhattJ[m];
854 }
855 }
856 }
857}
858
859/**
860 * @brief Method to get the integral contribution over [taus(i+1)
861 * taus(i)]. Stored in m_uInt.
862 */
864 const size_t timeStep, const size_t tauml, const Instance &instance)
865{
866 // Assume y has already been updated to time level m
867 for (size_t q = 0; q < m_nQuadPts; ++q)
868 {
869 m_expFactor[q] =
870 exp(instance.z[q] * m_deltaT * NekDouble(timeStep - tauml)) *
871 pow(instance.z[q], -m_alpha) * instance.w[q];
872 }
873
874 for (size_t i = 0; i < m_nvars; ++i)
875 {
876 for (size_t j = 0; j < m_npoints; ++j)
877 {
878 m_uInt[i][j] = 0;
879
880 for (size_t q = 0; q < m_nQuadPts; ++q)
881 {
882 m_uInt[i][j] += instance.stage_y[i][j][q] * m_expFactor[q];
883 }
884
885 if (m_uInt[i][j].real() < 1e8)
886 {
887 m_uInt[i][j] = m_uInt[i][j].real();
888 }
889 }
890 }
891}
892
893/**
894 * @brief Method to get the solution to y' = z*y + f(u), using an
895 * exponential integrator with implicit order (m_order + 1)
896 * interpolation of the f(u) term.
897 */
899 Instance &instance,
901{
902 size_t order;
903
904 // Try automated high-order method.
905 if (timeStep <= m_order)
906 {
907 // Not enough history. For now, demote to lower-order method.
908 // TODO: use multi-stage method.
909 order = timeStep;
910
911 // Prep for the time step.
912 for (size_t m = 0; m <= order; ++m)
913 {
914 for (size_t q = 0; q < m_nQuadPts; ++q)
915 {
916 instance.AtEh[m][q] = 0;
917
918 for (size_t i = 0; i <= order; ++i)
919 {
920 instance.AtEh[m][q] +=
921 instance.As[order + 1][m][i] * instance.Eh[i][q];
922 }
923 }
924 }
925 }
926 else
927 {
928 order = m_order;
929 }
930
931 // y = y * instance.E + F * instance.AtEh;
932 for (size_t m = 0; m <= order; ++m)
933 {
934 m_op.DoOdeRhs(m_u[m], m_F, m_deltaT * (timeStep - m));
935
936 for (size_t i = 0; i < m_nvars; ++i)
937 {
938 for (size_t j = 0; j < m_npoints; ++j)
939 {
940 for (size_t q = 0; q < m_nQuadPts; ++q)
941 {
942 // y * instance.E
943 if (m == 0)
944 {
945 y[i][j][q] *= instance.E[q];
946 }
947
948 // F * instance.AtEh
949 y[i][j][q] += m_F[i][j] * instance.AtEh[m][q];
950 }
951 }
952 }
953 }
954}
955
956/**
957 * @brief Method to update sandboxes to the current time.
958 *
959 * (1) advances ceiling sandbox
960 * (2) moves ceiling sandbox ---> stash
961 * (3) activates floor sandboxing
962 * (4) advances floor sandbox
963 * (5) moves floor sandbox ---> stash
964 */
966 Instance &instance)
967{
968 // (1)
969 // update(instance.csandbox.y)
970 timeAdvance(timeStep, instance, instance.csandbox_y);
971 instance.csandbox_ind.second = timeStep;
972
973 // (2)
974 // Determine if ceiling stashing is necessary
975 instance.cstash_counter =
976 modIncrement(instance.cstash_counter, instance.cstash_base);
977
978 if (timeStep % instance.cstash_base == 0)
979 {
980 // Then need to stash
981 // instance.cstash_y = instance.csandbox_y;
982 instance.cstash_ind = instance.csandbox_ind;
983
984 // Stash the c sandbox value. This step has to be a deep copy
985 // because the values in the sandbox are still needed for the
986 // time advance.
987 for (size_t i = 0; i < m_nvars; ++i)
988 {
989 for (size_t j = 0; j < m_npoints; ++j)
990 {
991 for (size_t q = 0; q < m_nQuadPts; ++q)
992 {
993 instance.cstash_y[i][j][q] = instance.csandbox_y[i][j][q];
994 }
995 }
996 }
997 }
998
999 if (instance.fsandbox_active)
1000 {
1001 // (4)
1002 timeAdvance(timeStep, instance, instance.fsandbox_y);
1003
1004 instance.fsandbox_ind.second = timeStep;
1005
1006 // (5) Move floor sandbox to stash
1007 if ((instance.fsandbox_ind.second - instance.fsandbox_ind.first) %
1008 instance.fsandbox_stashincrement ==
1009 0)
1010 {
1011 // instance.fstash_y = instance.fsandbox_y;
1012 instance.fstash_ind = instance.fsandbox_ind;
1013
1014 // Stash the f sandbox values. This step has to be a deep
1015 // copy because the values in the sandbox are still needed
1016 // for the time advance.
1017 for (size_t i = 0; i < m_nvars; ++i)
1018 {
1019 for (size_t j = 0; j < m_npoints; ++j)
1020 {
1021 for (size_t q = 0; q < m_nQuadPts; ++q)
1022 {
1023 instance.fstash_y[i][j][q] =
1024 instance.fsandbox_y[i][j][q];
1025 }
1026 }
1027 }
1028 }
1029 }
1030 else // Determine if advancing floor sandbox is necessary at next time
1031 {
1032 // (3)
1033 if (timeStep % instance.fsandbox_activebase == 0)
1034 {
1035 instance.fsandbox_active = true;
1036 instance.fsandbox_ind =
1037 std::pair<size_t, size_t>(timeStep, timeStep);
1038 }
1039 }
1040}
1041
1042/**
1043 * @brief Worker method to print details on the integration scheme
1044 */
1046{
1047 os << "Time Integration Scheme: " << GetFullName() << std::endl
1048 << " Alpha " << m_alpha << std::endl
1049 << " Base " << m_base << std::endl
1050 << " Number of instances " << m_Lmax << std::endl
1051 << " Number of quadature points " << m_nQuadPts << std::endl
1052 << " Talbot Parameter: sigma " << m_sigma << std::endl
1053 << " Talbot Parameter: mu0 " << m_mu0 << std::endl
1054 << " Talbot Parameter: nu " << m_nu << std::endl;
1055}
1056
1058{
1059 os << "Time Integration Scheme: " << GetFullName() << std::endl
1060 << " Alpha " << m_alpha << std::endl
1061 << " Base " << m_base << std::endl
1062 << " Number of instances " << m_Lmax << std::endl
1063 << " Number of quadature points " << m_nQuadPts << std::endl
1064 << " Talbot Parameter: sigma " << m_sigma << std::endl
1065 << " Talbot Parameter: mu0 " << m_mu0 << std::endl
1066 << " Talbot Parameter: nu " << m_nu << std::endl;
1067}
1068
1069// Friend Operators
1070std::ostream &operator<<(std::ostream &os,
1072{
1073 rhs.print(os);
1074
1075 return os;
1076}
1077
1078std::ostream &operator<<(std::ostream &os,
1080{
1081 os << *rhs.get();
1082
1083 return os;
1084}
1085
1086} // namespace Nektar::LibUtilities
NekDouble L
#define ASSERTL1(condition, msg)
Assert Level 1 – Debugging which is used whether in FULLDEBUG or DEBUG compilation mode....
Definition: ErrorUtil.hpp:242
size_t computeTaus(const size_t base, const size_t m)
Method to compute the demarcation interval marker tau_{m, ell}.
LUE void v_InitializeScheme(const NekDouble deltaT, ConstDoubleArray &y_0, const NekDouble time, const TimeIntegrationSchemeOperators &op) override
Worker method to initialize the integration scheme.
void timeAdvance(const size_t timeStep, Instance &instance, ComplexTripleArray &y)
Method to get the solution to y' = z*y + f(u), using an exponential integrator with implicit order (m...
void updateStage(const size_t timeStep, Instance &instance)
Method to rearrange of staging/stashing for current time.
void advanceSandbox(const size_t timeStep, Instance &instance)
Method to update sandboxes to the current time.
void integralClassInitialize(const size_t index, Instance &instance) const
Method to initialize the integral class.
LUE void v_print(std::ostream &os) const override
Worker method to print details on the integration scheme.
size_t modIncrement(const size_t counter, const size_t base) const
Method that increments the counter then performs mod calculation.
void finalIncrement(const size_t timeStep)
Method to approximate the integral.
void integralContribution(const size_t timeStep, const size_t tauml, const Instance &instance)
Method to get the integral contribution over [taus(i+1) taus(i)]. Stored in m_uInt.
void talbotQuadrature(const size_t nQuadPts, const NekDouble mu, const NekDouble nu, const NekDouble sigma, ComplexSingleArray &lamb, ComplexSingleArray &w) const
Method to compute the quadrature rule over Tablot contour.
size_t computeL(const size_t base, const size_t m) const
Method to compute the smallest integer L such that base < 2.
LUE ConstDoubleArray & v_TimeIntegrate(const size_t timestep, const NekDouble delta_t) override
Worker method that performs the time integration.
FractionalInTimeIntegrationScheme(std::string variant, size_t order, std::vector< NekDouble > freeParams)
Constructor.
size_t computeQML(const size_t base, const size_t m)
Method to compute the demarcation integers q_{m, ell}.
Base class for time integration schemes.
Binds a set of functions for use by time integration schemes.
void DoOdeRhs(InArrayType &inarray, OutArrayType &outarray, const NekDouble time) const
AT< AT< NekDouble > > DoubleArray
AT< AT< std::complex< NekDouble > > > ComplexDoubleArray
AT< std::complex< NekDouble > > ComplexSingleArray
AT< AT< AT< NekDouble > > > TripleArray
std::shared_ptr< FractionalInTimeIntegrationScheme > FractionalInTimeIntegrationSchemeSharedPtr
AT< AT< AT< std::complex< NekDouble > > > > ComplexTripleArray
std::ostream & operator<<(std::ostream &os, const BasisKey &rhs)
std::vector< double > w(NPUPPER)
std::vector< double > q(NPUPPER *NPUPPER)
double NekDouble
scalarT< T > log(scalarT< T > in)
Definition: scalar.hpp:303