-
Notifications
You must be signed in to change notification settings - Fork 27
/
index.html
3 lines (3 loc) · 26.7 KB
/
index.html
1
2
3
<!DOCTYPE html>
<html lang="en"><head><meta charset="UTF-8"/><meta name="viewport" content="width=device-width, initial-scale=1.0"/><title>Utilities · 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.15.0/css/fontawesome.min.css" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/font-awesome/5.15.0/css/solid.min.css" rel="stylesheet" type="text/css"/><link href="https://cdnjs.cloudflare.com/ajax/libs/font-awesome/5.15.0/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" data-theme-primary-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><a class="tocitem" href="../representation/">Representation</a></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 class="is-active"><a class="tocitem" href>Utilities</a><ul class="internal"><li><a class="tocitem" href="#Operations"><span>Operations</span></a></li><li><a class="tocitem" href="#Volume"><span>Volume</span></a></li><li><a class="tocitem" href="#Largest-inscribed-ball-with-center"><span>Largest inscribed ball with center</span></a></li><li><a class="tocitem" href="#Chebyshev-center"><span>Chebyshev center</span></a></li><li><a class="tocitem" href="#Defining-new-representation"><span>Defining new representation</span></a></li></ul></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><li><a class="tocitem" href="../generated/Convex hull of a set of points/">Convex hull of a set of points</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>Utilities</a></li></ul><ul class="is-hidden-tablet"><li class="is-active"><a href>Utilities</a></li></ul></nav><div class="docs-right"><a class="docs-edit-link" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/master/docs/src/utilities.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="Utilities"><a class="docs-heading-anchor" href="#Utilities">Utilities</a><a id="Utilities-1"></a><a class="docs-heading-anchor-permalink" href="#Utilities" title="Permalink"></a></h1><h2 id="Operations"><a class="docs-heading-anchor" href="#Operations">Operations</a><a id="Operations-1"></a><a class="docs-heading-anchor-permalink" href="#Operations" title="Permalink"></a></h2><article class="docstring"><header><a class="docstring-binding" id="Base.:+" href="#Base.:+"><code>Base.:+</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">+(p1::VRep, p2::VRep)</code></pre><p>Minkowski sum between <code>p1</code> and <code>p2</code> using the V-representation. If the V-representation is not computed for <code>p1</code> or <code>p2</code>, it is computed.</p><pre><code class="language-none">+(p::Rep, el::Union{Line, Ray})
+(el::Union{Line, Ray}, p::Rep)</code></pre><p>Same as <code>p + vrep([el])</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/repop.jl#L123-L133">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Base.:*" href="#Base.:*"><code>Base.:*</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">*(p1::Rep, p2::Rep)</code></pre><p>Cartesian product between the polyhedra <code>p1</code> and <code>p2</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/repop.jl#L181-L185">source</a></section><section><div><pre><code class="language-none">*(P::Union{AbstractMatrix, UniformScaling}, p::VRep)</code></pre><p>Transform the polyhedron represented by <span>$p$</span> into <span>$P p$</span> by transforming each element of the V-representation (points, symmetric points, rays and lines) <code>x</code> into <span>$P x$</span>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/repop.jl#L224-L228">source</a></section><section><div><pre><code class="language-none">*(α::Number, p::Rep)</code></pre><p>Transform the polyhedron represented by <span>$p$</span> into <span>$\alpha p$</span> by transforming each element of the V-representation (points, symmetric points, rays and lines) <code>x</code> into <span>$\alpha x$</span>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/repop.jl#L240-L246">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Base.:\\" href="#Base.:\\"><code>Base.:\</code></a> — <span class="docstring-category">Function</span></header><section><div><p>(P::Union{AbstractMatrix, UniformScaling}, p::HRep)</p><p>Transform the polyhedron represented by <span>$p$</span> into <span>$P^{-1} p$</span> by transforming each halfspace <span>$\langle a, x \rangle \le \beta$</span> into <span>$\langle P^\top a, x \rangle \le \beta$</span> and each hyperplane <span>$\langle a, x \rangle = \beta$</span> into <span>$\langle P^\top a, x \rangle = \beta$</span>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/repop.jl#L188-L192">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Base.:/" href="#Base.:/"><code>Base.:/</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">/(p::HRep, P::Union{AbstractMatrix, UniformScaling})</code></pre><p>Transform the polyhedron represented by <span>$p$</span> into <span>$P^{-T} p$</span> by transforming each halfspace <span>$\langle a, x \rangle \le \beta$</span> into <span>$\langle P a, x \rangle \le \beta$</span> and each hyperplane <span>$\langle a, x \rangle = \beta$</span> into <span>$\langle P a, x \rangle = \beta$</span>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/repop.jl#L202-L206">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Base.intersect" href="#Base.intersect"><code>Base.intersect</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">intersect(P1::HRep, P2::HRep)</code></pre><p>Takes the intersection of <code>P1</code> and <code>P2</code> <span>$\{\, x : x \in P_1, x \in P_2 \,\}$</span>. It is very efficient between two H-representations or between two polyhedron for which the H-representation has already been computed. However, if <code>P1</code> (resp. <code>P2</code>) is a polyhedron for which the H-representation has not been computed yet, it will trigger a representation conversion which is costly. See the <a href="http://www.cs.mcgill.ca/~fukuda/soft/polyfaq/node25.html">Polyhedral Computation FAQ</a> for a discussion on this operation.</p><p>The type of the result will be chosen closer to the type of <code>P1</code>. For instance, if <code>P1</code> is a polyhedron (resp. H-representation) and <code>P2</code> is a H-representation (resp. polyhedron), <code>intersect(P1, P2)</code> will be a polyhedron (resp. H-representation). If <code>P1</code> and <code>P2</code> are both polyhedra (resp. H-representation), the resulting polyhedron type (resp. H-representation type) will be computed according to the type of <code>P1</code>. The coefficient type however, will be promoted as required taking both the coefficient type of <code>P1</code> and <code>P2</code> into account.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/repop.jl#L6-L17">source</a></section><section><div><pre><code class="language-none">intersect(v::VRepresentation{T}, h::HRepElement)</code></pre><p>Compute the intersection of <code>v</code> with an halfspace or hyperplane <code>h</code>. The method used by default is to keep the V-representation element of <code>v</code> that are in <code>h</code> and add new ones generated as the intersection between the hyperplane defining <code>h</code> and the segment between two adjacent V-representation elements of <code>v</code> that are in either sides of the hyperplane. See Lemma 3 of [FP96] for more detail on the method.</p><p>[FP96] Fukuda, K. and Prodon, A. <strong>Double description method revisited</strong> <em>Combinatorics and computer science</em>, <em>Springer</em>, <strong>1996</strong>, 91-111</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/repelemop.jl#L241-L254">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Base.intersect!" href="#Base.intersect!"><code>Base.intersect!</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">intersect!(p::HRep, h::Union{HRepresentation, HRepElement})</code></pre><p>Same as <a href="#Base.intersect"><code>intersect</code></a> except that <code>p</code> is modified to be equal to the intersection.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/repop.jl#L44-L48">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.convexhull" href="#Polyhedra.convexhull"><code>Polyhedra.convexhull</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">convexhull(P1::VRep, P2::VRep)</code></pre><p>Takes the convex hull of <code>P1</code> and <code>P2</code> <span>$\{\, \lambda x + (1-\lambda) y : x \in P_1, y \in P_2 \,\}$</span>. It is very efficient between two V-representations or between two polyhedron for which the V-representation has already been computed. However, if <code>P1</code> (resp. <code>P2</code>) is a polyhedron for which the V-representation has not been computed yet, it will trigger a representation conversion which is costly.</p><p>The type of the result will be chosen closer to the type of <code>P1</code>. For instance, if <code>P1</code> is a polyhedron (resp. V-representation) and <code>P2</code> is a V-representation (resp. polyhedron), <code>convexhull(P1, P2)</code> will be a polyhedron (resp. V-representation). If <code>P1</code> and <code>P2</code> are both polyhedra (resp. V-representation), the resulting polyhedron type (resp. V-representation type) will be computed according to the type of <code>P1</code>. The coefficient type however, will be promoted as required taking both the coefficient type of <code>P1</code> and <code>P2</code> into account.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/repop.jl#L56-L66">source</a></section><section><div><pre><code class="language-none">convexhull(p1::HRepresentation, p2::HRepresentation)</code></pre><p>Returns the Balas [Theorem 3.3, B85] extended H-representation of the convex hull of <code>p1</code> and <code>p2</code>.</p><p>[B85] Balas, E., 1985. <em>Disjunctive programming and a hierarchy of relaxations for discrete optimization problems</em>. SIAM Journal on Algebraic Discrete Methods, 6(3), pp.466-486.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/extended.jl#L74-L83">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.convexhull!" href="#Polyhedra.convexhull!"><code>Polyhedra.convexhull!</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">convexhull!(p1::VRep, p2::VRep)</code></pre><p>Same as <a href="#Polyhedra.convexhull"><code>convexhull</code></a> except that <code>p1</code> is modified to be equal to the convex hull.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/repop.jl#L95-L99">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.translate" href="#Polyhedra.translate"><code>Polyhedra.translate</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">translate(p::Polyhedra.Rep, v::AbstractVector)</code></pre><p>Computes translation of the polyhedron <code>p</code> with the vector <code>v</code>. That is, computes</p><p class="math-container">\[\{\, x + v \mid x \in p \,\}.\]</p><p>By default, if the H-representation, it simply translates every hyperplanes and halfspace, otherwise, it translates every points of the V-representation. That is, this operation can be achieved both in the H-representation and V-representation hence does not trigger any representation conversion.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/repelemop.jl#L32-L44">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.polar" href="#Polyhedra.polar"><code>Polyhedra.polar</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">polar(rep::Representation)</code></pre><p>Return the polar of the polyhedron <code>rep</code> assumed to contain the origin. The polar of a convex set <code>S</code> is defined as the set of <code>y</code> such that <code>⟨x, y⟩ ≤ 1</code> for all <code>x in S</code>. Note that the polar of a V-representation is a H-representation and vice versa.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/repop.jl#L261-L270">source</a></section></article><h2 id="Volume"><a class="docs-heading-anchor" href="#Volume">Volume</a><a id="Volume-1"></a><a class="docs-heading-anchor-permalink" href="#Volume" title="Permalink"></a></h2><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.volume" href="#Polyhedra.volume"><code>Polyhedra.volume</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">volume(p::Polyhedron{T}) where {T}</code></pre><p>Returns the <code>fulldim(p)</code>-dimensional hyper-volume of the polyhedron <code>p</code>. Returns <code>Inf</code> or <code>-one(T)</code> if it is infinite depending on whether the type <code>T</code> has an infinite value.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/polyhedron.jl#L41-L46">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.surface" href="#Polyhedra.surface"><code>Polyhedra.surface</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">surface(p::Polyhedron{T}) where {T}</code></pre><p>Returns the <code>fulldim(p)-1</code>-dimensional hyper-volume of the surface of the polyhedron <code>p</code>. Returns <code>Inf</code> or <code>-one(T)</code> if it is infinite depending on whether the type <code>T</code> has an infinite value.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/polyhedron.jl#L87-L92">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.center_of_mass" href="#Polyhedra.center_of_mass"><code>Polyhedra.center_of_mass</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">center_of_mass(p::Polyhedron{T}) where {T}</code></pre><p>Returns the center of mass of <code>p</code>, represented as a <code>Vector{T}</code> of length <code>fulldim(p)</code>. Throws an error if <code>p</code> is degenerate.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/polyhedron.jl#L64-L69">source</a></section></article><h2 id="Largest-inscribed-ball-with-center"><a class="docs-heading-anchor" href="#Largest-inscribed-ball-with-center">Largest inscribed ball with center</a><a id="Largest-inscribed-ball-with-center-1"></a><a class="docs-heading-anchor-permalink" href="#Largest-inscribed-ball-with-center" title="Permalink"></a></h2><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.maximum_radius_with_center" href="#Polyhedra.maximum_radius_with_center"><code>Polyhedra.maximum_radius_with_center</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">maximum_radius_with_center(h::HRep, center)</code></pre><p>Return the maximum radius <code>r</code> such that the Euclidean ball of center <code>center</code> and radius <code>r</code> is included in the polyhedron <code>h</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/center.jl#L4-L9">source</a></section></article><h2 id="Chebyshev-center"><a class="docs-heading-anchor" href="#Chebyshev-center">Chebyshev center</a><a id="Chebyshev-center-1"></a><a class="docs-heading-anchor-permalink" href="#Chebyshev-center" title="Permalink"></a></h2><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.chebyshevcenter" href="#Polyhedra.chebyshevcenter"><code>Polyhedra.chebyshevcenter</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">chebyshevcenter(p::Rep[, solver])</code></pre><p>If <code>p</code> is a H-representation or is a polyhedron for which the H-representation has already been computed, calls <code>hchebyshevcenter</code>, otherwise, call <code>vchebyshevcenter</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/center.jl#L139-L143">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.hchebyshevcenter" href="#Polyhedra.hchebyshevcenter"><code>Polyhedra.hchebyshevcenter</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">hchebyshevcenter(p::HRep[, solver]; linearity_detected=false, proper=true)</code></pre><p>Return a tuple with the center and radius of the largest euclidean ball contained in the polyhedron <code>p</code>. Throws an error if the polyhedron is empty or if the radius is infinite. Linearity is detected first except if <code>linearity_detected</code>.</p><p>Note that a polytope may have several Chebyshev center. In general, the set of Chebyshev center of a polytope <code>p</code> is a polytope which has a lower dimension than <code>p</code> if <code>p</code> has a positive dimension. For instance, if <code>p</code> is the rectangle <code>[-2, 2] x [-1, 1]</code>, the Chebyshev radius of <code>p</code> is 1 and the set of Chebyshev centers is <code>[-1, 1] x {0}</code>. The <em>proper</em> Chebyshev center is <code>(0, 0)</code>, the Chebyshev center of <code>[-1, 1] x {0}</code>. If <code>!proper</code> then any Chebyshev center is returned (the one returned depends on the solver). Otherwise the proper Chebyshev center is computed. The proper Chebyshev center is defined by induction on the dimension of <code>p</code>. If <code>p</code> has dimension 0 then it is a singleton and its proper Chebyshev center is the only element of <code>p</code>. Otherwise, the dimension of the set <code>q</code> of Chebyshev centers of <code>p</code> is smaller than the dimension of <code>p</code> and the proper Chebyshev center of <code>p</code> is the proper Chebyshev center of <code>q</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/center.jl#L50-L69">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.vchebyshevcenter" href="#Polyhedra.vchebyshevcenter"><code>Polyhedra.vchebyshevcenter</code></a> — <span class="docstring-category">Function</span></header><section><div><pre><code class="language-julia">vchebyshevcenter(p::VRep[, solver])</code></pre><p>Return a tuple with the center and radius of the smallest euclidean ball containing the polyhedron <code>p</code>. Throws an error if the polyhedron is empty or if the radius is infinite (i.e. <code>p</code> is not a polytope, it contains rays).</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/center.jl#L129-L134">source</a></section></article><h2 id="Defining-new-representation"><a class="docs-heading-anchor" href="#Defining-new-representation">Defining new representation</a><a id="Defining-new-representation-1"></a><a class="docs-heading-anchor-permalink" href="#Defining-new-representation" title="Permalink"></a></h2><p>The following macros make it easy to define new representations:</p><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.@subrepelem" href="#Polyhedra.@subrepelem"><code>Polyhedra.@subrepelem</code></a> — <span class="docstring-category">Macro</span></header><section><div><p>The representation <code>rep</code> contain the elements <code>elem</code> inside a representation in the field <code>field</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/indices.jl#L165-L167">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.@norepelem" href="#Polyhedra.@norepelem"><code>Polyhedra.@norepelem</code></a> — <span class="docstring-category">Macro</span></header><section><div><p>The representation <code>rep</code> does not contain any <code>elem</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/indices.jl#L90-L92">source</a></section></article><article class="docstring"><header><a class="docstring-binding" id="Polyhedra.@vecrepelem" href="#Polyhedra.@vecrepelem"><code>Polyhedra.@vecrepelem</code></a> — <span class="docstring-category">Macro</span></header><section><div><p>The representation <code>rep</code> contain the elements <code>elem</code> inside a vector in the field <code>field</code>.</p></div><a class="docs-sourcelink" target="_blank" href="https://github.com/JuliaPolyhedra/Polyhedra.jl/blob/0288f2758e1e84f2896122b08a010febd0b6a2a2/src/indices.jl#L138-L140">source</a></section></article></article><nav class="docs-footer"><a class="docs-footer-prevpage" href="../optimization/">« Optimization</a><a class="docs-footer-nextpage" href="../generated/Convex hull and intersection/">Convex hull and intersection »</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="Thursday 11 November 2021 14:13">Thursday 11 November 2021</span>. Using Julia version 1.6.3.</p></section><footer class="modal-card-foot"></footer></div></div></div></body></html>