Tensors

Type hierarchy

The type hierarchy of tensors is as follows:

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.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
TensorKit.BraidingTensorType
struct BraidingTensor{S<:IndexSpace} <: AbstractTensorMap{S, 2, 2}
BraidingTensor(V1::S, V2::S, adjoint::Bool=false) where {S<:IndexSpace}

Specific subtype of AbstractTensorMap for representing the braiding tensor that braids the first input over the second input; its inverse can be obtained as the adjoint.

It holds that domain(BraidingTensor(V1, V2)) == V1 ⊗ V2 and codomain(BraidingTensor(V1, V2)) == V2 ⊗ V1.

source

Some aliases are provided for convenience:

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.TensorType
Tensor{S, N, I, A, F₁, F₂} = TensorMap{S, N, 0, I, A, F₁, F₂}

Specific subtype of AbstractTensor for representing tensors whose data is stored in blocks of some subtype of DenseMatrix.

A Tensor{S, N, I, A, F₁, F₂} is actually a special case TensorMap{S, N, 0, I, A, F₁, F₂}, i.e. a tensor map with only a non-trivial output space.

source
TensorKit.TrivialTensorMapType
TrivialTensorMap{S<:IndexSpace, N₁, N₂, A<:DenseMatrix} = TensorMap{S, N₁, N₂, Trivial, 
                                                                    A, Nothing, Nothing}

A special case of TensorMap for representing tensor maps with trivial symmetry, i.e., whose sectortype is Trivial.

source
TensorKit.TrivialTensorType
TrivialTensor{S, N, A} = TrivialTensorMap{S, N, 0, A}

A special case of Tensor for representing tensors with trivial symmetry, i.e., whose sectortype is Trivial.

source

TensorMap constructors

General constructors

A general TensorMap can be constructed by specifying its data, codmain and domain in one of the following ways:

TensorKit.TensorMapMethod
TensorMap(data::AbstractDict{<:Sector,<:DenseMatrix}, codomain::ProductSpace{S,N₁},
            domain::ProductSpace{S,N₂}) where {S<:ElementarySpace,N₁,N₂}
TensorMap(data, codomain ← domain)
TensorMap(data, domain → codomain)

Construct a TensorMap by explicitly specifying its block data.

Arguments

  • data::AbstractDict{<:Sector,<:DenseMatrix}: dictionary containing the block data for each coupled sector c as a DenseMatrix of size (blockdim(codomain, c), blockdim(domain, c)).
  • codomain::ProductSpace{S,N₁}: the codomain as a ProductSpace of N₁ spaces of type S<:ElementarySpace.
  • domain::ProductSpace{S,N₂}: the domain as a ProductSpace of N₂ spaces of type S<:ElementarySpace.

Alternatively, the domain and codomain can be specified by passing a HomSpace using the syntax codomain ← domain or domain → codomain.

source
TensorKit.TensorMapMethod
TensorMap([f, eltype,] codomain::ProductSpace{S,N₁}, domain::ProductSpace{S,N₂})
            where {S<:ElementarySpace,N₁,N₂}
TensorMap([f, eltype,], codomain ← domain)
TensorMap([f, eltype,], domain → codomain)

Construct a TensorMap from a general callable that produces block data for each coupled sector.

Arguments

  • f: callable object that returns a DenseMatrix, or UndefInitializer.
  • eltype::Type{<:Number}: element type of the data.
  • codomain::ProductSpace{S,N₁}: the codomain as a ProductSpace of N₁ spaces of type S<:ElementarySpace.
  • domain::ProductSpace{S,N₂}: the domain as a ProductSpace of N₂ spaces of type S<:ElementarySpace.

If eltype is left unspecified, f should support the calling syntax f(::Tuple{Int,Int}) such that f((m, n)) returns a DenseMatrix with size(f((m, n))) == (m, n). If eltype is specified, f is instead called as f(eltype, (m, n)). In the case where f is left unspecified or undef is passed explicitly, a TensorMap with uninitialized data is generated.

Alternatively, the domain and codomain can be specified by passing a HomSpace using the syntax codomain ← domain or domain → codomain.

source
TensorKit.TensorMapMethod
TensorMap(data::DenseArray, codomain::ProductSpace{S,N₁}, domain::ProductSpace{S,N₂};
                tol=sqrt(eps(real(float(eltype(data)))))) where {S<:ElementarySpace,N₁,N₂}
TensorMap(data, codomain ← domain)
TensorMap(data, domain → codomain)

Construct a TensorMap from a plain multidimensional array.

Arguments

  • data::DenseArray: tensor data as a plain array.
  • codomain::ProductSpace{S,N₁}: the codomain as a ProductSpace of N₁ spaces of type S<:ElementarySpace.
  • domain::ProductSpace{S,N₂}: the domain as a ProductSpace of N₂ spaces of type S<:ElementarySpace.
  • tol=sqrt(eps(real(float(eltype(data)))))::Float64:

Here, data can be specified in two ways. It can either be a DenseArray of rank N₁ + N₂ whose size matches that of the domain and codomain spaces, size(data) == (dims(codomain)..., dims(domain)...), or a DenseMatrix where size(data) == (dim(codomain), dim(domain)). The TensorMap constructor will then reconstruct the tensor data such that the resulting tensor t satisfies data == convert(Array, t). For the case where sectortype(S) == Trivial, the data array is simply reshaped into matrix form and referred to as such in the resulting TensorMap instance. When S<:GradedSpace, the data array has to be compatible with how how each sector in every space V is assigned to an index range within 1:dim(V).

Alternatively, the domain and codomain can be specified by passing a HomSpace using the syntax codomain ← domain or domain → codomain.

Note

This constructor only works for sectortype values for which conversion to a plain array is possible, and only in the case where the data actually respects the specified symmetry structure.

source

Several special-purpose methods exist to generate data according to specific distributions:

TensorKit.randisometryFunction
randisometry([::Type{T}=Float64], dims::Dims{2}) -> Array{T,2}
randhaar([::Type{T}=Float64], dims::Dims{2}) -> Array{T,2}

Create a random isometry of size dims, uniformly distributed according to the Haar measure.

See also randuniform and randnormal.

source

Specific constructors

Additionally, several special-purpose constructors exist to generate data according to specific distributions:

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 InnerProductStyle(cod) === EuclideanProduct().

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 have an inner product. 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)) has an Euclidean inner product. 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

Accessing properties and data

The following methods exist to obtain type information:

TensorKit.spacetypeFunction
spacetype(::Union{T,Type{T}}) where {T<:AbstractTensorMap} -> Type{S<:IndexSpace}

Return the type of the elementary space S of a tensor.

source
TensorKit.sectortypeMethod
sectortype(::Union{T,Type{T}}) where {T<:AbstractTensorMap} -> Type{I<:Sector}

Return the type of sector I of a tensor.

source
TensorKit.storagetypeFunction
storagetype(::Union{T,Type{T}}) where {T<:TensorMap} -> Type{A<:DenseMatrix}

Return the type of the storage A of the tensor map.

source
TensorKit.tensormaptypeFunction
tensormaptype(::Type{S}, N₁::Int, N₂::Int, [::Type{T}]) where {S<:IndexSpace,T} -> ::Type{<:TensorMap}

Return the fully specified type of a tensor map with elementary space S, N₁ output spaces and N₂ input spaces, either with scalar type T or with storage type T.

source

To obtain information about the indices, you can use:

TensorKit.domainFunction
domain(t::AbstractTensorMap{S,N₁,N₂}, [i::Int]) -> ProductSpace{S,N₂}

Return the domain of a tensor, i.e. the product space of the input spaces. If i is specified, return the i-th input space. Implementations should provide domain(t).

See also codomain and space.

source
TensorKit.codomainFunction
codomain(t::AbstractTensorMap{S,N₁,N₂}, [i::Int]) -> ProductSpace{S,N₁}

Return the codomain of a tensor, i.e. the product space of the output spaces. If i is specified, return the i-th output space. Implementations should provide codomain(t).

See also domain and space.

source
TensorKit.spaceMethod
space(t::AbstractTensorMap{S,N₁,N₂}, [i::Int]) -> HomSpace{S,N₁,N₂}

The index information of a tensor, i.e. the HomSpace of its domain and codomain. If i is specified, return the i-th index space.

source
TensorKit.numinFunction
numin(::Union{T,Type{T}}) where {T<:AbstractTensorMap} -> Int

Return the number of input spaces of a tensor. This is equivalent to the number of spaces in the domain of that tensor.

See also numout and numind.

source
TensorKit.numoutFunction
numout(::Union{T,Type{T}}) where {T<:AbstractTensorMap} -> Int

Return the number of output spaces of a tensor. This is equivalent to the number of spaces in the codomain of that tensor.

See also numin and numind.

source
TensorKit.numindFunction
numind(::Union{T,Type{T}}) where {T<:AbstractTensorMap} -> Int

Return the total number of input and output spaces of a tensor. This is equivalent to the total number of spaces in the domain and codomain of that tensor.

See also numout and numin.

source
TensorKit.allindFunction
allind(::Union{T,Type{T}}) where {T<:AbstractTensorMap} -> Tuple{Int}

Return all indices of a tensor, i.e. the indices of its domain and codomain.

See also codomainind and domainind.

source

To obtain information about the data, the following methods exist:

TensorKit.blockdimMethod
blockdim(t::AbstractTensorMap, c::Sector) -> Base.Dims

Return the dimensions of the block of a tensor corresponding to a coupled sector c.

source
TensorKit.fusiontreesMethod
fusiontrees(t::AbstractTensorMap)

Return an iterator over all splitting - fusion tree pairs of a tensor.

source
TensorKit.hasblockFunction
hasblock(t::AbstractTensorMap, c::Sector) -> Bool

Verify whether a tensor has a block corresponding to a coupled sector c.

source

For TensorMaps with Trivial sectortype, the data can be directly accessed and manipulated in a straightforward way:

Base.getindexMethod
Base.getindex(t::TrivialTensorMap)
t[]

Return a view into the data of t as a StridedViews.StridedView of size (dims(codomain(t))..., dims(domain(t))...).

source
Base.getindexMethod
Base.getindex(t::TrivialTensorMap, indices::Vararg{Int})
t[indices]

Return a view into the data slice of t corresponding to indices, by slicing the StridedViews.StridedView into the full data array.

source
Base.setindex!Method
Base.setindex!(t::TrivialTensorMap, v, indices::Vararg{Int})
t[indices] = v

Assigns v to the data slice of t corresponding to indices.

source

For general TensorMaps, this can be done using custom getindex and setindex! methods:

Base.getindexMethod
Base.getindex(t::TensorMap{<:IndexSpace,N₁,N₂,I},
              sectors::NTuple{N₁+N₂,I}) where {N₁,N₂,I<:Sector} 
    -> StridedViews.StridedView
t[sectors]

Return a view into the data slice of t corresponding to the splitting - fusion tree pair with combined uncoupled charges sectors. In particular, if sectors == (s1..., s2...) where s1 and s2 correspond to the coupled charges in the codomain and domain respectively, then a StridedViews.StridedView of size (dims(codomain(t), s1)..., dims(domain(t), s2)) is returned.

This method is only available for the case where FusionStyle(I) isa UniqueFusion, since it assumes a uniquely defined coupled charge.

source
Base.getindexMethod
Base.getindex(t::TensorMap{<:IndexSpace,N₁,N₂,I},
              f₁::FusionTree{I,N₁},
              f₂::FusionTree{I,N₂}) where {N₁,N₂,I<:Sector}
    -> StridedViews.StridedView
t[f₁, f₂]

Return a view into the data slice of t corresponding to the splitting - fusion tree pair (f₁, f₂). In particular, if f₁.coupled == f₂.coupled == c, then a StridedViews.StridedView of size (dims(codomain(t), f₁.uncoupled)..., dims(domain(t), f₂.uncoupled)) is returned which represents the slice of block(t, c) whose row indices correspond to f₁.uncoupled and column indices correspond to f₂.uncoupled.

source
Base.setindex!Method
Base.setindex!(t::TensorMap{<:IndexSpace,N₁,N₂,I},
               v,
               f₁::FusionTree{I,N₁},
               f₂::FusionTree{I,N₂}) where {N₁,N₂,I<:Sector}
t[f₁, f₂] = v

Copies v into the data slice of t corresponding to the splitting - fusion tree pair (f₁, f₂). Here, v can be any object that can be copied into a StridedViews.StridedView of size (dims(codomain(t), f₁.uncoupled)..., dims(domain(t), f₂.uncoupled)) using Base.copy!.

See also Base.getindex(::TensorMap{<:IndexSpace,N₁,N₂,I<:Sector}, ::FusionTree{I<:Sector,N₁}, ::FusionTree{I<:Sector,N₂})

source

TensorMap operations

The operations that can be performed on a TensorMap can be organized into index manipulations, (planar) traces and (planar) contractions.

Index manipulations

A general index manipulation of a TensorMap object can be built up by considering some transformation of the fusion trees, along with a permutation of the stored data. They come in three flavours, which are either of the type transform(!) which are exported, or of the type add_transform!, for additional expert-mode options that allows for addition and scaling, as well as the selection of a custom backend.

TensorKit.permuteMethod
permute(tsrc::AbstractTensorMap{S}, (p₁, p₂)::Tuple{IndexTuple{N₁},IndexTuple{N₂}};
        copy::Bool=false) where {S,N₁,N₂}
    -> tdst::TensorMap{S,N₁,N₂}

Return tensor tdst obtained by permuting the indices of tsrc. The codomain and domain of tdst correspond to the indices in p₁ and p₂ of tsrc respectively.

If copy=false, tdst might share data with tsrc whenever possible. Otherwise, a copy is always made.

To permute into an existing destination, see permute! and add_permute!

source
TensorKit.braidMethod
braid(tsrc::AbstractTensorMap{S}, (p₁, p₂)::Tuple{IndexTuple{N₁},IndexTuple{N₂}}, levels::Tuple;
      copy::Bool = false) where {S,N₁,N₂}
    -> tdst::TensorMap{S,N₁,N₂}

Return tensor tdst obtained by braiding the indices of tsrc. The codomain and domain of tdst correspond to the indices in p₁ and p₂ of tsrc respectively. Here, levels is a tuple of length numind(tsrc) that assigns a level or height to the indices of tsrc, which determines whether they will braid over or under any other index with which they have to change places.

If copy=false, tdst might share data with tsrc whenever possible. Otherwise, a copy is always made.

To braid into an existing destination, see braid! and add_braid!

source
Base.transposeMethod
transpose(tsrc::AbstractTensorMap{S}, (p₁, p₂)::Tuple{IndexTuple{N₁},IndexTuple{N₂}};
          copy::Bool=false) where {S,N₁,N₂}
    -> tdst::TensorMap{S,N₁,N₂}

Return tensor tdst obtained by transposing the indices of tsrc. The codomain and domain of tdst correspond to the indices in p₁ and p₂ of tsrc respectively. The new index positions should be attainable without any indices crossing each other, i.e., the permutation (p₁..., reverse(p₂)...) should constitute a cyclic permutation of (codomainind(tsrc)..., reverse(domainind(tsrc))...).

If copy=false, tdst might share data with tsrc whenever possible. Otherwise, a copy is always made.

To permute into an existing destination, see permute! and add_permute!

source
TensorKit.twistMethod
twist(t::AbstractTensorMap, i::Int; inv::Bool=false)
    -> t

Apply a twist to the ith index of t and return the result as a new tensor. If inv=true, use the inverse twist.

See twist! for storing the result in place.

source
Base.permute!Method
permute!(tdst::AbstractTensorMap{S,N₁,N₂}, tsrc::AbstractTensorMap{S},
         (p₁, p₂)::Tuple{IndexTuple{N₁},IndexTuple{N₂}}) where {S,N₁,N₂}
    -> tdst

Write into tdst the result of permuting the indices of tsrc. The codomain and domain of tdst correspond to the indices in p₁ and p₂ of tsrc respectively.

See permute for creating a new tensor and add_permute! for a more general version.

source
TensorKit.braid!Function
braid!(tdst::AbstractTensorMap{S,N₁,N₂}, tsrc::AbstractTensorMap{S},
       (p₁, p₂)::Tuple{IndexTuple{N₁},IndexTuple{N₂}}, levels::Tuple) where {S,N₁,N₂}
    -> tdst

Write into tdst the result of braiding the indices of tsrc. The codomain and domain of tdst correspond to the indices in p₁ and p₂ of tsrc respectively. Here, levels is a tuple of length numind(tsrc) that assigns a level or height to the indices of tsrc, which determines whether they will braid over or under any other index with which they have to change places.

See braid for creating a new tensor and add_braid! for a more general version.

source
LinearAlgebra.transpose!Function
transpose!(tdst::AbstractTensorMap{S,N₁,N₂}, tsrc::AbstractTensorMap{S},
           (p₁, p₂)::Tuple{IndexTuple{N₁},IndexTuple{N₂}}) where {S,N₁,N₂}
    -> tdst

Write into tdst the result of transposing the indices of tsrc. The codomain and domain of tdst correspond to the indices in p₁ and p₂ of tsrc respectively. The new index positions should be attainable without any indices crossing each other, i.e., the permutation (p₁..., reverse(p₂)...) should constitute a cyclic permutation of (codomainind(tsrc)..., reverse(domainind(tsrc))...).

See transpose for creating a new tensor and add_transpose! for a more general version.

source
TensorKit.twist!Function
twist!(t::AbstractTensorMap, i::Int; inv::Bool=false)
    -> t

Apply a twist to the ith index of t, storing the result in t. If inv=true, use the inverse twist.

See twist for creating a new tensor.

source
Missing docstring.

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

Missing docstring.

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

Missing docstring.

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

Missing docstring.

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

Traces and contractions

Missing docstring.

Missing docstring for trace_permute!. 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, rightind)::Index2Tuple;
            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 unique (no residual freedom) so that they always return the same result for the same input tensor t.

Orthogonality requires InnerProductStyle(t) <: HasInnerProduct, and leftorth(!) is currently only implemented for InnerProductStyle(t) === EuclideanProduct().

source
TensorKit.rightorthFunction
rightorth(t::AbstractTensorMap, (leftind, rightind)::Index2Tuple;
            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 unique (no residual freedom) so that they always return the same result for the same input tensor t.

Orthogonality requires InnerProductStyle(t) <: HasInnerProduct, and rightorth(!) is currently only implemented for InnerProductStyle(t) === EuclideanProduct().

source
TensorKit.leftnullFunction
leftnull(t::AbstractTensor, (leftind, rightind)::Index2Tuple;
            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 InnerProductStyle(t) <: HasInnerProduct, and leftnull(!) is currently only implemented for InnerProductStyle(t) === EuclideanProduct().

source
TensorKit.rightnullFunction
rightnull(t::AbstractTensor, (leftind, rightind)::Index2Tuple;
            alg::OrthogonalFactorizationAlgorithm = LQ(),
            atol::Real = 0.0,
            rtol::Real = eps(real(float(one(scalartype(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 InnerProductStyle(t) <: HasInnerProduct, and rightnull(!) is currently only implemented for InnerProductStyle(t) === EuclideanProduct().

source
TensorKit.tsvdFunction
tsvd(t::AbstractTensorMap, (leftind, rightind)::Index2Tuple;
    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.
  • truncbelow(χ::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 InnerProductStyle(t) <: HasInnerProduct, and tsvd(!) is currently only implemented for InnerProductStyle(t) === EuclideanProduct().

source
TensorKit.eighFunction
eigh(t::AbstractTensorMap, (leftind, rightind)::Index2Tuple) -> 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 scalartype 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 InnerProductStyle(t) === EuclideanProduct().

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, rightind)::Index2Tuple; 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 and permute as eigen of dense matrices. See the corresponding documentation for more information.

See also eigen and eigh.

source
LinearAlgebra.eigenFunction
eigen(t::AbstractTensor, (leftind, rightind)::Index2Tuple; 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 and permute as eigen of dense matrices. See the corresponding documentation for more information.

See also eig and eigh

source