BDF¶
- class pathsim.solvers.bdf.BDF(*solver_args, **solver_kwargs)[source]¶
Bases:
ImplicitSolverBase class for the backward differentiation formula (BDF) integrators.
Notes
This solver class is not intended to be used directly
- x¶
internal ‘working’ state
- Type:
numeric, array[numeric]
- opt¶
optimizer instance to solve the implicit update equation
- Type:
NewtonAnderson, Anderson, etc.
- history¶
internal history of past results
- Type:
deque[numeric]
- startup¶
internal solver instance for startup (building history) of multistep methods (using ‘DIRK3’ for ‘BDF’ methods)
- Type:
- classmethod cast(other, parent, **solver_kwargs)[source]¶
cast to this solver needs special handling of startup method
- classmethod create(initial_value, parent=None, from_engine=None, **solver_kwargs)[source]¶
Create a new BDF solver, properly initializing the startup solver.
- Parameters:
- Returns:
engine – new BDF solver instance
- Return type:
- stages(t, dt)[source]¶
Generator that yields the intermediate evaluation time during the timestep ‘t + ratio * dt’.
- buffer(dt)[source]¶
buffer the state for the multistep method
- Parameters:
dt (float) – integration timestep
- step(f, dt)[source]¶
Performs the explicit timestep for (t+dt) based on the state and input at (t).
Note
This is only required for the startup solver.
- Parameters:
f (numeric, array[numeric]) – evaluation of rhs function
dt (float) – integration timestep
- Returns:
success (bool) – True if the timestep was successful
error (float) – estimated error of the internal error controller
scale (float) – estimated timestep rescale factor for error control
- class pathsim.solvers.bdf.BDF2(*solver_args, **solver_kwargs)[source]¶
Bases:
BDFFixed-step 2nd order BDF method. A-stable.
\[x_{n+1} = \tfrac{4}{3}\,x_n - \tfrac{1}{3}\,x_{n-1} + \tfrac{2}{3}\,h\,f(x_{n+1}, t_{n+1})\]Uses
DIRK3as startup solver for the first step.Characteristics¶
Order: 2
Implicit linear multistep, fixed timestep
A-stable
Note
The workhorse fixed-step stiff solver. A-stability means no eigenvalue in the left half-plane causes instability, regardless of the timestep. Well suited for block diagrams with a fixed simulation clock and moderately-to-very stiff dynamics. For adaptive stepping, use
GEAR21orESDIRK43.References
- class pathsim.solvers.bdf.BDF3(*solver_args, **solver_kwargs)[source]¶
Bases:
BDFFixed-step 3rd order BDF method. \(A(\alpha)\)-stable with \(\alpha \approx 86°\).
Uses
DIRK3as startup solver for the first two steps.Characteristics¶
Order: 3
Implicit linear multistep, fixed timestep
\(A(\alpha)\)-stable, \(\alpha \approx 86°\)
Note
Higher accuracy than
BDF2with only a slight reduction in the stability wedge. The \(86°\) sector covers nearly the entire left half-plane, so most stiff block diagrams remain well-handled. For adaptive stepping, useGEAR32.References
[1] Gear, C. W. (1971). “Numerical Initial Value Problems in Ordinary Differential Equations”. Prentice-Hall.
[2] Hairer, E., & Wanner, G. (1996). “Solving Ordinary Differential Equations II: Stiff and Differential-Algebraic Problems”. Springer Series in Computational Mathematics, Vol. 14. :doi:`10.1007/978-3-642-05221-7`
- class pathsim.solvers.bdf.BDF4(*solver_args, **solver_kwargs)[source]¶
Bases:
BDFFixed-step 4th order BDF method. \(A(\alpha)\)-stable with \(\alpha \approx 73°\).
Uses
DIRK3as startup solver for the first three steps.Characteristics¶
Order: 4
Implicit linear multistep, fixed timestep
\(A(\alpha)\)-stable, \(\alpha \approx 73°\)
Note
The narrower stability wedge compared to
BDF2/BDF3means eigenvalues close to the imaginary axis may be poorly damped. Safe for block diagrams whose stiff modes are strongly dissipative (well inside the left half-plane). For problems with near-imaginary eigenvalues (e.g. lightly damped oscillators), prefer lower-order BDF or an ESDIRK solver.References
[1] Gear, C. W. (1971). “Numerical Initial Value Problems in Ordinary Differential Equations”. Prentice-Hall.
[2] Hairer, E., & Wanner, G. (1996). “Solving Ordinary Differential Equations II: Stiff and Differential-Algebraic Problems”. Springer Series in Computational Mathematics, Vol. 14. :doi:`10.1007/978-3-642-05221-7`
- class pathsim.solvers.bdf.BDF5(*solver_args, **solver_kwargs)[source]¶
Bases:
BDFFixed-step 5th order BDF method. \(A(\alpha)\)-stable with \(\alpha \approx 51°\).
Uses
DIRK3as startup solver for the first four steps.Characteristics¶
Order: 5
Implicit linear multistep, fixed timestep
\(A(\alpha)\)-stable, \(\alpha \approx 51°\)
Note
The stability wedge is noticeably smaller than
BDF3orBDF4. Only appropriate when the stiff eigenvalues of the block diagram are concentrated well inside the left half-plane and high accuracy per step is essential. For most stiff systems,BDF2orBDF3with a smaller timestep is more robust.References
[1] Gear, C. W. (1971). “Numerical Initial Value Problems in Ordinary Differential Equations”. Prentice-Hall.
[2] Hairer, E., & Wanner, G. (1996). “Solving Ordinary Differential Equations II: Stiff and Differential-Algebraic Problems”. Springer Series in Computational Mathematics, Vol. 14. :doi:`10.1007/978-3-642-05221-7`
- class pathsim.solvers.bdf.BDF6(*solver_args, **solver_kwargs)[source]¶
Bases:
BDFFixed-step 6th order BDF method. Not A-stable (\(\alpha \approx 18°\)).
Uses
DIRK3as startup solver for the first five steps.Characteristics¶
Order: 6
Implicit linear multistep, fixed timestep
\(A(\alpha)\)-stable, \(\alpha \approx 18°\) (not A-stable)
Note
The very narrow stability wedge means that most stiff problems will be unstable at practical timestep sizes. Provided mainly for completeness. For 6th order accuracy on non-stiff systems, the explicit
RKV65is cheaper. For stiff systems,BDF2–BDF4with a smaller timestep or an ESDIRK solver are far more robust.References
[1] Gear, C. W. (1971). “Numerical Initial Value Problems in Ordinary Differential Equations”. Prentice-Hall.
[2] Hairer, E., & Wanner, G. (1996). “Solving Ordinary Differential Equations II: Stiff and Differential-Algebraic Problems”. Springer Series in Computational Mathematics, Vol. 14. :doi:`10.1007/978-3-642-05221-7`