/
index.html
2 lines (2 loc) · 12.9 KB
/
index.html
1
2
<!DOCTYPE html>
<html lang="en"><head><meta charset="UTF-8"/><meta name="viewport" content="width=device-width, initial-scale=1.0"/><title>Containment/Redundancy · Polyhedra</title><link href="https://cdnjs.cloudflare.com/ajax/libs/normalize/4.2.0/normalize.min.css" rel="stylesheet" type="text/css"/><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/4.6.3/css/font-awesome.min.css" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/highlight.js/9.12.0/styles/default.min.css" rel="stylesheet" type="text/css"/><script>documenterBaseURL=".."</script><script src="https://cdnjs.cloudflare.com/ajax/libs/require.js/2.2.0/require.min.js" data-main="../assets/documenter.js"></script><script src="../siteinfo.js"></script><script src="../../versions.js"></script><link href="../assets/documenter.css" rel="stylesheet" type="text/css"/></head><body><nav class="toc"><h1>Polyhedra</h1><select id="version-selector" onChange="window.location.href=this.value" style="visibility: hidden"></select><form class="search" id="search-form" action="../search/"><input id="search-query" name="q" type="text" placeholder="Search docs"/></form><ul><li><a class="toctext" href="../">Index</a></li><li><a class="toctext" href="../installation/">Installation</a></li><li><a class="toctext" href="../representation/">Representation</a></li><li><a class="toctext" href="../polyhedron/">Polyhedron</a></li><li><a class="toctext" href="../plot/">Plot</a></li><li class="current"><a class="toctext" href>Containment/Redundancy</a><ul class="internal"><li><a class="toctext" href="#Containment-1">Containment</a></li><li><a class="toctext" href="#Linearity-1">Linearity</a></li><li><a class="toctext" href="#Duplicates-1">Duplicates</a></li><li><a class="toctext" href="#Redundancy-1">Redundancy</a></li></ul></li><li><a class="toctext" href="../projection/">Projection/Elimination</a></li><li><a class="toctext" href="../optimization/">Optimization</a></li><li><a class="toctext" href="../utilities/">Utilities</a></li></ul></nav><article id="docs"><header><nav><ul><li><a href>Containment/Redundancy</a></li></ul><a class="edit-page" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/master/docs/src/redundancy.md"><span class="fa"></span> Edit on GitHub</a></nav><hr/><div id="topbar"><span>Containment/Redundancy</span><a class="fa fa-bars" href="#"></a></div></header><h1><a class="nav-anchor" id="Containment/Redundancy-1" href="#Containment/Redundancy-1">Containment/Redundancy</a></h1><h2><a class="nav-anchor" id="Containment-1" href="#Containment-1">Containment</a></h2><section class="docstring"><div class="docstring-header"><a class="docstring-binding" id="Base.in" href="#Base.in"><code>Base.in</code></a> — <span class="docstring-category">Function</span>.</div><div><div><pre><code class="language-none">in(p::VRepElement, h::HRepElement)</code></pre><p>Returns whether <code>p</code> is in <code>h</code>. If <code>h</code> is an hyperplane, it returns whether <span>$\langle a, x \rangle \approx \beta$</span>. If <code>h</code> is an halfspace, it returns whether <span>$\langle a, x \rangle \le \beta$</span>.</p><pre><code class="language-none">in(p::VRepElement, h::HRep)</code></pre><p>Returns whether <code>p</code> is in <code>h</code>, e.g. in all the hyperplanes and halfspaces supporting <code>h</code>.</p></div></div></section><section class="docstring"><div class="docstring-header"><a class="docstring-binding" id="Base.issubset" href="#Base.issubset"><code>Base.issubset</code></a> — <span class="docstring-category">Function</span>.</div><div><div><pre><code class="language-none">issubset(p::Rep, h::HRepElement)</code></pre><p>Returns whether <code>p</code> is a subset of <code>h</code>, i.e. whether <code>h</code> supports the polyhedron <code>p</code>.</p></div></div></section><section class="docstring"><div class="docstring-header"><a class="docstring-binding" id="Polyhedra.ininterior" href="#Polyhedra.ininterior"><code>Polyhedra.ininterior</code></a> — <span class="docstring-category">Function</span>.</div><div><div><pre><code class="language-none">ininterior(p::VRepElement, h::HRepElement)</code></pre><p>Returns whether <code>p</code> is in the interior of <code>h</code>. If <code>h</code> is an hyperplane, it always returns <code>false</code>. If <code>h</code> is an halfspace <span>$\langle a, x \rangle \leq \beta$</span>, it returns whether <code>p</code> is in the open halfspace <span>$\langle a, x \rangle < \beta$</span></p><pre><code class="language-none">ininterior(p::VRepElement, h::HRep)</code></pre><p>Returns whether <code>p</code> is in the interior of <code>h</code>, e.g. in the interior of all the hyperplanes and halfspaces supporting <code>h</code>.</p></div></div></section><section class="docstring"><div class="docstring-header"><a class="docstring-binding" id="Polyhedra.inrelativeinterior" href="#Polyhedra.inrelativeinterior"><code>Polyhedra.inrelativeinterior</code></a> — <span class="docstring-category">Function</span>.</div><div><div><pre><code class="language-none">inrelativeinterior(p::VRepElement, h::HRepElement)</code></pre><p>Returns whether <code>p</code> is in the relative interior of <code>h</code>. If <code>h</code> is an hyperplane, it is equivalent to <code>p in h</code> since the relative interior of an hyperplane is itself. If <code>h</code> is an halfspace, it is equivalent to <code>ininterior(p, h)</code>.</p><pre><code class="language-none">inrelativeinterior(p::VRepElement, h::HRep)</code></pre><p>Returns whether <code>p</code> is in the relative interior of <code>h</code>, e.g. in the relative interior of all the hyperplanes and halfspaces supporting <code>h</code>.</p></div></div></section><section class="docstring"><div class="docstring-header"><a class="docstring-binding" id="Polyhedra.support_function" href="#Polyhedra.support_function"><code>Polyhedra.support_function</code></a> — <span class="docstring-category">Function</span>.</div><div><div><p>support<em>function(h::AbstractVector, rep::Rep, solver=Polyhedra.linear</em>objective_solver(p))</p><p>Return the value of the support function of <code>rep</code> at <code>h</code>. See [Section 13, R15] for more details.</p><p>[R15] Rockafellar, R.T. <em>Convex analysis</em>. Princeton university press, <strong>2015</strong>.</p></div></div></section><h2><a class="nav-anchor" id="Linearity-1" href="#Linearity-1">Linearity</a></h2><section class="docstring"><div class="docstring-header"><a class="docstring-binding" id="Polyhedra.detecthlinearity!" href="#Polyhedra.detecthlinearity!"><code>Polyhedra.detecthlinearity!</code></a> — <span class="docstring-category">Function</span>.</div><div><div><pre><code class="language-none">detecthlinearity!(p::VRep)</code></pre><p>Detects all the hyperplanes contained in the H-representation and remove all redundant hyperplanes.</p><p><strong>Examples</strong></p><p>The representation</p><pre><code class="language-julia">h = HalfSpace([1, 1], 1]) ∩ HalfSpace([-1, -1], -1)</code></pre><p>contains the hyperplane <code>HyperPlane([1, 1], 1)</code>.</p></div></div></section><section class="docstring"><div class="docstring-header"><a class="docstring-binding" id="Polyhedra.detectvlinearity!" href="#Polyhedra.detectvlinearity!"><code>Polyhedra.detectvlinearity!</code></a> — <span class="docstring-category">Function</span>.</div><div><div><pre><code class="language-none">detectvlinearity!(p::VRep)</code></pre><p>Detects all the lines contained in the V-representation and remove all redundant lines.</p><p><strong>Examples</strong></p><p>The representation</p><pre><code class="language-julia">v = conichull([1, 1], [-1, -1])</code></pre><p>contains the line <code>Line([1, 1])</code>.</p></div></div></section><section class="docstring"><div class="docstring-header"><a class="docstring-binding" id="Polyhedra.dim" href="#Polyhedra.dim"><code>Polyhedra.dim</code></a> — <span class="docstring-category">Function</span>.</div><div><div><pre><code class="language-none">dim(h::HRep, current=false)</code></pre><p>Returns the dimension of the affine hull of the polyhedron. That is the number of non-redundant hyperplanes that define it. If <code>current</code> is <code>true</code> 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 <a href="#Polyhedra.detecthlinearity!"><code>detecthlinearity!</code></a>.</p></div></div></section><h2><a class="nav-anchor" id="Duplicates-1" href="#Duplicates-1">Duplicates</a></h2><section class="docstring"><div class="docstring-header"><a class="docstring-binding" id="Polyhedra.removeduplicates" href="#Polyhedra.removeduplicates"><code>Polyhedra.removeduplicates</code></a> — <span class="docstring-category">Function</span>.</div><div><div><pre><code class="language-none">removeduplicates(rep::Representation)</code></pre><p>Removes the duplicates in the Representation.</p><ul><li>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, <code>HalfSpace([1, 1], 1)</code> is equal to <code>HalfSpace([1, 0], 0)</code> in the affine hull generated by <code>HyperPlane([0, 1], 1])</code>.</li><li>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 <code>Line([0, 1])</code>, <code>[1, 1]</code> is equal to <code>[1, 0]</code> and <code>Ray([2, 2])</code> is equal to <code>Ray([1, 0])</code>.</li></ul></div></div></section><h2><a class="nav-anchor" id="Redundancy-1" href="#Redundancy-1">Redundancy</a></h2><section class="docstring"><div class="docstring-header"><a class="docstring-binding" id="Polyhedra.isredundant" href="#Polyhedra.isredundant"><code>Polyhedra.isredundant</code></a> — <span class="docstring-category">Function</span>.</div><div><div><pre><code class="language-none">isredundant(p::Rep, idx::Index; strongly=false)</code></pre><p>Return a <code>Bool</code> indicating whether the element with index <code>idx</code> can be removed without changing the polyhedron represented by <code>p</code>. If <code>strongly</code> is <code>true</code>,</p><ul><li>if <code>idx</code> is an H-representation element <code>h</code>, it returns <code>true</code> only if no V-representation element of <code>p</code> is in the hyperplane of <code>h</code>.</li><li>if <code>idx</code> is a V-representation element <code>v</code>, it returns <code>true</code> only if <code>v</code> is in the relative interior of <code>p</code>.</li></ul></div></div></section><section class="docstring"><div class="docstring-header"><a class="docstring-binding" id="Polyhedra.removehredundancy!" href="#Polyhedra.removehredundancy!"><code>Polyhedra.removehredundancy!</code></a> — <span class="docstring-category">Function</span>.</div><div><div><pre><code class="language-none">removehredundancy!(p::HRep)</code></pre><p>Removes the elements of the H-representation of <code>p</code> that can be removed without changing the polyhedron represented by <code>p</code>. That is, it only keeps the halfspaces corresponding to facets of the polyhedron.</p></div></div></section><section class="docstring"><div class="docstring-header"><a class="docstring-binding" id="Polyhedra.removevredundancy" href="#Polyhedra.removevredundancy"><code>Polyhedra.removevredundancy</code></a> — <span class="docstring-category">Function</span>.</div><div><div><pre><code class="language-none">removevredundancy(vr::VRepresentation)</code></pre><p>Return a V-representation of the polyhedron represented by <code>vr</code> all the elements of <code>vr</code> except the redundant ones, i.e. the elements that can be expressed as convex combination of other ones.</p></div></div></section><section class="docstring"><div class="docstring-header"><a class="docstring-binding" id="Polyhedra.removevredundancy!" href="#Polyhedra.removevredundancy!"><code>Polyhedra.removevredundancy!</code></a> — <span class="docstring-category">Function</span>.</div><div><div><pre><code class="language-none">removevredundancy!(p::VRep; strongly=false)</code></pre><p>Removes the elements of the V-representation of <code>p</code> that can be removed without changing the polyhedron represented by <code>p</code>. 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 <code>strongly=true</code>, weakly redundant points, i.e., points that are not extreme but are not in the relative interior either, may be kept.</p></div></div></section><footer><hr/><a class="previous" href="../plot/"><span class="direction">Previous</span><span class="title">Plot</span></a><a class="next" href="../projection/"><span class="direction">Next</span><span class="title">Projection/Elimination</span></a></footer></article></body></html>