API Documentation#
Constants#
- __version__ = '0.6.1'#
PyCI version string.
Numerical Types#
- c_long = dtype('int64')#
Signed integer C++ dtype.
- c_ulong = dtype('uint64')#
Unsigned integer C++ dtype.
- c_double = dtype('float64')#
Floating point C++ dtype.
Classes#
Second-quantized operator class#
- class secondquant_op#
Second-quantized operator class.
For arbitrary-seniority systems:
\[O = \sum_{pq}{t_{pq} a^\dagger_p a_q} + \sum_{pqrs}{g_{pqrs} a^\dagger_p a^\dagger_q a_s a_r}\]For seniority-zero systems:
\[\]O = sum_{p}{h_p N_p} + sum_{p neq q}{v_{pq} P^dagger_p P_q} + sum_{pq}{w_{pq} N_p N_q}
where
\[h_{p} = \left<p|T|p\right> = t_{pp}\]\[v_{pq} = \left<pp|V|qq\right> = g_{ppqq}\]\[w_{pq} = 2 \left<pq|V|pq\right> - \left<pq|V|qp\right> = 2 * g_{pqpq} - g_{pqqp}\]- property ecore#
Constant (or "zero-particle") integral.
- Returns:
ecore – Constant (or "zero-particle") integral.
- Return type:
float
- property h#
Seniority-zero one-particle molecular integral array.
- Returns:
h – Seniority-zero one-particle molecular integral array.
- Return type:
numpy.ndarray
- property nbasis#
Number of spatial orbital basis functions.
- Returns:
nbasis – Number of spatial orbital basis functions.
- Return type:
int
- property one_mo#
One-particle molecular integral array.
- Returns:
one_mo – One-particle molecular integral array.
- Return type:
numpy.ndarray
- to_file()#
Write this operator to an FCIDUMP file.
- Parameters:
filename (TextIO) – Name of FCIDUMP file to write.
nelec (int, default=0) – Number of electrons to write.
ms2 (int, default=0) – Spin number to write.
tol (float, default=0.0) – Only write integrals with magnitude larger than this value.
- property two_mo#
Two-particle molecular integral array.
- Returns:
two_mo – Two-particle molecular integral array.
- Return type:
numpy.ndarray
- property v#
Seniority-zero two-particle molecular integral array.
- Returns:
v – Seniority-zero two-particle molecular integral array.
- Return type:
numpy.ndarray
- property w#
Seniority-two two-particle molecular integral array.
- Returns:
w – Seniority-two two-particle molecular integral array.
- Return type:
numpy.ndarray
- hamiltonian#
alias of
secondquant_op
Base wave function classes#
- class wavefunction#
Wave function base class.
- property nbasis#
Number of spatial orbital functions.
- Returns:
nbasis – Number of spatial orbital functions.
- Return type:
int
- property nocc#
Number of occupied spin-orbitals.
- Returns:
nocc – Number of occupied spin-orbitals.
- Return type:
int
- property nocc_dn#
Number of occupied spin-down orbitals.
- Returns:
nocc_dn – Number of occupied spin-down orbitals.
- Return type:
int
- property nocc_up#
Number of occupied spin-up orbitals.
- Returns:
nocc_up – Number of occupied spin-up orbitals.
- Return type:
int
- property nvir#
Number of virtual spin-orbitals.
- Returns:
nvir – Number of virtual spin-orbitals.
- Return type:
int
- property nvir_dn#
Number of virtual spin-down orbitals.
- Returns:
nvir_dn – Number of virtual spin-down orbitals.
- Return type:
int
- property nvir_up#
Number of virtual spin-up orbitals.
- Returns:
nvir_up – Number of virtual spin-up orbitals.
- Return type:
int
- squeeze()#
Free any unused memory allocated to this object.
- class one_spin_wfn#
One-spin wave function base class.
- add_all_dets()#
Add all determinants to the wave function.
- Parameters:
nthread (int) – Number of threads to use.
- add_det()#
Add determinant
det
to the wave function.- Parameters:
det (numpy.ndarray) – Determinant.
- Returns:
index – Index of added determinant in wave function, or -1 if the determinant was not added.
- Return type:
int
- add_dets_from_wfn()#
Add the determinants from another wave function.
- Parameters:
wfn (pyci.one_spin_wfn) – Wave function.
- add_excited_dets()#
Add excited determinants to the wave function.
- Parameters:
exc (int) – Excitation order.
ref (numpy.ndarray, default=None) – Reference determinant. Default is the Hartree-Fock determinant.
- add_hartreefock_det()#
Add the Hartree-Fock determinant to the wave function.
- add_occs()#
Add occupation vector
occs
to the wave function.- Parameters:
occs (numpy.ndarray) – Occupation vector.
- Returns:
index – Index of added determinant in wave function, or -1 if the determinant was not added.
- Return type:
int
- index_det()#
Return the index of determinant
det
in the wave function.If the determinant is not in the wave function, this function returns -1.
- Parameters:
det (numpy.ndarray) – Determinant.
- Returns:
index – Index of determinant or -1.
- Return type:
int
- index_det_from_rank()#
Return the index of determinant with rank
rank
in the wave function.If the determinant is not in the wave function, this function returns -1.
- Parameters:
rank (int) – Rank of determinant.
- Returns:
index – Index of determinant or -1.
- Return type:
int
- rank_det()#
Return the rank of determinant
det
.- Parameters:
det (numpy.ndarray) – Determinant.
- Returns:
rank – Rank of determinant.
- Return type:
int
- reserve()#
Reserve space in memory for
n
determinants in the wave function object.- Parameters:
n (int) – Number of determinants for which to reserve space.
- to_det_array()#
Return a section of the wave function as a numpy.ndarray of determinants.
Arguments
low
andhigh
behave like the first two arguments torange
. By default, if neither are specified, the whole determinant is returned.- Parameters:
low (int, default=-1) – Beginning of section.
end (int, default=-1) – End of section.
- Returns:
array – Array of determinants.
- Return type:
numpy.ndarray
- to_file()#
Write the wave function to a binary file.
- Parameters:
filename (TextIO) – Name of the file to write.
- to_occ_array()#
Return a section of the wave function as a numpy.ndarray of occupation vectors.
Arguments
low
andhigh
behave like the first two arguments torange()
. By default, if neither are specified, the whole determinant is returned.- Parameters:
low (int, default=-1) – Beginning of section.
end (int, default=-1) – End of section.
- Returns:
array – Array of occupation vectors.
- Return type:
numpy.ndarray
- class two_spin_wfn#
Two-spin wave function base class.
- add_all_dets()#
Add all determinants to the wave function.
- Parameters:
nthread (int) – Number of threads to use.
- add_det()#
Add determinant
det
to the wave function.- Parameters:
det (numpy.ndarray) – Determinant.
- Returns:
index – Index of added determinant in wave function, or -1 if the determinant was not added.
- Return type:
int
- add_dets_from_wfn()#
Add the determinants from another wave function.
- Parameters:
wfn (pyci.two_spin_wfn) – Wave function.
- add_excited_dets()#
Add excited determinants to the wave function.
- Parameters:
exc (int) – Excitation order.
ref (numpy.ndarray, default=None) – Reference determinant. Default is the Hartree-Fock determinant.
- add_hartreefock_det()#
Add the Hartree-Fock determinant to the wave function.
- add_occs()#
Add occupation vector
occs
to the wave function.- Parameters:
occs (numpy.ndarray) – Occupation vector.
- Returns:
index – Index of added determinant in wave function, or -1 if the determinant was not added.
- Return type:
int
- index_det()#
Return the index of determinant
det
in the wave function.If the determinant is not in the wave function, this function returns -1.
- Parameters:
det (numpy.ndarray) – Determinant.
- Returns:
index – Index of determinant or -1.
- Return type:
int
- index_det_from_rank()#
Return the index of determinant with rank
rank
in the wave function.If the determinant is not in the wave function, this function returns -1.
- Parameters:
rank (int) – Rank of determinant.
- Returns:
index – Index of determinant or -1.
- Return type:
int
- rank_det()#
Return the rank of determinant
det
.- Parameters:
det (numpy.ndarray) – Determinant.
- Returns:
rank – Rank of determinant.
- Return type:
int
- reserve()#
Reserve space in memory for
n
determinants in the wave function object.- Parameters:
n (int) – Number of determinants for which to reserve space.
- to_det_array()#
Return a section of the wave function as a numpy.ndarray of determinants.
Arguments
low
andhigh
behave like the first two arguments torange
. By default, if neither are specified, the whole determinant is returned.- Parameters:
low (int, default=-1) – Beginning of section.
end (int, default=-1) – End of section.
- Returns:
array – Array of determinants.
- Return type:
numpy.ndarray
- to_file()#
Write the wave function to a binary file.
- Parameters:
filename (TextIO) – Name of the file to write.
- to_occ_array()#
Return a section of the wave function as a numpy.ndarray of occupation vectors.
Arguments
low
andhigh
behave like the first two arguments torange()
. By default, if neither are specified, the whole determinant is returned.- Parameters:
low (int, default=-1) – Beginning of section.
end (int, default=-1) – End of section.
- Returns:
array – Array of occupation vectors.
- Return type:
numpy.ndarray
DOCI wave function#
- class doci_wfn#
DOCI wave function base class.
FullCI wave function#
- class fullci_wfn#
FullCI wave function base class.
Generalized CI wave function#
- class genci_wfn#
Generalized CI wave function base class.
Sparse matrix operator#
- class sparse_op#
Sparse matrix operator class.
- data()#
Return CSR matrix data vector
- property dtype#
Data type of matrix.
- Returns:
dtype – Data type of matrix.
- Return type:
numpy.dtype
- property ecore#
Constant (or “zero-particle”) integral.
- Returns:
ecore – Constant (or “zero-particle”) integral.
- Return type:
float
- get_element()#
Return the \(\left(i, j\right)\)-th element of the sparse matrix operator.
- indices()#
Return CSR matrix indices vector
- indptr()#
Return CSR matrix index pointer vector
- matvec()#
Compute the matrix vector product of the sparse matrix operator with vector
x
.\[A \mathbf{x} = \mathbf{y}\]- Parameters:
x (numpy.ndarray) – Vector to which the operator will be applied.
out (numpy.ndarray, default=None) – Array in which to store the result. One will be created if this is not specified.
- Returns:
y – Result vector.
- Return type:
numpy.ndarray
- reserve()#
Reserve space in memory for
n
nonzero elements in the sparse matrix operator.- Parameters:
n (int) – Number of elements for which to reserve space.
- property shape#
Shape of the matrix.
- Returns:
nrow (int) – Number of rows.
ncol (int) – Number of columns.
- property size#
Number of non-zero matrix elements.
- Returns:
size – Number of non-zero matrix elements.
- Return type:
int
- solve()#
Solve a CI eigenproblem.
- Parameters:
n (int, default=1) – Number of lowest eigenpairs to find.
c0 (np.ndarray, default=[1,0,...,0]) – Initial guess for lowest eigenvector.
ncv (int, default=min(nrow, max(2 * n + 1, 20))) – Number of Lanczos vectors to use.
maxiter (int, default=nrow * n * 10) – Maximum number of iterations to perform.
tol (float, default=1.0e-12) – Convergence tolerance.
- Returns:
es (np.ndarray) – Energies.
cs (np.ndarray) – Coefficient vectors.
- squeeze()#
Free any unused memory allocated to this object.
- property symmetric#
Whether the sparse matrix operator is symmetric/Hermitian.
- Returns:
symmetric – Whether the sparse matrix operator is symmetric/Hermitian.
- Return type:
bool
- update()#
Update a sparse matrix operator for the HCI algorithm.
- Parameters:
ham (pyci.secondquant_op) – Hamiltonian.
wfn (pyci.wavefunction) – Wave function.
Notes
The user is responsible for using the same
ham
andwfn
. This method doesn’t work with rectangular operators. Those must be re-initialized from the wave function object.
Functions#
Threading#
The number of threads to use is decided as follows, in decreasing priority:
The value of a keyword argument
nthread=
to a particular PyCI function or method callThe value of
n
from the most recent invokation ofpyci.set_num_threads(n)
The value of environment variable
PYCI_NUM_THREADS
when Python is startedThe number of threads supported by the hardware (
std::thread::hardware_concurrency()
)
- get_num_threads()#
Return the default number of threads to use.
- Returns:
nthread – Number of threads.
- Return type:
int
- set_num_threads()#
Set the default number of threads to use.
- Parameters:
nthread (int) – Number of threads.
Selected CI routines#
- add_hci()#
Add determinants to a wave function by running an iteration of Heat-Bath CI [HCI1].
[HCI1]Holmes, Adam A., Norm M. Tubman, and C. J. Umrigar. “Heat-bath configuration interaction: An efficient selected configuration interaction algorithm inspired by heat-bath sampling.” Journal of chemical theory and computation 12.8 (2016): 3674-3680.
- Parameters:
ham (pyci.secondquant_op) – Hamiltonian.
wfn (pyci.wavefunction) – Wave function.
coeffs (numpy.ndarray) – Coefficient vector.
eps (float, default=1.0e-5) – \(\epsilon\) value for Heat-Bath CI routine.
nthread (int) – Number of threads to use.
- Returns:
ndet – Number of determinants added.
- Return type:
int
- add_excitations(wfn, *excitations, ref=None)#
Add excited determinants to a wave function.
- Parameters:
wfn (pyci.wavefunction) – Wave function.
excitations (Sequence[int]) – List of excitation levels of determinants to add.
ref (numpy.ndarray, optional) – Reference determinant by which to determine excitation levels. Default is the Hartree-Fock determinant.
- add_seniorities(wfn, *seniorities)#
Add determinants of the specified seniority/ies to the wave function.
- Parameters:
wfn (pyci.fullci_wfn) – FullCI wave function.
seniorities (Sequence[int]) – List of seniorities of determinants to add.
- add_gkci(wfn, t=-0.5, p=1.0, mode='cntsp', dim=3, energies=None, width=None)#
Add determinants to the wave function according to the odometer algorithm (Griebel-Knapeck CI) [GKCI1].
[GKCI1]Anderson, James SM, Farnaz Heidar-Zadeh, and Paul W. Ayers. “Breaking the curse of dimension for the electronic Schrödinger equation with functional analysis.” Computational and Theoretical Chemistry 1142 (2018): 66-77.
- Parameters:
wfn (pyci.wavefunction) – Wave function.
t (float, default=-0.5) – Smoothness factor.
p (float, default=1.0) – Cost factor.
mode (Sequence[int] or ('cntsp' | 'gamma' | 'interval'), default='cntsp') – Node pattern.
dim (int, default=3) – Number of nodes (for ‘gamma’ mode).
energies (np.ndarray, optional) – Orbital energies (required for ‘interval’ mode).
width (float, optional) – Width of one interval (required for ‘interval’ mode).
Post-CI routines#
- compute_rdms()#
Compute the one- and two- particle reduced density matrices (RDMs) of a wave function.
- Parameters:
wfn (pyci.wavefunction) – Wave function.
coeffs (numpy.ndarray) – Coefficient vector.
- Returns:
d1 (numpy.ndarray) – One-particle RDM matrix.
d2 (numpy.ndarray) – Two-particle RDM matrix.
Notes
For DOCI wave functions, this method returns two nbasis-by-nbasis matrices, which include the unique seniority-zero and seniority-two terms from the full 2-RDMs:
\[D_0 = \left<pp|qq\right>\]\[D_2 = \left<pq|pq\right>\]The diagonal elements of \(D_0\) are equal to the 1-RDM elements \(\left<p|p\right>\).
For FullCI wave functions, the leading dimension of
rdm1
has length 2 and specifies the spin-block 0) “up-up” or 1) “down-down”, and the leading dimensions ofrdm2
has length 3 and specifies the spin-block 0) “up-up-up-up”, 1) “down-down-down-down’, or 2) “up-down-up-down”.For Generalized CI wave functions,
rdm1
andrdm2
are the full 1-RDM and 2-RDM, respectively.
- spinize_rdms(d1, d2)#
Convert the DOCI matrices or FullCI RDM spin-blocks to full, generalized RDMs.
- Parameters:
d1 (numpy.ndarray) – \(D_0\) matrix or FullCI 1-RDM spin-blocks.
d2 (numpy.ndarray) – \(D_2\) matrix or FullCI 2-RDM spin-blocks.
- Returns:
rdm1 (numpy.ndarray) – Generalized one-particle RDM.
rdm2 (numpy.ndarray) – Generalized two-particle RDM.
- compute_enpt2()#
Compute the second-order multi-reference Epstein-Nesbet (ENPT2) energy for a wave function.
- Parameters:
ham (pyci.secondquant_op) – Hamiltonian.
wfn (pyci.wavefunction) – Wave function.
coeffs (numpy.ndarray) – Coefficient vector.
energy (float) – Variational CI energy for this wave function and Hamiltonian.
eps (float, default=1.0e-5) – \(\epsilon\) value for ENPT2 routine.
nthread (int) – Number of threads to use.
- Returns:
pt_energy – ENPT2 energy.
- Return type:
float
- compute_overlap()#
Compute the overlap \(\left<\Psi_1|\Psi_2\right>\) of two wave functions.
- Parameters:
wfn1 (pyci.wavefunction) – First wave function.
wfn2 (pyci.wavefunction) – Second wave function.
coeffs1 (numpy.ndarray) – First coefficient vector.
coeffs2 (numpy.ndarray) – Second coefficient vector.
- Returns:
olp – Overlap.
- Return type:
float
Integral transformations#
- make_senzero_integrals(one_mo, two_mo)#
Return the non-zero chunks for seniority-zero of the full one- and two- particle integrals.
- Parameters:
one_mo (numpy.ndarray) – Full one-particle integral array.
two_mo (numpy.ndarray) – Full two-particle integral array.
- Returns:
h (numpy.ndarray) – Seniority-zero one-particle integrals.
v (numpy.ndarray) – Seniority-zero two-particle integrals.
w (numpy.ndarray) – Seniority-two two-particle integrals.
- reduce_senzero_integrals(h, v, w, nocc)#
Reduce the non-zero chunks for seniority-zero of the one- and two- particle integrals.
- Parameters:
h (numpy.ndarray) – Seniority-zero one-particle integrals.
v (numpy.ndarray) – Seniority-zero two-particle integrals.
w (numpy.ndarray) – Seniority-two two-particle integrals.
nocc (int) – Number of pair-occupied orbitals.
- Returns:
rv (numpy.ndarray) – Reduced seniority-zero two-particle integrals.
rw (numpy.ndarray) – Reduced seniority-two two-particle integrals.
FanCI classes#
- class AP1roG(ham: secondquant_op, nocc: int, nproj: int | None = None, wfn: doci_wfn | None = None, **kwargs: Any)#
AP1roG FanCI class.
The AP1roG wave function is the Antisymmetrized Product of 1-reference-orbital Geminals [AP1roG1], which has the form
\[\left|\Psi\right> = \prod_{p=1}^{N_\text{elec}/2}{\left( a^\dagger_{\alpha;p} a^\dagger_{\beta;p} + \sum_{k=N/2+1}^K{ c_{pk} a^\dagger_{\alpha;k} a^\dagger_{\beta;k} } \right)} \left|\psi_0\right> \,.\][AP1roG1]Limacher, Peter A., et al. “A new mean-field method suitable for strongly correlated electrons: Computationally facile antisymmetric products of nonorthogonal geminals.” Journal of chemical theory and computation 9.3 (2013): 1394-1401.
- compute_jacobian(x: ndarray) ndarray #
Compute the Jacobian of the FanCI objective function.
j : x[k] -> y[n, k]
- Parameters:
x (np.ndarray) – Parameter array, [p_0, p_1, …, p_n, E].
- Returns:
jac – Jacobian matrix.
- Return type:
np.ndarray
- compute_objective(x: ndarray) ndarray #
Compute the FanCI objective function.
f : x[k] -> y[n]
- Parameters:
x (np.ndarray) – Parameter array, [p_0, p_1, …, p_n, E].
- Returns:
obj – Objective vector.
- Return type:
np.ndarray
- compute_overlap(x: ndarray) ndarray #
Compute the FanCI overlap vector.
- Parameters:
x (np.ndarray) – Parameter array, [p_0, p_1, …, p_n].
- Returns:
ovlp – Overlap array.
- Return type:
np.ndarray
- compute_overlap_deriv(x: ndarray) ndarray #
Compute the FanCI overlap derivative matrix.
- Parameters:
x (np.ndarray) – Parameter array, [p_0, p_1, …, p_n].
- Returns:
ovlp – Overlap derivative array.
- Return type:
np.ndarray
- class APIG(ham: secondquant_op, nocc: int, nproj: int | None = None, wfn: doci_wfn | None = None, **kwargs: Any)#
APIG FanCI class.
The APIG wave function is the Antisymmetrized Product of Interacting Geminals [APIG1], which has the form
\[\left|\Psi\right> = \prod_{p=1}^{N_\text{elec}/2}{ \sum_{k=1}^K{ \left( c_{pk} a^\dagger_{\alpha;k} a^\dagger_{\beta;k} \right) } } \left|\psi_0\right> \,.\][APIG1]Silver, David M. “Natural orbital expansion of interacting geminals.” The Journal of Chemical Physics 50.12 (1969): 5108-5116.
- compute_jacobian(x: ndarray) ndarray #
Compute the Jacobian of the FanCI objective function.
j : x[k] -> y[n, k]
- Parameters:
x (np.ndarray) – Parameter array, [p_0, p_1, …, p_n, E].
- Returns:
jac – Jacobian matrix.
- Return type:
np.ndarray
- compute_objective(x: ndarray) ndarray #
Compute the FanCI objective function.
f : x[k] -> y[n]
- Parameters:
x (np.ndarray) – Parameter array, [p_0, p_1, …, p_n, E].
- Returns:
obj – Objective vector.
- Return type:
np.ndarray
- compute_overlap(x: ndarray) ndarray #
Compute the FanCI overlap vector.
- Parameters:
x (np.ndarray) – Parameter array, [p_0, p_1, …, p_n].
- Returns:
ovlp – Overlap array.
- Return type:
np.ndarray
- compute_overlap_deriv(x: ndarray) ndarray #
Compute the FanCI overlap derivative matrix.
- Parameters:
x (np.ndarray) – Parameter array, [p_0, p_1, …, p_n].
- Returns:
ovlp – Overlap derivative array.
- Return type:
np.ndarray
- class DetRatio(ham: secondquant_op, nocc: int, numerator: int, denominator: int, nproj: int | None = None, wfn: doci_wfn | None = None, **kwargs: Any)#
Determinant ratio FanCI class. This wave function has the form
\[\left|\Psi\right> = \sum_k{c_k\prod_{i=1}^N{a^\dagger_{k_i}}\left|\theta\right>}\]where
\[c_k = \frac{ \prod_{n=1}^{n_\text{numer.}}{\text{det}\left[\vec{x_n;}_{k_1}\mid\cdots\mid\vec{x_n;}_{k_N}\right]} }{ \prod_{d=1}^{n_\text{denom.}}{\text{det}\left[\vec{y_d;}_{k_1}\mid\cdots\mid\vec{y_d;}_{k_N}\right]} } \,,\]with the vectors \(\left\{x_n\right\}\) and \(\left\{y_d\right\}\) being the parameters of the wave function which are concatenated into matrices according to the occupied indices of each \(\psi_k\) prior to their determinants being evaluated. For seniority-zero systems, \(n_\text{numer.} + n_\text{denom.}\) is even, and otherwise it is odd.
- compute_overlap(x: ndarray, occs_array: ndarray | str) ndarray #
Compute the FanCI overlap vector.
- Parameters:
x (np.ndarray) – Parameter array, [p_0, p_1, …, p_n].
occs_array ((np.ndarray | 'P' | 'S')) – Array of determinant occupations for which to compute overlap. A string “P” or “S” can be passed instead that indicates whether
occs_array
corresponds to the “P” space or “S” space, so that a more efficient, specialized computation can be done for these.
- Returns:
ovlp – Overlap array.
- Return type:
np.ndarray
- compute_overlap_deriv(x: ndarray, occs_array: ndarray | str) ndarray #
Compute the FanCI overlap derivative matrix.
- Parameters:
x (np.ndarray) – Parameter array, [p_0, p_1, …, p_n].
occs_array ((np.ndarray | 'P' | 'S')) – Array of determinant occupations for which to compute overlap. A string “P” or “S” can be passed instead that indicates whether
occs_array
corresponds to the “P” space or “S” space, so that a more efficient, specialized computation can be done for these.
- Returns:
ovlp – Overlap derivative array.
- Return type:
np.ndarray
- class pCCDS(ham: secondquant_op, nocc_up: int, nocc_dn: int, nproj: int | None = None, wfn: fullci_wfn | None = None, **kwargs: Any)#
Pair coupled cluster doubles + singles (or “AP1roGSDspin”) FanCI class.
\[\left| \Psi_\text{pCCD+S} \right> = \prod_{i=1}^{N/2} \left( 1 + \sum_{a=N/2+1}^K t_{i\bar{i};a\bar{a}} \hat{\tau}^{a\bar{a}}_{i\bar{i}} \right) \prod_{i=1}^{N/2} \left( 1 + \sum_{a=N/2+1}^K t_{i;a\bar{a}} \hat{\tau}^{a}_{i} \hat{n}_{\bar{i}} \right) \left| \psi_0 \right> \,.\]This is AP1roG supplemented with single excitations that preserve the spin number, i.e., no alpha to beta excitation is allowed, or vice-versa. Furthermore, the single excitations added are those that break an occupied pair of spin orbitals, increasing the seniority number. This determines the sen-o part of the name (from every occupied pair only one element of the pair is excited) [pCCDS1].
[pCCDS1]Gaikwad, Pratiksha B., et al. “Coupled Cluster-Inspired Geminal Wavefunctions.” arXiv preprint arXiv:2310.01764 (2023).
- compute_overlap(x: ndarray, occs_array: ndarray | str) ndarray #
Compute the pCCDS overlap vector.
The occupation vector is described in terms of the occupied/virtual indexes of the excitation that generates it from the reference. The indexes are further identified as corresponding to excitations of pairs or single spin-orbitals. The overlap is computed as the product of the permanents of the matrices formed by the pair and single excitations indexes.
- Parameters:
x (np.ndarray) – Parameter array, [p_0, p_1, …, p_n].
occs_array ((np.ndarray | 'P' | 'S')) – Array of determinant occupations for which to compute overlap. A string “P” or “S” can be passed instead that indicates whether
occs_array
corresponds to the “P” space or “S” space, so that a more efficient, specialized computation can be done for these.
- Returns:
ovlp – Overlap array.
- Return type:
np.ndarray
- compute_overlap_deriv(x: ndarray, occs_array: ndarray | str) ndarray #
Compute the FanCI overlap derivative matrix.
- Parameters:
x (np.ndarray) – Parameter array, [p_0, p_1, …, p_n].
occs_array ((np.ndarray | 'P' | 'S')) – Array of determinant occupations for which to compute overlap. A string “P” or “S” can be passed instead that indicates whether
occs_array
corresponds to the “P” space or “S” space, so that a more efficient, specialized computation can be done for these.
- Returns:
ovlp – Overlap derivative array.
- Return type:
np.ndarray