Fermion (qiskit_addon_sqd.fermion)

Functions for the study of fermionic systems.

class SCIState(amplitudes, ci_strs_a, ci_strs_b, norb, nelec)[source]

Bases: object

The amplitudes and determinants describing a quantum state.

Parameters:
amplitudes: ndarray

An \(M \times N\) array where \(M =\) len(ci_strs_a) and \(N\) = len(ci_strs_b). amplitudes[i][j] is the amplitude of the determinant pair (ci_strs_a[i], ci_strs_b[j]).

ci_strs_a: ndarray

The alpha determinants.

ci_strs_b: ndarray

The beta determinants.

classmethod load(filename)[source]

Load an SCIState object from an .npz file.

nelec: tuple[int, int]

The numbers of alpha and beta electrons.

norb: int

The number of spatial orbitals.

orbital_occupancies()[source]

Average orbital occupancies.

Return type:

tuple[ndarray, ndarray]

rdm(rank=1, spin_summed=False)[source]

Compute reduced density matrix.

Parameters:
Return type:

ndarray

save(filename)[source]

Save the SCIState object to an .npz file.

bitstring_matrix_to_ci_strs(bitstring_matrix, open_shell=False)[source]

Convert bitstrings (rows) in a bitstring_matrix into integer representations of determinants.

This function separates each bitstring in bitstring_matrix in half, flips the bits and translates them into integer representations, and finally appends them to their respective (spin-up or spin-down) lists. Those lists are sorted and output from this function.

Parameters:
  • bitstring_matrix (ndarray) – A 2D array of bool representations of bit values such that each row represents a single bitstring

  • open_shell (bool) – A flag specifying whether unique configurations from the left and right halves of the bitstrings should be kept separate. If False, configurations from the left and right halves of the bitstrings are combined into a single set of unique configurations. That combined set will be returned for both the left and right bitstrings.

Returns:

A length-2 tuple of determinant lists representing the right (spin-up) and left (spin-down) halves of the bitstrings, respectively.

Return type:

tuple[ndarray, ndarray]

diagonalize_fermionic_hamiltonian(one_body_tensor, two_body_tensor, bit_array, samples_per_batch, norb, nelec, *, num_batches=1, energy_tol=1e-08, occupancies_tol=1e-05, max_iterations=100, sci_solver=None, symmetrize_spin=False, include_configurations=None, initial_occupancies=None, carryover_threshold=0.0001, callback=None, seed=None)[source]

Run the sample-based quantum diagonalization (SQD) algorithm.

Parameters:
  • one_body_tensor (ndarray) – The one-body tensor of the Hamiltonian.

  • two_body_tensor (ndarray) – The two-body tensor of the Hamiltonian.

  • bit_array (BitArray) – Array of sampled bitstrings. Each bitstring should have both the alpha part and beta part concatenated together, with the alpha part concatenated on the right-hand side, like this: [b_N, ..., b_0, a_N, ..., a_0].

  • samples_per_batch (int) – The number of bitstrings to include in each subsampled batch of bitstrings.

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

  • nelec (tuple[int, int]) – The numbers of alpha and beta electrons.

  • num_batches (int) – The number of batches to subsample in each configuration recovery iteration. This argument indirectly controls the dimensions of the diagonalization subspaces. A higher value will yield larger subspace dimensions.

  • energy_tol (float) – Numerical tolerance for convergence of the energy. If the change in energy between iterations is smaller than this value, then the configuration recovery loop will exit, if the occupancies have also converged (see the occupancies_tol argument).

  • occupancies_tol (float) – Numerical tolerance for convergence of the average orbital occupancies. If the maximum change in absolute value of the average occupancy of an orbital between iterations is smaller than this value, then the configuration recovery loop will exit, if the energy has also converged (see the energy_tol argument).

  • max_iterations (int) – Limit on the number of configuration recovery iterations.

  • sci_solver (Callable[[list[tuple[ndarray, ndarray]], ndarray, ndarray, int, tuple[int, int]], list[SCIResult]] | None) –

    Selected configuration interaction solver function.

    Inputs:

    • List of pairs (strings_a, strings_b) of arrays of spin-alpha CI strings and spin-beta CI strings whose Cartesian product give the basis of the subspace in which to perform a diagonalization. A list is passed to allow the solver function to perform the diagonalizations in parallel.

    • One-body tensor of the Hamiltonian.

    • Two-body tensor of the Hamiltonian.

    • The number of spatial orbitals.

    • A pair (n_alpha, n_beta) indicating the numbers of alpha and beta electrons.

    Output: List of (energy, sci_state, occupancies) triplets, where each triplet contains the result of the corresponding diagonalization.

  • symmetrize_spin (bool) – Whether to always merge spin-alpha and spin-beta CI strings into a single list, so that the diagonalization subspace is invariant with respect to the exchange of spin alpha with spin beta.

  • include_configurations (list[int] | tuple[list[int], list[int]] | None) – Configurations to always include in the diagonalization subspace. You can specify either a single list of single-spin strings to use for both spin sectors, or a pair (alpha_strings, beta_strings) of lists of single-spin strings, one for each spin.

  • initial_occupancies (tuple[ndarray, ndarray] | None) – Initial guess for the average occupancies of the orbitals.

  • carryover_threshold (float) – Threshold for carrying over bitstrings with large CI weight from one iteration of configuration recovery to the next. All single-spin CI strings associated with configurations whose coefficient has absolute value greater than this threshold will be included in the diagonalization subspace for the next iteration. A smaller threshold will retain more configurations, leading to a larger subspace and hence a more costly diagonalization.

  • callback (Callable[[list[SCIResult]], None] | None) – A callback function to be called after each configuration recovery iteration. The function will be passed the output of the sci_solver function, which is a list of (energy, sci_state, occupancies) triplets, where each triplet contains the result of a diagonalization.

  • seed (Generator | int | None) – A seed for the pseudorandom number generator.

Returns:

The estimate of the energy and the SCI state with that energy.

Return type:

SCIResult

enlarge_batch_from_transitions(bitstring_matrix, transition_operators)[source]

Apply the set of transition operators to the configurations represented in bitstring_matrix.

Parameters:
  • bitstring_matrix (ndarray) – A 2D array of bool representations of bit values such that each row represents a single bitstring.

  • transition_operators (ndarray) – A 1D or 2D array I, +, -, and n strings representing the action of the identity, creation, annihilation, or number operators. Each row represents a transition operator.

Returns:

Bitstring matrix representing the augmented set of electronic configurations after applying the excitation operators.

Return type:

ndarray

solve_fermion(bitstring_matrix, /, hcore, eri, *, open_shell=False, spin_sq=None, shift=0.1, **kwargs)[source]

Approximate the ground state given molecular integrals and a set of electronic configurations.

Parameters:
  • bitstring_matrix (tuple[ndarray, ndarray] | ndarray) –

    A set of configurations defining the subspace onto which the Hamiltonian will be projected and diagonalized.

    This may be specified in two ways:

    • Bitstring matrix: A 2D numpy.ndarray of bool values, where each row represents a bitstring. The spin-up configurations should occupy column indices (N, N/2], and the spin-down configurations should occupy column indices (N/2, 0], where N is the number of qubits.

    • CI strings: A tuple of two sequences containing integer representations of spin-up and spin-down determinants, respectively. The expected format is ([a_str_0, ..., a_str_N], [b_str_0, ..., b_str_M]).

  • hcore (ndarray) – Core Hamiltonian matrix representing single-electron integrals

  • eri (ndarray) – Electronic repulsion integrals representing two-electron integrals

  • open_shell (bool) – A flag specifying whether configurations from the left and right halves of the bitstrings should be kept separate. If False, CI strings from the left and right halves of the bitstrings are combined into a single set of unique configurations and used for both the alpha and beta subspaces.

  • spin_sq (float | None) – Target value for the total spin squared for the ground state, \(S^2 = s(s + 1)\). If None, no spin will be imposed.

  • shift (float) – Level shift for states which have different spin. \((H + shift * S^2)|ψ> = E|ψ>\)

  • **kwargs – Keyword arguments to pass to pyscf.fci.selected_ci.kernel_fixed_space

Returns:

  • Minimum energy from SCI calculation

  • The SCI ground state

  • Tuple containing orbital occupancies for spin-up and spin-down orbitals. Formatted as: (array([occ_a_0, ..., occ_a_N]), array([occ_b_0, ..., occ_b_N]))

  • Expectation value of spin-squared

Return type:

tuple[float, SCIState, tuple[ndarray, ndarray], float]

solve_sci(ci_strings, one_body_tensor, two_body_tensor, norb, nelec, *, spin_sq=None, **kwargs)[source]

Diagonalize Hamiltonian in subspace defined by CI strings.

Parameters:
  • ci_strings (tuple[ndarray, ndarray]) – Pair (strings_a, strings_b) of arrays of spin-alpha CI strings and spin-beta CI strings whose Cartesian product give the basis of the subspace in which to perform a diagonalization.

  • one_body_tensor (ndarray) – The one-body tensor of the Hamiltonian.

  • two_body_tensor (ndarray) – The two-body tensor of the Hamiltonian.

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

  • nelec (tuple[int, int]) – The numbers of alpha and beta electrons.

  • spin_sq (float | None) – Target value for the total spin squared for the ground state. If None, no spin will be imposed.

  • **kwargs

    Keyword arguments to pass to pyscf.fci.selected_ci.kernel_fixed_space

Returns:

The diagonalization result.

Return type:

SCIResult

solve_sci_batch(ci_strings, one_body_tensor, two_body_tensor, norb, nelec, *, spin_sq=None, **kwargs)[source]

Diagonalize Hamiltonian in subspaces.

Parameters:
  • ci_strings (list[tuple[ndarray, ndarray]]) – List of pairs (strings_a, strings_b) of arrays of spin-alpha CI strings and spin-beta CI strings whose Cartesian product give the basis of the subspace in which to perform a diagonalization.

  • one_body_tensor (ndarray) – The one-body tensor of the Hamiltonian.

  • two_body_tensor (ndarray) – The two-body tensor of the Hamiltonian.

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

  • nelec (tuple[int, int]) – The numbers of alpha and beta electrons.

  • spin_sq (float | None) – Target value for the total spin squared for the ground state. If None, no spin will be imposed.

  • **kwargs

    Keyword arguments to pass to pyscf.fci.selected_ci.kernel_fixed_space

Returns:

The results of the diagonalizations in the subspaces given by ci_strings.

Return type:

list[SCIResult]

optimize_orbitals(bitstring_matrix, /, hcore, eri, k_flat, *, open_shell=False, spin_sq=0.0, num_iters=10, num_steps_grad=10000, learning_rate=0.01, **kwargs)[source]

Optimize orbitals to produce a minimal ground state.

The process involves iterating over 3 steps:

For num_iters iterations:
  • Rotate the integrals with respect to the parameters, k_flat

  • Diagonalize and approximate the groundstate energy and wavefunction amplitudes

  • Optimize k_flat using gradient descent and the wavefunction amplitudes found in Step 2

Refer to Sec. II A 4 for more detailed discussion on this orbital optimization technique.

Parameters:
  • bitstring_matrix (tuple[ndarray, ndarray] | ndarray) –

    A set of configurations defining the subspace onto which the Hamiltonian will be projected and diagonalized.

    This may be specified in two ways:

    • Bitstring matrix: A 2D numpy.ndarray of bool values, where each row represents a bitstring. The spin-up configurations should occupy column indices (N, N/2], and the spin-down configurations should occupy column indices (N/2, 0], where N is the number of qubits.

    • CI strings: A tuple of two sequences containing integer representations of spin-up and spin-down determinants, respectively. The expected format is ([a_str_0, ..., a_str_N], [b_str_0, ..., b_str_M]).

  • hcore (ndarray) – Core Hamiltonian matrix representing single-electron integrals

  • eri (ndarray) – Electronic repulsion integrals representing two-electron integrals

  • k_flat (ndarray) – 1D array defining the orbital transform, K. The array should specify the upper triangle of the anti-symmetric transform operator in row-major order, excluding the diagonal.

  • open_shell (bool) – A flag specifying whether configurations from the left and right halves of the bitstrings should be kept separate. If False, CI strings from the left and right halves of the bitstrings are combined into a single set of unique configurations and used for both the alpha and beta subspaces.

  • spin_sq (float) – Target value for the total spin squared for the ground state

  • num_iters (int) – The number of iterations of orbital optimization to perform

  • num_steps_grad (int) – The number of steps of gradient descent to perform during each optimization iteration

  • learning_rate (float) – The learning rate to use during gradient descent

  • **kwargs

    Keyword arguments to pass to pyscf.fci.selected_ci.kernel_fixed_space

Returns:

  • The groundstate energy found during the last optimization iteration

  • An optimized 1D array defining the orbital transform

  • Tuple containing orbital occupancies for spin-up and spin-down orbitals. Formatted as: (array([occ_a_0, ..., occ_a_N]), array([occ_b_0, ..., occ_b_N]))

Return type:

tuple[float, ndarray, tuple[ndarray, ndarray]]

rotate_integrals(hcore, eri, k_flat)[source]

Perform a similarity transform on the integrals.

The transformation is described as:

\[\hat{\widetilde{H}} = \hat{U^{\dagger}}(k)\hat{H}\hat{U}(k)\]

For more information on how \(\hat{U}\) and \(\hat{U^{\dagger}}\) are generated from k_flat and applied to the one- and two-body integrals, refer to Sec. II A 4.

Parameters:
  • hcore (ndarray) – Core Hamiltonian matrix representing single-electron integrals

  • eri (ndarray) – Electronic repulsion integrals representing two-electron integrals

  • k_flat (ndarray) – 1D array defining the orbital transform, K. The array should specify the upper triangle of the anti-symmetric transform operator in row-major order, excluding the diagonal.

Returns:

  • The rotated core Hamiltonian matrix

  • The rotated ERI matrix

Return type:

tuple[ndarray, ndarray]