MFEM v4.9.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
31std::string ODESolver::IMEXTypes =
32 "\n\tIMEX solver: \n\t"
33 " (L-Stab): 61 - Forward Backward Euler, 62 - IMEXRK2(2,2,2),\n\t"
34 " 63 - IMEXRK2(2,3,2), 64 - IMEX_DIRK_RK3\n";
35
39
40std::unique_ptr<ODESolver> ODESolver::Select(int ode_solver_type)
41{
42 if (ode_solver_type < 20)
43 {
44 return SelectExplicit(ode_solver_type);
45 }
46 else
47 {
48 return SelectImplicit(ode_solver_type);
49 }
50}
51
52std::unique_ptr<ODESolver> ODESolver::SelectExplicit(int ode_solver_type)
53{
54 using ode_ptr = std::unique_ptr<ODESolver>;
55 switch (ode_solver_type)
56 {
57 // Explicit RK methods
58 case 1: return ode_ptr(new ForwardEulerSolver);
59 case 2: return ode_ptr(new RK2Solver(0.5)); // midpoint method
60 case 3: return ode_ptr(new RK3SSPSolver);
61 case 4: return ode_ptr(new RK4Solver);
62 case 6: return ode_ptr(new RK6Solver);
63
64 // Explicit AB methods
65 case 11: return ode_ptr(new AB1Solver);
66 case 12: return ode_ptr(new AB2Solver);
67 case 13: return ode_ptr(new AB3Solver);
68 case 14: return ode_ptr(new AB4Solver);
69 case 15: return ode_ptr(new AB5Solver);
70
71 default:
72 MFEM_ABORT("Unknown ODE solver type: " << ode_solver_type);
73 }
74}
75
76std::unique_ptr<ODESolver> ODESolver::SelectImplicit(int ode_solver_type)
77{
78 using ode_ptr = std::unique_ptr<ODESolver>;
79 switch (ode_solver_type)
80 {
81 // Implicit L-stable methods
82 case 21: return ode_ptr(new BackwardEulerSolver);
83 case 22: return ode_ptr(new SDIRK23Solver(2));
84 case 23: return ode_ptr(new SDIRK33Solver);
85
86 // Implicit A-stable methods (not L-stable)
87 case 32: return ode_ptr(new ImplicitMidpointSolver);
88 case 33: return ode_ptr(new SDIRK23Solver);
89 case 34: return ode_ptr(new SDIRK34Solver);
90
91 // Implicit generalized alpha
92 case 40: return ode_ptr(new GeneralizedAlphaSolver(0.0));
93 case 41: return ode_ptr(new GeneralizedAlphaSolver(0.1));
94 case 42: return ode_ptr(new GeneralizedAlphaSolver(0.2));
95 case 43: return ode_ptr(new GeneralizedAlphaSolver(0.3));
96 case 44: return ode_ptr(new GeneralizedAlphaSolver(0.4));
97 case 45: return ode_ptr(new GeneralizedAlphaSolver(0.5));
98 case 46: return ode_ptr(new GeneralizedAlphaSolver(0.6));
99 case 47: return ode_ptr(new GeneralizedAlphaSolver(0.7));
100 case 48: return ode_ptr(new GeneralizedAlphaSolver(0.8));
101 case 49: return ode_ptr(new GeneralizedAlphaSolver(0.9));
102 case 50: return ode_ptr(new GeneralizedAlphaSolver(1.0));
103
104 // Implicit AM methods
105 case 51: return ode_ptr(new AM1Solver);
106 case 52: return ode_ptr(new AM2Solver);
107 case 53: return ode_ptr(new AM3Solver);
108 case 54: return ode_ptr(new AM4Solver);
109
110 default:
111 MFEM_ABORT("Unknown ODE solver type: " << ode_solver_type );
112 }
113}
114
115std::unique_ptr<ODESolver> ODESolver::SelectIMEX(const int ode_solver_type)
116{
117 using ode_ptr = std::unique_ptr<ODESolver>;
118 switch (ode_solver_type)
119 {
120 // L-stable IMEX methods
121 case 61: return ode_ptr(new IMEXExpImplEuler);
122 case 62: return ode_ptr(new IMEXRK2);
123 case 63: return ode_ptr(new IMEXRK2_3StageExplicit);
124 case 64: return ode_ptr(new IMEX_DIRK_RK3);
125
126 default: MFEM_ABORT("Unknown ODE solver type: " << ode_solver_type );
127 }
128}
129
131{
132 mem_type = m_t;
133 for (int i = 0; i < smax; i++)
134 {
135 idx[i] = smax - i - 1;
136 data[i].SetSize(vsize, mem_type);
137 }
138
139 ss = 0;
140}
141
143{
144 MFEM_ASSERT_INDEX_IN_RANGE(i,0,ss);
145 return data[idx[i]];
146}
147
149{
150 MFEM_ASSERT_INDEX_IN_RANGE(i,0,ss);
151 return data[idx[i]];
152}
153
154void ODEStateDataVector::Get(int i, Vector &vec) const
155{
156 MFEM_ASSERT_INDEX_IN_RANGE(i,0,ss);
157 vec = data[idx[i]];
158}
159
161{
162 MFEM_ASSERT_INDEX_IN_RANGE(i,0,smax);
163 data[idx[i]] = state;
164}
165
167{
168 ShiftStages();
169 data[idx[0]] = state;
170 Increment();
171}
172
173void ODEStateDataVector::Print(std::ostream &os) const
174{
175 os << ss <<"/" <<smax<<std::endl;
176 idx.Print(os);
177 for (int i = 0; i < ss; i++) { data[idx[i]].Print(os); }
178}
179
180
182{
183 this->f = &f_;
185}
186
192
194{
195 f->SetTime(t);
196 f->Mult(x, dxdt);
197 x.Add(dt, dxdt);
198 t += dt;
199}
200
201
203{
204 ODESolver::Init(f_);
205 int n = f->Width();
206 dxdt.SetSize(n, mem_type);
207 x1.SetSize(n, mem_type);
208}
209
211{
212 // 0 |
213 // a | a
214 // ---+--------
215 // | 1-b b b = 1/(2a)
216
217 const real_t b = 0.5/a;
218
219 f->SetTime(t);
220 f->Mult(x, dxdt);
221 add(x, (1. - b)*dt, dxdt, x1);
222 x.Add(a*dt, dxdt);
223
224 f->SetTime(t + a*dt);
225 f->Mult(x, dxdt);
226 add(x1, b*dt, dxdt, x);
227 t += dt;
228}
229
230
232{
233 ODESolver::Init(f_);
234 int n = f->Width();
235 y.SetSize(n, mem_type);
236 k.SetSize(n, mem_type);
237}
238
240{
241 // x0 = x, t0 = t, k0 = dt*f(t0, x0)
242 f->SetTime(t);
243 f->Mult(x, k);
244
245 // x1 = x + k0, t1 = t + dt, k1 = dt*f(t1, x1)
246 add(x, dt, k, y);
247 f->SetTime(t + dt);
248 f->Mult(y, k);
249
250 // x2 = 3/4*x + 1/4*(x1 + k1), t2 = t + 1/2*dt, k2 = dt*f(t2, x2)
251 y.Add(dt, k);
252 add(3./4, x, 1./4, y, y);
253 f->SetTime(t + dt/2);
254 f->Mult(y, k);
255
256 // x3 = 1/3*x + 2/3*(x2 + k2), t3 = t + dt
257 y.Add(dt, k);
258 add(1./3, x, 2./3, y, x);
259 t += dt;
260}
261
262
264{
265 ODESolver::Init(f_);
266 int n = f->Width();
267 y.SetSize(n, mem_type);
268 k.SetSize(n, mem_type);
269 z.SetSize(n, mem_type);
270}
271
273{
274 // 0 |
275 // 1/2 | 1/2
276 // 1/2 | 0 1/2
277 // 1 | 0 0 1
278 // -----+-------------------
279 // | 1/6 1/3 1/3 1/6
280
281 f->SetTime(t);
282 f->Mult(x, k); // k1
283 add(x, dt/2, k, y);
284 add(x, dt/6, k, z);
285
286 f->SetTime(t + dt/2);
287 f->Mult(y, k); // k2
288 add(x, dt/2, k, y);
289 z.Add(dt/3, k);
290
291 f->Mult(y, k); // k3
292 add(x, dt, k, y);
293 z.Add(dt/3, k);
294
295 f->SetTime(t + dt);
296 f->Mult(y, k); // k4
297 add(z, dt/6, k, x);
298 t += dt;
299}
300
302 const real_t *c_)
303{
304 s = s_;
305 a = a_;
306 b = b_;
307 c = c_;
308 k = new Vector[s];
309}
310
312{
313 ODESolver::Init(f_);
314 int n = f->Width();
315 y.SetSize(n, mem_type);
316 for (int i = 0; i < s; i++)
317 {
318 k[i].SetSize(n, mem_type);
319 }
320}
321
323{
324 // 0 |
325 // c[0] | a[0]
326 // c[1] | a[1] a[2]
327 // ... | ...
328 // c[s-2] | ... a[s(s-1)/2-1]
329 // --------+---------------------
330 // | b[0] b[1] ... b[s-1]
331
332 f->SetTime(t);
333 f->Mult(x, k[0]);
334 for (int l = 0, i = 1; i < s; i++)
335 {
336 add(x, a[l++]*dt, k[0], y);
337 for (int j = 1; j < i; j++)
338 {
339 y.Add(a[l++]*dt, k[j]);
340 }
341
342 f->SetTime(t + c[i-1]*dt);
343 f->Mult(y, k[i]);
344 }
345 for (int i = 0; i < s; i++)
346 {
347 x.Add(b[i]*dt, k[i]);
348 }
349 t += dt;
350}
351
353{
354 delete [] k;
355}
356
357const real_t RK6Solver::a[] =
358{
359 .6e-1,
360 .1923996296296296296296296296296296296296e-1,
361 .7669337037037037037037037037037037037037e-1,
362 .35975e-1,
363 0.,
364 .107925,
365 1.318683415233148260919747276431735612861,
366 0.,
367 -5.042058063628562225427761634715637693344,
368 4.220674648395413964508014358283902080483,
369 -41.87259166432751461803757780644346812905,
370 0.,
371 159.4325621631374917700365669070346830453,
372 -122.1192135650100309202516203389242140663,
373 5.531743066200053768252631238332999150076,
374 -54.43015693531650433250642051294142461271,
375 0.,
376 207.0672513650184644273657173866509835987,
377 -158.6108137845899991828742424365058599469,
378 6.991816585950242321992597280791793907096,
379 -.1859723106220323397765171799549294623692e-1,
380 -54.66374178728197680241215648050386959351,
381 0.,
382 207.9528062553893734515824816699834244238,
383 -159.2889574744995071508959805871426654216,
384 7.018743740796944434698170760964252490817,
385 -.1833878590504572306472782005141738268361e-1,
386 -.5119484997882099077875432497245168395840e-3
387};
388const real_t RK6Solver::b[] =
389{
390 .3438957868357036009278820124728322386520e-1,
391 0.,
392 0.,
393 .2582624555633503404659558098586120858767,
394 .4209371189673537150642551514069801967032,
395 4.405396469669310170148836816197095664891,
396 -176.4831190242986576151740942499002125029,
397 172.3641334014150730294022582711902413315
398};
399const real_t RK6Solver::c[] =
400{
401 .6e-1,
402 .9593333333333333333333333333333333333333e-1,
403 .1439,
404 .4973,
405 .9725,
406 .9995,
407 1.,
408};
409
410const real_t RK8Solver::a[] =
411{
412 .5e-1,
413 -.69931640625e-2,
414 .1135556640625,
415 .399609375e-1,
416 0.,
417 .1198828125,
418 .3613975628004575124052940721184028345129,
419 0.,
420 -1.341524066700492771819987788202715834917,
421 1.370126503900035259414693716084313000404,
422 .490472027972027972027972027972027972028e-1,
423 0.,
424 0.,
425 .2350972042214404739862988335493427143122,
426 .180855592981356728810903963653454488485,
427 .6169289044289044289044289044289044289044e-1,
428 0.,
429 0.,
430 .1123656831464027662262557035130015442303,
431 -.3885046071451366767049048108111244567456e-1,
432 .1979188712522045855379188712522045855379e-1,
433 -1.767630240222326875735597119572145586714,
434 0.,
435 0.,
436 -62.5,
437 -6.061889377376669100821361459659331999758,
438 5.650823198222763138561298030600840174201,
439 65.62169641937623283799566054863063741227,
440 -1.180945066554970799825116282628297957882,
441 0.,
442 0.,
443 -41.50473441114320841606641502701994225874,
444 -4.434438319103725011225169229846100211776,
445 4.260408188586133024812193710744693240761,
446 43.75364022446171584987676829438379303004,
447 .787142548991231068744647504422630755086e-2,
448 -1.281405999441488405459510291182054246266,
449 0.,
450 0.,
451 -45.04713996013986630220754257136007322267,
452 -4.731362069449576477311464265491282810943,
453 4.514967016593807841185851584597240996214,
454 47.44909557172985134869022392235929015114,
455 .1059228297111661135687393955516542875228e-1,
456 -.5746842263844616254432318478286296232021e-2,
457 -1.724470134262485191756709817484481861731,
458 0.,
459 0.,
460 -60.92349008483054016518434619253765246063,
461 -5.95151837622239245520283276706185486829,
462 5.556523730698456235979791650843592496839,
463 63.98301198033305336837536378635995939281,
464 .1464202825041496159275921391759452676003e-1,
465 .6460408772358203603621865144977650714892e-1,
466 -.7930323169008878984024452548693373291447e-1,
467 -3.301622667747079016353994789790983625569,
468 0.,
469 0.,
470 -118.011272359752508566692330395789886851,
471 -10.14142238845611248642783916034510897595,
472 9.139311332232057923544012273556827000619,
473 123.3759428284042683684847180986501894364,
474 4.623244378874580474839807625067630924792,
475 -3.383277738068201923652550971536811240814,
476 4.527592100324618189451265339351129035325,
477 -5.828495485811622963193088019162985703755
478};
479const real_t RK8Solver::b[] =
480{
481 .4427989419007951074716746668098518862111e-1,
482 0.,
483 0.,
484 0.,
485 0.,
486 .3541049391724448744815552028733568354121,
487 .2479692154956437828667629415370663023884,
488 -15.69420203883808405099207034271191213468,
489 25.08406496555856261343930031237186278518,
490 -31.73836778626027646833156112007297739997,
491 22.93828327398878395231483560344797018313,
492 -.2361324633071542145259900641263517600737
493};
494const real_t RK8Solver::c[] =
495{
496 .5e-1,
497 .1065625,
498 .15984375,
499 .39,
500 .465,
501 .155,
502 .943,
503 .901802041735856958259707940678372149956,
504 .909,
505 .94,
506 1.,
507};
508
509
511 stages(s_), state(s_)
512{
513 a = a_;
514}
515
517{
518 ODESolver::Init(f_);
519 if (RKsolver) { RKsolver->Init(f_); }
520 state.SetSize(f->Width(), mem_type);
521 dt_ = -1.0;
522}
523
525{
526 CheckTimestep(dt);
527
528 if (state.Size() >= stages -1)
529 {
530 f->SetTime(t);
531 f->Mult(x, state[0]);
532 state.Increment();
533 for (int i = 0; i < stages; i++)
534 {
535 x.Add(a[i]*dt, state[i]);
536 }
537 t += dt;
538 }
539 else
540 {
541 f->Mult(x,state[0]);
542 RKsolver->Step(x,t,dt);
543 state.Increment();
544 }
545
546 state.ShiftStages();
547}
548
550{
551 if (dt_ < 0.0)
552 {
553 dt_ = dt;
554 return;
555 }
556 else if (fabs(dt-dt_) >10*std::numeric_limits<real_t>::epsilon())
557 {
558 state.Reset();
559 dt_ = dt;
560
561 if (print())
562 {
563 mfem::out << "WARNING:" << std::endl;
564 mfem::out << " - Time step changed" << std::endl;
565 mfem::out << " - Purging time stepping history" << std::endl;
566 mfem::out << " - Will run Runge-Kutta to rebuild history" << std::endl;
567 }
568 }
569}
570
571const real_t AB1Solver::a[] =
572{1.0};
573const real_t AB2Solver::a[] =
574{1.5,-0.5};
575const real_t AB3Solver::a[] =
576{23.0/12.0,-4.0/3.0, 5.0/12.0};
577const real_t AB4Solver::a[] =
578{55.0/24.0,-59.0/24.0, 37.0/24.0,-9.0/24.0};
579const real_t AB5Solver::a[] =
580{1901.0/720.0,-2774.0/720.0, 2616.0/720.0,-1274.0/720.0, 251.0/720.0};
581
582
584 stages(s_), state(s_)
585{
586 a = a_;
587}
588
590{
591 ODESolver::Init(f_);
592 if (RKsolver) { RKsolver->Init(f_); }
593 state.SetSize(f->Width(), mem_type);
594 dt_ = -1.0;
595}
596
598{
599 if (dt_ < 0.0)
600 {
601 dt_ = dt;
602 }
603 else if (fabs(dt-dt_) > 10*std::numeric_limits<real_t>::epsilon())
604 {
605 state.Reset();
606 dt_ = dt;
607
608 if (print())
609 {
610 mfem::out << "WARNING:" << std::endl;
611 mfem::out << " - Time step changed" << std::endl;
612 mfem::out << " - Purging time stepping history" << std::endl;
613 mfem::out << " - Will run Runge-Kutta to rebuild history" << std::endl;
614 }
615 }
616
617 if ((state.Size() == 0)&&(stages>1))
618 {
619 f->Mult(x,state[0]);
620 state.Increment();
621 }
622
623 if (state.Size() >= stages )
624 {
625 f->SetTime(t);
626 for (int i = 0; i < stages; i++)
627 {
628 x.Add(a[i+1]*dt, state[i]);
629 }
630 state.ShiftStages();
631 f->ImplicitSolve(a[0]*dt, x, state[0]);
632 x.Add(a[0]*dt, state[0]);
633 t += dt;
634 }
635 else
636 {
637 state.ShiftStages();
638 RKsolver->Step(x,t,dt);
639 f->Mult(x,state[0]);
640 state.Increment();
641 }
642}
643
644const real_t AM1Solver::a[] =
645{0.5, 0.5};
646const real_t AM2Solver::a[] =
647{5.0/12.0, 2.0/3.0, -1.0/12.0};
648const real_t AM3Solver::a[] =
649{3.0/8.0, 19.0/24.0,-5.0/24.0, 1.0/24.0};
650const real_t AM4Solver::a[] =
651{251.0/720.0,646.0/720.0,-264.0/720.0, 106.0/720.0, -19.0/720.0};
652
653
659
661{
662 f->SetTime(t + dt);
663 f->ImplicitSolve(dt, x, k); // solve for k: k = f(x + dt*k, t + dt)
664 x.Add(dt, k);
665 t += dt;
666}
667
668
674
676{
677 f->SetTime(t + dt/2);
678 f->ImplicitSolve(dt/2, x, k);
679 x.Add(dt, k);
680 t += dt;
681}
682
683
685{
686 if (gamma_opt == 0)
687 {
688 gamma = (3. - sqrt(3.))/6.; // not A-stable, order 3
689 }
690 else if (gamma_opt == 2)
691 {
692 gamma = (2. - sqrt(2.))/2.; // L-stable, order 2
693 }
694 else if (gamma_opt == 3)
695 {
696 gamma = (2. + sqrt(2.))/2.; // L-stable, order 2
697 }
698 else
699 {
700 gamma = (3. + sqrt(3.))/6.; // A-stable, order 3
701 }
702}
703
710
712{
713 // with a = gamma:
714 // a | a
715 // 1-a | 1-2a a
716 // ------+-----------
717 // | 1/2 1/2
718 // note: with gamma_opt=3, both solve are outside [t,t+dt] since a>1
719 f->SetTime(t + gamma*dt);
720 f->ImplicitSolve(gamma*dt, x, k);
721 add(x, (1.-2.*gamma)*dt, k, y); // y = x + (1-2*gamma)*dt*k
722 x.Add(dt/2, k);
723
724 f->SetTime(t + (1.-gamma)*dt);
725 f->ImplicitSolve(gamma*dt, y, k);
726 x.Add(dt/2, k);
727 t += dt;
728}
729
730
738
740{
741 // a | a
742 // 1/2 | 1/2-a a
743 // 1-a | 2a 1-4a a
744 // ------+--------------------
745 // | b 1-2b b
746 // note: two solves are outside [t,t+dt] since c1=a>1, c3=1-a<0
747 const real_t a = 1./sqrt(3.)*cos(M_PI/18.) + 0.5;
748 const real_t b = 1./(6.*(2.*a-1.)*(2.*a-1.));
749
750 f->SetTime(t + a*dt);
751 f->ImplicitSolve(a*dt, x, k);
752 add(x, (0.5-a)*dt, k, y);
753 add(x, (2.*a)*dt, k, z);
754 x.Add(b*dt, k);
755
756 f->SetTime(t + dt/2);
757 f->ImplicitSolve(a*dt, y, k);
758 z.Add((1.-4.*a)*dt, k);
759 x.Add((1.-2.*b)*dt, k);
760
761 f->SetTime(t + (1.-a)*dt);
762 f->ImplicitSolve(a*dt, z, k);
763 x.Add(b*dt, k);
764 t += dt;
765}
766
767
774
776{
777 // a | a
778 // c | c-a a
779 // 1 | b 1-a-b a
780 // -----+----------------
781 // | b 1-a-b a
782 const real_t a = 0.435866521508458999416019;
783 const real_t b = 1.20849664917601007033648;
784 const real_t c = 0.717933260754229499708010;
785
786 f->SetTime(t + a*dt);
787 f->ImplicitSolve(a*dt, x, k);
788 add(x, (c-a)*dt, k, y);
789 x.Add(b*dt, k);
790
791 f->SetTime(t + c*dt);
792 f->ImplicitSolve(a*dt, y, k);
793 x.Add((1.0-a-b)*dt, k);
794
795 f->SetTime(t + dt);
796 f->ImplicitSolve(a*dt, x, k);
797 x.Add(a*dt, k);
798 t += dt;
799}
800
807
809{
810 // 0 | 0 0
811 // 1 | 1/2 1/2
812 // ------+-----------
813 // | 1/2 1/2
814 f->SetTime(t);
815 f->Mult(x,k);
816 add(x, dt/2.0, k, y);
817 x.Add(dt/2.0, k);
818
819 f->SetTime(t + dt);
820 f->ImplicitSolve(dt/2.0, y, k);
821 x.Add(dt/2.0, k);
822 t += dt;
823}
824
832
834{
835 // 0 | 0 0 0
836 // 2a | a a 0
837 // 1 | 1-b-a b a
838 // ------+--------------------
839 // | 1-b-a b a
840 const real_t a = (2.0 - sqrt(2.0)) / 2.0;
841 const real_t b = (1.0 - 2.0*a) / (4.0*a);
842
843 f->SetTime(t);
844 f->Mult(x,k);
845 add(x, a*dt, k, y);
846 add(x, (1.0-b-a)*dt, k, z);
847 x.Add((1.0-b-a)*dt, k);
848
849 f->SetTime(t + (2.0*a)*dt);
850 f->ImplicitSolve(a*dt, y, k);
851 z.Add(b*dt, k);
852 x.Add(b*dt, k);
853
854 f->SetTime(t + dt);
855 f->ImplicitSolve(a*dt, z, k);
856 x.Add(a*dt, k);
857 t += dt;
858}
859
867
869{
870 // 0 | 0 0 0
871 // 2a | a a 0
872 // 1 | 1-b-a b a
873 // ------+----------------------------
874 // | 1-b_2-b_3 b_2 b_3
875 const real_t a = (3.0 + sqrt(3.0)) / 6.0;
876 const real_t b = (1.0 - 2.0*a) / (4.0*a);
877 const real_t b_2 = 1.0 / ( 12.0*a*(1.0 - 2.0*a) );
878 const real_t b_3 = (1.0 - 3.0*a) / ( 3.0*(1.0 - 2.0*a) );
879
880 f->SetTime(t);
881 f->Mult(x,k);
882 add(x, a*dt, k, y);
883 add(x, (1.0-b-a)*dt, k, z);
884 x.Add((1.0-b_2-b_3)*dt, k);
885
886 f->SetTime(t + (2.0*a)*dt);
887 f->ImplicitSolve(a*dt, y, k);
888 z.Add(b*dt, k);
889 x.Add(b_2*dt, k);
890
891 f->SetTime(t + dt);
892 f->ImplicitSolve(a*dt, z, k);
893 x.Add(b_3*dt, k);
894 t += dt;
895}
896
904
906{
907 rho_inf = (rho_inf > 1.0) ? 1.0 : rho_inf;
908 rho_inf = (rho_inf < 0.0) ? 0.0 : rho_inf;
909
910 // According to Jansen
911 alpha_m = 0.5*(3.0 - rho_inf)/(1.0 + rho_inf);
912 alpha_f = 1.0/(1.0 + rho_inf);
913 gamma = 0.5 + alpha_m - alpha_f;
914}
915
917{
918 os << "Generalized alpha time integrator:" << std::endl;
919 os << "alpha_m = " << alpha_m << std::endl;
920 os << "alpha_f = " << alpha_f << std::endl;
921 os << "gamma = " << gamma << std::endl;
922
923 if (gamma == 0.5 + alpha_m - alpha_f)
924 {
925 os<<"Second order"<<" and ";
926 }
927 else
928 {
929 os<<"First order"<<" and ";
930 }
931
932 if ((alpha_m >= alpha_f)&&(alpha_f >= 0.5))
933 {
934 os<<"Stable"<<std::endl;
935 }
936 else
937 {
938 os<<"Unstable"<<std::endl;
939 }
940}
941
942// This routine state[0] represents xdot
944{
945 if (state.Size() == 0)
946 {
947 f->Mult(x,state[0]);
948 state.Increment();
949 }
950
951 // Set y = x + alpha_f*(1.0 - (gamma/alpha_m))*dt*xdot
952 add(x, alpha_f*(1.0 - (gamma/alpha_m))*dt, state[0], y);
953
954 // Solve k = f(y + dt_eff*k)
955 real_t dt_eff = (gamma*alpha_f/alpha_m)*dt;
956 f->SetTime(t + alpha_f*dt);
957 f->ImplicitSolve(dt_eff, y, k);
958
959 // Update x and xdot
960 x.Add((1.0 - (gamma/alpha_m))*dt, state[0]);
961 x.Add( (gamma/alpha_m) *dt, k);
962
963 state[0] *= (1.0-(1.0/alpha_m));
964 state[0].Add((1.0/alpha_m),k);
965
966 t += dt;
967}
968
969
970void
972{
973 P_ = &P; F_ = &F;
974
975 dp_.SetSize(F_->Height());
976 dq_.SetSize(P_->Height());
977}
978
979void
981{
982 F_->SetTime(t);
983 F_->Mult(q,dp_);
984 p.Add(dt,dp_);
985
986 P_->Mult(p,dq_);
987 q.Add(dt,dq_);
988
989 t += dt;
990}
991
992void
994{
995 P_->Mult(p,dq_);
996 q.Add(0.5*dt,dq_);
997
998 F_->SetTime(t+0.5*dt);
999 F_->Mult(q,dp_);
1000 p.Add(dt,dp_);
1001
1002 P_->Mult(p,dq_);
1003 q.Add(0.5*dt,dq_);
1004
1005 t += dt;
1006}
1007
1009 : order_(order)
1010{
1011 a_.SetSize(order);
1012 b_.SetSize(order);
1013
1014 switch (order_)
1015 {
1016 case 1:
1017 a_[0] = 1.0;
1018 b_[0] = 1.0;
1019 break;
1020 case 2:
1021 a_[0] = 0.5;
1022 a_[1] = 0.5;
1023 b_[0] = 0.0;
1024 b_[1] = 1.0;
1025 break;
1026 case 3:
1027 a_[0] = 2.0/3.0;
1028 a_[1] = -2.0/3.0;
1029 a_[2] = 1.0;
1030 b_[0] = 7.0/24.0;
1031 b_[1] = 0.75;
1032 b_[2] = -1.0/24.0;
1033 break;
1034 case 4:
1035 a_[0] = (2.0+pow(2.0,1.0/3.0)+pow(2.0,-1.0/3.0))/6.0;
1036 a_[1] = (1.0-pow(2.0,1.0/3.0)-pow(2.0,-1.0/3.0))/6.0;
1037 a_[2] = a_[1];
1038 a_[3] = a_[0];
1039 b_[0] = 0.0;
1040 b_[1] = 1.0/(2.0-pow(2.0,1.0/3.0));
1041 b_[2] = 1.0/(1.0-pow(2.0,2.0/3.0));
1042 b_[3] = b_[1];
1043 break;
1044 default:
1045 MFEM_ASSERT(false, "Unsupported order in SIAVSolver");
1046 };
1047}
1048
1049void
1051{
1052 for (int i=0; i<order_; i++)
1053 {
1054 if ( b_[i] != 0.0 )
1055 {
1056 F_->SetTime(t);
1057 if ( F_->isExplicit() )
1058 {
1059 F_->Mult(q, dp_);
1060 }
1061 else
1062 {
1063 F_->ImplicitSolve(b_[i] * dt, q, dp_);
1064 }
1065 p.Add(b_[i] * dt, dp_);
1066 }
1067
1068 P_->Mult(p, dq_);
1069 q.Add(a_[i] * dt, dq_);
1070
1071 t += a_[i] * dt;
1072 }
1073}
1074
1075std::string SecondOrderODESolver::Types =
1076 "ODE solver: \n\t"
1077 " [0--10] - GeneralizedAlpha(0.1 * s),\n\t"
1078 " 11 - Average Acceleration, 12 - Linear Acceleration\n\t"
1079 " 13 - CentralDifference, 14 - FoxGoodwin";
1080
1082{
1083 SecondOrderODESolver* ode_solver = NULL;
1084 switch (ode_solver_type)
1085 {
1086 // Implicit methods
1087 case 0: ode_solver = new GeneralizedAlpha2Solver(0.0); break;
1088 case 1: ode_solver = new GeneralizedAlpha2Solver(0.1); break;
1089 case 2: ode_solver = new GeneralizedAlpha2Solver(0.2); break;
1090 case 3: ode_solver = new GeneralizedAlpha2Solver(0.3); break;
1091 case 4: ode_solver = new GeneralizedAlpha2Solver(0.4); break;
1092 case 5: ode_solver = new GeneralizedAlpha2Solver(0.5); break;
1093 case 6: ode_solver = new GeneralizedAlpha2Solver(0.6); break;
1094 case 7: ode_solver = new GeneralizedAlpha2Solver(0.7); break;
1095 case 8: ode_solver = new GeneralizedAlpha2Solver(0.8); break;
1096 case 9: ode_solver = new GeneralizedAlpha2Solver(0.9); break;
1097 case 10: ode_solver = new GeneralizedAlpha2Solver(1.0); break;
1098
1099 case 11: ode_solver = new AverageAccelerationSolver(); break;
1100 case 12: ode_solver = new LinearAccelerationSolver(); break;
1101 case 13: ode_solver = new CentralDifferenceSolver(); break;
1102 case 14: ode_solver = new FoxGoodwinSolver(); break;
1103
1104 default:
1105 MFEM_ABORT("Unknown ODE solver type: " << ode_solver_type);
1106 }
1107 return ode_solver;
1108}
1109
1110// In this routine state[0] represents d2xdt2
1112 real_t &dt)
1113{
1114 x.Add(dt, dxdt);
1115
1116 f->SetTime(t + dt);
1117 f->ImplicitSolve(0.5*dt*dt, dt, x, dxdt, state[0]);
1118
1119 x .Add(0.5*dt*dt, state[0]);
1120 dxdt.Add(dt, state[0]);
1121 t += dt;
1122}
1123
1124// In this routine state[0] represents d2xdt2
1126 real_t &dt)
1127{
1128 x.Add(0.5*dt, dxdt);
1129
1130 f->SetTime(t + dt);
1131 f->ImplicitSolve(0.25*dt*dt, 0.5*dt, x, dxdt, state[0]);
1132
1133 x.Add(0.5*dt, dxdt);
1134 x.Add(0.5*dt*dt, state[0]);
1135 dxdt.Add(dt, state[0]);
1136 t += dt;
1137}
1138
1145
1147{
1148 os << "Newmark time integrator:" << std::endl;
1149 os << "beta = " << beta << std::endl;
1150 os << "gamma = " << gamma << std::endl;
1151
1152 if (gamma == 0.5)
1153 {
1154 os<<"Second order"<<" and ";
1155 }
1156 else
1157 {
1158 os<<"First order"<<" and ";
1159 }
1160
1161 if ((gamma >= 0.5) && (beta >= (gamma + 0.5)*(gamma + 0.5)/4))
1162 {
1163 os<<"A-Stable"<<std::endl;
1164 }
1165 else if ((gamma >= 0.5) && (beta >= 0.5*gamma))
1166 {
1167 os<<"Conditionally stable"<<std::endl;
1168 }
1169 else
1170 {
1171 os<<"Unstable"<<std::endl;
1172 }
1173}
1174
1175// In this routine state[0] represents d2xdt2
1177{
1178 real_t fac0 = 0.5 - beta;
1179 real_t fac2 = 1.0 - gamma;
1180 real_t fac3 = beta;
1181 real_t fac4 = gamma;
1182
1183 // In the first pass compute d2xdt2 directly from operator.
1184 if (state.Size() == 0)
1185 {
1186 if (no_mult)
1187 {
1188 MidPointStep(x, dxdt, t, dt);
1189 return;
1190 }
1191 else
1192 {
1193 f->Mult(x, dxdt, state[0]);
1194 }
1195 state.Increment();
1196 }
1197 f->SetTime(t + dt);
1198
1199 x.Add(dt, dxdt);
1200 x.Add(fac0*dt*dt, state[0]);
1201 dxdt.Add(fac2*dt, state[0]);
1202
1203 f->SetTime(t + dt);
1204 f->ImplicitSolve(fac3*dt*dt, fac4*dt, x, dxdt, state[0]);
1205
1206 x .Add(fac3*dt*dt, state[0]);
1207 dxdt.Add(fac4*dt, state[0]);
1208 t += dt;
1209}
1210
1218
1220{
1221 os << "Generalized alpha time integrator:" << std::endl;
1222 os << "alpha_m = " << alpha_m << std::endl;
1223 os << "alpha_f = " << alpha_f << std::endl;
1224 os << "beta = " << beta << std::endl;
1225 os << "gamma = " << gamma << std::endl;
1226
1227 if (gamma == 0.5 + alpha_m - alpha_f)
1228 {
1229 os<<"Second order"<<" and ";
1230 }
1231 else
1232 {
1233 os<<"First order"<<" and ";
1234 }
1235
1236 if ((alpha_m >= alpha_f)&&
1237 (alpha_f >= 0.5) &&
1238 (beta >= 0.25 + 0.5*(alpha_m - alpha_f)))
1239 {
1240 os<<"Stable"<<std::endl;
1241 }
1242 else
1243 {
1244 os<<"Unstable"<<std::endl;
1245 }
1246}
1247
1248// In this routine state[0] represents d2xdt2
1250 real_t &t, real_t &dt)
1251{
1252 real_t fac0 = (0.5 - (beta/alpha_m));
1253 real_t fac1 = alpha_f;
1254 real_t fac2 = alpha_f*(1.0 - (gamma/alpha_m));
1255 real_t fac3 = beta*alpha_f/alpha_m;
1256 real_t fac4 = gamma*alpha_f/alpha_m;
1257 real_t fac5 = alpha_m;
1258
1259 // In the first pass compute d2xdt2 directly from operator.
1260 if (state.Size() == 0)
1261 {
1262 if (no_mult)
1263 {
1264 MidPointStep(x, dxdt, t, dt);
1265 return;
1266 }
1267 else
1268 {
1269 f->Mult(x, dxdt, state[0]);
1270 }
1271 state.Increment();
1272 }
1273
1274 // Predict alpha levels
1275 add(dxdt, fac0*dt, state[0], va);
1276 add(x, fac1*dt, va, xa);
1277 add(dxdt, fac2*dt, state[0], va);
1278
1279 // Solve alpha levels
1280 f->SetTime(t + dt);
1281 f->ImplicitSolve(fac3*dt*dt, fac4*dt, xa, va, aa);
1282
1283 // Correct alpha levels
1284 xa.Add(fac3*dt*dt, aa);
1285 va.Add(fac4*dt, aa);
1286
1287 // Extrapolate
1288 x *= 1.0 - 1.0/fac1;
1289 x.Add (1.0/fac1, xa);
1290
1291 dxdt *= 1.0 - 1.0/fac1;
1292 dxdt.Add (1.0/fac1, va);
1293
1294 state[0] *= 1.0 - 1.0/fac5;
1295 state[0].Add (1.0/fac5, aa);
1296
1297 t += dt;
1298}
1299
1301{
1302 ODESolver::Init(f_);
1303 int n = f->Width();
1304 k1.SetSize(n, mem_type);
1305 k2.SetSize(n, mem_type);
1306}
1307
1309{
1310 f->SetTime(t);
1312 f->Mult(x, k1);
1313
1314 f->SetTime(t+dt);
1316 f->ImplicitSolve(dt, x, k2);
1317
1318 f->SetTime(t);
1319 x.Add(dt, k1);
1320 x.Add(dt, k2);
1321 t += dt;
1322}
1323
1325{
1326 ODESolver::Init(f_);
1327 int n = f->Width();
1328 k1_exp.SetSize(n, mem_type);
1329 k2_exp.SetSize(n, mem_type);
1330 k_imp.SetSize(n, mem_type);
1331 y.SetSize(n, mem_type);
1332}
1333
1335{
1336 double gamma = 1 - sqrt(2)/2;
1337 double delta = 1 - 1/(2*gamma);
1338
1339 f->SetTime(t);
1340
1341 //K1 exp is just f_1(t, x)
1343 f->Mult(x, k1_exp);
1344
1345 //K2 exp is f_1(t + gamma dt, x + dt gamma K1)
1346 f->SetTime(t + gamma*dt);
1347 add(x, dt*gamma, k1_exp, y);
1349 f->Mult(y, k2_exp);
1350
1351 //K2_imp = f_2(t + gamma dt, x + dt gamma K2_imp)
1353 f->ImplicitSolve(dt*gamma, x, k_imp);
1354 //reuse k_imp to avoid extra vector
1355
1356 //K3_imp = f_2(t+dt,x + dt(1-gamma)K2_imp + dt gamma K3_imp)
1357 f -> SetTime(t + dt);
1358 //add(x, dt*(1-gamma), k2_imp, z);
1359 //optimization to avoid extra vector
1360 x.Add(dt*(1-gamma), k_imp);
1362 //f->ImplicitSolve(dt*gamma, z, k3_imp);
1363 //reuse k_imp to avoid extra vector
1364 f->ImplicitSolve(dt*gamma, x, k_imp);
1365
1366 //add it all up
1367 x.Add(dt*delta, k1_exp);
1368 x.Add(dt*(1-delta), k2_exp);
1369 //x.Add(dt*(1-gamma), k2_imp); it is already added to x above
1370 x.Add(dt*gamma, k_imp);
1371 t += dt;
1372}
1373
1375{
1376 ODESolver::Init(f_);
1377 int n = f->Width();
1378 k1_exp.SetSize(n, mem_type);
1379 k2_exp.SetSize(n, mem_type);
1380 k3_exp.SetSize(n, mem_type);
1381 k_imp.SetSize(n, mem_type);
1382 y.SetSize(n, mem_type);
1383}
1384
1386{
1387 double gamma = 1 - sqrt(2)/2;
1388 double delta = -2*sqrt(2)/3;
1389
1390 f->SetTime(t);
1391
1392 //K1 exp is just f_1(t, x)
1394 f->Mult(x, k1_exp);
1395
1396 //K2 exp is f_1(t + gamma dt, x + dt gamma K1)
1397 f->SetTime(t + gamma*dt);
1398 add(x, dt*gamma, k1_exp, y);
1400 f->Mult(y, k2_exp);
1401
1402 //K3 Exp is f_1(t + dt, x + dt delta K1_exp + dt (1-delta) K2_exp)
1403 f->SetTime(t + dt);
1404 add(x, dt*delta, k1_exp, y);
1405 //add(y, dt*(1-delta), k2_exp, w);
1406 //optimization to avoid extra vector
1407 y.Add(dt*(1-delta), k2_exp);
1409 f->Mult(y, k3_exp);
1410
1411 //K2_imp = f_2(t + gamma dt, x + dt gamma K2_imp)
1412 f->SetTime(t + gamma*dt);
1414 f->ImplicitSolve(dt*gamma, x, k_imp);
1415
1416 //K3_imp = f_2(t+dt,x + dt(1-gamma)K2_imp + dt gamma K3_imp)
1417 f -> SetTime(t + dt);
1418 //add(x, dt*(1-gamma), k2_imp, z);
1419 x.Add(dt*(1-gamma), k_imp);
1421 f->ImplicitSolve(dt*gamma, x, k_imp);
1422
1423 //add it all up
1424 x.Add(dt*delta, k2_exp);
1425 x.Add(dt*(1-delta), k3_exp);
1426 //x.Add(dt*(1-gamma), k2_imp); // it is already added to x above
1427 x.Add(dt*gamma, k_imp);
1428 t += dt;
1429}
1430
1432{
1433 ODESolver::Init(f_);
1434 int n = f->Width();
1435 k1_exp.SetSize(n, mem_type);
1436 k2_exp.SetSize(n, mem_type);
1437 k3_exp.SetSize(n, mem_type);
1438 k4_exp.SetSize(n, mem_type);
1439 k2_imp.SetSize(n, mem_type);
1440 k3_imp.SetSize(n, mem_type);
1441 y.SetSize(n, mem_type);
1442}
1443
1445{
1446 double gamma = 0.4358665215;
1447 double b1 = 1.208496649;
1448 double b2 = -0.644363171;
1449 double a_31 = 0.3212788860;
1450 double a_32 = 0.3966543747;
1451 double a_41 = -0.105858296;
1452 double a_42 = 0.5529291479;
1453 double a_43 = 0.5529291479;
1454
1455 //K1_exp
1456 f->SetTime(t);
1458 f->Mult(x, k1_exp);
1459
1460 //K2_imp, K2_exp
1461 f->SetTime(t + gamma*dt);
1462 add(x, dt*gamma, k1_exp, y);
1464 f->Mult(y, k2_exp);
1466 f->ImplicitSolve(dt*gamma, x, k2_imp);
1467
1468 //K3_imp, K3_exp
1469 f->SetTime(t + (1+gamma)/2*dt);
1470 add(x, dt*a_31, k1_exp, y);
1471 //add(y, dt*a_32, k2_exp, w);
1472 //optimization to avoid extra vector
1473 y.Add(dt*a_32, k2_exp);
1475 f->Mult(y, k3_exp);
1476 add(x, dt*(1-gamma)/2, k2_imp, y);
1478 f->ImplicitSolve(dt*gamma, y, k3_imp);
1479
1480 //K4_imp, K4_exp
1481 f->SetTime(t+dt);
1482 add(x, dt*a_41, k1_exp, y);
1483 //add(y, dt*a_42, k2_exp, v);
1484 y.Add(dt*a_42, k2_exp);
1485 //add(w, dt*a_43, k3_exp, v);
1486 y.Add(dt*a_43, k3_exp);
1488 f->Mult(y, k4_exp);
1489 //add(x, dt*b1, k2_imp, z);
1490 //add(z, dt*b2, k3_imp, u);
1491 //optimization to avoid extra vector
1492 x.Add(dt*b1, k2_imp);
1493 x.Add(dt*b2, k3_imp);
1495 f->ImplicitSolve(dt*gamma, x, k3_imp);
1496
1497 //add it all together
1498 x.Add(dt*b1, k2_exp);
1499 x.Add(dt*b2, k3_exp);
1500 x.Add(dt*gamma, k4_exp);
1501 //x.Add(dt*b1, k2_imp); //already added above
1502 //x.Add(dt*b2, k3_imp); //already added above
1503 x.Add(dt*gamma, k3_imp);
1504 t += dt;
1505
1506}
1507
1508}
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:516
void CheckTimestep(real_t dt)
Definition ode.cpp:549
AdamsBashforthSolver(int s_, const real_t *a_)
Definition ode.cpp:510
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:524
std::unique_ptr< ODESolver > RKsolver
Definition ode.hpp:506
std::unique_ptr< ODESolver > RKsolver
Definition ode.hpp:589
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:597
AdamsMoultonSolver(int s_, const real_t *a_)
Definition ode.cpp:583
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:589
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:840
void Print(std::ostream &out=mfem::out, int width=4) const
Prints array to stream with width elements per row.
Definition array.cpp:24
The classical midpoint method.
Definition ode.hpp:899
Backward Euler ODE solver. L-stable.
Definition ode.hpp:356
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:654
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:660
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:825
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:833
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:860
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:868
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:322
virtual ~ExplicitRKSolver()
Definition ode.cpp:352
ExplicitRKSolver(int s_, const real_t *a_, const real_t *b_, const real_t *c_)
Definition ode.cpp:301
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:311
The classical forward Euler method.
Definition ode.hpp:245
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:193
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:187
void PrintProperties(std::ostream &os=mfem::out)
Definition ode.cpp:1219
void Init(SecondOrderTimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:1211
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:1249
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:897
void SetRhoInf(real_t rho_inf)
Definition ode.cpp:905
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:943
void PrintProperties(std::ostream &os=mfem::out)
Definition ode.cpp:916
Forward-backward Euler method.
Definition ode.hpp:955
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:1300
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:1308
Second order, 2/3-stage implicit-explicit (IMEX) Runge-Kutta (RK) method.
Definition ode.hpp:987
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:1374
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:1385
Second order, two-stage implicit-explicit (IMEX) Runge-Kutta (RK) method.
Definition ode.hpp:971
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:1324
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:1334
Third order, 3/4-stage implicit-explicit (IMEX) Runge-Kutta (RK) method.
Definition ode.hpp:1004
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:1431
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:1444
Implicit midpoint method. A-stable, not L-stable.
Definition ode.hpp:369
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:675
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:669
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:1176
void PrintProperties(std::ostream &os=mfem::out)
Definition ode.cpp:1146
TimeDependentOperator * f
Pointer to the associated TimeDependentOperator.
Definition ode.hpp:124
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:181
static MFEM_EXPORT std::unique_ptr< ODESolver > SelectImplicit(const int ode_solver_type)
Definition ode.cpp:76
static MFEM_EXPORT std::string Types
Definition ode.hpp:199
static MFEM_EXPORT std::unique_ptr< ODESolver > Select(const int ode_solver_type)
Definition ode.cpp:40
MemoryType mem_type
Definition ode.hpp:125
static MFEM_EXPORT std::string ImplicitTypes
Definition ode.hpp:197
static MFEM_EXPORT std::unique_ptr< ODESolver > SelectIMEX(const int ode_solver_type)
Definition ode.cpp:115
static MFEM_EXPORT std::string ExplicitTypes
Definition ode.hpp:196
static MFEM_EXPORT std::unique_ptr< ODESolver > SelectExplicit(const int ode_solver_type)
Definition ode.cpp:52
static MFEM_EXPORT std::string IMEXTypes
Definition ode.hpp:198
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:142
void Set(int i, Vector &state) override
Set the ith state vector.
Definition ode.cpp:160
void SetSize(int vsize, MemoryType mem_type)
Set the number of stages and the size of the vectors.
Definition ode.cpp:130
void Append(Vector &state) override
Add state vector and increment state size.
Definition ode.cpp:166
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:173
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:202
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:210
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition ode.hpp:278
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:231
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:239
The classical explicit forth-order Runge-Kutta method, RK4.
Definition ode.hpp:291
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:263
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:272
SDIRK23Solver(int gamma_opt=1)
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:711
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:704
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:775
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:768
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:739
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:731
void Step(Vector &q, Vector &p, real_t &t, real_t &dt) override
Definition ode.cpp:980
void Step(Vector &q, Vector &p, real_t &t, real_t &dt) override
Definition ode.cpp:993
TimeDependentOperator * F_
Definition ode.hpp:683
Vector dq_
Definition ode.hpp:687
Vector dp_
Definition ode.hpp:686
virtual void Init(Operator &P, TimeDependentOperator &F)
Definition ode.cpp:971
Operator * P_
Definition ode.hpp:684
void Step(Vector &q, Vector &p, real_t &t, real_t &dt) override
Definition ode.cpp:1050
SIAVSolver(int order)
Definition ode.cpp:1008
Abstract class for solving systems of ODEs: d2x/dt2 = f(x,dx/dt,t)
Definition ode.hpp:724
SecondOrderTimeDependentOperator * f
Pointer to the associated TimeDependentOperator.
Definition ode.hpp:727
void MidPointStep(Vector &x, Vector &dxdt, real_t &t, real_t &dt)
Definition ode.cpp:1125
ODEStateDataVector state
Definition ode.hpp:729
static MFEM_EXPORT std::string Types
Help info for SecondOrderODESolver options.
Definition ode.hpp:815
void EulerStep(Vector &x, Vector &dxdt, real_t &t, real_t &dt)
Definition ode.cpp:1111
static MFEM_EXPORT SecondOrderODESolver * Select(const int ode_solver_type)
Function selecting the desired SecondOrderODESolver.
Definition ode.cpp:1081
virtual void Init(SecondOrderTimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:1139
Base abstract class for second order time dependent operators.
Definition operator.hpp:744
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:344
bool isExplicit() const
True if type is EXPLICIT.
Definition operator.hpp:409
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 SetEvalMode(const EvalMode new_eval_mode)
Set the evaluation mode of the time-dependent operator.
Definition operator.hpp:429
virtual void SetTime(const real_t t_)
Set the current time.
Definition operator.hpp:406
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:808
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:801
Vector data type.
Definition vector.hpp:82
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:584
Vector & Add(const real_t a, const Vector &Va)
(*this) += a * Va
Definition vector.cpp:326
real_t b
Definition lissajous.cpp:42
real_t delta
Definition lissajous.cpp:43
real_t a
Definition lissajous.cpp:41
mfem::real_t real_t
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:414
MemoryType GetMemoryType(MemoryClass mc)
Return a suitable MemoryType for a given MemoryClass.
float real_t
Definition config.hpp:46
MemoryType
Memory types supported by MFEM.
real_t p(const Vector &x, real_t t)