Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagesmc
Path: blob/master/src/sage/homology/chain_complex.py
8818 views
1
r"""
2
Chain complexes
3
4
AUTHORS:
5
6
- John H. Palmieri (2009-04)
7
8
This module implements bounded chain complexes of free `R`-modules,
9
for any commutative ring `R` (although the interesting things, like
10
homology, only work if `R` is the integers or a field).
11
12
Fix a ring `R`. A chain complex over `R` is a collection of
13
`R`-modules `\{C_n\}` indexed by the integers, with `R`-module maps
14
`d_n : C_n \rightarrow C_{n+1}` such that `d_{n+1} \circ d_n = 0` for
15
all `n`. The maps `d_n` are called *differentials*.
16
17
One can vary this somewhat: the differentials may decrease degree by
18
one instead of increasing it: sometimes a chain complex is defined
19
with `d_n : C_n \rightarrow C_{n-1}` for each `n`. Indeed, the
20
differentials may change dimension by any fixed integer.
21
22
Also, the modules may be indexed over an abelian group other than the
23
integers, e.g., `\ZZ^{m}` for some integer `m \geq 1`, in which case
24
the differentials may change the grading by any element of that
25
grading group. The elements of the grading group are generally called
26
degrees, so `C_n` is the module in degree `n` and so on.
27
28
In this implementation, the ring `R` must be commutative and the
29
modules `C_n` must be free `R`-modules. As noted above, homology
30
calculations will only work if the ring `R` is either `\ZZ` or a
31
field. The modules may be indexed by any free abelian group. The
32
differentials may increase degree by 1 or decrease it, or indeed
33
change it by any fixed amount: this is controlled by the
34
``degree_of_differential`` parameter used in defining the chain
35
complex.
36
"""
37
38
39
########################################################################
40
# Copyright (C) 2013 John H. Palmieri <[email protected]>
41
# Volker Braun <[email protected]>
42
#
43
# Distributed under the terms of the GNU General Public License (GPL)
44
# as published by the Free Software Foundation; either version 2 of
45
# the License, or (at your option) any later version.
46
#
47
# http://www.gnu.org/licenses/
48
########################################################################
49
50
from copy import copy
51
52
from sage.structure.parent import Parent
53
from sage.structure.element import ModuleElement
54
from sage.structure.element import is_Vector
55
from sage.misc.cachefunc import cached_method
56
57
from sage.rings.integer_ring import ZZ
58
from sage.rings.rational_field import QQ
59
from sage.modules.free_module import FreeModule
60
from sage.modules.free_module_element import vector
61
from sage.matrix.matrix0 import Matrix
62
from sage.matrix.constructor import matrix, prepare_dict
63
from sage.misc.latex import latex
64
from sage.rings.all import GF, prime_range
65
from sage.misc.decorators import rename_keyword
66
from sage.homology.homology_group import HomologyGroup
67
68
69
def _latex_module(R, m):
70
"""
71
LaTeX string representing a free module over ``R`` of rank ``m``.
72
73
INPUT:
74
75
- ``R`` -- a commutative ring
76
- ``m`` -- non-negative integer
77
78
This is used by the ``_latex_`` method for chain complexes.
79
80
EXAMPLES::
81
82
sage: from sage.homology.chain_complex import _latex_module
83
sage: _latex_module(ZZ, 3)
84
'\\Bold{Z}^{3}'
85
sage: _latex_module(ZZ, 0)
86
'0'
87
sage: _latex_module(GF(3), 1)
88
'\\Bold{F}_{3}^{1}'
89
"""
90
if m == 0:
91
return str(latex(0))
92
return str(latex(FreeModule(R, m)))
93
94
95
@rename_keyword(deprecation=15151, check_products='check', check_diffs='check')
96
def ChainComplex(data=None, **kwds):
97
r"""
98
Define a chain complex.
99
100
INPUT:
101
102
- ``data`` -- the data defining the chain complex; see below for
103
more details.
104
105
The following keyword arguments are supported:
106
107
- ``base_ring`` -- a commutative ring (optional), the ring over
108
which the chain complex is defined. If this is not specified,
109
it is determined by the data defining the chain complex.
110
111
- ``grading_group`` -- a additive free abelian group (optional,
112
default ``ZZ``), the group over which the chain complex is
113
indexed.
114
115
- ``degree_of_differential`` -- element of grading_group
116
(optional, default ``1``). The degree of the differential.
117
118
- ``degree`` -- alias for ``degree_of_differential``.
119
120
- ``check`` -- boolean (optional, default ``True``). If ``True``,
121
check that each consecutive pair of differentials are
122
composable and have composite equal to zero.
123
124
OUTPUT:
125
126
A chain complex.
127
128
.. WARNING::
129
130
Right now, homology calculations will only work if the base
131
ring is either `\ZZ` or a field, so please take this into account
132
when defining a chain complex.
133
134
Use data to define the chain complex. This may be in any of the
135
following forms.
136
137
1. a dictionary with integers (or more generally, elements of
138
grading_group) for keys, and with ``data[n]`` a matrix representing
139
(via left multiplication) the differential coming from degree
140
`n`. (Note that the shape of the matrix then determines the
141
rank of the free modules `C_n` and `C_{n+d}`.)
142
143
2. a list/tuple/iterable of the form `[C_0, d_0, C_1, d_1, C_2,
144
d_2, ...]`, where each `C_i` is a free module and each `d_i` is
145
a matrix, as above. This only makes sense if ``grading_group``
146
is `\ZZ` and ``degree`` is 1.
147
148
3. a list/tuple/iterable of the form `[r_0, d_0, r_1, d_1, r_2,
149
d_2, \ldots]`, where `r_i` is the rank of the free module `C_i`
150
and each `d_i` is a matrix, as above. This only makes sense if
151
``grading_group`` is `\ZZ` and ``degree`` is 1.
152
153
4. a list/tuple/iterable of the form `[d_0, d_1, d_2, \ldots]` where
154
each `d_i` is a matrix, as above. This only makes sense if
155
``grading_group`` is `\ZZ` and ``degree`` is 1.
156
157
.. NOTE::
158
159
In fact, the free modules `C_i` in case 2 and the ranks `r_i`
160
in case 3 are ignored: only the matrices are kept, and from
161
their shapes, the ranks of the modules are determined.
162
(Indeed, if ``data`` is a list or tuple, then any element which
163
is not a matrix is discarded; thus the list may have any number
164
of different things in it, and all of the non-matrices will be
165
ignored.) No error checking is done to make sure, for
166
instance, that the given modules have the appropriate ranks for
167
the given matrices. However, as long as ``check`` is True, the
168
code checks to see if the matrices are composable and that each
169
appropriate composite is zero.
170
171
If the base ring is not specified, then the matrices are examined
172
to determine a ring over which they are all naturally defined, and
173
this becomes the base ring for the complex. If no such ring can
174
be found, an error is raised. If the base ring is specified, then
175
the matrices are converted automatically to this ring when
176
defining the chain complex. If some matrix cannot be converted,
177
then an error is raised.
178
179
EXAMPLES::
180
181
sage: ChainComplex()
182
Trivial chain complex over Integer Ring
183
184
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])})
185
sage: C
186
Chain complex with at most 2 nonzero terms over Integer Ring
187
188
sage: m = matrix(ZZ, 2, 2, [0, 1, 0, 0])
189
sage: D = ChainComplex([m, m], base_ring=GF(2)); D
190
Chain complex with at most 3 nonzero terms over Finite Field of size 2
191
sage: D == loads(dumps(D))
192
True
193
sage: D.differential(0)==m, m.is_immutable(), D.differential(0).is_immutable()
194
(True, False, True)
195
196
Note that when a chain complex is defined in Sage, new
197
differentials may be created: every nonzero module in the chain
198
complex must have a differential coming from it, even if that
199
differential is zero::
200
201
sage: IZ = ChainComplex({0: identity_matrix(ZZ, 1)})
202
sage: IZ.differential() # the differentials in the chain complex
203
{0: [1], 1: [], -1: []}
204
sage: IZ.differential(1).parent()
205
Full MatrixSpace of 0 by 1 dense matrices over Integer Ring
206
sage: mat = ChainComplex({0: matrix(ZZ, 3, 4)}).differential(1)
207
sage: mat.nrows(), mat.ncols()
208
(0, 3)
209
210
Defining the base ring implicitly::
211
212
sage: ChainComplex([matrix(QQ, 3, 1), matrix(ZZ, 4, 3)])
213
Chain complex with at most 3 nonzero terms over Rational Field
214
sage: ChainComplex([matrix(GF(125, 'a'), 3, 1), matrix(ZZ, 4, 3)])
215
Chain complex with at most 3 nonzero terms over Finite Field in a of size 5^3
216
217
If the matrices are defined over incompatible rings, an error results::
218
219
sage: ChainComplex([matrix(GF(125, 'a'), 3, 1), matrix(QQ, 4, 3)])
220
Traceback (most recent call last):
221
...
222
TypeError: unable to find a common ring for all elements
223
224
If the base ring is given explicitly but is not compatible with
225
the matrices, an error results::
226
227
sage: ChainComplex([matrix(GF(125, 'a'), 3, 1)], base_ring=QQ)
228
Traceback (most recent call last):
229
...
230
TypeError: Unable to coerce 0 (<type
231
'sage.rings.finite_rings.element_givaro.FiniteField_givaroElement'>) to Rational
232
"""
233
check = kwds.get('check', True)
234
base_ring = kwds.get('base_ring', None)
235
grading_group = kwds.get('grading_group', ZZ)
236
degree = kwds.get('degree_of_differential', kwds.get('degree', 1))
237
try:
238
degree = grading_group(degree)
239
except Exception:
240
raise ValueError('degree is not an element of the grading group')
241
242
# transform data into data_dict
243
if data is None or (isinstance(data, (list, tuple)) and len(data) == 0):
244
# the zero chain complex
245
try:
246
zero = grading_group.identity()
247
except AttributeError:
248
zero = grading_group.zero_element()
249
if base_ring is None:
250
base_ring = ZZ
251
data_dict = dict()
252
elif isinstance(data, dict): # data is dictionary
253
data_dict = data
254
else: # data is list/tuple/iterable
255
data_matrices = filter(lambda x: isinstance(x, Matrix), data)
256
if degree != 1:
257
raise ValueError('degree must be +1 if the data argument is a list or tuple')
258
if grading_group != ZZ:
259
raise ValueError('grading_group must be ZZ if the data argument is a list or tuple')
260
data_dict = dict((grading_group(i), m) for i,m in enumerate(data_matrices))
261
262
if base_ring is None:
263
_, base_ring = prepare_dict(dict([n, data_dict[n].base_ring()(0)] for n in data_dict))
264
265
# make sure values in data_dict are appropriate matrices
266
for n in data_dict.keys():
267
if not n in grading_group:
268
raise ValueError('one of the dictionary keys is not an element of the grading group')
269
mat = data_dict[n]
270
if not isinstance(mat, Matrix):
271
raise TypeError('one of the differentials in the data is not a matrix')
272
if mat.base_ring() is base_ring:
273
if not mat.is_immutable():
274
mat = copy(mat) # do not make any arguments passed immutable
275
mat.set_immutable()
276
else:
277
mat = mat.change_ring(base_ring)
278
mat.set_immutable()
279
data_dict[n] = mat
280
281
# include any "obvious" zero matrices that are not 0x0
282
for n in data_dict.keys(): # note: data_dict will be mutated in this loop
283
mat1 = data_dict[n]
284
if (mat1.nrows(), mat1.ncols()) == (0, 0):
285
del data_dict[n]
286
if (mat1.nrows() != 0) and (n+degree not in data_dict):
287
if n+2*degree in data_dict:
288
mat2 = matrix(base_ring, data_dict[n+2*degree].ncols(), mat1.nrows())
289
else:
290
mat2 = matrix(base_ring, 0, mat1.nrows())
291
mat2.set_immutable()
292
data_dict[n+degree] = mat2
293
if (mat1.ncols() != 0) and (n-degree not in data_dict):
294
if n-2*degree in data_dict:
295
mat0 = matrix(base_ring, mat1.ncols(), data_dict[n-2*degree].nrows())
296
else:
297
mat0 = matrix(base_ring, mat1.ncols(), 0)
298
mat0.set_immutable()
299
data_dict[n-degree] = mat0
300
301
# check that this is a complex: going twice is zero
302
if check:
303
for n in data_dict.keys():
304
mat0 = data_dict[n]
305
try:
306
mat1 = data_dict[n+degree]
307
except KeyError:
308
continue
309
try:
310
prod = mat1 * mat0
311
except TypeError:
312
raise TypeError('the differentials d_{{{}}} and d_{{{}}} are not compatible: '
313
'their product is not defined'.format(n, n+degree))
314
if not prod.is_zero():
315
raise ValueError('the differentials d_{{{}}} and d_{{{}}} are not compatible: '
316
'their composition is not zero.'.format(n, n+degree))
317
318
return ChainComplex_class(grading_group, degree, base_ring, data_dict)
319
320
321
class Chain_class(ModuleElement):
322
323
def __init__(self, parent, vectors, check=True):
324
r"""
325
A Chain in a Chain Complex
326
327
A chain is collection of module elements for each module `C_n`
328
of the chain complex `(C_n, d_n)`. There is no restriction on
329
how the differentials `d_n` act on the elements of the chain.
330
331
.. NOTE::
332
333
You must use the chain complex to construct chains.
334
335
EXAMPLES::
336
337
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])}, base_ring=GF(7))
338
sage: C.category()
339
Category of chain complexes over Finite Field of size 7
340
341
TESTS::
342
343
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])})
344
sage: c = C({0:vector([0, 1, 2]), 1:vector([3, 4])})
345
sage: TestSuite(c).run()
346
"""
347
# only nonzero vectors shall be stored, ensuring this is the
348
# job of the _element constructor_
349
assert all(v.is_immutable() and not v.is_zero()
350
and v.base_ring() is parent.base_ring()
351
for v in vectors.values())
352
self._vec = vectors
353
super(Chain_class, self).__init__(parent)
354
355
def vector(self, degree):
356
"""
357
Return the free module element in ``degree``.
358
359
EXAMPLES::
360
361
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])})
362
sage: c = C({0:vector([1, 2, 3]), 1:vector([4, 5])})
363
sage: c.vector(0)
364
(1, 2, 3)
365
sage: c.vector(1)
366
(4, 5)
367
sage: c.vector(2)
368
()
369
"""
370
try:
371
return self._vec[degree]
372
except KeyError:
373
return self.parent().free_module(degree).zero()
374
375
def _repr_(self):
376
"""
377
Print representation.
378
379
EXAMPLES::
380
381
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])})
382
sage: C()
383
Trivial chain
384
sage: C({0:vector([1, 2, 3])})
385
Chain(0:(1, 2, 3))
386
sage: c = C({0:vector([1, 2, 3]), 1:vector([4, 5])}); c
387
Chain with 2 nonzero terms over Integer Ring
388
sage: c._repr_()
389
'Chain with 2 nonzero terms over Integer Ring'
390
"""
391
n = len(self._vec)
392
if n == 0:
393
return 'Trivial chain'
394
395
if n == 1:
396
deg, vec = self._vec.iteritems().next()
397
return 'Chain({0}:{1})'.format(deg, vec)
398
399
return 'Chain with {0} nonzero terms over {1}'.format(
400
n, self.parent().base_ring())
401
402
def _ascii_art_(self):
403
"""
404
Return an ascii art representation.
405
406
Note that arrows go to the left so that composition of
407
differentials is the usual matrix multiplication.
408
409
EXAMPLES::
410
411
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0]), 1:zero_matrix(1,2)})
412
sage: c = C({0:vector([1, 2, 3]), 1:vector([4, 5])})
413
sage: ascii_art(c)
414
d_2 d_1 d_0 [1] d_-1
415
0 <---- [0] <---- [4] <---- [2] <----- 0
416
[5] [3]
417
"""
418
from sage.misc.ascii_art import AsciiArt
419
420
def arrow_art(d):
421
d_str = [' d_{0} '.format(d)]
422
arrow = ' <' + '-'*(len(d_str[0])-3) + ' '
423
d_str.append(arrow)
424
return AsciiArt(d_str, baseline=0)
425
426
def vector_art(d):
427
v = self.vector(d)
428
if v.degree() == 0:
429
return AsciiArt(['0'])
430
v = str(v.column()).splitlines()
431
return AsciiArt(v, baseline=len(v)/2)
432
433
result = []
434
chain_complex = self.parent()
435
for ordered in chain_complex.ordered_degrees():
436
ordered = list(reversed(ordered))
437
if len(ordered) == 0:
438
return AsciiArt(['0'])
439
result_ordered = vector_art(ordered[0] + chain_complex.degree_of_differential())
440
for n in ordered:
441
result_ordered += arrow_art(n) + vector_art(n)
442
result = [result_ordered] + result
443
concatenated = result[0]
444
for r in result[1:]:
445
concatenated += AsciiArt([' ... ']) + r
446
return concatenated
447
448
def is_cycle(self):
449
"""
450
Return whether the chain is a cycle.
451
452
OUTPUT:
453
454
Boolean. Whether the elements of the chain are in the kernel
455
of the differentials.
456
457
EXAMPLES::
458
459
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])})
460
sage: c = C({0:vector([0, 1, 2]), 1:vector([3, 4])})
461
sage: c.is_cycle()
462
True
463
"""
464
chain_complex = self.parent()
465
for d, v in self._vec.iteritems():
466
dv = chain_complex.differential(d) * v
467
if not dv.is_zero():
468
return False
469
return True
470
471
def is_boundary(self):
472
"""
473
Return whether the chain is a boundary.
474
475
OUTPUT:
476
477
Boolean. Whether the elements of the chain are in the image of
478
the differentials.
479
480
EXAMPLES::
481
482
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])})
483
sage: c = C({0:vector([0, 1, 2]), 1:vector([3, 4])})
484
sage: c.is_boundary()
485
False
486
sage: z3 = C({1:(1, 0)})
487
sage: z3.is_cycle()
488
True
489
sage: (2*z3).is_boundary()
490
False
491
sage: (3*z3).is_boundary()
492
True
493
"""
494
chain_complex = self.parent()
495
for d, v in self._vec.iteritems():
496
d = chain_complex.differential(d - chain_complex.degree_of_differential()).transpose()
497
if v not in d.image():
498
return False
499
return True
500
501
def _add_(self, other):
502
"""
503
Module addition
504
505
EXAMPLES::
506
507
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])})
508
sage: c = C({0:vector([0, 1, 2]), 1:vector([3, 4])})
509
sage: c + c
510
Chain with 2 nonzero terms over Integer Ring
511
sage: ascii_art(c + c)
512
d_1 d_0 [0] d_-1
513
0 <---- [6] <---- [2] <----- 0
514
[8] [4]
515
"""
516
vectors = dict()
517
for d in set(self._vec.keys() + other._vec.keys()):
518
v = self.vector(d) + other.vector(d)
519
if not v.is_zero():
520
v.set_immutable()
521
vectors[d] = v
522
parent = self.parent()
523
return parent.element_class(parent, vectors)
524
525
def _rmul_(self, scalar):
526
"""
527
Scalar multiplication
528
529
EXAMPLES::
530
531
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])})
532
sage: c = C({0:vector([0, 1, 2]), 1:vector([3, 4])})
533
sage: 2 * c
534
Chain with 2 nonzero terms over Integer Ring
535
sage: 2 * c == c + c == c * 2
536
True
537
"""
538
vectors = dict()
539
for d, v in self._vec.iteritems():
540
v = scalar * v
541
if not v.is_zero():
542
v.set_immutable()
543
vectors[d] = v
544
parent = self.parent()
545
return parent.element_class(parent, vectors)
546
547
def __cmp__(self, other):
548
"""
549
Compare two chains
550
551
EXAMPLES::
552
553
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])})
554
sage: c = C({0:vector([0, 1, 2]), 1:vector([3, 4])})
555
sage: c == c
556
True
557
sage: c == C(0)
558
False
559
"""
560
c = cmp(type(self), type(other))
561
if c != 0:
562
return c
563
c = cmp(self.parent(), other.parent())
564
if c != 0:
565
return c
566
return cmp(self._vec, other._vec)
567
568
569
class ChainComplex_class(Parent):
570
r"""
571
See :func:`ChainComplex` for full documentation.
572
573
The differentials are required to be in the following canonical form:
574
575
* All differentials that are not `0 \times 0` must be specified
576
(even if they have zero rows or zero columns), and
577
578
* Differentials that are `0 \times 0` must not be specified.
579
580
* Immutable matrices over the ``base_ring``
581
582
This and more is ensured by the assertions in the
583
constructor. The :func:`ChainComplex` factory function must
584
ensure that only valid input is passed.
585
586
EXAMPLES::
587
588
sage: C = ChainComplex(); C
589
Trivial chain complex over Integer Ring
590
591
sage: D = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])})
592
sage: D
593
Chain complex with at most 2 nonzero terms over Integer Ring
594
"""
595
def __init__(self, grading_group, degree_of_differential, base_ring, differentials):
596
"""
597
Initialize ``self``.
598
599
TESTS::
600
601
sage: ChainComplex().base_ring()
602
Integer Ring
603
604
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])})
605
sage: TestSuite(C).run()
606
"""
607
if any(d.base_ring() != base_ring or not d.is_immutable() or
608
(d.ncols(), d.nrows()) == (0, 0)
609
for d in differentials.values()):
610
raise ValueError('invalid differentials')
611
if degree_of_differential.parent() is not grading_group:
612
raise ValueError('the degree_of_differential.parent() must be grading_group')
613
if grading_group is not ZZ and grading_group.is_multiplicative():
614
raise ValueError('grading_group must be either ZZ or multiplicative')
615
# all differentials (excluding the 0x0 ones) must be specified to the constructor
616
if any(dim+degree_of_differential not in differentials and d.nrows() != 0
617
for dim, d in differentials.iteritems()):
618
raise ValueError('invalid differentials')
619
if any(dim-degree_of_differential not in differentials and d.ncols() != 0
620
for dim, d in differentials.iteritems()):
621
raise ValueError('invalid differentials')
622
self._grading_group = grading_group
623
self._degree_of_differential = degree_of_differential
624
self._diff = differentials
625
626
from sage.categories.all import ChainComplexes
627
category = ChainComplexes(base_ring)
628
super(ChainComplex_class, self).__init__(base=base_ring, category=category)
629
630
Element = Chain_class
631
632
def _element_constructor_(self, vectors, check=True):
633
"""
634
The element constructor.
635
636
This is part of the Parent/Element framework. Calling the
637
parent uses this method to construct elements.
638
639
TESTS::
640
641
sage: D = ChainComplex({0: matrix(ZZ, 2, 2, [1,0,0,2])})
642
sage: D._element_constructor_(0)
643
Trivial chain
644
sage: D({0:[2, 3]})
645
Chain(0:(2, 3))
646
"""
647
if not vectors: # special case: the zero chain
648
return self.element_class(self, {})
649
if isinstance(vectors, Chain_class):
650
vectors = vectors._vec
651
data = dict()
652
for degree, vec in vectors.iteritems():
653
if not is_Vector(vec):
654
vec = vector(self.base_ring(), vec)
655
vec.set_immutable()
656
if check and vec.degree() != self.free_module_rank(degree):
657
raise ValueError('vector dimension does not match module dimension')
658
if vec.is_zero():
659
continue
660
if vec.base_ring() != self.base_ring():
661
vec = vec.change_ring(self.base_ring())
662
if not vec.is_immutable():
663
vec = copy(vec)
664
vec.set_immutable()
665
data[degree] = vec
666
return self.element_class(self, data)
667
668
def random_element(self):
669
"""
670
Return a random element.
671
672
EXAMPLES::
673
674
sage: D = ChainComplex({0: matrix(ZZ, 2, 2, [1,0,0,2])})
675
sage: D.random_element() # random output
676
Chain with 1 nonzero terms over Integer Ring
677
"""
678
vec = dict()
679
for d in self.nonzero_degrees():
680
vec[d] = self.free_module(d).random_element()
681
return self(vec)
682
683
_an_element_ = random_element
684
685
@cached_method
686
def rank(self, degree, ring=None):
687
r"""
688
Return the rank of a differential
689
690
INPUT:
691
692
- ``degree`` -- an element `\delta` of the grading
693
group. Which differential `d_{\delta}` we want to know the
694
rank of
695
696
- ``ring`` -- (optional) a commutative ring `S`;
697
if specified, the rank is computed after changing to this ring
698
699
OUTPUT:
700
701
The rank of the differential `d_{\delta} \otimes_R S`, where
702
`R` is the base ring of the chain complex.
703
704
EXAMPLES::
705
706
sage: C = ChainComplex({0:matrix(ZZ, [[2]])})
707
sage: C.differential(0)
708
[2]
709
sage: C.rank(0)
710
1
711
sage: C.rank(0, ring=GF(2))
712
0
713
"""
714
degree = self.grading_group()(degree)
715
try:
716
d = self._diff[degree]
717
except IndexError:
718
return ZZ.zero()
719
if d.nrows() == 0 or d.ncols() == 0:
720
return ZZ.zero()
721
if ring is None:
722
return d.rank()
723
return d.change_ring(ring).rank()
724
725
def grading_group(self):
726
r"""
727
Return the grading group.
728
729
OUTPUT:
730
731
The discrete abelian group that indexes the individual modules
732
of the complex. Usually `\ZZ`.
733
734
EXAMPLES::
735
736
sage: G = AdditiveAbelianGroup([0, 3])
737
sage: C = ChainComplex(grading_group=G, degree=G([1,2]))
738
sage: C.grading_group()
739
Additive abelian group isomorphic to Z/3 + Z
740
sage: C.degree_of_differential()
741
(2, 1)
742
"""
743
return self._grading_group
744
745
@cached_method
746
def nonzero_degrees(self):
747
r"""
748
Return the degrees in which the module is non-trivial.
749
750
See also :meth:`ordered_degrees`.
751
752
OUTPUT:
753
754
The tuple containing all degrees `n` (grading group elements)
755
such that the module `C_n` of the chain is non-trivial.
756
757
EXAMPLES::
758
759
sage: one = matrix(ZZ, [[1]])
760
sage: D = ChainComplex({0: one, 2: one, 6:one})
761
sage: ascii_art(D)
762
[1] [1] [0] [1]
763
0 <-- C_7 <---- C_6 <-- 0 ... 0 <-- C_3 <---- C_2 <---- C_1 <---- C_0 <-- 0
764
sage: D.nonzero_degrees()
765
(0, 1, 2, 3, 6, 7)
766
"""
767
return tuple(sorted(n for n,d in self._diff.iteritems() if d.ncols() > 0))
768
769
@cached_method
770
def ordered_degrees(self, start=None, exclude_first=False):
771
r"""
772
Sort the degrees in the order determined by the differential
773
774
INPUT:
775
776
- ``start`` -- (default: ``None``) a degree (element of the grading
777
group) or ``None``
778
779
- ``exclude_first`` -- boolean (optional; default:
780
``False``); whether to exclude the lowest degree -- this is a
781
handy way to just get the degrees of the non-zero modules,
782
as the domain of the first differential is zero.
783
784
OUTPUT:
785
786
If ``start`` has been specified, the longest tuple of degrees
787
788
* containing ``start`` (unless ``start`` would be the first
789
and ``exclude_first=True``),
790
791
* in ascending order relative to :meth:`degree_of_differential`, and
792
793
* such that none of the corresponding differentials are `0\times 0`.
794
795
If ``start`` has not been specified, a tuple of such tuples of
796
degrees. One for each sequence of non-zero differentials. They
797
are returned in sort order.
798
799
EXAMPLES::
800
801
sage: one = matrix(ZZ, [[1]])
802
sage: D = ChainComplex({0: one, 2: one, 6:one})
803
sage: ascii_art(D)
804
[1] [1] [0] [1]
805
0 <-- C_7 <---- C_6 <-- 0 ... 0 <-- C_3 <---- C_2 <---- C_1 <---- C_0 <-- 0
806
sage: D.ordered_degrees()
807
((-1, 0, 1, 2, 3), (5, 6, 7))
808
sage: D.ordered_degrees(exclude_first=True)
809
((0, 1, 2, 3), (6, 7))
810
sage: D.ordered_degrees(6)
811
(5, 6, 7)
812
sage: D.ordered_degrees(5, exclude_first=True)
813
(6, 7)
814
"""
815
if start is None:
816
result = []
817
degrees = set(self._diff.keys())
818
while len(degrees) > 0:
819
ordered = self.ordered_degrees(degrees.pop())
820
degrees.difference_update(ordered)
821
if exclude_first:
822
ordered = tuple(ordered[1:])
823
result.append(ordered)
824
result.sort()
825
return tuple(result)
826
827
import collections
828
deg = start
829
result = collections.deque()
830
result.append(start)
831
832
next_deg = start + self.degree_of_differential()
833
while next_deg in self._diff:
834
result.append(next_deg)
835
next_deg += self.degree_of_differential()
836
837
prev_deg = start - self.degree_of_differential()
838
while prev_deg in self._diff:
839
result.appendleft(prev_deg)
840
prev_deg -= self.degree_of_differential()
841
842
if exclude_first:
843
result.popleft()
844
return tuple(result)
845
846
def degree_of_differential(self):
847
"""
848
Return the degree of the differentials of the complex
849
850
OUTPUT:
851
852
An element of the grading group.
853
854
EXAMPLES::
855
856
sage: D = ChainComplex({0: matrix(ZZ, 2, 2, [1,0,0,2])})
857
sage: D.degree_of_differential()
858
1
859
"""
860
return self._degree_of_differential
861
862
def differential(self, dim=None):
863
"""
864
The differentials which make up the chain complex.
865
866
INPUT:
867
868
- ``dim`` -- element of the grading group (optional, default
869
``None``); if this is ``None``, return a dictionary of all
870
of the differentials, or if this is a single element, return
871
the differential starting in that dimension
872
873
OUTPUT:
874
875
Either a dictionary of all of the differentials or a single
876
differential (i.e., a matrix).
877
878
EXAMPLES::
879
880
sage: D = ChainComplex({0: matrix(ZZ, 2, 2, [1,0,0,2])})
881
sage: D.differential()
882
{0: [1 0]
883
[0 2], 1: [], -1: []}
884
sage: D.differential(0)
885
[1 0]
886
[0 2]
887
sage: C = ChainComplex({0: identity_matrix(ZZ, 40)})
888
sage: C.differential()
889
{0: 40 x 40 dense matrix over Integer Ring, 1: [],
890
-1: 40 x 0 dense matrix over Integer Ring}
891
"""
892
if dim is None:
893
return copy(self._diff)
894
dim = self.grading_group()(dim)
895
try:
896
return self._diff[dim]
897
except KeyError:
898
pass
899
# all differentials that are not 0x0 are in self._diff
900
return matrix(self.base_ring(), 0, 0)
901
902
def dual(self):
903
"""
904
The dual chain complex to ``self``.
905
906
Since all modules in ``self`` are free of finite rank, the
907
dual in dimension `n` is isomorphic to the original chain
908
complex in dimension `n`, and the corresponding boundary
909
matrix is the transpose of the matrix in the original complex.
910
This converts a chain complex to a cochain complex and vice versa.
911
912
EXAMPLES::
913
914
sage: C = ChainComplex({2: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])})
915
sage: C.degree_of_differential()
916
1
917
sage: C.differential(2)
918
[3 0 0]
919
[0 0 0]
920
sage: C.dual().degree_of_differential()
921
-1
922
sage: C.dual().differential(3)
923
[3 0]
924
[0 0]
925
[0 0]
926
"""
927
data = {}
928
deg = self.degree_of_differential()
929
for d in self.differential():
930
data[(d+deg)] = self.differential()[d].transpose()
931
return ChainComplex(data, degree=-deg)
932
933
def free_module_rank(self, degree):
934
r"""
935
Return the rank of the free module at the given ``degree``.
936
937
INPUT:
938
939
- ``degree`` -- an element of the grading group
940
941
OUTPUT:
942
943
Integer. The rank of the free module `C_n` at the given degree
944
`n`.
945
946
EXAMPLES::
947
948
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0]), 1: matrix(ZZ, [[0, 1]])})
949
sage: [C.free_module_rank(i) for i in range(-2, 5)]
950
[0, 0, 3, 2, 1, 0, 0]
951
"""
952
try:
953
return self._diff[degree].ncols()
954
except KeyError:
955
return ZZ.zero()
956
957
def free_module(self, degree=None):
958
r"""
959
Return the free module at fixed ``degree``, or their sum.
960
961
INPUT:
962
963
- ``degree`` -- an element of the grading group or ``None`` (default).
964
965
OUTPUT:
966
967
The free module `C_n` at the given degree `n`. If the degree
968
is not specified, the sum `\bigoplus C_n` is returned.
969
970
EXAMPLES::
971
972
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0]), 1: matrix(ZZ, [[0, 1]])})
973
sage: C.free_module()
974
Ambient free module of rank 6 over the principal ideal domain Integer Ring
975
sage: C.free_module(0)
976
Ambient free module of rank 3 over the principal ideal domain Integer Ring
977
sage: C.free_module(1)
978
Ambient free module of rank 2 over the principal ideal domain Integer Ring
979
sage: C.free_module(2)
980
Ambient free module of rank 1 over the principal ideal domain Integer Ring
981
"""
982
if degree is None:
983
rank = sum([mat.ncols() for mat in self.differential().values()])
984
else:
985
rank = self.free_module_rank(degree)
986
return FreeModule(self.base_ring(), rank)
987
988
def __cmp__(self, other):
989
"""
990
Return ``True`` iff this chain complex is the same as other: that
991
is, if the base rings and the matrices of the two are the
992
same.
993
994
EXAMPLES::
995
996
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])}, base_ring=GF(2))
997
sage: D = ChainComplex({0: matrix(GF(2), 2, 3, [1, 0, 0, 0, 0, 0]), 1: matrix(ZZ, 0, 2), 3: matrix(ZZ, 0, 0)}) # base_ring determined from the matrices
998
sage: C == D
999
True
1000
"""
1001
if not isinstance(other, ChainComplex_class):
1002
return cmp(type(other), ChainComplex_class)
1003
c = cmp(self.base_ring(), other.base_ring())
1004
if c != 0:
1005
return c
1006
R = self.base_ring()
1007
equal = True
1008
for d,mat in self.differential().iteritems():
1009
if d not in other.differential():
1010
equal = equal and mat.ncols() == 0 and mat.nrows() == 0
1011
else:
1012
equal = (equal and
1013
other.differential()[d].change_ring(R) == mat.change_ring(R))
1014
for d,mat in other.differential().iteritems():
1015
if d not in self.differential():
1016
equal = equal and mat.ncols() == 0 and mat.nrows() == 0
1017
if equal:
1018
return 0
1019
return -1
1020
1021
def _homology_chomp(self, deg, base_ring, verbose, generators):
1022
"""
1023
Helper function for :meth:`homology`.
1024
1025
EXAMPLES::
1026
1027
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])}, base_ring=GF(2))
1028
sage: C._homology_chomp(None, GF(2), False, False) # optional - CHomP
1029
{0: Vector space of dimension 2 over Finite Field of size 2, 1: Vector space of dimension 1 over Finite Field of size 2}
1030
"""
1031
from sage.interfaces.chomp import homchain
1032
H = homchain(self, base_ring=base_ring, verbose=verbose, generators=generators)
1033
if H is None:
1034
raise RuntimeError('ran CHomP, but no output')
1035
if deg is None:
1036
return H
1037
try:
1038
return H[deg]
1039
except KeyError:
1040
return HomologyGroup(0, base_ring)
1041
1042
1043
@rename_keyword(deprecation=15151, dim='deg')
1044
def homology(self, deg=None, **kwds):
1045
r"""
1046
The homology of the chain complex.
1047
1048
INPUT:
1049
1050
- ``deg`` -- an element of the grading group for the chain
1051
complex (default: ``None``); the degree in which
1052
to compute homology -- if this is ``None``, return the
1053
homology in every degree in which the chain complex is
1054
possibly nonzero
1055
1056
- ``base_ring`` -- a commutative ring (optional, default is the
1057
base ring for the chain complex); must be either the
1058
integers `\ZZ` or a field
1059
1060
- ``generators`` -- boolean (optional, default ``False``); if
1061
``True``, return generators for the homology groups along with
1062
the groups. See :trac:`6100`
1063
1064
- ``verbose`` - boolean (optional, default ``False``); if
1065
``True``, print some messages as the homology is computed
1066
1067
- ``algorithm`` - string (optional, default ``'auto'``); the
1068
options are:
1069
1070
* ``'auto'``
1071
* ``'chomp'``
1072
* ``'dhsw'``
1073
* ``'pari'``
1074
* ``'no_chomp'``
1075
1076
see below for descriptions
1077
1078
OUTPUT:
1079
1080
If the degree is specified, the homology in degree ``deg``.
1081
Otherwise, the homology in every dimension as a dictionary
1082
indexed by dimension.
1083
1084
ALGORITHM:
1085
1086
If ``algorithm`` is set to ``'auto'`` (the default), then use
1087
CHomP if available. CHomP is available at the web page
1088
http://chomp.rutgers.edu/. It is also an experimental package
1089
for Sage. If ``algorithm`` is ``chomp``, always use chomp.
1090
1091
CHomP computes homology, not cohomology, and only works over
1092
the integers or finite prime fields. Therefore if any of
1093
these conditions fails, or if CHomP is not present, or if
1094
``algorithm`` is set to 'no_chomp', go to plan B: if ``self``
1095
has a ``_homology`` method -- each simplicial complex has
1096
this, for example -- then call that. Such a method implements
1097
specialized algorithms for the particular type of cell
1098
complex.
1099
1100
Otherwise, move on to plan C: compute the chain complex of
1101
``self`` and compute its homology groups. To do this: over a
1102
field, just compute ranks and nullities, thus obtaining
1103
dimensions of the homology groups as vector spaces. Over the
1104
integers, compute Smith normal form of the boundary matrices
1105
defining the chain complex according to the value of
1106
``algorithm``. If ``algorithm`` is ``'auto'`` or ``'no_chomp'``,
1107
then for each relatively small matrix, use the standard Sage
1108
method, which calls the Pari package. For any large matrix,
1109
reduce it using the Dumas, Heckenbach, Saunders, and Welker
1110
elimination algorithm [DHSW]_: see
1111
:func:`~sage.homology.matrix_utils.dhsw_snf` for details.
1112
1113
Finally, ``algorithm`` may also be ``'pari'`` or ``'dhsw'``, which
1114
forces the named algorithm to be used regardless of the size
1115
of the matrices and regardless of whether CHomP is available.
1116
1117
As of this writing, CHomP is by far the fastest option,
1118
followed by the ``'auto'`` or ``'no_chomp'`` setting of using the
1119
Dumas, Heckenbach, Saunders, and Welker elimination algorithm
1120
[DHSW]_ for large matrices and Pari for small ones.
1121
1122
.. WARNING::
1123
1124
This only works if the base ring is the integers or a
1125
field. Other values will return an error.
1126
1127
EXAMPLES::
1128
1129
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])})
1130
sage: C.homology()
1131
{0: Z x Z, 1: Z x C3}
1132
sage: C.homology(deg=1, base_ring = GF(3))
1133
Vector space of dimension 2 over Finite Field of size 3
1134
sage: D = ChainComplex({0: identity_matrix(ZZ, 4), 4: identity_matrix(ZZ, 30)})
1135
sage: D.homology()
1136
{0: 0, 1: 0, 4: 0, 5: 0}
1137
1138
Generators: generators are given as
1139
a list of cycles, each of which is an element in the
1140
appropriate free module, and hence is represented as a vector::
1141
1142
sage: C.homology(1, generators=True) # optional - CHomP
1143
(Z x C3, [(0, 1), (1, 0)])
1144
1145
Tests for :trac:`6100`, the Klein bottle with generators::
1146
1147
sage: d0 = matrix(ZZ, 0,1)
1148
sage: d1 = matrix(ZZ, 1,3, [[0,0,0]])
1149
sage: d2 = matrix(ZZ, 3,2, [[1,1], [1,-1], [-1,1]])
1150
sage: C_k = ChainComplex({0:d0, 1:d1, 2:d2}, degree=-1)
1151
sage: C_k.homology(generators=true) # optional - CHomP
1152
{0: (Z, [(1)]), 1: (Z x C2, [(0, 0, 1), (0, 1, -1)])}
1153
1154
From a torus using a field::
1155
1156
sage: T = simplicial_complexes.Torus()
1157
sage: C_t = T.chain_complex()
1158
sage: C_t.homology(base_ring=QQ, generators=True)
1159
{0: [(Vector space of dimension 1 over Rational Field, Chain(0:(0, 0, 0, 0, 0, 0, 1)))],
1160
1: [(Vector space of dimension 1 over Rational Field,
1161
Chain(1:(0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, -1, -1, 0, 0, 0, -1, 0, 1, 0))),
1162
(Vector space of dimension 1 over Rational Field,
1163
Chain(1:(0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 1, 0, -1, -1)))],
1164
2: [(Vector space of dimension 1 over Rational Field,
1165
Chain(2:(1, -1, -1, -1, 1, -1, -1, 1, 1, 1, 1, 1, -1, -1)))]}
1166
"""
1167
from sage.interfaces.chomp import have_chomp, homchain
1168
1169
if deg is not None and deg not in self.grading_group():
1170
raise ValueError('degree is not an element of the grading group')
1171
1172
verbose = kwds.get('verbose', False)
1173
generators = kwds.get('generators', False)
1174
base_ring = kwds.get('base_ring', self.base_ring())
1175
if not (base_ring.is_field() or base_ring is ZZ):
1176
raise NotImplementedError('can only compute homology if the base ring is the integers or a field')
1177
1178
algorithm = kwds.get('algorithm', 'auto')
1179
if algorithm not in ['dhsw', 'pari', 'auto', 'no_chomp', 'chomp']:
1180
raise NotImplementedError('algorithm not recognized')
1181
if algorithm == 'auto' \
1182
and (base_ring == ZZ or (base_ring.is_prime_field() and base_ring != QQ)) \
1183
and have_chomp('homchain'):
1184
algorithm = 'chomp'
1185
if algorithm == 'chomp':
1186
return self._homology_chomp(deg, base_ring, verbose, generators)
1187
1188
if deg is None:
1189
deg = self.nonzero_degrees()
1190
if isinstance(deg, (list, tuple)):
1191
answer = {}
1192
for deg in self.nonzero_degrees():
1193
answer[deg] = self._homology_in_degree(deg, base_ring, verbose, generators, algorithm)
1194
return answer
1195
else:
1196
return self._homology_in_degree(deg, base_ring, verbose, generators, algorithm)
1197
1198
def _homology_in_degree(self, deg, base_ring, verbose, generators, algorithm):
1199
"""
1200
Helper method for :meth:`homology`.
1201
1202
EXAMPLES::
1203
1204
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])})
1205
sage: C.homology(1) == C._homology_in_degree(1, ZZ, False, False, 'auto')
1206
True
1207
"""
1208
if deg not in self.nonzero_degrees():
1209
zero_homology = HomologyGroup(0, base_ring)
1210
if generators:
1211
return (zero_homology, vector(base_ring, []))
1212
else:
1213
return zero_homology
1214
if verbose:
1215
print('Computing homology of the chain complex in dimension %s...' % n)
1216
1217
fraction_field = base_ring.fraction_field()
1218
def change_ring(X):
1219
if X.base_ring() is base_ring:
1220
return X
1221
return X.change_ring(base_ring)
1222
1223
# d_out is the differential going out of degree deg,
1224
# d_in is the differential entering degree deg
1225
differential = self.degree_of_differential()
1226
d_in = change_ring(self.differential(deg - differential))
1227
d_out = change_ring(self.differential(deg))
1228
d_out_rank = self.rank(deg, ring=fraction_field)
1229
d_out_nullity = d_out.ncols() - d_out_rank
1230
1231
if d_in.is_zero():
1232
if generators: #Include the generators of the nullspace
1233
return [(HomologyGroup(1, base_ring), self({deg:gen}))
1234
for gen in d_out.right_kernel().basis()]
1235
else:
1236
return HomologyGroup(d_out_nullity, base_ring)
1237
1238
if generators:
1239
orders, gens = self._homology_generators_snf(d_in, d_out, d_out_rank)
1240
answer = [(HomologyGroup(1, base_ring, [order]), self({deg:gen}))
1241
for order, gen in zip(orders, gens)]
1242
else:
1243
if base_ring.is_field():
1244
d_in_rank = self.rank(deg-differential, ring=base_ring)
1245
answer = HomologyGroup(d_out_nullity - d_in_rank, base_ring)
1246
elif base_ring == ZZ:
1247
if d_in.ncols() == 0:
1248
all_divs = [0] * d_out_nullity
1249
else:
1250
if algorithm in ['auto', 'no_chomp']:
1251
if ((d_in.ncols() > 300 and d_in.nrows() > 300)
1252
or (min(d_in.ncols(), d_in.nrows()) > 100 and
1253
d_in.ncols() + d_in.nrows() > 600)):
1254
algorithm = 'dhsw'
1255
else:
1256
algorithm = 'pari'
1257
if algorithm == 'dhsw':
1258
from sage.homology.matrix_utils import dhsw_snf
1259
all_divs = dhsw_snf(d_in, verbose=verbose)
1260
elif algorithm == 'pari':
1261
all_divs = d_in.elementary_divisors(algorithm)
1262
else:
1263
raise ValueError('unsupported algorithm')
1264
all_divs = all_divs[:d_out_nullity]
1265
# divisors equal to 1 produce trivial
1266
# summands, so filter them out
1267
divisors = filter(lambda x: x != 1, all_divs)
1268
answer = HomologyGroup(len(divisors), base_ring, divisors)
1269
else:
1270
raise NotImplementedError('only base rings ZZ and fields are supported')
1271
return answer
1272
1273
def _homology_generators_snf(self, d_in, d_out, d_out_rank):
1274
"""
1275
Compute the homology generators using the Smith normal form.
1276
1277
EXAMPLES::
1278
1279
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])})
1280
sage: C.homology(1)
1281
Z x C3
1282
sage: C._homology_generators_snf(C.differential(0), C.differential(1), 0)
1283
([3, 0], [(1, 0), (0, 1)])
1284
"""
1285
# Find the kernel of the out-going differential.
1286
K = d_out.right_kernel().matrix().transpose().change_ring(d_out.base_ring())
1287
1288
# Compute the induced map to the kernel
1289
S = K.augment(d_in).hermite_form()
1290
d_in_induced = S.submatrix(row=0, nrows=d_in.nrows()-d_out_rank,
1291
col=d_in.nrows()-d_out_rank, ncols=d_in.ncols())
1292
1293
# Find the SNF of the induced matrix and appropriate generators
1294
(N, P, Q) = d_in_induced.smith_form()
1295
all_divs = [0]*N.nrows()
1296
non_triv = 0
1297
for i in range(0, N.nrows()):
1298
if i >= N.ncols():
1299
break
1300
all_divs[i] = N[i][i]
1301
if N[i][i] == 1:
1302
non_triv = non_triv + 1
1303
divisors = filter(lambda x: x != 1, all_divs)
1304
gens = (K * P.inverse().submatrix(col=non_triv)).columns()
1305
return divisors, gens
1306
1307
def betti(self, deg=None, base_ring=None):
1308
"""
1309
The Betti number the chain complex.
1310
1311
That is, write the homology in this degree as a direct sum
1312
of a free module and a torsion module; the Betti number is the
1313
rank of the free summand.
1314
1315
INPUT:
1316
1317
- ``deg`` -- an element of the grading group for the chain
1318
complex or None (default ``None``); if ``None``,
1319
then return every Betti number, as a dictionary indexed by
1320
degree, or if an element of the grading group, then return
1321
the Betti number in that degree
1322
1323
- ``base_ring`` -- a commutative ring (optional, default is the
1324
base ring for the chain complex); compute homology with
1325
these coefficients -- must be either the integers or a
1326
field
1327
1328
OUTPUT:
1329
1330
The Betti number in degree ``deg`` -- the rank of the free
1331
part of the homology module in this degree.
1332
1333
EXAMPLES::
1334
1335
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])})
1336
sage: C.betti(0)
1337
2
1338
sage: [C.betti(n) for n in range(5)]
1339
[2, 1, 0, 0, 0]
1340
sage: C.betti()
1341
{0: 2, 1: 1}
1342
"""
1343
if base_ring is None:
1344
base_ring = QQ
1345
try:
1346
base_ring = base_ring.fraction_field()
1347
except AttributeError:
1348
raise NotImplementedError('only implemented if the base ring is ZZ or a field')
1349
H = self.homology(deg, base_ring=base_ring)
1350
if isinstance(H, dict):
1351
return dict([deg, homology_group.dimension()] for deg, homology_group in H.iteritems())
1352
else:
1353
return H.dimension()
1354
1355
def torsion_list(self, max_prime, min_prime=2):
1356
r"""
1357
Look for torsion in this chain complex by computing its mod `p`
1358
homology for a range of primes `p`.
1359
1360
INPUT:
1361
1362
- ``max_prime`` -- prime number; search for torsion mod `p` for
1363
all `p` strictly less than this number
1364
1365
- ``min_prime`` -- prime (optional, default 2); search for
1366
torsion mod `p` for primes at least as big as this
1367
1368
Return a list of pairs `(p, d)` where `p` is a prime at which
1369
there is torsion and `d` is a list of dimensions in which this
1370
torsion occurs.
1371
1372
The base ring for the chain complex must be the integers; if
1373
not, an error is raised.
1374
1375
ALGORITHM:
1376
1377
let `C` denote the chain complex. Let `P` equal
1378
``max_prime``. Compute the mod `P` homology of `C`, and use
1379
this as the base-line computation: the assumption is that this
1380
is isomorphic to the integral homology tensored with
1381
`\GF{P}`. Then compute the mod `p` homology for a range of
1382
primes `p`, and record whenever the answer differs from the
1383
base-line answer.
1384
1385
EXAMPLES::
1386
1387
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])})
1388
sage: C.homology()
1389
{0: Z x Z, 1: Z x C3}
1390
sage: C.torsion_list(11)
1391
[(3, [1])]
1392
sage: C = ChainComplex([matrix(ZZ, 1, 1, [2]), matrix(ZZ, 1, 1), matrix(1, 1, [3])])
1393
sage: C.homology(1)
1394
C2
1395
sage: C.homology(3)
1396
C3
1397
sage: C.torsion_list(5)
1398
[(2, [1]), (3, [3])]
1399
"""
1400
if self.base_ring() != ZZ:
1401
raise NotImplementedError('only implemented for base ring the integers')
1402
answer = []
1403
torsion_free = self.betti(base_ring=GF(max_prime))
1404
for p in prime_range(min_prime, max_prime):
1405
mod_p_betti = self.betti(base_ring=GF(p))
1406
if mod_p_betti != torsion_free:
1407
diff_dict = {}
1408
temp_diff = {}
1409
D = self.degree_of_differential()
1410
for i in torsion_free:
1411
temp_diff[i] = mod_p_betti.get(i, 0) - torsion_free[i]
1412
for i in temp_diff:
1413
if temp_diff[i] > 0:
1414
if i+D in diff_dict:
1415
lower = diff_dict[i+D]
1416
else:
1417
lower = 0
1418
current = temp_diff[i]
1419
if current > lower:
1420
diff_dict[i] = current - lower
1421
if i-D in diff_dict:
1422
diff_dict[i-D] -= current - lower
1423
differences = []
1424
for i in diff_dict:
1425
if diff_dict[i] != 0:
1426
differences.append(i)
1427
answer.append((p,differences))
1428
return answer
1429
1430
def _Hom_(self, other, category=None):
1431
"""
1432
Return the set of chain maps between chain complexes ``self``
1433
and ``other``.
1434
1435
EXAMPLES::
1436
1437
sage: S = simplicial_complexes.Sphere(2)
1438
sage: T = simplicial_complexes.Torus()
1439
sage: C = S.chain_complex(augmented=True,cochain=True)
1440
sage: D = T.chain_complex(augmented=True,cochain=True)
1441
sage: Hom(C,D) # indirect doctest
1442
Set of Morphisms from Chain complex with at most 4 nonzero terms over
1443
Integer Ring to Chain complex with at most 4 nonzero terms over Integer
1444
Ring in Category of chain complexes over Integer Ring
1445
"""
1446
from sage.homology.chain_complex_homspace import ChainComplexHomspace
1447
return ChainComplexHomspace(self, other)
1448
1449
def _flip_(self):
1450
"""
1451
Flip chain complex upside down (degree `n` gets changed to
1452
degree `-n`), thus turning a chain complex into a cochain complex
1453
without changing the homology (except for flipping it, too).
1454
1455
EXAMPLES::
1456
1457
sage: C = ChainComplex({2: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])})
1458
sage: C.degree_of_differential()
1459
1
1460
sage: C.differential(2)
1461
[3 0 0]
1462
[0 0 0]
1463
sage: C._flip_().degree_of_differential()
1464
-1
1465
sage: C._flip_().differential(-2)
1466
[3 0 0]
1467
[0 0 0]
1468
"""
1469
data = {}
1470
deg = self.degree_of_differential()
1471
for d in self.differential():
1472
data[-d] = self.differential()[d]
1473
return ChainComplex(data, degree=-deg)
1474
1475
def _chomp_repr_(self):
1476
r"""
1477
String representation of ``self`` suitable for use by the CHomP
1478
program.
1479
1480
Since CHomP can only handle chain complexes, not cochain
1481
complexes, and since it likes its complexes to start in degree
1482
0, flip the complex over if necessary, and shift it to start
1483
in degree 0. Note also that CHomP only works over the
1484
integers or a finite prime field.
1485
1486
EXAMPLES::
1487
1488
sage: C = ChainComplex({-2: matrix(ZZ, 1, 3, [3, 0, 0])}, degree=-1)
1489
sage: C._chomp_repr_()
1490
'chain complex\n\nmax dimension = 1\n\ndimension 0\n boundary a1 = 0\n\ndimension 1\n boundary a1 = + 3 * a1 \n boundary a2 = 0\n boundary a3 = 0\n\n'
1491
sage: C = ChainComplex({-2: matrix(ZZ, 1, 3, [3, 0, 0])}, degree=1)
1492
sage: C._chomp_repr_()
1493
'chain complex\n\nmax dimension = 1\n\ndimension 0\n boundary a1 = 0\n\ndimension 1\n boundary a1 = + 3 * a1 \n boundary a2 = 0\n boundary a3 = 0\n\n'
1494
"""
1495
deg = self.degree_of_differential()
1496
if (self.grading_group() != ZZ or
1497
(deg != 1 and deg != -1)):
1498
raise ValueError('CHomP only works on Z-graded chain complexes with '
1499
'differential of degree 1 or -1')
1500
base_ring = self.base_ring()
1501
if (base_ring == QQ) or (base_ring != ZZ and not (base_ring.is_prime_field())):
1502
raise ValueError('CHomP doesn\'t compute over the rationals, only over Z or F_p')
1503
if deg == -1:
1504
diffs = self.differential()
1505
else:
1506
diffs = self._flip_().differential()
1507
1508
if len(diffs) == 0:
1509
diffs = {0: matrix(ZZ, 0,0)}
1510
1511
maxdim = max(diffs)
1512
mindim = min(diffs)
1513
# will shift chain complex by subtracting mindim from
1514
# dimensions, so its bottom dimension is zero.
1515
s = "chain complex\n\nmax dimension = %s\n\n" % (maxdim - mindim - 1,)
1516
1517
for i in range(0, maxdim - mindim):
1518
s += "dimension %s\n" % i
1519
mat = diffs.get(i + mindim, matrix(base_ring, 0, 0))
1520
for idx in range(mat.ncols()):
1521
s += " boundary a%s = " % (idx + 1)
1522
# construct list of bdries
1523
col = mat.column(idx)
1524
if col.nonzero_positions():
1525
for j in col.nonzero_positions():
1526
entry = col[j]
1527
if entry > 0:
1528
sgn = "+"
1529
else:
1530
sgn = "-"
1531
entry = -entry
1532
s += "%s %s * a%s " % (sgn, entry, j+1)
1533
else:
1534
s += "0"
1535
s += "\n"
1536
s += "\n"
1537
return s
1538
1539
def _repr_(self):
1540
"""
1541
Print representation.
1542
1543
EXAMPLES::
1544
1545
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])})
1546
sage: C
1547
Chain complex with at most 2 nonzero terms over Integer Ring
1548
"""
1549
diffs = filter(lambda mat: mat.nrows() + mat.ncols() > 0,
1550
self._diff.values())
1551
if len(diffs) == 0:
1552
s = 'Trivial chain complex'
1553
else:
1554
s = 'Chain complex with at most {0} nonzero terms'.format(len(diffs)-1)
1555
s += ' over {0}'.format(self.base_ring())
1556
return s
1557
1558
def _ascii_art_(self):
1559
"""
1560
Return an ascii art representation.
1561
1562
Note that arrows go to the left so that composition of
1563
differentials is the usual matrix multiplication.
1564
1565
EXAMPLES::
1566
1567
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0]), 1:zero_matrix(1,2)})
1568
sage: ascii_art(C)
1569
[3 0 0]
1570
[0 0] [0 0 0]
1571
0 <-- C_2 <------ C_1 <-------- C_0 <-- 0
1572
1573
sage: one = matrix(ZZ, [[1]])
1574
sage: D = ChainComplex({0: one, 2: one, 6:one})
1575
sage: ascii_art(D)
1576
[1] [1] [0] [1]
1577
0 <-- C_7 <---- C_6 <-- 0 ... 0 <-- C_3 <---- C_2 <---- C_1 <---- C_0 <-- 0
1578
"""
1579
from sage.misc.ascii_art import AsciiArt
1580
1581
def arrow_art(n):
1582
d_n = self.differential(n)
1583
if d_n.nrows() == 0 or d_n.ncols() == 0:
1584
return AsciiArt(['<--'])
1585
d_str = [' '+line+' ' for line in str(d_n).splitlines()]
1586
arrow = '<' + '-'*(len(d_str[0])-1)
1587
d_str.append(arrow)
1588
return AsciiArt(d_str)
1589
1590
def module_art(n):
1591
C_n = self.free_module(n)
1592
if C_n.rank() == 0:
1593
return AsciiArt([' 0 '])
1594
else:
1595
return AsciiArt([' C_{0} '.format(n)])
1596
1597
result = []
1598
for ordered in self.ordered_degrees():
1599
ordered = list(reversed(ordered))
1600
if len(ordered) == 0:
1601
return AsciiArt(['0'])
1602
result_ordered = module_art(ordered[0] + self.degree_of_differential())
1603
for n in ordered:
1604
result_ordered += arrow_art(n) + module_art(n)
1605
result = [result_ordered] + result
1606
concatenated = result[0]
1607
for r in result[1:]:
1608
concatenated += AsciiArt([' ... ']) + r
1609
return concatenated
1610
1611
def _latex_(self):
1612
"""
1613
LaTeX print representation.
1614
1615
EXAMPLES::
1616
1617
sage: C = ChainComplex({0: matrix(ZZ, 2, 3, [3, 0, 0, 0, 0, 0])})
1618
sage: C._latex_()
1619
'\\Bold{Z}^{3} \\xrightarrow{d_{0}} \\Bold{Z}^{2}'
1620
"""
1621
# Warning: this is likely to screw up if, for example, the
1622
# degree of the differential is 2 and there are nonzero terms
1623
# in consecutive dimensions (e.g., in dimensions 0 and 1). In
1624
# such cases, the representation might show a differential
1625
# connecting these terms, although the differential goes from
1626
# dimension 0 to dimension 2, and from dimension 1 to
1627
# dimension 3, etc. I don't know how much effort should be
1628
# put into trying to fix this.
1629
string = ""
1630
dict = self._diff
1631
deg = self.degree_of_differential()
1632
ring = self.base_ring()
1633
if self.grading_group() != ZZ:
1634
guess = dict.keys()[0]
1635
if guess-deg in dict:
1636
string += "\\dots \\xrightarrow{d_{%s}} " % latex(guess-deg)
1637
string += _latex_module(ring, mat.ncols())
1638
string += " \\xrightarrow{d_{%s}} \\dots" % latex(guess)
1639
else:
1640
backwards = (deg < 0)
1641
sorted_list = sorted(dict.keys(), reverse=backwards)
1642
if len(dict) <= 6:
1643
for n in sorted_list[1:-1]:
1644
mat = dict[n]
1645
string += _latex_module(ring, mat.ncols())
1646
string += " \\xrightarrow{d_{%s}} " % latex(n)
1647
mat = dict[sorted_list[-1]]
1648
string += _latex_module(ring, mat.ncols())
1649
else:
1650
for n in sorted_list[:2]:
1651
mat = dict[n]
1652
string += _latex_module(ring, mat.ncols())
1653
string += " \\xrightarrow{d_{%s}} " % latex(n)
1654
string += "\\dots "
1655
n = sorted_list[-2]
1656
string += "\\xrightarrow{d_{%s}} " % latex(n)
1657
mat = dict[sorted_list[-1]]
1658
string += _latex_module(ring, mat.ncols())
1659
return string
1660
1661
1662
from sage.structure.sage_object import register_unpickle_override
1663
register_unpickle_override('sage.homology.chain_complex', 'ChainComplex', ChainComplex_class)
1664
1665
1666