pyscf.pbc.ci package#
Submodules#
pyscf.pbc.ci.cisd module#
pyscf.pbc.ci.kcis_rhf module#
- class pyscf.pbc.ci.kcis_rhf.KCIS(mf, frozen=None, mo_coeff=None, mo_occ=None)[source]#
Bases:
StreamObject
- get_diag(kshift, eris=None)#
Diagonal elements of CIS Hamiltonian.
- Arguments:
cis {KCIS} – A KCIS instance kshift {int} – k-shift index. A k-shift vector is an exciton momentum.
Available k-shift indices depend on the k-point mesh. For example, a 2 by 2 by 2 k-point mesh allows at most 8 k-shift values, which can be targeted by 0, 1, 2, 3, 4, 5, 6, or 7.
- Keyword Arguments:
- eris {_CIS_ERIS} – Depending on cis.direct, eris may
contain 4-center (cis.direct=False) or 3-center (cis.direct=True) electron repulsion integrals (default: {None})
- Returns:
1D array – an array formed by diagonal elements of CIS Hamiltonian
- get_frozen_mask()#
Boolean mask for orbitals in k-point post-HF method.
Creates a boolean mask to remove frozen orbitals and keep other orbitals for post-HF calculations.
- Args:
mp (
MP2
): An instantiation of an SCF or post-Hartree-Fock object.- Returns:
moidx (list of
ndarray
of bool): Boolean mask of orbitals to include.
- get_kconserv_r(kshift)[source]#
Get the momentum conservation array for a set of k-points.
Given k-point index m the array kconserv_r1[m] returns the index n that satisfies momentum conservation,
(k(m) - k(n) - kshift) dot a = 2npi
This is used for symmetry of 1p-1h excitation operator vector R_{m k_m}^{n k_n} is zero unless n satisfies the above.
Note that this method is adapted from kpts_helper.get_kconserv().
- Arguments:
kshift {int} – index of momentum vector. It can be chosen as any of the available k-point index based on the specified kpt mesh. E.g. int from 0 to 7 can be chosen for a [2,2,2] grid.
- Returns:
list – a list of k(n) corresponding to k(m) that ranges from 0 to max_k_index
- get_nmo(per_kpoint=False)#
Number of orbitals for k-point calculations.
Number of orbitals for use in a calculation with k-points, taking into account frozen orbitals.
- Note:
If per_kpoint is False, then the number of orbitals here is equal to max(nocc) + max(nvir), where each max is done over all k-points. Otherwise the number of orbitals is returned as a list of number of orbitals at each k-point.
- Args:
mp (
MP2
): An instantiation of an SCF or post-Hartree-Fock object. per_kpoint (bool, optional): True returns the number of orbitals at each k-point.For a description of False, see Note.
- Returns:
- nmo (int, list of int): Number of orbitals. For return type, see description of arg
per_kpoint.
- get_nocc(per_kpoint=False)#
Number of occupied orbitals for k-point calculations.
Number of occupied orbitals for use in a calculation with k-points, taking into account frozen orbitals.
- Args:
mp (
MP2
): An instantiation of an SCF or post-Hartree-Fock object. per_kpoint (bool, optional): True returns the number of occupiedorbitals at each k-point. False gives the max of this list.
- Returns:
- nocc (int, list of int): Number of occupied orbitals. For return type, see description of arg
per_kpoint.
- kernel(nroots=1, eris=None, kptlist=None, **kargs)#
CIS excitation energy with k-point sampling.
- Arguments:
cis {KCIS} – A KCIS instance
- Keyword Arguments:
nroots {int} – Number of requested excitation energies (default: {1}) eris {_CIS_ERIS} – Depending on cis.direct, eris may
contain 4-center (cis.direct=False) or 3-center (cis.direct=True) electron repulsion integrals (default: {None})
- kptlist {list} – A list of indices for k-shift, i.e. the exciton momentum.
Available k-shift indices depend on the k-point mesh. For example, a 2 by 2 by 2 k-point mesh allows at most 8 k-shift values, which can be targeted by [0, 1, 2, 3, 4, 5, 6, 7]. When kptlist=None, all k-shift will be computed. (default: {None})
- Returns:
tuple – A tuple of excitation energies and corresponding eigenvectors
- matvec(vector, kshift, eris=None)#
Compute matrix-vector product of the Hamiltonian matrix and a CIS c oefficient vector, in the space of single excitation.
- Arguments:
cis {KCIS} – A KCIS instance vector {1D array} – CIS coefficient vector kshift {int} – k-shift index. A k-shift vector is an exciton momentum.
Available k-shift indices depend on the k-point mesh. For example, a 2 by 2 by 2 k-point mesh allows at most 8 k-shift values, which can be targeted by 0, 1, 2, 3, 4, 5, 6, or 7.
- Keyword Arguments:
- eris {_CIS_ERIS} – Depending on cis.direct, eris may
contain 4-center (cis.direct=False) or 3-center (cis.direct=True) electron repulsion integrals (default: {None})
- Returns:
- 1D array – matrix-vector product of the Hamiltonian matrix and the
input vector.
- property nkpts#
- property nmo#
- property nocc#
- pyscf.pbc.ci.kcis_rhf.cis_H(cis, kshift, eris=None)[source]#
Build full Hamiltonian matrix in the space of single excitation, i.e. CIS Hamiltonian.
- Arguments:
cis {KCIS} – A KCIS instance kshift {int} – k-shift index. A k-shift vector is an exciton momentum.
Available k-shift indices depend on the k-point mesh. For example, a 2 by 2 by 2 k-point mesh allows at most 8 k-shift values, which can be targeted by 0, 1, 2, 3, 4, 5, 6, or 7.
- Keyword Arguments:
- eris {_CIS_ERIS} – Depending on cis.direct, eris may
contain 4-center (cis.direct=False) or 3-center (cis.direct=True) electron repulsion integrals (default: {None})
- Raises:
- MemoryError: MemoryError will be raise if there is not enough space to
store the full Hamiltonian matrix, which scales as Nk^2 O^2 V^2
- Returns:
2D array – the Hamiltonian matrix reshaped into (ki,i,a) by (kj,j,b)
- pyscf.pbc.ci.kcis_rhf.cis_diag(cis, kshift, eris=None)[source]#
Diagonal elements of CIS Hamiltonian.
- Arguments:
cis {KCIS} – A KCIS instance kshift {int} – k-shift index. A k-shift vector is an exciton momentum.
Available k-shift indices depend on the k-point mesh. For example, a 2 by 2 by 2 k-point mesh allows at most 8 k-shift values, which can be targeted by 0, 1, 2, 3, 4, 5, 6, or 7.
- Keyword Arguments:
- eris {_CIS_ERIS} – Depending on cis.direct, eris may
contain 4-center (cis.direct=False) or 3-center (cis.direct=True) electron repulsion integrals (default: {None})
- Returns:
1D array – an array formed by diagonal elements of CIS Hamiltonian
- pyscf.pbc.ci.kcis_rhf.cis_matvec_singlet(cis, vector, kshift, eris=None)[source]#
Compute matrix-vector product of the Hamiltonian matrix and a CIS c oefficient vector, in the space of single excitation.
- Arguments:
cis {KCIS} – A KCIS instance vector {1D array} – CIS coefficient vector kshift {int} – k-shift index. A k-shift vector is an exciton momentum.
Available k-shift indices depend on the k-point mesh. For example, a 2 by 2 by 2 k-point mesh allows at most 8 k-shift values, which can be targeted by 0, 1, 2, 3, 4, 5, 6, or 7.
- Keyword Arguments:
- eris {_CIS_ERIS} – Depending on cis.direct, eris may
contain 4-center (cis.direct=False) or 3-center (cis.direct=True) electron repulsion integrals (default: {None})
- Returns:
- 1D array – matrix-vector product of the Hamiltonian matrix and the
input vector.
- pyscf.pbc.ci.kcis_rhf.kernel(cis, nroots=1, eris=None, kptlist=None, **kargs)[source]#
CIS excitation energy with k-point sampling.
- Arguments:
cis {KCIS} – A KCIS instance
- Keyword Arguments:
nroots {int} – Number of requested excitation energies (default: {1}) eris {_CIS_ERIS} – Depending on cis.direct, eris may
contain 4-center (cis.direct=False) or 3-center (cis.direct=True) electron repulsion integrals (default: {None})
- kptlist {list} – A list of indices for k-shift, i.e. the exciton momentum.
Available k-shift indices depend on the k-point mesh. For example, a 2 by 2 by 2 k-point mesh allows at most 8 k-shift values, which can be targeted by [0, 1, 2, 3, 4, 5, 6, 7]. When kptlist=None, all k-shift will be computed. (default: {None})
- Returns:
tuple – A tuple of excitation energies and corresponding eigenvectors
Module contents#
- pyscf.pbc.ci.CIS(mf, frozen=None, mo_coeff=None, mo_occ=None)#
- pyscf.pbc.ci.CISD(mf, frozen=None, mo_coeff=None, mo_occ=None)#