quimb.tensor.array_ops

Backend agnostic array operations.

Classes

PArray

Simple array-like object that lazily generates the actual array by

Functions

asarray(array)

Maybe convert data for a tensor to use. If array already has a

calc_fuse_perm_and_shape(shape, axes_groups)

fuse(x, *axes_groups[, backend])

Fuse the give group or groups of axes. The new fused axes will be

ndim(array)

The number of dimensions of an array.

multiply_diagonal(x, v, axis[, backend])

Multiply v into x as if contracting in a diagonal matrix.

align_axes(*arrays, axes[, backend])

Prepare a set of arrays that should be contractible along axes.

iscomplex(x)

Does x have a complex dtype?

norm_fro(x)

The frobenius norm of an array.

sensibly_scale(x)

Take an array and scale it very roughly such that random tensor

_numba_find_diag_axes(x[, atol])

Numba-compiled array diagonal axis finder.

find_diag_axes(x[, atol])

Try and find a pair of axes of x in which it is diagonal.

_numba_find_antidiag_axes(x[, atol])

Numba-compiled array antidiagonal axis finder.

find_antidiag_axes(x[, atol])

Try and find a pair of axes of x in which it is anti-diagonal.

_numba_find_columns(x[, atol])

Numba-compiled single non-zero column axis finder.

find_columns(x[, atol])

Try and find columns of axes which are zero apart from a single index.

Module Contents

quimb.tensor.array_ops.asarray(array)[source]

Maybe convert data for a tensor to use. If array already has a .shape attribute, i.e. looks like an array, it is left as-is. Else the elements are inspected to see which libraries’ array constructor should be used, defaulting to numpy if everything is builtin or numpy numbers.

quimb.tensor.array_ops.calc_fuse_perm_and_shape(shape, axes_groups)[source]
quimb.tensor.array_ops.fuse(x, *axes_groups, backend=None)[source]

Fuse the give group or groups of axes. The new fused axes will be inserted at the minimum index of any fused axis (even if it is not in the first group). For example, fuse(x, [5, 3], [7, 2, 6]) will produce an array with axes like:

groups inserted at axis 2, removed beyond that.
        ......<--
(0, 1, g0, g1, 4, 8, ...)
        |   |
        |   g1=(7, 2, 6)
        g0=(5, 3)
Parameters:

axes_groups (sequence of sequences of int) – The axes to fuse. Each group of axes will be fused into a single axis.

quimb.tensor.array_ops.ndim(array)[source]

The number of dimensions of an array.

quimb.tensor.array_ops.multiply_diagonal(x, v, axis, backend=None)[source]

Multiply v into x as if contracting in a diagonal matrix.

quimb.tensor.array_ops.align_axes(*arrays, axes, backend=None)[source]

Prepare a set of arrays that should be contractible along axes.

For example, block symmetric arrays need to have aligned sectors prior to fusing.

quimb.tensor.array_ops.iscomplex(x)[source]

Does x have a complex dtype?

quimb.tensor.array_ops.norm_fro(x)[source]

The frobenius norm of an array.

quimb.tensor.array_ops.sensibly_scale(x)[source]

Take an array and scale it very roughly such that random tensor networks consisting of such arrays do not have gigantic norms.

quimb.tensor.array_ops._numba_find_diag_axes(x, atol=1e-12)[source]

Numba-compiled array diagonal axis finder.

Parameters:
  • x (numpy.ndarray) – The array to search for diagonal axes.

  • atol (float) – The tolerance with which to compare to zero.

Returns:

diag_axes – The set of pairs of axes which are diagonal.

Return type:

set[tuple[int]]

quimb.tensor.array_ops.find_diag_axes(x, atol=1e-12)[source]

Try and find a pair of axes of x in which it is diagonal.

Parameters:
  • x (array-like) – The array to search.

  • atol (float, optional) – Tolerance with which to compare to zero.

Returns:

The two axes if found else None.

Return type:

tuple[int] or None

Examples

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

Which means we can reduce x without loss of information to:

>>> np.einsum('abb->ab', x)
array([[1, 2],
       [3, 4]])
quimb.tensor.array_ops._numba_find_antidiag_axes(x, atol=1e-12)[source]

Numba-compiled array antidiagonal axis finder.

Parameters:
  • x (numpy.ndarray) – The array to search for anti-diagonal axes.

  • atol (float) – The tolerance with which to compare to zero.

Returns:

antidiag_axes – The set of pairs of axes which are anti-diagonal.

Return type:

set[tuple[int]]

quimb.tensor.array_ops.find_antidiag_axes(x, atol=1e-12)[source]

Try and find a pair of axes of x in which it is anti-diagonal.

Parameters:
  • x (array-like) – The array to search.

  • atol (float, optional) – Tolerance with which to compare to zero.

Returns:

The two axes if found else None.

Return type:

tuple[int] or None

Examples

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

Which means we can reduce x without loss of information to:

>>> np.einsum('aba->ab', x[::-1, :, :])
array([[3, 4],
       [1, 2]])

as long as we flip the order of dimensions on other tensors corresponding to the the same index.

quimb.tensor.array_ops._numba_find_columns(x, atol=1e-12)[source]

Numba-compiled single non-zero column axis finder.

Parameters:
  • x (array) – The array to search.

  • atol (float, optional) – Absolute tolerance to compare to zero with.

Returns:

Set of pairs (axis, index) defining lone non-zero columns.

Return type:

set[tuple[int]]

quimb.tensor.array_ops.find_columns(x, atol=1e-12)[source]

Try and find columns of axes which are zero apart from a single index.

Parameters:
  • x (array-like) – The array to search.

  • atol (float, optional) – Tolerance with which to compare to zero.

Returns:

If found, the first integer is which axis, and the second is which column of that axis, else None.

Return type:

tuple[int] or None

Examples

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

Which means we can happily slice x without loss of information to:

>>> x[:, :, 1]
array([[1, 2],
       [3, 4]])
class quimb.tensor.array_ops.PArray(fn, params, shape=None)[source]

Simple array-like object that lazily generates the actual array by calling a function with a set of parameters.

Parameters:
  • fn (callable) – The function that generates the tensor data from params.

  • params (sequence of numbers) – The initial parameters supplied to the generating function like fn(params).

See also

PTensor

__slots__ = ('_fn', '_params', '_data', '_shape', '_shape_fn_id')
property fn
property params
_shape
_shape_fn_id
copy()[source]
property data
property shape
property ndim
add_function(g)[source]

Chain the new function g on top of current function f like g(f(params)).