In this tutorial we describe the process of making the molecular data files necessary to run the HFVQE code. We focus on how to use the OpenFermion plugin modules to generate molecular files with canonical Hartree-Fock and generate integrals in a given atomic orbital basis set. We also provide helper functions to run variational Hartree-Fock simulating the experiment and generating initial parameters.
This tutorial will follow the code in recirq/hfvqe/molecular_data/
for constructing MolecularData
objects and getting atomic orbital integrals.
In addition to the standard requirement of ReCirq and its dependencies, this notebook uses OpenFermion-pyscf (and pyscf) to compute some quantities. We install it below if you don't already have it.
View on QuantumAI | Run in Google Colab | View source on GitHub | Download notebook |
Setup
try:
import recirq
except ImportError:
!pip install --quiet git+https://github.com/quantumlib/ReCirq
try:
import openfermionpyscf
except ImportError:
!pip install --quiet openfermionpyscf~=0.5.0
Now we can import the packages required for this notebook.
import os
import numpy as np
import scipy
from recirq.hfvqe.molecular_data.molecular_data_construction import (
h6_linear_molecule, h8_linear_molecule,
h10_linear_molecule, h12_linear_molecule,
get_ao_integrals)
from recirq.hfvqe.gradient_hf import rhf_minimization, rhf_func_generator
from recirq.hfvqe.objective import RestrictedHartreeFockObjective, generate_hamiltonian
def make_rhf_objective(molecule):
S, Hcore, TEI = get_ao_integrals(molecule)
_, X = scipy.linalg.eigh(Hcore, S)
molecular_hamiltonian = generate_hamiltonian(
general_basis_change(Hcore, X, (1, 0)),
numpy.einsum('psqr', general_basis_change(TEI, X, (1, 0, 1, 0)),
molecule.nuclear_repulsion))
rhf_objective = RestrictedHartreeFockObjective(molecular_hamiltonian,
molecule.n_electrons)
return rhf_objective, S, Hcore, TEI
Hydrogen Chain MolecularData
We provide helper functions in the hfvqe
module to generate the Hydrogen chain data. Each chain is constructed using OpenFermion and Psi4 via the OpenFermion-Psi4 plugin. We will use creating H6 with equal spacing between Hydrogen atoms as an example.
from openfermion import MolecularData, general_basis_change
from openfermionpyscf import run_pyscf
n_hydrogens = 6
bond_distance = 1.3 # in Angstroms
molecule = MolecularData(
geometry=[('H', (0, 0, i * bond_distance)) for i in range(n_hydrogens)],
charge=0,
basis='6-31g',
multiplicity=1,
description=f"linear_r-{bond_distance}")
The previous lines set up the MolecularData file. We can now use pyscf to either run a full self-consistent-field Hartree-Fock calculation or get atomic integrals. Via Openfermion-Pyscf we provide an interface to running Hartree-Fock, coupled-cluster, second order perturbation theory, configuration-interaction singles-doubles, and full configuration interaction. Many of these methods depend on parameters such as convergence criteria or initial vectors in the subspace expansion. run_pyscf
assumes common defaults which are appropriate for most systems. Below we will run default Hartree-Fock and CISD.
molecule = run_pyscf(molecule, run_scf=True, run_cisd=True)
print("Hartree-Fock energy:", molecule.hf_energy,
"\nCISD energy:", molecule.cisd_energy)
Hartree-Fock energy: -3.1025910099208867 CISD energy: -3.2221034851022132
The MolecularData
file holds almost all information that is required for post-Hartree-Fock correlated calculations. For example, we provide access to integrals as attributes of MolecularData
.
print("Overlap Ints")
print(molecule.overlap_integrals)
print()
print("One-electron integrals")
print(molecule.one_body_integrals)
Overlap Ints [[1.00000000e+00 6.58292049e-01 1.07450262e-01 2.97950190e-01 2.92936786e-04 2.79155456e-02 1.86658314e-08 5.52805483e-04 2.50532789e-14 2.34306448e-06 7.06186034e-22 2.12361769e-09] [6.58292049e-01 1.00000000e+00 2.97950190e-01 6.14673658e-01 2.79155456e-02 1.42750753e-01 5.52805483e-04 1.25256830e-02 2.34306448e-06 4.15253811e-04 2.12361769e-09 5.20133758e-06] [1.07450262e-01 2.97950190e-01 1.00000000e+00 6.58292049e-01 1.07450262e-01 2.97950190e-01 2.92936786e-04 2.79155456e-02 1.86658314e-08 5.52805483e-04 2.50532789e-14 2.34306448e-06] [2.97950190e-01 6.14673658e-01 6.58292049e-01 1.00000000e+00 2.97950190e-01 6.14673658e-01 2.79155456e-02 1.42750753e-01 5.52805483e-04 1.25256830e-02 2.34306448e-06 4.15253811e-04] [2.92936786e-04 2.79155456e-02 1.07450262e-01 2.97950190e-01 1.00000000e+00 6.58292049e-01 1.07450262e-01 2.97950190e-01 2.92936786e-04 2.79155456e-02 1.86658314e-08 5.52805483e-04] [2.79155456e-02 1.42750753e-01 2.97950190e-01 6.14673658e-01 6.58292049e-01 1.00000000e+00 2.97950190e-01 6.14673658e-01 2.79155456e-02 1.42750753e-01 5.52805483e-04 1.25256830e-02] [1.86658314e-08 5.52805483e-04 2.92936786e-04 2.79155456e-02 1.07450262e-01 2.97950190e-01 1.00000000e+00 6.58292049e-01 1.07450262e-01 2.97950190e-01 2.92936786e-04 2.79155456e-02] [5.52805483e-04 1.25256830e-02 2.79155456e-02 1.42750753e-01 2.97950190e-01 6.14673658e-01 6.58292049e-01 1.00000000e+00 2.97950190e-01 6.14673658e-01 2.79155456e-02 1.42750753e-01] [2.50532789e-14 2.34306448e-06 1.86658314e-08 5.52805483e-04 2.92936786e-04 2.79155456e-02 1.07450262e-01 2.97950190e-01 1.00000000e+00 6.58292049e-01 1.07450262e-01 2.97950190e-01] [2.34306448e-06 4.15253811e-04 5.52805483e-04 1.25256830e-02 2.79155456e-02 1.42750753e-01 2.97950190e-01 6.14673658e-01 6.58292049e-01 1.00000000e+00 2.97950190e-01 6.14673658e-01] [7.06186034e-22 2.12361769e-09 2.50532789e-14 2.34306448e-06 1.86658314e-08 5.52805483e-04 2.92936786e-04 2.79155456e-02 1.07450262e-01 2.97950190e-01 1.00000000e+00 6.58292049e-01] [2.12361769e-09 5.20133758e-06 2.34306448e-06 4.15253811e-04 5.52805483e-04 1.25256830e-02 2.79155456e-02 1.42750753e-01 2.97950190e-01 6.14673658e-01 6.58292049e-01 1.00000000e+00]] One-electron integrals [[-1.87432470e+00 8.30474957e-16 -1.16409641e-01 2.40115605e-16 -6.15072095e-02 -1.04502880e-17 -1.25730762e-01 -7.62100301e-15 -5.89306300e-02 -1.89510298e-15 2.71517033e-03 9.69295149e-16] [ 7.65303957e-16 -1.69304399e+00 7.17422320e-16 1.65551295e-01 -4.80081553e-16 -4.41326519e-02 -5.85602351e-15 1.00494584e-01 -3.51191974e-15 8.00969374e-02 -6.02009053e-16 6.99057737e-03] [-1.16409641e-01 7.85633889e-16 -1.61481536e+00 6.58186268e-16 1.42976433e-01 -1.94644870e-15 1.99568166e-02 8.97271573e-17 -9.93300277e-02 -4.28908794e-15 -8.38901334e-02 7.50068949e-16] [ 1.16002929e-16 1.65551295e-01 6.80250560e-16 -1.47182180e+00 -1.04088434e-15 -1.55554498e-01 -2.90217185e-15 6.81490634e-02 2.98734315e-15 -1.47404734e-01 3.01474057e-15 4.83499815e-02] [-6.15072095e-02 -1.48052550e-16 1.42976433e-01 -6.42123324e-16 -1.29158710e+00 1.73874001e-16 5.71957393e-02 6.42512014e-16 -1.32455830e-01 -7.37015749e-16 1.56350662e-01 -3.08589299e-16] [-7.81003040e-16 -4.41326519e-02 -1.98643038e-15 -1.55554498e-01 6.10304892e-16 -1.24806947e+00 6.19169910e-15 -1.25035099e-01 -1.09022263e-16 1.02022935e-01 9.71194681e-16 2.06045942e-01] [-1.25730762e-01 -5.86275040e-15 1.99568166e-02 -3.03540303e-15 5.71957393e-02 6.27433178e-15 -5.80091614e-01 -1.11597849e-14 -1.60790676e-01 3.74232901e-16 -2.49740453e-02 -2.53549601e-15] [-7.45751485e-15 1.00494584e-01 6.11535851e-17 6.81490634e-02 6.14781730e-16 -1.25035099e-01 -1.14836135e-14 -4.13964503e-01 -7.94898605e-15 -9.28605333e-02 4.01370699e-16 4.59666718e-02] [-5.89306300e-02 -3.52082184e-15 -9.93300277e-02 2.84289510e-15 -1.32455830e-01 -1.26247702e-16 -1.60790676e-01 -7.85377204e-15 -4.68340636e-01 2.80584737e-15 -3.67279868e-02 -1.26066754e-15] [-2.17172175e-15 8.00969374e-02 -4.22628982e-15 -1.47404734e-01 -3.82562175e-16 1.02022935e-01 4.43626024e-16 -9.28605333e-02 2.84974480e-15 -5.21989869e-01 -1.15881513e-15 3.96033476e-02] [ 2.71517033e-03 -7.06811865e-16 -8.38901334e-02 3.07019658e-15 1.56350662e-01 7.31457785e-16 -2.49740453e-02 5.56288697e-16 -3.67279868e-02 -1.14688228e-15 -5.04870445e-01 -2.44815045e-16] [ 6.24473839e-16 6.99057737e-03 1.14463887e-15 4.83499815e-02 5.05984493e-17 2.06045942e-01 -2.59556507e-15 4.59666718e-02 -1.31941582e-15 3.96033476e-02 1.68309057e-16 -4.15306461e-01]]
For the Hartree-Fock experiment we will need to get the atomic basis integrals from the molecular integrals. We can use the identity \(C^{\dagger}SC = I\) to reverse the transformation on the one and two electron integrals.
oei_mo, tei_mo = molecule.one_body_integrals, molecule.two_body_integrals
C = molecule.canonical_orbitals
S = molecule.overlap_integrals
oei_ao = general_basis_change(oei_mo, C.conj().T @ S, key=(1, 0))
print(oei_ao)
[[-1.18698197e+00 -1.09398249e+00 -2.90382070e-01 -5.73478152e-01 -1.91588194e-03 -7.08618226e-02 -2.34028697e-07 -1.97458486e-03 -5.37021242e-13 -1.18105831e-05 -6.67512531e-16 -1.46391011e-08] [-1.09398249e+00 -1.30863711e+00 -6.31325886e-01 -9.69723130e-01 -7.72940958e-02 -2.75295274e-01 -2.08668247e-03 -2.98813560e-02 -1.21151858e-05 -1.23417943e-03 -1.46391013e-08 -1.91822048e-05] [-2.90382070e-01 -6.31325886e-01 -1.51260669e+00 -1.30790890e+00 -3.09964749e-01 -6.50637018e-01 -1.94571033e-03 -7.66467737e-02 -2.34897054e-07 -2.06566594e-03 -5.37947355e-13 -1.21151858e-05] [-5.73478152e-01 -9.69723130e-01 -1.30790890e+00 -1.61455088e+00 -6.69270617e-01 -1.08041337e+00 -7.82781452e-02 -2.89819977e-01 -2.06566594e-03 -3.04640626e-02 -1.18105831e-05 -1.23417943e-03] [-1.91588194e-03 -7.72940958e-02 -3.09964749e-01 -6.69270617e-01 -1.61437154e+00 -1.37489989e+00 -3.14983314e-01 -6.74320764e-01 -1.94571033e-03 -7.82781452e-02 -2.34028698e-07 -2.08668247e-03] [-7.08618226e-02 -2.75295274e-01 -6.50637018e-01 -1.08041337e+00 -1.37489989e+00 -1.71629957e+00 -6.74320764e-01 -1.10900856e+00 -7.66467737e-02 -2.89819977e-01 -1.97458486e-03 -2.98813560e-02] [-2.34028697e-07 -2.08668247e-03 -1.94571033e-03 -7.82781452e-02 -3.14983314e-01 -6.74320764e-01 -1.61437154e+00 -1.37489989e+00 -3.09964749e-01 -6.69270617e-01 -1.91588194e-03 -7.72940958e-02] [-1.97458486e-03 -2.98813560e-02 -7.66467737e-02 -2.89819977e-01 -6.74320764e-01 -1.10900856e+00 -1.37489989e+00 -1.71629957e+00 -6.50637018e-01 -1.08041337e+00 -7.08618226e-02 -2.75295274e-01] [-5.37188257e-13 -1.21151858e-05 -2.34897054e-07 -2.06566594e-03 -1.94571033e-03 -7.66467737e-02 -3.09964749e-01 -6.50637018e-01 -1.51260669e+00 -1.30790890e+00 -2.90382070e-01 -6.31325886e-01] [-1.18105831e-05 -1.23417943e-03 -2.06566594e-03 -3.04640626e-02 -7.82781452e-02 -2.89819977e-01 -6.69270617e-01 -1.08041337e+00 -1.30790890e+00 -1.61455088e+00 -5.73478152e-01 -9.69723130e-01] [-2.74367917e-16 -1.46391010e-08 -5.37802203e-13 -1.18105831e-05 -2.34028697e-07 -1.97458486e-03 -1.91588194e-03 -7.08618226e-02 -2.90382070e-01 -5.73478152e-01 -1.18698197e+00 -1.09398249e+00] [-1.46391010e-08 -1.91822048e-05 -1.21151858e-05 -1.23417943e-03 -2.08668247e-03 -2.98813560e-02 -7.72940958e-02 -2.75295274e-01 -6.31325886e-01 -9.69723130e-01 -1.09398249e+00 -1.30863711e+00]]
print(oei_mo)
[[-1.87432470e+00 8.30474957e-16 -1.16409641e-01 2.40115605e-16 -6.15072095e-02 -1.04502880e-17 -1.25730762e-01 -7.62100301e-15 -5.89306300e-02 -1.89510298e-15 2.71517033e-03 9.69295149e-16] [ 7.65303957e-16 -1.69304399e+00 7.17422320e-16 1.65551295e-01 -4.80081553e-16 -4.41326519e-02 -5.85602351e-15 1.00494584e-01 -3.51191974e-15 8.00969374e-02 -6.02009053e-16 6.99057737e-03] [-1.16409641e-01 7.85633889e-16 -1.61481536e+00 6.58186268e-16 1.42976433e-01 -1.94644870e-15 1.99568166e-02 8.97271573e-17 -9.93300277e-02 -4.28908794e-15 -8.38901334e-02 7.50068949e-16] [ 1.16002929e-16 1.65551295e-01 6.80250560e-16 -1.47182180e+00 -1.04088434e-15 -1.55554498e-01 -2.90217185e-15 6.81490634e-02 2.98734315e-15 -1.47404734e-01 3.01474057e-15 4.83499815e-02] [-6.15072095e-02 -1.48052550e-16 1.42976433e-01 -6.42123324e-16 -1.29158710e+00 1.73874001e-16 5.71957393e-02 6.42512014e-16 -1.32455830e-01 -7.37015749e-16 1.56350662e-01 -3.08589299e-16] [-7.81003040e-16 -4.41326519e-02 -1.98643038e-15 -1.55554498e-01 6.10304892e-16 -1.24806947e+00 6.19169910e-15 -1.25035099e-01 -1.09022263e-16 1.02022935e-01 9.71194681e-16 2.06045942e-01] [-1.25730762e-01 -5.86275040e-15 1.99568166e-02 -3.03540303e-15 5.71957393e-02 6.27433178e-15 -5.80091614e-01 -1.11597849e-14 -1.60790676e-01 3.74232901e-16 -2.49740453e-02 -2.53549601e-15] [-7.45751485e-15 1.00494584e-01 6.11535851e-17 6.81490634e-02 6.14781730e-16 -1.25035099e-01 -1.14836135e-14 -4.13964503e-01 -7.94898605e-15 -9.28605333e-02 4.01370699e-16 4.59666718e-02] [-5.89306300e-02 -3.52082184e-15 -9.93300277e-02 2.84289510e-15 -1.32455830e-01 -1.26247702e-16 -1.60790676e-01 -7.85377204e-15 -4.68340636e-01 2.80584737e-15 -3.67279868e-02 -1.26066754e-15] [-2.17172175e-15 8.00969374e-02 -4.22628982e-15 -1.47404734e-01 -3.82562175e-16 1.02022935e-01 4.43626024e-16 -9.28605333e-02 2.84974480e-15 -5.21989869e-01 -1.15881513e-15 3.96033476e-02] [ 2.71517033e-03 -7.06811865e-16 -8.38901334e-02 3.07019658e-15 1.56350662e-01 7.31457785e-16 -2.49740453e-02 5.56288697e-16 -3.67279868e-02 -1.14688228e-15 -5.04870445e-01 -2.44815045e-16] [ 6.24473839e-16 6.99057737e-03 1.14463887e-15 4.83499815e-02 5.05984493e-17 2.06045942e-01 -2.59556507e-15 4.59666718e-02 -1.31941582e-15 3.96033476e-02 1.68309057e-16 -4.15306461e-01]]
# Use pyscf to get atomic integrals and compare to transformed integrals from above
pyscf_mol = molecule._pyscf_data['mol']
t = pyscf_mol.intor('int1e_kin')
v = pyscf_mol.intor('int1e_nuc')
hcore = t + v
assert np.allclose(hcore, oei_ao)
The two-electron integrals can also be acquired from the pyscf_mol
object or we can use the general_basis_change
to transform the two-electron integrals back into the AO basis.
# Transform the two-electron integrals
tei_ao = general_basis_change(tei_mo, C.conj().T @ S, key=(1, 1, 0, 0))
# re-ordering for chem->physics storage of the integrals
eri_ao = np.einsum('ijkl->iklj', pyscf_mol.intor('int2e', aosym='s1'))
assert np.allclose(tei_ao, eri_ao)
assert not np.allclose(tei_ao, tei_mo)
We also provide a function in recirq.hfvqe.molecular_data.molecular_data_construction
that uses psi4 to generate atomic integrals for the HF-VQE study. Once the atomic orbital integrals are obtained we can perform the first step in setting up the HF-VQE experiment. This involes transforming the integrals to the core-orbital basis and building an RestrictedHartreeFockObjective
.
# diagonalize the AO-core Hamiltonian
_, X = scipy.linalg.eigh(oei_ao, S)
obi = general_basis_change(oei_ao, X, (1, 0))
tbi = np.einsum('psqr', general_basis_change(pyscf_mol.intor('int2e', aosym='s1'), X, (1, 0, 1, 0)))
molecular_hamiltonian = generate_hamiltonian(obi, tbi,
molecule.nuclear_repulsion)
rhf_objective = RestrictedHartreeFockObjective(molecular_hamiltonian,
molecule.n_electrons)
To get initial parameters we can simulate the Hartree-Fock experiment which is performing variational-Hartree-Fock theory.
from recirq.hfvqe.gradient_hf import rhf_minimization, rhf_func_generator
# Uses conjugate gradient to solve non-linear hartree-fock functional:
scipy_result = rhf_minimization(rhf_objective, verbose=True)
Optimization terminated successfully. Current function value: -3.102591 Iterations: 12 Function evaluations: 25 Gradient evaluations: 25
print(molecule.hf_energy, scipy_result.fun)
-3.1025910099208867 -3.1025910099026532
print("Initial Parameters for HF-VQE Study ", scipy_result.x)
Initial Parameters for HF-VQE Study [ 1.58563595e-19 4.73763008e-01 1.95920191e-19 -1.36102259e-02 -7.81094141e-21 4.44832962e-01 5.19837860e-20 -2.59144562e-01 -1.67779731e-19 -1.17062375e-01 1.48974667e-20 -7.61172516e-03 -8.15871275e-20 -9.47442448e-02 6.66883010e-20 -2.71237049e-04 6.23737351e-20 -7.30836704e-02 -4.75683713e-19 4.10271133e-02 1.88267517e-19 -4.52955414e-03 -2.14967064e-20 -8.54284651e-02 -1.54063058e-19 -4.71103578e-02 -3.04830409e-21]
# Explicitly build kappa matrix
from recirq.hfvqe.circuits import rhf_params_to_matrix
import matplotlib.pyplot as plt
kappa = rhf_params_to_matrix(scipy_result.x, len(rhf_objective.occ) + len(rhf_objective.virt), rhf_objective.occ,
rhf_objective.virt)
plt.imshow(kappa)
plt.colorbar()
<matplotlib.colorbar.Colorbar at 0x7f46e0827a60>