Doc-Strings
CTMRGRuntime{LT}

a struct to hold the tensors during the ctmrg algorithm, containing

  • D × D × D × Dbulk tensor
  • χ × χcorner tensor
  • χ × D × χedge tensor

and LT is a AbstractLattice to define the lattice type.

Heisenberg(Jz::T,Jx::T,Jy::T) where {T<:Real}

return a struct representing the heisenberg model with magnetisation fields Jz, Jx and Jy..

IPEPS{LT<:AbstractLattice, T, N}

Infinite projected entangled pair of states. LT is the type of lattice, T and N are bulk tensor element type and order.

SquareCTMRGRuntime(bulk::AbstractArray{T,4}, env::Val, χ::Int)

create a SquareCTMRGRuntime with bulk-tensor bulk. The corner and edge tensors are initialized according to env. If env = Val(:random), the corner is initialized as a random χ×χ tensor and the edge is initialized as a random χ×D×χ tensor where D = size(bulk,1). If env = Val(:raw), corner- and edge-tensor are initialized by summing over one or two indices of bulk respectively and embedding the result in zeros-tensors of the appropriate size, truncating if necessary.

example

julia> rt = SquareCTMRGRuntime(randn(2,2,2,2), Val(:raw), 4);

julia> rt.corner[1:2,1:2] ≈ dropdims(sum(rt.bulk, dims = (3,4)), dims = (3,4))
true

julia> rt.edge[1:2,1:2,1:2] ≈ dropdims(sum(rt.bulk, dims = 4), dims = 4)
true
TFIsing(hx::Real)

return a struct representing the transverse field ising model with magnetisation hx.

ctmrg(rt::CTMRGRuntime; tol, maxit)

return a CTMRGRuntime with an environment consisting of corner and edge tensor that have either been iterated for maxit iterations or converged according to tol. Convergence is tested by looking at the sum of the absolut differences in the corner singular values. If it is less than tol, convergence is reached.

example

julia> a = model_tensor(Ising(),β);

julia> rt = SquareCTMRGRuntime(a, Val(:random), χ);

julia> env = ctmrg(rt; tol=1e-6, maxit=100);

for the environment of an isingmodel at inverse temperature β on an infinite square lattice.

hamiltonian(model<:HamiltonianModel)

return the hamiltonian of the model as a two-site tensor operator.

hamiltonian(model::Heisenberg)

return the heisenberg hamiltonian for the model as a two-site operator.

hamiltonian(model::TFIsing)

return the transverse field ising hamiltonian for the provided model as a two-site operator.

mag_tensor(::Ising,β)

return the operator for the magnetisation at inverse temperature β at a site in the two-dimensional ising model on a square lattice in tensor-network form.

model_tensor(::Ising,β)

return the isingtensor at inverse temperature β for a two-dimensional square lattice tensor-network.

num_grad(f, K::AbstractArray; [δ = 1e-5])

return the numerical gradient of f for each element of K.

example

julia> TensorNetworkAD.num_grad(tr, rand(2,2)) ≈ I
true
num_grad(f, K::Real; [δ = 1e-5])

return the numerical gradient of f at K calculated with (f(K+δ/2) - f(K-δ/2))/δ

example

julia> TensorNetworkAD.num_grad(x -> x * x, 3) ≈ 6
true
optimiseipeps(ipeps, h; χ, tol, maxit, optimargs = (), optimmethod = LBFGS(m = 20))

return the tensor bulk' that describes an ipeps that minimises the energy of the two-site hamiltonian h. The minimization is done using Optim with default-method LBFGS. Alternative methods can be specified by loading LineSearches and providing optimmethod. Other options to optim can be passed with optimargs. The energy is calculated using ctmrg with parameters χ, tol and maxit.

trg(a, χ, niter)

return the partition-function of a two-dimensional system of size 2^niter described by the tensor a calculated via the tensor renormalization group algorithm. a is a rank-4 tensor with the following indices:

    |1
4--[a]--2
   3|
(st::StopFunction)(state)

stopfunction for ctmrg, returning true if singular values are converged or the maximum number of iterations is reached. Implemented as a closure since it needs to remember the last singular values it saw for comparison.

ctmrgstep(rt,vals)

evaluate one step of the ctmrg-algorithm, returning a tuple of an updated CTMRGRuntime with updated corner and edge tensor and a vector of singular values to test convergence with.

diaglocalhamiltonian(diag::Vector)

return the 2-site Hamiltonian with single-body terms given by the diagonal diag.

energy(h, ipeps; χ, tol, maxit)

return the energy of the ipeps 2-site hamiltonian h and calculated via a ctmrg with parameters χ, tol and maxit.

expectationvalue(h, ap, rt)

return the expectationvalue of a two-site operator h with the sites described by rank-6 tensor ap each and an environment described by a SquareCTMRGRuntimert.

fixedpoint(f, guess, stopfun)

return the result of applying guess = f(guess) until convergence. Convergence is decided by applying stopfun(guess) which returns a Boolean.

indexperm_symmetrize(ipeps::SquareIPEPS)

return a SquareIPEPS based on ipeps that is symmetric under permutation of its virtual indices.

magnetisation(model<:HamiltonianModel, β, χ)

return the magnetisation of the model as a function of the inverse temperature β and the environment bonddimension χ as calculated with ctmrg. Requires that mag_tensor and model_tensor are defined for model.

magofβ(::Ising,β)

return the analytical result for the magnetisation at inverse temperature β for the 2d classical ising model.

tensorfromclassical(h::Matrix)

given a classical 2-body hamiltonian h, return the corresponding tensor for use in e.g. trg for a two-dimensional square-lattice.

Example

julia> model_tensor(Ising(),β) ≈ tensorfromclassical([β -β; -β β])

true