scipy.spatial.distance.squareform#

scipy.spatial.distance.squareform(X, force='no', checks=True)[source]#

Convert a vector-form distance vector to a square-form distance matrix, and vice-versa.

Parameters:
Xarray_like

Either a condensed or redundant distance matrix.

forcestr, optional

As with MATLAB(TM), if force is equal to 'tovector' or 'tomatrix', the input will be treated as a distance matrix or distance vector respectively.

checksbool, optional

If set to False, no checks will be made for matrix symmetry nor zero diagonals. This is useful if it is known that X - X.T1 is small and diag(X) is close to zero. These values are ignored any way so they do not disrupt the squareform transformation.

Returns:
Yndarray

If a condensed distance matrix is passed, a redundant one is returned, or if a redundant one is passed, a condensed distance matrix is returned.

Notes

  1. v = squareform(X)

    Given a square n-by-n symmetric distance matrix X, v = squareform(X) returns a n * (n-1) / 2 (i.e. binomial coefficient n choose 2) sized vector v where \(v[{n \choose 2} - {n-i \choose 2} + (j-i-1)]\) is the distance between distinct points i and j. If X is non-square or asymmetric, an error is raised.

  2. X = squareform(v)

    Given a n * (n-1) / 2 sized vector v for some integer n >= 1 encoding distances as described, X = squareform(v) returns a n-by-n distance matrix X. The X[i, j] and X[j, i] values are set to \(v[{n \choose 2} - {n-i \choose 2} + (j-i-1)]\) and all diagonal elements are zero.

In SciPy 0.19.0, squareform stopped casting all input types to float64, and started returning arrays of the same dtype as the input.

Examples

>>> import numpy as np
>>> from scipy.spatial.distance import pdist, squareform

x is an array of five points in three-dimensional space.

>>> x = np.array([[2, 0, 2], [2, 2, 3], [-2, 4, 5], [0, 1, 9], [2, 2, 4]])

pdist(x) computes the Euclidean distances between each pair of points in x. The distances are returned in a one-dimensional array with length 5*(5 - 1)/2 = 10.

>>> distvec = pdist(x)
>>> distvec
array([2.23606798, 6.40312424, 7.34846923, 2.82842712, 4.89897949,
       6.40312424, 1.        , 5.38516481, 4.58257569, 5.47722558])

squareform(distvec) returns the 5x5 distance matrix.

>>> m = squareform(distvec)
>>> m
array([[0.        , 2.23606798, 6.40312424, 7.34846923, 2.82842712],
       [2.23606798, 0.        , 4.89897949, 6.40312424, 1.        ],
       [6.40312424, 4.89897949, 0.        , 5.38516481, 4.58257569],
       [7.34846923, 6.40312424, 5.38516481, 0.        , 5.47722558],
       [2.82842712, 1.        , 4.58257569, 5.47722558, 0.        ]])

When given a square distance matrix m, squareform(m) returns the one-dimensional condensed distance vector associated with the matrix. In this case, we recover distvec.

>>> squareform(m)
array([2.23606798, 6.40312424, 7.34846923, 2.82842712, 4.89897949,
       6.40312424, 1.        , 5.38516481, 4.58257569, 5.47722558])