Polytope Class
This class handles all computations relating to lattice polytopes, such as the computation of lattice points and faces. When using reflexive polytopes, it also allows the computation of topological properties of the arising Calabi-Yau hypersurfaces that only depend on the polytope.
Constructor
cytools.polytope.Polytope
Description:
Constructs a Polytope
object describing a lattice polytope. This is
handled by the hidden __init__
function.
- CYTools only supports lattice polytopes, so any floating point numbers will be truncated to integers.
- The Polytope class is also imported to the root of the CYTools package,
so it can be imported from
cytools.polytope
or fromcytools
.
Arguments:
points
: A list of lattice points defining the polytope as their convex hull.labels
: A list of labels to specify the points. I.e., points[i] is labelled/accessed as labels[i]. If no labels are provided, then the points are given semi-arbitrary default labels.backend
: A string that specifies the backend used to construct the convex hull. The available options are "ppl", "qhull", or "palp". When not specified, it uses PPL for dimensions up to four, and palp otherwise.
Example
We construct two polytopes from lists of points.
from cytools import Polytope
p1 = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
print(p1)
# A 4-dimensional reflexive lattice polytope in ZZ^4
p2 = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[-1,-1,-1,0]])
print(p2)
# A 3-dimensional lattice polytope in ZZ^4
Functions
all_triangulations
Description: Computes all triangulations of the polytope using TOPCOM. There is the option to only compute fine, regular or fine triangulations.
Polytopes with more than around 15 points usually have too many triangulations, so this function may take too long or run out of memory.
Arguments:
only_fine
: Restricts to only fine triangulations.only_regular
: Restricts to only regular triangulations.only_star
: Restricts to only star triangulations. When not specified it defaults to True for reflexive polytopes and False otherwise.star_origin
: The index of the point that will be used as the star origin. If the polytope is reflexive this is set to 0, but otherwise it must be specified.include_points_interior_to_facets
: Whether to include points interior to facets from the triangulation.points
: List of point labels that will be used. Note that if this option is used then the parameterinclude_points_interior_to_facets
is ignored.backend
: The optimizer used to check regularity computation. The available options are the backends of theis_solid
function of theCone
class. If not specified, it will be picked automatically. Note that TOPCOM is not used to check regularity since it is much slower.as_list
: By default this function returns a generator object, which is usually desired for efficiency. However, this flag can be set to True so that it returns the full list of triangulations at once.raw_output
: Return the triangulations as lists of simplices instead of as Triangulation objects.
Returns:
A generator of Triangulation
objects, or a list of
Triangulation
objects if as_list
is set to True.
Example
We construct a polytope and find all of its triangulations. We try picking different restrictions and see how the number of triangulations changes.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-2,-1,-1],[-2,-1,-1,-1]])
g = p.all_triangulations()
next(g) # Takes some time while TOPCOM finds all the triangulations
# A fine, regular, star triangulation of a 4-dimensional point
# configuration with 7 points in ZZ^4
next(g) # Produces the next triangulation immediately
# A fine, regular, star triangulation of a 4-dimensional point
# configuration with 7 points in ZZ^4
len(p.all_triangulations(as_list=True)) # Number of fine, regular, star triangulations
# 2
len(p.all_triangulations(only_regular=False, only_star=False, only_fine=False, as_list=True) )# Number of triangularions, no matter if fine, regular, or star
# 6
ambient_dimension
Description: Returns the dimension of the ambient lattice.
Arguments: None.
Returns: The dimension of the ambient lattice.
Aliases:
ambient_dim
.
Example
We construct a polytope and check the dimension of the ambient lattice.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[-1,-1,-1,0]])
p.ambient_dimension()
# 4
automorphisms
Description: Returns the matrices that leave the polytope invariant. These matrices act on the points by multiplication on the right.
Arguments:
square_to_one
: Flag that restricts to only matrices that square to the identity.action
: Flag that specifies whether the returned matrices act on the left or the right. This option is ignored whenas_dictionary
is set to True.as_dictionary
: Return each automphism as a dictionary that describes the action on the indices of the points.
Returns: A list of automorphism matrices or dictionaries.
Example
We construct a polytope, and find its automorphisms. We also check that one of the non-trivial automorphisms is indeed an automorphism by checking that it acts as a permutation on the vertices. We also show how to get matrices that act on the left, which are simply the transpose matrices, and we show how to get dictionaries that describe how the indices of the points transform.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
autos = p.automorphisms()
a = autos[1]
print(a)
# [[ 1, 0, 0, 0],
# [-1, -1, -6, -9],
# [ 0, 0, 1, 0],
# [ 0, 0, 0, 1]])
print(f"{p.vertices()}\n{p.vertices().dot(a)}") # Print vertices before and after applying the automorphism
# [[ 1 0 0 0]
# [ 0 1 0 0]
# [ 0 0 1 0]
# [ 0 0 0 1]
# [-1 -1 -6 -9]]
# [[ 1 0 0 0]
# [-1 -1 -6 -9]
# [ 0 0 1 0]
# [ 0 0 0 1]
# [ 0 1 0 0]]
autos2 = p.automorphisms(square_to_one=True)
a2 = autos2[1]
print(f"{a2}\n{a2.dot(a2)}") # Print the automorphism and its square
# [[ 1 0 0 0]
# [-1 -1 -6 -9]
# [ 0 0 1 0]
# [ 0 0 0 1]]
# [[1 0 0 0]
# [0 1 0 0]
# [0 0 1 0]
# [0 0 0 1]]
autos_left = p.automorphisms(square_to_one=True, action="left")
print(autos_left[1].dot(p.vertices().T)) # The vertices are now columns
# [[ 1 -1 0 0 0]
# [ 0 -1 0 0 1]
# [ 0 -6 1 0 0]
# [ 0 -9 0 1 0]]
autos_dict = p.automorphisms(as_dictionary=True)
print(autos_dict)
# [{0: 0, 1: 1, 2: 2, 3: 3, 4: 4, 5: 5, 6: 6, 7: 7, 8: 8, 9: 9},
# {0: 0, 1: 4, 2: 2, 3: 3, 4: 1, 5: 5, 6: 6, 7: 7, 8: 8, 9: 9},
# {0: 0, 1: 1, 2: 2, 3: 3, 4: 5, 5: 4, 6: 6, 7: 7, 8: 8, 9: 9},
# {0: 0, 1: 4, 2: 2, 3: 3, 4: 5, 5: 1, 6: 6, 7: 7, 8: 8, 9: 9},
# {0: 0, 1: 5, 2: 2, 3: 3, 4: 1, 5: 4, 6: 6, 7: 7, 8: 8, 9: 9},
# {0: 0, 1: 5, 2: 2, 3: 3, 4: 4, 5: 1, 6: 6, 7: 7, 8: 8, 9: 9}]
backend
Description: Returns the backend.
Arguments: None.
Returns: The computational backend
chi
Description: Computes the Euler characteristic of the Calabi-Yau obtained as the anticanonical hypersurface in the toric variety given by a desingularization of the face or normal fan of the polytope when the lattice is specified as "N" or "M", respectively.
Only reflexive polytopes of dimension 2-5 are currently supported.
Arguments:
lattice
: Specifies the lattice on which the polytope is defined. Options are "N" and "M".
Returns: The Euler characteristic of the arising Calabi-Yau manifold.
Example
We construct a polytope and compute the Euler characteristic of the associated hypersurfaces.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
p.chi(lattice="N")
# -540
p.chi(lattice="M")
# 540
clear_cache
Description: Clears the cached results of any previous computation.
Arguments: None.
Returns: Nothing.
Example
We compute the lattice points of a large polytope.
p = Polytope([[-1,-1,-1,-1,-1],[3611,-1,-1,-1,-1],[-1,42,-1,-1,-1],[-1,-1,6,-1,-1],[-1,-1,-1,2,-1],[-1,-1,-1,-1,1]])
pts = p.points() # Takes a few seconds
pts = p.points() # It runs instantly because the result is cached
p.clear_cache() # Clears the results of any previos computation
pts = p.points() # Again it takes a few seconds since the cache was cleared
dimension
Description: Returns the dimension of the polytope.
Arguments: None.
Returns: The dimension of the polytope.
Aliases:
dim
.
Example
We construct a polytope and check its dimension.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[-1,-1,-1,0]])
p.dimension()
# 3
dual_polytope
Description: Returns the dual polytope (also called polar polytope). Only lattice polytopes are currently supported, so only duals of reflexive polytopes can be computed.
If is a lattice polytope, the dual polytope of is . A lattice polytope is reflexive if its dual is also a lattice polytope.
Arguments: None.
Returns: The dual polytope.
Aliases:
dual
, polar_polytope
, polar
.
Example
We construct a reflexive polytope and find its dual. We then verify that the dual of the dual is the original polytope.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
p_dual = p.dual_polytope()
print(p_dual)
# A 4-dimensional reflexive lattice polytope in ZZ^4
p_dual_dual = p_dual.dual_polytope()
p_dual_dual is p
# True
faces
Description: Computes the faces of a polytope.
When the polytope is 4-dimensional it calls the slightly more optimized
_faces4d()
function.
Arguments:
d
: Optional parameter that specifies the dimension of the desired faces.
Returns:
A tuple of PolytopeFace
objects of dimension d, if
specified. Otherwise, a tuple of tuples of
PolytopeFace
objects organized in ascending
dimension.
Example
We show that this function returns a tuple of 2-faces if d
is set to
2. Otherwise, the function returns all faces in tuples organized in
ascending dimension. We verify that the first element in the tuple of
2-faces is the same as the first element in the corresponding subtuple
in the tuple of all faces.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
faces2d = p.faces(2)
allfaces = p.faces()
print(faces2d[0]) # Print first face in tuple of 2-faces
# A 2-dimensional face of a 4-dimensional polytope in ZZ^4
faces2d[0] is allfaces[2][0]
# True
facets
Description: Returns the facets (codimension-1 faces) of the polytope.
Arguments: None.
Returns:
A list of PolytopeFace
objects of codimension 1.
Example
We construct a polytope and find its facets.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
facets = p.facets()
find_2d_reflexive_subpolytopes
Description: Use the algorithm by Huang and Taylor described in 1907.09482 to find 2D reflexive subpolytopes in 4D polytopes.
Arguments: None.
Returns: The list of 2D reflexive subpolytopes.
Example
We construct a polytope and find its 2D reflexive subpolytopes.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
p.find_2d_reflexive_subpolytopes()
# [A 2-dimensional lattice polytope in ZZ^4]
glsm_basis
Description: Computes a basis of columns of the GLSM charge matrix.
Arguments:
include_origin
: Indicates whether to use the origin in the calculation. This corresponds to the inclusion of the canonical divisor.include_points_interior_to_facets
: By default only boundary points not interior to facets are used. If this flag is set to true then points interior to facets are also used.points
: The list of indices of the points that will be used. Note that if this option is used then the parametersinclude_origin
andinclude_points_interior_to_facets
are ignored. Also, note that the indices returned here will be the indices of the sorted list of points.integral
: Indicates whether to find an integral basis for the columns of the GLSM charge matrix. (i.e. so that remaining columns can be written as an integer linear combination of the basis elements.)
Returns: A list of column indices that form a basis.
Example
We construct a polytope, find its GLSM charge matrix and a basis of columns.
import numpy as np
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
p.glsm_basis()
# array([1, 6])
glsm = p.glsm_charge_matrix()
np.linalg.matrix_rank(glsm) == np.linalg.matrix_rank(glsm[:,p.glsm_basis()]) # This shows that the columns form a basis
# True
glsm_charge_matrix
Description: Computes the GLSM charge matrix of the theory resulting from this polytope.
Arguments:
include_origin
: Indicates whether to use the origin in the calculation. This corresponds to the inclusion of the canonical divisor.include_points_interior_to_facets
: By default only boundary points not interior to facets are used. If this flag is set to true then points interior to facets are also used.points
: The list of indices of the points that will be used. Note that if this option is used then the parametersinclude_origin
andinclude_points_interior_to_facets
are ignored.integral
: Indicates whether to find an integral basis for the columns of the GLSM charge matrix. (i.e. so that remaining columns can be written as an integer linear combination of the basis elements.)
Returns: The GLSM charge matrix.
Example
We construct a polytope and find the GLSM charge matrix with different parameters.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
p.glsm_charge_matrix()
# array([[-18, 1, 9, 6, 1, 1, 0],
# [ -6, 0, 3, 2, 0, 0, 1]])
p.glsm_charge_matrix().dot(p.points_not_interior_to_facets()) # By definition this product must be zero
# array([[0, 0, 0, 0],
# [0, 0, 0, 0]])
p.glsm_charge_matrix(include_origin=False) # Excludes the canonical divisor
# array([[1, 9, 6, 1, 1, 0],
# [0, 3, 2, 0, 0, 1]])
p.glsm_charge_matrix(include_points_interior_to_facets=True) # Includes points interior to facets
# array([[-18, 1, 9, 6, 1, 1, 0, 0, 0, 0],
# [ -6, 0, 3, 2, 0, 0, 1, 0, 0, 0],
# [ -4, 0, 2, 1, 0, 0, 0, 1, 0, 0],
# [ -3, 0, 1, 1, 0, 0, 0, 0, 1, 0],
# [ -2, 0, 1, 0, 0, 0, 0, 0, 0, 1]])
glsm_linear_relations
Description: Computes the linear relations of the GLSM charge matrix.
Arguments:
include_origin
: Indicates whether to use the origin in the calculation. This corresponds to the inclusion of the canonical divisor.include_points_interior_to_facets
: By default only boundary points not interior to facets are used. If this flag is set to true then points interior to facets are also used.points
: The list of indices of the points that will be used. Note that if this option is used then the parametersinclude_origin
andinclude_points_interior_to_facets
are ignored.integral
: Indicates whether to find an integral basis for the columns of the GLSM charge matrix. (i.e. so that remaining columns can be written as an integer linear combination of the basis elements.)
Returns: A matrix of linear relations of the columns of the GLSM charge matrix.
Example
We construct a polytope and find its GLSM charge matrix and linear relations with different parameters.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
p.glsm_linear_relations()
# array([[ 1, 1, 1, 1, 1, 1, 1],
# [ 0, 9, -1, 0, 0, 0, 3],
# [ 0, 6, 0, -1, 0, 0, 2],
# [ 0, 1, 0, 0, -1, 0, 0],
# [ 0, 1, 0, 0, 0, -1, 0]])
p.glsm_linear_relations().dot(p.glsm_charge_matrix().T) # By definition this product must be zero
# array([[0, 0],
# [0, 0],
# [0, 0],
# [0, 0],
# [0, 0]])
p.glsm_linear_relations(include_origin=False) # Excludes the canonical divisor
# array([[ 9, -1, 0, 0, 0, 3],
# [ 6, 0, -1, 0, 0, 2],
# [ 1, 0, 0, -1, 0, 0],
# [ 1, 0, 0, 0, -1, 0]])
p.glsm_linear_relations(include_points_interior_to_facets=True) # Includes points interior to facets
# array([[ 1, 1, 1, 1, 1, 1, 1, 1, 1, 1],
# [ 0, 9, -1, 0, 0, 0, 3, 2, 1, 1],
# [ 0, 6, 0, -1, 0, 0, 2, 1, 1, 0],
# [ 0, 1, 0, 0, -1, 0, 0, 0, 0, 0],
# [ 0, 1, 0, 0, 0, -1, 0, 0, 0, 0]])
hpq
Description: Returns the Hodge number of the Calabi-Yau obtained as the anticanonical hypersurface in the toric variety given by a desingularization of the face or normal fan of the polytope when the lattice is specified as "N" or "M", respectively.
- Only reflexive polytopes of dimension 2-5 are currently supported.
Arguments:
p
: The holomorphic index of the Dolbeault cohomology of interest.q
: The anti-holomorphic index of the Dolbeault cohomology of interest.lattice
: Specifies the lattice on which the polytope is defined. Options are "N" and "M".
Returns: The Hodge number of the arising Calabi-Yau manifold.
Example
We construct a polytope and check some Hodge numbers of the associated hypersurfaces.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
p.hpq(0,0,lattice="N")
# 1
p.hpq(0,1,lattice="N")
# 0
p.hpq(1,1,lattice="N")
# 2
p.hpq(1,2,lattice="N")
# 272
inequalities
Description: Returns the inequalities giving the hyperplane representation of the polytope. The inequalities are given in the form
.
Note, however, that equalities are not included.
Arguments: None.
Returns: The inequalities defining the polytope.
Example
We construct a polytope and find the defining inequalities.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
p.inequalities()
# array([[ 4, -1, -1, -1, 1],
# [-1, 4, -1, -1, 1],
# [-1, -1, 4, -1, 1],
# [-1, -1, -1, 4, 1],
# [-1, -1, -1, -1, 1]])
is_affinely_equivalent
Description: Returns True if the polytopes can be transformed into each other by an integral affine transformation.
Arguments:
other
: The other polytope being compared.backend
: Selects which backend to use to compute the normal form. Options are "native", which uses native python code, or "palp", which uses PALP for the computation.
Returns: The truth value of the polytopes being affinely equivalent.
Example
We construct two polytopes and check if they are affinely equivalent.
p1 = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
p2 = Polytope([[1,0,0,1],[0,1,0,1],[0,0,1,1],[0,0,0,2],[-1,-1,-1,0]])
p1.is_affinely_equivalent(p2)
# True
is_favorable
Description: Returns True if the Calabi-Yau hypersurface arising from this polytope is favorable (i.e. all Kahler forms descend from Kahler forms on the ambient toric variety) and False otherwise.
Only reflexive polytopes of dimension 2-5 are currently supported.
Arguments:
lattice
: Specifies the lattice on which the polytope is defined. Options are "N" and "M".
The truth value of the polytope being favorable.
Example
We construct two reflexive polytopes and find whether they are favorable when considered in the N lattice.
p1 = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
p2 = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-3,-6]])
p1.is_favorable(lattice="N")
# True
p2.is_favorable(lattice="N")
# False
is_linearly_equivalent
Description: Returns True if the polytopes can be transformed into each other by an transformation.
Arguments:
other
: The other polytope being compared.backend
: Selects which backend to use to compute the normal form. Options are "native", which uses native python code, or "palp", which uses PALP for the computation.
Returns: The truth value of the polytopes being linearly equivalent.
Example
We construct two polytopes and check if they are linearly equivalent.
p1 = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
p2 = Polytope([[-1,0,0,0],[0,-1,0,0],[0,0,-1,0],[0,0,0,-1],[1,1,1,1]])
p1.is_linearly_equivalent(p2)
# True
is_reflexive
Description: Returns True if the polytope is reflexive and False otherwise.
Arguments:
allow_translations
: Whether to allow the polytope to be translated.
Returns: The truth value of the polytope being reflexive.
Example
We construct a polytope and check if it is reflexive.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
p.is_reflexive()
# True
is_solid
Description: Returns True if the polytope is solid (i.e. full-dimensional) and False otherwise.
Arguments: None.
Returns: The truth value of the polytope being full-dimensional.
Example
We construct a polytope and check if it is solid.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[-1,-1,-1,0]])
p.is_solid()
# False
labels
Description: Returns the point labels, in order.
minkowski_sum
Description: Returns the Minkowski sum of the two polytopes.
Arguments:
other
: The other polytope used for the Minkowski sum.
Returns: The Minkowski sum.
Example
We construct two polytops and compute their Minkowski sum.
p1 = Polytope([[1,0,0],[0,1,0],[-1,-1,0]])
p2 = Polytope([[0,0,1],[0,0,-1]])
p1.minkowski_sum(p2)
# A 3-dimensional reflexive lattice polytope in ZZ^3
nef_partitions
Description: Computes the nef partitions of the polytope using PALP.
This is currently an experimental feature and may change significantly in future versions.
Arguments:
keep_symmetric
: Keep symmetric partitions related by lattice automorphisms.keep_products
: Keep product partitions corresponding to complete intersections being direct products.keep_projections
: Keep projection partitions, i.e. partitions where one of the parts consists of a single vertex.codim
: The number of parts in the partition or, equivalently, the codimension of the complete intersection Calabi-Yau.compute_hodge_numbers
: Indicates whether Hodge numbers of the CICY are computed.return_hodge_numbers
: Indicates whether to return the Hodge numbers along with the nef partitions. They are returned in a separate tuple and they are ordered as in the Hodge diamond from top to bottom and left to right.
Returns: The nef partitions of the polytope. If return_hodge_numbers is set to True then two tuples are returned, one with the nef partitions and one with the corresponding Hodge numbers.
Example
We construct a tesseract and find the 2- and 3-part nef partitions.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,0,0,0],[0,-1,0,0],[0,0,-1,0],[0,0,0,-1]])
nef_part_2 = p.nef_partitions() # Default codimension is 2
print(nef_part_2[0]) # Print the first of the nef partitions
# ((5, 2, 3, 4), (8, 7, 6, 1))
nef_part_3 = p.nef_partitions(codim=3) # Codimension 3
print(nef_part_3[0]) # Print the first of the nef partitions
# ((6, 5, 3), (2, 4), (8, 7, 1))
normal_form
Description: Returns the normal form of the polytope as defined by Kreuzer-Skarke.
Arguments:
affine_transform
: Flag that determines whether to only use transformations or also allow translations.backend
: Selects which backend to use. Options are "native", which uses native python code, or "palp", which uses PALP for the computation. There is a different convention for affine normal forms between the native algorithm and PALP, and PALP generally works better.
Returns: The list of vertices in normal form.
Example
We construct a polytope, and find its linear and affine normal forms.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
p.normal_form()
# array([[ 1, 0, 0, 0],
# [ 0, 1, 0, 0],
# [ 0, 0, 1, 0],
# [ 0, 0, 0, 1],
# [-9, -6, -1, -1]])
p.normal_form(affine_transform=True)
# array([[ 1, 0, 0, 0],
# [ 0, 1, 0, 0],
# [ 0, 0, 1, 0],
# [12, 17, 17, 18],
# [ 0, 0, 0, 0]])
points
Description: Returns the lattice points of the polytope.
Points are sorted so that interior points are first, and then the rest are arranged by decreasing number of saturated inequalities and lexicographically. For reflexive polytopes this is useful since the origin will be at index 0 and boundary points interior to facets will be last.
Arguments:
which
: Which points to return. Specified by a (list of) labels. NOT INDICES!!!optimal
: Whether to return the points in their optimal coordinates.as_indices
: Return the points as indices of the full list of points of the polytope.
Returns: The list of lattice points of the polytope.
Aliases:
pts
.
Example
We construct a polytope and compute the lattice points. One can verify that the first point is the only interior point, and the last three points are the ones interior to facets. Thus it follows the aforementioned ordering.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
p.points()
# array([[ 0, 0, 0, 0],
# [-1, -1, -6, -9],
# [ 0, 0, 0, 1],
# [ 0, 0, 1, 0],
# [ 0, 1, 0, 0],
# [ 1, 0, 0, 0],
# [ 0, 0, -2, -3],
# [ 0, 0, -1, -2],
# [ 0, 0, -1, -1],
# [ 0, 0, 0, -1]])
points_to_indices
Description: Returns the list of indices corresponding to the given points. It also accepts a single point, in which case it returns the corresponding index.
Arguments:
points
: A point or a list of points.is_optimal
: Whether the points argument represents points in the optimal (True) or input (False) basis
Returns: The list of indices corresponding to the given points, or the index of the point if only one is given.
Example
We construct a polytope and find the indices of some of its points.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
p.points_to_indices([-1,-1,-6,-9]) # We input a single point, so a single index is returned
# 1
p.points_to_indices([[-1,-1,-6,-9],[0,0,0,0],[0,0,1,0]]) # We input a list of points, so a list of indices is returned
# array([1, 0, 3])
points_to_labels
Description: Returns the list of labels corresponding to the given points. It also accepts a single point, in which case it returns the corresponding label.
Arguments:
points
: A point or a list of points.is_optimal
: Whether the points argument represents points in the optimal (True) or input (False) basis
Returns: The list of labels corresponding to the given points, or the label of the point if only one is given.
random_triangulations_fair
Description: Constructs pseudorandom regular (optionally star) triangulations of a given point set. Implements Algorithm #3 from the paper Bounding the Kreuzer-Skarke Landscape by Mehmet Demirtas, Liam McAllister, and Andres Rios-Tascon. arXiv:2008.01730
This is a Markov chain Monte Carlo algorithm that involves taking random walks inside the subset of the secondary fan corresponding to fine triangulations and performing random flips. For details, please see Section 4.1 in the paper.
- By default this function tries to guess reasonable parameters for the polytope that is being used. However, it may take some trial-and-error to find optimal parameters that produce a good sampling and prevent the algorithm from getting stuck.
- This function is designed mainly for large polytopes where sampling triangulations is challenging. When small polytopes are used it is likely to get stuck or not produce an optimal set.
Arguments:
N
: Number of desired unique triangulations. If not specified, it will generate as many triangulations as it can find until it has to retry more thanmax_retries
times to obtain a new triangulation. This parameter is required when settingas_list
to True.n_walk
: Number of hit-and-run steps per triangulation. Default is n_points//10+10.n_flip
: Number of random flips performed per triangulation. Default is n_points//10+10.initial_walk_steps
: Number of hit-and-run steps to take before starting to record triangulations. Small values may result in a bias towards Delaunay-like triangulations. Default is 2*n_pts//10+10.walk_step_size
: Determines the size of random steps taken in the secondary fan. The algorithm may stall if too small.max_steps_to_wall
: Maximum number of steps to take towards a wall of the subset of the secondary fan that correspond to fine triangulations. If a wall is not found, a new random direction is selected. Setting this to be very large (>100) reduces performance. If this is set to be too low, the algorithm may stall.fine_tune_steps
: Number of steps to determine the location of a wall. Decreasing improves performance, but might result in biased samples.max_retries
: Maximum number of attempts to obtain a new triangulation before the process is terminated.make_star
: Converts the obtained triangulations into star triangulations. If not specified, defaults to True for reflexive polytopes, and False for other polytopes.include_points_interior_to_facets
: Whether to include points interior to facets from the triangulation. If not specified, it is set to False for reflexive polytopes and True otherwise.points
: List of point labels that will be used. Note that if this option is used then the parameterinclude_points_interior_to_facets
is ignored.backend
: Specifies the backend used to compute the triangulation. The available options are "cgal" and "qhull".as_list
: By default this function returns a generator object, which is usually desired for efficiency. However, this flag can be set to True so that it returns the full list of triangulations at once.progress_bar
: Shows number of triangulations obtained and progress bar. Note that this option is only available when returning a list instead of a generator.seed
: A seed for the random number generator. This can be used to obtain reproducible results.
Returns:
A generator of Triangulation
objects, or a list of
Triangulation
objects if as_list
is set to True.
Example
We construct a polytope and find some random triangulations. The
computation takes considerable time, but they should be a fair sample
from the full set of triangulations (if the parameters are chosen
correctly). For (some) machine learning purposes or when the fairness
of the sample is not crucial, the
random_triangulations_fast
function
should be used instead.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]]).dual()
g = p.random_triangulations_fast()
next(g) # Takes a long time (around a minute)
# A fine, regular, star triangulation of a 4-dimensional point
# configuration with 106 points in ZZ^4
next(g) # Takes slightly shorter (still around a minute)
# A fine, regular, star triangulation of a 4-dimensional point
# configuration with 106 points in ZZ^4
rand_triangs = p.random_triangulations_fair(N=10, as_list=True) # Produces the list of 10 triangulations, but takes a long time (around 10 minutes)
It is worth noting that the time it takes to obtain each triangulation varies very significantly on the parameters used. The function tries to guess reasonable parameters, but it is better to adjust them to your desired balance between speed and fairness of the sampling.
random_triangulations_fast
Description: Constructs pseudorandom regular (optionally fine and star) triangulations of a given point set. This is done by picking random heights around the Delaunay heights from a Gaussian distribution.
This function produces random triangulations very quickly, but it does
not produce a fair sample. When a fair sampling is required the
random_triangulations_fair
function should be used.
Arguments:
N
: Number of desired unique triangulations. If not specified, it will generate as many triangulations as it can find until it has to retry more thanmax_retries
times to obtain a new triangulation. This parameter is required when settingas_list
to True.c
: A contant used as the standard deviation of the Gaussian distribution used to pick the heights. A largerc
results in a wider range of possible triangulations, but with a larger fraction of them being non-fine, which slows down the process whenonly_fine
is set to True.max_retries
: Maximum number of attempts to obtain a new triangulation before the process is terminated.make_star
: Converts the obtained triangulations into star triangulations. If not specified, defaults to True for reflexive polytopes, and False for other polytopes.only_fine
: Restricts to fine triangulations.include_points_interior_to_facets
: Whether to include points interior to facets from the triangulation. If not specified, it is set to False for reflexive polytopes and True otherwise.points
: List of point labels that will be used. Note that if this option is used then the parameterinclude_points_interior_to_facets
is ignored.backend
: Specifies the backend used to compute the triangulation. The available options are "cgal" and "qhull".as_list
: By default this function returns a generator object, which is usually desired for efficiency. However, this flag can be set to True so that it returns the full list of triangulations at once.progress_bar
: Shows the number of triangulations obtained and progress bar. Note that this option is only available when returning a list instead of a generator.seed
: A seed for the random number generator. This can be used to obtain reproducible results.
Returns:
A generator of Triangulation
objects, or a list of
Triangulation
objects if as_list
is set to True.
Example
We construct a polytope and find some random triangulations. The
triangulations are obtained very quicly, but they are not a fair sample
of the space of triangulations. For a fair sample, the
random_triangulations_fair
function
should be used.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]]).dual()
g = p.random_triangulations_fast()
next(g) # Runs very quickly
# A fine, regular, star triangulation of a 4-dimensional point configuration with 106 points in ZZ^4
next(g) # Keeps producing triangulations until it has trouble finding more
# A fine, regular, star triangulation of a 4-dimensional point configuration with 106 points in ZZ^4
rand_triangs = p.random_triangulations_fast(N=10, as_list=True) # Produces the list of 10 triangulations very quickly
triangulate
Description: Returns a single regular triangulation of the polytope.
When reflexive polytopes are used, it defaults to returning a fine, regular, star triangulation.
Arguments:
heights
: A list of heights specifying the regular triangulation. When not specified, it will return the Delaunay triangulation when using CGAL, a triangulation obtained from random heights near the Delaunay when using QHull, or the placing triangulation when using TOPCOM. Heights can only be specified when using CGAL or QHull as the backend.make_star
: Indicates whether to turn the triangulation into a star triangulation by deleting internal lines and connecting all points to the origin, or equivalently by decreasing the height of the origin to be much lower than the rest. By default, this flag is set to true if the polytope is reflexive and neither heights or simplices are inputted. Otherwise, it is set to False.include_points_interior_to_facets
: Whether to include points interior to facets from the triangulation. If not specified, it is set to False for reflexive polytopes and True otherwise.points
: List of point labels that will be used. Note that if this option is used then the parameterinclude_points_interior_to_facets
is ignored.simplices
: A list of simplices specifying the triangulation. This is useful when a triangulation was previously computed and it needs to be used again. Note that the order of the points needs to be consistent with the order that thePolytope
class uses.check_input_simplices
: Flag that specifies whether to check if the input simplices define a valid triangulation.backend
: Specifies the backend used to compute the triangulation. The available options are "qhull", "cgal", and "topcom". CGAL is the default one as it is very fast and robust.
Returns:
A Triangulation
object describing a triangulation
of the polytope.
Example
We construct a triangulation of a reflexive polytope and check that by default it is a fine, regular, star triangulation. We also try constructing triangulations with heights, input simplices, and using the other backends.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-2,-1,-1],[-2,-1,-1,-1]])
p.triangulate()
# A fine, regular, star triangulation of a 4-dimensional polytope in
# ZZ^4
p.triangulate(heights=[-30,5,5,24,-19,-14,29])
# A fine, regular, star triangulation of a 4-dimensional polytope in
# ZZ^4
p.triangulate(simplices=[[0,1,2,3,4],[0,1,2,3,5],[0,1,2,4,6],[0,1,2,5,6],[0,1,3,4,5],[0,1,4,5,6],[0,2,3,4,5],[0,2,4,5,6]])
# A fine, regular, star triangulation of a 4-dimensional polytope in
# ZZ^4
p.triangulate(backend="qhull")
# A fine, regular, star triangulation of a 4-dimensional polytope in
# ZZ^4
p.triangulate(backend="topcom")
# A fine, regular, star triangulation of a 4-dimensional polytope in
# ZZ^4
vertices
Description: Returns the vertices of the polytope.
Arguments:
optimal
: Whether to return the points in their optimal coordinates.as_indices
: Return the points as indices of the full list of points of the polytope.
Returns: The list of vertices of the polytope.
Example
We construct a polytope and find its vertices. We can see that they match the points that we used to construct the polytope.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
p.vertices()
# array([[ 1, 0, 0, 0],
# [ 0, 1, 0, 0],
# [ 0, 0, 1, 0],
# [ 0, 0, 0, 1],
# [-1, -1, -1, -1]])
volume
Description: Returns the volume of the polytope.
By convention, the standard simplex has unit volume. To get the more typical Euclidean volume it must be multiplied by .
Arguments: None.
Returns: The volume of the polytope.
Example
We construct a standard simplex and a cube, and find their volumes.
p1 = Polytope([[1,0,0],[0,1,0],[0,0,1],[0,0,0]])
p2 = Polytope([[1,0,0],[0,1,0],[0,0,1],[0,0,0],[0,1,1],[1,0,1],[1,1,0],[1,1,1]])
p1.volume()
# 1
p2.volume()
# 6
Hidden Functions
__add__
Description:
Implements addition of polytopes with the
minkowski_sum
function.
Arguments:
other
: The other polytope used for the Minkowski sum.
Returns: The Minkowski sum.
Example
We construct two polytops and compute their Minkowski sum.
p1 = Polytope([[1,0,0],[0,1,0],[-1,-1,0]])
p2 = Polytope([[0,0,1],[0,0,-1]])
p1 + p2
# A 3-dimensional reflexive lattice polytope in ZZ^3
__eq__
Description: Implements comparison of polytopes with ==.
Arguments:
other
: The other polytope that is being compared.
Returns: The truth value of the polytopes being equal.
Example
We construct two polytopes and compare them.
p1 = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
p2 = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
p1 == p2
# True
__getstate__
Description: Gets the state of the class instance, for pickling.
Arguments: None
Returns: Nothing.
__hash__
Description: Implements the ability to obtain hash values from polytopes.
Arguments: None.
Returns: The hash value of the polytope.
Example
We compute the hash value of a polytope. Also, we construct a set and a dictionary with a polytope, which make use of the hash function.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
h = hash(p) # Obtain hash value
d = {p: 1} # Create dictionary with polytope keys
s = {p} # Create a set of polytopes
__init__
Description:
Initializes a Polytope
object describing a lattice polytope.
CYTools only supports lattice polytopes, so any floating point numbers will be truncated to integers.
Arguments:
points
: A list of lattice points defining the polytope as their convex hull.labels
: A list of labels to specify the points. I.e., points[i] is labelled/accessed as labels[i]. If no labels are provided, then the points are given semi-arbitrary default labels.backend
: A string that specifies the backend used to construct the convex hull. The available options are "ppl", "qhull", or "palp". When not specified, it uses PPL for dimensions up to four, and palp otherwise.
Returns: Nothing.
Example
This is the function that is called when creating a new Polytope
object. Thus, it is used in the following example.
from cytools import Polytope
p1 = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
print(p1)
# A 4-dimensional reflexive lattice polytope in ZZ^4
p2 = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[-1,-1,-1,0]])
print(p2)
# A 3-dimensional lattice polytope in ZZ^4
__ne__
Description: Implements comparison of polytopes with !=.
Arguments:
other
: The other polytope that is being compared.
Returns: The truth value of the polytopes being different.
Example
We construct two polytopes and compare them.
p1 = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
p2 = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
p1 != p2
# False
__repr__
Description: Returns an umabiguous string describing the polytope.
Arguments: None.
Returns: A string describing the polytope.
Example
This function can be used to convert the polytope to a string or to print information about the polytope.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
print(repr(p)) # Prints polytope info
# A 4-dimensional reflexive lattice polytope in ZZ^4
__setstate__
Description: Gets the state of the class instance, for pickling.
Arguments:
state
: The dictionary of the instance state, read from pickle.
Returns: Nothing.
__str__
Description: Returns a human-readable string describing the polytope.
Arguments: None.
Returns: A string describing the polytope.
Example
This function can be used to convert the polytope to a string or to print information about the polytope.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
poly_info = str(p) # Converts to string
print(p) # Prints polytope info
# A 4-dimensional reflexive lattice polytope in ZZ^4
_dump
Description: Get every class variable.
No copying is done, so be careful with mutability!
Arguments: None.
Returns: Dictionary mapping variable name to value.
_faces4d
Description: Computes the faces of a 4D polytope.
Arguments: None.
Returns: A tuple of tuples of faces organized in ascending dimension.
Example
We construct a 4D polytope and compute its faces. Since this function
generally should not be directly used, we do this with the
faces
function.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
allfaces = p.faces() # the _faces4d function is used since it is a 4d polytope
print(allfaces[1][0]) # Print the first face in the tuple of 1-dimensional faces
# A 1-dimensional face of a 4-dimensional polytope in ZZ^4
_optimal_to_input
Description: We internally store the points in an 'optimal' representation (translated, LLL-reduced, ...). This eases computations. We will always want to return answers in original representation. This function performs the mapping.
This is a kind of costly map so, whenever possible, use _optimalpts2labels and _labels2inputPts
Arguments:
pts_opt
: The points in the optimal representation.
Returns: The points in the original representation.
_process_points
Description: Internal function for processing input points. Should only be called once (in the initializer). Abstracted here to clarify logic.
Sets: self._transl_vector, self._transf_mat_inv, self._poly_optimal, self._ineqs_optimal, self._ineqs_input, self._labels2optPts, self._labels2inputPts, self._pts_saturating, self._pts_order, and self._pts_indices
Arguments:
pts_input
: The points input from the user.labels
: The point labels input from the user.
Returns: Nothing.
_triang_labels
Description: Constructs the list of point labels of the points that will be used in a triangulation.
Typically this function should not be called by the user. Instead, it is called by various other functions in the Polytope class.
Arguments:
include_points_interior_to_facets
: Whether to include points interior to facets from the triangulation. If not specified, it is set to False for reflexive polytopes and True otherwise.
Returns: A tuple of the indices of the points that will be included in a triangulation
Example
We construct triangulations in various ways. We use the
triangulate
function instead of using this function
directly.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
t1 = p.triangulate()
print(t1)
# A fine, regular, star triangulation of a 4-dimensional point
# configuration with 7 points in ZZ^4
t2 = p.triangulate(include_points_interior_to_facets=True)
print(t2)
# A fine, regular, star triangulation of a 4-dimensional point
# configuration with 10 points in ZZ^4
t3 = p.triangulate(points=[1,2,3,4,5])
print(t3)
# A fine, regular, non-star triangulation of a 4-dimensional point
# configuration with 5 points in ZZ^4