Using the Exact Geodesic Transport with Conjugate Gradients(EGT-CG) Optimizer#

The Exact Geodesic Transport with Conjugate Gradients (EGT-CG) optimizer is a curvature-aware Riemannian optimizer designed specifically for variational circuits based on the Hamming-weight encoder (HWE) ansatze.

It updates parameters along exact geodesic paths on the hyperspherical manifold defined by the HWE, combining analytic metric computation, conjugate-gradient memory, and dynamic learning rates for fast, globally convergent optimization.

For more details, see:

    1. Ferreira‑Martins et al., Variational quantum algorithms with exact geodesic transport,

arXiv:2506.17395 (2025).

The current configuration focuses on ground state estimation problems where given a variational state \(\ket{\psi(\boldsymbol{\theta})}\) parameterized by hyperspherical angles \(\boldsymbol{\theta}\),

the loss function minimized is the energy expectation value of the given hamiltonian:

\[\mathcal{L}(\boldsymbol{\theta}) = \bra{\psi(\boldsymbol{\theta})} \, H \, \ket{\psi(\boldsymbol{\theta})} \, ,\]

But it is important to note that this framework generalizes to arbitrary loss functions.

Defining and running the optimizer#

In qiboml, the EGT-CG optimizer is implemented in the class qiboml.models.optimizers.ExactGeodesicTransportCG.

class qiboml.models.optimizers.ExactGeodesicTransportCG(nqubits: int, weight: int, hamiltonian, angles, backtrack_rate: Optional[float] = None, geometric_gradient: bool = False, multiplicative_factor: float = 1.0, c1: float = 0.0001, c2: float = 0.9, backend=None)[source]#

Bases: object

Exact Geodesic Transport with Conjugate Gradients Optimizer.

Implements the Exact Geodesic Transport with Conjugate Gradients (EGT-CG) optimizer, a curvature-aware Riemannian optimizer designed specifically for variational circuits based on the Hamming-weight encoder (HWE) ansatze. It updates parameters along exact geodesic paths on the hyperspherical manifold defined by the HWE, combining analytic metric computation, conjugate-gradient memory, and dynamic learning rates for fast, globally convergent optimization.

Parameters:
  • nqubits (int) – Number of qubits in the quantum circuit.

  • weight (int) – Hamming weight to encode.

  • ( (hamiltonian) – class:’qibo.hamiltonians.Hamiltonian’): Hamiltonian whose expectation value defines the loss to minimize.

  • angles (ndarray) – Initial hyperspherical angles parameterizing the amplitudes.

  • backtrack_rate (float, optional) – Backtracking rate for Wolfe condition line search. If None, defaults to \(0.5\). Defaults to None.

  • geometric_gradient (bool, optional) – If True, uses the geometric gradient estimator instead of numerical finite differences. Defaults to False.

  • multiplicative_factor (float, optional) – Scaling factor applied to the initial learning rate during optimization. Defaults to \(1\).

  • c1 (float, optional) – Constant for Armijo condition (sufficient decrease) in Wolfe line search. Defaults to \(10^{-3}\).

  • c2 (float, optional) – Constant for curvature condition in Wolfe line search. It should satisfy c1 < c2 < 1. Defaults to \(0.9\).

  • backend (qibo.backends.abstract.Backend, optional) – backend to be used in the execution. If None, it uses the current backend. Defaults to None.

Returns:

Instantiated optimizer object.

Return type:

ExactGeodesicTransportCG

References

A. J. Ferreira‑Martins, R. M. S. Farias, G. Camilo, T. O. Maciel, A. Tosta, R. Lin, A. Alhajri, T. Haug, and L. Aolita, Variational quantum algorithms with exact geodesic transport, arXiv:2506.17395 (2025).

initialize_cg_state()[source]#

Initialize CG state.

Sets up the internal variables x, u, v, and initial step size eta based on the current angles.

angles_to_amplitudes(angles)[source]#

Convert angles to amplitudes.

Parameters:

angles (ndarray) – Angles in hyperspherical coordinates.

Returns:

Amplitudes calculated from the hyperspherical coordinates.

Return type:

ndarray

encoder()[source]#

Build and return the Hamming-weight encoder circuit for the given amplitudes.

Returns:

Circuit prepared with current amplitudes.

Return type:

circuit (qibo.models.hamming_weight_encoder)

state(initial_state=None, nshots=1000)[source]#

Return the statevector after encoding.

Parameters:
  • initial_state (ndarray, optional) – Initial statevector. Defaults to None.

  • nshots (int, optional) – Number of measurement shots. Defaults to 1000.

Returns:

Statevector of the encoded quantum state.

Return type:

statevector (ndarray)

loss()[source]#

Loss function to be minimized.

Given a quantum state \(\ket{\psi}\) and a Hamiltonian \(H\), the loss function is defined as

\[\mathcal{L} = \bra{\psi} \, H \, \ket{\psi} \, .\]
Returns:

Expectation value of hamiltonain.

Return type:

float

gradient(epsilon=1e-08)[source]#

Numerically compute gradient of loss wrt angles.

Returns:

Gradient of loss w.r.t. angles.

Return type:

ndarray

amplitudes_to_full_state(amps)[source]#

Convert amplitudes to the full quantum statevector.

Parameters:

amps (ndarray) – Amplitude vector.

Returns:

Statevector corresponding to the given amplitudes.

Return type:

ndarray

geom_gradient()[source]#

Compute geometric gradient using the diagonal metric tensor and Jacobian.

Returns:

Geometric gradient vector.

Return type:

ndarray

jacobian()[source]#

Compute Jacobian of amplitudes wrt angles.

Returns:

Jacobian matrix.

Return type:

ndarray

metric_tensor()[source]#

Compute the diagonal metric tensor in hyperspherical coordinates.

Returns:

Diagonal elements of the metric tensor.

Return type:

ndarray

tangent_vector()[source]#

Compute the Riemannian gradient (tangent vector) at the current point on the hypersphere.

Returns:

Tangent vector in the tangent space of the hypersphere.

Return type:

ndarray

optimize_step_size(x_prev, u_prev, v_prev, loss_prev)[source]#

Perform Wolfe line search to determine optimal step size eta.

Parameters:
  • x_prev (ndarray) – Previous position on the sphere.

  • u_prev (ndarray) – Previous search direction.

  • v_prev (ndarray) – Previous gradient vector.

  • loss_prev (float) – Loss at previous position.

Returns:

Respectively, updated position, angles, gradient, and step size.

Return type:

tuple

exponential_map_with_direction(direction, eta=None)[source]#

Exponential map from current point along specified direction.

Parameters:
  • direction (ndarray) – Tangent vector direction.

  • eta (float, optional) – Step size. Defaults to current eta.

Returns:

New point on the hypersphere.

Return type:

ndarray

amplitudes_to_angles(x)[source]#

Convert amplitude vector back to hyperspherical angles.

Parameters:

x (ndarray) – Amplitude vector.

Returns:

Corresponding angles.

Return type:

ndarray

parallel_transport(u, v, a, eta=None)[source]#

Parallel transport a tangent vector u along geodesic defined by v.

Parameters:
  • u (ndarray) – Vector to transport.

  • v (ndarray) – Direction of geodesic.

  • a (ndarray) – Starting point on sphere.

  • eta (float, optional) – Step size. If None, defaults to current eta. Defaults to None.

Returns:

Transported vector.

Return type:

ndarray

sphere_inner_product(u, v, x)[source]#

Compute inner product on tangent space at x on the sphere.

Parameters:
  • u (ndarray) – First tangent vector.

  • v (ndarray) – Second tangent vector.

  • x (ndarray) – Base point on the sphere.

Returns:

Inner product value.

Return type:

float

beta_dy(v_next, x_next, transported_u, st)[source]#

Compute Dai and Yuan Beta.

Parameters:
  • v_next (ndarray) – Next gradient.

  • x_next (ndarray) – Next point.

  • transported_u (ndarray) – Parallel-transported u.

  • st (float) – Scaling factor.

Returns:

Dai-Yuan beta value.

Return type:

float

beta_hs(v_next, x_next, transported_u, transported_v, lt, st)[source]#

Compute Hestenes-Stiefel conjugate gradient beta.

Parameters:
  • v_next (ndarray) – Next gradient.

  • x_next (ndarray) – Next point.

  • transported_u (ndarray) – Parallel-transported u.

  • transported_v (ndarray) – Parallel-transported v.

  • lt (float) – Scaling factor.

  • st (float) – Scaling factor.

Returns:

Hestenes-Stiefel beta value.

Return type:

float

run_egt_cg(steps: int = 10, tolerance: float = 1e-08)[source]#

Run the EGT-CG optimizer for a specified number of steps.

Parameters:
  • steps (int, optional) – Number of optimization iterations. Defaults to \(10\).

  • tolerance (float, optional) – Maximum tolerance for the residue of the gradient update. Defaults to \(10^{-8}\).

Returns:

(final_loss, losses, final_parameters) final_loss (float): Final loss value. losses (list): Loss at each iteration. final_parameters (ndarray): Final angles.

Return type:

tuple

It requires as inputs the number of qubits, the Hamming weight, a Hamiltonian, and your initial parameters.

Example usage:

import numpy as np
from qibo import hamiltonians
from scipy.special import comb
from qiboml.models.optimizers import ExactGeodesicTransportCG

nqubits = 4
weight = 2
dim = int(comb(nqubits, weight))
np.random.seed(42)

# Initial angles in hyperspherical coordinates
theta_init = np.random.uniform(low=0, high=np.pi, size=dim - 1)

# Define the Hamiltonian
hamiltonian = hamiltonians.XXZ(nqubits=nqubits)

# Instantiate the optimizer
optimizer = ExactGeodesicTransportCG(
    nqubits=nqubits,
    weight=weight,
    hamiltonian=hamiltonian,
    angles=theta_init,
)

# Run optimization for 20 steps
final_loss, losses, final_params = optimizer(steps=20)

Available Arguments#

When constructing an ExactGeodesicTransportCG, you may specify:

  • nqubits: Number of qubits in the circuit.

  • weight: Hamming weight of the state.

  • hamiltonian: qibo.hamiltonians.Hamiltonian instance.

  • angles: Initial hyperspherical angles (numpy array).

  • backtrack_rate: Backtracking factor for line search (default: 0.5).

  • geometric_gradient: Whether to use geometric gradient or finite difference (default: False).

  • multiplicative_factor: Scaling of initial learning rate (default: 1.0).

  • c1, c2: Wolfe line search constants (defaults: 1e-3, 0.9).

  • backend: Optional qibo backend (default: global backend)

Implementation Details#

The optimizer constructs the state at each iteration using the qiboml.models.encodings.hamming_weight_encoder() circuit, computes the loss as the Hamiltonian expectation value, and updates the angles by:

  1. Computing the Riemannian (natural) gradient.

  2. Performing a conjugate-gradient step along the geodesic direction, with the step size determined by a Wolfe condition line search

At the end of the run, it returns:

  • final_loss: Loss at optimal parameters.

  • losses: List of losses per iteration.

  • final_params: Optimal angles.