Polyhedra¶
In this module, a polyhedron is a convex (possibly unbounded) set in Euclidean space cut out by a finite set of linear inequalities and linear equations. Note that the dimension of the polyhedron can be less than the dimension of the ambient space. There are two complementary representations of the same data:
- H(alf-space/Hyperplane)-representation
This describes a polyhedron as the common solution set of a finite number of
linear inequalities \(A \vec{x} + b \geq 0\), and
linear equations \(C \vec{x} + d = 0\).
- V(ertex)-representation
The other representation is as the convex hull of vertices (and rays and lines to all for unbounded polyhedra) as generators. The polyhedron is then the Minkowski sum
\[P = \text{conv}\{v_1,\dots,v_k\} + \sum_{i=1}^m \RR_+ r_i + \sum_{j=1}^n \RR \ell_j\]where
vertices \(v_1,\dots,v_k\) are a finite number of points. Each vertex is specified by an arbitrary vector, and two points are equal if and only if the vector is the same.
rays \(r_1,\dots,r_m\) are a finite number of directions (directions of infinity). Each ray is specified by a nonzero vector, and two rays are equal if and only if the vectors are the same up to rescaling with a positive constant.
lines \(\ell_1,\dots,\ell_n\) are a finite number of unoriented directions. In other words, a line is equivalent to the set \(\{r, -r\}\) for a ray \(r\). Each line is specified by a nonzero vector, and two lines are equivalent if and only if the vectors are the same up to rescaling with a nonzero (possibly negative) constant.
When specifying a polyhedron, you can input a non-minimal set of inequalities/equations or generating vertices/rays/lines. The non-minimal generators are usually called points, non-extremal rays, and non-extremal lines, but for our purposes it is more convenient to always talk about vertices/rays/lines. Sage will remove any superfluous representation objects and always return a minimal representation. For example, \((0,0)\) is a superfluous vertex here:
sage: triangle = Polyhedron(vertices=[(0,2), (-1,0), (1,0), (0,0)])
sage: triangle.vertices()
(A vertex at (-1, 0), A vertex at (1, 0), A vertex at (0, 2))
>>> from sage.all import *
>>> triangle = Polyhedron(vertices=[(Integer(0),Integer(2)), (-Integer(1),Integer(0)), (Integer(1),Integer(0)), (Integer(0),Integer(0))])
>>> triangle.vertices()
(A vertex at (-1, 0), A vertex at (1, 0), A vertex at (0, 2))
See also
If one only needs to keep track of a system of linear system of inequalities, one should also consider the class for mixed integer linear programming.
Unbounded Polyhedra¶
A polytope is defined as a bounded polyhedron. In this case, the minimal representation is unique and a vertex of the minimal representation is equivalent to a 0-dimensional face of the polytope. This is why one generally does not distinguish vertices and 0-dimensional faces. But for non-bounded polyhedra we have to allow for a more general notion of “vertex” in order to make sense of the Minkowski sum presentation:
sage: half_plane = Polyhedron(ieqs=[(0,1,0)])
sage: half_plane.Hrepresentation()
(An inequality (1, 0) x + 0 >= 0,)
sage: half_plane.Vrepresentation()
(A line in the direction (0, 1), A ray in the direction (1, 0), A vertex at (0, 0))
>>> from sage.all import *
>>> half_plane = Polyhedron(ieqs=[(Integer(0),Integer(1),Integer(0))])
>>> half_plane.Hrepresentation()
(An inequality (1, 0) x + 0 >= 0,)
>>> half_plane.Vrepresentation()
(A line in the direction (0, 1), A ray in the direction (1, 0), A vertex at (0, 0))
Note how we need a point in the above example to anchor the ray and line. But any point on the boundary of the half-plane would serve the purpose just as well. Sage picked the origin here, but this choice is not unique. Similarly, the choice of ray is arbitrary but necessary to generate the half-plane.
Finally, note that while rays and lines generate unbounded edges of the polyhedron they are not in a one-to-one correspondence with them. For example, the infinite strip has two infinite edges (1-faces) but only one generating line:
sage: strip = Polyhedron(vertices=[(1,0),(-1,0)], lines=[(0,1)])
sage: strip.Hrepresentation()
(An inequality (1, 0) x + 1 >= 0, An inequality (-1, 0) x + 1 >= 0)
sage: strip.lines()
(A line in the direction (0, 1),)
sage: [f.ambient_V_indices() for f in strip.faces(1)]
[(0, 2), (0, 1)]
sage: for face in strip.faces(1):
....: print(face.ambient_V_indices())
(0, 2)
(0, 1)
sage: for face in strip.faces(1):
....: print("{} = {}".format(face.ambient_V_indices(), face.as_polyhedron().Vrepresentation()))
(0, 2) = (A line in the direction (0, 1), A vertex at (1, 0))
(0, 1) = (A line in the direction (0, 1), A vertex at (-1, 0))
>>> from sage.all import *
>>> strip = Polyhedron(vertices=[(Integer(1),Integer(0)),(-Integer(1),Integer(0))], lines=[(Integer(0),Integer(1))])
>>> strip.Hrepresentation()
(An inequality (1, 0) x + 1 >= 0, An inequality (-1, 0) x + 1 >= 0)
>>> strip.lines()
(A line in the direction (0, 1),)
>>> [f.ambient_V_indices() for f in strip.faces(Integer(1))]
[(0, 2), (0, 1)]
>>> for face in strip.faces(Integer(1)):
... print(face.ambient_V_indices())
(0, 2)
(0, 1)
>>> for face in strip.faces(Integer(1)):
... print("{} = {}".format(face.ambient_V_indices(), face.as_polyhedron().Vrepresentation()))
(0, 2) = (A line in the direction (0, 1), A vertex at (1, 0))
(0, 1) = (A line in the direction (0, 1), A vertex at (-1, 0))
EXAMPLES:
sage: trunc_quadr = Polyhedron(vertices=[[1,0],[0,1]], rays=[[1,0],[0,1]])
sage: trunc_quadr
A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 2 vertices and 2 rays
sage: v = next(trunc_quadr.vertex_generator()) # the first vertex in the internal enumeration
sage: v
A vertex at (0, 1)
sage: v.vector()
(0, 1)
sage: list(v)
[0, 1]
sage: len(v)
2
sage: v[0] + v[1]
1
sage: v.is_vertex()
True
sage: type(v)
<class 'sage.geometry.polyhedron.representation.Vertex'>
sage: type( v() )
<class 'sage.modules.vector_rational_dense.Vector_rational_dense'>
sage: v.polyhedron()
A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 2 vertices and 2 rays
sage: r = next(trunc_quadr.ray_generator())
sage: r
A ray in the direction (0, 1)
sage: r.vector()
(0, 1)
sage: list( v.neighbors() )
[A ray in the direction (0, 1), A vertex at (1, 0)]
>>> from sage.all import *
>>> trunc_quadr = Polyhedron(vertices=[[Integer(1),Integer(0)],[Integer(0),Integer(1)]], rays=[[Integer(1),Integer(0)],[Integer(0),Integer(1)]])
>>> trunc_quadr
A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 2 vertices and 2 rays
>>> v = next(trunc_quadr.vertex_generator()) # the first vertex in the internal enumeration
>>> v
A vertex at (0, 1)
>>> v.vector()
(0, 1)
>>> list(v)
[0, 1]
>>> len(v)
2
>>> v[Integer(0)] + v[Integer(1)]
1
>>> v.is_vertex()
True
>>> type(v)
<class 'sage.geometry.polyhedron.representation.Vertex'>
>>> type( v() )
<class 'sage.modules.vector_rational_dense.Vector_rational_dense'>
>>> v.polyhedron()
A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 2 vertices and 2 rays
>>> r = next(trunc_quadr.ray_generator())
>>> r
A ray in the direction (0, 1)
>>> r.vector()
(0, 1)
>>> list( v.neighbors() )
[A ray in the direction (0, 1), A vertex at (1, 0)]
Inequalities \(A \vec{x} + b \geq 0\) (and, similarly, equations) are
specified by a list [b, A]
:
sage: Polyhedron(ieqs=[(0,1,0),(0,0,1),(1,-1,-1)]).Hrepresentation()
(An inequality (-1, -1) x + 1 >= 0,
An inequality (1, 0) x + 0 >= 0,
An inequality (0, 1) x + 0 >= 0)
>>> from sage.all import *
>>> Polyhedron(ieqs=[(Integer(0),Integer(1),Integer(0)),(Integer(0),Integer(0),Integer(1)),(Integer(1),-Integer(1),-Integer(1))]).Hrepresentation()
(An inequality (-1, -1) x + 1 >= 0,
An inequality (1, 0) x + 0 >= 0,
An inequality (0, 1) x + 0 >= 0)
See Polyhedron()
for a detailed description of all possible ways
to construct a polyhedron.
Base Rings¶
The base ring of the polyhedron can be specified by the base_ring
optional keyword argument. If not specified, a suitable common base
ring for all coordinates and coefficients will be chosen
automatically. Important cases are:
base_ring=QQ
uses a fast implementation for exact rational numbers.base_ring=ZZ
is similar toQQ
, but the resulting polyhedron object will have extra methods for lattice polyhedra.base_ring=RDF
uses floating point numbers, this is fast but susceptible to numerical errors.
Polyhedra with symmetries often are defined over some algebraic field
extension of the rationals. As a simple example, consider the
equilateral triangle whose vertex coordinates involve \(\sqrt{3}\). An
exact way to work with roots in Sage is the
Algebraic Real Field
sage: triangle = Polyhedron([(0,0), (1,0), (1/2, sqrt(3)/2)], base_ring=AA) # needs sage.rings.number_field sage.symbolic
sage: triangle.Hrepresentation() # needs sage.rings.number_field sage.symbolic
(An inequality (-1, -0.5773502691896258?) x + 1 >= 0,
An inequality (1, -0.5773502691896258?) x + 0 >= 0,
An inequality (0, 1.154700538379252?) x + 0 >= 0)
>>> from sage.all import *
>>> triangle = Polyhedron([(Integer(0),Integer(0)), (Integer(1),Integer(0)), (Integer(1)/Integer(2), sqrt(Integer(3))/Integer(2))], base_ring=AA) # needs sage.rings.number_field sage.symbolic
>>> triangle.Hrepresentation() # needs sage.rings.number_field sage.symbolic
(An inequality (-1, -0.5773502691896258?) x + 1 >= 0,
An inequality (1, -0.5773502691896258?) x + 0 >= 0,
An inequality (0, 1.154700538379252?) x + 0 >= 0)
Without specifying the base_ring
, the sqrt(3)
would be a
symbolic ring element and, therefore, the polyhedron defined over the
symbolic ring. This is currently not supported as SR is not exact:
sage: Polyhedron([(0,0), (1,0), (1/2, sqrt(3)/2)]) # needs sage.symbolic
Traceback (most recent call last):
...
ValueError: no default backend for computations with Symbolic Ring
sage: SR.is_exact() # needs sage.symbolic
False
>>> from sage.all import *
>>> Polyhedron([(Integer(0),Integer(0)), (Integer(1),Integer(0)), (Integer(1)/Integer(2), sqrt(Integer(3))/Integer(2))]) # needs sage.symbolic
Traceback (most recent call last):
...
ValueError: no default backend for computations with Symbolic Ring
>>> SR.is_exact() # needs sage.symbolic
False
Even faster than all algebraic real numbers (the field AA
) is
to take the smallest extension field. For the equilateral
triangle, that would be:
sage: x = polygen(ZZ, 'x')
sage: K.<sqrt3> = NumberField(x^2 - 3, embedding=AA(3)**(1/2)) # needs sage.rings.number_field
sage: Polyhedron([(0,0), (1,0), (1/2, sqrt3/2)]) # needs sage.rings.number_field
A 2-dimensional polyhedron in
(Number Field in sqrt3 with defining polynomial x^2 - 3 with sqrt3 = 1.732050807568878?)^2
defined as the convex hull of 3 vertices
>>> from sage.all import *
>>> x = polygen(ZZ, 'x')
>>> K = NumberField(x**Integer(2) - Integer(3), embedding=AA(Integer(3))**(Integer(1)/Integer(2)), names=('sqrt3',)); (sqrt3,) = K._first_ngens(1)# needs sage.rings.number_field
>>> Polyhedron([(Integer(0),Integer(0)), (Integer(1),Integer(0)), (Integer(1)/Integer(2), sqrt3/Integer(2))]) # needs sage.rings.number_field
A 2-dimensional polyhedron in
(Number Field in sqrt3 with defining polynomial x^2 - 3 with sqrt3 = 1.732050807568878?)^2
defined as the convex hull of 3 vertices
Warning
Be careful when you construct polyhedra with floating point numbers. The only
available backend for such computation is cdd
which uses machine floating
point numbers which have limited precision. If the input consists of
floating point numbers and the base_ring
is not specified, the base ring is
set to be the RealField
with the precision given by the minimal bit precision
of the input. Then, if the obtained minimum is 53 bits of precision, the
constructor converts automatically the base ring to RDF
. Otherwise,
it returns an error:
sage: Polyhedron(vertices = [[1.12345678901234, 2.12345678901234]])
A 0-dimensional polyhedron in RDF^2 defined as the convex hull of 1 vertex
sage: Polyhedron(vertices = [[1.12345678901234, 2.123456789012345]])
A 0-dimensional polyhedron in RDF^2 defined as the convex hull of 1 vertex
sage: Polyhedron(vertices = [[1.123456789012345, 2.123456789012345]]) # needs sage.rings.real_mpfr
Traceback (most recent call last):
...
ValueError: the only allowed inexact ring is 'RDF' with backend 'cdd'
>>> from sage.all import *
>>> Polyhedron(vertices = [[RealNumber('1.12345678901234'), RealNumber('2.12345678901234')]])
A 0-dimensional polyhedron in RDF^2 defined as the convex hull of 1 vertex
>>> Polyhedron(vertices = [[RealNumber('1.12345678901234'), RealNumber('2.123456789012345')]])
A 0-dimensional polyhedron in RDF^2 defined as the convex hull of 1 vertex
>>> Polyhedron(vertices = [[RealNumber('1.123456789012345'), RealNumber('2.123456789012345')]]) # needs sage.rings.real_mpfr
Traceback (most recent call last):
...
ValueError: the only allowed inexact ring is 'RDF' with backend 'cdd'
The strongly suggested method to input floating point numbers is to specify the
base_ring
to be RDF
:
sage: Polyhedron(vertices = [[1.123456789012345, 2.123456789012345]], base_ring=RDF)
A 0-dimensional polyhedron in RDF^2 defined as the convex hull of 1 vertex
>>> from sage.all import *
>>> Polyhedron(vertices = [[RealNumber('1.123456789012345'), RealNumber('2.123456789012345')]], base_ring=RDF)
A 0-dimensional polyhedron in RDF^2 defined as the convex hull of 1 vertex
See also
Base classes¶
Depending on the chosen base ring, a specific class is used to represent the polyhedron object.
See also
The most important base class is Base class for polyhedra from which other base classes and backends inherit.
Backends¶
There are different backends available to deal with polyhedron objects.
See also
Note
Depending on the backend used, it may occur that different methods are available or not.
Appendix¶
REFERENCES:
Komei Fukuda’s FAQ in Polyhedral Computation
AUTHORS:
Marshall Hampton: first version, bug fixes, and various improvements, 2008 and 2009
Arnaud Bergeron: improvements to triangulation and rendering, 2008
Sebastien Barthelemy: documentation improvements, 2008
Volker Braun: refactoring, handle non-compact case, 2009 and 2010
Andrey Novoseltsev: added lattice_from_incidences, 2010
Volker Braun: rewrite to use PPL instead of cddlib, 2011
Volker Braun: Add support for arbitrary subfields of the reals
- sage.geometry.polyhedron.constructor.Polyhedron(vertices=None, rays=None, lines=None, ieqs=None, eqns=None, ambient_dim=None, base_ring=None, minimize=True, verbose=False, backend=None, mutable=False)[source]¶
Construct a polyhedron object.
You may either define it with vertex/ray/line or inequalities/equations data, but not both. Redundant data will automatically be removed (unless
minimize=False
), and the complementary representation will be computed.INPUT:
vertices
– iterable of points. Each point can be specified as any iterable container ofbase_ring
elements. Ifrays
orlines
are specified but novertices
, the origin is taken to be the single vertex.Instead of vertices, the first argument can also be an object that can be converted to a
Polyhedron()
via anas_polyhedron()
orpolyhedron()
method. In this case, the following 5 arguments cannot be provided.rays
– list of rays; each ray can be specified as any iterable container ofbase_ring
elementslines
– list of lines; each line can be specified as any iterable container ofbase_ring
elementsieqs
– list of inequalities; each line can be specified as any iterable container ofbase_ring
elements. An entry equal to[-1,7,3,4]
represents the inequality \(7x_1+3x_2+4x_3\geq 1\).eqns
– list of equalities; each line can be specified as any iterable container ofbase_ring
elements. An entry equal to[-1,7,3,4]
represents the equality \(7x_1+3x_2+4x_3= 1\).ambient_dim
– integer; the ambient space dimension. Usually can be figured out automatically from the H/Vrepresentation dimensions.base_ring
– a sub-field of the reals implemented in Sage. The field over which the polyhedron will be defined. ForQQ
and algebraic extensions, exact arithmetic will be used. ForRDF
, floating point numbers will be used. Floating point arithmetic is faster but might give the wrong result for degenerate input.backend
– string orNone
(default). The backend to use. Valid choices are'cdd'
: use cdd (backend_cdd
) with \(\QQ\) or \(\RDF\) coefficients depending onbase_ring
'normaliz'
: use normaliz (backend_normaliz
) with \(\ZZ\) or \(\QQ\) coefficients depending onbase_ring
'polymake'
: use polymake (backend_polymake
) with \(\QQ\), \(\RDF\) orQuadraticField
coefficients depending onbase_ring
'ppl'
: use ppl (backend_ppl
) with \(\ZZ\) or \(\QQ\) coefficients depending onbase_ring
'field'
: use python implementation (backend_field
) for any field
Some backends support further optional arguments:
minimize
– boolean (default:True
); whether to immediately remove redundant H/V-representation data; currently not used.verbose
– boolean (default:False
); whether to print verbose output for debugging purposes; only supported by the cdd and normaliz backendsmutable
– boolean (default:False
); whether the polyhedron is mutable
OUTPUT: the polyhedron defined by the input data
EXAMPLES:
Construct some polyhedra:
sage: square_from_vertices = Polyhedron(vertices = [[1, 1], [1, -1], [-1, 1], [-1, -1]]) sage: square_from_ieqs = Polyhedron(ieqs = [[1, 0, 1], [1, 1, 0], [1, 0, -1], [1, -1, 0]]) sage: list(square_from_ieqs.vertex_generator()) [A vertex at (1, -1), A vertex at (1, 1), A vertex at (-1, 1), A vertex at (-1, -1)] sage: list(square_from_vertices.inequality_generator()) [An inequality (1, 0) x + 1 >= 0, An inequality (0, 1) x + 1 >= 0, An inequality (-1, 0) x + 1 >= 0, An inequality (0, -1) x + 1 >= 0] sage: p = Polyhedron(vertices = [[1.1, 2.2], [3.3, 4.4]], base_ring=RDF) sage: p.n_inequalities() 2
>>> from sage.all import * >>> square_from_vertices = Polyhedron(vertices = [[Integer(1), Integer(1)], [Integer(1), -Integer(1)], [-Integer(1), Integer(1)], [-Integer(1), -Integer(1)]]) >>> square_from_ieqs = Polyhedron(ieqs = [[Integer(1), Integer(0), Integer(1)], [Integer(1), Integer(1), Integer(0)], [Integer(1), Integer(0), -Integer(1)], [Integer(1), -Integer(1), Integer(0)]]) >>> list(square_from_ieqs.vertex_generator()) [A vertex at (1, -1), A vertex at (1, 1), A vertex at (-1, 1), A vertex at (-1, -1)] >>> list(square_from_vertices.inequality_generator()) [An inequality (1, 0) x + 1 >= 0, An inequality (0, 1) x + 1 >= 0, An inequality (-1, 0) x + 1 >= 0, An inequality (0, -1) x + 1 >= 0] >>> p = Polyhedron(vertices = [[RealNumber('1.1'), RealNumber('2.2')], [RealNumber('3.3'), RealNumber('4.4')]], base_ring=RDF) >>> p.n_inequalities() 2
The same polyhedron given in two ways:
sage: p = Polyhedron(ieqs = [[0,1,0,0],[0,0,1,0]]) sage: p.Vrepresentation() (A line in the direction (0, 0, 1), A ray in the direction (1, 0, 0), A ray in the direction (0, 1, 0), A vertex at (0, 0, 0)) sage: q = Polyhedron(vertices=[[0,0,0]], rays=[[1,0,0],[0,1,0]], lines=[[0,0,1]]) sage: q.Hrepresentation() (An inequality (1, 0, 0) x + 0 >= 0, An inequality (0, 1, 0) x + 0 >= 0)
>>> from sage.all import * >>> p = Polyhedron(ieqs = [[Integer(0),Integer(1),Integer(0),Integer(0)],[Integer(0),Integer(0),Integer(1),Integer(0)]]) >>> p.Vrepresentation() (A line in the direction (0, 0, 1), A ray in the direction (1, 0, 0), A ray in the direction (0, 1, 0), A vertex at (0, 0, 0)) >>> q = Polyhedron(vertices=[[Integer(0),Integer(0),Integer(0)]], rays=[[Integer(1),Integer(0),Integer(0)],[Integer(0),Integer(1),Integer(0)]], lines=[[Integer(0),Integer(0),Integer(1)]]) >>> q.Hrepresentation() (An inequality (1, 0, 0) x + 0 >= 0, An inequality (0, 1, 0) x + 0 >= 0)
Finally, a more complicated example. Take \(\mathbb{R}_{\geq 0}^6\) with coordinates \(a, b, \dots, f\) and
The inequality \(e+b \geq c+d\)
The inequality \(e+c \geq b+d\)
The equation \(a+b+c+d+e+f = 31\)
sage: positive_coords = Polyhedron(ieqs=[ ....: [0, 1, 0, 0, 0, 0, 0], [0, 0, 1, 0, 0, 0, 0], [0, 0, 0, 1, 0, 0, 0], ....: [0, 0, 0, 0, 1, 0, 0], [0, 0, 0, 0, 0, 1, 0], [0, 0, 0, 0, 0, 0, 1]]) sage: P = Polyhedron(ieqs=positive_coords.inequalities() + ( ....: [0,0,1,-1,-1,1,0], [0,0,-1,1,-1,1,0]), eqns=[[-31,1,1,1,1,1,1]]) sage: P A 5-dimensional polyhedron in QQ^6 defined as the convex hull of 7 vertices sage: P.dim() 5 sage: P.Vrepresentation() (A vertex at (31, 0, 0, 0, 0, 0), A vertex at (0, 0, 0, 0, 0, 31), A vertex at (0, 0, 0, 0, 31, 0), A vertex at (0, 0, 31/2, 0, 31/2, 0), A vertex at (0, 31/2, 31/2, 0, 0, 0), A vertex at (0, 31/2, 0, 0, 31/2, 0), A vertex at (0, 0, 0, 31/2, 31/2, 0))
>>> from sage.all import * >>> positive_coords = Polyhedron(ieqs=[ ... [Integer(0), Integer(1), Integer(0), Integer(0), Integer(0), Integer(0), Integer(0)], [Integer(0), Integer(0), Integer(1), Integer(0), Integer(0), Integer(0), Integer(0)], [Integer(0), Integer(0), Integer(0), Integer(1), Integer(0), Integer(0), Integer(0)], ... [Integer(0), Integer(0), Integer(0), Integer(0), Integer(1), Integer(0), Integer(0)], [Integer(0), Integer(0), Integer(0), Integer(0), Integer(0), Integer(1), Integer(0)], [Integer(0), Integer(0), Integer(0), Integer(0), Integer(0), Integer(0), Integer(1)]]) >>> P = Polyhedron(ieqs=positive_coords.inequalities() + ( ... [Integer(0),Integer(0),Integer(1),-Integer(1),-Integer(1),Integer(1),Integer(0)], [Integer(0),Integer(0),-Integer(1),Integer(1),-Integer(1),Integer(1),Integer(0)]), eqns=[[-Integer(31),Integer(1),Integer(1),Integer(1),Integer(1),Integer(1),Integer(1)]]) >>> P A 5-dimensional polyhedron in QQ^6 defined as the convex hull of 7 vertices >>> P.dim() 5 >>> P.Vrepresentation() (A vertex at (31, 0, 0, 0, 0, 0), A vertex at (0, 0, 0, 0, 0, 31), A vertex at (0, 0, 0, 0, 31, 0), A vertex at (0, 0, 31/2, 0, 31/2, 0), A vertex at (0, 31/2, 31/2, 0, 0, 0), A vertex at (0, 31/2, 0, 0, 31/2, 0), A vertex at (0, 0, 0, 31/2, 31/2, 0))
Regular icosahedron, centered at \(0\) with edge length \(2\), with vertices given by the cyclic shifts of \((0, \pm 1, \pm (1+\sqrt(5))/2)\), cf. Wikipedia article Regular_icosahedron. It needs a number field:
sage: # needs sage.rings.number_field sage: R0.<r0> = QQ[] sage: R1.<r1> = NumberField(r0^2-5, embedding=AA(5)**(1/2)) sage: gold = (1+r1)/2 sage: v = [[0, 1, gold], [0, 1, -gold], [0, -1, gold], [0, -1, -gold]] sage: pp = Permutation((1, 2, 3)) sage: icosah = Polyhedron( # needs sage.combinat ....: [(pp^2).action(w) for w in v] + [pp.action(w) for w in v] + v, ....: base_ring=R1) sage: len(icosah.faces(2)) # needs sage.combinat 20
>>> from sage.all import * >>> # needs sage.rings.number_field >>> R0 = QQ['r0']; (r0,) = R0._first_ngens(1) >>> R1 = NumberField(r0**Integer(2)-Integer(5), embedding=AA(Integer(5))**(Integer(1)/Integer(2)), names=('r1',)); (r1,) = R1._first_ngens(1) >>> gold = (Integer(1)+r1)/Integer(2) >>> v = [[Integer(0), Integer(1), gold], [Integer(0), Integer(1), -gold], [Integer(0), -Integer(1), gold], [Integer(0), -Integer(1), -gold]] >>> pp = Permutation((Integer(1), Integer(2), Integer(3))) >>> icosah = Polyhedron( # needs sage.combinat ... [(pp**Integer(2)).action(w) for w in v] + [pp.action(w) for w in v] + v, ... base_ring=R1) >>> len(icosah.faces(Integer(2))) # needs sage.combinat 20
When the input contains elements of a Number Field, they require an embedding:
sage: # needs sage.rings.number_field sage: x = polygen(ZZ, 'x') sage: K = NumberField(x^2 - 2,'s') sage: s = K.0 sage: L = NumberField(x^3 - 2,'t') sage: t = L.0 sage: P = Polyhedron(vertices=[[0,s], [t,0]]) Traceback (most recent call last): ... ValueError: invalid base ring
>>> from sage.all import * >>> # needs sage.rings.number_field >>> x = polygen(ZZ, 'x') >>> K = NumberField(x**Integer(2) - Integer(2),'s') >>> s = K.gen(0) >>> L = NumberField(x**Integer(3) - Integer(2),'t') >>> t = L.gen(0) >>> P = Polyhedron(vertices=[[Integer(0),s], [t,Integer(0)]]) Traceback (most recent call last): ... ValueError: invalid base ring
Converting from a given polyhedron:
sage: cb = polytopes.cube(); cb A 3-dimensional polyhedron in ZZ^3 defined as the convex hull of 8 vertices sage: Polyhedron(cb, base_ring=QQ) A 3-dimensional polyhedron in QQ^3 defined as the convex hull of 8 vertices
>>> from sage.all import * >>> cb = polytopes.cube(); cb A 3-dimensional polyhedron in ZZ^3 defined as the convex hull of 8 vertices >>> Polyhedron(cb, base_ring=QQ) A 3-dimensional polyhedron in QQ^3 defined as the convex hull of 8 vertices
Converting from other objects to a polyhedron:
sage: quadrant = Cone([(1,0), (0,1)]) sage: Polyhedron(quadrant) A 2-dimensional polyhedron in ZZ^2 defined as the convex hull of 1 vertex and 2 rays sage: Polyhedron(quadrant, base_ring=QQ) A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 1 vertex and 2 rays sage: o = lattice_polytope.cross_polytope(2) sage: Polyhedron(o) A 2-dimensional polyhedron in ZZ^2 defined as the convex hull of 4 vertices sage: Polyhedron(o, base_ring=QQ) A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 4 vertices sage: p = MixedIntegerLinearProgram(solver='PPL') sage: x, y = p['x'], p['y'] sage: p.add_constraint(x <= 1) sage: p.add_constraint(x >= -1) sage: p.add_constraint(y <= 1) sage: p.add_constraint(y >= -1) sage: Polyhedron(p, base_ring=ZZ) A 2-dimensional polyhedron in ZZ^2 defined as the convex hull of 4 vertices sage: Polyhedron(p) A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 4 vertices sage: # needs sage.combinat sage: H.<x,y> = HyperplaneArrangements(QQ) sage: h = x + y - 1; h Hyperplane x + y - 1 sage: Polyhedron(h, base_ring=ZZ) A 1-dimensional polyhedron in ZZ^2 defined as the convex hull of 1 vertex and 1 line sage: Polyhedron(h) A 1-dimensional polyhedron in QQ^2 defined as the convex hull of 1 vertex and 1 line
>>> from sage.all import * >>> quadrant = Cone([(Integer(1),Integer(0)), (Integer(0),Integer(1))]) >>> Polyhedron(quadrant) A 2-dimensional polyhedron in ZZ^2 defined as the convex hull of 1 vertex and 2 rays >>> Polyhedron(quadrant, base_ring=QQ) A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 1 vertex and 2 rays >>> o = lattice_polytope.cross_polytope(Integer(2)) >>> Polyhedron(o) A 2-dimensional polyhedron in ZZ^2 defined as the convex hull of 4 vertices >>> Polyhedron(o, base_ring=QQ) A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 4 vertices >>> p = MixedIntegerLinearProgram(solver='PPL') >>> x, y = p['x'], p['y'] >>> p.add_constraint(x <= Integer(1)) >>> p.add_constraint(x >= -Integer(1)) >>> p.add_constraint(y <= Integer(1)) >>> p.add_constraint(y >= -Integer(1)) >>> Polyhedron(p, base_ring=ZZ) A 2-dimensional polyhedron in ZZ^2 defined as the convex hull of 4 vertices >>> Polyhedron(p) A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 4 vertices >>> # needs sage.combinat >>> H = HyperplaneArrangements(QQ, names=('x', 'y',)); (x, y,) = H._first_ngens(2) >>> h = x + y - Integer(1); h Hyperplane x + y - 1 >>> Polyhedron(h, base_ring=ZZ) A 1-dimensional polyhedron in ZZ^2 defined as the convex hull of 1 vertex and 1 line >>> Polyhedron(h) A 1-dimensional polyhedron in QQ^2 defined as the convex hull of 1 vertex and 1 line
Note
Once constructed, a
Polyhedron
object is immutable.Although the option
base_ring=RDF
allows numerical data to be used, it might not give the right answer for degenerate input data - the results can depend upon the tolerance setting of cdd.
See also