Documentation ¶
Overview ¶
package ode provides interfaces and structures to solve multivariable ODE initial value problems.
Example (Solve) ¶
package main import ( "fmt" "log" "gonum.org/v1/exp/ode" "gonum.org/v1/gonum/mat" ) func main() { const ( g = -10.0 // local gravity field [m.s^-2] ) // we declare our physical model. First argument is initial time, which is 0 seconds. // Next is the initial state vector, which corresponds to 100 meters above the ground // with 0 m/s velocity. ballModel, err := ode.NewIVP(0, mat.NewVecDense(2, []float64{100., 0.}), func(yvec *mat.VecDense, _ float64, xvec mat.Vector) { // this anonymous function defines the physics. // The first variable xvec[0] corresponds to position // second variable xvec[1] is velocity. Dx := xvec.AtVec(1) // yvec represents change in xvec, or derivative with respect to domain // Change in position will be equal to velocity, which is the second variable: // thus yvec[0] = xvec[1], which is the same as saying "change in xvec[0]" is equal to xvec[1] yvec.SetVec(0, Dx) // change in velocity is acceleration. We suppose our ball is on earth accelerating at `g` yvec.SetVec(1, g) }) if err != nil { log.Fatal(err) } // Here we choose our algorithm. Runge-Kutta 4th order is used var solver ode.Integrator = ode.NewDormandPrince5(ode.DefaultParam) // Solve function makes it easy to integrate a problem without having // to implement the `for` loop. This example integrates the IVP with a step size // of 0.1 over a domain of 10. arbitrary units, in this case, 10 seconds. results, err := ode.SolveIVP(ballModel, solver, 0.1, 10.) fmt.Println(results) }
Output:
Index ¶
Examples ¶
Constants ¶
This section is empty.
Variables ¶
var DefaultParam = Parameters{}
Functions ¶
This section is empty.
Types ¶
type DoPri5 ¶
type DoPri5 struct {
// contains filtered or unexported fields
}
func NewDormandPrince5 ¶
func NewDormandPrince5(cfg Parameters) *DoPri5
NewDormandPrince5 returns a adaptive-ready Dormand Prince solver of order 5
To enable step size adaptivity minimum step size must be set and absolute tolerance must be set. i.e:
NewDormandPrince5(ConfigScalarTolerance(0, 0.1), ConfigStepLimits(1, 1e-3))
If a invalid configuration is passed the function panics.
type IVP ¶
type IVP struct { // Initial values for the state vector Y0 mat.Vector // Independent variable point at which Y0 is evaluated T0 float64 // Func are the differential equations f(t,y(t)) such that // dst = y'(t) = Func(t, y(t)) Func func(dst *mat.VecDense, t float64, y mat.Vector) }
IVP defines a multivariable, initial value problem represented by a system of ordinary differential equations.
These problems have the form
y'(t) = f(t, y(t)) y(t_0) = y_0
Where: t is a scalar representing the integration domain, which is time for most physical problems. y is the state vector. y' is the vector of first derivatives of the state vector y. f are the differential equations represented by Func. An initial value problem is characterized by the initial conditions imposed on the state vector y at the beginning of the integration domain. These initial conditions are returned by the IV() method for the state vector as y0.
The term "state vector" and "state variables" are used interchangeably throughout the code and refer to y vector of independent variables.
type IVP2 ¶
type IVP2 struct { Y0 mat.Vector DY0 mat.Vector T0 float64 // Func are the second derivatives of the solution such that // dst = y”(t) = Func(t, y(t)) Func func(dst *mat.VecDense, t float64, y mat.Vector) }
IVP2 defines a multivariable, initial value problem represented by a second order system of ordinary differential equations.
These problems have the form
y''(t) = f(t, y(t)) y(t_0) = y_0 y'(t_0) = y'_0
type Integrator ¶
type Integrator interface { // Init initializes the integrator and sets the initial condition. Init(IVP) // Step advances the current state by taking at most the given step. It returns a proposed step size // for the next step and an error indicating whether the step was successful. Step(step float64) (stepNext float64, err error) // State stores the current state of the integrator in-place in dst. State(dst *State) }
Integrator can integrate an initial-value problem (IVP) for a first-order system of ordinary differential equations (ODEs).
type Parameters ¶
type Parameters struct { // Permissible tolerance given an adaptive method AbsTolerance float64 // Minimum/Maximum step allowed for a single iteration MinStep, MaxStep float64 }
Parameters represents the most common configuration parameters
type RKN1210 ¶
type RKN1210 struct {
// contains filtered or unexported fields
}
RKN1210 Runge-Kutta-Nyström integrator of order 12(10). This method was sourced from J. R. DORMAND, M. E. A. EL-MIKKAWY, P. J. PRINCE, High-Order Embedded Runge-Kutta-Nystrom Formulae, IMA Journal of Numerical Analysis, Volume 7, Issue 4, October 1987, Pages 423–430, https://doi.org/10.1093/imanum/7.4.423
func NewRKN1210 ¶
func NewRKN1210(cfg Parameters) *RKN1210
NewRKN1210 configures a new RKN1210 instance.
func (*RKN1210) Init ¶
Init initializes the Runge-Kutta-Nyström integrator with a second order ODE initial value problem.