ffsim.DoubleFactorizedHamiltonian¶
- class ffsim.DoubleFactorizedHamiltonian(one_body_tensor, diag_coulomb_mats, orbital_rotations, constant=0.0, z_representation=False)[source]¶
Bases:
SupportsApproximateEquality,SupportsDiagonal,SupportsFermionOperator,SupportsLinearOperatorA Hamiltonian in the double-factorized representation.
The double-factorized form of the molecular Hamiltonian is
\[\begin{split}H = \sum_{\substack{pq \\ \sigma}} \kappa_{pq} a^\dagger_{p\sigma} a_{q\sigma} + \frac12 \sum_t \sum_{\substack{ij \\ \sigma\tau}} J^{(t)}_{ij} n^{(t)}_{i\sigma} n^{(t)}_{j\tau} + \text{constant}.\end{split}\]where
\[n^{(t)}_{i\sigma} = \sum_{pq} U^{(t)}_{pi} a^\dagger_{p\sigma} a_{q\sigma} U^{(t)*}_{qi}.\]Here each \(U^{(t)}\) is a unitary matrix and each \(J^{(t)}\) is a real symmetric matrix.
“Z” representation
The “Z” representation of the double factorization is an alternative representation that sometimes yields simpler quantum circuits.
Under the Jordan-Wigner transformation, the number operators take the form
\[n^{(t)}_{i\sigma} = \frac{(1 - z^{(t)}_{i\sigma})}{2}\]where \(z^{(t)}_{i\sigma}\) is the Pauli Z operator in the rotated basis. The “Z” representation is obtained by rewriting the two-body part in terms of these Pauli Z operators and updating the one-body term as appropriate:
\[\begin{split}H = \sum_{\substack{pq \\ \sigma}} \kappa'_{pq} a^\dagger_{p\sigma} a_{q\sigma} + \frac18 \sum_t \sum_{\substack{ij \\ \sigma\tau}}^* J^{(t)}_{ij} z^{(t)}_{i\sigma} z^{(t)}_{j\tau} + \text{constant}''\end{split}\]where the asterisk denotes summation over indices \(ij, \sigma\tau\) where \(\sigma \neq \tau\) or \(i \neq j\).
References
Attributes
The constant.
The number of spatial orbitals.
Whether the Hamiltonian is in the "Z" representation rather than the "number" representation.
The one-body tensor \(\kappa\).
The diagonal Coulomb matrices.
The orbital rotations.
- norb¶
The number of spatial orbitals.
-
z_representation:
bool= False¶ Whether the Hamiltonian is in the “Z” representation rather than the “number” representation.
Methods
from_molecular_hamiltonian(hamiltonian, *[, ...])Initialize a DoubleFactorizedHamiltonian from a MolecularHamiltonian.
Convert the DoubleFactorizedHamiltonian to a MolecularHamiltonian.
Return the Hamiltonian in the "number" representation.
Return the Hamiltonian in the "Z" representation.
- static from_molecular_hamiltonian(hamiltonian, *, z_representation=False, tol=1e-08, max_vecs=None, optimize=False, method='L-BFGS-B', callback=None, options=None, diag_coulomb_indices=None, cholesky=True)[source]¶
Initialize a DoubleFactorizedHamiltonian from a MolecularHamiltonian.
This function takes as input a
MolecularHamiltonian, which stores a one-body tensor, two-body tensor, and constant. It performs a double-factorized decomposition of the two-body tensor and computes a new one-body tensor and constant, and returns aDoubleFactorizedHamiltonianstoring the results.See
DoubleFactorizedHamiltonianfor a description of thez_representationargument. Seeffsim.linalg.double_factorized()for a description of the rest of the arguments.- Parameters:
hamiltonian (
MolecularHamiltonian) – The Hamiltonian whose double-factorized representation to compute.z_representation (
bool) – Whether to use the “Z” representation of the decomposition.tol (
float) – Tolerance for error in the decomposition. The error is defined as the maximum absolute difference between an element of the original tensor and the corresponding element of the reconstructed tensor.max_vecs (
int|None) – An optional limit on the number of terms to keep in the decomposition of the two-body tensor. This argument overridestol.optimize (
bool) – Whether to optimize the tensors returned by the decomposition.method (
str) – The optimization method. See the documentation of scipy.optimize.minimize for possible values.callback – Callback function for the optimization. See the documentation of scipy.optimize.minimize for usage.
options (
dict|None) – Options for the optimization. See the documentation of scipy.optimize.minimize for usage.diag_coulomb_indices (
list[tuple[int,int]] |None) – Allowed indices for nonzero values of the diagonal Coulomb matrices. Matrix entries corresponding to indices not in this list will be set to zero. This list should contain only upper triangular indices, i.e., pairs \((i, j)\) where \(i \leq j\). Passing a list with lower triangular indices will raise an error. This parameter is only used ifoptimizeis set to True.cholesky (
bool) – Whether to perform the factorization using a modified Cholesky decomposition. If False, a full eigenvalue decomposition is used instead, which can be much more expensive. This argument is ignored ifoptimizeis set to True.
- Return type:
- Returns:
The double-factorized Hamiltonian.
- to_molecular_hamiltonian()[source]¶
Convert the DoubleFactorizedHamiltonian to a MolecularHamiltonian.
- to_number_representation()[source]¶
Return the Hamiltonian in the “number” representation.
- Return type: