SciPy

scipy.linalg.fiedler

scipy.linalg.fiedler(a)[source]

Returns a symmetric Fiedler matrix

Given an sequence of numbers a, Fiedler matrices have the structure F[i, j] = np.abs(a[i] - a[j]), and hence zero diagonals and nonnegative entries. A Fiedler matrix has a dominant positive eigenvalue and other eigenvalues are negative. Although not valid generally, for certain inputs, the inverse and the determinant can be derived explicitly as given in [1].

Parameters
a(n,) array_like

coefficient array

Returns
F(n, n) ndarray

See also

circulant, toeplitz

Notes

New in version 1.3.0.

References

1(1,2)

J. Todd, “Basic Numerical Mathematics: Vol.2 : Numerical Algebra”, 1977, Birkhauser, DOI:10.1007/978-3-0348-7286-7

Examples

>>> from scipy.linalg import det, inv, fiedler
>>> a = [1, 4, 12, 45, 77]
>>> n = len(a)
>>> A = fiedler(a)
>>> A
array([[ 0,  3, 11, 44, 76],
       [ 3,  0,  8, 41, 73],
       [11,  8,  0, 33, 65],
       [44, 41, 33,  0, 32],
       [76, 73, 65, 32,  0]])

The explicit formulas for determinant and inverse seem to hold only for monotonically increasing/decreasing arrays. Note the tridiagonal structure and the corners.

>>> Ai = inv(A)
>>> Ai[np.abs(Ai) < 1e-12] = 0.  # cleanup the numerical noise for display
>>> Ai
array([[-0.16008772,  0.16666667,  0.        ,  0.        ,  0.00657895],
       [ 0.16666667, -0.22916667,  0.0625    ,  0.        ,  0.        ],
       [ 0.        ,  0.0625    , -0.07765152,  0.01515152,  0.        ],
       [ 0.        ,  0.        ,  0.01515152, -0.03077652,  0.015625  ],
       [ 0.00657895,  0.        ,  0.        ,  0.015625  , -0.00904605]])
>>> det(A)
15409151.999999998
>>> (-1)**(n-1) * 2**(n-2) * np.diff(a).prod() * (a[-1] - a[0])
15409152

Previous topic

scipy.linalg.dft

Next topic

scipy.linalg.fiedler_companion