pyscf.pbc.tdscf package

Submodules

pyscf.pbc.tdscf.kproxy module

pyscf.pbc.tdscf.kproxy_supercell 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:

  • pyscf.tdscf.proxy: the molecular implementation;

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

  • (this module) 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.pbc.tdscf.kproxy_supercell.PhysERI(model, proxy, x, mf_constructor, frozen=None, **kwargs)[source]

Bases: PeriodicMFMixin, TDProxyMatrixBlocks

proxy_choices = {'dft': <function KTDDFT>, 'hf': <function KTDHF>}
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.pbc.tdscf.kproxy_supercell.TDProxy(mf, proxy, x, mf_constructor, frozen=None, **kwargs)[source]

Bases: TDProxy

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 (tuple): numbers of occupied orbitals; nmo (tuple): the total numbers of AOs per k-point;

Returns:

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

pyscf.pbc.tdscf.kproxy_supercell.assert_scf_converged(model, threshold=1e-07)[source]

Tests if scf is converged. Args:

model: a mean-field model to test; threshold (float): threshold for eigenvalue comparison;

Returns:

True if model is converged and False otherwise.

pyscf.pbc.tdscf.kproxy_supercell.get_sparse_ov_transform(oo, vv)[source]

Retrieves a sparse ovov transform out of sparse oo and vv transforms. Args:

oo (ndarray): the transformation in the occupied space; vv (ndarray): the transformation in the virtual space;

Returns:

The resulting matrix representing the sparse transform in the ov space.

pyscf.pbc.tdscf.kproxy_supercell.k2s(model, grid_spec, mf_constructor, threshold=None, degeneracy_threshold=None, imaginary_threshold=None)[source]

Converts k-point model into a supercell with real orbitals. Args:

model: a mean-field pbc model; grid_spec (Iterable): integer dimensions of the k-grid in the mean-field model; mf_constructor (Callable): a function constructing the mean-field object; threshold (float): a threshold for determining the negative k-point index; degeneracy_threshold (float): a threshold for assuming degeneracy when composing real-valued orbitals; imaginary_threshold (float): a threshold for asserting real-valued supercell orbitals;

Returns:

The same class where the Cell object was replaced by the supercell and all fields were adjusted accordingly.

pyscf.pbc.tdscf.kproxy_supercell.ko_mask(nocc, nmo)[source]

Prepares a mask of an occupied space. Args:

nocc (Iterable): occupation numbers per k-point; nmo (Iterable): numbers of orbitals per k-point;

Returns:

The mask where True denotes occupied orbitals. Basis order: [k, orb=o+v]

pyscf.pbc.tdscf.kproxy_supercell.minus_k(model, threshold=None, degeneracy_threshold=None)[source]

Retrieves an array of indexes of negative k. Args:

model: a mean-field pbc model; threshold (float): a threshold for determining the negative; degeneracy_threshold (float): a threshold for assuming degeneracy;

Returns:

A list of integers with indexes of the corresponding k-points.

pyscf.pbc.tdscf.kproxy_supercell.orb2ov(space, nocc, nmo)[source]

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

space (ndarray): the obital space. Basis order: [k, orb=o+v]; nocc (Iterable): the numbers of occupied orbitals per k-point; nmo (Iterable): the total numbers of orbitals per k-point;

Returns:

The ov space specification. Basis order: [k_o, o, k_v, v].

pyscf.pbc.tdscf.kproxy_supercell.ov2orb(space, nocc, nmo)[source]

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

space (ndarray): the ov space. Basis order: [k_o, o, k_v, v]; nocc (Iterable): the numbers of occupied orbitals per k-point; nmo (Iterable): the total numbers of orbitals per k-point;

Returns:

The orbital space specification. Basis order: [k, orb=o+v].

pyscf.pbc.tdscf.kproxy_supercell.sparse_transform(m, *args)[source]

Performs a sparse transform of a dense tensor. Args:

m (ndarray): a tensor to transform; *args: alternating indexes and bases to transform into;

Returns:

The transformed tensor.

pyscf.pbc.tdscf.kproxy_supercell.split_transform(transform, nocc, nmo, tolerance=1e-14)[source]

Splits the transform into oo and vv blocks. Args:

transform (numpy.ndarray): the original transform. The basis order for the transform is [real orb=o+v; k, orb=o+v]; nocc (Iterable): occupation numbers per k-point; nmo (Iterable): the number of orbitals per k-point; tolerance (float): tolerance to check zeros at the ov block;

Returns:

oo and vv blocks of the transform.

pyscf.pbc.tdscf.kproxy_supercell.supercell_response(vind, space, nocc, nmo, double, rot_bloch, log_dest)[source]

Retrieves a raw response matrix. Args:

vind (Callable): a pyscf matvec routine; space (ndarray): the active space: either for both rows and columns (1D array) or for rows and columns separately (2D array). Basis order: [k, orb=o+v]; nocc (ndarray): the numbers of occupied orbitals (frozen and active) per k-point; nmo (ndarray): the total number of orbitals per k-point; double (bool): set to True if vind returns the double-sized (i.e. full) matrix; rot_bloch (ndarray): a matrix specifying the rotation from real orbitals returned from pyscf to Bloch functions; log_dest (object): pyscf logging;

Returns:

The TD matrix.

pyscf.pbc.tdscf.kproxy_supercell.supercell_response_ov(vind, space_ov, nocc, nmo, double, rot_bloch, 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). Basis order: [k_o, o, k_v, v]; nocc (ndarray): the numbers of occupied orbitals (frozen and active) per k-point; nmo (ndarray): the total number of orbitals per k-point; double (bool): set to True if vind returns the double-sized (i.e. full) matrix; rot_bloch (ndarray): a matrix specifying the rotation from real orbitals returned from pyscf to Bloch functions; log_dest (object): pyscf logging;

Returns:

The TD matrix.

pyscf.pbc.tdscf.kproxy_supercell.supercell_space_required(transform_oo, transform_vv, final_space)[source]

For a given orbital transformation and a given ov mask in the transformed space, calculates a minimal ov mask in the original space required to achieve this transform. Args:

transform_oo (ndarray): the transformation in the occupied space; transform_vv (ndarray): the transformation in the virtual space; final_space (ndarray): the final ov space. Basis order: [k_o, o, k_v, v];

Returns:

The initial active space. Basis order: [k_o, o, k_v, v].

pyscf.pbc.tdscf.krhf module

pyscf.pbc.tdscf.krhf.CIS

alias of TDA

pyscf.pbc.tdscf.krhf.KTDA

alias of TDA

class pyscf.pbc.tdscf.krhf.KTDBase(mf, kshift_lst=None)[source]

Bases: TDBase

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.

dump_flags(verbose=None)[source]
get_nto(*args, **kwargs)

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.pbc.tdscf.krhf.KTDHF

alias of TDHF

pyscf.pbc.tdscf.krhf.RPA

alias of TDHF

class pyscf.pbc.tdscf.krhf.TDA(mf, kshift_lst=None)[source]

Bases: KTDBase

conv_tol = 1e-06
gen_vind(mf, kshift)[source]
init_guess(mf, kshift, nstates=None)[source]
kernel(x0=None)[source]

TDA diagonalization solver

Args:
x0: list of init guess arrays for each k-shift specified in self.kshift_lst

[x0_1, x0_2, …, x0_nshift]

x0_i ~ (nstates, nkpts*nocc*nvir)

class pyscf.pbc.tdscf.krhf.TDHF(mf, kshift_lst=None)[source]

Bases: TDA

gen_vind(mf, kshift)[source]

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

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

TDHF diagonalization with non-Hermitian eigenvalue solver

pyscf.pbc.tdscf.krhf.purify_krlyov_heff(precision, hermi, log)[source]

pyscf.pbc.tdscf.krhf_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:

  • 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;

  • (this module) 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.pbc.tdscf.krhf_slow.PhysERI(model, frozen=None)[source]

Bases: PhysERI

eri_mknj(item, pair_row, pair_column)[source]

Retrieves the merged ERI block using ‘mknj’ notation with pairs of k-indexes (k1, k1, k2, k2). Args:

item (str): a 4-character string of ‘mknj’ letters; pair_row (Iterable): a k-point pair k2 = pair_row[k1] for each k1 (row indexes in the final matrix); pair_column (Iterable): a k-point pair k4 = pair_row[k3] for each k3 (column indexes in the final matrix);

Returns:

The corresponding block of ERI (phys notation).

get_k_ix(item, like)[source]

Retrieves block indexes: row and column. Args:

item (str): a string of ‘mknj’ letters; like (tuple): a 2-tuple with sample pair of k-points;

Returns:

Row and column indexes of a sub-block with conserving momentum.

primary_driver = 'full'
tdhf_diag(block)[source]

Retrieves the merged diagonal block only with specific pairs of k-indexes (k, block[k]). Args:

block (Iterable): a k-point pair k2 = pair[k1] for each k1;

Returns:

The diagonal block.

tdhf_primary_form(k)[source]

A primary form of TDHF matrixes (full). Args:

k (tuple, int): momentum transfer: either a pair of k-point indexes specifying the momentum transfer vector or a single integer with the second index assuming the first index being zero;

Returns:

Output type: “full”, and the corresponding matrix.

class pyscf.pbc.tdscf.krhf_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.pbc.tdscf.krhf_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.pbc.tdscf.krhf_slow.TDRHF(mf, frozen=None)[source]

Bases: TDRHF

eri4

alias of PhysERI4

eri8

alias of PhysERI8

kernel(k=None)[source]

Calculates eigenstates and eigenvalues of the TDHF problem. Args:

k (tuple, int): momentum transfer: either an index specifying the momentum transfer or a list of such indexes;

Returns:

Positive eigenvalues and eigenvectors.

static v2a(vectors, nocc, nmo)

Transforms (reshapes) and normalizes vectors into amplitudes. Args:

vectors (numpy.ndarray): raw eigenvectors to transform; nocc (tuple): numbers of occupied orbitals; nmo (int): the total number of orbitals per k-point;

Returns:

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

pyscf.pbc.tdscf.krhf_slow.get_block_k_ix(eri, k)[source]

Retrieves k indexes of the block with a specific momentum transfer. Args:

eri (TDDFTMatrixBlocks): ERI of the problem; k (tuple, int): momentum transfer: either a pair of k-point indexes specifying the momentum transfer vector or a single integer with the second index assuming the first index being zero;

Returns:

4 arrays: r1, r2, c1, c2 specifying k-indexes of the ERI matrix block.

k34=0,c1[0]

k34=1,c1[1]

k34=nk-1,c1[-1]

k34=0,c2[0]

k34=1,c2[1]

k34=nk-1,c2[-1]

k12=0,r1[0]

Block r1, c1

Block r1, c2

k12=1,r1[1]

k12=nk-1,r1[-1]

k12=0,r2[0]

Block r2, c1

Block r2, c2

k12=1,r2[1]

k12=nk-1,r2[-1]

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

Transforms (reshapes) and normalizes vectors into amplitudes. Args:

vectors (numpy.ndarray): raw eigenvectors to transform; nocc (tuple): numbers of occupied orbitals; nmo (int): the total number of orbitals per k-point;

Returns:

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

pyscf.pbc.tdscf.krhf_slow_gamma 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:

  • 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;

  • (this module) 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.pbc.tdscf.krhf_slow_gamma.PhysERI(model, frozen=None)[source]

Bases: PhysERI

eri_mknj(item)[source]

Retrieves the merged ERI block using ‘mknj’ notation with pairs of k-indexes (k1, k1, k2, k2). Args:

item (str): a 4-character string of ‘mknj’ letters;

Returns:

The corresponding block of ERI (phys notation).

tdhf_diag()[source]

Retrieves the merged diagonal block with equal pairs of k-indexes (k, k) only.

Returns:

The diagonal block.

class pyscf.pbc.tdscf.krhf_slow_gamma.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.pbc.tdscf.krhf_slow_gamma.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.pbc.tdscf.krhf_slow_gamma.TDRHF(mf, frozen=None)[source]

Bases: TDRHF

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 (tuple): numbers of occupied orbitals; nmo (int): the total number of orbitals per k-point;

Returns:

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

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

Transforms (reshapes) and normalizes vectors into amplitudes. Args:

vectors (numpy.ndarray): raw eigenvectors to transform; nocc (tuple): numbers of occupied orbitals; nmo (int): the total number of orbitals per k-point;

Returns:

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

pyscf.pbc.tdscf.krhf_slow_supercell 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:

  • 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;

  • (this module) 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.pbc.tdscf.krhf_slow_supercell.PhysERI(model, frozen=None)[source]

Bases: PeriodicMFMixin, TDERIMatrixBlocks

ao2mo_k(coeff, k)[source]

Phys ERI in MO basis. Args:

coeff (Iterable): MO orbitals; k (Iterable): the 4 k-points MOs correspond to;

Returns:

ERI in MO basis.

eri_mknj(item, pairs_row=None, pairs_column=None)[source]

Retrieves the merged ERI block using ‘mknj’ notation with all k-indexes. Args:

item (str): a 4-character string of ‘mknj’ letters; pairs_row (Iterable): iterator for pairs of row k-points (first index in the output matrix); pairs_column (Iterable): iterator for pairs of column k-points (second index in the output matrix);

Returns:

The corresponding block of ERI (phys notation).

eri_mknj_k(item, k)[source]

Retrieves ERI block using ‘mknj’ notation. Args:

item (str): a 4-character string of ‘mknj’ letters; k (Iterable): k indexes;

Returns:

The corresponding block of ERI (phys notation).

tdhf_diag(pairs=None)[source]

Retrieves the merged diagonal block with specified or all possible k-index pairs. Args:

pairs (Iterable): pairs of k-points to assmble;

Returns:

The diagonal block.

tdhf_diag_k(k1, k2)[source]

Retrieves the diagonal block. Args:

k1 (int): first k-index (row); k2 (int): second k-index (column);

Returns:

The diagonal block.

class pyscf.pbc.tdscf.krhf_slow_supercell.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.pbc.tdscf.krhf_slow_supercell.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.pbc.tdscf.krhf_slow_supercell.TDRHF(mf, frozen=None)[source]

Bases: TDRHF

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 (tuple): numbers of occupied orbitals; nmo (tuple): the total numbers of AOs per k-point;

Returns:

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

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

Transforms (reshapes) and normalizes vectors into amplitudes. Args:

vectors (numpy.ndarray): raw eigenvectors to transform; nocc (tuple): numbers of occupied orbitals; nmo (tuple): the total numbers of AOs per k-point;

Returns:

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

pyscf.pbc.tdscf.krks module

pyscf.pbc.tdscf.krks.KTDA

alias of TDA

pyscf.pbc.tdscf.krks.KTDDFT

alias of TDDFT

pyscf.pbc.tdscf.krks.RPA

alias of TDDFT

class pyscf.pbc.tdscf.krks.TDA(mf, kshift_lst=None)[source]

Bases: TDA

kernel(x0=None)[source]

TDA diagonalization solver

Args:
x0: list of init guess arrays for each k-shift specified in self.kshift_lst

[x0_1, x0_2, …, x0_nshift]

x0_i ~ (nstates, nkpts*nocc*nvir)

class pyscf.pbc.tdscf.krks.TDDFT(mf, kshift_lst=None)[source]

Bases: TDHF

kernel(x0=None)[source]

TDHF diagonalization with non-Hermitian eigenvalue solver

pyscf.pbc.tdscf.kuhf module

pyscf.pbc.tdscf.kuhf.CIS

alias of TDA

pyscf.pbc.tdscf.kuhf.KTDA

alias of TDA

pyscf.pbc.tdscf.kuhf.KTDHF

alias of TDHF

pyscf.pbc.tdscf.kuhf.RPA

alias of TDHF

class pyscf.pbc.tdscf.kuhf.TDA(mf, kshift_lst=None)[source]

Bases: KTDBase

conv_tol = 1e-06
gen_vind(mf, kshift)[source]

Compute Ax

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

TDA diagonalization solver

class pyscf.pbc.tdscf.kuhf.TDHF(mf, kshift_lst=None)[source]

Bases: TDA

gen_vind(mf, kshift)[source]

Compute Ax

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

TDHF diagonalization with non-Hermitian eigenvalue solver

pyscf.pbc.tdscf.kuks module

pyscf.pbc.tdscf.kuks.KTDA

alias of TDA

pyscf.pbc.tdscf.kuks.KTDDFT

alias of TDDFT

pyscf.pbc.tdscf.kuks.RPA

alias of TDDFT

class pyscf.pbc.tdscf.kuks.TDA(mf, kshift_lst=None)[source]

Bases: TDA

kernel(x0=None)[source]

TDA diagonalization solver

class pyscf.pbc.tdscf.kuks.TDDFT(mf, kshift_lst=None)[source]

Bases: TDHF

kernel(x0=None)[source]

TDHF diagonalization with non-Hermitian eigenvalue solver

pyscf.pbc.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:

  • pyscf.tdscf.proxy: the molecular implementation;

  • (this module) 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.pbc.tdscf.proxy.PhysERI(model, proxy, frozen=None)[source]

Bases: PhysERI

proxy_choices = {'dft': <function KTDDFT>, 'hf': <function KTDHF>}
class pyscf.pbc.tdscf.proxy.TDProxy(mf, proxy, frozen=None)[source]

Bases: TDProxy

proxy_eri

alias of PhysERI

pyscf.pbc.tdscf.rhf module

pyscf.pbc.tdscf.rhf.CIS

alias of TDA

pyscf.pbc.tdscf.rhf.RPA

alias of TDHF

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

Bases: TDBase

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

TDA diagonalization solver

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

Bases: TDBase

analyze(*args, **kwargs)
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]
oscillator_strength(*args, **kwargs)
transition_dipole(*args, **kwargs)

Transition dipole moments in the length gauge

transition_magnetic_dipole(*args, **kwargs)

Transition magnetic dipole moments (imaginary part only)

transition_magnetic_quadrupole(*args, **kwargs)

Transition magnetic quadrupole moments (imaginary part only)

transition_octupole(*args, **kwargs)

Transition octupole moments in the length gauge

transition_quadrupole(*args, **kwargs)

Transition quadrupole moments in the length gauge

transition_velocity_dipole(*args, **kwargs)

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

transition_velocity_octupole(*args, **kwargs)

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

transition_velocity_quadrupole(*args, **kwargs)

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

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

Bases: TDA

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

TDHF diagonalization with non-Hermitian eigenvalue solver

pyscf.pbc.tdscf.rhf.TDRHF

alias of TDHF

pyscf.pbc.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:

  • pyscf.tdscf.rhf_slow: the molecular implementation;

  • (this module) 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).

pyscf.pbc.tdscf.rks module

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

Bases: TDA

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

TDDFT diagonalization solver

pyscf.pbc.tdscf.rks.TDDFTNoHybrid

alias of CasidaTDDFT

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

Driver to create TDDFT or CasidaTDDFT object

pyscf.pbc.tdscf.uhf module

pyscf.pbc.tdscf.uhf.CIS

alias of TDA

pyscf.pbc.tdscf.uhf.RPA

alias of TDHF

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

Bases: TDBase

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

TDA diagonalization solver

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

Bases: TDA

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

TDHF diagonalization with non-Hermitian eigenvalue solver

singlet = None
pyscf.pbc.tdscf.uhf.TDUHF

alias of TDHF

pyscf.pbc.tdscf.uks module

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

Bases: TDA

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

TDDFT diagonalization solver

pyscf.pbc.tdscf.uks.TDDFTNoHybrid

alias of CasidaTDDFT

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

Driver to create TDDFT or CasidaTDDFT object

Module contents

pyscf.pbc.tdscf.KTD(mf)
pyscf.pbc.tdscf.KTDA(mf)[source]
pyscf.pbc.tdscf.KTDDFT(mf)[source]
pyscf.pbc.tdscf.KTDHF(mf)[source]
pyscf.pbc.tdscf.TDA(mf)[source]
pyscf.pbc.tdscf.TDDFT(mf)[source]
pyscf.pbc.tdscf.TDHF(mf)[source]