Skip to content

abstractsimplicialcomplex

AbstractSimplicialComplex

Bases: UserList

Embodies the concept of abstract simplicial complex.

Abstract Simplicial Complexes, ASC, are closed set of elements that are topologically related. No geometrical properties are considered here.

Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
class AbstractSimplicialComplex(UserList):
    """Embodies the concept of abstract simplicial complex.

    Abstract Simplicial Complexes, ASC, are closed set of elements
    that are topologically related. No geometrical properties are
    considered here.
    """

    def __init__(self, indices: list[list[int]], 
                 name: Optional[str]=None) -> None:
        """Initializes an `AbstractSimplicialComplex` object.

        Parameters
        ----------
        indices : list of list of int
            The list of vertex indices forming the highest degree simplices.
        name : str, optional
            A name for the complex.
        """
        self.data = []
        self._top_simplex_orientations = None
        self._0_simplex_orientations = None
        self._chain_complex = []
        self._name = name
        self._adjacency_tensors = {}


        if _top_simplices_are_well_defined(indices):
            self.build_topology(indices)

    def __str__(self) -> str:
        """Returns the name of the complex.

        Returns
        -------
        str
            The name of the complex.
        """
        return self.name

    @property
    def name(self) -> str:
        """Gets the name of the complex.

        Returns
        -------
        str
            The name of the complex.
        """
        if self._name is None:
            return (
        f'{self.dim}D Abstract Simplicial Complex of shape {self.shape}.')
        else: 
            return self._name

    @name.setter
    def name(self, new_name: str) -> None:
        """Sets a custom name for the complex.

        Parameters
        ----------
        new_name : str
            The new name for the complex.
        """
        self._name = new_name

    def __contains__(self, other: AbstractSimplicialComplex) -> bool:
        """Checks if one abstract simplicial complex contains another one.

        Parameters
        ----------
        other : AbstractSimplicialComplex
            The other abstract simplicial complex to check.

        Returns
        -------
        bool
            True if the other abstract simplicial complex is contained, False otherwise.

        Examples
        --------
        >>> asc0 = AbstractSimplicialComplex([[1,2,3,4]])
        >>> asc1 = AbstractSimplicialComplex([[1,2,4]])
        >>> asc1 in asc0
        True
        """
        if self.dim < other.dim:
            return False
        else:
            other_splcs = np.array(other.data[-1])
            self_splcs = np.array(self.data[other.dim])
            for other_splx in other_splcs:
                if (other_splx == self_splcs).all(1).any()  == False:
                    return False
            return True

    @property
    def isclosed(self) -> bool:
        """Tests if the n-complex is homeomorphic to a n-ball.

        Returns
        -------
        bool
            True if the complex is homeomorphic to a n-ball, False otherwise.

        Notes
        -----
        * The name of this property is confusing for it tests if the complex 
          is homeomorphic to a ball. For instance, a torus would return False.
        * Maybe a better definition would be based on Betty numbers?
        * The Euler characteristics of a n-dimensional ball is 1 + (-1)**n.

        See also
        --------
        * The wikipedia page on euler characteristics,  
          [its a start...](https://en.wikipedia.org/wiki/Euler_characteristic)
        """
        n = self.dim 
        return self.euler_characteristic == 1 + (-1) ** n

    @property
    def ispure(self) -> bool:
        """Checks if the considered simplicial complex is pure.

        Returns
        -------
        bool
            True if the complex is pure, False otherwise.

        Notes
        -----
        * A simplicial n-complex is said pure iff
            every k-simplices within it is a face of at least one n-simplex.
        * For more insight, see:
          [course by Keenan Crane](https://brickisland.net/DDGSpring2019/).
        """
        if self.dim == 0:
            return True

        for mtrx in self._chain_complex[1:-1]:
            vct = abs(mtrx).sum(axis=1)

            if np.any(vct == 0):
                return False

        return True

    @property
    def dim(self) -> int:
        """Gets the topologic dimension of the complex.

        Returns
        -------
        int
            The topologic dimension of the complex.
        """
        return len(self) - 1

    @property
    def shape(self) -> tuple[int]:
        """Gets the numbers of simplices for each dimension.

        Returns
        -------
        tuple of int
            The numbers of simplices for each dimension.
        """
        return tuple([mdl.size for mdl in self])

    @property
    def euler_characteristic(self) -> int:
        """Computes the Euler characteristic of the complex.

        Returns
        -------
        int
            The Euler characteristic of the complex.

        Notes
        -----
        * The Euler characteristic (chi) of a CW-complex is given by:
                         chi = sum_i((-1)^i * #splx_i)
        """
        return sum([(-1) ** i * nbr for i, nbr in enumerate(self.shape)])


    def get_indices(self, splcs: list[list[int]]) -> list[int]:
        """Gets the indices of the given simplices.

        Parameters
        ----------
        splcs : list of list of int
            The list of k-simplices to get the indices of, 
            all the simplices should have the same dimensions

        Returns
        -------
        list of int
            The indices of the given simplices
        """
        idxs = []
        dim = len(splcs[0])-1
        splcs = np.sort(splcs)
        for i, splx_self in enumerate(self.data[dim]._vertex_indices):
            if len(np.where((splcs == splx_self).all(axis=1))[0])!=0:
                idxs.append(i)
        return idxs


    def orientation(self, dim:Optional[int]=None) -> Optional[np.ndarray[int]]:
        """Returns the orientation of highest or lowest degree simplices.

        Parameters
        ----------
        dim : int, optional
            The topological dimension of the simplices we want the orientation of. 
            Should be either `self.dim` or 0. If None, the orientation of the top simplices is returned.

        Returns
        -------
        numpy.ndarray of int, optional
            An (N,)-shape array filled with 1 and -1.
        """

        if (dim is None) or (dim == self.dim):
            return self._top_simplex_orientations
        elif dim == 0:
            return self._0_simplex_orientations
        else:
            logger.warning('Orientation only defined on highest simplices.')
            return None


    def build_topology(self, indices: list[list[int]]) -> None:
        """Computes all faces and incidence matrices from the given simplices.

        Parameters
        ----------
        indices : list of list of int
            List of the vertex indices forming the highest degree simplices.

        Notes
        -----
        * This method is the heart of the `AbstractSimplicialComplex` class.
        * It computes the simplices of all degrees from the list of 
          the highest degree ones.
        * It also computes a coherent orientation for the top simplices.
        """

        simplices = []
        highest_splcs, *smlr_splcs = _format_simplices(indices)

        n = len(highest_splcs[0]) - 1

        logger.info(f'Building a {n}-Abstract Simplicial Complex with '
                    + f'{len(highest_splcs)} {n}-simplices')

        simplices.append(highest_splcs)

        k = n
        while k >= 1:
            faces = _compute_faces(simplices[0])

            if (len(smlr_splcs)!=0) and (len(smlr_splcs[0][0])==k):
                faces = add_lower_level_simplices_to_faces(faces, smlr_splcs)

            mtrx = _compute_incidence_matrix(faces)

            if k == n:
                self._top_simplex_orientations = _compute_orientation(mtrx)
                # We multiply each column of the incidence matrix by the
                # orientation of the top-simplices so for each row (i.e. face)
                # we have one +1 entry and one -1 entry.
                mtrx = mtrx.multiply(self._top_simplex_orientations).tocsr()


            self._chain_complex.insert(0, mtrx)
            if k == 1:
                self._0_simplex_orientations = _compute_orientation(mtrx.T)

            trimed_faces = remove_duplicated_simplices(faces)
            simplices.insert(0, trimed_faces)

            k -= 1

        add_null_incidence_matrices(self._chain_complex, simplices)

        self.data = _format_complex(simplices, self._chain_complex)


    def star(self, simplices:dict[int, list[int]]|list[int]|int,
             dim:Optional[int]=None) -> dict[int, list[int]]:
        """Gets the list of all coboundaries of a given simplex.

        Parameters
        ----------
        simplices : dict of int to list of int or list of int or int
            The list of simplex indices forming the chain we want the star of.
        dim : int, optional
            The topological dimension of the simplices we want the star of.

        Returns
        -------
        dict of int to list of int
            The keys are topological dimensions k.
            The values are lists of the k-simplex indices belonging to 
            the seeked star ensemble.

        Notes
        -----
        * See the `.math.topology` module for details.
        """

        if isinstance(simplices, dict):
            return star(self._chain_complex, simplices)
        elif dimension_is_specified(dim):
            return star(self._chain_complex, {dim: simplices})
        else:
            return None


    def closure(self, simplices:dict[int, list[int]]|list[int]|int,
                dim:Optional[int]=None) -> dict[int, list[int]]:
        """Gets the smallest simplicial complex containing the given simplices.

        Parameters
        ----------
        simplices : dict of int to list of int or list of int or int
            The list of simplex indices forming the chain we want the closure of.
        dim : int, optional
            The topological dimension of the simplices we want the closure of.

        Returns
        -------
        dict of int to list of int
            The keys are topological dimensions k.
            The values are lists of the k-simplex indices belonging to 
            the seeked closure ensemble.

        Notes
        -----
        * See the `.math.topology` module for details.
        """

        if isinstance(simplices, dict):
            return closure(self._chain_complex, simplices)
        elif dimension_is_specified(dim):
            return closure(self._chain_complex, {dim: simplices})
        else:
            return None


    def link(self, simplices:dict[int, list[int]]|list[int]|int,
            dim:Optional[int]=None) -> dict[int, list[int]]:
        """Gets the topological sphere surrounding the given simplices.

        Parameters
        ----------
        simplices : dict of int to list of int or list of int or int
            The list of simplex indices forming the chain we want the link of.
        dim : int, optional
            The topological dimension of the simplices we want the link of.

        Returns
        -------
        dict of int to list of int
            The keys are topological dimensions k.
            The values are lists of the k-simplex indices belonging to 
            the seeked link ensemble.

        Notes
        -----
        * See the `.math.topology` module for details.
        """

        if isinstance(simplices, dict):
            return link(self._chain_complex, simplices)
        elif dimension_is_specified(dim):
            return link(self._chain_complex, {dim: simplices})
        else:
            return None


    def border(self, simplices:Optional[dict[int, list[int]]]=None
               ) -> Optional[dict[int, list[int]]]:
        """Gets the subcomplex boundary of a pure complex.

        Parameters
        ----------
        simplices : dict of int to list of int, optional
            The list of simplex indices forming the chain we want the border of.
            If None, the whole ASC is considered.

        Returns
        -------
        dict of int to list of int, optional
            The keys are topological dimensions k.
            The values are lists of the k-simplex indices belonging to 
            the seeked border ensemble.

        Notes
        -----
        * If no simplices are specified, the whole ASC is considered.
        * To simplicify we only compute borders for pure complexes.
        * The provided indices should specify a pure (sub-)complex.
        * 0-simplices are considered to be their own borders.
        * If the provided complex is closed, the output is None.
        """

        if simplices is None:
            if self.isclosed:
                return None
            else:
                simplices = {self.dim: [i for i,_ in enumerate(self[-1])]}

        return border(self._chain_complex, simplices)


    def interior(self, simplices:Optional[dict[int, list[int]]]=None
                 ) -> dict[int, list[int]]:
        """Computes the interior of a subset of simplices.

        Parameters
        ----------
        simplices : dict of int to list of int, optional
            The list of simplex indices forming the chain we want the interior of.
            If None, the whole ASC is considered.

        Returns
        -------
        dict of int to list of int
            The keys are topological dimensions k.
            The values are lists of the k-simplex indices belonging to 
            the seeked interior ensemble.

        Notes
        -----
        * If no simplices are specified, the whole ASC is considered.
        * The provided indices should specify a pure (sub-)complex.
        * 0-simplices have no interior. In this case the output is set to None.
        """

        if simplices is None:
            simplices = {self.dim: [i for i,_ in enumerate(self[-1])]}

        return interior(self._chain_complex, simplices)


    def faces(self, simplex_dim:int,
              face_dim:Optional[int]=None) -> list[list[int]]:
        """Returns the list of face indices for all k-simplices.

        Parameters
        ----------
        simplex_dim : int
            The topological dimension of the considered simplices.
        face_dim : int, optional
            The topological dimension of the desired faces. 
            If None, the method consider `face_dim = simplex_dim-1`.

        Returns
        -------
        list of list of int
            The list of face indices for all simplices in the considered 
            k-Module.

        Notes
        -----
        * The parameters must verify: 0 <= face_dim <= simplex_dim <=self.dim.
        * If face_dim == simplex_dim, the list of simplex indices is returned.
        * This method may seem redundant with the `closure` method.
          Actualy, it is meant to be applied globally on the whole complex
          while the `closure` method is restricted to small subsets of
          simplices.
        * Its design makes it much faster and therefore suitable for processing
          the whole structure at once, especially during instantiation.
        """

        face_dim = face_dim if face_dim is not None else simplex_dim-1

        if not _dimensions_are_well_ordered(face_dim, simplex_dim, self.dim):
            return None

        if face_dim == simplex_dim:
            N = self[simplex_dim].size
            return np.arange(N).reshape(N, 1)

        nbr_splx = self.shape[simplex_dim]
        mtrx = sp.identity(nbr_splx, dtype=int)

        k = simplex_dim
        while k > face_dim:
            mtrx = abs(self[k].boundary) @ mtrx / (k+1)
            k -= 1

        faces = [[] for _ in range(mtrx.shape[1])]
        for fid, sid  in zip(mtrx.nonzero()[0], mtrx.nonzero()[1]):
            faces[sid].append(fid)

        if not isinstance(faces[0], Iterable): faces = [[fid] for fid in faces]

        return faces


    def cofaces(self, simplex_dim:int, 
                coface_dim:Optional[int]=None, 
                ordered:bool=False) -> list[list[int]]:
        """Returns the list of k'-coface indices for all k-simplices (k'>k).

        Parameters
        ----------
        simplex_dim : int
            The topological dimension of the considered simplices.
        coface_dim : int, optional
            The topological dimension of the desired cofaces.
            Default is simplex_dim+1. 
        ordered : bool, optional
            If True, for each simplex, the coface indices 
            are returned in the following order so we can loop around from one 
            coface to the neighboring one. Default is False.

        Returns
        -------
        list of list of int
            The list of coface indices for all simplices
            in the considered k-Module.

        Notes
        -----
        * The parameters must verify: 0<=simplex_dim<=coface_dim<=self.dim
        * This method may seem redundant with the `star` method.
          Actually, it is meant to be applied globally on the whole complex
          while the `star` method is designed to be applied on small subsets 
          of simplices.
        * Its design makes it much faster and therefore suitable for 
          processing the whole structure at once, especially during 
          instantiation.
        """

        coface_dim = coface_dim if coface_dim is not None else simplex_dim+1

        if not _dimensions_are_well_ordered(simplex_dim, coface_dim, self.dim):
            return None

        nbr_splx = self.shape[simplex_dim]
        mtrx = sp.identity(nbr_splx, dtype=int)

        k = simplex_dim
        while k < coface_dim:
            mtrx = abs(self[k].coboundary) @ mtrx / (k+1)
            k += 1

        cofaces = [[] for _ in range(mtrx.shape[1])]
        for cfid, sid  in zip(mtrx.nonzero()[0], mtrx.nonzero()[1]):
            cofaces[sid].append(cfid)

        if ordered:
            return ordered_splx_ids_loops(cofaces, self[coface_dim].adjacency)
        else:
            return cofaces


    def incidence_matrix(self, top_dim:int, low_dim:int
                         ) -> Optional[sp.csr_matrix[int]]:
        """Computes the incidence matrix between k-simplices and l-simplices.

        Parameters
        ----------
        top_dim : int
            The highest topological dimension (k) of the simplices to consider.
        low_dim : int
            The smallest topological dimension (l) of the simplices to consider.

        Returns
        -------
        scipy.sparse.csr_matrix of int, optional
            A sparse matrix of shape (Nl, Nk), filled with 0s and 1s, where Nl 
            and Nk are respectively the numbers of l-simplices and k-simplices.

        Notes
        -----
        * The (i,j)-element is 1 if the 0-simplex of index i is a l-face of 
          the k-simplex of index j.
        * If the 2 dimension are egal, the identity matrix is returned. 
        """
        try:
            assert low_dim <= top_dim, 'low_dim must be smaller than top_dim.'
        except AssertionError as msg:
            logger.warning(msg)
            return None

        if low_dim == top_dim:
            return sp.identity(self[top_dim].size, dtype=int).tocsr()

        mtrx = abs(self[top_dim].boundary)
        k = top_dim - 1
        while k > low_dim:
            mtrx = abs(self[k].boundary) @ mtrx / (k+1)
            k -= 1

        return mtrx.astype(bool).astype(int)


    def adjacency_tensor(self, k:int, l:int) -> sprs.COO[int]:
        """Returns the vertex-based adjacency between k, l & k+l simplices.

        Parameters
        ----------
        k : int
            The first topological dimension to consider.
        l : int
            The second topological dimension to consider.

        Returns
        -------
        sparse.COO of int
            A 3rd-order tensor filled with 0, 1 & -1

        Notes
        -----
        * This method is needed to compute the wedge product between `Cochain` 
          instances.
        * This method uses the `COO` class from the `sparse` library. Because 
          implementing such a tensor in its dense version in `numpy` would have 
          been too costing for large structures and higher-order sparse arrays 
          at not yet available in `scpipy.sparse`.
        """
        key = k, l

        if self._adjacency_tensors.get(key) is None:
            self._adjacency_tensors[key] = _compute_adjacency_tensor(self, k, l)

            if k != l :
                self._adjacency_tensors[key[::-1]] = \
                    self._adjacency_tensors[key].transpose((0, 2, 1))

        return self._adjacency_tensors[key]

dim property

Gets the topologic dimension of the complex.

Returns:

Type Description
int

The topologic dimension of the complex.

euler_characteristic property

Computes the Euler characteristic of the complex.

Returns:

Type Description
int

The Euler characteristic of the complex.

Notes
  • The Euler characteristic (chi) of a CW-complex is given by: chi = sum_i((-1)^i * #splx_i)

isclosed property

Tests if the n-complex is homeomorphic to a n-ball.

Returns:

Type Description
bool

True if the complex is homeomorphic to a n-ball, False otherwise.

Notes
  • The name of this property is confusing for it tests if the complex is homeomorphic to a ball. For instance, a torus would return False.
  • Maybe a better definition would be based on Betty numbers?
  • The Euler characteristics of a n-dimensional ball is 1 + (-1)**n.
See also

ispure property

Checks if the considered simplicial complex is pure.

Returns:

Type Description
bool

True if the complex is pure, False otherwise.

Notes
  • A simplicial n-complex is said pure iff every k-simplices within it is a face of at least one n-simplex.
  • For more insight, see: course by Keenan Crane.

name property writable

Gets the name of the complex.

Returns:

Type Description
str

The name of the complex.

shape property

Gets the numbers of simplices for each dimension.

Returns:

Type Description
tuple of int

The numbers of simplices for each dimension.

__contains__(other)

Checks if one abstract simplicial complex contains another one.

Parameters:

Name Type Description Default
other AbstractSimplicialComplex

The other abstract simplicial complex to check.

required

Returns:

Type Description
bool

True if the other abstract simplicial complex is contained, False otherwise.

Examples:

>>> asc0 = AbstractSimplicialComplex([[1,2,3,4]])
>>> asc1 = AbstractSimplicialComplex([[1,2,4]])
>>> asc1 in asc0
True
Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
def __contains__(self, other: AbstractSimplicialComplex) -> bool:
    """Checks if one abstract simplicial complex contains another one.

    Parameters
    ----------
    other : AbstractSimplicialComplex
        The other abstract simplicial complex to check.

    Returns
    -------
    bool
        True if the other abstract simplicial complex is contained, False otherwise.

    Examples
    --------
    >>> asc0 = AbstractSimplicialComplex([[1,2,3,4]])
    >>> asc1 = AbstractSimplicialComplex([[1,2,4]])
    >>> asc1 in asc0
    True
    """
    if self.dim < other.dim:
        return False
    else:
        other_splcs = np.array(other.data[-1])
        self_splcs = np.array(self.data[other.dim])
        for other_splx in other_splcs:
            if (other_splx == self_splcs).all(1).any()  == False:
                return False
        return True

__init__(indices, name=None)

Initializes an AbstractSimplicialComplex object.

Parameters:

Name Type Description Default
indices list of list of int

The list of vertex indices forming the highest degree simplices.

required
name str

A name for the complex.

None
Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
def __init__(self, indices: list[list[int]], 
             name: Optional[str]=None) -> None:
    """Initializes an `AbstractSimplicialComplex` object.

    Parameters
    ----------
    indices : list of list of int
        The list of vertex indices forming the highest degree simplices.
    name : str, optional
        A name for the complex.
    """
    self.data = []
    self._top_simplex_orientations = None
    self._0_simplex_orientations = None
    self._chain_complex = []
    self._name = name
    self._adjacency_tensors = {}


    if _top_simplices_are_well_defined(indices):
        self.build_topology(indices)

__str__()

Returns the name of the complex.

Returns:

Type Description
str

The name of the complex.

Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
72
73
74
75
76
77
78
79
80
def __str__(self) -> str:
    """Returns the name of the complex.

    Returns
    -------
    str
        The name of the complex.
    """
    return self.name

adjacency_tensor(k, l)

Returns the vertex-based adjacency between k, l & k+l simplices.

Parameters:

Name Type Description Default
k int

The first topological dimension to consider.

required
l int

The second topological dimension to consider.

required

Returns:

Type Description
sparse.COO of int

A 3rd-order tensor filled with 0, 1 & -1

Notes
  • This method is needed to compute the wedge product between Cochain instances.
  • This method uses the COO class from the sparse library. Because implementing such a tensor in its dense version in numpy would have been too costing for large structures and higher-order sparse arrays at not yet available in scpipy.sparse.
Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
def adjacency_tensor(self, k:int, l:int) -> sprs.COO[int]:
    """Returns the vertex-based adjacency between k, l & k+l simplices.

    Parameters
    ----------
    k : int
        The first topological dimension to consider.
    l : int
        The second topological dimension to consider.

    Returns
    -------
    sparse.COO of int
        A 3rd-order tensor filled with 0, 1 & -1

    Notes
    -----
    * This method is needed to compute the wedge product between `Cochain` 
      instances.
    * This method uses the `COO` class from the `sparse` library. Because 
      implementing such a tensor in its dense version in `numpy` would have 
      been too costing for large structures and higher-order sparse arrays 
      at not yet available in `scpipy.sparse`.
    """
    key = k, l

    if self._adjacency_tensors.get(key) is None:
        self._adjacency_tensors[key] = _compute_adjacency_tensor(self, k, l)

        if k != l :
            self._adjacency_tensors[key[::-1]] = \
                self._adjacency_tensors[key].transpose((0, 2, 1))

    return self._adjacency_tensors[key]

border(simplices=None)

Gets the subcomplex boundary of a pure complex.

Parameters:

Name Type Description Default
simplices dict of int to list of int

The list of simplex indices forming the chain we want the border of. If None, the whole ASC is considered.

None

Returns:

Type Description
dict of int to list of int, optional

The keys are topological dimensions k. The values are lists of the k-simplex indices belonging to the seeked border ensemble.

Notes
  • If no simplices are specified, the whole ASC is considered.
  • To simplicify we only compute borders for pure complexes.
  • The provided indices should specify a pure (sub-)complex.
  • 0-simplices are considered to be their own borders.
  • If the provided complex is closed, the output is None.
Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
def border(self, simplices:Optional[dict[int, list[int]]]=None
           ) -> Optional[dict[int, list[int]]]:
    """Gets the subcomplex boundary of a pure complex.

    Parameters
    ----------
    simplices : dict of int to list of int, optional
        The list of simplex indices forming the chain we want the border of.
        If None, the whole ASC is considered.

    Returns
    -------
    dict of int to list of int, optional
        The keys are topological dimensions k.
        The values are lists of the k-simplex indices belonging to 
        the seeked border ensemble.

    Notes
    -----
    * If no simplices are specified, the whole ASC is considered.
    * To simplicify we only compute borders for pure complexes.
    * The provided indices should specify a pure (sub-)complex.
    * 0-simplices are considered to be their own borders.
    * If the provided complex is closed, the output is None.
    """

    if simplices is None:
        if self.isclosed:
            return None
        else:
            simplices = {self.dim: [i for i,_ in enumerate(self[-1])]}

    return border(self._chain_complex, simplices)

build_topology(indices)

Computes all faces and incidence matrices from the given simplices.

Parameters:

Name Type Description Default
indices list of list of int

List of the vertex indices forming the highest degree simplices.

required
Notes
  • This method is the heart of the AbstractSimplicialComplex class.
  • It computes the simplices of all degrees from the list of the highest degree ones.
  • It also computes a coherent orientation for the top simplices.
Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
def build_topology(self, indices: list[list[int]]) -> None:
    """Computes all faces and incidence matrices from the given simplices.

    Parameters
    ----------
    indices : list of list of int
        List of the vertex indices forming the highest degree simplices.

    Notes
    -----
    * This method is the heart of the `AbstractSimplicialComplex` class.
    * It computes the simplices of all degrees from the list of 
      the highest degree ones.
    * It also computes a coherent orientation for the top simplices.
    """

    simplices = []
    highest_splcs, *smlr_splcs = _format_simplices(indices)

    n = len(highest_splcs[0]) - 1

    logger.info(f'Building a {n}-Abstract Simplicial Complex with '
                + f'{len(highest_splcs)} {n}-simplices')

    simplices.append(highest_splcs)

    k = n
    while k >= 1:
        faces = _compute_faces(simplices[0])

        if (len(smlr_splcs)!=0) and (len(smlr_splcs[0][0])==k):
            faces = add_lower_level_simplices_to_faces(faces, smlr_splcs)

        mtrx = _compute_incidence_matrix(faces)

        if k == n:
            self._top_simplex_orientations = _compute_orientation(mtrx)
            # We multiply each column of the incidence matrix by the
            # orientation of the top-simplices so for each row (i.e. face)
            # we have one +1 entry and one -1 entry.
            mtrx = mtrx.multiply(self._top_simplex_orientations).tocsr()


        self._chain_complex.insert(0, mtrx)
        if k == 1:
            self._0_simplex_orientations = _compute_orientation(mtrx.T)

        trimed_faces = remove_duplicated_simplices(faces)
        simplices.insert(0, trimed_faces)

        k -= 1

    add_null_incidence_matrices(self._chain_complex, simplices)

    self.data = _format_complex(simplices, self._chain_complex)

closure(simplices, dim=None)

Gets the smallest simplicial complex containing the given simplices.

Parameters:

Name Type Description Default
simplices dict of int to list of int or list of int or int

The list of simplex indices forming the chain we want the closure of.

required
dim int

The topological dimension of the simplices we want the closure of.

None

Returns:

Type Description
dict of int to list of int

The keys are topological dimensions k. The values are lists of the k-simplex indices belonging to the seeked closure ensemble.

Notes
  • See the .math.topology module for details.
Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
def closure(self, simplices:dict[int, list[int]]|list[int]|int,
            dim:Optional[int]=None) -> dict[int, list[int]]:
    """Gets the smallest simplicial complex containing the given simplices.

    Parameters
    ----------
    simplices : dict of int to list of int or list of int or int
        The list of simplex indices forming the chain we want the closure of.
    dim : int, optional
        The topological dimension of the simplices we want the closure of.

    Returns
    -------
    dict of int to list of int
        The keys are topological dimensions k.
        The values are lists of the k-simplex indices belonging to 
        the seeked closure ensemble.

    Notes
    -----
    * See the `.math.topology` module for details.
    """

    if isinstance(simplices, dict):
        return closure(self._chain_complex, simplices)
    elif dimension_is_specified(dim):
        return closure(self._chain_complex, {dim: simplices})
    else:
        return None

cofaces(simplex_dim, coface_dim=None, ordered=False)

Returns the list of k'-coface indices for all k-simplices (k'>k).

Parameters:

Name Type Description Default
simplex_dim int

The topological dimension of the considered simplices.

required
coface_dim int

The topological dimension of the desired cofaces. Default is simplex_dim+1.

None
ordered bool

If True, for each simplex, the coface indices are returned in the following order so we can loop around from one coface to the neighboring one. Default is False.

False

Returns:

Type Description
list of list of int

The list of coface indices for all simplices in the considered k-Module.

Notes
  • The parameters must verify: 0<=simplex_dim<=coface_dim<=self.dim
  • This method may seem redundant with the star method. Actually, it is meant to be applied globally on the whole complex while the star method is designed to be applied on small subsets of simplices.
  • Its design makes it much faster and therefore suitable for processing the whole structure at once, especially during instantiation.
Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
def cofaces(self, simplex_dim:int, 
            coface_dim:Optional[int]=None, 
            ordered:bool=False) -> list[list[int]]:
    """Returns the list of k'-coface indices for all k-simplices (k'>k).

    Parameters
    ----------
    simplex_dim : int
        The topological dimension of the considered simplices.
    coface_dim : int, optional
        The topological dimension of the desired cofaces.
        Default is simplex_dim+1. 
    ordered : bool, optional
        If True, for each simplex, the coface indices 
        are returned in the following order so we can loop around from one 
        coface to the neighboring one. Default is False.

    Returns
    -------
    list of list of int
        The list of coface indices for all simplices
        in the considered k-Module.

    Notes
    -----
    * The parameters must verify: 0<=simplex_dim<=coface_dim<=self.dim
    * This method may seem redundant with the `star` method.
      Actually, it is meant to be applied globally on the whole complex
      while the `star` method is designed to be applied on small subsets 
      of simplices.
    * Its design makes it much faster and therefore suitable for 
      processing the whole structure at once, especially during 
      instantiation.
    """

    coface_dim = coface_dim if coface_dim is not None else simplex_dim+1

    if not _dimensions_are_well_ordered(simplex_dim, coface_dim, self.dim):
        return None

    nbr_splx = self.shape[simplex_dim]
    mtrx = sp.identity(nbr_splx, dtype=int)

    k = simplex_dim
    while k < coface_dim:
        mtrx = abs(self[k].coboundary) @ mtrx / (k+1)
        k += 1

    cofaces = [[] for _ in range(mtrx.shape[1])]
    for cfid, sid  in zip(mtrx.nonzero()[0], mtrx.nonzero()[1]):
        cofaces[sid].append(cfid)

    if ordered:
        return ordered_splx_ids_loops(cofaces, self[coface_dim].adjacency)
    else:
        return cofaces

faces(simplex_dim, face_dim=None)

Returns the list of face indices for all k-simplices.

Parameters:

Name Type Description Default
simplex_dim int

The topological dimension of the considered simplices.

required
face_dim int

The topological dimension of the desired faces. If None, the method consider face_dim = simplex_dim-1.

None

Returns:

Type Description
list of list of int

The list of face indices for all simplices in the considered k-Module.

Notes
  • The parameters must verify: 0 <= face_dim <= simplex_dim <=self.dim.
  • If face_dim == simplex_dim, the list of simplex indices is returned.
  • This method may seem redundant with the closure method. Actualy, it is meant to be applied globally on the whole complex while the closure method is restricted to small subsets of simplices.
  • Its design makes it much faster and therefore suitable for processing the whole structure at once, especially during instantiation.
Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
def faces(self, simplex_dim:int,
          face_dim:Optional[int]=None) -> list[list[int]]:
    """Returns the list of face indices for all k-simplices.

    Parameters
    ----------
    simplex_dim : int
        The topological dimension of the considered simplices.
    face_dim : int, optional
        The topological dimension of the desired faces. 
        If None, the method consider `face_dim = simplex_dim-1`.

    Returns
    -------
    list of list of int
        The list of face indices for all simplices in the considered 
        k-Module.

    Notes
    -----
    * The parameters must verify: 0 <= face_dim <= simplex_dim <=self.dim.
    * If face_dim == simplex_dim, the list of simplex indices is returned.
    * This method may seem redundant with the `closure` method.
      Actualy, it is meant to be applied globally on the whole complex
      while the `closure` method is restricted to small subsets of
      simplices.
    * Its design makes it much faster and therefore suitable for processing
      the whole structure at once, especially during instantiation.
    """

    face_dim = face_dim if face_dim is not None else simplex_dim-1

    if not _dimensions_are_well_ordered(face_dim, simplex_dim, self.dim):
        return None

    if face_dim == simplex_dim:
        N = self[simplex_dim].size
        return np.arange(N).reshape(N, 1)

    nbr_splx = self.shape[simplex_dim]
    mtrx = sp.identity(nbr_splx, dtype=int)

    k = simplex_dim
    while k > face_dim:
        mtrx = abs(self[k].boundary) @ mtrx / (k+1)
        k -= 1

    faces = [[] for _ in range(mtrx.shape[1])]
    for fid, sid  in zip(mtrx.nonzero()[0], mtrx.nonzero()[1]):
        faces[sid].append(fid)

    if not isinstance(faces[0], Iterable): faces = [[fid] for fid in faces]

    return faces

get_indices(splcs)

Gets the indices of the given simplices.

Parameters:

Name Type Description Default
splcs list of list of int

The list of k-simplices to get the indices of, all the simplices should have the same dimensions

required

Returns:

Type Description
list of int

The indices of the given simplices

Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
def get_indices(self, splcs: list[list[int]]) -> list[int]:
    """Gets the indices of the given simplices.

    Parameters
    ----------
    splcs : list of list of int
        The list of k-simplices to get the indices of, 
        all the simplices should have the same dimensions

    Returns
    -------
    list of int
        The indices of the given simplices
    """
    idxs = []
    dim = len(splcs[0])-1
    splcs = np.sort(splcs)
    for i, splx_self in enumerate(self.data[dim]._vertex_indices):
        if len(np.where((splcs == splx_self).all(axis=1))[0])!=0:
            idxs.append(i)
    return idxs

incidence_matrix(top_dim, low_dim)

Computes the incidence matrix between k-simplices and l-simplices.

Parameters:

Name Type Description Default
top_dim int

The highest topological dimension (k) of the simplices to consider.

required
low_dim int

The smallest topological dimension (l) of the simplices to consider.

required

Returns:

Type Description
scipy.sparse.csr_matrix of int, optional

A sparse matrix of shape (Nl, Nk), filled with 0s and 1s, where Nl and Nk are respectively the numbers of l-simplices and k-simplices.

Notes
  • The (i,j)-element is 1 if the 0-simplex of index i is a l-face of the k-simplex of index j.
  • If the 2 dimension are egal, the identity matrix is returned.
Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
def incidence_matrix(self, top_dim:int, low_dim:int
                     ) -> Optional[sp.csr_matrix[int]]:
    """Computes the incidence matrix between k-simplices and l-simplices.

    Parameters
    ----------
    top_dim : int
        The highest topological dimension (k) of the simplices to consider.
    low_dim : int
        The smallest topological dimension (l) of the simplices to consider.

    Returns
    -------
    scipy.sparse.csr_matrix of int, optional
        A sparse matrix of shape (Nl, Nk), filled with 0s and 1s, where Nl 
        and Nk are respectively the numbers of l-simplices and k-simplices.

    Notes
    -----
    * The (i,j)-element is 1 if the 0-simplex of index i is a l-face of 
      the k-simplex of index j.
    * If the 2 dimension are egal, the identity matrix is returned. 
    """
    try:
        assert low_dim <= top_dim, 'low_dim must be smaller than top_dim.'
    except AssertionError as msg:
        logger.warning(msg)
        return None

    if low_dim == top_dim:
        return sp.identity(self[top_dim].size, dtype=int).tocsr()

    mtrx = abs(self[top_dim].boundary)
    k = top_dim - 1
    while k > low_dim:
        mtrx = abs(self[k].boundary) @ mtrx / (k+1)
        k -= 1

    return mtrx.astype(bool).astype(int)

interior(simplices=None)

Computes the interior of a subset of simplices.

Parameters:

Name Type Description Default
simplices dict of int to list of int

The list of simplex indices forming the chain we want the interior of. If None, the whole ASC is considered.

None

Returns:

Type Description
dict of int to list of int

The keys are topological dimensions k. The values are lists of the k-simplex indices belonging to the seeked interior ensemble.

Notes
  • If no simplices are specified, the whole ASC is considered.
  • The provided indices should specify a pure (sub-)complex.
  • 0-simplices have no interior. In this case the output is set to None.
Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
def interior(self, simplices:Optional[dict[int, list[int]]]=None
             ) -> dict[int, list[int]]:
    """Computes the interior of a subset of simplices.

    Parameters
    ----------
    simplices : dict of int to list of int, optional
        The list of simplex indices forming the chain we want the interior of.
        If None, the whole ASC is considered.

    Returns
    -------
    dict of int to list of int
        The keys are topological dimensions k.
        The values are lists of the k-simplex indices belonging to 
        the seeked interior ensemble.

    Notes
    -----
    * If no simplices are specified, the whole ASC is considered.
    * The provided indices should specify a pure (sub-)complex.
    * 0-simplices have no interior. In this case the output is set to None.
    """

    if simplices is None:
        simplices = {self.dim: [i for i,_ in enumerate(self[-1])]}

    return interior(self._chain_complex, simplices)

Gets the topological sphere surrounding the given simplices.

Parameters:

Name Type Description Default
simplices dict of int to list of int or list of int or int

The list of simplex indices forming the chain we want the link of.

required
dim int

The topological dimension of the simplices we want the link of.

None

Returns:

Type Description
dict of int to list of int

The keys are topological dimensions k. The values are lists of the k-simplex indices belonging to the seeked link ensemble.

Notes
  • See the .math.topology module for details.
Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
def link(self, simplices:dict[int, list[int]]|list[int]|int,
        dim:Optional[int]=None) -> dict[int, list[int]]:
    """Gets the topological sphere surrounding the given simplices.

    Parameters
    ----------
    simplices : dict of int to list of int or list of int or int
        The list of simplex indices forming the chain we want the link of.
    dim : int, optional
        The topological dimension of the simplices we want the link of.

    Returns
    -------
    dict of int to list of int
        The keys are topological dimensions k.
        The values are lists of the k-simplex indices belonging to 
        the seeked link ensemble.

    Notes
    -----
    * See the `.math.topology` module for details.
    """

    if isinstance(simplices, dict):
        return link(self._chain_complex, simplices)
    elif dimension_is_specified(dim):
        return link(self._chain_complex, {dim: simplices})
    else:
        return None

orientation(dim=None)

Returns the orientation of highest or lowest degree simplices.

Parameters:

Name Type Description Default
dim int

The topological dimension of the simplices we want the orientation of. Should be either self.dim or 0. If None, the orientation of the top simplices is returned.

None

Returns:

Type Description
numpy.ndarray of int, optional

An (N,)-shape array filled with 1 and -1.

Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
def orientation(self, dim:Optional[int]=None) -> Optional[np.ndarray[int]]:
    """Returns the orientation of highest or lowest degree simplices.

    Parameters
    ----------
    dim : int, optional
        The topological dimension of the simplices we want the orientation of. 
        Should be either `self.dim` or 0. If None, the orientation of the top simplices is returned.

    Returns
    -------
    numpy.ndarray of int, optional
        An (N,)-shape array filled with 1 and -1.
    """

    if (dim is None) or (dim == self.dim):
        return self._top_simplex_orientations
    elif dim == 0:
        return self._0_simplex_orientations
    else:
        logger.warning('Orientation only defined on highest simplices.')
        return None

star(simplices, dim=None)

Gets the list of all coboundaries of a given simplex.

Parameters:

Name Type Description Default
simplices dict of int to list of int or list of int or int

The list of simplex indices forming the chain we want the star of.

required
dim int

The topological dimension of the simplices we want the star of.

None

Returns:

Type Description
dict of int to list of int

The keys are topological dimensions k. The values are lists of the k-simplex indices belonging to the seeked star ensemble.

Notes
  • See the .math.topology module for details.
Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
def star(self, simplices:dict[int, list[int]]|list[int]|int,
         dim:Optional[int]=None) -> dict[int, list[int]]:
    """Gets the list of all coboundaries of a given simplex.

    Parameters
    ----------
    simplices : dict of int to list of int or list of int or int
        The list of simplex indices forming the chain we want the star of.
    dim : int, optional
        The topological dimension of the simplices we want the star of.

    Returns
    -------
    dict of int to list of int
        The keys are topological dimensions k.
        The values are lists of the k-simplex indices belonging to 
        the seeked star ensemble.

    Notes
    -----
    * See the `.math.topology` module for details.
    """

    if isinstance(simplices, dict):
        return star(self._chain_complex, simplices)
    elif dimension_is_specified(dim):
        return star(self._chain_complex, {dim: simplices})
    else:
        return None

add_lower_level_simplices_to_faces(faces, smlr_splcs)

Add lower degree simplices provided as input to the complex.

Parameters:

Name Type Description Default
faces numpy.ndarray of int

The array of faces to which lower degree simplices will be added.

required
smlr_splcs list of int

The list of lower degree simplices to add.

required

Returns:

Type Description
numpy.ndarray of int

The updated array of faces with lower degree simplices added.

Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
def add_lower_level_simplices_to_faces(faces:np.ndarray[int], 
                              smlr_splcs:list[int]) -> np.ndarray[int]:
    """Add lower degree simplices provided as input to the complex.

    Parameters
    ----------
    faces : numpy.ndarray of int
        The array of faces to which lower degree simplices will be added.
    smlr_splcs : list of int
        The list of lower degree simplices to add.

    Returns
    -------
    numpy.ndarray of int
        The updated array of faces with lower degree simplices added.
    """
    splcs = smlr_splcs.pop(0)
    splx_nbr = splcs.shape[0]

    splcs = np.hstack((splcs, np.zeros((splx_nbr, 2),
                                        dtype=int)))

    return np.vstack((faces, splcs))

add_null_incidence_matrices(chain_complex, simplices)

Add null matrices on both sides of the chain complex.

Parameters:

Name Type Description Default
chain_complex list of scipy.sparse.csr_matrix of int

The list of incidence matrices forming the chain complex.

required
simplices numpy.ndarray of int

The array of simplices in the complex.

required
Notes

In order to compute properly the differential operators, we need to be able to send 0-simplices to 0 with the boundary function and the n-simplices as well with the coboundary one.

Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
def add_null_incidence_matrices(chain_complex: list[sp.csr_matrix[int]],
                                simplices: np.array[int]) -> None:
    """Add null matrices on both sides of the chain complex.

    Parameters
    ----------
    chain_complex : list of scipy.sparse.csr_matrix of int
        The list of incidence matrices forming the chain complex.
    simplices : numpy.ndarray of int
        The array of simplices in the complex.

    Notes
    -----
    In order to compute properly the differential operators, we need to be able
    to send 0-simplices to 0 with the `boundary` function and the n-simplices as 
    well with the `coboundary` one.
    """

    for k in [0, -1]:
        Nk = simplices[k].shape[0]
        mtrx = sp.csr_matrix((Nk, Nk), dtype=int)

        if k == 0:
            chain_complex.insert(0, mtrx)
        else:
            chain_complex.append(mtrx.T)

dimension_is_specified(dim)

Checks that the parameter dim is not None.

Parameters:

Name Type Description Default
dim int

The dimension to check.

required

Returns:

Type Description
bool

True if the dimension is specified, False otherwise.

Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
def dimension_is_specified(dim: Optional[int]) -> bool:
    """Checks that the parameter `dim` is not None.

    Parameters
    ----------
    dim : int, optional
        The dimension to check.

    Returns
    -------
    bool
        True if the dimension is specified, False otherwise.
    """
    try:
        assert dim is not None, 'A dimension must be specified.'
        return True

    except AssertionError as msg:
        logger.warning(msg)
        return False

order_index_loop(indices, rows, cols)

Sorts a list of simplex indices forming a loop in a following order.

Parameters:

Name Type Description Default
indices list of int

List of indices to order.

required
rows numpy.ndarray of int

Rows of nonzero elements in the simplex adjacency matrix.

required
cols numpy.ndarray of int

Columns of nonzero elements in the simplex adjacency matrix.

required

Returns:

Type Description
list of int

The seeked sorted list.

Notes
  • Prior to running this function, an adjacency matrix must be computed.
  • The attributes rows and cols are extracted from such a matrix: cols, rows = adjacency.nonzero()
Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
def order_index_loop(indices: list[int],
                     rows: np.ndarray[int],
                     cols: np.ndarray[int]) -> list[int]:
    """Sorts a list of simplex indices forming a loop in a following order.

    Parameters
    ----------
    indices : list of int
        List of indices to order.
    rows : numpy.ndarray of int
        Rows of nonzero elements in the simplex adjacency matrix.
    cols : numpy.ndarray of int
        Columns of nonzero elements in the simplex adjacency matrix.

    Returns
    -------
    list of int
        The seeked sorted list.

    Notes
    -----
    * Prior to running this function, an adjacency matrix must be computed.
    * The attributes rows and cols are extracted from such a matrix:
        `cols, rows = adjacency.nonzero()`
    """

    ordered_indices = [indices.pop(0)]

    while len(indices) > 1:
        next_cfid = list(set(rows[cols == ordered_indices[-1]]) &
                         set(indices))[0]
        ordered_indices.append(next_cfid)
        indices.remove(next_cfid)

    ordered_indices += indices

    return ordered_indices

ordered_splx_ids_loops(indices, adjacency)

Organizes face/coface indices in following order.

Parameters:

Name Type Description Default
indices list of list of int

All the lists of the k-simplex indices to sort.

required
adjacency scipy.sparse.csr_matrix of int

The adjacency matrix of the corresponding k-simplices.

required

Returns:

Type Description
list of list of int

Same as the input but each element is now in following order.

Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
def ordered_splx_ids_loops(indices:list[list[int]], 
                           adjacency:sp.csr_matrix[int]) -> list[list[int]]:
    """Organizes face/coface indices in following order.

    Parameters
    ----------
    indices : list of list of int
        All the lists of the k-simplex indices to sort.
    adjacency : scipy.sparse.csr_matrix of int
        The adjacency matrix of the corresponding k-simplices.

    Returns
    -------
    list of list of int
        Same as the input but each element is now in following order.
    """
    rows, cols = adjacency.nonzero()

    return [order_index_loop(ids, rows, cols) for ids in indices]

remove_duplicated_simplices(faces)

Remove simplices that could have been generated twice.

Parameters:

Name Type Description Default
faces numpy.ndarray of int

The array of faces from which duplicated simplices will be removed.

required

Returns:

Type Description
numpy.ndarray of int

The updated array of faces with duplicated simplices removed.

Source code in src/dxtr/complexes/abstractsimplicialcomplex.py
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
def remove_duplicated_simplices(faces: np.ndarray[int]) -> np.ndarray[int]:
    """Remove simplices that could have been generated twice.

    Parameters
    ----------
    faces : numpy.ndarray of int
        The array of faces from which duplicated simplices will be removed.

    Returns
    -------
    numpy.ndarray of int
        The updated array of faces with duplicated simplices removed.
    """
    _, ids = np.unique(faces[:, :-2], axis=0, return_index=True)

    return faces[ids, :-2]