pycvxset.Ellipsoid (API details)

pycvxset.Ellipsoid(**kwargs)

Ellipsoid class.

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

Bases: object

Ellipsoid class.

We can define a bounded, non-empty ellipsoid \(\mathcal{P}\) using one of the following combinations:

  1. \((c, Q)\) for a full-dimensional ellipsoid in the quadratic form \(\mathcal{P}=\{x \in \mathbb{R}^n\ |\ (x - c)^T Q^{-1} (x - c) \leq 1\}\) with a n-dimensional positive-definite matrix \(Q\) and a n-dimensional vector \(c\). Here, pycvxset computes a n-dimensional lower-triangular, square matrix \(G\) that satisfies \(GG^T=Q\).

  2. \((c, G)\) for a full-dimensional or a degenerate ellipsoid as an affine transformation of a unit-ball \(\mathcal{P} = \{x \in \mathbb{R}^n\ |\ \exists u\in\mathbb{R}^N,\ x = c + G u,\ {\|u\|}_2 \leq 1\}\) with a n x N matrix \(G\). Here, pycvxset computes \(Q=GG^T\).

Parameters:
  • c (array_like) – Center of the ellipsoid c. Vector of length (self.dim,)

  • Q (array_like, optional) – Shape matrix of the ellipsoid Q. Q must be a positive definite matrix (self.dim times self.dim).

  • G (array_like, optional) – Square root of the shape matrix of the ellipsoid G that satisfies \(GG^T=Q\). Need not be a square matrix, but must have self.dim rows.

  • r (scalar, optional) – Non-negative scalar that provides the radius of the self.dim-dimensional ball.

Raises:
  • ValueError – When more than one of Q, G, r was provided

  • ValueError – When c or Q or G or r does not satisfy implicit properties

Notes

  1. Empty ellipsoids are not permitted (c is a required keyword argument).

  2. When provided G is such that \(Q=GG^T\) is positive definite, we overwrite \(G\) with a lower-triangular, square, n-dimensional matrix for consistency. Here, \(G\) has strictly positive diagonal elements, and its determinant is the product of its diagonal elements (see volume computation).

  3. We use the eigenvalues of \(Q\) to determine the radii of the maximum volume inscribing ball (Chebyshev radius \(R^-\)) and the minimum volume circumscribing ball \(R^+\geq R^-\).

    1. The ellipsoid represents a singleton when \(R^+\) is negligible.

    2. The ellipsoid is full-dimensional when \(R^-\) is non-trivial.

__init__(**kwargs)[source]

Constructor for Ellipsoid

affine_hull()[source]

Compute the left null space of self.G to identify the affine hull. Affine hull is the entire self.dim-dimensional space when self is full-dimensional.

affine_map(M)

Compute the affine transformation of the given ellipsoid based on a given scalar/matrix.

Parameters:

M (int | float | array_like) – Scalar or matrix (m times self.dim) for the affine map.

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

  • ValueError – When M has columns not equal to self.dim

Returns:

Affine transformation of the given ellipsoid \(\mathcal{R}=M\mathcal{P}=\{Mx|x\in\mathcal{P}\}\)

Return type:

Ellipsoid

chebyshev_centering()[source]

Compute the Chebyshev center and radius of the ellipsoid.

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 an ellipsoid.

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.

Returns:

A tuple with two items:

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

  2. xi (cvxpy.Variable): CVXPY variable representing the latent dimension variable with length G.shape[1]

Return type:

tuple

contains(Q)

Check containment of a set or a collection of points in an ellipsoid.

Parameters:

Q (array_like | Polytope | Ellipsoid) – Polytope/Ellipsoid or a collection of points (each row is a point) to be tested for containment. When providing a collection of points, Q is a matrix (N times self.dim) with each row is a point.

Raises:
  • ValueError – Test point(s) are NOT of the same dimension

  • ValueError – Test point(s) can not be converted into a 2D numpy array of floats

  • ValueError – Q is a constrained zonotope

  • NotImplementedError – Unable to perform ellipsoidal containment check using CVXPY

Returns:

An element of the array is True if the point is in the ellipsoid, with as many

elements as the number of rows in test_points.

Return type:

bool or numpy.ndarray[bool]

Notes

  • Containment of a polytope: This function requires the polytope Q to be in V-Rep. If Q is in H-Rep, a vertex enumeration is performed. This function then checks if all vertices of Q are in the given ellipsoid, which occurs if and only if the polytope is contained within the ellipsoid [BV04].

  • Containment of an ellipsoid: This function solves a semi-definite program (S-procedure) [BV04].

  • Containment of points: For each point \(v\), the function checks if there is a \(u\in\mathbb{R}^{\mathcal{P}.\text{dim}}\) such that \(\|u\|_2 \leq 1\) and \(Gu + c=v\) [BV04]. This can be efficiently done via numpy.linalg.lstsq.

copy()[source]

Create a copy of the ellipsoid. Copy (c, G) to preserve degenerate ellipsoids.

classmethod deflate(set_to_be_centered)

Compute the minimum volume ellipsoid that covers the given set (also known as Lowner-John Ellipsoid).

Parameters:

set_to_be_centered (Polytope | ConstrainedZonotope) – Set to be circumscribed.

Returns:

Minimum volume circumscribing ellipsoid

Return type:

Ellipsoid

Notes

This function is a wrapper for minimum_volume_circumscribing_ellipsoid() of the set set_to_be_centered. Please check that function for more details including raising exceptions. [EllipsoidalTbx-Min_verticesolEll]

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.

classmethod inflate(set_to_be_centered)

Compute the maximum volume ellipsoid that fits within the given set.

Parameters:

set_to_be_centered (Polytope | ConstrainedZonotope) – Set to be inscribed.

Returns:

Maximum volume inscribing ellipsoid

Return type:

Ellipsoid

Notes

This function is a wrapper for maximum_volume_inscribing_ellipsoid() of the set set_to_expand_within. Please check that function for more details including raising exceptions. [EllipsoidalTbx-MinVolEll]

classmethod inflate_ball(set_to_be_centered)

Compute the largest ball (Chebyshev ball) of a given set.

Parameters:

set_to_be_centered (Polytope | ConstrainedZonotope) – Set to compute Chebyshev ball for.

Returns:

Maximum volume inscribing ellipsoid

Return type:

Ellipsoid

Notes

This function is a wrapper for chebyshev_center() of attr:set_to_be_centered. Please check that function for more details including raising exceptions.

interior_point()[source]

Compute an interior point to the Ellipsoid

intersection_with_affine_set(Ae, be)

Compute the intersection of an ellipsoid 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 an ellipsoid with the affine set.

Return type:

Ellipsoid

Notes

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

inverse_affine_map_under_invertible_matrix(M)

Compute the inverse affine transformation of an ellipsoid based on a given scalar/matrix.

Parameters:

M (int | float | array_like) – Scalar or invertible square matrix for the affine map

Raises:
  • TypeError – When M is not convertible into a 2D square numpy matrix

  • TypeError – When M is not invertible

Returns:

Inverse affine transformation of the given ellipsoid \(\mathcal{R}=\mathcal{P}M=\{x|Mx\in\mathcal{P}\}\)

Return type:

Ellipsoid

Notes

  1. Since M is invertible, \(\mathcal{R}\) is also a bounded ellipsoid.

maximum_volume_inscribing_ellipsoid()[source]

Compute the maximum volume inscribing ellipsoid for a given 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:
  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_ball()[source]

Compute the radius of the ball that circumscribes the ellipsoid and has the minimum volume.

minimum_volume_circumscribing_ellipsoid()[source]

Compute the minimum volume circumscribing ellipsoid for a given ellipsoid.

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}\).

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(point)

Add a point to an ellipsoid

Parameters:

point (array_like) – Vector (self.dim,) that describes the point to be added.

Raises:
  • ValueError – point is a set (ConstrainedZonotope or Ellipsoid or Polytope)

  • TypeError – point can not be converted into a numpy array of float

  • ValueError – point can not be converted into a 1D numpy array of float

  • ValueError – Mismatch in dimension

Returns:

Sum of the ellipsoid and the point.

Return type:

Ellipsoid

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

For a point \(y\in\mathbb{R}^{\mathcal{P}.\text{dim}}\) and an ellipsoid \(\mathcal{P}=\{G u + c\ |\ \| u \|_2 \leq 1 \}\) with \(GG^T=Q\), this function solves a convex program with decision variables \(x,u\in\mathbb{R}^{\mathcal{P}.\text{dim}}\),

\[\begin{split}\text{minimize} &\quad \|x - y\|_p\\ \text{subject to} &\quad x = G u + c\\ &\quad {\| u \|}_2 \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:

Ellipsoid

quadratic_form_as_a_symmetric_matrix()[source]

Define a (self.dim + 1)-dimensional symmetric matrix M where self = {x | [x, 1] @ M @ [x, 1] <= 0}. Here, when Q is not positive definite, we use pseudo-inverse of Q.

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:

Ellipsoid that has been sliced at the specified dimensions.

Return type:

Ellipsoid

Notes

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:

Ellipsoid

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

Using duality, the support function and vector of an ellipsoid has a closed-form expressions. For a support direction \(\eta\in\mathbb{R}^{\mathcal{P}.\text{dim}}\) and an ellipsoid \(\mathcal{P}=\{G u + c | \| u \|_2 \leq 1 \}\) with \(GG^T=Q\),

\[\begin{split}\rho_{\mathcal{P}}(\eta) &= \eta^\top c + \sqrt{\eta^\top Q \eta} = \eta^\top c + \|G^T \eta\|_2\\ \nu_{\mathcal{P}}(\eta) &= c + \frac{G G^\top \eta}{\|G^T \eta\|_2} = c + \frac{Q \eta}{\|G^T \eta\|_2}\end{split}\]

For degenerate (not full-dimensional) ellipsoids and \(\eta\) not in the low-dimensional affine hull containing the ellipsoid,

\[\begin{split}\rho_{\mathcal{P}}(\eta) &= \eta^\top c \\ \nu_{\mathcal{P}}(\eta) &= c\end{split}\]
volume()[source]

Compute the volume of the ellipsoid.

Notes

Volume of the ellipsoid is zero if it is not full-dimensional. For full-dimensional ellipsoid, we used the following observations:

  1. from [BV04] , the volume of an ellipsoid is proportional to \(det(G)\).

  2. Square-root of the determinant of the shape matrix coincides with the determinant of G

  3. Since G is lower-triangular, its determinant is the product of its diagonal elements.

property G

Affine transformation matrix \(G\) that satisfies \(GG^T=Q\)

property Q

Shape matrix of the ellipsoid \(Q\)

property c

Center of the ellipsoid \(c\)

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 ellipsoid \(dim\)

property is_empty

Check if the ellipsoid is empty. Always False by construction.

property is_full_dimensional

Check if the ellipsoid is full-dimensional, i.e., sqrt of all eigenvalues of Q are above PYCVXSET_ZERO

property is_singleton

Check if the ellipsoid is a singleton, i.e., sqrt of all eigenvalues of Q are below PYCVXSET_ZERO

property latent_dim

Latent dimension of the ellipsoid \(dim\)