Skip to content

decomp

High-level tensor decomposition with multiple output modes.

decomp

decomp(
    T: Tensor,
    axes: Union[int, str, Sequence[Union[int, str]]],
    mode: str = "SVD",
    flow: str = "><",
    itag: Optional[Union[str, Tuple[str, str]]] = None,
    trunc: Optional[Dict[str, Union[int, float]]] = None,
) -> Union[Tuple[Tensor, Tensor], Tuple[Tensor, Tensor, Tensor]]

Perform tensor decomposition with flexible output modes.

Parameters:

Name Type Description Default
T Tensor

Tensor to be decomposed.

required
axes Union[int, str, Sequence[Union[int, str]]]

Index or indices to separate from all others. Can be:

  • Single integer position or string tag
  • Sequence of integer positions or string tags (merges multiple axes first)

required
mode str

Decomposition mode:

  • "SVD": Returns (U, S, Vh) where S is diagonal matrix tensor (full SVD)
  • "UR": Returns (U, R) where R = S*Vh (singular values multiplied into Vh)
  • "LV": Returns (L, V) where L = U*S (singular values multiplied into U)
  • "QR": Returns (Q, R) where Q is orthogonal and R is upper triangular

'SVD'
flow str

Arrow direction control. Default is "><" (both arrows incoming).

  • For SVD mode: Controls S matrix arrow directions ("><", ">>", or "<<")
  • For UR mode: Both ">>" and "><" normalize to ">>" (outward bonds); "<<" is also accepted
  • For LV mode: Both "<<" and "><" normalize to "<<" (inward bonds); ">>" is also accepted
  • For QR mode: Controls bond arrow directions
Note: The underlying svd naturally produces ">>" or "<<" depending on left_index.direction. This parameter uses capcup() to adjust from the natural flow to the desired flow.

'><'
itag Optional[Union[str, Tuple[str, str]]]

Index tag(s) for the bond dimension(s). Can be:

  • None: Use default tags "_bond_L" and "_bond_R"
  • str: Use same tag for both left and right bonds
  • tuple[str, str]: Use (left_tag, right_tag) for left and right bonds respectively

None
trunc Optional[Dict[str, Union[int, float]]]

Truncation specification as a dict. If None, no truncation. Supported keys:

  • "nkeep": Keep at most n singular values globally
  • "thresh": Keep singular values >= t per block
Both can be specified together: thresh is applied first, then nkeep.

None

Returns:

Type Description
tuple[Tensor, Tensor] or tuple[Tensor, Tensor, Tensor]
  • "SVD" mode: (U, S, Vh) with S as diagonal matrix tensor
  • "UR" mode: (U, R) where R incorporates singular values
  • "LV" mode: (L, V) where L incorporates singular values
  • "QR" mode: (Q, R) where Q is orthogonal and R is upper triangular

Raises:

Type Description
ValueError

If mode is not one of "SVD", "UR", "LV", or "QR", or if flow is not ">>", "<<", or "><"

Examples:

>>> from nicole import Tensor, decomp, U1Group, Direction, Index, Sector
>>> 
>>> # Create a sample tensor
>>> group = U1Group()
>>> idx1 = Index(Direction.OUT, group, sectors=(Sector(0, 2), Sector(1, 3)))
>>> idx2 = Index(Direction.IN, group, sectors=(Sector(0, 3),))
>>> T = Tensor.random([idx1, idx2], itags=["a", "b"], seed=1)
>>> 
>>> # UR mode: Get U and R=S*Vh (most efficient for reconstruction)
>>> U, R = decomp(T, axes=0, mode="UR")
>>> 
>>> # SVD mode: Get full SVD with diagonal S
>>> U, S, Vh = decomp(T, axes=0, mode="SVD")
>>> 
>>> # LV mode: Get L=U*S and V
>>> L, V = decomp(T, axes=0, mode="LV")
Notes
  • UR and LV modes are more memory and computationally efficient than SVD mode
  • SVD mode constructs a full diagonal matrix tensor for S
  • All modes produce mathematically equivalent decompositions
  • When multiple axes are specified, they are first merged using an n-to-1 isometry, decomposed, and then the U tensor is unmerged back to the original axes

Description

Performs symmetry-preserving decomposition with four output modes:

Modes

  • "SVD": Returns (U, S, Vh) with S as diagonal tensor
  • "UR": Returns (U, R) where R = S @ Vh (most efficient for reconstruction)
  • "LV": Returns (L, V) where L = U @ S
  • "QR": Returns (Q, R) where Q is orthogonal and R is upper triangular (no truncation)

Each symmetry sector is decomposed independently.

Parameters

axes

Index or indices to separate from all others. Can be:

  • Single integer position or string tag
  • Sequence of integer positions or string tags (merges multiple axes first)

mode

Decomposition mode: "SVD", "UR", "LV", or "QR" (default: "SVD")

  • "SVD", "UR", "LV": SVD-based decomposition with optional truncation
  • "QR": QR decomposition (orthogonal factorization, no truncation)

flow

Arrow direction control for bond indices (default: "><"):

  • "><": Both bond arrows incoming — arrows converge from U and Vh into S; S has (IN, IN) (default)
  • ">>": Arrow chain flows left to right — S has (IN, OUT)
  • "<<": Arrow chain flows right to left — S has (OUT, IN)

itag

Index tag(s) for the bond dimension(s):

  • None: Use default tags "_bond_L" and "_bond_R"
  • str: Use same tag for both left and right bonds
  • tuple[str, str]: Use (left_tag, right_tag) for left and right bonds

trunc

Optional truncation with trunc parameter (dict) - only for SVD-based modes:

  • "nkeep": Keep n largest singular values globally
  • "thresh": Keep singular values ≥ t per block
  • Both can be specified together (thresh applied first, then nkeep)
  • Not applicable for "QR" mode (QR decomposition doesn't support truncation)

Usage Examples

from nicole import decomp

# UR mode (most efficient)
U, R = decomp(T, axes=0, mode="UR")

# SVD mode with custom bond tags
U, S, Vh = decomp(T, axes=0, mode="SVD", itag=("left", "right"))

# LV mode with truncation
L, V = decomp(T, axes=0, mode="LV", trunc={"nkeep": 100})

# QR mode (orthogonal factorization, no truncation)
Q, R = decomp(T, axes=0, mode="QR")

# Decompose merging multiple axes
U, R = decomp(T, axes=[0, 1, 2], mode="UR")

See Also

Notes

  • UR and LV modes are more efficient than SVD mode for reconstruction
  • QR mode provides orthogonal factorization without truncation
  • Bond dimension after truncation may be smaller (SVD-based modes only)
  • Charge sectors with zero singular values are automatically removed (SVD-based modes)
  • When multiple axes specified, they are merged first using n-to-1 isometry
  • QR mode is useful for obtaining canonical forms without compression