DormandPrince.Checks
— Type@enum Checks
Enum used to represent status of checks on user options to the solver in the Report
type.
Values
INPUT_CHECKS_SUCCESSFUL
: All checks on user options passed.MAX_ALLOWED_STEPS_NEGATIVE
: The maximum allowed steps value is negative.UNSUPPORTED_UROUND
: The uround value is either too small (<= 1e-35) or too large (>= 1.0).CURIOUS_BETA
: The Beta value for stabilized step control is greater than 0.2.CURIOUS_SAFETY_FACTOR
: The safety factor is either too large (>= 1.0) or too small (<= 1e-4).
DormandPrince.DP5Solver
— Typestruct DP5Solver
A 5th Order Dormand-Prince solver object that contains:
- the ODE problem of interest
- The solution to the ODE at the last integrated-to time point
- intermediate arrays used in the Runge-Kutta method
- constants and variables used by the solver
- user-provided options for the solver
The storage of such intermediate information allows for efficient integration from a previously integrated-to time point and a future time point.
DormandPrince.DP5Solver
— MethodDP5Solver(f, x, y; kw...)
Create a 5th Order Dormand-Prince solver object to solve the ODE problem y' = f(x, y)
.
Examples
julia> fcn(x, y, f)
f[1] = 0.85y[1]
end
julia> solver = DP5Solver(fcn, 0.0, [19.0]; atol = 1e-10, rtol = 1e-10)
Arguments
f
: The function representing the ODE, should be in the formf(x, y, f)
.x
: The starting time point of the ODE problem.y
: The initial value of the ODE in vector form
Keyword Arguments
atol
: Absolute Tolerance. Default is 1e-10rtol
: Relative Tolerance. Default is 1e-10uround
: Rounding unit, default is eps(T) where T <: Real.safety_factor
: Safety factor in step size prediction, default is 0.9.- Step Size Selection parameters, with the new step size being subject to the constraint:
step_size_selection_one <= new_step/old_step <= step_size_selection_two
step_size_selection_one
: Defaut is 0.2step_size_selection_two
: Default is 10.0
beta
: Beta for stabilized step size control. Large values of beta (<= 0.1) make step size control more stable.maximal_step_size
: The largest the step size can be. Default is 0.0, which internally translates toxend - x
.initial_step_size
: Initial step size, default is 0.0. An initial guess is computed internally.maximum_allowed_steps
: Maximum number of allowed steps, default is 100000.print_error_messages
: Whether to print error messages, default is true.stiffness_test_activation_step
: After integer multiples of this number of steps, perform stiffness detection. Default is 1000.
DormandPrince.DP8Solver
— Typestruct DP8Solver
An 8th Order Dormand-Prince solver object that contains:
- the ODE problem of interest
- The solution to the ODE at the last integrated-to time point
- intermediate arrays used in the Runge-Kutta method
- constants and variables used by the solver
- user-provided options for the solver
The storage of such intermediate information allows for efficient integration from a previously integrated-to time point and a future time point.
DormandPrince.DP8Solver
— MethodDP8Solver(f, x, y; kw...)
Create an 8th Order Dormand-Prince solver object to solve the ODE problem y' = f(x, y)
.
Examples
julia> fcn(x, y, f)
f[1] = 0.85y[1]
end
julia> solver = DP8Solver(fcn, 0.0, [19.0]; atol = 1e-10, rtol = 1e-10)
Arguments
f
: The function representing the ODE, should be in the formf(x, y, f)
.x
: The starting time point of the ODE problem.y
: The initial value of the ODE in vector form
Keyword Arguments
atol
: Absolute Tolerance. Default is 1e-10rtol
: Relative Tolerance. Default is 1e-10uround
: Rounding unit, default is eps(T) where T <: Real.safety_factor
: Safety factor in step size prediction, default is 0.9.- Step Size Selection parameters, with the new step size being subject to the constraint:
step_size_selection_one <= new_step/old_step <= step_size_selection_two
step_size_selection_one
: Defaut is 0.2step_size_selection_two
: Default is 10.0
beta
: Beta for stabilized step size control. Large values of beta (<= 0.1) make step size control more stable.maximal_step_size
: The largest the step size can be. Default is 0.0, which internally translates toxend - x
.initial_step_size
: Initial step size, default is 0.0. An initial guess is computed internally.maximum_allowed_steps
: Maximum number of allowed steps, default is 100000.print_error_messages
: Whether to print error messages, default is true.stiffness_test_activation_step
: After integer multiples of this number of steps, perform stiffness detection. Default is 1000.
DormandPrince.Idid
— Type@enum Idid
Enum used to represent status of the integration process in the Report
type.
Values
COMPUTATION_SUCCESSFUL
: Integration completed successfully.INPUT_NOT_CONSISTENT
: The options given to the solver violate acceptable ranges (see the
checks
field of the Report
type for more information).
LARGER_NMAX_NEEDED
: The maximum number of allowed steps is too small.STEP_SIZE_BECOMES_TOO_SMALL
: The step size becomes too small.STEP_SIZE_BECOMES_NAN
: The step size becomes NaN.
DormandPrince.Options
— Typestruct Options{T <: Real}
Holds the options used in the integration process by a solver object.
Fields
atol
: Absolute Tolerance. Default is 1e-10rtol
: Relative Tolerance. Default is 1e-10uround
: Rounding unit, default is eps(T) where T <: Real.safety_factor
: Safety factor in step size prediction, default is 0.9.- Step Size Selection parameters, with the new step size being subject to the constraint:
step_size_selection_one <= new_step/old_step <= step_size_selection_two
step_size_selection_one
: Defaut is 0.2step_size_selection_two
: Default is 10.0
beta
: Beta for stabilized step size control. Large values of beta (<= 0.1) make step size control more stable.maximal_step_size
: The largest the step size can be. Default is 0.0, which internally translates toxend - x
.initial_step_size
: Initial step size, default is 0.0. An initial guess is computed internally.maximum_allowed_steps
: Maximum number of allowed steps, default is 100000.print_error_messages
: Whether to print error messages, default is true.stiffness_test_activation_step
: After integer multiples of this number of steps, perform stiffness detection. Default is 1000.
DormandPrince.Report
— Typestruct Report{T <: Real}
Contains data on the integration process including after integration has completed and if an error was detected within options provided for integration.
Fields
x_final
: The final time point integration reached. Should be equal to the time point provided by the user with successful usage.checks
: Status of checks performed on the options provided by the user, represented by aChecks
element. Should beINPUT_CHECKS_SUCCESSFUL
with successful usage.idid
: Status of the integration process, represented by anIdid
element. Should beCOMPUTATION_SUCCESSFUL
with successful usage.num_func_evals
: Number of function evaluations performed during integration.num_computed_steps
: Number of computed steps during integration.num_accpeted_steps
: Number of accepted steps during integration.num_rejected_steps
: Number of rejected steps during integration.
DormandPrince.SolverIterator
— Typestruct SolverIterator{T <: Real}
SolverIterator(solver, times)
Given a solver and a vector of times, this iterator will return the state of the solver at each time.
The solver will be mutated to hold the solution from the last time given.
Examples
julia> solver = DP5Solver(fcn, 0.0, [0.0])
julia> times = [1.0, 2.0, 3.0]
julia> intermediate_vals = []
julia> for (time, value) in SolverIterator(solver, times)
push!(intermediate_values, value[])
end
DormandPrince.get_current_state
— Methodget_current_state(solver::DP5Solver)
Gets the current state/solution held by the DP5Solver
object.
DormandPrince.get_current_state
— Methodget_current_state(solver::DP8Solver)
Gets the current state/solution held by the DP8Solver
object.
DormandPrince.integrate!
— Methodintegrate!(solver, time)
Integrate the ODE problem that is part of solver to the end time time
. solver
can be a DP5Solver
or a DP8Solver
type.
The solver
is mutated to hold the solution and solver
state at the time time
, allowing for efficient integration to future end times of interest through subsequent calls to integrate!
with later times.
Examples
julia> function fcn(x, y, f)
f[1] = 0.85y[1]
end
julia> solver = DP5Solver(fcn, 0.0, [19.0])
# Integrate from initial time of 0.0 to 1.0
julia> integrate!(solver, 1.0)
# integrate from the last time of 1.0 to 2.0
julia> integrate!(solver, 2.0)
# Solver object now holds solution and state at 2.0
julia> get_current_state(solver)
DormandPrince.integrate!
— Methodintegrate!(callback, solver, times)
Integrate the ODE problem that is part of solver
to the end times times
and apply the callback on the time and solution at each time. solver
can be a DP5Solver
or a DP8Solver
type.
At the end of all the times the solver holds the solution and solver state at the last time in times
.
Examples
julia> function fcn(x, y, f)
f[1] = 0.85y[1]
end
julia> times = [1.0, 1.1, 1.9, 2.4]
julia> solver = DP5Solver(fcn, 0.0, [19.0])
julia> intermediate_values = []
julia> integrate!(solver, times) do time, val
push!(intermediate_values, val[])
end