ITensor

Description

ITensors.ITensorType
ITensor{N}

An ITensor is a tensor whose interface is independent of its memory layout. Therefore it is not necessary to know the ordering of an ITensor's indices, only which indices an ITensor has. Operations like contraction and addition of ITensors automatically handle any memory permutations.

Examples

julia> i = Index(2, "i")
(dim=2|id=287|"i")

julia> A = randomITensor(i', i)
ITensor ord=2 (dim=2|id=287|"i")' (dim=2|id=287|"i")
NDTensors.Dense{Float64,Array{Float64,1}}

julia> @show A;
A = ITensor ord=2
Dim 1: (dim=2|id=287|"i")'
Dim 2: (dim=2|id=287|"i")
NDTensors.Dense{Float64,Array{Float64,1}}
 2×2
 0.28358594718392427   1.4342219756446355
 1.6620103556283987   -0.40952231269251566

julia> @show inds(A);
inds(A) = IndexSet{2} (dim=2|id=287|"i")' (dim=2|id=287|"i") 

julia> A[i => 1, i' => 2] = 1;

julia> @show A;
A = ITensor ord=2
Dim 1: (dim=2|id=287|"i")'
Dim 2: (dim=2|id=287|"i")
NDTensors.Dense{Float64,Array{Float64,1}}
 2×2
 0.28358594718392427   1.4342219756446355
 1.0                  -0.40952231269251566

julia> @show store(A);
store(A) = [0.28358594718392427, 1.0, 1.4342219756446355, -0.40952231269251566]

julia> B = randomITensor(i, i');

julia> @show B;
B = ITensor ord=2
Dim 1: (dim=2|id=287|"i")
Dim 2: (dim=2|id=287|"i")'
NDTensors.Dense{Float64,Array{Float64,1}}
 2×2
 -0.6510816500352691   0.2579101497658179
  0.256266641521826   -0.9464735926768166

julia> @show A + B;
A + B = ITensor ord=2
Dim 1: (dim=2|id=287|"i")'
Dim 2: (dim=2|id=287|"i")
NDTensors.Dense{Float64,Array{Float64,1}}
 2×2
 -0.3674957028513448   1.6904886171664615
  1.2579101497658178  -1.3559959053693322
source

Dense Constructors

ITensors.ITensorMethod
ITensor([::Type{ElT} = Float64, ]inds)
ITensor([::Type{ElT} = Float64, ]inds::Index...)

Construct an ITensor filled with zeros having indices inds and element type ElT. If the element type is not specified, it defaults to Float64.

The storage will have NDTensors.Dense type.

source
ITensors.ITensorMethod
ITensor([::Type{ElT} = Float64, ]::UndefInitializer, inds)
ITensor([::Type{ElT} = Float64, ]::UndefInitializer, inds::Index...)

Construct an ITensor filled with undefined elements having indices inds and element type ElT. If the element type is not specified, it defaults to Float64.

The storage will have NDTensors.Dense type.

source
ITensors.randomITensorMethod
randomITensor([::Type{ElT <: Number} = Float64, ]inds)
randomITensor([::Type{ElT <: Number} = Float64, ]inds::Index...)

Construct an ITensor with type ElT and indices inds, whose elements are normally distributed random numbers. If the element type is not specified, it defaults to Float64.

source
ITensors.seteltMethod
setelt(ivs...)

Create an ITensor with all zeros except the specified value, which is set to 1.

Examples

i = Index(2,"i")
A = setelt(i=>2)
# A[i=>2] == 1, all other elements zero

j = Index(3,"j")
B = setelt(i=>1,j=>3)
# B[i=>1,j=>3] == 1, all other element zero
source

QN BlockSparse Constructors

ITensors.ITensorMethod
ITensor([::Type{ElT} = Float64, ][flux::QN = QN(), ]inds)
ITensor([::Type{ElT} = Float64, ][flux::QN = QN(), ]inds::Index...)

Construct an ITensor with BlockSparse storage filled with zero(ElT) where the nonzero blocks are determined by flux.

If ElT is not specified it defaults to Float64.

source
ITensors.ITensorMethod
ITensor(::Array, ::IndexSet; tol = 0)

ITensor(::Array, ::Index...; tol = 0)

Create a block sparse ITensor from the input Array, and collection of QN indices. Zeros are dropped and nonzero blocks are determined from the zero values of the array.

Optionally, you can set a tolerance such that elements less than or equal to the tolerance are dropped.

Examples

julia> i = Index([QN(0)=>1, QN(1)=>2], "i");

julia> A = [1e-9 0.0 0.0;
            0.0 2.0 3.0;
            0.0 1e-10 4.0];

julia> @show ITensor(A, i', dag(i); tol = 1e-8);
ITensor(A, i', dag(i); tol = 1.0e-8) = ITensor ord=2
Dim 1: (dim=3|id=468|"i")' <Out>
 1: QN(0) => 1
 2: QN(1) => 2
Dim 2: (dim=3|id=468|"i") <In>
 1: QN(0) => 1
 2: QN(1) => 2
NDTensors.BlockSparse{Float64,Array{Float64,1},2}
 3×3
Block: (2, 2)
 [2:3, 2:3]
 2.0  3.0
 0.0  4.0
source

Empty Constructors

ITensors.emptyITensorMethod
emptyITensor([::Type{ElT} = Float64, ]inds)
emptyITensor([::Type{ElT} = Float64, ]inds::Index...)

Construct an ITensor with storage type NDTensors.Empty, indices inds, and element type ElT. If the element type is not specified, it defaults to Float64.

source

QN Empty Constructors

ITensors.emptyITensorMethod
emptyITensor([::Type{ElT} = Float64, ]inds)
emptyITensor([::Type{ElT} = Float64, ]inds::Index...)

Construct an ITensor with storage type NDTensors.Empty, indices inds, and element type ElT. If the element type is not specified, it defaults to Float64.

source
emptyITensor([::Type{ElT} = Float64, ]inds)
emptyITensor([::Type{ElT} = Float64, ]inds::QNIndex...)

Construct an ITensor with NDTensors.BlockSparse storage of element type ElT with the no blocks.

If ElT is not specified it defaults to Float64.

In the future, this will use the storage NDTensors.EmptyBlockSparse.

source

Diagonal constructors

ITensors.diagITensorMethod
diagITensor([::Type{ElT} = Float64, ]inds)
diagITensor([::Type{ElT} = Float64, ]inds::Index...)

Make a sparse ITensor of element type ElT with only elements along the diagonal stored. Defaults to having zero(T) along the diagonal.

The storage will have NDTensors.Diag type.

source
ITensors.diagITensorMethod
diagITensor(v::Vector{T}, inds)
diagITensor(v::Vector{T}, inds::Index...)

Make a sparse ITensor with non-zero elements only along the diagonal. The diagonal elements will be set to the values stored in v and the ITensor will have element type float(T). The storage will have type NDTensors.Diag.

source
ITensors.diagITensorMethod
diagITensor(x::Number, inds)
diagITensor(x::Number, inds::Index...)

Make a sparse ITensor with non-zero elements only along the diagonal. The diagonal elements will be set to the value float(x) and the ITensor will have element type float(eltype(x)). The storage will have NDTensors.Diag type.

source
ITensors.deltaMethod
delta([::Type{ElT} = Float64, ]inds)
delta([::Type{ElT} = Float64, ]inds::Index...)

Make a uniform diagonal ITensor with all diagonal elements one(ElT). Only a single diagonal element is stored.

This function has an alias δ.

source

QN Diagonal constructors

ITensors.diagITensorMethod
diagITensor([::Type{ElT} = Float64, ][flux::QN = QN(), ]is)
diagITensor([::Type{ElT} = Float64, ][flux::QN = QN(), ]is::Index...)

Make an ITensor with storage type NDTensors.DiagBlockSparse with elements zero(ElT). The ITensor only has diagonal blocks consistent with the specified flux.

If the element type is not specified, it defaults to Float64. If theflux is not specified, it defaults to QN().

source
ITensors.deltaMethod
delta([::Type{ElT} = Float64, ][flux::QN = QN(), ]is)
delta([::Type{ElT} = Float64, ][flux::QN = QN(), ]is::Index...)

Make an ITensor with storage type NDTensors.DiagBlockSparse with uniform elements one(ElT). The ITensor only has diagonal blocks consistent with the specified flux.

If the element type is not specified, it defaults to Float64. If theflux is not specified, it defaults to QN().

source

Convert to Array

Core.ArrayMethod
Array{ElT}(T::ITensor, i:Index...)
Array(T::ITensor, i:Index...)

Matrix{ElT}(T::ITensor, row_i:Index, col_i::Index)
Matrix(T::ITensor, row_i:Index, col_i::Index)

Vector{ElT}(T::ITensor)
Vector(T::ITensor)

Given an ITensor T with indices i..., returns an Array with a copy of the ITensor's elements. The order in which the indices are provided indicates the order of the data in the resulting Array.

source
NDTensors.arrayMethod
array(T::ITensor)

Given an ITensor T, returns an Array with a copy of the ITensor's elements, or a view in the case the the ITensor's storage is Dense. The ordering of the elements in the Array, in terms of which Index is treated as the row versus column, depends on the internal layout of the ITensor. Therefore this method is intended for developer use only and not recommended for use in ITensor applications.

source
NDTensors.matrixMethod
matrix(T::ITensor)

Given an ITensor T with two indices, returns a Matrix with a copy of the ITensor's elements, or a view in the case the ITensor's storage is Dense. The ordering of the elements in the Matrix, in terms of which Index is treated as the row versus column, depends on the internal layout of the ITensor. Therefore this method is intended for developer use only and not recommended for use in ITensor applications.

source
NDTensors.vectorMethod
vector(T::ITensor)

Given an ITensor T with one index, returns a Vector with a copy of the ITensor's elements, or a view in the case the ITensor's storage is Dense.

source

Getting and setting elements

Base.getindexMethod
getindex(T::ITensor, ivs...)

Get the specified element of the ITensor, using a list of IndexVals or Pair{<:Index, Int}.

Example

i = Index(2; tags = "i")
A = ITensor(2.0, i, i')
A[i => 1, i' => 2] # 2.0, same as: A[i' => 2, i => 1]
source
Base.getindexMethod
getindex(T::ITensor, I::Int...)

Get the specified element of the ITensor, using internal Index ordering of the ITensor.

Example

i = Index(2; tags = "i")
A = ITensor(2.0, i, i')
A[1, 2] # 2.0, same as: A[i => 1, i' => 2]
source
Base.setindex!Method
setindex!(T::ITensor, x::Number, ivs...)

Set the specified element of the ITensor using a list of IndexVals or Pair{<:Index, Int}.

Example

i = Index(2; tags = "i")
A = ITensor(i, i')
A[i => 1, i' => 2] = 1.0 # same as: A[i' => 2, i => 1] = 1.0
A[i => 2, i' => :] = [2.0 3.0]
source
Base.setindex!Method
setindex!(T::ITensor, x::Number, I::Int...)

setindex!(T::ITensor, x::Number, I::CartesianIndex)

Set the specified element of the ITensor, using internal Index ordering of the ITensor.

Example

i = Index(2; tags = "i")
A = ITensor(i, i')
A[1, 2] = 1.0 # same as: A[i => 1, i' => 2] = 1.0
A[2, :] = [2.0 3.0]
source

Properties

NDTensors.indMethod
ind(T::ITensor, i::Int)

Get the Index of the ITensor along dimension i.

source

Priming and tagging

ITensors.primeMethod
prime[!](A::ITensor, plinc::Int = 1; <keyword arguments>) -> ITensor

prime(is::IndexSet, plinc::Int = 1; <keyword arguments>) -> IndexSet

Increase the prime level of the indices of an ITensor or IndexSet.

Optionally, only modify the indices with the specified keyword arguments.

Arguments

  • tags = nothing: if specified, only modify Index i if hastags(i, tags) == true.
  • plev = nothing: if specified, only modify Index i if hasplev(i, plev) == true.

The ITensor functions come in two versions, f and f!. The latter modifies the ITensor in-place. In both versions, the ITensor storage is not modified or copied (so it returns an ITensor with a view of the original storage).

source
ITensors.setprimeMethod
setprime[!](A::ITensor, plev::Int; <keyword arguments>) -> ITensor

setprime(is::IndexSet, plev::Int; <keyword arguments>) -> IndexSet

Set the prime level of the indices of an ITensor or IndexSet.

Optionally, only modify the indices with the specified keyword arguments.

Arguments

  • tags = nothing: if specified, only modify Index i if hastags(i, tags) == true.
  • plev = nothing: if specified, only modify Index i if hasplev(i, plev) == true.

The ITensor functions come in two versions, f and f!. The latter modifies the ITensor in-place. In both versions, the ITensor storage is not modified or copied (so it returns an ITensor with a view of the original storage).

source
ITensors.noprimeMethod
noprime[!](A::ITensor; <keyword arguments>) -> ITensor

noprime(is::IndexSet; <keyword arguments>) -> IndexSet

Set the prime level of the indices of an ITensor or IndexSet to zero.

Optionally, only modify the indices with the specified keyword arguments.

Arguments

  • tags = nothing: if specified, only modify Index i if hastags(i, tags) == true.
  • plev = nothing: if specified, only modify Index i if hasplev(i, plev) == true.

The ITensor functions come in two versions, f and f!. The latter modifies the ITensor in-place. In both versions, the ITensor storage is not modified or copied (so it returns an ITensor with a view of the original storage).

source
ITensors.mapprimeMethod
replaceprime[!](A::ITensor, plold::Int, plnew::Int; <keyword arguments>) -> ITensor
replaceprime[!](A::ITensor, plold => plnew; <keyword arguments>) -> ITensor
mapprime[!](A::ITensor, <arguments>; <keyword arguments>) -> ITensor

replaceprime(is::IndexSet, plold::Int, plnew::Int; <keyword arguments>) -> IndexSet
replaceprime(is::IndexSet, plold => plnew; <keyword arguments>) -> IndexSet
mapprime(is::IndexSet, <arguments>; <keyword arguments>) -> IndexSet

Set the prime level of the indices of an ITensor or IndexSet with prime level plold to plnew.

Optionally, only modify the indices with the specified keyword arguments.

Arguments

  • tags = nothing: if specified, only modify Index i if hastags(i, tags) == true.
  • plev = nothing: if specified, only modify Index i if hasplev(i, plev) == true.

The ITensor functions come in two versions, f and f!. The latter modifies the ITensor in-place. In both versions, the ITensor storage is not modified or copied (so it returns an ITensor with a view of the original storage).

source
ITensors.swapprimeMethod
swapprime[!](A::ITensor, pl1::Int, pl2::Int; <keyword arguments>) -> ITensor
swapprime[!](A::ITensor, pl1 => pl2; <keyword arguments>) -> ITensor

swapprime(is::ITensor, pl1::Int, pl2::Int; <keyword arguments>) -> IndexSet
swapprime(is::ITensor, pl1 => pl2; <keyword arguments>) -> IndexSet

Set the prime level of the indices of an ITensor or IndexSetwith prime level pl1 to pl2, and those with prime level pl2 to pl1.

Optionally, only modify the indices with the specified keyword arguments.

Arguments

  • tags = nothing: if specified, only modify Index i if hastags(i, tags) == true.
  • plev = nothing: if specified, only modify Index i if hasplev(i, plev) == true.

The ITensor functions come in two versions, f and f!. The latter modifies the ITensor in-place. In both versions, the ITensor storage is not modified or copied (so it returns an ITensor with a view of the original storage).

source
ITensors.addtagsMethod
addtags[!](A::ITensor, ts::String; <keyword arguments>) -> ITensor

addtags(is::IndexSet, ts::String; <keyword arguments>) -> IndexSet

Add the tags ts to the indices of an ITensor or IndexSet.

Optionally, only modify the indices with the specified keyword arguments.

Arguments

  • tags = nothing: if specified, only modify Index i if hastags(i, tags) == true.
  • plev = nothing: if specified, only modify Index i if hasplev(i, plev) == true.

The ITensor functions come in two versions, f and f!. The latter modifies the ITensor in-place. In both versions, the ITensor storage is not modified or copied (so it returns an ITensor with a view of the original storage).

source
ITensors.removetagsMethod
removetags[!](A::ITensor, ts::String; <keyword arguments>) -> ITensor

removetags(is::IndexSet, ts::String; <keyword arguments>) -> IndexSet

Remove the tags ts from the indices of an ITensor or IndexSet.

Optionally, only modify the indices with the specified keyword arguments.

Arguments

  • tags = nothing: if specified, only modify Index i if hastags(i, tags) == true.
  • plev = nothing: if specified, only modify Index i if hasplev(i, plev) == true.

The ITensor functions come in two versions, f and f!. The latter modifies the ITensor in-place. In both versions, the ITensor storage is not modified or copied (so it returns an ITensor with a view of the original storage).

source
ITensors.replacetagsMethod
replacetags[!](A::ITensor, tsold::String, tsnew::String; <keyword arguments>) -> ITensor

replacetags(is::IndexSet, tsold::String, tsnew::String; <keyword arguments>) -> IndexSet

Replace the tags tsold with tsnew for the indices of an ITensor.

Optionally, only modify the indices with the specified keyword arguments.

Arguments

  • tags = nothing: if specified, only modify Index i if hastags(i, tags) == true.
  • plev = nothing: if specified, only modify Index i if hasplev(i, plev) == true.

The ITensor functions come in two versions, f and f!. The latter modifies the ITensor in-place. In both versions, the ITensor storage is not modified or copied (so it returns an ITensor with a view of the original storage).

source
ITensors.settagsMethod
settags[!](A::ITensor, ts::String; <keyword arguments>) -> ITensor

settags(is::IndexSet, ts::String; <keyword arguments>) -> IndexSet

Set the tags of the indices of an ITensor or IndexSet to ts.

Optionally, only modify the indices with the specified keyword arguments.

Arguments

  • tags = nothing: if specified, only modify Index i if hastags(i, tags) == true.
  • plev = nothing: if specified, only modify Index i if hasplev(i, plev) == true.

The ITensor functions come in two versions, f and f!. The latter modifies the ITensor in-place. In both versions, the ITensor storage is not modified or copied (so it returns an ITensor with a view of the original storage).

source
ITensors.swaptagsMethod
swaptags[!](A::ITensor, ts1::String, ts2::String; <keyword arguments>) -> ITensor

swaptags(is::IndexSet, ts1::String, ts2::String; <keyword arguments>) -> IndexSet

Swap the tags ts1 with ts2 for the indices of an ITensor.

Optionally, only modify the indices with the specified keyword arguments.

Arguments

  • tags = nothing: if specified, only modify Index i if hastags(i, tags) == true.
  • plev = nothing: if specified, only modify Index i if hasplev(i, plev) == true.

The ITensor functions come in two versions, f and f!. The latter modifies the ITensor in-place. In both versions, the ITensor storage is not modified or copied (so it returns an ITensor with a view of the original storage).

source

IndexSet set operations

ITensors.commonindsFunction
commoninds(A, B; kwargs...)
commoninds(::Order{N}, A, B; kwargs...)

Return an IndexSet with indices that are common between the indices of A and B (the set intersection, similar to Base.intersect).

Optionally, specify the desired number of indices as Order(N), which adds a check and can be a bit more efficient.

source
ITensors.uniqueindsFunction
uniqueinds(A, B; kwargs...)
uniqueinds(::Order{N}, A, B; kwargs...)

Return an IndexSet with indices that are unique to the set of indices of A and not in B (the set difference, similar to Base.setdiff).

Optionally, specify the desired number of indices as Order(N), which adds a check and can be a bit more efficient.

source
ITensors.noncommonindsFunction
noncommoninds(A, B; kwargs...)
noncommoninds(::Order{N}, A, B; kwargs...)

Return an IndexSet with indices that are not common between the indices of A and B (the symmetric set difference, similar to Base.symdiff).

Optionally, specify the desired number of indices as Order(N), which adds a check and can be a bit more efficient.

source
ITensors.unionindsFunction
unioninds(A, B; kwargs...)
unioninds(::Order{N}, A, B; kwargs...)

Return an IndexSet with indices that are the union of the indices of A and B (the set union, similar to Base.union).

Optionally, specify the desired number of indices as Order(N), which adds a check and can be a bit more efficient.

source
ITensors.hascommonindsFunction
hascommoninds(A, B; kwargs...)

hascommoninds(B; kwargs...) -> f::Function

Check if the ITensors or sets of indices A and B have common indices.

If only one ITensor or set of indices B is passed, return a function f such that f(A) = hascommoninds(A, B; kwargs...)

source

Index Manipulations

ITensors.replaceindMethod
replaceind[!](A::ITensor, i1::Index, i2::Index) -> ITensor

Replace the Index i1 with the Index i2 in the ITensor.

The indices must have the same space (i.e. the same dimension and QNs, if applicable).

source
ITensors.replaceindsMethod
replaceinds(A::ITensor, inds1, inds2) -> ITensor

replaceinds!(A::ITensor, inds1, inds2)

Replace the Index inds1[n] with the Index inds2[n] in the ITensor, where n runs from 1 to length(inds1) == length(inds2).

The indices must have the same space (i.e. the same dimension and QNs, if applicable).

The storage of the ITensor is not modified or copied (the output ITensor is a view of the input ITensor).

source
ITensors.swapindMethod
swapind(A::ITensor, i1::Index, i2::Index) -> ITensor

swapind!(A::ITensor, i1::Index, i2::Index)

Swap the Index i1 with the Index i2 in the ITensor.

The indices must have the same space (i.e. the same dimension and QNs, if applicable).

source
ITensors.swapindsMethod
swapinds(A::ITensor, inds1, inds2) -> ITensor

swapinds!(A::ITensor, inds1, inds2)

Swap the Index inds1[n] with the Index inds2[n] in the ITensor, where n runs from 1 to length(inds1) == length(inds2).

The indices must have the same space (i.e. the same dimension and QNs, if applicable).

The storage of the ITensor is not modified or copied (the output ITensor is a view of the input ITensor).

source

Math operations

Base.:*Method
A::ITensor * B::ITensor

Contract ITensors A and B to obtain a new ITensor. This contraction * operator finds all matching indices common to A and B and sums over them, such that the result will have only the unique indices of A and B. To prevent indices from matching, their prime level or tags can be modified such that they no longer compare equal - for more information see the documentation on Index objects.

source
Base.expMethod
exp(A::ITensor, Linds=Rinds', Rinds=inds(A,plev=0); ishermitian = false)

Compute the exponential of the tensor A by treating it as a matrix $A_{lr}$ with the left index l running over all indices in Linds and r running over all indices in Rinds.

Only accepts index lists Linds,Rinds such that: (1) length(Linds) + length(Rinds) == length(inds(A)) (2) length(Linds) == length(Rinds) (3) For each pair of indices (Linds[n],Rinds[n]), Linds[n] and Rinds[n] represent the same Hilbert space (the same QN structure in the QN case, or just the same length in the dense case), and appear in A with opposite directions.

When ishermitian=true the exponential of Hermitian(A_{lr}) is computed internally.

source

Decompositions

LinearAlgebra.svdMethod
svd(A::ITensor, inds::Index...; <keyword arguments>)

Singular value decomposition (SVD) of an ITensor A, computed by treating the "left indices" provided collectively as a row index, and the remaining "right indices" as a column index (matricization of a tensor).

The first three return arguments are U, S, and V, such that A ≈ U * S * V.

Whether or not the SVD performs a trunction depends on the keyword arguments provided.

Arguments

  • maxdim::Int: the maximum number of singular values to keep.
  • mindim::Int: the minimum number of singular values to keep.
  • cutoff::Float64: set the desired truncation error of the SVD, by default defined as the sum of the squares of the smallest singular values.
  • lefttags::String = "Link,u": set the tags of the Index shared by U and S.
  • righttags::String = "Link,v": set the tags of the Index shared by S and V.
  • alg::String = "divide_and_conquer". Options:
    • "divide_and_conquer" - A divide-and-conquer algorithm. LAPACK's gesdd. Fast, but may lead to some innacurate singular values for very ill-conditioned matrices. Also may sometimes fail to converge, leading to errors (in which case "qr_iteration" or "recursive" can be tried).
    • "qr_iteration" - Typically slower but more accurate for very ill-conditioned matrices compared to "divide_and_conquer". LAPACK's gesvd.
    • "recursive" - ITensor's custom svd. Very reliable, but may be slow if high precision is needed. To get an svd of a matrix A, an eigendecomposition of $A^{\dagger} A$ is used to compute U and then a qr of $A^{\dagger} U$ is used to compute V. This is performed recursively to compute small singular values.
  • use_absolute_cutoff::Bool = false: set if all probability weights below the cutoff value should be discarded, rather than the sum of discarded weights.
  • use_relative_cutoff::Bool = true: set if the singular values should be normalized for the sake of truncation.

See also: factorize

source
LinearAlgebra.factorizeMethod
factorize(A::ITensor, Linds::Index...; <keyword arguments>)

Perform a factorization of A into ITensors L and R such that A ≈ L * R.

Arguments

  • ortho::String = "left": Choose orthogonality properties of the factorization.
    • "left": the left factor L is an orthogonal basis such that L * dag(prime(L, commonind(L,R))) ≈ I.
    • "right": the right factor R forms an orthogonal basis.
    • "none", neither of the factors form an orthogonal basis, and in general are made as symmetrically as possible (depending on the decomposition used).
  • which_decomp::Union{String, Nothing} = nothing: choose what kind of decomposition is used.
    • nothing: choose the decomposition automatically based on the other arguments. For example, when nothing is chosen and ortho = "left" or "right", and a cutoff is provided, svd or eigen is used depending on the provided cutoff (eigen is only used when the cutoff is greater than 1e-12, since it has a lower precision). When no truncation is requested qr is used for dense ITensors and svd for block-sparse ITensors (in the future qr will be used also for block-sparse ITensors in this case).
    • "svd": L = U and R = S * V for ortho = "left", L = U * S and R = V for ortho = "right", and L = U * sqrt.(S) and R = sqrt.(S) * V for ortho = "none". To control which svd algorithm is choose, use the svd_alg keyword argument. See the documentation for svd for the supported algorithms, which are the same as those accepted by the alg keyword argument.
    • "eigen": L = U and $R = U^{\dagger} A$ where U is determined from the eigendecompositon $A A^{\dagger} = U D U^{\dagger}$ for ortho = "left" (and vice versa for ortho = "right"). "eigen" is not supported for ortho = "none".
    • "qr": L=Q and R an upper-triangular matrix when ortho = "left", and R = Q and L a lower-triangular matrix when ortho = "right" (currently supported for dense ITensors only).

In the future, other decompositions like QR (for block-sparse ITensors), polar, cholesky, LU, etc. are expected to be supported.

For truncation arguments, see: svd

source

Operations

ITensors.permuteMethod
permute(T::ITensors, inds; always_copy::Bool = false)
permute(T::ITensors, inds::Index...; always_copy::Bool = false)

Return a new ITensor T with indices permuted according to the input indices inds. The storage of the ITensor is permuted accordingly.

If always_copy = false, it avoids copying data if possible. Therefore, it may return a view. Use always_copy = true if you never want it to return an ITensor with a view of the original ITensor.

source