The polymake backend for polyhedral computations¶
Note
This backend requires polymake.
To install it, type sage -i polymake
in the terminal.
AUTHORS:
Matthias Köppe (2017-03): initial version
- class sage.geometry.polyhedron.backend_polymake.Polyhedron_QQ_polymake(parent, Vrep, Hrep, polymake_polytope=None, **kwds)¶
Bases:
sage.geometry.polyhedron.backend_polymake.Polyhedron_polymake
,sage.geometry.polyhedron.base_QQ.Polyhedron_QQ
Polyhedra over \(\QQ\) with polymake.
INPUT:
Vrep
– a list[vertices, rays, lines]
orNone
Hrep
– a list[ieqs, eqns]
orNone
EXAMPLES:
sage: p = Polyhedron(vertices=[(0,0),(1,0),(0,1)], # optional - polymake ....: rays=[(1,1)], lines=[], ....: backend='polymake', base_ring=QQ) sage: TestSuite(p).run() # optional - polymake
- class sage.geometry.polyhedron.backend_polymake.Polyhedron_ZZ_polymake(parent, Vrep, Hrep, polymake_polytope=None, **kwds)¶
Bases:
sage.geometry.polyhedron.backend_polymake.Polyhedron_polymake
,sage.geometry.polyhedron.base_ZZ.Polyhedron_ZZ
Polyhedra over \(\ZZ\) with polymake.
INPUT:
Vrep
– a list[vertices, rays, lines]
orNone
Hrep
– a list[ieqs, eqns]
orNone
EXAMPLES:
sage: p = Polyhedron(vertices=[(0,0),(1,0),(0,1)], # optional - polymake ....: rays=[(1,1)], lines=[], ....: backend='polymake', base_ring=ZZ) sage: TestSuite(p).run() # optional - polymake
- class sage.geometry.polyhedron.backend_polymake.Polyhedron_polymake(parent, Vrep, Hrep, polymake_polytope=None, **kwds)¶
Bases:
sage.geometry.polyhedron.base.Polyhedron_base
Polyhedra with polymake
INPUT:
parent
–Polyhedra
the parentVrep
– a list[vertices, rays, lines]
orNone
; the V-representation of the polyhedron; ifNone
, the polyhedron is determined by the H-representationHrep
– a list[ieqs, eqns]
orNone
; the H-representation of the polyhedron; ifNone
, the polyhedron is determined by the V-representationpolymake_polytope
– a polymake polytope object
Only one of
Vrep
,Hrep
, orpolymake_polytope
can be different fromNone
.EXAMPLES:
sage: p = Polyhedron(vertices=[(0,0),(1,0),(0,1)], rays=[(1,1)], # optional - polymake ....: lines=[], backend='polymake') sage: TestSuite(p).run() # optional - polymake
A lower-dimensional affine cone; we test that there are no mysterious inequalities coming in from the homogenization:
sage: P = Polyhedron(vertices=[(1, 1)], rays=[(0, 1)], # optional - polymake ....: backend='polymake') sage: P.n_inequalities() # optional - polymake 1 sage: P.equations() # optional - polymake (An equation (1, 0) x - 1 == 0,)
The empty polyhedron:
sage: Polyhedron(eqns=[[1, 0, 0]], backend='polymake') # optional - polymake The empty polyhedron in QQ^2
It can also be obtained differently:
sage: P=Polyhedron(ieqs=[[-2, 1, 1], [-3, -1, -1], [-4, 1, -2]], # optional - polymake ....: backend='polymake') sage: P # optional - polymake The empty polyhedron in QQ^2 sage: P.Vrepresentation() # optional - polymake () sage: P.Hrepresentation() # optional - polymake (An equation -1 == 0,)
The full polyhedron:
sage: Polyhedron(eqns=[[0, 0, 0]], backend='polymake') # optional - polymake A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 1 vertex and 2 lines sage: Polyhedron(ieqs=[[0, 0, 0]], backend='polymake') # optional - polymake A 2-dimensional polyhedron in QQ^2 defined as the convex hull of 1 vertex and 2 lines
Quadratic fields work:
sage: V = polytopes.dodecahedron().vertices_list() sage: Polyhedron(vertices=V, backend='polymake') # optional - polymake A 3-dimensional polyhedron in (Number Field in sqrt5 with defining polynomial x^2 - 5 with sqrt5 = 2.236067977499790?)^3 defined as the convex hull of 20 vertices