MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
ode.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
13#include "operator.hpp"
14#include "ode.hpp"
15
16namespace mfem
17{
18
19std::string ODESolver::ExplicitTypes =
20 "\n\tExplicit solver: \n\t"
21 " RK : 1 - Forward Euler, 2 - RK2(0.5), 3 - RK3 SSP, 4 - RK4, 6 - RK6,\n\t"
22 " AB : 11 - AB1, 12 - AB2, 13 - AB3, 14 - AB4, 15 - AB5\n";
23
24std::string ODESolver::ImplicitTypes =
25 "\n\tImplicit solver: \n\t"
26 " (L-Stab): 21 - Backward Euler, 22 - SDIRK23(2), 23 - SDIRK33,\n\t"
27 " (A-Stab): 32 - Implicit Midpoint, 33 - SDIRK23, 34 - SDIRK34,\n\t"
28 " GA : 40 -- 50 - Generalized-alpha,\n\t"
29 " AM : 51 - AM1, 52 - AM2, 53 - AM3, 54 - AM4\n";
30
33
34std::unique_ptr<ODESolver> ODESolver::Select(int ode_solver_type)
35{
36 if (ode_solver_type < 20)
37 {
38 return SelectExplicit(ode_solver_type);
39 }
40 else
41 {
42 return SelectImplicit(ode_solver_type);
43 }
44}
45
46std::unique_ptr<ODESolver> ODESolver::SelectExplicit(int ode_solver_type)
47{
48 using ode_ptr = std::unique_ptr<ODESolver>;
49 switch (ode_solver_type)
50 {
51 // Explicit RK methods
52 case 1: return ode_ptr(new ForwardEulerSolver);
53 case 2: return ode_ptr(new RK2Solver(0.5)); // midpoint method
54 case 3: return ode_ptr(new RK3SSPSolver);
55 case 4: return ode_ptr(new RK4Solver);
56 case 6: return ode_ptr(new RK6Solver);
57
58 // Explicit AB methods
59 case 11: return ode_ptr(new AB1Solver);
60 case 12: return ode_ptr(new AB2Solver);
61 case 13: return ode_ptr(new AB3Solver);
62 case 14: return ode_ptr(new AB4Solver);
63 case 15: return ode_ptr(new AB5Solver);
64
65 default:
66 MFEM_ABORT("Unknown ODE solver type: " << ode_solver_type);
67 }
68}
69
70std::unique_ptr<ODESolver> ODESolver::SelectImplicit(int ode_solver_type)
71{
72 using ode_ptr = std::unique_ptr<ODESolver>;
73 switch (ode_solver_type)
74 {
75 // Implicit L-stable methods
76 case 21: return ode_ptr(new BackwardEulerSolver);
77 case 22: return ode_ptr(new SDIRK23Solver(2));
78 case 23: return ode_ptr(new SDIRK33Solver);
79
80 // Implicit A-stable methods (not L-stable)
81 case 32: return ode_ptr(new ImplicitMidpointSolver);
82 case 33: return ode_ptr(new SDIRK23Solver);
83 case 34: return ode_ptr(new SDIRK34Solver);
84
85 // Implicit generalized alpha
86 case 40: return ode_ptr(new GeneralizedAlphaSolver(0.0));
87 case 41: return ode_ptr(new GeneralizedAlphaSolver(0.1));
88 case 42: return ode_ptr(new GeneralizedAlphaSolver(0.2));
89 case 43: return ode_ptr(new GeneralizedAlphaSolver(0.3));
90 case 44: return ode_ptr(new GeneralizedAlphaSolver(0.4));
91 case 45: return ode_ptr(new GeneralizedAlphaSolver(0.5));
92 case 46: return ode_ptr(new GeneralizedAlphaSolver(0.6));
93 case 47: return ode_ptr(new GeneralizedAlphaSolver(0.7));
94 case 48: return ode_ptr(new GeneralizedAlphaSolver(0.8));
95 case 49: return ode_ptr(new GeneralizedAlphaSolver(0.9));
96 case 50: return ode_ptr(new GeneralizedAlphaSolver(1.0));
97
98 // Implicit AM methods
99 case 51: return ode_ptr(new AM1Solver);
100 case 52: return ode_ptr(new AM2Solver);
101 case 53: return ode_ptr(new AM3Solver);
102 case 54: return ode_ptr(new AM4Solver);
103
104 default:
105 MFEM_ABORT("Unknown ODE solver type: " << ode_solver_type );
106 }
107}
108
109
111{
112 mem_type = m_t;
113 for (int i = 0; i < smax; i++)
114 {
115 idx[i] = smax - i - 1;
116 data[i].SetSize(vsize, mem_type);
117 }
118
119 ss = 0;
120}
121
123{
124 MFEM_ASSERT_INDEX_IN_RANGE(i,0,ss);
125 return data[idx[i]];
126}
127
129{
130 MFEM_ASSERT_INDEX_IN_RANGE(i,0,ss);
131 return data[idx[i]];
132}
133
134void ODEStateDataVector::Get(int i, Vector &vec) const
135{
136 MFEM_ASSERT_INDEX_IN_RANGE(i,0,ss);
137 vec = data[idx[i]];
138}
139
141{
142 MFEM_ASSERT_INDEX_IN_RANGE(i,0,smax);
143 data[idx[i]] = state;
144}
145
147{
148 ShiftStages();
149 data[idx[0]] = state;
150 Increment();
151}
152
153void ODEStateDataVector::Print(std::ostream &os) const
154{
155 os << ss <<"/" <<smax<<std::endl;
156 idx.Print(os);
157 for (int i = 0; i < ss; i++) { data[idx[i]].Print(os); }
158}
159
160
162{
163 this->f = &f_;
165}
166
172
174{
175 f->SetTime(t);
176 f->Mult(x, dxdt);
177 x.Add(dt, dxdt);
178 t += dt;
179}
180
181
183{
184 ODESolver::Init(f_);
185 int n = f->Width();
186 dxdt.SetSize(n, mem_type);
187 x1.SetSize(n, mem_type);
188}
189
191{
192 // 0 |
193 // a | a
194 // ---+--------
195 // | 1-b b b = 1/(2a)
196
197 const real_t b = 0.5/a;
198
199 f->SetTime(t);
200 f->Mult(x, dxdt);
201 add(x, (1. - b)*dt, dxdt, x1);
202 x.Add(a*dt, dxdt);
203
204 f->SetTime(t + a*dt);
205 f->Mult(x, dxdt);
206 add(x1, b*dt, dxdt, x);
207 t += dt;
208}
209
210
212{
213 ODESolver::Init(f_);
214 int n = f->Width();
215 y.SetSize(n, mem_type);
216 k.SetSize(n, mem_type);
217}
218
220{
221 // x0 = x, t0 = t, k0 = dt*f(t0, x0)
222 f->SetTime(t);
223 f->Mult(x, k);
224
225 // x1 = x + k0, t1 = t + dt, k1 = dt*f(t1, x1)
226 add(x, dt, k, y);
227 f->SetTime(t + dt);
228 f->Mult(y, k);
229
230 // x2 = 3/4*x + 1/4*(x1 + k1), t2 = t + 1/2*dt, k2 = dt*f(t2, x2)
231 y.Add(dt, k);
232 add(3./4, x, 1./4, y, y);
233 f->SetTime(t + dt/2);
234 f->Mult(y, k);
235
236 // x3 = 1/3*x + 2/3*(x2 + k2), t3 = t + dt
237 y.Add(dt, k);
238 add(1./3, x, 2./3, y, x);
239 t += dt;
240}
241
242
244{
245 ODESolver::Init(f_);
246 int n = f->Width();
247 y.SetSize(n, mem_type);
248 k.SetSize(n, mem_type);
249 z.SetSize(n, mem_type);
250}
251
253{
254 // 0 |
255 // 1/2 | 1/2
256 // 1/2 | 0 1/2
257 // 1 | 0 0 1
258 // -----+-------------------
259 // | 1/6 1/3 1/3 1/6
260
261 f->SetTime(t);
262 f->Mult(x, k); // k1
263 add(x, dt/2, k, y);
264 add(x, dt/6, k, z);
265
266 f->SetTime(t + dt/2);
267 f->Mult(y, k); // k2
268 add(x, dt/2, k, y);
269 z.Add(dt/3, k);
270
271 f->Mult(y, k); // k3
272 add(x, dt, k, y);
273 z.Add(dt/3, k);
274
275 f->SetTime(t + dt);
276 f->Mult(y, k); // k4
277 add(z, dt/6, k, x);
278 t += dt;
279}
280
282 const real_t *c_)
283{
284 s = s_;
285 a = a_;
286 b = b_;
287 c = c_;
288 k = new Vector[s];
289}
290
292{
293 ODESolver::Init(f_);
294 int n = f->Width();
295 y.SetSize(n, mem_type);
296 for (int i = 0; i < s; i++)
297 {
298 k[i].SetSize(n, mem_type);
299 }
300}
301
303{
304 // 0 |
305 // c[0] | a[0]
306 // c[1] | a[1] a[2]
307 // ... | ...
308 // c[s-2] | ... a[s(s-1)/2-1]
309 // --------+---------------------
310 // | b[0] b[1] ... b[s-1]
311
312 f->SetTime(t);
313 f->Mult(x, k[0]);
314 for (int l = 0, i = 1; i < s; i++)
315 {
316 add(x, a[l++]*dt, k[0], y);
317 for (int j = 1; j < i; j++)
318 {
319 y.Add(a[l++]*dt, k[j]);
320 }
321
322 f->SetTime(t + c[i-1]*dt);
323 f->Mult(y, k[i]);
324 }
325 for (int i = 0; i < s; i++)
326 {
327 x.Add(b[i]*dt, k[i]);
328 }
329 t += dt;
330}
331
333{
334 delete [] k;
335}
336
337const real_t RK6Solver::a[] =
338{
339 .6e-1,
340 .1923996296296296296296296296296296296296e-1,
341 .7669337037037037037037037037037037037037e-1,
342 .35975e-1,
343 0.,
344 .107925,
345 1.318683415233148260919747276431735612861,
346 0.,
347 -5.042058063628562225427761634715637693344,
348 4.220674648395413964508014358283902080483,
349 -41.87259166432751461803757780644346812905,
350 0.,
351 159.4325621631374917700365669070346830453,
352 -122.1192135650100309202516203389242140663,
353 5.531743066200053768252631238332999150076,
354 -54.43015693531650433250642051294142461271,
355 0.,
356 207.0672513650184644273657173866509835987,
357 -158.6108137845899991828742424365058599469,
358 6.991816585950242321992597280791793907096,
359 -.1859723106220323397765171799549294623692e-1,
360 -54.66374178728197680241215648050386959351,
361 0.,
362 207.9528062553893734515824816699834244238,
363 -159.2889574744995071508959805871426654216,
364 7.018743740796944434698170760964252490817,
365 -.1833878590504572306472782005141738268361e-1,
366 -.5119484997882099077875432497245168395840e-3
367};
368const real_t RK6Solver::b[] =
369{
370 .3438957868357036009278820124728322386520e-1,
371 0.,
372 0.,
373 .2582624555633503404659558098586120858767,
374 .4209371189673537150642551514069801967032,
375 4.405396469669310170148836816197095664891,
376 -176.4831190242986576151740942499002125029,
377 172.3641334014150730294022582711902413315
378};
379const real_t RK6Solver::c[] =
380{
381 .6e-1,
382 .9593333333333333333333333333333333333333e-1,
383 .1439,
384 .4973,
385 .9725,
386 .9995,
387 1.,
388};
389
390const real_t RK8Solver::a[] =
391{
392 .5e-1,
393 -.69931640625e-2,
394 .1135556640625,
395 .399609375e-1,
396 0.,
397 .1198828125,
398 .3613975628004575124052940721184028345129,
399 0.,
400 -1.341524066700492771819987788202715834917,
401 1.370126503900035259414693716084313000404,
402 .490472027972027972027972027972027972028e-1,
403 0.,
404 0.,
405 .2350972042214404739862988335493427143122,
406 .180855592981356728810903963653454488485,
407 .6169289044289044289044289044289044289044e-1,
408 0.,
409 0.,
410 .1123656831464027662262557035130015442303,
411 -.3885046071451366767049048108111244567456e-1,
412 .1979188712522045855379188712522045855379e-1,
413 -1.767630240222326875735597119572145586714,
414 0.,
415 0.,
416 -62.5,
417 -6.061889377376669100821361459659331999758,
418 5.650823198222763138561298030600840174201,
419 65.62169641937623283799566054863063741227,
420 -1.180945066554970799825116282628297957882,
421 0.,
422 0.,
423 -41.50473441114320841606641502701994225874,
424 -4.434438319103725011225169229846100211776,
425 4.260408188586133024812193710744693240761,
426 43.75364022446171584987676829438379303004,
427 .787142548991231068744647504422630755086e-2,
428 -1.281405999441488405459510291182054246266,
429 0.,
430 0.,
431 -45.04713996013986630220754257136007322267,
432 -4.731362069449576477311464265491282810943,
433 4.514967016593807841185851584597240996214,
434 47.44909557172985134869022392235929015114,
435 .1059228297111661135687393955516542875228e-1,
436 -.5746842263844616254432318478286296232021e-2,
437 -1.724470134262485191756709817484481861731,
438 0.,
439 0.,
440 -60.92349008483054016518434619253765246063,
441 -5.95151837622239245520283276706185486829,
442 5.556523730698456235979791650843592496839,
443 63.98301198033305336837536378635995939281,
444 .1464202825041496159275921391759452676003e-1,
445 .6460408772358203603621865144977650714892e-1,
446 -.7930323169008878984024452548693373291447e-1,
447 -3.301622667747079016353994789790983625569,
448 0.,
449 0.,
450 -118.011272359752508566692330395789886851,
451 -10.14142238845611248642783916034510897595,
452 9.139311332232057923544012273556827000619,
453 123.3759428284042683684847180986501894364,
454 4.623244378874580474839807625067630924792,
455 -3.383277738068201923652550971536811240814,
456 4.527592100324618189451265339351129035325,
457 -5.828495485811622963193088019162985703755
458};
459const real_t RK8Solver::b[] =
460{
461 .4427989419007951074716746668098518862111e-1,
462 0.,
463 0.,
464 0.,
465 0.,
466 .3541049391724448744815552028733568354121,
467 .2479692154956437828667629415370663023884,
468 -15.69420203883808405099207034271191213468,
469 25.08406496555856261343930031237186278518,
470 -31.73836778626027646833156112007297739997,
471 22.93828327398878395231483560344797018313,
472 -.2361324633071542145259900641263517600737
473};
474const real_t RK8Solver::c[] =
475{
476 .5e-1,
477 .1065625,
478 .15984375,
479 .39,
480 .465,
481 .155,
482 .943,
483 .901802041735856958259707940678372149956,
484 .909,
485 .94,
486 1.,
487};
488
489
491 stages(s_), state(s_)
492{
493 a = a_;
494}
495
497{
498 ODESolver::Init(f_);
499 if (RKsolver) { RKsolver->Init(f_); }
500 state.SetSize(f->Width(), mem_type);
501 dt_ = -1.0;
502}
503
505{
506 CheckTimestep(dt);
507
508 if (state.Size() >= stages -1)
509 {
510 f->SetTime(t);
511 f->Mult(x, state[0]);
512 state.Increment();
513 for (int i = 0; i < stages; i++)
514 {
515 x.Add(a[i]*dt, state[i]);
516 }
517 t += dt;
518 }
519 else
520 {
521 f->Mult(x,state[0]);
522 RKsolver->Step(x,t,dt);
523 state.Increment();
524 }
525
526 state.ShiftStages();
527}
528
530{
531 if (dt_ < 0.0)
532 {
533 dt_ = dt;
534 return;
535 }
536 else if (fabs(dt-dt_) >10*std::numeric_limits<real_t>::epsilon())
537 {
538 state.Reset();
539 dt_ = dt;
540
541 if (print())
542 {
543 mfem::out << "WARNING:" << std::endl;
544 mfem::out << " - Time step changed" << std::endl;
545 mfem::out << " - Purging time stepping history" << std::endl;
546 mfem::out << " - Will run Runge-Kutta to rebuild history" << std::endl;
547 }
548 }
549}
550
551const real_t AB1Solver::a[] =
552{1.0};
553const real_t AB2Solver::a[] =
554{1.5,-0.5};
555const real_t AB3Solver::a[] =
556{23.0/12.0,-4.0/3.0, 5.0/12.0};
557const real_t AB4Solver::a[] =
558{55.0/24.0,-59.0/24.0, 37.0/24.0,-9.0/24.0};
559const real_t AB5Solver::a[] =
560{1901.0/720.0,-2774.0/720.0, 2616.0/720.0,-1274.0/720.0, 251.0/720.0};
561
562
564 stages(s_), state(s_)
565{
566 a = a_;
567}
568
570{
571 ODESolver::Init(f_);
572 if (RKsolver) { RKsolver->Init(f_); }
573 state.SetSize(f->Width(), mem_type);
574 dt_ = -1.0;
575}
576
578{
579 if (dt_ < 0.0)
580 {
581 dt_ = dt;
582 }
583 else if (fabs(dt-dt_) > 10*std::numeric_limits<real_t>::epsilon())
584 {
585 state.Reset();
586 dt_ = dt;
587
588 if (print())
589 {
590 mfem::out << "WARNING:" << std::endl;
591 mfem::out << " - Time step changed" << std::endl;
592 mfem::out << " - Purging time stepping history" << std::endl;
593 mfem::out << " - Will run Runge-Kutta to rebuild history" << std::endl;
594 }
595 }
596
597 if ((state.Size() == 0)&&(stages>1))
598 {
599 f->Mult(x,state[0]);
600 state.Increment();
601 }
602
603 if (state.Size() >= stages )
604 {
605 f->SetTime(t);
606 for (int i = 0; i < stages; i++)
607 {
608 x.Add(a[i+1]*dt, state[i]);
609 }
610 state.ShiftStages();
611 f->ImplicitSolve(a[0]*dt, x, state[0]);
612 x.Add(a[0]*dt, state[0]);
613 t += dt;
614 }
615 else
616 {
617 state.ShiftStages();
618 RKsolver->Step(x,t,dt);
619 f->Mult(x,state[0]);
620 state.Increment();
621 }
622}
623
624const real_t AM1Solver::a[] =
625{0.5, 0.5};
626const real_t AM2Solver::a[] =
627{5.0/12.0, 2.0/3.0, -1.0/12.0};
628const real_t AM3Solver::a[] =
629{3.0/8.0, 19.0/24.0,-5.0/24.0, 1.0/24.0};
630const real_t AM4Solver::a[] =
631{251.0/720.0,646.0/720.0,-264.0/720.0, 106.0/720.0, -19.0/720.0};
632
633
639
641{
642 f->SetTime(t + dt);
643 f->ImplicitSolve(dt, x, k); // solve for k: k = f(x + dt*k, t + dt)
644 x.Add(dt, k);
645 t += dt;
646}
647
648
654
656{
657 f->SetTime(t + dt/2);
658 f->ImplicitSolve(dt/2, x, k);
659 x.Add(dt, k);
660 t += dt;
661}
662
663
665{
666 if (gamma_opt == 0)
667 {
668 gamma = (3. - sqrt(3.))/6.; // not A-stable, order 3
669 }
670 else if (gamma_opt == 2)
671 {
672 gamma = (2. - sqrt(2.))/2.; // L-stable, order 2
673 }
674 else if (gamma_opt == 3)
675 {
676 gamma = (2. + sqrt(2.))/2.; // L-stable, order 2
677 }
678 else
679 {
680 gamma = (3. + sqrt(3.))/6.; // A-stable, order 3
681 }
682}
683
690
692{
693 // with a = gamma:
694 // a | a
695 // 1-a | 1-2a a
696 // ------+-----------
697 // | 1/2 1/2
698 // note: with gamma_opt=3, both solve are outside [t,t+dt] since a>1
699 f->SetTime(t + gamma*dt);
700 f->ImplicitSolve(gamma*dt, x, k);
701 add(x, (1.-2.*gamma)*dt, k, y); // y = x + (1-2*gamma)*dt*k
702 x.Add(dt/2, k);
703
704 f->SetTime(t + (1.-gamma)*dt);
705 f->ImplicitSolve(gamma*dt, y, k);
706 x.Add(dt/2, k);
707 t += dt;
708}
709
710
718
720{
721 // a | a
722 // 1/2 | 1/2-a a
723 // 1-a | 2a 1-4a a
724 // ------+--------------------
725 // | b 1-2b b
726 // note: two solves are outside [t,t+dt] since c1=a>1, c3=1-a<0
727 const real_t a = 1./sqrt(3.)*cos(M_PI/18.) + 0.5;
728 const real_t b = 1./(6.*(2.*a-1.)*(2.*a-1.));
729
730 f->SetTime(t + a*dt);
731 f->ImplicitSolve(a*dt, x, k);
732 add(x, (0.5-a)*dt, k, y);
733 add(x, (2.*a)*dt, k, z);
734 x.Add(b*dt, k);
735
736 f->SetTime(t + dt/2);
737 f->ImplicitSolve(a*dt, y, k);
738 z.Add((1.-4.*a)*dt, k);
739 x.Add((1.-2.*b)*dt, k);
740
741 f->SetTime(t + (1.-a)*dt);
742 f->ImplicitSolve(a*dt, z, k);
743 x.Add(b*dt, k);
744 t += dt;
745}
746
747
754
756{
757 // a | a
758 // c | c-a a
759 // 1 | b 1-a-b a
760 // -----+----------------
761 // | b 1-a-b a
762 const real_t a = 0.435866521508458999416019;
763 const real_t b = 1.20849664917601007033648;
764 const real_t c = 0.717933260754229499708010;
765
766 f->SetTime(t + a*dt);
767 f->ImplicitSolve(a*dt, x, k);
768 add(x, (c-a)*dt, k, y);
769 x.Add(b*dt, k);
770
771 f->SetTime(t + c*dt);
772 f->ImplicitSolve(a*dt, y, k);
773 x.Add((1.0-a-b)*dt, k);
774
775 f->SetTime(t + dt);
776 f->ImplicitSolve(a*dt, x, k);
777 x.Add(a*dt, k);
778 t += dt;
779}
780
787
789{
790 // 0 | 0 0
791 // 1 | 1/2 1/2
792 // ------+-----------
793 // | 1/2 1/2
794 f->SetTime(t);
795 f->Mult(x,k);
796 add(x, dt/2.0, k, y);
797 x.Add(dt/2.0, k);
798
799 f->SetTime(t + dt);
800 f->ImplicitSolve(dt/2.0, y, k);
801 x.Add(dt/2.0, k);
802 t += dt;
803}
804
812
814{
815 // 0 | 0 0 0
816 // 2a | a a 0
817 // 1 | 1-b-a b a
818 // ------+--------------------
819 // | 1-b-a b a
820 const real_t a = (2.0 - sqrt(2.0)) / 2.0;
821 const real_t b = (1.0 - 2.0*a) / (4.0*a);
822
823 f->SetTime(t);
824 f->Mult(x,k);
825 add(x, a*dt, k, y);
826 add(x, (1.0-b-a)*dt, k, z);
827 x.Add((1.0-b-a)*dt, k);
828
829 f->SetTime(t + (2.0*a)*dt);
830 f->ImplicitSolve(a*dt, y, k);
831 z.Add(b*dt, k);
832 x.Add(b*dt, k);
833
834 f->SetTime(t + dt);
835 f->ImplicitSolve(a*dt, z, k);
836 x.Add(a*dt, k);
837 t += dt;
838}
839
847
849{
850 // 0 | 0 0 0
851 // 2a | a a 0
852 // 1 | 1-b-a b a
853 // ------+----------------------------
854 // | 1-b_2-b_3 b_2 b_3
855 const real_t a = (3.0 + sqrt(3.0)) / 6.0;
856 const real_t b = (1.0 - 2.0*a) / (4.0*a);
857 const real_t b_2 = 1.0 / ( 12.0*a*(1.0 - 2.0*a) );
858 const real_t b_3 = (1.0 - 3.0*a) / ( 3.0*(1.0 - 2.0*a) );
859
860 f->SetTime(t);
861 f->Mult(x,k);
862 add(x, a*dt, k, y);
863 add(x, (1.0-b-a)*dt, k, z);
864 x.Add((1.0-b_2-b_3)*dt, k);
865
866 f->SetTime(t + (2.0*a)*dt);
867 f->ImplicitSolve(a*dt, y, k);
868 z.Add(b*dt, k);
869 x.Add(b_2*dt, k);
870
871 f->SetTime(t + dt);
872 f->ImplicitSolve(a*dt, z, k);
873 x.Add(b_3*dt, k);
874 t += dt;
875}
876
884
886{
887 rho_inf = (rho_inf > 1.0) ? 1.0 : rho_inf;
888 rho_inf = (rho_inf < 0.0) ? 0.0 : rho_inf;
889
890 // According to Jansen
891 alpha_m = 0.5*(3.0 - rho_inf)/(1.0 + rho_inf);
892 alpha_f = 1.0/(1.0 + rho_inf);
893 gamma = 0.5 + alpha_m - alpha_f;
894}
895
897{
898 os << "Generalized alpha time integrator:" << std::endl;
899 os << "alpha_m = " << alpha_m << std::endl;
900 os << "alpha_f = " << alpha_f << std::endl;
901 os << "gamma = " << gamma << std::endl;
902
903 if (gamma == 0.5 + alpha_m - alpha_f)
904 {
905 os<<"Second order"<<" and ";
906 }
907 else
908 {
909 os<<"First order"<<" and ";
910 }
911
912 if ((alpha_m >= alpha_f)&&(alpha_f >= 0.5))
913 {
914 os<<"Stable"<<std::endl;
915 }
916 else
917 {
918 os<<"Unstable"<<std::endl;
919 }
920}
921
922// This routine state[0] represents xdot
924{
925 if (state.Size() == 0)
926 {
927 f->Mult(x,state[0]);
928 state.Increment();
929 }
930
931 // Set y = x + alpha_f*(1.0 - (gamma/alpha_m))*dt*xdot
932 add(x, alpha_f*(1.0 - (gamma/alpha_m))*dt, state[0], y);
933
934 // Solve k = f(y + dt_eff*k)
935 real_t dt_eff = (gamma*alpha_f/alpha_m)*dt;
936 f->SetTime(t + alpha_f*dt);
937 f->ImplicitSolve(dt_eff, y, k);
938
939 // Update x and xdot
940 x.Add((1.0 - (gamma/alpha_m))*dt, state[0]);
941 x.Add( (gamma/alpha_m) *dt, k);
942
943 state[0] *= (1.0-(1.0/alpha_m));
944 state[0].Add((1.0/alpha_m),k);
945
946 t += dt;
947}
948
949
950void
952{
953 P_ = &P; F_ = &F;
954
955 dp_.SetSize(F_->Height());
956 dq_.SetSize(P_->Height());
957}
958
959void
961{
962 F_->SetTime(t);
963 F_->Mult(q,dp_);
964 p.Add(dt,dp_);
965
966 P_->Mult(p,dq_);
967 q.Add(dt,dq_);
968
969 t += dt;
970}
971
972void
974{
975 P_->Mult(p,dq_);
976 q.Add(0.5*dt,dq_);
977
978 F_->SetTime(t+0.5*dt);
979 F_->Mult(q,dp_);
980 p.Add(dt,dp_);
981
982 P_->Mult(p,dq_);
983 q.Add(0.5*dt,dq_);
984
985 t += dt;
986}
987
989 : order_(order)
990{
991 a_.SetSize(order);
992 b_.SetSize(order);
993
994 switch (order_)
995 {
996 case 1:
997 a_[0] = 1.0;
998 b_[0] = 1.0;
999 break;
1000 case 2:
1001 a_[0] = 0.5;
1002 a_[1] = 0.5;
1003 b_[0] = 0.0;
1004 b_[1] = 1.0;
1005 break;
1006 case 3:
1007 a_[0] = 2.0/3.0;
1008 a_[1] = -2.0/3.0;
1009 a_[2] = 1.0;
1010 b_[0] = 7.0/24.0;
1011 b_[1] = 0.75;
1012 b_[2] = -1.0/24.0;
1013 break;
1014 case 4:
1015 a_[0] = (2.0+pow(2.0,1.0/3.0)+pow(2.0,-1.0/3.0))/6.0;
1016 a_[1] = (1.0-pow(2.0,1.0/3.0)-pow(2.0,-1.0/3.0))/6.0;
1017 a_[2] = a_[1];
1018 a_[3] = a_[0];
1019 b_[0] = 0.0;
1020 b_[1] = 1.0/(2.0-pow(2.0,1.0/3.0));
1021 b_[2] = 1.0/(1.0-pow(2.0,2.0/3.0));
1022 b_[3] = b_[1];
1023 break;
1024 default:
1025 MFEM_ASSERT(false, "Unsupported order in SIAVSolver");
1026 };
1027}
1028
1029void
1031{
1032 for (int i=0; i<order_; i++)
1033 {
1034 if ( b_[i] != 0.0 )
1035 {
1036 F_->SetTime(t);
1037 if ( F_->isExplicit() )
1038 {
1039 F_->Mult(q, dp_);
1040 }
1041 else
1042 {
1043 F_->ImplicitSolve(b_[i] * dt, q, dp_);
1044 }
1045 p.Add(b_[i] * dt, dp_);
1046 }
1047
1048 P_->Mult(p, dq_);
1049 q.Add(a_[i] * dt, dq_);
1050
1051 t += a_[i] * dt;
1052 }
1053}
1054
1055std::string SecondOrderODESolver::Types =
1056 "ODE solver: \n\t"
1057 " [0--10] - GeneralizedAlpha(0.1 * s),\n\t"
1058 " 11 - Average Acceleration, 12 - Linear Acceleration\n\t"
1059 " 13 - CentralDifference, 14 - FoxGoodwin";
1060
1062{
1063 SecondOrderODESolver* ode_solver = NULL;
1064 switch (ode_solver_type)
1065 {
1066 // Implicit methods
1067 case 0: ode_solver = new GeneralizedAlpha2Solver(0.0); break;
1068 case 1: ode_solver = new GeneralizedAlpha2Solver(0.1); break;
1069 case 2: ode_solver = new GeneralizedAlpha2Solver(0.2); break;
1070 case 3: ode_solver = new GeneralizedAlpha2Solver(0.3); break;
1071 case 4: ode_solver = new GeneralizedAlpha2Solver(0.4); break;
1072 case 5: ode_solver = new GeneralizedAlpha2Solver(0.5); break;
1073 case 6: ode_solver = new GeneralizedAlpha2Solver(0.6); break;
1074 case 7: ode_solver = new GeneralizedAlpha2Solver(0.7); break;
1075 case 8: ode_solver = new GeneralizedAlpha2Solver(0.8); break;
1076 case 9: ode_solver = new GeneralizedAlpha2Solver(0.9); break;
1077 case 10: ode_solver = new GeneralizedAlpha2Solver(1.0); break;
1078
1079 case 11: ode_solver = new AverageAccelerationSolver(); break;
1080 case 12: ode_solver = new LinearAccelerationSolver(); break;
1081 case 13: ode_solver = new CentralDifferenceSolver(); break;
1082 case 14: ode_solver = new FoxGoodwinSolver(); break;
1083
1084 default:
1085 MFEM_ABORT("Unknown ODE solver type: " << ode_solver_type);
1086 }
1087 return ode_solver;
1088}
1089
1090// In this routine state[0] represents d2xdt2
1092 real_t &dt)
1093{
1094 x.Add(dt, dxdt);
1095
1096 f->SetTime(t + dt);
1097 f->ImplicitSolve(0.5*dt*dt, dt, x, dxdt, state[0]);
1098
1099 x .Add(0.5*dt*dt, state[0]);
1100 dxdt.Add(dt, state[0]);
1101 t += dt;
1102}
1103
1104// In this routine state[0] represents d2xdt2
1106 real_t &dt)
1107{
1108 x.Add(0.5*dt, dxdt);
1109
1110 f->SetTime(t + dt);
1111 f->ImplicitSolve(0.25*dt*dt, 0.5*dt, x, dxdt, state[0]);
1112
1113 x.Add(0.5*dt, dxdt);
1114 x.Add(0.5*dt*dt, state[0]);
1115 dxdt.Add(dt, state[0]);
1116 t += dt;
1117}
1118
1125
1127{
1128 os << "Newmark time integrator:" << std::endl;
1129 os << "beta = " << beta << std::endl;
1130 os << "gamma = " << gamma << std::endl;
1131
1132 if (gamma == 0.5)
1133 {
1134 os<<"Second order"<<" and ";
1135 }
1136 else
1137 {
1138 os<<"First order"<<" and ";
1139 }
1140
1141 if ((gamma >= 0.5) && (beta >= (gamma + 0.5)*(gamma + 0.5)/4))
1142 {
1143 os<<"A-Stable"<<std::endl;
1144 }
1145 else if ((gamma >= 0.5) && (beta >= 0.5*gamma))
1146 {
1147 os<<"Conditionally stable"<<std::endl;
1148 }
1149 else
1150 {
1151 os<<"Unstable"<<std::endl;
1152 }
1153}
1154
1155// In this routine state[0] represents d2xdt2
1157{
1158 real_t fac0 = 0.5 - beta;
1159 real_t fac2 = 1.0 - gamma;
1160 real_t fac3 = beta;
1161 real_t fac4 = gamma;
1162
1163 // In the first pass compute d2xdt2 directly from operator.
1164 if (state.Size() == 0)
1165 {
1166 if (no_mult)
1167 {
1168 MidPointStep(x, dxdt, t, dt);
1169 return;
1170 }
1171 else
1172 {
1173 f->Mult(x, dxdt, state[0]);
1174 }
1175 state.Increment();
1176 }
1177 f->SetTime(t + dt);
1178
1179 x.Add(dt, dxdt);
1180 x.Add(fac0*dt*dt, state[0]);
1181 dxdt.Add(fac2*dt, state[0]);
1182
1183 f->SetTime(t + dt);
1184 f->ImplicitSolve(fac3*dt*dt, fac4*dt, x, dxdt, state[0]);
1185
1186 x .Add(fac3*dt*dt, state[0]);
1187 dxdt.Add(fac4*dt, state[0]);
1188 t += dt;
1189}
1190
1198
1200{
1201 os << "Generalized alpha time integrator:" << std::endl;
1202 os << "alpha_m = " << alpha_m << std::endl;
1203 os << "alpha_f = " << alpha_f << std::endl;
1204 os << "beta = " << beta << std::endl;
1205 os << "gamma = " << gamma << std::endl;
1206
1207 if (gamma == 0.5 + alpha_m - alpha_f)
1208 {
1209 os<<"Second order"<<" and ";
1210 }
1211 else
1212 {
1213 os<<"First order"<<" and ";
1214 }
1215
1216 if ((alpha_m >= alpha_f)&&
1217 (alpha_f >= 0.5) &&
1218 (beta >= 0.25 + 0.5*(alpha_m - alpha_f)))
1219 {
1220 os<<"Stable"<<std::endl;
1221 }
1222 else
1223 {
1224 os<<"Unstable"<<std::endl;
1225 }
1226}
1227
1228// In this routine state[0] represents d2xdt2
1230 real_t &t, real_t &dt)
1231{
1232 real_t fac0 = (0.5 - (beta/alpha_m));
1233 real_t fac1 = alpha_f;
1234 real_t fac2 = alpha_f*(1.0 - (gamma/alpha_m));
1235 real_t fac3 = beta*alpha_f/alpha_m;
1236 real_t fac4 = gamma*alpha_f/alpha_m;
1237 real_t fac5 = alpha_m;
1238
1239 // In the first pass compute d2xdt2 directly from operator.
1240 if (state.Size() == 0)
1241 {
1242 if (no_mult)
1243 {
1244 MidPointStep(x, dxdt, t, dt);
1245 return;
1246 }
1247 else
1248 {
1249 f->Mult(x, dxdt, state[0]);
1250 }
1251 state.Increment();
1252 }
1253
1254 // Predict alpha levels
1255 add(dxdt, fac0*dt, state[0], va);
1256 add(x, fac1*dt, va, xa);
1257 add(dxdt, fac2*dt, state[0], va);
1258
1259 // Solve alpha levels
1260 f->SetTime(t + dt);
1261 f->ImplicitSolve(fac3*dt*dt, fac4*dt, xa, va, aa);
1262
1263 // Correct alpha levels
1264 xa.Add(fac3*dt*dt, aa);
1265 va.Add(fac4*dt, aa);
1266
1267 // Extrapolate
1268 x *= 1.0 - 1.0/fac1;
1269 x.Add (1.0/fac1, xa);
1270
1271 dxdt *= 1.0 - 1.0/fac1;
1272 dxdt.Add (1.0/fac1, va);
1273
1274 state[0] *= 1.0 - 1.0/fac5;
1275 state[0].Add (1.0/fac5, aa);
1276
1277 t += dt;
1278}
1279
1280}
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:496
void CheckTimestep(real_t dt)
Definition ode.cpp:529
AdamsBashforthSolver(int s_, const real_t *a_)
Definition ode.cpp:490
void Step(Vector &x, real_t &t, real_t &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Definition ode.cpp:504
std::unique_ptr< ODESolver > RKsolver
Definition ode.hpp:487
std::unique_ptr< ODESolver > RKsolver
Definition ode.hpp:570
void Step(Vector &x, real_t &t, real_t &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Definition ode.cpp:577
AdamsMoultonSolver(int s_, const real_t *a_)
Definition ode.cpp:563
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:569
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:758
void Print(std::ostream &out=mfem::out, int width=4) const
Prints array to stream with width elements per row.
Definition array.cpp:23
The classical midpoint method.
Definition ode.hpp:880
Backward Euler ODE solver. L-stable.
Definition ode.hpp:337
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:634
void Step(Vector &x, real_t &t, real_t &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Definition ode.cpp:640
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:805
void Step(Vector &x, real_t &t, real_t &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Definition ode.cpp:813
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:840
void Step(Vector &x, real_t &t, real_t &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Definition ode.cpp:848
void Step(Vector &x, real_t &t, real_t &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Definition ode.cpp:302
virtual ~ExplicitRKSolver()
Definition ode.cpp:332
ExplicitRKSolver(int s_, const real_t *a_, const real_t *b_, const real_t *c_)
Definition ode.cpp:281
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:291
The classical forward Euler method.
Definition ode.hpp:226
void Step(Vector &x, real_t &t, real_t &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Definition ode.cpp:173
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:167
void PrintProperties(std::ostream &os=mfem::out)
Definition ode.cpp:1199
void Init(SecondOrderTimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:1191
void Step(Vector &x, Vector &dxdt, real_t &t, real_t &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Definition ode.cpp:1229
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:877
void SetRhoInf(real_t rho_inf)
Definition ode.cpp:885
void Step(Vector &x, real_t &t, real_t &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Definition ode.cpp:923
void PrintProperties(std::ostream &os=mfem::out)
Definition ode.cpp:896
Implicit midpoint method. A-stable, not L-stable.
Definition ode.hpp:350
void Step(Vector &x, real_t &t, real_t &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Definition ode.cpp:655
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:649
void Step(Vector &x, Vector &dxdt, real_t &t, real_t &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Definition ode.cpp:1156
void PrintProperties(std::ostream &os=mfem::out)
Definition ode.cpp:1126
TimeDependentOperator * f
Pointer to the associated TimeDependentOperator.
Definition ode.hpp:113
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:161
static MFEM_EXPORT std::unique_ptr< ODESolver > SelectImplicit(const int ode_solver_type)
Definition ode.cpp:70
static MFEM_EXPORT std::string Types
Definition ode.hpp:187
static MFEM_EXPORT std::unique_ptr< ODESolver > Select(const int ode_solver_type)
Definition ode.cpp:34
MemoryType mem_type
Definition ode.hpp:114
static MFEM_EXPORT std::string ImplicitTypes
Definition ode.hpp:186
static MFEM_EXPORT std::string ExplicitTypes
Definition ode.hpp:185
static MFEM_EXPORT std::unique_ptr< ODESolver > SelectExplicit(const int ode_solver_type)
Definition ode.cpp:46
void Increment()
Increment the stage counter.
Definition ode.hpp:80
const Vector & Get(int i) const override
Get the ith state vector.
Definition ode.cpp:122
void Set(int i, Vector &state) override
Set the ith state vector.
Definition ode.cpp:140
void SetSize(int vsize, MemoryType mem_type)
Set the number of stages and the size of the vectors.
Definition ode.cpp:110
void Append(Vector &state) override
Add state vector and increment state size.
Definition ode.cpp:146
void ShiftStages()
Shift the stage counter for the next timestep.
Definition ode.hpp:74
void Print(std::ostream &os=mfem::out) const
Print state data.
Definition ode.cpp:153
int Size() const override
Get the current number of stored stages.
Definition ode.hpp:96
void Reset()
Reset the stage counter.
Definition ode.hpp:83
Abstract operator.
Definition operator.hpp:25
virtual MemoryClass GetMemoryClass() const
Return the MemoryClass preferred by the Operator.
Definition operator.hpp:86
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition operator.hpp:72
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:182
void Step(Vector &x, real_t &t, real_t &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Definition ode.cpp:190
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition ode.hpp:259
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:211
void Step(Vector &x, real_t &t, real_t &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Definition ode.cpp:219
The classical explicit forth-order Runge-Kutta method, RK4.
Definition ode.hpp:272
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:243
void Step(Vector &x, real_t &t, real_t &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Definition ode.cpp:252
SDIRK23Solver(int gamma_opt=1)
Definition ode.cpp:664
void Step(Vector &x, real_t &t, real_t &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Definition ode.cpp:691
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:684
void Step(Vector &x, real_t &t, real_t &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Definition ode.cpp:755
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:748
void Step(Vector &x, real_t &t, real_t &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Definition ode.cpp:719
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:711
void Step(Vector &q, Vector &p, real_t &t, real_t &dt) override
Definition ode.cpp:960
void Step(Vector &q, Vector &p, real_t &t, real_t &dt) override
Definition ode.cpp:973
TimeDependentOperator * F_
Definition ode.hpp:664
Vector dq_
Definition ode.hpp:668
Vector dp_
Definition ode.hpp:667
virtual void Init(Operator &P, TimeDependentOperator &F)
Definition ode.cpp:951
Operator * P_
Definition ode.hpp:665
void Step(Vector &q, Vector &p, real_t &t, real_t &dt) override
Definition ode.cpp:1030
SIAVSolver(int order)
Definition ode.cpp:988
Abstract class for solving systems of ODEs: d2x/dt2 = f(x,dx/dt,t)
Definition ode.hpp:705
SecondOrderTimeDependentOperator * f
Pointer to the associated TimeDependentOperator.
Definition ode.hpp:708
void MidPointStep(Vector &x, Vector &dxdt, real_t &t, real_t &dt)
Definition ode.cpp:1105
ODEStateDataVector state
Definition ode.hpp:710
static MFEM_EXPORT std::string Types
Help info for SecondOrderODESolver options.
Definition ode.hpp:796
void EulerStep(Vector &x, Vector &dxdt, real_t &t, real_t &dt)
Definition ode.cpp:1091
static MFEM_EXPORT SecondOrderODESolver * Select(const int ode_solver_type)
Function selecting the desired SecondOrderODESolver.
Definition ode.cpp:1061
virtual void Init(SecondOrderTimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:1119
Base abstract class for second order time dependent operators.
Definition operator.hpp:732
virtual void Mult(const Vector &x, const Vector &dxdt, Vector &y) const
Perform the action of the operator: y = k = f(x,@ dxdt, t), where k solves the algebraic equation F(x...
Definition operator.cpp:352
virtual void ImplicitSolve(const real_t fac0, const real_t fac1, const Vector &x, const Vector &dxdt, Vector &k)
Solve the equation: k = f(x + fac0 k, dxdt + fac1 k, t), for the unknown k at the current time t.
Definition operator.cpp:359
Base abstract class for first order time dependent operators.
Definition operator.hpp:332
bool isExplicit() const
True if type is EXPLICIT.
Definition operator.hpp:397
virtual void ImplicitSolve(const real_t gamma, const Vector &u, Vector &k)
Solve for the unknown k, at the current time t, the following equation: F(u + gamma k,...
Definition operator.cpp:298
void Mult(const Vector &u, Vector &k) const override
Perform the action of the operator (u,t) -> k(u,t) where t is the current time set by SetTime() and k...
Definition operator.cpp:293
virtual void SetTime(const real_t t_)
Set the current time.
Definition operator.hpp:394
void Step(Vector &x, real_t &t, real_t &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Definition ode.cpp:788
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:781
Vector data type.
Definition vector.hpp:82
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:558
Vector & Add(const real_t a, const Vector &Va)
(*this) += a * Va
Definition vector.cpp:322
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition globals.hpp:66
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition vector.cpp:391
MemoryType GetMemoryType(MemoryClass mc)
Return a suitable MemoryType for a given MemoryClass.
float real_t
Definition config.hpp:43
MemoryType
Memory types supported by MFEM.
real_t p(const Vector &x, real_t t)