ffsim.UCJOpSpinBalanced

class ffsim.UCJOpSpinBalanced(diag_coulomb_mats, orbital_rotations, final_orbital_rotation=None, validate=True, rtol=1e-05, atol=1e-08)[source]

Bases: SupportsApplyUnitary, SupportsApproximateEquality

A spin-balanced unitary cluster Jastrow operator.

A unitary cluster Jastrow (UCJ) operator has the form

\[\prod_{k = 1}^L \mathcal{U}_k e^{i \mathcal{J}_k} \mathcal{U}_k^\dagger\]

where each \(\mathcal{U_k}\) is an orbital rotation and each \(\mathcal{J}\) is a diagonal Coulomb operator of the form

\[\begin{split}\mathcal{J} = \frac12\sum_{\substack{ij \\ \sigma \tau}} \mathbf{J}^{(\sigma \tau)}_{ij} n_{i\sigma} n_{j\tau}.\end{split}\]

For the spin-balanced operator, we require that \(\mathbf{J}^{(\alpha \alpha)} = \mathbf{J}^{(\beta \beta)}\) and \(\mathbf{J}^{(\alpha \beta)} = \mathbf{J}^{(\beta \alpha)}\). Therefore, each diagonal Coulomb operator is described by 2 matrices, \(\mathbf{J}^{(\alpha \alpha)}\) and \(\mathbf{J}^{(\alpha \beta)}\), and both of these matrices are symmetric. Furthermore, each orbital rotation is described by a single matrix because the same orbital rotation is applied to both spin alpha and spin beta. The number of terms \(L\) is referred to as the number of ansatz repetitions and is accessible via the n_reps attribute.

To support variational optimization of the orbital basis, an optional final orbital rotation can be included in the operator, to be performed at the end.

diag_coulomb_mats

The diagonal Coulomb matrices, as a Numpy array of shape (n_reps, 2, norb, norb) The last two axes index the rows and columns of the matrices, and the third from last axis, which has 2 dimensions, indexes the spin interaction type of the matrix: alpha-alpha, and then alpha-beta. The first axis indexes the ansatz repetitions.

Type:

np.ndarray

orbital_rotations

The orbital rotations, as a Numpy array of shape (n_reps, norb, norb).

Type:

np.ndarray

final_orbital_rotation

The optional final orbital rotation, as a Numpy array of shape (norb, norb).

Type:

np.ndarray | None

Parameters:
  • validate (InitVar) – Whether to validate the operator attributes. Setting this to False skips validation, which is useful if you need to create many instances of this class and are confident that the attributes are valid.

  • rtol (InitVar) – Relative numerical tolerance for validation checks.

  • atol (InitVar) – Absolute numerical tolerance for validation checks.

Methods

from_cisd_vec(cisd_vec, *, norb, nocc[, ...])

Initialize the UCJ operator from a CISD vector.

from_parameters(params, *, norb, n_reps[, ...])

Initialize the UCJ operator from a real-valued parameter vector.

from_t_amplitudes(t2, *[, t1, n_reps, ...])

Initialize the UCJ operator from t2 (and optionally t1) amplitudes.

n_params(norb, n_reps, *[, ...])

Return the number of parameters of an ansatz with given settings.

to_parameters(*[, interaction_pairs])

Convert the UCJ operator to a real-valued parameter vector.

static from_cisd_vec(cisd_vec, *, norb, nocc, c0_threshold=1e-12, n_reps=None, interaction_pairs=None, tol=1e-08, optimize=False, method='L-BFGS-B', callback=None, options=None, regularization=0, multi_stage_start=None, multi_stage_step=None)[source]

Initialize the UCJ operator from a CISD vector.

The CISD amplitudes are converted to CC-style amplitudes using

\[t_{1,ia} = \frac{c^{(1)}_{ia}}{c_0}, \qquad t_{2,ijab} = \frac{c^{(2)}_{ijab}}{c_0} - \frac{t_{1,ia} t_{1,jb}}{2},\]

and the resulting amplitudes are then passed to from_t_amplitudes().

Parameters:
  • cisd_vec (ndarray) – The CISD vector. This is a one-dimensional array containing the reference coefficient \(c_0\) in the first entry, followed by the singles and then doubles coefficients.

  • norb (int) – The number of spatial orbitals.

  • nocc (int) – The number of occupied orbitals.

  • c0_threshold (float) – Absolute value threshold for the reference coefficient. An error is raised if the absolute value of \(c_0\) is smaller than this threshold.

  • n_reps (int | None) – The number of ansatz repetitions. If not specified, then it is set to the number of terms resulting from the double factorization of the t2 amplitudes. If the specified number of repetitions is larger than the number of terms resulting from the double factorization, then the ansatz is padded with additional identity operators up to the specified number of repetitions.

  • interaction_pairs (tuple[list[tuple[int, int]] | None, list[tuple[int, int]] | None] | None) – Optional restrictions on allowed orbital interactions for the diagonal Coulomb operators. If specified, interaction_pairs should be a pair of lists, for alpha-alpha and alpha-beta interactions, in that order. Either list can be substituted with None to indicate no restrictions on interactions. Each list should contain pairs of integers representing the orbitals that are allowed to interact. These pairs can also be interpreted as indices of diagonal Coulomb matrix entries that are allowed to be nonzero. Each integer pair must be upper triangular, that is, of the form \((i, j)\) where \(i \leq j\).

  • tol (float) – Tolerance for error in the double-factorized decomposition of the converted t2 amplitudes. The error is defined as the maximum absolute difference between an element of the original tensor and the corresponding element of the reconstructed tensor.

  • optimize (bool) – Whether to optimize the tensors returned by the decomposition to to minimize the error in the factorization.

  • method (str) – The optimization method. See the documentation of scipy.optimize.minimize for possible values. This argument is ignored if optimize is set to False.

  • callback – Callback function for the optimization. See the documentation of scipy.optimize.minimize for usage. This argument is ignored if optimize is set to False.

  • options (dict | None) – Options for the optimization. See the documentation of scipy.optimize.minimize for usage. This argument is ignored if optimize is set to False.

  • regularization (float) – See ffsim.linalg.double_factorized_t2() for a description of this argument. This argument is ignored if optimize is set to False.

  • multi_stage_start (int | None) – See ffsim.linalg.double_factorized_t2() for a description of this argument. This argument is ignored if optimize is set to False.

  • multi_stage_step (int | None) – See ffsim.linalg.double_factorized_t2() for a description of this argument. This argument is ignored if optimize is set to False.

Return type:

UCJOpSpinBalanced

Returns:

The UCJ operator with parameters initialized from the CISD amplitudes.

Raises:

ValueError – The CISD reference coefficient is smaller than the specified threshold c0_tol.

static from_parameters(params, *, norb, n_reps, interaction_pairs=None, with_final_orbital_rotation=False)[source]

Initialize the UCJ operator from a real-valued parameter vector.

Parameters:
  • params (ndarray) – The real-valued parameter vector.

  • norb (int) – The number of spatial orbitals.

  • n_reps (int) – The number of ansatz repetitions.

  • interaction_pairs (tuple[list[tuple[int, int]] | None, list[tuple[int, int]] | None] | None) – Optional restrictions on allowed orbital interactions for the diagonal Coulomb operators. If specified, interaction_pairs should be a pair of lists, for alpha-alpha and alpha-beta interactions, in that order. Either list can be substituted with None to indicate no restrictions on interactions. Each list should contain pairs of integers representing the orbitals that are allowed to interact. These pairs can also be interpreted as indices of diagonal Coulomb matrix entries that are allowed to be nonzero. Each integer pair must be upper triangular, that is, of the form \((i, j)\) where \(i \leq j\).

  • with_final_orbital_rotation (bool) – Whether to include a final orbital rotation in the operator.

Return type:

UCJOpSpinBalanced

Returns:

The UCJ operator constructed from the given parameters.

Raises:
  • ValueError – The number of parameters passed did not match the number expected based on the function inputs.

  • ValueError – Interaction pairs list contained duplicate interactions.

  • ValueError – Interaction pairs list contained lower triangular pairs.

static from_t_amplitudes(t2, *, t1=None, n_reps=None, interaction_pairs=None, tol=1e-08, optimize=False, method='L-BFGS-B', callback=None, options=None, regularization=0, multi_stage_start=None, multi_stage_step=None)[source]

Initialize the UCJ operator from t2 (and optionally t1) amplitudes.

Performs a double factorization of the t2 amplitudes and constructs the ansatz repetitions from the terms of the decomposition, up to an optionally specified number of ansatz repetitions.

The default behavior of this routine is to perform a straightforward “exact” factorization of the t2 amplitudes tensor based on a nested eigenvalue decomposition, and then truncate the terms based on the values of tol and n_reps. If optimize is set to True, then the entries of the resulting tensors (the diagonal Coulomb matrices and orbital rotations) are further optimized with scipy.optimize.minimize to reduce the error in the factorization. See ffsim.linalg.double_factorized_t2() for details.

Parameters:
  • t2 (ndarray) – The t2 amplitudes.

  • t1 (ndarray | None) – The t1 amplitudes.

  • n_reps (int | None) – The number of ansatz repetitions. If not specified, then it is set to the number of terms resulting from the double factorization of the t2 amplitudes. If the specified number of repetitions is larger than the number of terms resulting from the double factorization, then the ansatz is padded with additional identity operators up to the specified number of repetitions.

  • interaction_pairs (tuple[list[tuple[int, int]] | None, list[tuple[int, int]] | None] | None) – Optional restrictions on allowed orbital interactions for the diagonal Coulomb operators. If specified, interaction_pairs should be a pair of lists, for alpha-alpha and alpha-beta interactions, in that order. Either list can be substituted with None to indicate no restrictions on interactions. Each list should contain pairs of integers representing the orbitals that are allowed to interact. These pairs can also be interpreted as indices of diagonal Coulomb matrix entries that are allowed to be nonzero. Each integer pair must be upper triangular, that is, of the form \((i, j)\) where \(i \leq j\).

  • tol (float) – Tolerance for error in the double-factorized decomposition of the t2 amplitudes. The error is defined as the maximum absolute difference between an element of the original tensor and the corresponding element of the reconstructed tensor.

  • optimize (bool) – Whether to optimize the tensors returned by the decomposition to to minimize the error in the factorization.

  • method (str) – The optimization method. See the documentation of scipy.optimize.minimize for possible values. This argument is ignored if optimize is set to False.

  • callback – Callback function for the optimization. See the documentation of scipy.optimize.minimize for usage. This argument is ignored if optimize is set to False.

  • options (dict | None) – Options for the optimization. See the documentation of scipy.optimize.minimize for usage. This argument is ignored if optimize is set to False.

  • regularization (float) – See ffsim.linalg.double_factorized_t2() for a description of this argument. This argument is ignored if optimize is set to False.

  • multi_stage_start (int | None) – See ffsim.linalg.double_factorized_t2() for a description of this argument. This argument is ignored if optimize is set to False.

  • multi_stage_step (int | None) – See ffsim.linalg.double_factorized_t2() for a description of this argument. This argument is ignored if optimize is set to False.

Return type:

UCJOpSpinBalanced

Returns:

The UCJ operator with parameters initialized from the t2 amplitudes.

Raises:
  • ValueError – Interaction pairs list contained duplicate interactions.

  • ValueError – Interaction pairs list contained lower triangular pairs.

static n_params(norb, n_reps, *, interaction_pairs=None, with_final_orbital_rotation=False)[source]

Return the number of parameters of an ansatz with given settings.

Parameters:
  • norb (int) – The number of spatial orbitals.

  • n_reps (int) – The number of ansatz repetitions.

  • interaction_pairs (tuple[list[tuple[int, int]] | None, list[tuple[int, int]] | None] | None) – Optional restrictions on allowed orbital interactions for the diagonal Coulomb operators. If specified, interaction_pairs should be a pair of lists, for alpha-alpha and alpha-beta interactions, in that order. Either list can be substituted with None to indicate no restrictions on interactions. Each list should contain pairs of integers representing the orbitals that are allowed to interact. These pairs can also be interpreted as indices of diagonal Coulomb matrix entries that are allowed to be nonzero. Each integer pair must be upper triangular, that is, of the form \((i, j)\) where \(i \leq j\).

  • with_final_orbital_rotation (bool) – Whether to include a final orbital rotation in the operator.

Return type:

int

Returns:

The number of parameters of the ansatz.

Raises:
  • ValueError – Interaction pairs list contained duplicate interactions.

  • ValueError – Interaction pairs list contained lower triangular pairs.

to_parameters(*, interaction_pairs=None)[source]

Convert the UCJ operator to a real-valued parameter vector.

Note

If interaction_pairs is specified, the returned parameter vector will incorporate only the diagonal Coulomb matrix entries corresponding to the specified interactions, so the original operator will not be recoverable from the parameter vector.

Parameters:

interaction_pairs (tuple[list[tuple[int, int]] | None, list[tuple[int, int]] | None] | None) – Optional restrictions on allowed orbital interactions for the diagonal Coulomb operators. If specified, interaction_pairs should be a pair of lists, for alpha-alpha and alpha-beta interactions, in that order. Either list can be substituted with None to indicate no restrictions on interactions. Each list should contain pairs of integers representing the orbitals that are allowed to interact. These pairs can also be interpreted as indices of diagonal Coulomb matrix entries that are allowed to be nonzero. Each integer pair must be upper triangular, that is, of the form \((i, j)\) where \(i \leq j\).

Return type:

ndarray

Returns:

The real-valued parameter vector.

Raises:
  • ValueError – Interaction pairs list contained duplicate interactions.

  • ValueError – Interaction pairs list contained lower triangular pairs.

Attributes

atol

final_orbital_rotation

n_reps

The number of ansatz repetitions.

norb

The number of spatial orbitals.

rtol

validate

diag_coulomb_mats

orbital_rotations