Base classes for triangulations#

We provide (fast) cython implementations here.

AUTHORS:

  • Volker Braun (2010-09-14): initial version.

class sage.geometry.triangulation.base.ConnectedTriangulationsIterator[source]#

Bases: SageObject

A Python shim for the C++-class ‘triangulations’

INPUT:

  • point_configuration – a PointConfiguration.

  • seed – a regular triangulation or None (default). In the latter case, a suitable triangulation is generated automatically. Otherwise, you can explicitly specify the seed triangulation as

    • A Triangulation object, or

    • an iterable of iterables specifying the vertices of the simplices, or

    • an iterable of integers, which are then considered the enumerated simplices (see simplex_to_int().

  • star – either None (default) or an integer. If an integer is passed, all returned triangulations will be star with respect to the

  • fine – boolean (default: False). Whether to return only fine triangulations, that is, simplicial decompositions that make use of all the points of the configuration.

OUTPUT:

An iterator. The generated values are tuples of integers, which encode simplices of the triangulation. The output is a suitable input to Triangulation.

EXAMPLES:

sage: p = PointConfiguration([[0,0],[0,1],[1,0],[1,1],[-1,-1]])
sage: from sage.geometry.triangulation.base import ConnectedTriangulationsIterator
sage: ci = ConnectedTriangulationsIterator(p)
sage: next(ci)
(9, 10)
sage: next(ci)
(2, 3, 4, 5)
sage: next(ci)
(7, 8)
sage: next(ci)
(1, 3, 5, 7)
sage: next(ci)
Traceback (most recent call last):
...
StopIteration
>>> from sage.all import *
>>> p = PointConfiguration([[Integer(0),Integer(0)],[Integer(0),Integer(1)],[Integer(1),Integer(0)],[Integer(1),Integer(1)],[-Integer(1),-Integer(1)]])
>>> from sage.geometry.triangulation.base import ConnectedTriangulationsIterator
>>> ci = ConnectedTriangulationsIterator(p)
>>> next(ci)
(9, 10)
>>> next(ci)
(2, 3, 4, 5)
>>> next(ci)
(7, 8)
>>> next(ci)
(1, 3, 5, 7)
>>> next(ci)
Traceback (most recent call last):
...
StopIteration

You can reconstruct the triangulation from the compressed output via:

sage: from sage.geometry.triangulation.element import Triangulation
sage: Triangulation((2, 3, 4, 5), p)
(<0,1,3>, <0,1,4>, <0,2,3>, <0,2,4>)
>>> from sage.all import *
>>> from sage.geometry.triangulation.element import Triangulation
>>> Triangulation((Integer(2), Integer(3), Integer(4), Integer(5)), p)
(<0,1,3>, <0,1,4>, <0,2,3>, <0,2,4>)

How to use the restrictions:

sage: ci = ConnectedTriangulationsIterator(p, fine=True)
sage: list(ci)
[(2, 3, 4, 5), (1, 3, 5, 7)]
sage: ci = ConnectedTriangulationsIterator(p, star=1)
sage: list(ci)
[(7, 8)]
sage: ci = ConnectedTriangulationsIterator(p, star=1, fine=True)
sage: list(ci)
[]
>>> from sage.all import *
>>> ci = ConnectedTriangulationsIterator(p, fine=True)
>>> list(ci)
[(2, 3, 4, 5), (1, 3, 5, 7)]
>>> ci = ConnectedTriangulationsIterator(p, star=Integer(1))
>>> list(ci)
[(7, 8)]
>>> ci = ConnectedTriangulationsIterator(p, star=Integer(1), fine=True)
>>> list(ci)
[]
class sage.geometry.triangulation.base.Point[source]#

Bases: SageObject

A point of a point configuration.

Note that the coordinates of the points of a point configuration are somewhat arbitrary. What counts are the abstract linear relations between the points, for example encoded by the circuits().

Warning

You should not create Point objects manually. The constructor of PointConfiguration_base takes care of this for you.

INPUT:

  • point_configurationPointConfiguration_base. The point configuration to which the point belongs.

  • i – integer. The index of the point in the point configuration.

  • projective – the projective coordinates of the point.

  • affine – the affine coordinates of the point.

  • reduced – the reduced (with linearities removed) coordinates of the point.

EXAMPLES:

sage: pc = PointConfiguration([(0,0)])
sage: from sage.geometry.triangulation.base import Point
sage: Point(pc, 123, (0,0,1), (0,0), ())
P(0, 0)
>>> from sage.all import *
>>> pc = PointConfiguration([(Integer(0),Integer(0))])
>>> from sage.geometry.triangulation.base import Point
>>> Point(pc, Integer(123), (Integer(0),Integer(0),Integer(1)), (Integer(0),Integer(0)), ())
P(0, 0)
affine()[source]#

Return the affine coordinates of the point in the ambient space.

OUTPUT:

A tuple containing the coordinates.

EXAMPLES:

sage: pc = PointConfiguration([[10, 0, 1], [10, 0, 0], [10, 2, 3]])
sage: p = pc.point(2); p
P(10, 2, 3)
sage: p.affine()
(10, 2, 3)
sage: p.projective()
(10, 2, 3, 1)
sage: p.reduced_affine()
(2, 2)
sage: p.reduced_projective()
(2, 2, 1)
sage: p.reduced_affine_vector()
(2, 2)
>>> from sage.all import *
>>> pc = PointConfiguration([[Integer(10), Integer(0), Integer(1)], [Integer(10), Integer(0), Integer(0)], [Integer(10), Integer(2), Integer(3)]])
>>> p = pc.point(Integer(2)); p
P(10, 2, 3)
>>> p.affine()
(10, 2, 3)
>>> p.projective()
(10, 2, 3, 1)
>>> p.reduced_affine()
(2, 2)
>>> p.reduced_projective()
(2, 2, 1)
>>> p.reduced_affine_vector()
(2, 2)
index()[source]#

Return the index of the point in the point configuration.

EXAMPLES:

sage: pc = PointConfiguration([[0, 1], [0, 0], [1, 0]])
sage: p = pc.point(2); p
P(1, 0)
sage: p.index()
2
>>> from sage.all import *
>>> pc = PointConfiguration([[Integer(0), Integer(1)], [Integer(0), Integer(0)], [Integer(1), Integer(0)]])
>>> p = pc.point(Integer(2)); p
P(1, 0)
>>> p.index()
2
point_configuration()[source]#

Return the point configuration to which the point belongs.

OUTPUT:

A PointConfiguration.

EXAMPLES:

sage: pc = PointConfiguration([ (0,0), (1,0), (0,1) ])
sage: p = pc.point(0)
sage: p is pc.point(0)
True
sage: p.point_configuration() is pc
True
>>> from sage.all import *
>>> pc = PointConfiguration([ (Integer(0),Integer(0)), (Integer(1),Integer(0)), (Integer(0),Integer(1)) ])
>>> p = pc.point(Integer(0))
>>> p is pc.point(Integer(0))
True
>>> p.point_configuration() is pc
True
projective()[source]#

Return the projective coordinates of the point in the ambient space.

OUTPUT:

A tuple containing the coordinates.

EXAMPLES:

sage: pc = PointConfiguration([[10, 0, 1], [10, 0, 0], [10, 2, 3]])
sage: p = pc.point(2); p
P(10, 2, 3)
sage: p.affine()
(10, 2, 3)
sage: p.projective()
(10, 2, 3, 1)
sage: p.reduced_affine()
(2, 2)
sage: p.reduced_projective()
(2, 2, 1)
sage: p.reduced_affine_vector()
(2, 2)
>>> from sage.all import *
>>> pc = PointConfiguration([[Integer(10), Integer(0), Integer(1)], [Integer(10), Integer(0), Integer(0)], [Integer(10), Integer(2), Integer(3)]])
>>> p = pc.point(Integer(2)); p
P(10, 2, 3)
>>> p.affine()
(10, 2, 3)
>>> p.projective()
(10, 2, 3, 1)
>>> p.reduced_affine()
(2, 2)
>>> p.reduced_projective()
(2, 2, 1)
>>> p.reduced_affine_vector()
(2, 2)
reduced_affine()[source]#

Return the affine coordinates of the point on the hyperplane spanned by the point configuration.

OUTPUT:

A tuple containing the coordinates.

EXAMPLES:

sage: pc = PointConfiguration([[10, 0, 1], [10, 0, 0], [10, 2, 3]])
sage: p = pc.point(2); p
P(10, 2, 3)
sage: p.affine()
(10, 2, 3)
sage: p.projective()
(10, 2, 3, 1)
sage: p.reduced_affine()
(2, 2)
sage: p.reduced_projective()
(2, 2, 1)
sage: p.reduced_affine_vector()
(2, 2)
>>> from sage.all import *
>>> pc = PointConfiguration([[Integer(10), Integer(0), Integer(1)], [Integer(10), Integer(0), Integer(0)], [Integer(10), Integer(2), Integer(3)]])
>>> p = pc.point(Integer(2)); p
P(10, 2, 3)
>>> p.affine()
(10, 2, 3)
>>> p.projective()
(10, 2, 3, 1)
>>> p.reduced_affine()
(2, 2)
>>> p.reduced_projective()
(2, 2, 1)
>>> p.reduced_affine_vector()
(2, 2)
reduced_affine_vector()[source]#

Return the affine coordinates of the point on the hyperplane spanned by the point configuration.

OUTPUT:

A tuple containing the coordinates.

EXAMPLES:

sage: pc = PointConfiguration([[10, 0, 1], [10, 0, 0], [10, 2, 3]])
sage: p = pc.point(2); p
P(10, 2, 3)
sage: p.affine()
(10, 2, 3)
sage: p.projective()
(10, 2, 3, 1)
sage: p.reduced_affine()
(2, 2)
sage: p.reduced_projective()
(2, 2, 1)
sage: p.reduced_affine_vector()
(2, 2)
>>> from sage.all import *
>>> pc = PointConfiguration([[Integer(10), Integer(0), Integer(1)], [Integer(10), Integer(0), Integer(0)], [Integer(10), Integer(2), Integer(3)]])
>>> p = pc.point(Integer(2)); p
P(10, 2, 3)
>>> p.affine()
(10, 2, 3)
>>> p.projective()
(10, 2, 3, 1)
>>> p.reduced_affine()
(2, 2)
>>> p.reduced_projective()
(2, 2, 1)
>>> p.reduced_affine_vector()
(2, 2)
reduced_projective()[source]#

Return the projective coordinates of the point on the hyperplane spanned by the point configuration.

OUTPUT:

A tuple containing the coordinates.

EXAMPLES:

sage: pc = PointConfiguration([[10, 0, 1], [10, 0, 0], [10, 2, 3]])
sage: p = pc.point(2); p
P(10, 2, 3)
sage: p.affine()
(10, 2, 3)
sage: p.projective()
(10, 2, 3, 1)
sage: p.reduced_affine()
(2, 2)
sage: p.reduced_projective()
(2, 2, 1)
sage: p.reduced_affine_vector()
(2, 2)
>>> from sage.all import *
>>> pc = PointConfiguration([[Integer(10), Integer(0), Integer(1)], [Integer(10), Integer(0), Integer(0)], [Integer(10), Integer(2), Integer(3)]])
>>> p = pc.point(Integer(2)); p
P(10, 2, 3)
>>> p.affine()
(10, 2, 3)
>>> p.projective()
(10, 2, 3, 1)
>>> p.reduced_affine()
(2, 2)
>>> p.reduced_projective()
(2, 2, 1)
>>> p.reduced_affine_vector()
(2, 2)
reduced_projective_vector()[source]#

Return the affine coordinates of the point on the hyperplane spanned by the point configuration.

OUTPUT:

A tuple containing the coordinates.

EXAMPLES:

sage: pc = PointConfiguration([[10, 0, 1], [10, 0, 0], [10, 2, 3]])
sage: p = pc.point(2); p
P(10, 2, 3)
sage: p.affine()
(10, 2, 3)
sage: p.projective()
(10, 2, 3, 1)
sage: p.reduced_affine()
(2, 2)
sage: p.reduced_projective()
(2, 2, 1)
sage: p.reduced_affine_vector()
(2, 2)
sage: type(p.reduced_affine_vector())
<class 'sage.modules.vector_rational_dense.Vector_rational_dense'>
>>> from sage.all import *
>>> pc = PointConfiguration([[Integer(10), Integer(0), Integer(1)], [Integer(10), Integer(0), Integer(0)], [Integer(10), Integer(2), Integer(3)]])
>>> p = pc.point(Integer(2)); p
P(10, 2, 3)
>>> p.affine()
(10, 2, 3)
>>> p.projective()
(10, 2, 3, 1)
>>> p.reduced_affine()
(2, 2)
>>> p.reduced_projective()
(2, 2, 1)
>>> p.reduced_affine_vector()
(2, 2)
>>> type(p.reduced_affine_vector())
<class 'sage.modules.vector_rational_dense.Vector_rational_dense'>
class sage.geometry.triangulation.base.PointConfiguration_base[source]#

Bases: Parent

The cython abstract base class for PointConfiguration.

Warning

You should not instantiate this base class, but only its derived class PointConfiguration.

ambient_dim()[source]#

Return the dimension of the ambient space of the point configuration.

See also dimension()

EXAMPLES:

sage: p = PointConfiguration([[0,0,0]])
sage: p.ambient_dim()
3
sage: p.dim()
0
>>> from sage.all import *
>>> p = PointConfiguration([[Integer(0),Integer(0),Integer(0)]])
>>> p.ambient_dim()
3
>>> p.dim()
0
base_ring()[source]#

Return the base ring, that is, the ring containing the coordinates of the points.

OUTPUT:

A ring.

EXAMPLES:

sage: p = PointConfiguration([(0,0)])
sage: p.base_ring()
Integer Ring

sage: p = PointConfiguration([(1/2,3)])
sage: p.base_ring()
Rational Field

sage: p = PointConfiguration([(0.2, 5)])
sage: p.base_ring()
Real Field with 53 bits of precision
>>> from sage.all import *
>>> p = PointConfiguration([(Integer(0),Integer(0))])
>>> p.base_ring()
Integer Ring

>>> p = PointConfiguration([(Integer(1)/Integer(2),Integer(3))])
>>> p.base_ring()
Rational Field

>>> p = PointConfiguration([(RealNumber('0.2'), Integer(5))])
>>> p.base_ring()
Real Field with 53 bits of precision
dim()[source]#

Return the actual dimension of the point configuration.

See also ambient_dim()

EXAMPLES:

sage: p = PointConfiguration([[0,0,0]])
sage: p.ambient_dim()
3
sage: p.dim()
0
>>> from sage.all import *
>>> p = PointConfiguration([[Integer(0),Integer(0),Integer(0)]])
>>> p.ambient_dim()
3
>>> p.dim()
0
int_to_simplex(s)[source]#

Reverse the enumeration of possible simplices in simplex_to_int().

The enumeration is compatible with [PUNTOS].

INPUT:

  • s – int. An integer that uniquely specifies a simplex.

OUTPUT:

An ordered tuple consisting of the indices of the vertices of the simplex.

EXAMPLES:

sage: U=matrix([
....:    [ 0, 0, 0, 0, 0, 2, 4,-1, 1, 1, 0, 0, 1, 0],
....:    [ 0, 0, 0, 1, 0, 0,-1, 0, 0, 0, 0, 0, 0, 0],
....:    [ 0, 2, 0, 0, 0, 0,-1, 0, 1, 0, 1, 0, 0, 1],
....:    [ 0, 1, 1, 0, 0, 1, 0,-2, 1, 0, 0,-1, 1, 1],
....:    [ 0, 0, 0, 0, 1, 0,-1, 0, 0, 0, 0, 0, 0, 0]
....: ])
sage: pc = PointConfiguration(U.columns())
sage: pc.simplex_to_int([1,3,4,7,10,13])
1678
sage: pc.int_to_simplex(1678)
(1, 3, 4, 7, 10, 13)
>>> from sage.all import *
>>> U=matrix([
...    [ Integer(0), Integer(0), Integer(0), Integer(0), Integer(0), Integer(2), Integer(4),-Integer(1), Integer(1), Integer(1), Integer(0), Integer(0), Integer(1), Integer(0)],
...    [ Integer(0), Integer(0), Integer(0), Integer(1), Integer(0), Integer(0),-Integer(1), Integer(0), Integer(0), Integer(0), Integer(0), Integer(0), Integer(0), Integer(0)],
...    [ Integer(0), Integer(2), Integer(0), Integer(0), Integer(0), Integer(0),-Integer(1), Integer(0), Integer(1), Integer(0), Integer(1), Integer(0), Integer(0), Integer(1)],
...    [ Integer(0), Integer(1), Integer(1), Integer(0), Integer(0), Integer(1), Integer(0),-Integer(2), Integer(1), Integer(0), Integer(0),-Integer(1), Integer(1), Integer(1)],
...    [ Integer(0), Integer(0), Integer(0), Integer(0), Integer(1), Integer(0),-Integer(1), Integer(0), Integer(0), Integer(0), Integer(0), Integer(0), Integer(0), Integer(0)]
... ])
>>> pc = PointConfiguration(U.columns())
>>> pc.simplex_to_int([Integer(1),Integer(3),Integer(4),Integer(7),Integer(10),Integer(13)])
1678
>>> pc.int_to_simplex(Integer(1678))
(1, 3, 4, 7, 10, 13)
is_affine()[source]#

Return whether the configuration is defined by affine points.

OUTPUT:

Boolean. If true, the homogeneous coordinates all have \(1\) as their last entry.

EXAMPLES:

sage: p = PointConfiguration([(0.2, 5), (3, 0.1)])
sage: p.is_affine()
True

sage: p = PointConfiguration([(0.2, 5, 1), (3, 0.1, 1)], projective=True)
sage: p.is_affine()
False
>>> from sage.all import *
>>> p = PointConfiguration([(RealNumber('0.2'), Integer(5)), (Integer(3), RealNumber('0.1'))])
>>> p.is_affine()
True

>>> p = PointConfiguration([(RealNumber('0.2'), Integer(5), Integer(1)), (Integer(3), RealNumber('0.1'), Integer(1))], projective=True)
>>> p.is_affine()
False
n_points()[source]#

Return the number of points.

Same as len(self).

EXAMPLES:

sage: p = PointConfiguration([[0,0],[0,1],[1,0],[1,1],[-1,-1]])
sage: p
A point configuration in affine 2-space over Integer Ring
consisting of 5 points. The triangulations of this point
configuration are assumed to be connected, not necessarily
fine, not necessarily regular.
sage: len(p)
5
sage: p.n_points()
5
>>> from sage.all import *
>>> p = PointConfiguration([[Integer(0),Integer(0)],[Integer(0),Integer(1)],[Integer(1),Integer(0)],[Integer(1),Integer(1)],[-Integer(1),-Integer(1)]])
>>> p
A point configuration in affine 2-space over Integer Ring
consisting of 5 points. The triangulations of this point
configuration are assumed to be connected, not necessarily
fine, not necessarily regular.
>>> len(p)
5
>>> p.n_points()
5
point(i)[source]#

Return the i-th point of the configuration.

Same as __getitem__()

INPUT:

  • i – integer.

OUTPUT:

A point of the point configuration.

EXAMPLES:

sage: pconfig = PointConfiguration([[0,0],[0,1],[1,0],[1,1],[-1,-1]])
sage: list(pconfig)
[P(0, 0), P(0, 1), P(1, 0), P(1, 1), P(-1, -1)]
sage: [ p for p in pconfig.points() ]
[P(0, 0), P(0, 1), P(1, 0), P(1, 1), P(-1, -1)]
sage: pconfig.point(0)
P(0, 0)
sage: pconfig[0]
P(0, 0)
sage: pconfig.point(1)
P(0, 1)
sage: pconfig.point( pconfig.n_points()-1 )
P(-1, -1)
>>> from sage.all import *
>>> pconfig = PointConfiguration([[Integer(0),Integer(0)],[Integer(0),Integer(1)],[Integer(1),Integer(0)],[Integer(1),Integer(1)],[-Integer(1),-Integer(1)]])
>>> list(pconfig)
[P(0, 0), P(0, 1), P(1, 0), P(1, 1), P(-1, -1)]
>>> [ p for p in pconfig.points() ]
[P(0, 0), P(0, 1), P(1, 0), P(1, 1), P(-1, -1)]
>>> pconfig.point(Integer(0))
P(0, 0)
>>> pconfig[Integer(0)]
P(0, 0)
>>> pconfig.point(Integer(1))
P(0, 1)
>>> pconfig.point( pconfig.n_points()-Integer(1) )
P(-1, -1)
points()[source]#

Return a list of the points.

OUTPUT:

A list of the points. See also the __iter__() method, which returns the corresponding generator.

EXAMPLES:

sage: pconfig = PointConfiguration([[0,0],[0,1],[1,0],[1,1],[-1,-1]])
sage: list(pconfig)
[P(0, 0), P(0, 1), P(1, 0), P(1, 1), P(-1, -1)]
sage: [ p for p in pconfig.points() ]
[P(0, 0), P(0, 1), P(1, 0), P(1, 1), P(-1, -1)]
sage: pconfig.point(0)
P(0, 0)
sage: pconfig.point(1)
P(0, 1)
sage: pconfig.point( pconfig.n_points()-1 )
P(-1, -1)
>>> from sage.all import *
>>> pconfig = PointConfiguration([[Integer(0),Integer(0)],[Integer(0),Integer(1)],[Integer(1),Integer(0)],[Integer(1),Integer(1)],[-Integer(1),-Integer(1)]])
>>> list(pconfig)
[P(0, 0), P(0, 1), P(1, 0), P(1, 1), P(-1, -1)]
>>> [ p for p in pconfig.points() ]
[P(0, 0), P(0, 1), P(1, 0), P(1, 1), P(-1, -1)]
>>> pconfig.point(Integer(0))
P(0, 0)
>>> pconfig.point(Integer(1))
P(0, 1)
>>> pconfig.point( pconfig.n_points()-Integer(1) )
P(-1, -1)
reduced_affine_vector_space()[source]#

Return the vector space that contains the affine points.

OUTPUT:

A vector space over the fraction field of base_ring().

EXAMPLES:

sage: p = PointConfiguration([[0,0,0], [1,2,3]])
sage: p.base_ring()
Integer Ring
sage: p.reduced_affine_vector_space()
Vector space of dimension 1 over Rational Field
sage: p.reduced_projective_vector_space()
Vector space of dimension 2 over Rational Field
>>> from sage.all import *
>>> p = PointConfiguration([[Integer(0),Integer(0),Integer(0)], [Integer(1),Integer(2),Integer(3)]])
>>> p.base_ring()
Integer Ring
>>> p.reduced_affine_vector_space()
Vector space of dimension 1 over Rational Field
>>> p.reduced_projective_vector_space()
Vector space of dimension 2 over Rational Field
reduced_projective_vector_space()[source]#

Return the vector space that is spanned by the homogeneous coordinates.

OUTPUT:

A vector space over the fraction field of base_ring().

EXAMPLES:

sage: p = PointConfiguration([[0,0,0], [1,2,3]])
sage: p.base_ring()
Integer Ring
sage: p.reduced_affine_vector_space()
Vector space of dimension 1 over Rational Field
sage: p.reduced_projective_vector_space()
Vector space of dimension 2 over Rational Field
>>> from sage.all import *
>>> p = PointConfiguration([[Integer(0),Integer(0),Integer(0)], [Integer(1),Integer(2),Integer(3)]])
>>> p.base_ring()
Integer Ring
>>> p.reduced_affine_vector_space()
Vector space of dimension 1 over Rational Field
>>> p.reduced_projective_vector_space()
Vector space of dimension 2 over Rational Field
simplex_to_int(simplex)[source]#

Return an integer that uniquely identifies the given simplex.

See also the inverse method int_to_simplex().

The enumeration is compatible with [PUNTOS].

INPUT:

  • simplex – iterable, for example a list. The elements are the vertex indices of the simplex.

OUTPUT:

An integer that uniquely specifies the simplex.

EXAMPLES:

sage: U=matrix([
....:    [ 0, 0, 0, 0, 0, 2, 4,-1, 1, 1, 0, 0, 1, 0],
....:    [ 0, 0, 0, 1, 0, 0,-1, 0, 0, 0, 0, 0, 0, 0],
....:    [ 0, 2, 0, 0, 0, 0,-1, 0, 1, 0, 1, 0, 0, 1],
....:    [ 0, 1, 1, 0, 0, 1, 0,-2, 1, 0, 0,-1, 1, 1],
....:    [ 0, 0, 0, 0, 1, 0,-1, 0, 0, 0, 0, 0, 0, 0]
....: ])
sage: pc = PointConfiguration(U.columns())
sage: pc.simplex_to_int([1,3,4,7,10,13])
1678
sage: pc.int_to_simplex(1678)
(1, 3, 4, 7, 10, 13)
>>> from sage.all import *
>>> U=matrix([
...    [ Integer(0), Integer(0), Integer(0), Integer(0), Integer(0), Integer(2), Integer(4),-Integer(1), Integer(1), Integer(1), Integer(0), Integer(0), Integer(1), Integer(0)],
...    [ Integer(0), Integer(0), Integer(0), Integer(1), Integer(0), Integer(0),-Integer(1), Integer(0), Integer(0), Integer(0), Integer(0), Integer(0), Integer(0), Integer(0)],
...    [ Integer(0), Integer(2), Integer(0), Integer(0), Integer(0), Integer(0),-Integer(1), Integer(0), Integer(1), Integer(0), Integer(1), Integer(0), Integer(0), Integer(1)],
...    [ Integer(0), Integer(1), Integer(1), Integer(0), Integer(0), Integer(1), Integer(0),-Integer(2), Integer(1), Integer(0), Integer(0),-Integer(1), Integer(1), Integer(1)],
...    [ Integer(0), Integer(0), Integer(0), Integer(0), Integer(1), Integer(0),-Integer(1), Integer(0), Integer(0), Integer(0), Integer(0), Integer(0), Integer(0), Integer(0)]
... ])
>>> pc = PointConfiguration(U.columns())
>>> pc.simplex_to_int([Integer(1),Integer(3),Integer(4),Integer(7),Integer(10),Integer(13)])
1678
>>> pc.int_to_simplex(Integer(1678))
(1, 3, 4, 7, 10, 13)