ffsim.UCJOpSpinUnbalanced

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

Bases: SupportsApplyUnitary, SupportsApproximateEquality

A spin-unbalanced 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}\]

Since \(\mathbf{J}^{(\sigma \tau)}_{ij} = \mathbf{J}^{(\tau \sigma)}_{ji}\), each diagonal Coulomb operator requires 3 matrices for its description: \(\mathbf{J}^{(\alpha \alpha)}\), \(\mathbf{J}^{(\alpha \beta)}\), and \(\mathbf{J}^{(\beta \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.

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.

Attributes

atol

final_orbital_rotation

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

n_reps

The number of ansatz repetitions.

norb

The number of spatial orbitals.

rtol

validate

diag_coulomb_mats

The diagonal Coulomb matrices, as a Numpy array of shape (n_reps, 3, norb, norb).

orbital_rotations

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

atol: InitVar = 1e-08
final_orbital_rotation: ndarray | None = None

The optional final orbital rotation, as a Numpy array of shape (2, norb, norb). This can be viewed as a list of two orbital rotations, the first one for spin alpha and the second one for spin beta.

n_reps

The number of ansatz repetitions.

norb

The number of spatial orbitals.

rtol: InitVar = 1e-05
validate: InitVar = True
diag_coulomb_mats: ndarray

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

orbital_rotations: ndarray

The orbital rotations, as a Numpy array of shape (n_reps, 2, norb, norb). The last two axes index the rows and columns of the orbital rotations, and the third from last axis, which has 2 dimensions, indexes the spin sector of the orbital rotation: first alpha, then beta. The first axis indexes the ansatz repetitions.

Methods

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_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, list[tuple[int, int]] | None] | None) – Optional restrictions on allowed orbital interactions for the diagonal Coulomb operators. If specified, interaction_pairs should be a tuple of 3 lists, for alpha-alpha, alpha-beta, and beta-beta interactions, in that order. Any 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. For the alpha-alpha and beta-beta interactions, 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:

UCJOpSpinUnbalanced

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 for alpha-alpha or beta-beta interactions 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 repetitions. Terms are included in decreasing order of the magnitude of the corresponding singular value in the factorization.

Parameters:
  • t2 (tuple[ndarray, ndarray, ndarray]) – The t2 amplitudes. This should be a tuple of 3 Numpy arrays, (t2aa, t2ab, t2bb), containing the alpha-alpha, alpha-beta, and beta-beta t2 amplitudes.

  • t1 (tuple[ndarray, ndarray] | None) – The t1 amplitudes. This should be a pair of Numpy arrays, (t1a, t1b), containing the alpha and beta t1 amplitudes.

  • n_reps (int | tuple[int, int] | None) – The number of ansatz repetitions. You can pass a single integer or a pair of integers. If a single integer, terms from the alpha-beta t2 amplitudes are used before including any terms from the alpha-alpha and beta-beta t2 amplitudes. If a pair of integers, then the first integer specifies the number of terms to use from the alpha-beta t2 amplitudes, and the second integer specifies the number of terms to use from the alpha-alpha and beta-beta t2 amplitudes. 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, list[tuple[int, int]] | None] | None) – Optional restrictions on allowed orbital interactions for the diagonal Coulomb operators. If specified, interaction_pairs should be a tuple of 3 lists, for alpha-alpha, alpha-beta, and beta-beta interactions, in that order. Any 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. For the alpha-alpha and beta-beta interactions, 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:

UCJOpSpinUnbalanced

Returns:

The UCJ operator with parameters initialized from the t2 amplitudes.

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

  • ValueError – Interaction pairs list for alpha-alpha or beta-beta interactions 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, list[tuple[int, int]] | None] | None) – Optional restrictions on allowed orbital interactions for the diagonal Coulomb operators. If specified, interaction_pairs should be a tuple of 3 lists, for alpha-alpha, alpha-beta, and beta-beta interactions, in that order. Any 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. For the alpha-alpha and beta-beta interactions, 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 for alpha-alpha or beta-beta interactions 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, list[tuple[int, int]] | None] | None) – Optional restrictions on allowed orbital interactions for the diagonal Coulomb operators. If specified, interaction_pairs should be a tuple of 3 lists, for alpha-alpha, alpha-beta, and beta-beta interactions, in that order. Any 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. For the alpha-alpha and beta-beta interactions, 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 for alpha-alpha or beta-beta interactions contained lower triangular pairs.