cuquantum.cutensornet.tensor.decompose

cuquantum.cutensornet.tensor.decompose(subscripts, operand, *, method=None, options=None, stream=None, return_info=False)[source]

Perform the tensor decomposition of the operand based on the expression described by subscripts.

The expression adopts a similar notation as Einstein summation (einsum) expression, but in the reversed order. Meanwhile the input is now a single operand while the output contains two or three tensors. Unlike einsum expression, the mode labels for all input and output operands must be specified explicitly.

See the notes and examples for clarification. See also Tensor Decomposition.

Parameters
  • subscripts – The mode labels (subscripts) defining the decomposition as a comma-separated sequence of characters. Unicode characters are allowed in the expression thereby expanding the rank (dimensionality) of the tensor that can be specified.

  • operand – A ndarray-like tensor object. The currently supported types are numpy.ndarray, cupy.ndarray, and torch.Tensor.

  • method – Specify decomposition method as a cuquantum.cutensornet.tensor.QRMethod or a cuquantum.cutensornet.tensor.SVDMethod object. Alternatively, a dict containing the parameters for the QRMethod or SVDMethod constructor can also be provided. If not specified, the value will be set to the default-constructed QRMethod. Note that both SVD and QR method operate in reduced fashion, similar to full_matrices=False for numpy.linalg.svd and reduced=True for numpy.linalg.qr.

  • options – Specify the computational options for the decomposition as a cuquantum.cutensornet.tensor.DecompositionOptions object. Alternatively, a dict containing the parameters for the DecompositionOptions constructor can also be provided. If not specified, the value will be set to the default-constructed DecompositionOptions.

  • stream – Provide the CUDA stream to use for the decomposition. Acceptable inputs include cudaStream_t (as Python int), cupy.cuda.Stream, and torch.cuda.Stream. If a stream is not provided, the current stream will be used.

  • return_info – If true, information about the decomposition will be returned via a cuquantum.cutensornet.tensor.SVDInfo object. Currently this option is only supported for SVD decomposition (which is specified via method).

Returns

  • For QR decomposition (default), the output tensors Q and R (ndarray-like objects) of the same type and on the same device as the input operand are returned as the result of the decomposition.

  • For SVD decomposition, if return_info is False, a 3-tuple of output tensors U, S and V (ndarray-like objects) of the same type as the input operand are returned as the result of the decomposition. If return_info is True, a 4-tuple of output tensors U, S, V and a dict object that contains information about the decomposition will be returned. Note, depending on the choice of cuquantum.cutensornet.tensor.SVDMethod.partition, the returned S operand may be None. Also see partition.

Return type

Depending on the decomposition method specified in method, the results returned may vary

The decomposition expression adopts a similar notation as einsum expression. The subscripts string is a list of subscript labels where each label refers to a mode of the corresponding operand. The subscript labels are separated by either comma or identifier ->. The subscript labels before the identifier -> are viewed as input, and the ones after are viewed as outputs, respectively. The requirements on the subscripts for SVD and QR decomposition are summarized below:

  • For SVD and QR decomposition, the subscripts string is expected to contain exactly one input and two outputs (the modes for s is not needed in the case of SVD).

  • One and only one identical mode is expected to exist in the two output mode labels.

  • When inverted, the decomposition subscript yields a valid einsum subscript that can specify the contraction of the outputs to reproduce the input (modes for s excluded for SVD).

Examples

>>> # equivalent:
>>> # q, r = numpy.linalg.qr(a)
>>> q, r = tensor.decompose('ij->ik,kj', a)
>>> # equivalent:
>>> # u, s, v = numpy.linalg.svd(a, full_matrices=False)
>>> u, s, v = tensor.decompose('ij->ik,kj', a, method=tensor.SVDMethod())

For generalization to multi-dimensional tensors (here a is a rank-4 tensor):

>>> u, s, v = tensor.decompose('ijab->ixb,jax', a, method=tensor.SVDMethod())
>>> # u is unitary
>>> identity = cuquantum.contract('ixb,iyb->xy', u, u.conj())
>>> # re-construct the tensor a by inverting the expression
>>> a_reconstructed = cuquantum.contract('ixb,x,jax->ijab', u, s, v)

Broadcasting is supported for certain cases via ellipsis notation. One may add an ellipsis in the input to represent all the modes that are not explicitly specified in the labels. The ellipsis must also appear in one of the outputs to indicate which output the represented modes will all be partitioned onto.

Example

Given a rank-6 operand a:

>>> # equivalent:
>>> # q, r = tensor.decompose('ijabcd->ixj,abcdx', a)
>>> q, r = tensor.decompose('ij...->ixj,...x', a)

Note

It is encouraged for users to maintain the library handle themselves so as to reduce the context initialization time:

from cuquantum import cutensornet as cutn
from cuquantum.cutensornet.tensor import decompose, QRMethod

handle = cutn.create()
q, r = decompose(..., method=QRMethod(), options={"handle": handle}, ...)
# ... the same handle can be reused for further calls ...
# when it's done, remember to destroy the handle
cutn.destroy(handle)

Below we give more pedagogical examples.

Examples

Use NumPy operands:

>>> from cuquantum.cutensornet.tensor import decompose, SVDMethod
>>> import numpy as np
>>> T = np.random.random((4,4,6,6))

Perform tensor QR decomposition such that T[i,j,a,b] = Q[i,k,a] R[k,j,b]. The results q and r are NumPy ndarrays (with the computation performed on the GPU):

>>> q, r = decompose('ijab->ika,kjb', T)

Perform exact tensor SVD decomposition such that T[i,j,a,b] = U[i,k,a] S[k] V[k,j,b]:

>>> u, s, v = decompose('ijab->ika,kjb', T, method=SVDMethod())

Perform exact tensor SVD decomposition such that T[i,j,a,b] = US[i,k,a] V[k,j,b] where US[i,k,a] represents the product of U[i,k,a] and S[k]:

>>> us, _, v = decompose('ijab->ika,kjb', T, method=SVDMethod(partition="U"))

Perform exact tensor SVD decomposition such that T[i,j,a,b] = U[i,k,a] S[k] V[k,j,b] then normalize the L2 norm of output singular values to 1:

>>> u, s_normalized, v = decompose('ijab->ika,kjb', T, method=SVDMethod(normalization="L2") )
>>> print(np.linalg.norm(s_normalized)) # 1.0

Perform truncated SVD decomposition to keep at most 8 singular values. Meanwhile, request the runtime information from the SVD truncation with return_info=True.

>>> u, s, v, info = decompose('ijab->ika,kjb', T, method=SVDMethod(max_extent=8), return_info=True)
>>> print(s.shape) # (8,)
>>> print(info)

We can also perform truncated SVD decomposition with all requirements below:

  • the number of remaining singular values shall not exceed 8.

  • remaining singular values are all larger than 0.01

  • remaining singular values are all larger than 0.1 * largest singular values

  • the L1 norm of the remaining singular values are normalized to 1.

  • the remaining singular values (after truncation and normalization) are equally partitioned onto U and V

>>> method = {"max_extent": 8,
...           "abs_cutoff": 0.01,
...           "rel_cutoff": 0.1,
...           "normalization": "L1",
...           "partition": "UV"}
>>> svd_method = SVDMethod(**method)
>>> us, _, sv, info = decompose('ijab->ika,kjb', T, method=svd_method, return_info=True)

Alternatively, the options can be provided as a dict object:

>>> us, _, sv, info = tensor_decompose('ijab->ika,kjb', T, method=method, return_info=True)

Use CuPy operands. The results q and r are CuPy ndarrays on the same device as the input operand, and dev is any valid device ID on your system that you wish to use to store the tensors and perform the decomposition:

>>> import cupy
>>> dev = 0
>>> with cupy.cuda.Device(dev):
...     T = cupy.random.random((4,4,6,6))
>>> q, r = decompose('ijab->ika,kjb', T)

Use PyTorch operands. The results q and r are PyTorch tensors on the same device (dev) as the input operand:

>>> import torch
>>> dev = 0
>>> T = torch.rand(4,4,6,6, device=f'cuda:{dev}')
>>> q, r = decompose('ijab->ika,kjb', T)