Skip to content

Commit

Permalink
Ready for progress
Browse files Browse the repository at this point in the history
  • Loading branch information
leburgel committed Sep 5, 2023
1 parent f87847f commit 7fb9ead
Show file tree
Hide file tree
Showing 6 changed files with 49 additions and 34 deletions.
20 changes: 9 additions & 11 deletions examples/heisenberg.jl
Original file line number Diff line number Diff line change
Expand Up @@ -29,24 +29,23 @@ function iCtmGsEh(ψ::InfinitePEPS, env::PEPSKit.CTMRGEnv, H::AbstractTensorMap{
E = 0.0
for r in 1:size(ψ,1), c in 1:size(ψ,2)
ρ = two_site_rho(r, c, ψ, env)
nn = @tensor ρ[1,2;1,2]
Eh = @tensor H[1,2;3,4]*ρ[1,2;3,4]
@tensor nn = ρ[1 2; 1 2]
@tensor Eh = H[1 2; 3 4] * ρ[1 2; 3 4]
Eh = Eh / nn
E = E + Eh
E = E + Eh
#@diffset Es[r,c] = Eh;
end
return real(E)
end

function H_expectation_value::InfinitePEPS, env::PEPSKit.CTMRGEnv, H::AbstractTensorMap{S,2,2}) where S
Eh = iCtmGsEh(ψ, env, H)

ψ1 = rotl90(ψ)
env1 = PEPSKit.rotate_north(env,EAST);
env1 = PEPSKit.rotate_north(env, EAST)
Ev = iCtmGsEh(ψ1, env1, H)
E = real(Eh + Ev)
return E
end
end

function SqLatHeisenberg()
Sx,Sy,Sz,_ = spinmatrices(1//2)
Expand All @@ -72,13 +71,12 @@ function cfun(x)

function fun(peps)
env = leading_boundary(peps, alg_ctm, env)
x = H_expectation_value(peps, env, H)
x = H_expectation_value(peps, env, H)
return x
end

∂E = fun'(ψ)
env = leading_boundary(ψ, alg_ctm, env)
E = H_expectation_value(ψ, env, H)
∂E = fun'(ψ)

@assert !isnan(norm(∂E))
return E,∂E
Expand Down Expand Up @@ -124,7 +122,7 @@ end


alg_ctm = CTMRG(
verbose=10000,
verbose=1,
tol=1e-10,
trscheme=truncdim(10),
miniter=4,
Expand All @@ -137,7 +135,7 @@ function main(;d=2,D=2,Lx=1,Ly=1)
optimize(
cfun,
(ψ,env),
ConjugateGradient(verbosity=3);
ConjugateGradient(verbosity=2);
inner=my_inner,
retract=my_retract,
scale! = my_scale!,
Expand Down
12 changes: 6 additions & 6 deletions src/algorithms/ctmrg.jl
Original file line number Diff line number Diff line change
Expand Up @@ -29,7 +29,7 @@ function MPSKit.leading_boundary(peps_above::InfinitePEPS,peps_below::InfinitePE
while (err>alg.tol&&iter<=alg.maxiter) || iter<=alg.miniter
ϵ = 0.0
for i in 1:4
envs,ϵ₀ = left_move(peps_above, peps_below,alg,envs);
envs,ϵ₀ = left_move(peps_above, peps_below, alg, envs);
ϵ = max(ϵ,ϵ₀)
envs = rotate_north(envs,EAST);
peps_above = envs.peps_above;
Expand All @@ -40,7 +40,7 @@ function MPSKit.leading_boundary(peps_above::InfinitePEPS,peps_below::InfinitePE

err = abs(old_norm-new_norm)
= abs((ϵ₁-ϵ)/ϵ₁)
@ignore_derivatives alg.verbose > 1 && @printf("%4d %.2e %.10e %.2e %.2e\n",
@ignore_derivatives alg.verbose > 1 && @printf("CTMRG: \titeration: %4d\t\terror: %.2e\t\tnorm: %.10e\t\tϵ: %.2e\t\tdϵ: %.2e\n",
iter,err,abs(new_norm),ϵ,dϵ)

old_norm = new_norm
Expand All @@ -58,17 +58,17 @@ function left_move(peps_above::InfinitePEPS{PType},peps_below::InfinitePEPS{PTyp
corners::typeof(envs.corners) = copy(envs.corners);
edges::typeof(envs.edges) = copy(envs.edges);

above_projector_type = tensormaptype(spacetype(PType),1,3,storagetype(PType));
below_projector_type = tensormaptype(spacetype(PType),3,1,storagetype(PType));
above_projector_type = tensormaptype(spacetype(PType),1,3,storagetype(PType))
below_projector_type = tensormaptype(spacetype(PType),3,1,storagetype(PType))
ϵ = 0.0
n0 = 1.0
n1 = 1.0
for col in 1:size(peps_above,2)
cop = mod1(col+1,size(peps_above,2))
com = mod1(col-1,size(peps_above,2))

above_projs = Vector{above_projector_type}(undef,size(peps_above,1));
below_projs = Vector{below_projector_type}(undef,size(peps_above,1));
above_projs = Vector{above_projector_type}(undef,size(peps_above,1))
below_projs = Vector{below_projector_type}(undef,size(peps_above,1))

# find all projectors
for row in 1:size(peps_above,1)
Expand Down
26 changes: 13 additions & 13 deletions src/environments/ctmrgenv.jl
Original file line number Diff line number Diff line change
Expand Up @@ -9,36 +9,36 @@ end
CTMRGEnv(peps::InfinitePEPS) = CTMRGEnv(peps,peps);

function CTMRGEnv(peps_above::InfinitePEPS{P},peps_below::InfinitePEPS{P}) where P
ou = oneunit(space(peps_above,1,1)); # the bogus space
ou = oneunit(space(peps_above,1,1)) # the bogus space

C_type = tensormaptype(spacetype(P),1,1,storagetype(P));
T_type = tensormaptype(spacetype(P),3,1,storagetype(P)); # debatable how we should do the legs?
C_type = tensormaptype(spacetype(P),1,1,storagetype(P))
T_type = tensormaptype(spacetype(P),3,1,storagetype(P)) # debatable how we should do the legs?

#first index is de
corners = Array{C_type}(undef,4,size(peps_above)...);
edges = Array{T_type}(undef,4,size(peps_above)...);
corners = Array{C_type}(undef,4,size(peps_above)...)
edges = Array{T_type}(undef,4,size(peps_above)...)

for dir in 1:4, i in 1:size(peps_above,1),j in 1:size(peps_above,2)
@diffset corners[dir,i,j] = TensorMap(randn,eltype(P),ou,ou)
@diffset edges[dir,i,j] = TensorMap(randn,eltype(P),ou*space(peps_above[i,j],dir+1)'*space(peps_below[i,j],dir+1),ou)
end

@diffset corners[:,:,:]./=norm.(corners[:,:,:]);
@diffset edges[:,:,:]./=norm.(edges[:,:,:]);
@diffset corners[:,:,:]./=norm.(corners[:,:,:])
@diffset edges[:,:,:]./=norm.(edges[:,:,:])

CTMRGEnv(peps_above,peps_below,corners,edges)
end

function Base.rotl90(envs::CTMRGEnv{P,C,T}) where {P,C,T}
n_peps_above = rotl90(envs.peps_above);
n_peps_below = rotl90(envs.peps_below);
n_corners = Array{C,3}(undef,size(envs.corners)...);
n_edges = Array{T,3}(undef,size(envs.edges)...);
n_peps_above = rotl90(envs.peps_above)
n_peps_below = rotl90(envs.peps_below)
n_corners = Array{C,3}(undef,size(envs.corners)...)
n_edges = Array{T,3}(undef,size(envs.edges)...)

for dir in 1:4
dirm = mod1(dir-1,4)
@diffset n_corners[dirm,:,:] .= rotl90(envs.corners[dir,:,:]);
@diffset n_edges[dirm,:,:] .= rotl90(envs.edges[dir,:,:]);
@diffset n_corners[dirm,:,:] .= rotl90(envs.corners[dir,:,:])
@diffset n_edges[dirm,:,:] .= rotl90(envs.edges[dir,:,:])
end

CTMRGEnv(n_peps_above,n_peps_below,n_corners,n_edges)
Expand Down
19 changes: 18 additions & 1 deletion src/operators/transferpeps.jl
Original file line number Diff line number Diff line change
Expand Up @@ -61,4 +61,21 @@ end
function MPSKit.transfer_right(GR::GenericMPSTensor{S,3}, O::NTuple{2, PEPSTensor}, A::GenericMPSTensor{S,3}, Ā::GenericMPSTensor{S,3}) where S
return @tensor GR′[-1 -2 -3; -4] := GR[7 4 2; 1] * conj(Ā[-4 6 3; 1]) *
O[1][5; 9 4 6 -2] * conj(O[2][5; 8 2 3 -3]) * A[-1 9 8 7]
end
end

function MPSKit.expectation_value(st::MPSMultiline, O::TransferPEPSMultiline)

end

function MPSKit.expectation_value(st::MPSMultiline, ca::MPSKit.PerMPOInfEnv{H,V,S,A}) where {H<:TransferPEPSMultiline,V,S,A}
opp = ca.opp
retval = PeriodicArray{eltype(st.AC[1, 1]),2}(undef, size(st, 1), size(st, 2))
for (i, j) in product(1:size(st, 1), 1:size(st, 2))
retval[i, j] = @plansor leftenv(ca, i, j, st)[1 2; 3] *
opp[i, j][2 4; 6 5] *
st.AC[i, j][3 6; 7] *
rightenv(ca, i, j, st)[7 5; 8] *
conj(st.AC[i+1, j][1 4; 8])
end
return retval
end
4 changes: 2 additions & 2 deletions src/states/abstractpeps.jl
Original file line number Diff line number Diff line change
Expand Up @@ -40,5 +40,5 @@ Abstract supertype for a 2D projected entangled pairs operator.
abstract type AbstractPEPO end


Base.rotl90(t::PEPSTensor) = permute(t,(1,),(3,4,5,2));
Base.rotl90(t::PEPOTensor) = permute(t,(1,2),(4,5,6,3));
Base.rotl90(t::PEPSTensor) = permute(t, ((1,), (3, 4, 5, 2)))
Base.rotl90(t::PEPOTensor) = permute(t, ((1, 2), (4, 5, 6, 3)));
2 changes: 1 addition & 1 deletion src/utility/util.jl
Original file line number Diff line number Diff line change
Expand Up @@ -47,7 +47,7 @@ function _setindex(a::AbstractArray,v,args...)
b
end

function ChainRulesCore.rrule(::typeof(_setindex),a::AbstractArray,tv,args...)
function ChainRulesCore.rrule(::typeof(_setindex), a::AbstractArray, tv, args...)
t = _setindex(a,tv,args...);

function toret(v)
Expand Down

0 comments on commit 7fb9ead

Please sign in to comment.