Triangulation Class
This class handles triangulations of lattice polytopes. It can compute various properties of the triangulation, as well as construct a ToricVariety or CalabiYau object if the triangulation is suitable.
Generally, objects of this class should not be constructed directly by the
end user. Instead, they should be created by various functions of the
Polytope
class.
Constructor
cytools.triangulation.Triangulation
Description:
Constructs a Triangulation
object describing a triangulation of a lattice
polytope. This is handled by the hidden __init__
function.
If you construct a triangulation object directly by inputting a list of
points they may be reordered to match the ordering of the points from the
Polytope
class. This is to ensure that computations of
toric varieties and Calabi-Yau manifolds are correct. To avoid this
subtlety we discourage users from constructing Triangulation objects
directly, and instead use the triangulation functions in the
Polytope
class.
Arguments:
poly
: The ambient polytope of the points to be triangulated. If not specified, it's constructed as the convex hull of the given points.pts
: The list of points to be triangulated. Specified by labels.heights
: The heights specifying the regular triangulation. When not specified, construct based off of the backend:Heights can only be specified when using CGAL or QHull as the backend.- (CGAL) Delaunay triangulation,
- (QHULL) triangulation from random heights near Delaunay, or
- (TOPCOM) placing triangulation.make_star
: 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 until it is much lower than all other heights.simplices
: Array-like of simplices specifying the triangulation. Each simplex is a list of point labels. This is useful when a triangulation was previously computed and it needs to be used again. Note that the ordering of the points needs to be consistent.check_input_simplices
: Whether to check if the input simplices define a valid triangulation.backend
: The backend used to compute the triangulation. Options are "qhull", "cgal", and "topcom". CGAL is the default as it is very fast and robust.verbosity
: The verbosity level.
Example
We construct a triangulation of a polytope. Since this class is not
intended to by initialized by the end user, we create it via the
triangulate
function of the
Polytope
class. In this example the polytope is reflexive,
so by default the triangulation is fine, regular, and star. Also, since the
polytope is reflexive then by default only the lattice points not interior
to facets are included in the triangulation.
from cytools import Polytope
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
t = p.triangulate()
print(t)
# A fine, regular, star triangulation of a 4-dimensional point
# configuration with 7 points in ZZ^4
Functions
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 triangulation and find its ambient dimension.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
t = p.triangulate()
t.ambient_dim()
# 4
automorphism_orbit
Description: Returns all of the triangulations of the polytope that can be obtained by applying one or more polytope automorphisms to the triangulation.
It also has the option of restricting the simplices to faces of the polytope of a particular dimension or codimension. This restriction is useful for checking CY equivalences from different triangulations.
Depending on how the point configuration was constructed, it may be the case that the automorphism group of the point configuration is larger or smaller than the one from the polytope. This function only uses the subset of automorphisms of the polytope that are also automorphisms of the point configuration.
Arguments:
automorphism
: The index or list of indices of the polytope automorphisms to use. If not specified it uses all automorphisms.on_faces_dim
: Restrict the simplices to faces of the polytope of a given dimension.on_faces_codim
: Restrict the simplices to faces of the polytope of a given codimension.
Returns: The list of triangulations obtained by performing automorphisms transformations.
Example
We construct a triangulation and find some of its automorphism orbits.
p = Polytope([[-1,0,0,0],[-1,1,0,0],[-1,0,1,0],[2,-1,0,-1],[2,0,-1,-1],[2,-1,-1,-1],[-1,0,0,1],[-1,1,0,1],[-1,0,1,1]])
t = p.triangulate()
orbit_all_autos = t.automorphism_orbit()
print(len(orbit_all_autos))
# 36
orbit_all_autos_2faces = t.automorphism_orbit(on_faces_dim=2)
print(len(orbit_all_autos_2faces))
# 36
orbit_sixth_auto = t.automorphism_orbit(automorphism=5)
print(len(orbit_sixth_auto))
# 2
orbit_list_autos = t.automorphism_orbit(automorphism=[5,6,9])
print(len(orbit_list_autos))
# 12
check_heights
Description: Check if the heights uniquely define a triangulation. That is, if they do not lie on a wall of the generated secondary cone.
If heights don't uniquely correspond to a trinagulation, delete the heights so they are re-calculated (keep triangulation, though...)
Arguments:
verbosity
: The verbosity level.
Returns: Nothing.
clear_cache
Description: Clears the cached results of any previous computation.
Arguments:
recursive
: Whether to also clear the cache of the ambient polytope.
Returns: Nothing.
Example
We construct a triangulation, compute its GKZ vector, clear the cache and then compute it again.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
t = p.triangulate()
gkz_phi = t.gkz_phi() # Computes the GKZ vector
t.clear_cache()
gkz_phi = t.gkz_phi() # The GKZ vector is recomputed
dimension
Description: Returns the dimension of the triangulated point configuration.
Arguments: None.
Returns: The dimension of the triangulation.
Aliases:
dim
.
Example
We construct a triangulation and find its dimension.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
t = p.triangulate()
t.dimension()
# 4
get_cy
Description: Returns a CalabiYau object corresponding to the anti-canonical hypersurface on the toric variety defined by the fine, star, regular triangulation. If a nef-partition is specified then it returns the complete intersection Calabi-Yau that it specifies.
Only Calabi-Yau 3-fold hypersurfaces are fully supported. Other dimensions and CICYs require enabling the experimental features of CYTools. See experimental features for more details.
Arguments:
nef_partition
: A list of tuples of indices specifying a nef-partition of the polytope, which correspondingly defines a complete intersection Calabi-Yau.
Returns: The Calabi-Yau arising from the triangulation.
Example
We construct a triangulation and obtain the Calabi-Yau hypersurface in the resulting toric variety.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
t = p.triangulate()
t.get_cy()
# A Calabi-Yau 3-fold hypersurface with h11=2 and h21=272 in a 4-dimensional toric variety
get_toric_variety
Description: Returns a ToricVariety object corresponding to the fan defined by the triangulation.
Arguments: None.
Returns: The toric variety arising from the triangulation.
Example
We construct a triangulation and obtain the resulting toric variety.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
t = p.triangulate()
t.get_toric_variety()
# A simplicial compact 4-dimensional toric variety with 9 affine patches
gkz_phi
Description: Returns the GKZ phi vector of the triangulation. The -th component is defined to be the sum of the volumes of the simplices that have the -th point as a vertex.
Arguments: None.
Returns: The GKZ phi vector of the triangulation.
Example
We construct a triangulation and find its GKZ vector.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
t = p.triangulate()
t.gkz_phi()
# array([18, 12, 9, 12, 12, 12, 15])
heights
Description: Returns the a height vector if the triangulation is regular. An exception is raised if a height vector could not be found either because the optimizer failed or because the triangulation is not regular.
Arguments:
integral
: Whether to find an integral height vector.backend
: The optimizer used for the computation. The available options are the backends of thefind_interior_point
function of theCone
class. If not specified, it will be picked automatically.
Returns: A height vector giving rise to the triangulation.
Example
We construct a triangulation and find a height vector that generates it.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
t = p.triangulate()
t.heights()
# array([0., 0., 0., 0., 0., 1.])
is_equivalent
Description: Compares two triangulations with or without allowing automorphism transformations and with the option of restricting to faces of a particular dimension.
Arguments:
other
: The other triangulation that is being compared.use_automorphisms
: Whether to check the equivalence using automorphism transformations. This flag is ignored for point configurations that are not full dimensional.on_faces_dim
: Restrict the simplices to faces of the polytope of a given dimension.on_faces_codim
: Restrict the simplices to faces of the polytope of a given codimension.
Returns: The truth value of the triangulations being equivalent under the specified parameters.
Example
We construct two triangulations and check whether they are equivalent under various conditions.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[-1,1,1,0],[0,-1,-1,0],[0,0,0,1],[1,-2,1,1],[-2,2,1,-1],[1,1,-1,-1]])
triangs_gen = p.all_triangulations()
t1 = next(triangs_gen)
t2 = next(triangs_gen)
t1.is_equivalent(t2)
# False
t1.is_equivalent(t2, on_faces_dim=2)
# True
t1.is_equivalent(t2, on_faces_dim=2, use_automorphisms=False)
# True
is_fine
Description: Returns True if the triangulation is fine (all the points are used), and False otherwise. Note that this only checks if it is fine with respect to the point configuration, not with respect to the full set of lattice points of the polytope.
Arguments: None.
Returns: The truth value of the triangulation being fine.
Example
We construct a triangulation and check if it is fine.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
t = p.triangulate()
t.is_fine()
# True
is_regular
Description: Returns True if the triangulation is regular and False otherwise.
Arguments:
backend
: The optimizer used for the computation. The available options are the backends of theis_solid
function of theCone
class. If not specified, it will be picked automatically.
Returns: The truth value of the triangulation being regular.
Example
We construct a triangulation and check if it is regular.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
t = p.triangulate()
t.is_regular()
# True
is_star
Description:
Returns True if the triangulation is star and False otherwise. The star
origin is assumed to be the origin, so for polytopes that don't contain
the origin this function will always return False unless star_origin
is specified.
Arguments:
star_origin
: The index of the origin of the star triangulation
Returns: The truth value of the triangulation being star.
Example
We construct a triangulation and check if it is star.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
t = p.triangulate()
t.is_star()
# True
is_valid
Description: Returns True if the presumed triangulation meets all requirements to be a triangulation. The simplices must cover the full volume of the convex hull, and they cannot intersect at full-dimensional regions.
Arguments:
backend
: The optimizer used for the computation. The available options are the backends of theis_solid
function of theCone
class. If not specified, it will be picked automatically.verbosity
: The verbosity level.
Returns: The truth value of the triangulation being valid.
Example
This function is useful when constructing a triangulation from a given set of simplices. We show this in this example.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
t = p.triangulate(simplices=[[0,1,2,3,4],[0,1,2,3,5],[0,1,2,4,5],[0,1,3,4,5],[0,2,3,4,5]]) # It is already used here unless using check_input_simplices=False
t.is_valid()
# True
labels
Description: Returns the labels of lattice points in the triangulation.
Arguments: None.
Returns: The labels of lattice points in the triangulation.
neighbor_triangulations
Description: Returns the list of triangulations that differ by one bistellar flip from the current triangulation. The computation is performed with a modified version of TOPCOM. There is the option of limiting the flips to fine, regular, and star triangulations. An additional backend is used to check regularity, as checking this with TOPCOM is very slow for large polytopes.
Arguments:
only_fine
: Restricts to fine triangulations.only_regular
: Restricts the to regular triangulations.only_star
: Restricts to star triangulations.backend
: The backend used to check regularity. The options are any backend available for theis_solid
function of theCone
class. If not specified, it will be picked automatically.verbose
: Whether to print extra info from the TOPCOM command.
Returns: The list of triangulations that differ by one bistellar flip from the current triangulation.
Example
We construct a triangulation and find its neighbor triangulations.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]]).dual()
t = p.triangulate()
triangs = t.neighbor_triangulations()
len(triangs) # Print how many triangulations it found
# 263
points
Description: Returns the points of the triangulation. Note that these are not necessarily equal to the lattice points of the polytope they define.
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 points of the triangulation.
Aliases:
pts
.
Example
We construct a triangulation and print the points in the point configuration. Note that since the polytope is reflexive, then by default only the lattice points not interior to facets were used.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
t = p.triangulate()
t.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]])
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.
poly
Description: Returns the polytope being triangulated.
Arguments: None.
Returns: The ambient polytope.
random_flips
Description: Returns a triangulation obtained by performing N random bistellar flips. The computation is performed with a modified version of TOPCOM. There is the option of limiting the flips to fine, regular, and star triangulations. An additional backend is used to check regularity, as checking this with TOPCOM is very slow for large polytopes.
Arguments:
N
: The number of bistellar flips to perform.only_fine
: Restricts to flips to fine triangulations. If not specified, it is set to True if the triangulation is fine, and False otherwise.only_regular
: Restricts the flips to regular triangulations. If not specified, it is set to True if the triangulation is regular, and False otherwise.only_star
: Restricts the flips to star triangulations. If not specified, it is set to True if the triangulation is star, and False otherwise.backend
: The backend used to check regularity. The options are any backend available for theis_solid
function of theCone
class. If not specified, it will be picked automatically.seed
: A seed for the random number generator. This can be used to obtain reproducible results.
Returns: A new triangulation obtained by performing N random flips.
Example
We construct a triangulation and perform 5 random flips to find a new triangulation.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]]).dual()
t = p.triangulate()
t.random_flips(5) # Takes a few seconds
# A fine, star triangulation of a 4-dimensional point configuration
# with 106 points in ZZ^4
restrict
Description: Restrict the triangulation to some face(s).
Arguments:
restrict_to
: The face(s) to restrict to. If none provided, gives restriction to each dim-face (see next argument).restrict_dim
: If restrict_to is None, sets the dimension of the faces to restrict to. Only used if restrict_to == None.as_poly
: Construct the formal Triangulation objects.
Returns: The restrictions
secondary_cone
Description: Computes the (hyperplanes defining the) secondary cone of the triangulation. This cone is also called the 'cone of strictly convex piecewise linear functions'.
The triangulation is regular if and only if this cone is solid (i.e. full-dimensional), in which case the points in the strict interior correspond to heights that give rise to the triangulation.
Also, allow calculation of the 'secondary cone of the N-skeleton'. This
cone is defined as
1) the intersection of the secondary cones of all N-dim faces of
this triangulation or, equivalently,
2) the union of all secondary cones arising from triangulations
with the same N-face restrictions.
Any point in the strict interior of this cone correspond to heights
which generate a triangulation with the imposed N-face restrictions.
Set N via argument on_faces_dim
.
Some cases of interest: 1) N=self.dim() -> return the secondary cone of the triangulation 2) N=2 -> return the union of all secondary cones arising from triangulations with the same 2-face restrictions. Akin to Kcup.
Arguments:
backend
: The backend to use. Options are "native", which uses a native implementation of an algorithm by Berglund, Katz and Klemm, or "topcom" which uses differences of GKZ vectors.computation.include_points_not_in_triangulation
: This flag allows the exclusion of points that are not part of the triangulation. This can be done to check regularity faster, but this cannot be used if the actual cone in the secondary fan is needed.as_cone
: Return a cone or just the defining hyperplanes.on_faces_dim
: Compute the secondary cone for each face with this dimension and then take their intersection. This has the interpretation of enforcing the triangulations on each of these faces, but being agnostic to the rest of the structure.use_cache
: Whether to use cached values of the secondary cone.
Returns: The secondary cone.
Aliases:
cpl_cone
.
Example
We construct a triangulation and find its secondary cone.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
t = p.triangulate()
t.secondary_cone()
# A rational polyhedral cone in RR^7 defined by 3 hyperplanes normals
simplices
Description: Returns the simplices of the triangulation. It also has the option of restricting the simplices to faces of the polytope of a particular dimension or codimension. This restriction is useful for checking CY equivalences from different triangulations.
Arguments:
- `on_faces_dim: Restrict the simplices to faces of the polytope of a given dimension.
on_faces_codim
: Restrict the simplices to faces of the polytope of a given codimension.split_by_face
: Return the simplices for each face. Don't merge the collection.as_np_array
: Return the simplices as a numpy array. Otherwise, they are returned as a set of frozensets.as_indices
: Whether to map the simplices from labels to point indices (in the triangulations).
Returns: The simplices of the triangulation.
Example
We construct a triangulation and find its simplices. We also find the simplices the lie on 2-faces of the polytope.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
t = p.triangulate()
t.simplices()
# array([[0, 1, 2, 3, 4],
# [0, 1, 2, 3, 5],
# [0, 1, 2, 4, 5],
# [0, 1, 3, 4, 5],
# [0, 2, 3, 4, 5]])
t.simplices(on_faces_dim=2)
# [[1 2 3]
# [1 2 4]
# [1 2 5]
# [1 3 4]
# [1 3 5]
# [1 4 5]
# [2 3 4]
# [2 3 5]
# [2 4 5]
# [3 4 5]]
sr_ideal
Description: Returns the Stanley-Reisner ideal if the triangulation is star.
That is, find all sets of points which are not subsets of any simplex. The SR-ideal is generated by such subsets.
N.B.: This function returns the generators of the ideal. I.e., we don't return multiples of the generators.
E.g., take the simplicial complex [[1,2,3],[1,2,4],[2,3,4],[1,2,5]] It has the following sets of points not appearing together (3,5), (4,5), (1,3,4),(1,3,5),(1,4,5), (2,3,5), (2,4,5), (3,4,5) So, the SR ideal is generated by x3x5, x4x5, x1x3x4 Don't include things like x2x3x5 since it is a multiple of x3x5.
Arguments: None.
Returns: The Stanley-Reisner ideal of the triangulation.
Example
We construct a triangulation and find its SR ideal.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
t = p.triangulate()
t.sr_ideal()
# array([[1, 4, 5],
# [2, 3, 6]])
triangulation_to_polytope_indices
Description: Takes a list of indices of points of the triangulation and it returns the corresponding indices of the polytope. It also accepts a single entry, in which case it returns the corresponding index.
Arguments:
points
: A list of indices of points.
Returns: The list of indices corresponding to the given points. Or the index of the point if only one is given.
Hidden Functions
__eq__
Description: Implements comparison of triangulations with ==.
Arguments:
other
: The other triangulation that is being compared.
Returns: The truth value of the triangulations being equal.
Example
We construct two triangulations and compare them.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
t1 = p.triangulate(backend="qhull")
t2 = p.triangulate(backend="topcom")
t1 == t2
# True
__hash__
Description: Implements the ability to obtain hash values from triangulations.
Arguments: None.
Returns: The hash value of the triangulation.
Example
We compute the hash value of a triangulation. Also, we construct a set and a dictionary with a triangulation, 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]])
t = p.triangulate()
h = hash(t) # Obtain hash value
d = {t: 1} # Create dictionary with triangulation keys
s = {t} # Create a set of triangulations
__init__
Description:
Initializes a Triangulation
object.
Arguments:
poly
: The ambient polytope of the points to be triangulated.pts
: The list of points to be triangulated. Specified by labels.make_star
: 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 until it is much lower than all other heights.simplices
: Array-like of simplices specifying the triangulation. Each simplex is a list of point labels. This is useful when a triangulation was previously computed and it needs to be used again. Note that the ordering of the points needs to be consistent.check_input_simplices
: Whether to check if the input simplices define a valid triangulation.heights
: The heights specifying the regular triangulation. When not specified, construct based off of the backend:Heights can only be specified when using CGAL or QHull as the backend.- (CGAL) a Delaunay triangulation,
- (QHULL) triangulation from random heights near Delaunay, or
- (TOPCOM) placing triangulation.check_heights
: Whether to check if the input/default heights define a valid/unique triangulation.backend
: The backend used to compute the triangulation. Options are "qhull", "cgal", and "topcom". CGAL is the default as it is very fast and robust.verbosity
: The verbosity level.
Returns: Nothing.
Example
This is the function that is called when creating a new
Triangulation
object. In this example the polytope is reflexive, so
by default the triangulation is fine, regular, and star. Also, since
the polytope is reflexive then by default only the lattice points not
interior to facets are included in the triangulation.
from cytools import Polytope
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
t = p.triangulate()
print(t)
# A fine, regular, star triangulation of a 4-dimensional point
# configuration with 7 points in ZZ^4
__ne__
Description: Implements comparison of triangulations with !=.
Arguments:
other
: The other triangulation that is being compared.
Returns: The truth value of the triangulations being different.
Example
We construct two triangulations and compare them.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-1,-1]])
t1 = p.triangulate(backend="qhull")
t2 = p.triangulate(backend="topcom")
t1 != t2
# False
__repr__
Description: Returns a string describing the triangulation.
Arguments: None.
Returns: A string describing the triangulation.
Example
This function can be used to convert the triangulation to a string or to print information about the triangulation.
p = Polytope([[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1],[-1,-1,-6,-9]])
t = p.triangulate()
poly_info = str(t) # Converts to string
print(t) # Prints triangulation info
# A fine, regular, star triangulation of a 4-dimensional point configuration with 7 points 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.
_fine_neighbors_2d
Description: An optimized variant of neighbor_triangulations for triangulations that are: 1) 2D (in dimension... ambient dimension doesn't matter) 2) fine 3) fine neighbors are desired In this case, _fine_neighbors_2d runs much quicker than the corresponding TOPCOM calculation
Arguments:
only_regular
: Restricts the to regular triangulations.only_star
: Restricts to star triangulations.backend
: The backend used to check regularity. The options are any backend available for theis_solid
function of theCone
class. If not specified, it will be picked automatically.
Returns: The list of triangulations that differ by one diagonal flip from the current triangulation.