Skip to content

Commit c2a32f2

Browse files
Add full update bond environment and gauge fixing (#203)
* Add full update bond environment and gauge fixing * There's still room to improve gauge fixing by avoiding `inv` of `R`, `L`
1 parent 32ef92d commit c2a32f2

File tree

6 files changed

+289
-0
lines changed

6 files changed

+289
-0
lines changed

src/PEPSKit.jl

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -45,7 +45,9 @@ include("algorithms/contractions/ctmrg_contractions.jl")
4545
include("algorithms/contractions/localoperator.jl")
4646
include("algorithms/contractions/vumps_contractions.jl")
4747
include("algorithms/contractions/bondenv/benv_tools.jl")
48+
include("algorithms/contractions/bondenv/gaugefix.jl")
4849
include("algorithms/contractions/bondenv/als_solve.jl")
50+
include("algorithms/contractions/bondenv/benv_ctm.jl")
4951

5052
include("algorithms/ctmrg/sparse_environments.jl")
5153
include("algorithms/ctmrg/ctmrg.jl")
Lines changed: 70 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,70 @@
1+
"""
2+
Construct the environment (norm) tensor
3+
```
4+
C1---T1---------T1---C2 r-1
5+
| ‖ ‖ |
6+
T4===XX== ==YY===T2 r
7+
| ‖ ‖ |
8+
C4---T3---------T3---C3 r+1
9+
c-1 c c+1 c+2
10+
```
11+
where `XX = X' X` and `YY = Y' Y` (stacked together).
12+
13+
Axis order: `[DX1 DY1; DX0 DY0]`, as in
14+
```
15+
┌---------------------┐
16+
| ┌----┐ |
17+
└-| |---DX0 DY0---┘
18+
|benv|
19+
┌-| |---DX1 DY1---┐
20+
| └----┘ |
21+
└---------------------┘
22+
```
23+
"""
24+
function bondenv_fu(row::Int, col::Int, X::PEPSOrth, Y::PEPSOrth, env::CTMRGEnv)
25+
Nr, Nc = size(env.corners)[[2, 3]]
26+
cm1 = _prev(col, Nc)
27+
cp1 = _next(col, Nc)
28+
cp2 = _next(cp1, Nc)
29+
rm1 = _prev(row, Nr)
30+
rp1 = _next(row, Nr)
31+
c1 = env.corners[1, rm1, cm1]
32+
c2 = env.corners[2, rm1, cp2]
33+
c3 = env.corners[3, rp1, cp2]
34+
c4 = env.corners[4, rp1, cm1]
35+
t1X, t1Y = env.edges[1, rm1, col], env.edges[1, rm1, cp1]
36+
t2 = env.edges[2, row, cp2]
37+
t3X, t3Y = env.edges[3, rp1, col], env.edges[3, rp1, cp1]
38+
t4 = env.edges[4, row, cm1]
39+
#= index labels
40+
41+
C1--χ4--T1X---------χ6---------T1Y--χ8---C2 r-1
42+
| ‖ ‖ |
43+
χ2 DNX DNY χ10
44+
| ‖ ‖ |
45+
T4==DWX==XX===DX== ==DY===YY==DEY==T2 r
46+
| ‖ ‖ |
47+
χ1 DSX DSY χ9
48+
| ‖ ‖ |
49+
C4--χ3--T3X---------χ5---------T3Y--χ7---C3 r+1
50+
c-1 c c+1 c+2
51+
=#
52+
@autoopt @tensor benv[DX1 DY1; DX0 DY0] := (
53+
c4[χ3 χ1] *
54+
t4[χ1 DWX0 DWX1 χ2] *
55+
c1[χ2 χ4] *
56+
t3X[χ5 DSX0 DSX1 χ3] *
57+
X[DNX0 DX0 DSX0 DWX0] *
58+
conj(X[DNX1 DX1 DSX1 DWX1]) *
59+
t1X[χ4 DNX0 DNX1 χ6] *
60+
c3[χ9 χ7] *
61+
t2[χ10 DEY0 DEY1 χ9] *
62+
c2[χ8 χ10] *
63+
t3Y[χ7 DSY0 DSY1 χ5] *
64+
Y[DNY0 DEY0 DSY0 DY0] *
65+
conj(Y[DNY1 DEY1 DSY1 DY1]) *
66+
t1Y[χ6 DNY0 DNY1 χ8]
67+
)
68+
normalize!(benv, Inf)
69+
return benv
70+
end
Lines changed: 120 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,120 @@
1+
"""
2+
Replace bond environment `benv` by its positive approximant `Z† Z`
3+
(returns the "half environment" `Z`)
4+
```
5+
┌-----------------┐ ┌---------------┐
6+
| ┌----┐ | | |
7+
└-| |-- 3 4 --┘ └-- Z -- 3 4 --┘
8+
|benv| = ↓
9+
┌-| |-- 1 2 --┐ ┌-- Z†-- 1 2 --┐
10+
| └----┘ | | |
11+
└-----------------┘ └---------------┘
12+
```
13+
"""
14+
function positive_approx(benv::BondEnv)
15+
# eigen-decomposition: benv = U D U'
16+
D, U = eigh((benv + benv') / 2)
17+
# determine if `env` is (mostly) positive or negative
18+
sgn = sign(sum(D.data))
19+
# When optimizing the truncation of a bond,
20+
# its environment can always be multiplied by a number.
21+
# If `benv` is negative (e.g. obtained approximately from CTMRG),
22+
# we can multiply it by (-1).
23+
data = D.data
24+
@inbounds for i in eachindex(data)
25+
d = (sgn == -1) ? -data[i] : data[i]
26+
data[i] = (d > 0) ? sqrt(d) : zero(d)
27+
end
28+
Z = D * U'
29+
return Z
30+
end
31+
32+
"""
33+
Use QR decomposition to fix gauge of the half bond environment `Z`.
34+
The reduced bond tensors `a`, `b` and `Z` are arranged as
35+
```
36+
┌---------------┐
37+
| |
38+
└---Z---a---b---┘
39+
| ↓ ↓
40+
41+
```
42+
Reference:
43+
- Physical Review B 90, 064425 (2014)
44+
- Physical Review B 92, 035142 (2015)
45+
"""
46+
function fixgauge_benv(
47+
Z::AbstractTensorMap{T,S,1,2},
48+
a::AbstractTensorMap{T,S,1,2},
49+
b::AbstractTensorMap{T,S,2,1},
50+
) where {T<:Number,S<:ElementarySpace}
51+
@assert !isdual(space(Z, 1))
52+
@assert !isdual(space(a, 2))
53+
@assert !isdual(space(b, 2))
54+
#= QR/LQ decomposition of Z
55+
56+
3 - Z - 2 = 2 - L - 1 3 - QL - 1
57+
↓ ↓
58+
1 2
59+
60+
= 1 - QR - 3 1 - R - 2
61+
62+
2
63+
=#
64+
QL, L = leftorth(Z, ((2, 1), (3,)))
65+
QR, R = leftorth(Z, ((3, 1), (2,)))
66+
@debug "cond(L) = $(LinearAlgebra.cond(L)); cond(R) = $(LinearAlgebra.cond(R))"
67+
# TODO: find a better way to fix gauge that avoids `inv`
68+
Linv, Rinv = inv(L), inv(R)
69+
#= fix gauge of Z, a, b
70+
┌---------------------------------------┐
71+
| |
72+
└---Z---Rinv)---(R--a)--(b--L)---(Linv--┘
73+
| ↓ ↓
74+
75+
76+
-1 - R - 1 - a - -3 -1 - b - 1 - L - -3
77+
↓ ↓
78+
-2 -2
79+
80+
┌-----------------------------------------┐
81+
| |
82+
└---Z-- 1 --Rinv-- -2 -3 --Linv-- 2 -┘
83+
84+
-1
85+
=#
86+
@plansor a[-1; -2 -3] := R[-1; 1] * a[1; -2 -3]
87+
@plansor b[-1 -2; -3] := b[-1 -2; 1] * L[-3; 1]
88+
@plansor Z[-1; -2 -3] := Z[-1; 1 2] * Rinv[1; -2] * Linv[2; -3]
89+
(isdual(space(R, 1)) == isdual(space(R, 2))) && twist!(a, 1)
90+
(isdual(space(L, 1)) == isdual(space(L, 2))) && twist!(b, 3)
91+
return Z, a, b, (Linv, Rinv)
92+
end
93+
94+
"""
95+
When the (half) bond environment `Z` consists of two `PEPSOrth` tensors `X`, `Y` as
96+
```
97+
┌---------------┐ ┌-------------------┐
98+
| | = | | ,
99+
└---Z-- --┘ └--Z0---X-- --Y--┘
100+
↓ ↓
101+
```
102+
apply the gauge transformation `Linv`, `Rinv` for `Z` to `X`, `Y`:
103+
```
104+
-1 -1
105+
| |
106+
-4 - X - 1 - Rinv - -2 -4 - Linv - 1 - Y - -2
107+
| |
108+
-3 -3
109+
```
110+
"""
111+
function _fixgauge_benvXY(
112+
X::PEPSOrth{T,S},
113+
Y::PEPSOrth{T,S},
114+
Linv::AbstractTensorMap{T,S,1,1},
115+
Rinv::AbstractTensorMap{T,S,1,1},
116+
) where {T<:Number,S<:ElementarySpace}
117+
@plansor X[-1 -2 -3 -4] := X[-1 1 -3 -4] * Rinv[1; -2]
118+
@plansor Y[-1 -2 -3 -4] := Y[-1 -2 -3 1] * Linv[1; -4]
119+
return X, Y
120+
end

test/bondenv/benv_fu.jl

Lines changed: 48 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,48 @@
1+
using Test
2+
using TensorKit
3+
using PEPSKit
4+
using LinearAlgebra
5+
using KrylovKit
6+
using Random
7+
8+
Random.seed!(100)
9+
Nr, Nc = 2, 2
10+
# create Hubbard iPEPS using simple update
11+
function get_hubbard_state(t::Float64=1.0, U::Float64=8.0)
12+
H = hubbard_model(ComplexF64, Trivial, U1Irrep, InfiniteSquare(Nr, Nc); t, U, mu=U / 2)
13+
Vphy = Vect[FermionParity U1Irrep]((0, 0) => 2, (1, 1//2) => 1, (1, -1//2) => 1)
14+
wpeps = InfiniteWeightPEPS(rand, ComplexF64, Vphy, Vphy; unitcell=(Nr, Nc))
15+
alg = SimpleUpdate(1e-2, 1e-8, 10000, truncerr(1e-10) & truncdim(4))
16+
wpeps, = simpleupdate(wpeps, H, alg; bipartite=false, check_interval=2000)
17+
peps = InfinitePEPS(wpeps)
18+
normalize!.(peps.A, Inf)
19+
return peps
20+
end
21+
22+
peps = get_hubbard_state()
23+
# calculate CTMRG environment
24+
Envspace = Vect[FermionParity U1Irrep](
25+
(0, 0) => 4, (1, 1//2) => 1, (1, -1//2) => 1, (0, 1) => 1, (0, -1) => 1
26+
)
27+
ctm_alg = SequentialCTMRG(; tol=1e-10, verbosity=2, trscheme=truncerr(1e-10) & truncdim(8))
28+
env, = leading_boundary(CTMRGEnv(rand, ComplexF64, peps, Envspace), peps, ctm_alg)
29+
for row in 1:Nr, col in 1:Nc
30+
cp1 = PEPSKit._next(col, Nc)
31+
A, B = peps.A[row, col], peps.A[row, cp1]
32+
X, a, b, Y = PEPSKit._qr_bond(A, B)
33+
benv = PEPSKit.bondenv_fu(row, col, X, Y, env)
34+
@assert [isdual(space(benv, ax)) for ax in 1:numind(benv)] == [0, 0, 1, 1]
35+
Z = PEPSKit.positive_approx(benv)
36+
# verify that gauge fixing can greatly reduce
37+
# condition number for physical state bond envs
38+
cond1 = cond(Z' * Z)
39+
Z2, a2, b2, (Linv, Rinv) = PEPSKit.fixgauge_benv(Z, a, b)
40+
benv2 = Z2' * Z2
41+
cond2 = cond(benv2)
42+
@test 1 <= cond2 < cond1
43+
@info "benv cond number: (gauge-fixed) $(cond2)$(cond1) (initial)"
44+
# verify gauge fixing is done correctly
45+
@tensor half[:] := Z[-1; 1 3] * a[1; -2 2] * b[2 -3; 3]
46+
@tensor half2[:] := Z2[-1; 1 3] * a2[1; -2 2] * b2[2 -3; 3]
47+
@test half half2
48+
end

test/bondenv/benv_gaugefix.jl

Lines changed: 43 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,43 @@
1+
using Test
2+
using TensorKit
3+
using PEPSKit
4+
using LinearAlgebra
5+
using KrylovKit
6+
using Random
7+
8+
Vphy = Vect[FermionParity U1Irrep]((0, 0) => 1, (1, 1) => 1, (1, -1) => 2)
9+
Vin = Vect[FermionParity U1Irrep]((0, 0) => 1, (1, 1) => 3, (1, -1) => 2)
10+
V = Vect[FermionParity U1Irrep]((0, 0) => 1, (1, 1) => 2, (1, -1) => 3)
11+
Vs = (V, V')
12+
for V1 in Vs, V2 in Vs, V3 in Vs
13+
#=
14+
┌---┬---------------┬---┐
15+
| | | | ┌--------------┐
16+
├---X--- -2 -3 ---Y---┤ = | |
17+
| | | | └--Z-- -2 -3 -┘
18+
└---┴-------Z0------┴---┘ ↓
19+
↓ -1
20+
-1
21+
=#
22+
X = rand(ComplexF64, Vin V1' Vin' Vin)
23+
Y = rand(ComplexF64, Vin Vin Vin' V3)
24+
Z0 = randn(ComplexF64, Vphy Vin Vin' Vin Vin Vin Vin')
25+
@tensor Z[p; Xe Yw] := Z0[p; Xn Xs Xw Yn Ye Ys] * X[Xn Xe Xs Xw] * Y[Yn Ye Ys Yw]
26+
#=
27+
┌---------------------------┐
28+
| |
29+
└---Z-- 1 --a-- 2 --b-- 3 --┘
30+
↓ ↓ ↓
31+
-1 -2 -3
32+
=#
33+
a = randn(ComplexF64, V1 Vphy' V2)
34+
b = randn(ComplexF64, V2 Vphy V3)
35+
@tensor half[:] := Z[-1; 1 3] * a[1; -2 2] * b[2 -3; 3]
36+
Z2, a2, b2, (Linv, Rinv) = PEPSKit.fixgauge_benv(Z, a, b)
37+
@tensor half2[:] := Z2[-1; 1 3] * a2[1; -2 2] * b2[2 -3; 3]
38+
@test half half2
39+
# test gauge transformation of X, Y
40+
X2, Y2 = PEPSKit._fixgauge_benvXY(X, Y, Linv, Rinv)
41+
@tensor Z2_[p; Xe Yw] := Z0[p; Xn Xs Xw Yn Ye Ys] * X2[Xn Xe Xs Xw] * Y2[Yn Ye Ys Yw]
42+
@test Z2 Z2_
43+
end

test/runtests.jl

Lines changed: 6 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -50,6 +50,12 @@ end
5050
@time @safetestset "Iterative optimization after truncation" begin
5151
include("bondenv/bond_truncate.jl")
5252
end
53+
@time @safetestset "Gauge fixing" begin
54+
include("bondenv/benv_gaugefix.jl")
55+
end
56+
@time @safetestset "Full update bond environment" begin
57+
include("bondenv/benv_fu.jl")
58+
end
5359
end
5460
if GROUP == "ALL" || GROUP == "TIMEEVOL"
5561
@time @safetestset "Cluster truncation with projectors" begin

0 commit comments

Comments
 (0)