qml.qchem.factorize¶
- factorize(two_electron, tol_factor=1e-05, tol_eigval=1e-05, cholesky=False, compressed=False, regularization=None, **compression_kwargs)[source]¶
Return the double-factorized form of a two-electron integral tensor in spatial basis.
The two-electron tensor \(V\), in the chemist notation, can be decomposed in terms of orthonormal matrices \(U\) (leaf tensors) and symmetric matrices \(Z\) (core tensors) such that \(V_{ijkl} \approx \sum_r^R \sum_{pq} U_{ip}^{(r)} U_{jp}^{(r)} Z_{pq}^{(r)} U_{kq}^{(r)} U_{lq}^{(r)}\), where the rank \(R\) is determined by a threshold error.
For explicit double factorization, i.e., when
compressed=False
, the above decomposition is done using an eigenvalue or Cholesky decomposition to obtain symmetric matrices \(L^{(r)}\) such that \(V_{ijkl} = \sum_r^R L_{ij}^{(r)} L_{kl}^{(r) T}\), where core and leaf tensors are obtained by further diagonalizing each matrix \(L^{(r)}\) and truncating its eigenvalues (and the corresponding eigenvectors) at a threshold error. See theory section for more details.For compressed double factorization (CDF), i.e., when
compressed=True
, the above decomposition is done by optimizing the following cost function \(\mathcal{L}\) in a greedy layered-wise manner:\[\mathcal{L}(U, Z) = \frac{1}{2} \bigg|V_{ijkl} - \sum_r^R \sum_{pq} U_{ip}^{(r)} U_{jp}^{(r)} Z_{pq}^{(r)} U_{kq}^{(r)} U_{lq}^{(r)}\bigg|_{\text{F}} + \rho \sum_r^R \sum_{pq} \bigg|Z_{pq}^{(r)}\bigg|^{\gamma},\]where leaf tensors \(U\) are defined by the antisymmetric orbital rotations \(X\) such that \(U^{(r)} = \exp{(X^{(r)})}\), \(|\cdot|_{\text{F}}\) computes the Frobenius norm, \(\rho\) is a constant scaling factor, and \(|\cdot|^\gamma\) specifies the optional L1 and L2 regularization. See references arXiv:2104.08957 and arxiv:2212.07957 for more details.
Note
Packages JAX and Optax are required when performing CDF with
compressed=True
. Install them usingpip install jax optax
.- Parameters
two_electron (array[array[float]]) – Two-electron integral tensor in the molecular orbital basis arranged in chemist notation.
tol_factor (float) – Threshold error value for discarding the negligible factors. This will be used only when
compressed=False
.tol_eigval (float) – Threshold error value for discarding the negligible factor eigenvalues. This will be used only when
compressed=False
.cholesky (bool) – Use Cholesky decomposition for obtaining the symmetric matrices \(L^{(r)}\) instead of eigendecomposition. Default is
False
.compressed (bool) – Use compressed double factorization to optimize the factors returned in the decomposition. Look at the keyword arguments (
compression_kwargs
) for the available options which must be provided only whencompressed=True
.regularization (string | None) – Type of regularization (
"L1"
or"L2"
) to be used for optimizing the factors. Default is to not include any regularization term.
- Keyword Arguments
num_factors (int) – Maximum number of factors that should be optimized for compressed double factorization. Default is \(2\times N\), where N is the number of dimensions of two-electron tensor.
num_steps (int) – Maximum number of epochs for optimizing each factor. Default is
1000
.optimizer (optax.optimizer) – An optax optimizer instance. If not provided, Adam is used with
0.001
learning rate.init_params (dict[str, TensorLike] | None) – Intial values of the orbital rotations (\(X\)) and core tensors (\(Z\)) of shape
(num_factors, N, N)
given as a dictionary with keys"X"
and"Z"
, where N is the number of dimension of two-electron tensor. If not given, zero matrices will be used ifcholesky=False
and the core and leaf tensors corresponding to the firstnum_factors
will be used ifcholesky=True
.norm_prefactor (float) – Prefactor for scaling the regularization term. Default is
1e-5
.
- Returns
Tuple containing symmetric matrices (factors) approximating the two-electron integral tensor and core tensors and leaf tensors of the generated factors. In the explicit case where the core and leaf tensors could be truncated, they will be returned as a list.
- Return type
tuple(TensorLike, TensorLike, TensorLike)
- Raises
ValueError – If the specified regularization type is not supported.
ImportError – If the specified packages are not installed when
compressed=True
.
Example
>>> symbols = ['H', 'H'] >>> geometry = np.array([[0.0, 0.0, 0.0], ... [1.398397361, 0.0, 0.0]], requires_grad=False) >>> mol = qml.qchem.Molecule(symbols, geometry) >>> core, one, two = qml.qchem.electron_integrals(mol)() >>> two = np.swapaxes(two, 1, 3) # convert to chemist notation >>> factors, cores, leaves = qml.qchem.factorize(two, 1e-5, 1e-5) >>> print(factors) [[[-1.06723440e-01 6.42958741e-15] [ 7.71977824e-15 1.04898533e-01]] [[ 1.71099288e-13 -4.25688222e-01] [-4.25688222e-01 2.31561666e-13]] [[-8.14472856e-01 -3.89054708e-13] [-3.88994463e-13 -8.28642140e-01]]]
Theory
The second quantized electronic Hamiltonian is constructed in terms of fermionic creation, \(a^{\dagger}\) , and annihilation, \(a\), operators as [arXiv:1902.02134]
\[H = \sum_{\alpha \in \{\uparrow, \downarrow \} } \sum_{pq} h_{pq} a_{p,\alpha}^{\dagger} a_{q, \alpha} + \frac{1}{2} \sum_{\alpha, \beta \in \{\uparrow, \downarrow \} } \sum_{pqrs} h_{pqrs} a_{p, \alpha}^{\dagger} a_{q, \beta}^{\dagger} a_{r, \beta} a_{s, \alpha},\]where \(h_{pq}\) and \(h_{pqrs}\) are the one- and two-electron integrals computed as
\[h_{pq} = \int \phi_p(r)^* \left ( -\frac{\nabla_r^2}{2} - \sum_i \frac{Z_i}{|r-R_i|} \right) \phi_q(r) dr,\]and
\[h_{pqrs} = \int \frac{\phi_p(r_1)^* \phi_q(r_2)^* \phi_r(r_2) \phi_s(r_1)}{|r_1 - r_2|} dr_1 dr_2.\]The two-electron integrals can be rearranged in the so-called chemist notation which gives
\[V_{pqrs} = \int \frac{\phi_p(r_1)^* \phi_q(r_1)^* \phi_r(r_2) \phi_s(r_2)}{|r_1 - r_2|} dr_1 dr_2,\]and the molecular Hamiltonian can be rewritten as
\[H = \sum_{\alpha \in \{\uparrow, \downarrow \} } \sum_{pq} T_{pq} a_{p,\alpha}^{\dagger} a_{q, \alpha} + \frac{1}{2} \sum_{\alpha, \beta \in \{\uparrow, \downarrow \} } \sum_{pqrs} V_{pqrs} a_{p, \alpha}^{\dagger} a_{q, \alpha} a_{r, \beta}^{\dagger} a_{s, \beta},\]with
\[T_{pq} = h_{pq} - \frac{1}{2} \sum_s h_{pssq}.\]This notation allows a low-rank factorization of the two-electron integral. The objective of the factorization is to find a set of symmetric matrices, \(L^{(r)}\), such that
\[V_{ijkl} = \sum_r^R L_{ij}^{(r)} L_{kl}^{(r) T},\]with the rank \(R \leq n^2\) where \(n\) is the number of molecular orbitals. The matrices \(L^{(r)}\) are diagonalized and for each matrix the eigenvalues that are smaller than a given threshold (and their corresponding eigenvectors) are discarded. These can be used to further decompose \(V_{ijkl}\) in terms of orthonormal matrices \(U\) (leaf tensors) and symmetric matrices \(Z\) (core tensors), such that
\[V_{ijkl} = \sum_r^R \sum_{pq} U_{ip}^{(r)} U_{jp}^{(r)} Z_{pq}^{(r)} U_{kq}^{(r)} U_{lq}^{(r)},\]where \(U^{(r)}\) are the eigenvectors of \(L^{(r)}\) and \(Z^{(r)}\) are the outer product of the eigenvalues of \(L^{(r)}\).
The factorization algorithm has the following steps [arXiv:1902.02134]:
Reshape the \(n \times n \times n \times n\) two-electron tensor to a \(n^2 \times n^2\) matrix where \(n\) is the number of orbitals.
Decompose the resulting matrix either via eigendecomposition or Cholesky decomposition.
For the eigendecomposition, keep the \(r\) eigenvectors with corresponding eigenvalues larger than the threshold. Multiply these eigenvectors by the square root of the eigenvalues and reshape them to \(r \times n \times n\) matrices to obtain \(L^{(r)}\).
Whereas for the Cholesky decomposition, keep the first \(r\) Cholesky vectors that result in an residual error below the threshold and reshape them to \(r \times n \times n\) matrices to obtain \(L^{(r)}\).
Diagonalize the \(L^{(r)}\) (\(n \times n\)) matrices and for each matrix keep the eigenvalues (and their corresponding eigenvectors) that are larger than a threshold.
Compute the orthonormal matrices \(U\) and the symmetric matrices \(Z\) from the retained eigenvalues and eigenvectors to get the core and leaf tensors.