ITensorMPOCompression.jl

Block respecting compression of MPOs and iMPOs
Author JanReimers
Popularity
0 Stars
Updated Last
6 Months Ago
Started In
February 2022

ITensorMPOCompression

A package to enable block respecting compression of Matrix Product Operators (MPOs) for use with ITensors.jl.

Reference: Daniel E. Parker, Xiangyu Cao, and Michael P. Zaletel Phys. Rev. B 102, 035147

ITensors has built in support for generating MPOs from local operators. This framework is called AutoMPO. MPOs generated by AutoMPO are already compressed, so you only need this package if you:

1. Are using MPOs that are created by hand, outside the `AutoMPO` framework.  
2. Want your MPO to be in orthogonal/canonical form.
3. Need to see the bond spectrums throughout the lattice. 

Support for compression of infinite lattice MPOs (iMPOs) is in the ITensorInfiniteMPS package which can be found here.

Technical notes on what this package does can be founs here

Installation

You can install this package through the Julia package manager:

julia> ] add ITensorMPOCompression

Examples

Here are is an example of orthogonalizing an hand generated (not using AutoMPO) MPO

julia> using ITensors

julia> using ITensorMPOCompression

julia> include("../test/hamiltonians/hamiltonians.jl");

julia> N = 10; #10 sites

julia> NNN = 7; #Include up to 7th nearest neighbour interactions

julia> sites = siteinds("S=1/2", N);

julia> H = transIsing_MPO(sites, NNN);

julia> is_regular_form(H,lower) == true
true

julia> @show get_Dw(H); #Show bond dimensions
get_Dw(H) = [30, 30, 30, 30, 30, 30, 30, 30, 30]

julia> pprint(H[2]) #schematic view of operators at site #2
I 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 
S 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 
S 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 
0 0 I 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 
S 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 
0 0 0 0 I 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 
0 0 0 0 0 I 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 
S 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 
0 0 0 0 0 0 0 I 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 
0 0 0 0 0 0 0 0 I 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 
0 0 0 0 0 0 0 0 0 I 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 
S 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 
0 0 0 0 0 0 0 0 0 0 0 I 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 
0 0 0 0 0 0 0 0 0 0 0 0 I 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 
0 0 0 0 0 0 0 0 0 0 0 0 0 I 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 
0 0 0 0 0 0 0 0 0 0 0 0 0 0 I 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 
S 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 
0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 I 0 0 0 0 0 0 0 0 0 0 0 0 0 
0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 I 0 0 0 0 0 0 0 0 0 0 0 0 
0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 I 0 0 0 0 0 0 0 0 0 0 0 
0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 I 0 0 0 0 0 0 0 0 0 0 
0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 I 0 0 0 0 0 0 0 0 0 
S 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 
0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 I 0 0 0 0 0 0 0 
0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 I 0 0 0 0 0 0 
0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 I 0 0 0 0 0 
0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 I 0 0 0 0 
0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 I 0 0 0 
0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 I 0 0 
0 S 0 S 0 0 S 0 0 0 S 0 0 0 0 S 0 0 0 0 0 S 0 0 0 0 0 0 S I 


julia> orthogonalize!(H,left); #Orthogonalize into left canonical form.  Also does rank reduction.

julia> orthogonalize!(H,right); #Orthogonalize into right canoncial form.

julia> pprint(H[2]) #Show vastly reduced matrix of operators at site #2
I 0 0 0 
S S S 0 
0 S 0 I 


julia> @show get_Dw(H); #Show reduced bond dimensions
get_Dw(H) = [3, 4, 5, 6, 7, 6, 5, 4, 3]

julia> is_regular_form(H,lower) == true
true

julia> isortho(H, right) == true #looks at cached ortho center limits
true

julia> check_ortho(H,right) == true #Does the more expensive V*V_dagger==Id contraction and test
true

julia> pprint(H) #High level view of what is in the MPO.
  n    Dw1  Dw2   d   Reg.  Orth.
                      Form  Form 
   1    1    3    2    L      M
   2    3    4    2    L      R
   3    4    5    2    L      R
   4    5    6    2    L      R
   5    6    7    2    L      R
   6    7    6    2    L      R
   7    6    5    2    L      R
   8    5    4    2    L      R
   9    4    3    2    L      R
  10    3    1    2    L      B

Here are is an example of truncating an hand generated (not using AutoMPO) MPO

julia> using ITensors

julia> using ITensorMPOCompression

julia> include("../test/hamiltonians/hamiltonians.jl");

julia> N = 10; #10 sites

julia> NNN = 7; #Include up to 7th nearest neighbour interactions

julia> sites = siteinds("S=1/2", N);

julia> H = transIsing_MPO(sites, NNN);

julia> is_regular_form(H,lower) == true
true

julia> spectrums = truncate!(H,left);

julia> pprint(H[5])
I 0 0 0 0 0 0 
S S S S S S 0 
S S S S S S 0 
S S S S S S 0 
S S S S S S 0 
0 S S S S S I 


julia> @show get_Dw(H);
get_Dw(H) = [3, 4, 5, 6, 7, 6, 5, 4, 3]

julia> @show spectrums;
spectrums = 
Bond  Ns   max(s)     min(s)    Entropy  Tr. Error
   1    1  0.30739   3.07e-01   0.22292  0.00e+00
   2    2  0.35392   3.49e-02   0.26838  0.00e+00
   3    3  0.37473   2.06e-02   0.29133  0.00e+00
   4    4  0.38473   1.77e-02   0.30255  0.00e+00
   5    5  0.38773   7.25e-04   0.30588  0.00e+00
   6    4  0.38473   1.77e-02   0.30255  0.00e+00
   7    3  0.37473   2.06e-02   0.29133  0.00e+00
   8    2  0.35392   3.49e-02   0.26838  0.00e+00
   9    1  0.30739   3.07e-01   0.22292  0.00e+00


julia> is_regular_form(H,lower) == true
true

julia> isortho(H, left) == true
true

julia> check_ortho(H, left) == true
true

Generating this README

This file was generated with weave.jl with the following commands:

using ITensorMPOCompression, Weave
weave(
  joinpath(pkgdir(ITensorMPOCompression), "examples", "README.jl");
  doctype="github",
  out_path=pkgdir(ITensorMPOCompression),
)

Used By Packages

No packages found.