Containment/Redundancy
Containment
Base.in — Functionin(p::VRepElement, h::HRepElement)Returns whether p is in h. If h is an hyperplane, it returns whether $\langle a, x \rangle \approx \beta$. If h is an halfspace, it returns whether $\langle a, x \rangle \le \beta$.
in(p::VRepElement, h::HRep)Returns whether p is in h, e.g. in all the hyperplanes and halfspaces supporting h.
Base.issubset — Functionissubset(p::Rep, h::HRepElement)Returns whether p is a subset of h, i.e. whether h supports the polyhedron p.
Polyhedra.ininterior — Functionininterior(p::VRepElement, h::HRepElement)Returns whether p is in the interior of h. If h is an hyperplane, it always returns false. If h is an halfspace $\langle a, x \rangle \leq \beta$, it returns whether p is in the open halfspace $\langle a, x \rangle < \beta$
ininterior(p::VRepElement, h::HRep)Returns whether p is in the interior of h, e.g. in the interior of all the hyperplanes and halfspaces supporting h.
Polyhedra.inrelativeinterior — Functioninrelativeinterior(p::VRepElement, h::HRepElement; tol)Returns whether p is in the relative interior of h. If h is an hyperplane, it is equivalent to p in h since the relative interior of an hyperplane is itself. If h is an halfspace, it is equivalent to ininterior(p, h).
inrelativeinterior(p::VRepElement, h::HRep; tol)Returns whether p is in the relative interior of h, e.g. in the relative interior of all the hyperplanes and halfspaces supporting h.
Polyhedra.support_function — Functionsupportfunction(h::AbstractVector, rep::Rep, solver=Polyhedra.linearobjective_solver(p))
Return the value of the support function of rep at h. See [Section 13, R15] for more details.
[R15] Rockafellar, R.T. Convex analysis. Princeton university press, 2015.
Linearity
Polyhedra.detecthlinearity! — Functiondetecthlinearity!(p::HRep; kws...)Detects all the hyperplanes contained in the H-representation and remove all redundant hyperplanes.
The remaining keyword arguments kws are passed to detecthlinearity.
Examples
The representation
h = HalfSpace([1, 1], 1]) ∩ HalfSpace([-1, -1], -1)contains the hyperplane HyperPlane([1, 1], 1).
Polyhedra.detecthlinearity — Functiondetecthlinearity(hr::HRepresentation, solver; kws...)Return a new H-representation with linearity detected using solver.
The remaining keyword arguments kws are passed to detect_new_linearities.
Polyhedra.detectvlinearity! — Functiondetectvlinearity!(p::VRep, solver=default_solver(p); kws...)Detects all the lines contained in the V-representation and remove all redundant lines.
The remaining keyword arguments kws are passed to detectvlinearity.
Examples
The representation
v = conichull([1, 1], [-1, -1])contains the line Line([1, 1]).
Polyhedra.detectvlinearity — Functiondetectvlinearity(vr::VRepresentation, solver; kws...)Return a new V-representation with linearity detected using solver.
The remaining keyword arguments kws are passed to detect_new_linearities.
Polyhedra.detect_new_linearities — Functiondetect_new_linearities(rep::HRepresentation{T}, solver; verbose=0, tol=Base.rtoldefault(T)) where {T}Given a polyhedron with H-representation rep, detect whether a new hyperplane can be generated from the halfspaces in halfspaces using an linear program solved by solver. The method is similar to the method used for lines described as follows. This function is automatically called by removehredundancy if a solver is provided.
detect_new_linearities(rep::VRepresentation{T}, solver; verbose=0, tol=Base.rtoldefault(T)) where {T}Given a cone defined by the V-representation rep (ignoring the points in the representation if any), detect whether a new line can be generated from the rays in rays using an linear program solved by solver. The method is as follows (suppose lines is empty for simplicity). This function is automatically called by removevredundancy if a solver is provided.
The keyword argument tol is used as a tolerance to decide whether a number is zero.
If there was a line l in the cone, it would mean that there exist μ >= 0 and ν >= 0 such that Σ μ_i r_i = l and Σ ν_i r_i = -l. We deduce from this that Σ λ_i r_i = 0 where λ = μ + ν.
Conversely, if there are λ >= 0 such that Σ λ_i r_i = 0 then let j be the index of λ with largest magnitude (to make sure it is nonzero). We have Σ_{i != j} λ_i/λ_j r_i = -r_j. As both r_j and -r_j are in the cone, r_j generates a line in the cone. However, this means that we now have Σ_{i != j} λ_i/λ_j r_i ≡ 0 (mod r_j) so if there is another λ_i with nonzero value, we can transform it to a line as well. In summary, we have a line r_i for each i such that λ_i != 0.
The dual program is:
max z
s.t. r_i'x ≥ zWhen the primal is feasible, the dual program may still be feasible. We know that z = 0 by strong duality as the objective value needs to be equal to the objective of the primal which is zero. So the constraints are r_i'x ≥ 0. If we have r_i'x > 0 for some i, it means that -r_i does not belong to the cone hence r_i can be dropped for the purpose of searching for lines.
Note
In CDDLib, the dual program is solved, if the objective value is zero then linearity are found by, for each i such that r_i'x = 0, solve an LP to find whether -r_i belongs to the cone. CDDLib ignores the primal results provided in λ which directly gives linearity without the need to solve an LP for each ray. The method implemented in Polyhedra is therefore significantly more efficient as its complexity is O(dimension of linespace) which is upper bounded by O(fulldim) while the method of CDDLib is O(number of rays).
Polyhedra.dim — Functiondim(h::HRep, current=false)Returns the dimension of the affine hull of the polyhedron. That is the number of non-redundant hyperplanes that define it. If current is true then it simply returns the dimension according the current number of hyperplanes, assuming that the H-linearity has already been detected. Otherwise, it first calls detecthlinearity!.
Duplicates
Polyhedra.removeduplicates — Functionremoveduplicates(rep::Representation)Removes the duplicates in the Representation.
- In an H-representation, it removes the redundant hyperplanes and it remove an halfspace when it is equal to another halfspace in the affine hull. For instance,
HalfSpace([1, 1], 1)is equal toHalfSpace([1, 0], 0)in the affine hull generated byHyperPlane([0, 1], 1]). - In a V-representation, it removes the redundant lines and it remove a point (resp. ray) when it is equal to another point (resp. ray) in the line hull. For instance, in the line hull generated by
Line([0, 1]),[1, 1]is equal to[1, 0]andRay([2, 2])is equal toRay([1, 0]).
Redundancy
Polyhedra.Redundancy — Type@enum Redundancy UNKNOWN_REDUNDANCY LINEARITY_DETECTED NO_REDUNDANCYRedundancy of a H-representation or V-representation.
UNKNOWN_REDUNDANCY: It is unknown whether there are any undetected linearity or redundancy.LINEARITY_DETECTED: There are no undetected linearity.NO_REDUNDANCY: There are no undetected linearity not redundancy.
An undetected linearity for a V-representation is a Line that is implicit in a conic hull of Rays. For instance, conichull([1, 1], [-1, -1]) has undetected linearity because it contains the Line Line([1, 1]). Some undetected linearity is less obvious, e.g., conichull([1, 0, -1], [0, 1, 1], [-1, -1, 0]) contains the Line Line([1, 1, 0]) as the sum of the first two Rays is Ray([1, 1, 0]).
An undetected linearity for a H-representation is a HyperPlane that is implicit in an intersection of HalfSpaces. For instance, HalfSpace([1, 1], 1) ∩ HalfSpace([-1, -1], 1) has undetected linearity because it contains the HyperPlane HyperPlane([1, 1], 1). Some undetected linearity is less obvious, e.g., HalfSpace([1, 0], -1) ∩ HalfSpace([0, 1], 1) ∩ HalfSpace([-1, -1], 0) contains the HyperPlane HyperPlane([1, 1], 0) as the sum of the first two HalfSpaces is HalfSpace([1, 1], 0).
Polyhedra.hredundancy — Functionhredundancy(p::Polyhedron)Return the Redundancy of hrep(p).
Polyhedra.vredundancy — Functionvredundancy(p::Polyhedron)Return the Redundancy of vrep(p).
Polyhedra.isredundant — Functionisredundant(p::Rep, idx::Index; strongly=false)Return a Bool indicating whether the element with index idx can be removed without changing the polyhedron represented by p. If strongly is true,
- if
idxis an H-representation elementh, it returnstrueonly if no V-representation element ofpis in the hyperplane ofh. - if
idxis a V-representation elementv, it returnstrueonly ifvis in the relative interior ofp.
Polyhedra.removehredundancy — Functionremovehredundancy(hr::HRepresentation, solver)Return a H-representation of the polyhedron represented by hr with all the elements of hr except the redundant ones, i.e. the elements that can be expressed as convex combination of other ones.
Polyhedra.removehredundancy! — Functionremovehredundancy!(p::HRep; kws...)Removes the elements of the H-representation of p that can be removed without changing the polyhedron represented by p. That is, it only keeps the halfspaces corresponding to facets of the polyhedron.
The remaining keyword arguments kws are passed to detectvlinearity! and detecthlinearity!.
Polyhedra.removevredundancy — Functionremovevredundancy(vr::VRepresentation, solver)Return a V-representation of the polyhedron represented by vr with all the elements of vr except the redundant ones, i.e. the elements that can be expressed as convex combination of other ones.
Polyhedra.removevredundancy! — Functionremovevredundancy!(p::VRep; strongly=false, planar=true, kws...)Removes the elements of the V-representation of p that can be removed without changing the polyhedron represented by p. That is, it only keeps the extreme points and rays. This operation is often called "convex hull" as the remaining points are the extreme points of the convex hull of the initial set of points. If strongly=true, weakly redundant points, i.e., points that are not extreme but are not in the relative interior either, may be kept. If fulldim(p) is 2, strongly is false and planar is true, a planar convex hull algorithm is used.
The remaining keyword arguments kws are passed to detectvlinearity! and detecthlinearity!.