How to simulate the local unitary cluster Jastrow (LUCJ) ansatz¶
In this guide, we show how to use ffsim to simulate the local unitary cluster Jastrow (LUCJ) ansatz. We’ll use it to calculate an approximation to the ground state energy of an ethene molecule.
First, let’s build the molecule.
[1]:
import pyscf
import pyscf.mcscf
import ffsim
# Build an ethene molecule
bond_distance = 1.339
a = 0.5 * bond_distance
b = a + 0.5626
c = 0.9289
mol = pyscf.gto.Mole()
mol.build(
atom=[
["C", (0, 0, a)],
["C", (0, 0, -a)],
["H", (0, c, b)],
["H", (0, -c, b)],
["H", (0, c, -b)],
["H", (0, -c, -b)],
],
basis="sto-6g",
symmetry="d2h",
)
# Define active space
active_space = range(mol.nelectron // 2 - 2, mol.nelectron // 2 + 2)
# Get molecular data and molecular Hamiltonian (one- and two-body tensors)
scf = pyscf.scf.RHF(mol).run()
mol_data = ffsim.MolecularData.from_scf(scf, active_space=active_space)
norb = mol_data.norb
nelec = mol_data.nelec
mol_hamiltonian = mol_data.hamiltonian
# Compute FCI energy
mol_data.run_fci()
print(f"norb = {norb}")
print(f"nelec = {nelec}")
converged SCF energy = -77.8266321248744
Parsing /tmp/tmpvc94rcvf
converged SCF energy = -77.8266321248744
CASCI E = -77.8742165643862 E(CI) = -4.02122442107773 S^2 = 0.0000000
norb = 4
nelec = (2, 2)
Overwritten attributes get_hcore get_ovlp of <class 'pyscf.scf.hf.RHF'>
/home/runner/work/ffsim/ffsim/.tox/docs/lib/python3.12/site-packages/pyscf/gto/mole.py:1294: UserWarning: Function mol.dumps drops attribute energy_nuc because it is not JSON-serializable
warnings.warn(msg)
/home/runner/work/ffsim/ffsim/.tox/docs/lib/python3.12/site-packages/pyscf/gto/mole.py:1294: UserWarning: Function mol.dumps drops attribute intor_symmetric because it is not JSON-serializable
warnings.warn(msg)
General UCJ ansatz¶
Since our molecule has a closed-shell Hartree-Fock state, we’ll use the spin-balanced variant of the UCJ ansatz, UCJOpSpinBalanced. We’ll initialize the ansatz from t2 amplitudes obtained from a CCSD calculations. We’ll first demonstrate the general UCJ ansatz, without adding the locality constraints of the LUCJ ansatz just yet.
The following code cell initializes the ansatz operator, applies it to the Hartree-Fock state, and computes the energy of the resulting state.
[2]:
import numpy as np
from pyscf import cc
# Get CCSD t2 amplitudes for initializing the ansatz
ccsd = cc.CCSD(
scf, frozen=[i for i in range(mol.nao_nr()) if i not in active_space]
).run()
# Construct UCJ operator
n_reps = 2
operator = ffsim.UCJOpSpinBalanced.from_t_amplitudes(ccsd.t2, n_reps=n_reps)
# Construct the Hartree-Fock state to use as the reference state
reference_state = ffsim.hartree_fock_state(norb, nelec)
# Apply the operator to the reference state
ansatz_state = ffsim.apply_unitary(reference_state, operator, norb=norb, nelec=nelec)
# Compute the energy ⟨ψ|H|ψ⟩ of the ansatz state
hamiltonian = ffsim.linear_operator(mol_hamiltonian, norb=norb, nelec=nelec)
energy = np.real(np.vdot(ansatz_state, hamiltonian @ ansatz_state))
print(f"Energy at initialization: {energy}")
E(CCSD) = -77.87421536374032 E_corr = -0.04758323886585756
Energy at initialization: -77.87160024816279
<class 'pyscf.cc.ccsd.CCSD'> does not have attributes converged
To variationally optimize the ansatz, we’ll take advantage of methods for conversion to and from real-valued parameter vectors. In the following code cell, we define an objective function that takes a parameter vector as input and outputs the energy of the associated ansatz state. We then optimize this objective function using scipy.optimize.minimize
, with an initial guess obtained from the operator we initialized previously from t2 amplitudes.
[3]:
import scipy.optimize
def fun(x):
# Initialize the ansatz operator from the parameter vector
operator = ffsim.UCJOpSpinBalanced.from_parameters(x, norb=norb, n_reps=n_reps)
# Apply the ansatz operator to the reference state
final_state = ffsim.apply_unitary(reference_state, operator, norb=norb, nelec=nelec)
# Return the energy ⟨ψ|H|ψ⟩ of the ansatz state
return np.real(np.vdot(final_state, hamiltonian @ final_state))
result = scipy.optimize.minimize(
fun, x0=operator.to_parameters(), method="L-BFGS-B", options=dict(maxiter=10)
)
print(f"Number of parameters: {len(result.x)}")
print(result)
Number of parameters: 72
message: STOP: TOTAL NO. of ITERATIONS REACHED LIMIT
success: False
status: 1
fun: -77.87387391547978
x: [-1.152e+00 -4.714e-04 ... 2.006e-04 1.287e-01]
nit: 10
jac: [ 3.553e-05 1.990e-05 ... -1.421e-06 7.105e-06]
nfev: 949
njev: 13
hess_inv: <72x72 LbfgsInvHessProduct with dtype=float64>
LUCJ ansatz¶
Now, let’s add locality constraints to simulate the LUCJ ansatz. We’ll restrict same-spin interactions to a line topology, and opposite-spin interactions to those within the same spatial orbital. As explained in The local unitary cluster Jastrow (LUCJ) ansatz, these constraints allow the ansatz to be simulated directly on a square lattice.
In the following code cell, we demonstrate the optimization of the ansatz with these constraints imposed. Notice that with the same number of ansatz repetitions, the number of parameters of the ansatz has decreased from 72 to 46.
[4]:
pairs_aa = [(p, p + 1) for p in range(norb - 1)]
pairs_ab = [(p, p) for p in range(norb)]
interaction_pairs = (pairs_aa, pairs_ab)
def fun(x):
operator = ffsim.UCJOpSpinBalanced.from_parameters(
x, norb=norb, n_reps=n_reps, interaction_pairs=interaction_pairs
)
final_state = ffsim.apply_unitary(reference_state, operator, norb=norb, nelec=nelec)
return np.real(np.vdot(final_state, hamiltonian @ final_state))
result = scipy.optimize.minimize(
fun,
x0=operator.to_parameters(interaction_pairs=interaction_pairs),
method="L-BFGS-B",
options=dict(maxiter=10),
)
print(f"Number of parameters: {len(result.x)}")
print(result)
Number of parameters: 46
message: CONVERGENCE: REL_REDUCTION_OF_F_<=_FACTR*EPSMCH
success: True
status: 0
fun: -77.87363426573044
x: [-1.152e+00 -2.894e-05 ... 3.519e-02 2.561e-01]
nit: 5
jac: [ 3.268e-05 1.137e-05 ... 2.842e-06 -5.684e-06]
nfev: 329
njev: 7
hess_inv: <46x46 LbfgsInvHessProduct with dtype=float64>
Optimize with the linear method¶
ffsim includes an implementation of the “linear method” for optimization of a variational wavefunction. The linear method often converges faster than a standard optimization algorithm like L-BFGS-B. The interface is similar to that of scipy.optimize.minimize
, the main difference being that instead of passing a callable that directly returns the function value to be optimized, you pass two objects: a callable that returns the wavefunction, and the
Hamiltonian representing the energy to be optimized as a LinearOperator
. The code cell below shows how to use the linear method to optimize the LUCJ ansatz from the previous example. It also shows how you can use an optional callback function to save intermediate results of the optimization.
[5]:
from collections import defaultdict
from ffsim.optimize import minimize_linear_method
# Define function that converts a list of parameters to the corresponding state vector
def params_to_vec(x: np.ndarray) -> np.ndarray:
operator = ffsim.UCJOpSpinBalanced.from_parameters(
x, norb=norb, n_reps=n_reps, interaction_pairs=interaction_pairs
)
return ffsim.apply_unitary(reference_state, operator, norb=norb, nelec=nelec)
# Define a callback function used to save optimization information (this is optional)
info = defaultdict(list)
def callback(intermediate_result: scipy.optimize.OptimizeResult):
# The callback function is called after each iteration. It accepts
# an OptimizeResult object storing the parameters and function value at
# the current iteration, and possibly other information
info["x"].append(intermediate_result.x)
info["fun"].append(intermediate_result.fun)
if hasattr(intermediate_result, "jac"):
info["jac"].append(intermediate_result.jac)
if hasattr(intermediate_result, "regularization"):
info["regularization"].append(intermediate_result.regularization)
if hasattr(intermediate_result, "variation"):
info["variation"].append(intermediate_result.variation)
# Optimize with the linear method
result = minimize_linear_method(
params_to_vec,
hamiltonian,
x0=operator.to_parameters(interaction_pairs=interaction_pairs),
maxiter=10,
callback=callback,
)
# Print some information
print(f"Number of parameters: {len(result.x)}")
print(result)
print()
for i, (fun, jac, regularization, variation) in enumerate(
zip(info["fun"], info["jac"], info["regularization"], info["variation"])
):
print(f"Iteration {i + 1}")
print(f" Energy: {fun}")
print(f" Norm of gradient: {np.linalg.norm(jac)}")
print(f" Regularization hyperparameter: {np.linalg.norm(regularization)}")
print(f" Variation hyperparameter: {np.linalg.norm(variation)}")
Number of parameters: 46
message: Convergence: Relative reduction of objective function <= ftol.
success: True
fun: -77.87363432754856
x: [-1.152e+00 -3.593e-04 ... 3.469e-02 2.561e-01]
nit: 3
jac: [-1.009e-07 -2.803e-07 ... -2.121e-07 -3.081e-07]
nfev: 523
njev: 4
nlinop: 339
Iteration 1
Energy: -77.87362955402136
Norm of gradient: 0.0017165491164075387
Regularization hyperparameter: 0.0014208280242812012
Variation hyperparameter: 0.9875543513636535
Iteration 2
Energy: -77.8736343016204
Norm of gradient: 7.317103416655098e-05
Regularization hyperparameter: 0.00362634756942963
Variation hyperparameter: 0.9878164041461586
Iteration 3
Energy: -77.87363432754856
Norm of gradient: 8.879352563944313e-06
Regularization hyperparameter: 0.003626667521704698
Variation hyperparameter: 0.9878164640925233