linear_multistep_method

Examples:

>>> import nodepy.linear_multistep_method as lm
>>> ab3=lm.Adams_Bashforth(3)
>>> ab3.order()
3
>>> bdf2=lm.backward_difference_formula(2)
>>> bdf2.order()
2
>>> bdf2.is_zero_stable()
True
>>> bdf7=lm.backward_difference_formula(7)
>>> bdf7.is_zero_stable()
False
>>> bdf3=lm.backward_difference_formula(3)
>>> bdf3.A_alpha_stability()
86
>>> ssp32=lm.elm_ssp2(3)
>>> ssp32.order()
2
>>> ssp32.ssp_coefficient()
1/2
>>> ssp32.plot_stability_region() 
<Figure...
class nodepy.linear_multistep_method.LinearMultistepMethod(alpha, beta, name='Linear multistep method', shortname='LMM', description='')[source]

Implementation of linear multistep methods in the form:

$$ \alpha_k y_{n+k} + \alpha_{k-1} y_{n+k-1} + \cdots + \alpha_0 y_n = h ( \beta_k f_{n+k} + \cdots + \beta_0 f_n ) $$

Methods are automatically normalized so that \(\alpha_k=1\).

Notes: Representation follows Hairer & Wanner p. 368, NOT Butcher.

References:
characteristic_polynomials()[source]

Returns the characteristic polynomials (also known as generating polynomials) of a linear multistep method. They are:

$$ \rho(z) = \sum_{j=0}^k \alpha_k z^k $$

$$ \sigma(z) = \sum_{j=0}^k \beta_k z^k $$

Examples:

>>> from nodepy import lm
>>> ab5 = lm.Adams_Bashforth(5)
>>> rho,sigma = ab5.characteristic_polynomials()
>>> print(rho)
   5     4
1 x - 1 x

>>> print(sigma)
      4         3         2
2.64 x - 3.853 x + 3.633 x - 1.769 x + 0.3486

Reference: [HNorW93] p. 370, eq. 2.4

order(tol=1e-10)[source]

Return the order of the local truncation error of a linear multistep method.

Examples:

>>> from nodepy import lm
>>> am3=lm.Adams_Moulton(3)
>>> am3.order()
4
absolute_monotonicity_radius()[source]
property p
latex()[source]

Print a LaTeX representation of a linear multistep formula.

Example:

>>> from nodepy import lm
>>> print(lm.Adams_Bashforth(2).latex())
\begin{align} y_{n + 2} - y_{n + 1} = \frac{3}{2}h f(y_{n + 1}) - \frac{1}{2}h f(y_{n})\end{align}
ssp_coefficient()[source]

Return the SSP coefficient of the method.

The SSP coefficient is given by

$$ \min_{0 \le j < k} - \alpha_k/ \beta_k $$

if \(\alpha_j<0\) and \(\beta_j>0\) for all \(j\), and is equal to zero otherwise.

Examples:

>>> from nodepy import lm
>>> ssp32=lm.elm_ssp2(3)
>>> ssp32.ssp_coefficient()
1/2

>>> bdf2=lm.backward_difference_formula(2)
>>> bdf2.ssp_coefficient()
0
plot_stability_region(N=100, bounds=None, color='r', filled=True, alpha=1.0, to_file=False, longtitle=False)[source]

The region of absolute stability of a linear multistep method is the set

$$ \{ z \in C \mid \rho(\zeta) - z \sigma(\zeta) \text{ satisfies the root condition} \} $$

where \(\rho(\zeta)\) and \(\sigma(\zeta)\) are the characteristic functions of the method.

Also plots the boundary locus, which is given by the set of points z:

$$ \{z \mid z = \rho(\exp(\imath \theta)) / \sigma(\exp(\imath\theta)), 0 \le \theta \le 2\pi \} $$

Here \(\rho\) and \(\sigma\) are the characteristic polynomials of the method.

Reference: [LeV07] section 7.6.1

Input: (all optional)
  • N – Number of gridpoints to use in each direction

  • bounds – limits of plotting region

  • color – color to use for this plot

  • filled – if true, stability region is filled in (solid); otherwise it is outlined

plot_boundary_locus(N=1000, figsize=None)[source]

Plot the boundary locus, which is given by the set of points

$$ \{z \mid z = \rho(\exp(\imath\theta)) / \sigma(\exp(\imath\theta)), 0\le \theta \le 2\pi \} $$

where \(\rho\) and \(\sigma\) are the characteristic polynomials of the method.

Reference: [LeV07] section 7.6.1

A_alpha_stability(N=1000, tol=1e-14)[source]

Angle of \(A(\alpha)\)-stability.

The result is given in degrees. The result is only accurate to about 1 degree, so we round down.

Examples:

>>> from nodepy import lm
>>> bdf5 = lm.backward_difference_formula(5)
>>> bdf5.A_alpha_stability()
51
is_explicit()[source]
is_zero_stable(tol=1e-13)[source]

True if the method is zero-stable.

Examples:

>>> from nodepy import lm
>>> bdf5=lm.backward_difference_formula(5)
>>> bdf5.is_zero_stable()
True
class nodepy.linear_multistep_method.AdditiveLinearMultistepMethod(alpha, beta, gamma, name='Additive linear multistep method')[source]

Method for solving equations of the form

$$ y’(t) = f(y) + g(y) $$

The method takes the form

$$ \alpha_k y_{n+k} + \alpha_{k-1} y_{n+k-1} + \cdots + \alpha_0 y_n = h ( \beta_k f_{n+k} + \cdots + \beta_0 f_n + \gamma_k f_{n+k} + \cdots + \gamma_0 f_n ) $$

Methods are automatically normalized so that \(\alpha_k=1\).

The usual reference for these is Ascher, Ruuth, and Whetton. But we follow a different notation (as just described).

order(tol=1e-10)[source]

Return the order of the local truncation error of an additive linear multistep method. The output is the minimum of the order of the component methods.

plot_imex_stability_region(both_real=False, N=100, color='r', filled=True, alpha=1.0, fignum=None, bounds=[-10, 1, -5, 5])[source]
Input: (all optional)
  • N – Number of gridpoints to use in each direction

  • bounds – limits of plotting region

  • color – color to use for this plot

  • filled – if true, stability region is filled in (solid); otherwise it is outlined

stiff_damping_factor(epsilon=1e-07)[source]

Return the magnitude of the largest root at \(z = -\infty\). This routine just computes a numerical approximation to the true value (with absolute accuracy epsilon).

nodepy.linear_multistep_method.Adams_Bashforth(k)[source]

Construct the k-step, Adams-Bashforth method. The methods are explicit and have order k. They have the form:

$$ y_{n+1} = y_n + h \sum_{j=0}^{k-1} \beta_j f(y_{n-k+j+1}) $$

They are generated using equations (1.5) and (1.7) from [HNorW93] III.1, along with the binomial expansion.

Examples:

>>> import nodepy.linear_multistep_method as lm
>>> ab3=lm.Adams_Bashforth(3)
>>> print(ab3)
3-step Adams-Bashforth
Explicit 3-step method of order 3
alpha = [ 0      0      -1     1     ]
beta  = [ 5/12   -4/3   23/12  0     ]

>>> ab3.order()
3

Reference: [HNorW93]

nodepy.linear_multistep_method.Nystrom(k)[source]

Construct the \(k\)-step explicit Nystrom linear multistep method. The methods are explicit and have order \(k\).

They have the form:

$$y_{n+1} = y_{n-1} + h \sum_{j=0}^{k-1} \beta_j f(y_{n-k+j+1})$$

They are generated using equations (1.13) and (1.7) from [HNorW93] III.1, along with the binomial expansion and the relation in exercise 4 on p. 367.

Note that the term “Nystrom method” is also commonly used to refer to a class of methods for second-order ODEs; those are NOT the methods generated by this function.

Examples:

>>> import nodepy.linear_multistep_method as lm
>>> nys3=lm.Nystrom(6)
>>> nys3.order()
6

Reference: [HNorW93]

nodepy.linear_multistep_method.Adams_Moulton(k)[source]

Construct the \(k\)-step, Adams-Moulton method. The methods are implicit and have order \(k+1\). They have the form:

$$ y_{n+1} = y_n + h \sum_{j=0}^{k} \beta_j f(y_{n-k+j+1}) $$

They are generated using equation (1.9) and the equation in Exercise 3 from Hairer & Wanner III.1, along with the binomial expansion.

Examples:

>>> import nodepy.linear_multistep_method as lm
>>> am3=lm.Adams_Moulton(3)
>>> am3.order()
4

Reference: [HNorW93]

nodepy.linear_multistep_method.Milne_Simpson(k)[source]

Construct the \(k\)-step, Milne-Simpson method. The methods are implicit and (for \(k \ge 3\)) have order \(k+1\). They have the form:

$$ y_{n+1} = y_{n-1} + h \sum_{j=0}^{k} \beta_j f(y_{n-k+j+1}) $$

They are generated using equation (1.15), the equation in Exercise 3, and the relation in exercise 4, all from Hairer & Wanner III.1, along with the binomial expansion.

Examples:

>>> import nodepy.linear_multistep_method as lm
>>> ms3=lm.Milne_Simpson(3)
>>> ms3.order()
4

Reference: [HNorW93]

nodepy.linear_multistep_method.backward_difference_formula(k)[source]

Construct the \(k\)-step backward differentiation method. The methods are implicit and have order \(k\).

They have the form:

$$ \sum_{j=0}^{k} \alpha_j y_{n+k-j+1} = h \beta_j f(y_{n+1}) $$

They are generated using equation (1.22’) from Hairer & Wanner III.1, along with the binomial expansion.

Examples:

>>> import nodepy.linear_multistep_method as lm
>>> bdf4=lm.backward_difference_formula(4)
>>> bdf4.A_alpha_stability()
73

Reference: [HNorW93] pp. 364-365

nodepy.linear_multistep_method.elm_ssp2(k)[source]

Returns the optimal SSP \(k\)-step linear multistep method of order 2.

Examples:

>>> import nodepy.linear_multistep_method as lm
>>> lm10=lm.elm_ssp2(10)
>>> lm10.ssp_coefficient()
8/9
nodepy.linear_multistep_method.sand_cc(s)[source]

Construct Sand’s circle-contractive method of order \(p=2(s+1)\) that uses \(2^s + 1\) steps.

Examples:

>>> import nodepy.linear_multistep_method as lm
>>> cc4 = lm.sand_cc(4)
>>> cc4.order()
10
>>> cc4.ssp_coefficient()
1/8

Reference: [San86]

nodepy.linear_multistep_method.arw2(gam, c)[source]

Returns the second order IMEX additive multistep method based on the parametrization in Section 3.2 of Ascher, Ruuth, & Whetton. The parameters are gam and c. Known methods are obtained with the following values:

(1/2,0): CNAB (1/2,1/8): MCNAB (0,1): CNLF (1,0): SBDF

Examples:

>>> from nodepy import lm
>>> import sympy
>>> CNLF = lm.arw2(0,sympy.Rational(1))
>>> CNLF.order()
2
>>> CNLF.method1.ssp_coefficient()
1
>>> CNLF.method2.ssp_coefficient()
0
>>> print(CNLF.stiff_damping_factor()) 
0.999...
nodepy.linear_multistep_method.arw3(gam, theta, c)[source]

Returns the third order IMEX additive multistep method based on the parametrization in Section 3.3 of Ascher, Ruuth, & Whetton. The parameters are gam, theta, and c. Known methods are obtained with the following values:

(1,0,0): SBDF3

Note that there is one sign error in the ARW paper; it is corrected here.

nodepy.linear_multistep_method.loadLMM(which='All')[source]

Load a set of standard linear multistep methods for testing.

Examples:

>>> from nodepy import lm
>>> ebdf5 = lm.loadLMM('eBDF5')
>>> print(ebdf5)
eBDF 5
Explicit 5-step method of order 5
alpha = [ -12/137   75/137    -200/137  300/137   -300/137  1        ]
beta  = [ 60/137    -300/137  600/137   -600/137  300/137   0        ]

>>> ebdf5.is_zero_stable()
True