EditURL = "../../../../examples/bose_hubbard/main.jl"
This example demonstrates the simulation of the two-dimensional Bose-Hubbard model. In particular, the point will be to showcase the use of internal symmetries and finite particle densities in PEPS ground state searches. As we will see, incorporating symmetries into the simulation consists of initializing a symmetric Hamiltonian, PEPS state and CTM environment - made possible through TensorKit.
But first let's seed the RNG and import the required modules:
using Random
using TensorKit, PEPSKit
using MPSKit: add_physical_charge
Random.seed!(2928528935);We will construct the Bose-Hubbard model Hamiltonian through the
bose_hubbard_model,
function from MPSKitModels as reexported by PEPSKit. We'll simulate the model in its
Mott-insulating phase where the ratio
t = 1.0
U = 30.0
cutoff = 2
mu = 0.0
lattice = InfiniteSquare(1, 1);Next, we impose an explicit global symmetry argument of the
Hamiltonian constructor to U1Irrep and passing one as the particle number density
keyword argument n:
symmetry = U1Irrep
n = 1
H = bose_hubbard_model(ComplexF64, symmetry, lattice; cutoff, t, U, n);Before we continue, it might be interesting to inspect the corresponding lattice physical
spaces (which is here just a
physical_spaces = physicalspace(H)1×1 Matrix{TensorKit.GradedSpace{TensorKitSectors.U1Irrep, TensorKit.SortedVectorDict{TensorKitSectors.U1Irrep, Int64}}}:
(0 => 1, 1 => 1, -1 => 1)
Note that the physical space contains
When running PEPS simulations with explicit internal symmetries, specifying the structure of
the virtual spaces of the PEPS and its environment becomes a bit more involved. For the
environment, one could in principle allow the virtual space to be chosen dynamically during
the boundary contraction using CTMRG by using a truncation scheme that allows for this
(e.g. using alg=:truncrank or alg=:trunctol to truncate to a fixed total bond dimension
or singular value cutoff respectively). For the PEPS virtual space however, the structure
has to be specified before the optimization.
While there are a host of techniques to do this in an informed way (e.g. starting from a
simple update result), here we just specify the virtual space manually. Since we're dealing
with a model at unit filling our physical space only contains integer
V_peps = U1Space(0 => 2, 1 => 1, -1 => 1)
V_env = U1Space(0 => 6, 1 => 4, -1 => 4, 2 => 2, -2 => 2);Having defined our Hamiltonian and spaces, it is just a matter of plugging this into the optimization framework in the usual way to find the ground state. So, we first specify all algorithms and their tolerances:
boundary_alg = (; tol = 1.0e-8, alg = :simultaneous, trunc = (; alg = :fixedspace))
gradient_alg = (; tol = 1.0e-6, maxiter = 10, alg = :linsolver, iterscheme = :fixed)
optimizer_alg = (; tol = 1.0e-4, alg = :lbfgs, maxiter = 150, ls_maxiter = 2, ls_maxfg = 2);!!! note
Taking CTMRG gradients and optimizing symmetric tensors tends to be more problematic
than with dense tensors. In particular, this means that one frequently needs to tweak
the boundary_alg, gradient_alg and optimizer_alg settings. There rarely is a
general-purpose set of settings which will always work, so instead one has to adjust
the simulation settings for each specific application. For example, it might help to
switch between the CTMRG flavors alg=:simultaneous and alg=:sequential to
improve convergence. The evaluation of the CTMRG gradient can be instable, so there it
is advised to try the different iterscheme=:diffgauge and iterscheme=:fixed schemes
as well as different alg keywords. Of course the tolerances of the algorithms and
their subalgorithms also have to be compatible. For more details on the available
options, see the fixedpoint docstring.
Keep in mind that the PEPS is constructed from a unit cell of spaces, so we have to make a
matrix of V_peps spaces:
virtual_spaces = fill(V_peps, size(lattice)...)
peps₀ = InfinitePEPS(randn, ComplexF64, physical_spaces, virtual_spaces)
env₀, = leading_boundary(CTMRGEnv(peps₀, V_env), peps₀; boundary_alg...);[ Info: CTMRG init: obj = +1.693461429863e+00 +8.390974048721e-02im err = 1.0000e+00
[ Info: CTMRG conv 19: obj = +1.181834754305e+01 -1.515735205612e-11im err = 3.6943029805e-09 time = 0.24 sec
And at last, we optimize (which might take a bit):
peps, env, E, info = fixedpoint(
H, peps₀, env₀; boundary_alg, gradient_alg, optimizer_alg, verbosity = 3
)
@show E;[ Info: LBFGS: initializing with f = 9.360531870693e+00, ‖∇f‖ = 1.6944e+01
[ Info: LBFGS: iter 1, Δt 8.22 s: f = 1.243263804654e-01, ‖∇f‖ = 6.2838e+00, α = 1.56e+02, m = 0, nfg = 7
[ Info: LBFGS: iter 2, Δt 3.63 s: f = 6.592923309831e-02, ‖∇f‖ = 7.6343e+00, α = 5.34e-01, m = 1, nfg = 2
[ Info: LBFGS: iter 3, Δt 722.6 ms: f = -4.434585356894e-02, ‖∇f‖ = 1.6162e+00, α = 1.00e+00, m = 2, nfg = 1
[ Info: LBFGS: iter 4, Δt 669.8 ms: f = -7.581641308235e-02, ‖∇f‖ = 1.4777e+00, α = 1.00e+00, m = 3, nfg = 1
[ Info: LBFGS: iter 5, Δt 2.22 s: f = -1.270142972187e-01, ‖∇f‖ = 3.1212e+00, α = 5.23e-01, m = 4, nfg = 3
[ Info: LBFGS: iter 6, Δt 652.8 ms: f = -1.633483775361e-01, ‖∇f‖ = 1.2563e+00, α = 1.00e+00, m = 5, nfg = 1
[ Info: LBFGS: iter 7, Δt 667.0 ms: f = -1.937189397513e-01, ‖∇f‖ = 9.6787e-01, α = 1.00e+00, m = 6, nfg = 1
[ Info: LBFGS: iter 8, Δt 1.45 s: f = -2.086095056473e-01, ‖∇f‖ = 7.3708e-01, α = 1.49e-01, m = 7, nfg = 2
[ Info: LBFGS: iter 9, Δt 1.56 s: f = -2.213320160969e-01, ‖∇f‖ = 4.5556e-01, α = 3.29e-01, m = 8, nfg = 2
[ Info: LBFGS: iter 10, Δt 623.9 ms: f = -2.287763259808e-01, ‖∇f‖ = 6.3644e-01, α = 1.00e+00, m = 9, nfg = 1
[ Info: LBFGS: iter 11, Δt 601.4 ms: f = -2.354779361407e-01, ‖∇f‖ = 4.5587e-01, α = 1.00e+00, m = 10, nfg = 1
[ Info: LBFGS: iter 12, Δt 567.9 ms: f = -2.452786754607e-01, ‖∇f‖ = 3.7293e-01, α = 1.00e+00, m = 11, nfg = 1
[ Info: LBFGS: iter 13, Δt 547.8 ms: f = -2.511737944315e-01, ‖∇f‖ = 3.5345e-01, α = 1.00e+00, m = 12, nfg = 1
[ Info: LBFGS: iter 14, Δt 483.4 ms: f = -2.567831916068e-01, ‖∇f‖ = 2.7925e-01, α = 1.00e+00, m = 13, nfg = 1
[ Info: LBFGS: iter 15, Δt 486.6 ms: f = -2.644210470846e-01, ‖∇f‖ = 2.5999e-01, α = 1.00e+00, m = 14, nfg = 1
[ Info: LBFGS: iter 16, Δt 464.1 ms: f = -2.667624388742e-01, ‖∇f‖ = 3.8750e-01, α = 1.00e+00, m = 15, nfg = 1
[ Info: LBFGS: iter 17, Δt 456.4 ms: f = -2.691854381316e-01, ‖∇f‖ = 1.1926e-01, α = 1.00e+00, m = 16, nfg = 1
[ Info: LBFGS: iter 18, Δt 450.9 ms: f = -2.697593549258e-01, ‖∇f‖ = 8.8060e-02, α = 1.00e+00, m = 17, nfg = 1
[ Info: LBFGS: iter 19, Δt 993.8 ms: f = -2.704149848038e-01, ‖∇f‖ = 8.0340e-02, α = 1.00e+00, m = 18, nfg = 1
[ Info: LBFGS: iter 20, Δt 364.3 ms: f = -2.709871859088e-01, ‖∇f‖ = 5.6227e-02, α = 1.00e+00, m = 19, nfg = 1
[ Info: LBFGS: iter 21, Δt 434.9 ms: f = -2.713280444404e-01, ‖∇f‖ = 8.7455e-02, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 22, Δt 413.8 ms: f = -2.716153824055e-01, ‖∇f‖ = 4.4395e-02, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 23, Δt 441.7 ms: f = -2.717554344943e-01, ‖∇f‖ = 3.7465e-02, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 24, Δt 457.9 ms: f = -2.720212512862e-01, ‖∇f‖ = 4.9598e-02, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 25, Δt 463.5 ms: f = -2.721510510078e-01, ‖∇f‖ = 7.7025e-02, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 26, Δt 466.9 ms: f = -2.723111439950e-01, ‖∇f‖ = 3.3781e-02, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 27, Δt 472.3 ms: f = -2.724110171666e-01, ‖∇f‖ = 1.9197e-02, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 28, Δt 462.7 ms: f = -2.724644440119e-01, ‖∇f‖ = 1.9745e-02, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 29, Δt 440.9 ms: f = -2.726363658069e-01, ‖∇f‖ = 3.1087e-02, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 30, Δt 917.0 ms: f = -2.727060629591e-01, ‖∇f‖ = 2.3834e-02, α = 5.02e-01, m = 20, nfg = 2
[ Info: LBFGS: iter 31, Δt 445.2 ms: f = -2.727581806789e-01, ‖∇f‖ = 1.2723e-02, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 32, Δt 573.6 ms: f = -2.727851137084e-01, ‖∇f‖ = 1.3596e-02, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 33, Δt 806.6 ms: f = -2.728202281569e-01, ‖∇f‖ = 1.3577e-02, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 34, Δt 423.9 ms: f = -2.728822332116e-01, ‖∇f‖ = 1.7053e-02, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 35, Δt 928.6 ms: f = -2.729248153492e-01, ‖∇f‖ = 3.1384e-02, α = 4.87e-01, m = 20, nfg = 2
[ Info: LBFGS: iter 36, Δt 477.7 ms: f = -2.729812649791e-01, ‖∇f‖ = 1.4232e-02, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 37, Δt 458.5 ms: f = -2.730056542779e-01, ‖∇f‖ = 9.5005e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 38, Δt 469.2 ms: f = -2.730253098317e-01, ‖∇f‖ = 1.1766e-02, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 39, Δt 493.4 ms: f = -2.730452014699e-01, ‖∇f‖ = 1.3500e-02, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 40, Δt 427.7 ms: f = -2.730601701290e-01, ‖∇f‖ = 6.9352e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 41, Δt 463.8 ms: f = -2.730697005888e-01, ‖∇f‖ = 5.7942e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 42, Δt 464.8 ms: f = -2.730721704336e-01, ‖∇f‖ = 1.0410e-02, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 43, Δt 432.4 ms: f = -2.730768282437e-01, ‖∇f‖ = 5.8856e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 44, Δt 607.3 ms: f = -2.730873408022e-01, ‖∇f‖ = 7.8023e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 45, Δt 795.2 ms: f = -2.730968572102e-01, ‖∇f‖ = 1.0220e-02, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 46, Δt 422.7 ms: f = -2.731119154643e-01, ‖∇f‖ = 1.0811e-02, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 47, Δt 808.0 ms: f = -2.731222065641e-01, ‖∇f‖ = 1.4362e-02, α = 4.12e-01, m = 20, nfg = 2
[ Info: LBFGS: iter 48, Δt 430.5 ms: f = -2.731365310275e-01, ‖∇f‖ = 6.3727e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 49, Δt 447.2 ms: f = -2.731445175618e-01, ‖∇f‖ = 8.2368e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 50, Δt 454.2 ms: f = -2.731516662222e-01, ‖∇f‖ = 9.6146e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 51, Δt 423.9 ms: f = -2.731573199914e-01, ‖∇f‖ = 8.4020e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 52, Δt 445.3 ms: f = -2.731616184908e-01, ‖∇f‖ = 3.3090e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 53, Δt 423.3 ms: f = -2.731636905345e-01, ‖∇f‖ = 3.6654e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 54, Δt 444.9 ms: f = -2.731657365792e-01, ‖∇f‖ = 4.0111e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 55, Δt 876.4 ms: f = -2.731675789298e-01, ‖∇f‖ = 6.3344e-03, α = 5.06e-01, m = 20, nfg = 2
[ Info: LBFGS: iter 56, Δt 967.8 ms: f = -2.731704113587e-01, ‖∇f‖ = 3.5012e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 57, Δt 357.5 ms: f = -2.731731471721e-01, ‖∇f‖ = 2.9705e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 58, Δt 419.5 ms: f = -2.731756296229e-01, ‖∇f‖ = 4.3506e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 59, Δt 404.0 ms: f = -2.731771933453e-01, ‖∇f‖ = 8.4929e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 60, Δt 424.1 ms: f = -2.731802824513e-01, ‖∇f‖ = 3.5372e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 61, Δt 425.9 ms: f = -2.731825100077e-01, ‖∇f‖ = 3.1938e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 62, Δt 446.7 ms: f = -2.731849162403e-01, ‖∇f‖ = 3.8723e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 63, Δt 423.9 ms: f = -2.731890276955e-01, ‖∇f‖ = 6.0233e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 64, Δt 463.8 ms: f = -2.731906762173e-01, ‖∇f‖ = 7.5548e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 65, Δt 422.5 ms: f = -2.731931250629e-01, ‖∇f‖ = 2.3343e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 66, Δt 448.3 ms: f = -2.731937854638e-01, ‖∇f‖ = 2.2001e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 67, Δt 422.6 ms: f = -2.731951572123e-01, ‖∇f‖ = 3.1454e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 68, Δt 451.6 ms: f = -2.731982612016e-01, ‖∇f‖ = 4.5688e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 69, Δt 1.34 s: f = -2.732004425120e-01, ‖∇f‖ = 7.4161e-03, α = 4.97e-01, m = 20, nfg = 2
[ Info: LBFGS: iter 70, Δt 415.0 ms: f = -2.732034978258e-01, ‖∇f‖ = 4.3965e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 71, Δt 400.2 ms: f = -2.732060317353e-01, ‖∇f‖ = 2.5959e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 72, Δt 420.2 ms: f = -2.732067653412e-01, ‖∇f‖ = 4.0226e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 73, Δt 418.0 ms: f = -2.732075286193e-01, ‖∇f‖ = 2.6208e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 74, Δt 449.4 ms: f = -2.732086394826e-01, ‖∇f‖ = 2.8144e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 75, Δt 428.7 ms: f = -2.732098613586e-01, ‖∇f‖ = 2.9348e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 76, Δt 441.6 ms: f = -2.732104250744e-01, ‖∇f‖ = 6.1063e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 77, Δt 427.8 ms: f = -2.732117991912e-01, ‖∇f‖ = 2.1887e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 78, Δt 458.8 ms: f = -2.732127290127e-01, ‖∇f‖ = 2.1978e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 79, Δt 432.6 ms: f = -2.732142378452e-01, ‖∇f‖ = 4.1889e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 80, Δt 452.3 ms: f = -2.732161076841e-01, ‖∇f‖ = 4.9611e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 81, Δt 1.04 s: f = -2.732168467935e-01, ‖∇f‖ = 1.1651e-02, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 82, Δt 342.3 ms: f = -2.732207559069e-01, ‖∇f‖ = 4.3107e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 83, Δt 421.8 ms: f = -2.732224453896e-01, ‖∇f‖ = 1.5951e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 84, Δt 413.5 ms: f = -2.732236351981e-01, ‖∇f‖ = 2.3806e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 85, Δt 385.1 ms: f = -2.732245493598e-01, ‖∇f‖ = 2.1850e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 86, Δt 843.3 ms: f = -2.732250468594e-01, ‖∇f‖ = 3.5908e-03, α = 3.99e-01, m = 20, nfg = 2
[ Info: LBFGS: iter 87, Δt 439.7 ms: f = -2.732258864564e-01, ‖∇f‖ = 1.7634e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 88, Δt 443.9 ms: f = -2.732264191346e-01, ‖∇f‖ = 1.5207e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 89, Δt 435.0 ms: f = -2.732270911530e-01, ‖∇f‖ = 2.0339e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 90, Δt 456.3 ms: f = -2.732275332266e-01, ‖∇f‖ = 3.4911e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 91, Δt 422.8 ms: f = -2.732281198242e-01, ‖∇f‖ = 1.8175e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 92, Δt 452.7 ms: f = -2.732286156101e-01, ‖∇f‖ = 1.7374e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 93, Δt 984.3 ms: f = -2.732290620520e-01, ‖∇f‖ = 2.4746e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 94, Δt 389.9 ms: f = -2.732298053289e-01, ‖∇f‖ = 2.7607e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 95, Δt 846.5 ms: f = -2.732303397137e-01, ‖∇f‖ = 4.7337e-03, α = 3.65e-01, m = 20, nfg = 2
[ Info: LBFGS: iter 96, Δt 402.4 ms: f = -2.732315115338e-01, ‖∇f‖ = 2.6853e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 97, Δt 434.7 ms: f = -2.732327488682e-01, ‖∇f‖ = 2.1740e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 98, Δt 447.3 ms: f = -2.732338965315e-01, ‖∇f‖ = 2.8814e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 99, Δt 422.9 ms: f = -2.732348853203e-01, ‖∇f‖ = 2.6726e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 100, Δt 455.0 ms: f = -2.732360755628e-01, ‖∇f‖ = 3.3694e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 101, Δt 462.7 ms: f = -2.732374157251e-01, ‖∇f‖ = 2.2761e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 102, Δt 425.7 ms: f = -2.732384441348e-01, ‖∇f‖ = 1.8476e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 103, Δt 435.0 ms: f = -2.732391060163e-01, ‖∇f‖ = 2.1827e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 104, Δt 880.6 ms: f = -2.732393121228e-01, ‖∇f‖ = 1.1870e-03, α = 4.35e-01, m = 20, nfg = 2
[ Info: LBFGS: iter 105, Δt 988.0 ms: f = -2.732394415105e-01, ‖∇f‖ = 9.6213e-04, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 106, Δt 352.8 ms: f = -2.732396967055e-01, ‖∇f‖ = 1.1000e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 107, Δt 413.8 ms: f = -2.732401458599e-01, ‖∇f‖ = 2.8831e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 108, Δt 386.8 ms: f = -2.732406523075e-01, ‖∇f‖ = 2.0902e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 109, Δt 383.9 ms: f = -2.732410323389e-01, ‖∇f‖ = 1.1560e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 110, Δt 440.7 ms: f = -2.732412518839e-01, ‖∇f‖ = 1.0759e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 111, Δt 424.2 ms: f = -2.732414186268e-01, ‖∇f‖ = 1.1730e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 112, Δt 448.3 ms: f = -2.732418417966e-01, ‖∇f‖ = 2.1278e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 113, Δt 422.0 ms: f = -2.732427145624e-01, ‖∇f‖ = 3.0607e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 114, Δt 451.6 ms: f = -2.732436618205e-01, ‖∇f‖ = 2.4125e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 115, Δt 419.3 ms: f = -2.732443955919e-01, ‖∇f‖ = 3.9514e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 116, Δt 453.2 ms: f = -2.732454002143e-01, ‖∇f‖ = 1.8433e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 117, Δt 411.4 ms: f = -2.732458974944e-01, ‖∇f‖ = 1.4613e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 118, Δt 419.5 ms: f = -2.732463785883e-01, ‖∇f‖ = 1.7104e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 119, Δt 1.40 s: f = -2.732466044311e-01, ‖∇f‖ = 2.1166e-03, α = 5.19e-01, m = 20, nfg = 2
[ Info: LBFGS: iter 120, Δt 356.3 ms: f = -2.732468613336e-01, ‖∇f‖ = 1.2965e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 121, Δt 420.7 ms: f = -2.732472718617e-01, ‖∇f‖ = 1.6179e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 122, Δt 453.3 ms: f = -2.732475792491e-01, ‖∇f‖ = 2.3507e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 123, Δt 395.6 ms: f = -2.732482303394e-01, ‖∇f‖ = 2.6828e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 124, Δt 903.9 ms: f = -2.732485536506e-01, ‖∇f‖ = 2.7793e-03, α = 4.74e-01, m = 20, nfg = 2
[ Info: LBFGS: iter 125, Δt 419.1 ms: f = -2.732489011759e-01, ‖∇f‖ = 1.3803e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 126, Δt 458.5 ms: f = -2.732492455012e-01, ‖∇f‖ = 1.5179e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 127, Δt 428.4 ms: f = -2.732494608995e-01, ‖∇f‖ = 2.1760e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 128, Δt 449.3 ms: f = -2.732500173039e-01, ‖∇f‖ = 2.8077e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 129, Δt 423.0 ms: f = -2.732505577698e-01, ‖∇f‖ = 2.7372e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 130, Δt 448.1 ms: f = -2.732510414354e-01, ‖∇f‖ = 1.4427e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 131, Δt 429.0 ms: f = -2.732513479891e-01, ‖∇f‖ = 1.2583e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 132, Δt 429.8 ms: f = -2.732516964375e-01, ‖∇f‖ = 2.2233e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 133, Δt 995.9 ms: f = -2.732520411325e-01, ‖∇f‖ = 2.4544e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 134, Δt 365.6 ms: f = -2.732522275980e-01, ‖∇f‖ = 4.5553e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 135, Δt 421.6 ms: f = -2.732529533185e-01, ‖∇f‖ = 1.3649e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 136, Δt 406.0 ms: f = -2.732531854735e-01, ‖∇f‖ = 1.0602e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 137, Δt 418.9 ms: f = -2.732534308208e-01, ‖∇f‖ = 1.6110e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 138, Δt 479.5 ms: f = -2.732537475106e-01, ‖∇f‖ = 1.9066e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 139, Δt 432.7 ms: f = -2.732541958553e-01, ‖∇f‖ = 2.7130e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 140, Δt 458.5 ms: f = -2.732543692518e-01, ‖∇f‖ = 3.4246e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 141, Δt 447.5 ms: f = -2.732549342383e-01, ‖∇f‖ = 1.0514e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 142, Δt 452.3 ms: f = -2.732550885399e-01, ‖∇f‖ = 9.8876e-04, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 143, Δt 442.1 ms: f = -2.732553016391e-01, ‖∇f‖ = 1.7312e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 144, Δt 497.9 ms: f = -2.732555304767e-01, ‖∇f‖ = 1.9487e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 145, Δt 463.5 ms: f = -2.732557706175e-01, ‖∇f‖ = 1.2634e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 146, Δt 618.1 ms: f = -2.732559583620e-01, ‖∇f‖ = 9.4509e-04, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 147, Δt 814.0 ms: f = -2.732560645489e-01, ‖∇f‖ = 1.1628e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 148, Δt 396.5 ms: f = -2.732563117685e-01, ‖∇f‖ = 1.9114e-03, α = 1.00e+00, m = 20, nfg = 1
[ Info: LBFGS: iter 149, Δt 421.7 ms: f = -2.732568091557e-01, ‖∇f‖ = 2.6653e-03, α = 1.00e+00, m = 20, nfg = 1
┌ Warning: LBFGS: not converged to requested tol after 150 iterations and time 1.96 m: f = -2.732575225334e-01, ‖∇f‖ = 2.5121e-03
└ @ OptimKit ~/.julia/packages/OptimKit/OEwMx/src/lbfgs.jl:199
E = -0.2732575225334225
We can compare our PEPS result to the energy obtained using a cylinder-MPS calculation
using a cylinder circumference of
E_ref = -0.273284888
@show (E - E_ref) / E_ref;(E - E_ref) / E_ref = -0.0001001353085337828
This page was generated using Literate.jl.