Skip to content

Commit 899c5a5

Browse files
Add iPEPS SUWeight to CTMRGEnv conversion (#275)
Co-authored-by: leburgel <lander.burgelman@gmail.com>
1 parent 40d909b commit 899c5a5

3 files changed

Lines changed: 89 additions & 0 deletions

File tree

src/environments/suweight.jl

Lines changed: 46 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -370,3 +370,49 @@ function Base.rot180(wts::SUWeight)
370370
wts_y = circshift(rot180(wts[2, :, :]), (1, 0))
371371
return SUWeight(wts_x, wts_y)
372372
end
373+
374+
function _CTMRGEnv(wts::SUWeight, flips::Array{Bool, 3})
375+
@assert size(wts) == size(flips)
376+
_, Nr, Nc = size(wts)
377+
edges = map(Iterators.product(1:4, 1:Nr, 1:Nc)) do (d, r, c)
378+
wt_idx = if d == NORTH
379+
CartesianIndex(2, _next(r, Nr), c)
380+
elseif d == EAST
381+
CartesianIndex(1, r, _prev(c, Nc))
382+
elseif d == SOUTH
383+
CartesianIndex(2, r, c)
384+
else # WEST
385+
CartesianIndex(1, r, c)
386+
end
387+
wt = deepcopy(wts[wt_idx])
388+
if (flips[wt_idx] && d in (NORTH, EAST)) || (!flips[wt_idx] && d in (SOUTH, WEST))
389+
wt = permute(wt, ((2,), (1,)))
390+
end
391+
twistdual!(wt, 1)
392+
wt = insertrightunit(insertleftunit(wt, 1))
393+
return permute(wt, ((1, 2, 3), (4,)))
394+
end
395+
corners = map(CartesianIndices(edges)) do idx
396+
return TensorKit.id(Float64, codomain(edges[idx], 1))
397+
end
398+
return CTMRGEnv(corners, edges)
399+
end
400+
401+
"""
402+
CTMRGEnv(wts::SUWeight, peps::InfinitePEPS)
403+
404+
Construct a CTMRG environment with bond dimension χ = 1
405+
for an InfinitePEPS `peps` from the accompanying SUWeight `wts`.
406+
The scalartype of the returned environment is always `Float64`.
407+
"""
408+
function CTMRGEnv(wts::SUWeight, peps::InfinitePEPS)
409+
_, Nr, Nc = size(wts)
410+
flips = map(Iterators.product(1:2, 1:Nr, 1:Nc)) do (d, r, c)
411+
return if d == 1
412+
isdual(domain(peps[r, c], EAST))
413+
else
414+
isdual(domain(peps[r, c], NORTH))
415+
end
416+
end
417+
return _CTMRGEnv(wts, flips)
418+
end

test/ctmrg/suweight.jl

Lines changed: 40 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,40 @@
1+
using Test
2+
using Random
3+
using TensorKit
4+
using PEPSKit
5+
using PEPSKit: str
6+
7+
function rand_wts(peps::InfinitePEPS)
8+
wts = SUWeight(peps)
9+
for idx in CartesianIndices(wts.data)
10+
n = length(wts.data[idx].data)
11+
wts.data[idx].data[:] = sort(rand(Float64, n); lt = !isless)
12+
end
13+
return wts
14+
end
15+
16+
function su_rdm_1x1(row::Int, col::Int, peps::InfinitePEPS, wts::SUWeight)
17+
Nr, Nc = size(peps)
18+
@assert 1 <= row <= Nr && 1 <= col <= Nc
19+
t = peps.A[row, col]
20+
t = absorb_weight(t, wts, row, col, Tuple(1:4))
21+
ρ = t * t'
22+
return ρ / str(ρ)
23+
end
24+
25+
Vps = [Vect[U1Irrep](0 => 2, 1 => 2, -1 => 1), Vect[FermionParity](0 => 1, 1 => 2)]
26+
Vvs = [Vect[U1Irrep](0 => 3, 1 => 1, -1 => 2), Vect[FermionParity](0 => 2, 1 => 2)]
27+
28+
for (Vp, Vv) in zip(Vps, Vvs)
29+
Nspaces = [Vv Vv' Vv; Vv' Vv Vv']
30+
Espaces = [Vv Vv Vv'; Vv Vv' Vv']
31+
Pspaces = fill(Vp, size(Nspaces))
32+
peps = InfinitePEPS(randn, ComplexF64, Pspaces, Nspaces, Espaces)
33+
wts = rand_wts(peps)
34+
env = CTMRGEnv(wts, peps)
35+
for (r, c) in Tuple.(CartesianIndices(peps.A))
36+
ρ1 = su_rdm_1x1(r, c, peps, wts)
37+
ρ2 = reduced_densitymatrix(((r, c),), peps, env)
38+
@test ρ1 ρ2
39+
end
40+
end

test/runtests.jl

Lines changed: 3 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -20,6 +20,9 @@ end
2020
@time @safetestset ":fixed CTMRG iteration scheme" begin
2121
include("ctmrg/fixed_iterscheme.jl")
2222
end
23+
@time @safetestset "SUWeight conversion" begin
24+
include("ctmrg/suweight.jl")
25+
end
2326
@time @safetestset "Flavors" begin
2427
include("ctmrg/flavors.jl")
2528
end

0 commit comments

Comments
 (0)