ITensor—Intelligent Tensor—is a C++ library for implementing tensor product wavefunction calculations. It is efficient and flexible enough to be used for research-grade simulations.
Features include:

ITensors have an Einstein summation interface making them nearly as easy to multiply as scalars: tensors indices have unique identities and matching indices automatically contract when two ITensors are multiplied. This type of interface makes it simple to transcribe tensor network diagrams into correct, efficient code.

For example, the diagram below (resembling the overlap of matrix product states) can be converted to code as

Installing ITensor:

  1. Make sure you have an up-to-date C++11 compiler and LAPACK installed. On UNIX systems, use your package manager; on Mac OS install the free Xcode app from the app store; for Windows install cygwin.
  2. Clone the latest version of ITensor:
    git clone itensor
    (Or download the zip file if you do not have git.)
    Cloning with git allows you to track changes to ITensor and is the preferred method; for more see our git quickstart guide.
  3. Create the file: cp Follow the instructions in this file to customize for your machine.
  4. Type make to build ITensor.

For more details, read the full installation instructions.
Browse the documentation to learn more about ITensor.

We are grateful for generous support from the Simons Foundation.

Code Samples

Perform a DMRG Calculation

//Define Hilbert space of N spin-one sites
int N = 100;
auto sites = SpinOne(N);

//Create 1d Heisenberg Hamiltonian
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;
auto H = MPO(ampo);

//Set up random initial wavefunction
auto psi = MPS(sites);

//Perform 5 sweeps of DMRG
auto sweeps = Sweeps(5);
//Specify max number of states kept each sweep
sweeps.maxm() = 50, 50, 100, 100, 200;

//Run the DMRG algorithm

//Continue to analyze wavefunction afterward 
Real energy = overlap(psi,H,psi); //<psi|H|psi>
for(int j = 1; j <= N; ++j)
    //Make site j the MPS "orthogonality center"
    //Measure magnetization
    Real Szj = (psi.A(j)
               * sites.op("Sz",j)
               * dag(prime(psi.A(j),Site))).real();
    println("Sz_",j," = ",Szj);

Contract Two ITensors

Index a("a",2), 
ITensor Z(a,b), 


//the * operator finds and //contracts common index 'b' //regardless of index order:
ITensor R = Z * X; Print( R.real(a(2),c(1)) ); //output: R.real(a(1),c(2)) = -1