sf_ip_ea package¶
Submodules¶
sf_ip_ea.bloch module¶
Bloch effective Hamiltonian solver.
This module constructs a Bloch effective Hamiltonian for a given SF-IP/EA wavefunction. Note that this module is currently dependent on Psi4.
-
sf_ip_ea.bloch.do_bloch(wfn, n_sites, site_list=None, site_list_orbs=None, molden_file='orbs.molden', skip_localization=False, neutral=False)¶ Bloch effective Hamiltonian solver.
Solves the Bloch effective Hamiltonian and returns a matrix containing J coupling information. Sites are designated using
site_listorsite_list_orbs; if neither of these keywords are specified, each orbital in the CAS/RAS2 space is assumed to be its own site. A Molden file containing the orbitals is written toorbs.molden(or whichever file is specified by the user). The J coupling values are also written to standard output.- Parameters
- sf_wfnsf_wfn
SF-IP-EA wfn object containing info about the calculation.
- n_sitesint
The number of sites.
- site_listlist
List of which atoms are “sites”. If this is not given, the program assumes one orbital per site. Atomic center ordering starts at zero. Optional.
- site_list_orbslist
A list of sites, using lists of orbitals rather than atomic centers. (It’s a list of lists. For example, for a two-site case where MOs 55, 56, and 59 are on one site and the remaining orbitals are on the other, use
[[55,56,59],[57,58,60]].) Note that ordering starts at 1, not zero, so it follows the same indexing as the MO printing in the Psi4 output files. Optional.- molden_filestring
Molden filename to which orbitals are written. Optional. Defaults to
orbs.molden.- skip_localizationbool
Whether to skip orbital localization. If true, the user should localize the orbitals in wfn.wfn beforehand! Optional. Defaults to False.
- neutralbool
Calculate the J couplings using neutral determinants only. Optional. Defaults to False.
- Returns
- numpy.ndarray
NumPy matrix of J coupling values
-
sf_ip_ea.bloch.lowdin_orth(A)¶ Performs Lowdin orthonormalization on a given matrix A.
Orthonormalizes a given NumPy matrix A based on Lowdin’s approach (i.e. based on the SVD of A).
- Parameters
- Anumpy.ndarray
NumPy array to orthogonalize.
- Returns
- numpy.ndarray
An orthogonalized version of A.
sf_ip_ea.f module¶
Handling of the Fock matrix.
This module is responsible for handling things related to the Fock matrix.
-
sf_ip_ea.f.get_F(wfn)¶ Gets alpha and beta Fock matrices.
Given a Psi4 Wavefunction object, this returns NumPy representations of the alpha and beta Fock matrices in the spatial orbital MO basis.
- Parameters
- wfnpsi4.core.Wavefunction
Psi4 wavefunction object from which to obtain alpha and beta Fock matrices.
- Returns
- tuple
A tuple (Fa, Fb) containing NumPy representations of the alpha and beta Fock matrices, respectively.
sf_ip_ea.linop module¶
LinOp module responsible for performing Hamiltonian-vector multiplication.
This module contains a derived class of NumPy’s LinearOperator. This
class contains all of the information necessary to do a Hamiltonian-vector
multiply (the most expensive step in Davidson) using NumPy’s einsum
method.
-
class
sf_ip_ea.linop.LinOpH(shape_in, offset_in, ras1_in, ras2_in, ras3_in, Fa_in, Fb_in, tei_in, n_SF_in, delta_ec_in, conf_space_in)¶ Bases:
scipy.sparse.linalg.interface.LinearOperatorLinear operator for Hamiltonian-vector multiplication.
This is a derived class of NumPy’s LinearOperator class. It contains information on how to perform tensor-contraction-based multiplication of the Hamiltonian with a given vector or set of vectors.
- Attributes
- offsetdouble
Diagonal offset to Hamiltonian (usually reference energy).
- ras1int
The RAS1 space (doubly occupied in reference).
- ras2int
The RAS2 space (singly occupied in reference).
- ras3int
The RAS3 space (virtual in reference).
- n_SFint
Number of spin-flips to perform.
- delta_ecint
Change in electron count (indicates IP/EA).
- conf_spacestring
Excitations to include (hole, particle, etc).
- num_detsint
Number of determinants.
- Fanumpy.ndarray
Alpha Fock matrix.
- Fbnumpy.ndarray
Beta Fock matrix.
- teiTEI
Two-electron integrals.
Methods
__call__(self, x)Call self as a function.
adjoint(self)Hermitian adjoint.
diag(self)Returns approximate diagonal of Hamiltonian.
do_cas_1sf(self, v, Fa, Fb, tei, offset_v, …)Do CAS-1SF.
do_cas_1sf_ea(self, v, Fa, Fb, tei, …)Do CAS-1SF-EA.
do_cas_1sf_ip(self, v, Fa, Fb, tei, …)Do CAS-1SF-IP.
do_cas_1sf_neutral(self, v, Fa, Fb, tei, …)Do CAS-1SF with neutral determinants only.
do_cas_2sf(self, v, Fa, Fb, tei, offset_v, …)Do CAS-2SF.
do_cas_ea(self, v, Fa, Fb, tei, offset_v, …)Do CAS-EA.
do_cas_ip(self, v, Fa, Fb, tei, offset_v, …)Do CAS-IP.
do_h_1sf(self, v, Fa, Fb, tei, offset_v, …)Do RAS(h)-1SF.
do_h_1sf_ip(self, v, Fa, Fb, tei, offset_v, …)Do RAS(h)-1SF-IP.
do_h_ea(self, v, Fa, Fb, tei, offset_v, …)Do RAS(h)-EA.
do_h_ip(self, v, Fa, Fb, tei, offset_v, …)Do RAS(h)-IP.
do_hp_1sf(self, v, Fa, Fb, tei, offset_v, …)Do RAS(h,p)-1SF.
do_p_1sf(self, v, Fa, Fb, tei, offset_v, …)Do RAS(p)-1SF.
do_p_1sf_ea(self, v, Fa, Fb, tei, offset_v, …)Do RAS(p)-1SF-EA.
do_p_ea(self, v, Fa, Fb, tei, offset_v, …)Do RAS(p)-EA.
do_p_ip(self, v, Fa, Fb, tei, offset_v, …)Do RAS(p)-IP.
dot(self, x)Matrix-matrix or matrix-vector multiplication.
matmat(self, X)Matrix-matrix multiplication.
matvec(self, x)Matrix-vector multiplication.
rmatvec(self, x)Adjoint matrix-vector multiplication.
transpose(self)Transpose this linear operator.
-
diag(self)¶ Returns approximate diagonal of Hamiltonian.
This returns the approximate diagonal of the Hamiltonian for the Davidson code (used in the preconditioner step). This is calculated for each determinant by summing over the occupied orbital energies.
- Returns
- np.ndarray
1-D NumPy representation of the diagonal of the Hamiltonian.
-
do_cas_1sf(self, v, Fa, Fb, tei, offset_v, ras1, ras2, ras3)¶ Do CAS-1SF.
Defines H*v for a CAS-1SF calculation. Blocks are defined as:
block1 = v(ia:ab)
- Parameters
- vnumpy.ndarray
Guess vector to multiply.
- Fanumpy.ndarray
Alpha Fock matrix.
- Fbnumpy.ndarray
Beta Fock matrix.
- teiTEI
Two-electron integral object.
- offset_vnumpy.ndarray
Vector of offset values for guess vector (usually ROHF energy)
- ras1int
Number of RAS1 orbitals.
- ras2int
Number of RAS2 orbitals.
- ras3int
Number of RAS3 orbitals.
- Returns
- The result of H*v.
-
do_cas_1sf_ea(self, v, Fa, Fb, tei, offset_v, ras1, ras2, ras3)¶ Do CAS-1SF-EA.
Defines H*v for a CAS-1SF-EA calculation. Blocks are defined as:
block1 = v(ija:aab)
- Parameters
- vnumpy.ndarray
Guess vector to multiply.
- Fanumpy.ndarray
Alpha Fock matrix.
- Fbnumpy.ndarray
Beta Fock matrix.
- teiTEI
Two-electron integral object.
- offset_vnumpy.ndarray
Vector of offset values for guess vector (usually ROHF energy)
- ras1int
Number of RAS1 orbitals.
- ras2int
Number of RAS2 orbitals.
- ras3int
Number of RAS3 orbitals.
- Returns
- The result of H*v.
-
do_cas_1sf_ip(self, v, Fa, Fb, tei, offset_v, ras1, ras2, ras3)¶ Do CAS-1SF-IP.
Defines H*v for a CAS-1SF-IP calculation. Blocks are defined as:
block1 = v(ija:aab)
- Parameters
- vnumpy.ndarray
Guess vector to multiply.
- Fanumpy.ndarray
Alpha Fock matrix.
- Fbnumpy.ndarray
Beta Fock matrix.
- teiTEI
Two-electron integral object.
- offset_vnumpy.ndarray
Vector of offset values for guess vector (usually ROHF energy)
- ras1int
Number of RAS1 orbitals.
- ras2int
Number of RAS2 orbitals.
- ras3int
Number of RAS3 orbitals.
- Returns
- The result of H*v.
-
do_cas_1sf_neutral(self, v, Fa, Fb, tei, offset_v, ras1, ras2, ras3)¶ Do CAS-1SF with neutral determinants only.
Defines H*v for a CAS-1SF calculation, including only neutral determinants in the Fock space. Blocks are defined as:
block1 = v(ia:ab) where i==a
- Parameters
- vnumpy.ndarray
Guess vector to multiply.
- Fanumpy.ndarray
Alpha Fock matrix.
- Fbnumpy.ndarray
Beta Fock matrix.
- teiTEI
Two-electron integral object.
- offset_vnumpy.ndarray
Vector of offset values for guess vector (usually ROHF energy)
- ras1int
Number of RAS1 orbitals.
- ras2int
Number of RAS2 orbitals.
- ras3int
Number of RAS3 orbitals.
- Returns
- The result of H*v.
-
do_cas_2sf(self, v, Fa, Fb, tei, offset_v, ras1, ras2, ras3)¶ Do CAS-2SF.
Defines H*v for a CAS-2SF calculation. Blocks are defined as:
block1 = v(ijab:aabb)
- Parameters
- vnumpy.ndarray
Guess vector to multiply.
- Fanumpy.ndarray
Alpha Fock matrix.
- Fbnumpy.ndarray
Beta Fock matrix.
- teiTEI
Two-electron integral object.
- offset_vnumpy.ndarray
Vector of offset values for guess vector (usually ROHF energy)
- ras1int
Number of RAS1 orbitals.
- ras2int
Number of RAS2 orbitals.
- ras3int
Number of RAS3 orbitals.
- Returns
- The result of H*v.
-
do_cas_ea(self, v, Fa, Fb, tei, offset_v, ras1, ras2, ras3)¶ Do CAS-EA.
Defines H*v for a CAS-EA calculation. Blocks are defined as:
block1 = v(a:b)
- Parameters
- vnumpy.ndarray
Guess vector to multiply.
- Fanumpy.ndarray
Alpha Fock matrix.
- Fbnumpy.ndarray
Beta Fock matrix.
- teiTEI
Two-electron integral object.
- offset_vnumpy.ndarray
Vector of offset values for guess vector (usually ROHF energy)
- ras1int
Number of RAS1 orbitals.
- ras2int
Number of RAS2 orbitals.
- ras3int
Number of RAS3 orbitals.
- Returns
- The result of H*v.
-
do_cas_ip(self, v, Fa, Fb, tei, offset_v, ras1, ras2, ras3)¶ Do CAS-IP.
Defines H*v for a CAS-IP calculation. Blocks are defined as:
block1 = v(i:a)
- Parameters
- vnumpy.ndarray
Guess vector to multiply.
- Fanumpy.ndarray
Alpha Fock matrix.
- Fbnumpy.ndarray
Beta Fock matrix.
- teiTEI
Two-electron integral object.
- offset_vnumpy.ndarray
Vector of offset values for guess vector (usually ROHF energy)
- ras1int
Number of RAS1 orbitals.
- ras2int
Number of RAS2 orbitals.
- ras3int
Number of RAS3 orbitals.
- Returns
- The result of H*v.
-
do_h_1sf(self, v, Fa, Fb, tei, offset_v, ras1, ras2, ras3)¶ Do RAS(h)-1SF.
Defines H*v for a RAS(h)-1SF calculation. Blocks are defined as:
block1 = v(ai:ba) block2 = v(Ia:ab) block3 = v(Iiab:babb)
- Parameters
- vnumpy.ndarray
Guess vector to multiply.
- Fanumpy.ndarray
Alpha Fock matrix.
- Fbnumpy.ndarray
Beta Fock matrix.
- teiTEI
Two-electron integral object.
- offset_vnumpy.ndarray
Vector of offset values for guess vector (usually ROHF energy)
- ras1int
Number of RAS1 orbitals.
- ras2int
Number of RAS2 orbitals.
- ras3int
Number of RAS3 orbitals.
- Returns
- The result of H*v.
-
do_h_1sf_ip(self, v, Fa, Fb, tei, offset_v, ras1, ras2, ras3)¶ Do RAS(h)-1SF-IP.
Defines H*v for a RAS(h)-1SF-IP calculation. Blocks are defined as:
block1 = v(ija:aab) block2 = v(Iia:aab) block3 = v(Iijab:baabb)
- Parameters
- vnumpy.ndarray
Guess vector to multiply.
- Fanumpy.ndarray
Alpha Fock matrix.
- Fbnumpy.ndarray
Beta Fock matrix.
- teiTEI
Two-electron integral object.
- offset_vnumpy.ndarray
Vector of offset values for guess vector (usually ROHF energy)
- ras1int
Number of RAS1 orbitals.
- ras2int
Number of RAS2 orbitals.
- ras3int
Number of RAS3 orbitals.
- Returns
- The result of H*v.
-
do_h_ea(self, v, Fa, Fb, tei, offset_v, ras1, ras2, ras3)¶ Do RAS(h)-EA.
Defines H*v for a RAS(h)-EA calculation. Blocks are defined as:
block1 = v(a:b) block2 = v(Iab:bbb)
- Parameters
- vnumpy.ndarray
Guess vector to multiply.
- Fanumpy.ndarray
Alpha Fock matrix.
- Fbnumpy.ndarray
Beta Fock matrix.
- teiTEI
Two-electron integral object.
- offset_vnumpy.ndarray
Vector of offset values for guess vector (usually ROHF energy)
- ras1int
Number of RAS1 orbitals.
- ras2int
Number of RAS2 orbitals.
- ras3int
Number of RAS3 orbitals.
- Returns
- The result of H*v.
-
do_h_ip(self, v, Fa, Fb, tei, offset_v, ras1, ras2, ras3)¶ Do RAS(h)-IP.
Defines H*v for a RAS(h)-IP calculation. Blocks are defined as:
block1 = v(i:a) block2 = v(I:a) block3 = v(Iia:bab)
- Parameters
- vnumpy.ndarray
Guess vector to multiply.
- Fanumpy.ndarray
Alpha Fock matrix.
- Fbnumpy.ndarray
Beta Fock matrix.
- teiTEI
Two-electron integral object.
- offset_vnumpy.ndarray
Vector of offset values for guess vector (usually ROHF energy)
- ras1int
Number of RAS1 orbitals.
- ras2int
Number of RAS2 orbitals.
- ras3int
Number of RAS3 orbitals.
- Returns
- The result of H*v.
-
do_hp_1sf(self, v, Fa, Fb, tei, offset_v, ras1, ras2, ras3)¶ Do RAS(h,p)-1SF.
Defines H*v for a RAS(h,p)-1SF calculation. Blocks are defined as:
block1 = v(ia:ab) block2 = v(Ia:ab) block3 = v(iA:ab) block4 = v(Iiab:babb) block5 = v(ijAb:abaa)
- Parameters
- vnumpy.ndarray
Guess vector to multiply.
- Fanumpy.ndarray
Alpha Fock matrix.
- Fbnumpy.ndarray
Beta Fock matrix.
- teiTEI
Two-electron integral object.
- offset_vnumpy.ndarray
Vector of offset values for guess vector (usually ROHF energy)
- ras1int
Number of RAS1 orbitals.
- ras2int
Number of RAS2 orbitals.
- ras3int
Number of RAS3 orbitals.
- Returns
- The result of H*v.
-
do_p_1sf(self, v, Fa, Fb, tei, offset_v, ras1, ras2, ras3)¶ Do RAS(p)-1SF.
Defines H*v for a RAS(p)-1SF calculation. Blocks are defined as:
block1 = v(ai:ba) block2 = v(Ai:ba) block3 = v(Aaij:abaa)
- Parameters
- vnumpy.ndarray
Guess vector to multiply.
- Fanumpy.ndarray
Alpha Fock matrix.
- Fbnumpy.ndarray
Beta Fock matrix.
- teiTEI
Two-electron integral object.
- offset_vnumpy.ndarray
Vector of offset values for guess vector (usually ROHF energy)
- ras1int
Number of RAS1 orbitals.
- ras2int
Number of RAS2 orbitals.
- ras3int
Number of RAS3 orbitals.
- Returns
- The result of H*v.
-
do_p_1sf_ea(self, v, Fa, Fb, tei, offset_v, ras1, ras2, ras3)¶ Do RAS(p)-1SF-EA.
Defines H*v for a RAS(p)-1SF-EA calculation. Blocks are defined as:
block1 = v(iab:abb) block2 = v(Aia:bab) block3 = v(Aijab:baabb)
- Parameters
- vnumpy.ndarray
Guess vector to multiply.
- Fanumpy.ndarray
Alpha Fock matrix.
- Fbnumpy.ndarray
Beta Fock matrix.
- teiTEI
Two-electron integral object.
- offset_vnumpy.ndarray
Vector of offset values for guess vector (usually ROHF energy)
- ras1int
Number of RAS1 orbitals.
- ras2int
Number of RAS2 orbitals.
- ras3int
Number of RAS3 orbitals.
- Returns
- The result of H*v.
-
do_p_ea(self, v, Fa, Fb, tei, offset_v, ras1, ras2, ras3)¶ Do RAS(p)-EA.
Defines H*v for a RAS(p)-EA calculation. Blocks are defined as:
block1 = v(a:b) block2 = v(A:b) block3 = v(iAa:aab)
- Parameters
- vnumpy.ndarray
Guess vector to multiply.
- Fanumpy.ndarray
Alpha Fock matrix.
- Fbnumpy.ndarray
Beta Fock matrix.
- teiTEI
Two-electron integral object.
- offset_vnumpy.ndarray
Vector of offset values for guess vector (usually ROHF energy)
- ras1int
Number of RAS1 orbitals.
- ras2int
Number of RAS2 orbitals.
- ras3int
Number of RAS3 orbitals.
- Returns
- The result of H*v.
-
do_p_ip(self, v, Fa, Fb, tei, offset_v, ras1, ras2, ras3)¶ Do RAS(p)-IP.
Defines H*v for a RAS(p)-IP calculation. Blocks are defined as:
block1 = v(i:a) block2 = v(ijA:aaa)
- Parameters
- vnumpy.ndarray
Guess vector to multiply.
- Fanumpy.ndarray
Alpha Fock matrix.
- Fbnumpy.ndarray
Beta Fock matrix.
- teiTEI
Two-electron integral object.
- offset_vnumpy.ndarray
Vector of offset values for guess vector (usually ROHF energy)
- ras1int
Number of RAS1 orbitals.
- ras2int
Number of RAS2 orbitals.
- ras3int
Number of RAS3 orbitals.
- Returns
- The result of H*v.
sf_ip_ea.np_sf module¶
NumPy-based Fock-space CI.
Runs a RAS-nSF-IP/EA calculation, given the integrals in the form of NumPy arrays as input.
-
sf_ip_ea.np_sf.do_sf_np(delta_a, delta_b, ras1, ras2, ras3, Fa, Fb, tei_int, e, conf_space='', sf_opts={}, J_in=None, C_in=None)¶ Performs the RAS-SF-IP/EA calculation.
Performs a RAS-SF-IP/EA calculation. The number of spin flips and IP/EA is determined by the given parameters delta_a (number of alpha electrons removed) and delta_b (number of beta electrons added). One- and two-electron integrals are passed in as NumPy arrays. Whether to perform hole and/or particle excitations is specified using conf_space.
- Parameters
- delta_aint
Desired number of alpha electrons to remove
- delta_bint
Desired number of beta electrons to add
- ras1int
Size of RAS1 space
- ras2int
Size of RAS2 space
- ras3int
Size of RAS3 space
- Fanumpy.ndarray
Alpha Fock matrix
- Fbnumpy.ndarray
Beta Fock matrix
- tei_intnumpy.ndarray
Two-electron integrals
- efloat
ROHF energy
- conf_spacestr
Desired excitation scheme
- sf_optsdict
Additional options for spin-flip
- J_innumpy.ndarray
J matrix (for DF calculations)
- C_innumpy.ndarray
MO coefficients matrix (for DF calculations)
- Returns
- A FockWfn object with calculation information and results
sf_ip_ea.post_ci_analysis module¶
Handling post-CI analysis.
These are functions to handle post-CI analysis, primarily the S**2 values and information about the determinants.
-
sf_ip_ea.post_ci_analysis.calc_s2(vect, wfn)¶ Calculate S**2 for a given Fock CI eigenvector.
Calculates the S**2 expectation value for a given Fock CI eigenvector corresponding to a particular state.
- Parameters
- vectnumpy.ndarray
Eigenvector for which to compute S**2.
- wfnFockWfn
Reference Fock CI wavefunction object for the calculation.
- Returns
- S**2 value for the given eigenvector.
-
sf_ip_ea.post_ci_analysis.calc_sz(vect, wfn)¶ Calculate Sz for a given Fock CI eigenvector.
Calculates the Sz expectation value for a given Fock CI eigenvector corresponding to a particular state.
- Parameters
- vectnumpy.ndarray
Eigenvector for which to compute Sz.
- wfnFockWfn
Reference Fock CI wavefunction object for the calculation.
- Returns
- Sz value for the given eigenvector.
-
sf_ip_ea.post_ci_analysis.generate_dets(n_SF, delta_ec, conf_space, ras1, ras2, ras3)¶ Generates an ordered list of determinants.
This generates an ordered list of all of the determinants for the Fock space described by the input parameters.
Determinants are in the following form (indexing starts at zero):
det = [...] (det[0] is 0th determinant, det[1] is 1st, etc.) det[count] = [[[elim (alpha)], [elim (beta)]], [[add (a)], [add (b)]]]
- Parameters
- n_SFint
Number of spin-flips
- delta_ecint
Change in electron count
- conf_spacestr
Configuration space
- ras1int
Number of RAS1 orbitals
- ras2int
Number of RAS2 orbitals
- ras3int
Number of RAS3 orbitals
- Returns
- Ordered list of determinants
-
sf_ip_ea.post_ci_analysis.print_det_list(wfn)¶ Prints the full list of determinants to standard output.
This prints the full determinant list to standard output, in order. This does not print by default, but is useful for some post-CI analysis. Information is printed to standard output.
- Parameters
- wfnFockWfn
Fock CI wavefunction object with determinants to print.
-
sf_ip_ea.post_ci_analysis.print_dets(vect, wfn, dets_to_print=10)¶ Given a CI vector/state, prints information about the most important determinants.
- Parameters
- vectnumpy.ndarray
Eigenvector corresponding to the desired state.
- wfnFockWfn
Fock CI wavefunction object for the calculation.
- dets_to_printint
Number of determinants to print. Optional. Defaults to 10.
sf_ip_ea.psi4_sf module¶
Fock-space CI using Psi4’s interface.
Runs a RAS-nSF-IP/EA calculation using Psi4 to construct the integrals.
-
sf_ip_ea.psi4_sf.do_sf_psi4(delta_a, delta_b, mol, conf_space='', ref_opts={}, sf_opts={})¶ Runs RAS-nSF-IP/EA using Psi4.
This runs a RAS-nSF-IP/EA calculation using Psi4 to solve for the reference state ROHF orbitals, if needed, and construct the one- and two-electron integral objects to pass along to the NumPy-based solver.
- Parameters
- delta_aint
Number of alpha electrons to eliminate
- delta_bint
Number of beta electrons to create
- molPsi4.core.Molecule
Psi4 Molecule object on which to run the calculation
- conf_spacestr
Inclusion of holes/particles (Options: “”, “h”, “p”, “h,p”)
- ref_optsdict
Additional options for Psi4 (see Psi4 docs)
- sf_optsdict
Additional options for Fock CI
- Returns
- A FockWfn object containing calculation results
sf_ip_ea.solvers module¶
Davidson solver module.
This module holds the Davidson solver and related functions.
-
sf_ip_ea.solvers.davidson(A, vInit, e_conv=1e-08, r_conv=0.0001, vect_cutoff=1e-05, maxIter=200, collapse_per_root=20)¶ Solves for the eigenvalues and eigenvectors of a Hamiltonian.
Uses the Davidson method to solve for the eigenvalues and eigenvectors of a given Hermitian matrix A.
- Parameters
- Asf_ip_ea.linop.LinOpH
LinOp object defining our matrix-vector multiply.
- vInitnumpy.ndarray
Initial guess vector(s).
- e_convfloat
Cutoff for energy convergence. (Default: 1e-8)
- r_convfloat
Cutoff for residual squared convergence. (Default: 1e-4)
- vect_cutofffloat
Cutoff for adding vectors to Krylov space. (Default: 1e-5)
- maxIterint
Maximum number of iterations. (Default: 200)
- collapse_per_rootint
Number of vectors per root to save after collapsing the Krylov space. (Default: 20)
- Returns
- numpy.ndarray
Eigenvalues/eigenvectors for the Hamiltonian.
sf_ip_ea.tei module¶
Two-electron integral object handling.
This module handles the two-electron integrals. It generates and stores integrals in the form of NumPy arrays. The TEI subclasses handle full (TEIFullBase) and density-fitted integrals (TEIDFBase), and both of these have subclasses which handle the multiple ways of constructing the integrals (by using Psi4, by passing in pre-constructed NumPy arrays, and so on). When an interface to a new program is added, a corresponding TEI object specific to that program should be added here as a subclass of TEIFullBase or TEIDFBase.
Used Psi4NumPy Tutorials for reference for the density fitting.
-
class
sf_ip_ea.tei.TEI¶ Bases:
objectAbstract parent class for two-electron integral object handling.
This class handles the two-electron integrals. It generates and stores integrals in the form of NumPy arrays. Relevant sub-blocks can be easily accessed via the
get_subblockroutine.Methods
get_subblock(self, a, b, c, d)Returns a given subblock of the ERI matrix.
-
abstract
get_subblock(self, a, b, c, d)¶ Returns a given subblock of the ERI matrix.
-
abstract
-
class
sf_ip_ea.tei.TEIDFBase¶ Bases:
sf_ip_ea.tei.TEIBase class for constructing density-fitted two-electron integrals.
Methods
get_subblock(self, a, b, c, d)Returns a given subblock of the two-electron integrals (DF).
-
get_subblock(self, a, b, c, d)¶ Returns a given subblock of the two-electron integrals (DF).
Returns a given subblock of the DF two-electron integral object. The RAS space to return is given by
a,b,c, andd. This function performs the contraction of the given bra (B_ab) and ket (B_cd) matrices using einsum. The output is given in physicists’ notation with the form <ab|cd>.- Parameters
- aint
RAS space of index 1.
- bint
RAS space of index 2.
- cint
RAS space of index 3.
- dint
RAS space of index 4.
- Returns
- numpy.ndarray
NumPy representation of the desired subblock.
-
-
class
sf_ip_ea.tei.TEIDFNumPy(C, ras1, ras2, ras3, conf_space, np_tei, np_J)¶ Bases:
sf_ip_ea.tei.TEIDFBaseClass for constructing TEIs using NumPy arrays.
This allows for construction of DF integrals using NumPy arrays as input. B matrices are constructed in the initialization step and are subsequently multiplied to form the appropriate subblocks. B matrices are only formed for the necessary subblocks, given the excitation scheme.
- Attributes
- B11numpy.ndarray
RAS1/RAS1 B-matrix.
- B12numpy.ndarray
RAS1/RAS2 B-matrix.
- B13numpy.ndarray
RAS1/RAS3 B-matrix.
- B21numpy.ndarray
RAS2/RAS1 B-matrix.
- B22numpy.ndarray
RAS2/RAS2 B-matrix.
- B23numpy.ndarray
RAS2/RAS3 B-matrix.
- B31numpy.ndarray
RAS3/RAS1 B-matrix.
- B32numpy.ndarray
RAS3/RAS2 B-matrix.
- B33numpy.ndarray
RAS3/RAS3 B-matrix.
Methods
get_subblock(self, a, b, c, d)Returns a given subblock of the two-electron integrals (DF).
-
class
sf_ip_ea.tei.TEIDFPsi4(C, basis, aux, ras1, ras2, ras3, conf_space)¶ Bases:
sf_ip_ea.tei.TEIDFBaseMethods
get_subblock(self, a, b, c, d)Returns a given subblock of the two-electron integrals (DF).
-
class
sf_ip_ea.tei.TEIFullBase¶ Bases:
sf_ip_ea.tei.TEIBase class for constructing full two-electron integrals.
Methods
get_full(self)Returns the full set of two-electron integrals as a NumPy array.
get_subblock(self, a, b, c, d)Returns a given subblock of the two-electron integrals.
-
get_full(self)¶ Returns the full set of two-electron integrals as a NumPy array.
This returns the full two-electron integral (for the necessary RAS spaces, given the excitations) as a NumPy array. This is useful for storing the array for use in future calculations.
- Returns
- numpy.ndarray
NumPy representation of two-electron integral object.
-
get_subblock(self, a, b, c, d)¶ Returns a given subblock of the two-electron integrals.
This returns a NumPy representation of a given subblock of the full two-electron integrals. The subblock to return is given by parameters
a,b,c, andd, where each indicates a RAS space (1, 2, or 3). The returned integral is in physicists’ notation and has the form <ab|cd>.So to get the block with a and c in RAS1 and b and d in RAS2, one would use get_subblock(1, 2, 1, 2).
- Parameters
- aint
RAS space of index 1.
- bint
RAS space of index 2.
- cint
RAS space of index 3.
- dint
RAS space of index 4.
- Returns
- numpy.ndarray
NumPy representation of the desired subblock.
-
-
class
sf_ip_ea.tei.TEIFullNumPy(ras1, ras2, ras3, conf_space, np_tei)¶ Bases:
sf_ip_ea.tei.TEIFullBaseClass for constructing full TEI integrals from a NumPy array.
This class stores the integrals as a NumPy array, given a NumPy array. Note that the relevant subset of the array should be given. So, in the case of a RAS(h) calculation, one would give the TEI constructed in the basis of RAS1 and RAS2 orbitals only (not RAS3).
- Attributes
- erinumpy.ndarray
Numpy representation of the relevant two-electron integrals.
- indlist
A list of index ranges for the subblocks (RAS1, RAS2, RAS3).
Methods
get_full(self)Returns the full set of two-electron integrals as a NumPy array.
get_subblock(self, a, b, c, d)Returns a given subblock of the two-electron integrals.
-
class
sf_ip_ea.tei.TEIFullPsi4(C, basis, ras1, ras2, ras3, conf_space)¶ Bases:
sf_ip_ea.tei.TEIFullBaseClass for constructing full TEI integrals using Psi4.
This class constructs the full two-electron integrals using Psi4. It then stores the integrals as NumPy arrays.
- Attributes
- erinumpy.ndarray
Numpy representation of the relevant two-electron integrals.
- indlist
A list of index ranges for the subblocks (RAS1, RAS2, RAS3).
Methods
get_full(self)Returns the full set of two-electron integrals as a NumPy array.
get_subblock(self, a, b, c, d)Returns a given subblock of the two-electron integrals.
Module contents¶
A program for running RAS-SF-IP/EA calculations.
This program runs RAS-SF-IP/EA calculations using an efficient
tensor-contraction-based scheme. The contractions have been hand-derived and
use NumPy’s einsum for efficiency. The program is run primarily through
the main fock_ci function.
-
sf_ip_ea.fock_ci(delta_a, delta_b, mol, conf_space='', ref_opts={}, sf_opts={}, program='PSI4')¶ Performs Fock-space CI (SF-IP/EA).
This is the main function call for the program. Given the changes in alpha and beta electron counts, it performs the correct number of spin-flips and IP/EAs.
- Parameters
- delta_aint
Desired number of alpha electrons to remove.
- delta_bint
Desired number of beta electrons to add.
- molMolecule
The molecule object to run the calculation on. This should be built in whichever program you’ll use to run the reference, and should be handled properly by the reference program.
- conf_spacestring
- Desired configuration space/additional excitations.
""CAS"h"1 hole excitation"p"1 particle excitation"h,p"1 hole + 1 particle excitation
- ref_optsdict
Options for the reference program. See relevant reference program docs (ex. Psi4 docs) for details.
- sf_optsdict
- Additional options for stand-alone SF code.
sf_diag_method: Diagonalization method to use.RSPDirect (deprecated)LANCZOSUse NumPy’s LanczosDAVIDSONUse our Davidson
num_roots: Number of roots to solve for.guess_type: Type of guess vector to useCASDo CAS first and use that as an initial guess.RANDOMRandom orthonormal basisREADRead guess from a NumPy file (TODO)
integral_type: Which integrals to use (DF or FULL)FULLUse full integrals (no density fitting)DFUse density fit integrals
- Returns
- sf_wfn
Wavefunction object (sf_wfn) containing calculation data and results