ZonoOpt Module
Classes and tailored optimization routines for zonotopes, constrained zonotopes, and hybrid zonotopes.
See the github page for documentation: https://github.com/psu-PAC-Lab/ZonoOpt
More information about ZonoOpt can be found in the the following publication. Please cite this if you publish work based on ZonoOpt: Robbins, J.A., Siefert, J.A., and Pangborn, H.C., “Sparsity-Promoting Reachability Analysis and Optimization of Constrained Zonotopes,” 2025.**
- class zonoopt.Box(self: zonoopt._core.Box, x_lb: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]'], x_ub: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]'])
Bases:
pybind11_objectBox (i.e., interval vector) class
Constructor from intervals of lower and upper bounds
- Parameters:
x_lb (numpy.array) – vector of lower bounds
x_ub (numpy.array) – vector of upper bounds
- __add__(self: zonoopt._core.Box, other: zonoopt._core.Box) zonoopt._core.Box
Elementwise addition
- __getitem__(self: zonoopt._core.Box, i: SupportsInt) zonoopt._core.Interval
Get interval at index i
- Parameters:
i (int) – index
- Returns:
interval at index i in Box
- Return type:
- __mul__(*args, **kwargs)
Overloaded function.
__mul__(self: zonoopt._core.Box, other: zonoopt._core.Box) -> zonoopt._core.Box
Elementwise multiplication
- Args:
other (Box): rhs box
- Returns:
Box: self * other (elementwise)
__mul__(self: zonoopt._core.Box, alpha: typing.SupportsFloat) -> zonoopt._core.Box
Elementwise multiplication with scalar
- Args:
alpha (float): scalar multiplier
- Returns:
Box: alpha * self (elementwise)
- __setitem__(self: zonoopt._core.Box, i: SupportsInt, val: zonoopt._core.Interval) None
Set indexed interval in box to specified value
- Parameters:
i (int) – index
val (Interval) – new interval for index i in Box
- __sub__(self: zonoopt._core.Box, other: zonoopt._core.Box) zonoopt._core.Box
Elementwise subtraction
- __truediv__(self: zonoopt._core.Box, other: zonoopt._core.Box) zonoopt._core.Box
Elementwise division
- center(self: zonoopt._core.Box) Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']
Gets center of box (x_ub + x_lb) / 2
- Returns:
center of interval
- Return type:
numpy.array
- contract(self: zonoopt._core.Box, A: scipy.sparse.csr_matrix[numpy.float64], b: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]'], iter: SupportsInt) bool
Interval contractor.
Executes a forward-backward interval contractor for the equality constraint A*x=b. For points x in the box, this shrinks the box without removing any points x that satisfy A*x=b. If the contractor detects that the box does not intersect A*x=b, then this function will return false.
- Parameters:
A (scipy.sparse.csr_matrix) – constraint matrix
b (numpy.vector) – constraint vector
iter (int) – number of contractor iterations
- Returns:
flag indicating that the contractor did not detect that A*x=b and the box do not intersect
- Return type:
bool
- copy(self: zonoopt._core.Box) zonoopt._core.Box
Copies Box object
- Returns:
copy of object
- Return type:
- dot(self: zonoopt._core.Box, x: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]']) zonoopt._core.Interval
Linear map with vector
- Parameters:
x (numpy.array) – vector
- Returns:
result of linear map of box with vector
- Return type:
- linear_map(self: zonoopt._core.Box, A: scipy.sparse.csr_matrix[numpy.float64]) zonoopt._core.Box
Linear map of box based on interval arithmetic
- Parameters:
A (scipy.sparse.csr_matrix) – linear map matrix
- Returns:
linear mapped box
- Return type:
- lower(self: zonoopt._core.Box) Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']
Get reference to lower bounds
- Returns:
lower bounds
- Return type:
numpy.array
- project(self: zonoopt._core.Box, x: Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]', 'flags.writeable']) None
Projects vector onto the Box (in place)
- Parameters:
x (numpy.array) – vector to be projected
- size(self: zonoopt._core.Box) int
Get size of Box object
- Returns:
size of box
- Return type:
int
- upper(self: zonoopt._core.Box) Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']
Get reference to upper bounds
- Returns:
upper bounds
- Return type:
numpy.array
- width(self: zonoopt._core.Box) float
Get width of box.
Specifically, this returns the sum of the widths of each interval in the box
- Returns:
width of box
- Return type:
float
- class zonoopt.ConZono(self: zonoopt._core.ConZono, G: scipy.sparse.csc_matrix[numpy.float64], c: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]'], A: scipy.sparse.csc_matrix[numpy.float64], b: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]'], zero_one_form: bool = False)
Bases:
HybZonoConstrained zonotope class
A constrained zonotope is defined as: Z = {G xi + c | A xi = b, xi in [-1, 1]^nG}. Equivalently, the following shorthand can be used: Z = <G, c, A, b>. Optionally, in 0-1 form, the factors are xi in [0,1]. The set dimension is n, and the number of equality constraints is nC.
ConZono constructor
- Parameters:
G (scipy.sparse.csc_matrix) – generator matrix
c (numpy.array) – center
A (scipy.sparse.csc_matrix) – constraint matrix
b (numpy.array) – constraint vector
zero_one_form (bool, optional) – true if set is in 0-1 form
- constraint_reduction(self: zonoopt._core.ConZono) None
Execute constraint reduction algorithm from Scott et. al. 2016
Removes one constraint and one generator from the constrained zonotope. The resulting set is an over-approximation of the original set.
- set(self: zonoopt._core.ConZono, G: scipy.sparse.csc_matrix[numpy.float64], c: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]'], A: scipy.sparse.csc_matrix[numpy.float64], b: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]'], zero_one_form: bool = False) None
Reset constrained zonotope object with the given parameters.
- Parameters:
G (scipy.sparse.csc_matrix) – generator matrix
c (numpy.array) – center
A (scipy.sparse.csc_matrix) – constraint matrix
b (numpy.array) – constraint vector
zero_one_form (bool, optional) – true if set is in 0-1 form
- to_zono_approx(self: zonoopt._core.ConZono) ZonoOpt::Zono
Compute outer approximation of constrained zonotope as zonotope using SVD
- Returns:
Zonotope over-approximation
- Return type:
- class zonoopt.EmptySet(self: zonoopt._core.EmptySet, n: SupportsInt)
Bases:
ConZonoEmpty Set class
Used to facilitate set operations with trivial solutions when one of the sets is an empty set.
EmptySet constructor
- Parameters:
n (int) – dimension
- class zonoopt.HybZono(self: zonoopt._core.HybZono, Gc: scipy.sparse.csc_matrix[numpy.float64], Gb: scipy.sparse.csc_matrix[numpy.float64], c: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]'], Ac: scipy.sparse.csc_matrix[numpy.float64], Ab: scipy.sparse.csc_matrix[numpy.float64], b: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]'], zero_one_form: bool = False, sharp: bool = False)
Bases:
pybind11_objectHybrid zonotope class
A hybrid zonotope is defined as: Z = {Gc * xi_c + Gb * xi_b + c | Ac * xi_c + Ab * xi_b = b, xi_c in [-1, 1]^nGc, xi_b in {-1, 1}^nGb}. Equivalently, the following shorthand can be used: Z = <Gc, Gb, c, Ac, Ab, b>. Optionally, in 0-1 form, the factors are xi_c in [0, 1]^nGc, xi_b in {0, 1}^nGb. The set dimension is n, and the number of equality constraints is nC.
HybZono constructor
- Parameters:
Gc (scipy.sparse.csc_matrix) – continuous generator matrix
Gb (scipy.sparse.csc_matrix) – binary generator matrix
c (numpy.array) – center
Ac (scipy.sparse.csc_matrix) – continuous constraint matrix
Ab (scipy.sparse.csc_matrix) – binary constraint matrix
b (numpy.array) – constraint vector
zero_one_form (bool, optional) – true if set is in 0-1 form
sharp (bool, optional) – true if set is known to be sharp, i.e., convex relaxation = convex hull
- bounding_box(self: zonoopt._core.HybZono, settings: zonoopt._core.OptSettings = OptSettings structure: verbose: false verbosity_interval: 100 t_max: 1.79769e+308 k_max_admm: 5000 rho: 10 eps_dual: 0.01 eps_prim: 0.001 k_inf_check: 10 inf_norm_conv: true use_interval_contractor: true contractor_iter: 1 search_mode: 0 polish: 1 eps_dual_search: 0.1 eps_prim_search: 0.01 eps_r: 0.01 eps_a: 0.1 k_max_bnb: 100000 n_threads_bnb: 4 max_nodes: 100000 contractor_tree_search_depth: 10, solution: zonoopt._core.OptSolution = None) zonoopt._core.Box
Computes a bounding box of the set object as a Box object.
- Parameters:
settings (OptSettings, optional) – optimization settings structure
solution (OptSolution, optional) – optimization solution structure pointer, populated with result
- Returns:
bounding box of the set
- Return type:
In general, solves 2*n support optimizations where n is the set dimension to compute a bounding box.
- complement(self: zonoopt._core.HybZono, delta_m: typing.SupportsFloat = 100, remove_redundancy: bool = True, settings: zonoopt._core.OptSettings = OptSettings structure: verbose: false verbosity_interval: 100 t_max: 1.79769e+308 k_max_admm: 5000 rho: 10 eps_dual: 0.01 eps_prim: 0.001 k_inf_check: 10 inf_norm_conv: true use_interval_contractor: true contractor_iter: 1 search_mode: 0 polish: 1 eps_dual_search: 0.1 eps_prim_search: 0.01 eps_r: 0.01 eps_a: 0.1 k_max_bnb: 100000 n_threads_bnb: 4 max_nodes: 100000 contractor_tree_search_depth: 10, solution: zonoopt._core.OptSolution = None, n_leaves: typing.SupportsInt = 2147483647, contractor_iter: typing.SupportsInt = 100) zonoopt._core.HybZono
Computes the complement of the set Z.
- Parameters:
delta_m (float, optional) – parameter defining range of complement
remove_redundancy (bool, optional) – remove redundant constraints and unused generators in get_leaves function call
settings (OptSettings, optional) – optimization settings for get_leaves function call
solution (OptSolution, optional) – optimization solution for get_leaves function call
n_leaves (int, optional) – maximum number of leaves to return in get_leaves function call
contractor_iter (int, optional) – number of interval contractor iterations in remove_redundancy if using
- Returns:
Hybrid zonotope complement of the given set
- Return type:
Computes the complement according to the method of Bird and Jain: “Unions and Complements of Hybrid Zonotopes” delta_m is a parameter that defines the set over which the complement is defined. For a constrained zonotope, the complement is restricted to the set X = {G xi + c | A xi = b, xi in [-1-delta_m, 1+delta+m]^{nG}}.
- contains_point(self: zonoopt._core.HybZono, x: typing.Annotated[numpy.typing.ArrayLike, numpy.float64, "[m, 1]"], settings: zonoopt._core.OptSettings = OptSettings structure: verbose: false verbosity_interval: 100 t_max: 1.79769e+308 k_max_admm: 5000 rho: 10 eps_dual: 0.01 eps_prim: 0.001 k_inf_check: 10 inf_norm_conv: true use_interval_contractor: true contractor_iter: 1 search_mode: 0 polish: 1 eps_dual_search: 0.1 eps_prim_search: 0.01 eps_r: 0.01 eps_a: 0.1 k_max_bnb: 100000 n_threads_bnb: 4 max_nodes: 100000 contractor_tree_search_depth: 10, solution: zonoopt._core.OptSolution = None) bool
Checks whether the point x is contained in the set object.
- Parameters:
x (numpy.array) – point to be checked for set containment
settings (OptSettings, optional) – optimization settings structure
solution (OptSolution, optional) – optimization solution structure pointer, populated with result
- Returns:
true if set contains point, false otherwise
- Return type:
bool
False positives are possible; will return true if the optimization converges within the specified tolerances. Will return false only if an infeasibility certificate is found, i.e., false negatives are not possible.
- convert_form(self: zonoopt._core.HybZono) None
Converts the set representation between -1-1 and 0-1 forms.
This method converts the set representation between -1-1 and 0-1 forms. If the set is in -1-1 form, then xi_c in [-1,1] and xi_b in {-1,1}. If the set is in 0-1 form, then xi_c in [0,1] and xi_b in {0,1}.
- convex_relaxation(self: zonoopt._core.HybZono) ZonoOpt::ConZono
Computes the convex relaxation of the hybrid zonotope.
- Returns:
Constrained zonotope Z = <[Gc, Gb], c, [Ac, Ab,], b>
- Return type:
This method returns the convex relaxation of the hybrid zonotope. If the set is sharp, the convex relaxation is the convex hull.
- copy(self: zonoopt._core.HybZono) zonoopt._core.HybZono
Creates a copy of the hybrid zonotope object.
- Returns:
A copy of the hybrid zonotope object.
- Return type:
- get_A(self: zonoopt._core.HybZono) scipy.sparse.csc_matrix[numpy.float64]
Returns constraint matrix
- Returns:
A
- Return type:
scipy.sparse.csc_matrix
- get_Ab(self: zonoopt._core.HybZono) scipy.sparse.csc_matrix[numpy.float64]
Returns binary constraint matrix
- Returns:
Ab
- Return type:
scipy.sparse.csc_matrix
- get_Ac(self: zonoopt._core.HybZono) scipy.sparse.csc_matrix[numpy.float64]
Returns continuous constraint matrix
- Returns:
Ac
- Return type:
scipy.sparse.csc_matrix
- get_G(self: zonoopt._core.HybZono) scipy.sparse.csc_matrix[numpy.float64]
Returns generator matrix
- Returns:
G
- Return type:
scipy.sparse.csc_matrix
- get_Gb(self: zonoopt._core.HybZono) scipy.sparse.csc_matrix[numpy.float64]
Returns binary generator matrix
- Returns:
Gb
- Return type:
scipy.sparse.csc_matrix
- get_Gc(self: zonoopt._core.HybZono) scipy.sparse.csc_matrix[numpy.float64]
Returns continuous generator matrix
- Returns:
Gc
- Return type:
scipy.sparse.csc_matrix
- get_b(self: zonoopt._core.HybZono) Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']
Returns constraint vector
- Returns:
b
- Return type:
numpy.array
- get_c(self: zonoopt._core.HybZono) Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']
Returns center vector
- Returns:
c
- Return type:
numpy.array
- get_leaves(self: zonoopt._core.HybZono, remove_redundancy: bool = True, settings: zonoopt._core.OptSettings = OptSettings structure: verbose: false verbosity_interval: 100 t_max: 1.79769e+308 k_max_admm: 5000 rho: 10 eps_dual: 0.01 eps_prim: 0.001 k_inf_check: 10 inf_norm_conv: true use_interval_contractor: true contractor_iter: 1 search_mode: 0 polish: 1 eps_dual_search: 0.1 eps_prim_search: 0.01 eps_r: 0.01 eps_a: 0.1 k_max_bnb: 100000 n_threads_bnb: 4 max_nodes: 100000 contractor_tree_search_depth: 10, solution: zonoopt._core.OptSolution = None, n_leaves: typing.SupportsInt = 2147483647, contractor_iter: typing.SupportsInt = 100) list[ZonoOpt::ConZono]
Computes individual constrained zonotopes whose union is the hybrid zonotope object.
- Parameters:
remove_redundancy (bool, optional) – flag to make call to remove_redundancy for each identified leaf
settings (OptSettings, optional) – optimization settings structure
solution (OptSolution, optional) – optimization solution structure pointer, populated with result
n_leaves (int, optional) – max number of leaves to find
contractor_iter (int, optional) – number of interval contractor iterations to run if using remove_redundancy
- Returns:
vector of constrained zonotopes [Z0, Z1, …] such that Zi is a subset of the current set for all i
- Return type:
list[ConZono]
Searches for constrained zonotopes that correspond to feasible combinations of the hybrid zonotope binary variables. If the branch and bound converges (i.e., did not hit max time, max number of branch and bound iterations, or max nodes in queue) and the n_leaves argument does not stop the optimization before exhausting all possibilities, then the resulting vector of constrained zonotopes can be unioned to recover the original set. It is possible for a leaf to be the empty set if the optimization converges before detecting an infeasibility certificate.
- get_n(self: zonoopt._core.HybZono) int
Returns dimension of set
- Returns:
- Return type:
int
- get_nC(self: zonoopt._core.HybZono) int
Returns number of constraints in set definition
- Returns:
nC
- Return type:
int
- get_nG(self: zonoopt._core.HybZono) int
Returns number of generators in set definition
- Returns:
nG
- Return type:
int
- get_nGb(self: zonoopt._core.HybZono) int
Returns number of binary generators in set definition
- Returns:
nGb
- Return type:
int
- get_nGc(self: zonoopt._core.HybZono) int
Returns number of continuous generators in set definition
- Returns:
nGc
- Return type:
int
- is_0_1_form(self: zonoopt._core.HybZono) bool
Returns true if factors are in range [0,1], false if they are in range [-1,1].
- Returns:
zero_one_form flag
- Return type:
bool
- is_conzono(self: zonoopt._core.HybZono) bool
Polymorphic type checking
- Returns:
true if set is a constrained zonotope
- Return type:
bool
- is_empty(self: zonoopt._core.HybZono, settings: zonoopt._core.OptSettings = OptSettings structure: verbose: false verbosity_interval: 100 t_max: 1.79769e+308 k_max_admm: 5000 rho: 10 eps_dual: 0.01 eps_prim: 0.001 k_inf_check: 10 inf_norm_conv: true use_interval_contractor: true contractor_iter: 1 search_mode: 0 polish: 1 eps_dual_search: 0.1 eps_prim_search: 0.01 eps_r: 0.01 eps_a: 0.1 k_max_bnb: 100000 n_threads_bnb: 4 max_nodes: 100000 contractor_tree_search_depth: 10, solution: zonoopt._core.OptSolution = None) bool
Returns true if the set is provably empty, false otherwise.
- Parameters:
settings (OptSettings, optional) – optimization settings structure
solution (OptSolution, optional) – optimization solution structure pointer, populated with result
- Returns:
flag indicating whether set is provably empty
- Return type:
bool
- is_empty_set(self: zonoopt._core.HybZono) bool
Polymorphic type checking
- Returns:
true if set is a empty set object
- Return type:
bool
- is_hybzono(self: zonoopt._core.HybZono) bool
Polymorphic type checking
- Returns:
true if set is a hybrid zonotope
- Return type:
bool
- is_point(self: zonoopt._core.HybZono) bool
Polymorphic type checking
- Returns:
true if set is a point
- Return type:
bool
- is_sharp(self: zonoopt._core.HybZono) bool
Returns true if set is known to be sharp
- Returns:
sharp flag
- Return type:
bool
- is_zono(self: zonoopt._core.HybZono) bool
Polymorphic type checking
- Returns:
true if set is a zonotope
- Return type:
bool
- optimize_over(self: zonoopt._core.HybZono, P: scipy.sparse.csc_matrix[numpy.float64], q: typing.Annotated[numpy.typing.ArrayLike, numpy.float64, "[m, 1]"], c: typing.SupportsFloat = 0, settings: zonoopt._core.OptSettings = OptSettings structure: verbose: false verbosity_interval: 100 t_max: 1.79769e+308 k_max_admm: 5000 rho: 10 eps_dual: 0.01 eps_prim: 0.001 k_inf_check: 10 inf_norm_conv: true use_interval_contractor: true contractor_iter: 1 search_mode: 0 polish: 1 eps_dual_search: 0.1 eps_prim_search: 0.01 eps_r: 0.01 eps_a: 0.1 k_max_bnb: 100000 n_threads_bnb: 4 max_nodes: 100000 contractor_tree_search_depth: 10, solution: zonoopt._core.OptSolution = None) Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']
Solves optimization problem with quadratic objective over the current set
- Parameters:
P (scipy.sparse.csc_matrix) – quadratic objective matrix
q (numpy.array) – linear objective vector
c (float, optional) – constant term in objective function
settings (OptSettings, optional) – optimization settings structure
solution (OptSolution, optional) – optimization solution structure pointer, populated with result
- Returns:
point z in the current set
- Return type:
numpy.array
Solves optimization problem of the form min 0.5*z^T*P*z + q^T*z + c where z is a vector in the current set
- project_point(self: zonoopt._core.HybZono, x: typing.Annotated[numpy.typing.ArrayLike, numpy.float64, "[m, 1]"], settings: zonoopt._core.OptSettings = OptSettings structure: verbose: false verbosity_interval: 100 t_max: 1.79769e+308 k_max_admm: 5000 rho: 10 eps_dual: 0.01 eps_prim: 0.001 k_inf_check: 10 inf_norm_conv: true use_interval_contractor: true contractor_iter: 1 search_mode: 0 polish: 1 eps_dual_search: 0.1 eps_prim_search: 0.01 eps_r: 0.01 eps_a: 0.1 k_max_bnb: 100000 n_threads_bnb: 4 max_nodes: 100000 contractor_tree_search_depth: 10, solution: zonoopt._core.OptSolution = None) Annotated[numpy.typing.NDArray[numpy.float64], '[m, 1]']
Returns the projection of the point x onto the set object.
- Parameters:
x (numpy.array) – point to be projected
settings (OptSettings, optional) – optimization settings structure
solution (OptSolution, optional) – optimization solution structure pointer, populated with result
- Returns:
point z in the current set
- Return type:
numpy.array
- remove_redundancy(self: zonoopt._core.HybZono, contractor_iter: SupportsInt = 100) None
Removes redundant constraints and any unused generators
This method uses an interval contractor to detect generators that can be removed. Additionally, any linearly dependent rows of the constraint matrix A are removed. If the linearly dependent constraints are not consistent (e.g., if A = [1, 0.1; 1, 0.1] and b = [1; 0.8]), the returned set is not equivalent to the original set. Unused factors are also removed.
- Parameters:
contractor_iter (int) – number of interval contractor iterations to run
- set(self: zonoopt._core.HybZono, Gc: scipy.sparse.csc_matrix[numpy.float64], Gb: scipy.sparse.csc_matrix[numpy.float64], c: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]'], Ac: scipy.sparse.csc_matrix[numpy.float64], Ab: scipy.sparse.csc_matrix[numpy.float64], b: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]'], zero_one_form: bool = False, sharp: bool = False) None
Reset hybrid zonotope object with the given parameters.
- Parameters:
Gc (scipy.sparse.csc_matrix) – continuous generator matrix
Gb (scipy.sparse.csc_matrix) – binary generator matrix
c (numpy.array) – center
Ac (scipy.sparse.csc_matrix) – continuous constraint matrix
Ab (scipy.sparse.csc_matrix) – binary constraint matrix
b (numpy.array) – constraint vector
zero_one_form (bool) – true if set is in 0-1 form
sharp (bool) – true if set is known to be sharp, i.e., convex relaxation = convex hull
- support(self: zonoopt._core.HybZono, d: typing.Annotated[numpy.typing.ArrayLike, numpy.float64, "[m, 1]"], settings: zonoopt._core.OptSettings = OptSettings structure: verbose: false verbosity_interval: 100 t_max: 1.79769e+308 k_max_admm: 5000 rho: 10 eps_dual: 0.01 eps_prim: 0.001 k_inf_check: 10 inf_norm_conv: true use_interval_contractor: true contractor_iter: 1 search_mode: 0 polish: 1 eps_dual_search: 0.1 eps_prim_search: 0.01 eps_r: 0.01 eps_a: 0.1 k_max_bnb: 100000 n_threads_bnb: 4 max_nodes: 100000 contractor_tree_search_depth: 10, solution: zonoopt._core.OptSolution = None) float
Computes support function of the set in the direction d.
- Parameters:
d (numpy.array) – vector defining direction for support function
settings (OptSettings, optional) – optimization settings structure
solution (OptSolution, optional) – optimization solution structure pointer, populated with result
- Returns:
support value
- Return type:
float
Solves max_{z in Z} <z, d> where <., .> is the inner product
- class zonoopt.IneqTerm(self: zonoopt._core.IneqTerm, idx: SupportsInt, coeff: SupportsFloat)
Bases:
pybind11_objectStructure containing term in 0-1 inequality.
IneqTerm constructor
- property coeff
coefficient of variable
- property idx
index of variable
- class zonoopt.IneqType(self: zonoopt._core.IneqType, value: SupportsInt)
Bases:
pybind11_objectEnumeration to select inequality direction / use equality.
Members:
LESS : Strictly less than
LESS_OR_EQUAL : Less than or equal to
EQUAL : Equal to
GREATER_OR_EQUAL : Greater than or equal to
GREATER : Strictly greater than
- EQUAL = <IneqType.EQUAL: 0>
- GREATER = <IneqType.GREATER: 2>
- GREATER_OR_EQUAL = <IneqType.GREATER_OR_EQUAL: 1>
- LESS = <IneqType.LESS: -2>
- LESS_OR_EQUAL = <IneqType.LESS_OR_EQUAL: -1>
- property name
- property value
- class zonoopt.Inequality(self: zonoopt._core.Inequality, n_dims: SupportsInt)
Bases:
pybind11_objectInequality class
Constructs inequality, must specify number of dimensions.
- Parameters:
n_dims (int) – number of dimensions in the inequality
- add_term(self: zonoopt._core.Inequality, idx: SupportsInt, coeff: SupportsFloat) None
Adds a term to the inequality.
- Parameters:
idx (int) – index of variable in the inequality
coeff (float) – coefficient of the variable
- get_ineq_type(self: zonoopt._core.Inequality) zonoopt._core.IneqType
Get inequality type / direction.
- Returns:
inequality type (type member)
- Return type:
- get_n_dims(self: zonoopt._core.Inequality) int
Get number of dimensions for inequality.
- Returns:
number of dimensions (n_dims member)
- Return type:
int
- get_rhs(self: zonoopt._core.Inequality) float
Get right-hand side of the inequality.
- Returns:
right-hand side value (rhs member)
- Return type:
float
- get_terms(self: zonoopt._core.Inequality) list[zonoopt._core.IneqTerm]
Get reference to terms in the inequality.
- Returns:
reference to terms member
- Return type:
list[IneqTerm]
- set_ineq_type(self: zonoopt._core.Inequality, type: zonoopt._core.IneqType) None
Sets the direction of the inequality or sets it to be an equality
- Parameters:
type (IneqType) – inequality type (e.g., less than or equal, greater than or equal, or equal)
- set_rhs(self: zonoopt._core.Inequality, rhs: SupportsFloat) None
Sets right-hand side of the inequality.
- Parameters:
rhs (float) – right-hand side value
- class zonoopt.Interval(self: zonoopt._core.Interval, y_min: SupportsFloat, y_max: SupportsFloat)
Bases:
pybind11_objectInterval class
Implements interface from IntervalBase. This class owns its lower and upper bounds.
Interval constructor.
- Parameters:
y_min (float) – lower bound
y_max (float) – upper bound
- __add__(self: zonoopt._core.Interval, other: zonoopt._core.Interval) zonoopt._core.Interval
Interval addition
- __mul__(*args, **kwargs)
Overloaded function.
__mul__(self: zonoopt._core.Interval, other: zonoopt._core.Interval) -> zonoopt._core.Interval
Interval multiplication
- Args:
other (Interval): rhs interval
- Returns:
Interval: self * other
__mul__(self: zonoopt._core.Interval, alpha: typing.SupportsFloat) -> zonoopt._core.Interval
Interval multiplication with scalar
- Args:
alpha (float): scalar multiplier
- Returns:
Interval: alpha * self
- __sub__(self: zonoopt._core.Interval, other: zonoopt._core.Interval) zonoopt._core.Interval
Interval subtraction
- __truediv__(self: zonoopt._core.Interval, other: zonoopt._core.Interval) zonoopt._core.Interval
Interval division
- arccos(self: zonoopt._core.Interval) zonoopt._core.Interval
Compute interval containing arccos(x) for all x in interval
- Returns:
interval containing arccos(x)
- Return type:
- arcsin(self: zonoopt._core.Interval) zonoopt._core.Interval
Compute interval containing arcsin(x) for all x in interval
- Returns:
interval containing arcsin(x)
- Return type:
- arctan(self: zonoopt._core.Interval) zonoopt._core.Interval
Compute interval containing arctan(x) for all x in interval
- Returns:
interval containing arctan(x)
- Return type:
- center(self: zonoopt._core.Interval) float
Gets center of interval (ub + lb) / 2
- Returns:
center of interval
- Return type:
float
- contains(self: zonoopt._core.Interval, y: SupportsFloat) bool
Checks whether interval contains a value
- Parameters:
y (float) – scalar value
- Returns:
flag indicating if interval contains y
- Return type:
bool
- cos(self: zonoopt._core.Interval) zonoopt._core.Interval
Compute interval containing cos(x) for all x in interval
- Returns:
interval containing cos(x)
- Return type:
- exp(self: zonoopt._core.Interval) zonoopt._core.Interval
Compute interval containing exp(x) for all x in interval
- Returns:
interval containing exp(x)
- Return type:
- intersect(self: zonoopt._core.Interval, other: zonoopt._core.Interval) zonoopt._core.Interval
Interval intersection
- inv(self: zonoopt._core.Interval) zonoopt._core.Interval
Interval inverse
- Returns:
inverse of self
- Return type:
- is_empty(self: zonoopt._core.Interval) bool
Checks whether interval is empty
- Returns:
flag indicating whether interval is empty
- Return type:
bool
- is_single_valued(self: zonoopt._core.Interval) bool
Checks whether interval is single-valued (i.e., width is 0 within numerical tolerance)
- Returns:
flag indicating if interval is single-value
- Return type:
bool
- property lb
lower bound
- sin(self: zonoopt._core.Interval) zonoopt._core.Interval
Compute interval containing sin(x) for all x in interval
- Returns:
interval containing sin(x)
- Return type:
- tan(self: zonoopt._core.Interval) zonoopt._core.Interval
Compute interval containing tan(x) for all x in interval
- Returns:
interval containing tan(x)
- Return type:
- property ub
upper bound
- width(self: zonoopt._core.Interval) float
Gets width of interval (ub - lb)
- Returns:
width of interval
- Return type:
float
- class zonoopt.OptSettings(self: zonoopt._core.OptSettings)
Bases:
pybind11_objectSettings for optimization routines in ZonoOpt library.
- property contractor_iter
number of interval contractor iterations
- property contractor_tree_search_depth
when applying interval contractor in branch and bound, this is how deep to search the constraint tree for affected variables
- property eps_a
absolute convergence tolerance
- property eps_dual
dual convergence tolerance
- property eps_dual_search
dual residual convergence tolerance during branch and bound and search
- property eps_prim
primal convergence tolerance
- property eps_prim_search
primal residual convergence tolerance during branch and bound and search
- property eps_r
relative convergence tolerance
- property inf_norm_conv
use infinity norm for convergence check (if false, scaled 2-norm is used)
- property k_inf_check
check infeasibility every k_inf_check iterations
- property k_max_admm
max admm iterations
- property k_max_bnb
max number of branch-and-bound iterations
- property max_nodes
terminate if more than this many nodes are in branch and bound queue
- property n_threads_bnb
max threads for branch and bound
- property polish
flag to perform solution polishing
- property rho
admm penalty parameter, higher prioritizes feasibility during iterations, lower prioritizes optimality
- property search_mode
0-> best first, 1-> best dive
- settings_valid(self: zonoopt._core.OptSettings) bool
check whether settings struct is valid
- property t_max
max time for optimization
- property use_interval_contractor
flag to use interval contractor for constraint tightening / implication
- property verbose
display optimization progress
- property verbosity_interval
print every verbose_interval iterations
- class zonoopt.OptSolution(self: zonoopt._core.OptSolution)
Bases:
pybind11_objectSolution data structure for optimization routines in ZonoOpt library.
- property J
objective
- property converged
true if optimization has converged
- property dual_residual
dual residual, corresponds to optimality
- property infeasible
true if optimization problem is provably infeasible
- property iter
number of iterations
- property primal_residual
primal residual, corresponds to feasibility
- property run_time
time to compute solution
- property startup_time
time to factorize matrices and run interval contractors
- property u
ADMM dual variable
- property x
ADMM primal variable, approximately equal to z when converged
- property z
solution vector
- class zonoopt.Point(self: zonoopt._core.Point, c: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]'])
Bases:
ZonoPoint class
A point is defined entirely by the center vector c.
Point constructor
- Parameters:
c (numpy.array) – center vector
- set(self: zonoopt._core.Point, c: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]']) None
Reset point object with the given parameters.
- Parameters:
c (numpy.array) – center vector
- class zonoopt.Zono(self: zonoopt._core.Zono, G: scipy.sparse.csc_matrix[numpy.float64], c: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]'], zero_one_form: bool = False)
Bases:
ConZonoZonotope class
A zonotope is defined as: Z = {G xi + c | xi in [-1, 1]^nG}. Equivalently, the following shorthand can be used: Z = <G, c>. Optionally, in 0-1 form, the factors are xi in [0,1]. The set dimension is n, and the number of generators is nG.
Zono constructor
- Parameters:
G (scipy.sparse.csc_matrix) – generator matrix
c (numpy.array) – center
zero_one_form (bool, optional) – true if set is in 0-1 form
- get_volume(self: zonoopt._core.Zono) float
Get volume of zonotope.
Reference: Gover and Krikorian 2010, “Determinants and the volumes of parallelotopes and zonotopes” Requires nG choose n determinant computations.
- Returns:
volume of zonotope
- Return type:
float
- reduce_order(self: zonoopt._core.Zono, n_o: SupportsInt) zonoopt._core.Zono
Perform zonotope order reduction.
- Parameters:
n_o (int) – desired order, must be greater than or equal to the dimension of the set
- Returns:
zonotope with order n_o
- Return type:
- set(self: zonoopt._core.Zono, G: scipy.sparse.csc_matrix[numpy.float64], c: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]'], zero_one_form: bool = False) None
Reset zonotope object with the given parameters.
- Parameters:
G (scipy.sparse.csc_matrix) – generator matrix
c (numpy.array) – center
zero_one_form (bool, optional) – true if set is in 0-1 form
- zonoopt.affine_map(Z: zonoopt._core.HybZono, R: scipy.sparse.csc_matrix[numpy.float64], s: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, 1]'] = array([], dtype=float64)) zonoopt._core.HybZono
Returns affine map R*Z + s of set Z
- zonoopt.cartesian_product(Z1: zonoopt._core.HybZono, Z2: zonoopt._core.HybZono) zonoopt._core.HybZono
Computes the Cartesian product of two sets Z1 and Z2.
- zonoopt.constrain(Z: zonoopt._core.HybZono, ineqs: collections.abc.Sequence[zonoopt._core.Inequality], R: scipy.sparse.csc_matrix[numpy.float64] = <Compressed Sparse Column sparse matrix of dtype 'float64' with 0 stored elements and shape (0, 0)>) zonoopt._core.HybZono
Applies inequalities to set.
- Parameters:
Z (HybZono) – Set for inequalities to be applied to
ineqs (list[Inequality]) – list of inequalities
R (scipy.sparse.csc_matrix, optional) – For generalized intersection with the inequalities. Defaults to identity.
- Returns:
zonotopic set
- Return type:
Constrains set Z by applying the given inequalities to the set. For example, given a set Z with states z0, z1, z2, the constraint z0 + z1 - z2 <= 2 could be added via an inequality object. R is used for generalized intersection-like operations. For instance, when all the inequalities are <= inequalities, this function returns Z int_R (Hx<=f) where H is the halfspace represented by the inequalities.
- zonoopt.get_vertices(Z, t_max=60.0)
Get vertices of zonotopic set using scipy linprog.
- Parameters:
Z (HybZono) – Zonotopic set.
t_max (float, optional) – Maximum time to spend on finding vertices. Defaults to 60.0 seconds.
- Returns:
Vertices of the zonotopic set. If Z is a point, returns its coordinates.
- Return type:
numpy.ndarray
- zonoopt.halfspace_intersection(Z: zonoopt._core.HybZono, H: scipy.sparse.csc_matrix[numpy.float64], f: typing.Annotated[numpy.typing.ArrayLike, numpy.float64, "[m, 1]"], R: scipy.sparse.csc_matrix[numpy.float64] = <Compressed Sparse Column sparse matrix of dtype 'float64' with 0 stored elements and shape (0, 0)>) zonoopt._core.HybZono
Computes the intersection generalized intersection of set Z with halfspace H*x <= f over matrix R.
- zonoopt.intersection(Z1: zonoopt._core.HybZono, Z2: zonoopt._core.HybZono, R: scipy.sparse.csc_matrix[numpy.float64] = <Compressed Sparse Column sparse matrix of dtype 'float64' with 0 stored elements and shape (0, 0)>) zonoopt._core.HybZono
Computes the generalized intersection of sets Z1 and Z2 over the matrix R.
- zonoopt.intersection_over_dims(Z1: zonoopt._core.HybZono, Z2: zonoopt._core.HybZono, dims: collections.abc.Sequence[SupportsInt]) zonoopt._core.HybZono
Computes the intersection of sets Z1 and Z2 over the specified dimensions.
- zonoopt.interval_2_zono(box: zonoopt._core.Box) zonoopt._core.Zono
Builds a zonotope from a Box object.
- zonoopt.make_regular_zono_2D(radius: SupportsFloat, n_sides: SupportsInt, outer_approx: bool = False, c: Annotated[numpy.typing.ArrayLike, numpy.float64, '[2, 1]'] = array([0., 0.])) zonoopt._core.Zono
Builds a 2D regular zonotope with a given radius and number of sides.
- Parameters:
radius (float) – radius of the zonotope
n_sides (int) – number of sides (must be an even number >= 4)
outer_approx (bool, optional) – flag to do an outer approximation instead of an inner approximation
c (numpy.array, optional) – center vector
- Returns:
zonotope
- Return type:
- zonoopt.minkowski_sum(Z1: zonoopt._core.HybZono, Z2: zonoopt._core.HybZono) zonoopt._core.HybZono
Computes Minkowski sum of two sets Z1 and Z2.
- zonoopt.plot(Z, ax=None, settings=OptSettings structure: verbose: false verbosity_interval: 100 t_max: 1.79769e+308 k_max_admm: 5000 rho: 10 eps_dual: 0.01 eps_prim: 0.001 k_inf_check: 10 inf_norm_conv: true use_interval_contractor: true contractor_iter: 1 search_mode: 0 polish: 1 eps_dual_search: 0.1 eps_prim_search: 0.01 eps_r: 0.01 eps_a: 0.1 k_max_bnb: 100000 n_threads_bnb: 4 max_nodes: 100000 contractor_tree_search_depth: 10, t_max=60.0, **kwargs)
Plots zonotopic set using matplotlib.
- Parameters:
Z (HybZono) – zonotopic set to be plotted
ax (matplotlib.axes.Axes, optional) – Axes to plot on. If None, current axes are used.
settings (OptSettings, optional) – Settings for the optimization. Defaults to OptSettings().
t_max (float, optional) – Maximum time to spend on finding vertices. Defaults to 60.0 seconds.
**kwargs – Additional keyword arguments passed to the plotting function (e.g., color, alpha).
- Returns:
List of matplotlib objects representing the plotted zonotope.
- Return type:
list
- zonoopt.pontry_diff(Z1: zonoopt._core.HybZono, Z2: zonoopt._core.HybZono, exact: bool = False) zonoopt._core.HybZono
Computes the Pontryagin difference Z1 - Z2.
For inner approximations (exact=false), the algorithm from Vinod et. al. 2025 is used. Note that this algorithm is exact when the minuend is a constrained zonotope and the matrix [G;A] is invertible. Exact Pontryagin difference can only be computed when the subtrahend is a zonotope. If subtrahend is a constrained zonotope, it will first be over-approximated as a zonotope. If subtrahend is a hybrid zonotope, a get_leaves operation will first be performed to produce a union of constrained zonotopes. If the minuend is a hybrid zonotope and exact is false, a get_leaves operation will be performed for Z1 to reduce the number of leaves in the resulting set.
- zonoopt.project_onto_dims(Z: zonoopt._core.HybZono, dims: collections.abc.Sequence[SupportsInt]) zonoopt._core.HybZono
Projects set Z onto the dimensions specified in dims.
- zonoopt.set_diff(Z1: zonoopt._core.HybZono, Z2: zonoopt._core.HybZono, delta_m: typing.SupportsFloat = 100, remove_redundancy: bool = True, settings: zonoopt._core.OptSettings = OptSettings structure: verbose: false verbosity_interval: 100 t_max: 1.79769e+308 k_max_admm: 5000 rho: 10 eps_dual: 0.01 eps_prim: 0.001 k_inf_check: 10 inf_norm_conv: true use_interval_contractor: true contractor_iter: 1 search_mode: 0 polish: 1 eps_dual_search: 0.1 eps_prim_search: 0.01 eps_r: 0.01 eps_a: 0.1 k_max_bnb: 100000 n_threads_bnb: 4 max_nodes: 100000 contractor_tree_search_depth: 10, solution: zonoopt._core.OptSolution = None, n_leaves: typing.SupportsInt = 2147483647, contractor_iter: typing.SupportsInt = 100) zonoopt._core.HybZono
Set difference Z1 \ Z2
- Parameters:
Z1 (HybZono) – zonotopic set
Z2 (HybZono) – zonotopic set
delta_m (float, optional) – parameter defining range of complement
remove_redundancy (bool, optional) – remove redundant constraints and unused generators in get_leaves function call
settings (OptSettings, optional) – optimization settings for get_leaves function call
solution (OptSolution, optional) – optimization solution for get_leaves function call
n_leaves (int, optional) – maximum number of leaves to return in get_leaves function call
contractor_iter (int, optional) – number of interval contractor iterations if using remove_redundancy
- Returns:
zonotopic set
- Return type:
- zonoopt.union_of_many(Z_list: collections.abc.Sequence[zonoopt._core.HybZono], preserve_sharpness: bool = False, expose_indicators: bool = False) zonoopt._core.HybZono
Computes the union of several sets
- Parameters:
Z_list (list[HybZono]) – sets to be unioned
preserve_sharpness (bool, optional) – flag to preserve sharpness of the union at expense of complexity.
expose_indicators (bool, optional) – flag to append indicator set to the union.
- Returns:
zonotopic set
- Return type:
Computes union of sets {Z0, Z1, …, Zn}. If expose_indicators is true, returns union({Z0, …, Zn}) x I where I is the indicator set for the union. Specifically, each dimension of I corresponds to one of the Zi in the union. So for union_of_many({Z0, Z1, Z2}, true) with Z0, Z1, Z2 not intersecting, if a vector [z, i] is in union({Z0, Z1, Z2}) x I, then i = [1, 0, 0] if z is in Z0, etc.
- zonoopt.vrep_2_conzono(Vpoly: Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, n]']) zonoopt._core.ConZono
Builds a constrained zonotope from a vertex representation polytope.
- Parameters:
Vpoly (numpy.array) – vertices of V-rep polytope
- Returns:
constrained zonotope
- Return type:
Vpoly is a matrix where each row is a vertex of the polytope.
- zonoopt.vrep_2_hybzono(Vpolys: collections.abc.Sequence[Annotated[numpy.typing.ArrayLike, numpy.float64, '[m, n]']], expose_indicators: bool = False) zonoopt._core.HybZono
Computes a hybrid zonotope from a union of vertex representation polytopes.
- Parameters:
Vpolys (list[numpy.array]) – V-rep polytopes to be unioned
expose_indicators (bool, optional) – flag to append indicator set to the union.
- Returns:
hybrid zonotope
- Return type:
Vpolys is a vector of matrices, where each matrix represents a polytope in vertex representation. Each row in each polytope matrix is a vertex of the polytope, and each column corresponds to a dimension. The function constructs a hybrid zonotope in [0,1] form that represents the union of these polytopes. This function computes union of sets {V0, V1, …, Vn}. If expose_indicators is true, returns union({V0, …, Vn}) x I where I is the indicator set for the union. Specifically, each dimension of I corresponds to one of the Vi in the union. So for vrep_2_hybzono({V0, V1, V2}, true) with V0, V1, V2 not intersecting, if a vector [z, i] is in union({V0, V1, V2}) x I, then i = [1, 0, 0] if z is in V0, etc.
- zonoopt.zono_union_2_hybzono(Zs: collections.abc.Sequence[zonoopt._core.Zono], expose_indicators: bool = False) zonoopt._core.HybZono
Computes a hybrid zonotope from a union of zonotopes.
- Parameters:
Zs (list[Zono]) – zonotopes to be unioned
expose_indicators (bool, optional) – flag to append indicator set to the union.
- Returns:
hybrid zonotope
- Return type:
This function computes union of sets {Z0, Z1, …, Zn}. This can be more efficient than union_of_many if all sets are zonotopes because generators can be reused. If expose_indicators is true, returns union({Z0, …, Zn}) x I where I is the indicator set for the union. Specifically, each dimension of I corresponds to one of the Zi in the union. So for zono_union_2_hybzono({Z0, Z1, Z2}, true) with Z0, Z1, VZ2 not intersecting, if a vector [z, i] is in union({Z0, Z1, Z2}) x I, then i = [1, 0, 0] if z is in Z0, etc.