3.2. Define Your Model#

3.2.1. Overview#

In this notebook we will introduce 3 basic components for Renormalizer: Op, BasisSet and Model. These components are essential for MPS/MPO construction/manipulation for a given physical model.

[1]:
# add renormalizer to path
import sys
sys.path.append("..")

3.2.2. Op#

Op is the abbreviation for “operators”. It offers a convenient interface to represent operators of your interest. Op is the bridge between symbolic math equations and numerical tensors in MPS/MPO.

[2]:
from renormalizer import Op
2022-10-28 13:29:10,474[INFO] Use NumPy as backend
2022-10-28 13:29:10,475[INFO] use 64 bits

3.2.2.1. Basics#

There are 3 basic attributes for an operator:

  • The operation or its symbol, such as \(\hat X\), \(\hat a^\dagger\), \(\hat p^2\)

  • The DOF at wich the operation is operated on, such as “spin 0”, “exiton 1”. In mathematical expression this is usually the subscript for an operator, i.e., \(\hat X_0\) or \(\hat a^\dagger_1\)

  • The factor or coefficient for the operator.

There is an additional useful attributes that is sometimes quite useful: the quantum number.

Op are constructed from the above 4 attributes. The first 2 are necessary and the last 2 are optional.

[3]:
Op("X", 0, factor=0.5, qn=0)
[3]:
Op('X', [0], 0.5, [0])

In renormalizer, you may use anything that can be hashed and compared to denote a DOF. Common examples include int, str, and tuple of them.

[4]:
Op("X", ("spin 0", "cool stuff", 2077))
[4]:
Op('X', [('spin 0', 'cool stuff', 2077)], 1.0)

You may wonder what are the allowed symbols for operators. For that please refer to the BasisSet section.

Products of operators can be constructed intuitively

[5]:
Op(r"a^\dagger a", [0, 1])
[5]:
Op('a^\\dagger a', [0, 1], 1.0)

Here, DOFs for the operators are specified through a list, where the first element is the DOF for \(a^\dagger\) and the second element is the DOF for \(a\).

Note that using tuple to specify the DOF has a totally different meaning. Renormalizer will recognize the tuple as a single DOF and set all DOFs in the operator to that DOF.

[6]:
Op(r"a^\dagger a", (0, 1))
[6]:
Op('a^\\dagger a', [(0, 1), (0, 1)], 1.0)

Op also has a lot of useful functions and attributes, please refer to the API document for more details.

3.2.2.2. Operations and OpSum#

Common and simple operations between operators are supported

[7]:
op1 = Op("X", 0, 0.5)
op2 = Op("Z", 1, 0.5)
op1 * op2
[7]:
Op('X Z', [0, 1], 0.25)

Addition between Op will result in an OpSum instance that is a subclass of list.

[8]:
op1 + op2, type(op1 + op2)
[8]:
([Op('X', [0], 0.5), Op('Z', [1], 0.5)], renormalizer.model.op.OpSum)

OpSum supports simple operator algebra such as multiplication/addition, etc.

[9]:
# multiplication
opsum = op1 + op2
opsum * op1
[9]:
[Op('X X', [0, 0], 0.25), Op('Z X', [1, 0], 0.25)]
[10]:
opsum * opsum
[10]:
[Op('X X', [0, 0], 0.25),
 Op('X Z', [0, 1], 0.25),
 Op('Z X', [1, 0], 0.25),
 Op('Z Z', [1, 1], 0.25)]
[11]:
# addition/subtraction
opsum -= op2
opsum
[11]:
[Op('X', [0], 0.5), Op('Z', [1], 0.5), Op('Z', [1], -0.5)]
[12]:
opsum.simplify()
[12]:
[Op('X', [0], 0.5)]

However, in general the functionalities are limited and the performance is not optimized. We recommand using SymPy for advanced symbolic mathematics and convert the final result to Renormalizer Op.

3.2.3. BasisSet#

An essential step for converting symbolic operators into numeric tensors is to specify a set of basis. Renormalizer includes a zoo of basis set with an emphasis for electron-phonon or vibronic models.

[13]:
import renormalizer
[s for s in dir(renormalizer) if s.startswith("Basis")]
[13]:
['BasisHalfSpin',
 'BasisHopsBoson',
 'BasisMultiElectron',
 'BasisMultiElectronVac',
 'BasisSHO',
 'BasisSimpleElectron',
 'BasisSineDVR']

Each BasisSet is associated with one or more DOFs. For example, we can setup the spin-1/2 basis set for a spin DOF denoted as ("spin", 0)

[14]:
from renormalizer import BasisHalfSpin
b = BasisHalfSpin(("spin", 0))
b
[14]:
(dof: ('spin', 0), nbas: 2, qn: [0, 0])

Each basis set supports a variety of operations. The spin-1/2 basis set naturally supports Pauli operators such as "X", "sigma_y".

[15]:
b.op_mat("X"), b.op_mat("sigma_y")
[15]:
(array([[0., 1.],
        [1., 0.]]),
 array([[0.+0.j, 0.-1.j],
        [0.+1.j, 0.+0.j]]))

For phonon DOF, it is necessary to specify the quanta truncated in the basis set

[16]:
# SHO represent simple harmonic oscillator
from renormalizer import BasisSHO
# omega is the vibration frequency, nbas is the number of basis
b = BasisSHO(dof=0, omega=50, nbas=4)
# note the different value because of omega
b.op_mat("b^\dagger+b"), b.op_mat("x")
[16]:
(array([[0.        , 1.        , 0.        , 0.        ],
        [1.        , 0.        , 1.41421356, 0.        ],
        [0.        , 1.41421356, 0.        , 1.73205081],
        [0.        , 0.        , 1.73205081, 0.        ]]),
 array([[0.        , 0.1       , 0.        , 0.        ],
        [0.1       , 0.        , 0.14142136, 0.        ],
        [0.        , 0.14142136, 0.        , 0.17320508],
        [0.        , 0.        , 0.17320508, 0.        ]]))

Note that BasisSHO supports both first and second-quantized operators such as

\[\hat x = \sqrt{\frac{1}{2\omega}} (b^\dagger + b)\]

here \(x\) is the mass-weighted coordinate.

The number of possible basis sets is infinite and highly customized basis set for a particular problem is rairly common. You may subclass the BasisSet parent class in renormalizer.model.basis.BasisSet to customize the basis set and the numerical behavior of the operators.

3.2.4. Model#

A Model is basically made up with a list of BasisSet and the Hamiltonian in OpSum/list form. A Model instance is necessary for any MPS/MPO construction.

[17]:
from renormalizer import Model
model = Model([BasisHalfSpin(0)], [Op("X", 0)])

For example, the corresponding MPO can be constructed from the above model conveniently

[18]:
from renormalizer import Mpo
mpo = Mpo(model)
mpo
2022-10-28 13:29:10,610[DEBUG] # of operator terms: 1
[18]:
<class 'renormalizer.mps.mpo.Mpo'> with 1 sites
[19]:
mpo.todense()
[19]:
array([[0., 1.],
       [1., 0.]])

You may wonder

  • What are the limitations for the model, in particular the operator terms, to construct MPO? For example, is nearest-neighbour interaction or at most two-body operator enforced?

  • If any numerical compression is involved in the construction. If so, how to control the accuracy?

  • What are the time cost of constructing MPO for an arbitrary model?

The answer is: renormalizer offers a unique algorithm for exact and efficient MPO construction of arbitrary operators. It means

  • There’s NO limitation for the model. Long-range interaction, three-body operators, anything you like.

  • Numerical compression is not involved and the construction is exact.

  • The construction time cost is neglegible for most models

Additionally, we also guarantee that the constructed MPO is optimal in terms of bond dimension.

For details of our algorithm, see A general automatic method for optimal construction of matrix product operators using bipartite graph theory.

Renormalizer also offers a set of built-in models for fast construction of common models. Please refer to the API document for more details.

[20]:
import renormalizer
[s for s in dir(renormalizer) if s.endswith("Model")]
[20]:
['HolsteinModel', 'Model', 'SpinBosonModel', 'TI1DModel']
[ ]: