Types
TensorTrains.AbstractTensorTrain — TypeAbstractTensorTrain{F<:Number, N}An abstract type representing a Tensor Train.
TensorTrains.AbstractPeriodicTensorTrain — TypeAbstractPeriodicTensorTrain{F<:Number, N} <: AbstractTensorTrain{F,N}An abstract type representing a Tensor Train with periodic boundary conditions.
TensorTrains.TensorTrain — TypeTensorTrain{F<:Number, N} <: AbstractTensorTrain{F,N}A type for representing a Tensor Train
Fis the type of the matrix entriesNis the number of indices of each tensor (2 virtual ones +N-2physical ones)
TensorTrains.PeriodicTensorTrain — TypePeriodicTensorTrain{F<:Number, N} <: AbstractPeriodicTensorTrain{F,N}A type for representing a Tensor Train with periodic boundary conditions
Fis the type of the matrix entriesNis the number of indices of each tensor (2 virtual ones +N-2physical ones)
Functions
TensorTrains.normalize_eachmatrix! — Functionnormalize_eachmatrix!(A::AbstractTensorTrain)Divide each matrix by its maximum (absolute) element and return the sum of the logs of the individual normalizations. This is used to keep the entries from exploding during computations
TensorTrains.flat_tt — Functionflat_tt([T = Float64], bondsizes::AbstractVector{<:Integer}, q...)
flat_tt([T = Float64], d::Integer, L::Integer, q...)Construct a (normalized) Tensor Train filled with one(T), by specifying either:
bondsizes: the size of each bondda fixed size for all bonds,Lthe length
and
qa Tuple/Vector specifying the number of values taken by each variable on a single site
TensorTrains.rand_tt — Functionrand_tt([rng=default_rng()], [T = Float64], bondsizes::AbstractVector{<:Integer}, q...)
rand_tt([rng=default_rng()], [T = Float64], d::Integer, L::Integer, q...)Construct a Tensor Train with rand(T) entries, by specifying either:
bondsizes: the size of each bondda fixed size for all bonds,Lthe length
and
qa Tuple/Vector specifying the number of values taken by each variable on a single site
TensorTrains.flat_periodic_tt — Functionflat_periodic_tt(bondsizes::AbstractVector{<:Integer}, q...)
flat_periodic_tt(d::Integer, L::Integer, q...)Construct a (normalized) Tensor Train with periodic boundary conditions filled with a constant, by specifying either:
bondsizes: the size of each bondda fixed size for all bonds,Lthe length
and
qa Tuple/Vector specifying the number of values taken by each variable on a single site
TensorTrains.rand_periodic_tt — Functionrand_periodic_tt(bondsizes::AbstractVector{<:Integer}, q...)
rand_periodic_tt(d::Integer, L::Integer, q...)Construct a Tensor Train with periodic boundary conditions with entries random in [0,1], by specifying either:
bondsizes: the size of each bondda fixed size for all bonds,Lthe length
and
qa Tuple/Vector specifying the number of values taken by each variable on a single site
TensorTrains.bond_dims — Functionbond_dims(A::AbstractTensorTrain)Return a vector with the dimensions of the virtual bonds
TensorTrains.evaluate — Functionevaluate(p::PMS, X...)Return the value of p for input X. If ψ is the tensor train wrapped by p, then the output is $\lvert\psi (X)\rvert^2$
evaluate(A::AbstractTensorTrain, X...)Evaluate the Tensor Train A at input X
Example:
L = 3
q = (2, 3)
A = rand_tt(4, L, q...)
X = [[rand(1:qi) for qi in q] for l in 1:L]
evaluate(A, X)TensorTrains.nparams — Functionnparams(A::AbstractTensorTrain)Return the number of parameters of the tensor train. One complex number counts as two parameters.
TensorTrains.is_in_domain — Functionis_in_domain(A::AbstractTensorTrain, X...)Return true if data X is in the domain of A. After this check it should be safe to use @inbounds when indexing A with X
TensorTrains.marginals — Functionmarginals(p::MPS; l, r)Compute the marginal distributions $p(x^l)$ at each site
Optional arguments
l = accumulate_L(A)[1],r = accumulate_R(A)[1]pre-computed partial normalizations
marginals(A::AbstractTensorTrain; l, r)Compute the marginal distributions $p(x^l)$ at each site
Optional arguments
l = accumulate_L(A)[1],r = accumulate_R(A)[1]pre-computed partial normalizations
TensorTrains.twovar_marginals — Functiontwovar_marginals(p::MPS; l, r, M, maxdist)Compute the marginal distributions for each pair of sites $p(x^l, x^m)$
Optional arguments
l = accumulate_L(p)[1],r = accumulate_R(p)[1],M = accumulate_M(p)pre-computed partial normalizationsmaxdist = length(p): compute marginals only at distancemaxdist: $|l-m|\le maxdist$
twovar_marginals(A::AbstractTensorTrain; l, r, M, Δlmax)Compute the marginal distributions for each pair of sites $p(x^l, x^m)$
Optional arguments
l = accumulate_L(A)[1],r = accumulate_R(A)[1],M = accumulate_M(A)pre-computed partial normalizationsmaxdist = length(A): compute marginals only at distancemaxdist: $|l-m|\le maxdist$
TensorTrains.lognormalization — Functionlognormalization(A::AbstractTensorTrain)Compute the natural logarithm of the normalization $\log Z=\log \sum_{x^1,\ldots,x^L} A^1(x^1)\cdots A^L(x^L)$. Throws an error if the normalization is negative.
TensorTrains.normalization — Functionnormalization(p::MPS)Compute the normalization $Z=\sum_{x^1,\ldots,x^L} \lvert A^1(x^1)\cdots A^L(x^L)\rvert ^2$ and return it as a Logarithmic, which stores the sign and the logarithm of the absolute value (see the docs of LogarithmicNumbers.jl https://github.com/cjdoris/LogarithmicNumbers.jl?tab=readme-ov-file#documentation)
normalization(A::AbstractTensorTrain)Compute the normalization $Z=\sum_{x^1,\ldots,x^L} A^1(x^1)\cdots A^L(x^L)$ and return it as a Logarithmic, which stores the sign and the logarithm of the absolute value (see the docs of LogarithmicNumbers.jl https://github.com/cjdoris/LogarithmicNumbers.jl?tab=readme-ov-file#documentation)
LinearAlgebra.normalize! — Functionnormalize!(p::MPS)Compute the normalization of $Z=\sum_{x^1,\ldots,x^L} \lvert A^1(x^1)\cdots A^L(x^L)\rvert ^2$ (see normalization) and rescale the tensors in p such that, after this call, $|Z|^2=1$. Return the natural logarithm of the absolute normalization $\log|Z|$
LinearAlgebra.normalize!(A::AbstractTensorTrain)Compute the normalization of $Z=\sum_{x^1,\ldots,x^L} A^1(x^1)\cdots A^L(x^L)$ (see normalization) and rescale the tensors in A such that, after this call, $|Z|=1$. Return the natural logarithm of the absolute normalization $\log|Z|$
Base.:+ — FunctionBase.:(+)(A::AbstracTensorTrain, B::AbstracTensorTrain)Compute the sum of two Tensor Trains. Matrix sizes are doubled
Base.:- — FunctionBase.:(-)(A::AbstracTensorTrain, B::AbstracTensorTrain)Compute the difference of two Tensor Trains. Matrix sizes are doubled
LinearAlgebra.dot — FunctionLinearAlgebra.dot(A::AbstractTensorTrain, B::AbstractTensorTrain)Compute the inner product between tensor trains A and B
\[A\cdot B = \sum_{x^1,x^2,\ldots,x^L}A^1(x^1)A^2(x^2)\cdots A^L(x^L)B^1(x^1)B^2(x^2)\cdots B^L(x^L)\]
dot(q::InfiniteUniformTensorTrain, p::InfiniteUniformTensorTrain)Return the "dot product" between the infinite tensor trains q and p.
Since two infinite tensor trains are either equal or orthogonal (orthogonality catastrophe), what is actually returned here is the leading eigenvalue of the transfer operator obtained from the matrices of q and p
LinearAlgebra.norm — FunctionLinearAlgebra.norm(A::AbstractTensorTrain)Compute the 2-norm (Frobenius norm) of tensor train A
\[\lVert A\rVert_2 = \sqrt{\sum_{x^1,x^2,\ldots,x^L}\left[A^1(x^1)A^2(x^2)\cdots A^L(x^L)\right]^2} = \sqrt{A\cdot A}\]
TensorTrains.norm2m — Functionnorm2m(A::AbstractTensorTrain, B::AbstractTensorTrain)Given two tensor trains A,B, compute norm(A - B)^2 as
\[\lVert A-B\rVert_2^2 = \lVert A \rVert_2^2 + \lVert B \rVert_2^2 - 2A\cdot B\]
StatsBase.sample! — Functionsample!([rng], x, p::MPS; r)Draw an exact sample from p and store the result in x.
Optionally specify a random number generator rng as the first argument (defaults to Random.default_rng()) and provide a pre-computed r = accumulate_R(p).
The output is x,q, the sampled sequence and its probability
sample!([rng], x, A::AbstractTensorTrain; r)Draw an exact sample from A interpreted as a probability distribution and store the result in x. A doesn't need to be normalized, however error will be raised if it is found to take negative values.
Optionally specify a random number generator rng as the first argument (defaults to Random.default_rng()) and provide a pre-computed rz = accumulate_R(A).
The output is x, p, the sampled sequence and its probability
StatsBase.sample — Functionsample([rng], p::MPS; r)Draw an exact sample from p.
Optionally specify a random number generator rng as the first argument (defaults to Random.default_rng()) and provide a pre-computed r = accumulate_R(p).
The output is x,q, the sampled sequence and its probability
sample([rng], A::AbstractTensorTrain; r)Draw an exact sample from A interpreted as a probability distribution. A doesn't need to be normalized, however error will be raised if it is found to take negative values.
Optionally specify a random number generator rng as the first argument (defaults to Random.default_rng()) and provide a pre-computed rz = accumulate_R(A).
The output is x, p, the sampled sequence and its probability
TensorTrains.orthogonalize_right! — Functionorthogonalize_right!(A::AbstractTensorTrain; svd_trunc::SVDTrunc, indices)Bring A to right-orthogonal form by means of SVD decompositions.
Optionally perform truncations by passing a SVDTrunc. Optionally pass a range of indices to perform the orthogonalizations only on those.
TensorTrains.orthogonalize_left! — Functionorthogonalize_left!(A::AbstractTensorTrain; svd_trunc::SVDTrunc, indices::UnitRange)Bring A to left-orthogonal form by means of SVD decompositions.
Optionally perform truncations by passing a SVDTrunc. Optionally pass a range of indices to perform the orthogonalizations only on those.
TensorTrains.orthogonalize_two_site_center! — Functionorthogonalize_two_site_center!(C::TensorTrain, k::Integer; svd_trunc=TruncThresh(1e-6))Orthogonalize the tensor train for a two-site DMRG update at positions k and k+1. This puts sites 1:k-1 in left-canonical form, sites k+2:N in right-canonical form, and leaves sites k and k+1 as the non-orthogonal center for merging.
TensorTrains.is_two_site_canonical — Functionis_two_site_canonical(A, k)Check if the tensor train is in canonical form for a two-site update at positions k and k+1. This means sites 1:k-1 are left-canonical, sites k+2:N are right-canonical, and sites k and k+1 are non-canonical (ready for merging).
TensorTrains.compress! — Functioncompress!(A::AbstractTensorTrain; svd_trunc::SVDTrunc)Compress A by means of SVD decompositions + truncations
Uniform Tensor Trains
TensorTrains.AbstractUniformTensorTrain — TypeAbstractUniformTensorTrain{F,N} <: AbstractPeriodicTensorTrain{F,N}An abstract type for representing tensor trains with periodic boundary conditions and all matrices equal
Fis the type of the matrix entriesNis the number of indices of each tensor (2 virtual ones +N-2physical ones)
TensorTrains.UniformTensorTrain — TypeUniformTensorTrain{F<:Number, N} <: AbstractUniformTensorTrain{F,N}A type for representing a tensor train with periodic boundary conditions, all matrices equal and a fixed length L.
Fis the type of the matrix entriesNis the number of indices of each tensor (2 virtual ones +N-2physical ones)
FIELDS
tensoronly one is storedLthe length of the tensor trainza re-scaling constant
TensorTrains.InfiniteUniformTensorTrain — TypeInfiniteUniformTensorTrain{F<:Number, N} <: AbstractUniformTensorTrain{F,N}A type for representing an infinite tensor train with periodic boundary conditions and all matrices equal.
Fis the type of the matrix entriesNis the number of indices of each tensor (2 virtual ones +N-2physical ones)
FIELDS
tensoronly one is storedza re-scaling constant
TensorTrains.symmetrized_uniform_tensor_train — Functionsymmetrized_uniform_tensor_train(A::AbstractTensorTrain)Convert a tensor train $f(x^1,x^2,\ldots,x^L)$ into a UniformTensorTrain $g$ such that
\[g(\mathcal P[x^1,x^2,\ldots,x^L]) = g(x^1,x^2,\ldots,x^L) = \sum_{l=1}^L f(x^l, x^{l+1},\ldots,x^L,x^1,\ldots,x^{l-1})\]
for any cyclic permutation $\mathcal P$
TensorTrains.periodic_tensor_train — Functionperiodic_tensor_train(A::UniformTensorTrain)Produce a PeriodicTensorTrain corresponding to A, with the matrix concretely repeated length(A) times
Infinite Tensor Trains
TensorTrains.TruncVUMPS — TypeTruncVUMPS{TI, TF} <: SVDTruncA type used to perform truncations of an InfiniteUniformTensorTrain to a target bond size d. It uses the Variational Uniform Matrix Product States (VUMPS) algorithm from MPSKit.jl.
FIELDS
d: target bond dimensionmaxiter = 100: max number of iterations for the VUMPS algorithmtol = 1e-14: tolerance for the VUMPS algorithm
p = rand_infinite_uniform_tt(10, 2, 2)
compress!(p, TruncVUMPS(5))Truncators
TensorTrains.SVDTrunc — Typeabstract type SVDTruncSVD truncator. Can be threshold-based or bond size-based
TensorTrains.TruncThresh — TypeTruncThresh{T} <: SVDTruncA type used to perform SVD-based truncations based on a threshold ε. Given a vector $\lambda$ of $m$ singular values, those below $\varepsilon\sqrt{\sum_{k=1}^m \lambda_k^2}$ are truncated to zero.
FIELDS
ε: threshold.
svd_trunc = TruncThresh(1e-5)
M = rand(5,6)
M_trunc = svd_trunc(M)TensorTrains.TruncBond — TypeTruncBond{T} <: SVDTruncA type used to perform SVD-based truncations based on bond size m'. Given a vector $\lambda$ of $m$ singular values, only the $m'$ largest are kept, the others are truncated to zero.
FIELDS
mprime: number of singular values to retain
svd_trunc = TruncBond(3)
M = rand(5,6)
M_trunc = svd_trunc(M)TensorTrains.TruncBondMax — TypeTruncBondMax{T} <: SVDTruncSimilar to TruncBond, but also stores the maximum error $\sqrt{\frac{\sum_{k=m'+1}^m\lambda_k^2}{\sum_{k=1}^m\lambda_k^2}}$ made since the creation of the object
FIELDS
mprime: number of singular values to retainmaxerr: a 1-element vector storing the maximum error
TensorTrains.TruncBondThresh — TypeTruncBondThresh{T} <: SVDTruncA mixture of TruncBond and TruncThresh, truncates to the most stringent criterion.
Gradients
TensorTrains.grad_squareloss_two_site — Functiongrad_squareloss_two_site(ψ::TensorTrain, k::Integer, X, Y) -> grad_sl, slCompute the gradient of the square loss from fitting data (X,Y) with the tensor train ψ with respect to the merged tensors Aᵏ and Aᵏ⁺¹. Return also the loss, which is a byproduct of the computation.
DMRG-like optimization
TensorTrains.two_site_dmrg! — Functiontwo_site_dmrg!(p::MPS, X, nsweeps; weights, kw...)Fit a MPS to data X using a MPS ansatz and the 2site-DMRG-like gradient descent.
two_site_dmrg!(p::TensorTrain, X, Y, nsweeps; kw...)Perform regression on data (X, Y) using a tensor train ansatz and the 2site-DMRG-like gradient descent.
Matrix Product States
Can be accessed, for now, via
using TensorTrains.MatrixProductStatesTensorTrains.MatrixProductStates.MPS — TypeMPS{T<:AbstractTensorTrain}With a little abuse of nomenclature, a type for representing a probability distribution whose value is the square module of the value of the contained tensor train.
FIELDS
ψ: anAbstractTensorTrain
Example:
L = 3
q = (2, 3)
ψ = rand_mps(4, L, q...)
p = MPS(ψ)
X = [[rand(1:qi) for qi in q] for l in 1:L]
evaluate(p, X), abs2(evaluate(ψ, X)) # are the sameTensorTrains.MatrixProductStates.rand_mps — Functionrand_mps([rng=default_rng()], [T = Float64], bondsizes::AbstractVector{<:Integer}, q...)
rand_mps([rng=default_rng()], [T = Float64], d::Integer, L::Integer, q...)Construct a Matrix Product States with rand(T) entries, by specifying either:
bondsizes: the size of each bondda fixed size for all bonds,Lthe length
and
qa Tuple/Vector specifying the number of values taken by each variable on a single site
StatsAPI.loglikelihood — FunctionStatsBase.loglikelihood(p::MPS, X)Compute the loglikelihood of the data X under the MPS distribution p.
TensorTrains.MatrixProductStates.grad_loglikelihood_two_site — Functiongrad_loglikelihood_two_site(p::MPS, k::Integer, X) -> grad_ll, llCompute the gradient of the loglikelihood of data X under the MPS distribution p with respect to the merged tensors Aᵏ and Aᵏ⁺¹. Return also the loglikelihood, which is a byproduct of the computation.
TensorTrains.MatrixProductStates.grad_normalization_two_site_canonical — Functiongrad_normalization_two_site_canonical(p::MPS, k::Integer) -> gradz, zCompute the gradient of the normalization of p with respect to the merged tensors Aᵏ and Aᵏ⁺¹. Return also the normalization, which is a byproduct of the computation.