Integrals

Exported functions

Two-index integrals

ElemCo.Integrals.eri_2e2idx!Method
eri_2e2idx!(out, basis::BasisSet)

Compute the two-electron 2-index electron-repulsion integral matrix. The result is stored in out.

source

Three-index integrals

ElemCo.Integrals.eri_2e3idx!Method
eri_2e3idx!(out, buffer, batch::BasisBatch)

Compute the two-electron three-index electron-repulsion integral in batch mode (see BasisBatcher). The result is stored in out. buffer is a preallocated buffer [MAlloc]Buffer{Cdouble} of size buffer_size_3idx(batch.bb).

source
ElemCo.Integrals.eri_2e3idx!Method
eri_2e3idx!(out, ao_basis::BasisSet, fit_basis::BasisSet)

Compute the two-electron three-index electron-repulsion integral. The result is stored in out.

source

Calculation in batches

ElemCo.Integrals.BasisBatchType
BasisBatch

A structure to represent a batch of basis functions. The structure is used to loop over basis functions in a batched manner. The structure is used in conjunction with BasisBatcher.

source
ElemCo.Integrals.BasisBatcherType
BasisBatcher

A structure to loop of basis functions in a batched manner. This is useful for computing integrals in batches. The structure is initialized with two basis sets and provides a method to loop over all basis functions of the second basis in a batched manner:

function BasisBatcher(basis1::BasisSet, basis2::BasisSet, target_length::Int=100)

The batch length is determined by the target_length parameter. The actual batch length depends on the number of basis functions in the shells of the second basis set and can be smaller or larger than the target_length. One can use the max_batch_length function to determine the maximum batch length.

Example

bbatches = BasisBatcher(basis1, basis2)
for batch in bbatches
  range = range(batch) # range of basis functions in the batch
  eri_2e3idx!(@view(pqP[:,:,range]), batch)
end
source
ElemCo.Integrals.buffer_size_3idxMethod
buffer_size_3idx(bb::BasisBatcher)

Return the buffer size needed in the 3-index integral calculation.

The buffer has to be of type [MAlloc]Buffer{Cdouble}(lenbuf) or [MAlloc]ThreadsBuffer{Cdouble}(lenbuf).

source

Internal functions

Base.rangeMethod
range(bb::BasisBatch)

Return the range of basis functions in the batch.

source
ElemCo.Integrals.eri_2e2idx_cart!Method
eri_2e2idx!(out, ash1::AngularShell, ash2::AngularShell, basis::BasisSet)

Compute the two-electron 2-index electron-repulsion integral between two angular shells. The result is stored in out.

source
ElemCo.Integrals.eri_2e2idx_cart!Method
eri_2e2idx_cart!(out, i::Int, j::Int, basis::BasisSet)

Compute the two-electron 2-index electron-repulsion integral between two angular shells. The result is stored in out.

source
ElemCo.Integrals.eri_2e2idx_cartMethod
eri_2e2idx_cart(ash1::AngularShell, ash2::AngularShell, basis::BasisSet)

Compute the two-electron 2-index electron-repulsion integral between two angular shells.

source
ElemCo.Integrals.eri_2e2idx_sph!Method
eri_2e2idx!(out, ash1::AngularShell, ash2::AngularShell, basis::BasisSet)

Compute the two-electron 2-index electron-repulsion integral between two angular shells. The result is stored in out.

source
ElemCo.Integrals.eri_2e2idx_sph!Method
eri_2e2idx_sph!(out, i::Int, j::Int, basis::BasisSet)

Compute the two-electron 2-index electron-repulsion integral between two angular shells. The result is stored in out.

source
ElemCo.Integrals.eri_2e2idx_sphMethod
eri_2e2idx_sph(ash1::AngularShell, ash2::AngularShell, basis::BasisSet)

Compute the two-electron 2-index electron-repulsion integral between two angular shells.

source
ElemCo.Integrals.eri_2e3idx_cart!Method
eri_2e3idx!(out, ash1ao::AngularShell, ash2ao::AngularShell, ashfit::AngularShell, basis::BasisSet)

Compute the two-electron three-index electron-repulsion integral $v_{a_1}^{a_2 P}$ for given angular shells. basis has to contain ao and fit bases. The result is stored in out.

source
ElemCo.Integrals.eri_2e3idx_cart!Method
eri_2e3idx_cart!(out, i::Int, j::Int, P::Int, basis::BasisSet)

Compute the two-electron three-index electron-repulsion integral $v_i^{j P}$ for given angular shells. basis has to contain ao and fit bases. The result is stored in out.

source
ElemCo.Integrals.eri_2e3idx_cartMethod
eri_2e3idx_cart(ash1ao::AngularShell, ash2ao::AngularShell, ashfit::AngularShell, basis::BasisSet)

Compute the two-electron three-index electron-repulsion integral $v_{a_1}^{a_2 P}$ for given angular shells. basis has to contain ao and fit bases.

source
ElemCo.Integrals.eri_2e3idx_sph!Method
eri_2e3idx!(out, ash1ao::AngularShell, ash2ao::AngularShell, ashfit::AngularShell, basis::BasisSet)

Compute the two-electron three-index electron-repulsion integral $v_{a_1}^{a_2 P}$ for given angular shells. basis has to contain ao and fit bases. The result is stored in out.

source
ElemCo.Integrals.eri_2e3idx_sph!Method
eri_2e3idx_sph!(out, i::Int, j::Int, P::Int, basis::BasisSet)

Compute the two-electron three-index electron-repulsion integral $v_i^{j P}$ for given angular shells. basis has to contain ao and fit bases. The result is stored in out.

source
ElemCo.Integrals.eri_2e3idx_sphMethod
eri_2e3idx_sph(ash1ao::AngularShell, ash2ao::AngularShell, ashfit::AngularShell, basis::BasisSet)

Compute the two-electron three-index electron-repulsion integral $v_{a_1}^{a_2 P}$ for given angular shells. basis has to contain ao and fit bases.

source
ElemCo.Integrals.kinetic_cart!Method
kinetic!(out, ash1::AngularShell, ash2::AngularShell, basis::BasisSet)

Compute the kinetic integral between two angular shells. The result is stored in out.

source
ElemCo.Integrals.kinetic_cart!Method
kinetic_cart!(out, i::Int, j::Int, basis::BasisSet)

Compute the kinetic integral between two angular shells. The result is stored in out.

source
ElemCo.Integrals.kinetic_sph!Method
kinetic!(out, ash1::AngularShell, ash2::AngularShell, basis::BasisSet)

Compute the kinetic integral between two angular shells. The result is stored in out.

source
ElemCo.Integrals.kinetic_sph!Method
kinetic_sph!(out, i::Int, j::Int, basis::BasisSet)

Compute the kinetic integral between two angular shells. The result is stored in out.

source
ElemCo.Integrals.nuclear_cart!Method
nuclear!(out, ash1::AngularShell, ash2::AngularShell, basis::BasisSet)

Compute the nuclear integral between two angular shells. The result is stored in out.

source
ElemCo.Integrals.nuclear_cart!Method
nuclear_cart!(out, i::Int, j::Int, basis::BasisSet)

Compute the nuclear integral between two angular shells. The result is stored in out.

source
ElemCo.Integrals.nuclear_sph!Method
nuclear!(out, ash1::AngularShell, ash2::AngularShell, basis::BasisSet)

Compute the nuclear integral between two angular shells. The result is stored in out.

source
ElemCo.Integrals.nuclear_sph!Method
nuclear_sph!(out, i::Int, j::Int, basis::BasisSet)

Compute the nuclear integral between two angular shells. The result is stored in out.

source
ElemCo.Integrals.overlap_cart!Method
overlap!(out, ash1::AngularShell, ash2::AngularShell, basis::BasisSet)

Compute the overlap integral between two angular shells. The result is stored in out.

source
ElemCo.Integrals.overlap_cart!Method
overlap_cart!(out, i::Int, j::Int, basis::BasisSet)

Compute the overlap integral between two angular shells. The result is stored in out.

source
ElemCo.Integrals.overlap_sph!Method
overlap!(out, ash1::AngularShell, ash2::AngularShell, basis::BasisSet)

Compute the overlap integral between two angular shells. The result is stored in out.

source
ElemCo.Integrals.overlap_sph!Method
overlap_sph!(out, i::Int, j::Int, basis::BasisSet)

Compute the overlap integral between two angular shells. The result is stored in out.

source

Integral libraries

ElemCo.LibcintModule
Libcint

Minimal wrap around the integral library libcint. This module exposes libcint functions to the Julia interface.

(adapted from GaussianBasis.jl)

source