Path: blob/master/src/sage/homology/simplicial_complex.py
8818 views
r"""1Finite simplicial complexes23AUTHORS:45- John H. Palmieri (2009-04)67- D. Benjamin Antieau (2009-06) - added is_connected, generated_subcomplex,8remove_facet, and is_flag_complex methods;9cached the output of the graph() method.1011- Travis Scrimshaw (2012-08-17): Made :class:`SimplicialComplex` have an12immutable option, and added ``__hash__()`` function which checks to make13sure it is immutable. Made :meth:`SimplicialComplex.remove_face()` into a14mutator. Deprecated the ``vertex_set`` parameter.1516- Christian Stump (2011-06) - implementation of is_cohen_macaulay1718- Travis Scrimshaw (2013-02-16): Allowed :class:`SimplicialComplex` to make19mutable copies.2021This module implements the basic structure of finite simplicial22complexes. Given a set `V` of "vertices", a simplicial complex on `V`23is a collection `K` of subsets of `V` satisfying the condition that if24`S` is one of the subsets in `K`, then so is every subset of `S`. The25subsets `S` are called the 'simplices' of `K`.2627A simplicial complex `K` can be viewed as a purely combinatorial28object, as described above, but it also gives rise to a topological29space `|K|` (its *geometric realization*) as follows: first, the30points of `V` should be in general position in euclidean space. Next,31if `\{v\}` is in `K`, then the vertex `v` is in `|K|`. If `\{v, w\}`32is in `K`, then the line segment from `v` to `w` is in `|K|`. If `\{u,33v, w\}` is in `K`, then the triangle with vertices `u`, `v`, and `w`34is in `|K|`. In general, `|K|` is the union of the convex hulls of35simplices of `K`. Frequently, one abuses notation and uses `K` to36denote both the simplicial complex and the associated topological37space.3839.. image:: ../../media/simplices.png4041For any simplicial complex `K` and any commutative ring `R` there is42an associated chain complex, with differential of degree `-1`. The43`n^{th}` term is the free `R`-module with basis given by the44`n`-simplices of `K`. The differential is determined by its value on45any simplex: on the `n`-simplex with vertices `(v_0, v_1, ..., v_n)`,46the differential is the alternating sum with `i^{th}` summand `(-1)^i`47multiplied by the `(n-1)`-simplex obtained by omitting vertex `v_i`.4849In the implementation here, the vertex set must be finite. To define a50simplicial complex, specify its vertex set: this should be a list,51tuple, or set, or it can be a non-negative integer `n`, in which case52the vertex set is `(0, ..., n)`. Also specify the facets: the maximal53faces.5455.. NOTE::5657The elements of the vertex set are not automatically contained in58the simplicial complex: each one is only included if and only if it59is a vertex of at least one of the specified facets.6061.. NOTE::6263This class derives from64:class:`~sage.homology.cell_complex.GenericCellComplex`, and so65inherits its methods. Some of those methods are not listed here;66see the :mod:`Generic Cell Complex <sage.homology.cell_complex>`67page instead.6869EXAMPLES::7071sage: SimplicialComplex([[1], [3, 7]])72Simplicial complex with vertex set (1, 3, 7) and facets {(3, 7), (1,)}73sage: SimplicialComplex() # the empty simplicial complex74Simplicial complex with vertex set () and facets {()}75sage: X = SimplicialComplex([[0,1], [1,2], [2,3], [3,0]])76sage: X77Simplicial complex with vertex set (0, 1, 2, 3) and facets {(1, 2), (2, 3), (0, 3), (0, 1)}78sage: X.stanley_reisner_ring()79Quotient of Multivariate Polynomial Ring in x0, x1, x2, x3 over Integer Ring by the ideal (x1*x3, x0*x2)80sage: X.is_pure()81True8283Sage can perform a number of operations on simplicial complexes, such84as the join and the product, and it can also compute homology::8586sage: S = SimplicialComplex([[0,1], [1,2], [0,2]]) # circle87sage: T = S.product(S) # torus88sage: T89Simplicial complex with 9 vertices and 18 facets90sage: T.homology() # this computes reduced homology91{0: 0, 1: Z x Z, 2: Z}92sage: T.euler_characteristic()9309495Sage knows about some basic combinatorial data associated to a96simplicial complex::9798sage: X = SimplicialComplex([[0,1], [1,2], [2,3], [0,3]])99sage: X.f_vector()100[1, 4, 4]101sage: X.face_poset()102Finite poset containing 8 elements103sage: X.stanley_reisner_ring()104Quotient of Multivariate Polynomial Ring in x0, x1, x2, x3 over Integer Ring by the ideal (x1*x3, x0*x2)105106Mutability (see :trac:`12587`)::107108sage: S = SimplicialComplex([[1,4], [2,4]])109sage: S.add_face([1,3])110sage: S.remove_face([1,3]); S111Simplicial complex with vertex set (1, 2, 3, 4) and facets {(2, 4), (1, 4), (3,)}112sage: hash(S)113Traceback (most recent call last):114...115ValueError: This simplicial complex must be immutable. Call set_immutable().116sage: S = SimplicialComplex([[1,4], [2,4]])117sage: S.set_immutable()118sage: S.add_face([1,3])119Traceback (most recent call last):120...121ValueError: This simplicial complex is not mutable122sage: S.remove_face([1,3])123Traceback (most recent call last):124...125ValueError: This simplicial complex is not mutable126sage: hash(S) == hash(S)127True128129sage: S2 = SimplicialComplex([[1,4], [2,4]], is_mutable=False)130sage: hash(S2) == hash(S)131True132133We can also make mutable copies of an immutable simplicial complex134(see :trac:`14142`)::135136sage: S = SimplicialComplex([[1,4], [2,4]])137sage: S.set_immutable()138sage: T = copy(S)139sage: T.is_mutable()140True141sage: S == T142True143"""144145# possible future directions for SimplicialComplex:146#147# make compatible with GAP (see http://linalg.org/gap.html)148# compare to and make compatible with polymake149# see Macaulay: http://www.math.uiuc.edu/Macaulay2/doc/Macaulay2-1.1/share/doc/Macaulay2/SimplicialComplexes/html/___Simplicial__Complex.html; compare performance and make compatible150# should + have any meaning?151# cohomology: compute cup products (and Massey products?)152153from copy import copy154from sage.homology.cell_complex import GenericCellComplex155from sage.structure.sage_object import SageObject156from sage.rings.integer import Integer157from sage.rings.polynomial.polynomial_ring_constructor import PolynomialRing158from sage.sets.set import Set159from sage.rings.integer_ring import ZZ160from sage.structure.parent_gens import normalize_names161from sage.misc.latex import latex162from sage.matrix.constructor import matrix163from sage.homology.chain_complex import ChainComplex164from sage.graphs.graph import Graph165166def lattice_paths(t1, t2, length=None):167"""168Given lists (or tuples or ...) ``t1`` and ``t2``, think of them as169labelings for vertices: ``t1`` labeling points on the x-axis,170``t2`` labeling points on the y-axis, both increasing. Return the171list of rectilinear paths along the grid defined by these points172in the plane, starting from ``(t1[0], t2[0])``, ending at173``(t1[last], t2[last])``, and at each grid point, going either174right or up. See the examples.175176:param t1: labeling for vertices177:param t2: labeling for vertices178:param length: if not ``None``, then an integer, the length of the desired179path.180:type length: integer or ``None``; optional, default ``None``181:type t1: tuple, list, other iterable182:type t2: tuple, list, other iterable183:return: list of lists of vertices making up the paths as described above184:rtype: list of lists185186This is used when triangulating the product of simplices. The187optional argument ``length`` is used for `\Delta`-complexes, to188specify all simplices in a product: in the triangulation of a189product of two simplices, there is a `d`-simplex for every path of190length `d+1` in the lattice. The path must start at the bottom191left and end at the upper right, and it must use at least one192point in each row and in each column, so if ``length`` is too193small, there will be no paths.194195EXAMPLES::196197sage: from sage.homology.simplicial_complex import lattice_paths198sage: lattice_paths([0,1,2], [0,1,2])199[[(0, 0), (0, 1), (0, 2), (1, 2), (2, 2)],200[(0, 0), (0, 1), (1, 1), (1, 2), (2, 2)],201[(0, 0), (1, 0), (1, 1), (1, 2), (2, 2)],202[(0, 0), (0, 1), (1, 1), (2, 1), (2, 2)],203[(0, 0), (1, 0), (1, 1), (2, 1), (2, 2)],204[(0, 0), (1, 0), (2, 0), (2, 1), (2, 2)]]205sage: lattice_paths(('a', 'b', 'c'), (0, 3, 5))206[[('a', 0), ('a', 3), ('a', 5), ('b', 5), ('c', 5)],207[('a', 0), ('a', 3), ('b', 3), ('b', 5), ('c', 5)],208[('a', 0), ('b', 0), ('b', 3), ('b', 5), ('c', 5)],209[('a', 0), ('a', 3), ('b', 3), ('c', 3), ('c', 5)],210[('a', 0), ('b', 0), ('b', 3), ('c', 3), ('c', 5)],211[('a', 0), ('b', 0), ('c', 0), ('c', 3), ('c', 5)]]212sage: lattice_paths(range(3), range(3), length=2)213[]214sage: lattice_paths(range(3), range(3), length=3)215[[(0, 0), (1, 1), (2, 2)]]216sage: lattice_paths(range(3), range(3), length=4)217[[(0, 0), (1, 1), (1, 2), (2, 2)],218[(0, 0), (0, 1), (1, 2), (2, 2)],219[(0, 0), (1, 1), (2, 1), (2, 2)],220[(0, 0), (1, 0), (2, 1), (2, 2)],221[(0, 0), (0, 1), (1, 1), (2, 2)],222[(0, 0), (1, 0), (1, 1), (2, 2)]]223"""224if length is None:225# 0 x n (or k x 0) rectangle:226if len(t1) == 0 or len(t2) == 0:227return [[]]228# 1 x n (or k x 1) rectangle:229elif len(t1) == 1:230return [[(t1[0], w) for w in t2]]231elif len(t2) == 1:232return [[(v, t2[0]) for v in t1]]233else:234# recursive: paths in rectangle with either one fewer row235# or column, plus the upper right corner236return ([path + [(t1[-1], t2[-1])] for path237in lattice_paths(t1[:-1], t2)] +238[path + [(t1[-1], t2[-1])] for path239in lattice_paths(t1, t2[:-1])])240else:241if length > len(t1) + len(t2) - 1:242return []243# as above, except make sure that lengths are correct. if244# not, return an empty list.245#246# 0 x n (or k x 0) rectangle:247elif len(t1) == 0 or len(t2) == 0:248if length == 0:249return [[]]250else:251return []252# 1 x n (or k x 1) rectangle:253elif len(t1) == 1:254if length == len(t2):255return [[(t1[0], w) for w in t2]]256else:257return []258elif len(t2) == 1:259if length == len(t1):260return [[(v, t2[0]) for v in t1]]261else:262return []263else:264# recursive: paths of length one fewer in rectangle with265# either one fewer row, one fewer column, or one fewer of266# each, and then plus the upper right corner267return ([path + [(t1[-1], t2[-1])] for path268in lattice_paths(t1[:-1], t2, length=length-1)] +269[path + [(t1[-1], t2[-1])] for path270in lattice_paths(t1, t2[:-1], length=length-1)] +271[path + [(t1[-1], t2[-1])] for path272in lattice_paths(t1[:-1], t2[:-1], length=length-1)])273274def rename_vertex(n, keep, left = True):275"""276Rename a vertex: the vertices from the list ``keep`` get277relabeled 0, 1, 2, ..., in order. Any other vertex (e.g. 4) gets278renamed to by prepending an 'L' or an 'R' (thus to either 'L4' or279'R4'), depending on whether the argument left is ``True`` or ``False``.280281:param n: a 'vertex': either an integer or a string282:param keep: a list of three vertices283:param left: if ``True``, rename for use in left factor284:type left: boolean; optional, default ``True``285286This is used by the :meth:`~SimplicialComplex.connected_sum` method for287simplicial complexes.288289EXAMPLES::290291sage: from sage.homology.simplicial_complex import rename_vertex292sage: rename_vertex(6, [5, 6, 7])2931294sage: rename_vertex(3, [5, 6, 7, 8, 9])295'L3'296sage: rename_vertex(3, [5, 6, 7], left=False)297'R3'298"""299lookup = dict(zip(keep, range(len(keep))))300try:301return lookup[n]302except KeyError:303if left:304return "L" + str(n)305else:306return "R" + str(n)307308class Simplex(SageObject):309"""310Define a simplex.311312Topologically, a simplex is the convex hull of a collection of313vertices in general position. Combinatorially, it is defined just314by specifying a set of vertices. It is represented in Sage by the315tuple of the vertices.316317:param X: set of vertices318:type X: integer or list, tuple, or other iterable319:return: simplex with those vertices320321``X`` may be a non-negative integer `n`, in which case the322simplicial complex will have `n+1` vertices `(0, 1, ..., n)`, or323it may be anything which may be converted to a tuple, in which324case the vertices will be that tuple. In the second case, each325vertex must be hashable, so it should be a number, a string, or a326tuple, for instance, but not a list.327328.. WARNING::329330The vertices should be distinct, and no error checking is done331to make sure this is the case.332333EXAMPLES::334335sage: Simplex(4)336(0, 1, 2, 3, 4)337sage: Simplex([3, 4, 1])338(3, 4, 1)339sage: X = Simplex((3, 'a', 'vertex')); X340(3, 'a', 'vertex')341sage: X == loads(dumps(X))342True343344Vertices may be tuples but not lists::345346sage: Simplex([(1,2), (3,4)])347((1, 2), (3, 4))348sage: Simplex([[1,2], [3,4]])349Traceback (most recent call last):350...351TypeError: unhashable type: 'list'352"""353354def __init__(self, X):355"""356Define a simplex. See :class:`Simplex` for full documentation.357358EXAMPLES::359360sage: Simplex(2)361(0, 1, 2)362sage: Simplex(('a', 'b', 'c'))363('a', 'b', 'c')364"""365try:366N = int(X) + 1367self.__tuple = tuple(range(N))368except TypeError:369self.__tuple = tuple(X)370self.__set = frozenset(self.__tuple)371372def tuple(self):373"""374The tuple attached to this simplex.375376EXAMPLES::377378sage: Simplex(3).tuple()379(0, 1, 2, 3)380381Although simplices are printed as if they were tuples, they382are not the same type::383384sage: type(Simplex(3).tuple())385<type 'tuple'>386sage: type(Simplex(3))387<class 'sage.homology.simplicial_complex.Simplex'>388"""389return self.__tuple390391def set(self):392"""393The frozenset attached to this simplex.394395EXAMPLES::396397sage: Simplex(3).set()398frozenset([0, 1, 2, 3])399"""400return self.__set401402def is_face(self, other):403"""404Return ``True`` iff this simplex is a face of other.405406EXAMPLES::407408sage: Simplex(3).is_face(Simplex(5))409True410sage: Simplex(5).is_face(Simplex(2))411False412sage: Simplex(['a', 'b', 'c']).is_face(Simplex(8))413False414"""415return self.__set.issubset(other.__set)416417def __contains__(self, x):418"""419Return ``True`` iff ``x`` is a vertex of this simplex.420421EXAMPLES::422423sage: 3 in Simplex(5)424True425sage: 3 in Simplex(2)426False427"""428return self.__set.__contains__(x)429430def __getitem__(self, n):431"""432Return the `n`-th vertex in this simplex.433434EXAMPLES::435436sage: Simplex(5)[2]4372438sage: Simplex(['a', 'b', 'c'])[1]439'b'440"""441return self.__tuple.__getitem__(n)442443def __iter__(self):444"""445Iterator for the vertices of this simplex.446447EXAMPLES::448449sage: [v**2 for v in Simplex(3)]450[0, 1, 4, 9]451"""452return self.__tuple.__iter__()453454def __add__(self, other):455"""456Simplex obtained by concatenating the underlying tuples of the457two arguments.458459:param other: another simplex460461EXAMPLES::462463sage: Simplex((1,2,3)) + Simplex((5,6))464(1, 2, 3, 5, 6)465"""466return Simplex(self.__tuple.__add__(other.__tuple))467468def face(self, n):469"""470The `n`-th face of this simplex.471472:param n: an integer between 0 and the dimension of this simplex473:type n: integer474:return: the simplex obtained by removing the `n`-th vertex from this475simplex476477EXAMPLES::478479sage: S = Simplex(4)480sage: S.face(0)481(1, 2, 3, 4)482sage: S.face(3)483(0, 1, 2, 4)484"""485if n >= 0 and n <= self.dimension():486return Simplex(self.__tuple[:n] + self.__tuple[n+1:])487else:488raise IndexError, "%s does not have an nth face for n=%s." % (self, n)489490def faces(self):491"""492The list of faces (of codimension 1) of this simplex.493494EXAMPLES::495496sage: S = Simplex(4)497sage: S.faces()498[(1, 2, 3, 4), (0, 2, 3, 4), (0, 1, 3, 4), (0, 1, 2, 4), (0, 1, 2, 3)]499sage: len(Simplex(10).faces())50011501"""502return [self.face(i) for i in range(self.dimension()+1)]503504def dimension(self):505"""506The dimension of this simplex.507508The dimension of a simplex is the number of vertices minus 1.509510EXAMPLES::511512sage: Simplex(5).dimension() == 5513True514sage: Simplex(5).face(1).dimension()5154516"""517return len(self.__tuple) - 1518519def is_empty(self):520"""521Return ``True`` iff this simplex is the empty simplex.522523EXAMPLES::524525sage: [Simplex(n).is_empty() for n in range(-1,4)]526[True, False, False, False, False]527"""528return self.dimension() < 0529530def join(self, right, rename_vertices=True):531"""532The join of this simplex with another one.533534The join of two simplices `[v_0, ..., v_k]` and `[w_0, ...,535w_n]` is the simplex `[v_0, ..., v_k, w_0, ..., w_n]`.536537:param right: the other simplex (the right-hand factor)538539:param rename_vertices: If this is ``True``, the vertices in the540join will be renamed by this formula: vertex "v" in the541left-hand factor --> vertex "Lv" in the join, vertex "w"542in the right-hand factor --> vertex "Rw" in the join. If543this is false, this tries to construct the join without544renaming the vertices; this may cause problems if the two545factors have any vertices with names in common.546547:type rename_vertices: boolean; optional, default ``True``548549EXAMPLES::550551sage: Simplex(2).join(Simplex(3))552('L0', 'L1', 'L2', 'R0', 'R1', 'R2', 'R3')553sage: Simplex(['a', 'b']).join(Simplex(['x', 'y', 'z']))554('La', 'Lb', 'Rx', 'Ry', 'Rz')555sage: Simplex(['a', 'b']).join(Simplex(['x', 'y', 'z']), rename_vertices=False)556('a', 'b', 'x', 'y', 'z')557"""558if rename_vertices:559vertex_set = (["L" + str(v) for v in self]560+ ["R" + str(w) for w in right])561else:562vertex_set = self.__tuple + right.__tuple563return Simplex(vertex_set)564565def product(self, other, rename_vertices=True):566r"""567The product of this simplex with another one, as a list of simplices.568569:param other: the other simplex570571:param rename_vertices: If this is ``False``, then the vertices in572the product are the set of ordered pairs `(v,w)` where `v`573is a vertex in the left-hand factor (``self``) and `w` is574a vertex in the right-hand factor (``other``). If this is575``True``, then the vertices are renamed as "LvRw" (e.g., the576vertex (1,2) would become "L1R2"). This is useful if you577want to define the Stanley-Reisner ring of the complex:578vertex names like (0,1) are not suitable for that, while579vertex names like "L0R1" are.580581:type rename_vertices: boolean; optional, default ``True``582583Algorithm: see Hatcher, p. 277-278 [Hat]_ (who in turn refers to584Eilenberg-Steenrod, p. 68): given ``S = Simplex(m)`` and585``T = Simplex(n)``, then `S \times T` can be586triangulated as follows: for each path `f` from `(0,0)` to587`(m,n)` along the integer grid in the plane, going up or right588at each lattice point, associate an `(m+n)`-simplex with589vertices `v_0`, `v_1`, ..., where `v_k` is the `k^{th}` vertex590in the path `f`.591592Note that there are `m+n` choose `n` such paths. Note also593that each vertex in the product is a pair of vertices `(v,w)`594where `v` is a vertex in the left-hand factor and `w`595is a vertex in the right-hand factor.596597.. NOTE::598599This produces a list of simplices -- not a :class:`Simplex`, not600a :class:`SimplicialComplex`.601602EXAMPLES::603604sage: len(Simplex(2).product(Simplex(2)))6056606sage: Simplex(1).product(Simplex(1))607[('L0R0', 'L0R1', 'L1R1'), ('L0R0', 'L1R0', 'L1R1')]608sage: Simplex(1).product(Simplex(1), rename_vertices=False)609[((0, 0), (0, 1), (1, 1)), ((0, 0), (1, 0), (1, 1))]610"""611if not rename_vertices:612return [Simplex(x) for x in lattice_paths(self.tuple(), other.tuple())]613614answer = []615for x in lattice_paths(self.tuple(), other.tuple()):616new = tuple(["L" + str(v) + "R" + str(w) for (v,w) in x])617answer.append(Simplex(new))618return answer619620def __cmp__(self, other):621"""622Return ``True`` iff this simplex is the same as ``other``: that623is, if the vertices of the two are the same, even with a624different ordering625626:param other: the other simplex627628EXAMPLES::629630sage: Simplex([0,1,2]) == Simplex([0,2,1])631True632sage: Simplex([0,1,2]) == Simplex(['a','b','c'])633False634sage: Simplex([1]) < Simplex([2])635True636sage: Simplex([1]) > Simplex([2])637False638"""639if not isinstance(other, Simplex):640return -1641if self.__set == other.__set:642return 0643return cmp(sorted(tuple(self.__set)), sorted(tuple(other.__set)))644645def __hash__(self):646"""647Hash value for this simplex. This computes the hash value of648the Python frozenset of the underlying tuple, since this is649what's important when testing equality.650651EXAMPLES::652653sage: Simplex([1,2,0]).__hash__() == Simplex(2).__hash__()654True655sage: Simplex([1,2,0,1,1,2]).__hash__() == Simplex(2).__hash__()656True657"""658return hash(self.__set)659660def _repr_(self):661"""662Print representation.663664EXAMPLES::665666sage: S = Simplex(5)667sage: S._repr_()668'(0, 1, 2, 3, 4, 5)'669"""670return self.__tuple.__repr__()671672def _latex_(self):673r"""674LaTeX representation.675676EXAMPLES::677678sage: Simplex(3)._latex_()679\left(0,6801,6812,6823\right)683"""684return latex(self.__tuple)685686class SimplicialComplex(GenericCellComplex):687r"""688Define a simplicial complex.689690:param maximal_faces: set of maximal faces691:param maximality_check: see below692:type maximality_check: boolean; optional, default ``True``693:param sort_facets: see below694:type sort_facets: boolean; optional, default ``True``695:param name_check: see below696:type name_check: boolean; optional, default ``False``697:param is_mutable: Set to ``False`` to make this immutable698:type is_mutable: boolean; optional, default ``True``699:return: a simplicial complex700701``maximal_faces`` should be a list or tuple or set (indeed,702anything which may be converted to a set) whose elements are lists703(or tuples, etc.) of vertices. Maximal faces are also known as704'facets'.705706If ``maximality_check`` is ``True``, check that each maximal face is,707in fact, maximal. In this case, when producing the internal708representation of the simplicial complex, omit those that are not.709It is highly recommended that this be ``True``; various methods for710this class may fail if faces which are claimed to be maximal are711in fact not.712713If ``sort_facets`` is ``True``, sort the vertices in each facet. If714the vertices in different facets are not ordered compatibly (e.g.,715if you have facets ``(1, 3, 5)`` and ``(5, 3, 8)``), then homology716calculations may have unpredictable results.717718If ``name_check`` is ``True``, check the names of the vertices to see719if they can be easily converted to generators of a polynomial ring720-- use this if you plan to use the Stanley-Reisner ring for the721simplicial complex.722723.. WARNING::724725Earlier versions of Sage supported a ``vertex_set`` argument726to specify the vertices. This is now deprecated -- see727:trac:`12587` -- the set of vertices is determined from the728maximal faces.729730EXAMPLES::731732sage: SimplicialComplex([[1,2], [1,4]])733Simplicial complex with vertex set (1, 2, 4) and facets {(1, 2), (1, 4)}734sage: SimplicialComplex([[0,2], [0,3], [0]])735Simplicial complex with vertex set (0, 2, 3) and facets {(0, 2), (0, 3)}736sage: SimplicialComplex([[0,2], [0,3], [0]], maximality_check=False)737Simplicial complex with vertex set (0, 2, 3) and facets {(0, 2), (0, 3), (0,)}738sage: S = SimplicialComplex((('a', 'b'), ['a', 'c'], ('b', 'c')))739sage: S740Simplicial complex with vertex set ('a', 'b', 'c') and facets {('b', 'c'), ('a', 'c'), ('a', 'b')}741742Finally, if there is only one argument and it is a743simplicial complex, return that complex. If it is an object with744a built-in conversion to simplicial complexes (via a745``_simplicial_`` method), then the resulting simplicial complex is746returned::747748sage: S = SimplicialComplex([[0,2], [0,3], [0,6]])749sage: SimplicialComplex(S) == S750True751sage: Tc = cubical_complexes.Torus(); Tc752Cubical complex with 16 vertices and 64 cubes753sage: Ts = SimplicialComplex(Tc); Ts754Simplicial complex with 16 vertices and 32 facets755sage: Ts.homology()756{0: 0, 1: Z x Z, 2: Z}757758TESTS:759760Check that we can make mutable copies (see :trac:`14142`)::761762sage: S = SimplicialComplex([[0,2], [0,3]], is_mutable=False)763sage: S.is_mutable()764False765sage: C = copy(S)766sage: C.is_mutable()767True768sage: SimplicialComplex(S, is_mutable=True).is_mutable()769True770sage: SimplicialComplex(S, is_immutable=False).is_mutable()771True772"""773774def __init__(self, vertex_set=None, maximal_faces=None, **kwds):775"""776Define a simplicial complex. See ``SimplicialComplex`` for more777documentation.778779.. WARNING::780781We are deprecating the option ``vertex_set`` in :trac:`12587`.782783EXAMPLES::784785sage: SimplicialComplex([[0,2], [0,3], [0]])786Simplicial complex with vertex set (0, 2, 3) and facets {(0, 2), (0, 3)}787sage: SimplicialComplex((('a', 'b'), ('a', 'c'), ('b', 'c')))788Simplicial complex with vertex set ('a', 'b', 'c') and facets {('b', 'c'), ('a', 'c'), ('a', 'b')}789790TESTS::791792sage: S = SimplicialComplex([[1,4], [2,4]])793sage: S2 = SimplicialComplex([[1,4], [2,4]], is_mutable=False)794sage: S == S2795True796sage: S3 = SimplicialComplex(maximal_faces=[[1,4], [2,4]])797sage: S == S3798True799800sage: S = SimplicialComplex((('a', 'b'), ('a', 'c'), ('b', 'c')))801sage: S == loads(dumps(S))802True803804sage: Y = SimplicialComplex([1,2,3,4], [[1,2], [2,3], [3,4]])805doctest:1: DeprecationWarning: vertex_set is deprecated.806See http://trac.sagemath.org/12587 for details.807sage: Y = SimplicialComplex([1,2,3,4], [[1,2], [2,3], [3,4]], vertex_check=False)808doctest:1: DeprecationWarning: vertex_check is deprecated.809See http://trac.sagemath.org/12587 for details.810"""811from sage.misc.misc import union812# process kwds813sort_facets = kwds.get('sort_facets', True)814maximality_check = kwds.get('maximality_check', True)815name_check = kwds.get('name_check', False)816# done with kwds except mutability817818# For deprecation #12587819if maximal_faces is None:820maximal_faces = vertex_set821elif vertex_set is not None:822# We've passed in both vertex_set and maximal_faces823from sage.misc.superseded import deprecation824deprecation(12587, "vertex_set is deprecated.")825826if 'vertex_check' in kwds:827from sage.misc.superseded import deprecation828deprecation(12587, "vertex_check is deprecated.")829830C = None831if maximal_faces is None:832vertex_set = []833maximal_faces = []834elif isinstance(maximal_faces, SimplicialComplex):835C = maximal_faces836else:837try:838C = maximal_faces._simplicial_()839except AttributeError:840if not isinstance(maximal_faces, (list, tuple, Simplex)):841# Convert it into a list (in case it is an iterable)842maximal_faces = list(maximal_faces)843if len(maximal_faces) != 0:844vertex_set = reduce(union, maximal_faces)845if C is not None:846self._vertex_set = copy(C.vertices())847self._facets = list(C.facets())848self._faces = copy(C._faces)849self._gen_dict = copy(C._gen_dict)850self._complex = copy(C._complex)851self.__contractible = copy(C.__contractible)852self.__enlarged = copy(C.__enlarged)853self._graph = copy(C._graph)854self._is_mutable = True855return856857if sort_facets:858try: # vertex_set is an iterable859vertices = Simplex(sorted(vertex_set))860except TypeError: # vertex_set is an integer861vertices = Simplex(vertex_set)862else:863vertices = Simplex(vertex_set)864gen_dict = {}865for v in vertices:866if name_check:867try:868if int(v) < 0:869raise ValueError("The vertex %s does not have an appropriate name."%v)870except ValueError: # v is not an integer871try:872normalize_names(1, v)873except ValueError:874raise ValueError("The vertex %s does not have an appropriate name."%v)875# build dictionary of generator names876try:877gen_dict[v] = 'x%s'%int(v)878except Exception:879gen_dict[v] = v880# build set of facets881good_faces = []882maximal_simplices = [Simplex(f) for f in maximal_faces]883for face in maximal_simplices:884# check whether each given face is actually maximal885face_is_maximal = True886if maximality_check:887faces_to_be_removed = []888for other in good_faces:889if other.is_face(face):890faces_to_be_removed.append(other)891elif face_is_maximal:892face_is_maximal = not face.is_face(other)893for x in faces_to_be_removed:894good_faces.remove(x)895if sort_facets:896face = Simplex(sorted(face.tuple()))897if face_is_maximal:898good_faces += [face]899# if no maximal faces, add the empty face as a facet900if len(maximal_simplices) == 0:901good_faces.append(Simplex(-1))902# now record the attributes for self903# self._vertex_set: the Simplex formed by the vertices904self._vertex_set = vertices905# self._facets: list of facets906self._facets = good_faces907# self._sorted: True if the vertex set should be sorted. This908# gets used by the add_faces method.909self._sorted = sort_facets910# self._faces: dictionary of dictionaries of faces. The main911# dictionary is keyed by subcomplexes, and each value is a912# dictionary keyed by dimension. This should be empty until913# needed -- that is, until the faces method is called914self._faces = {}915# self._gen_dict: dictionary of names for the polynomial916# generators of the Stanley-Reisner ring917self._gen_dict = gen_dict918# self._complex: dictionary indexed by dimension d, subcomplex,919# etc.: differential from dim d to dim d-1 in the associated920# chain complex. thus to get the differential in the cochain921# complex from dim d-1 to dim d, take the transpose of this922# one.923self._complex = {}924# self.__contractible: if not None, a contractible subcomplex925# of self, as found by the _contractible_subcomplex method.926self.__contractible = None927# self.__enlarged: dictionary of enlarged subcomplexes,928# indexed by subcomplexes. For use in the _enlarge_subcomplex929# method.930self.__enlarged = {}931# initialize self._graph to None.932self._graph = None933934# Handle mutability keywords935self._is_mutable = True936if not kwds.get('is_mutable', True) or kwds.get('is_immutable', False):937self.set_immutable()938939def __hash__(self):940"""941Compute the hash value of ``self``.942943If this simplicial complex is immutable, it computes the hash value944based upon the facets. Otherwise it raises a ``ValueError``.945946EXAMPLES::947948sage: S = SimplicialComplex([[1,4], [2,4]])949sage: hash(S)950Traceback (most recent call last):951...952ValueError: This simplicial complex must be immutable. Call set_immutable().953sage: S.set_immutable()954sage: hash(S) == hash(S)955True956sage: S2 = SimplicialComplex([[1,4], [2,4]], is_mutable=False)957sage: S == S2958True959sage: hash(S) == hash(S2)960True961"""962if self._is_mutable:963raise ValueError("This simplicial complex must be immutable. Call set_immutable().")964return hash(self._facets)965966def __cmp__(self,right):967"""968Two simplicial complexes are equal iff their vertex sets are969equal and their sets of facets are equal.970971EXAMPLES::972973sage: SimplicialComplex([[1,2], [2,3], [4]]) == SimplicialComplex([[4], [2,3], [3], [2,1]])974True975sage: X = SimplicialComplex()976sage: X.add_face([1,3])977sage: X == SimplicialComplex([[1,3]])978True979"""980if set(self._facets) == set(right._facets):981return 0982else:983return -1984985def __copy__(self):986"""987Return a mutable copy of ``self``.988989EXAMPLES::990991sage: S = SimplicialComplex([[0,2], [0,3]], is_mutable=False)992sage: S.is_mutable()993False994sage: C = copy(S)995sage: C.is_mutable()996True997sage: C == S998True999sage: S.is_mutable()1000False1001sage: T = copy(C)1002sage: T == C1003True1004"""1005return SimplicialComplex(self, is_mutable=True)10061007def vertices(self):1008"""1009The vertex set of this simplicial complex.10101011EXAMPLES::10121013sage: S = SimplicialComplex([[i] for i in range(16)] + [[0,1], [1,2]])1014sage: S1015Simplicial complex with 16 vertices and 15 facets1016sage: S.vertices()1017(0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15)10181019Note that this actually returns a simplex::10201021sage: type(S.vertices())1022<class 'sage.homology.simplicial_complex.Simplex'>1023"""1024return self._vertex_set10251026def maximal_faces(self):1027"""1028The maximal faces (a.k.a. facets) of this simplicial complex.10291030This just returns the set of facets used in defining the1031simplicial complex, so if the simplicial complex was defined1032with no maximality checking, none is done here, either.10331034EXAMPLES::10351036sage: Y = SimplicialComplex([[0,2], [1,4]])1037sage: Y.maximal_faces()1038{(1, 4), (0, 2)}10391040``facets`` is a synonym for ``maximal_faces``::10411042sage: S = SimplicialComplex([[0,1], [0,1,2]])1043sage: S.facets()1044{(0, 1, 2)}1045"""1046return Set(self._facets)10471048facets = maximal_faces10491050def faces(self, subcomplex=None):1051"""1052The faces of this simplicial complex, in the form of a1053dictionary of sets keyed by dimension. If the optional1054argument ``subcomplex`` is present, then return only the1055faces which are *not* in the subcomplex.10561057:param subcomplex: a subcomplex of this simplicial complex.1058Return faces which are not in this subcomplex.10591060:type subcomplex: optional, default ``None``10611062EXAMPLES::10631064sage: Y = SimplicialComplex([[1,2], [1,4]])1065sage: Y.faces()1066{0: set([(4,), (2,), (1,)]), 1: set([(1, 2), (1, 4)]), -1: set([()])}1067sage: L = SimplicialComplex([[1,2]])1068sage: Y.faces(subcomplex=L)1069{0: set([(4,)]), 1: set([(1, 4)]), -1: set([])}1070"""1071# Make the subcomplex immutable if it is not1072if subcomplex is not None and subcomplex._is_mutable:1073subcomplex = SimplicialComplex(subcomplex._facets, maximality_check=False,1074sort_facets=False, is_mutable=False)10751076if subcomplex not in self._faces:1077# Faces is the dictionary of faces in self but not in1078# subcomplex, indexed by dimension1079Faces = {}1080# sub_facets is the dictionary of facets in the subcomplex1081sub_facets = {}1082dimension = max([face.dimension() for face in self._facets])1083for i in range(-1,dimension+1):1084Faces[i] = set([])1085sub_facets[i] = set([])1086for f in self._facets:1087dim = f.dimension()1088Faces[dim].add(f)1089if subcomplex is not None:1090for g in subcomplex._facets:1091dim = g.dimension()1092Faces[dim].discard(g)1093sub_facets[dim].add(g)1094# bad_faces is the set of faces in the subcomplex in the1095# current dimension1096bad_faces = sub_facets[dimension]1097for dim in range(dimension, -1, -1):1098# bad_bdries = boundaries of bad_faces: things to be1099# discarded in dim-11100bad_bdries = sub_facets[dim-1]1101for f in bad_faces:1102bad_bdries.update(f.faces())1103for f in Faces[dim]:1104Faces[dim-1].update(set(f.faces()).difference(bad_bdries))1105bad_faces = bad_bdries1106self._faces[subcomplex] = Faces1107return self._faces[subcomplex]11081109def face_iterator(self, increasing=True):1110"""1111An iterator for the faces in this simplicial complex.11121113INPUTS:11141115- ``increasing`` -- (optional, default ``True``) if ``True``, return1116faces in increasing order of dimension, thus starting with1117the empty face. Otherwise it returns faces in decreasing order of1118dimension.11191120EXAMPLES::11211122sage: S1 = simplicial_complexes.Sphere(1)1123sage: [f for f in S1.face_iterator()]1124[(), (2,), (0,), (1,), (1, 2), (0, 2), (0, 1)]1125"""1126Fs = self.faces()1127dim_index = xrange(-1,self.dimension()+1)1128if not increasing:1129dim_index = reversed(dim_index)1130for i in dim_index:1131for F in Fs[i]:1132yield F11331134cells = faces11351136def n_faces(self, n, subcomplex=None):1137"""1138The set of simplices of dimension ``n`` of this simplicial complex.1139If the optional argument ``subcomplex`` is present, then1140return the ``n``-dimensional faces which are *not* in the1141subcomplex.11421143:param n: non-negative integer1144:param subcomplex: a subcomplex of this simplicial complex.1145Return ``n``-dimensional faces which are not in this1146subcomplex.1147:type subcomplex: optional, default ``None``11481149EXAMPLES::11501151sage: S = Set(range(1,5))1152sage: Z = SimplicialComplex(S.subsets())1153sage: Z1154Simplicial complex with vertex set (1, 2, 3, 4) and facets {(1, 2, 3, 4)}1155sage: Z.n_faces(2)1156set([(1, 3, 4), (1, 2, 3), (2, 3, 4), (1, 2, 4)])1157sage: K = SimplicialComplex([[1,2,3], [2,3,4]])1158sage: Z.n_faces(2, subcomplex=K)1159set([(1, 3, 4), (1, 2, 4)])1160"""1161if n in self.faces(subcomplex):1162return self.faces(subcomplex)[n]1163else:1164return set([])11651166def is_pure(self):1167"""1168Return ``True`` iff this simplicial complex is pure.11691170A simplicial complex is pure if and only if all of its maximal faces1171have the same dimension.11721173.. WARNING::11741175This may give the wrong answer if the simplicial complex1176was constructed with ``maximality_check`` set to ``False``.11771178EXAMPLES::11791180sage: U = SimplicialComplex([[1,2], [1, 3, 4]])1181sage: U.is_pure()1182False1183sage: X = SimplicialComplex([[0,1], [0,2], [1,2]])1184sage: X.is_pure()1185True11861187Demonstration of the warning::11881189sage: S = SimplicialComplex([[0,1], [0]], maximality_check=False)1190sage: S.is_pure()1191False1192"""1193dims = [face.dimension() for face in self._facets]1194return max(dims) == min(dims)11951196def h_vector(self):1197r"""1198The `h`-vector of this simplicial complex.11991200If the complex has dimension `d` and `(f_{-1}, f_0, f_1, ...,1201f_d)` is its `f`-vector (with `f_{-1} = 1`, representing the1202empy simplex), then the `h`-vector `(h_0, h_1, ..., h_d,1203h_{d+1})` is defined by12041205.. MATH::12061207\sum_{i=0}^{d+1} h_i x^{d+1-i} = \sum_{i=0}^{d+1} f_{i-1} (x-1)^{d+1-i}.12081209Alternatively,12101211.. MATH::12121213h_j = \sum_{i=-1}^{j-1} (-1)^{j-i-1} \binom{d-i}{j-i-1} f_i.12141215EXAMPLES:12161217The `f`- and `h`-vectors of the boundary of an octahedron are1218computed in Wikipedia's page on simplicial complexes,1219http://en.wikipedia.org/wiki/Simplicial_complex::12201221sage: square = SimplicialComplex([[0,1], [1,2], [2,3], [0,3]])1222sage: S0 = SimplicialComplex([[0], [1]])1223sage: octa = square.join(S0) # boundary of an octahedron1224sage: octa.f_vector()1225[1, 6, 12, 8]1226sage: octa.h_vector()1227[1, 3, 3, 1]1228"""1229from sage.rings.arith import binomial1230d = self.dimension()1231f = self.f_vector() # indexed starting at 0, since it's a Python list1232h = []1233for j in range(0, d+2):1234s = 01235for i in range(-1, j):1236s += (-1)**(j-i-1) * binomial(d-i, j-i-1) * f[i+1]1237h.append(s)1238return h12391240def g_vector(self):1241r"""1242The `g`-vector of this simplicial complex.12431244If the `h`-vector of the complex is `(h_0, h_1, ..., h_d,1245h_{d+1})` -- see :meth:`h_vector` -- then its `g`-vector1246`(g_0, g_1, ..., g_{[(d+1)/2]})` is defined by `g_0 = 1` and1247`g_i = h_i - h_{i-1}` for `i > 0`.12481249EXAMPLES::12501251sage: S3 = simplicial_complexes.Sphere(3).barycentric_subdivision()1252sage: S3.f_vector()1253[1, 30, 150, 240, 120]1254sage: S3.h_vector()1255[1, 26, 66, 26, 1]1256sage: S3.g_vector()1257[1, 25, 40]1258"""1259from sage.functions.other import floor1260d = self.dimension()1261h = self.h_vector()1262g = [1]1263for i in range(1, floor((d+1)/2) + 1):1264g.append(h[i] - h[i-1])1265return g12661267def flip_graph(self):1268"""1269If ``self`` is pure, then it returns the the flip graph of ``self``,1270otherwise, it returns ``None``.12711272The flip graph of a pure simplicial complex is the (undirected) graph1273with vertices being the facets, such that two facets are joined by1274an edge if they meet in a codimension `1` face.12751276The flip graph is used to detect if ``self`` is a pseudomanifold.12771278EXAMPLES::12791280sage: S0 = simplicial_complexes.Sphere(0)1281sage: G = S0.flip_graph()1282sage: G.vertices(); G.edges(labels=False)1283[(0,), (1,)]1284[((0,), (1,))]12851286sage: G = (S0.wedge(S0)).flip_graph()1287sage: G.vertices(); G.edges(labels=False)1288[(0,), ('L1',), ('R1',)]1289[((0,), ('L1',)), ((0,), ('R1',)), (('L1',), ('R1',))]12901291sage: S1 = simplicial_complexes.Sphere(1)1292sage: S2 = simplicial_complexes.Sphere(2)1293sage: G = (S1.wedge(S1)).flip_graph()1294sage: G.vertices(); G.edges(labels=False)1295[(0, 'L1'), (0, 'L2'), (0, 'R1'), (0, 'R2'), ('L1', 'L2'), ('R1', 'R2')]1296[((0, 'L1'), (0, 'L2')),1297((0, 'L1'), (0, 'R1')),1298((0, 'L1'), (0, 'R2')),1299((0, 'L1'), ('L1', 'L2')),1300((0, 'L2'), (0, 'R1')),1301((0, 'L2'), (0, 'R2')),1302((0, 'L2'), ('L1', 'L2')),1303((0, 'R1'), (0, 'R2')),1304((0, 'R1'), ('R1', 'R2')),1305((0, 'R2'), ('R1', 'R2'))]13061307sage: (S1.wedge(S2)).flip_graph() is None1308True13091310sage: G = S2.flip_graph()1311sage: G.vertices(); G.edges(labels=False)1312[(0, 1, 2), (0, 1, 3), (0, 2, 3), (1, 2, 3)]1313[((0, 1, 2), (0, 1, 3)),1314((0, 1, 2), (0, 2, 3)),1315((0, 1, 2), (1, 2, 3)),1316((0, 1, 3), (0, 2, 3)),1317((0, 1, 3), (1, 2, 3)),1318((0, 2, 3), (1, 2, 3))]13191320sage: T = simplicial_complexes.Torus()1321sage: G = T.suspension(4).flip_graph()1322sage: len(G.vertices()); len(G.edges(labels=False))13234613241611325"""1326from collections import defaultdict1327if not self.is_pure():1328return None1329d = self.dimension()1330Fs = self.facets()1331flipG = Graph()1332flipG.add_vertices(Fs)1333edges = defaultdict(list)1334# go through all codim 1 faces to build the edge1335for F in Fs:1336F_tuple = sorted(F._Simplex__set)1337for i in range(d+1):1338coF = tuple(F_tuple[:i]+F_tuple[i+1:])1339if coF in edges:1340for G in edges[coF]:1341flipG.add_edge((F,G))1342edges[coF].append(F)1343return flipG13441345def is_pseudomanifold(self):1346"""1347Return True if self is a pseudomanifold.13481349A pseudomanifold is a simplicial complex with the following properties:13501351- it is pure of some dimension `d` (all of its facets are `d`-dimensional)1352- every `(d-1)`-dimensional simplex is the face of exactly two facets1353- for every two facets `S` and `T`, there is a sequence of1354facets13551356.. math::13571358S = f_0, f_1, ..., f_n = T13591360such that for each `i`, `f_i` and `f_{i-1}` intersect in a1361`(d-1)`-simplex.13621363By convention, `S^0` is the only 0-dimensional pseudomanifold.13641365EXAMPLES::13661367sage: S0 = simplicial_complexes.Sphere(0)1368sage: S0.is_pseudomanifold()1369True1370sage: (S0.wedge(S0)).is_pseudomanifold()1371False1372sage: S1 = simplicial_complexes.Sphere(1)1373sage: S2 = simplicial_complexes.Sphere(2)1374sage: (S1.wedge(S1)).is_pseudomanifold()1375False1376sage: (S1.wedge(S2)).is_pseudomanifold()1377False1378sage: S2.is_pseudomanifold()1379True1380sage: T = simplicial_complexes.Torus()1381sage: T.suspension(4).is_pseudomanifold()1382True1383"""1384if not self.is_pure():1385return False1386d = self.dimension()1387if d == 0:1388return len(self.facets()) == 21389F = self.facets()1390X = self.n_faces(d-1)1391# is each (d-1)-simplex is the face of exactly two facets?1392for s in X:1393if len([a for a in [s.is_face(f) for f in F] if a]) != 2:1394return False1395# construct a graph with one vertex for each facet, one edge1396# when two facets intersect in a (d-1)-simplex, and see1397# whether that graph is connected.1398V = [f.set() for f in self.facets()]1399E = (lambda a,b: len(a.intersection(b)) == d)1400g = Graph([V,E])1401return g.is_connected()14021403def product(self, right, rename_vertices=True, is_mutable=True):1404"""1405The product of this simplicial complex with another one.14061407:param right: the other simplicial complex (the right-hand1408factor)14091410:param rename_vertices: If this is False, then the vertices in1411the product are the set of ordered pairs `(v,w)` where `v`1412is a vertex in ``self`` and `w` is a vertex in1413``right``. If this is ``True``, then the vertices are renamed1414as "LvRw" (e.g., the vertex (1,2) would become "L1R2").1415This is useful if you want to define the Stanley-Reisner1416ring of the complex: vertex names like (0,1) are not1417suitable for that, while vertex names like "L0R1" are.14181419:type rename_vertices: boolean; optional, default ``True``14201421:param is_mutable: Determines if the output is mutable1422:type is_mutable: boolean; optional, default ``True``14231424The vertices in the product will be the set of ordered pairs1425`(v,w)` where `v` is a vertex in self and `w` is a vertex in1426right.14271428.. WARNING::14291430If ``X`` and ``Y`` are simplicial complexes, then ``X*Y``1431returns their join, not their product.14321433EXAMPLES::14341435sage: S = SimplicialComplex([[0,1], [1,2], [0,2]]) # circle1436sage: K = SimplicialComplex([[0,1]]) # edge1437sage: S.product(K).vertices() # cylinder1438('L0R0', 'L0R1', 'L1R0', 'L1R1', 'L2R0', 'L2R1')1439sage: S.product(K, rename_vertices=False).vertices()1440((0, 0), (0, 1), (1, 0), (1, 1), (2, 0), (2, 1))1441sage: T = S.product(S) # torus1442sage: T1443Simplicial complex with 9 vertices and 18 facets1444sage: T.homology()1445{0: 0, 1: Z x Z, 2: Z}14461447These can get large pretty quickly::14481449sage: T = simplicial_complexes.Torus(); T1450Simplicial complex with vertex set (0, 1, 2, 3, 4, 5, 6) and 14 facets1451sage: K = simplicial_complexes.KleinBottle(); K1452Simplicial complex with vertex set (0, 1, 2, 3, 4, 5, 6, 7) and 16 facets1453sage: T.product(K) # long time: 5 or 6 seconds1454Simplicial complex with 56 vertices and 1344 facets1455"""1456facets = []1457for f in self._facets:1458for g in right._facets:1459facets.extend(f.product(g, rename_vertices))1460return SimplicialComplex(facets, is_mutable=is_mutable)14611462def join(self, right, rename_vertices=True, is_mutable=True):1463"""1464The join of this simplicial complex with another one.14651466The join of two simplicial complexes `S` and `T` is the1467simplicial complex `S*T` with simplices of the form `[v_0,1468..., v_k, w_0, ..., w_n]` for all simplices `[v_0, ..., v_k]` in1469`S` and `[w_0, ..., w_n]` in `T`.14701471:param right: the other simplicial complex (the right-hand factor)14721473:param rename_vertices: If this is True, the vertices in the1474join will be renamed by the formula: vertex "v" in the1475left-hand factor --> vertex "Lv" in the join, vertex "w" in1476the right-hand factor --> vertex "Rw" in the join. If this1477is false, this tries to construct the join without renaming1478the vertices; this will cause problems if the two factors1479have any vertices with names in common.14801481:type rename_vertices: boolean; optional, default ``True``14821483:param is_mutable: Determines if the output is mutable1484:type is_mutable: boolean; optional, default ``True``14851486EXAMPLES::14871488sage: S = SimplicialComplex([[0], [1]])1489sage: T = SimplicialComplex([[2], [3]])1490sage: S.join(T)1491Simplicial complex with vertex set ('L0', 'L1', 'R2', 'R3') and 4 facets1492sage: S.join(T, rename_vertices=False)1493Simplicial complex with vertex set (0, 1, 2, 3) and facets {(1, 3), (1, 2), (0, 2), (0, 3)}14941495The notation '*' may be used, as well::14961497sage: S * S1498Simplicial complex with vertex set ('L0', 'L1', 'R0', 'R1') and 4 facets1499sage: S * S * S * S * S * S * S * S1500Simplicial complex with 16 vertices and 256 facets1501"""1502facets = []1503for f in self._facets:1504for g in right._facets:1505facets.append(f.join(g, rename_vertices))1506return SimplicialComplex(facets, is_mutable=is_mutable)15071508# Use * to mean 'join':1509__mul__ = join15101511def cone(self, is_mutable=True):1512"""1513The cone on this simplicial complex.15141515:param is_mutable: Determines if the output is mutable1516:type is_mutable: boolean; optional, default ``True``15171518The cone is the simplicial complex formed by adding a new1519vertex `C` and simplices of the form `[C, v_0, ..., v_k]` for1520every simplex `[v_0, ..., v_k]` in the original simplicial1521complex. That is, the cone is the join of the original1522complex with a one-point simplicial complex.15231524EXAMPLES::15251526sage: S = SimplicialComplex([[0], [1]])1527sage: S.cone()1528Simplicial complex with vertex set ('L0', 'L1', 'R0') and facets {('L0', 'R0'), ('L1', 'R0')}1529"""1530return self.join(SimplicialComplex([["0"]], is_mutable=is_mutable),1531rename_vertices = True)15321533def suspension(self, n=1, is_mutable=True):1534r"""1535The suspension of this simplicial complex.15361537:param n: positive integer -- suspend this many times.15381539:type n: optional, default 115401541:param is_mutable: Determines if the output is mutable1542:type is_mutable: boolean; optional, default ``True``15431544The suspension is the simplicial complex formed by adding two1545new vertices `S_0` and `S_1` and simplices of the form `[S_0,1546v_0, ..., v_k]` and `[S_1, v_0, ..., v_k]` for every simplex1547`[v_0, ..., v_k]` in the original simplicial complex. That1548is, the suspension is the join of the original complex with a1549two-point simplicial complex.15501551If the simplicial complex `M` happens to be a pseudomanifold1552(see :meth:`is_pseudomanifold`), then this instead constructs1553Datta's one-point suspension (see p. 434 in the cited1554article): choose a vertex `u` in `M` and choose a new vertex1555`w` to add. Denote the join of simplices by "`*`". The1556facets in the one-point suspension are of the two forms15571558- `u * \alpha` where `\alpha` is a facet of `M` not containing1559`u`15601561- `w * \beta` where `\beta` is any facet of `M`.15621563REFERENCES:15641565- Basudeb Datta, "Minimal triangulations of manifolds",1566J. Indian Inst. Sci. 87 (2007), no. 4, 429-449.15671568EXAMPLES::15691570sage: S0 = SimplicialComplex([[0], [1]])1571sage: S0.suspension() == simplicial_complexes.Sphere(1)1572True1573sage: S3 = S0.suspension(3) # the 3-sphere1574sage: S3.homology()1575{0: 0, 1: 0, 2: 0, 3: Z}15761577For pseudomanifolds, the complex constructed here will be1578smaller than that obtained by taking the join with the15790-sphere: the join adds two vertices, while this construction1580only adds one. ::15811582sage: T = simplicial_complexes.Torus()1583sage: T.join(S0).vertices() # 9 vertices1584('L0', 'L1', 'L2', 'L3', 'L4', 'L5', 'L6', 'R0', 'R1')1585sage: T.suspension().vertices() # 8 vertices1586(0, 1, 2, 3, 4, 5, 6, 7)1587"""1588if n<0:1589raise ValueError, "n must be non-negative."1590if n==0:1591return self1592if n==1:1593if self.is_pseudomanifold():1594# Use one-point compactification of Datta. The1595# construction is a bit slower, but the resulting1596# complex is smaller.1597V = self.vertices()1598u = V[0]1599w = 01600while w in V:1601w += 11602w = Simplex([w])1603new_facets = []1604for f in self.facets():1605if u not in f:1606new_facets.append(f.join(Simplex([u]), rename_vertices=False))1607new_facets.append(f.join(w, rename_vertices=False))1608return SimplicialComplex(new_facets)1609else:1610return self.join(SimplicialComplex([["0"], ["1"]], is_mutable=is_mutable),1611rename_vertices = True)1612return self.suspension(1, is_mutable).suspension(int(n-1), is_mutable)16131614def disjoint_union(self, right, rename_vertices=True, is_mutable=True):1615"""1616The disjoint union of this simplicial complex with another one.16171618:param right: the other simplicial complex (the right-hand factor)16191620:param rename_vertices: If this is True, the vertices in the1621disjoint union will be renamed by the formula: vertex "v"1622in the left-hand factor --> vertex "Lv" in the disjoint1623union, vertex "w" in the right-hand factor --> vertex "Rw"1624in the disjoint union. If this is false, this tries to1625construct the disjoint union without renaming the vertices;1626this will cause problems if the two factors have any1627vertices with names in common.16281629:type rename_vertices: boolean; optional, default True16301631EXAMPLES::16321633sage: S1 = simplicial_complexes.Sphere(1)1634sage: S2 = simplicial_complexes.Sphere(2)1635sage: S1.disjoint_union(S2).homology()1636{0: Z, 1: Z, 2: Z}1637"""1638facets = []1639for f in self._facets:1640facets.append(tuple(["L" + str(v) for v in f]))1641for f in right._facets:1642facets.append(tuple(["R" + str(v) for v in f]))1643return SimplicialComplex(facets, is_mutable=is_mutable)16441645def wedge(self, right, rename_vertices=True, is_mutable=True):1646"""1647The wedge (one-point union) of this simplicial complex with1648another one.16491650:param right: the other simplicial complex (the right-hand factor)16511652:param rename_vertices: If this is ``True``, the vertices in the1653wedge will be renamed by the formula: first vertex in each1654are glued together and called "0". Otherwise, each vertex1655"v" in the left-hand factor --> vertex "Lv" in the wedge,1656vertex "w" in the right-hand factor --> vertex "Rw" in the1657wedge. If this is ``False``, this tries to construct the wedge1658without renaming the vertices; this will cause problems if1659the two factors have any vertices with names in common.16601661:type rename_vertices: boolean; optional, default ``True``16621663:param is_mutable: Determines if the output is mutable1664:type is_mutable: boolean; optional, default ``True``16651666.. NOTE::16671668This operation is not well-defined if ``self`` or1669``other`` is not path-connected.16701671EXAMPLES::16721673sage: S1 = simplicial_complexes.Sphere(1)1674sage: S2 = simplicial_complexes.Sphere(2)1675sage: S1.wedge(S2).homology()1676{0: 0, 1: Z, 2: Z}1677"""1678left_vertices = list(self.vertices())1679left_0 = left_vertices.pop(0)1680right_vertices = list(right.vertices())1681right_0 = right_vertices.pop(0)1682left_dict = {left_0: 0}1683right_dict = {right_0: 0}1684if rename_vertices:1685facets = []1686for v in left_vertices:1687left_dict[v] = "L" + str(v)1688for v in right_vertices:1689right_dict[v] = "R" + str(v)16901691for f in self._facets:1692facets.append(tuple([left_dict[v] for v in f]))1693for f in right._facets:1694facets.append(tuple([right_dict[v] for v in f]))1695else:1696facets = self._facets + right._facets1697return SimplicialComplex(facets, is_mutable=is_mutable)16981699def chain_complex(self, **kwds):1700"""1701The chain complex associated to this simplicial complex.17021703:param dimensions: if ``None``, compute the chain complex in all1704dimensions. If a list or tuple of integers, compute the1705chain complex in those dimensions, setting the chain groups1706in all other dimensions to zero.1707:param base_ring: commutative ring1708:type base_ring: optional, default ``ZZ``1709:param subcomplex: a subcomplex of this simplicial complex.1710Compute the chain complex relative to this subcomplex.1711:type subcomplex: optional, default empty1712:param augmented: If ``True``, return the augmented chain complex1713(that is, include a class in dimension `-1` corresponding1714to the empty cell). This is ignored if ``dimensions`` is1715specified.1716:type augmented: boolean; optional, default ``False``1717:param cochain: If ``True``, return the cochain complex (that is,1718the dual of the chain complex).1719:type cochain: boolean; optional, default ``False``1720:param verbose: If ``True``, print some messages as the chain1721complex is computed.1722:type verbose: boolean; optional, default ``False``1723:param check_diffs: If ``True``, make sure that the chain complex1724is actually a chain complex: the differentials are1725composable and their product is zero.1726:type check_diffs: boolean; optional, default ``False``17271728.. NOTE::17291730If subcomplex is nonempty, then the argument ``augmented``1731has no effect: the chain complex relative to a nonempty1732subcomplex is zero in dimension `-1`.17331734EXAMPLES::17351736sage: circle = SimplicialComplex([[0,1], [1,2], [0, 2]])1737sage: circle.chain_complex()1738Chain complex with at most 2 nonzero terms over Integer Ring1739sage: circle.chain_complex()._latex_()1740'\\Bold{Z}^{3} \\xrightarrow{d_{1}} \\Bold{Z}^{3}'1741sage: circle.chain_complex(base_ring=QQ, augmented=True)1742Chain complex with at most 3 nonzero terms over Rational Field1743"""1744augmented = kwds.get('augmented', False)1745cochain = kwds.get('cochain', False)1746verbose = kwds.get('verbose', False)1747check_diffs = kwds.get('check_diffs', False)1748base_ring = kwds.get('base_ring', ZZ)1749dimensions = kwds.get('dimensions', None)1750subcomplex = kwds.get('subcomplex', None)17511752# initialize subcomplex1753if subcomplex is None:1754subcomplex = SimplicialComplex(is_mutable=False)1755else:1756# subcomplex is not empty, so don't augment the chain complex1757augmented = False1758# Use an immutable copy of the subcomplex1759if not subcomplex._is_mutable:1760subcomplex = SimplicialComplex(subcomplex._facets, maximality_check=False,1761sort_facets=False, is_mutable=False)1762# now construct the range of dimensions in which to compute1763if dimensions is None:1764dimensions = range(0, self.dimension()+1)1765first = 01766else:1767augmented = False1768first = dimensions[0]1769differentials = {}1770# in the chain complex, compute the first dimension by hand,1771# and don't cache it: it may be differ from situation to1772# situation because of boundary effects.1773current = None1774current_dim = None1775if augmented: # then first == 01776current = list(self.n_faces(0, subcomplex=subcomplex))1777current_dim = 01778if cochain:1779differentials[-1] = matrix(base_ring, len(current), 1,1780[1]*len(current))1781else:1782differentials[0] = matrix(base_ring, 1, len(current),1783[1]*len(current))1784elif first == 0 and not augmented:1785current = list(self.n_faces(0, subcomplex=subcomplex))1786current_dim = 01787if not cochain:1788differentials[0] = matrix(base_ring, 0, len(current))1789else: # first > 01790current = list(self.n_faces(first, subcomplex=subcomplex))1791current_dim = first1792if not cochain:1793differentials[first] = matrix(base_ring, 0, len(current))1794for n in dimensions[1:]:1795if verbose:1796print " starting dimension %s" % n1797if (n, subcomplex) in self._complex:1798if cochain:1799differentials[n-1] = self._complex[(n, subcomplex)].transpose().change_ring(base_ring)1800mat = differentials[n-1]1801else:1802differentials[n] = self._complex[(n, subcomplex)].change_ring(base_ring)1803mat = differentials[n]1804if verbose:1805print " boundary matrix (cached): it's %s by %s." % (mat.nrows(), mat.ncols())1806else:1807# 'current' is the list of faces in dimension n1808#1809# 'old' is a dictionary, with keys the faces in the1810# previous dimension (dim n-1 for the chain complex,1811# n+1 for the cochain complex), values the integers 0,1812# 1, 2, ... (the index of the face). finding an entry1813# in a dictionary seems to be faster than finding the1814# index of an entry in a list.1815if current_dim == n-1:1816old = dict(zip(current, range(len(current))))1817else:1818set_of_faces = list(self.n_faces(n-1, subcomplex=subcomplex))1819old = dict(zip(set_of_faces, range(len(set_of_faces))))1820current = list(self.n_faces(n, subcomplex=subcomplex))1821current_dim = n1822# construct matrix. it is easiest to construct it as1823# a sparse matrix, specifying which entries are1824# nonzero via a dictionary.1825matrix_data = {}1826col = 01827if len(old) and len(current):1828for simplex in current:1829for i in range(n+1):1830face_i = simplex.face(i)1831try:1832matrix_data[(old[face_i], col)] = (-1)**i1833except KeyError:1834pass1835col += 11836mat = matrix(ZZ, len(old), len(current), matrix_data)1837if cochain:1838self._complex[(n, subcomplex)] = mat1839differentials[n-1] = mat.transpose().change_ring(base_ring)1840else:1841self._complex[(n, subcomplex)] = mat1842differentials[n] = mat.change_ring(base_ring)1843if verbose:1844print " boundary matrix computed: it's %s by %s." % (mat.nrows(), mat.ncols())1845# now for the cochain complex, compute the last dimension by1846# hand, and don't cache it.1847if cochain:1848n = dimensions[-1] + 11849if current_dim != n-1:1850current = list(self.n_faces(n-1, subcomplex=subcomplex))1851differentials[n-1] = matrix(base_ring, 0, len(current))1852# finally, return the chain complex1853if cochain:1854return ChainComplex(data=differentials, degree=1, **kwds)1855else:1856return ChainComplex(data=differentials, degree=-1, **kwds)18571858def _homology_(self, dim=None, **kwds):1859"""1860The reduced homology of this simplicial complex.18611862:param dim: If ``None``, then return the homology in every1863dimension. If ``dim`` is an integer or list, return the1864homology in the given dimensions. (Actually, if ``dim`` is1865a list, return the homology in the range from ``min(dim)``1866to ``max(dim)``.)18671868:type dim: integer or list of integers or ``None``; optional,1869default ``None``18701871:param base_ring: commutative ring. Must be ``ZZ`` or a field.18721873:type base_ring: optional, default ``ZZ``18741875:param subcomplex: a subcomplex of this simplicial complex.1876Compute homology relative to this subcomplex.18771878:type subcomplex: optional, default ``None``18791880:param cohomology: If ``True``, compute cohomology rather than1881homology.18821883:type cohomology: boolean; optional, default ``False``18841885:param enlarge: If ``True``, find a new subcomplex homotopy1886equivalent to, and probably larger than, the given one.18871888:type enlarge: boolean; optional, default ``True``18891890:param algorithm: The options are ``'auto'``, ``'dhsw'``,1891``'pari'`` or ``'no_chomp'``. If ``'auto'``, first try CHomP,1892then use the Dumas, Heckenbach, Saunders, and Welker elimination1893algorithm for large matrices, Pari for small ones. If1894``'no_chomp'``, then don't try CHomP, but behave the same1895otherwise. If ``'pari'``, then compute elementary divisors1896using Pari. If ``'dhsw'``, then use the DHSW algorithm to1897compute elementary divisors. (As of this writing, CHomP is1898by far the fastest option, followed by the ``'auto'`` or1899``'no_chomp'`` setting of using DHSW for large matrices and1900Pari for small ones.)19011902:type algorithm: string; optional, default ``'auto'``19031904:param verbose: If ``True``, print some messages as the homology1905is computed.19061907:type verbose: boolean; optional, default ``False``19081909Algorithm: if ``subcomplex`` is ``None``, replace it with a facet1910-- a contractible subcomplex of the original complex. Then no1911matter what ``subcomplex`` is, replace it with a subcomplex1912`L` which is homotopy equivalent and as large as possible.1913Compute the homology of the original complex relative to `L`:1914if `L` is large, then the relative chain complex will be small1915enough to speed up computations considerably.19161917EXAMPLES::19181919sage: circle = SimplicialComplex([[0,1], [1,2], [0, 2]])1920sage: circle._homology_()1921{0: 0, 1: Z}1922sage: sphere = SimplicialComplex([[0,1,2,3]])1923sage: sphere.remove_face([0,1,2,3])1924sage: sphere1925Simplicial complex with vertex set (0, 1, 2, 3) and facets {(0, 2, 3), (0, 1, 2), (1, 2, 3), (0, 1, 3)}1926sage: sphere._homology_()1927{0: 0, 1: 0, 2: Z}19281929Another way to get a two-sphere: take a two-point space and take its1930three-fold join with itself::19311932sage: S = SimplicialComplex([[0], [1]])1933sage: (S*S*S)._homology_(dim=2, cohomology=True)1934Z19351936Relative homology::19371938sage: T = SimplicialComplex([[0,1,2]])1939sage: U = SimplicialComplex([[0,1], [1,2], [0,2]])1940sage: T._homology_(subcomplex=U)1941{0: 0, 1: 0, 2: Z}1942"""1943from sage.modules.all import VectorSpace1944from sage.homology.homology_group import HomologyGroup19451946base_ring = kwds.get('base_ring', ZZ)1947cohomology = kwds.get('cohomology', False)1948enlarge = kwds.get('enlarge', True)1949verbose = kwds.get('verbose', False)1950subcomplex = kwds.get('subcomplex', None)19511952if dim is not None:1953if isinstance(dim, (list, tuple)):1954low = min(dim) - 11955high = max(dim) + 21956else:1957low = dim - 11958high = dim + 21959dims = range(low, high)1960else:1961dims = None19621963if verbose:1964print "starting calculation of the homology of this",1965print "%s-dimensional simplicial complex" % self.dimension()1966if subcomplex is None:1967if enlarge:1968if verbose:1969print "Constructing contractible subcomplex..."1970L = self._contractible_subcomplex(verbose=verbose)1971if verbose:1972print "Done finding contractible subcomplex."1973vec = [len(self.n_faces(n-1, subcomplex=L)) for n in range(self.dimension()+2)]1974print "The difference between the f-vectors is:"1975print " %s" % vec1976else:1977L = SimplicialComplex([[self.vertices().tuple()[0]]])1978else:1979if enlarge:1980if verbose:1981print "Enlarging subcomplex..."1982L = self._enlarge_subcomplex(subcomplex, verbose=verbose)1983if verbose:1984print "Done enlarging subcomplex:"1985else:1986L = subcomplex1987L.set_immutable()19881989if verbose:1990print "Computing the chain complex..."1991kwds['subcomplex']=L1992C = self.chain_complex(dimensions=dims, augmented=True,1993cochain=cohomology, **kwds)1994if verbose:1995print " Done computing the chain complex. "1996print "Now computing homology..."1997if 'subcomplex' in kwds:1998del kwds['subcomplex']1999answer = C.homology(**kwds)20002001if dim is None:2002dim = range(self.dimension()+1)2003zero = HomologyGroup(0, base_ring)2004if isinstance(dim, (list, tuple)):2005return dict([d, answer.get(d, zero)] for d in dim)2006return answer.get(dim, zero)20072008def add_face(self, face):2009"""2010Add a face to this simplicial complex20112012:param face: a subset of the vertex set20132014This *changes* the simplicial complex, adding a new face and all2015of its subfaces.20162017EXAMPLES::20182019sage: X = SimplicialComplex([[0,1], [0,2]])2020sage: X.add_face([0,1,2,]); X2021Simplicial complex with vertex set (0, 1, 2) and facets {(0, 1, 2)}2022sage: Y = SimplicialComplex(); Y2023Simplicial complex with vertex set () and facets {()}2024sage: Y.add_face([0,1])2025sage: Y.add_face([1,2,3])2026sage: Y2027Simplicial complex with vertex set (0, 1, 2, 3) and facets {(1, 2, 3), (0, 1)}20282029If you add a face which is already present, there is no effect::20302031sage: Y.add_face([1,3]); Y2032Simplicial complex with vertex set (0, 1, 2, 3) and facets {(1, 2, 3), (0, 1)}20332034Check that the bug reported at :trac:`14354` has been fixed::20352036sage: T = SimplicialComplex([range(1,5)]).n_skeleton(1)2037sage: T.homology()2038{0: 0, 1: Z x Z x Z}2039sage: T.add_face([1,2,3])2040sage: T.homology()2041{0: 0, 1: Z x Z, 2: 0}20422043Check we've fixed the bug reported at :trac:`14578`::20442045sage: t0 = SimplicialComplex()2046sage: t0.add_face(('a', 'b'))2047sage: t0.add_face(('c', 'd', 'e'))2048sage: t0.add_face(('e', 'f', 'c'))2049sage: t0.homology()2050{0: Z, 1: 0, 2: 0}2051"""2052if not self._is_mutable:2053raise ValueError("This simplicial complex is not mutable")20542055if self._sorted:2056new_face = Simplex(sorted(face))2057else:2058new_face = Simplex(face)20592060face_is_maximal = True2061for other in self._facets:2062if face_is_maximal:2063face_is_maximal = not new_face.is_face(other)2064if face_is_maximal:2065# remove any old facets which are no longer maximal2066Facets = list(self._facets)2067for old_face in self._facets:2068if old_face.is_face(new_face):2069Facets.remove(old_face)2070# add new_face to facet list2071Facets.append(new_face)2072self._facets = Facets20732074# Update the vertex set2075from sage.misc.misc import union20762077if self._sorted:2078self._vertex_set = Simplex(sorted(reduce(union, [self._vertex_set, new_face])))2079else:2080self._vertex_set = Simplex(reduce(union, [self._vertex_set, new_face]))20812082# update self._faces if necessary2083if None in self._faces:2084all_new_faces = SimplicialComplex([new_face]).faces()2085for dim in range(0, new_face.dimension()+1):2086if dim in self._faces[None]:2087self._faces[None][dim] = self._faces[None][dim].union(all_new_faces[dim])2088else:2089self._faces[None][dim] = all_new_faces[dim]2090# update self._graph if necessary2091if self._graph is not None:2092d = new_face.dimension()+12093for i in range(d):2094for j in range(i+1,d):2095self._graph.add_edge(new_face[i],new_face[j])2096self._complex = {}2097self.__contractible = None2098self.__enlarged = {}20992100def remove_face(self, face):2101"""2102Remove a face from this simplicial complex and return the2103resulting simplicial complex.21042105:param face: a face of the simplicial complex21062107This *changes* the simplicial complex.21082109ALGORITHM:21102111The facets of the new simplicial complex are2112the facets of the original complex not containing ``face``,2113together with those of ``link(face)*boundary(face)``.21142115EXAMPLES::21162117sage: S = range(1,5)2118sage: Z = SimplicialComplex([S]); Z2119Simplicial complex with vertex set (1, 2, 3, 4) and facets {(1, 2, 3, 4)}2120sage: Z.remove_face([1,2])2121sage: Z2122Simplicial complex with vertex set (1, 2, 3, 4) and facets {(1, 3, 4), (2, 3, 4)}21232124sage: S = SimplicialComplex([[0,1,2],[2,3]])2125sage: S2126Simplicial complex with vertex set (0, 1, 2, 3) and facets {(0, 1, 2), (2, 3)}2127sage: S.remove_face([0,1,2])2128sage: S2129Simplicial complex with vertex set (0, 1, 2, 3) and facets {(1, 2), (2, 3), (0, 2), (0, 1)}2130"""2131if not self._is_mutable:2132raise ValueError("This simplicial complex is not mutable")21332134simplex = Simplex(face)2135facets = self.facets()2136if all([not simplex.is_face(F) for F in facets]):2137# face is not in self: nothing to remove2138return self2139link = self.link(simplex)2140join_facets = []2141for f in simplex.faces():2142for g in link.facets():2143join_facets.append(f.join(g, rename_vertices=False))2144# join_facets is the list of facets in the join bdry(face) * link(face)2145remaining = join_facets + [elem for elem in facets if not simplex.is_face(elem)]21462147# Check to see if there are any non-maximial faces2148# build set of facets2149self._facets = []2150for f in remaining:2151face = Simplex(f)2152face_is_maximal = True2153faces_to_be_removed = []2154for other in self._facets:2155if other.is_face(face):2156faces_to_be_removed.append(other)2157elif face_is_maximal:2158face_is_maximal = not face.is_face(other)2159for x in faces_to_be_removed:2160self._facets.remove(x)2161face = Simplex(sorted(face.tuple()))2162if face_is_maximal:2163self._facets.append(face)2164# if no maximal faces, add the empty face as a facet2165if len(remaining) == 0:2166self._facets.append(Simplex(-1))21672168# Recreate the vertex set2169from sage.misc.misc import union2170if self._sorted:2171self._vertex_set = Simplex(sorted(reduce(union, self._facets)))2172else:2173self._vertex_set = Simplex(reduce(union, self._facets))21742175# Update self._faces and self._graph if necessary2176if None in self._faces:2177self._faces = {}2178self.faces()2179if self._graph is not None:2180# Only if removing a 1 or 2 dim face will the graph be affected2181if len(face) == 1:2182self._graph.delete_vertex(face[0])2183self._graph.add_vertex(face[0])2184elif len(face) == 2:2185self._graph.delete_edge(face[0], face[1])2186self._complex = {}2187self.__contractible = None2188self.__enlarged = {}21892190def connected_sum(self, other, is_mutable=True):2191"""2192The connected sum of this simplicial complex with another one.21932194:param other: another simplicial complex2195:param is_mutable: Determines if the output is mutable2196:type is_mutable: boolean; optional, default ``True``2197:return: the connected sum ``self # other``21982199.. WARNING::22002201This does not check that ``self`` and ``other`` are manifolds,2202only that their facets all have the same dimension. Since a2203(more or less) random facet is chosen from each complex and2204then glued together, this method may return random2205results if applied to non-manifolds, depending on which2206facet is chosen.22072208Algorithm: a facet is chosen from each surface, and removed.2209The vertices of these two facets are relabeled to2210``(0,1,...,dim)``. Of the remaining vertices, the ones from2211the left-hand factor are renamed by prepending an "L", and2212similarly the remaining vertices in the right-hand factor are2213renamed by prepending an "R".22142215EXAMPLES::22162217sage: S1 = simplicial_complexes.Sphere(1)2218sage: S1.connected_sum(S1.connected_sum(S1)).homology()2219{0: 0, 1: Z}2220sage: P = simplicial_complexes.RealProjectivePlane(); P2221Simplicial complex with vertex set (0, 1, 2, 3, 4, 5) and 10 facets2222sage: P.connected_sum(P) # the Klein bottle2223Simplicial complex with 9 vertices and 18 facets22242225The notation '+' may be used for connected sum, also::22262227sage: P + P # the Klein bottle2228Simplicial complex with 9 vertices and 18 facets2229sage: (P + P).homology()[1]2230Z x C22231"""2232if not (self.is_pure() and other.is_pure() and2233self.dimension() == other.dimension()):2234raise ValueError, "Complexes are not pure of the same dimension."2235# first find a top-dimensional simplex to remove from each surface2236keep_left = self._facets[0]2237keep_right = other._facets[0]2238# construct the set of vertices:2239left = set(self.vertices()).difference(set(keep_left))2240right = set(other.vertices()).difference(set(keep_right))2241# construct the set of facets:2242left = set(self._facets).difference(set([keep_left]))2243right = set(other._facets).difference(set([keep_right]))2244facet_set = ([[rename_vertex(v, keep=list(keep_left))2245for v in face] for face in left]2246+ [[rename_vertex(v, keep=list(keep_right), left=False)2247for v in face] for face in right])2248# return the new surface2249return SimplicialComplex(facet_set, is_mutable=is_mutable)22502251__add__ = connected_sum22522253def link(self, simplex, is_mutable=True):2254"""2255The link of a simplex in this simplicial complex.22562257The link of a simplex `F` is the simplicial complex formed by2258all simplices `G` which are disjoint from `F` but for which `F2259\cup G` is a simplex.22602261:param simplex: a simplex in this simplicial complex.2262:param is_mutable: Determines if the output is mutable2263:type is_mutable: boolean; optional, default ``True``22642265EXAMPLES::22662267sage: X = SimplicialComplex([[0,1,2], [1,2,3]])2268sage: X.link(Simplex([0]))2269Simplicial complex with vertex set (1, 2) and facets {(1, 2)}2270sage: X.link([1,2])2271Simplicial complex with vertex set (0, 3) and facets {(3,), (0,)}2272sage: Y = SimplicialComplex([[0,1,2,3]])2273sage: Y.link([1])2274Simplicial complex with vertex set (0, 2, 3) and facets {(0, 2, 3)}2275"""2276faces = []2277s = Simplex(simplex)2278for f in self._facets:2279if s.is_face(f):2280faces.append(Simplex(list(f.set().difference(s.set()))))2281return SimplicialComplex(faces, is_mutable=is_mutable)22822283def is_cohen_macaulay(self, ncpus=0):2284r"""2285Returns True if ``self`` is Cohen-Macaulay, i.e., if2286`\tilde{H}_i(\operatorname{lk}_\Delta(F);\ZZ) = 0` for all2287`F \in \Delta` and `i < \operatorname{dim}\operatorname{lk}_\Delta(F)`.2288Here, `\Delta` is ``self``, and `\operatorname{lk}` denotes the2289link operator on ``self``.22902291INPUT:22922293- ``ncpus`` -- (default: 0) number of cpus used for the2294computation. If this is 0, determine the number of cpus2295automatically based on the hardware being used.22962297For finite simplicial complexes, this is equivalent to the2298statement that the Stanley-Reisner ring of ``self`` is2299Cohen-Macaulay.23002301EXAMPLES:23022303Spheres are Cohen-Macaulay::23042305sage: S = SimplicialComplex([[1,2],[2,3],[3,1]])2306sage: S.is_cohen_macaulay(ncpus=3)2307True23082309The following example is taken from Bruns, Herzog - Cohen-Macaulay2310rings, Figure 5.3::23112312sage: S = SimplicialComplex([[1,2,3],[1,4,5]])2313sage: S.is_cohen_macaulay(ncpus=3)2314...2315False2316"""2317from sage.parallel.decorate import parallel2318from sage.rings.rational_field import QQ23192320if ncpus == 0:2321import os2322try:2323ncpus = int(os.environ['SAGE_NUM_THREADS'])2324except KeyError:2325ncpus = 123262327facs = [ x for x in self.face_iterator() ]2328n = len(facs)2329facs_divided = [ [] for i in range(ncpus) ]2330for i in range(n):2331facs_divided[i%ncpus].append(facs[i])23322333def all_homologies_vanish(F):2334S = self.link(F)2335H = S.homology(base_ring=QQ)2336return all( H[j].dimension() == 0 for j in xrange(S.dimension()) )23372338@parallel(ncpus=ncpus)2339def all_homologies_in_list_vanish(Fs):2340return all( all_homologies_vanish(F) for F in Fs )23412342return all( answer[1] for answer in all_homologies_in_list_vanish(facs_divided) )23432344def effective_vertices(self):2345"""2346The set of vertices belonging to some face. Returns the list of2347vertices.23482349.. WARNING::23502351This method is deprecated. See :trac:`12587`.23522353EXAMPLES::23542355sage: S = SimplicialComplex([[0,1,2,3],[6,7]])2356sage: S2357Simplicial complex with vertex set (0, 1, 2, 3, 6, 7) and facets {(6, 7), (0, 1, 2, 3)}2358sage: S.effective_vertices()2359doctest:1: DeprecationWarning: effective_vertices is deprecated. Use vertices instead2360See http://trac.sagemath.org/12587 for details.2361(0, 1, 2, 3, 6, 7)2362"""2363from sage.misc.superseded import deprecation2364deprecation(12587, "effective_vertices is deprecated. Use vertices instead")2365return self._vertex_set23662367def generated_subcomplex(self,sub_vertex_set, is_mutable=True):2368"""2369Returns the largest sub-simplicial complex of ``self`` containing2370exactly ``sub_vertex_set`` as vertices.23712372:param sub_vertex_set: The sub-vertex set.2373:param is_mutable: Determines if the output is mutable2374:type is_mutable: boolean; optional, default ``True``23752376EXAMPLES::23772378sage: S = simplicial_complexes.Sphere(2)2379sage: S2380Simplicial complex with vertex set (0, 1, 2, 3) and facets {(0, 2, 3), (0, 1, 2), (1, 2, 3), (0, 1, 3)}2381sage: S.generated_subcomplex([0,1,2])2382Simplicial complex with vertex set (0, 1, 2) and facets {(0, 1, 2)}23832384"""2385if not self.vertices().set().issuperset(sub_vertex_set):2386raise ValueError, "input must be a subset of the vertex set."2387faces = []2388for i in range(self.dimension()+1):2389for j in self.faces()[i]:2390if j.set().issubset(sub_vertex_set):2391faces.append(j)2392return SimplicialComplex(faces, maximality_check=True,2393is_mutable=is_mutable)23942395def _complement(self, simplex):2396"""2397Return the complement of a simplex in the vertex set of this2398simplicial complex.23992400:param simplex: a simplex (need not be in the simplicial complex)24012402OUTPUT: its complement: the simplex formed by the vertices not2403contained in ``simplex``.24042405Note that this only depends on the vertex set of the2406simplicial complex, not on its simplices.24072408EXAMPLES::24092410sage: X = SimplicialComplex([[0,1,2,3,4,5]])2411sage: X._complement([1,2,3])2412(0, 4, 5)2413sage: X._complement([0,1,3,4])2414(2, 5)2415sage: X._complement([0,4,1,3])2416(2, 5)2417"""2418return Simplex(set(self.vertices()).difference(simplex))24192420def _transpose_simplices(self, *simplices):2421"""2422Given tuple ``L`` of simplices, returns new list, where each2423simplex is formed by taking a vertex from each simplex from2424``L``.24252426:param simplices: a bunch of simplices24272428If ``simplices`` consists of `(f_0, f_1, f_2, ...)`, then the2429output consists of all possible simplices of the form `(v_0,2430v_1, v_2, ...)`, where `v_i` is a vertex of `f_i`. If a2431vertex appears more than once in such a simplex, remove all2432but one of its appearances. If such a simplex contains others2433already produced, then ignore that larger simplex -- the2434output should be a list of minimal simplices constructed in2435this way.24362437This is used in computing the minimal nonfaces and hence the2438Stanley-Reisner ring.24392440Note that this only depends on the vertex set of the2441simplicial complex, not on its simplices.24422443I don't know if there is a standard name for this, but it2444looked sort of like the transpose of a matrix; hence the name2445for this method.24462447EXAMPLES::24482449sage: X = SimplicialComplex()2450sage: X._transpose_simplices([1,2])2451[(1,), (2,)]2452sage: X._transpose_simplices([1,2], [3,4])2453[(1, 3), (1, 4), (2, 3), (2, 4)]24542455In the following example, one can construct the simplices2456``(1,2)`` and ``(1,3)``, but you can also construct ``(1,1) = (1,)``,2457which is a face of both of the others. So the answer omits2458``(1,2)`` and ``(1,3)``::24592460sage: X._transpose_simplices([1,2], [1,3])2461[(1,), (2, 3)]2462"""2463answer = []2464if len(simplices) == 1:2465answer = [Simplex((v,)) for v in simplices[0]]2466elif len(simplices) > 1:2467face = simplices[0]2468rest = simplices[1:]2469for v in face:2470for partial in self._transpose_simplices(*rest):2471if v not in partial:2472L = [v] + list(partial)2473L.sort()2474simplex = Simplex(L)2475else:2476simplex = partial2477add_simplex = True2478simplices_to_delete = []2479for already in answer:2480if add_simplex:2481if already.is_face(simplex):2482add_simplex = False2483if add_simplex and simplex.is_face(already):2484simplices_to_delete.append(already)2485if add_simplex:2486answer.append(simplex)2487for x in simplices_to_delete:2488answer.remove(x)2489return answer24902491def minimal_nonfaces(self):2492"""2493Set consisting of the minimal subsets of the vertex set of2494this simplicial complex which do not form faces.24952496Algorithm: first take the complement (within the vertex set)2497of each facet, obtaining a set `(f_1, f_2, ...)` of simplices.2498Now form the set of all simplices of the form `(v_1, v_2,2499...)` where vertex `v_i` is in face `f_i`. This set will2500contain the minimal nonfaces and may contain some non-minimal2501nonfaces also, so loop through the set to find the minimal2502ones. (The last two steps are taken care of by the2503``_transpose_simplices`` routine.)25042505This is used in computing the2506:meth:`Stanley-Reisner ring<stanley_reisner_ring>` and the2507:meth:`Alexander dual<alexander_dual>`.25082509EXAMPLES::25102511sage: X = SimplicialComplex([[1,3],[1,2]])2512sage: X.minimal_nonfaces()2513{(2, 3)}2514sage: Y = SimplicialComplex([[0,1], [1,2], [2,3], [3,0]])2515sage: Y.minimal_nonfaces()2516{(1, 3), (0, 2)}2517"""2518complements = [self._complement(facet) for facet in self._facets]2519return Set(self._transpose_simplices(*complements))25202521def _stanley_reisner_base_ring(self, base_ring=ZZ):2522"""2523The polynomial algebra of which the Stanley-Reisner ring is a2524quotient.25252526:param base_ring: a commutative ring2527:type base_ring: optional, default ``ZZ``2528:return: a polynomial algebra with coefficients in base_ring,2529with one generator for each vertex in the simplicial complex.25302531See the documentation for :meth:`stanley_reisner_ring` for a2532warning about the names of the vertices.25332534EXAMPLES::25352536sage: X = SimplicialComplex([[1,2], [0], [3]])2537sage: X._stanley_reisner_base_ring()2538Multivariate Polynomial Ring in x0, x1, x2, x3 over Integer Ring2539sage: Y = SimplicialComplex([['a', 'b', 'c']])2540sage: Y._stanley_reisner_base_ring(base_ring=QQ)2541Multivariate Polynomial Ring in a, c, b over Rational Field2542"""2543return PolynomialRing(base_ring, self._gen_dict.values())25442545def stanley_reisner_ring(self, base_ring=ZZ):2546"""2547The Stanley-Reisner ring of this simplicial complex.25482549:param base_ring: a commutative ring2550:type base_ring: optional, default ``ZZ``2551:return: a quotient of a polynomial algebra with coefficients2552in ``base_ring``, with one generator for each vertex in the2553simplicial complex, by the ideal generated by the products2554of those vertices which do not form faces in it.25552556Thus the ideal is generated by the products corresponding to2557the minimal nonfaces of the simplicial complex.25582559.. WARNING::25602561This may be quite slow!25622563Also, this may behave badly if the vertices have the2564'wrong' names. To avoid this, define the simplicial complex2565at the start with the flag ``name_check`` set to ``True``.25662567More precisely, this is a quotient of a polynomial ring2568with one generator for each vertex. If the name of a2569vertex is a non-negative integer, then the corresponding2570polynomial generator is named ``'x'`` followed by that integer2571(e.g., ``'x2'``, ``'x3'``, ``'x5'``, ...). Otherwise, the2572polynomial generators are given the same names as the vertices.2573Thus if the vertex set is ``(2, 'x2')``, there will be problems.25742575EXAMPLES::25762577sage: X = SimplicialComplex([[0,1], [1,2], [2,3], [0,3]])2578sage: X.stanley_reisner_ring()2579Quotient of Multivariate Polynomial Ring in x0, x1, x2, x3 over Integer Ring by the ideal (x1*x3, x0*x2)2580sage: Y = SimplicialComplex([[0,1,2,3,4]]); Y2581Simplicial complex with vertex set (0, 1, 2, 3, 4) and facets {(0, 1, 2, 3, 4)}2582sage: Y.add_face([0,1,2,3,4])2583sage: Y.stanley_reisner_ring(base_ring=QQ)2584Multivariate Polynomial Ring in x0, x1, x2, x3, x4 over Rational Field2585"""2586R = self._stanley_reisner_base_ring(base_ring)2587products = []2588for f in self.minimal_nonfaces():2589prod = 12590for v in f:2591prod *= R(self._gen_dict[v])2592products.append(prod)2593return R.quotient(products)25942595def alexander_dual(self, is_mutable=True):2596"""2597The Alexander dual of this simplicial complex: according to2598the Macaulay2 documentation, this is the simplicial complex2599whose faces are the complements of its nonfaces.26002601Thus find the minimal nonfaces and take their complements to2602find the facets in the Alexander dual.26032604:param is_mutable: Determines if the output is mutable2605:type is_mutable: boolean; optional, default ``True``26062607EXAMPLES::26082609sage: Y = SimplicialComplex([[i] for i in range(5)]); Y2610Simplicial complex with vertex set (0, 1, 2, 3, 4) and facets {(4,), (2,), (3,), (0,), (1,)}2611sage: Y.alexander_dual()2612Simplicial complex with vertex set (0, 1, 2, 3, 4) and 10 facets2613sage: X = SimplicialComplex([[0,1], [1,2], [2,3], [3,0]])2614sage: X.alexander_dual()2615Simplicial complex with vertex set (0, 1, 2, 3) and facets {(1, 3), (0, 2)}2616"""2617nonfaces = self.minimal_nonfaces()2618return SimplicialComplex([self._complement(f) for f in nonfaces], is_mutable=is_mutable)26192620def barycentric_subdivision(self):2621"""2622The barycentric subdivision of this simplicial complex.26232624See http://en.wikipedia.org/wiki/Barycentric_subdivision for a2625definition.26262627EXAMPLES::26282629sage: triangle = SimplicialComplex([[0,1], [1,2], [0, 2]])2630sage: hexagon = triangle.barycentric_subdivision()2631sage: hexagon2632Simplicial complex with 6 vertices and 6 facets2633sage: hexagon.homology(1) == triangle.homology(1)2634True26352636Barycentric subdivisions can get quite large, since each2637`n`-dimensional facet in the original complex produces2638`(n+1)!` facets in the subdivision::26392640sage: S4 = simplicial_complexes.Sphere(4)2641sage: S42642Simplicial complex with vertex set (0, 1, 2, 3, 4, 5) and 6 facets2643sage: S4.barycentric_subdivision()2644Simplicial complex with 62 vertices and 720 facets2645"""2646return self.face_poset().order_complex()26472648def graph(self):2649"""2650The 1-skeleton of this simplicial complex, as a graph.26512652.. WARNING::26532654This may give the wrong answer if the simplicial complex2655was constructed with ``maximality_check`` set to ``False``.26562657EXAMPLES::26582659sage: S = SimplicialComplex([[0,1,2,3]])2660sage: G = S.graph(); G2661Graph on 4 vertices2662sage: G.edges()2663[(0, 1, None), (0, 2, None), (0, 3, None), (1, 2, None), (1, 3, None), (2, 3, None)]2664"""2665if self._graph is None:2666edges = self.n_faces(1)2667vertices = map(min, filter(lambda f: f.dimension() == 0, self._facets))2668used_vertices = [] # vertices which are in an edge2669d = {}2670for e in edges:2671v = min(e)2672if v in d:2673d[v].append(max(e))2674else:2675d[v] = [max(e)]2676used_vertices.extend(list(e))2677for v in vertices:2678if v not in used_vertices:2679d[v] = []2680self._graph = Graph(d)2681return self._graph26822683def delta_complex(self, sort_simplices=False):2684r"""2685Returns ``self`` as a `\Delta`-complex. The `\Delta`-complex2686is essentially identical to the simplicial complex: it has2687same simplices with the same boundaries.26882689:param sort_simplices: if ``True``, sort the list of simplices in2690each dimension2691:type sort_simplices: boolean; optional, default ``False``26922693EXAMPLES::26942695sage: T = simplicial_complexes.Torus()2696sage: Td = T.delta_complex()2697sage: Td2698Delta complex with 7 vertices and 43 simplices2699sage: T.homology() == Td.homology()2700True2701"""2702from delta_complex import DeltaComplex2703data = {}2704dim = self.dimension()2705n_cells = self.n_cells(dim)2706if sort_simplices:2707n_cells.sort()2708for n in range(dim, -1, -1):2709bdries = self.n_cells(n-1)2710if sort_simplices:2711bdries.sort()2712data[n] = []2713for f in n_cells:2714data[n].append([bdries.index(f.face(i)) for i in range(n+1)])2715n_cells = bdries2716return DeltaComplex(data)27172718def is_flag_complex(self):2719"""2720Returns ``True`` if and only if ``self`` is a flag complex.27212722A flag complex is a simplicial complex that is the largest simplicial2723complex on its 1-skeleton. Thus a flag complex is the clique complex2724of its graph.27252726EXAMPLES::27272728sage: h = Graph({0:[1,2,3,4],1:[2,3,4],2:[3]})2729sage: x = h.clique_complex()2730sage: x2731Simplicial complex with vertex set (0, 1, 2, 3, 4) and facets {(0, 1, 4), (0, 1, 2, 3)}2732sage: x.is_flag_complex()2733True27342735sage: X = simplicial_complexes.ChessboardComplex(3,3)2736sage: X.is_flag_complex()2737True2738"""2739return self == self.graph().clique_complex()27402741def is_connected(self):2742"""2743Returns ``True`` if and only if ``self`` is connected.27442745.. WARNING::27462747This may give the wrong answer if the simplicial complex2748was constructed with ``maximality_check`` set to ``False``.27492750EXAMPLES::27512752sage: V = SimplicialComplex([[0,1,2],[3]])2753sage: V2754Simplicial complex with vertex set (0, 1, 2, 3) and facets {(0, 1, 2), (3,)}2755sage: V.is_connected()2756False27572758sage: X = SimplicialComplex([[0,1,2]])2759sage: X.is_connected()2760True27612762sage: U = simplicial_complexes.ChessboardComplex(3,3)2763sage: U.is_connected()2764True27652766sage: W = simplicial_complexes.Sphere(3)2767sage: W.is_connected()2768True27692770sage: S = SimplicialComplex([[0,1],[2,3]])2771sage: S.is_connected()2772False2773"""2774return self.graph().is_connected()27752776def n_skeleton(self, n):2777"""2778The `n`-skeleton of this simplicial complex.27792780The `n`-skeleton of a simplicial complex is obtained by discarding2781all of the simplices in dimensions larger than `n`.27822783:param n: non-negative integer27842785EXAMPLES::27862787sage: X = SimplicialComplex([[0,1], [1,2,3], [0,2,3]])2788sage: X.n_skeleton(1)2789Simplicial complex with vertex set (0, 1, 2, 3) and facets {(2, 3), (0, 2), (1, 3), (1, 2), (0, 3), (0, 1)}2790sage: X.set_immutable()2791sage: X.n_skeleton(2)2792Simplicial complex with vertex set (0, 1, 2, 3) and facets {(0, 2, 3), (1, 2, 3), (0, 1)}2793"""2794# make sure it's a list (it will be a tuple if immutable)2795facets = list(filter(lambda f: f.dimension()<n, self._facets))2796facets.extend(self.n_faces(n))2797return SimplicialComplex(facets, is_mutable=self._is_mutable)27982799def _contractible_subcomplex(self, verbose=False):2800"""2801Find a contractible subcomplex `L` of this simplicial complex,2802preferably one which is as large as possible.28032804:param verbose: If ``True``, print some messages as the simplicial2805complex is computed.2806:type verbose: boolean; optional, default ``False``28072808Motivation: if `K` is the original complex and if `L` is2809contractible, then the relative homology `H_*(K,L)` is2810isomorphic to the reduced homology of `K`. If `L` is large,2811then the relative chain complex will be a good deal smaller2812than the augmented chain complex for `K`, and this leads to a2813speed improvement for computing the homology of `K`.28142815This just passes an immutable subcomplex consisting of a facet to the2816method ``_enlarge_subcomplex``.28172818.. NOTE::28192820Thus when the simplicial complex is empty, so is the2821resulting 'contractible subcomplex', which is therefore not2822technically contractible. In this case, that doesn't2823matter because the homology is computed correctly anyway.28242825EXAMPLES::28262827sage: sphere = SimplicialComplex([[0,1,2,3]])2828sage: sphere.remove_face([0,1,2,3])2829sage: sphere2830Simplicial complex with vertex set (0, 1, 2, 3) and facets {(0, 2, 3), (0, 1, 2), (1, 2, 3), (0, 1, 3)}2831sage: L = sphere._contractible_subcomplex(); L2832Simplicial complex with vertex set (0, 1, 2, 3) and facets {(0, 2, 3), (1, 2, 3), (0, 1, 3)}2833sage: L.homology()2834{0: 0, 1: 0, 2: 0}2835"""2836facets = [self._facets[0]]2837return self._enlarge_subcomplex(SimplicialComplex(facets, is_mutable=False), verbose=verbose)28382839def _enlarge_subcomplex(self, subcomplex, verbose=False):2840"""2841Given a subcomplex `S` of this simplicial complex `K`, find a2842subcomplex `L`, as large as possible, containing `S` which is2843homotopy equivalent to `S` (so that `H_{*}(K,S)` is isomorphic2844to `H_{*}(K,L)`). This way, the chain complex for computing2845`H_{*}(K,L)` will be smaller than that for computing2846`H_{*}(K,S)`, so the computations should be faster.28472848:param subcomplex: a subcomplex of this simplicial complex2849:param verbose: If ``True``, print some messages as the simplicial2850complex is computed.2851:type verbose: boolean; optional, default ``False``2852:return: a complex `L` containing ``subcomplex`` and contained2853in ``self``, homotopy equivalent to ``subcomplex``.28542855Algorithm: start with the subcomplex `S` and loop through the2856facets of `K` which are not in `S`. For each one, see whether2857its intersection with `S` is contractible, and if so, add it.2858This is recursive: testing for contractibility calls this2859routine again, via ``_contractible_subcomplex``.28602861EXAMPLES::28622863sage: T = simplicial_complexes.Torus(); T2864Simplicial complex with vertex set (0, 1, 2, 3, 4, 5, 6) and 14 facets28652866Inside the torus, define a subcomplex consisting of a loop::28672868sage: S = SimplicialComplex([[0,1], [1,2], [0,2]], is_mutable=False)2869sage: S.homology()2870{0: 0, 1: Z}2871sage: L = T._enlarge_subcomplex(S)2872sage: L2873Simplicial complex with vertex set (0, 1, 2, 3, 4, 5, 6) and 8 facets2874sage: L.facets()2875{(0, 1, 5), (1, 3, 6), (1, 2), (1, 2, 4), (1, 3, 4), (0, 2), (1, 5, 6), (0, 1)}2876sage: L.homology()[1]2877Z2878"""2879# Make the subcomplex immutable if not2880if subcomplex is not None and subcomplex._is_mutable:2881subcomplex = SimplicialComplex(subcomplex._facets,2882maximality_check=False,2883sort_facets=False,2884is_mutable=False)28852886if subcomplex in self.__enlarged:2887return self.__enlarged[subcomplex]2888faces = filter(lambda x: x not in subcomplex._facets, list(self._facets))2889done = False2890new_facets = list(subcomplex._facets)2891while not done:2892done = True2893remove_these = []2894if verbose:2895print " looping through %s facets" % len(faces)2896for f in faces:2897f_set = f.set()2898int_facets = set( a.set().intersection(f_set) for a in new_facets )2899intersection = SimplicialComplex(int_facets)2900if not intersection._facets[0].is_empty():2901if (len(intersection._facets) == 1 or2902intersection == intersection._contractible_subcomplex()):2903new_facets.append(f)2904remove_these.append(f)2905done = False2906if verbose and not done:2907print " added %s facets" % len(remove_these)2908for f in remove_these:2909faces.remove(f)2910if verbose:2911print " now constructing a simplicial complex with %s vertices and %s facets" % (self.vertices().dimension()+1, len(new_facets))2912L = SimplicialComplex(new_facets, maximality_check=False,2913sort_facets=False, is_mutable=self._is_mutable)2914self.__enlarged[subcomplex] = L2915return L29162917def _cubical_(self):2918r"""2919Cubical complex constructed from ``self``.29202921ALGORITHM:29222923The algorithm comes from a paper by Shtan'ko and Shtogrin, as2924reported by Bukhshtaber and Panov. Let `I^m` denote the unit2925`m`-cube, viewed as a cubical complex. Let `[m] = \{1, 2,2926..., m\}`; then each face of `I^m` has the following form, for2927subsets `I \subset J \subset [m]`:29282929.. MATH::29302931F_{I \subset J} = \{ (y_1,...,y_m) \in I^m \,:\, y_i =0 \text{2932for } i \in I, y_j = 1 \text{ for } j \not \in J\}.29332934If `K` is a simplicial complex on vertex set `[m]` and if `I2935\subset [m]`, write `I \in K` if `I` is a simplex of `K`.2936Then we associate to `K` the cubical subcomplex of `I^m` with2937faces29382939.. MATH::29402941\{F_{I \subset J} \,:\, J \in K, I \neq \emptyset \}29422943The geometric realization of this cubical complex is2944homeomorphic to the geometric realization of the original2945simplicial complex.29462947REFERENCES:29482949.. [BP2000] V. M. Bukhshtaber and T. E. Panov, "Moment-angle complexes2950and combinatorics of simplicial manifolds," *Uspekhi2951Mat. Nauk* 55 (2000), 171--172.29522953.. [SS1992] M. A. Shtan'ko and and M. I. Shtogrin, "Embedding cubic2954manifolds and complexes into a cubic lattice", *Uspekhi2955Mat. Nauk* 47 (1992), 219-220.29562957EXAMPLES::29582959sage: T = simplicial_complexes.Torus()2960sage: T.homology()2961{0: 0, 1: Z x Z, 2: Z}2962sage: Tc = T._cubical_()2963sage: Tc2964Cubical complex with 42 vertices and 168 cubes2965sage: Tc.homology()2966{0: 0, 1: Z x Z, 2: Z}2967"""2968from sage.homology.cubical_complex import CubicalComplex2969V = self.vertices()2970embed = V.dimension() + 12971# dictionary to translate vertices to the numbers 1, ..., embed2972vd = dict(zip(V, range(1, embed + 1)))2973cubes = []2974for JJ in self.facets():2975J = [vd[i] for i in JJ]2976for i in J:2977# loop over indices from 1 to embed. if equal to i,2978# set to 0. if not in J, set to 1. Otherwise, range2979# from 0 to 12980cube = []2981for n in range(1, embed+1):2982if n == i:2983cube.append([0,])2984elif n not in J:2985cube.append([1,])2986else:2987cube.append([0,1])2988cubes.append(cube)2989return CubicalComplex(cubes)29902991def connected_component(self, simplex=None):2992"""2993Return the connected component of this simplicial complex2994containing ``simplex``. If ``simplex`` is omitted, then return2995the connected component containing the zeroth vertex in the2996vertex list. (If the simplicial complex is empty, raise an2997error.)29982999EXAMPLES::30003001sage: S1 = simplicial_complexes.Sphere(1)3002sage: S1 == S1.connected_component()3003True3004sage: X = S1.disjoint_union(S1)3005sage: X == X.connected_component()3006False3007sage: v0 = X.vertices()[0]3008sage: v1 = X.vertices()[-1]3009sage: X.connected_component(Simplex([v0])) == X.connected_component(Simplex([v1]))3010False30113012sage: S0 = simplicial_complexes.Sphere(0)3013sage: S0.vertices()3014(0, 1)3015sage: S0.connected_component()3016Simplicial complex with vertex set (0,) and facets {(0,)}3017sage: S0.connected_component(Simplex((1,)))3018Simplicial complex with vertex set (1,) and facets {(1,)}30193020sage: SimplicialComplex([[]]).connected_component()3021Traceback (most recent call last):3022...3023ValueError: the empty simplicial complex has no connected components.3024"""3025if self.dimension() == -1:3026raise ValueError("the empty simplicial complex has no connected components.")3027if simplex is None:3028v = self.vertices()[0]3029else:3030v = simplex[0]3031vertices = self.graph().connected_component_containing_vertex(v)3032facets = [f for f in self.facets() if f.is_face(Simplex(vertices))]3033return SimplicialComplex(facets)30343035def fundamental_group(self, base_point=None, simplify=True):3036r"""3037Return the fundamental group of this simplicial complex.30383039INPUT:30403041- ``base_point`` (optional, default None) -- if this complex is3042not path-connected, then specify a vertex; the fundamental3043group is computed with that vertex as a base point. If the3044complex is path-connected, then you may specify a vertex or3045leave this as its default setting of ``None``. (If this3046complex is path-connected, then this argument is ignored.)30473048- ``simplify`` (bool, optional True) -- if False, then return a3049presentation of the group in terms of generators and3050relations. If True, the default, simplify as much as GAP is3051able to.30523053Algorithm: we compute the edge-path group -- see3054:wikipedia:`Fundamental_group`. Choose a spanning tree for the30551-skeleton, and then the group's generators are given by the3056edges in the 1-skeleton; there are two types of relations:3057`e=1` if `e` is in the spanning tree, and for every 2-simplex,3058if its edges are `e_0`, `e_1`, and `e_2`, then we impose the3059relation `e_0 e_1^{-1} e_2 = 1`.30603061EXAMPLES::30623063sage: S1 = simplicial_complexes.Sphere(1)3064sage: S1.fundamental_group()3065Finitely presented group < e | >30663067If we pass the argument ``simplify=False``, we get generators and3068relations in a form which is not usually very helpful. Here is the3069cyclic group of order 2, for instance::30703071sage: RP2 = simplicial_complexes.RealProjectiveSpace(2)3072sage: C2 = RP2.fundamental_group(simplify=False)3073sage: C23074Finitely presented group < e0, e1, e2, e3, e4, e5, e6, e7, e8, e9 | e6, e5, e3, e9, e4*e7^-1*e6, e9*e7^-1*e0, e0*e1^-1*e2, e5*e1^-1*e8, e4*e3^-1*e8, e2 >3075sage: C2.simplified()3076Finitely presented group < e0 | e0^2 >30773078This is the same answer given if the argument ``simplify`` is True3079(the default)::30803081sage: RP2.fundamental_group()3082Finitely presented group < e0 | e0^2 >30833084You must specify a base point to compute the fundamental group3085of a non-connected complex::30863087sage: K = S1.disjoint_union(RP2)3088sage: K.fundamental_group()3089Traceback (most recent call last):3090...3091ValueError: this complex is not connected, so you must specify a base point.3092sage: v0 = list(K.vertices())[0]3093sage: K.fundamental_group(base_point=v0)3094Finitely presented group < e | >3095sage: v1 = list(K.vertices())[-1]3096sage: K.fundamental_group(base_point=v1)3097Finitely presented group < e0 | e0^2 >30983099Some other examples::31003101sage: S1.wedge(S1).fundamental_group()3102Finitely presented group < e0, e1 | >3103sage: simplicial_complexes.Torus().fundamental_group()3104Finitely presented group < e0, e3 | e0*e3^-1*e0^-1*e3 >3105sage: simplicial_complexes.MooreSpace(5).fundamental_group()3106Finitely presented group < e1 | e1^5 >3107"""3108if not self.is_connected():3109if base_point is None:3110raise ValueError("this complex is not connected, so you must specify a base point.")3111return self.connected_component(Simplex([base_point])).fundamental_group(simplify=simplify)31123113from sage.groups.free_group import FreeGroup3114from sage.interfaces.gap import gap3115spanning_tree = [e[:2] for e in self.graph().min_spanning_tree()]3116gens = [tuple(e) for e in self.n_cells(1) if tuple(e) not in spanning_tree]31173118if len(gens) == 0:3119return gap.TrivialGroup()31203121gens_dict = dict(zip(gens, range(len(gens))))3122FG = FreeGroup(len(gens), 'e')3123rels = []3124for f in self.n_cells(2):3125bdry = [tuple(e) for e in f.faces()]3126z = dict()3127for i in range(3):3128if bdry[i] in spanning_tree:3129z[i] = FG.one()3130else:3131z[i] = FG.gen(gens_dict[bdry[i]])3132rels.append(z[0]*z[1].inverse()*z[2])3133if simplify:3134return FG.quotient(rels).simplified()3135else:3136return FG.quotient(rels)31373138def category(self):3139"""3140Return the category to which this simplicial complex belongs: the3141category of all simplicial complexes.31423143EXAMPLES::31443145sage: SimplicialComplex([[0,1], [1,2,3,4,5]]).category()3146Category of simplicial complexes3147"""3148import sage.categories.all3149return sage.categories.all.SimplicialComplexes()31503151def is_isomorphic(self,other, certify = False):3152r"""3153Checks whether two simplicial complexes are isomorphic31543155INPUT:31563157- ``certify`` - if ``True``, then output is ``(a,b)``, where ``a``3158is a boolean and ``b`` is either a map or ``None``.31593160This is done by creating two graphs and checking whether they3161are isomorphic.31623163EXAMPLES::31643165sage: Z1 = SimplicialComplex([[0,1],[1,2],[2,3,4],[4,5]])3166sage: Z2 = SimplicialComplex([['a','b'],['b','c'],['c','d','e'],['e','f']])3167sage: Z3 = SimplicialComplex([[1,2,3]])3168sage: Z1.is_isomorphic(Z2)3169True3170sage: Z1.is_isomorphic(Z2, certify=True)3171(True, {0: 'a', 1: 'b', 2: 'c', 3: 'd', 4: 'e', 5: 'f'})3172sage: Z3.is_isomorphic(Z2)3173False3174"""3175g1 = Graph()3176g2 = Graph()3177g1.add_edges((v,f) for f in self.facets() for v in f)3178g2.add_edges((v,f) for f in other.facets() for v in f)3179g1.add_edges(("fake_vertex",v,"special_edge") for v in self.vertices())3180g2.add_edges(("fake_vertex",v,"special_edge") for v in other.vertices())3181if not certify:3182return g1.is_isomorphic(g2)3183isisom, tr = g1.is_isomorphic(g2, certify = True)31843185if isisom:3186for f in self.facets():3187tr.pop(f)3188tr.pop("fake_vertex")31893190return isisom,tr31913192def automorphism_group(self):3193r"""3194Returns the automorphism group of the simplicial complex31953196This is done by creating a bipartite graph, whose vertices are3197vertices and facets of the simplicial complex, and computing3198its automorphism group.31993200.. WARNING::32013202Since :trac:`14319` the domain of the automorphism group is equal to3203the graph's vertex set, and the ``translation`` argument has become3204useless.32053206EXAMPLES::32073208sage: S = simplicial_complexes.Simplex(3)3209sage: S.automorphism_group().is_isomorphic(SymmetricGroup(4))3210True32113212sage: P = simplicial_complexes.RealProjectivePlane()3213sage: P.automorphism_group().is_isomorphic(AlternatingGroup(5))3214True32153216sage: Z = SimplicialComplex([['1','2'],['2','3','a']])3217sage: Z.automorphism_group().is_isomorphic(CyclicPermutationGroup(2))3218True3219sage: group = Z.automorphism_group()3220sage: group.domain()3221{'1', '2', '3', 'a'}3222"""3223from sage.groups.perm_gps.permgroup import PermutationGroup32243225G = Graph()3226G.add_vertices(self.vertices())3227G.add_edges((f.tuple(),v) for f in self.facets() for v in f)3228groupe = G.automorphism_group(partition=[list(self.vertices()),3229[f.tuple() for f in self.facets()]])323032313232gens = [ [tuple(c) for c in g.cycle_tuples() if not isinstance(c[0],tuple)]3233for g in groupe.gens()]32343235permgroup = PermutationGroup(gens = gens, domain = self.vertices())32363237return permgroup32383239def _Hom_(self, other, category=None):3240"""3241Return the set of simplicial maps between simplicial complexes3242``self`` and ``other``.32433244EXAMPLES::32453246sage: S = simplicial_complexes.Sphere(1)3247sage: T = simplicial_complexes.Sphere(2)3248sage: H = Hom(S,T) # indirect doctest3249sage: H3250Set of Morphisms from Simplicial complex with vertex set (0, 1, 2) and facets {(1, 2), (0, 2), (0, 1)} to Simplicial complex with vertex set (0, 1, 2, 3) and facets {(0, 2, 3), (0, 1, 2), (1, 2, 3), (0, 1, 3)} in Category of simplicial complexes3251sage: f = {0:0,1:1,2:3}3252sage: x = H(f)3253sage: x3254Simplicial complex morphism {0: 0, 1: 1, 2: 3} from Simplicial complex with vertex set (0, 1, 2) and facets {(1, 2), (0, 2), (0, 1)} to Simplicial complex with vertex set (0, 1, 2, 3) and facets {(0, 2, 3), (0, 1, 2), (1, 2, 3), (0, 1, 3)}3255"""3256from sage.homology.simplicial_complex_homset import SimplicialComplexHomset3257return SimplicialComplexHomset(self, other)32583259# @cached_method when we switch to immutable SimplicialComplex3260def _is_numeric(self):3261"""3262Test whether all vertices are labeled by integers32633264OUTPUT:32653266Boolean. Whether all vertices are labeled by (not necessarily3267consecutive) integers.32683269EXAMPLES::32703271sage: s = SimplicialComplex()3272sage: s._is_numeric()3273True3274sage: s.add_face(['a', 'b', 123])3275sage: s._is_numeric()3276False3277"""3278return all([isinstance(v, (int, Integer, long)) for v in self._vertex_set])32793280# @cached_method when we switch to immutable SimplicialComplex3281def _translation_to_numeric(self):3282"""3283Return a dictionary enumerating the vertices32843285See also :meth:`_translation_from_numeric`, which returns the3286inverse map.32873288OUTPUT:32893290A dictionary. The keys are the vertices, and the associated3291values are integers from 0 to number of vertices - 1.32923293EXAMPLES::32943295sage: s = SimplicialComplex()3296sage: s._translation_to_numeric()3297{}3298sage: s.add_face(['a', 'b', 123])3299sage: s._translation_to_numeric() # random output3300{'a': 1, 123: 0, 'b': 2}3301sage: set(s._translation_to_numeric().keys()) == set(['a', 'b', 123])3302True3303sage: sorted(s._translation_to_numeric().values())3304[0, 1, 2]3305"""3306return dict((vertex, i) for i, vertex in enumerate(self._vertex_set))33073308# @cached_method when we switch to immutable SimplicialComplex3309def _translation_from_numeric(self):3310"""3311Return a dictionary mapping vertex indices to vertices33123313See also :meth:`_translation_to_numeric`, which returns the3314inverse map.33153316OUTPUT:33173318A dictionary. The keys are integers from 0 to the number of3319vertices - 1. The associated values are the vertices.33203321EXAMPLES::33223323sage: s = SimplicialComplex()3324sage: s._translation_from_numeric()3325{}3326sage: s.add_face(['a', 'b', 123])3327sage: s._translation_from_numeric() # random output3328{0: 123, 1: 'a', 2: 'b'}3329sage: sorted(s._translation_from_numeric().keys())3330[0, 1, 2]3331sage: set(s._translation_from_numeric().values()) == set(['a', 'b', 123])3332True3333"""3334return dict(enumerate(self._vertex_set))33353336def _chomp_repr_(self):3337r"""3338String representation of ``self`` suitable for use by the CHomP3339program. This lists each facet on its own line, and makes3340sure vertices are listed as numbers.33413342EXAMPLES::33433344sage: S = SimplicialComplex([(0,1,2), (2,3,5)])3345sage: print S._chomp_repr_()3346(2, 3, 5)3347(0, 1, 2)33483349A simplicial complex whose vertices are tuples, not integers::33503351sage: S = SimplicialComplex([[(0,1), (1,2), (3,4)]])3352sage: S._chomp_repr_()3353'(0, 1, 2)\n'3354"""3355s = ""3356numeric = self._is_numeric()3357if not numeric:3358d = self._translation_to_numeric()3359for f in self.facets():3360if numeric:3361s += str(f)3362else:3363s += '(' + ', '.join(str(d[a]) for a in f) + ')'3364s += '\n'3365return s33663367# this function overrides the standard one for GenericCellComplex,3368# because it lists the maximal faces, not the total number of faces.3369def _repr_(self):3370"""3371Print representation.33723373If there are only a few vertices or faces, they are listed. If3374there are lots, the number is given.33753376EXAMPLES::33773378sage: X = SimplicialComplex([[0,1], [1,2]])3379sage: X._repr_()3380'Simplicial complex with vertex set (0, 1, 2) and facets {(1, 2), (0, 1)}'3381sage: SimplicialComplex([[i for i in range(16)]])3382Simplicial complex with 16 vertices and 1 facets3383"""3384vertex_limit = 453385facet_limit = 553386vertices = self.vertices()3387facets = Set(self._facets)3388vertex_string = "with vertex set %s" % vertices3389if len(vertex_string) > vertex_limit:3390vertex_string = "with %s vertices" % str(1+vertices.dimension())3391facet_string = "facets %s" % facets3392if len(facet_string) > facet_limit:3393facet_string = "%s facets" % len(facets)3394return "Simplicial complex " + vertex_string + " and " + facet_string33953396def set_immutable(self):3397"""3398Make this simplicial complex immutable.33993400EXAMPLES::34013402sage: S = SimplicialComplex([[1,4], [2,4]])3403sage: S.is_mutable()3404True3405sage: S.set_immutable()3406sage: S.is_mutable()3407False3408"""3409self._is_mutable = False3410self._facets = tuple(self._facets)34113412def is_mutable(self):3413"""3414Return ``True`` if mutable.34153416EXAMPLES::34173418sage: S = SimplicialComplex([[1,4], [2,4]])3419sage: S.is_mutable()3420True3421sage: S.set_immutable()3422sage: S.is_mutable()3423False3424sage: S2 = SimplicialComplex([[1,4], [2,4]], is_mutable=False)3425sage: S2.is_mutable()3426False3427sage: S3 = SimplicialComplex([[1,4], [2,4]], is_mutable=False)3428sage: S3.is_mutable()3429False3430"""3431return self._is_mutable34323433def is_immutable(self):3434"""3435Return ``True`` if immutable.34363437EXAMPLES::34383439sage: S = SimplicialComplex([[1,4], [2,4]])3440sage: S.is_immutable()3441False3442sage: S.set_immutable()3443sage: S.is_immutable()3444True3445"""3446return not self._is_mutable3447344834493450