Skip to content

module

Module

Bases: UserList

A container for all simplices of the same degree.

Notes
  • This k-Module object differs a bit from the exact mathematical definition: In theory, the k-Module of a simplicial complex is composed by all the k'-simplices with k' <= k. Here, our k-Module only contains k-simplices.
Source code in src/dxtr/complexes/module.py
 38
 39
 40
 41
 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
class Module(UserList):
    """A container for all simplices of the same degree.

    Notes
    -----
    * This k-Module object differs a bit from
      the exact mathematical definition: In theory,
      the k-Module of a simplicial complex is composed
      by all the k'-simplices  with k' <= k. Here, 
      our k-Module only contains k-simplices.
    """

    def __init__(self, simplices:np.ndarray[int],
                 chain_complex:list[sp.csr_matrix[int]]) -> None:
        """Instanciates a k-Module"""
        super().__init__([Simplex(tuple(splx)) for splx in simplices])

        splx_nbr, splx_dim = simplices.shape

        self._vertex_indices = simplices.reshape(splx_nbr) if splx_dim==1 else simplices
        self._chain_complex = chain_complex
        self._vertices = None
        self._volumes = None
        self._covolumes = None
        self._circumcenters = None
        self._deficit_angles = None
        self._dihedral_angles = None
        self._well_centeredness = None

    @property
    def dim(self) -> int:
        """Topological dimension (k) of the considered k-Module.
        """
        return len(self[0]) - 1

    @property
    def size(self) -> int:
        """Number of simplices within the k-Module.
        """
        return len(self)

    @property
    def boundary(self) -> sp.csr_matrix[int]:
        """Incidence matrix between these k-simplices and their (k-1)-faces.
        """
        return self._chain_complex[self.dim]

    @property
    def coboundary(self) -> sp.csc_matrix[int]:
        """Transpose of the incidence matrix between the (k+1)-cofaces 
           and these k-simplices.
        """
        return self._chain_complex[self.dim+1].T.tocsr(copy=False)

    @property
    def adjacency(self) -> Optional[sp.csc_matrix[int]]:
        """Adjacency matrix, computed through the faces.
        """
        try:
            bnd = self.boundary
            assert bnd.data.shape[0] != 0, ('Cannot compute adjacency '
                + 'for 0-simplices, try the `coadjacency` method instead.')

            adj = bnd.T @ bnd
            adj -= sp.diags(adj.diagonal(), dtype=int)
            return adj

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

    @property
    def coadjacency(self) -> Optional[sp.csr_matrix[int]]:
        """Coadjacency matrix, i,e. adjacency computed through the cofaces.
        """
        try:
            cbnd = self.coboundary
            assert cbnd.data.shape[0] != 0, ('Cannot compute co-adjacency '
                + 'for n-simplices, try the `adjacency` method instead.')

            cadj = cbnd.T @ cbnd
            cadj -= sp.diags(cadj.diagonal(), dtype=int)
            return cadj

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

    @property
    def vertex_indices(self) -> np.ndarray[int]:
        """Indices of all the k-simplices within the Module."""
        return self._vertex_indices

    @property
    def vertices(self) -> Optional[np.ndarray[float]]:
        """Position vectors of all the vertices around each k-simplex.

        Note
        -----
        * The returned array is of shape N*(k+1)*D where:
            - N is the number of k-simplices.
            - k+1 is the number of vertices defining a k-simplex
            - D is the dimension of the embedding euclidean space.
        * These positions are computed by the `build_geometry` method 
          within the `SimplicialComplex` Class.
        * Returns None, in the case of an `AbstractSimplicialComplex`.
        """
        if self._vertices is None:
            return None
        else:
            return self._vertices[self._vertex_indices]

    @property
    def volumes(self) -> Optional[np.ndarray[float]]:
        """Volumes of all k-simplices.

        Note
        -----
        * These volumes are computed by the `build_geometry` method
          within the `SimplicialComplex` class.
        """
        return self._volumes

    @property
    def covolumes(self) -> Optional[np.ndarray[float]]:
        """Covolumes of all k-simplices.

        Note
        -----
        * These covolumes are computed by the `build_dual_geometry` method 
          within the `SimplicialComplex` class.
        """
        return self._covolumes

    @property
    def circumcenters(self) -> Optional[np.ndarray[float]]:
        """Circumcenter position vectors for all k-simplex.

        Note
        -----
        * The returned array is of shape N*D where:
            - N is the number of k-simplices.
            - D is the dimension of the embedding euclidean space.
        * These circumcenters are computed by the `build_dual_geometry` method 
          within the `SimplicialManifold` Class.
        """
        return self._circumcenters

    @property
    def deficit_angles(self) -> Optional[np.ndarray[float]]:
        """Gaussian curvature for all k-simplices.

        Notes
        -----
          * For now, gaussian curvature can only be computed for 
            0- & 1-simplices.
        """
        return self._deficit_angles

    @property
    def dihedral_angles(self) -> Optional[np.ndarray[float]]:
        """Gaussian curvature for all k-simplices.

        Notes
        -----
          * For now, gaussian curvature can only be computed for 
            0- & 1-simplices.
        """
        return self._dihedral_angles

    @property
    def well_centeredness(self) -> np.ndarray[bool]:
        """Well-centeredness for all k-simplices.

        Notes
        -----
          * The 'well-centeredness' of a k-simplex corresponds to the fact 
            that its circumcenter lies within its volume.
          * It is a measure of the *quality* of the complex.
          * By construction, 0- and 1-simplices are always well-centered.
          * This property is important to compute the covolumes of the faces of 
            the simplices.
        """
        return self._well_centeredness

    def closest_simplex_to(self, target:Iterable, 
                           index_only:bool=False) -> Optional[int|Simplex]:
        """Gets the k-simplex closest to a given position.

        Parameters
        ----------
        target
            The desired position.
        index_only 
            Optional (default is False)
            If True, only the index of the simplex within the provided Module 
            is returned otherwise, the Simplex instance is returned.

        Returns
        -------
            The seeked Simplex or its index.
        """

        target = np.asarray(target)
        positions = self.circumcenters
        n = positions.shape[-1]

        if not _has_proper_shape(target, n): 
            return None

        distance_to_target = lng.norm(positions - target, axis=1)

        sidx = distance_to_target.argmin()

        if index_only:
            return sidx
        else:
            return self[sidx]

    def set(self, property_name:str, values:np.ndarray[float]) -> None:
        """Sets the values of the simplex geometrical properties.

        Notes
        -----
          * Property_name should either be `vertices`, `volumes`, 
            `covolumes`, `circumcenters`, `deficit angles`, `dihedral angles` 
            or 'well-centeredness'.
        """
        try:
            assert property_name in ['vertex indices', 'vertices', 'volumes',
                                     'covolumes', 'circumcenters', 
                                     'deficit angles', 'dihedral angles', 
                                     'well-centeredness'], (
                f'The provided property ({property_name}) is not supported.')

            if property_name == 'vertex indices':
                self._vertex_indices = values

            if property_name == 'vertices':
                self._vertices = values

                for splx in self:
                    splx._vertices = values

            elif property_name == 'volumes':
                self._volumes = values

            elif property_name == 'covolumes':
                self._covolumes = values

            elif property_name == 'circumcenters':
                self._circumcenters = values

            elif property_name == 'deficit angles':
                self._deficit_angles = values

            elif property_name == 'dihedral angles':
                self._dihedral_angles = values

            elif property_name == 'well-centeredness':
                self._well_centeredness = values     

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

adjacency property

Adjacency matrix, computed through the faces.

boundary property

Incidence matrix between these k-simplices and their (k-1)-faces.

circumcenters property

Circumcenter position vectors for all k-simplex.

Note
  • The returned array is of shape N*D where:
    • N is the number of k-simplices.
    • D is the dimension of the embedding euclidean space.
  • These circumcenters are computed by the build_dual_geometry method within the SimplicialManifold Class.

coadjacency property

Coadjacency matrix, i,e. adjacency computed through the cofaces.

coboundary property

Transpose of the incidence matrix between the (k+1)-cofaces and these k-simplices.

covolumes property

Covolumes of all k-simplices.

Note
  • These covolumes are computed by the build_dual_geometry method within the SimplicialComplex class.

deficit_angles property

Gaussian curvature for all k-simplices.

Notes
  • For now, gaussian curvature can only be computed for 0- & 1-simplices.

dihedral_angles property

Gaussian curvature for all k-simplices.

Notes
  • For now, gaussian curvature can only be computed for 0- & 1-simplices.

dim property

Topological dimension (k) of the considered k-Module.

size property

Number of simplices within the k-Module.

vertex_indices property

Indices of all the k-simplices within the Module.

vertices property

Position vectors of all the vertices around each k-simplex.

Note
  • The returned array is of shape N(k+1)D where:
    • N is the number of k-simplices.
    • k+1 is the number of vertices defining a k-simplex
    • D is the dimension of the embedding euclidean space.
  • These positions are computed by the build_geometry method within the SimplicialComplex Class.
  • Returns None, in the case of an AbstractSimplicialComplex.

volumes property

Volumes of all k-simplices.

Note
  • These volumes are computed by the build_geometry method within the SimplicialComplex class.

well_centeredness property

Well-centeredness for all k-simplices.

Notes
  • The 'well-centeredness' of a k-simplex corresponds to the fact that its circumcenter lies within its volume.
  • It is a measure of the quality of the complex.
  • By construction, 0- and 1-simplices are always well-centered.
  • This property is important to compute the covolumes of the faces of the simplices.

__init__(simplices, chain_complex)

Instanciates a k-Module

Source code in src/dxtr/complexes/module.py
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
def __init__(self, simplices:np.ndarray[int],
             chain_complex:list[sp.csr_matrix[int]]) -> None:
    """Instanciates a k-Module"""
    super().__init__([Simplex(tuple(splx)) for splx in simplices])

    splx_nbr, splx_dim = simplices.shape

    self._vertex_indices = simplices.reshape(splx_nbr) if splx_dim==1 else simplices
    self._chain_complex = chain_complex
    self._vertices = None
    self._volumes = None
    self._covolumes = None
    self._circumcenters = None
    self._deficit_angles = None
    self._dihedral_angles = None
    self._well_centeredness = None

closest_simplex_to(target, index_only=False)

Gets the k-simplex closest to a given position.

Parameters:

Name Type Description Default
target Iterable

The desired position.

required
index_only bool

Optional (default is False) If True, only the index of the simplex within the provided Module is returned otherwise, the Simplex instance is returned.

False

Returns:

Type Description
The seeked Simplex or its index.
Source code in src/dxtr/complexes/module.py
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
def closest_simplex_to(self, target:Iterable, 
                       index_only:bool=False) -> Optional[int|Simplex]:
    """Gets the k-simplex closest to a given position.

    Parameters
    ----------
    target
        The desired position.
    index_only 
        Optional (default is False)
        If True, only the index of the simplex within the provided Module 
        is returned otherwise, the Simplex instance is returned.

    Returns
    -------
        The seeked Simplex or its index.
    """

    target = np.asarray(target)
    positions = self.circumcenters
    n = positions.shape[-1]

    if not _has_proper_shape(target, n): 
        return None

    distance_to_target = lng.norm(positions - target, axis=1)

    sidx = distance_to_target.argmin()

    if index_only:
        return sidx
    else:
        return self[sidx]

set(property_name, values)

Sets the values of the simplex geometrical properties.

Notes
  • Property_name should either be vertices, volumes, covolumes, circumcenters, deficit angles, dihedral angles or 'well-centeredness'.
Source code in src/dxtr/complexes/module.py
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
def set(self, property_name:str, values:np.ndarray[float]) -> None:
    """Sets the values of the simplex geometrical properties.

    Notes
    -----
      * Property_name should either be `vertices`, `volumes`, 
        `covolumes`, `circumcenters`, `deficit angles`, `dihedral angles` 
        or 'well-centeredness'.
    """
    try:
        assert property_name in ['vertex indices', 'vertices', 'volumes',
                                 'covolumes', 'circumcenters', 
                                 'deficit angles', 'dihedral angles', 
                                 'well-centeredness'], (
            f'The provided property ({property_name}) is not supported.')

        if property_name == 'vertex indices':
            self._vertex_indices = values

        if property_name == 'vertices':
            self._vertices = values

            for splx in self:
                splx._vertices = values

        elif property_name == 'volumes':
            self._volumes = values

        elif property_name == 'covolumes':
            self._covolumes = values

        elif property_name == 'circumcenters':
            self._circumcenters = values

        elif property_name == 'deficit angles':
            self._deficit_angles = values

        elif property_name == 'dihedral angles':
            self._dihedral_angles = values

        elif property_name == 'well-centeredness':
            self._well_centeredness = values     

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