API
Main interface
IntervalRootFinding.BreadthFirstSearch
— TypeBreadthFirstSearch <: BreadthFirstBBSearch
Type implementing the BreadthFirstBBSearch
interface for interval roots finding.
Fields:
initial
: region (as aRoot
object) in which roots are searched.contractor
: contractor to use (Bisection
,Newton
orKrawczyk
)tol
: tolerance of the search
IntervalRootFinding.DepthFirstSearch
— TypeDepthFirstSearch <: DepthFirstBBSearch
Type implementing the DepthFirstBBSearch
interface for interval roots finding.
Fields:
initial
: region (as aRoot
object) in which roots are searched.contractor
: contractor to use (Bisection
,Newton
orKrawczyk
)tol
: tolerance of the search
IntervalRootFinding.branch_and_prune
— Methodbranch_and_prune(X, contractor, strategy, tol)
Generic branch and prune routine for finding isolated roots using the given contractor to determine the status of a given box X
.
See the documentation of the roots
function for explanation of the other arguments.
Polynomials.roots
— Methodroots(f, X, contractor=Newton, strategy=BreadthFirstSearch, tol=1e-15)
roots(f, deriv, X, contractor=Newton, strategy=BreadthFirstSearch, tol=1e-15)
roots(f, X, contractor, tol)
roots(f, deriv, X, contractor, tol)
Uses a generic branch and prune routine to find in principle all isolated roots of a function f:R^n → R^n
in a region X
, if the number of roots is finite.
Inputs:
f
: function whose roots will be foundX
:Interval
orIntervalBox
in which roots are searchedcontractor
: function that, when applied to the functionf
, determines the status of a given boxX
. It returns the new box and a symbol indicating the status. Current possible values areBisection
,Newton
andKrawczyk
deriv
: explicit derivative off
forNewton
andKrawczyk
strategy
:SearchStrategy
determining the order in which regions are processed.tol
: Absolute tolerance. If a region has a diameter smaller thantol
, it is returned with status:unknown
.
Root
object
IntervalRootFinding.Root
— TypeRoot
Object representing a possible root inside a given region. The field status
is either :unknown
or :unique
. If status
is :unique
then we know that there is a unique root of the function in question inside the given region.
Internally the status may also be :empty
for region guaranteed to contain no root, however such Root
s are discarded by default and thus never returned by the roots
function.
Fields
interval
: a region (eitherInterval
ofIntervalBox
) searched for roots.status
: the status of the region, valid values are:empty
,unkown
and:unique
.
IntervalRootFinding.isunique
— Methodisunique(rt)
Return whether a Root
is unique.
IntervalRootFinding.root_region
— Methodroot_region(rt)
Return the region associated to a Root
.
IntervalRootFinding.root_status
— Methodroot_status(rt)
Return the status of a Root
.
Contractors
IntervalRootFinding.Bisection
— TypeBisection{F} <: Contractor{F}
Contractor type for the bisection method.
IntervalRootFinding.Contractor
— TypeContractor{F}
Abstract type for contractors.
IntervalRootFinding.Krawczyk
— TypeKrawczyk{F, FP} <: Contractor{F}
Contractor type for the Krawczyk method.
Fields
f::F
: function whose roots are searchedf::FP
: derivative or jacobian off
(C::Krawczyk)(X, tol, α=where_bisect)
Contract an interval X
using Krawczyk operator and return the contracted interval together with its status.
Inputs
R
: Root object containing the interval to contract.tol
: Precision to which unique solutions are refined.α
: Point of bisection of intervals.
IntervalRootFinding.Newton
— TypeNewton{F, FP} <: Contractor{F}
Contractor type for the Newton method.
Fields
f::F
: function whose roots are searchedf::FP
: derivative or jacobian off
(C::Newton)(X, tol, α=where_bisect)
Contract an interval X
using Newton operator and return the contracted interval together with its status.
Inputs
R
: Root object containing the interval to contract.tol
: Precision to which unique solutions are refined.α
: Point of bisection of intervals.
IntervalRootFinding.determine_region_status
— Methoddetermine_region_status(contract, f, R)
Contraction operation for contractors using the first derivative of the function.
Currently Newton
and Krawczyk
contractors use this.
IntervalRootFinding.refine
— Methodrefine(op, X::Tuple{Symbol, Region}, tol)
Wrap the refine method to leave unchanged intervals that are not guaranteed to contain an unique solution.
IntervalRootFinding.refine
— Methodrefine(op, X::Region, tol)
Generic refine operation for Krawczyk and Newton. This function assumes that it is already known that X
contains a unique root.
IntervalRootFinding.safe_isempty
— Methodsafe_isempty(X)
Similar to isempty
function for IntervalBox
, but also works for SVector
of Interval
.
IntervalRootFinding.𝒦
— Method𝒦(f, jacobian, X, α)
Multi-variable Krawczyk operator.
IntervalRootFinding.𝒦
— Method𝒦(f, f′, X, α)
Single-variable Krawczyk operator.
The symbol for the operator is accessed with \scrK<tab>
.
IntervalRootFinding.𝒩
— Method𝒩(f, jacobian, X, α)
Multi-variable Newton operator.
IntervalRootFinding.𝒩
— Method𝒩(f, f′, X, α)
Single-variable Newton operator.
The symbol for the operator is accessed with \scrN<tab>
.
Branch-and-bound search interface
IntervalRootFinding.BBLeaf
— TypeBBLeaf{DATA} <: AbstractBBLeaf
Leaf node of a BBTree
that contains some data. Its status is either
:working
: the leaf will be further processed.:final
: the leaf won't be touched anymore.
IntervalRootFinding.BBNode
— TypeBBNode <: AbstractBBNode
Intermediate node of a BBTree
. Does not contain any data by itself, only redirect toward its children.
IntervalRootFinding.BBSearch
— TypeBBSearch{DATA}
Branch and bound search interface in element of type DATA.
This interface provide an iterable that perform the search.
There is currently three types of search supported BreadFirstBBSearch
, DepthFirstBBSearch
and KeyBBSearch
, each one processing the element of the tree in a different order. When subtyping one of these, the following methods must be implemented:
root_element(::BBSearch)
: return the element with which the search is startedprocess(::BBSearch, elem::DATA)
: return a symbol representing the action to perform with the elementelem
and an object of typeDATA
reprensenting the state of the element after processing (may returnelem
unchanged).bisect(::BBSearch, elem::DATA)
: return two elements of typeDATA
build by bisectingelem
Subtyping BBSearch
directly allows to have control over the order in which the elements are process. To do this the following methods must be implemented:
root_element(::BBSearch)
: return the first element to be processed. Use to build the initial tree.get_leaf_id!(::BBSearch, wt::BBTree)
: return the id of the next leaf that will be processed and remove it from the list of working leaves ofwt
.insert_leaf!(::BBSearch, wt::BBTree, leaf::BBLeaf)
: insert a leaf in the list of working leaves.
Valid symbols returned by the process function
:store
: the element is considered as final and is stored, it will not be further processed:bisect
: the element is bisected and each of the two resulting part will be processed:discard
: the element is discarded from the tree, allowing to free memory
IntervalRootFinding.BBTree
— TypeBBTree{DATA}
Tree storing the data used and produced by a branch and bound search in a structured way.
Nodes and leaves can be accessed using their index using the bracket syntax wt[node_id]
. However this is slow, as nodes and leaves are stored separately.
Support the iterator interface. The element yielded by the iteration are tuples (node_id, lvl)
where lvl
is the depth of the node in the tree.
IntervalRootFinding.KeyBBSearch
— TypeKeyBBSearch{DATA} <: BBSearch{DATA}
Interface to a branch and bound search that use a key function to decide which element to process first. The search process first the element with the largest key as computed by keyfunc(ks::KeyBBSearch, elem)
.
Untested.
IntervalRootFinding.data
— Methoddata(leaf::BBLeaf)
Return the data stored in the leaf.
IntervalRootFinding.data
— Methoddata(wt::BBTree)
Return all the data stored in a BBTree
as a list. The ordering of the elements is arbitrary.
IntervalRootFinding.discard_leaf!
— Methoddiscard_leaf!(wt::BBTree, id::Int)
Delete the BBLeaf
with index id
and all its ancestors to which it is the last descendant.
IntervalRootFinding.get_leaf_id!
— Methodget_leaf_id!(::BBSearch, wt::BBTree)
Return the id of the next leaf that will be processed and remove it from the list of working leaves.
Must be define for custom searches that are direct subtype of BBSearch
.
IntervalRootFinding.insert_leaf!
— Methodinsert_leaf!(::BBSearch, wt::BBTree, leaf::BBLeaf)
Insert the id of a new leaf that has been produced by bisecting an older leaf into the list of working leaves.
Must be define for custom searches that are direct subtype of BBSearch
.
IntervalRootFinding.nnodes
— Methodnnodes(wt::BBTree)
Number of nodes (including leaves) in a BBTree
.
IntervalRootFinding.root_element
— Methodroot_element(search::BBSearch)
Return the initial element of the search. The BBTree
will be build around it.
Can be define for custom searches that are direct subtype of BBSearch
, default behavior is to fetch the field initial
of the search.
Others
IntervalRootFinding.Compl
— TypeComplex numbers as 2-vectors, enough for polynomials.
IntervalRootFinding.realify
— MethodTakes a complex (polynomial) function f and returns a function g:R^2 -> R^2 that implements it.
IntervalRootFinding.realify_derivative
— MethodTakes the derivative of a complex function and returns the real jacobian that implements it.
IntervalRootFinding.gauss_elimination_interval!
— MethodSolves the system of linear equations using Gaussian Elimination. Preconditioning is used when the precondition
keyword argument is true
.
REF: Luc Jaulin et al., Applied Interval Analysis, pg. 72
IntervalRootFinding.gauss_elimination_interval1!
— MethodUsing Base.
`
IntervalRootFinding.gauss_seidel_interval!
— MethodIteratively solves the system of interval linear equations and returns the solution set. Uses the Gauss-Seidel method (Hansen-Sengupta version) to solve the system. Keyword precondition
to turn preconditioning off. Eldon Hansen and G. William Walster : Global Optimization Using Interval Analysis - Chapter 5 - Page 115
IntervalRootFinding.preconditioner
— MethodPreconditions the matrix A and b with the inverse of mid(A)
IntervalRootFinding.newton1d
— Methodnewton1d
performs the interval Newton method on the given function f
with its derivative f′
and initial interval x
. Optional keyword arguments give the tolerances reltol
and abstol
. reltol
is the tolerance on the relative error whereas abstol
is the tolerance on |f(X)|, and a debug
boolean argument that prints out diagnostic information.
IntervalRootFinding.newton1d
— Methodnewton1d
performs the interval Newton method on the given function f
and initial interval x
. Optional keyword arguments give the tolerances reltol
and abstol
. reltol
is the tolerance on the relative error whereas abstol
is the tolerance on |f(X)|, and a debug
boolean argument that prints out diagnostic information.
IntervalRootFinding.quadratic_helper!
— MethodHelper function for quadratic_interval that computes roots of a real quadratic using interval arithmetic to bound rounding errors.
IntervalRootFinding.quadratic_roots
— MethodFunction to solve a quadratic equation where the coefficients are intervals. Returns an array of intervals of the roots. Arguments a
, b
and c
are interval coefficients of x²
, x
and 1
respectively. The interval case differs from the non-interval case in that there might be three disjoint interval roots. In the third case, one interval root extends to −∞ and another extends to +∞. This algorithm finds the set of points where F.lo(x) ≥ 0
and the set of points where F.hi(x) ≤ 0
and takes the intersection of these two sets. Eldon Hansen and G. William Walster : Global Optimization Using Interval Analysis - Chapter 8