Tensors

Type hierarchy

TensorKit.AbstractTensorMapType
abstract type AbstractTensorMap{S<:IndexSpace, N₁, N₂} end

Abstract supertype of all tensor maps, i.e. linear maps between tensor products of vector spaces of type S<:IndexSpace. An AbstractTensorMap maps from an input space of type ProductSpace{S, N₂} to an output space of type ProductSpace{S, N₁}.

source
TensorKit.AbstractTensorType
AbstractTensor{S<:IndexSpace, N} = AbstractTensorMap{S, N, 0}

Abstract supertype of all tensors, i.e. elements in the tensor product space of type ProductSpace{S, N}, built from elementary spaces of type S<:IndexSpace.

An AbstractTensor{S, N} is actually a special case AbstractTensorMap{S, N, 0}, i.e. a tensor map with only a non-trivial output space.

source
TensorKit.TensorMapType
struct TensorMap{S<:IndexSpace, N₁, N₂, ...} <: AbstractTensorMap{S, N₁, N₂}

Specific subtype of AbstractTensorMap for representing tensor maps (morphisms in a tensor category) whose data is stored in blocks of some subtype of DenseMatrix.

source

Specific TensorMap constructors

TensorKit.idFunction
id([A::Type{<:DenseMatrix} = Matrix{Float64},] space::VectorSpace) -> TensorMap

Construct the identity endomorphism on space space, i.e. return a t::TensorMap with domain(t) == codomain(t) == V, where storagetype(t) = A can be specified.

source
TensorKit.isomorphismFunction
isomorphism([A::Type{<:DenseMatrix} = Matrix{Float64},]
                cod::VectorSpace, dom::VectorSpace)
-> TensorMap

Return a t::TensorMap that implements a specific isomorphism between the codomain cod and the domain dom, and for which storagetype(t) can optionally be chosen to be of type A. If the two spaces do not allow for such an isomorphism, and are thus not isomorphic, and error will be thrown. When they are isomorphic, there is no canonical choice for a specific isomorphism, but the current choice is such that isomorphism(cod, dom) == inv(isomorphism(dom, cod)).

See also unitary when spacetype(cod) isa EuclideanSpace.

source
TensorKit.unitaryFunction
unitary([A::Type{<:DenseMatrix} = Matrix{Float64},] cod::VectorSpace, dom::VectorSpace)
-> TensorMap

Return a t::TensorMap that implements a specific unitary isomorphism between the codomain cod and the domain dom, for which spacetype(dom) (== spacetype(cod)) must be a subtype of EuclideanSpace. Furthermore, storagetype(t) can optionally be chosen to be of type A. If the two spaces do not allow for such an isomorphism, and are thus not isomorphic, and error will be thrown. When they are isomorphic, there is no canonical choice for a specific isomorphism, but the current choice is such that unitary(cod, dom) == inv(unitary(dom, cod)) = adjoint(unitary(dom, cod)).

source
TensorKit.isometryFunction
isometry([A::Type{<:DenseMatrix} = Matrix{Float64},] cod::VectorSpace, dom::VectorSpace)
-> TensorMap

Return a t::TensorMap that implements a specific isometry that embeds the domain dom into the codomain cod, and which requires that spacetype(dom) (== spacetype(cod)) is a subtype of EuclideanSpace. An isometry t is such that its adjoint t' is the left inverse of t, i.e. t'*t = id(dom), while t*t' is some idempotent endomorphism of cod, i.e. it squares to itself. When dom and cod do not allow for such an isometric inclusion, an error will be thrown.

source

TensorMap operations

TensorKit.permuteMethod
permute(tsrc::AbstractTensorMap{S}, p1::NTuple{N₁, Int}, p2::NTuple{N₂, Int} = ())
    -> tdst::TensorMap{S, N₁, N₂}

Permute the indices of tsrc::AbstractTensorMap{S} such that a new tensor tdst::TensorMap{S, N₁, N₂} is obtained, with indices in p1 playing the role of the codomain or range of the map, and indices in p2 indicating the domain.

To permute into an existing tdst, see add!

source
Base.permute!Function
permute!(tdst::AbstractTensorMap{S, N₁, N₂}, tsrc::AbstractTensorMap{S},
    p1::IndexTuple{N₁}, p2::IndexTuple{N₂}=()
) -> tdst

Permute the indices of tsrc and write the result into tdst, with indices in p1 playing the role of the codomain or range of the map tdst, and indices in p2 indicating the domain of tdst.

For more details see add!, of which this is a special case.

source
TensorKit.braidFunction
braid(f::FusionTree{<:Sector, N}, levels::NTuple{N, Int}, p::NTuple{N, Int})
-> <:AbstractDict{typeof(t), <:Number}

Perform a braiding of the uncoupled indices of the fusion tree f and return the result as a <:AbstractDict of output trees and corresponding coefficients. The braiding is determined by specifying that the new sector at position k corresponds to the sector that was originally at the position i = p[k], and assigning to every index i of the original fusion tree a distinct level or depth levels[i]. This permutation is then decomposed into elementary swaps between neighbouring indices, where the swaps are applied as braids such that if i and j cross, $τ_{i,j}$ is applied if levels[i] < levels[j] and $τ_{j,i}^{-1}$ if levels[i] > levels[j]. This does not allow to encode the most general braid, but a general braid can be obtained by combining such operations.

source
braid(f1::FusionTree{I}, f2::FusionTree{I},
        levels1::IndexTuple, levels2::IndexTuple,
        p1::IndexTuple{N₁}, p2::IndexTuple{N₂}) where {I<:Sector, N₁, N₂}
-> <:AbstractDict{Tuple{FusionTree{I, N₁}, FusionTree{I, N₂}}, <:Number}

Input is a fusion-splitting tree pair that describes the fusion of a set of incoming uncoupled sectors to a set of outgoing uncoupled sectors, represented using the splitting tree f1 and fusion tree f2, such that the incoming sectors f2.uncoupled are fused to f1.coupled == f2.coupled and then to the outgoing sectors f1.uncoupled. Compute new trees and corresponding coefficients obtained from repartitioning and braiding the tree such that sectors p1 become outgoing and sectors p2 become incoming. The uncoupled indices in splitting tree f1 and fusion tree f2 have levels (or depths) levels1 and levels2 respectively, which determines how indices braid. In particular, if i and j cross, $τ_{i,j}$ is applied if levels[i] < levels[j] and $τ_{j,i}^{-1}$ if levels[i] > levels[j]. This does not allow to encode the most general braid, but a general braid can be obtained by combining such operations.

source
Missing docstring.

Missing docstring for braid!. Check Documenter's build log for details.

Missing docstring.

Missing docstring for twist. Check Documenter's build log for details.

Missing docstring.

Missing docstring for twist!. Check Documenter's build log for details.

Missing docstring.

Missing docstring for add!. Check Documenter's build log for details.

Missing docstring.

Missing docstring for trace!. Check Documenter's build log for details.

Missing docstring.

Missing docstring for contract!. Check Documenter's build log for details.

TensorMap factorizations

TensorKit.leftorthFunction
leftorth(t::AbstractTensorMap, leftind::Tuple, rightind::Tuple;
            alg::OrthogonalFactorizationAlgorithm = QRpos()) -> Q, R

Create orthonormal basis Q for indices in leftind, and remainder R such that permute(t, leftind, rightind) = Q*R.

If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using leftorth!(t, alg = QRpos()).

Different algorithms are available, namely QR(), QRpos(), SVD() and Polar(). QR() and QRpos() use a standard QR decomposition, producing an upper triangular matrix R. Polar() produces a Hermitian and positive semidefinite R. QRpos() corrects the standard QR decomposition such that the diagonal elements of R are positive. Only QRpos() and Polar() are uniqe (no residual freedom) so that they always return the same result for the same input tensor t.

Orthogonality requires spacetype(t)<:InnerProductSpace, and leftorth(!) is currently only implemented for spacetype(t)<:EuclideanSpace.

source
TensorKit.rightorthFunction
rightorth(t::AbstractTensorMap, leftind::Tuple, rightind::Tuple;
            alg::OrthogonalFactorizationAlgorithm = LQpos()) -> L, Q

Create orthonormal basis Q for indices in rightind, and remainder L such that permute(t, leftind, rightind) = L*Q.

If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using rightorth!(t, alg = LQpos()).

Different algorithms are available, namely LQ(), LQpos(), RQ(), RQpos(), SVD() and Polar(). LQ() and LQpos() produce a lower triangular matrix L and are computed using a QR decomposition of the transpose. RQ() and RQpos() produce an upper triangular remainder L and only works if the total left dimension is smaller than or equal to the total right dimension. LQpos() and RQpos() add an additional correction such that the diagonal elements of L are positive. Polar() produces a Hermitian and positive semidefinite L. Only LQpos(), RQpos() and Polar() are uniqe (no residual freedom) so that they always return the same result for the same input tensor t.

Orthogonality requires spacetype(t)<:InnerProductSpace, and rightorth(!) is currently only implemented for spacetype(t)<:EuclideanSpace.

source
TensorKit.leftnullFunction
leftnull(t::AbstractTensor, leftind::Tuple, rightind::Tuple;
            alg::OrthogonalFactorizationAlgorithm = QRpos()) -> N

Create orthonormal basis for the orthogonal complement of the support of the indices in leftind, such that N' * permute(t, leftind, rightind) = 0.

If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using leftnull!(t, alg = QRpos()).

Different algorithms are available, namely QR() (or equivalently, QRpos()), SVD() and SDD(). The first assumes that the matrix is full rank and requires iszero(atol) and iszero(rtol). With SVD() and SDD(), rightnull will use the corresponding singular value decomposition, and one can specify an absolute or relative tolerance for which singular values are to be considered zero, where max(atol, norm(t)*rtol) is used as upper bound.

Orthogonality requires spacetype(t)<:InnerProductSpace, and leftnull(!) is currently only implemented for spacetype(t)<:EuclideanSpace.

source
TensorKit.rightnullFunction
rightnull(t::AbstractTensor, leftind::Tuple, rightind::Tuple;
            alg::OrthogonalFactorizationAlgorithm = LQ(),
            atol::Real = 0.0,
            rtol::Real = eps(real(float(one(eltype(t)))))*iszero(atol)) -> N

Create orthonormal basis for the orthogonal complement of the support of the indices in rightind, such that permute(t, leftind, rightind)*N' = 0.

If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using rightnull!(t, alg = LQpos()).

Different algorithms are available, namely LQ() (or equivalently, LQpos), SVD() and SDD(). The first assumes that the matrix is full rank and requires iszero(atol) and iszero(rtol). With SVD() and SDD(), rightnull will use the corresponding singular value decomposition, and one can specify an absolute or relative tolerance for which singular values are to be considered zero, where max(atol, norm(t)*rtol) is used as upper bound.

Orthogonality requires spacetype(t)<:InnerProductSpace, and rightnull(!) is currently only implemented for spacetype(t)<:EuclideanSpace.

source
TensorKit.tsvdFunction
tsvd(t::AbstractTensorMap, leftind::Tuple, rightind::Tuple;
    trunc::TruncationScheme = notrunc(), p::Real = 2, alg::Union{SVD, SDD} = SDD())
    -> U, S, V, ϵ

Compute the (possibly truncated)) singular value decomposition such that norm(permute(t, leftind, rightind) - U * S *V) ≈ ϵ, where ϵ thus represents the truncation error.

If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using tsvd!(t, trunc = notrunc(), p = 2).

A truncation parameter trunc can be specified for the new internal dimension, in which case a truncated singular value decomposition will be computed. Choices are:

  • notrunc(): no truncation (default);
  • truncerr(η::Real): truncates such that the p-norm of the truncated singular values is smaller than η times the p-norm of all singular values;
  • truncdim(χ::Int): truncates such that the equivalent total dimension of the internal vector space is no larger than χ;
  • truncspace(V): truncates such that the dimension of the internal vector space is smaller than that of V in any sector.
  • trunbelow(χ::Real): truncates such that every singular value is larger then χ ;

The method tsvd also returns the truncation error ϵ, computed as the p norm of the singular values that were truncated.

THe keyword alg can be equal to SVD() or SDD(), corresponding to the underlying LAPACK algorithm that computes the decomposition (_gesvd or _gesdd).

Orthogonality requires spacetype(t)<:InnerProductSpace, and svd(!) is currently only implemented for spacetype(t)<:EuclideanSpace.

source
TensorKit.eighFunction
eigh(t::AbstractTensorMap{<:EuclideanSpace}, leftind::Tuple, rightind::Tuple) -> D, V

Compute eigenvalue factorization of tensor t as linear map from rightind to leftind. The function eigh assumes that the linear map is hermitian and D and V tensors with the same eltype as t. See eig and eigen for non-hermitian tensors. Hermiticity requires that the tensor acts on inner product spaces, and the current implementation requires spacetyp(t) <: EuclideanSpace.

If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using eigh!(t). Note that the permuted tensor on which eigh! is called should have equal domain and codomain, as otherwise the eigenvalue decomposition is meaningless and cannot satisfy

permute(t, leftind, rightind) * V = V * D

See also eigen and eig.

source
TensorKit.eigFunction
eig(t::AbstractTensor, leftind::Tuple, rightind::Tuple; kwargs...) -> D, V

Compute eigenvalue factorization of tensor t as linear map from rightind to leftind. The function eig assumes that the linear map is not hermitian and returns type stable complex valued D and V tensors for both real and complex valued t. See eigh for hermitian linear maps

If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using eig!(t). Note that the permuted tensor on which eig! is called should have equal domain and codomain, as otherwise the eigenvalue decomposition is meaningless and cannot satisfy

permute(t, leftind, rightind) * V = V * D

Accepts the same keyword arguments scale, permute and sortby as eigen of dense matrices. See the corresponding documentation for more information.

See also eigen and eigh.

source
LinearAlgebra.eigenFunction
eigen(t::AbstractTensor, leftind::Tuple, rightind::Tuple; kwargs...) -> D, V

Compute eigenvalue factorization of tensor t as linear map from rightind to leftind.

If leftind and rightind are not specified, the current partition of left and right indices of t is used. In that case, less memory is allocated if one allows the data in t to be destroyed/overwritten, by using eigen!(t). Note that the permuted tensor on which eigen! is called should have equal domain and codomain, as otherwise the eigenvalue decomposition is meaningless and cannot satisfy

permute(t, leftind, rightind) * V = V * D

Accepts the same keyword arguments scale, permute and sortby as eigen of dense matrices. See the corresponding documentation for more information.

See also eig and eigh

source