Utilities

Utilities

Operations

Base.:*Function.
*(p1::Rep, p2::Rep)

Cartesian product between the polyhedra p1 and p2.

*(P::AbstractMatrix, p::VRep)

Transform the polyhedron represented by $p$ into $P p$ by transforming each element of the V-representation (points, symmetric points, rays and lines) x into $P x$.

Base.:\Function.
\(P::AbstractMatrix, p::HRep)

Transform the polyhedron represented by $p$ into $P^{-1} p$ by transforming each halfspace $\langle a, x \rangle \le \beta$ into $\langle P^\top a, x \rangle \le \beta$ and each hyperplane $\langle a, x \rangle = \beta$ into $\langle P^\top a, x \rangle = \beta$.

Base.:/Function.
/(p::HRep, P::AbstractMatrix)

Transform the polyhedron represented by $p$ into $P^{-T} p$ by transforming each halfspace $\langle a, x \rangle \le \beta$ into $\langle P a, x \rangle \le \beta$ and each hyperplane $\langle a, x \rangle = \beta$ into $\langle P a, x \rangle = \beta$.

Base.intersectFunction.
intersect(P1::HRep, P2::HRep)

Takes the intersection of P1 and P2 $\{\, x : x \in P_1, x \in P_2 \,\}$. It is very efficient between two H-representations or between two polyhedron for which the H-representation has already been computed. However, if P1 (resp. P2) is a polyhedron for which the H-representation has not been computed yet, it will trigger a representation conversion which is costly. See the Polyhedral Computation FAQ for a discussion on this operation.

The type of the result will be chosen closer to the type of P1. For instance, if P1 is a polyhedron (resp. H-representation) and P2 is a H-representation (resp. polyhedron), intersect(P1, P2) will be a polyhedron (resp. H-representation). If P1 and P2 are both polyhedra (resp. H-representation), the resulting polyhedron type (resp. H-representation type) will be computed according to the type of P1. The coefficient type however, will be promoted as required taking both the coefficient type of P1 and P2 into account.

intersect(v::VRepresentation{T}, h::HRepElement)

Compute the intersection of v with an halfspace or hyperplane h. The method used by default is to keep the V-representation element of v that are in h and add new ones generated as the intersection between the hyperplane defining h and the segment between two adjacent V-representation elements of v that are in either sides of the hyperplane. See Lemma 3 of [FP96] for more detail on the method.

[FP96] Fukuda, K. and Prodon, A. Double description method revisited Combinatorics and computer science, Springer, 1996, 91-111

Base.intersect!Function.
intersect!(p::HRep, h::Union{HRepresentation, HRepElement})

Same as intersect except that p is modified to be equal to the intersection.

Polyhedra.convexhullFunction.
convexhull(P1::VRep, P2::VRep)

Takes the convex hull of P1 and P2 $\{\, \lambda x + (1-\lambda) y : x \in P_1, y \in P_2 \,\}$. It is very efficient between two V-representations or between two polyhedron for which the V-representation has already been computed. However, if P1 (resp. P2) is a polyhedron for which the V-representation has not been computed yet, it will trigger a representation conversion which is costly.

The type of the result will be chosen closer to the type of P1. For instance, if P1 is a polyhedron (resp. V-representation) and P2 is a V-representation (resp. polyhedron), convexhull(P1, P2) will be a polyhedron (resp. V-representation). If P1 and P2 are both polyhedra (resp. V-representation), the resulting polyhedron type (resp. V-representation type) will be computed according to the type of P1. The coefficient type however, will be promoted as required taking both the coefficient type of P1 and P2 into account.

Polyhedra.convexhull!Function.
convexhull!(p1::VRep, p2::VRep)

Same as convexhull except that p1 is modified to be equal to the convex hull.

Volume

Polyhedra.volumeFunction.
volume(p::Polyhedron{T}) where {T}

Returns the fulldim(p)-dimensional hyper-volume of the polyhedron p. Returns Inf or -one(T) if it is infinite depending on whether the type T has an infinite value.

Polyhedra.surfaceFunction.
surface(p::Polyhedron{T}) where {T}

Returns the fulldim(p)-1-dimensional hyper-volume of the surface of the polyhedron p. Returns Inf or -one(T) if it is infinite depending on whether the type T has an infinite value.

Chebyshev center

chebyshevcenter(p::Rep[, solver])

If p is a H-representation or is a polyhedron for which the H-representation has already been computed, calls hchebyshevcenter, otherwise, call vchebyshevcenter.

hchebyshevcenter(p::HRep[, solver])

Return a tuple with the center and radius of the largest euclidean ball contained in the polyhedron p. Throws an error if the polyhedron is empty or if the radius is infinite.

vchebyshevcenter(p::VRep[, solver])

Return a tuple with the center and radius of the smallest euclidean ball containing the polyhedron p. Throws an error if the polyhedron is empty or if the radius is infinite (i.e. p is not a polytope, it contains rays).

Defining new representation

The following macros make it easy to define new representations:

The representation rep contain the elements elem inside a representation in the field field.

The representation rep does not contain any elem.

The representation rep contain the elements elem inside a vector in the field field.