VQE Electronic structure problem#

Problem description#

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

\[H\ket{\Psi}=E\ket{\Psi}\]

The Hamiltonian \(H\) for a system of \(M\) nuclei and \(N\) electrons with coordinates \(R_A\in\mathbb R^3\) and \(r_i\in\mathbb R^3\), respectively, is defined as

\[H= -\sum_{i=1}^N\frac12\nabla_i^2-\sum_{i=1}^M\frac{1}{2M_A}\nabla_A^2-\sum_{i=1}^N\sum_{A=1}^M\frac{Z_A}{r_{iA}}+\sum_{i=1}^N\sum_{j>i}^N\frac{1}{r_{ij}}+\sum_{A=1}^M\sum_{B>A}^M\frac{Z_AZ_B}{R_{AB}}\]

where \(r_{iA}=|r_i-R_A|\), \(r_{ij}=|r_i-r_j|\), and \(R_{AB}=|R_A-R_B|\), are the distances between electrons and nuclei, \(M_A\) is the ratio of the mass of nucleus \(A\) to an electron, and \(Z_A\) 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

\[H_{\text{el}} = -\sum_{i=1}^N\frac12\nabla_i^2-\sum_{i=1}^N\sum_{A=1}^M\frac{Z_A}{r_{iA}}+\sum_{i=1}^N\sum_{j>i}^N\frac{1}{r_{ij}}\]

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

\[H_{\text{el}}\ket{\Psi_{\text{el}}}=E_{\text{el}}\ket{\Psi_{\text{el}}}\]

where \(\ket{\Psi_{\text{el}}}\) is the electronic wave function.

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

\[E_0 = \braket{\Psi_0|H_{\text{el}}|\Psi_0}\]

where \(\ket{\Psi_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 \(\{\psi_i\}\) of \(K\) spacial orbitals that correspond to \(M=2K\) spin orbitals \(\chi_{2i}=\psi_{i}\alpha\) and \(\chi_{2i+1}=\psi_{i}\beta\), for \(i=0,\dotsc,K-1\), 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:

\[\hat H_{\text{el}} = \sum\limits_{i,j=0}^{M-1}h_{ij}a_i^{\dagger}a_j+\frac{1}{2}\sum\limits_{i,j,k,l=0}^{M-1}h_{i,j,k,l}a_i^{\dagger}a_j^{\dagger}a_ka_l\]

where \(a_i^{\dagger}\) and \(a_i\) are the fermionic creation and annihilation operators, respectively.

Here, the coefficients \(h_{ij}\) and \(h_{ijkl}\) are one- and two-electron integrals which can be computed classically:

  • One-electron integrals:

    \[h_{i,j}=\int\mathrm dx \chi_i^*(x)\chi_j(x)\]
  • Two-electron integrals:

    \[h_{i,j,k,l} = \int\mathrm dx_1\mathrm dx_2\chi_i^*(x_1)\chi_j^*(x_2)\frac{1}{|r_1-r_2|}\chi_k(x_1)\chi_l(x_2)\]

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

The Hartree-Fock state

\[\ket{\Psi_{\text{HF}}} = \ket{1_0,\dotsc,1_{N-1},0_N,\dotsc,0_{M-1}}\]

where the first (i.e., the lowest energy) \(N\) orbitals are occupied is the best approximation of the ground state \(\ket{\hat\Psi_0}\) of the Hamiltonian \(\hat H_{\text{el}}\) in this form.

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

  • single electron excitation states:

    \[\ket{\Psi_i^r}=\ket{1_0,\dotsc,0_i,\dotsc,1_{N-1}0_{N},\dotsc,1_r,\dotsc,0_{M-1}}\]

    One electron is moved from the occupied orbital \(\chi_i\) to the unoccupied (virtual) orbital \(\chi_r\).

  • double electron excitation states:

    \[\ket{\Psi_{ij}^{rs}}=\ket{1_0,\dotsc,0_i,\dotsc,0_j,\dotsc,1_{N-1},0_{N},\dotsc,1_r,\dotsc,1_s,\dotsc,0_{M-1}}\]

    Two electrons are moved from the occupied orbitals \(\chi_i, \chi_j\) to the unoccupied (virtual) orbitals \(\chi_r, \chi_s\).

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

That is, a ground state can be expressed as

\[\ket{\hat\Psi_{0}}=c_0\ket{\Psi_{\text{HF}}}+\sum_{i,r}c_i^r\ket{\Psi_i^r}+\sum_{i<j,r<s}c_{ij}^{rs}\ket{\Psi_i^r}+\dotsb\]

Solving the electronic structure problem, i.e., finding a ground state \(\ket{\hat\Psi_0}\) of the electronic Hamiltonian \(\hat H_{\text{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 = \sum\limits_{i,j=0}^{M-1}h_{i,j}a^{\dagger}_ia_j + \sum\limits_{i,j,k,l=0}^{M-1}h_{i,j,k,l}a^{\dagger}_ia^{\dagger}_ja_ka_l\]

for one-electron integrals:

\[h_{i,j}=\int\mathrm dx \chi_i^*(x)\chi_j(x)\]

and two-electron integrals:

\[h_{i,j,k,l} = \int\mathrm dx_1\mathrm dx_2\chi_i^*(x_1)\chi_j^*(x_2)\frac{1}{|r_1-r_2|}\chi_k(x_1)\chi_l(x_2)\]
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:

\[ \begin{align}\begin{aligned}-&0.812170607248714 - 0.0453026155037992 X_0X_1Y_2Y_3 + 0.0453026155037992 X_0Y_1Y_2X_3 - 0.0453026155037992 Y_0Y_1X_2X_3\\&+0.171412826447769 Z_0 + 0.168688981703612 Z_0Z_1 + 0.120625234833904 Z_0Z_2 + 0.165927850337703 Z_0Z_3 + 0.171412826447769 Z_1\\&+0.165927850337703 Z_1Z_2 + 0.120625234833904 Z_1Z_3 - 0.223431536908133 Z_2 + 0.174412876122615Z Z_2Z_3 - 0.223431536908133 Z_3\end{aligned}\end{align} \]

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

\[\ket{\Psi_{\text{HF}}}=\ket{1_0,\dotsc,1_N,0_{N+1},\dotsc,0_M}\]

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

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

\[\begin{split}U_i^r(\theta) = \begin{pmatrix} 1&0&0&0\\ 0&\cos(\theta)&-\sin(\theta)&0\\ 0&\sin(\theta)&\cos(\theta)&0\\ 0&0&0&1 \end{pmatrix}\end{split}\]

Similarly, the double (D) excitation unitaries \(U_{ij}^{rs}(\theta)\) 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 \(\ket{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.