DUNE PDELab (git)

onestepparameter.hh
1// -*- tab-width: 2; indent-tabs-mode: nil -*-
2// vi: set et ts=2 sw=2 sts=2:
3
4#ifndef DUNE_PDELAB_INSTATIONARY_ONESTEPPARAMETER_HH
5#define DUNE_PDELAB_INSTATIONARY_ONESTEPPARAMETER_HH
6
10
11namespace Dune {
12 namespace PDELab {
13
19
42 template<class R>
44 {
45 public:
46 typedef R RealType;
47
50 virtual bool implicit () const = 0;
51
54 virtual unsigned s () const = 0;
55
59 virtual R a (int r, int i) const = 0;
60
64 virtual R b (int r, int i) const = 0;
65
69 virtual R d (int r) const = 0;
70
73 virtual std::string name () const = 0;
74
77 };
78
79
88 template<class R>
90 {
93
94 public:
97 : theta(theta_)
98 {
99 D[0] = 0.0; D[1] = 1.0;
100 A[0][0] = -1.0; A[0][1] = 1.0;
101 B[0][0] = 1.0-theta; B[0][1] = theta;
102 }
103
106 virtual bool implicit () const override
107 {
108 return (theta > 0.0);
109 }
110
113 virtual unsigned s () const override
114 {
115 return 1;
116 }
117
121 virtual R a (int r, int i) const override
122 {
123 return A[r-1][i];
124 }
125
129 virtual R b (int r, int i) const override
130 {
131 return B[r-1][i];
132 }
133
137 virtual R d (int i) const override
138 {
139 return D[i];
140 }
141
144 virtual std::string name () const override
145 {
146 return std::string("one step theta");
147 }
148
149 private:
150 R theta;
154 };
155
156
163 template<class R>
165 {
166 public:
167
170 {}
171
174 virtual std::string name () const override
175 {
176 return std::string("explicit Euler");
177 }
178
179 };
180
181
188 template<class R>
190 {
191 public:
192
195 {}
196
199 virtual std::string name () const override
200 {
201 return std::string("implicit Euler");
202 }
203
204 };
205
206
213 template<class R>
215 {
216 public:
217
219 {
220 D[0] = 0.0; D[1] = 1.0; D[2] = 1.0;
221
222 A[0][0] = -1.0; A[0][1] = 1.0; A[0][2] = 0.0;
223 A[1][0] = -0.5; A[1][1] = -0.5; A[1][2] = 1.0;
224
225 B[0][0] = 1.0; B[0][1] = 0.0; B[0][2] = 0.0;
226 B[1][0] = 0.0; B[1][1] = 0.5; B[1][2] = 0.0;
227 }
228
231 virtual bool implicit () const override
232 {
233 return false;
234 }
235
238 virtual unsigned s () const override
239 {
240 return 2;
241 }
242
246 virtual R a (int r, int i) const override
247 {
248 return A[r-1][i];
249 }
250
254 virtual R b (int r, int i) const override
255 {
256 return B[r-1][i];
257 }
258
262 virtual R d (int i) const override
263 {
264 return D[i];
265 }
266
269 virtual std::string name () const override
270 {
271 return std::string("Heun");
272 }
273
274 private:
278 };
279
286 template<class R>
288 {
289 public:
290
292 {
293 D[0] = 0.0; D[1] = 1.0; D[2] = 0.5; D[3] = 1.0;
294
295 A[0][0] = -1.0; A[0][1] = 1.0; A[0][2] = 0.0; A[0][3] = 0.0;
296 A[1][0] = -0.75; A[1][1] = -0.25; A[1][2] = 1.0; A[1][3] = 0.0;
297 A[2][0] = -1.0/3.0; A[2][1] = 0.0; A[2][2] = -2.0/3.0; A[2][3] = 1.0;
298
299 B[0][0] = 1.0; B[0][1] = 0.0; B[0][2] = 0.0; B[0][3] = 0.0;
300 B[1][0] = 0.0; B[1][1] = 0.25; B[1][2] = 0.0; B[1][3] = 0.0;
301 B[2][0] = 0.0; B[2][1] = 0.0; B[2][2] = 2.0/3.0; B[2][3] = 0.0;
302
303 }
304
307 virtual bool implicit () const override
308 {
309 return false;
310 }
311
314 virtual unsigned s () const override
315 {
316 return 3;
317 }
318
322 virtual R a (int r, int i) const override
323 {
324 return A[r-1][i];
325 }
326
330 virtual R b (int r, int i) const override
331 {
332 return B[r-1][i];
333 }
334
338 virtual R d (int i) const override
339 {
340 return D[i];
341 }
342
345 virtual std::string name () const override
346 {
347 return std::string("Shu's third order method");
348 }
349
350 private:
354 };
355
356
363 template<class R>
365 {
366 public:
367
368 RK4Parameter ()
369 {
370 D[0] = 0.0; D[1] = 0.5; D[2] = 0.5; D[3] = 1.0; D[4] = 1.0;
371
372 A[0][0] = -1.0; A[0][1] = 1.0; A[0][2] = 0.0; A[0][3] = 0.0; A[0][4] = 0.0;
373 A[1][0] = -1.0; A[1][1] = 0.0; A[1][2] = 1.0; A[1][3] = 0.0; A[1][4] = 0.0;
374 A[2][0] = -1.0; A[2][1] = 0.0; A[2][2] = 0.0; A[2][3] = 1.0; A[2][4] = 0.0;
375 A[3][0] = -1.0; A[3][1] = 0.0; A[3][2] = 0.0; A[3][3] = 0.0; A[3][4] = 1.0;
376
377 B[0][0] = 0.5; B[0][1] = 0.0; B[0][2] = 0.0; B[0][3] = 0.0; B[0][4] = 0.0;
378 B[1][0] = 0.0; B[1][1] = 0.5; B[1][2] = 0.0; B[1][3] = 0.0; B[1][4] = 0.0;
379 B[2][0] = 0.0; B[2][1] = 0.0; B[2][2] = 1.0; B[2][3] = 0.0; B[2][4] = 0.0;
380 B[3][0] = 1.0/6.0; B[3][1] = 1.0/3.0; B[3][2] = 1.0/3.0; B[3][3] = 1.0/6.0; B[3][4] = 0.0;
381
382 }
383
386 virtual bool implicit () const override
387 {
388 return false;
389 }
390
393 virtual unsigned s () const override
394 {
395 return 4;
396 }
397
401 virtual R a (int r, int i) const override
402 {
403 return A[r-1][i];
404 }
405
409 virtual R b (int r, int i) const override
410 {
411 return B[r-1][i];
412 }
413
417 virtual R d (int i) const override
418 {
419 return D[i];
420 }
421
424 virtual std::string name () const override
425 {
426 return std::string("RK4");
427 }
428
429 private:
433 };
434
435
436
437
444 template<class R>
446 {
447 public:
448
450 {
451 using std::sqrt;
452 alpha = 1.0 - 0.5*sqrt(2.0);
453
454 D[0] = 0.0; D[1] = alpha; D[2] = 1.0;
455
456 A[0][0] = -1.0; A[0][1] = 1.0; A[0][2] = 0.0;
457 A[1][0] = -1.0; A[1][1] = 0.0; A[1][2] = 1.0;
458
459 B[0][0] = 0.0; B[0][1] = alpha; B[0][2] = 0.0;
460 B[1][0] = 0.0; B[1][1] = 1.0-alpha; B[1][2] = alpha;
461 }
462
465 virtual bool implicit () const override
466 {
467 return true;
468 }
469
472 virtual unsigned s () const override
473 {
474 return 2;
475 }
476
480 virtual R a (int r, int i) const override
481 {
482 return A[r-1][i];
483 }
484
488 virtual R b (int r, int i) const override
489 {
490 return B[r-1][i];
491 }
492
496 virtual R d (int i) const override
497 {
498 return D[i];
499 }
500
503 virtual std::string name () const override
504 {
505 return std::string("Alexander (order 2)");
506 }
507
508 private:
509 R alpha;
513 };
514
521 template<class R>
523 {
524 public:
525
527 {
528 using std::sqrt;
529 R alpha, theta, thetap, beta;
530 theta = 1.0 - 0.5*sqrt(2.0);
531 thetap = 1.0-2.0*theta;
532 alpha = 2.0-sqrt(2.0);
533 beta = 1.0-alpha;
534
535 D[0] = 0.0; D[1] = theta; D[2] = 1.0-theta; D[3] = 1.0;
536
537 A[0][0] = -1.0; A[0][1] = 1.0; A[0][2] = 0.0; A[0][3] = 0.0;
538 A[1][0] = 0.0; A[1][1] = -1.0; A[1][2] = 1.0; A[1][3] = 0.0;
539 A[2][0] = 0.0; A[2][1] = 0.0; A[2][2] = -1.0; A[2][3] = 1.0;
540
541 B[0][0] = beta*theta; B[0][1] = alpha*theta; B[0][2] = 0.0; B[0][3] = 0.0;
542 B[1][0] = 0.0; B[1][1] = alpha*thetap; B[1][2] = alpha*theta; B[1][3] = 0.0;
543 B[2][0] = 0.0; B[2][1] = 0.0; B[2][2] = beta*theta; B[2][3] = alpha*theta;
544 }
545
548 virtual bool implicit () const override
549 {
550 return true;
551 }
552
555 virtual unsigned s () const override
556 {
557 return 3;
558 }
559
563 virtual R a (int r, int i) const override
564 {
565 return A[r-1][i];
566 }
567
571 virtual R b (int r, int i) const override
572 {
573 return B[r-1][i];
574 }
575
579 virtual R d (int i) const override
580 {
581 return D[i];
582 }
583
586 virtual std::string name () const override
587 {
588 return std::string("Fractional step theta");
589 }
590
591 private:
595 };
596
604 template<class R>
606 {
607 public:
608
610 {
611 R alpha = 0.4358665215;
612
613 // Newton iteration for alpha
614 for (int i=1; i<=10; i++)
615 {
616 alpha = alpha - (alpha*(alpha*alpha-3.0*(alpha-0.5))-1.0/6.0)/(3.0*alpha*(alpha-2.0)+1.5);
617 // std::cout.precision(16);
618 // std::cout << "alpha "
619 // << std::setw(8) << i << " "
620 // << std::scientific << alpha << std::endl;
621 }
622
623 R tau2 = (1.0+alpha)*0.5;
624 R b1 = -(6.0*alpha*alpha -16.0*alpha + 1.0)*0.25;
625 R b2 = (6*alpha*alpha - 20.0*alpha + 5.0)*0.25;
626
627 // std::cout.precision(16);
628 // std::cout << "alpha " << std::scientific << alpha << std::endl;
629 // std::cout << "tau2 " << std::scientific << tau2 << std::endl;
630 // std::cout << "b1 " << std::scientific << b1 << std::endl;
631 // std::cout << "b2 " << std::scientific << b2 << std::endl;
632
633 D[0] = 0.0; D[1] = alpha; D[2] = tau2; D[3] = 1.0;
634
635 A[0][0] = -1.0; A[0][1] = 1.0; A[0][2] = 0.0; A[0][3] = 0.0;
636 A[1][0] = -1.0; A[1][1] = 0.0; A[1][2] = 1.0; A[1][3] = 0.0;
637 A[2][0] = -1.0; A[2][1] = 0.0; A[2][2] = 0.0; A[2][3] = 1.0;
638
639 B[0][0] = 0.0; B[0][1] = alpha; B[0][2] = 0.0; B[0][3] = 0.0;
640 B[1][0] = 0.0; B[1][1] = tau2-alpha; B[1][2] = alpha; B[1][3] = 0.0;
641 B[2][0] = 0.0; B[2][1] = b1; B[2][2] = b2; B[2][3] = alpha;
642 }
643
646 virtual bool implicit () const override
647 {
648 return true;
649 }
650
653 virtual unsigned s () const override
654 {
655 return 3;
656 }
657
661 virtual R a (int r, int i) const override
662 {
663 return A[r-1][i];
664 }
665
669 virtual R b (int r, int i) const override
670 {
671 return B[r-1][i];
672 }
673
677 virtual R d (int i) const override
678 {
679 return D[i];
680 }
681
684 virtual std::string name () const override
685 {
686 return std::string("Alexander (claims order 3)");
687 }
688
689 private:
690 R alpha, theta, thetap, beta;
694 };
695
696 } // end namespace PDELab
697} // end namespace Dune
698#endif // DUNE_PDELAB_INSTATIONARY_ONESTEPPARAMETER_HH
Parameters to turn the OneStepMethod into an Alexander scheme.
Definition: onestepparameter.hh:446
Parameters to turn the OneStepMethod into an Alexander3 scheme.
Definition: onestepparameter.hh:606
Parameters to turn the ExplicitOneStepMethod into an explicit Euler method.
Definition: onestepparameter.hh:165
Parameters to turn the OneStepMethod into a fractional step theta scheme.
Definition: onestepparameter.hh:523
Parameters to turn the ExplicitOneStepMethod into a Heun scheme.
Definition: onestepparameter.hh:215
Parameters to turn the OneStepMethod into an implicit Euler method.
Definition: onestepparameter.hh:190
Parameters to turn the OneStepMethod into an one step theta method.
Definition: onestepparameter.hh:90
Parameters to turn the ExplicitOneStepMethod into a classical fourth order Runge-Kutta method.
Definition: onestepparameter.hh:365
Parameters to turn the ExplicitOneStepMethod into a third order strong stability preserving (SSP) sch...
Definition: onestepparameter.hh:288
Base parameter class for time stepping scheme parameters.
Definition: onestepparameter.hh:44
A few common exception classes.
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
virtual bool implicit() const override
Return true if method is implicit.
Definition: onestepparameter.hh:231
virtual R d(int i) const override
Return entries of the d Vector.
Definition: onestepparameter.hh:417
OneStepThetaParameter(R theta_)
construct OneStepThetaParameter class
Definition: onestepparameter.hh:96
virtual R d(int i) const override
Return entries of the d Vector.
Definition: onestepparameter.hh:338
virtual std::string name() const override
Return name of the scheme.
Definition: onestepparameter.hh:345
virtual R b(int r, int i) const override
Return entries of the B matrix.
Definition: onestepparameter.hh:571
virtual R a(int r, int i) const override
Return entries of the A matrix.
Definition: onestepparameter.hh:563
virtual R b(int r, int i) const override
Return entries of the B matrix.
Definition: onestepparameter.hh:129
virtual R b(int r, int i) const override
Return entries of the B matrix.
Definition: onestepparameter.hh:254
virtual unsigned s() const override
Return number of stages s of the method.
Definition: onestepparameter.hh:653
virtual R a(int r, int i) const override
Return entries of the A matrix.
Definition: onestepparameter.hh:480
virtual bool implicit() const override
Return true if method is implicit.
Definition: onestepparameter.hh:548
virtual bool implicit() const override
Return true if method is implicit.
Definition: onestepparameter.hh:646
virtual std::string name() const override
Return name of the scheme.
Definition: onestepparameter.hh:199
virtual unsigned s() const override
Return number of stages s of the method.
Definition: onestepparameter.hh:393
virtual std::string name() const override
Return name of the scheme.
Definition: onestepparameter.hh:503
virtual bool implicit() const override
Return true if method is implicit.
Definition: onestepparameter.hh:386
virtual R d(int i) const override
Return entries of the d Vector.
Definition: onestepparameter.hh:137
virtual unsigned s() const override
Return number of stages s of the method.
Definition: onestepparameter.hh:314
virtual R a(int r, int i) const override
Return entries of the A matrix.
Definition: onestepparameter.hh:322
virtual std::string name() const override
Return name of the scheme.
Definition: onestepparameter.hh:269
virtual R b(int r, int i) const override
Return entries of the B matrix.
Definition: onestepparameter.hh:488
virtual std::string name() const override
Return name of the scheme.
Definition: onestepparameter.hh:174
virtual std::string name() const override
Return name of the scheme.
Definition: onestepparameter.hh:424
virtual std::string name() const override
Return name of the scheme.
Definition: onestepparameter.hh:586
virtual R d(int i) const override
Return entries of the d Vector.
Definition: onestepparameter.hh:677
virtual bool implicit() const =0
Return true if method is implicit.
virtual R b(int r, int i) const override
Return entries of the B matrix.
Definition: onestepparameter.hh:330
virtual std::string name() const override
Return name of the scheme.
Definition: onestepparameter.hh:144
virtual R b(int r, int i) const override
Return entries of the B matrix.
Definition: onestepparameter.hh:669
virtual unsigned s() const override
Return number of stages s of the method.
Definition: onestepparameter.hh:472
virtual std::string name() const =0
Return name of the scheme.
virtual unsigned s() const override
Return number of stages s of the method.
Definition: onestepparameter.hh:238
virtual bool implicit() const override
Return true if method is implicit.
Definition: onestepparameter.hh:465
virtual std::string name() const override
Return name of the scheme.
Definition: onestepparameter.hh:684
virtual R d(int i) const override
Return entries of the d Vector.
Definition: onestepparameter.hh:262
virtual R a(int r, int i) const override
Return entries of the A matrix.
Definition: onestepparameter.hh:661
virtual ~TimeSteppingParameterInterface()
every abstract base class has a virtual destructor
Definition: onestepparameter.hh:76
virtual R d(int r) const =0
Return entries of the d Vector.
virtual R b(int r, int i) const =0
Return entries of the B matrix.
virtual unsigned s() const =0
Return number of stages of the method.
virtual R a(int r, int i) const =0
Return entries of the A matrix.
virtual unsigned s() const override
Return number of stages s of the method.
Definition: onestepparameter.hh:113
virtual R a(int r, int i) const override
Return entries of the A matrix.
Definition: onestepparameter.hh:121
virtual R b(int r, int i) const override
Return entries of the B matrix.
Definition: onestepparameter.hh:409
virtual R a(int r, int i) const override
Return entries of the A matrix.
Definition: onestepparameter.hh:401
virtual R d(int i) const override
Return entries of the d Vector.
Definition: onestepparameter.hh:496
virtual R a(int r, int i) const override
Return entries of the A matrix.
Definition: onestepparameter.hh:246
virtual bool implicit() const override
Return true if method is implicit.
Definition: onestepparameter.hh:106
virtual bool implicit() const override
Return true if method is implicit.
Definition: onestepparameter.hh:307
virtual R d(int i) const override
Return entries of the d Vector.
Definition: onestepparameter.hh:579
virtual unsigned s() const override
Return number of stages s of the method.
Definition: onestepparameter.hh:555
Dune namespace.
Definition: alignedallocator.hh:13
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)