pycvxset.ConstrainedZonotope (API details)

pycvxset.ConstrainedZonotope(**kwargs)

Constrained zonotope class

class pycvxset.ConstrainedZonotope.ConstrainedZonotope(**kwargs)[source]

Bases: object

Constrained zonotope class

Constrained zonotope defines a polytope in the working dimension \(\mathbb{R}^n\) as an affine transformation of a polytope defined in latent space \(B_\infty(A_e, b_e)\subset \mathbb{R}^{N_C}\). Here, \(B_\infty(A_e, b_e)\) is defined as the intersection of a unit \(\ell_\infty\)-norm ball and a collection of \(M_C\) linear constraints \(\{\xi\in\mathbb{R}^{N_C}|A_e \xi = b_e\}.\)

Formally, a constrained zonotope is defined as follows,

\[\mathcal{C} = \{G \xi + c\ |\ \xi \in B_\infty(A_e, b_e)\} \subset \mathbb{R}^n,\]

where

\[B_\infty(A_e, b_e)= \{\xi\ |\ \| \xi \|_\infty \leq 1, A_e \xi = b_e\} \subset \mathbb{R}^{N_C},\]

with \(G\in\mathbb{R}^{n\times N_C}\), \(c\in\mathbb{R}^{n}\), \(A_e\in\mathbb{R}^{M_C\times N_C}\), and \(b\in\mathbb{R}^{M_C}\).

A constrained zonotope provide an alternative and equivalent representation of any convex and compact polytope. Furthermore, a constrained zonotope admits closed-form expressions for several set manipulations that can often be accomplished without invoking any optimization solvers. See [SDGR16] [RK22] [VWD24] for more details.

A zonotope is a special class of constrained zonotopes, and are defined as

\[\mathcal{Z} = \{G \xi + c\ |\ \|\xi\|_\infty \leq 1\} \subset \mathbb{R}^n.\]

In other words, a zonotope is a constrained zonotope with no equality constraints in the latent dimension space. In ConstrainedZonotope, we model zonotopes by having (Ae,be) be empty (n_equalities is zero).

Constrained zonotope object construction admits one of the following combinations (as keyword arguments):

  1. dim for an empty constrained zonotope of dimension dim,

  2. (G, c, Ae, be) for a constrained zonotope,

  3. (G, c) for a zonotope,

  4. (lb, ub) for a zonotope equivalent to an axis-aligned cuboid with appropriate bounds \(\{x\ |\ lb\leq x \leq ub\}\), and

  5. (c, h) for a zonotope equivalent to an axis-aligned cuboid centered at c with specified scalar/vector half-sides \(h\), \(\{x\ |\ \forall i\in\{1,2,...,n\}, |x_i - c_i| \leq h_i\}\).

  6. (c=p, G=None) for a zonotope equivalent to a single point p,

  7. P for a constrained zonotope equivalent to the pycvxset.Polytope.Polytope object P,

Parameters:
  • dim (int, optional) – Dimension of the empty constrained zonotope. If NOTHING is provided, dim=0 is assumed.

  • c (array_like, optional) – Affine transformation translation vector. Must be 1D array, and the constrained zonotope dimension is determined by number of elements in c. When c is provided, either (G) or (G, Ae, be) or (h) must be provided additionally. When h is provided, c is the centroid of the resulting zonotope.

  • G (array_like) – Affine transformation matrix. The vectors are stacked vertically with matching number of rows as c. When G is provided, (c, Ae, be) OR (c) must also be provided. To define a constrained zonotope with a single point, set c to the point AND G to None (do not set (Ae, be) or set them to (None, None)).

  • Ae (array_like) – Equality coefficient vectors. The vectors are stacked vertically with matching number of columns as G. When Ae is provided, (G, c, be) must also be provided.

  • be (array_like) – Equality coefficient constants. The constants are expected to be in a 1D numpy array. When be is provided, (G, c, Ae) must also be provided.

  • lb (array_like, optional) – Lower bounds of the axis-aligned cuboid. Must be 1D array, and the constrained zonotope dimension is determined by number of elements in lb. When lb is provided, ub must also be provided.

  • ub (array_like, optional) – Upper bounds of the axis-aligned cuboid. Must be 1D array of length as same as lb. When ub is provided, lb must also be provided.

  • h (array_like, optional) – Half-side length of the axis-aligned cuboid. Can be a scalar or a vector of length as same as c. When h is provided, c must also be provided.

  • polytope (Polytope, optional) – Polytope to use to construct constrained zonotope.

Raises:
  • ValueError – (G, c) is not compatible.

  • ValueError – (G, c, Ae, be) is not compatible.

  • ValueError – (lb, ub) is not valid

  • ValueError – (c, h) is not valid

  • ValueError – Provided polytope is not bounded.

  • UserWarning – When a row with all zeros in Ae and be.

__init__(**kwargs)[source]

Constructor for ConstrainedZonotope class

affine_map(M)

Multiply a matrix or a scalar with a constrained zonotope

Parameters:

M (array_like) – Matrix (N times self.dim) or a scalar to be multiplied with a constrained zonotope

Returns:

Constrained zonotope which is the product of M and self. Specifically, given a constrained zonotope \(\mathcal{P}\), and a matrix \(M\in\mathbb{R}^{m\times P.\text{dim}}\) or scalar \(M\), then this function returns a constrained zonotope \(\mathcal{R}=\{Mx|x\in\mathcal{P}\}\).

Return type:

ConstrainedZonotope

Notes

This function implements (11) of [SDGR16].

approximate_pontryagin_difference(norm_type, G_S, c_S, method='inner-least-squares')

Approximate Pontryagin difference between a constrained zonotope and an affine transformation of a unit-norm ball.

Specifically, we approximate the Pontryagin difference \(\mathcal{P}\ominus\mathcal{S}=\{x\in\mathbb{R}^n\ |\ x + \mathcal{S}\subseteq\mathcal{P}\}\) between a constrained zonotope \(\mathcal{P}\) and a subtrahend \(\mathcal{S}\). The subtrahend must be specifically of the form \(\mathcal{S}=\{G_S \xi + c_S | \|\xi\|_p\leq 1\}\) for norm_type \(p\in\{1,2,\infty\}\).

  • For p=1, the subtrahend is a convex hull of intervals specified by the columns of \(G_S\) with each interval is symmetric about \(c_S\).

  • For p=2, the subtrahend is an ellipsoid is characterized by affine transformation \(G_S\) of a unit Euclidean norm ball which is then shifted to \(c_S\). Here, \(G_S\) is given by pycvxset.Ellipsoid.Ellipsoid.G. See pycvxset.Ellipsoid.Ellipsoid() to obtain \(G_S\) and \(c_S\) for broader classes of ellipsoids.

  • For p=’inf, the subtrahend is a zonotope characterized by (\(G_S\), \(c_S\)).

Parameters:
  • norm_type (int | str) – Norm type.

  • G_S (array_like) – Affine transformation matrix (self.dim times N) for scaling the ball

  • c_S (array_like) – Affine transformation vector (self.dim,) for translating the ball

  • method (str, optional) – Approximation method to use. Can be one of [‘inner-least-squares’]. Defaults to ‘inner-least-squares’.

Raises:
  • ValueError – Norm type is not in [1, 2, ‘inf]

  • ValueError – Mismatch in dimension (no. of columns of G_S, length of c_S, and self.dim should match)

Returns:

Approximated Pontryagin difference between self and the affine transformation of the unit-norm ball.

Return type:

ConstrainedZonotope

Notes

For method ‘inner-least-squares’, this function implements Table 1 of [VWD24].

cartesian_product(sequence_Q)

Generate the Cartesian product of a set \(\mathcal{Q}\) (or a list of mathcal{Q}) with \(\mathcal{P}\).

Parameters:

sequence_Q (list | tuple | Polytope | ConstrainedZonotope) – List of sets to take Cartesian product with

Returns:

Cartesian product of self and all sets in sequence_Q

Return type:

ConstrainedZonotope

chebyshev_centering()

Computes a ball with the largest radius that fits within the constrained zonotope. The ball’s center is known as the Chebyshev center, and its radius is the Chebyshev radius.

Raises:
  • ValueError – When the constrained zonotope is not full-dimensional

  • ValueError – When the constrained zonotope is empty

  • NotImplementedError – Unable to solve the linear program using CVXPY

Returns:

A tuple with two items
  1. center (numpy.ndarray): Approximate Chebyshev radius of the constrained zonotope

  2. radius (float): Approximate Chebyshev radius of the constrained zonotope

Return type:

tuple

Notes

Unlike pycvxset.Polytope.Polytope.chebyshev_centering(), this function computes an approximate Chebyshev center and radius. Specifically, it guarantees that a ball of the computed radius, centered at the computed center is contained in the constrained zonotope. However, it does not guarantee that the computed radius is the radius of the largest possible ball contained in the given constrained zonotope. For more details, see [VWS24].

closest_point(points, p=2)

Wrapper for project() to compute the point in the convex set closest to the given point.

Parameters:
  • points (array_like) – Points to project. Matrix (N times self.dim), where each row is a point.

  • p (str | int) – Norm-type. It can be 1, 2, or ‘inf’. Defaults to 2.

Returns:

Projection of points to the set as a 2D numpy.ndarray. These arrays have as many rows as points.

Return type:

numpy.ndarray

Notes

For more detailed description, see documentation for project() function.

containment_constraints(x, flatten_order='F')[source]

Get CVXPY constraints for containment of x (a cvxpy.Variable) in a constrained zonotope.

Parameters:
  • x (cvxpy.Variable) – CVXPY variable to be optimized

  • flatten_order (char) – Order to use for flatten (choose between “F”, “C”). Defaults to “F”, which implements column-major flatten. In 2D, column-major flatten results in stacking rows horizontally to achieve a single horizontal row.

Raises:

ValueError – When constrained zonotope is empty

Returns:

A tuple with two items:

  1. constraint_list (list): CVXPY constraints for the containment of x in the constrained zonotope.

  2. xi (cvxpy.Variable | None): CVXPY variable representing the latent dimension variable. It is None, when the constrained zonotope is a single point.

Return type:

tuple

contains(Q, verbose=False, time_limit=60)

Check containment of a set \(\mathcal{Q}\) (could be a polytope or a constrained zonotope), or a collection of points \(Q \in \mathbb{R}^{n_Q \times \mathcal{P}.\text{dim}}\) in the given constrained zonotope.

Parameters:
  • Q (array_like | Polytope | ConstrainedZonotope) – Polytope object, ConstrainedZonotope object, or a collection of points to be tested for containment within the constrained zonotope. When providing a collection of points, Q is a matrix (N times self.dim) with each row is a point.

  • verbose (bool, optional) – Verbosity flag to provide cvxpy when solving the MIQP related to checking containment of a constrained zonotope within another constrained zonotope. Defaults to False.

  • time_limit (int, optional) – Time limit in seconds for GUROBI solver when solving the MIQP related to checking containment of a constrained zonotope within another constrained zonotope. Set time_limit to np.inf if no limit is desired. Defaults to 60s (see constants.py).

Raises:
  • ValueError – Dimension mismatch between Q and the constrained zonotope

  • ValueError – Q is Ellipsoid

  • ValueError – Unable to perform containment check between two constrained zonotopes (including time_limit issues)

  • NotImplementedError – GUROBI is not installed when checking containment of two constrained zonotopes

  • NotImplementedError – Failed to check containment between two constrained zonotopes, due to an unhandled status

Returns:

Boolean corresponding to \(\mathcal{Q}\subseteq\mathcal{P}\) or \(\mathcal{Q}\in\mathcal{P}\) .

Return type:

bool | numpy.ndarray([bool])

Notes

  • We use \(\mathcal{P}\) to denote the constrained zonotope characterized by self.

  • When Q is a constrained zonotope, a bool is returned which is True if and only if \(\mathcal{Q}\subseteq \mathcal{P}\). This function uses the non-convex programming capabilities of GUROBI (via CVXPY) to check containment between two constrained zonotopes.

  • When Q is a polytope, a bool is returned which is True if and only if \(\mathcal{Q}\subseteq \mathcal{P}\).

    • When Q is in V-Rep, the containment problem simplifies to checking if all vertices of Q are contained \(\mathcal{P}\).

    • When Q is in H-Rep, this function converts Q into a constrained zonotope, and then checks for containment.

  • When Q is a single n-dimensional point, a bool is returned which is True if and only if \(Q\in \mathcal{P}\). This function uses distance() to check containment of a point in a constrained zonotope.

  • When Q is a collection of n-dimensional points (Q is a matrix with each row describing a point), an array of is returned where each element is True if and only if \(Q_i\in \mathcal{P}\). This function uses distance() to check containment of a point in a constrained zonotope.

Warning

This function requires gurobipy when checking if the ConstrainedZonotope object represented by self contains a ConstrainedZonotope object Q or pycvxset.Polytope.Polytope object Q in H-Rep.

copy()[source]

Get a copy of the constrained zonotope

distance(points, p=2)

Wrapper for project() to compute distance of a point to a convex set.

Parameters:
  • points (array_like) – Points to project. Matrix (N times self.dim), where each row is a point.

  • p (int | str) – Norm-type. It can be 1, 2, or ‘inf’. Defaults to 2.

Returns:

Distance of points to the set as a 1D numpy.ndarray. These arrays have as many rows as points.

Return type:

numpy.ndarray

Notes

For more detailed description, see documentation for project() function.

extreme(eta)

Wrapper for support() to compute the extreme point.

Parameters:

eta (array_like) – Support directions. Matrix (N times self.dim), where each row is a support direction.

Returns:

Support vector evaluation(s) as a 2D numpy.ndarray. The array has as many rows as eta.

Return type:

numpy.ndarray

Notes

For more detailed description, see documentation for support() function.

interior_point()

Compute an interior point of the constrained zonotope.

Returns:

A point that lies in the (relative) interior of the constrained zonotope.

Return type:

numpy.ndarray

Notes

This function returns \(G v + c\), where \(v\) is the Chebyshev center of the polytope \(B_\infty(A_e, b_e)= \{\xi\ |\ \| \xi \|_\infty \leq 1, A_e \xi = b_e\} \subset \mathbb{R}^{N_C}\). See pycvxset.Polytope.Polytope.chebyshev_centering() for more details on Chebyshev centering.

intersection(Q)

Compute the intersection of constrained zonotope with another constrained zonotope or polytope.

Parameters:

Q (ConstrainedZonotope | Polytope) – Set to intersect with

Raises:

TypeError – When Q is neither a ConstrainedZonotope object or a Polytope object

Returns:

Intersection of self with Q

Return type:

ConstrainedZonotope

Notes

intersection_under_inverse_affine_map(Y, R)

Compute the intersection of constrained zonotope with another constrained zonotope under an inverse affine map

Parameters:
  • Y (ConstrainedZonotope) – Set to intersect with

  • R (array_like) – Matrix (Y.dim times self.dim) for the inverse-affine map.

Raises:
  • ValueError – When Y is not a ConstrainedZonotope

  • ValueError – When R is not of correct dimension

Returns:

Intersection of self with Y under an inverse affine map R. Specifically, given constrained zonotopes \(\mathcal{P}\) (self) and \(\mathcal{Y}\), and a matrix \(R\), we compute the set \(\mathcal{P}\cap_R\mathcal{Y}=\{x \in \mathcal{P}| Rx \in \mathcal{Y}\}\). When \(R\) is invertible, \(\mathcal{P}\cap_R\mathcal{Y}=(R^{-1}\mathcal{P})\cap\mathcal{Y}\).

Return type:

ConstrainedZonotope

Notes

This function implements (13) of [SDGR16]. This function does not require R to be invertible.

intersection_with_affine_set(Ae, be)

Compute the intersection of a constrained zonotope with an affine set.

Parameters:
  • Ae (array_like) – Equality coefficient matrix (N times self.dim) that define the affine set \(\{x|A_ex = b_e\}\).

  • be (array_like) – Equality constant vector (N,) that define the affine set \(\{x| A_ex = b_e\}\).

Raises:

ValueError – When the number of columns in Ae is different from self.dim

Returns:

The intersection of a constrained zonotope with the affine set.

Return type:

ConstrainedZonotope

Notes

This function implements imposes the constraints \(\{A_ex = b_e\}\) as constraints in the latent dimension of the constrained zonotope — \(A_e (G \xi + c) = b_e\) for every feasible \(\xi\).

intersection_with_halfspaces(A, b)

Compute the intersection of a constrained zonotope with a collection of halfspaces (polyhedron).

Parameters:
  • A (array_like) – Inequality coefficient matrix (N times self.dim) that define the polyhedron \(\{x|Ax \leq b\}\).

  • b (array_like) – Inequality constant vector (N,) that define the polyhedron \(\{x| Ax \leq b\}\).

Returns:

The intersection of a constrained zonotope with a collection of halfspaces.

Return type:

ConstrainedZonotope

Notes

This function implements (10) of [RK22] for each halfspace. We skip redundant inequalities when encountered and we return empty set when intersection yields an empty set.

inverse_affine_map_under_invertible_matrix(M)

Compute the inverse affine map of a constrained zonotope for a given matrix M.

Parameters:

M (array_like) – Square invertible matrix of dimension self.dim

Raises:
  • ValueError – M is not convertible into a 2D float array

  • ValueError – M is not square

  • ValueError – M is not invertible

Returns:

Inverse affine map of self under M. Specifically, Given a constrained zonotope \(\mathcal{P}\) and an invertible self.dim-dimensional matrix \(M\in\mathbb{R}^{\mathcal{P}.dim\times \mathcal{P}.\text{dim}}\), this function computes the constrained zonotope \(\mathcal{R}=M^{-1}\mathcal{P}\).

Return type:

ConstrainedZonotope

Notes

This function is a wrapper for affine_map(). We require M to be invertible in order to ensure that the resulting set is representable as a constrained zonotope.

maximum_volume_inscribing_ellipsoid()

Compute the maximum volume ellipsoid that fits within the given constrained zonotope.

Raises:
  • ValueError – When the constrained zonotope is not full-dimensional

  • ValueError – When the constrained zonotope is empty

  • NotImplementedError – Unable to solve the convex program using CVXPY

Returns:

A tuple with three items:
  1. center (numpy.ndarray): Approximate maximum volume inscribed ellipsoid’s center

  2. shape_matrix (numpy.ndarray): Approximate maximum volume inscribed ellipsoid’s shape matrix

  3. sqrt_shape_matrix (numpy.ndarray): Approximate maximum volume inscribed ellipsoid’s square root of shape matrix.

Return type:

tuple

Notes

Unlike pycvxset.Polytope.Polytope.maximum_volume_inscribing_ellipsoid(), this function computes an approximate maximum volume inscribed ellipsoid. Specifically, it guarantees that the computed ellipsoid is contained in the constrained zonotope. However, it does not guarantee that the computed ellipsoid is the largest possible ellipsoid (in terms of volume) contained in the given constrained zonotope. For more details, see [VWS24].

minimize(x, objective_to_minimize, cvxpy_args, task_str='')

Solve a convex program with CVXPY objective subject to containment constraints.

Parameters:
  • x (cvxpy.Variable) – CVXPY variable to be optimized

  • objective_to_minimize (cvxpy.Expression) – CVXPY expression to be minimized

  • cvxpy_args (dict) – CVXPY arguments to be passed to the solver

  • task_str (str, optional) – Task string to be used in error messages. Defaults to ‘’.

Raises:

NotImplementedError – Unable to solve problem using CVXPY

Returns:

A tuple with three items:
  1. x.value (numpy.ndarray): Optimal value of x. np.nan * np.ones((self.dim,)) if the problem is not solved.

  2. problem.value (float): Optimal value of the convex program. np.inf if the problem is infeasible, -np.inf if problem is unbounded, and finite otherwise.

  3. problem_status (str): Status of the problem

Return type:

tuple

Notes

This function uses containment_constraints() to obtain the list of CVXPY expressions that form the containment constraints on x.

Warning

Please pay attention to the NotImplementedError generated by this function. It may be possible to get CVXPY to solve the same problem by switching the solver. For example, consider the following code block.

from pycvxset import Polytope
P = Polytope(A=[[1, 1], [-1, -1]], b=[1, 1])
P.cvxpy_args_lp = {'solver': 'CLARABEL'}    # Default solver used in pycvxset
try:
    print('Is polytope bounded?', P.is_bounded)
except NotImplementedError as err:
    print(str(err))
P.cvxpy_args_lp = {'solver': 'OSQP'}
print('Is polytope bounded?', P.is_bounded)

This code block produces the following output:

Unable to solve the task (support function evaluation of the set at eta = [-0. -1.]). CVXPY returned error:
Solver 'CLARABEL' failed. Try another solver, or solve with verbose=True for more information.

Is polytope bounded? False
minimum_volume_circumscribing_rectangle()

Compute the minimum volume circumscribing rectangle for a set.

Raises:
  • ValueError – When set is empty

  • ValueError – When set is unbounded OR solver error!

Returns:

A tuple of two elements
  • lb (numpy.ndarray): Lower bound \(l\) on the set, \(\mathcal{P}\subseteq\{l\}\oplus\mathbb{R}_{\geq 0}\).

  • ub (numpy.ndarray): Upper bound \(u\) on the set, \(\mathcal{P}\subseteq\{u\}\oplus(-\mathbb{R}_{\geq 0})\).

Return type:

tuple

Notes

This function computes the lower/upper bound by an element-wise support computation (2n linear programs), where n is attr:self.dim. To reuse the support() function for the lower bound computation, we solve the optimization for each \(i\in\{1,2,...,n\}\),

\[\inf_{x\in\mathcal{P}} e_i^\top x=-\sup_{x\in\mathcal{P}} -e_i^\top x=-\rho_{\mathcal{P}}(-e_i),\]

where \(e_i\in\mathbb{R}^n\) denotes the standard coordinate vector, and \(\rho_{\mathcal{P}}\) is the support function of \(\mathcal{P}\).

minus(Q)

Implement - operation: Pontryagin difference with a constrained zonotope minuend

When Q is an ellipsoid or a zonotope, minus returns an inner-approximation of the set corresponding to the Pontryagin difference of a constrained zonotope and Q. When Q is a point or a singleton set, an exact set corresponding to the translation by -Q or -Q.c is returned.

Parameters:

Q (array_like | Ellipsoid | ConstrainedZonotope) – Point/set to use as subtrahend in the Pontryagin difference.

Raises:
  • TypeError – When Q is not one of the following — convertible into a 1D numpy array, or an ellipsoid, or a zonotope.

  • ValueError – When Q has a dimension mismatch with self.

  • UserWarning – When using with Q that is a set, warn that an inner-approximation is returned.

Returns:

Pontryagin difference of self and Q

Return type:

ConstrainedZonotope

Notes

This function uses approximate_pontryagin_difference() when Q is a set and uses plus() when Q is a point.

plot(method='inner', ax=None, direction_vectors=None, n_vertices=None, n_halfspaces=None, patch_args=None, vertex_args=None, center_args=None, autoscale_enable=True, decimal_precision=3)

Plot a polytopic approximation of the set.

Parameters:
  • method (str, optional) – Type of polytopic approximation to use. Can be [“inner” or “outer”]. Defaults to “inner”.

  • ax (axis object, optional) – Axis on which the patch is to be plotted

  • direction_vectors (array_like, optional) – Directions to use when performing ray shooting. Matrix (N times self.dim) for some N >= 1. Defaults to None, in which case we use pycvxset.common.spread_points_on_a_unit_sphere() to compute the direction vectors.

  • n_vertices (int, optional) – Number of vertices to use when computing the polytopic inner-approximation. Ignored if method is “outer” or when direction_vectors are provided. More than n_vertices may be used in some cases (see notes). Defaults to None.

  • n_halfspaces (int, optional) – Number of halfspaces to use when computing the polytopic outer-approximation. Ignored if method is “outer” or when direction_vectors are provided. More than n_halfspaces may be used in some cases (see notes). Defaults to None.

  • patch_args (dict, optional) – Arguments to pass for plotting faces and edges. See pycvxset.Polytope.Polytope.plot() for more details. Defaults to None.

  • vertex_args (dict, optional) – Arguments to pass for plotting vertices. See pycvxset.Polytope.Polytope.plot() for more details. Defaults to None.

  • center_args (dict, optional) – For ellipsoidal set, arguments to pass to scatter plot for the center. If a label is desired, pass it in center_args.

  • autoscale_enable (bool, optional) – When set to True, matplotlib adjusts axes to view full polytope. See pycvxset.Polytope.Polytope.plot() for more details. Defaults to True.

  • decimal_precision (int, optional) – When plotting a 3D polytope that is in V-Rep and not in H-Rep, we round vertex to the specified precision to avoid numerical issues. Defaults to PLOTTING_DECIMAL_PRECISION_CDD specified in pycvxset.common.constants.

Returns:

See pycvxset.Polytope.Polytope.plot() for details.

Return type:

tuple

Notes

This function is a wrapper for polytopic_inner_approximation() and polytopic_outer_approximation() for more details in polytope construction.

plus(Q)

Add a point or a set Q to a constrained zonotope (Minkowski sum).

Parameters:

Q (array_like | Polytope | ConstrainedZonotope) – The point or set to add

Raises:
  • TypeError – When Q is not one of the following — convertible into a 1D numpy array or a constrained zonotope.

  • ValueError – When Q has a dimension mismatch with self.

  • UserWarning – When using with Q that is a set, warn that an inner-approximation is returned.

Returns:

Minkowski sum of self and Q.

Return type:

ConstrainedZonotope

Notes

Given a constrained zonotope \(\mathcal{P}\) and a set \(Q\), this function computes the Minkowski sum of Q and the constrained zonotope, defined as \(\mathcal{R}=\{x + q|x\in\mathcal{P}, q\in\mathcal{Q}\}\). On the other hand, when Q is a point, this function computes the constrained zonotope \(\mathcal{R}=\{x + Q|x\in\mathcal{P}\}\).

This function implements (12) of [SDGR16] when Q is a constrained zonotope. When Q is a Polytope, this function converts it into a constrained zonotope, and then uses (12) of [SDGR16].

polytopic_inner_approximation(direction_vectors=None, n_vertices=None, verbose=False)

Compute a polytopic inner-approximation of a given set via ray shooting.

Parameters:
  • direction_vectors (array_like, optional) – Directions to use when performing ray shooting. Matrix (N times self.dim) for some \(N \geq 1\). Defaults to None.

  • n_vertices (int, optional) – Number of vertices to be used for the inner-approximation. n_vertices is overridden whenever direction_vectors are provided. Defaults to None.

  • verbose (bool, optional) – If true, pycvxset.common.spread_points_on_a_unit_sphere() is passed with verbose. Defaults to False.

Returns:

Polytopic inner-approximation in V-Rep of a given set with n_vertices no smaller than user-provided n_vertices.

Return type:

Polytope

Notes

We compute the polytope using extreme() evaluated along the direction vectors computed by pycvxset.common.spread_points_on_a_unit_sphere(). When direction_vectors is None and n_vertices is None, we select \(\text{n\_vertices} = 2 \text{self.dim} + 2^\text{self.dim} \text{SPOAUS\_DIRECTIONS\_PER\_QUADRANT}\) (as defined in pycvxset.common.constants()). [BV04]

polytopic_outer_approximation(direction_vectors=None, n_halfspaces=None, verbose=False)

Compute a polytopic outer-approximation of a given set via ray shooting.

Parameters:
  • direction_vectors (array_like, optional) – Directions to use when performing ray shooting. Matrix (N times self.dim) for some \(N \geq 1\). Defaults to None.

  • n_halfspaces (int, optional) – Number of halfspaces to be used for the inner-approximation. n_vertices is overridden whenever direction_vectors are provided. Defaults to None.

  • verbose (bool, optional) – If true, pycvxset.common.spread_points_on_a_unit_sphere() is passed with verbose. Defaults to False.

Returns:

Polytopic outer-approximation in H-Rep of a given set with n_halfspaces no smaller than user-provided n_vertices.

Return type:

Polytope

Notes

We compute the polytope using support() evaluated along the direction vectors computed by pycvxset.common.spread_points_on_a_unit_sphere(). When direction_vectors is None and n_halfspaces is None, we select \(\text{n\_halfspaces} = 2 \text{self.dim} + 2^\text{self.dim} \text{SPOAUS\_DIRECTIONS\_PER\_QUADRANT}\) (as defined in pycvxset.common.constants). [BV04]

project(x, p=2)[source]

Project a point or a collection of points on to a set.

Given a set \(\mathcal{P}\) and a test point \(y\in\mathbb{R}^{\mathcal{P}.\text{dim}}\), this function solves a convex program,

\[\begin{split}\text{minimize} &\quad \|x - y\|_p\\ \text{subject to} &\quad x \in \mathcal{P}\\\end{split}\]
Parameters:
  • points (array_like) – Points to project. Matrix (N times self.dim), where each row is a point.

  • p (str | int) – Norm-type. It can be 1, 2, or ‘inf’. Defaults to 2, which is the Euclidean norm.

Raises:
  • ValueError – Set is empty

  • ValueError – Dimension mismatch — no. of columns in points is different from self.dim.

  • ValueError – Points is not convertible into a 2D array

  • NotImplementedError – Unable to solve problem using CVXPY

Returns:

A tuple with two items:
  1. projected_point (numpy.ndarray): Projection point(s) as a 2D numpy.ndarray. Matrix (N times self.dim), where each row is a projection of the point in points to the set \(\mathcal{P}\).

  2. distance (numpy.ndarray): Distance(s) as a 1D numpy.ndarray. Vector (N,), where each row is a projection of the point in points to the set \(\mathcal{P}\).

Return type:

tuple

Notes

Given a constrained zonotope \(\mathcal{P}\) and a test point \(y\in\mathbb{R}^{\mathcal{P}.\text{dim}}\), this function solves a convex program with decision variables \(x\in\mathbb{R}^{\mathcal{P}.\text{dim}}\) and \(\xi\in\mathbb{R}^{\mathcal{P}.\text{latent\_dim}}\),

\[\begin{split}\text{minimize} &\quad \|x - y\|_p\\ \text{subject to} &\quad x = G_P \xi + c_P\\ &\quad A_P \xi = b_P\\ &\quad \|\xi\|_\infty \leq 1\end{split}\]
projection(project_away_dims)[source]

Orthogonal projection of a set \(\mathcal{P}\) after removing some user-specified dimensions.

\[\mathcal{R} = \{r \in \mathbb{R}^{m}\ |\ \exists v \in \mathbb{R}^{n - m},\ \text{Lift}(r,v)\in \mathcal{P}\}\]

Here, \(m = \mathcal{P}.\text{dim} - \text{length}(\text{project\_away\_dim})\), and \(\text{Lift}(r,v)\) lifts (“undo”s the projection) using the appropriate components of v. This function uses affine_map() to implement the projection by designing an appropriate affine map \(M \in \{0,1\}^{m\times\mathcal{P}.\text{dim}}\) with each row of \(M\) corresponding to some standard axis vector \(e_i\in\mathbb{R}^m\).

Parameters:

project_away_dims (array_like) – Dimensions to projected away in integer interval [0, 1, …, n - 1].

Raises:

ValueError – When project_away_dims are not in the integer interval | All dimensions are projected away

Returns:

m-dimensional set obtained via projection.

Return type:

ConstrainedZonotope

remove_redundancies()

Remove any redundancies in the equality system using pycddlib and other geometric properties

slice(dims, constants)[source]

Slice a set restricting certain dimensions to constants.

This function uses intersection_with_affine_set() to implement the slicing by designing an appropriate affine set from dims and constants.

Parameters:
  • dims (array_like) – List of dims to restrict to a constant in the integer interval [0, 1, …, n - 1].

  • constants (array_like) – List of constants

Raises:
  • ValueError – dims has entries beyond n

  • ValueError – dims and constants are not 1D arrays of same size

Returns:

Constrained zonotope that has been sliced at the specified dimensions.

Return type:

ConstrainedZonotope

slice_then_projection(dims, constants)[source]

Wrapper for slice() and projection().

The function first restricts a set at certain dimensions to constants, and then projects away those dimensions. Useful for visual inspection of higher dimensional sets.

Parameters:
  • dims (array_like) – List of dims to restrict to a constant in the integer interval [0, 1, …, dim - 1], and then project away.

  • constants (array_like) – List of constants

Raises:
  • ValueError – dims has entries beyond n

  • ValueError – dims and constants are not 1D arrays of same size

  • ValueError – When dims are not in the integer interval | All dimensions are projected away

Returns:

m-dimensional set obtained via projection after slicing.

Return type:

ConstrainedZonotope

support(eta)[source]

Evaluates the support function and support vector of a set.

The support function of a set \(\mathcal{P}\) is defined as \(\rho_{\mathcal{P}}(\eta) = \max_{x\in\mathcal{P}} \eta^\top x\). The support vector of a set \(\mathcal{P}\) is defined as \(\nu_{\mathcal{P}}(\eta) = \arg\max_{x\in\mathcal{P}} \eta^\top x\).

Parameters:

eta (array_like) – Support directions. Matrix (N times self.dim), where each row is a support direction.

Raises:
  • ValueError – Set is empty

  • ValueError – Mismatch in eta dimension

  • ValueError – eta is not convertible into a 2D array

  • NotImplementedError – Unable to solve problem using CVXPY

Returns:

A tuple with two items:
  1. support_function_evaluations (numpy.ndarray): Support function evaluation(s) as a 2D numpy.ndarray. Vector (N,) with as many rows as eta.

  2. support_vectors (numpy.ndarray): Support vectors as a 2D numpy.ndarray. Matrix N x self.dim with as many rows as eta.

Return type:

tuple

Notes

Given a constrained zonotope \(\mathcal{P}\) and a support direction \(\eta\in\mathbb{R}^{\mathcal{P}.\text{dim}}\), this function solves a convex program with decision variables \(x\in\mathbb{R}^{\mathcal{P}.\text{dim}}\) and \(\xi\in\mathbb{R}^{\mathcal{P}.\text{latent\_dim}}\),

\[\begin{split}\text{maximize} &\quad \eta^\top x\\ \text{subject to} &\quad x = G_P \xi + c_P\\ &\quad A_P \xi = b_P\\ &\quad \|\xi\|_\infty \leq 1\end{split}\]
property Ae

Equality coefficient vectors Ae for the constrained zonotope.

Returns:

Equality coefficient vectors Ae for the constrained zonotope. Ae is np.empty((0, self.latent_dim)) for a zonotope.

Return type:

numpy.ndarray

property G

Affine transformation matrix G for the constrained zonotope.

Returns:

Affine transformation matrix G.

Return type:

numpy.ndarray

property He

Equality constraints He=[Ae, be] for the constrained zonotope.

Returns:

H-Rep in [Ae, be]. He is np.empty((0, self.latent_dim + 1)) for a zonotope.

Return type:

numpy.ndarray

property be

Equality constants be for the constrained zonotope.

Returns:

Equality constants be for the constrained zonotope. be is np.empty((0,)) for a zonotope.

Return type:

numpy.ndarray

property c

Affine transformation vector c for the constrained zonotope.

Returns:

Affine transformation vector c.

Return type:

numpy.ndarray

property cvxpy_args_lp

CVXPY arguments in use when solving a linear program

Returns:

CVXPY arguments in use when solving a linear program. Defaults to dictionary in pycvxset.common.DEFAULT_CVXPY_ARGS_LP.

Return type:

dict

property cvxpy_args_socp

CVXPY arguments in use when solving a second-order cone program

Returns:

CVXPY arguments in use when solving a second-order cone program. Defaults to dictionary in pycvxset.common.DEFAULT_CVXPY_ARGS_SOCP.

Return type:

dict

property dim

Dimension of the constrained zonotope.

Returns:

Dimension of the constrained zonotope.

Return type:

int

Notes

We determine dimension from G, since c is set to None in case of empty (constrained) zonotope.

property is_bounded

Check if the constrained zonotope is bounded (which is always True)

property is_empty

Check if the constrained zonotope is empty

Raises:

NotImplementedError – Unable to solve the feasibility problem using CVXPY

Returns:

When True, the polytope is empty

Return type:

bool

property is_full_dimensional

Check if the affine dimension of the constrained zonotope is the same as the constrained zonotope dimension

Returns:

True when the affine hull containing the constrained zonotope has the dimension self.dim

Return type:

bool

Notes

An empty polytope is full dimensional if dim=0, otherwise it is not full-dimensional. See Sec. 2.1.3 of [BV04] for discussion on affine dimension.

A non-empty zonotope is full-dimensional if and only if G has full row rank.

A non-empty constrained zonotope is full-dimensional if and only if [G; A] has full row rank.

property is_singleton

Check if the constrained zonotope is a singleton

property is_zonotope

Check if the constrained zonotope is a zonotope

property latent_dim

Latent dimension of the constrained zonotope.

Returns:

Latent dimension of the constrained zonotope.

Return type:

int

property n_equalities

Number of equality constraints used when defining the constrained zonotope.

Returns:

Number of equality constraints used when defining the constrained zonotope

Return type:

int