Python API

Interval

class pyibex.Interval

Bases: pybind11_builtins.pybind11_object

An Interval represents a closed sub set of R

The docs string is taken from ibex_Interval.h source file For more information read doc from http://www.ibex-lib.org/

Examples

>>> a = Interval(1,2)
>>> a = Interval(2)
>>> a = Interval.ALL_REALS
>>> b = Interval(a)

Warning

For C++ user, recall that the assigment operator in python copies only references and not the object. a and b are referencing the same object. To make a copy used:

For instance:
>>> a = Interval(4)
>>> b = a
>>> b = a.copy()
>>> b = Interval(a)

It is posible to acess ub and lower bound with lb() and ub() but also with an array like syntax.

Example

>>> a = Interval(1,4)
>>> a.lb()
1
>>> a[0]
1

Note

Because intervals are immutable object, it is not possible to set directly the upper / lower bounds of the Interval object

assign(self: pyibex.Interval, itv: pyibex.Interval) → None

assign the value of itv to this

bisect(self: pyibex.Interval, ratio: float = 0.5) → Tuple[pyibex.Interval, pyibex.Interval]

Bisect self into two subintervals.

The Interval musts be bisectable i.e. is_bisectable() must be true.

Parameters:ratio (float) – says where to split (in [0, 1] , 0.5=middle)
Returns:two sub intervals
Return type:list (itv1, itv2)
complementary(self: pyibex.Interval) → tuple
Returns:The complementary of x.
contains(self: pyibex.Interval, x: float) → bool

True iff self contains a d.

Note: d can also be an “open bound”, i.e., infinity.
So this function is not restricted to a set-membership interpretation.
copy(self: pyibex.Interval) → pyibex.Interval

return a new objec which is the copy of self

diam(self: pyibex.Interval) → float

Return the diameter of self. By convention, 0 if self is empty.

diff(self: pyibex.Interval, y: pyibex.Interval) → tuple
Returns:xy
inflate(*args, **kwargs)

Overloaded function.

  1. inflate(self: pyibex.pyibex.Interval, arg0: float) -> pyibex.pyibex.Interval

Add [-rad,+rad] to self. :returns: Return self.

  1. inflate(self: pyibex.pyibex.Interval, arg0: float, arg1: float) -> pyibex.pyibex.Interval

Add [-rad,+rad] to self. :returns: Return self.

interior_contains(self: pyibex.Interval, x: float) → bool

True iff the interior of self contains a d.

intersects(self: pyibex.Interval, x: pyibex.Interval) → bool

True iff self and a x intersect.

is_bisectable(self: pyibex.Interval) → bool

True iff self can be bisected into two non-degenerated intervals.

Examples of non bisectable intervals are [0,next_float(0)] or [DBL_MAX,+oo).

is_degenerated(self: pyibex.Interval) → bool

True iff self is degenerated.

An interval is degenerated if it is of the form [a, a]

Note

An empty interval is considered here as degenerated.

is_disjoint(self: pyibex.Interval, x: pyibex.Interval) → bool

True iff self and a x do not intersect.

is_empty(self: pyibex.Interval) → bool

True iff self is empty.

is_interior_subset(self: pyibex.Interval, x: pyibex.Interval) → bool

Returns: bool: True iff this interval is in the interior of x.

Note

In particular, (-oo,oo) is in the interior of (-oo,oo) and the empty set is in the interior of the empty set.

Note

Always return true if self is empty.

is_strict_interior_subset(self: pyibex.Interval, x: pyibex.Interval) → bool
Returns:True iff this interval is in the interior of a x and different from x.
Return type:bool

Note

In particular, (-oo,oo) is not “strictly” in the interior of (-oo,oo) and the empty set is not “strictly” in the interior of the empty set.

is_strict_subset(self: pyibex.Interval, x: pyibex.Interval) → bool

True iff this interval is a subset of x and not x itself.

Note

In particular, (-oo,oo) is not a strict subset of (-oo,oo) and the empty set is not a strict subset of the empty set although in both cases, the first is inside the interior of the second.

is_strict_superset(self: pyibex.Interval, x: pyibex.Interval) → bool

True iff this interval is a superset of a x different from x. See is_strict_subset(const Interval&) const.

is_subset(self: pyibex.Interval, x: pyibex.Interval) → bool

True iff this interval is a subset of a x.

Note

Always return true if self is empty.

is_superset(self: pyibex.Interval, x: pyibex.Interval) → bool

True iff this interval is a superset of a x. .. note:: Always return true if x is empty.

is_unbounded(self: pyibex.Interval) → bool
Returns:True if one bound of self is infinite.

Note

An empty interval is always bounded.

lb(self: pyibex.Interval) → float

return the upper bound

mag(self: pyibex.Interval) → float

Returns the magnitude of self

mag(self)=max(|lower bound|, |upper bound|).
mid(self: pyibex.Interval) → float

Returns the midpoint of self.

The return point is guaranteed to be included in self but not necessarily to be the closest floating point from the real midpoint. Cases are:

  • emptyset -> Quiet NaN
  • [-oo, +oo] -> midP = 0.0
  • [-oo, b] -> midP = -MAXREAL
  • [a, +oo] -> midP = MAXREAL
  • [a, b] -> midP ~ a + .5*(b-a)
mig(self: pyibex.Interval) → float

Returns the mignitude of self:

  • +(lower bound) if self > 0
  • -(upper bound) if self < 0
  • 0 otherwise
overlaps(self: pyibex.Interval, x: pyibex.Interval) → bool
True iff self and a x intersect and the intersection
has a non-null volume. Equivalently, some interior points (of this or x) must belong to the intersection.
rad(self: pyibex.Interval) → float

Return the diameter of self. By convention, 0 if self is empty.

rel_distance(self: pyibex.Interval, x: pyibex.Interval) → float

Relative Hausdorff distance between self and x.

The relative distance is basically distance(x)/diam(self).

See:
#ibex::distance (const ibex::Interval &x1, const ibex::Interval &x2).
set_empty(self: pyibex.Interval) → None

set self to empty set

ub(self: pyibex.Interval) → float

return the lower bound

ALL_REALS = [ ENTIRE ]
EMPTY_SET = [ empty ]
HALF_PI = [1.5708, 1.5708]
NEG_REALS = [-inf, 0]
ONE = [1, 1]
PI = [3.14159, 3.14159]
POS_REALS = [0, inf]
TWO_PI = [6.28319, 6.28319]
ZERO = [0, 0]

IntervalVector

class pyibex.IntervalVector

Bases: pybind11_builtins.pybind11_object

An IntervalVector is a cartesian product of Intervals The docs string is taken from ibex_IntervalVector.h source file For more information read doc from http://www.ibex-lib.org/

An IntervalVector can be initialized with:

  • a int (i.e. size) only ( all elements are initialized to [-oo,oo])
  • a int (i.e. size) and a value for all elements
  • another IntervalVector (make a copy)
  • a list of bounds [[lb_1, ub_1], […], [lb_n, ub_n]]
  • a tuple of intervals

Example

>>> a = IntervalVector(2)
([ ENTIRE ] ; [ ENTIRE ])
>>> a = IntervalVector(2, Interval(1,5))
([1, 5] ; [1, 5])
>>> a = IntervalVector(2, [1,5])
([1, 5] ; [1, 5])
>>> a = IntervalVector([[1,2], [3,4]])
([1, 2] ; [3, 4])
>>> a = IntervalVector( (Interval(2,4), Interval(5,6) ) )
([2, 4] ; [5, 6])
assign(self: pyibex.IntervalVector, arg0: pyibex.IntervalVector) → None
bisect(self: pyibex.IntervalVector, i: int, ratio: float = 0.5) → Tuple[pyibex.IntervalVector, pyibex.IntervalVector]

Bisect the box The box is bisected along the dimension a i and with a ratio a ratio. If (self)[i] is the interval [a,a+d]:

  • The first box of the result is (self)[0]x…x(self)[i-1]x[a+ratio*d]x…
  • The second box is (self)[0]x…x(self)[i-1]x[a+ratio*d,a+d]x…
Parameters:
  • i (int) – dimension to besect along
  • radio (float) – ratio Default is 0.5. Need to belong to [0, 1]
Returns:

two sub intervalvectors.

Return type:

list<IntervalVector>

clear(self: pyibex.IntervalVector) → None

Set all the elements to 0 (even if empty). .. note:: Emptiness is “overridden”.

complementary(self: pyibex.IntervalVector) → List[pyibex.IntervalVector]

Return the complementary of self. Store the complementary under the form of a union of non-overlapping IntervalVectors, into a result, and return the size of the union.

If (self) is the empty set with n components, the complementary of (self) is the n-dimensional box (-oo,oo)x…(-oo,oo).

If the complementary is empty, a result is an array of one element set to the empty box. It is not a zero-sized array containing no element.

Returns:complementaries boxes.
Return type:list<IntervalVector>

Example

>>> a = IntervalVector([[2,3], [-2, 4]])
([2, 3] ; [-2, 4])
>>> a.complementary()
[([-inf, 2] ; [ ENTIRE ]),
 ([3, inf] ; [ ENTIRE ]),
 ([2, 3] ; [-inf, -2]),
 ([2, 3] ; [4, inf])]
contains(self: pyibex.IntervalVector, x: Vector) → bool

True iff this interval vector contains a x. Dimension of a x must be equal to the dimension of (self). :param x: to be compared :type x: IntervalVector

See also

pyibex.pyibex.Interval.contains.

copy(self: pyibex.IntervalVector) → pyibex.IntervalVector
diam(self: pyibex.IntervalVector) → Vector

Returns the vector of diameters.

diff(self: pyibex.IntervalVector, y: pyibex.IntervalVector, compactness: bool = True) → List[pyibex.IntervalVector]

Return self \ y (set difference).

Returns:the difference under the form of a union of non-overlapping IntervalVectors. The lenght of the list is <= 2*n
Return type:list<IntervalVector>

Note

If the difference is empty, result is an array of one element set to the empty box. It is <b>not</b> a zero-sized array containing no element.

Parameters:y (IntervalVector) – to be compared
static empty(n: int) → pyibex.IntervalVector

Create [empty; …; empty]

Create an empty IntervalVector of dimension a n (all the components being empty Intervals)

Args:
n (int): dimension > 0
extr_diam_index(self: pyibex.IntervalVector, min: bool) → int

Return the index of a component with minimal/maximal diameter. :param min: true => minimal diameter :type min: bool

Raises:InvalidIntervalVectorOp if the IntervalVector is empty.
inflate(*args, **kwargs)

Overloaded function.

  1. inflate(self: pyibex.pyibex.IntervalVector, rad: float) -> pyibex.pyibex.IntervalVector

Add [-rad,+rad] to all the components of self.

Parameters:rad (float) – radius to use
Returns:self
  1. inflate(self: pyibex.pyibex.IntervalVector, rad: float, chi: float) -> pyibex.pyibex.IntervalVector

Add [-rad,+rad] to all the components of self.

Parameters:rad (float) – radius to use
Returns:self
init(self: pyibex.IntervalVector, x: pyibex.Interval) → None

Set all the elements to x (even if empty). :param x: value for elements :type x: Interval

Note

Emptiness is “overridden”.

interior_contains(self: pyibex.IntervalVector, x: Vector) → bool

True iff a x is in the interior of this interval vector Dimension of a x must be equal to the dimension of (self). :param x: to be compared :type x: IntervalVector

See also

pyibex.pyibex.Interval.strictly_contain : same for Interval

intersects(self: pyibex.IntervalVector, x: pyibex.IntervalVector) → bool

True iff this interval vector intersects a x. Dimension of a x must be equal to the dimension of (self). :param x: to be compared :type x: IntervalVector

See also

Interval.intersects.

is_bisectable(self: pyibex.IntervalVector) → bool
Returns:True iff self can be bisected along one dimension.
Return type:bool

See also

Interval.is_bisectable.

is_disjoint(self: pyibex.IntervalVector, x: pyibex.IntervalVector) → bool

True iff this interval vector does not intersect a x. Dimension of a x must be equal to the dimension of (self). :param x: to be compared :type x: IntervalVector

See also

Interval.is_disjoint.

is_empty(self: pyibex.IntervalVector) → bool

Return true iff this IntervalVector is empty

is_flat(self: pyibex.IntervalVector) → bool

Return true iff this IntervalVector is flat.

An IntervalVector is “flat” if the radius is 0 on at least one dimension An empty interval vector is considered as flat.
is_interior_subset(self: pyibex.IntervalVector, x: pyibex.IntervalVector) → bool

Return True iff this interval vector is a subset of the interior of a x. Dimension of a x must be equal to the dimension of this vector. .. note:: Always return true if this interval vector is empty.

See also

Interval.is_interior_subse

is_strict_interior_subset(self: pyibex.IntervalVector, x: pyibex.IntervalVector) → bool

Return True iff this interval vector is a subset of the interior of a x and different from x. Dimension of a x must be equal to the dimension of this vector. :param x: to be compared :type x: IntervalVector

Note

Always return true if this interval vector is empty.

Sea Also:
Interval.is_interior_subset.
is_strict_subset(self: pyibex.IntervalVector, x: pyibex.IntervalVector) → bool

Return Return True iff this interval vector is inside the interior of a x. Dimension of a x must be equal to the dimension of this vector. .. note:: return true if this interval vector is empty and a x not.

See also

Interval.is_strict_subset.

is_strict_superset(self: pyibex.IntervalVector, x: pyibex.IntervalVector) → bool

Return True iff a x is inside the interior of (self). Dimension of a x must be equal to the dimension of this vector. :param x: to be compared :type x: IntervalVector

Note

return true if x is empty and not (self).

See also

Interval.is_strict_superset.

is_subset(self: pyibex.IntervalVector, x: pyibex.IntervalVector) → bool

Return True iff this interval vector is a subset of a x. Dimension of a x must be equal to the dimension of this vector. .. note:: Always return true if this interval vector is empty.

See also

#Interval.is_subset.

is_superset(self: pyibex.IntervalVector, x: pyibex.IntervalVector) → bool

Return True iff this interval vector is a superset of a x. Dimension of a x must be equal to the dimension of this vector. :param x: to be compared :type x: IntervalVector

Note

Always return true if a x is empty.

See also

Interval.is_superset.

is_unbounded(self: pyibex.IntervalVector) → bool

Return True iff this interval vector contains an infinite bound.

Note

An empty interval vector is always bounded.

is_zero(self: pyibex.IntervalVector) → bool
Returns:True iff self is a vector of zeros.
Return type:bool
lb(self: pyibex.IntervalVector) → Vector

Return the lower bound vector .. note:: (self) must be nonempty

mag(self: pyibex.IntervalVector) → Vector

Return the magnitude vector. .. note:: (self) must be nonempty

max_diam(self: pyibex.IntervalVector) → float

Return the maximal diameter among all the components. :raises: InvalidIntervalVectorOp if the IntervalVector is empty.

mid(self: pyibex.IntervalVector) → Vector

Return the midpoint .. note:: (self) must be nonempty

mig(self: pyibex.IntervalVector) → Vector

Return the mignitude vector. .. note:: (self) must be nonempty

min_diam(self: pyibex.IntervalVector) → float

Return the minimal diameter among all the components. :raises: InvalidIntervalVectorOp if the IntervalVector is empty

overlaps(self: pyibex.IntervalVector, x: pyibex.IntervalVector) → bool

True iff this interval vector intersects a x and the intersection has a non-null volume Dimension of a x must be equal to the dimension of (self). :param x: to be compared :type x: IntervalVector

See also

Interval.strictly_intersects.

perimeter(self: pyibex.IntervalVector) → float

Return the perimeter of this interval vector. .. note:: Return c POS_INFINITY if unbounded.

put(self: pyibex.IntervalVector, arg0: int, arg1: pyibex.IntervalVector) → None

Put a subvector into self at a given position. :param start_index - the position where the subvector: :param subvec: the subvector :type subvec: IntervalVector

Note

(self) must not be empty

rad(self: pyibex.IntervalVector) → Vector
Returns:Vector of radii.
Return type:list
rel_distance(self: pyibex.IntervalVector, x: pyibex.IntervalVector) → float

Return the relative distance with x. :param x: :type x: IntervalVector

Returns:f$displaystyle max_{i=1..n} rel_distance([this]_i, x_i)/diam([this]_i)f$.
Sea also:
distance Interval.rel_distance
resize(self: pyibex.IntervalVector, n: int) → None

Resize this IntervalVector. If the size is increased, the existing components are not modified and the new ones are set to (-inf,+inf), even if (self) is the empty Interval (however, in this case, the status of (self) remains “empty”).

Parameters:n (int) – size of the new vector
set_empty(self: pyibex.IntervalVector) → None
Set this IntervalVector to the empty IntervalVector
The dimension remains the same.
size(*args, **kwargs)

Overloaded function.

  1. size(self: pyibex.pyibex.IntervalVector) -> int
  2. size(self: pyibex.pyibex.IntervalVector) -> int
Returns:The dimension (number of components)
sort_indices(self: pyibex.IntervalVector, min: bool) → List[int]

Return the indices of all the components, sorted by increasing/decreasing diameter. :returns: list of indices sorted by increasing/decreasing diameter. :rtype: list

Parameters:min (bool) – true => first indice corresponds to smallest diameter
subvector(self: pyibex.IntervalVector, start_index: int, end_index: int) → pyibex.IntervalVector

Return a subvector :returns: [ (self)[start_index]; …; (self)[end_index] ].

Parameters:
  • start_index (int) – first index
  • end_index (int) – last index

Note

self must not be empty

tolist(self: pyibex.IntervalVector) → List[float]

Convert an IntervalVector into a list

Returns:[x1_lb, x1_ub, x2_lb, x2_ub, …., ]
Return type:list<double>
ub(self: pyibex.IntervalVector) → Vector

Return the upper bound vector .. note:: (self) must be nonempty

/
volume(self: pyibex.IntervalVector) → float

Return the volume of this interval vector. .. note:: Return c POS_INFINITY if the vector is unbounded and not flat.

Note

Return 0 if the vector is flat and not unbounded.

Warning

If the interval vector is both flat and unbounded, the result is undefined.

See also

flat() unbounded()

width(self: pyibex.IntervalVector) → float

Function

class pyibex.Function

Bases: pybind11_builtins.pybind11_object

backward(*args, **kwargs)

Overloaded function.

  1. backward(self: pyibex.pyibex.Function, arg0: pyibex.pyibex.Interval, arg1: pyibex.pyibex.IntervalVector) -> bool
  2. backward(self: pyibex.pyibex.Function, arg0: pyibex.pyibex.IntervalVector, arg1: pyibex.pyibex.IntervalVector) -> bool
  3. backward(self: pyibex.pyibex.Function, arg0: ibex::IntervalMatrix, arg1: pyibex.pyibex.IntervalVector) -> bool
diff(self: pyibex.Function) → pyibex.Function
eval(*args, **kwargs)

Overloaded function.

  1. eval(self: pyibex.pyibex.Function, arg0: pyibex.pyibex.IntervalVector) -> pyibex.pyibex.Interval
  2. eval(self: pyibex.pyibex.Function, arg0: pyibex.pyibex.Interval) -> pyibex.pyibex.IntervalVector
eval_matrix(self: pyibex.Function, arg0: pyibex.IntervalVector) → ibex::IntervalMatrix
eval_vector(self: pyibex.Function, arg0: pyibex.IntervalVector) → pyibex.IntervalVector
nb_arg(self: pyibex.Function) → int

Contractors (Ctc)

class pyibex.Ctc

Bases: pybind11_builtins.pybind11_object

Contractor abscract class All contractor need to inherit from this class and define the method contract and initilize the Ctc class with the dimension of the contractor i.e. the dimension of the box passed to contract function.

In python you can define yout own contractor .. rubric:: Example

>>> #  contractor on a 2 dimensionnal boxes
>>> class myCtc(Ctc):
>>> def __init__(self):
>>>     Ctc.__init__(self, 2) # init Ctc abstract class
>>>
>>> def contract(self, box):
>>>     pass

You can also make the intersection / union of contractors using & and |

Example

>>> ctc1 = CtcXXX
>>> ctc2 = CtcXXX
>>> # intersection
>>> ctc3 = ctc1 & ctc2 # intersection
>>> ctc4 = ctc1 | ctc3 & ctc2 # combinaison
nb_var

Return the number of variables this contractor works with.

CtcUnion

class pyibex.CtcUnion

Bases: pyibex.pyibex.Ctc

Union of contractors For a box [x] the union of {c_0,…c_n} performs \(c_1([x]) \cup ... \cup c_n([x])\)

Parameters:list of contractors (list<Ctc>) –
Returns:
Return type:CtcUnion
contract(self: pyibex.CtcUnion, arg0: pyibex.IntervalVector) → None

CtcCompo

class pyibex.CtcCompo

Bases: pyibex.pyibex.Ctc

Intersection (composition) of contractors For a box [x] the composition of {c_0,…c_n} performs c_n(…(c_1(c_0([x])))).

Parameters:list of contractors (list<Ctc>) –
Returns:
Return type:CtcCompo
contract(self: pyibex.CtcCompo, arg0: pyibex.IntervalVector) → None

CtcFwdBwd

class pyibex.CtcFwdBwd

Bases: pyibex.pyibex.Ctc

Forward-backward contractor (HC4Revise).

Can be initialized with:
  • a function ( constraint f(x) cmp 0) default cmp = EQ
  • a function and a interval/box f(x) in itv
Parameters:
  • f (Function) – Function used for the Forward-backward
  • op (CmpOp) – comparaison operator
  • itv_y (Interval) – for constraint f(x) in itv_y
  • box_y (IntervalVector) – for constraint f(x) in box_y
Returns:

a Forward-backward contractor

Return type:

CtcFwdBwd

Example

>>> f = Function("x","y", "(sin(x)^2 + 3*y)")
>>> # Constaint f(x) <= O
>>> ctc = CtcFwdBwd(f, LEQ)
>>> # Constaint f(x) <= [0, 1]
>>> ctc = CtcFwdBwd(f, Interval(0, 1))
contract(self: pyibex.CtcFwdBwd, arg0: pyibex.IntervalVector) → None

CtcInverse

class pyibex.CtcInverse

Bases: pyibex.pyibex.Ctc

Image of a contractor by a function in a forward-backaward manner [x] = [x] cap f^{-1}( ctc( f([x]) ) )

Parameters:
  • ctc (Ctc) – contrctor to use
  • f (Function) – function to use
Returns:

contractor

Return type:

CtcInverse

contract(self: pyibex.CtcInverse, arg0: pyibex.IntervalVector) → None

CtcNotIn

class pyibex.CtcNotIn

Bases: pyibex.pyibex.Ctc

Contractor for f(x) not-in y wrt to f, y can be an Interval or a box. .. note:: Only Interval and IntervalVector types are supported for now.

contract(self: pyibex.CtcNotIn, arg0: pyibex.IntervalVector) → None

CtcFixPoint

class pyibex.CtcFixPoint

Bases: pyibex.pyibex.Ctc

Fix point of a contractor When the Hausdorff distance between two iterations is less than ratio*diameter the fix-point is considered to be reached. :param ctc: a contractor :type ctc: Ctc :param ratio: stop ratio between two iterationl :type ratio: double

contract(self: pyibex.CtcFixPoint, arg0: pyibex.IntervalVector) → None

CtcQInter

class pyibex.CtcQInter

Bases: pyibex.pyibex.Ctc

Exact relaxed intersection of a list of contractors :param list: list of contractors :type list: list<Ctc> :param q: The number of contractors which have to intersect the result :type q: int

contract(self: pyibex.CtcQInter, arg0: pyibex.IntervalVector) → None

Separators (Sep)

class pyibex.Sep

Bases: pybind11_builtins.pybind11_object

Separator interface. A separator is an operator that performs two independent and complementary contractions. The separator is associated with a set (noted S) and the first contraction (called “inner”) removes points inside S. The second contraction (called “outer”) removes points outside S. In summary: Given a box [x], the separator produces two sub-boxes [x_in] and [x_out] that verify:

  • ([x] \cap [x_in]) \subset S
  • ([x] \cap [x_out]) \cap S = \emptyset

For efficiency reasons, the separate(…) function takes only two input-output arguments, x_in and x_out, each containing initially a copy of the box [x]. A separator can also be viewed a as pair of contractors.

See also

pyibex.SepCtcPair.

See also

L. Jaulin and B. Desrochers (2014). “Introduction to the Algebra of Separators with Application to Path Planning”. Engineering Applications of Artificial Intelligence volume 33, pp. 141-147.

All separator need to inherit from this class and define the method separate and initilize the Sep class with the dimension of the separator i.e. the dimension of the box passed to separate function.

In python you can define yout own separator

Example

>>> #  separator on a 2 dimensionnal boxes
>>> class mySep(Sep):
>>>     def __init__(self):
>>>         Sep.__init__(self, 2) # init Sep abstract class
>>>     def separate(self, box_in, box_out):
>>>         pass

You can also make the intersection / union / restriction of separators using &, | and ~

Example

>>> sep1 = SepXXX
>>> sep2 = SepXXX
>>> # intersection
>>> sep3 = sep1 & sep2 # intersection
>>> sep4 = sep1 | sep3 & sep2 # combinaison
>>> sep5 =   ~sep4 | spe1
separate(self: pyibex.Sep, x_in: pyibex.IntervalVector, x_out: pyibex.IntervalVector) → None

Separate a box in two sub-boxes.

Parameters:
  • x_in (IntervalVector) – As input: the initial box. As output: result of the first (“inner”) contraction
  • x_out (IntervalVector) – As input: the initial box. As output: the result of the second (“outer”) contraction

Note

Precondition: x_in and x_out must be the same boxes.

nb_var

The number of variables this separator works with.

SepUnion

class pyibex.SepUnion

Bases: pyibex.pyibex.Sep

Union of separators

Parameters:list<Sep> – list of separators
separate(self: pyibex.SepUnion, arg0: pyibex.IntervalVector, arg1: pyibex.IntervalVector) → None

SepInter

class pyibex.SepInter

Bases: pyibex.pyibex.Sep

Intersection of separators

Parameters:list<Sep> – list of separators
separate(self: pyibex.SepInter, arg0: pyibex.IntervalVector, arg1: pyibex.IntervalVector) → None

SepCtcPair

class pyibex.SepCtcPair

Bases: pyibex.pyibex.Sep

Build a separator with a pair of complementary contractors
Parameters:
  • ctc_in (Ctc) – inner contractor
  • ctc_out (Ctc) – outer contractor

Example

>>> f = Function("x[2]",  "x[0] + 2*x[1]")
>>> ctc_out = CtcFwdBwd(f, LEQ)
>>> ctc_in = CtcFwdBwd(f, GEQ)
>>> sep = SepCtcPair(ctc_in, ctc_out)
ctc_in(self: pyibex.SepCtcPair) → pyibex.Ctc
ctc_out(self: pyibex.SepCtcPair) → pyibex.Ctc
separate(self: pyibex.SepCtcPair, arg0: pyibex.IntervalVector, arg1: pyibex.IntervalVector) → None

SepFwdBwd

class pyibex.SepFwdBwd

Bases: pyibex.pyibex.SepCtcPair

This separator applies inner and outer forward-backward contractors to separate a box w.r.t. a constraint \(f(x) in [y]\).

Parameters:
  • f (Function) – Function used for the Forward-backward
  • op (CmpOp) – comparaison operator
  • itv_y (Interval) – for constraint f(x) in itv_y
  • box_y (IntervalVector) – for constraint f(x) in box_y

Example

>>> # disk of raduis 2
>>> f = Function("x", "y", "x^2+y^2 - 4")
>>> sep = SepFwdBwd(f, LEQ)
separate(self: pyibex.SepFwdBwd, arg0: pyibex.IntervalVector, arg1: pyibex.IntervalVector) → None

SepNot

class pyibex.SepNot

Bases: pyibex.pyibex.Sep

Separator for \(f(x) \notin [y]\)

See also

pyibex.CtCNotIn

separate(self: pyibex.SepNot, arg0: pyibex.IntervalVector, arg1: pyibex.IntervalVector) → None

SepInverse

class pyibex.SepInverse

Bases: pyibex.pyibex.Sep

Image of a separator by a function in a forward-backward manner

Parameters:
  • sep (Sep) – Separator to use
  • f (Function) – Function to use
See Also
pyibex.CtcInverse
separate(self: pyibex.SepInverse, arg0: pyibex.IntervalVector, arg1: pyibex.IntervalVector) → None

SepQInter

class pyibex.SepQInter

Bases: pyibex.pyibex.Sep

Exact relaxed intersection of a list of separators

The complexity of this algorithm is exponential wrt the dimension of the input box.

Parameters:
  • list (list<Ctc>) – list of separators
  • q (int) – The number of separators which have to intersect the result
separate(self: pyibex.SepQInter, arg0: pyibex.IntervalVector, arg1: pyibex.IntervalVector) → None

Geometry