List of faces

This module provides a class to store faces of a polyhedron in Bit-representation.

This class allocates memory to store the faces in. A face will be stored as vertex-incidences, where each Bit represents an incidence. In conversions there a methods to actually convert facets of a polyhedron to bit-representations of vertices stored in ListOfFaces.

Moreover, ListOfFaces calculates the dimension of a polyhedron, assuming the faces are the facets of this polyhedron.

Each face is stored over-aligned according to the chunktype.

EXAMPLES:

Provide enough space to store \(20\) faces as incidences to \(60\) vertices:

sage: from sage.geometry.polyhedron.combinatorial_polyhedron.list_of_faces \
....: import ListOfFaces
sage: face_list = ListOfFaces(20, 60, 20)
sage: face_list.matrix().is_zero()
True
>>> from sage.all import *
>>> from sage.geometry.polyhedron.combinatorial_polyhedron.list_of_faces import ListOfFaces
>>> face_list = ListOfFaces(Integer(20), Integer(60), Integer(20))
>>> face_list.matrix().is_zero()
True

Obtain the facets of a polyhedron:

sage: from sage.geometry.polyhedron.combinatorial_polyhedron.conversions \
....:         import incidence_matrix_to_bit_rep_of_facets
sage: P = polytopes.cube()
sage: face_list = incidence_matrix_to_bit_rep_of_facets(P.incidence_matrix())
sage: face_list = incidence_matrix_to_bit_rep_of_facets(P.incidence_matrix())
sage: face_list.compute_dimension()
3
>>> from sage.all import *
>>> from sage.geometry.polyhedron.combinatorial_polyhedron.conversions         import incidence_matrix_to_bit_rep_of_facets
>>> P = polytopes.cube()
>>> face_list = incidence_matrix_to_bit_rep_of_facets(P.incidence_matrix())
>>> face_list = incidence_matrix_to_bit_rep_of_facets(P.incidence_matrix())
>>> face_list.compute_dimension()
3

Obtain the Vrepresentation of a polyhedron as facet-incidences:

sage: # needs sage.combinat
sage: from sage.geometry.polyhedron.combinatorial_polyhedron.conversions \
....:         import incidence_matrix_to_bit_rep_of_Vrep
sage: P = polytopes.associahedron(['A',3])
sage: face_list = incidence_matrix_to_bit_rep_of_Vrep(P.incidence_matrix())
sage: face_list.compute_dimension()
3
>>> from sage.all import *
>>> # needs sage.combinat
>>> from sage.geometry.polyhedron.combinatorial_polyhedron.conversions         import incidence_matrix_to_bit_rep_of_Vrep
>>> P = polytopes.associahedron(['A',Integer(3)])
>>> face_list = incidence_matrix_to_bit_rep_of_Vrep(P.incidence_matrix())
>>> face_list.compute_dimension()
3

Obtain the facets of a polyhedron as ListOfFaces from a facet list:

sage: from sage.geometry.polyhedron.combinatorial_polyhedron.conversions \
....:         import facets_tuple_to_bit_rep_of_facets
sage: facets = ((0,1,2), (0,1,3), (0,2,3), (1,2,3))
sage: face_list = facets_tuple_to_bit_rep_of_facets(facets, 4)
>>> from sage.all import *
>>> from sage.geometry.polyhedron.combinatorial_polyhedron.conversions         import facets_tuple_to_bit_rep_of_facets
>>> facets = ((Integer(0),Integer(1),Integer(2)), (Integer(0),Integer(1),Integer(3)), (Integer(0),Integer(2),Integer(3)), (Integer(1),Integer(2),Integer(3)))
>>> face_list = facets_tuple_to_bit_rep_of_facets(facets, Integer(4))

Likewise for the Vrepresentatives as facet-incidences:

sage: from sage.geometry.polyhedron.combinatorial_polyhedron.conversions \
....:         import facets_tuple_to_bit_rep_of_Vrep
sage: facets = ((0,1,2), (0,1,3), (0,2,3), (1,2,3))
sage: face_list = facets_tuple_to_bit_rep_of_Vrep(facets, 4)
>>> from sage.all import *
>>> from sage.geometry.polyhedron.combinatorial_polyhedron.conversions         import facets_tuple_to_bit_rep_of_Vrep
>>> facets = ((Integer(0),Integer(1),Integer(2)), (Integer(0),Integer(1),Integer(3)), (Integer(0),Integer(2),Integer(3)), (Integer(1),Integer(2),Integer(3)))
>>> face_list = facets_tuple_to_bit_rep_of_Vrep(facets, Integer(4))

Obtain the matrix of a list of faces:

sage: face_list.matrix()
[1 1 1 0]
[1 1 0 1]
[1 0 1 1]
[0 1 1 1]
>>> from sage.all import *
>>> face_list.matrix()
[1 1 1 0]
[1 1 0 1]
[1 0 1 1]
[0 1 1 1]

See also

base, face_iterator, conversions, polyhedron_faces_lattice.

AUTHOR:

  • Jonathan Kliem (2019-04)

class sage.geometry.polyhedron.combinatorial_polyhedron.list_of_faces.ListOfFaces[source]

Bases: object

A class to store the Bit-representation of faces in.

This class will allocate the memory for the faces.

INPUT:

  • n_faces – the number of faces to be stored

  • n_atoms – the total number of atoms the faces contain

  • n_coatoms – the total number of coatoms of the polyhedron

See also

incidence_matrix_to_bit_rep_of_facets(), incidence_matrix_to_bit_rep_of_Vrep(), facets_tuple_to_bit_rep_of_facets(), facets_tuple_to_bit_rep_of_Vrep(), FaceIterator, CombinatorialPolyhedron.

EXAMPLES:

sage: from sage.geometry.polyhedron.combinatorial_polyhedron.list_of_faces \
....:     import ListOfFaces
sage: facets = ListOfFaces(5, 13, 5)
sage: facets.matrix().dimensions()
(5, 13)
>>> from sage.all import *
>>> from sage.geometry.polyhedron.combinatorial_polyhedron.list_of_faces     import ListOfFaces
>>> facets = ListOfFaces(Integer(5), Integer(13), Integer(5))
>>> facets.matrix().dimensions()
(5, 13)
compute_dimension()[source]

Compute the dimension of a polyhedron by its facets.

This assumes that self is the list of facets of a polyhedron.

EXAMPLES:

sage: from sage.geometry.polyhedron.combinatorial_polyhedron.conversions \
....:     import facets_tuple_to_bit_rep_of_facets, \
....:            facets_tuple_to_bit_rep_of_Vrep
sage: bi_pyr = ((0,1,4), (1,2,4), (2,3,4), (3,0,4),
....:           (0,1,5), (1,2,5), (2,3,5), (3,0,5))
sage: facets = facets_tuple_to_bit_rep_of_facets(bi_pyr, 6)
sage: Vrep = facets_tuple_to_bit_rep_of_Vrep(bi_pyr, 6)
sage: facets.compute_dimension()
3
sage: Vrep.compute_dimension()
3
>>> from sage.all import *
>>> from sage.geometry.polyhedron.combinatorial_polyhedron.conversions     import facets_tuple_to_bit_rep_of_facets,            facets_tuple_to_bit_rep_of_Vrep
>>> bi_pyr = ((Integer(0),Integer(1),Integer(4)), (Integer(1),Integer(2),Integer(4)), (Integer(2),Integer(3),Integer(4)), (Integer(3),Integer(0),Integer(4)),
...           (Integer(0),Integer(1),Integer(5)), (Integer(1),Integer(2),Integer(5)), (Integer(2),Integer(3),Integer(5)), (Integer(3),Integer(0),Integer(5)))
>>> facets = facets_tuple_to_bit_rep_of_facets(bi_pyr, Integer(6))
>>> Vrep = facets_tuple_to_bit_rep_of_Vrep(bi_pyr, Integer(6))
>>> facets.compute_dimension()
3
>>> Vrep.compute_dimension()
3

ALGORITHM:

This is done by iteration:

Computes the facets of one of the facets (i.e. the ridges contained in one of the facets). Then computes the dimension of the facet, by considering its facets.

Repeats until a face has only one facet. Usually this is a vertex.

However, in the unbounded case, this might be different. The face with only one facet might be a ray or a line. So the correct dimension of a polyhedron with one facet is the number of [lines, rays, vertices] that the facet contains.

Hence, we know the dimension of a face, which has only one facet and iteratively we know the dimension of entire polyhedron we started from.

matrix()[source]

Obtain the matrix of self.

Each row represents a face and each column an atom.

EXAMPLES:

sage: from sage.geometry.polyhedron.combinatorial_polyhedron.conversions \
....:     import facets_tuple_to_bit_rep_of_facets, \
....:     facets_tuple_to_bit_rep_of_Vrep
sage: bi_pyr = ((0,1,4), (1,2,4), (2,3,4), (3,0,4), (0,1,5), (1,2,5), (2,3,5), (3,0,5))
sage: facets = facets_tuple_to_bit_rep_of_facets(bi_pyr, 6)
sage: Vrep = facets_tuple_to_bit_rep_of_Vrep(bi_pyr, 6)
sage: facets.matrix()
[1 1 0 0 1 0]
[0 1 1 0 1 0]
[0 0 1 1 1 0]
[1 0 0 1 1 0]
[1 1 0 0 0 1]
[0 1 1 0 0 1]
[0 0 1 1 0 1]
[1 0 0 1 0 1]
sage: facets.matrix().transpose() == Vrep.matrix()
True
>>> from sage.all import *
>>> from sage.geometry.polyhedron.combinatorial_polyhedron.conversions     import facets_tuple_to_bit_rep_of_facets,     facets_tuple_to_bit_rep_of_Vrep
>>> bi_pyr = ((Integer(0),Integer(1),Integer(4)), (Integer(1),Integer(2),Integer(4)), (Integer(2),Integer(3),Integer(4)), (Integer(3),Integer(0),Integer(4)), (Integer(0),Integer(1),Integer(5)), (Integer(1),Integer(2),Integer(5)), (Integer(2),Integer(3),Integer(5)), (Integer(3),Integer(0),Integer(5)))
>>> facets = facets_tuple_to_bit_rep_of_facets(bi_pyr, Integer(6))
>>> Vrep = facets_tuple_to_bit_rep_of_Vrep(bi_pyr, Integer(6))
>>> facets.matrix()
[1 1 0 0 1 0]
[0 1 1 0 1 0]
[0 0 1 1 1 0]
[1 0 0 1 1 0]
[1 1 0 0 0 1]
[0 1 1 0 0 1]
[0 0 1 1 0 1]
[1 0 0 1 0 1]
>>> facets.matrix().transpose() == Vrep.matrix()
True
pyramid()[source]

Return the list of faces of the pyramid.

EXAMPLES:

sage: from sage.geometry.polyhedron.combinatorial_polyhedron.conversions \
....:         import facets_tuple_to_bit_rep_of_facets
sage: facets = ((0,1,2), (0,1,3), (0,2,3), (1,2,3))
sage: face_list = facets_tuple_to_bit_rep_of_facets(facets, 4)
sage: face_list.matrix()
[1 1 1 0]
[1 1 0 1]
[1 0 1 1]
[0 1 1 1]
sage: face_list.pyramid().matrix()
[1 1 1 0 1]
[1 1 0 1 1]
[1 0 1 1 1]
[0 1 1 1 1]
[1 1 1 1 0]
>>> from sage.all import *
>>> from sage.geometry.polyhedron.combinatorial_polyhedron.conversions         import facets_tuple_to_bit_rep_of_facets
>>> facets = ((Integer(0),Integer(1),Integer(2)), (Integer(0),Integer(1),Integer(3)), (Integer(0),Integer(2),Integer(3)), (Integer(1),Integer(2),Integer(3)))
>>> face_list = facets_tuple_to_bit_rep_of_facets(facets, Integer(4))
>>> face_list.matrix()
[1 1 1 0]
[1 1 0 1]
[1 0 1 1]
[0 1 1 1]
>>> face_list.pyramid().matrix()
[1 1 1 0 1]
[1 1 0 1 1]
[1 0 1 1 1]
[0 1 1 1 1]
[1 1 1 1 0]

Incorrect facets that illustrate how this method works:

sage: facets = ((0,1,2,3), (0,1,2,3), (0,1,2,3), (0,1,2,3))
sage: face_list = facets_tuple_to_bit_rep_of_facets(facets, 4)
sage: face_list.matrix()
[1 1 1 1]
[1 1 1 1]
[1 1 1 1]
[1 1 1 1]
sage: face_list.pyramid().matrix()
[1 1 1 1 1]
[1 1 1 1 1]
[1 1 1 1 1]
[1 1 1 1 1]
[1 1 1 1 0]
>>> from sage.all import *
>>> facets = ((Integer(0),Integer(1),Integer(2),Integer(3)), (Integer(0),Integer(1),Integer(2),Integer(3)), (Integer(0),Integer(1),Integer(2),Integer(3)), (Integer(0),Integer(1),Integer(2),Integer(3)))
>>> face_list = facets_tuple_to_bit_rep_of_facets(facets, Integer(4))
>>> face_list.matrix()
[1 1 1 1]
[1 1 1 1]
[1 1 1 1]
[1 1 1 1]
>>> face_list.pyramid().matrix()
[1 1 1 1 1]
[1 1 1 1 1]
[1 1 1 1 1]
[1 1 1 1 1]
[1 1 1 1 0]

sage: facets = ((), (), (), ())
sage: face_list = facets_tuple_to_bit_rep_of_facets(facets, 4)
sage: face_list.matrix()
[0 0 0 0]
[0 0 0 0]
[0 0 0 0]
[0 0 0 0]
sage: face_list.pyramid().matrix()
[0 0 0 0 1]
[0 0 0 0 1]
[0 0 0 0 1]
[0 0 0 0 1]
[1 1 1 1 0]
>>> from sage.all import *
>>> facets = ((), (), (), ())
>>> face_list = facets_tuple_to_bit_rep_of_facets(facets, Integer(4))
>>> face_list.matrix()
[0 0 0 0]
[0 0 0 0]
[0 0 0 0]
[0 0 0 0]
>>> face_list.pyramid().matrix()
[0 0 0 0 1]
[0 0 0 0 1]
[0 0 0 0 1]
[0 0 0 0 1]
[1 1 1 1 0]