Efficiently Simulating the VQE Algorithm with the Quantum Chemistry Toolbox
1. Introduction
mqchem
is an integrated toolbox designed specifically for solving quantum chemistry problems. Its core advantage lies in its simulation method based on Configuration Interaction (CI). Unlike general-purpose simulators that operate in the full \(2^{N}\)-dimensional Hilbert space, mqchem
’s simulator evolves within a subspace defined by a fixed number of electrons (n_electrons
). This approach significantly reduces runtime, enabling us to simulate chemical systems of 20 qubits
or more in a very short amount of time.
The mqchem
module primarily consists of three core components that work in synergy:
``CIHamiltonian``: For efficiently representing the Hamiltonian of a chemical problem.
``UCCExcitationGate``: The fundamental building block for variational ansatze like UCCSD.
``MQChemSimulator``: A high-performance simulator that performs quantum circuit evolution and measurements within the CI space.
In this tutorial, you will learn to follow a natural workflow using the components of the mqchem
module to solve a chemistry problem:
Understand the core components of
mqchem
and their respective responsibilities.Learn how to combine these components to apply a quantum circuit and calculate expectation values.
Utilize the high-level function
prepare_uccsd_vqe
to quickly set up a complete VQE workflow.Work through a hands-on example to calculate the ground state energy of a linear H₆ molecule.
Prerequisites: Before you begin, please ensure you have MindQuantum installed, along with openfermion
and openfermionpyscf
for chemistry calculations.
pip install mindquantum openfermion openfermionpyscf
2. The Core Components of mqchem
When solving a quantum chemistry problem, we typically need to define the problem itself (the Hamiltonian), construct a solution strategy (the ansatz circuit), and execute it in a simulation environment. The three core components of mqchem
correspond precisely to these steps.
2.1 CIHamiltonian
: Defining the Problem
CIHamiltonian
is a wrapper class for FermionOperator
that converts a fermionic Hamiltonian into a format suitable for efficient computation within the CI space. It serves as the objective function for our VQE calculations.
[1]:
from mindquantum.core.operators import FermionOperator
from mindquantum.simulator import mqchem
# Any FermionOperator can be wrapped.
# For example, let's define a simple number operator H = n_0 + n_1
ham_op = FermionOperator('0^ 0') + FermionOperator('1^ 1')
# Wrap it as a CIHamiltonian
ci_ham = mqchem.CIHamiltonian(ham_op)
print("CI Hamiltonian:")
print(ci_ham)
CI Hamiltonian:
1 [0^ 0] +
1 [1^ 1]
2.2 UCCExcitationGate
: Constructing the Solution
The UCCExcitationGate
is the fundamental building block for constructing a UCC (Unitary Coupled-Cluster) ansatz. It represents the unitary operator \(U(\theta) = e^{\theta(T - T^\dagger)}\), where \(T\) is a fermionic excitation operator (e.g., \(a_p^\dagger a_q\)) and \(\theta\) is a trainable parameter.
The simulators in the mqchem
module are specially optimized to apply this type of gate with high efficiency.
[2]:
# Create an excitation operator T = a_2^\dagger a_0 that excites from orbital 0 to 2.
# 'theta' is the name of a trainable parameter.
t_op = FermionOperator('2^ 0', 'theta')
# Create a UCC excitation gate from this operator
ucc_gate = mqchem.UCCExcitationGate(t_op)
print("Created UCC Gate:")
print(ucc_gate)
Created UCC Gate:
exp{(theta)([2^ 0] + [2 0^])}
2.3 MQChemSimulator
: Providing the Execution Environment
The MQChemSimulator
is the engine that executes the simulation. It maintains a quantum state on the CI basis and provides methods to apply UCCExcitationGate
s and calculate the expectation value of a CIHamiltonian
.
When initializing, you need to specify the system’s n_qubits
(total number of spin-orbitals) and n_electrons
(total number of electrons). By default, the simulator starts from the Hartree-Fock (HF) state, which is the most common reference state in quantum chemistry calculations.
[3]:
from mindquantum.core.circuit import Circuit
# Initialize a simulator for a 4-qubit, 2-electron system
sim = mqchem.MQChemSimulator(n_qubits=4, n_electrons=2, seed=42)
print(f"Simulator n_qubits: {sim.n_qubits}")
print(f"Simulator n_electrons: {sim.n_electrons}")
# The default HF state is |0011> (bits are numbered 0, 1, 2, 3 from right to left)
print("\nInitial Hartree-Fock state:")
print(sim.get_qs(ket=True))
Simulator n_qubits: 4
Simulator n_electrons: 2
Initial Hartree-Fock state:
1¦0011⟩
3. Putting It All Together: Simulation and Calculation
Now that we understand the three main components, let’s see how to combine them to execute a complete simulation workflow.
3.1 Applying Gates and Circuits
We can use the apply_gate
or apply_circuit
methods to act with a UCCExcitationGate
on the quantum state maintained by the MQChemSimulator
.
[4]:
# Reset the simulator to the HF state
sim.reset()
# Apply the previously created ucc_gate and assign the value 0.1 to the parameter 'theta'
sim.apply_gate(ucc_gate, pr={'theta': 0.1})
print("State after applying a single gate:")
print(sim.get_qs(ket=True))
# You can also apply a circuit containing multiple gates
sim.reset()
t_op2 = FermionOperator('3^ 1', 0.2) # Use a constant as the parameter
ucc_gate2 = mqchem.UCCExcitationGate(t_op2)
circuit = Circuit([ucc_gate, ucc_gate2])
# Apply the circuit and assign a value to the parameter 'theta'
sim.apply_circuit(circuit, pr={'theta': 0.1})
print("\nState after applying the circuit:")
print(sim.get_qs(ket=True))
State after applying a single gate:
0.99500417¦0011⟩
-0.09983342¦0110⟩
State after applying the circuit:
0.97517033¦0011⟩
-0.0978434¦0110⟩
0.19767681¦1001⟩
0.01983384¦1100⟩
3.2 Calculating Expectation Values
The get_expectation
method is used to calculate the expectation value \(\langle\psi|H|\psi\rangle\) of a CIHamiltonian
with respect to the current quantum state.
For example, let’s verify the expectation value of the Hamiltonian \(H = n_0 + n_1\) in the HF state \(|0011\rangle\). Since orbitals 0 and 1 are both occupied, we expect the result to be \(1+1=2\).
[5]:
sim.reset()
# Use the previously created ci_ham
expectation = sim.get_expectation(ci_ham)
print(f"Expectation value of <HF| (n_0 + n_1) |HF> is: {expectation}")
Expectation value of <HF| (n_0 + n_1) |HF> is: 2.0
4. VQE in Action: Calculating the Ground State Energy of an H₆ Molecule
Now, we will apply what we’ve learned to calculate the ground state energy of a linear H₆ molecule using the VQE algorithm. This is a classic example that showcases the power of the mqchem
module.
4.1 Preparation: Defining the Molecule
We use openfermion
and openfermionpyscf
to obtain the data for the H₆ molecule. run_pyscf
will perform classical quantum chemistry calculations (like HF and CCSD), and these results will serve as input for our quantum algorithm.
Key Point: We must set run_ccsd=True
because the prepare_uccsd_vqe
function requires the results of the CCSD calculation to provide high-quality initial parameters for the VQE.
[6]:
from openfermion import MolecularData
from openfermionpyscf import run_pyscf
# Define an H₆ molecule with a bond length of 1 Angstrom
geometry = [("H", (i, 0, 0)) for i in range(6)]
basis = "sto-3g"
multiplicity = 1 # Spin multiplicity (2S+1), singlet state S=0
charge = 0
# Create the molecular data object and run PySCF
# Note: run_ccsd=True is crucial for obtaining initial amplitudes
molecule = MolecularData(geometry, basis, multiplicity, charge)
mol = run_pyscf(molecule, run_ccsd=True, run_fci=True)
print(f"Molecule: H₆ linear chain")
print(f"Number of qubits (spin-orbitals): {mol.n_qubits}")
print(f"Number of electrons: {mol.n_electrons}")
print(f"Exact ground state energy (FCI): {mol.fci_energy:.8f} Hartrees")
Molecule: H₆ linear chain
Number of qubits (spin-orbitals): 12
Number of electrons: 6
Exact ground state energy (FCI): -3.23606628 Hartrees
4.2 Streamlined VQE Preparation with prepare_uccsd_vqe
prepare_uccsd_vqe
is a major highlight of the mqchem
module. This high-level factory function automates the most tedious steps of VQE preparation, encapsulating the component creation process we learned earlier:
Generates the fermionic Hamiltonian from the molecular data and wraps it as a
CIHamiltonian
.Constructs a parameterized
Circuit
based on the UCCSD (Unitary Coupled-Cluster Singles and Doubles) theory, containing a series ofUCCExcitationGate
s.Extracts the coupled-cluster amplitudes from the classical CCSD calculation to use as initial parameters for the VQE optimization. This is often an excellent starting point.
[7]:
# Generate the Hamiltonian, UCCSD ansatz circuit, and initial parameters in one step
# threshold=1e-6 is used to filter out excitation terms with small CCSD amplitudes, simplifying the ansatz
hamiltonian, ansatz_circuit, initial_amplitudes = mqchem.prepare_uccsd_vqe(mol, threshold=1e-6)
print("CI Hamiltonian (first few terms):")
print(str(hamiltonian.fermion_hamiltonian)[:200] + "...")
print(f"\nGenerated Ansatz circuit contains {len(ansatz_circuit)} gates.")
print("Circuit parameter names:", ansatz_circuit.params_name)
print(f"\nInitial amplitudes extracted from CCSD (total {len(initial_amplitudes)}):")
print(initial_amplitudes)
CI Hamiltonian (first few terms):
4.6038 [] +
-2.2818 [0^ 0] +
-0.1459 [4 0^] +
0.0642 [8 0^] +
-2.2818 [1^ 1] +
-0.1459 [5 1^] +
0.0642 [9 1^] +
-2.0409 [2^ 2] +
-0.2111 [6 2^] +
0.0417 [10 2^] +
-2.0409 [3^ 3] +
-0.2111 [7 3^] +...
Generated Ansatz circuit contains 59 gates.
Circuit parameter names: ['d1_0', 's_1', 'd1_1', 'd1_2', 's_3', 'd1_3', 'd1_4', 's_5', 'd1_5', 'd1_6', 's_7', 'd1_7', 'd1_8', 'd2_1', 'd2_3', 'd2_9', 'd2_5', 'd2_18', 'd2_7', 'd2_11', 'd2_16', 'd2_13', 'd2_20', 'd2_22', 'd2_24', 'd2_27', 'd2_29', 'd2_31', 'd2_34']
Initial amplitudes extracted from CCSD (total 29):
[-0.02172895 -0.00966027 -0.03288994 -0.10380602 0.0028648 -0.0143801
-0.03201595 -0.00681583 -0.02420889 -0.01695636 -0.00239016 -0.01001419
-0.01221332 0.00432926 0.01629623 0.00904958 -0.00799714 0.01818314
0.03079327 0.04406509 0.02507472 0.01007596 -0.00356389 -0.00956368
-0.01592823 -0.01108764 -0.01197809 -0.00586536 0.00582737]
4.3 Getting the Expectation and Gradient Function
The core of VQE is to minimize the energy expectation value. Efficient optimization algorithms (like L-BFGS-B) require the gradient of the energy with respect to the parameters. The get_expectation_with_grad
method is designed precisely for this. It returns a function that takes parameter values and returns both the energy expectation value and the gradient simultaneously.
[8]:
# Initialize a new simulator for the VQE calculation
vqe_sim = mqchem.MQChemSimulator(mol.n_qubits, mol.n_electrons)
# Get the function that can calculate both expectation and gradient
# This function is compatible with the interface of scipy.optimize.minimize
grad_ops = vqe_sim.get_expectation_with_grad(hamiltonian, ansatz_circuit)
# We can test this function
energy, gradient = grad_ops(initial_amplitudes)
print(f"Energy calculated with initial parameters: {energy:.8f}")
print(f"Corresponding gradient (first 5 elements): {gradient[:5]}")
Energy calculated with initial parameters: -3.21777544
Corresponding gradient (first 5 elements): [-0.03673877 0.02322724 -0.06124962 -0.10907011 -0.00220752]
4.4 Running the Optimizer
Finally, we use a classical optimizer, scipy.optimize.minimize
, to find the optimal parameters and thus obtain the variational approximation of the ground state energy.
[9]:
from scipy.optimize import minimize
from time import time
t0 = time()
# Perform optimization using the L-BFGS-B algorithm
# jac=True indicates that our grad_ops function returns both the value and the Jacobian (gradient)
result = minimize(grad_ops, initial_amplitudes, method="L-BFGS-B", jac=True)
vqe_energy = result.fun
print("\n--- VQE Optimization Complete ---")
print(f"VQE calculated ground state energy: {vqe_energy:.8f} Hartrees")
print(f"Exact FCI ground state energy: {mol.fci_energy:.8f} Hartrees")
print(f"Error: {abs(vqe_energy - mol.fci_energy):.2e} Hartrees")
print(f"Runtime: {time() - t0:.2f} seconds")
--- VQE Optimization Complete ---
VQE calculated ground state energy: -3.23544945 Hartrees
Exact FCI ground state energy: -3.23606628 Hartrees
Error: 6.17e-04 Hartrees
Runtime: 0.02 seconds
As you can see, the VQE result is very close to the exact FCI energy, and the entire optimization process completed in under 0.02 seconds! This demonstrates the correctness and high performance of the mqchem
components working together.
5. Summary
In this tutorial, we systematically introduced the mqchem
module. We started with its three core components—CIHamiltonian
, UCCExcitationGate
, and MQChemSimulator
—to understand their respective roles in the quantum chemistry simulation workflow.
We not only learned how to use these components individually but, more importantly, saw how they seamlessly integrate to solve practical problems. Finally, through a VQE example of calculating the ground state energy of an H₆ molecule, we demonstrated how to use convenient tools like prepare_uccsd_vqe
to efficiently build and run a complete quantum chemistry simulation pipeline.
The mqchem
module provides powerful support for high-performance, high-precision quantum chemistry research on the MindQuantum platform. We hope this tutorial helps you get started quickly and encourages you to explore more interesting chemistry problems