ADAPT-VQE

ADAPT-VQE is an advanced quantum algorithm designed to improve upon the standard Variational Quantum Eigensolver (VQE) approach for solving quantum chemistry problems. It addresses key challenges faced by traditional VQE methods by dynamically constructing a problem-specific ansatz, offering several advantages:

  • Faster convergence: Adaptively selects the most impactful operators, potentially achieving convergence more quickly than fixed-ansatz VQE methods.

  • Enhanced efficiency: Builds a compact ansatz tailored to the specific problem, potentially reducing overall circuit depth.

  • Increased accuracy: Has demonstrated the ability to outperform standard VQE approaches in terms of accuracy for certain molecular systems.

  • Adaptability: Automatically adjusts to different molecular systems without requiring significant user intervention or prior knowledge of the system’s electronic structure.

The ADAPT-VQE algorithm works by iteratively growing the quantum circuit ansatz, selecting operators from a predefined pool based on their gradient magnitudes. This adaptive approach allows the algorithm to focus computational resources on the most relevant aspects of the problem, potentially leading to more efficient and accurate simulations of molecular systems on quantum computers.

Here we demonstrate how to use the CUDA-Q Solvers library to execute the ADAPT-VQE algorithm.

# ============================================================================ #
# 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

# Run this script with
# python3 adapt_h2.py
#
# In order to leverage CUDA-Q MQPU and distribute the work across
# multiple QPUs (thereby observing a speed-up), set the target and
# use MPI:
#
# cudaq.set_target('nvidia', mqpu=True)
# cudaq.mpi.initialize()
#
# run with
#
# mpiexec -np N and vary N to see the speedup...
# e.g. mpiexec -np 2 python3 adapt_h2_mqpu.py
#
# End the script with
# cudaq.mpi.finalize()

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

# Create the ADAPT operator pool
operators = solvers.get_operator_pool("spin_complement_gsd",
                                      num_orbitals=molecule.n_orbitals)

# Get the number of electrons so we can
# capture it in the initial state kernel
numElectrons = molecule.n_electrons


# Define the initial Hartree Fock state
@cudaq.kernel
def initState(q: cudaq.qview):
    for i in range(numElectrons):
        x(q[i])


# Run ADAPT-VQE
energy, thetas, ops = solvers.adapt_vqe(initState, molecule.hamiltonian,
                                        operators)

# Print the result.
print("<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/adapt.h"
#include "cudaq/solvers/operators.h"

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

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;

  // Create the operator pool
  std::vector<cudaq::spin_op> opPool = cudaq::solvers::get_operator_pool(
      "spin_complement_gsd", {{"num-orbitals", h.num_qubits() / 2}});

  // Run ADAPT
  auto [energy, thetas, ops] = cudaq::solvers::adapt_vqe(
      [](cudaq::qvector<> &q) __qpu__ {
        x(q[0]);
        x(q[1]);
      },
      h, opPool, {{"grad_norm_tolerance", 1e-3}});

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

Compile and run with

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