pyscf.agf2 package

Submodules

pyscf.agf2.aux_space module

Auxiliary space class and helper functions.

class pyscf.agf2.aux_space.AuxiliarySpace(energy, coupling, chempot=0.0)[source]

Bases: object

Simple container to hold the energies, couplings and chemical

potential associated with an auxiliary space.

Attributes:
energy1D array

Energies of the poles

coupling2D array

Coupling vector of the poles to each physical state

chempotfloat

Chemical potential associated with the energies

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

See subclasses.

copy()[source]

Returns a copy of the current object.

Returns:

AuxiliarySpace

dot(phys, vec, out=None, chempot=0.0)[source]

Returns the dot product of get_array() with a vector.

Args:
phys2D array

Physical (1p + 1h) part of the matrix

vecndarray

Vector to compute dot product with

Kwargs:
out2D array

If provided, use to store output

chempotfloat

Scale energies (by default, chempot is not used and energies retain their values). Default 0.0

Returns:

ndarray with shape of vec

eig(phys, out=None, chempot=0.0)[source]
Computes the eigenvalues and eigenvectors of the array

returned by get_array().

Args:
phys2D array

Physical (1p + 1h) part of the matrix

Kwargs:
out2D array

If provided, use to store output

chempotfloat

Scale energies (by default, chempot is not used and energies retain their values). Default 0.0

Returns:

tuple of ndarrays (eigenvalues, eigenvectors)

get_array(phys, out=None, chempot=0.0)[source]
Expresses the auxiliaries as an array, i.e. the extended

Fock matrix in AGF2 or Hamiltonian of ADC(2).

Args:
phys2D array

Physical (1p + 1h) part of the matrix

Kwargs:
out2D array

If provided, use to store output

chempotfloat

Scale energies (by default, chempot is not used and energies retain their values). Default 0.0

Returns:

Array representing the coupling of the auxiliary space to the physical space

get_occupied()[source]
Returns a copy of the current AuxiliarySpace object

containing only the poles with energy less than the chemical potential. The object should already be sorted.

Returns:

AuxiliarySpace with only the occupied auxiliaries

get_virtual()[source]
Returns a copy of the current AuxiliarySpace object

containing only the poles with energy greater than the chemical potential. The object should already be sorted.

Returns:

AuxiliarySpace with only the virtual auxiliaries

classmethod load(chkfile, key=None)[source]

Loads the auxiliaries from a chkfile

Args:
chkfilestr

Name of chkfile

keystr

Key to be used in h5py object. It can contain “/” to represent the path in the HDF5 storage structure.

moment(n, squeeze=True)[source]

Builds the nth moment of the spectral distribution.

Args:
nint or list of int

Moment(s) to compute

Kwargs:
squeezebool

If True, use np.squeeze() on output so that in the case of n being an int, a 2D array is returned. If False, output is always 3D. Default True.

Returns:

ndarray of moments

property naux
property nphys
real_freq_spectrum(*args, **kwargs)[source]

See subclasses.

remove_uncoupled(tol)[source]
Removes poles with very low spectral weight (uncoupled

to the physical space) in-place.

Args:
tolfloat

Threshold for the spectral weight (squared norm)

save(chkfile, key=None)[source]

Saves the auxiliaries in chkfile

Args:
chkfilestr

Name of chkfile

keystr

Key to be used in h5py object. It can contain “/” to represent the path in the HDF5 storage structure.

sort()[source]

Sort in-place via the energies to make slicing easier.

class pyscf.agf2.aux_space.GreensFunction(energy, coupling, chempot=0.0)[source]

Bases: AuxiliarySpace

Defines a Green’s function represented as a AuxiliarySpace object.

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

See subclasses.

make_rdm1(chempot=None, occupancy=2)[source]
Returns the first-order reduced density matrix associated

with the Green’s function.

Kwargs:
chempotfloat

If provided, use instead of self.chempot

occupancyint

Occupancy of the states, i.e. 2 for RHF and 1 for UHF

real_freq_spectrum(grid, eta=0.02)[source]
Express the auxiliaries as a spectral function on the real

frequency axis.

Args:
grid1D array

Real frequency grid

Kwargs:
etafloat

Peak broadening factor in Hartrees. Default is 0.02.

Returns:

ndarray of the spectrum, with the first index being the frequency

class pyscf.agf2.aux_space.SelfEnergy(energy, coupling, chempot=0.0)[source]

Bases: AuxiliarySpace

Defines a self-energy represented as a AuxiliarySpace object.

compress(phys=None, n=(None, 0), tol=1e-12)[source]
Compress the auxiliaries via moments of the particle and

hole Green’s function and self-energy. Resulting naux depends on the chosen n.

Kwargs:
phys2D array or None

Physical space (1p + 1h), typically the Fock matrix. Only required if n[0] is not None.

ntuple of int

Compression level of the Green’s function and self-energy, respectively.

tolfloat

Linear dependency tolerance. Default value is 1e-12

Returns:

SelfEnergy with reduced auxiliary dimension.

Raises:

MemoryError if the compression according to Green’s function moments will exceed the maximum allowed memory.

get_greens_function(phys)[source]
Returns a GreensFunction by solving the Dyson

equation.

Args:
phys2D array

Physical space (1p + 1h), typically the Fock matrix

Returns:

GreensFunction

make_rdm1(phys, chempot=None, occupancy=2)[source]
Returns the first-order reduced density matrix associated

with the self-energy via the GreensFunction.

Args:
phys2D array

Physical space (1p + 1h), typically the Fock matrix

Kwargs:
chempotfloat

If provided, use instead of self.chempot

occupancyint

Occupancy of the states, i.e. 2 for RHF and 1 for UHF

real_freq_spectrum(grid, eta=0.02)[source]

See subclasses.

pyscf.agf2.aux_space.combine(*auxspcs)[source]

Combine a set of AuxiliarySpace objects. attr:chempot is inherited from the first element.

pyscf.agf2.aux_space.compress_via_gf(se, phys, n=0, tol=1e-12)[source]
Compress the auxiliaries of the separate occupied and

virtual parts of the self-energy according to consistency in the moments of the Green’s function

Args:
seSelfEnergy

Auxiliaries of the self-energy

phys2D array

Physical space (1p + 1h), typically the Fock matrix

Kwargs:
nint

Truncation parameter, conserves the separate particle and hole moments to order 2*n+1.

tolfloat

Linear dependency tolerance. Default value is 1e-12

Returns:

SelfEnergy with reduced auxiliary dimension

pyscf.agf2.aux_space.compress_via_se(se, n=0)[source]
Compress the auxiliaries of the separate occupied and

virtual parts of the self-energy according to consistency in its moments.

Args:
seSelfEnergy

Auxiliaries of the self-energy

Kwargs:
nint

Truncation parameter, conserves the separate particle and hole moments to order 2*n+1.

Returns:

SelfEnergy with reduced auxiliary dimension

Ref:
[1] H. Muther, T. Taigel and T.T.S. Kuo, Nucl. Phys., 482,

1988, pp. 601-616.

[2] D. Van Neck, K. Piers and M. Waroquier, J. Chem. Phys.,

115, 2001, pp. 15-25.

[3] H. Muther and L.D. Skouras, Nucl. Phys., 55, 1993,

pp. 541-562.

[4] Y. Dewulf, D. Van Neck, L. Van Daele and M. Waroquier,

Phys. Lett. B, 396, 1997, pp. 7-14.

pyscf.agf2.aux_space.davidson(auxspc, phys, chempot=None, nroots=1, which='SM', tol=1e-14, maxiter=None, ntrial=None)[source]
Diagonalise the result of AuxiliarySpace.get_array() using

the sparse AuxiliarySpace.dot() method, with the Davidson algorithm.

This algorithm may perform poorly for IPs or EAs if they are not extremal eigenvalues, which they are not in standard AGF2.

Args:
auxspcAuxiliarySpace or subclass

Auxiliary space object to solve for

phys2D array

Physical space (1p + 1h), typically the Fock matrix

Kwargs:
chempotfloat

If provided, use instead of self.chempot

nrootsint

Number of roots to solve for. Default 1.

whichstr
Which eigenvalues to solve for. Options are:

LM : Largest (in magnitude) eigenvalues. SM : Smallest (in magnitude) eigenvalues. LA : Largest (algebraic) eigenvalues. SA : Smallest (algebraic) eigenvalues.

Default ‘SM’.

tolfloat

Convergence threshold

maxiterint

Maximum number of iterations. Default 10*dim

ntrialint

Maximum number of trial vectors. Default min(dim, max(2*nroots+1, 20))

Returns:

tuple of ndarrays (eigenvalues, eigenvectors)

pyscf.agf2.chempot module

Functions for tuning the chemical potential.

pyscf.agf2.chempot.binsearch_chempot(fock, nphys, nelec, occupancy=2)[source]
Finds a chemical potential which best agrees with the number

of physical electrons and abides by the Aufbau principal via a binary search.

Args:
fock2D array or tuple of arrays

Fock matrix to diagonalise, may be the physical Fock matrix or extended Fock matrix. Can also be the output of np.linalg.eigh() for this matrix, i.e. a tuple of the eigenvalues and eigenvectors.

nphysint

Number of physical degrees of freedom

nelecint

Number of physical electrons

Kwargs:
occupancyint

Occupancy of the states, i.e. 2 for RHF and 1 for UHF. Default 2.

Returns:

chemical potential, and the error in the number of electrons

pyscf.agf2.chempot.minimize_chempot(se, fock, nelec, occupancy=2, x0=0.0, tol=1e-06, maxiter=200, jac=True)[source]
Finds a set of auxiliary energies and chemical potential on

the physical space which best satisfy the number of electrons.

Args:
seAuxiliarySpace

Auxiliary space

fock2D array
phys2D array

Physical space (1p + 1h), typically the Fock matrix

nelecint

Number of physical electrons

Kwargs:
occupancyint

Occupancy of the states, i.e. 2 for RHF and 1 for UHF. Default 2.

x0float

Initial guess for chempot. Default 0.0

tolfloat

Convergence threshold (units are the same as nelec). Default 1e-6.

maxiterint

Maximum number of iterations. Default 200.

jacbool

If True, use gradient. Default True.

Returns:

AuxiliarySpace object with altered energy and chempot, and the SciPy OptimizeResult object.

pyscf.agf2.chkfile module

Functions to support chkfiles with MPI

pyscf.agf2.chkfile.dump_agf2(agf2, chkfile=None, key='agf2', gf=None, se=None, frozen=None, nmom=None, mo_energy=None, mo_coeff=None, mo_occ=None)[source]

Save the AGF2 calculatuion to a chkfile.

pyscf.agf2.chkfile.load(chkfile, key)[source]

Load array(s) from chkfile.

See pyscf.lib.chkfile

pyscf.agf2.chkfile.load_agf2(chkfile)[source]

Load the AGF2 data from the chkfile.

pyscf.agf2.chkfile.load_mol(chkfile)[source]

Load the mol from the chkfile.

See pyscf.lib.chkfile

pyscf.agf2.dfragf2 module

Auxiliary second-order Green’s function perturbation theory with density fitting

class pyscf.agf2.dfragf2.DF(mol, auxbasis=None)[source]

Bases: DF

Replaces the DF.prange function with one which natively supports MPI, if used.

prange(start=None, stop=None, step=None)[source]
class pyscf.agf2.dfragf2.DFRAGF2(mf, frozen=None, mo_energy=None, mo_coeff=None, mo_occ=None)[source]

Bases: RAGF2

Restricted AGF2 with canonical HF reference with density fitting

Attributes:
verboseint

Print level. Default value equals to Mole.verbose

max_memoryfloat or int

Allowed memory in MB. Default value equals to Mole.max_memory

incore_completebool

Avoid all I/O. Default is False.

allow_lowmem_buildbool or str

Allow the self-energy build to switch to a serially slower code with lower thread-private memory overhead if needed. One of True, False or ‘force’. Default value is True.

conv_tolfloat

Convergence threshold for AGF2 energy. Default value is 1e-7

conv_tol_rdm1float

Convergence threshold for first-order reduced density matrix. Default value is 1e-8.

conv_tol_nelecfloat

Convergence threshold for the number of electrons. Default value is 1e-6.

max_cycleint

Maximum number of AGF2 iterations. Default value is 50.

max_cycle_outerint

Maximum number of outer Fock loop iterations. Default value is 20.

max_cycle_innerint

Maximum number of inner Fock loop iterations. Default value is 50.

weight_tolfloat

Threshold in spectral weight of auxiliaries to be considered zero. Default 1e-11.

diisbool or lib.diis.DIIS

Whether to use DIIS, can also be a lib.diis.DIIS object. Default value is True.

diis_spaceint

DIIS space size. Default value is 8.

diis_min_spaceint

Minimum space of DIIS. Default value is 1.

fock_diis_spaceint

DIIS space size for Fock loop iterations. Default value is 6.

fock_diis_min_space :

Minimum space of DIIS. Default value is 1.

os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

dampingfloat

Damping factor for the self-energy. Default value is 0.0

Saved results

e_corrfloat

AGF2 correlation energy

e_totfloat

Total energy (HF + correlation)

e_1bfloat

One-body part of e_tot

e_2bfloat

Two-body part of e_tot

e_initfloat

Initial correlation energy (truncated MP2)

convergedbool

Whether convergence was successful

seSelfEnergy

Auxiliaries of the self-energy

gfGreensFunction

Auxiliaries of the Green’s function

ao2mo(mo_coeff=None)[source]

Get the density-fitted electronic repulsion integrals in MO basis.

build_se_part(eri, gf_occ, gf_vir, os_factor=1.0, ss_factor=1.0)
Builds either the auxiliaries of the occupied self-energy,

or virtual if gf_occ and gf_vir are swapped.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gf_occGreensFunction

Occupied Green’s function

gf_virGreensFunction

Virtual Green’s function

Kwargs:
os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

Returns:

SelfEnergy

get_jk(eri, rdm1, with_j=True, with_k=True)

Get the J/K matrices.

Args:
erindarray or H5 dataset

Electronic repulsion integrals (NOT as _ChemistsERIs). In the case of no bra/ket symmetry, a tuple can be passed.

rdm12D array

Reduced density matrix

Kwargs:
with_jbool

Whether to compute J. Default value is True

with_kbool

Whether to compute K. Default value is True

Returns:

tuple of ndarrays corresponding to J and K, if either are not requested then they are set to None.

reset(mol=None)[source]
property with_df
pyscf.agf2.dfragf2.build_se_part(agf2, eri, gf_occ, gf_vir, os_factor=1.0, ss_factor=1.0)[source]
Builds either the auxiliaries of the occupied self-energy,

or virtual if gf_occ and gf_vir are swapped.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gf_occGreensFunction

Occupied Green’s function

gf_virGreensFunction

Virtual Green’s function

Kwargs:
os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

Returns:

SelfEnergy

pyscf.agf2.dfragf2.get_jk(agf2, eri, rdm1, with_j=True, with_k=True)[source]

Get the J/K matrices.

Args:
erindarray or H5 dataset

Electronic repulsion integrals (NOT as _ChemistsERIs). In the case of no bra/ket symmetry, a tuple can be passed.

rdm12D array

Reduced density matrix

Kwargs:
with_jbool

Whether to compute J. Default value is True

with_kbool

Whether to compute K. Default value is True

Returns:

tuple of ndarrays corresponding to J and K, if either are not requested then they are set to None.

pyscf.agf2.dfuagf2 module

Auxiliary second-order Green’s function perturbation theory for unrestricted references with density fitting

class pyscf.agf2.dfuagf2.DFUAGF2(mf, frozen=None, mo_energy=None, mo_coeff=None, mo_occ=None)[source]

Bases: UAGF2

Unrestricted AGF2 with canonical HF reference with density fitting

Attributes:
verboseint

Print level. Default value equals to Mole.verbose

max_memoryfloat or int

Allowed memory in MB. Default value equals to Mole.max_memory

incore_completebool

Avoid all I/O. Default is False.

allow_lowmem_buildbool

Allow the self-energy build to switch to a serially slower code with lower thread-private memory overhead if needed. One of True, False or ‘force’. Default value is True.

conv_tolfloat

Convergence threshold for AGF2 energy. Default value is 1e-7

conv_tol_rdm1float

Convergence threshold for first-order reduced density matrix. Default value is 1e-8.

conv_tol_nelecfloat

Convergence threshold for the number of electrons. Default value is 1e-6.

max_cycleint

Maximum number of AGF2 iterations. Default value is 50.

max_cycle_outerint

Maximum number of outer Fock loop iterations. Default value is 20.

max_cycle_innerint

Maximum number of inner Fock loop iterations. Default value is 50.

weight_tolfloat

Threshold in spectral weight of auxiliaries to be considered zero. Default 1e-11.

diisbool or lib.diis.DIIS

Whether to use DIIS, can also be a lib.diis.DIIS object. Default value is True.

diis_spaceint

DIIS space size. Default value is 8.

diis_min_spaceint

Minimum space of DIIS. Default value is 1.

fock_diis_spaceint

DIIS space size for Fock loop iterations. Default value is 6.

fock_diis_min_space :

Minimum space of DIIS. Default value is 1.

os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

dampingfloat

Damping factor for the self-energy. Default value is 0.0

Saved results

e_corrfloat

AGF2 correlation energy

e_totfloat

Total energy (HF + correlation)

e_1bfloat

One-body part of e_tot

e_2bfloat

Two-body part of e_tot

e_initfloat

Initial correlation energy (truncated MP2)

convergedbool

Whether convergence was successful

setuple of SelfEnergy

Auxiliaries of the self-energy for each spin

gftuple of GreensFunction

Auxiliaries of the Green’s function for each spin

ao2mo(mo_coeff=None)[source]

Get the density-fitted electronic repulsion integrals in MO basis.

build_se_part(eri, gf_occ, gf_vir, os_factor=1.0, ss_factor=1.0)
Builds either the auxiliaries of the occupied self-energy,

or virtual if gf_occ and gf_vir are swapped.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gf_occGreensFunction

Occupied Green’s function

gf_virGreensFunction

Virtual Green’s function

Kwargs:
os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

Returns:

SelfEnergy

get_jk(eri, rdm1, with_j=True, with_k=True)

Get the J/K matrices.

Args:
erindarray or H5 dataset

Electronic repulsion integrals (NOT as _ChemistsERIs). In the case of no bra/ket symmetry, a tuple can be passed.

rdm12D array

Reduced density matrix

Kwargs:
with_jbool

Whether to compute J. Default value is True

with_kbool

Whether to compute K. Default value is True

Returns:

tuple of ndarrays corresponding to J and K, if either are not requested then they are set to None.

reset(mol=None)[source]
property with_df
pyscf.agf2.dfuagf2.build_se_part(agf2, eri, gf_occ, gf_vir, os_factor=1.0, ss_factor=1.0)[source]
Builds either the auxiliaries of the occupied self-energy,

or virtual if gf_occ and gf_vir are swapped.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gf_occGreensFunction

Occupied Green’s function

gf_virGreensFunction

Virtual Green’s function

Kwargs:
os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

Returns:

SelfEnergy

pyscf.agf2.mpi_helper module

MPI helper functions using mpi4py

pyscf.agf2.mpi_helper.allreduce(sendbuf, root=0, op=None)[source]
pyscf.agf2.mpi_helper.allreduce_safe_inplace(array)[source]
pyscf.agf2.mpi_helper.barrier()[source]
pyscf.agf2.mpi_helper.bcast(buf, root=0)[source]
pyscf.agf2.mpi_helper.bcast_dict(buf, root=0)[source]
pyscf.agf2.mpi_helper.nrange(start, stop=None, step=1)[source]
pyscf.agf2.mpi_helper.prange(start, stop, step)[source]

lib.prange() distributed over MPI processes. Returns the range for a single MPI rank.

pyscf.agf2.mpi_helper.reduce(sendbuf, root=0, op=None)[source]

pyscf.agf2.ragf2 module

Auxiliary second-order Green’s function perturbation theory

class pyscf.agf2.ragf2.RAGF2(mf, frozen=None, mo_energy=None, mo_coeff=None, mo_occ=None)[source]

Bases: StreamObject

Restricted AGF2 with canonical HF reference

Attributes:
verboseint

Print level. Default value equals to Mole.verbose

max_memoryfloat or int

Allowed memory in MB. Default value equals to Mole.max_memory

incore_completebool

Avoid all I/O. Default is False.

conv_tolfloat

Convergence threshold for AGF2 energy. Default value is 1e-7

conv_tol_rdm1float

Convergence threshold for first-order reduced density matrix. Default value is 1e-8.

conv_tol_nelecfloat

Convergence threshold for the number of electrons. Default value is 1e-6.

max_cycleint

Maximum number of AGF2 iterations. Default value is 50.

max_cycle_outerint

Maximum number of outer Fock loop iterations. Default value is 20.

max_cycle_innerint

Maximum number of inner Fock loop iterations. Default value is 50.

weight_tolfloat

Threshold in spectral weight of auxiliaries to be considered zero. Default 1e-11.

diisbool or lib.diis.DIIS

Whether to use DIIS, can also be a lib.diis.DIIS object. Default value is True.

diis_spaceint

DIIS space size. Default value is 8.

diis_min_spaceint

Minimum space of DIIS. Default value is 1.

fock_diis_spaceint

DIIS space size for Fock loop iterations. Default value is 6.

fock_diis_min_spaceint

Minimum space of DIIS. Default value is 1.

os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

dampingfloat

Damping factor for the self-energy. Default value is 0.0

Saved results

e_corrfloat

AGF2 correlation energy

e_totfloat

Total energy (HF + correlation)

e_1bfloat

One-body part of e_tot

e_2bfloat

Two-body part of e_tot

e_initfloat

Initial correlation energy (truncated MP2)

convergedbool

Whether convergence was successful

seSelfEnergy

Auxiliaries of the self-energy

gfGreensFunction

Auxiliaries of the Green’s function

ao2mo(mo_coeff=None)[source]

Get the electronic repulsion integrals in MO basis.

async_io = True
build_gf(eri=None, gf=None, se=None)[source]
Builds the auxiliaries of the Green’s function by solving

the Dyson equation.

Kwargs:
eri_ChemistsERIs

Electronic repulsion integrals

gfGreensFunction

Auxiliaries of the Green’s function

seSelfEnergy

Auxiliaries of the self-energy

Returns:

GreensFunction

build_se(eri=None, gf=None, os_factor=None, ss_factor=None, se_prev=None)[source]

Builds the auxiliaries of the self-energy.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gfGreensFunction

Auxiliaries of the Green’s function

Kwargs:
os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

se_prevSelfEnergy

Previous self-energy for damping. Default value is None

Returns:

SelfEnergy

build_se_part(eri, gf_occ, gf_vir, os_factor=1.0, ss_factor=1.0)
Builds either the auxiliaries of the occupied self-energy,

or virtual if gf_occ and gf_vir are swapped.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gf_occGreensFunction

Occupied Green’s function

gf_virGreensFunction

Virtual Green’s function

Kwargs:
os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

Returns:

SelfEnergy

density_fit(auxbasis=None, with_df=None)[source]
dump_chk(chkfile=None, key='agf2', gf=None, se=None, frozen=None, nmom=None, mo_energy=None, mo_coeff=None, mo_occ=None)[source]
dump_flags(verbose=None)[source]
property e_corr
property e_tot
eaagf2(nroots=5)[source]
Find the (N+1)-electron charged excitations, corresponding

to the smallest nroots poles of the virtual Green’s function.

Kwargs:

See ipagf2()

energy_1body(eri, gf)

Calculates the one-body energy according to the RHF form.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gfGreensFunction

Auxiliaries of Green’s function

Returns:

One-body energy

energy_2body(gf, se)
Calculates the two-body energy using analytically integrated

Galitskii-Migdal formula. The formula is symmetric and only one side needs to be calculated.

Args:
gfGreensFunction

Auxiliaries of the Green’s function

seSelfEnergy

Auxiliaries of the self-energy

Returns

Two-body energy

energy_mp2(mo_energy=None, se=None)[source]
energy_nuc()[source]
fock_loop(eri, gf, se)
Self-consistent loop for the density matrix via the HF self-

consistent field.

Args:
eri_ChemistERIs

Electronic repulsion integrals

gfGreensFunction

Auxiliaries of the Green’s function

seSelfEnergy

Auxiliaries of the self-energy

Returns:

SelfEnergy, GreensFunction and a boolean indicating whether convergence was successful.

get_ea(gf, nroots=5)[source]
get_fock(eri=None, gf=None, rdm1=None)[source]

Computes the physical space Fock matrix in MO basis.

get_ip(gf, nroots=5)[source]
get_jk(eri, rdm1, with_j=True, with_k=True)

Get the J/K matrices.

Args:
erindarray or H5 dataset

Electronic repulsion integrals (NOT as _ChemistsERIs)

rdm12D array

Reduced density matrix

Kwargs:
with_jbool

Whether to compute J. Default value is True

with_kbool

Whether to compute K. Default value is True

Returns:

tuple of ndarrays corresponding to J and K, if either are not requested then they are set to None.

incore_complete = False
init_gf(frozen=False)[source]

Builds the Hartree-Fock Green’s function.

Returns:

GreensFunction, SelfEnergy

ipagf2(nroots=5)[source]
Find the (N-1)-electron charged excitations, corresponding

to the largest nroots poles of the occupied Green’s function.

Kwargs:
nrootsint

Number of roots (poles) requested. Default 1.

Returns:

IP and transition moment (float, 1D array) if nroots = 1, or array of IPs and moments (1D array, 2D array) if nroots > 1.

kernel(eri=None, gf=None, se=None, dump_chk=True)[source]

Kernel function is the main driver of a method. Every method should define the kernel function as the entry of the calculation. Note the return value of kernel function is not strictly defined. It can be anything related to the method (such as the energy, the wave-function, the DFT mesh grids etc.).

make_rdm1(gf=None)[source]

Computes the one-body reduced density matrix in MO basis.

Kwargs:
gfGreensFunction

Auxiliaries of the Green’s function

Returns:

ndarray of density matrix

property nmo
property nocc
property qmo_coeff

Gives the couplings in AO basis

property qmo_energy
property qmo_occ
reset(mol=None)[source]
run_diis(se, diis=None)[source]
Runs the direct inversion of the iterative subspace for the

self-energy.

Args:
seSelfEnergy

Auxiliaries of the self-energy

diislib.diis.DIIS

DIIS object

Returns:

SelfEnergy

update(chkfile=None, key='agf2')
update_from_chk(chkfile=None, key='agf2')
update_from_chk_(chkfile=None, key='agf2')[source]
pyscf.agf2.ragf2.build_se_part(agf2, eri, gf_occ, gf_vir, os_factor=1.0, ss_factor=1.0)[source]
Builds either the auxiliaries of the occupied self-energy,

or virtual if gf_occ and gf_vir are swapped.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gf_occGreensFunction

Occupied Green’s function

gf_virGreensFunction

Virtual Green’s function

Kwargs:
os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

Returns:

SelfEnergy

pyscf.agf2.ragf2.energy_1body(agf2, eri, gf)[source]

Calculates the one-body energy according to the RHF form.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gfGreensFunction

Auxiliaries of Green’s function

Returns:

One-body energy

pyscf.agf2.ragf2.energy_2body(agf2, gf, se)[source]
Calculates the two-body energy using analytically integrated

Galitskii-Migdal formula. The formula is symmetric and only one side needs to be calculated.

Args:
gfGreensFunction

Auxiliaries of the Green’s function

seSelfEnergy

Auxiliaries of the self-energy

Returns

Two-body energy

pyscf.agf2.ragf2.energy_mp2(agf2, mo_energy, se)[source]
Calculates the two-body energy using analytically integrated

Galitskii-Migdal formula for an MP2 self-energy. Per the definition of one- and two-body partitioning in the Dyson equation, this result is half of energy_2body().

Args:
gfGreensFunction

Auxiliaries of the Green’s function

seSelfEnergy

Auxiliaries of the self-energy

Returns

MP2 energy

pyscf.agf2.ragf2.fock_loop(agf2, eri, gf, se)[source]
Self-consistent loop for the density matrix via the HF self-

consistent field.

Args:
eri_ChemistERIs

Electronic repulsion integrals

gfGreensFunction

Auxiliaries of the Green’s function

seSelfEnergy

Auxiliaries of the self-energy

Returns:

SelfEnergy, GreensFunction and a boolean indicating whether convergence was successful.

pyscf.agf2.ragf2.get_fock(agf2, eri, gf=None, rdm1=None)[source]
Computes the physical space Fock matrix in MO basis. If rdm1

is not supplied, it is built from gf, which defaults to the mean-field Green’s function.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

Kwargs:
gfGreensfunction

Auxiliaries of the Green’s function

rdm12D array

Reduced density matrix.

Returns:

ndarray of physical space Fock matrix

pyscf.agf2.ragf2.get_frozen_mask(agf2)[source]
pyscf.agf2.ragf2.get_jk(agf2, eri, rdm1, with_j=True, with_k=True)[source]

Get the J/K matrices.

Args:
erindarray or H5 dataset

Electronic repulsion integrals (NOT as _ChemistsERIs)

rdm12D array

Reduced density matrix

Kwargs:
with_jbool

Whether to compute J. Default value is True

with_kbool

Whether to compute K. Default value is True

Returns:

tuple of ndarrays corresponding to J and K, if either are not requested then they are set to None.

pyscf.agf2.ragf2.kernel(agf2, eri=None, gf=None, se=None, verbose=None, dump_chk=True)[source]

pyscf.agf2.ragf2_slow module

Auxiliary second-order Green’s function perturbation theory for arbitrary moment consistency

class pyscf.agf2.ragf2_slow.RAGF2(mf, nmom=(None, 0), frozen=None, mo_energy=None, mo_coeff=None, mo_occ=None)[source]

Bases: RAGF2

Restricted AGF2 with canonical HF reference for arbitrary

moment consistency

Attributes:
nmomtuple of int

Compression level of the Green’s function and self-energy, respectively

verboseint

Print level. Default value equals to Mole.verbose

max_memoryfloat or int

Allowed memory in MB. Default value equals to Mole.max_memory

conv_tolfloat

Convergence threshold for AGF2 energy. Default value is 1e-7

conv_tol_rdm1float

Convergence threshold for first-order reduced density matrix. Default value is 1e-8.

conv_tol_nelecfloat

Convergence threshold for the number of electrons. Default value is 1e-6.

max_cycleint

Maximum number of AGF2 iterations. Default value is 50.

max_cycle_outerint

Maximum number of outer Fock loop iterations. Default value is 20.

max_cycle_innerint

Maximum number of inner Fock loop iterations. Default value is 50.

weight_tolfloat

Threshold in spectral weight of auxiliaries to be considered zero. Default 1e-11.

fock_diis_spaceint

DIIS space size for Fock loop iterations. Default value is 6.

fock_diis_min_space :

Minimum space of DIIS. Default value is 1.

os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

dampingfloat

Damping factor for the self-energy. Default value is 0.0

Saved results

e_corrfloat

AGF2 correlation energy

e_totfloat

Total energy (HF + correlation)

e_1bfloat

One-body part of e_tot

e_2bfloat

Two-body part of e_tot

e_initfloat

Initial correlation energy (truncated MP2)

convergedbool

Whether convergence was successful

seSelfEnergy

Auxiliaries of the self-energy

gfGreensFunction

Auxiliaries of the Green’s function

build_se(eri=None, gf=None, os_factor=None, ss_factor=None, se_prev=None)[source]

Builds the auxiliaries of the self-energy.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gfGreensFunction

Auxiliaries of the Green’s function

Kwargs:
os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

se_prevSelfEnergy

Previous self-energy for damping. Default value is None

Returns:

SelfEnergy

build_se_part(eri, gf_occ, gf_vir, os_factor=1.0, ss_factor=1.0)
Builds either the auxiliaries of the occupied self-energy,

or virtual if gf_occ and gf_vir are swapped.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gf_occGreensFunction

Occupied Green’s function

gf_virGreensFunction

Virtual Green’s function

Kwargs:
os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

Returns:

SelfEnergy

dump_flags(verbose=None)[source]
run_diis(se, diis=None)[source]
Runs the direct inversion of the iterative subspace for the

self-energy.

Args:
seSelfEnergy

Auxiliaries of the self-energy

diislib.diis.DIIS

DIIS object

Returns:

SelfEnergy

pyscf.agf2.ragf2_slow.build_se_part(agf2, eri, gf_occ, gf_vir, os_factor=1.0, ss_factor=1.0)[source]
Builds either the auxiliaries of the occupied self-energy,

or virtual if gf_occ and gf_vir are swapped.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gf_occGreensFunction

Occupied Green’s function

gf_virGreensFunction

Virtual Green’s function

Kwargs:
os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

Returns:

SelfEnergy

pyscf.agf2.uagf2 module

Auxiliary second-order Green’s function perturbation theory for unrestricted references

class pyscf.agf2.uagf2.UAGF2(mf, frozen=None, mo_energy=None, mo_coeff=None, mo_occ=None)[source]

Bases: RAGF2

Unrestricted AGF2 with canonical HF reference

Attributes:
verboseint

Print level. Default value equals to Mole.verbose

max_memoryfloat or int

Allowed memory in MB. Default value equals to Mole.max_memory

incore_completebool

Avoid all I/O. Default is False.

conv_tolfloat

Convergence threshold for AGF2 energy. Default value is 1e-7

conv_tol_rdm1float

Convergence threshold for first-order reduced density matrix. Default value is 1e-8.

conv_tol_nelecfloat

Convergence threshold for the number of electrons. Default value is 1e-6.

max_cycleint

Maximum number of AGF2 iterations. Default value is 50.

max_cycle_outerint

Maximum number of outer Fock loop iterations. Default value is 20.

max_cycle_innerint

Maximum number of inner Fock loop iterations. Default value is 50.

weight_tolfloat

Threshold in spectral weight of auxiliaries to be considered zero. Default 1e-11.

diisbool or lib.diis.DIIS

Whether to use DIIS, can also be a lib.diis.DIIS object. Default value is True.

diis_spaceint

DIIS space size. Default value is 8.

diis_min_spaceint

Minimum space of DIIS. Default value is 1.

fock_diis_spaceint

DIIS space size for Fock loop iterations. Default value is 6.

fock_diis_min_space :

Minimum space of DIIS. Default value is 1.

os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

dampingfloat

Damping factor for the self-energy. Default value is 0.0

Saved results

e_corrfloat

AGF2 correlation energy

e_totfloat

Total energy (HF + correlation)

e_1bfloat

One-body part of e_tot

e_2bfloat

Two-body part of e_tot

e_initfloat

Initial correlation energy (truncated MP2)

convergedbool

Whether convergence was successful

setuple of SelfEnergy

Auxiliaries of the self-energy for each spin

gftuple of GreensFunction

Auxiliaries of the Green’s function for each spin

ao2mo(mo_coeff=None)[source]

Get the electronic repulsion integrals in MO basis.

build_gf(eri=None, gf=None, se=None)[source]
Builds the auxiliaries of the Green’s functions by solving

the Dyson equation for each spin.

Kwargs:
eri_ChemistsERIs

Electronic repulsion integrals

gftuple of GreensFunction

Auxiliaries of the Green’s function for each spin

setuple of SelfEnergy

Auxiliaries of the self-energy for each spin

Returns:

tuple of GreensFunction

build_se(eri=None, gf=None, os_factor=None, ss_factor=None, se_prev=None)[source]

Builds the auxiliaries of the self-energy.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gftuple of GreensFunction

Auxiliaries of the Green’s function

Kwargs:
os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

se_prevSelfEnergy

Previous self-energy for damping. Default value is None

Returns

tuple of SelfEnergy

build_se_part(eri, gf_occ, gf_vir, os_factor=1.0, ss_factor=1.0)
Builds either the auxiliaries of the occupied self-energy,

or virtual if gf_occ and gf_vir are swapped, for a single spin.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gf_occtuple of GreensFunction

Occupied Green’s function for each spin

gf_virtuple of GreensFunction

Virtual Green’s function for each spin

Kwargs:
os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

Returns:

SelfEnergy

density_fit(auxbasis=None, with_df=None)[source]
eaagf2(nroots=5)[source]
Find the (N+1)-electron charged excitations, corresponding

to the smallest nroots poles of the virtual Green’s function.

Kwargs:

See ipagf2()

energy_1body(eri, gf)

Calculates the one-body energy according to the UHF form.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gftuple of GreensFunction

Auxiliaries of the Green’s function for each spin

Returns:

One-body energy

energy_2body(gf, se)
Calculates the two-body energy using analytically integrated

Galitskii-Migdal formula. The formula is symmetric and only one side needs to be calculated.

Args:
gftuple of GreensFunction

Auxiliaries of the Green’s function for each spin

setuple of SelfEnergy

Auxiliaries of the self-energy for each spin

Returns:

Two-body energy

energy_mp2(mo_energy=None, se=None)[source]
fock_loop(eri, gf, se)
Self-consistent loop for the density matrix via the HF self-

consistent field.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gftuple of GreensFunction

Auxiliaries of the Green’s function for each spin

setuple of SelfEnergy

Auxiliaries of the self-energy for each spin

Returns:

SelfEnergy, GreensFunction and a boolean indicating whether convergence was successful.

get_ea(gf, nroots=5)[source]
get_fock(eri=None, gf=None, rdm1=None)[source]

Computes the physical space Fock matrix in MO basis.

get_ip(gf, nroots=5)[source]
init_gf(frozen=False)[source]

Builds the Hartree-Fock Green’s function.

Returns:

tuple of GreensFunction, tuple of SelfEnergy

ipagf2(nroots=5)[source]
Find the (N-1)-electron charged excitations, corresponding

to the largest nroots poles of the occupied Green’s function.

Kwargs:
nrootsint

Number of roots (poles) requested. Default 1.

Returns:

IP and transition moment (float, 1D array) if nroots = 1, or array of IPs and moments (1D array, 2D array) if nroots > 1.

make_rdm1(gf=None)[source]

Compute the one-body reduced density matrix in MO basis.

Kwargs:
gftuple of GreensFunction

Auxiliaries of the Green’s functions for each spin

Returns:

tuple of ndarray of density matrices

property nmo
property nocc
property qmo_coeff

Gives the couplings in AO basis

property qmo_energy
property qmo_occ
run_diis(se, diis=None)[source]
Runs the direct inversion of the iterative subspace for the

self-energy.

Args:
seSelfEnergy

Auxiliaries of the self-energy

diislib.diis.DIIS

DIIS object

Returns:

tuple of SelfEnergy

pyscf.agf2.uagf2.build_se_part(agf2, eri, gf_occ, gf_vir, os_factor=1.0, ss_factor=1.0)[source]
Builds either the auxiliaries of the occupied self-energy,

or virtual if gf_occ and gf_vir are swapped, for a single spin.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gf_occtuple of GreensFunction

Occupied Green’s function for each spin

gf_virtuple of GreensFunction

Virtual Green’s function for each spin

Kwargs:
os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

Returns:

SelfEnergy

pyscf.agf2.uagf2.energy_1body(agf2, eri, gf)[source]

Calculates the one-body energy according to the UHF form.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gftuple of GreensFunction

Auxiliaries of the Green’s function for each spin

Returns:

One-body energy

pyscf.agf2.uagf2.energy_2body(agf2, gf, se)[source]
Calculates the two-body energy using analytically integrated

Galitskii-Migdal formula. The formula is symmetric and only one side needs to be calculated.

Args:
gftuple of GreensFunction

Auxiliaries of the Green’s function for each spin

setuple of SelfEnergy

Auxiliaries of the self-energy for each spin

Returns:

Two-body energy

pyscf.agf2.uagf2.energy_mp2(agf2, gf, se)[source]
Calculates the two-body energy using analytically integrated

Galitskii-Migdal formula for an MP2 self-energy. Per the definition of one- and two-body partitioning in the Dyson equation, this result is half of energy_2body().

Args:
gftuple of GreensFunction

Auxiliaries of the Green’s function for each spin

setuple of SelfEnergy

Auxiliaries of the self-energy for each spin

Returns:

MP2 energy

pyscf.agf2.uagf2.fock_loop(agf2, eri, gf, se)[source]
Self-consistent loop for the density matrix via the HF self-

consistent field.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gftuple of GreensFunction

Auxiliaries of the Green’s function for each spin

setuple of SelfEnergy

Auxiliaries of the self-energy for each spin

Returns:

SelfEnergy, GreensFunction and a boolean indicating whether convergence was successful.

pyscf.agf2.uagf2.get_fock(agf2, eri, gf=None, rdm1=None)[source]
Computes the physical space Fock matrix in MO basis. If rdm1

is not supplied, it is built from gf, which defaults to the mean-field Green’s function.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

Kwargs:
gfGreensFunction

Auxiliaries of the Green’s function

rdm12D array

Reduced density matrix

Returns:

ndarray of physical space Fock matrix

pyscf.agf2.uagf2.get_frozen_mask(agf2)[source]

pyscf.agf2.uagf2_slow module

Auxiliary second-order Green’s function perturbation theory for unrestricted references for arbitrary moment consistency

class pyscf.agf2.uagf2_slow.UAGF2(mf, nmom=(None, 0), frozen=None, mo_energy=None, mo_coeff=None, mo_occ=None)[source]

Bases: UAGF2

Unrestricted AGF2 with canonical HF reference for arbitrary

moment consistency

Attributes:
nmomtuple of int

Compression level of the Green’s function and self-energy, respectively

verboseint

Print level. Default value equals to Mole.verbose

max_memoryfloat or int

Allowed memory in MB. Default value equals to Mole.max_memory

conv_tolfloat

Convergence threshold for AGF2 energy. Default value is 1e-7

conv_tol_rdm1float

Convergence threshold for first-order reduced density matrix. Default value is 1e-8.

conv_tol_nelecfloat

Convergence threshold for the number of electrons. Default value is 1e-6.

max_cycleint

Maximum number of AGF2 iterations. Default value is 50.

max_cycle_outerint

Maximum number of outer Fock loop iterations. Default value is 20.

max_cycle_innerint

Maximum number of inner Fock loop iterations. Default value is 50.

weight_tolfloat

Threshold in spectral weight of auxiliaries to be considered zero. Default 1e-11.

fock_diis_spaceint

DIIS space size for Fock loop iterations. Default value is 6.

fock_diis_min_space :

Minimum space of DIIS. Default value is 1.

os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

dampingfloat

Damping factor for the self-energy. Default value is 0.0

Saved results

e_corrfloat

AGF2 correlation energy

e_totfloat

Total energy (HF + correlation)

e_1bfloat

One-body part of e_tot

e_2bfloat

Two-body part of e_tot

e_initfloat

Initial correlation energy (truncated MP2)

convergedbool

Whether convergence was successful

setuple of SelfEnergy

Auxiliaries of the self-energy for each spin

gftuple of GreensFunction

Auxiliaries of the Green’s function for each spin

build_se(eri=None, gf=None, os_factor=None, ss_factor=None, se_prev=None)[source]

Builds the auxiliaries of the self-energy.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gftuple of GreensFunction

Auxiliaries of the Green’s function

Kwargs:
os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

se_prevSelfEnergy

Previous self-energy for damping. Default value is None

Returns

SelfEnergy

build_se_part(eri, gf_occ, gf_vir, os_factor=1.0, ss_factor=1.0)
Builds either the auxiliaries of the occupied self-energy,

or virtual if gf_occ and gf_vir are swapped, for a single spin.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gf_occtuple of GreensFunction

Occupied Green’s function for each spin

gf_virtuple of GreensFunction

Virtual Green’s function for each spin

Kwargs:
os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

Returns:

SelfEnergy

dump_flags(verbose=None)[source]
run_diis(se, diis=None)[source]
Runs the direct inversion of the iterative subspace for the

self-energy.

Args:
seSelfEnergy

Auxiliaries of the self-energy

diislib.diis.DIIS

DIIS object

Returns:

tuple of SelfEnergy

pyscf.agf2.uagf2_slow.build_se_part(agf2, eri, gf_occ, gf_vir, os_factor=1.0, ss_factor=1.0)[source]
Builds either the auxiliaries of the occupied self-energy,

or virtual if gf_occ and gf_vir are swapped, for a single spin.

Args:
eri_ChemistsERIs

Electronic repulsion integrals

gf_occtuple of GreensFunction

Occupied Green’s function for each spin

gf_virtuple of GreensFunction

Virtual Green’s function for each spin

Kwargs:
os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

Returns:

SelfEnergy

Module contents

Auxiliary second-order Green’s function perturbation therory

The AGF2 method permits the computation of quasiparticle excitations and ground-state properties at the AGF2(None,0) level.

When using results of this code for publications, please cite the follow papers: “Wave function perspective and efficient truncation of renormalized second-order perturbation theory”, O. J. Backhouse, M. Nusspickel and G. H. Booth, J. Chem. Theory Comput., 16, 1090 (2020). “Efficient excitations and spectra within a perturbative renormalization approach”, O. J. Backhouse and G. H. Booth, J. Chem. Theory Comput., 16, 6294 (2020).

Simple usage:

>>> from pyscf import gto, scf, agf2
>>> mol = gto.M(atom='H 0 0 0; H 0 0 1')
>>> mf = scf.RHF(mol).run()
>>> gf2 = agf2.AGF2(mf).run()
>>> gf2.ipagf2()

agf2.AGF2() returns an instance of AGF2 class. Valid parameters to control the AGF2 method are:

verboseint

Print level. Default value equals to Mole.verbose

max_memoryfloat or int

Allowed memory in MB. Default value equals to Mole.max_memory

conv_tolfloat

Convergence threshold for AGF2 energy. Default value is 1e-7

conv_tol_rdm1float

Convergence threshold for first-order reduced density matrix. Default value is 1e-6.

conv_tol_nelecfloat

Convergence threshold for the number of electrons. Default value is 1e-6.

max_cycleint

Maximum number of AGF2 iterations. Default value is 50.

max_cycle_outerint

Maximum number of outer Fock loop iterations. Default value is 20.

max_cycle_innerint

Maximum number of inner Fock loop iterations. Default value is 50.

weight_tolfloat

Threshold in spectral weight of auxiliaries to be considered zero. Default 1e-11.

diis_spaceint

DIIS space size for Fock loop iterations. Default value is 6.

diis_min_space

Minimum space of DIIS. Default value is 1.

Saved result

e_corrfloat

AGF2 correlation energy

e_totfloat

Total energy (HF + correlation)

e_1bfloat

One-body part of e_tot

e_2bfloat

Two-body part of e_tot

e_initfloat

Initial correlation energy (truncated MP2)

convergedbool

Whether convergence was successful

seSelfEnergy

Auxiliaries of the self-energy

gfGreensFunction

Auxiliaries of the Green’s function

pyscf.agf2.AGF2(mf, nmom=(None, 0), frozen=None, mo_energy=None, mo_coeff=None, mo_occ=None)[source]

Restricted AGF2 with canonical HF reference

Attributes:
verboseint

Print level. Default value equals to Mole.verbose

max_memoryfloat or int

Allowed memory in MB. Default value equals to Mole.max_memory

incore_completebool

Avoid all I/O. Default is False.

conv_tolfloat

Convergence threshold for AGF2 energy. Default value is 1e-7

conv_tol_rdm1float

Convergence threshold for first-order reduced density matrix. Default value is 1e-8.

conv_tol_nelecfloat

Convergence threshold for the number of electrons. Default value is 1e-6.

max_cycleint

Maximum number of AGF2 iterations. Default value is 50.

max_cycle_outerint

Maximum number of outer Fock loop iterations. Default value is 20.

max_cycle_innerint

Maximum number of inner Fock loop iterations. Default value is 50.

weight_tolfloat

Threshold in spectral weight of auxiliaries to be considered zero. Default 1e-11.

diisbool or lib.diis.DIIS

Whether to use DIIS, can also be a lib.diis.DIIS object. Default value is True.

diis_spaceint

DIIS space size. Default value is 8.

diis_min_spaceint

Minimum space of DIIS. Default value is 1.

fock_diis_spaceint

DIIS space size for Fock loop iterations. Default value is 6.

fock_diis_min_spaceint

Minimum space of DIIS. Default value is 1.

os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

dampingfloat

Damping factor for the self-energy. Default value is 0.0

Saved results

e_corrfloat

AGF2 correlation energy

e_totfloat

Total energy (HF + correlation)

e_1bfloat

One-body part of e_tot

e_2bfloat

Two-body part of e_tot

e_initfloat

Initial correlation energy (truncated MP2)

convergedbool

Whether convergence was successful

seSelfEnergy

Auxiliaries of the self-energy

gfGreensFunction

Auxiliaries of the Green’s function

pyscf.agf2.RAGF2(mf, nmom=(None, 0), frozen=None, mo_energy=None, mo_coeff=None, mo_occ=None)[source]

Restricted AGF2 with canonical HF reference

Attributes:
verboseint

Print level. Default value equals to Mole.verbose

max_memoryfloat or int

Allowed memory in MB. Default value equals to Mole.max_memory

incore_completebool

Avoid all I/O. Default is False.

conv_tolfloat

Convergence threshold for AGF2 energy. Default value is 1e-7

conv_tol_rdm1float

Convergence threshold for first-order reduced density matrix. Default value is 1e-8.

conv_tol_nelecfloat

Convergence threshold for the number of electrons. Default value is 1e-6.

max_cycleint

Maximum number of AGF2 iterations. Default value is 50.

max_cycle_outerint

Maximum number of outer Fock loop iterations. Default value is 20.

max_cycle_innerint

Maximum number of inner Fock loop iterations. Default value is 50.

weight_tolfloat

Threshold in spectral weight of auxiliaries to be considered zero. Default 1e-11.

diisbool or lib.diis.DIIS

Whether to use DIIS, can also be a lib.diis.DIIS object. Default value is True.

diis_spaceint

DIIS space size. Default value is 8.

diis_min_spaceint

Minimum space of DIIS. Default value is 1.

fock_diis_spaceint

DIIS space size for Fock loop iterations. Default value is 6.

fock_diis_min_spaceint

Minimum space of DIIS. Default value is 1.

os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

dampingfloat

Damping factor for the self-energy. Default value is 0.0

Saved results

e_corrfloat

AGF2 correlation energy

e_totfloat

Total energy (HF + correlation)

e_1bfloat

One-body part of e_tot

e_2bfloat

Two-body part of e_tot

e_initfloat

Initial correlation energy (truncated MP2)

convergedbool

Whether convergence was successful

seSelfEnergy

Auxiliaries of the self-energy

gfGreensFunction

Auxiliaries of the Green’s function

pyscf.agf2.UAGF2(mf, nmom=(None, 0), frozen=None, mo_energy=None, mo_coeff=None, mo_occ=None)[source]

Unrestricted AGF2 with canonical HF reference

Attributes:
verboseint

Print level. Default value equals to Mole.verbose

max_memoryfloat or int

Allowed memory in MB. Default value equals to Mole.max_memory

incore_completebool

Avoid all I/O. Default is False.

conv_tolfloat

Convergence threshold for AGF2 energy. Default value is 1e-7

conv_tol_rdm1float

Convergence threshold for first-order reduced density matrix. Default value is 1e-8.

conv_tol_nelecfloat

Convergence threshold for the number of electrons. Default value is 1e-6.

max_cycleint

Maximum number of AGF2 iterations. Default value is 50.

max_cycle_outerint

Maximum number of outer Fock loop iterations. Default value is 20.

max_cycle_innerint

Maximum number of inner Fock loop iterations. Default value is 50.

weight_tolfloat

Threshold in spectral weight of auxiliaries to be considered zero. Default 1e-11.

diisbool or lib.diis.DIIS

Whether to use DIIS, can also be a lib.diis.DIIS object. Default value is True.

diis_spaceint

DIIS space size. Default value is 8.

diis_min_spaceint

Minimum space of DIIS. Default value is 1.

fock_diis_spaceint

DIIS space size for Fock loop iterations. Default value is 6.

fock_diis_min_space :

Minimum space of DIIS. Default value is 1.

os_factorfloat

Opposite-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

ss_factorfloat

Same-spin factor for spin-component-scaled (SCS) calculations. Default 1.0

dampingfloat

Damping factor for the self-energy. Default value is 0.0

Saved results

e_corrfloat

AGF2 correlation energy

e_totfloat

Total energy (HF + correlation)

e_1bfloat

One-body part of e_tot

e_2bfloat

Two-body part of e_tot

e_initfloat

Initial correlation energy (truncated MP2)

convergedbool

Whether convergence was successful

setuple of SelfEnergy

Auxiliaries of the self-energy for each spin

gftuple of GreensFunction

Auxiliaries of the Green’s function for each spin