MERA.jl reference
MERA.jl provides Julia implementations of some basic Multiscale Entaglement Renormalization Ansatz algorithms. For usage instructions, see the GitHub page. Below you can find the reference documentation, listing all the types and functions.
MERA types
MERA.GenericMERA
— TypeGenericMERA{N, LT <: Layer, OT}
A GenericMERA
is a collection of Layer
s. The type of these layers then determines whether the MERA is binary, ternary, etc.
On conventions and terminology:
- The physical indices of the MERA are at the "bottom", the scale invariant part at the "top".
- The counting of layers starts from the bottom, so the layer with physical indices is layer 1. The last layer is the scale invariant one, that then repeats upwards to infinity.
- Each layer is thought of as a linear map from its top, or input space to its bottom, or output space.
The type parameters are N
: The number of distinct layers (N-1 transition layers and one scale invariant one). LT
: Layer type. OT
: Operator type. The type of ascended and descended operators for this MERA. Determined from LT
. Typically a TensorMap
with input and output indices matching the causal cone width.
GenericMERA
is immutable, and the layers can not be changed after construction. All functions that modify a GenericMERA
return new objects.
MERA.TernaryMERA
— TypeTernaryMERA{N}
A ternary MERA is a MERA consisting of TernaryLayer
s.
MERA.BinaryMERA
— TypeBinaryMERA{N}
A binary MERA is a MERA consisting of BinaryLayer
s.
MERA.ModifiedBinaryMERA
— TypeModifiedBinaryMERA{N}
A modified binary MERA is a MERA consisting of ModifiedBinaryLayer
s.
Layer types
MERA.Layer
— TypeAbstract supertype of all layer types, e.g. BinaryLayer
and TernaryLayer
.
MERA.SimpleLayer
— TypeSimpleLayer <: Layer
A SimpleLayer
is a MERA layer that consists of a set of isometries and/or unitaries, and nothing else. This allows writing convenient generic versions of many methods, reducing code duplication for the concrete Layer
types. Every subtype of SimpleLayer
should implement the iteration and indexing interfaces to return the various tensors of the layer in the same order in which the constructor takes them in.
MERA.TernaryLayer
— TypeTernaryLayer{ST, ET, Tan} <: SimpleLayer
The type for layers of a ternary MERA.
Each layer consists of two tensors, a 2-to-2 disentangler, often called u
, and a 3-to-1 isometry, often called w
.
The type parameters are ST
for space type, e.g. ComplexSpace
or SU2Space
; ET
for element type, e.g. Complex{Float64}
; and Tan
for whether this layer is a tangent layer or not. If Tan = false
, the layer is question is an actual MERA layer. If Tan = true
it consists, instead of the actual tensors, of Stiefel/Grassmann tangent vectors of these tensors.
Index numbering convention is as follows, where the physical indices are at the bottom: Disentangler:
3| 4|
+------+
| u |
+------+
1| 2|
Isometry:
4|
+-------+
| w |
+-------+
1| 2| 3|
MERA.BinaryLayer
— TypeBinaryLayer{ST, ET, Tan} <: SimpleLayer
The type for layers of a binary MERA.
Each layer consists of two tensors, a 2-to-2 disentangler, often called u
, and a 2-to-1 isometry, often called w
.
The type parameters are ST
for space type, e.g. ComplexSpace
or SU2Space
; ET
for element type, e.g. Complex{Float64}
; and Tan
for whether this layer is a tangent layer or not. If Tan = false
, the layer is question is an actual MERA layer. If Tan = true
it consists, instead of the actual tensors, of Stiefel/Grassmann tangent vectors of these tensors.
Index numbering convention is as follows, where the physical indices are at the bottom: Disentangler:
3| 4|
+------+
| u |
+------+
1| 2|
Isometry:
3|
+------+
| w |
+------+
1| 2|
MERA.ModifiedBinaryLayer
— TypeModifiedBinaryLayer{ST, ET, Tan} <: SimpleLayer
The type for layers of a modified binary MERA.
Each layer consists of three tensors, a 2-to-2 disentangler, often called u
, and two 2-to-1 isometries, often called wl
and wr
, for left and right. Their relative locations are
| |
wl wr
| \ / |
| u |
| / \ |
The type parameters are ST
for space type, e.g. ComplexSpace
or SU2Space
; ET
for element type, e.g. Complex{Float64}
; and Tan
for whether this layer is a tangent layer or not. If Tan = false
, the layer is question is an actual MERA layer. If Tan = true
it consists, instead of the actual tensors, of Stiefel/Grassmann tangent vectors of these tensors.
Index numbering convention is as follows, where the physical indices are at the bottom: Disentangler:
3| 4|
+------+
| u |
+------+
1| 2|
Isometries:
3|
+------+
| w |
+------+
1| 2|
Utility functions
MERA.num_translayers
— Functionnum_translayers(m::GenericMERA)
Return the number of transition layers, i.e. layers below the scale invariant one, in the MERA.
MERA.layertype
— Functionlayertype(m::GenericMERA)
layertype(::Type{<:GenericMERA})
Return the type of the layers of the MERA.
MERA.operatortype
— Functionoperatortype(m::GenericMERA)
operatortype(::Type{<: GenericMERA})
Return the type of operator associate with this MERA or MERA type. That means the type of operator that fits in the causal cone, and is naturally emerges as one ascends local operators.
MERA.scalefactor
— Functionscalefactor(::Type{<: GenericMERA})
The ratio by which the number of sites changes when one descends by one layer, e.g. 2 for binary MERA, 3 for ternary.
MERA.causal_cone_width
— Functioncausal_cone_width(::Type{<: GenericMERA})
Return the width of the stable causal cone for this MERA type.
MERA.get_layer
— Functionget_layer(m::GenericMERA, depth)
Return the layer at the given depth. 1 is the lowest layer, i.e. the one with physical indices.
MERA.outputspace
— Functionoutputspace(layer::Layer)
Return the vector space of the downwards-pointing (towards the physical level) indices of `layer.
See also: inputspace
, internalspace
MERA.inputspace
— Functioninputspace(layer::Layer)
Return the vector space of the upwards-pointing (towards the scale invariance) indices of `layer.
See also: outputspace
, internalspace
MERA.internalspace
— Functioninternalspace(layer::Layer)
Return the internal vector space of `layer.
See also: outputspace
, inputspace
MERA.pseudoserialize
— Functionpseudoserialize(x)
Return a tuple of objects that can be used to reconstruct x
, and that are all of Julia base types.
The name refers to how this isn't quite serialization, since it doesn't break objects down to bit strings, but kinda serves the same purpose: pseudoserialised objects can easily be written and read from e.g. disk, without having to worry about quirks of the type system (like with JLD).
See also: depseudoserialize
MERA.depseudoserialize
— Functiondepseudoserialize(::Type{T}, args) where T <: GenericMERA
Reconstruct an object given the output of pseudoserialize
: x -> depseudoserialize(pseudoserialize(x)...)
should be an effective noop.
See also: pseudoserialize
MERA.remove_symmetry
— Functionremove_symmetry(V)
Strip a vector space of its symmetry structure, i.e. return the corresponding ℂ^n
or ℝ^n
.
remove_symmetry(t::AbstractTensorMap)
Strip an AbstractTensorMap
of its internal symmetries, and return the corresponding TensorMap
that operators on ComplexSpace
or CartesianSpace
.
remove_symmetry(m::GenericMERA)
Given a MERA which may possibly be built of symmetry preserving TensorMap
s, return another, equivalent MERA that has the symmetry structure stripped from it, and all tensors are dense.
TensorKitManifolds.projectisometric
— Functionprojectisometric(m::GenericMERA)
Project all the tensors of the MERA to respect the isometricity condition.
TensorKitManifolds.projectisometric!
— Functionprojectisometric!(m::GenericMERA)
Project all the tensors of the MERA to respect the isometricity condition, in place.
MERA.reset_storage
— Functionreset_storage(m::GenericMERA)
Reset cached operators, so that they will be recomputed when they are needed.
Generating and modifying MERAs
MERA.replace_layer
— Functionreplace_layer(c::MERACache, depth)
Create a new MERACache
that has all the stored pieces removed that are invalidated by changing the layer at depth
.
replace_layer(m::GenericMERA, layer, depth; check_invar=true)
Replace depth
layer of m
with layer
. If check_invar=true, check that the indices match afterwards.
MERA.release_transitionlayer
— Functionrelease_transitionlayer(m::GenericMERA)
Add one more transition layer at the top of the MERA, by taking the lowest of the scale invariant ones and releasing it to vary independently.
MERA.random_MERA
— Functionrandom_MERA(::Type{T <: GenericMERA}, ET, Vouts, Vints=Vouts; kwargs...)
Generate a random MERA of type T
.
ET
is the element type of the MERA, e.g. Float64
. Vouts
are the vector spaces between the various layers. Number of layers will be the length of Vouts
. Vouts[1]
will be the physical index space, and Vs[end]
will be the one at the scale invariant layer. Vints
can be a vector/tuple of the same length as Vouts
, and include vector spaces for the internal indices of each layer. Alternatively, it can hold any other extra parameters specific to each layer, as it is simply passed on to the function randomlayer
. Also passed to randomlayer
will be any additional keyword arguments, but these will all be the same for each layer.
See also: randomlayer
MERA.randomlayer
— Functionrandomlayer(::Type{T <: Layer}, T, Vin, Vout, Vint=Vout; random_disentangler=false)
Return a MERA layer with random tensors.
T
is the Layer
type, and Vin
and Vout
are the input and output spaces. Vint
is an internal vector space for the layer, connecting the disentanglers to the isometries. If random_disentangler=true
, the disentangler is also a random unitary, if false
(default), it is the identity or the product of two single-site isometries, depending on if the disentanler is supposed to be unitary or isometric.
Each subtype of Layer
should have its own method for this function.
MERA.expand_bonddim
— Functionexpand_bonddim(m::GenericMERA, depth, newdims; check_invar=true)
Expand the bond dimension of the MERA at the given depth.
The indices to expand are the input indices of the layer at depth
, i.e. depth=1
means the lowest virtual indices. The new bond dimension is given by newdims
, which for a non-symmetric MERA is just a number, and for a symmetric MERA is a dictionary of irrep => block dimension
. Not all irreps for a bond need to be listed, the ones left out are left untouched.
The expansion is done by padding tensors with zeros. Note that this breaks isometricity of the individual tensors. This is however of no consequence, since the MERA as a state remains exactly the same. A round of optimization on the MERA will restore isometricity of each tensor, or projectisometric
can be called to do so explicitly.
If check_invar = true
the function checks that the bond dimensions of various layers match after the expansion.
MERA.expand_internal_bonddim
— Functionexpand_internal_bonddim(m::GenericMERA, depth, newdims; check_invar=true)
Expand the bond dimension of the layer-internal indices of the MERA at the given depth.
The new bond dimension is given by newdims
, which for a non-symmetric MERA is just a number, and for a symmetric MERA is a dictionary of {irrep => block dimension}. Not all irreps for a bond need to be listed, the ones left out are left untouched.
The expansion is done by padding tensors with zeros. Note that this breaks isometricity of the individual tensors. This is however of no consequence, since the MERA as a state remains exactly the same. A round of optimization on the MERA will restore isometricity of each tensor, or projectisometric
can be called to do so explicitly.
Note that not all MERAs have an internal bond dimension, and some may have several, so this function will not make sense for all MERA types. Implementation relies on the function expand_internalspace
, defined for each Layer
type.
Measuring physical quantities
MERA.expect
— Functionexpect(op, m::GenericMERA, pars=(;), opscale=1, evalscale=1)
Return the expecation value of operator op
for the MERA m
.
The layer on which op
lives is set by opscale
, which by default is the physical one. evalscale
can be used to set whether the operator is ascended through the network or the density matrix is descended. pars
can hold additional options that are further down the line passed on to fixedpoint_densitymatrix
.
See also: fixedpoint_densitymatrix
MERA.scalingdimensions
— Functionscalingdimensions(m::GenericMERA, howmany=20)
Diagonalize the scale invariant ascending superoperator to compute the scaling dimensions of the underlying CFT.
The return value is a dictionary, the keys of which are symmetry sectors for a possible internal symmetry of the MERA (Trivial() if there is no internal symmetry), and values are scaling dimensions in this symmetry sector.
howmany
controls how many of lowest scaling dimensions are computed.
MERA.densitymatrix_entropies
— Functiondensitymatrix_entropies(m::GenericMERA)
Return a vector of entropies for the density matrices in the MERA. The first one is for the density matrix at the physical level, the last one is the scale invariant density matrix.
Ascending operators
MERA.ascend
— Functionascend(op, layer::Layer)
Ascend a local operator op
from the bottom of layer
to the top.
MERA.ascended_operator
— Functionascended_operator(m::GenericMERA, op, depth)
Return the operator op
ascended from the physical level to depth
.
This function utilises the cache, to avoid recomputation.
See also: scale_invariant_operator_sum
MERA.scale_invariant_operator_sum
— Functionscale_invariant_operator_sum(m::GenericMERA, op, pars)
Return the sum of the ascended versions of op
in the scale invariant part of the MERA.
To be more precise, this sum is of course infinite, and what we return is the component of it orthogonal to the dominant eigenoperator of the ascending superoperator (typically the identity). This component converges to a finite result like a geometric series, since all non-dominant eigenvalues of the ascending superoperator are smaller than 1.
To approximate the converging series, we use an iterative Krylov solver. The options for the solver should be in pars.scaleinvariant_krylovoptions
, they will be passed to KrylovKit.linsolve
.
Descending density matrices
MERA.descend
— Functiondescend(rho, layer::Layer)
Descend a local density matrix rho
from the top of layer
to the bottom.
MERA.densitymatrix
— Functiondensitymatrix(m::GenericMERA, depth, pars=(;))
Return the density matrix right below the layer at depth
.
pars
maybe a NamedTuple
of options, passed on to the function fixedpoint_densitymatrix
.
This function utilises the cache, to avoid recomputation.
See also: fixedpoint_densitymatrix
, densitymatrices
MERA.densitymatrices
— Functiondensitymatrices(m::GenericMERA, pars=(;))
Return all the distinct density matrices of the MERA, starting with the one at the physical level, and ending with the scale invariant one.
pars
maybe a NamedTuple
of options, passed on to the function fixedpoint_densitymatrix
.
See also: fixedpoint_densitymatrix
, densitymatrix
MERA.fixedpoint_densitymatrix
— Functionfixedpoint_densitymatrix(m::GenericMERA, pars::NamedTuple=(;))
Find the fixed point density matrix of the scale invariant part of the MERA.
To find the fixed point, we use an iterative Krylov solver. The options for the solver should be in pars.scaleinvariant_krylovoptions
, they will be passed to KrylovKit.eigsolve
.
Optimizing a MERA
MERA.minimize_expectation
— Functionminimize_expectation(m::GenericMERA, h, pars=(;);
finalize! = OptimKit._finalize!, vary_disentanglers=true,
kwargs...)
Return a MERA optimized to minimize the expectation value of operator h
, starting with m
as the initial guess.
pars
is a NamedTuple
of parameters for the optimisation. They are,
method
: ASymbol
that chooses which optimisation method to use. Options are:lbfgs
for L-BFGS (default),:ev
for Evenbly-Vidal,:cg
for conjugate gradient, and:gd
for gradient descent.:lbfgs
,:cg
, and:gd
are together known as the gradient methods.maxiter
: Maximum number of iterations to use. 2000 by default.gradient_delta
: Convergence threshold, as measured by the norm of the gradient.1e-14
by default.precondition
: Whether to apply preconditioning with the physical Hilbert space inner product.true
by default. See https://arxiv.org/abs/2007.03638 for more details.verbosity
: How much output to log. 2 by default.isometries_only_iters
: An integer for how many iterations should at first be done optimising only the isometries, and leaving the disentangler be. 0 by default.scaleinvariant_krylovoptions
: ANamedTuple
of keyword arguments passed toKrylovKit.linsolve
andKrylovKit.eigsolve
, when solving for the fixed-point density matrix and the scale invariant operator sum. The default is(tol = 1e-13, verbosity = 0, maxiter = 20)
.retraction
: Which retraction method to use. Options are:exp
for geodesics (default), and:cayley
for Cayley transforms. Only affects gradient methods.- transport: Which vector transport method to use. Currently each retraction
method only comes with a single compatible transport, so one should always use
transportto be the same as
retraction`. This may change. Only affects gradient methods. metric
: Which metric to use for Stiefel manifold. Options are:euclidean
(default) and:canonical
. Only affects gradient methods.ls_epsilon
: The ϵ parameter for the Hager-Zhang line search.1e-6
be default. Only affects gradient methods.lbfgs_m
: The rank of the approximation of the inverse Hessian in L-BFGS. 8 by default. Only affects the:lbfgs
method.cg_flavor
: The "flavor" of conjguate gradient to use.:HagerZhang
by default. Only affects the:cg:
method.ev_layer_iters
: How many times a single layer is optimised before moving to the next layer in the Evenbly-Vidal algorithm.1
by default. Only affects the:ev
method.
If any of these are specified in pars
, the specified values override the defaults.
finalize!
is a function that will be called at every iteration. It can be used to for instance log the development of some quantity during the optimisation, or modify the MERA in some custom (although undefined behavior may follow depending on how the state is changed). Its signature is finalize!(m, f, g, counter)
where m
is the current MERA, f
is its expectation value for h
, g
is the gradient MERA at m
, and counter
is the current iteration number. It should return the possibly modified m
, f
, and g
. If method = :ev
, then it should also be able to handle the case g = nothing
, since the Evenbly-Vidal algorithm does not use gradients.
vary_disentanglers
gives the option of running the optimisation but only for the isometries, leaving the disentanglers as they are.
MERA.environment
— Functionenvironment(layer::Layer, op, rho; vary_disentanglers=true)
Compute the environments with respect to op
of all the tensors in the layer, and return them as a Layer
. rho
is the local density matrix at the top indices of this layer.
If vary_disentanglers=false
, only compute the environments for the isometries, and set the environments for the disentanglers to zero.
MERA.gradient
— Functiongradient(h, m::GenericMERA, pars::NamedTuple; vary_disentanglers=true)
Compute the gradient of the expectation value of h
at the point m
.
pars
should have pars.metric
that specifies whether to use the :euclidean
or :canonical
metric for the Stiefel manifold, and pars.scaleinvariant_krylovoptions
that is passed on to environment
.
vary_disentanglers
allows computing the gradients only for the isometries, and setting the gradients for the disentanglers to zero.
The return value is a "tangent MERA": An object of a similar type as m
, but instead of regular layers with tensors that have isometricity constraints, instead each layer holds the corresponding gradients for each tensor.
TensorKitManifolds.retract
— Functionretract(m::GenericMERA, mtan::GenericMERA, alpha::Real; kwargs...)
Given a "tangent MERA" mtan
, at base point m
, retract in the direction of mtan
by distance alpha
. This is done tensor-by-tensor, i.e. each tensor is retracted along its respective Stiefel/Grassmann tangent.
The additional keyword argument are passed on to the respective TensorKitManifolds
function.
See TensorKitManifolds.retract
for more details.
See also: transport!
TensorKitManifolds.transport!
— Functiontransport!(mvec::GenericMERA, m::GenericMERA, mtan::GenericMERA, alpha::Real,
mend::GenericMERA; kwargs...)
Given a "tangent MERAs" mtan
and mvec
, at base point m
, transport mvec
in the direction of mtan
by distance alpha
. This is done tensor-by-tensor, i.e. each tensor is transported along its respective Stiefel/Grassmann tangent.
mend
is the endpoint on the manifold, i.e. the result of retracting by alpha
in the direction of mtan
. The additional keyword argument are passed on to the respective TensorKitManifolds
function.
See TensorKitManifolds.transport!
for more details.
See also: retract
TensorKitManifolds.inner
— Functioninner(m::GenericMERA, m1::GenericMERA, m2::GenericMERA; metric=:euclidean)
Given two tangent MERAs m1
and m2
, both at base point m
, compute their inner product. This means the sum of the inner products of the individual tensors.
See TensorKitManifolds.inner
for more details.
MERA.tensorwise_sum
— Functiontensorwise_sum(m1::T, m2::T) where T <: GenericMERA
Return a MERA for which each tensor is the sum of the corresponding tensors of m1
and m2
.
MERA.tensorwise_scale
— Functiontensorwise_scale(m::GenericMERA, alpha::Number)
Scale all the tensors of m
by alpha
.
Index
MERA.BinaryLayer
MERA.BinaryMERA
MERA.GenericMERA
MERA.Layer
MERA.ModifiedBinaryLayer
MERA.ModifiedBinaryMERA
MERA.SimpleLayer
MERA.TernaryLayer
MERA.TernaryMERA
MERA.ascend
MERA.ascended_operator
MERA.causal_cone_width
MERA.densitymatrices
MERA.densitymatrix
MERA.densitymatrix_entropies
MERA.depseudoserialize
MERA.descend
MERA.environment
MERA.expand_bonddim
MERA.expand_internal_bonddim
MERA.expect
MERA.fixedpoint_densitymatrix
MERA.get_layer
MERA.gradient
MERA.inputspace
MERA.internalspace
MERA.layertype
MERA.minimize_expectation
MERA.num_translayers
MERA.operatortype
MERA.outputspace
MERA.pseudoserialize
MERA.random_MERA
MERA.randomlayer
MERA.release_transitionlayer
MERA.remove_symmetry
MERA.replace_layer
MERA.reset_storage
MERA.scale_invariant_operator_sum
MERA.scalefactor
MERA.scalingdimensions
MERA.tensorwise_scale
MERA.tensorwise_sum
TensorKitManifolds.inner
TensorKitManifolds.projectisometric
TensorKitManifolds.projectisometric!
TensorKitManifolds.retract
TensorKitManifolds.transport!