scipy.sparse.linalg.

LaplacianNd#

class scipy.sparse.linalg.LaplacianNd(*args, **kwargs)[source]#

N 维度上的网格拉普拉斯算子和它的特征值/特征向量。

N 维度上的均匀矩形网格上构建拉普拉斯算子,并输出它的特征值和特征向量。拉普拉斯算子 L 是一个方形、负定的、实对称的数组,包含带符号的整数条目,其他位置为零。

参数::
grid_shape元组

一个长度为 N 的整数元组(对应于拉普拉斯算子的维度),其中每个条目给出该维度的尺寸。拉普拉斯矩阵的尺寸为 np.prod(grid_shape)

boundary_conditions{‘neumann’, ‘dirichlet’, ‘periodic’}, 可选

网格边界上的边界条件类型。有效值为 'dirichlet''neumann'``(默认) ``'periodic'

dtype数据类型

数组的数值类型。默认值为 np.int8

备注

与 MATLAB/Octave 中 1 维、2 维和 3 维拉普拉斯算子的实现相比 [1],此代码允许任意 N 维情况和无矩阵可调用选项,但目前仅限于纯狄利克雷、诺伊曼或周期边界条件。

矩形网格图的拉普拉斯矩阵 (scipy.sparse.csgraph.laplacian) 对应于具有诺伊曼条件的负拉普拉斯算子,即 boundary_conditions = 'neumann'

形状为 grid_shapeN 维规则网格的离散拉普拉斯算子的所有特征值和特征向量在分析上是已知的 [2]。

参考文献

[2]

“二阶导数的特征值和特征向量”,维基百科 https://en.wikipedia.org/wiki/Eigenvalues_and_eigenvectors_of_the_second_derivative

示例

>>> import numpy as np
>>> from scipy.sparse.linalg import LaplacianNd
>>> from scipy.sparse import diags, csgraph
>>> from scipy.linalg import eigvalsh

下面展示了在具有 n=6 个网格点的规则网格上,对纯诺伊曼边界条件进行的一维拉普拉斯算子演示,它与使用稀疏邻接矩阵 G 表示的具有 n 个顶点的无向线性图的负图拉普拉斯算子完全相同,该矩阵由著名的三对角矩阵表示。

>>> n = 6
>>> G = diags(np.ones(n - 1), 1, format='csr')
>>> Lf = csgraph.laplacian(G, symmetrized=True, form='function')
>>> grid_shape = (n, )
>>> lap = LaplacianNd(grid_shape, boundary_conditions='neumann')
>>> np.array_equal(lap.matmat(np.eye(n)), -Lf(np.eye(n)))
True

由于拉普拉斯算子的所有矩阵条目都是整数,因此 'int8' 是用于存储矩阵表示的默认数据类型。

>>> lap.tosparse()
<DIAgonal sparse array of dtype 'int8'
    with 16 stored elements (3 diagonals) and shape (6, 6)>
>>> lap.toarray()
array([[-1,  1,  0,  0,  0,  0],
       [ 1, -2,  1,  0,  0,  0],
       [ 0,  1, -2,  1,  0,  0],
       [ 0,  0,  1, -2,  1,  0],
       [ 0,  0,  0,  1, -2,  1],
       [ 0,  0,  0,  0,  1, -1]], dtype=int8)
>>> np.array_equal(lap.matmat(np.eye(n)), lap.toarray())
True
>>> np.array_equal(lap.tosparse().toarray(), lap.toarray())
True

可以计算任意数量的极值特征值和/或特征向量。

>>> lap = LaplacianNd(grid_shape, boundary_conditions='periodic')
>>> lap.eigenvalues()
array([-4., -3., -3., -1., -1.,  0.])
>>> lap.eigenvalues()[-2:]
array([-1.,  0.])
>>> lap.eigenvalues(2)
array([-1.,  0.])
>>> lap.eigenvectors(1)
array([[0.40824829],
       [0.40824829],
       [0.40824829],
       [0.40824829],
       [0.40824829],
       [0.40824829]])
>>> lap.eigenvectors(2)
array([[ 0.5       ,  0.40824829],
       [ 0.        ,  0.40824829],
       [-0.5       ,  0.40824829],
       [-0.5       ,  0.40824829],
       [ 0.        ,  0.40824829],
       [ 0.5       ,  0.40824829]])
>>> lap.eigenvectors()
array([[ 0.40824829,  0.28867513,  0.28867513,  0.5       ,  0.5       ,
         0.40824829],
       [-0.40824829, -0.57735027, -0.57735027,  0.        ,  0.        ,
         0.40824829],
       [ 0.40824829,  0.28867513,  0.28867513, -0.5       , -0.5       ,
         0.40824829],
       [-0.40824829,  0.28867513,  0.28867513, -0.5       , -0.5       ,
         0.40824829],
       [ 0.40824829, -0.57735027, -0.57735027,  0.        ,  0.        ,
         0.40824829],
       [-0.40824829,  0.28867513,  0.28867513,  0.5       ,  0.5       ,
         0.40824829]])

二维拉普拉斯算子在每个维度上具有 grid_shape = (2, 3) 个点的规则网格上进行了说明。

>>> grid_shape = (2, 3)
>>> n = np.prod(grid_shape)

网格点的编号如下

>>> np.arange(n).reshape(grid_shape + (-1,))
array([[[0],
        [1],
        [2]],

       [[3],
        [4],
        [5]]])

每个边界条件 'dirichlet''periodic''neumann' 都分别进行了说明;对于 'dirichlet'

>>> lap = LaplacianNd(grid_shape, boundary_conditions='dirichlet')
>>> lap.tosparse()
<Compressed Sparse Row sparse array of dtype 'int8'
    with 20 stored elements and shape (6, 6)>
>>> lap.toarray()
array([[-4,  1,  0,  1,  0,  0],
       [ 1, -4,  1,  0,  1,  0],
       [ 0,  1, -4,  0,  0,  1],
       [ 1,  0,  0, -4,  1,  0],
       [ 0,  1,  0,  1, -4,  1],
       [ 0,  0,  1,  0,  1, -4]], dtype=int8)
>>> np.array_equal(lap.matmat(np.eye(n)), lap.toarray())
True
>>> np.array_equal(lap.tosparse().toarray(), lap.toarray())
True
>>> lap.eigenvalues()
array([-6.41421356, -5.        , -4.41421356, -3.58578644, -3.        ,
       -1.58578644])
>>> eigvals = eigvalsh(lap.toarray().astype(np.float64))
>>> np.allclose(lap.eigenvalues(), eigvals)
True
>>> np.allclose(lap.toarray() @ lap.eigenvectors(),
...             lap.eigenvectors() @ np.diag(lap.eigenvalues()))
True

对于 'periodic'

>>> lap = LaplacianNd(grid_shape, boundary_conditions='periodic')
>>> lap.tosparse()
<Compressed Sparse Row sparse array of dtype 'int8'
    with 24 stored elements and shape (6, 6)>
>>> lap.toarray()
    array([[-4,  1,  1,  2,  0,  0],
           [ 1, -4,  1,  0,  2,  0],
           [ 1,  1, -4,  0,  0,  2],
           [ 2,  0,  0, -4,  1,  1],
           [ 0,  2,  0,  1, -4,  1],
           [ 0,  0,  2,  1,  1, -4]], dtype=int8)
>>> np.array_equal(lap.matmat(np.eye(n)), lap.toarray())
True
>>> np.array_equal(lap.tosparse().toarray(), lap.toarray())
True
>>> lap.eigenvalues()
array([-7., -7., -4., -3., -3.,  0.])
>>> eigvals = eigvalsh(lap.toarray().astype(np.float64))
>>> np.allclose(lap.eigenvalues(), eigvals)
True
>>> np.allclose(lap.toarray() @ lap.eigenvectors(),
...             lap.eigenvectors() @ np.diag(lap.eigenvalues()))
True

以及 'neumann'

>>> lap = LaplacianNd(grid_shape, boundary_conditions='neumann')
>>> lap.tosparse()
<Compressed Sparse Row sparse array of dtype 'int8'
    with 20 stored elements and shape (6, 6)>
>>> lap.toarray()
array([[-2,  1,  0,  1,  0,  0],
       [ 1, -3,  1,  0,  1,  0],
       [ 0,  1, -2,  0,  0,  1],
       [ 1,  0,  0, -2,  1,  0],
       [ 0,  1,  0,  1, -3,  1],
       [ 0,  0,  1,  0,  1, -2]])
>>> np.array_equal(lap.matmat(np.eye(n)), lap.toarray())
True
>>> np.array_equal(lap.tosparse().toarray(), lap.toarray())
True
>>> lap.eigenvalues()
array([-5., -3., -3., -2., -1.,  0.])
>>> eigvals = eigvalsh(lap.toarray().astype(np.float64))
>>> np.allclose(lap.eigenvalues(), eigvals)
True
>>> np.allclose(lap.toarray() @ lap.eigenvectors(),
...             lap.eigenvectors() @ np.diag(lap.eigenvalues()))
True
属性::
H

埃尔米特伴随矩阵。

T

转置此线性算子。

方法

toarray()

从拉普拉斯算子数据构造一个密集数组。

tosparse()

从拉普拉斯算子数据构造一个稀疏数组。

eigenvalues(m=None)

构造一个包含 m 个最大(按绝对值排列最小)拉普拉斯矩阵特征值的 1D 数组,以升序排列。

eigenvectors(m=None)

构造一个数组,其列由对应于 m 个排序特征值的 m 个特征向量 (float) 组成,这些特征向量是 Nd 拉普拉斯算子的特征向量。

.. versionadded:: 1.12.0