A class to keep information about faces of a polyhedron

This module gives you a tool to work with the faces of a polyhedron and their relative position. First, you need to find the faces. To get the faces in a particular dimension, use the face() method:

sage: P = polytopes.cross_polytope(3)
sage: P.faces(3)
(A 3-dimensional face of a Polyhedron in ZZ^3 defined as the convex hull of 6 vertices,)
sage: [f.ambient_V_indices() for f in P.facets()]
[(0, 1, 2),
 (0, 1, 3),
 (0, 2, 4),
 (0, 3, 4),
 (3, 4, 5),
 (2, 4, 5),
 (1, 3, 5),
 (1, 2, 5)]
sage: [f.ambient_V_indices() for f in P.faces(1)]
[(0, 1),
 (0, 2),
 (1, 2),
 (0, 3),
 (1, 3),
 (0, 4),
 (2, 4),
 (3, 4),
 (2, 5),
 (3, 5),
 (4, 5),
 (1, 5)]

or face_lattice() to get the whole face lattice as a poset:

sage: P.face_lattice()
Finite lattice containing 28 elements with distinguished linear extension

The faces are printed in shorthand notation where each integer is the index of a vertex/ray/line in the same order as the containing Polyhedron’s Vrepresentation()

sage: face = P.faces(1)[3];  face
A 1-dimensional face of a Polyhedron in ZZ^3 defined as the convex hull of 2 vertices
sage: face.ambient_V_indices()
(0, 3)
sage: P.Vrepresentation(0)
A vertex at (-1, 0, 0)
sage: P.Vrepresentation(3)
A vertex at (0, 0, 1)
sage: face.vertices()
(A vertex at (-1, 0, 0), A vertex at (0, 0, 1))

The face itself is not represented by Sage’s sage.geometry.polyhedron.constructor.Polyhedron() class, but by an auxiliary class to keep the information. You can get the face as a polyhedron with the PolyhedronFace.as_polyhedron() method:

sage: face.as_polyhedron()
A 1-dimensional polyhedron in ZZ^3 defined as the convex hull of 2 vertices
sage: _.equations()
(An equation (0, 1, 0) x + 0 == 0,
 An equation (1, 0, -1) x + 1 == 0)
class sage.geometry.polyhedron.face.PolyhedronFace(polyhedron, V_indices, H_indices)

Bases: sage.structure.sage_object.SageObject

A face of a polyhedron.

This class is for use in face_lattice().

INPUT:

No checking is performed whether the H/V-representation indices actually determine a face of the polyhedron. You should not manually create PolyhedronFace objects unless you know what you are doing.

OUTPUT:

A PolyhedronFace.

EXAMPLES:

sage: octahedron = polytopes.cross_polytope(3)
sage: inequality = octahedron.Hrepresentation(2)
sage: face_h = tuple([ inequality ])
sage: face_v = tuple( inequality.incident() )
sage: face_h_indices = [ h.index() for h in face_h ]
sage: face_v_indices = [ v.index() for v in face_v ]
sage: from sage.geometry.polyhedron.face import PolyhedronFace
sage: face = PolyhedronFace(octahedron, face_v_indices, face_h_indices)
sage: face
A 2-dimensional face of a Polyhedron in ZZ^3 defined as the convex hull of 3 vertices
sage: face.dim()
2
sage: face.ambient_V_indices()
(0, 1, 2)
sage: face.ambient_Hrepresentation()
(An inequality (1, 1, 1) x + 1 >= 0,)
sage: face.ambient_Vrepresentation()
(A vertex at (-1, 0, 0), A vertex at (0, -1, 0), A vertex at (0, 0, -1))
ambient_H_indices()

Return the indices of the H-representation objects of the ambient polyhedron that make up the H-representation of self.

See also ambient_Hrepresentation().

OUTPUT:

Tuple of indices

EXAMPLES:

sage: Q = polytopes.cross_polytope(3)
sage: F = Q.faces(1)
sage: [f.ambient_H_indices() for f in F]
[(1, 2),
 (2, 3),
 (2, 7),
 (0, 1),
 (1, 6),
 (0, 3),
 (3, 4),
 (0, 5),
 (4, 7),
 (5, 6),
 (4, 5),
 (6, 7)]
ambient_Hrepresentation(index=None)

Return the H-representation objects of the ambient polytope defining the face.

INPUT:

  • index – optional. Either an integer or None (default).

OUTPUT:

If the optional argument is not present, a tuple of H-representation objects. Each entry is either an inequality or an equation.

If the optional integer index is specified, the index-th element of the tuple is returned.

EXAMPLES:

sage: square = polytopes.hypercube(2)
sage: for face in square.face_lattice():
....:     print(face.ambient_Hrepresentation())
(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)
(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)
(An inequality (0, 1) x + 1 >= 0, An inequality (-1, 0) x + 1 >= 0)
(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,)
(An inequality (-1, 0) x + 1 >= 0,)
(An inequality (0, -1) x + 1 >= 0,)
()
ambient_V_indices()

Return the indices of the V-representation objects of the ambient polyhedron that make up the V-representation of self.

See also ambient_Vrepresentation().

OUTPUT:

Tuple of indices

EXAMPLES:

sage: P = polytopes.cube()
sage: F = P.faces(2)
sage: [f.ambient_V_indices() for f in F]
[(0, 1, 2, 3),
 (0, 1, 4, 5),
 (0, 2, 4, 6),
 (1, 3, 5, 7),
 (2, 3, 6, 7),
 (4, 5, 6, 7)]
ambient_Vrepresentation(index=None)

Return the V-representation objects of the ambient polytope defining the face.

INPUT:

  • index – optional. Either an integer or None (default).

OUTPUT:

If the optional argument is not present, a tuple of V-representation objects. Each entry is either a vertex, a ray, or a line.

If the optional integer index is specified, the index-th element of the tuple is returned.

EXAMPLES:

sage: square = polytopes.hypercube(2)
sage: for fl in square.face_lattice():
....:     print(fl.ambient_Vrepresentation())
()
(A vertex at (-1, -1),)
(A vertex at (-1, 1),)
(A vertex at (1, -1),)
(A vertex at (1, 1),)
(A vertex at (-1, -1), A vertex at (-1, 1))
(A vertex at (-1, -1), A vertex at (1, -1))
(A vertex at (1, -1), A vertex at (1, 1))
(A vertex at (-1, 1), A vertex at (1, 1))
(A vertex at (-1, -1), A vertex at (-1, 1),
 A vertex at (1, -1), A vertex at (1, 1))
ambient_dim()

Return the dimension of the containing polyhedron.

EXAMPLES:

sage: P = Polyhedron(vertices = [[1,0,0,0],[0,1,0,0]])
sage: face = P.faces(1)[0]
sage: face.ambient_dim()
4
as_polyhedron()

Return the face as an independent polyhedron.

OUTPUT:

A polyhedron.

EXAMPLES:

sage: P = polytopes.cross_polytope(3);  P
A 3-dimensional polyhedron in ZZ^3 defined as the convex hull of 6 vertices
sage: face = P.faces(2)[3]
sage: face
A 2-dimensional face of a Polyhedron in ZZ^3 defined as the convex hull of 3 vertices
sage: face.as_polyhedron()
A 2-dimensional polyhedron in ZZ^3 defined as the convex hull of 3 vertices

sage: P.intersection(face.as_polyhedron()) == face.as_polyhedron()
True
dim()

Return the dimension of the face.

OUTPUT:

Integer.

EXAMPLES:

sage: fl = polytopes.dodecahedron().face_lattice()
sage: [ x.dim() for x in fl ]
[-1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0,
  1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1, 1,
  1, 1, 1, 1, 1, 1, 1, 1, 1, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 2, 3]
line_generator()

Return a generator for the lines of the face.

EXAMPLES:

sage: pr = Polyhedron(rays = [[1,0],[-1,0],[0,1]], vertices = [[-1,-1]])
sage: face = pr.faces(1)[0]
sage: next(face.line_generator())
A line in the direction (1, 0)
lines()

Return all lines of the face.

OUTPUT:

A tuple of lines.

EXAMPLES:

sage: p = Polyhedron(rays = [[1,0],[-1,0],[0,1],[1,1]], vertices = [[-2,-2],[2,3]])
sage: p.lines()
(A line in the direction (1, 0),)
n_ambient_Hrepresentation()

Return the number of objects that make up the ambient H-representation of the polyhedron.

See also ambient_Hrepresentation().

OUTPUT:

Integer.

EXAMPLES:

sage: p = polytopes.cross_polytope(4)
sage: face = p.face_lattice()[10]
sage: face
A 1-dimensional face of a Polyhedron in ZZ^4 defined as the convex hull of 2 vertices
sage: face.ambient_Hrepresentation()
(An inequality (1, -1, 1, -1) x + 1 >= 0,
 An inequality (1, 1, 1, 1) x + 1 >= 0,
 An inequality (1, 1, 1, -1) x + 1 >= 0,
 An inequality (1, -1, 1, 1) x + 1 >= 0)
sage: face.n_ambient_Hrepresentation()
4
n_ambient_Vrepresentation()

Return the number of objects that make up the ambient V-representation of the polyhedron.

See also ambient_Vrepresentation().

OUTPUT:

Integer.

EXAMPLES:

sage: p = polytopes.cross_polytope(4)
sage: face = p.face_lattice()[10]
sage: face
A 1-dimensional face of a Polyhedron in ZZ^4 defined as the convex hull of 2 vertices
sage: face.ambient_Vrepresentation()
(A vertex at (-1, 0, 0, 0), A vertex at (0, 0, -1, 0))
sage: face.n_ambient_Vrepresentation()
2
n_lines()

Return the number of lines of the face.

OUTPUT:

Integer.

EXAMPLES:

sage: p = Polyhedron(rays = [[1,0],[-1,0],[0,1],[1,1]], vertices = [[-2,-2],[2,3]])
sage: p.n_lines()
1
n_rays()

Return the number of rays of the face.

OUTPUT:

Integer.

EXAMPLES:

sage: p = Polyhedron(ieqs = [[0,0,0,1],[0,0,1,0],[1,1,0,0]])
sage: face = p.faces(2)[0]
sage: face.n_rays()
2
n_vertices()

Return the number of vertices of the face.

OUTPUT:

Integer.

EXAMPLES:

sage: Q = polytopes.cross_polytope(3)
sage: face = Q.faces(2)[0]
sage: face.n_vertices()
3
polyhedron()

Return the containing polyhedron.

EXAMPLES:

sage: P = polytopes.cross_polytope(3); P
A 3-dimensional polyhedron in ZZ^3 defined as the convex hull of 6 vertices
sage: face = P.facets()[3]
sage: face
A 2-dimensional face of a Polyhedron in ZZ^3 defined as the convex hull of 3 vertices
sage: face.polyhedron()
A 3-dimensional polyhedron in ZZ^3 defined as the convex hull of 6 vertices
ray_generator()

Return a generator for the rays of the face.

EXAMPLES:

sage: pi = Polyhedron(ieqs = [[1,1,0],[1,0,1]])
sage: face = pi.faces(1)[0]
sage: next(face.ray_generator())
A ray in the direction (1, 0)
rays()

Return the rays of the face.

OUTPUT:

A tuple of rays.

EXAMPLES:

sage: p = Polyhedron(ieqs = [[0,0,0,1],[0,0,1,0],[1,1,0,0]])
sage: face = p.faces(2)[0]
sage: face.rays()
(A ray in the direction (1, 0, 0), A ray in the direction (0, 1, 0))
vertex_generator()

Return a generator for the vertices of the face.

EXAMPLES:

sage: triangle = Polyhedron(vertices=[[1,0],[0,1],[1,1]])
sage: face = triangle.facets()[0]
sage: for v in face.vertex_generator(): print(v)
A vertex at (0, 1)
A vertex at (1, 0)
sage: type(face.vertex_generator())
<... 'generator'>
vertices()

Return all vertices of the face.

OUTPUT:

A tuple of vertices.

EXAMPLES:

sage: triangle = Polyhedron(vertices=[[1,0],[0,1],[1,1]])
sage: face = triangle.faces(1)[0]
sage: face.vertices()
(A vertex at (0, 1), A vertex at (1, 0))