libode
Easy to compile, fast ODE integrators as C++ classes
Loading...
Searching...
No Matches
ode_ssp_3.cc
Go to the documentation of this file.
1
2
3#include "ode_ssp_3.h"
4
5namespace ode {
6
7OdeSsp3::OdeSsp3 (unsigned long neq) :
8 OdeAdaptive (neq, false),
9 OdeRK (neq, 3),
10 OdeERK (neq) {
11
12 method_ = "Ssp3";
13 //tableau of coefficients
14 c2 = 1.0; a21 = c2;
15 c3 = 1.0/2; a31 = 1.0/4; a32 = 1.0/4;
16 b1 = 1.0/6; b2 = 1.0/6; b3 = 2.0/3;
17}
18
19void OdeSsp3::step_ (double dt) {
20
21 unsigned long i;
22
23 //------------------------------------------------------------------
24 ode_fun_(sol_, k_[0]);
25
26 //------------------------------------------------------------------
27 for (i=0; i<neq_; i++) soltemp_[i] = sol_[i] + dt*a21*k_[0][i];
28 ode_fun_(soltemp_, k_[1]);
29
30 //------------------------------------------------------------------
31 for (i=0; i<neq_; i++) soltemp_[i] = sol_[i] + dt*(a31*k_[0][i] + a32*k_[1][i]);
32 ode_fun_(soltemp_, k_[2]);
33
34 //------------------------------------------------------------------
35 //store solution
36 for (i=0; i<neq_; i++) sol_[i] += dt*(b1*k_[0][i] + b2*k_[1][i] + b3*k_[2][i]);
37}
38
39} // namespace ode
Base class implementing solver functions with adaptive time steps.
unsigned long neq_
number of equations in the system of ODEs
Definition ode_base.h:327
std::string method_
the "name" of the solver/method, as in "Euler" or "RK4"
Definition ode_base.h:319
double * sol_
array for the solution, changing over time
Definition ode_base.h:333
void ode_fun_(double *solin, double *fout)
wrapper, calls ode_fun() and increases the neval counter by one
Definition ode_base.cc:87
Base class providing space for temporary solutions moving through RK stages.
Definition ode_erk.h:9
double * soltemp_
temporary solution vector
Definition ode_erk.h:22
Provides space for stage slope values, an array of arrays for k values.
Definition ode_rk.h:9
double ** k_
stage evaluations
Definition ode_rk.h:23
OdeSsp3(unsigned long neq)
constructs
Definition ode_ssp_3.cc:7