VQE Electronic structure problem#

Problem description#

Consider the problem of finding approximate solutions to the non-relativistic time-independent Schrödinger equation

H|Ψ=E|Ψ

The Hamiltonian H for a system of M nuclei and N electrons with coordinates RAR3 and riR3, respectively, is defined as

H=i=1N12i2i=1M12MAA2i=1NA=1MZAriA+i=1Nj>iN1rij+A=1MB>AMZAZBRAB

where riA=|riRA|, rij=|rirj|, and RAB=|RARB|, are the distances between electrons and nuclei, MA is the ratio of the mass of nucleus A to an electron, and ZA is the atomic number of nulecus A.

Within the Born-Oppenheimer approximation, one considers the electrons in a molecule to be moving a the field of fixed nuclei. Therefore, the kinetic energy of the nuclei is neglected, and the nuclear repulsion energy is constant.

The Hamiltonian describing the motion of N electrons in the field of M fixed nuclei is

Hel=i=1N12i2i=1NA=1MZAriA+i=1Nj>iN1rij

The electronic stucture problem consists of finding solutions to the Schrödinger equation for the electronic Hamiltonian

Hel|Ψel=Eel|Ψel

where |Ψel is the electronic wave function.

In the following, we focus on finding the ground state energy

E0=Ψ0|Hel|Ψ0

where |Ψ0 is a ground state wave function of the system.

A starting point for solving this problem is the (restricted) Hartree-Fock method which produces a set {ψi} of K spacial orbitals that correspond to M=2K spin orbitals χ2i=ψiα and χ2i+1=ψiβ, for i=0,,K1, for the N electrons of the molecule.

Utilizing the second quantization formalism, the electronic Hamiltonian is then expressed in the basis of the solutions (i.e., the spin orbitals) of the Hartree-Fock method:

H^el=i,j=0M1hijaiaj+12i,j,k,l=0M1hi,j,k,laiajakal

where ai and ai are the fermionic creation and annihilation operators, respectively.

Here, the coefficients hij and hijkl are one- and two-electron integrals which can be computed classically:

  • One-electron integrals:

    hi,j=dxχi(x)χj(x)
  • Two-electron integrals:

    hi,j,k,l=dx1dx2χi(x1)χj(x2)1|r1r2|χk(x1)χl(x2)

Note: There is a difference between the physicists’s notation (above) and the chemists’ notation for the two-electron integrals!

The Hartree-Fock state

|ΨHF=|10,,1N1,0N,,0M1

where the first (i.e., the lowest energy) N orbitals are occupied is the best approximation of the ground state |Ψ^0 of the Hamiltonian H^el in this form.

All feasible N-electron states are expressed as a superposition of the Hartree-Fock state and

  • single electron excitation states:

    |Ψir=|10,,0i,,1N10N,,1r,,0M1

    One electron is moved from the occupied orbital χi to the unoccupied (virtual) orbital χr.

  • double electron excitation states:

    |Ψijrs=|10,,0i,,0j,,1N1,0N,,1r,,1s,,0M1

    Two electrons are moved from the occupied orbitals χi,χj to the unoccupied (virtual) orbitals χr,χs.

  • higher order (triple, quadruple, ect.) excitation states.

That is, a ground state can be expressed as

|Ψ^0=c0|ΨHF+i,rcir|Ψir+i<j,r<scijrs|Ψir+

Solving the electronic structure problem, i.e., finding a ground state |Ψ^0 of the electronic Hamiltonian H^el within the feasible N-electron subspace, with a quantum computer requires transforming the fermionic representation into a qubit representation. This is achieved by, e.g., the Jordan-Wigner, Parity, or Bravyi-Kitaev transformation.

Electronic structure problem#

electronic_structure_problem(arg, active_orb=None, active_elec=None, ansatz_type='QCCSD', threshold=0.0001)[source]#

Creates a VQE problem instance for an electronic structure problem defined by the one-electron and two-electron integrals for the spin orbitals (in physicists’ notation).

The problem Hamiltonian is given by:

H=i,j=0M1hi,jaiaj+i,j,k,l=0M1hi,j,k,laiajakal

for one-electron integrals:

hi,j=dxχi(x)χj(x)

and two-electron integrals:

hi,j,k,l=dx1dx2χi(x1)χj(x2)1|r1r2|χk(x1)χl(x2)
Parameters:
argpyscf.gto.Mole or dict

A PySCF molecule or a dictionary specifying the electronic data for a molecule. The following data is required:

  • one_intnumpy.ndarray

    The one-electron integrals w.r.t. spin orbitals (in physicists’ notation).

  • two_intnumpy.ndarray

    The two-electron integrals w.r.t. spin orbitals (in physicists’ notation).

  • num_orbint

    The number of spin orbitals.

  • num_elecint

    The number of electrons.

active_orbint, optional

The number of active spin orbitals.

active_elecint, optional

The number of active electrons.

ansatz_typestring, optional

The ansatz type. Availabe is QCCSD. The default is QCCSD.

thresholdfloat, optional

The threshold for the absolute value of the coefficients of Pauli products in the quantum Hamiltonian. The default is 1e-4.

Returns:
VQEProblem

The VQE problem instance.

Examples

We calculate the electronic energy for the Hydrogen molecule at bond distance 0.74 angstroms:

from pyscf import gto
from qrisp import QuantumVariable
from qrisp.vqe.problems.electronic_structure import *

mol = gto.M(
    atom = '''H 0 0 0; H 0 0 0.74''',
    basis = 'sto-3g')

vqe = electronic_structure_problem(mol)
vqe.set_callback()

energy = vqe.run(QuantumVariable(4),depth=1,max_iter=50)
print(energy)
#Yields -1.8461290172512965

Helper functions#

electronic_data(mol)[source]#

A function that utilizes restricted Hartree-Fock (RHF) calculation in the PySCF quantum chemistry package to obtain the electronic data for defining an electronic structure problem.

Parameters:
molpyscf.gto.Mole

The molecule.

Returns:
datadict

A dictionary specifying the electronic data for a molecule. The following data is provided:

  • one_intnumpy.ndarray

    The one-electron integrals w.r.t. spin orbitals (in physicists’ notation).

  • two_intnumpy.ndarray

    The two-electron integrals w.r.t. spin orbitals (in physicists’ notation).

  • num_orbint

    The number of spin orbitals.

  • num_elecint

    The number of electrons.

  • energy_nuc

    The nuclear repulsion energy.

  • energy_hf

    The Hartree-Fock ground state energy.

Hamiltonian#

create_electronic_hamiltonian(arg, active_orb=None, active_elec=None)[source]#

Creates the qubit Hamiltonian for an electronic structure problem. If an Active Space (AS) is specified, the Hamiltonian is calculated following this paper.

Parameters:
argpyscf.gto.Mole or dict

A PySCF molecule or a dictionary specifying the electronic data for a molecule. The following data is required:

  • one_intnumpy.ndarray

    The one-electron integrals w.r.t. spin orbitals (in physicists’ notation).

  • two_intnumpy.ndarray

    The two-electron integrals w.r.t. spin orbitals (in physicists’ notation).

  • num_orbint

    The number of spin orbitals.

  • num_elecint

    The number of electrons.

active_orbint, optional

The number of active spin orbitals.

active_elecint, optional

The number of active electrons.

Returns:
HFermionicOperator

The fermionic Hamiltonian.

Examples

We calucalte the fermionic Hamiltonian for the Hydrogen molecule, and transform it to a Pauli Hamiltonian via Jordan-Wigner transform.

from pyscf import gto
from qrisp.vqe.problems.electronic_structure import *

mol = gto.M(
    atom = '''H 0 0 0; H 0 0 0.74''',
    basis = 'sto-3g')

H = create_electronic_hamiltonian(mol)
H.to_qubit_operator()   

Yields:

0.8121706072487140.0453026155037992X0X1Y2Y3+0.0453026155037992X0Y1Y2X30.0453026155037992Y0Y1X2X3+0.171412826447769Z0+0.168688981703612Z0Z1+0.120625234833904Z0Z2+0.165927850337703Z0Z3+0.171412826447769Z1+0.165927850337703Z1Z2+0.120625234833904Z1Z30.223431536908133Z2+0.174412876122615ZZ2Z30.223431536908133Z3

Ansatz#

create_QCCSD_ansatz(M, N)[source]#

This method creates a function for applying one layer of the QCCSD ansatz.

The chemistry-inspired Qubit Coupled Cluster Single Double (QCCSD) ansatz evolves the initial state, usually, the Hartree-Fock state

|ΨHF=|10,,1N,0N+1,,0M

under the action of parametrized (non-commuting) single and double excitation unitaries.

The single (S) excitation unitaries Uir implement a continuous swap for qubits i and r:

Uir(θ)=(10000cos(θ)sin(θ)00sin(θ)cos(θ)00001)

Similarly, the double (D) excitation unitaries Uijrs(θ) implement a continuous swap for qubit pairs i,j and r,s.

Parameters:
Mint

The number of (active) spin orbitals.

Nint

The number of (active) electrons.

Returns:
ansatzfunction

A function that can be applied to a QuantumVariable and a list of parameters.

num_paramsint

The number of parameters.

create_hartree_fock_init_function(M, N)[source]#

Creates the function that, when applied to a QuantumVariable, initializes the Hartee-Fock state: Consistent with the Jordan-Wigner mapping, the first N qubits are initialized in the |1 state.

Parameters:
Mint

The number of (active) spin orbitals.

Nint

The number of (active) electrons.

Returns:
init_functionfunction

A function that can be applied to a QuantumVariable.