Learn to Use ITensor

main / classes / autompo C++v3 | C++v2

AutoMPO

AutoMPO is a very powerful system for translating sums of local operators into an MPO (or IQMPO) tensor network. The notation for AutoMPO input is designed to be as close as possible to pencil-and-paper quantum mechanics notation.

The code for AutoMPO is located in the files "itensor/mps/autompo.h" and "itensor/mps/autompo.cc".

Synopsis

// Make a chain of N spin 1/2's with QN information
auto Nx = 12, Ny = 6;
auto N = Nx*Ny;
auto sites = SpinHalf(N,{"ConserveQNs=",true});

//
// Use AutoMPO to create the 
// next-neighbor Heisenberg model
//

auto ampo = AutoMPO(sites);
for(int j = 1; j < N; ++j)
    {
    ampo += 0.5,"S+",j,"S-",j+1;
    ampo += 0.5,"S-",j,"S+",j+1;
    ampo +=     "Sz",j,"Sz",j+1;
    }

  //Convert the AutoMPO object to an MPO
auto H = toMPO(ampo);

  //....

//
// Create a model with further-range interactions
  // capturing a 2D lattice (with a 1D ordering of sites)
//
auto ampo2D = AutoMPO(sites);
auto lattice = squareLattice(Nx,Ny);
for(auto b : lattice)
    {
    ampo2D += 0.5,"S+",b.s1,"S-",b.s2;
    ampo2D += 0.5,"S-",b.s1,"S+",b.s2;
    ampo2D +=     "Sz",b.s1,"Sz",b.s2;
    }
auto H2D = toMPO(ampo2D);

  //....

//
// Create the 1D Hubbard model
//
auto sitesElec = Electron(N);
auto ampoHub = AutoMPO(sitesElec);
for(int i = 1; i <= N; ++i)
    {
    ampoHub += U,"Nupdn",i;
    }
for(int b = 1; b < N; ++b)
    {
    ampoHub += -t,"Cdagup",b,"Cup",b+1;
    ampoHub += -t,"Cdagup",b+1,"Cup",b;
    ampoHub += -t,"Cdagdn",b,"Cdn",b+1;
    ampoHub += -t,"Cdagdn",b+1,"Cdn",b;
    }
auto H = toMPO(ampoHub);

//....

//
// Create a spin model with four-site terms
//
auto ampo4 = AutoMPO(sites);
for(auto i : range1(N-4))
    {
    ampo4 += "Sz",i,"Sz",i+1,"Sz",i+2,"Sz",i+4;
    }
for(auto i : range1(N-1))
    {
    ampo4 += 0.5,"S+",i,"S-",i+1;
    ampo4 += 0.5,"S-",i,"S+",i+1;
    }
auto H4 = toMPO(ampo4);

AutoMPO Interface

Converting an AutoMPO to an MPO

toMPO function

Call the toMPO function to create an MPO from an AutoMPO. You can pass various named arguments to control which backend is used to process the AutoMPO and to control the accuracy of this process.

Examples:

auto H1 = toMPO(ampo);
auto H2 = toMPO(ampo,{"Exact=",true});

Named arguments recognized:

toExpH function

The toExpH function converts an AutoMPO representing a sum of operators @@H@@ into an MPO which approximates the operator @@e^{-t H}@@ for a small time step t, making an error of order @@t^2@@ per time step. The time step t can be real or complex.

The method used to do the approximate exponentiation is based on the following article: Phys. Rev. B 91, 165112 (arxiv:1407.1832), and has the advantage that unlike naive approaches for exponentiating an MPO, the time-step error per site is independent of system size.

Note that the true amount of error per step, or the quality of the results can depend very highly on how short- or long-range the interactions are in the input Hamiltonian. It is recommended to test your results against a different method, such as Trotter gates, global Krylov, or TDVP to confirm that your results are accurate and controlled.

Examples:

auto ampo = AutoMPO(sites);
...
Real tau = 0.1;
//Real time evolution
auto expiH = toExpH(ampo,tau*Cplx_i);
//...
//Imaginary time evolution
auto expH = toExpH(ampo,tau);


This page current as of version 3.0.0


Back to Classes
Back to Main