Parents for Polyhedra#

sage.geometry.polyhedron.parent.Polyhedra(ambient_space_or_base_ring, ambient_dim, backend=None, ambient_space=None, base_ring=None)#

Construct a suitable parent class for polyhedra

INPUT:

  • base_ring – A ring. Currently there are backends for \(\ZZ\), \(\QQ\), and \(\RDF\).

  • ambient_dim – integer. The ambient space dimension.

  • ambient_space – A free module.

  • backend – string. The name of the backend for computations. There are

    several backends implemented:

    • backend="ppl" uses the Parma Polyhedra Library

    • backend="cdd" uses CDD

    • backend="normaliz" uses normaliz

    • backend="polymake" uses polymake

    • backend="field" a generic Sage implementation

OUTPUT:

A parent class for polyhedra over the given base ring if the backend supports it. If not, the parent base ring can be larger (for example, \(\QQ\) instead of \(\ZZ\)). If there is no implementation at all, a ValueError is raised.

EXAMPLES:

sage: from sage.geometry.polyhedron.parent import Polyhedra
sage: Polyhedra(AA, 3)                                                  # optional - sage.rings.number_field
Polyhedra in AA^3
sage: Polyhedra(ZZ, 3)
Polyhedra in ZZ^3
sage: type(_)
<class 'sage.geometry.polyhedron.parent.Polyhedra_ZZ_ppl_with_category'>
sage: Polyhedra(QQ, 3, backend='cdd')
Polyhedra in QQ^3
sage: type(_)
<class 'sage.geometry.polyhedron.parent.Polyhedra_QQ_cdd_with_category'>

CDD does not support integer polytopes directly:

sage: Polyhedra(ZZ, 3, backend='cdd')
Polyhedra in QQ^3

Using a more general form of the constructor:

sage: V = VectorSpace(QQ, 3)
sage: Polyhedra(V) is Polyhedra(QQ, 3)
True
sage: Polyhedra(V, backend='field') is Polyhedra(QQ, 3, 'field')
True
sage: Polyhedra(backend='field', ambient_space=V) is Polyhedra(QQ, 3, 'field')
True

sage: M = FreeModule(ZZ, 2)
sage: Polyhedra(M, backend='ppl') is Polyhedra(ZZ, 2, 'ppl')
True
class sage.geometry.polyhedron.parent.Polyhedra_QQ_cdd(base_ring, ambient_dim, backend)#

Bases: sage.geometry.polyhedron.parent.Polyhedra_base

Element#

alias of sage.geometry.polyhedron.backend_cdd.Polyhedron_QQ_cdd

class sage.geometry.polyhedron.parent.Polyhedra_QQ_normaliz(base_ring, ambient_dim, backend)#

Bases: sage.geometry.polyhedron.parent.Polyhedra_base

Element#

alias of sage.geometry.polyhedron.backend_normaliz.Polyhedron_QQ_normaliz

class sage.geometry.polyhedron.parent.Polyhedra_QQ_ppl(base_ring, ambient_dim, backend)#

Bases: sage.geometry.polyhedron.parent.Polyhedra_base

Element#

alias of sage.geometry.polyhedron.backend_ppl.Polyhedron_QQ_ppl

class sage.geometry.polyhedron.parent.Polyhedra_RDF_cdd(base_ring, ambient_dim, backend)#

Bases: sage.geometry.polyhedron.parent.Polyhedra_base

Element#

alias of sage.geometry.polyhedron.backend_cdd_rdf.Polyhedron_RDF_cdd

class sage.geometry.polyhedron.parent.Polyhedra_ZZ_normaliz(base_ring, ambient_dim, backend)#

Bases: sage.geometry.polyhedron.parent.Polyhedra_base

Element#

alias of sage.geometry.polyhedron.backend_normaliz.Polyhedron_ZZ_normaliz

class sage.geometry.polyhedron.parent.Polyhedra_ZZ_ppl(base_ring, ambient_dim, backend)#

Bases: sage.geometry.polyhedron.parent.Polyhedra_base

Element#

alias of sage.geometry.polyhedron.backend_ppl.Polyhedron_ZZ_ppl

class sage.geometry.polyhedron.parent.Polyhedra_base(base_ring, ambient_dim, backend)#

Bases: sage.structure.unique_representation.UniqueRepresentation, sage.structure.parent.Parent

Polyhedra in a fixed ambient space.

INPUT:

  • base_ring – either ZZ, QQ, or RDF. The base ring of the ambient module/vector space.

  • ambient_dim – integer. The ambient space dimension.

  • backend – string. The name of the backend for computations. There are

    several backends implemented:

    • backend="ppl" uses the Parma Polyhedra Library

    • backend="cdd" uses CDD

    • backend="normaliz" uses normaliz

    • backend="polymake" uses polymake

    • backend="field" a generic Sage implementation

EXAMPLES:

sage: from sage.geometry.polyhedron.parent import Polyhedra
sage: Polyhedra(ZZ, 3)
Polyhedra in ZZ^3
Hrepresentation_space()#

Return the linear space containing the H-representation vectors.

OUTPUT:

A free module over the base ring of dimension ambient_dim() + 1.

EXAMPLES:

sage: from sage.geometry.polyhedron.parent import Polyhedra
sage: Polyhedra(ZZ, 2).Hrepresentation_space()
Ambient free module of rank 3 over the principal ideal domain Integer Ring
Vrepresentation_space()#

Return the ambient vector space.

This is the vector space or module containing the Vrepresentation vectors.

OUTPUT:

A free module over the base ring of dimension ambient_dim().

EXAMPLES:

sage: from sage.geometry.polyhedron.parent import Polyhedra
sage: Polyhedra(QQ, 4).Vrepresentation_space()
Vector space of dimension 4 over Rational Field
sage: Polyhedra(QQ, 4).ambient_space()
Vector space of dimension 4 over Rational Field
ambient_dim()#

Return the dimension of the ambient space.

EXAMPLES:

sage: from sage.geometry.polyhedron.parent import Polyhedra
sage: Polyhedra(QQ, 3).ambient_dim()
3
ambient_space()#

Return the ambient vector space.

This is the vector space or module containing the Vrepresentation vectors.

OUTPUT:

A free module over the base ring of dimension ambient_dim().

EXAMPLES:

sage: from sage.geometry.polyhedron.parent import Polyhedra
sage: Polyhedra(QQ, 4).Vrepresentation_space()
Vector space of dimension 4 over Rational Field
sage: Polyhedra(QQ, 4).ambient_space()
Vector space of dimension 4 over Rational Field
an_element()#

Return a Polyhedron.

EXAMPLES:

sage: from sage.geometry.polyhedron.parent import Polyhedra
sage: Polyhedra(QQ, 4).an_element()
A 4-dimensional polyhedron in QQ^4 defined as the convex hull of 5 vertices
backend()#

Return the backend.

EXAMPLES:

sage: from sage.geometry.polyhedron.parent import Polyhedra
sage: Polyhedra(QQ, 3).backend()
'ppl'
base_extend(base_ring, backend=None, ambient_dim=None)#

Return the base extended parent.

INPUT:

  • base_ring, backend – see Polyhedron().

  • ambient_dim – if not None change ambient dimension accordingly.

EXAMPLES:

sage: from sage.geometry.polyhedron.parent import Polyhedra
sage: Polyhedra(ZZ,3).base_extend(QQ)
Polyhedra in QQ^3
sage: Polyhedra(ZZ,3).an_element().base_extend(QQ)
A 3-dimensional polyhedron in QQ^3 defined as the convex hull of 4 vertices
sage: Polyhedra(QQ, 2).base_extend(ZZ)
Polyhedra in QQ^2
change_ring(base_ring, backend=None, ambient_dim=None)#

Return the parent with the new base ring.

INPUT:

  • base_ring, backend – see Polyhedron().

  • ambient_dim – if not None change ambient dimension accordingly.

EXAMPLES:

sage: from sage.geometry.polyhedron.parent import Polyhedra
sage: Polyhedra(ZZ,3).change_ring(QQ)
Polyhedra in QQ^3
sage: Polyhedra(ZZ,3).an_element().change_ring(QQ)
A 3-dimensional polyhedron in QQ^3 defined as the convex hull of 4 vertices

sage: Polyhedra(RDF, 3).change_ring(QQ).backend()
'cdd'
sage: Polyhedra(QQ, 3).change_ring(ZZ, ambient_dim=4)
Polyhedra in ZZ^4
sage: Polyhedra(QQ, 3, backend='cdd').change_ring(QQ, ambient_dim=4).backend()
'cdd'
empty()#

Return the empty polyhedron.

EXAMPLES:

sage: from sage.geometry.polyhedron.parent import Polyhedra
sage: P = Polyhedra(QQ, 4)
sage: P.empty()
The empty polyhedron in QQ^4
sage: P.empty().is_empty()
True
list()#

Return the two polyhedra in ambient dimension 0, raise an error otherwise

EXAMPLES:

sage: from sage.geometry.polyhedron.parent import Polyhedra
sage: P = Polyhedra(QQ, 3)
sage: P.cardinality()
+Infinity

sage: P = Polyhedra(AA, 0)                                          # optional - sage.rings.number_field
sage: P.category()                                                  # optional - sage.rings.number_field
Category of finite enumerated polyhedral sets over Algebraic Real Field
sage: P.list()                                                      # optional - sage.rings.number_field
[The empty polyhedron in AA^0,
 A 0-dimensional polyhedron in AA^0 defined as the convex hull of 1 vertex]
sage: P.cardinality()                                               # optional - sage.rings.number_field
2
recycle(polyhedron)#

Recycle the H/V-representation objects of a polyhedron.

This speeds up creation of new polyhedra by reusing objects. After recycling a polyhedron object, it is not in a consistent state any more and neither the polyhedron nor its H/V-representation objects may be used any more.

INPUT:

  • polyhedron – a polyhedron whose parent is self.

EXAMPLES:

sage: p = Polyhedron([(0,0),(1,0),(0,1)])
sage: p.parent().recycle(p)
some_elements()#

Return a list of some elements of the semigroup.

EXAMPLES:

sage: from sage.geometry.polyhedron.parent import Polyhedra
sage: Polyhedra(QQ, 4).some_elements()
[A 3-dimensional polyhedron in QQ^4 defined as the convex hull of 4 vertices,
 A 4-dimensional polyhedron in QQ^4 defined as the convex hull of 1 vertex and 4 rays,
 A 2-dimensional polyhedron in QQ^4 defined as the convex hull of 2 vertices and 1 ray,
 The empty polyhedron in QQ^4]
sage: Polyhedra(ZZ,0).some_elements()
[The empty polyhedron in ZZ^0,
 A 0-dimensional polyhedron in ZZ^0 defined as the convex hull of 1 vertex]
universe()#

Return the entire ambient space as polyhedron.

EXAMPLES:

sage: from sage.geometry.polyhedron.parent import Polyhedra
sage: P = Polyhedra(QQ, 4)
sage: P.universe()
A 4-dimensional polyhedron in QQ^4 defined as the convex hull of 1 vertex and 4 lines
sage: P.universe().is_universe()
True
zero()#

Return the polyhedron consisting of the origin, which is the neutral element for Minkowski addition.

EXAMPLES:

sage: from sage.geometry.polyhedron.parent import Polyhedra
sage: p = Polyhedra(QQ, 4).zero();  p
A 0-dimensional polyhedron in QQ^4 defined as the convex hull of 1 vertex
sage: p+p == p
True
class sage.geometry.polyhedron.parent.Polyhedra_field(base_ring, ambient_dim, backend)#

Bases: sage.geometry.polyhedron.parent.Polyhedra_base

Element#

alias of sage.geometry.polyhedron.backend_field.Polyhedron_field

class sage.geometry.polyhedron.parent.Polyhedra_normaliz(base_ring, ambient_dim, backend)#

Bases: sage.geometry.polyhedron.parent.Polyhedra_base

Element#

alias of sage.geometry.polyhedron.backend_normaliz.Polyhedron_normaliz

class sage.geometry.polyhedron.parent.Polyhedra_number_field(base_ring, ambient_dim, backend)#

Bases: sage.geometry.polyhedron.parent.Polyhedra_base

Element#

alias of sage.geometry.polyhedron.backend_number_field.Polyhedron_number_field

class sage.geometry.polyhedron.parent.Polyhedra_polymake(base_ring, ambient_dim, backend)#

Bases: sage.geometry.polyhedron.parent.Polyhedra_base

Element#

alias of sage.geometry.polyhedron.backend_polymake.Polyhedron_polymake

sage.geometry.polyhedron.parent.does_backend_handle_base_ring(base_ring, backend)#

Return true, if backend can handle base_ring.

EXAMPLES:

sage: from sage.geometry.polyhedron.parent import does_backend_handle_base_ring
sage: does_backend_handle_base_ring(QQ, 'ppl')
True
sage: does_backend_handle_base_ring(QQ[sqrt(5)], 'ppl')                 # optional - sage.rings.number_field
False
sage: does_backend_handle_base_ring(QQ[sqrt(5)], 'field')               # optional - sage.rings.number_field
True