pyscf.tdscf package

Submodules

pyscf.tdscf.common_slow module

This and other _slow modules implement the time-dependent procedure. The primary performance drawback is that, unlike other ‘fast’ routines with an implicit construction of the eigenvalue problem, these modules construct TDHF matrices explicitly. As a result, regular numpy.linalg.eig can be used to retrieve TDHF roots in a reliable fashion without any issues related to the Davidson procedure.

This is a helper module defining basic interfaces.

class pyscf.tdscf.common_slow.MolecularMFMixin(model, frozen=None)[source]

Bases: object

property mo_coeff

MO coefficients.

property mo_coeff_full

MO coefficients.

property mo_energy

MO energies.

property mo_occ

MO occupation numbers.

property nmo

The total number of molecular orbitals.

property nmo_full

The true (including frozen degrees of freedom) total number of molecular orbitals.

property nocc

The number of occupied orbitals.

property nocc_full

The true (including frozen degrees of freedom) number of occupied orbitals.

squeeze(x)[source]

Squeezes quantities in the case of a PBC model.

class pyscf.tdscf.common_slow.PeriodicMFMixin(model, frozen=None)[source]

Bases: object

property mo_coeff

MO coefficients.

property mo_coeff_full

MO coefficients.

property mo_energy

MO energies.

property mo_occ

MO occupation numbers.

property nmo

The total number of molecular orbitals.

property nmo_full

The true (including frozen degrees of freedom) total number of molecular orbitals.

property nocc

The number of occupied orbitals.

property nocc_full

The true (including frozen degrees of freedom) number of occupied orbitals.

class pyscf.tdscf.common_slow.TDBase(mf, frozen=None)[source]

Bases: object

ao2mo()[source]

Picks ERI: either 4-fold or 8-fold symmetric.

Returns:

A suitable ERI.

kernel()[source]

Calculates eigenstates and eigenvalues of the TDHF problem.

Returns:

Positive eigenvalues and eigenvectors.

v2a = None
vector_to_amplitudes(vectors)[source]

Transforms (reshapes) and normalizes vectors into amplitudes.

Args:

vectors (numpy.ndarray): raw eigenvectors to transform;

Returns:

Amplitudes with the following shape: (# of roots, 2 (x or y), # of occupied orbitals, # of virtual orbitals).

class pyscf.tdscf.common_slow.TDERIMatrixBlocks[source]

Bases: TDMatrixBlocks

eri_mknj(item, *args)[source]

Retrieves ERI block using ‘mknj’ notation.

Args:

item (str): a 4-character string of ‘mknj’ letters; *args: other arguments passed to get_block_ov_notation;

Returns:

The corresponding block of ERI (matrix with paired dimensions).

eri_ov(item, *args)[source]

Retrieves ERI block using ‘ov’ notation.

Args:

item (str): a 4-character string of ‘o’ and ‘v’ letters; *args: other args passed to __calc_block__;

Returns:

The corresponding block of ERI (4-tensor, phys notation).

symmetries = [((0, 1, 2, 3), False)]
tdhf_diag(*args)[source]

Retrieves the diagonal block.

Args:

*args: args passed to __get_mo_energies__;

Returns:

The diagonal block.

tdhf_primary_form(*args, **kwargs)[source]

A primary form of TDHF matrixes (AB).

Returns:

Output type: “ab”, and the corresponding matrixes.

class pyscf.tdscf.common_slow.TDMatrixBlocks[source]

Bases: object

tdhf_ab_form(*args, **kwargs)[source]

The A-B form of the TD problem.

Returns:

A and B TD matrices.

tdhf_full_form(*args, **kwargs)[source]

The full form of the TD problem.

Returns:

The full TD matrix.

tdhf_mk_form(*args, **kwargs)[source]

The MK form of the TD problem.

Returns:

MK and K TD matrixes.

tdhf_primary_form(*args, **kwargs)[source]

A primary form of TDHF matrixes.

Returns:

Output type: “full”, “ab”, or “mk” and the corresponding matrix(es).

class pyscf.tdscf.common_slow.TDProxyMatrixBlocks(model)[source]

Bases: TDMatrixBlocks

tdhf_primary_form(*args, **kwargs)[source]

A primary form of TDHF matrixes.

Returns:

Output type: “full”, “ab”, or “mk” and the corresponding matrix(es).

class pyscf.tdscf.common_slow.VindTracker(vind)[source]

Bases: object

property elements_calc
property elements_total
property msize
property ncalls
property ratio
reset()[source]

Resets statistics.

text_stats()[source]
pyscf.tdscf.common_slow.ab2full(a, b)[source]

Transforms A and B TD matrices into a full matrix.

Args:

a (numpy.ndarray): TD A-matrix; b (numpy.ndarray): TD B-matrix;

Returns:

The full TD matrix.

pyscf.tdscf.common_slow.ab2mkk(a, b, tolerance=1e-12)[source]

Transforms A and B TD matrices into MK and K matrices.

Args:

a (numpy.ndarray): TD A-matrix; b (numpy.ndarray): TD B-matrix; tolerance (float): a tolerance for checking whether the input matrices are real;

Returns:

MK and K submatrices.

pyscf.tdscf.common_slow.eig(m, driver=None, nroots=None, half=True)[source]

Eigenvalue problem solver.

Args:

m (numpy.ndarray): the matrix to diagonalize; driver (str): one of the drivers; nroots (int): the number of roots ot calculate (ignored for driver == ‘eig’); half (bool): if True, implies spectrum symmetry and takes only a half of eigenvalues;

Returns:

pyscf.tdscf.common_slow.format_frozen_k(frozen, nmo, nk)[source]

Formats the argument into a mask array of bools where False values correspond to frozen orbitals for each k-point.

Args:

frozen (int, Iterable): the number of frozen valence orbitals or the list of frozen orbitals for all k-points or multiple lists of frozen orbitals for each k-point; nmo (int): the total number of molecular orbitals; nk (int): the total number of k-points;

Returns:

The mask array.

pyscf.tdscf.common_slow.format_frozen_mol(frozen, nmo)[source]

Formats the argument into a mask array of bools where False values correspond to frozen molecular orbitals.

Args:

frozen (int, Iterable): the number of frozen valence orbitals or the list of frozen orbitals; nmo (int): the total number of molecular orbitals;

Returns:

The mask array.

pyscf.tdscf.common_slow.format_mask(x)[source]

Formats a mask into a readable string.

Args:

x (ndarray): an array with the mask;

Returns:

A readable string with the mask.

pyscf.tdscf.common_slow.full2ab(full, tolerance=1e-12)[source]

Transforms a full TD matrix into A and B parts.

Args:

full (numpy.ndarray): the full TD matrix; tolerance (float): a tolerance for checking whether the full matrix is in the ABBA-form;

Returns:

A and B submatrices.

pyscf.tdscf.common_slow.full2mkk(full)[source]

Transforms a full TD matrix into MK and K parts.

Args:

full (numpy.ndarray): the full TD matrix;

Returns:

MK and K submatrices.

pyscf.tdscf.common_slow.k_nmo(model)[source]

Retrieves number of AOs per k-point.

Args:

model (RHF): the model;

Returns:

Numbers of AOs in the model.

pyscf.tdscf.common_slow.k_nocc(model)[source]

Retrieves occupation numbers.

Args:

model (RHF): the model;

Returns:

Numbers of occupied orbitals in the model.

pyscf.tdscf.common_slow.kernel(eri, driver=None, fast=True, nroots=None, **kwargs)[source]

Calculates eigenstates and eigenvalues of the TDHF problem.

Args:

eri (TDDFTMatrixBlocks): ERI; driver (str): one of the eigenvalue problem drivers; fast (bool): whether to run diagonalization on smaller matrixes; nroots (int): the number of roots to calculate; **kwargs: arguments to eri.tdhf_matrix;

Returns:

Positive eigenvalues and eigenvectors.

pyscf.tdscf.common_slow.mkk2ab(mk, k)[source]

Transforms MK and M TD matrices into A and B matrices.

Args:

mk (numpy.ndarray): TD MK-matrix; k (numpy.ndarray): TD K-matrix;

Returns:

A and B submatrices.

pyscf.tdscf.common_slow.mkk2full(mk, k)[source]

Transforms MK and M TD matrices into a full TD matrix.

Args:

mk (numpy.ndarray): TD MK-matrix; k (numpy.ndarray): TD K-matrix;

Returns:

The full TD matrix.

pyscf.tdscf.common_slow.mknj2i(item)[source]

Transforms “mknj” notation into tensor index order for the ERI.

Args:

item (str): an arbitrary transpose of “mknj” letters;

Returns:

4 indexes.

pyscf.tdscf.common_slow.msize(m)[source]

Checks whether the matrix is square and returns its size.

Args:

m (numpy.ndarray): the matrix to measure;

Returns:

An integer with the size.

pyscf.tdscf.dhf module

TDA and TDHF for no-pair DKS Hamiltonian

class pyscf.tdscf.dhf.TDA(mf)[source]

Bases: TDBase

Tamm-Dancoff approximation

Attributes:
conv_tolfloat

Diagonalization convergence tolerance. Default is 1e-9.

nstatesint

Number of TD states to be computed. Default is 3.

Saved results:

convergedbool

Diagonalization converged or not

e1D array

excitation energy for each excited state.

xyA list of two 2D arrays

The two 2D arrays are Excitation coefficients X (shape [nocc,nvir]) and de-excitation coefficients Y (shape [nocc,nvir]) for each excited state. (X,Y) are normalized to 1/2 in RHF/RKS methods and normalized to 1 for UHF/UKS methods. In the TDA calculation, Y = 0.

gen_vind(mf=None)[source]

Generate function to compute Ax

init_guess(mf, nstates=None, wfnsym=None)[source]
kernel(x0=None, nstates=None)[source]

TDA diagonalization solver

class pyscf.tdscf.dhf.TDBase(mf)[source]

Bases: TDBase

analyze(verbose=None)
get_ab(mf=None)[source]

A and B matrices for TDDFT response function.

A[i,a,j,b] = delta_{ab}delta_{ij}(E_a - E_i) + (ia||bj) B[i,a,j,b] = (ia||jb)

get_nto(state=1, threshold=0.3, verbose=None)

Natural transition orbital analysis.

The natural transition density matrix between ground state and excited state \(Tia = \langle \Psi_{ex} | i a^\dagger | \Psi_0 \rangle\) can be transformed to diagonal form through SVD \(T = O \sqrt{\lambda} V^\dagger\). O and V are occupied and virtual natural transition orbitals. The diagonal elements \(\lambda\) are the weights of the occupied-virtual orbital pair in the excitation.

Ref: Martin, R. L., JCP, 118, 4775-4777

Note in the TDHF/TDDFT calculations, the excitation part (X) is interpreted as the CIS coefficients and normalized to 1. The de-excitation part (Y) is ignored.

Args:

tdobj : TDA, or TDHF, or TDDFT object

stateint

Excited state ID. state = 1 means the first excited state. If state < 0, state ID is counted from the last excited state.

Kwargs:
thresholdfloat

Above which the NTO coefficients will be printed in the output.

Returns:

A list (weights, NTOs). NTOs are natural orbitals represented in AO basis. The first N_occ NTOs are occupied NTOs and the rest are virtual NTOs.

nuc_grad_method()[source]
class pyscf.tdscf.dhf.TDHF(mf)[source]

Bases: TDBase

gen_vind(mf=None)[source]

Generate function to compute

[ A B ][X] [-B* -A*][Y]

init_guess(mf, nstates=None, wfnsym=None)[source]
kernel(x0=None, nstates=None)[source]

TDHF diagonalization with non-Hermitian eigenvalue solver

pyscf.tdscf.dhf.analyze(tdobj, verbose=None)[source]
pyscf.tdscf.dhf.gen_tda_hop(mf, fock_ao=None)

A x

pyscf.tdscf.dhf.gen_tda_operation(mf, fock_ao=None)[source]

A x

pyscf.tdscf.dhf.gen_tdhf_operation(mf, fock_ao=None)[source]

Generate function to compute

[ A B ][X] [-B* -A*][Y]

pyscf.tdscf.dhf.get_ab(mf, mo_energy=None, mo_coeff=None, mo_occ=None)[source]

A and B matrices for TDDFT response function.

A[i,a,j,b] = delta_{ab}delta_{ij}(E_a - E_i) + (ia||bj) B[i,a,j,b] = (ia||jb)

pyscf.tdscf.dhf.get_nto(tdobj, state=1, threshold=0.3, verbose=None)[source]

pyscf.tdscf.dks module

pyscf.tdscf.dks.RPA

alias of TDDFT

class pyscf.tdscf.dks.TDA(mf)[source]

Bases: TDA

class pyscf.tdscf.dks.TDDFT(mf)[source]

Bases: TDHF

pyscf.tdscf.dks.TDDKS

alias of TDDFT

pyscf.tdscf.ghf module

pyscf.tdscf.ghf.CIS

alias of TDA

pyscf.tdscf.ghf.RPA

alias of TDHF

class pyscf.tdscf.ghf.TDA(mf)[source]

Bases: TDBase

Tamm-Dancoff approximation

Attributes:
conv_tolfloat

Diagonalization convergence tolerance. Default is 1e-9.

nstatesint

Number of TD states to be computed. Default is 3.

Saved results:

convergedbool

Diagonalization converged or not

e1D array

excitation energy for each excited state.

xyA list of two 2D arrays

The two 2D arrays are Excitation coefficients X (shape [nocc,nvir]) and de-excitation coefficients Y (shape [nocc,nvir]) for each excited state. (X,Y) are normalized to 1/2 in RHF/RKS methods and normalized to 1 for UHF/UKS methods. In the TDA calculation, Y = 0.

gen_vind(mf=None)[source]

Generate function to compute Ax

init_guess(mf, nstates=None, wfnsym=None)[source]
kernel(x0=None, nstates=None)[source]

TDA diagonalization solver

singlet = None
class pyscf.tdscf.ghf.TDBase(mf)[source]

Bases: TDBase

analyze(verbose=None)
get_ab(mf=None)[source]

A and B matrices for TDDFT response function.

A[i,a,j,b] = delta_{ab}delta_{ij}(E_a - E_i) + (ia||bj) B[i,a,j,b] = (ia||jb)

get_nto(state=1, threshold=0.3, verbose=None)

Natural transition orbital analysis.

The natural transition density matrix between ground state and excited state \(Tia = \langle \Psi_{ex} | i a^\dagger | \Psi_0 \rangle\) can be transformed to diagonal form through SVD \(T = O \sqrt{\lambda} V^\dagger\). O and V are occupied and virtual natural transition orbitals. The diagonal elements \(\lambda\) are the weights of the occupied-virtual orbital pair in the excitation.

Ref: Martin, R. L., JCP, 118, 4775-4777

Note in the TDHF/TDDFT calculations, the excitation part (X) is interpreted as the CIS coefficients and normalized to 1. The de-excitation part (Y) is ignored.

Args:

tdobj : TDA, or TDHF, or TDDFT object

stateint

Excited state ID. state = 1 means the first excited state. If state < 0, state ID is counted from the last excited state.

Kwargs:
thresholdfloat

Above which the NTO coefficients will be printed in the output.

Returns:

A list (weights, NTOs). NTOs are natural orbitals represented in AO basis. The first N_occ NTOs are occupied NTOs and the rest are virtual NTOs.

nuc_grad_method()[source]
pyscf.tdscf.ghf.TDGHF

alias of TDHF

class pyscf.tdscf.ghf.TDHF(mf)[source]

Bases: TDBase

gen_vind(mf=None)[source]

Generate function to compute

[ A B ][X] [-B* -A*][Y]

init_guess(mf, nstates=None, wfnsym=None)[source]
kernel(x0=None, nstates=None)[source]

TDHF diagonalization with non-Hermitian eigenvalue solver

singlet = None
pyscf.tdscf.ghf.analyze(tdobj, verbose=None)[source]
pyscf.tdscf.ghf.gen_tda_hop(mf, fock_ao=None, wfnsym=None)

A x

Kwargs:
wfnsymint or str

Point group symmetry irrep symbol or ID for excited CIS wavefunction.

pyscf.tdscf.ghf.gen_tda_operation(mf, fock_ao=None, wfnsym=None)[source]

A x

Kwargs:
wfnsymint or str

Point group symmetry irrep symbol or ID for excited CIS wavefunction.

pyscf.tdscf.ghf.gen_tdhf_operation(mf, fock_ao=None, wfnsym=None)[source]

Generate function to compute

[ A B ][X] [-B* -A*][Y]

pyscf.tdscf.ghf.get_ab(mf, mo_energy=None, mo_coeff=None, mo_occ=None)[source]

A and B matrices for TDDFT response function.

A[i,a,j,b] = delta_{ab}delta_{ij}(E_a - E_i) + (ia||bj) B[i,a,j,b] = (ia||jb)

pyscf.tdscf.ghf.get_nto(tdobj, state=1, threshold=0.3, verbose=None)[source]

pyscf.tdscf.gks module

class pyscf.tdscf.gks.CasidaTDDFT(mf)[source]

Bases: TDDFT, TDA

Solve the Casida TDDFT formula (A-B)(A+B)(X+Y) = (X+Y)w^2

gen_vind(mf=None)[source]

Generate function to compute

[ A B ][X] [-B* -A*][Y]

init_guess(mf, nstates=None, wfnsym=None)
kernel(x0=None, nstates=None)[source]

TDDFT diagonalization solver

nuc_grad_method()[source]
pyscf.tdscf.gks.RPA

alias of TDDFT

class pyscf.tdscf.gks.TDA(mf)[source]

Bases: TDA

class pyscf.tdscf.gks.TDDFT(mf)[source]

Bases: TDHF

pyscf.tdscf.gks.TDDFTNoHybrid

alias of CasidaTDDFT

pyscf.tdscf.gks.TDGKS

alias of TDDFT

pyscf.tdscf.gks.tddft(mf)[source]

Driver to create TDDFT or CasidaTDDFT object

pyscf.tdscf.proxy module

This and other proxy modules implement the time-dependent mean-field procedure using the existing pyscf implementations as a black box. The main purpose of these modules is to overcome the existing limitations in pyscf (i.e. real-only orbitals, davidson diagonalizer, incomplete Bloch space, etc). The primary performance drawback is that, unlike the original pyscf routines with an implicit construction of the eigenvalue problem, these modules construct TD matrices explicitly by proxying to pyscf density response routines with a O(N^4) complexity scaling. As a result, regular numpy.linalg.eig can be used to retrieve TD roots. Several variants of proxy-TD are available:

  • (this module) pyscf.tdscf.proxy: the molecular implementation;

  • pyscf.pbc.tdscf.proxy: PBC (periodic boundary condition) Gamma-point-only implementation;

  • pyscf.pbc.tdscf.kproxy_supercell: PBC implementation constructing supercells. Works with an arbitrary number of k-points but has an overhead due to ignoring the momentum conservation law. In addition, works only with time reversal invariant (TRI) models: i.e. the k-point grid has to be aligned and contain at least one TRI momentum.

  • pyscf.pbc.tdscf.kproxy: same as the above but respect the momentum conservation and, thus, diagonlizes smaller matrices (the performance gain is the total number of k-points in the model).

class pyscf.tdscf.proxy.PhysERI(model, proxy, frozen=None)[source]

Bases: MolecularMFMixin, TDProxyMatrixBlocks

proxy_choices = {'dft': <function TDDFT>, 'hf': <function TDHF>}
proxy_is_double()[source]

Determines if double-sized matrices are proxied.

Returns:

True if double-sized matrices are proxied.

proxy_response()[source]

A raw response matrix.

Returns:

A raw response matrix.

tdhf_primary_form(*args, **kwargs)[source]

A primary form of TD matrixes.

Returns:

Output type: “full”, “ab”, or “mk” and the corresponding matrix(es).

class pyscf.tdscf.proxy.TDProxy(mf, proxy, frozen=None)[source]

Bases: TDBase

ao2mo()[source]

Prepares ERI.

Returns:

A suitable ERI.

proxy_eri

alias of PhysERI

static v2a(vectors, nocc, nmo)

Transforms (reshapes) and normalizes vectors into amplitudes.

Args:

vectors (numpy.ndarray): raw eigenvectors to transform; nocc (int): number of occupied orbitals; nmo (int): the total number of orbitals;

Returns:

Amplitudes with the following shape: (# of roots, 2 (x or y), # of occupied orbitals, # of virtual orbitals).

pyscf.tdscf.proxy.mk_make_canonic(m, o, v, return_ov=False, space_ov=None)[source]

Makes the output of pyscf TDDFT matrix (MK form) to be canonic.

Args:

m (ndarray): the TDDFT matrix; o (ndarray): occupied orbital energies; v (ndarray): virtual orbital energies; return_ov (bool): if True, returns the K-matrix as well; space_ov (ndarray): an optional ov space;

Returns:

The rotated matrix as well as an optional K-matrix.

pyscf.tdscf.proxy.molecular_response(vind, space, nocc, nmo, double, log_dest)[source]

Retrieves a raw response matrix.

Args:

vind (Callable): a pyscf matvec routine; space (ndarray): the active orbital space mask: either the same mask for both rows and columns (1D array) or separate orbital masks for rows and columns (2D array); nocc (int): the number of occupied orbitals (frozen and active); nmo (int): the total number of orbitals; double (bool): set to True if vind returns the double-sized (i.e. full) matrix; log_dest (object): pyscf logging;

Returns:

The TD matrix.

pyscf.tdscf.proxy.molecular_response_ov(vind, space_ov, nocc, nmo, double, log_dest)[source]

Retrieves a raw response matrix.

Args:

vind (Callable): a pyscf matvec routine; space_ov (ndarray): the active ov space mask: either the same mask for both rows and columns (1D array) or separate ov masks for rows and columns (2D array); nocc (int): the number of occupied orbitals (frozen and active); nmo (int): the total number of orbitals; double (bool): set to True if vind returns the double-sized (i.e. full) matrix; log_dest (object): pyscf logging;

Returns:

The TD matrix.

Note:

The runtime scales with the size of the column mask space_ov[1] but not the row mask space_ov[0].

pyscf.tdscf.proxy.orb2ov(space, nocc)[source]

Converts orbital active space specification into ov-pairs space spec.

Args:

space (ndarray): the obital space; nocc (int): the number of occupied orbitals;

Returns:

The ov space specification.

pyscf.tdscf.rhf module

pyscf.tdscf.rhf.CIS

alias of TDA

pyscf.tdscf.rhf.RPA

alias of TDHF

class pyscf.tdscf.rhf.TDA(mf)[source]

Bases: TDBase

Tamm-Dancoff approximation

Attributes:
conv_tolfloat

Diagonalization convergence tolerance. Default is 1e-9.

nstatesint

Number of TD states to be computed. Default is 3.

Saved results:

convergedbool

Diagonalization converged or not

e1D array

excitation energy for each excited state.

xyA list of two 2D arrays

The two 2D arrays are Excitation coefficients X (shape [nocc,nvir]) and de-excitation coefficients Y (shape [nocc,nvir]) for each excited state. (X,Y) are normalized to 1/2 in RHF/RKS methods and normalized to 1 for UHF/UKS methods. In the TDA calculation, Y = 0.

Gradients(*args, **kwargs)
gen_vind(mf=None)[source]

Generate function to compute Ax

init_guess(mf, nstates=None, wfnsym=None)[source]
kernel(x0=None, nstates=None)[source]

TDA diagonalization solver

class pyscf.tdscf.rhf.TDBase(mf)[source]

Bases: StreamObject

DDCOSMO(solvent_obj=None, dm=None)

Add solvent model in TDDFT calculations.

Kwargs:
dmif given, solvent does not respond to the change of density

matrix. A frozen ddCOSMO potential is added to the results.

DDPCM(solvent_obj=None, dm=None)

Add solvent model in TDDFT calculations.

Kwargs:
dmif given, solvent does not respond to the change of density

matrix. A frozen ddCOSMO potential is added to the results.

PCM(solvent_obj=None, dm=None)

Add solvent model in TDDFT calculations.

Kwargs:
dmif given, solvent does not respond to the change of density

matrix. A frozen ddCOSMO potential is added to the results.

analyze(verbose=None)
as_scanner()

Generating a scanner/solver for TDA/TDHF/TDDFT PES.

The returned solver is a function. This function requires one argument “mol” as input and returns total TDA/TDHF/TDDFT energy.

The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the TDA/TDDFT and the underlying SCF objects (conv_tol, max_memory etc) are automatically applied in the solver.

Note scanner has side effects. It may change many underlying objects (_scf, with_df, with_x2c, …) during calculation.

Examples:

>>> from pyscf import gto, scf, tdscf
>>> mol = gto.M(atom='H 0 0 0; F 0 0 1')
>>> td_scanner = tdscf.TDHF(scf.RHF(mol)).as_scanner()
>>> de = td_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1'))
[ 0.34460866  0.34460866  0.7131453 ]
>>> de = td_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5'))
[ 0.14844013  0.14844013  0.47641829]
check_sanity()[source]

Check input of class/object attributes, check whether a class method is overwritten. It does not check the attributes which are prefixed with “_”. The return value of method set is the object itself. This allows a series of functions/methods to be executed in pipe.

conv_tol = 1e-09
ddCOSMO(solvent_obj=None, dm=None)

Add solvent model in TDDFT calculations.

Kwargs:
dmif given, solvent does not respond to the change of density

matrix. A frozen ddCOSMO potential is added to the results.

ddPCM(solvent_obj=None, dm=None)

Add solvent model in TDDFT calculations.

Kwargs:
dmif given, solvent does not respond to the change of density

matrix. A frozen ddCOSMO potential is added to the results.

deg_eia_thresh = 0.001
dump_flags(verbose=None)[source]
property e_tot

Excited state energies

gen_vind(mf=None)[source]
get_ab(mf=None)[source]

A and B matrices for TDDFT response function.

A[i,a,j,b] = delta_{ab}delta_{ij}(E_a - E_i) + (ia||bj) B[i,a,j,b] = (ia||jb)

get_nto(state=1, threshold=0.3, verbose=None)

Natural transition orbital analysis.

The natural transition density matrix between ground state and excited state \(Tia = \langle \Psi_{ex} | i a^\dagger | \Psi_0 \rangle\) can be transformed to diagonal form through SVD \(T = O \sqrt{\lambda} V^\dagger\). O and V are occupied and virtual natural transition orbitals. The diagonal elements \(\lambda\) are the weights of the occupied-virtual orbital pair in the excitation.

Ref: Martin, R. L., JCP, 118, 4775-4777

Note in the TDHF/TDDFT calculations, the excitation part (X) is interpreted as the CIS coefficients and normalized to 1. The de-excitation part (Y) is ignored.

Args:

tdobj : TDA, or TDHF, or TDDFT object

stateint

Excited state ID. state = 1 means the first excited state. If state < 0, state ID is counted from the last excited state.

Kwargs:
thresholdfloat

Above which the NTO coefficients will be printed in the output.

Returns:

A list (weights, NTOs). NTOs are natural orbitals represented in AO basis. The first N_occ NTOs are occupied NTOs and the rest are virtual NTOs.

get_precond(hdiag)[source]
level_shift = 0
lindep = 1e-12
max_cycle = 100
max_space = 50
property nroots
nstates = 3
nuc_grad_method()[source]
oscillator_strength(e=None, xy=None, gauge='length', order=0)
positive_eig_threshold = 0.001
reset(mol=None)[source]
singlet = True
transition_dipole(xy=None)

Transition dipole moments in the length gauge

transition_magnetic_dipole(xy=None)

Transition magnetic dipole moments (imaginary part only)

transition_magnetic_quadrupole(xy=None)

Transition magnetic quadrupole moments (imaginary part only)

transition_octupole(xy=None)

Transition octupole moments in the length gauge

transition_quadrupole(xy=None)

Transition quadrupole moments in the length gauge

transition_velocity_dipole(xy=None)

Transition dipole moments in the velocity gauge (imaginary part only)

transition_velocity_octupole(xy=None)

Transition octupole moments in the velocity gauge (imaginary part only)

transition_velocity_quadrupole(xy=None)

Transition quadrupole moments in the velocity gauge (imaginary part only)

class pyscf.tdscf.rhf.TDHF(mf)[source]

Bases: TDA

Time-dependent Hartree-Fock

Attributes:
conv_tolfloat

Diagonalization convergence tolerance. Default is 1e-9.

nstatesint

Number of TD states to be computed. Default is 3.

Saved results:

convergedbool

Diagonalization converged or not

e1D array

excitation energy for each excited state.

xyA list of two 2D arrays

The two 2D arrays are Excitation coefficients X (shape [nocc,nvir]) and de-excitation coefficients Y (shape [nocc,nvir]) for each excited state. (X,Y) are normalized to 1/2 in RHF/RKS methods and normalized to 1 for UHF/UKS methods. In the TDA calculation, Y = 0.

Gradients(*args, **kwargs)
gen_vind(mf=None)[source]

Generate function to compute

[ A B][X] [-B -A][Y]

init_guess(mf, nstates=None, wfnsym=None)[source]
kernel(x0=None, nstates=None)[source]

TDHF diagonalization with non-Hermitian eigenvalue solver

nuc_grad_method()[source]
pyscf.tdscf.rhf.TDRHF

alias of TDHF

class pyscf.tdscf.rhf.TD_Scanner(td)[source]

Bases: SinglePointScanner

pyscf.tdscf.rhf.analyze(tdobj, verbose=None)[source]
pyscf.tdscf.rhf.analyze_wfnsym(tdobj, x_sym, x)[source]

Guess the wfn symmetry of TDDFT X amplitude

pyscf.tdscf.rhf.as_scanner(td)[source]

Generating a scanner/solver for TDA/TDHF/TDDFT PES.

The returned solver is a function. This function requires one argument “mol” as input and returns total TDA/TDHF/TDDFT energy.

The solver will automatically use the results of last calculation as the initial guess of the new calculation. All parameters assigned in the TDA/TDDFT and the underlying SCF objects (conv_tol, max_memory etc) are automatically applied in the solver.

Note scanner has side effects. It may change many underlying objects (_scf, with_df, with_x2c, …) during calculation.

Examples:

>>> from pyscf import gto, scf, tdscf
>>> mol = gto.M(atom='H 0 0 0; F 0 0 1')
>>> td_scanner = tdscf.TDHF(scf.RHF(mol)).as_scanner()
>>> de = td_scanner(gto.M(atom='H 0 0 0; F 0 0 1.1'))
[ 0.34460866  0.34460866  0.7131453 ]
>>> de = td_scanner(gto.M(atom='H 0 0 0; F 0 0 1.5'))
[ 0.14844013  0.14844013  0.47641829]
pyscf.tdscf.rhf.gen_tda_hop(mf, fock_ao=None, singlet=True, wfnsym=None)

Generate function to compute A x

Kwargs:
wfnsymint or str

Point group symmetry irrep symbol or ID for excited CIS wavefunction.

pyscf.tdscf.rhf.gen_tda_operation(mf, fock_ao=None, singlet=True, wfnsym=None)[source]

Generate function to compute A x

Kwargs:
wfnsymint or str

Point group symmetry irrep symbol or ID for excited CIS wavefunction.

pyscf.tdscf.rhf.gen_tdhf_operation(mf, fock_ao=None, singlet=True, wfnsym=None)[source]

Generate function to compute

[ A B][X] [-B -A][Y]

pyscf.tdscf.rhf.get_ab(mf, mo_energy=None, mo_coeff=None, mo_occ=None)[source]

A and B matrices for TDDFT response function.

A[i,a,j,b] = delta_{ab}delta_{ij}(E_a - E_i) + (ia||bj) B[i,a,j,b] = (ia||jb)

pyscf.tdscf.rhf.get_nto(tdobj, state=1, threshold=0.3, verbose=None)[source]

Natural transition orbital analysis.

The natural transition density matrix between ground state and excited state \(Tia = \langle \Psi_{ex} | i a^\dagger | \Psi_0 \rangle\) can be transformed to diagonal form through SVD \(T = O \sqrt{\lambda} V^\dagger\). O and V are occupied and virtual natural transition orbitals. The diagonal elements \(\lambda\) are the weights of the occupied-virtual orbital pair in the excitation.

Ref: Martin, R. L., JCP, 118, 4775-4777

Note in the TDHF/TDDFT calculations, the excitation part (X) is interpreted as the CIS coefficients and normalized to 1. The de-excitation part (Y) is ignored.

Args:

tdobj : TDA, or TDHF, or TDDFT object

stateint

Excited state ID. state = 1 means the first excited state. If state < 0, state ID is counted from the last excited state.

Kwargs:
thresholdfloat

Above which the NTO coefficients will be printed in the output.

Returns:

A list (weights, NTOs). NTOs are natural orbitals represented in AO basis. The first N_occ NTOs are occupied NTOs and the rest are virtual NTOs.

pyscf.tdscf.rhf.oscillator_strength(tdobj, e=None, xy=None, gauge='length', order=0)[source]
pyscf.tdscf.rhf.transition_dipole(tdobj, xy=None)[source]

Transition dipole moments in the length gauge

pyscf.tdscf.rhf.transition_magnetic_dipole(tdobj, xy=None)[source]

Transition magnetic dipole moments (imaginary part only)

pyscf.tdscf.rhf.transition_magnetic_quadrupole(tdobj, xy=None)[source]

Transition magnetic quadrupole moments (imaginary part only)

pyscf.tdscf.rhf.transition_octupole(tdobj, xy=None)[source]

Transition octupole moments in the length gauge

pyscf.tdscf.rhf.transition_quadrupole(tdobj, xy=None)[source]

Transition quadrupole moments in the length gauge

pyscf.tdscf.rhf.transition_velocity_dipole(tdobj, xy=None)[source]

Transition dipole moments in the velocity gauge (imaginary part only)

pyscf.tdscf.rhf.transition_velocity_octupole(tdobj, xy=None)[source]

Transition octupole moments in the velocity gauge (imaginary part only)

pyscf.tdscf.rhf.transition_velocity_quadrupole(tdobj, xy=None)[source]

Transition quadrupole moments in the velocity gauge (imaginary part only)

pyscf.tdscf.rhf_slow module

This and other _slow modules implement the time-dependent Hartree-Fock procedure. The primary performance drawback is that, unlike other ‘fast’ routines with an implicit construction of the eigenvalue problem, these modules construct TDHF matrices explicitly via an AO-MO transformation, i.e. with a O(N^5) complexity scaling. As a result, regular numpy.linalg.eig can be used to retrieve TDHF roots in a reliable fashion without any issues related to the Davidson procedure. Several variants of TDHF are available:

  • (this module) pyscf.tdscf.rhf_slow: the molecular implementation;

  • pyscf.pbc.tdscf.rhf_slow: PBC (periodic boundary condition) implementation for RHF objects of pyscf.pbc.scf modules;

  • pyscf.pbc.tdscf.krhf_slow_supercell: PBC implementation for KRHF objects of pyscf.pbc.scf modules. Works with an arbitrary number of k-points but has a overhead due to an effective construction of a supercell.

  • pyscf.pbc.tdscf.krhf_slow_gamma: A Gamma-point calculation resembling the original pyscf.pbc.tdscf.krhf module. Despite its name, it accepts KRHF objects with an arbitrary number of k-points but finds only few TDHF roots corresponding to collective oscillations without momentum transfer;

  • pyscf.pbc.tdscf.krhf_slow: PBC implementation for KRHF objects of pyscf.pbc.scf modules. Works with an arbitrary number of k-points and employs k-point conservation (diagonalizes matrix blocks separately).

class pyscf.tdscf.rhf_slow.PhysERI(model, frozen=None)[source]

Bases: MolecularMFMixin, TDERIMatrixBlocks

ao2mo(coeff)[source]

Phys ERI in MO basis.

Args:

coeff (Iterable): MO orbitals;

Returns:

ERI in MO basis.

class pyscf.tdscf.rhf_slow.PhysERI4(model, frozen=None)[source]

Bases: PhysERI

symmetries = [((0, 1, 2, 3), False), ((1, 0, 3, 2), False), ((2, 3, 0, 1), True), ((3, 2, 1, 0), True)]
class pyscf.tdscf.rhf_slow.PhysERI8(model, frozen=None)[source]

Bases: PhysERI4

symmetries = [((0, 1, 2, 3), False), ((1, 0, 3, 2), False), ((2, 3, 0, 1), False), ((3, 2, 1, 0), False), ((2, 1, 0, 3), False), ((3, 0, 1, 2), False), ((0, 3, 2, 1), False), ((1, 2, 3, 0), False)]
class pyscf.tdscf.rhf_slow.TDRHF(mf, frozen=None)[source]

Bases: TDBase

ao2mo()[source]

Picks ERI: either 4-fold or 8-fold symmetric.

Returns:

A suitable ERI.

eri1

alias of PhysERI

eri4

alias of PhysERI4

eri8

alias of PhysERI8

static v2a(vectors, nocc, nmo)

Transforms (reshapes) and normalizes vectors into amplitudes.

Args:

vectors (numpy.ndarray): raw eigenvectors to transform; nocc (int): number of occupied orbitals; nmo (int): the total number of orbitals;

Returns:

Amplitudes with the following shape: (# of roots, 2 (x or y), # of occupied orbitals, # of virtual orbitals).

pyscf.tdscf.rhf_slow.vector_to_amplitudes(vectors, nocc, nmo)[source]

Transforms (reshapes) and normalizes vectors into amplitudes.

Args:

vectors (numpy.ndarray): raw eigenvectors to transform; nocc (int): number of occupied orbitals; nmo (int): the total number of orbitals;

Returns:

Amplitudes with the following shape: (# of roots, 2 (x or y), # of occupied orbitals, # of virtual orbitals).

pyscf.tdscf.rks module

class pyscf.tdscf.rks.CasidaTDDFT(mf)[source]

Bases: TDDFT, TDA

Solve the Casida TDDFT formula (A-B)(A+B)(X+Y) = (X+Y)w^2

gen_vind(mf=None)[source]

Generate function to compute

[ A B][X] [-B -A][Y]

init_guess(mf, nstates=None, wfnsym=None)
kernel(x0=None, nstates=None)[source]

TDDFT diagonalization solver

nuc_grad_method()[source]
pyscf.tdscf.rks.RPA

alias of TDDFT

class pyscf.tdscf.rks.TDA(mf)[source]

Bases: TDA

Gradients(*args, **kwargs)
nuc_grad_method()[source]
class pyscf.tdscf.rks.TDDFT(mf)[source]

Bases: TDHF

Gradients(*args, **kwargs)
nuc_grad_method()[source]
pyscf.tdscf.rks.TDDFTNoHybrid

alias of CasidaTDDFT

pyscf.tdscf.rks.TDH

alias of dRPA

pyscf.tdscf.rks.TDRKS

alias of TDDFT

class pyscf.tdscf.rks.dRPA(mf)[source]

Bases: CasidaTDDFT

class pyscf.tdscf.rks.dTDA(mf)[source]

Bases: TDA

pyscf.tdscf.rks.tddft(mf)[source]

Driver to create TDDFT or CasidaTDDFT object

pyscf.tdscf.uhf module

pyscf.tdscf.uhf.CIS

alias of TDA

pyscf.tdscf.uhf.RPA

alias of TDHF

class pyscf.tdscf.uhf.TDA(mf)[source]

Bases: TDBase

Tamm-Dancoff approximation

Attributes:
conv_tolfloat

Diagonalization convergence tolerance. Default is 1e-9.

nstatesint

Number of TD states to be computed. Default is 3.

Saved results:

convergedbool

Diagonalization converged or not

e1D array

excitation energy for each excited state.

xyA list of two 2D arrays

The two 2D arrays are Excitation coefficients X (shape [nocc,nvir]) and de-excitation coefficients Y (shape [nocc,nvir]) for each excited state. (X,Y) are normalized to 1/2 in RHF/RKS methods and normalized to 1 for UHF/UKS methods. In the TDA calculation, Y = 0.

Gradients(*args, **kwargs)
gen_vind(mf=None)[source]

Generate function to compute Ax

init_guess(mf, nstates=None, wfnsym=None)[source]
kernel(x0=None, nstates=None)[source]

TDA diagonalization solver

singlet = None
class pyscf.tdscf.uhf.TDBase(mf)[source]

Bases: TDBase

analyze(verbose=None)
get_ab(mf=None)[source]

A and B matrices for TDDFT response function.

A[i,a,j,b] = delta_{ab}delta_{ij}(E_a - E_i) + (ia||bj) B[i,a,j,b] = (ia||jb)

Spin symmetry is considered in the returned A, B lists. List A has three items: (A_aaaa, A_aabb, A_bbbb). A_bbaa = A_aabb.transpose(2,3,0,1). B has three items: (B_aaaa, B_aabb, B_bbbb). B_bbaa = B_aabb.transpose(2,3,0,1).

get_nto(state=1, threshold=0.3, verbose=None)

Natural transition orbital analysis.

The natural transition density matrix between ground state and excited state \(Tia = \langle \Psi_{ex} | i a^\dagger | \Psi_0 \rangle\) can be transformed to diagonal form through SVD \(T = O \sqrt{\lambda} V^\dagger\). O and V are occupied and virtual natural transition orbitals. The diagonal elements \(\lambda\) are the weights of the occupied-virtual orbital pair in the excitation.

Ref: Martin, R. L., JCP, 118, 4775-4777

Note in the TDHF/TDDFT calculations, the excitation part (X) is interpreted as the CIS coefficients and normalized to 1. The de-excitation part (Y) is ignored.

Args:
stateint

Excited state ID. state = 1 means the first excited state. If state < 0, state ID is counted from the last excited state.

Kwargs:
thresholdfloat

Above which the NTO coefficients will be printed in the output.

Returns:

A list (weights, NTOs). NTOs are natural orbitals represented in AO basis. The first N_occ NTOs are occupied NTOs and the rest are virtual NTOs.

nuc_grad_method()[source]
class pyscf.tdscf.uhf.TDHF(mf)[source]

Bases: TDA

Gradients(*args, **kwargs)
gen_vind(mf=None)[source]

Generate function to compute

[ A B][X] [-B -A][Y]

init_guess(mf, nstates=None, wfnsym=None)[source]
kernel(x0=None, nstates=None)[source]

TDHF diagonalization with non-Hermitian eigenvalue solver

singlet = None
pyscf.tdscf.uhf.TDUHF

alias of TDHF

pyscf.tdscf.uhf.analyze(tdobj, verbose=None)[source]
pyscf.tdscf.uhf.gen_tda_hop(mf, fock_ao=None, wfnsym=None)

A x

Kwargs:
wfnsymint or str

Point group symmetry irrep symbol or ID for excited CIS wavefunction.

pyscf.tdscf.uhf.gen_tda_operation(mf, fock_ao=None, wfnsym=None)[source]

A x

Kwargs:
wfnsymint or str

Point group symmetry irrep symbol or ID for excited CIS wavefunction.

pyscf.tdscf.uhf.gen_tdhf_operation(mf, fock_ao=None, singlet=True, wfnsym=None)[source]

Generate function to compute

[ A B][X] [-B -A][Y]

pyscf.tdscf.uhf.get_ab(mf, mo_energy=None, mo_coeff=None, mo_occ=None)[source]

A and B matrices for TDDFT response function.

A[i,a,j,b] = delta_{ab}delta_{ij}(E_a - E_i) + (ia||bj) B[i,a,j,b] = (ia||jb)

Spin symmetry is considered in the returned A, B lists. List A has three items: (A_aaaa, A_aabb, A_bbbb). A_bbaa = A_aabb.transpose(2,3,0,1). B has three items: (B_aaaa, B_aabb, B_bbbb). B_bbaa = B_aabb.transpose(2,3,0,1).

pyscf.tdscf.uhf.get_nto(tdobj, state=1, threshold=0.3, verbose=None)[source]

Natural transition orbital analysis.

The natural transition density matrix between ground state and excited state \(Tia = \langle \Psi_{ex} | i a^\dagger | \Psi_0 \rangle\) can be transformed to diagonal form through SVD \(T = O \sqrt{\lambda} V^\dagger\). O and V are occupied and virtual natural transition orbitals. The diagonal elements \(\lambda\) are the weights of the occupied-virtual orbital pair in the excitation.

Ref: Martin, R. L., JCP, 118, 4775-4777

Note in the TDHF/TDDFT calculations, the excitation part (X) is interpreted as the CIS coefficients and normalized to 1. The de-excitation part (Y) is ignored.

Args:
stateint

Excited state ID. state = 1 means the first excited state. If state < 0, state ID is counted from the last excited state.

Kwargs:
thresholdfloat

Above which the NTO coefficients will be printed in the output.

Returns:

A list (weights, NTOs). NTOs are natural orbitals represented in AO basis. The first N_occ NTOs are occupied NTOs and the rest are virtual NTOs.

pyscf.tdscf.uks module

class pyscf.tdscf.uks.CasidaTDDFT(mf)[source]

Bases: TDDFT, TDA

Solve the Casida TDDFT formula (A-B)(A+B)(X+Y) = (X+Y)w^2

gen_vind(mf=None)[source]

Generate function to compute

[ A B][X] [-B -A][Y]

init_guess(mf, nstates=None, wfnsym=None)
kernel(x0=None, nstates=None)[source]

TDDFT diagonalization solver

nuc_grad_method()[source]
pyscf.tdscf.uks.RPA

alias of TDDFT

class pyscf.tdscf.uks.TDA(mf)[source]

Bases: TDA

Gradients(*args, **kwargs)
nuc_grad_method()[source]
class pyscf.tdscf.uks.TDDFT(mf)[source]

Bases: TDHF

Gradients(*args, **kwargs)
nuc_grad_method()[source]
pyscf.tdscf.uks.TDDFTNoHybrid

alias of CasidaTDDFT

pyscf.tdscf.uks.TDH

alias of dRPA

pyscf.tdscf.uks.TDUKS

alias of TDDFT

class pyscf.tdscf.uks.dRPA(mf)[source]

Bases: CasidaTDDFT

class pyscf.tdscf.uks.dTDA(mf)[source]

Bases: TDA

pyscf.tdscf.uks.tddft(mf)[source]

Driver to create TDDFT or CasidaTDDFT object

Module contents

pyscf.tdscf.RPA(mf)[source]
pyscf.tdscf.TD(mf)
pyscf.tdscf.TDA(mf)[source]
pyscf.tdscf.TDDFT(mf)[source]
pyscf.tdscf.TDHF(mf)[source]
pyscf.tdscf.dRPA(mf)[source]
pyscf.tdscf.dTDA(mf)[source]