Base class for polyhedra: Methods for triangulation and volume computation

class sage.geometry.polyhedron.base7.Polyhedron_base7(parent, Vrep, Hrep, Vrep_minimal=None, Hrep_minimal=None, pref_rep=None, mutable=False, **kwds)[source]

Bases: Polyhedron_base6

Methods related to triangulation and volume.

centroid(engine='auto', **kwds)[source]

Return the center of the mass of the polytope.

The mass is taken with respect to the induced Lebesgue measure, see volume().

If the polyhedron is not compact, a NotImplementedError is raised.

INPUT:

  • engine – either ‘auto’ (default), ‘internal’, ‘TOPCOM’, or ‘normaliz’. The ‘internal’ and ‘TOPCOM’ instruct this package to always use its own triangulation algorithms or TOPCOM’s algorithms, respectively. By default (‘auto’), TOPCOM is used if it is available and internal routines otherwise.

  • **kwds – keyword arguments that are passed to the triangulation engine (see triangulate())

OUTPUT: the centroid as vector

ALGORITHM:

We triangulate the polytope and find the barycenter of the simplices. We add the individual barycenters weighted by the fraction of the total mass.

EXAMPLES:

sage: P = polytopes.hypercube(2).pyramid()
sage: P.centroid()
(1/4, 0, 0)

sage: P = polytopes.associahedron(['A', 2])                                 # needs sage.combinat
sage: P.centroid()                                                          # needs sage.combinat
(2/21, 2/21)

sage: P = polytopes.permutahedron(4, backend='normaliz')    # optional - pynormaliz
sage: P.centroid()                                          # optional - pynormaliz
(5/2, 5/2, 5/2, 5/2)
>>> from sage.all import *
>>> P = polytopes.hypercube(Integer(2)).pyramid()
>>> P.centroid()
(1/4, 0, 0)

>>> P = polytopes.associahedron(['A', Integer(2)])                                 # needs sage.combinat
>>> P.centroid()                                                          # needs sage.combinat
(2/21, 2/21)

>>> P = polytopes.permutahedron(Integer(4), backend='normaliz')    # optional - pynormaliz
>>> P.centroid()                                          # optional - pynormaliz
(5/2, 5/2, 5/2, 5/2)
P = polytopes.hypercube(2).pyramid()
P.centroid()
P = polytopes.associahedron(['A', 2])                                 # needs sage.combinat
P.centroid()                                                          # needs sage.combinat
P = polytopes.permutahedron(4, backend='normaliz')    # optional - pynormaliz
P.centroid()                                          # optional - pynormaliz

The method is not implemented for unbounded polyhedra:

sage: P = Polyhedron(vertices=[(0, 0)], rays=[(1, 0), (0, 1)])
sage: P.centroid()
Traceback (most recent call last):
...
NotImplementedError: the polyhedron is not compact
>>> from sage.all import *
>>> P = Polyhedron(vertices=[(Integer(0), Integer(0))], rays=[(Integer(1), Integer(0)), (Integer(0), Integer(1))])
>>> P.centroid()
Traceback (most recent call last):
...
NotImplementedError: the polyhedron is not compact
P = Polyhedron(vertices=[(0, 0)], rays=[(1, 0), (0, 1)])
P.centroid()

The centroid of an empty polyhedron is not defined:

sage: Polyhedron().centroid()
Traceback (most recent call last):
...
ZeroDivisionError: rational division by zero
>>> from sage.all import *
>>> Polyhedron().centroid()
Traceback (most recent call last):
...
ZeroDivisionError: rational division by zero
Polyhedron().centroid()
integrate(function, measure='ambient', **kwds)[source]

Return the integral of function over this polytope.

INPUT:

  • self – Polyhedron

  • function – a multivariate polynomial or a valid LattE description string for polynomials

  • measure – string, the measure to use

    Allowed values are:

    • ambient (default): Lebesgue measure of ambient space,

    • induced: Lebesgue measure of the affine hull,

    • induced_nonnormalized: Lebesgue measure of the affine hull without the normalization by \(\sqrt{\det(A^\top A)}\) (with \(A\) being the affine transformation matrix; see affine_hull()).

  • **kwds – additional keyword arguments that are passed to the engine

OUTPUT: the integral of the polynomial over the polytope

Note

The polytope triangulation algorithm is used. This function depends on LattE (i.e., the latte_int optional package).

EXAMPLES:

sage: P = polytopes.cube()
sage: x, y, z = polygens(QQ, 'x, y, z')
sage: P.integrate(x^2*y^2*z^2)                              # optional - latte_int
8/27
>>> from sage.all import *
>>> P = polytopes.cube()
>>> x, y, z = polygens(QQ, 'x, y, z')
>>> P.integrate(x**Integer(2)*y**Integer(2)*z**Integer(2))                              # optional - latte_int
8/27
P = polytopes.cube()
x, y, z = polygens(QQ, 'x, y, z')
P.integrate(x^2*y^2*z^2)                              # optional - latte_int

If the polyhedron has floating point coordinates, an inexact result can be obtained if we transform to rational coordinates:

sage: P = 1.4142*polytopes.cube()
sage: P_QQ = Polyhedron(vertices=[[QQ(vi) for vi in v] for v in P.vertex_generator()])
sage: RDF(P_QQ.integrate(x^2*y^2*z^2))                      # optional - latte_int
6.703841212195228
>>> from sage.all import *
>>> P = RealNumber('1.4142')*polytopes.cube()
>>> P_QQ = Polyhedron(vertices=[[QQ(vi) for vi in v] for v in P.vertex_generator()])
>>> RDF(P_QQ.integrate(x**Integer(2)*y**Integer(2)*z**Integer(2)))                      # optional - latte_int
6.703841212195228
P = 1.4142*polytopes.cube()
P_QQ = Polyhedron(vertices=[[QQ(vi) for vi in v] for v in P.vertex_generator()])
RDF(P_QQ.integrate(x^2*y^2*z^2))                      # optional - latte_int

Integral over a non full-dimensional polytope:

sage: x, y = polygens(QQ, 'x, y')
sage: P = Polyhedron(vertices=[[0,0], [1,1]])
sage: P.integrate(x*y)
0
sage: ixy = P.integrate(x*y, measure='induced'); ixy        # optional - latte_int
0.4714045207910317?
sage: ixy.parent()                                          # optional - latte_int
Algebraic Real Field
>>> from sage.all import *
>>> x, y = polygens(QQ, 'x, y')
>>> P = Polyhedron(vertices=[[Integer(0),Integer(0)], [Integer(1),Integer(1)]])
>>> P.integrate(x*y)
0
>>> ixy = P.integrate(x*y, measure='induced'); ixy        # optional - latte_int
0.4714045207910317?
>>> ixy.parent()                                          # optional - latte_int
Algebraic Real Field
x, y = polygens(QQ, 'x, y')
P = Polyhedron(vertices=[[0,0], [1,1]])
P.integrate(x*y)
ixy = P.integrate(x*y, measure='induced'); ixy        # optional - latte_int
ixy.parent()                                          # optional - latte_int

Convert to a symbolic expression:

sage: ixy.radical_expression()                              # optional - latte_int
1/3*sqrt(2)
>>> from sage.all import *
>>> ixy.radical_expression()                              # optional - latte_int
1/3*sqrt(2)
ixy.radical_expression()                              # optional - latte_int

Another non full-dimensional polytope integration:

sage: R.<x, y, z> = QQ[]
sage: P = polytopes.simplex(2)
sage: V = AA(P.volume(measure='induced'))                                   # needs sage.rings.number_field
sage: V.radical_expression()                                                # needs sage.rings.number_field sage.symbolic
1/2*sqrt(3)
sage: P.integrate(R(1), measure='induced') == V             # optional - latte_int, needs sage.rings.number_field sage.symbolic
True
>>> from sage.all import *
>>> R = QQ['x, y, z']; (x, y, z,) = R._first_ngens(3)
>>> P = polytopes.simplex(Integer(2))
>>> V = AA(P.volume(measure='induced'))                                   # needs sage.rings.number_field
>>> V.radical_expression()                                                # needs sage.rings.number_field sage.symbolic
1/2*sqrt(3)
>>> P.integrate(R(Integer(1)), measure='induced') == V             # optional - latte_int, needs sage.rings.number_field sage.symbolic
True
R.<x, y, z> = QQ[]
P = polytopes.simplex(2)
V = AA(P.volume(measure='induced'))                                   # needs sage.rings.number_field
V.radical_expression()                                                # needs sage.rings.number_field sage.symbolic
P.integrate(R(1), measure='induced') == V             # optional - latte_int, needs sage.rings.number_field sage.symbolic

Computing the mass center:

sage: (P.integrate(x, measure='induced')                    # optional - latte_int, needs sage.rings.number_field sage.symbolic
....:     / V).radical_expression()
1/3
sage: (P.integrate(y, measure='induced')                    # optional - latte_int, needs sage.rings.number_field sage.symbolic
....:     / V).radical_expression()
1/3
sage: (P.integrate(z, measure='induced')                    # optional - latte_int, needs sage.rings.number_field sage.symbolic
....:     / V).radical_expression()
1/3
>>> from sage.all import *
>>> (P.integrate(x, measure='induced')                    # optional - latte_int, needs sage.rings.number_field sage.symbolic
...     / V).radical_expression()
1/3
>>> (P.integrate(y, measure='induced')                    # optional - latte_int, needs sage.rings.number_field sage.symbolic
...     / V).radical_expression()
1/3
>>> (P.integrate(z, measure='induced')                    # optional - latte_int, needs sage.rings.number_field sage.symbolic
...     / V).radical_expression()
1/3
(P.integrate(x, measure='induced')                    # optional - latte_int, needs sage.rings.number_field sage.symbolic
    / V).radical_expression()
(P.integrate(y, measure='induced')                    # optional - latte_int, needs sage.rings.number_field sage.symbolic
    / V).radical_expression()
(P.integrate(z, measure='induced')                    # optional - latte_int, needs sage.rings.number_field sage.symbolic
    / V).radical_expression()
triangulate(engine='auto', connected=True, fine=False, regular=None, star=None)[source]

Return a triangulation of the polytope.

INPUT:

  • engine – either ‘auto’ (default), ‘internal’, ‘TOPCOM’, or ‘normaliz’. The ‘internal’ and ‘TOPCOM’ instruct this package to always use its own triangulation algorithms or TOPCOM’s algorithms, respectively. By default (‘auto’), TOPCOM is used if it is available and internal routines otherwise.

The remaining keyword parameters are passed through to the PointConfiguration constructor:

  • connected – boolean (default: True); whether the triangulations should be connected to the regular triangulations via bistellar flips. These are much easier to compute than all triangulations.

  • fine – boolean (default: False); whether the triangulations must be fine, that is, make use of all points of the configuration

  • regular – boolean or None (default: None); whether the triangulations must be regular. A regular triangulation is one that is induced by a piecewise-linear convex support function. In other words, the shadows of the faces of a polyhedron in one higher dimension.

    • True: Only regular triangulations.

    • False: Only non-regular triangulations.

    • None (default): Both kinds of triangulation.

  • star – either None (default) or a point. Whether the triangulations must be star. A triangulation is star if all maximal simplices contain a common point. The central point can be specified by its index (an integer) in the given points or by its coordinates (anything iterable.)

OUTPUT:

A triangulation of the convex hull of the vertices as a Triangulation. The indices in the triangulation correspond to the Vrepresentation() objects.

EXAMPLES:

sage: cube = polytopes.hypercube(3)
sage: triangulation = cube.triangulate(
....:    engine='internal')  # to make doctest independent of TOPCOM
sage: triangulation
(<0,1,2,7>, <0,1,5,7>, <0,2,3,7>, <0,3,4,7>, <0,4,5,7>, <1,5,6,7>)
sage: simplex_indices = triangulation[0]; simplex_indices
(0, 1, 2, 7)
sage: simplex_vertices = [cube.Vrepresentation(i) for i in simplex_indices]
sage: simplex_vertices
[A vertex at (1, -1, -1),
 A vertex at (1, 1, -1),
 A vertex at (1, 1, 1),
 A vertex at (-1, 1, 1)]
sage: Polyhedron(simplex_vertices)
A 3-dimensional polyhedron in ZZ^3 defined as the convex hull of 4 vertices
>>> from sage.all import *
>>> cube = polytopes.hypercube(Integer(3))
>>> triangulation = cube.triangulate(
...    engine='internal')  # to make doctest independent of TOPCOM
>>> triangulation
(<0,1,2,7>, <0,1,5,7>, <0,2,3,7>, <0,3,4,7>, <0,4,5,7>, <1,5,6,7>)
>>> simplex_indices = triangulation[Integer(0)]; simplex_indices
(0, 1, 2, 7)
>>> simplex_vertices = [cube.Vrepresentation(i) for i in simplex_indices]
>>> simplex_vertices
[A vertex at (1, -1, -1),
 A vertex at (1, 1, -1),
 A vertex at (1, 1, 1),
 A vertex at (-1, 1, 1)]
>>> Polyhedron(simplex_vertices)
A 3-dimensional polyhedron in ZZ^3 defined as the convex hull of 4 vertices
cube = polytopes.hypercube(3)
triangulation = cube.triangulate(
   engine='internal')  # to make doctest independent of TOPCOM
triangulation
simplex_indices = triangulation[0]; simplex_indices
simplex_vertices = [cube.Vrepresentation(i) for i in simplex_indices]
simplex_vertices
Polyhedron(simplex_vertices)

It is possible to use 'normaliz' as an engine. For this, the polyhedron should have the backend set to normaliz:

sage: P = Polyhedron(vertices=[[0,0,1], [1,0,1],            # optional - pynormaliz
....:                          [0,1,1], [1,1,1]],
....:                backend='normaliz')
sage: P.triangulate(engine='normaliz')                      # optional - pynormaliz
(<0,1,2>, <1,2,3>)

sage: P = Polyhedron(vertices=[[0,0,1], [1,0,1],
....:                          [0,1,1], [1,1,1]])
sage: P.triangulate(engine='normaliz')
Traceback (most recent call last):
...
TypeError: the polyhedron's backend should be 'normaliz'
>>> from sage.all import *
>>> P = Polyhedron(vertices=[[Integer(0),Integer(0),Integer(1)], [Integer(1),Integer(0),Integer(1)],            # optional - pynormaliz
...                          [Integer(0),Integer(1),Integer(1)], [Integer(1),Integer(1),Integer(1)]],
...                backend='normaliz')
>>> P.triangulate(engine='normaliz')                      # optional - pynormaliz
(<0,1,2>, <1,2,3>)

>>> P = Polyhedron(vertices=[[Integer(0),Integer(0),Integer(1)], [Integer(1),Integer(0),Integer(1)],
...                          [Integer(0),Integer(1),Integer(1)], [Integer(1),Integer(1),Integer(1)]])
>>> P.triangulate(engine='normaliz')
Traceback (most recent call last):
...
TypeError: the polyhedron's backend should be 'normaliz'
P = Polyhedron(vertices=[[0,0,1], [1,0,1],            # optional - pynormaliz
                         [0,1,1], [1,1,1]],
               backend='normaliz')
P.triangulate(engine='normaliz')                      # optional - pynormaliz
P = Polyhedron(vertices=[[0,0,1], [1,0,1],
                         [0,1,1], [1,1,1]])
P.triangulate(engine='normaliz')

The normaliz engine can triangulate pointed cones:

sage: # optional - pynormaliz
sage: C1 = Polyhedron(rays=[[0,0,1], [1,0,1],
....:                       [0,1,1], [1,1,1]],
....:                 backend='normaliz')
sage: C1.triangulate(engine='normaliz')
(<0,1,2>, <1,2,3>)
sage: C2 = Polyhedron(rays=[[1,0,1], [0,0,1],
....:                       [0,1,1], [1,1,10/9]],
....:                 backend='normaliz')
sage: C2.triangulate(engine='normaliz')
(<0,1,2>, <1,2,3>)
>>> from sage.all import *
>>> # optional - pynormaliz
>>> C1 = Polyhedron(rays=[[Integer(0),Integer(0),Integer(1)], [Integer(1),Integer(0),Integer(1)],
...                       [Integer(0),Integer(1),Integer(1)], [Integer(1),Integer(1),Integer(1)]],
...                 backend='normaliz')
>>> C1.triangulate(engine='normaliz')
(<0,1,2>, <1,2,3>)
>>> C2 = Polyhedron(rays=[[Integer(1),Integer(0),Integer(1)], [Integer(0),Integer(0),Integer(1)],
...                       [Integer(0),Integer(1),Integer(1)], [Integer(1),Integer(1),Integer(10)/Integer(9)]],
...                 backend='normaliz')
>>> C2.triangulate(engine='normaliz')
(<0,1,2>, <1,2,3>)
# optional - pynormaliz
C1 = Polyhedron(rays=[[0,0,1], [1,0,1],
                      [0,1,1], [1,1,1]],
                backend='normaliz')
C1.triangulate(engine='normaliz')
C2 = Polyhedron(rays=[[1,0,1], [0,0,1],
                      [0,1,1], [1,1,10/9]],
                backend='normaliz')
C2.triangulate(engine='normaliz')

They can also be affine cones:

sage: K = Polyhedron(vertices=[[1,1,1]],                    # optional - pynormaliz
....:                rays=[[1,0,0], [0,1,0], [1,1,-1], [1,1,1]],
....:                backend='normaliz')
sage: K.triangulate(engine='normaliz')                      # optional - pynormaliz
(<0,1,2>, <0,1,3>)
>>> from sage.all import *
>>> K = Polyhedron(vertices=[[Integer(1),Integer(1),Integer(1)]],                    # optional - pynormaliz
...                rays=[[Integer(1),Integer(0),Integer(0)], [Integer(0),Integer(1),Integer(0)], [Integer(1),Integer(1),-Integer(1)], [Integer(1),Integer(1),Integer(1)]],
...                backend='normaliz')
>>> K.triangulate(engine='normaliz')                      # optional - pynormaliz
(<0,1,2>, <0,1,3>)
K = Polyhedron(vertices=[[1,1,1]],                    # optional - pynormaliz
               rays=[[1,0,0], [0,1,0], [1,1,-1], [1,1,1]],
               backend='normaliz')
K.triangulate(engine='normaliz')                      # optional - pynormaliz
volume(measure='ambient', engine='auto', **kwds)[source]

Return the volume of the polytope.

INPUT:

  • measure – string. The measure to use. Allowed values are:

    • ambient (default): Lebesgue measure of ambient space (volume)

    • induced: Lebesgue measure of the affine hull (relative volume)

    • induced_rational: Scaling of the Lebesgue measure for rational polytopes, such that the unit hypercube has volume 1

    • induced_lattice: Scaling of the Lebesgue measure, such that the volume of the hypercube is factorial(n)

  • engine – string. The backend to use. Allowed values are:

    • 'auto' (default): choose engine according to measure

    • 'internal': see triangulate()

    • 'TOPCOM': see triangulate()

    • 'lrs': use David Avis’s lrs program (optional)

    • 'latte': use LattE integrale program (optional)

    • 'normaliz': use Normaliz program (optional)

  • **kwds – keyword arguments that are passed to the triangulation engine

OUTPUT: the volume of the polytope

EXAMPLES:

sage: polytopes.hypercube(3).volume()
8
sage: (polytopes.hypercube(3)*2).volume()
64
sage: polytopes.twenty_four_cell().volume()
2
>>> from sage.all import *
>>> polytopes.hypercube(Integer(3)).volume()
8
>>> (polytopes.hypercube(Integer(3))*Integer(2)).volume()
64
>>> polytopes.twenty_four_cell().volume()
2
polytopes.hypercube(3).volume()
(polytopes.hypercube(3)*2).volume()
polytopes.twenty_four_cell().volume()

Volume of the same polytopes, using the optional package lrslib (which requires a rational polytope):

sage: I3 = polytopes.hypercube(3)
sage: I3.volume(engine='lrs')                               # optional - lrslib
8
sage: C24 = polytopes.twenty_four_cell()
sage: C24.volume(engine='lrs')                              # optional - lrslib
2
>>> from sage.all import *
>>> I3 = polytopes.hypercube(Integer(3))
>>> I3.volume(engine='lrs')                               # optional - lrslib
8
>>> C24 = polytopes.twenty_four_cell()
>>> C24.volume(engine='lrs')                              # optional - lrslib
2
I3 = polytopes.hypercube(3)
I3.volume(engine='lrs')                               # optional - lrslib
C24 = polytopes.twenty_four_cell()
C24.volume(engine='lrs')                              # optional - lrslib

If the base ring is exact, the answer is exact:

sage: P5 = polytopes.regular_polygon(5)                                     # needs sage.rings.number_field
sage: P5.volume()                                                           # needs sage.rings.number_field
2.377641290737884?

sage: polytopes.icosahedron().volume()                                      # needs sage.groups sage.rings.number_field
5/12*sqrt5 + 5/4
sage: numerical_approx(_)  # abs tol 1e9                                    # needs sage.groups sage.rings.number_field
2.18169499062491
>>> from sage.all import *
>>> P5 = polytopes.regular_polygon(Integer(5))                                     # needs sage.rings.number_field
>>> P5.volume()                                                           # needs sage.rings.number_field
2.377641290737884?

>>> polytopes.icosahedron().volume()                                      # needs sage.groups sage.rings.number_field
5/12*sqrt5 + 5/4
>>> numerical_approx(_)  # abs tol 1e9                                    # needs sage.groups sage.rings.number_field
2.18169499062491
P5 = polytopes.regular_polygon(5)                                     # needs sage.rings.number_field
P5.volume()                                                           # needs sage.rings.number_field
polytopes.icosahedron().volume()                                      # needs sage.groups sage.rings.number_field
numerical_approx(_)  # abs tol 1e9                                    # needs sage.groups sage.rings.number_field

When considering lower-dimensional polytopes, we can ask for the ambient (full-dimensional), the induced measure (of the affine hull) or, in the case of lattice polytopes, for the induced rational measure. This is controlled by the parameter measure. Different engines may have different ideas on the definition of volume of a lower-dimensional object:

sage: P = Polyhedron([[0, 0], [1, 1]])
sage: P.volume()
0
sage: P.volume(measure='induced')                                           # needs sage.rings.number_field
1.414213562373095?
sage: P.volume(measure='induced_rational')                  # optional - latte_int
1

sage: # needs sage.rings.number_field
sage: S = polytopes.regular_polygon(6); S
A 2-dimensional polyhedron in AA^2 defined as the convex hull of 6 vertices
sage: edge = S.faces(1)[4].as_polyhedron()
sage: edge.vertices()
(A vertex at (0.866025403784439?, 1/2), A vertex at (0, 1))
sage: edge.volume()
0
sage: edge.volume(measure='induced')
1

sage: # optional - pynormaliz
sage: P = Polyhedron(backend='normaliz',
....:                vertices=[[1,0,0], [0,0,1],
....:                          [-1,1,1], [-1,2,0]])
sage: P.volume()
0
sage: P.volume(measure='induced')                                           # needs sage.rings.number_field
2.598076211353316?
sage: P.volume(measure='induced', engine='normaliz')
2.598076211353316
sage: P.volume(measure='induced_rational')                  # optional - latte_int
3/2
sage: P.volume(measure='induced_rational',
....:          engine='normaliz')
3/2
sage: P.volume(measure='induced_lattice')
3
>>> from sage.all import *
>>> P = Polyhedron([[Integer(0), Integer(0)], [Integer(1), Integer(1)]])
>>> P.volume()
0
>>> P.volume(measure='induced')                                           # needs sage.rings.number_field
1.414213562373095?
>>> P.volume(measure='induced_rational')                  # optional - latte_int
1

>>> # needs sage.rings.number_field
>>> S = polytopes.regular_polygon(Integer(6)); S
A 2-dimensional polyhedron in AA^2 defined as the convex hull of 6 vertices
>>> edge = S.faces(Integer(1))[Integer(4)].as_polyhedron()
>>> edge.vertices()
(A vertex at (0.866025403784439?, 1/2), A vertex at (0, 1))
>>> edge.volume()
0
>>> edge.volume(measure='induced')
1

>>> # optional - pynormaliz
>>> P = Polyhedron(backend='normaliz',
...                vertices=[[Integer(1),Integer(0),Integer(0)], [Integer(0),Integer(0),Integer(1)],
...                          [-Integer(1),Integer(1),Integer(1)], [-Integer(1),Integer(2),Integer(0)]])
>>> P.volume()
0
>>> P.volume(measure='induced')                                           # needs sage.rings.number_field
2.598076211353316?
>>> P.volume(measure='induced', engine='normaliz')
2.598076211353316
>>> P.volume(measure='induced_rational')                  # optional - latte_int
3/2
>>> P.volume(measure='induced_rational',
...          engine='normaliz')
3/2
>>> P.volume(measure='induced_lattice')
3
P = Polyhedron([[0, 0], [1, 1]])
P.volume()
P.volume(measure='induced')                                           # needs sage.rings.number_field
P.volume(measure='induced_rational')                  # optional - latte_int
# needs sage.rings.number_field
S = polytopes.regular_polygon(6); S
edge = S.faces(1)[4].as_polyhedron()
edge.vertices()
edge.volume()
edge.volume(measure='induced')
# optional - pynormaliz
P = Polyhedron(backend='normaliz',
               vertices=[[1,0,0], [0,0,1],
                         [-1,1,1], [-1,2,0]])
P.volume()
P.volume(measure='induced')                                           # needs sage.rings.number_field
P.volume(measure='induced', engine='normaliz')
P.volume(measure='induced_rational')                  # optional - latte_int
P.volume(measure='induced_rational',
         engine='normaliz')
P.volume(measure='induced_lattice')

The same polytope without normaliz backend:

sage: P = Polyhedron(vertices=[[1,0,0], [0,0,1], [-1,1,1], [-1,2,0]])
sage: P.volume(measure='induced_lattice', engine='latte')   # optional - latte_int
3

sage: # needs sage.groups sage.rings.number_field
sage: Dexact = polytopes.dodecahedron()
sage: F0 = Dexact.faces(2)[0].as_polyhedron()
sage: v = F0.volume(measure='induced', engine='internal'); v
1.53406271079097?
sage: F4 = Dexact.faces(2)[4].as_polyhedron()
sage: v = F4.volume(measure='induced', engine='internal'); v
1.53406271079097?
sage: RDF(v)    # abs tol 1e-9
1.53406271079044

sage: # needs sage.groups
sage: Dinexact = polytopes.dodecahedron(exact=False)
sage: F2 = Dinexact.faces(2)[2].as_polyhedron()
sage: w = F2.volume(measure='induced', engine='internal')
sage: RDF(w)    # abs tol 1e-9
1.5340627082974878

sage: all(polytopes.simplex(d).volume(measure='induced')                    # needs sage.rings.number_field sage.symbolic
....:        == sqrt(d+1)/factorial(d)
....:     for d in range(1,5))
True

sage: I = Polyhedron([[-3, 0], [0, 9]])
sage: I.volume(measure='induced')                                           # needs sage.rings.number_field
9.48683298050514?
sage: I.volume(measure='induced_rational')                  # optional - latte_int
3

sage: T = Polyhedron([[3, 0, 0], [0, 4, 0], [0, 0, 5]])
sage: T.volume(measure='induced')                                           # needs sage.rings.number_field
13.86542462386205?
sage: T.volume(measure='induced_rational')                  # optional - latte_int
1/2

sage: Q = Polyhedron(vertices=[(0, 0, 1, 1), (0, 1, 1, 0), (1, 1, 0, 0)])
sage: Q.volume(measure='induced')
1
sage: Q.volume(measure='induced_rational')                  # optional - latte_int
1/2
>>> from sage.all import *
>>> P = Polyhedron(vertices=[[Integer(1),Integer(0),Integer(0)], [Integer(0),Integer(0),Integer(1)], [-Integer(1),Integer(1),Integer(1)], [-Integer(1),Integer(2),Integer(0)]])
>>> P.volume(measure='induced_lattice', engine='latte')   # optional - latte_int
3

>>> # needs sage.groups sage.rings.number_field
>>> Dexact = polytopes.dodecahedron()
>>> F0 = Dexact.faces(Integer(2))[Integer(0)].as_polyhedron()
>>> v = F0.volume(measure='induced', engine='internal'); v
1.53406271079097?
>>> F4 = Dexact.faces(Integer(2))[Integer(4)].as_polyhedron()
>>> v = F4.volume(measure='induced', engine='internal'); v
1.53406271079097?
>>> RDF(v)    # abs tol 1e-9
1.53406271079044

>>> # needs sage.groups
>>> Dinexact = polytopes.dodecahedron(exact=False)
>>> F2 = Dinexact.faces(Integer(2))[Integer(2)].as_polyhedron()
>>> w = F2.volume(measure='induced', engine='internal')
>>> RDF(w)    # abs tol 1e-9
1.5340627082974878

>>> all(polytopes.simplex(d).volume(measure='induced')                    # needs sage.rings.number_field sage.symbolic
...        == sqrt(d+Integer(1))/factorial(d)
...     for d in range(Integer(1),Integer(5)))
True

>>> I = Polyhedron([[-Integer(3), Integer(0)], [Integer(0), Integer(9)]])
>>> I.volume(measure='induced')                                           # needs sage.rings.number_field
9.48683298050514?
>>> I.volume(measure='induced_rational')                  # optional - latte_int
3

>>> T = Polyhedron([[Integer(3), Integer(0), Integer(0)], [Integer(0), Integer(4), Integer(0)], [Integer(0), Integer(0), Integer(5)]])
>>> T.volume(measure='induced')                                           # needs sage.rings.number_field
13.86542462386205?
>>> T.volume(measure='induced_rational')                  # optional - latte_int
1/2

>>> Q = Polyhedron(vertices=[(Integer(0), Integer(0), Integer(1), Integer(1)), (Integer(0), Integer(1), Integer(1), Integer(0)), (Integer(1), Integer(1), Integer(0), Integer(0))])
>>> Q.volume(measure='induced')
1
>>> Q.volume(measure='induced_rational')                  # optional - latte_int
1/2
P = Polyhedron(vertices=[[1,0,0], [0,0,1], [-1,1,1], [-1,2,0]])
P.volume(measure='induced_lattice', engine='latte')   # optional - latte_int
# needs sage.groups sage.rings.number_field
Dexact = polytopes.dodecahedron()
F0 = Dexact.faces(2)[0].as_polyhedron()
v = F0.volume(measure='induced', engine='internal'); v
F4 = Dexact.faces(2)[4].as_polyhedron()
v = F4.volume(measure='induced', engine='internal'); v
RDF(v)    # abs tol 1e-9
# needs sage.groups
Dinexact = polytopes.dodecahedron(exact=False)
F2 = Dinexact.faces(2)[2].as_polyhedron()
w = F2.volume(measure='induced', engine='internal')
RDF(w)    # abs tol 1e-9
all(polytopes.simplex(d).volume(measure='induced')                    # needs sage.rings.number_field sage.symbolic
       == sqrt(d+1)/factorial(d)
    for d in range(1,5))
I = Polyhedron([[-3, 0], [0, 9]])
I.volume(measure='induced')                                           # needs sage.rings.number_field
I.volume(measure='induced_rational')                  # optional - latte_int
T = Polyhedron([[3, 0, 0], [0, 4, 0], [0, 0, 5]])
T.volume(measure='induced')                                           # needs sage.rings.number_field
T.volume(measure='induced_rational')                  # optional - latte_int
Q = Polyhedron(vertices=[(0, 0, 1, 1), (0, 1, 1, 0), (1, 1, 0, 0)])
Q.volume(measure='induced')
Q.volume(measure='induced_rational')                  # optional - latte_int

The volume of a full-dimensional unbounded polyhedron is infinity:

sage: P = Polyhedron(vertices=[[1, 0], [0, 1]], rays=[[1, 1]])
sage: P.volume()
+Infinity
>>> from sage.all import *
>>> P = Polyhedron(vertices=[[Integer(1), Integer(0)], [Integer(0), Integer(1)]], rays=[[Integer(1), Integer(1)]])
>>> P.volume()
+Infinity
P = Polyhedron(vertices=[[1, 0], [0, 1]], rays=[[1, 1]])
P.volume()

The volume of a non full-dimensional unbounded polyhedron depends on the measure used:

sage: P = Polyhedron(ieqs = [[1,1,1], [-1,-1,-1], [3,1,0]]); P
A 1-dimensional polyhedron in QQ^2 defined as the convex hull of 1 vertex and 1 ray
sage: P.volume()
0
sage: P.volume(measure='induced')
+Infinity
sage: P.volume(measure='ambient')
0
sage: P.volume(measure='induced_rational')                  # optional - pynormaliz
+Infinity
sage: P.volume(measure='induced_rational',engine='latte')
+Infinity
>>> from sage.all import *
>>> P = Polyhedron(ieqs = [[Integer(1),Integer(1),Integer(1)], [-Integer(1),-Integer(1),-Integer(1)], [Integer(3),Integer(1),Integer(0)]]); P
A 1-dimensional polyhedron in QQ^2 defined as the convex hull of 1 vertex and 1 ray
>>> P.volume()
0
>>> P.volume(measure='induced')
+Infinity
>>> P.volume(measure='ambient')
0
>>> P.volume(measure='induced_rational')                  # optional - pynormaliz
+Infinity
>>> P.volume(measure='induced_rational',engine='latte')
+Infinity
P = Polyhedron(ieqs = [[1,1,1], [-1,-1,-1], [3,1,0]]); P
P.volume()
P.volume(measure='induced')
P.volume(measure='ambient')
P.volume(measure='induced_rational')                  # optional - pynormaliz
P.volume(measure='induced_rational',engine='latte')

The volume in \(0\)-dimensional space is taken by counting measure:

sage: P = Polyhedron(vertices=[[]]); P
A 0-dimensional polyhedron in ZZ^0 defined as the convex hull of 1 vertex
sage: P.volume()
1
sage: P = Polyhedron(vertices=[]); P
The empty polyhedron in ZZ^0
sage: P.volume()
0
>>> from sage.all import *
>>> P = Polyhedron(vertices=[[]]); P
A 0-dimensional polyhedron in ZZ^0 defined as the convex hull of 1 vertex
>>> P.volume()
1
>>> P = Polyhedron(vertices=[]); P
The empty polyhedron in ZZ^0
>>> P.volume()
0
P = Polyhedron(vertices=[[]]); P
P.volume()
P = Polyhedron(vertices=[]); P
P.volume()