Getting Started
Since quantum mechanical systems are described by operators acting on state vectors living in Hilbert spaces, this library is mainly concerned with providing such operators and functions that manipulate state vectors. In this chapter, you will learn the essentials of using the HQS Quantum Solver software by considering an example computation.
Let us consider the tight-binding Hamiltonian defined on a chain lattice, i.e., the Hamiltonian is given by where is the number of lattice sites, the strength of the hopping term and the annihilation operator at site .
Imports and Parameters
First, we need to import the functions, classes and modules that are required for the script. If you want to follow along, copy the following imports into your script file. We will explain the classes and functions as we go along.
# Title : HQS Quantum Solver Getting Started Guide
# Filename : getting_started.py
# Python, NumPy and SciPy imports
from math import pi, sqrt
import numpy as np
from scipy.linalg import toeplitz
from scipy.sparse.linalg import eigsh
import matplotlib.pyplot as plt
# HQS Quantum Solver imports
from hqs_quantum_solver.operators import (
SpinlessFermionOperatorEigen as SpinlessFermionOperator
)
from hqs_quantum_solver.types import (
SpinlessFermionHamiltonianInput,
SpinlessFermionAnnihilationInput,
)
from hqs_quantum_solver.interface import scipy_operator_adapter
Next, we define the parameters that we will refer to in the script.
M = 16 # site count
N = 1 # particle count
E = 0 # energy level
With that out of the way, we can start setting up the tight-binding Hamiltonian.
Setting up the Hamiltonian
We want to look at a system given by the tight-binding Hamiltonian. Thus, the first step is to construct said Hamiltonian. Quantum Solver provides a generic interface for constructing Hamiltonians, as well as several convenience functions. Especially, Quantum Solver implements Hamiltonians of the form where is the, so called, hopping matrix. Note, that the additional terms default to zero unless we explicitly enable them.
If we define the matrix to be then the Hamiltonian in the previous formula becomes the tight-binding Hamiltonian defined at the beginning of this chapter.
This matrix can be easily constructed using the
toeplitz
function from SciPy.
t = 1
c = np.array([-t if j == 1 else 0 for j in range(M)])
h0 = toeplitz(c)
With this matrix at hand, we can construct the Hamiltonian.
Quantum Solver supports multiple linear algebra backends. For this reason,
the construction of quantum mechanical operators is a two-step process.
First, we create an "Input" object that describes the operator to be
constructed. Then, we pass this object to an operator class of the desired
backend. Note that we have imported the
SpinlessFermionOperatorEigen
as SpinlessFermionOperator
, which means that
we are using the Eigen backend, which implements its functionality using
the well-known Eigen library.
To create the "Input", we create a
SpinlessFermionHamiltonianInput
object. The constructor of this object expects
four requried arguments and a bunch of additional arguments. The additional
arguments define the additional terms mentioned above. The required arguments
are M (int)
, N (int)
, mod_N (int)
, and h0 (ndarray)
.
The parameter M
describes the number of sites in the system and
the parameter N
the number of particles in the system.
The parameter mod_N
sets the step size in which the number of particles
is allowed to change during the simulation. For example for
M = 10
, N = 5
and mod_N = 4
the system allows for states with
1, 5, and 9 particles. When mod_N = 0
the number of particles is conserved,
i.e., constant.
The parameter h0
is the "hopping matrix" given as a NumPy array.
With all this information, we can create the Hamiltonian by first creating the "Input" object and then handing it over to the linear algebra backend.
H_input = SpinlessFermionHamiltonianInput(M=M, N=N, mod_N=0, h0=h0)
H = SpinlessFermionOperator(H_input)
The Quantum Solver operators
primarily provide matrix-vector and matrix-matrix
routines.
Here, we will be using the
SpinlessFermionOperatorEigen
.
As a first example for using the
operator, we want to look at its spectrum. For that, we can use
the eigsh
function from SciPy. To use this function, however, we need to use an
adapter
(see scipy_operator_adapter
)
to make the Quantum Solver operator compatible with SciPy.
k = min(12, H.shape[0] - 1) # Number of eigenvalue to compute
evals, evecs = eigsh(scipy_operator_adapter(H), k=k, which="SA")
print(evals)
plt.figure("eigenvalues", clear=True)
plt.title(f"Eigenvalues (M = {M}, N = {N})")
plt.plot(evals, ".")
Now, we have enough code to run the script for the first time. Doing so should produce the plot shown below, which shows the ten smallest eigenvalues of the Hamiltonian, and hence, the lowest energy levels of the system.
Before continuing, modify the values of and , run the script, and look at the result a few times.
Note, for mod_N=0
the number of states Quantum Solver has to store is
. Hence, the number of states can grow rather quickly when
increasing N
while M
is large. This problem gets worse for mod_N != 0
.
Keep in mind, while
Quantum Solver provides you with highly accurate computations, you will
have to provide enough computational resources for it to shine.
Below, we have listed the 5 smallest eigenvalues for M = 16
and different
values of N
.
# | |||
---|---|---|---|
1 | -1.97 | -3.83 | -5.53 |
2 | -1.86 | -3.67 | -5.31 |
3 | -1.70 | -3.57 | -5.14 |
4 | -1.48 | -3.44 | -5.04 |
5 | -1.21 | -3.34 | -5.04 |
We can make our first observation: Adding the first two values of the column gives you the first value in the column, and adding the first three values of the column will give you the first value in the column. The reason for this is that in a two- and three-particle system, the lowest energy to obtain is the one of the particles occupying the lowest two and three energy levels, respectively. Since we are dealing with a fermionic system, the particles are not allowed to occupy the same state, which explains this behavior.
We shall later come back to this observation.
Note how few lines of code were necessary, to produce physically meaningful results.
Site Occupations
To take this example a little bit further, let us look at the computation of site occupation expectation values. More precisely, we want to compute for , where . Note that Hence, we just need to compute and then compute the square norm of this vector.
First, we have to construct the annihilation operator .
As before, we first have to create an input object that we then hand over
to the linear algebra backend.
This time we create an SpinlessFermionAnnihilationInput
object.
The constructor expects the parameters
M (int)
, N (int)
, and mod_N (int)
, which
describe the domain, i.e., the input, of the operator and have the same meaning as
before. Additionally, we need to provide an argument called
c_i (int | ndarray)
, which describes which annihilation operator is to be
constructed. When passing an integer j
, the annihilation operator for site
j
is constructed.
(See SpinlessFermionAnnihilationInput
.)
Note that the codomain of the operator has N - 1
particles, so the
operator maps to a different space.
Once we have constructed the operator, we can use the
dot
method, to apply it to a state vector.
Below, we define a function site_occupation
which computes the expectation
value for a given state vector and site. We then plot the expectation
value for every site and for the state vector corresponding to the
energy level E
.
def site_occupation(psi, j):
annihilation_input = SpinlessFermionAnnihilationInput(M=M, N=N, mod_N=0, c_i=j)
annihilation = SpinlessFermionOperator(annihilation_input)
c_psi = annihilation.dot(psi)
return c_psi.conj() @ c_psi
occupations = np.array([site_occupation(evecs[:, E], j) for j in range(M)])
plt.figure("site_occupation", clear=True)
plt.title(f"Site Occupation (M = {M}, N = {N}, E = {E}")
plt.plot(occupations, "x-")
plt.axis((0, M - 1, 0, 1))
plt.xlabel("j")
plt.ylabel("Expectation Value")
Energy occupations
While looking a site occupations is certainly interesting, for most applications it is more useful to look at energy occupations instead.
Consider the matrix given by The columns of this matrix are the eigenvectors of the matrix . It is a well-known fact that if we define we can write the Hamiltonian as where are the energy levels of the system. Consequently, we can measure expectation value of the occupation of the -th energy level, by computing
We can carry out this computation in a similar way to the computation of
the site occupation. All we have to do is to set the annihilation
variable
to a linear combination of annihilation operators, which is given by the
definition of .
Luckily, this is such a common operation that the Quantum Solver includes a
convenience method for exactly this case.
The c_i
parameter of the constructor of the
SpinlessFermionAnnihilationInput
also accepts an array of coefficients and creates the corresponding
linear combination of annihilation operators.
J, K = np.meshgrid(np.arange(M), np.arange(M), indexing='ij')
Q = sqrt(2) / sqrt(M+1) * np.sin((J+1) * (K+1) * pi / (M+1))
def energy_occupation(psi, k):
annihilation_input = SpinlessFermionAnnihilationInput(
M=M, N=N, mod_N=0, c_i=Q[:,k].conj()
)
annihilation = SpinlessFermionOperator(annihilation_input)
c_psi = annihilation.dot(psi)
return (c_psi.conj() @ c_psi).real
ks = np.arange(M)
expectation = np.array([energy_occupation(evecs[:, E], k) for k in ks])
plt.figure("energy", clear=True)
plt.title(f"Energy Occupation (M = {M}, N = {N}, E = {E})")
plt.plot(ks, expectation, "x-")
plt.xlabel("p")
plt.ylabel("Expectation Value")
Below, you can find the plots that this script produces for different parameter configurations.
In these figures, you can see how the fermions arrange in the different one-particle energy states to achieve the energy levels we have observed in the table in the Setting up the Hamiltonian section.
Note that if you do not know the eigenvectors of the matrix , you can
easily obtain them by using the
eigh
function of NumPy.
Summary
To conclude this introduction, let us sum up the main take-aways of this chapter.
- Quantum Solver works by providing classes for working with quantum mechanical operators.
- Operators are created by constructing an "Input" object and handing it to a linear algebra backend.
- The parameters
M
,N
,mod_N
describe the domain of the operator, whereM
is the number of sites,N
the number of particles, andmod_N
the step size by which the number of particles is allowed to change. - To convert Quantum Solver operators into SciPy linear operators, you can
use the
scipy_operator_adapter
function. - Operators provide the
dot
method, which applies the operator to a state vector.
Complete Code
# Title : HQS Quantum Solver Getting Started Guide
# Filename : getting_started.py
# Python, NumPy and SciPy imports
from math import pi, sqrt
import numpy as np
from scipy.linalg import toeplitz
from scipy.sparse.linalg import eigsh
import matplotlib.pyplot as plt
# HQS Quantum Solver imports
from hqs_quantum_solver.operators import (
SpinlessFermionOperatorEigen as SpinlessFermionOperator
)
from hqs_quantum_solver.types import (
SpinlessFermionHamiltonianInput,
SpinlessFermionAnnihilationInput,
)
from hqs_quantum_solver.interface import scipy_operator_adapter
# ===== Parameters =====
M = 16 # site count
N = 1 # particle count
E = 0 # energy level
# ===== Creation of the Hamiltonian =====
t = 1
c = np.array([-t if j == 1 else 0 for j in range(M)])
h0 = toeplitz(c)
H_input = SpinlessFermionHamiltonianInput(M=M, N=N, mod_N=0, h0=h0)
H = SpinlessFermionOperator(H_input)
# ===== Computing the Eigenvalues and Eigenvectors =====
k = min(12, H.shape[0] - 1) # Number of eigenvalue to compute
evals, evecs = eigsh(scipy_operator_adapter(H), k=k, which="SA")
print(evals)
plt.figure("eigenvalues", clear=True)
plt.title(f"Eigenvalues (M = {M}, N = {N})")
plt.plot(evals, ".")
# ===== Plotting the Site Occupations =====
def site_occupation(psi, j):
annihilation_input = SpinlessFermionAnnihilationInput(M=M, N=N, mod_N=0, c_i=j)
annihilation = SpinlessFermionOperator(annihilation_input)
c_psi = annihilation.dot(psi)
return c_psi.conj() @ c_psi
occupations = np.array([site_occupation(evecs[:, E], j) for j in range(M)])
plt.figure("site_occupation", clear=True)
plt.title(f"Site Occupation (M = {M}, N = {N}, E = {E}")
plt.plot(occupations, "x-")
plt.axis((0, M - 1, 0, 1))
plt.xlabel("j")
plt.ylabel("Expectation Value")
# ===== Plotting the Energy Occupations =====
J, K = np.meshgrid(np.arange(M), np.arange(M), indexing='ij')
Q = sqrt(2) / sqrt(M+1) * np.sin((J+1) * (K+1) * pi / (M+1))
def energy_occupation(psi, k):
annihilation_input = SpinlessFermionAnnihilationInput(
M=M, N=N, mod_N=0, c_i=Q[:,k].conj()
)
annihilation = SpinlessFermionOperator(annihilation_input)
c_psi = annihilation.dot(psi)
return (c_psi.conj() @ c_psi).real
ks = np.arange(M)
expectation = np.array([energy_occupation(evecs[:, E], k) for k in ks])
plt.figure("energy", clear=True)
plt.title(f"Energy Occupation (M = {M}, N = {N}, E = {E})")
plt.plot(ks, expectation, "x-")
plt.xlabel("p")
plt.ylabel("Expectation Value")