pycvxset.Polytope (API details)¶
|
Polytope class. |
- class pycvxset.Polytope.Polytope(**kwargs)[source]¶
Bases:
object
Polytope class.
Polytope object construction admits one of the following combinations (as keyword arguments):
(A, b) for a polytope in halfspace representation (H-rep) \(\{x\ |\ Ax \leq b\}\),
(A, b, Ae, be) for a polytope in halfspace representation (H-rep) with equality constraints \(\{x\ |\ Ax \leq b, A_e x = b_e\}\),
V for a polytope in vertex representation (V-rep) — \(\text{ConvexHull}(v_i)\) where \(v_i\) are rows of matrix V,
(lb, ub) for an axis-aligned cuboid with appropriate bounds \(\{x\ |\ lb\leq x \leq ub\}\), and
(c, h) for an axis-aligned cuboid centered at c with specified scalar/vector half-sides \(h\), \(\{x\ |\ \forall i\in\{1,\cdots,n\}, |x_i - c_i| \leq h_i\}.\)
dim for an empty Polytope of dimension dim (no argument generates a zero-dimensional empty Polytope),
- Parameters:
dim (int, optional) – Dimension of the empty polytope. If NOTHING is provided, dim=0 is assumed.
V (array_like, optional) – List of vertices of the polytope (V-Rep). The list must be 2-dimensional with vertices arranged row-wise and the polytope dimension determined by the column count.
A (array_like, optional) – Inequality coefficient vectors (H-Rep). The vectors are stacked vertically with the polytope dimension determined by the column count. When A is provided, b must also be provided.
b (array_like, optional) – Inequality constants (H-Rep). The constants are expected to be in a 1D numpy array. When b is provided, A must also be provided.
Ae (array_like) – Equality coefficient vectors (H-Rep). The vectors are stacked vertically with matching number of columns as A. When Ae is provided, A, b, and be must also be provided.
be (array_like) – Equality coefficient constants (H-Rep). The constants are expected to be in a 1D numpy array. When be is provided, A, b, and Ae must also be provided.
lb (array_like, optional) – Lower bounds of the axis-aligned cuboid. Must be 1D array, and the polytope 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.
c (array_like, optional) – Center of the axis-aligned cuboid. Must be 1D array, and the polytope dimension is determined by number of elements in c. When c is provided, h 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.
- Raises:
ValueError – When arguments provided is not one of [(A, b), (A, b, Ae, be), (lb, ub), (c, h), V, dim, NOTHING]
ValueError – Errors raised by issues with (A, b) and (Ae, be) — mismatch in dimensions, not convertible to appropriately-dimensioned numpy arrays, incompatible systems (A, b) and (Ae, be) etc.
ValueError – Errors raised by issues with V — not convertible to a 2-D numpy array, etc.
ValueError – Errors raised by issues with lb, ub — mismatch in dimensions, not convertible to 1D numpy arrays, etc.
ValueError – Errors raised by issues with c, h — mismatch in dimensions, not convertible to 1D numpy arrays, etc.
ValueError – Polytope is not bounded in any direction. We use a sufficient condition for speed, and therefore the detection may not be exhaustive.
ValueError – Polytope is not bounded in some directions. We use a sufficient condition for speed, and therefore the detection may not be exhaustive.
UserWarning – If some rows are removed from (A, b) due to all zeros in A or np.inf in b | (Ae, be) when all zeros in Ae and be.
- affine_map(M)¶
Compute the matrix times set.
- Parameters:
M (array_like) – A vector or array (0, 1, or 2-D numpy.ndarray-like object) with self.dim columns
- Raises:
ValueError – When M can not be converted EXACTLY into a 2D array of float
ValueError – When M does not have self.dim columns
- Returns:
The scaled polytope \(\mathcal{R} = M \mathcal{P} = \{M x: x\in \mathcal{P}\}\)
- Return type:
Notes
This function requires \(\mathcal{P}\) to be in V-Rep, and performs a vertex enumeration when \(\mathcal{P}\) is in H-Rep.
- chebyshev_centering()¶
Computes a ball with the largest radius that fits within the polytope. The ball’s center is known as the Chebyshev center, and its radius is the Chebyshev radius.
- Raises:
NotImplementedError – Unable to solve chebyshev centering problem using CVXPY
- Returns:
- A tuple with two items
center (numpy.ndarray): Chebyshev center of the polytope
radius (float): Chebyshev radius of the polytope
- Return type:
tuple
Notes
This function is called by the constructor for non-emptiness and boundedness check. It uses limited attributes (dim, A, b, Ae, be). This function requires H-Rep, and will perform a halfspace enumeration when a V-Rep polytope is provided.
We solve the LP (see Section 8.5.1 in [BV04]) for \(c\) (for center) and \(R\) (for radius),
\[\begin{split}\text{maximize} &\quad R \\ \text{subject to} &\quad A c + R ||A||_\text{row} \leq b,\\ &\quad A_e c = b_e, \\ &\quad R \geq 0,\end{split}\]where \(||A||_\text{row}\) is a vector of dimension
n_halfspaces
with each element as \(||a_i||_2\). When (Ae, be) is non-empty, then R is forced to zero post-solve. Consequently, chebyshev_centering also serves a method to find a feasible point in the relative interior of the polytope.- We can also infer the following from the Chebyshev radius R:
\(R=\infty\), the polytope is unbounded. Note that, this is just a sufficient condition, and an unbounded polytope can have a finite Chebyshev radius. For example, consider a 3-dimensional axis-aligned cuboid \([-1, 1] \times [-1, 1] \times \mathbb{R}\).
\(0 < R < \infty\), the polytope is nonempty and full-dimensional.
\(R=0\), the polytope is nonempty and but not full-dimensional.
\(R=-\infty\), the polytope is empty.
- 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 polytope.
- 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 polytope is empty
- Returns:
A tuple with two items:
constraint_list (list): CVXPY constraints for the containment of x in the polytope.
theta (cvxpy.Variable | None): CVXPY variable representing the convex combination coefficient when polytope is in V-Rep. It is None when the polytope is in H-Rep or empty.
- Return type:
tuple
- contains(Q)¶
Check containment of a set \(\mathcal{Q}\) (could be a polytope, an ellipsoid, or a constrained zonotope), or a collection of points \(Q \in \mathbb{R}^{n_Q \times \mathcal{P}.\text{dim}}\) in the polytope \(\mathcal{P}\).
- Parameters:
Q (array_like | Polytope | ConstrainedZonotope | Ellipsoid) – Set or a collection of points (each row is a point) to be tested for containment within \(\mathcal{P}\). When providing a collection of points, Q is a matrix (N times self.dim) with each row is a point.
- Raises:
ValueError – When two polytopes | the polytope and the test point(s) are NOT of the same dimension
- Returns:
When Q is a polytope, a bool is returned which is True if and only if \(Q\subseteq P\). On the other hand, if Q is array_like (Q is a point or a collection of points), then an numpy.ndarray of bool is returned, with as many elements as the number of rows in Q.
- Return type:
bool | numpy.ndarray[bool]
Notes
Containment of polytopes: This function accommodates the following combinations of representations for \(\mathcal{P}\) and \(\mathcal{Q}\). It eliminates the need for enumeration by comparing support function evaluations almost always. When P is H-Rep and Q is V-Rep, then we perform a quick check for containment of vertices of Q in P. Otherwise, * Q is empty: Return True * P is empty: Return False
Containment of points: This function accommodates \(\mathcal{P}\) to be in H-Rep or V-Rep. When testing if a point is in a polytope where only V-Rep is available for self, we solve a collection of second-order cone programs for each point using
project()
(check if distance between v and \(\mathcal{P}\) is nearly zero). Otherwise, we use the observation that for \(\mathcal{P}\) with (A, b, Ae, be), then \(v\in\mathcal{P}\) if and only if \(Av\leq b\) and \(A_ev=b_e\). For numerical precision considerations, we usenumpy.isclose()
.
- classmethod deflate_rectangle(set_to_be_centered)¶
Compute the rectangle with the smallest volume that contains the given set.
- Parameters:
set_to_be_centered (Polytope | ConstrainedZonotope | Ellipsoid) – Set to compute the.
- Returns:
Minimum volume circumscribing rectangle
- Return type:
Notes
This function is a wrapper for
minimum_volume_circumscribing_rectangle()
of attr:set_to_be_centered. Please check that function for more details including raising exceptions.
- determine_H_rep()¶
Determine the halfspace representation from a given vertex representation of the polytope.
- Raises:
ValueError – When H-rep computation fails OR Polytope is not bounded!
Notes
When the set is empty, we define an empty polytope.
Otherwise, we use cdd for the halfspace enumeration.
The computed vertex representation need not be minimal.
- determine_V_rep()¶
Determine the vertex representation from a given halfspace representation of the polytope.
- Raises:
ValueError – Vertex enumeration yields rays, which indicates that the polytope is unbounded. Numerical issues may also be a culprit.
ValueError – Polytope is not bounded!
Notes
We use cdd for the vertex enumeration. cdd uses the halfspace representation \([b, -A]\) where \(b - Ax \geq 0 \Leftrightarrow Ax \leq b\).
For a polyhedron described as
\[\mathcal{P} = \text{conv}(v_1, ..., v_n) + \text{nonneg}(r_1, ..., r_s),\]the V-representation matrix in cdd is [t V] where t is the column vector with n ones followed by s zeroes, and V is the stacked matrix of n vertex row vectors on top of s ray row vectors. pycvxset uses only bounded polyhedron, so we should never observe rays.
- 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(point_type='centroid')¶
Compute a point in the interior of the polytope. When the polytope is not full-dimensional, the point may lie on the boundary.
- Parameters:
point_type (str) – Type of interior point. Valid strings: {‘centroid’, ‘chebyshev’}. Defaults to ‘centroid’.
- Raises:
NotImplementedError – When an invalid point_type is provided.
ValueError – When the polytope provided is empty.
ValueError – When the polytope provided is not bounded.
- Returns:
A feasible point in the polytope in the interior as a 1D array
- Return type:
numpy.ndarray
Notes
point_type is ‘centroid’: Computes the average of the vertices. The function requires the polytope to be in V-Rep, and a vertex enumeration is performed if the polytope is in H-Rep.
point_type is ‘chebyshev’: Computes the Chebyshev center. The function requires the polytope to be in H-Rep, and a halfspace enumeration is performed if the polytope is in V-Rep.
- intersection(Q)¶
Intersect the polytope \(\mathcal{P}\) with another polytope \(\mathcal{Q}\).
- Parameters:
Q (Polytope) – Polytope to be intersected with self
- Raises:
ValueError – Q is not a polytope | Mismatch in dimensions
- Returns:
The intersection of P and Q
- Return type:
Notes
This function requires \(\mathcal{P}\) and Q to be of the same dimension and in H-Rep, and performs halfspace enumerations when P or Q are in V-Rep.
- intersection_under_inverse_affine_map(Q, R)¶
Compute the intersection of constrained zonotope under an inverse affine map
- Parameters:
Q (Polytope) – Set to intersect with
R (array_like) – Matrix of dimension Y.dim times self.dim
- Raises:
ValueError – When Q is not a Polytope
ValueError – When R is not of correct dimension
ValueError – When self is not bounded
- Returns:
The intersection of a polytope with another polytope under an inverse affine map. Specifically, given polytopes \(\mathcal{P}\) (self) and \(\mathcal{Q}\), and a matrix \(R\), we compute the set \(\{x \in \mathcal{P}| Rx \in \mathcal{Q}\}\).
- Return type:
Notes
This function requires both polytopes to be in H-Rep. Halfspace enumeration is performed when necessary.
Unlike
inverse_affine_map_under_invertible_matrix()
, this function does not require R to be invertible or square and also accommodates Q to be an unbounded polytope. However, self must be a bounded set.
- intersection_with_affine_set(Ae, be)¶
Intersect the polytope with an affine set.
- Parameters:
Ae (array_like) – Equality coefficient vectors (H-Rep). The vectors are stacked vertically.
be (array_like) – Equality constants (H-Rep). The constants are expected to be in a 1D numpy array.
- Raises:
ValueError – Mismatch in dimensions | (Ae, be) is not a valid collection of equality constraints.
- Returns:
The intersection of \(P\) and \(\{x:A_e x = b_e\}\)
- Return type:
Notes
This function requires \(\mathcal{P}.\text{dim}\) = Ae.shape[1] and \(\mathcal{P}\) should be in H-Rep, and performs halfspace enumeration if \(\mathcal{P}\) in V-Rep.
- intersection_with_halfspaces(A, b)¶
Intersect the polytope with a collection of halfspaces.
- Parameters:
A (array_like) – Inequality coefficient vectors (H-Rep). The vectors are stacked vertically.
b (array_like) – Inequality constants (H-Rep). The constants are expected to be in a 1D numpy array.
- Raises:
ValueError – Mismatch in dimensions | (A, b) is not a valid collection of halfspaces
- Returns:
The intersection of \(P\) and \(\{x\in\mathbb{R}^n\ |\ Ax\leq b\}\)
- Return type:
Notes
This function requires \(\mathcal{P}.\text{dim}\) = A.shape[1] and \(\mathcal{P}\) should be in H-Rep, and performs halfspace enumeration if \(\mathcal{P}\) in V-Rep.
- inverse_affine_map_under_invertible_matrix(M)¶
Compute the set times matrix, the inverse affine map of the set under an invertible matrix M.
- Parameters:
M (array_like) – A invertible array of size self.dim times self.dim
- Raises:
ValueError – When self is empty
TypeError – When M is not convertible into a 2D numpy array of float
ValueError – When M is not a square matrix
ValueError – When M is not invertible
- Returns:
The inverse-scaled polytope \(\mathcal{R} = \mathcal{P} \times M = \{x: M x\in \mathcal{P}\}\)
- Return type:
Notes
This function accommodates \(\mathcal{P}\) to be in H-Rep or in V-Rep. When \(\mathcal{P}\) is in H-Rep, \(\mathcal{P} M=\{x|Mx\in\mathcal{P}\}=\{x|AMx\leq b, A_eMx = b_e\}\). On the other hand, when \(\mathcal{P}\) is in V-Rep, \(\mathcal{P} M=ConvexHull(M^{-1}v_i)\).
We require M to be invertible in order to ensure that the resulting set is representable as a polytope.
- maximum_volume_inscribing_ellipsoid()¶
Compute the maximum volume ellipsoid that fits within the given polytope.
- Raises:
ValueError – When polytope is empty or has non-empty (Ae, be)
NotImplementedError – Unable to solve convex problem using CVXPY
- Returns:
- A tuple with three items:
center (numpy.ndarray): Maximum volume inscribed ellipsoid’s center
shape_matrix (numpy.ndarray): Maximum volume inscribed ellipsoid’s shape matrix
sqrt_shape_matrix (numpy.ndarray): Maximum volume inscribed ellipsoid’s square root of shape matrix. Returns None if the polytope is not full-dimensional.
- Return type:
tuple
Notes
This function requires H-Rep, and will perform a vertex enumeration when a V-Rep polytope is provided.
We solve the SDP for a positive definite matrix \(B\in\mathbb{S}^n_{++}\) and \(d\in\mathbb{R}^n\),
\[\begin{split}\text{maximize} &\quad \log \text{det} B \\ \text{subject to} &\quad \|B a_i\|_2 + a_i^T d \leq b_i,\\ &\quad A_e d = b_e,\\ &\quad B A_e^\top = 0,\end{split}\]where \((a_i,b_i)\) is the set of hyperplanes characterizing \(\mathcal{P}\), and the inscribing ellipsoid is given by \(\{Bu + d| {||u||}_2 \leq 1\}\). The center of the ellipsoid is given by \(c = d\), and the shape matrix is given by \(Q = (B B)^{-1}\). See [EllipsoidalTbx-Min_verticesolEll] and Section 8.4.2 in [BV04] for more details. The last two constraints arise from requiring the inscribing ellipsoid to lie in the affine set \(\{A_e x = b_e\}\).
When the polytope is full-dimensional, we can instead solve a second-order cone program (SOCP). Consider the full-dimensional ellipsoid \(\{Lu + c| {\|u\|}_2 \leq 1\}\), where \(G\) is a square, lower-triangular matrix with positive diagonal entries. Then, we solve the following (equivalent) optimization problem:
\[\begin{split}\text{minimize} &\quad \text{geomean}(L) \\ \text{subject to} &\quad \text{diag}(L) \geq 0\\ &\quad \|L^T a_i\|_2 + a_i^T d \leq b_i,\end{split}\]with decision variables \(G\) and \(c\). Here, we use the observation that \(\text{geomean}(L)\) is a monotone function of \(\log\det(GG^T)\) (the volume of the ellipsoid).
- 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:
x.value (numpy.ndarray): Optimal value of x. np.nan * np.ones((self.dim,)) if the problem is not solved.
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.
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
- minimize_H_rep()¶
Remove any redundant inequalities from the halfspace representation of the polytope using cdd.
- Raises:
ValueError – When minimal H-Rep computation fails OR polytope is not bounded!
- minimize_V_rep(prefer_qhull_over_cdd=True)¶
Remove any redundant vertices from the vertex representation of the polytope.
- Parameters:
prefer_qhull_over_cdd (bool, optional) – When True, minimize_V_rep uses qhull. Otherwise, we use cdd.
- Raises:
ValueError – When minimal V-Rep computation fails!
Notes
Use cdd or qhull for the reduction of vertices.
- minimum_volume_circumscribing_ellipsoid()¶
Compute the minimum volume ellipsoid that covers the given polytope (also known as Lowner-John Ellipsoid).
- Raises:
ValueError – When polytope is empty
NotImplementedError – Unable to solve convex problem using CVXPY
- Returns:
- A tuple with three items:
center (numpy.ndarray): Minimum volume circumscribed ellipsoid’s center
shape_matrix (numpy.ndarray): Minimum volume circumscribed ellipsoid’s shape matrix
sqrt_shape_matrix (numpy.ndarray): Minimum volume circumscribed ellipsoid’s square root of shape matrix. Returns None if the polytope is not full-dimensional.
- Return type:
tuple
Notes
This function requires V-Rep, and will perform a vertex enumeration when a H-Rep polytope is provided.
We solve the SDP for a positive definite matrix \(A\in\mathbb{S}^n_{++}\) and \(b\in\mathbb{R}^n\),
\[\begin{split}\text{maximize} &\quad \log \text{det} A^{-1} \\ \text{subject to} &\quad \|Av_i + b\|_2 \leq 1,\ \forall \text{ vertices } v_i \text{ of } \mathcal{P}\end{split}\]where the circumscribing ellipsoid is given by \(\{x| {||A x + b||}_2 \leq 1\}\), and we use the observation that \(\log \text{det} A^{-1}= -\log \text{det} A\). The center of the ellipsoid is given by \(c = -A^{-1}b\), and the shape matrix is given by \(Q = (A^T A)^{-1}\). See [EllipsoidalTbx-Min_verticesolEll] and Section 8.4.1 in [BV04] for more details.
When the polytope is full-dimensional, we can instead solve a second-order cone program (SOCP). Consider the full-dimensional ellipsoid \(\{Lu + c| {\|u\|}_2 \leq 1\}\), where \(G\) is a square, lower-triangular matrix with positive diagonal entries. Then, we solve the following (equivalent) optimization problem:
\[\begin{split}\text{minimize} &\quad \text{geomean}(L) \\ \text{subject to} &\quad \text{diag}(L) > 0\\ &\quad v_i = Lu_i + c,\ \forall \text{ vertices } v_i \text{ of } \mathcal{P},\\ &\quad u_i\in\mathbb{R}^,\ \|u_i\|_2 \leq 1,\end{split}\]with decision variables \(G\), \(c\), and \(u_i\), and \(\text{geomean}\) is the geometric mean of the diagonal elements of \(G\). Here, we use the observation that \(\text{geomean}(L)\) is a monotone function of \(\log\det(GG^T)\) (the volume of the ellipsoid). For the sake of convexity, we solve the following equivalent optimization problem after a change of variables \(L_\text{inv}=L^{-1}\) and \(c_{L_\text{inv}}=L^{-1}c\), and substituting for the variables \(u_i\):
\[\begin{split}\text{maximize} &\quad \text{geomean}(L_\text{inv}) \\ \text{subject to} &\quad \text{diag}(L_\text{inv}) > 0\\ &\quad \|L_\text{inv}v_i - c_{L_\text{inv}}\|_2 \leq 1,\ \forall \text{ vertices } v_i\end{split}\]The center of the ellipsoid is \(c = Lc_{L_\text{inv}}\), and the shape matrix is \(Q = GG^T\), where \(L=L_\text{inv}^{-1}\).
- minimum_volume_circumscribing_rectangle()¶
Compute the minimum volume circumscribing rectangle for a given polytope
- Raises:
ValueError – When polytope is empty
ValueError – When polytope is unbounded
- Returns:
- A tuple of two items
lb (numpy.ndarray): Lower bound \(l\) on the polytope, \(\mathcal{P}\subseteq\{l\}\oplus\mathbb{R}_{\geq 0}\).
ub (numpy.ndarray): Upper bound \(u\) on the polytope, \(\mathcal{P}\subseteq\{u\}\oplus(-\mathbb{R}_{\geq 0})\).
- Return type:
tuple
Notes
This function accommodates H-Rep and V-Rep. The lower/upper bound for V-Rep is given by an element-wise minimization/maximization operation, while the lower/upper bound for H-Rep is given by an element-wise support computation (2 * self.dim linear programs). For a H-Rep polytope, the function uses
support()
and is a wrapper forpycvxset.common.minimum_volume_circumscribing_rectangle()
.This function is called to check for boundedness of a polytope.
- minus(Q)¶
Implement - operation (Pontryagin difference when Q is a polytope, translation by -Q when Q is a point)
- Parameters:
Q (array_like | Polytope | ConstrainedZonotope | Ellipsoid) – Polytope to be subtracted in Pontryagin
polytope (difference sense from self or a vector for negative translation of the)
- Raises:
TypeError – When Q is neither a Polytope or a point
ValueError – When Q is not of the same dimension as self
- Returns:
The polytope \(\mathcal{R}\) that is the Pontryagin difference of P and Q or a negative translation of P by Q.
- Return type:
Notes
Subtraction with a point: This function accommodates \(\mathcal{P}\) in H-Rep and V-Rep. See
plus()
for more details.Subtraction with a polytope: This function requires \(\mathcal{P}\) and Q to be of the same dimension and :math:mathcal{P} in H-Rep, and performs halfspace enumerations when \(\mathcal{P}\) is in V-Rep. The H-rep of the polytope \(\mathcal{R}\) is, \(H_{\mathcal{R}} = [\mathcal{P}.A, \mathcal{P}.b - \rho_{ \mathcal{Q}}(\mathcal{P}.A)]\) where \(\rho_{\mathcal{Q}}\) is the support function (see
support()
). [KG98]
- normalize()¶
Normalize a H-Rep such that each row of A has unit \(\ell_2\)-norm.
Notes
This function requires P to be in H-Rep. If P is in V-Rep, a halfspace enumeration is performed.
- plot(ax=None, patch_args=None, vertex_args=None, autoscale_enable=True, decimal_precision=3)¶
Plot a 2D or 3D polytope.
- Parameters:
ax (axis object, optional) – Axis on which the patch is to be plotted
patch_args (dict, optional) – Arguments to pass for plotting faces and edges. See [Matplotlib-patch] for options for patch_args for 2D and [Matplotlib-Line3DCollection] for options for patch_args for 3D. Defaults to None, in which case we set edgecolor to black, and facecolor to skyblue.
vertex_args (dict, optional) – Arguments to pass for plotting vertices. See [Matplotlib-scatter] for options for vertex_args for 2D and [Matplotlib-Axes3D.scatter] for options for vertex_args for 3D. Defaults to None, in which case we skip plotting the vertices.
autoscale_enable (bool, optional) – When set to True, matplotlib adjusts axes to view full polytope. 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.
- Raises:
ValueError – When an incorrect axes is provided.
NotImplementedError – When a polytope.dim >= 4 or == 1 | autoscale_enable is set to False for 3D plotting.
UserWarning – When an empty polytope or an unbounded polytope is provided.
UserWarning – In 3D plotting, when all faces have less than 3 vertices
- Returns:
Tuple with axes containing the polygon, handle for plotting patch, handle for plotting vertices
- Return type:
(axes, handle, handle)
Notes
plot is a wrapper for
plot2d()
andplot3d()
functions. See their documentation for more details.Vertex-halfspace enumeration for 2D polytopes: If a 2D polytope is in H-Rep, vertex enumeration is performed to plot the polytope in 2D. No vertex enumeration is performed for a 2D polytope in V-Rep. This function calls minimize_V_rep to simplify plotting.
Vertex-halfspace enumeration for 3D polytopes: The 3D polytope needs to have both V-Rep and H-Rep. Consequently, at least a halfspace/vertex enumeration is performed for a 3D polytope in single representation.
Rounding vertices: This function calls
pycvxset.common.prune_and_round_vertices()
when a 3D polytope in V-Rep is given to be plotted.Handle returned for the patches: For 2D plot, axis handle of the single patch is returned. For 3D plotting, multiple patches are present, and plot returns the first patch’s axis handle. In this case, the order of patches are determined by the use of pycddlib.copy_input_incidence from the from the specified (A, b). Also, labeling the first patch is done using Poly3DCollection.
Disabling face colors: We can plot the polytope frames without filling it by setting patch_args[‘facecolor’] = None or patch_args[‘fill’] = False. If visual comparison of 3D polytopes is desired, it may be better to plot just the polytope frames instead by setting patch_args[‘facecolor’] = None.
Warning
This function may produce erroneous-looking plots when visually comparing multiple 3D sets. For more info on matplotlib limitations, see https://matplotlib.org/stable/api/toolkits/mplot3d/faq.html#my-3d-plot-doesn-t-look-right-at-certain-viewing-angles. If visual comparison is desired, it may be better to plot just the polytope frame instead by setting patch_args[‘facecolor’] = None.
- plot2d(ax=None, patch_args=None, vertex_args=None, autoscale_enable=True)¶
Plot a 2D polytope using matplotlib’s add_patch(Polygon()).
- Parameters:
ax (axis object, optional) – Axis on which the patch is to be plotted
patch_args (dict, optional) – Arguments to pass for plotting faces and edges. See [Matplotlib-patch] for options for patch_args. Defaults to None, in which case we set edgecolor to black, and facecolor to skyblue.
vertex_args (dict, optional) – Arguments to pass for plotting vertices. See [Matplotlib-scatter] for options for vertex_args. Defaults to None, in which case we skip plotting the vertices.
autoscale_enable (bool, optional) – When set to True (default), matplotlib adjusts axes to view full polytope.
- Raises:
ValueError – When a 2D polytope is provided
- Returns:
Tuple with axes containing the polygon, handle for plotting patch, handle for plotting vertices
- Return type:
(axes, handle, handle)
Notes
This function requires polytope in V-Rep, and performs a vertex enumeration if the polytope is H-Rep.
We sort the vertices in counter-clockwise direction with respect to the centroid, and then plot it using matplotlib’s Polygon. We can plot just the polytope frames without filling it by setting patch_args[‘facecolor’] = None or patch_args[‘fill’] = False.
- plot3d(ax=None, patch_args=None, vertex_args=None, autoscale_enable=True, decimal_precision=3)¶
Plot a 3D polytope using matplotlib’s Line3DCollection.
- Parameters:
ax (Axes, optional) – Axes to plot. Defaults to None, in which case a new axes is created. The function assumes that the provided axes was defined with projection=’3d’.
patch_args (dict, optional) – Arguments to pass for plotting faces and edges. See [Matplotlib-Line3DCollection] for options for patch_args. Defaults to None, in which case we set edgecolor to black, and facecolor to skyblue.
vertex_args (dict, optional) – Arguments to pass for plotting vertices. See [Matplotlib-Axes3D.scatter] for options for vertex_args. Defaults to None, in which case we skip plotting the vertices.
autoscale_enable (bool, optional) – When set to True (default), matplotlib adjusts axes to view full polytope.
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.
- Raises:
ValueError – When either a non-3D polytope is provided.
UserWarning – When all faces have less than 3 vertices
NotImplementedError – autoscale_enable needs to be currently enabled for 3D plotting
- Returns:
Tuple with axes containing the polygon, handle for plotting patch, handle for plotting vertices
- Return type:
(axes, handle, handle)
Notes
This function requires polytope in H-Rep and V-Rep, and uses CDD to compute incidences, and vertices. Consequently, at least one halfspace/vertex enumeration is performed.
Since 3D plotting involves multiple patches, the first patch’s axis handle is returned. In this case, the order of patches are determined by the use of pycddlib.copy_input_incidence from the specified (A, b).
This function calls
pycvxset.common.prune_and_round_vertices()
when a 3D polytope in V-Rep is given to be plotted.When label is passed in patch_args, the label is only applied to the first patch. Subsequent patches will not have a label. Moreover, the first patch is plotted using Poly3DCollection, while the subsequent patches are plotted using Line3DCollection in this case. Otherwise, when no label is provided, all patches are plotted using Line3DCollection.
We iterate over each halfspace, rotate the halfspace about its centroid so that halfspace is now parallel to XY plane, and then sort the vertices based on their XY coordinates in counter-clockwise direction with respect to the centroid, and then plot it using matplotlib’s Line3DCollection.
Warning
This function may produce erroneous-looking plots when visually comparing multiple 3D sets. For more info on matplotlib limitations, see https://matplotlib.org/stable/api/toolkits/mplot3d/faq.html#my-3d-plot-doesn-t-look-right-at-certain-viewing-angles. If visual comparison is desired, it may be better to plot just the polytope frame instead by setting patch_args[‘facecolor’] = None.
- plus(Q)¶
Add a point or a set to the polytope
- Parameters:
Q (array_like | Polytope) – Point or set to add to the polytope.
- Raises:
ValueError – When the point dimension does not match the polytope dimension.
- Returns:
Polytope which is the sum of self and Q.
- Return type:
Notes
Given a polytope \(\mathcal{P}\), and a set \(Q\), this function computes the Minkowski sum of Q and the polytope, 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 polytope \(\mathcal{R}=\{x + Q|x\in\mathcal{P}\}\).
Addition with a point: This function allows for \(\mathcal{P}\) to be in V-Rep or in H-Rep. For \(\mathcal{P}\) in V-rep, the polytope \(\mathcal{R}\) is the convex hull of the vertices of \(\mathcal{P}\) shifted by \(\text{point}\). Given \(\{v_i\}\) as the collection of vertices of \(\mathcal{P}\), \(\mathcal{R}=\mathrm{convexHull}(v_i + \text{point})\). For \(\mathcal{P}\) in H-rep, the polytope \(\mathcal{R}\) is defined by all points \(r = \text{point} + x\) with \(Ax\leq b, A_ex=b_e\). Thus, \(\mathcal{R}=\{x: A x \leq b + A \text{point}, A_e x = b_e + A_e \text{point}\}\).
Addition with a polytope: This function requires self and Q to be in V-Rep, and performs a vertex enumeration when self or Q are in H-Rep.. In vertex representation, \(\mathcal{R}\) is the convex hull of the pairwise sum of all combinations of points in \(\mathcal{P}\) and \(\mathcal{Q}\).
- 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:
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}\).
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
This function allows for \(\mathcal{P}\) to be in V-Rep or in H-Rep.
Given a polytope \(\mathcal{P}\) in V-Rep 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 \(\theta\in\mathbb{R}^{\mathcal{P}.\text{n\_vertices}}\),
\[\begin{split}\text{minimize} &\quad \|x - y\|_p\\ \text{subject to} &\quad x = \sum_i \theta_i v_i\\ &\quad \sum_i \theta_i = 1, \theta_i \geq 0\end{split}\]Given a polytope \(\mathcal{P}\) in H-Rep and a test point \(y\in\mathbb{R}^{\mathcal{P}.\text{dim}}\), this function solves a convex program with a decision variable \(x\in\mathbb{R}^{\mathcal{P}.\text{dim}}\),
\[\begin{split}\text{minimize} &\quad \|x - y\|_p\\ \text{subject to} &\quad A x \leq b\\ &\quad A_e x = b_e\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:
Notes
This function requires P to be in V-Rep, and performs a vertex enumeration when P is in H-Rep.
- 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:
Polytope that has been sliced at the specified dimensions.
- Return type:
Notes
This function requires \(\mathcal{P}\) to be in H-Rep, and performs a vertex halfspace when \(\mathcal{P}\) is in V-Rep.
- slice_then_projection(dims, constants)[source]¶
Wrapper for
slice()
andprojection()
.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:
Notes
This function requires \(\mathcal{P}\) to be in H-Rep, and performs a vertex halfspace when \(\mathcal{P}\) is in V-Rep.
- 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:
support_function_evaluations (numpy.ndarray): Support function evaluation(s) as a 2D numpy.ndarray. Vector (N,) with as many rows as eta.
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
This function allows for :math:mathcal{P} to be in V-Rep or in H-Rep.
Given a polytope \(\mathcal{P}\) in V-Rep 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 \(\theta\in\mathbb{R}^{\mathcal{P}.\text{n\_vertices}}\),
\[\begin{split}\text{maximize} &\quad \eta^\top x\\ \text{subject to} &\quad x = \sum_i \theta_i v_i\\ &\quad \sum_i \theta_i = 1, \theta_i \geq 0\end{split}\]Given a polytope \(\mathcal{P}\) in H-Rep and a support direction \(\eta\in\mathbb{R}^{\mathcal{P}.\text{dim}}\), this function solves a convex program with a decision variable \(x\in\mathbb{R}^{\mathcal{P}.\text{dim}}\),
\[\begin{split}\text{maximize} &\quad \eta^\top x\\ \text{subject to} &\quad A x \leq b\\ &\quad A_e x = b_e\end{split}\]
- volume()¶
Compute the volume of the polytope using QHull
- Returns:
Volume of the polytope
- Return type:
float
Notes
Works with V-representation: Yes
Works with H-representation: No
Performs a vertex-halfspace enumeration when H-rep is provided
Returns 0 when the polytope is empty
Requires polytope to be full-dimensional
- property A¶
Inequality coefficient vectors A for the polytope \(\{Ax \leq b, A_e x = b_e\}\).
- Returns:
Inequality coefficient vector (H-Rep). A is np.empty((0, self.dim)) for empty polytope.
- Return type:
numpy.ndarray
Notes
This function requires the polytope to be in H-Rep, and performs a halfspace enumeration if required.
- property Ae¶
Equality coefficient vectors Ae for the polytope \(\{Ax \leq b, A_e x = b_e\}\).
- Returns:
- Equality coefficient vector (H-Rep). Ae is np.empty((0, self.dim)) for empty or
full-dimensional polytope.
- Return type:
numpy.ndarray
Notes
This function requires the polytope to be in H-Rep, and performs a halfspace enumeration if required.
- property H¶
Inequality constraints in halfspace representation H=[A, b] for the polytope \(\{Ax \leq b, A_e x = b_e\}\).
- Returns:
H-Rep in [A, b]. H is np.empty((0, self.dim + 1)) for empty polytope.
- Return type:
numpy.ndarray
Notes
This function requires the polytope to be in H-Rep, and performs a halfspace enumeration if required.
- property He¶
Equality constraints in halfspace representation He=[Ae, be] for the polytope \(\{Ax \leq b, A_e x = b_e\}\).
- Returns:
H-Rep in [Ae, be]. He is np.empty((0, self.dim + 1)) for empty or full-dimensional polytope.
- Return type:
numpy.ndarray
Notes
This function requires the polytope to be in H-Rep, and performs a halfspace enumeration if required.
- property V¶
Vertex representation (V) where the polytope is given by \(\text{ConvexHull}(v_i)\) with \(v_i\) as the rows of \(V=[v_1;v_2;\ldots;v_{n_vertices}]\).
- Returns:
Vertices of the polytope, arranged row-wise. V is np.empty((0, self.dim)) if polytope is empty.
- Return type:
numpy.ndarray
Notes
This function requires the polytope to be in V-Rep, and performs a vertex enumeration if required.
- property b¶
Inequality constants b for the polytope \(\{Ax \leq b, A_e x = b_e\}\).
- Returns:
Inequality constants (H-Rep). b is np.empty((0,)) for empty polytope.
- Return type:
numpy.ndarray
Notes
This function requires the polytope to be in H-Rep, and performs a halfspace enumeration if required.
- property be¶
Equality constants be for the polytope \(\{Ax \leq b, A_e x = b_e\}\).
- Returns:
Equality constants (H-Rep). be is np.empty((0,)) for empty or full-dimensional polytope.
- Return type:
numpy.ndarray
Notes
This function requires the polytope to be in H-Rep, and performs a halfspace enumeration if required.
- 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_sdp¶
CVXPY arguments in use when solving a semi-definite program
- Returns:
CVXPY arguments in use when solving a semi-definite program. Defaults to dictionary in pycvxset.common.DEFAULT_CVXPY_ARGS_SDP.
- 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 polytope. In H-Rep polytope (A, b), this is the number of columns of A, while in V-Rep, this is the number of components of the vertices.
- Returns:
Dimension of the polytope
- Return type:
int
- property in_H_rep¶
Check if the polytope have a halfspace representation (H-Rep).
- Returns:
When True, the polytope has halfspace representation (H-Rep). Otherwise, False.
- Return type:
bool
- property in_V_rep¶
Check if the polytope have a vertex representation (V-Rep).
- Returns:
When True, the polytope has vertex representation (V-Rep). Otherwise, False.
- Return type:
bool
- property is_bounded¶
Check if the polytope is bounded.
- Returns:
True if the polytope is bounded, and False otherwise.
- Return type:
bool
Notes
This property is well-defined when the polytope is in V-Rep, but solves for a minimum volume circumscribing rectangle to check for boundedness.
- property is_empty¶
Check if the polytope is empty.
- Returns:
When True, the polytope is empty
- Return type:
bool
Notes
This property is well-defined when the polytope is in V-Rep, but may solve a Chebyshev centering problem when the polytope is in H-Rep. “may” because sometimes we can infer the emptiness of the polytope in H-Rep.
- property is_full_dimensional¶
Check if the affine dimension of the polytope is the same as the polytope dimension
- Returns:
True when the affine hull containing the polytope has the dimension self.dim
- Return type:
bool
Notes
This function can have self to be in V-Rep or H-Rep. See Sec. 2.1.3 of [BV04] for discussion on affine dimension.
An empty polytope is full dimensional if dim=0, otherwise it is not full-dimensional.
When the n-dimensional polytope is in V-rep, it is full-dimensional when its affine dimension is n. Recall that, the affine dimension is the dimension of the affine hull of the polytope is the linear subspace spanned by the vectors formed by subtracting the vertices with one of the vertices. Consequently, its dimension is given by the rank of the matrix P.V[1:] - P.V[0]. When there are fewer than self.dim + 1 vertices, we know it is coplanar without checking for matrix rank (simplex needs at least self.dim + 1 vertices and that is the polytope with the fewest vertices). For numerical stability, we zero-out all delta vertices below PYCVXSET_ZERO.
When the n-dimensional polytope is in H-Rep, it is full-dimensional if it can fit a n-dimensional ball of appropriate center and radius inside it (Chebyshev radius).
- property n_equalities¶
Number of linear equality constraints used to define the polytope \(\{Ax \leq b, A_e x = b_e\}\)
- Returns:
Number of linear equality constraints
- Return type:
int
Notes
A call to this property performs a halfspace enumeration if the polytope is in V-Rep.
- property n_halfspaces¶
Number of halfspaces used to define the polytope \(\{Ax \leq b, A_e x = b_e\}\)
- Returns:
Number of halfspaces
- Return type:
int
Notes
A call to this property performs a halfspace enumeration if the polytope is in V-Rep.
- property n_vertices¶
Number of vertices
- Returns:
Number of vertices
- Return type:
int
Notes
A call to this property performs a vertex enumeration if the polytope is in H-Rep.