pymolresponse package
Subpackages
Submodules
pymolresponse.ao2mo module
Tools for performing AO-to-MO transformations of two-electron integrals.
- class pymolresponse.ao2mo.AO2MO(C: ndarray, occupations: Sequence[int], verbose: int = 1, I: ndarray | None = None)[source]
Bases:
object
Interface for performing AO-to-MO tranformations of two-electron integrals.
- perform_rhf_full() None [source]
Perform the transformation \((\mu\nu|\lambda\sigma) \rightarrow (pq|rs)\).
- perform_rhf_partial() None [source]
Perform the transformation \((\mu\nu|\lambda\sigma) \rightarrow (ia|jb), (ij|ab)\).
- perform_uhf_full() None [source]
Perform the transformation \((\mu\nu|\lambda\sigma) \rightarrow (p^{\alpha}q^{\alpha}|r^{\alpha}s^{\alpha}), (p^{\alpha}q^{\alpha}|r^{\beta}s^{\beta}), (p^{\beta}q^{\beta}|r^{\alpha}s^{\alpha}), (p^{\beta}q^{\beta}|r^{\beta}s^{\beta})\).
- perform_uhf_partial() None [source]
Perform the transformation \((\mu\nu|\lambda\sigma) \rightarrow (i^{\alpha}a^{\alpha}|j^{\alpha}b^{\alpha}), (i^{\alpha}a^{\alpha}|j^{\beta}b^{\beta}), (i^{\beta}a^{\beta}|j^{\alpha}b^{\alpha}), (i^{\beta}a^{\beta}|j^{\beta}b^{\beta}), (i^{\alpha}j^{\alpha}|a^{\alpha}b^{\alpha}), (i^{\beta}j^{\beta}|a^{\beta}b^{\beta})\).
pymolresponse.constants module
Storage of common physical constants. Uses scipy.constants: https://docs.scipy.org/doc/scipy/reference/constants.html.
pymolresponse.core module
Core types.
- class pymolresponse.core.AO2MOTransformationType(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)[source]
Bases:
Enum
For routines that perform AO-to-MO transformations, specify the kind of transformation.
Starting from \((\mu\nu|\lambda\sigma)\),
a ‘partial’ transformation will produce \((ia|jb)\) and \((ij|ab)\) if needed
a ‘full’ transformation will produce \((pq|rs)\)
- full = 'full'
- partial = 'partial'
- class pymolresponse.core.Hamiltonian(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)[source]
Bases:
Enum
Specify which approximation for the orbital Hessian should be used.
RPA indicates the random phase approximation: the orbital Hessian is exact for the first-order polarization propagator.
TDA indicates the Tamm-Dancoff approximation: the B matrix is set to zero.
Currently this is only capable of handling the usual first-order polarization propagator for CPHF and TDHF.
- RPA = 'rpa'
- TDA = 'tda'
- class pymolresponse.core.Program(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)[source]
Bases:
Enum
Specify which program should be used for a particular step.
- DALTON = 'dalton'
- Psi4 = 'psi4'
- PySCF = 'pyscf'
- class pymolresponse.core.Spin(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)[source]
Bases:
Enum
Specify whether the calculation should conserve spin (singlet) or flip spin (triplet).
For response properties, this affects both the orbital Hessian and the perturbing operator. Note that these two specifications are currently independent!
For excitation energies, this only affects the orbital Hessian and determines whether singlet or triplet excitation energies are calculated.
- singlet = 0
- triplet = 1
pymolresponse.cphf module
Driver for solving the coupled perturbed Hartree-Fock (CPHF) equations.
- class pymolresponse.cphf.CPHF(solver: Solver)[source]
Bases:
Driver
Driver for solving the coupled perturbed Hartree-Fock (CPHF) equations.
- add_operator(operator: Operator) None [source]
Add an operator to the list of operators that will be used as the right-hand side perturbation.
- run(hamiltonian: Hamiltonian, spin: Spin, program: Program, program_obj) None [source]
pymolresponse.explicit_equations_full module
Explicit equations for orbital Hessian terms using fully-transformed MO-basis two-electron integrals, e.g., \((pq|rs)\).
- pymolresponse.explicit_equations_full.form_rpa_a_matrix_mo_singlet_full(E_MO: ndarray, TEI_MO: ndarray, nocc: int) ndarray [source]
Form the A (CIS) matrix in the MO basis. [singlet]
The equation for element \(\{ia,jb\}\) is \(\left<aj||ib\right> = \left<aj|ib\right> - \left<aj|bi\right> = [ai|jb] - [ab|ji] = 2(ai|jb) - (ab|ji)\). It also includes the virt-occ energy difference on the diagonal.
- pymolresponse.explicit_equations_full.form_rpa_a_matrix_mo_singlet_os_full(TEI_MO_xxyy: ndarray, nocc_x: int, nocc_y: int) ndarray [source]
Form the opposite-spin part of the A (CIS) matrix in the MO basis. [singlet]
The equation for element \(\{ia,jb\}\) is \(????\).
- pymolresponse.explicit_equations_full.form_rpa_a_matrix_mo_singlet_ss_full(E_MO: ndarray, TEI_MO: ndarray, nocc: int) ndarray [source]
Form the same-spin part of the A (CIS) matrix in the MO basis. [singlet]
The equation for element \(\{ia,jb\}\) is \(????\).
- pymolresponse.explicit_equations_full.form_rpa_a_matrix_mo_triplet_full(E_MO: ndarray, TEI_MO: ndarray, nocc: int) ndarray [source]
Form the A (CIS) matrix in the MO basis. [triplet]
The equation for element \(\{ia,jb\}\) is \(- \left<aj|bi\right> = - [ab|ji] = - (ab|ji)\). It also includes the virt-occ energy difference on the diagonal.
- pymolresponse.explicit_equations_full.form_rpa_b_matrix_mo_singlet_full(TEI_MO: ndarray, nocc: int) ndarray [source]
Form the B matrix for RPA in the MO basis. [singlet]
The equation for element \(\{ia,jb\}\) is \(\left<ab||ij\right> = \left<ab|ij\right> - \left<ab|ji\right> = [ai|bj] - [aj|bi] = 2(ai|bj) - (aj|bi)\).
- pymolresponse.explicit_equations_full.form_rpa_b_matrix_mo_singlet_os_full(TEI_MO_xxyy: ndarray, nocc_x: int, nocc_y: int) ndarray [source]
Form the opposite-spin part of the RPA B matrix in the MO basis. [singlet]
The equation for element \(\{ia,jb\}\) is \(????\).
pymolresponse.explicit_equations_partial module
Explicit equations for orbital Hessian terms using partially-transformed MO-basis two-electron integrals, e.g., \((ia|jb), (ij|ab)\).
- pymolresponse.explicit_equations_partial.form_rpa_a_matrix_mo_singlet_os_partial(TEI_MO_iajb_xxyy: ndarray) ndarray [source]
Form the opposite-spin part of the A (CIS) matrix in the MO basis. [singlet]
The equation for element \(\{ia,jb\}\) is \(????\).
- pymolresponse.explicit_equations_partial.form_rpa_a_matrix_mo_singlet_partial(E_MO: ndarray, TEI_MO_iajb: ndarray, TEI_MO_ijab: ndarray) ndarray [source]
Form the A (CIS) matrix in the MO basis. [singlet]
The equation for element \(\{ia,jb\}\) is \(\left<aj||ib\right> = \left<aj|ib\right> - \left<aj|bi\right> = [ai|jb] - [ab|ji] = 2(ai|jb) - (ab|ji)\). It also includes the virt-occ energy difference on the diagonal.
- pymolresponse.explicit_equations_partial.form_rpa_a_matrix_mo_singlet_ss_partial(E_MO: ndarray, TEI_MO_iajb: ndarray, TEI_MO_ijab: ndarray) ndarray [source]
Form the same-spin part of the A (CIS) matrix in the MO basis. [singlet]
The equation for element \(\{ia,jb\}\) is \(????\).
- pymolresponse.explicit_equations_partial.form_rpa_a_matrix_mo_triplet_partial(E_MO: ndarray, TEI_MO_ijab: ndarray) ndarray [source]
Form the A (CIS) matrix in the MO basis. [triplet]
The equation for element \(\{ia,jb\}\) is \(- \left<aj|bi\right> = - [ab|ji] = - (ab|ji)\). It also includes the virt-occ energy difference on the diagonal.
- pymolresponse.explicit_equations_partial.form_rpa_b_matrix_mo_singlet_os_partial(TEI_MO_iajb_xxyy: ndarray) ndarray [source]
Form the opposite-spin part of the RPA B matrix in the MO basis. [singlet]
The equation for element \(\{ia,jb\}\) is \(????\).
- pymolresponse.explicit_equations_partial.form_rpa_b_matrix_mo_singlet_partial(TEI_MO_iajb: ndarray) ndarray [source]
Form the B matrix for RPA in the MO basis. [singlet]
The equation for element \(\{ia,jb\}\) is \(\left<ab||ij\right> = \left<ab|ij\right> - \left<ab|ji\right> = [ai|bj] - [aj|bi] = 2(ai|bj) - (aj|bi)\).
pymolresponse.helpers module
Utility functions that are core to calculating physical values.
- pymolresponse.helpers.calc_center_of_nuclear_charge(coords: ndarray, charges: ndarray) ndarray [source]
- pymolresponse.helpers.get_uhf_values(mat_uhf_a: ndarray, mat_uhf_b: ndarray, pair_rohf: Tuple[int, int], nocc_a: int, nvirt_a: int, nocc_b: int, nvirt_b: int) List[float] [source]
For a pair ROHF 1-based indices, find the corresponing alpha- and beta-spin UHF values.
pymolresponse.integrals module
- class pymolresponse.integrals.IntegralLabel(label: str, comp: int | None = None, symmetry: IntegralSymmetry = IntegralSymmetry.SYMMETRIC)[source]
Bases:
object
- comp: int | None
- label: str
- symmetry: IntegralSymmetry
- class pymolresponse.integrals.IntegralSymmetry(value, names=None, *, module=None, qualname=None, type=None, start=1, boundary=None)[source]
Bases:
Enum
- ANTISYMMETRIC = 2
- SYMMETRIC = 1
- class pymolresponse.integrals.Integrals[source]
Bases:
ABC
- integrals(label: IntegralLabel) ndarray [source]
- class pymolresponse.integrals.JK[source]
Bases:
ABC
Interface for computing Coulomb and exchange integrals pre-contracted with the density.
pymolresponse.molecular_property module
- class pymolresponse.molecular_property.MolecularProperty(program: Program, program_obj, driver: Driver, mocoeffs: ndarray, moenergies: ndarray, occupations: ndarray)[source]
Bases:
ABC
A molecular property that is calculated from one or more operators.
The
- run(hamiltonian: Hamiltonian, spin: Spin) None [source]
- class pymolresponse.molecular_property.ResponseProperty(program: Program, program_obj, driver: CPHF, mocoeffs: ndarray, moenergies: ndarray, occupations: ndarray, *, frequencies: Sequence[float])[source]
Bases:
MolecularProperty
,ABC
A molecular property that is calculated as the response to one or more perturbations, each represented by an operator.
- class pymolresponse.molecular_property.TransitionProperty(program: Program, program_obj, driver: TDHF, mocoeffs: ndarray, moenergies: ndarray, occupations: ndarray, *, do_tda: bool = False)[source]
Bases:
MolecularProperty
,ABC
pymolresponse.operators module
- class pymolresponse.operators.Operator(label: str = '', is_imaginary: bool = False, is_spin_dependent: bool = False, triplet: bool = False, slice_idx: int = -1, ao_integrals: ndarray | None = None)[source]
Bases:
object
Handle property integrals, taking them from the AO basis to a representation of a right-hand side perturbation for CPHF or transition properties.
pymolresponse.solvers module
- class pymolresponse.solvers.EigSolver(mocoeffs: ndarray, moenergies: ndarray, occupations: ndarray)[source]
Bases:
Solver
,ABC
Base class for all solvers of the type $AX = 0$ (eigensolvers).
- class pymolresponse.solvers.EigSolverTDA(mocoeffs: ndarray, moenergies: ndarray, occupations: ndarray)[source]
Bases:
EigSolver
,ABC
Base class for eigensolvers that use the Tamm-Dancoff approximation (TDA) rather than the random phase approximation (RPA).
- class pymolresponse.solvers.ExactDiagonalizationSolver(mocoeffs: ndarray, moenergies: ndarray, occupations: ndarray)[source]
Bases:
EigSolver
Get the eigenvectors and eigenvalues of the RPA Hamiltonian by forming the electronic Hessian explicitly and diagonalizing it exactly.
- form_explicit_hessian(hamiltonian: Hamiltonian, spin: Spin, frequency: float) None [source]
- run(hamiltonian: Hamiltonian, spin: Spin, program: Program, program_obj) None [source]
Run the solver.
- class pymolresponse.solvers.ExactDiagonalizationSolverTDA(mocoeffs: ndarray, moenergies: ndarray, occupations: ndarray)[source]
Bases:
ExactDiagonalizationSolver
,EigSolverTDA
Get the eigenvectors and eigenvalues of the TDA Hamiltonian by forming the electronic Hessian explicitly and diagonalizing it exactly.
- form_explicit_hessian(hamiltonian: Hamiltonian, spin: Spin, frequency: float) None [source]
- class pymolresponse.solvers.ExactInv(mocoeffs: ndarray, moenergies: ndarray, occupations: ndarray, *, inv_func: Callable[[ndarray], ndarray] | None = None)[source]
Bases:
ExactLineqSolver
- class pymolresponse.solvers.ExactInvCholesky(mocoeffs: ndarray, moenergies: ndarray, occupations: ndarray)[source]
Bases:
ExactLineqSolver
- class pymolresponse.solvers.ExactLineqSolver(mocoeffs: ndarray, moenergies: ndarray, occupations: ndarray)[source]
Bases:
LineqSolver
,ABC
Base class for all solvers of the type $AX = B$, where $A$ is formed explicitly.
- form_explicit_hessian(hamiltonian: Hamiltonian, spin: Spin, frequency: float) None [source]
- abstract invert_explicit_hessian() None [source]
Invert the full explicitly-constructed electronic Hessian.
- run(hamiltonian: Hamiltonian, spin: Spin, program: Program, program_obj) None [source]
Run the solver.
- class pymolresponse.solvers.IterativeLinEqSolver(mocoeffs: ndarray, moenergies: ndarray, occupations: ndarray, jk_generator: JK, *, maxiter: int = 40, conv: float = 1e-09)[source]
Bases:
LineqSolver
- run(hamiltonian: Hamiltonian, spin: Spin, program: Program, program_obj) None [source]
Run the solver.
- class pymolresponse.solvers.LineqSolver(mocoeffs: ndarray, moenergies: ndarray, occupations: ndarray)[source]
Bases:
Solver
,ABC
Base class for all solvers of the type $AX = B$.
- class pymolresponse.solvers.Solver(mocoeffs: ndarray, moenergies: ndarray, occupations: ndarray)[source]
Bases:
ABC
A Solver does all of the TODO
There are two principal algorithm classes for forming response vectors:
Exact, where the full electronic Hessian is formed explicitly, diagonalized, and then contracted with property gradient vectors to response vectors, or
Iterative, where the electronic Hessian is never formed explicitly, TODO
- form_tei_mo(program: Program, program_obj, tei_mo_type: AO2MOTransformationType) None [source]
- abstract run(hamiltonian: Hamiltonian, spin: Spin, program: Program, program_obj) None [source]
Run the solver.
pymolresponse.td module
Drivers for solving the time-dependent Hartree-Fock (TDHF) equations.
- class pymolresponse.td.TDA(solver: EigSolverTDA)[source]
Bases:
TDHF
Driver for solving the time-dependent Hartree-Fock equations with the Tamm-Dancoff approximation (TDA), also called the configuration interaction with single excitation (CIS) equations when no XC contribution is present.
pymolresponse.utils module
Utility functions that are not core to calculating physical values.
- class pymolresponse.utils.Splitter(widths)[source]
Bases:
object
Split a line based not on a character, but a given number of field widths.
- pymolresponse.utils.fix_moenergies_shape(moenergies: Tuple[ndarray, ...] | ndarray) ndarray [source]
- pymolresponse.utils.flip_triangle_sign(A: ndarray, triangle: str = 'lower') ndarray [source]
Flip the sign of either the lower or upper triangle of a sqare matrix. Assume nothing about its symmetry.
- Parameters:
- Anp.ndarray
- triangle{‘lower’, ‘upper’}
- Returns:
- np.ndarray
- pymolresponse.utils.form_first_hyperpolarizability_averages(beta: ndarray) Tuple[ndarray, ndarray] [source]
- pymolresponse.utils.form_vec_energy_differences(moene_occ: ndarray, moene_virt: ndarray) ndarray [source]
- pymolresponse.utils.get_reference_value_from_file(filename: Path | str, hamiltonian: str, spin: str, frequency: str, label_1: str, label_2: str) float [source]
- pymolresponse.utils.matsym(amat: ndarray, thrzer: float = 1e-14) int [source]
Copied from
DALTON/gp/gphjj.F/MATSYM
.thrzer taken from
DALTON/include/thrzer.h
.
- Parameters:
- amatnp.ndarray
- thrzerfloat
Threshold below which a (floating point) number is considered zero.
- Returns:
- int
1 if the matrix is symmetric to threshold thrzer
2 if the matrix is antisymmetric to threshold thrzer
3 if all elements are below thrzer
0 otherwise (the matrix is unsymmetric about the diagonal)
- pymolresponse.utils.screen(mat: ndarray, thresh: float = 1e-16) ndarray [source]
Set all values smaller than the given threshold to zero (considering them as numerical noise).
- Parameters:
- matnp.ndarray
- threshfloat
Threshold below which all elements of mat smaller than thresh are set to zero.
- Returns:
- np.ndarray
Module contents
Molecular frequency-dependent response properties for arbitrary operators