Variational Quantum Eigensolver (VQE)

The Variational Quantum Eigensolver (VQE) is a hybrid quantum-classical algorithm designed to find the ground state energy of a quantum system. It combines quantum computation with classical optimization to iteratively improve an approximation of the ground state.

Key features of VQE:

  • Hybrid approach: Utilizes both quantum and classical resources efficiently.

  • Variational method: Uses a parameterized quantum circuit (ansatz) to prepare trial states.

  • Iterative optimization: Classical optimizer adjusts circuit parameters to minimize energy.

  • Flexibility: Can be applied to various problems in quantum chemistry and materials science.

VQE Algorithm Overview:

  1. Prepare an initial quantum state using a parameterized circuit (ansatz).

  2. Measure the expectation value of the Hamiltonian.

  3. Use a classical optimizer to adjust circuit parameters.

  4. Repeat steps 1-3 until convergence or a stopping criterion is met.

CUDA-Q Solvers Implementation

CUDA-Q Solvers provides a high-level interface for running VQE simulations. Here’s how to use it in both Python and C++:

# ============================================================================ #
# Copyright (c) 2024 NVIDIA Corporation & Affiliates.                          #
# All rights reserved.                                                         #
#                                                                              #
# This source code and the accompanying materials are made available under     #
# the terms of the Apache License 2.0 which accompanies this distribution.     #
# ============================================================================ #
import cudaq, cudaq_solvers as solvers
from scipy.optimize import minimize

# Create the molecular hamiltonian
geometry = [('H', (0., 0., 0.)), ('H', (0., 0., .7474))]
molecule = solvers.create_molecule(geometry, 'sto-3g', 0, 0, casci=True)

# Get the number of qubits and electrons
numQubits = molecule.n_orbitals * 2
numElectrons = molecule.n_electrons
spin = 0
initialX = [-.2] * solvers.stateprep.get_num_uccsd_parameters(
    numElectrons, numQubits)


# Define the UCCSD ansatz
@cudaq.kernel
def ansatz(thetas: list[float]):
    q = cudaq.qvector(numQubits)
    for i in range(numElectrons):
        x(q[i])
    solvers.stateprep.uccsd(q, thetas, numElectrons, spin)


# Run VQE
energy, params, all_data = solvers.vqe(ansatz,
                                       molecule.hamiltonian,
                                       initialX,
                                       optimizer=minimize,
                                       method='L-BFGS-B',
                                       jac='3-point',
                                       tol=1e-4,
                                       options={'disp': True})
print(f'Final <H> = {energy}')
/*******************************************************************************
 * Copyright (c) 2024 NVIDIA Corporation & Affiliates.                         *
 * All rights reserved.                                                        *
 *                                                                             *
 * This source code and the accompanying materials are made available under    *
 * the terms of the Apache License 2.0 which accompanies this distribution.    *
 ******************************************************************************/
#include "cudaq.h"
#include "cudaq/solvers/operators.h"
#include "cudaq/solvers/stateprep/uccsd.h"
#include "cudaq/solvers/vqe.h"

// Compile and run with
// nvq++ --enable-mlir -lcudaq-solvers uccsd_vqe.cpp -o uccsd_vqe
// ./uccsd_vqe

int main() {
  // Create the molecular hamiltonian
  cudaq::solvers::molecular_geometry geometry{{"H", {0., 0., 0.}},
                                              {"H", {0., 0., .7474}}};
  auto molecule = cudaq::solvers::create_molecule(
      geometry, "sto-3g", 0, 0, {.casci = true, .verbose = true});

  // Get the spin operator
  auto h = molecule.hamiltonian;

  // Get the number of electrons and qubits
  auto numElectrons = molecule.n_electrons;
  auto numQubits = molecule.n_orbitals * 2;

  // Create an initial set of parameters for the optimization
  auto numParams = cudaq::solvers::stateprep::get_num_uccsd_parameters(
      numElectrons, numQubits);
  std::vector<double> init(numParams, -2.);

  // Run VQE
  auto [energy, thetas, ops] = cudaq::solvers::vqe(
      [&](std::vector<double> params, std::size_t numQubits,
          std::size_t numElectrons) __qpu__ {
        cudaq::qvector q(numQubits);
        for (auto i : cudaq::range(numElectrons))
          x(q[i]);

        cudaq::solvers::stateprep::uccsd(q, params, numElectrons);
      },
      molecule.hamiltonian, init,
      [&](std::vector<double> x) {
        return std::make_tuple(x, numQubits, numElectrons);
      },
      {{"verbose", true}});

  printf("Final <H> = %.12lf\n", energy);
}

Compile and run with

nvq++ --enable-mlir -lcudaq-solvers uccsd_vqe.cpp -o uccsd_vqe
./uccsd_vqe

Code Explanation

  1. Molecule Creation: - Both examples start by defining the molecular geometry (H2 molecule). - The create_molecule function generates the molecular Hamiltonian.

  2. Ansatz Definition: - The UCCSD (Unitary Coupled Cluster Singles and Doubles) ansatz is used. - In Python, it’s defined as a cudaq.kernel. - In C++, it’s defined as a lambda function within the VQE call.

  3. VQE Execution: - The solvers.vqe function (Python) or solvers::vqe (C++) is called. - It takes the ansatz, Hamiltonian, initial parameters, and optimization settings.

  4. Optimization: - Python uses SciPy’s minimize function with L-BFGS-B method. - C++ uses CUDA-Q Solvers’ built-in optimizer. - Either language can make use of CUDA-QX builtin optimizers.

  5. Results: - Both versions print the final ground state energy.

The CUDA-Q Solvers implementation of VQE provides a high-level interface that handles the quantum-classical hybrid optimization loop, making it easy to apply VQE to molecular systems. Users can focus on defining the problem (molecule and ansatz) while CUDA-Q Solvers manages the complex interaction between quantum and classical resources.