-
Notifications
You must be signed in to change notification settings - Fork 27
/
index.html
36 lines (36 loc) · 39.1 KB
/
index.html
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
<!DOCTYPE html>
<html lang="en"><head><meta charset="UTF-8"/><meta name="viewport" content="width=device-width, initial-scale=1.0"/><title>Representation · Polyhedra</title><link href="https://fonts.googleapis.com/css?family=Lato|Roboto+Mono" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/font-awesome/5.11.2/css/fontawesome.min.css" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/font-awesome/5.11.2/css/solid.min.css" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/font-awesome/5.11.2/css/brands.min.css" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/KaTeX/0.11.1/katex.min.css" rel="stylesheet" type="text/css"/><script>documenterBaseURL=".."</script><script src="https://cdnjs.cloudflare.com/ajax/libs/require.js/2.3.6/require.min.js" data-main="../assets/documenter.js"></script><script src="../siteinfo.js"></script><script src="../../versions.js"></script><link class="docs-theme-link" rel="stylesheet" type="text/css" href="../assets/themes/documenter-dark.css" data-theme-name="documenter-dark"/><link class="docs-theme-link" rel="stylesheet" type="text/css" href="../assets/themes/documenter-light.css" data-theme-name="documenter-light" data-theme-primary/><script src="../assets/themeswap.js"></script></head><body><div id="documenter"><nav class="docs-sidebar"><div class="docs-package-name"><span class="docs-autofit">Polyhedra</span></div><form class="docs-search" action="../search/"><input class="docs-search-query" id="documenter-search-query" name="q" type="text" placeholder="Search docs"/></form><ul class="docs-menu"><li><a class="tocitem" href="../">Index</a></li><li><a class="tocitem" href="../installation/">Installation</a></li><li class="is-active"><a class="tocitem" href>Representation</a><ul class="internal"><li><a class="tocitem" href="#H-representation"><span>H-representation</span></a></li><li><a class="tocitem" href="#V-representation"><span>V-representation</span></a></li></ul></li><li><a class="tocitem" href="../polyhedron/">Polyhedron</a></li><li><a class="tocitem" href="../plot/">Plot</a></li><li><a class="tocitem" href="../redundancy/">Containment/Redundancy</a></li><li><a class="tocitem" href="../projection/">Projection/Elimination</a></li><li><a class="tocitem" href="../optimization/">Optimization</a></li><li><a class="tocitem" href="../utilities/">Utilities</a></li><li><span class="tocitem">Examples</span><ul><li><a class="tocitem" href="../generated/Convex hull and intersection/">Convex hull and intersection</a></li><li><a class="tocitem" href="../generated/Extended Formulation/">Extended Formulation</a></li><li><a class="tocitem" href="../generated/Minimal Robust Positively Invariant Set/">Minimal Robust Positively Invariant Set</a></li></ul></li></ul><div class="docs-version-selector field has-addons"><div class="control"><span class="docs-label button is-static is-size-7">Version</span></div><div class="docs-selector control is-expanded"><div class="select is-fullwidth is-size-7"><select id="documenter-version-selector"></select></div></div></div></nav><div class="docs-main"><header class="docs-navbar"><nav class="breadcrumb"><ul class="is-hidden-mobile"><li class="is-active"><a href>Representation</a></li></ul><ul class="is-hidden-tablet"><li class="is-active"><a href>Representation</a></li></ul></nav><div class="docs-right"><a class="docs-edit-link" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/master/docs/src/representation.md" title="Edit on GitHub"><span class="docs-icon fab"></span><span class="docs-label is-hidden-touch">Edit on GitHub</span></a><a class="docs-settings-button fas fa-cog" id="documenter-settings-button" href="#" title="Settings"></a><a class="docs-sidebar-button fa fa-bars is-hidden-desktop" id="documenter-sidebar-button" href="#"></a></div></header><article class="content" id="documenter-page"><h1 id="Representation"><a class="docs-heading-anchor" href="#Representation">Representation</a><a id="Representation-1"></a><a class="docs-heading-anchor-permalink" href="#Representation" title="Permalink"></a></h1><p>Polyhedra can be described in 2 different ways.</p><ol><li>H-representation: As the intersection of finitely many halfspaces given by its facets.</li><li>V-representation: As the convex hull of its vertices + the conic hull of its rays where '+' is the Minkowski sum.</li></ol><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.HRepresentation" href="#Polyhedra.HRepresentation"><code>Polyhedra.HRepresentation</code></a> — <span class="docstring-category">Type</span></header><section><div><pre><code class="language-julia">HRepresentation{T<:Real}</code></pre><p>Supertype for H-representations with coefficient type <code>T</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/representation.jl#L11-L15">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.VRepresentation" href="#Polyhedra.VRepresentation"><code>Polyhedra.VRepresentation</code></a> — <span class="docstring-category">Type</span></header><section><div><pre><code class="language-julia">VRepresentation{T<:Real}</code></pre><p>Supertype for V-representations coefficient type <code>T</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/representation.jl#L17-L21">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.Representation" href="#Polyhedra.Representation"><code>Polyhedra.Representation</code></a> — <span class="docstring-category">Type</span></header><section><div><pre><code class="language-julia">Representation{T<:Real}</code></pre><p>Supertype for H-(or V-)representations with coefficient type <code>T</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/representation.jl#L5-L9">source</a></section></article><p>In <code>Polyhedra.jl</code>, those representations are given the respective abstract types <code>HRepresentation</code> and <code>VRepresentation</code> which are themself subtypes of <code>Representation</code>.</p><p>These functions can be called on both H-representation and V-representation</p><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.fulldim" href="#Polyhedra.fulldim"><code>Polyhedra.fulldim</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">fulldim(rep::Rep)::Int</code></pre><p>Returns the dimension of the space in which polyhedron, representation, element or vector is defined. That is, a straight line in a 3D space has <code>fulldim</code> 3 even if its dimension is 1.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/dimension.jl#L40-L46">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.FullDim" href="#Polyhedra.FullDim"><code>Polyhedra.FullDim</code></a> — <span class="docstring-category">Constant</span></header><section><div><pre><code class="language-julia">FullDim(p)::FullDim</code></pre><p>Similar to <a href="#Polyhedra.fulldim"><code>fulldim</code></a> but used for type stability with the vector type. If the vector type is <code>StaticArrays.SVector</code> then it returns a <code>StaticArrays.Size</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/dimension.jl#L14-L20">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.coefficient_type" href="#Polyhedra.coefficient_type"><code>Polyhedra.coefficient_type</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">coefficient_type(rep::Rep)</code></pre><p>Returns the type of the coefficients used in the representation of <code>rep</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/representation.jl#L32-L36">source</a></section></article><h2 id="H-representation"><a class="docs-heading-anchor" href="#H-representation">H-representation</a><a id="H-representation-1"></a><a class="docs-heading-anchor-permalink" href="#H-representation" title="Permalink"></a></h2><p>The fundamental element of an H-representation is the halfspace</p><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.HalfSpace" href="#Polyhedra.HalfSpace"><code>Polyhedra.HalfSpace</code></a> — <span class="docstring-category">Type</span></header><section><div><pre><code class="language-julia">struct HalfSpace{T, AT} <: HRepElement{T, AT}
a::AT
β::T
end</code></pre><p>An halfspace defined by the set of points <span>$x$</span> such that <span>$\langle a, x \rangle \le \beta$</span>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/elements.jl#L14-L21">source</a></section></article><p>An H-representation can be created as the intersection of several halfspaces. For instance, the polytope</p><div>\[\begin{align*}
x_1 + x_2 &\leq 1 \\
x_1 - x_2 &\leq 0 \\
x_1 & \geq 0.
\end{align*}\]</div><p>can be created as follows:</p><pre><code class="language-julia">HalfSpace([1, 1], 1) ∩ HalfSpace([1, -1], 0) ∩ HalfSpace([-1, 0], 0)</code></pre><p>Even if <code>HalfSpace</code>s are enough to describe any polyhedron, it is sometimes important to represent the fact that the polyhedron is contained in an affine subspace. For instance, the simplex</p><div>\[\begin{align*}
x_1 + x_2 &= 1 \\
x_1 &\geq 0 \\
x_2 &\geq 0
\end{align*}\]</div><p>is 1-dimensional even if it is defined in a 2-dimensional space.</p><p>The fundamental element of an affine subspace is the hyperplane</p><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.HyperPlane" href="#Polyhedra.HyperPlane"><code>Polyhedra.HyperPlane</code></a> — <span class="docstring-category">Type</span></header><section><div><pre><code class="language-julia">struct HyperPlane{T, AT} <: HRepElement{T, AT}
a::AT
β::T
end</code></pre><p>An hyperplane defined by the set of points <span>$x$</span> such that <span>$\langle a, x \rangle = \beta$</span>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/elements.jl#L33-L40">source</a></section></article><p>An affine subspace can be created as the intersection of several hyperplanes. For instance</p><pre><code class="language-julia">HyperPlane([1, 1], 1) ∩ HyperPlane([1, 0], 0)</code></pre><p>represents the 0-dimensional affine subspace only containing the point <span>$(0, 1)$</span>.</p><p>To represent a polyhedron that is not full-dimensional, hyperplanes and halfspaces can be mixed in any order. For instance, the simplex defined above can be obtained as follows:</p><pre><code class="language-julia">HalfSpace([-1, 0], 0) ∩ HyperPlane([1, 1], 1) ∩ HalfSpace([0, -1], 0)</code></pre><p>In addition to being created incrementally with intersections, an H-representation can also be created using the <code>hrep</code> function</p><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.hrep" href="#Polyhedra.hrep"><code>Polyhedra.hrep</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">hrep(p::Polyhedron)</code></pre><p>Returns an H-representation for the polyhedron <code>p</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/polyhedron.jl#L20-L24">source</a></section><section><div><pre><code class="language-none">hrep(hyperplanes::HyperPlaneIt; d::FullDim)</code></pre><p>Creates an affine space of full dimension <code>d</code> from the list of hyperplanes <code>hyperplanes</code>.</p><p><strong>Examples</strong></p><pre><code class="language-julia">hrep([HyperPlane([0, 1, 0], 1), HyperPlane([0, 0, 1], 0)])</code></pre><p>creates the 1-dimensional affine subspace containing all the points <span>$(x_1, 0, 0)$</span>, i.e. the <span>$x_1$</span>-axis.</p><pre><code class="language-julia">hrep([HyperPlane([1, 1], 1), HyperPlane([1, 0], 0)])</code></pre><p>creates the 0-dimensional affine subspace only containing the point <span>$(0, 1)$</span>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/aff.jl#L20-L36">source</a></section><section><div><pre><code class="language-none">hrep(hyperplanes::HyperPlaneIt, halfspaces::HalfSpaceIt; d::FullDim)</code></pre><p>Creates an H-representation for the polyhedron of full dimension <code>d</code> equal to the intersection of the hyperplanes <code>hyperplanes</code> and halfspaces <code>halfspaces</code>.</p><p><strong>Examples</strong></p><p>For instance, the simplex</p><div>\[\begin{align*}
x_1 + x_2 &= 1 \\
x_1 &\geq 0 \\
x_2 &\geq 0
\end{align*}\]</div><p>can be created as follows:</p><pre><code class="language-julia-repl">julia> hrep([HyperPlane([1, 1], 1)], [HalfSpace([0, -1], 0), HalfSpace([-1, 0], 0)])
H-representation Polyhedra.Intersection{Int64,Array{Int64,1},Int64}:
1-element iterator of HyperPlane{Int64,Array{Int64,1}}:
HyperPlane([1, 1], 1),
2-element iterator of HalfSpace{Int64,Array{Int64,1}}:
HalfSpace([0, -1], 0)
HalfSpace([-1, 0], 0)</code></pre></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/vecrep.jl#L8-L33">source</a></section><section><div><pre><code class="language-none">hrep(halfspaces::HalfSpaceIt; d::FullDim)</code></pre><p>Creates an H-representation for the polyhedron of full dimension <code>d</code> equal to the intersection of the halfspaces <code>halfspaces</code>.</p><p><strong>Examples</strong></p><p>For instance, the polytope</p><div>\[\begin{align*}
x_1 + x_2 &\leq 1 \\
x_1 - x_2 &\leq 0 \\
x_1 & \geq 0.
\end{align*}\]</div><p>can be created as follows:</p><pre><code class="language-julia">hrep([HalfSpace([1, 1], 1), HalfSpace([1, -1], 0), HalfSpace([-1, 0], 0)])</code></pre></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/vecrep.jl#L39-L58">source</a></section><section><div><pre><code class="language-none">hrep(model::JuMP.Model)</code></pre><p>Builds an H-representation from the feasibility set of the JuMP model <code>model</code>. Note that if non-linear constraint are present in the model, they are ignored.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/jump.jl#L3-L8">source</a></section><section><div><pre><code class="language-none">hrep(A::AbstractMatrix, b::AbstractVector, linset::BitSet=BitSet())</code></pre><p>Creates an H-representation for the polyhedron defined by the inequalities <span>$\langle A_i, x \rangle = b_i$</span> if <code>i in linset</code> and <span>$\langle A_i, x \rangle \le b_i$</span> otherwise where <span>$A_i$</span> is the <span>$i$</span>th row of <code>A</code>, i.e. <code>A[i,:]</code> and <span>$b_i$</span> is <code>b[i]</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/matrep.jl#L4-L9">source</a></section></article><h3 id="Interface"><a class="docs-heading-anchor" href="#Interface">Interface</a><a id="Interface-1"></a><a class="docs-heading-anchor-permalink" href="#Interface" title="Permalink"></a></h3><p>An H-representation is represented as an intersection halfspaces and hyperplanes. The halfspaces can be obtained with <a href="#Polyhedra.halfspaces"><code>halfspaces</code></a> and the hyperplanes with <a href="#Polyhedra.hyperplanes"><code>hyperplanes</code></a>. As an hyperplane <span>$\langle a, x \rangle = \beta$</span> is the intersection of the two halfspaces <span>$\langle a, x \rangle \le \beta$</span> and <span>$\langle a, x \rangle \ge \beta$</span>, even if the H-representation contains hyperplanes, a list of halfspaces whose intersection is the polyhedron can be obtained with <a href="#Polyhedra.allhalfspaces"><code>allhalfspaces</code></a>, which has <code>nhalfspaces(p) + 2nhyperplanes(p)</code> elements for an H-representation <code>p</code> since each hyperplane is split in two halfspaces.</p><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.halfspaces" href="#Polyhedra.halfspaces"><code>Polyhedra.halfspaces</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">halfspaces(hrep::HRep)</code></pre><p>Returns an iterator over the halfspaces of the H-representation <code>hrep</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L140-L144">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.nhalfspaces" href="#Polyhedra.nhalfspaces"><code>Polyhedra.nhalfspaces</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">nhalfspaces(hrep::HRep)</code></pre><p>Returns the number of halfspaces of the H-representation <code>hrep</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L161-L165">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.hashalfspaces" href="#Polyhedra.hashalfspaces"><code>Polyhedra.hashalfspaces</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">hashalfspaces(hrep::HRep)</code></pre><p>Returns whether the H-representation <code>hrep</code> has any halfspace.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L168-L172">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.hyperplanes" href="#Polyhedra.hyperplanes"><code>Polyhedra.hyperplanes</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">hyperplanes(hrep::HRep)</code></pre><p>Returns an iterator over the hyperplanes of the H-representation <code>hrep</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L140-L144">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.nhyperplanes" href="#Polyhedra.nhyperplanes"><code>Polyhedra.nhyperplanes</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">nhyperplanes(hrep::HRep)</code></pre><p>Returns the number of hyperplanes of the H-representation <code>hrep</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L161-L165">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.hashyperplanes" href="#Polyhedra.hashyperplanes"><code>Polyhedra.hashyperplanes</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">hashyperplanes(hrep::HRep)</code></pre><p>Returns whether the H-representation <code>hrep</code> has any hyperplane.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L168-L172">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.allhalfspaces" href="#Polyhedra.allhalfspaces"><code>Polyhedra.allhalfspaces</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">allhalfspaces(hrep::HRep)</code></pre><p>Returns an iterator over the halfspaces and hyperplanes in the H-representation <code>hrep</code> splitting hyperplanes in two halfspaces.</p><p><strong>Examples</strong></p><pre><code class="language-julia">hrep = HyperPlane([1, 0], 1) ∩ HalfSpace([0, 1], 1)
collect(allhalfspaces(hrep)) # Returns [HalfSpace([1, 0]), HalfSpace([-1, 0]), HalfSpace([0, 1])]</code></pre></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L287-L298">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.nallhalfspaces" href="#Polyhedra.nallhalfspaces"><code>Polyhedra.nallhalfspaces</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">nallhalfspaces(hrep::HRep)</code></pre><p>Returns the number of halfspaces plus twice the number of hyperplanes in the H-representation <code>hrep</code>, i.e. <code>length(allhalfspaces(hrep))</code></p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L301-L305">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.hasallhalfspaces" href="#Polyhedra.hasallhalfspaces"><code>Polyhedra.hasallhalfspaces</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">hasallhalfspaces(hrep::HRep)</code></pre><p>Returns whether the H-representation <code>hrep</code> contains any halfspace or hyperplane.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L308-L312">source</a></section></article><h2 id="V-representation"><a class="docs-heading-anchor" href="#V-representation">V-representation</a><a id="V-representation-1"></a><a class="docs-heading-anchor-permalink" href="#V-representation" title="Permalink"></a></h2><p>The fundamental elements of an V-representation are the points (represented by <code>AbstractVector</code>s and and rays</p><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.Ray" href="#Polyhedra.Ray"><code>Polyhedra.Ray</code></a> — <span class="docstring-category">Type</span></header><section><div><pre><code class="language-julia">struct Ray{T, AT <: AbstractVector{T}}
a::AT
end</code></pre><p>The conic hull of <code>a</code>, i.e. the set of points <code>λa</code> where <code>λ</code> is any nonnegative real number.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/elements.jl#L109-L115">source</a></section></article><p>A V-representation can be created as the minkowski sum between a convex hull of points and a conic hull of rays. For instance, the positive orthant without the simplex defined in the H-representation section can be created as follows:</p><pre><code class="language-julia">convexhull([1, 0], [0, 1]) + conichull([1, 0], [0, 1])</code></pre><p>The V-representation represents the polyhedron as a minkowski sum of a polytope and a polyhedral cone. The polytope is represented using a <em>P-representation</em> : a convex hull of points. The polyhedral cone is represented using an <em>R-representation</em> : a conic hull of rays.</p><p>Even if rays are enough to describe any polyhedral cone, it is sometimes important to represent the fact that the polyhedron contains an affine subspace. For instance, the polyhedron created with</p><pre><code class="language-julia">convexhull([1, 0], [0, 1]) + conichull([1, 1], [-1, -1])</code></pre><p>contains the line <code>[1, 1]</code>.</p><p>The fundamental element of an affine subspace is the line</p><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.Line" href="#Polyhedra.Line"><code>Polyhedra.Line</code></a> — <span class="docstring-category">Type</span></header><section><div><pre><code class="language-julia">struct Line{T, AT <: AbstractVector{T}}
a::AT
end</code></pre><p>The conic hull of <code>a</code> and <code>-a</code>, i.e. the set of points <code>λa</code> where <code>λ</code> is any real number.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/elements.jl#L127-L133">source</a></section></article><p>An affine subspace can be created as the conic hull/minkownski sum of several lines. For instance</p><pre><code class="language-julia">conichull(Line([1, 0]), Line([0, 1]))</code></pre><p>represents the full space.</p><p>In addition to being created incrementally with convex hull and minkowsky addition, a V-representation can also be created using the <code>vrep</code> function</p><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.vrep" href="#Polyhedra.vrep"><code>Polyhedra.vrep</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">vrep(p::Polyhedron)</code></pre><p>Returns a V-representation for the polyhedron <code>p</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/polyhedron.jl#L34-L38">source</a></section><section><div><pre><code class="language-none">vrep(lines::LineIt; d::FullDim)</code></pre><p>Creates an affine space of full dimension <code>d</code> from the list of lines <code>lines</code>.</p><p><strong>Examples</strong></p><pre><code class="language-julia">vrep([Line([1, 0, 0]), Line([0, 1, 0])])</code></pre><p>creates the 2-dimensional affine subspace containing all the points <span>$(x_1, x_2, 0)$</span>, i.e. the <span>$x_1$</span><span>$x_2$</span>-plane.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/aff.jl#L138-L148">source</a></section><section><div><pre><code class="language-none">vrep(points::PointIt; d::FullDim)</code></pre><p>Creates a V-representation for the polytope of full dimension <code>d</code> equal to the convex hull of the points <code>points</code>.</p><p><strong>Examples</strong></p><p>The convex hull of <span>$(0, 0)$</span>, <span>$(0, 1)$</span> and <span>$(1/2, 1/2)$</span> can be created as follows using exact arithmetic</p><pre><code class="language-julia">vrep([[0, 0], [0, 1], [1//2, 1//2]])</code></pre><p>or as follows using floating point arithmetic</p><pre><code class="language-julia">vrep([[0, 0], [0, 1], [1/2, 1/2]])</code></pre></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/vecrep.jl#L144-L159">source</a></section><section><div><pre><code class="language-none">vrep(lines::LineIt, rays::RayIt; d::FullDim)</code></pre><p>Creates a V-representation for the polyhedral cone of full dimension <code>d</code> equal to the conic hull of the lines <code>lines</code> and rays <code>rays</code>.</p><p><strong>Examples</strong></p><pre><code class="language-julia">vrep([Line([0, 1])], [Ray([1, 0])])</code></pre><p>creates a V-representation for the halfspace <span>$x_1 \ge 0$</span>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/vecrep.jl#L181-L192">source</a></section><section><div><pre><code class="language-none">vrep(rays::RayIt)</code></pre><p>Creates a V-representation for the polyhedral cone of full dimension <code>d</code> equal to the conic hull of the rays <code>rays</code>.</p><p><strong>Examples</strong></p><pre><code class="language-julia">vrep([Ray([1, 0]), Ray([0, 1])])</code></pre><p>creates a V-representation for positive orthant.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/vecrep.jl#L197-L208">source</a></section><section><div><pre><code class="language-none">vrep(points::PointIt, lines::LineIt, rays::RayIt)</code></pre><p>Creates a V-representation for the polyhedron of full dimension <code>d</code> equal to the minkowski sum of the convex hull of <code>points</code> with the conic hull of <code>lines</code> and <code>rays</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/vecrep.jl#L233-L239">source</a></section><section><div><pre><code class="language-none">vrep(V::AbstractMatrix, R::AbstractMatrix, Rlinset::BitSet=BitSet())</code></pre><p>Creates a V-representation for the polyhedron defined by the points <span>$V_i$</span>, lines <span>$R_i$</span> if <code>i in Rlinset</code> and rays <span>$R_i$</span> otherwise where <span>$V_i$</span> (resp. <span>$R_i$</span>) is the <span>$i$</span>th row of <code>V</code> (resp. <code>R</code>), i.e. <code>V[i,:]</code> (resp. <code>R[i,:]</code>).</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/matrep.jl#L72-L78">source</a></section></article><h3 id="Interface-2"><a class="docs-heading-anchor" href="#Interface-2">Interface</a><a class="docs-heading-anchor-permalink" href="#Interface-2" title="Permalink"></a></h3><p>A P-representation is represented as a convex hull points. The points can be obtained with <a href="#Polyhedra.points"><code>points</code></a>.</p><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.points" href="#Polyhedra.points"><code>Polyhedra.points</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">points(vrep::VRep)</code></pre><p>Returns an iterator over the points of the V-representation <code>vrep</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L140-L144">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.npoints" href="#Polyhedra.npoints"><code>Polyhedra.npoints</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">npoints(vrep::VRep)</code></pre><p>Returns the number of points of the V-representation <code>vrep</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L161-L165">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.haspoints" href="#Polyhedra.haspoints"><code>Polyhedra.haspoints</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">haspoints(vrep::VRep)</code></pre><p>Returns whether the V-representation <code>vrep</code> has any point.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L168-L172">source</a></section></article><p>An R-representation is represented as a conic hull of lines and rays. The rays can be obtained with <a href="#Polyhedra.rays"><code>rays</code></a> and the lines with <a href="#Polyhedra.lines"><code>lines</code></a>. As a line <span>$r$</span> is the conic hull of of the two rays <span>$r$</span> and <span>$-r$</span>, even if the R-representation contains lines, a list of rays whose conic hull is the polyhedral cone can be obtained with <a href="#Polyhedra.allrays"><code>allrays</code></a>, which has <code>nrays(R) + 2nlines(R)</code> elements for an R-representation <code>R</code> since each line is split in two rays.</p><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.rays" href="#Polyhedra.rays"><code>Polyhedra.rays</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">rays(vrep::VRep)</code></pre><p>Returns an iterator over the rays of the V-representation <code>vrep</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L140-L144">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.nrays" href="#Polyhedra.nrays"><code>Polyhedra.nrays</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">nrays(vrep::VRep)</code></pre><p>Returns the number of rays of the V-representation <code>vrep</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L161-L165">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.hasrays" href="#Polyhedra.hasrays"><code>Polyhedra.hasrays</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">hasrays(vrep::VRep)</code></pre><p>Returns whether the V-representation <code>vrep</code> has any ray.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L168-L172">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.lines" href="#Polyhedra.lines"><code>Polyhedra.lines</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">lines(vrep::VRep)</code></pre><p>Returns an iterator over the lines of the V-representation <code>vrep</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L140-L144">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.nlines" href="#Polyhedra.nlines"><code>Polyhedra.nlines</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">nlines(vrep::VRep)</code></pre><p>Returns the number of lines of the V-representation <code>vrep</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L161-L165">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.haslines" href="#Polyhedra.haslines"><code>Polyhedra.haslines</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">haslines(vrep::VRep)</code></pre><p>Returns whether the V-representation <code>vrep</code> has any line.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L168-L172">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.allrays" href="#Polyhedra.allrays"><code>Polyhedra.allrays</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">allrays(vrep::VRep)</code></pre><p>Returns an iterator over the rays and lines in the V-representation <code>vrep</code> splitting lines in two rays.</p><p><strong>Examples</strong></p><pre><code class="language-julia">vrep = Line([1, 0]) + Ray([0, 1])
collect(allrays(vrep)) # Returns [Ray([1, 0]), Ray([-1, 0]), Ray([0, 1])]</code></pre></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L287-L298">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.nallrays" href="#Polyhedra.nallrays"><code>Polyhedra.nallrays</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">nallrays(vrep::VRep)</code></pre><p>Returns the number of rays plus twice the number of lines in the V-representation <code>vrep</code>, i.e. <code>length(allrays(vrep))</code></p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L301-L305">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.hasallrays" href="#Polyhedra.hasallrays"><code>Polyhedra.hasallrays</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">hasallrays(vrep::VRep)</code></pre><p>Returns whether the V-representation <code>vrep</code> contains any ray or line.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/a5b4772fef0aea563d16d577b43e43a10ab0b841/src/iterators.jl#L308-L312">source</a></section></article></article><nav class="docs-footer"><a class="docs-footer-prevpage" href="../installation/">« Installation</a><a class="docs-footer-nextpage" href="../polyhedron/">Polyhedron »</a><div class="flexbox-break"></div><p class="footer-message">Powered by <a href="https://github.com/JuliaDocs/Documenter.jl">Documenter.jl</a> and the <a href="https://julialang.org/">Julia Programming Language</a>.</p></nav></div><div class="modal" id="documenter-settings"><div class="modal-background"></div><div class="modal-card"><header class="modal-card-head"><p class="modal-card-title">Settings</p><button class="delete"></button></header><section class="modal-card-body"><p><label class="label">Theme</label><div class="select"><select id="documenter-themepicker"><option value="documenter-light">documenter-light</option><option value="documenter-dark">documenter-dark</option></select></div></p><hr/><p>This document was generated with <a href="https://github.com/JuliaDocs/Documenter.jl">Documenter.jl</a> on <span class="colophon-date" title="Friday 13 November 2020 20:41">Friday 13 November 2020</span>. Using Julia version 1.5.3.</p></section><footer class="modal-card-foot"></footer></div></div></div></body></html>