3.3.12. horton/gbasis/ints.h
– Evaluation of integrals of Gaussian basis functions¶
-
class
GB2Integral
¶ Inherits from GBCalculator
Subclassed by GB2AttractionIntegral, GB2KineticIntegral, GB2MomentIntegral, GB2OverlapIntegral
Public Functions
-
GB2Integral
(long max_shell_type)¶
-
void
reset
(long shell_type0, long shell_type1, const double *r0, const double *r1)¶
-
virtual void
add
(double coeff, double alpha0, double alpha1, const double *scales0, const double *scales1) = 0¶
-
void
cart_to_pure
()¶
-
const long
get_shell_type0
() const¶
-
const long
get_shell_type1
() const¶
-
-
class
GB2OverlapIntegral
¶ - #include <ints.h>
Compute the overlap integrals in a Gaussian orbital basis.
Inherits from GB2Integral
-
class
GB2KineticIntegral
¶ - #include <ints.h>
Compute the kinetic integrals in a Gaussian orbital basis.
Inherits from GB2Integral
-
class
GB2AttractionIntegral
¶ - #include <ints.h>
Compute the nuclear attraction integrals in a Gaussian orbital basis.
Inherits from GB2Integral
Subclassed by GB2ErfAttractionIntegral, GB2GaussAttractionIntegral, GB2NuclearAttractionIntegral
Public Functions
-
GB2AttractionIntegral
(long max_shell_type, double *charges, double *centers, long ncharge)¶ Initialize a GB2AttractionIntegral object.
- Parameters
max_shell_type
: Highest angular momentum index to be expected in the reset method.
-
~GB2AttractionIntegral
()¶
-
void
add
(double coeff, double alpha0, double alpha1, const double *scales0, const double *scales1)¶ Add results for a combination of Cartesian primitive shells to the work array.
- Parameters
coeff
: Product of the contraction coefficients of the four primitives.alpha0
: The exponent of primitive shell 0.alpha1
: The exponent of primitive shell 1.scales0
: The normalization prefactors for basis functions in primitive shell 0scales1
: The normalization prefactors for basis functions in primitive shell 1
-
virtual void
laplace_of_potential
(double gamma, double arg, long mmax, double *output) = 0¶ Evaluate the Laplace transform of the the potential applied to nuclear attraction terms.
For theoretical details and the precise definition of the Laplace transform, we refer to the following paper:
Ahlrichs, R. A simple algebraic derivation of the Obara-Saika scheme for general two-electron interaction potentials. Phys. Chem. Chem. Phys. 8, 3072–3077 (2006). 10.1039/B605188J
For the general definition of this transform, see Eq. (8) in the reference above. Section 5 contains solutions of the Laplace transform for several popular cases.
- Parameters
gamma
: Sum of the exponents of the two gaussian functions involved in the integral. Similar to the first term in Eq. (3) in Ahlrichs’ paper.arg
: Rescaled distance between the two centers obtained from the application of the Gaussian product theorem. Equivalent to Eq. (5) in Ahlrichs’ paper.mmax
: Maximum derivative of the Laplace transform to be considered.output
: Output array. The size must be at least mmax + 1.
Private Members
-
double *
charges
¶ Array with values of the nuclear charges.
-
double *
centers
¶ The centers where the charges are located.
-
long
ncharge
¶ Number of nuclear charges.
-
double *
work_g0
¶ Temporary array to store intermediate results.
-
double *
work_g1
¶ Temporary array to store intermediate results.
-
double *
work_g2
¶ Temporary array to store intermediate results.
-
double *
work_boys
¶ Temporary array to store the laplace of the interaction potential.
-
-
class
GB2NuclearAttractionIntegral
¶ - #include <ints.h>
Nuclear Electron Attraction two-center integrals.
The potential is 1/r.
Inherits from GB2AttractionIntegral
Public Functions
-
GB2NuclearAttractionIntegral
(long max_shell_type, double *charges, double *centers, long ncharge)¶ Initialize a GB2NuclearAttractionIntegral object.
- Parameters
max_shell_type
: Highest angular momentum index to be expected in the reset method.charges
: Array with values of the charges.centers
: The centers [[C1_x, C1_y, C1_z],[…],] around which the moment integrals are computed.ncharge
: Number of nuclear charges.
-
void
laplace_of_potential
(double gamma, double arg, long mmax, double *output)¶ Evaluate the Laplace transform of the ordinary Coulomb potential.
See Eq. (39) in Ahlrichs’ paper. This is basically a rescaled Boys function.
See base class for more details.
-
-
class
GB2ErfAttractionIntegral
¶ - #include <ints.h>
Short-range electron repulsion four-center integrals.
The potential is erf(mu*r)/r.
Inherits from GB2AttractionIntegral
Public Functions
-
GB2ErfAttractionIntegral
(long max_shell_type, double *charges, double *centers, long ncharge, double mu)¶ Initialize a GB2ErfAttractionIntegral object.
- Parameters
max_shell_type
: Highest angular momentum index to be expected in the reset method.charges
: Array with values of the charges.centers
: The centers [[C1_x, C1_y, C1_z],[…],] around which the moment integrals are computed.ncharge
: Number of nuclear charges.mu
: The range-separation parameter.
-
void
laplace_of_potential
(double gamma, double arg, long mmax, double *output)¶ Evaluate the Laplace transform of the long-range Coulomb potential. (The short-range part is damped away using an error function.) See (52) in Ahlrichs’ paper.
See base class for more details.
-
const double
get_mu
() const¶ The range-separation parameter.
Private Members
-
double
mu
¶ The range-separation parameter.
-
-
class
GB2GaussAttractionIntegral
¶ - #include <ints.h>
Gaussian nuclear electron attraction two-center integrals.
The potential is c exp(-alpha r^2).
Inherits from GB2AttractionIntegral
Public Functions
-
GB2GaussAttractionIntegral
(long max_shell_type, double *charges, double *centers, long ncharge, double c, double alpha)¶ Initialize a GB2GaussAttractionIntegral object.
- Parameters
max_shell_type
: Highest angular momentum index to be expected in the reset method.charges
: Array with values of the charges.centers
: The centers [[C1_x, C1_y, C1_z],[…],] around which the moment integrals are computed.ncharge
: Number of nuclear charges.c
: Coefficient of the gaussian.alpha
: Exponential parameter of the gaussian.
-
void
laplace_of_potential
(double gamma, double arg, long mmax, double *output)¶ Evaluate the Laplace transform of the Gaussian potential.
See Ahlrichs’ paper for details. This type of potential is used in the papers of P.M.W Gill et al. and J. Toulouse et al.:
Gill, P. M. W., & Adamson, R. D. (1996). A family of attenuated Coulomb operators. Chem. Phys. Lett., 261(1-2), 105–110. http://doi.org/10.1016/0009-2614(96)00931-1
Toulouse, J., Colonna, F., & Savin, A. (2004). Long-range-short-range separation of the electron-electron interaction in density-functional theory. Phys. Rev. A, 70, 62505. http://doi.org/10.1103/PhysRevA.70.062505
See base class for more details.
-
const double
get_c
() const¶ Coefficient of the gaussian.
-
const double
get_alpha
() const¶ Exponential parameter of the gaussian.
-
-
class
GB2MomentIntegral
¶ - #include <ints.h>
Compute the (multipole) moment integrals in a Gaussian orbital basis. < gto_a | (x - C_x)^l (y - C_y)^m (z - C_z)^n | gto_b >.
Inherits from GB2Integral
Public Functions
-
GB2MomentIntegral
(long max_shell_type, long *xyz, double *center)¶ Initialize Moment integral calculator.
- Parameters
max_shell_type
: The highest angular momentum index suportedxyz
: The powers of x,y,z in the integrals (l, m, n).center
: The center [C_x, C_y, C_z] around which the moment integrals arecomputed
-
void
add
(double coeff, double alpha0, double alpha1, const double *scales0, const double *scales1)¶ Add integrals for a pair of primite shells to the current contraction.
- Parameters
coeff
: The contraction coefficient for the current primitive.alpha0
: The exponent of the primitive shell 0.alpha1
: The exponent of the primitive shell 1.scales0
: The normalization constants for the basis functions in primitive shell 0.scales1
: The normalization constants for the basis functions in primitive shell 1.
-
-
class
GB4Integral
¶ - #include <ints.h>
Base class for four-center integrals.
Inherits from GBCalculator
Subclassed by GB4IntegralLibInt
Public Functions
-
GB4Integral
(long max_shell_type)¶ Initialize a GB4Integral object.
- Parameters
max_shell_type
: Highest angular momentum index to be expected in the reset method.
-
void
reset
(long shell_type0, long shell_type1, long shell_type2, long shell_type3, const double *r0, const double *r1, const double *r2, const double *r3)¶ Set internal parameters for a new group of four contractions.
- Parameters
shell_type0
: Angular momentum index for contraction 0.shell_type1
: Angular momentum index for contraction 1.shell_type2
: Angular momentum index for contraction 2.shell_type3
: Angular momentum index for contraction 3.r0
: Cartesian coordinates of center 0.r1
: Cartesian coordinates of center 1.r2
: Cartesian coordinates of center 2.r3
: Cartesian coordinates of center 3.
-
virtual void
add
(double coeff, double alpha0, double alpha1, double alpha2, double alpha3, const double *scales0, const double *scales1, const double *scales2, const double *scales3) = 0¶ Add results for a combination of Cartesian primitive shells to the work array.
- Parameters
coeff
: Product of the contraction coefficients of the four primitives.alpha0
: The exponent of primitive shell 0.alpha1
: The exponent of primitive shell 1.alpha2
: The exponent of primitive shell 2.alpha3
: The exponent of primitive shell 3.scales0
: The normalization prefactors for basis functions in primitive shell 0scales1
: The normalization prefactors for basis functions in primitive shell 1scales2
: The normalization prefactors for basis functions in primitive shell 2scales3
: The normalization prefactors for basis functions in primitive shell 3
-
void
cart_to_pure
()¶ Transform the results in the work array from Cartesian to pure functions where needed.
-
const long
get_shell_type0
() const¶ Shell type of contraction 0.
-
const long
get_shell_type1
() const¶ Shell type of contraction 1.
-
const long
get_shell_type2
() const¶ Shell type of contraction 2.
-
const long
get_shell_type3
() const¶ Shell type of contraction 3.
Protected Attributes
-
long
shell_type0
¶ Shell type of contraction 0.
-
long
shell_type1
¶ Shell type of contraction 1.
-
long
shell_type2
¶ Shell type of contraction 2.
-
long
shell_type3
¶ Shell type of contraction 3.
-
const double *
r0
¶ Center of contraction 0.
-
const double *
r1
¶ Center of contraction 1.
-
const double *
r2
¶ Center of contraction 2.
-
const double *
r3
¶ Center of contraction 3.
-
-
struct
libint_arg_t
¶ - #include <ints.h>
Arguments associated with one primitive shell in LibInt conventions.
-
class
GB4IntegralLibInt
¶ - #include <ints.h>
Base class for four-center integrals that use LibInt.
Inherits from GB4Integral
Subclassed by GB4ElectronRepulsionIntegralLibInt, GB4ErfIntegralLibInt, GB4GaussIntegralLibInt, GB4RAlphaIntegralLibInt
Public Functions
-
GB4IntegralLibInt
(long max_shell_type)¶ Initialize a GB4IntegralLibInt object.
- Parameters
max_shell_type
: Highest angular momentum index to be expected in the reset method.
-
~GB4IntegralLibInt
()¶
-
void
reset
(long shell_type0, long shell_type1, long shell_type2, long shell_type3, const double *r0, const double *r1, const double *r2, const double *r3)¶ Set internal parameters for a new group of four contractions.
See base class for details.
-
void
add
(double coeff, double alpha0, double alpha1, double alpha2, double alpha3, const double *scales0, const double *scales1, const double *scales2, const double *scales3)¶ Add results for a combination of Cartesian primitive shells to the work array.
See base class for details.
-
virtual void
laplace_of_potential
(double prefac, double rho, double t, long mmax, double *output) = 0¶ Evaluate the Laplace transform of the the potential.
For theoretical details and the precise definition of the Laplace transform, we refer to the following paper:
Ahlrichs, R. A simple algebraic derivation of the Obara-Saika scheme for general two-electron interaction potentials. Phys. Chem. Chem. Phys. 8, 3072–3077 (2006). 10.1039/B605188J
For the general definition of this transform, see Eq. (8) in the reference above. Section 5 contains solutions of the Laplace transform for several popular cases.
- Parameters
prefac
: Prefactor with which all results in the output array are multiplied.rho
: See Eq. (3) in Ahlrichs’ paper.t
: Rescaled distance between the two centers obtained from the application of the Gaussian product theorem. See Eq. (5) in Ahlrichs’ paper.mmax
: Maximum derivative of the Laplace transform to be considered.output
: Output array. The size must be at least mmax + 1.
Private Members
-
Libint_eri_t
erieval
¶ LibInt runtime object.
-
libint_arg_t
libint_args
[4]¶ Arguments (shell info) for libint.
-
long
order
[4]¶ Re-ordering of shells for compatibility with LibInt.
-
double
ab
[3]¶ Relative vector from shell 2 to 0 (LibInt order).
-
double
cd
[3]¶ Relative vector from shell 3 to 1 (LibInt order).
-
double
ab2
¶ Norm squared of ab.
-
double
cd2
¶ Norm squared of cd.
-
-
class
GB4ElectronRepulsionIntegralLibInt
¶ - #include <ints.h>
Electron repulsion four-center integrals.
The potential is 1/r.
Inherits from GB4IntegralLibInt
Public Functions
-
GB4ElectronRepulsionIntegralLibInt
(long max_shell_type)¶ Initialize a GB4ElectronRepulsionIntegralLibInt object.
- Parameters
max_shell_type
: Highest angular momentum index to be expected in the reset method.
-
void
laplace_of_potential
(double prefac, double rho, double t, long mmax, double *output)¶ Evaluate the Laplace transform of the ordinary Coulomb potential.
See Eq. (39) in Ahlrichs’ paper. This is basically a rescaled Boys function.
See base class for more details.
-
-
class
GB4ErfIntegralLibInt
¶ - #include <ints.h>
Short-range electron repulsion four-center integrals.
The potential is erf(mu*r)/r.
Inherits from GB4IntegralLibInt
Public Functions
-
GB4ErfIntegralLibInt
(long max_shell_type, double mu)¶ Initialize a GB4ErfIntegralLibInt object.
- Parameters
max_shell_type
: Highest angular momentum index to be expected in the reset method.mu
: The range-separation parameter
-
void
laplace_of_potential
(double prefac, double rho, double t, long mmax, double *output)¶ Evaluate the Laplace transform of the long-range Coulomb potential. (The short-range part is damped away using an error function.) See (52) in Ahlrichs’ paper.
See base class for more details.
-
const double
get_mu
() const¶ The range-separation parameter.
Private Members
-
double
mu
¶ The range-separation parameter.
-
-
class
GB4GaussIntegralLibInt
¶ - #include <ints.h>
Gaussian electron repulsion four-center integrals.
The potential is c exp(-alpha r^2).
Inherits from GB4IntegralLibInt
Public Functions
-
GB4GaussIntegralLibInt
(long max_shell_type, double c, double alpha)¶ Initialize a GB4GaussIntegralLibInt object.
- Parameters
max_shell_type
: Highest angular momentum index to be expected in the reset method.c
: Coefficient of the gaussian.alpha
: Exponential parameter of the gaussian.
-
void
laplace_of_potential
(double prefac, double rho, double t, long mmax, double *output)¶ Evaluate the Laplace transform of the Gaussian potential.
See Ahlrichs’ paper for details. This type of potential is used in the papers of P.M.W Gill et al. and J. Toulouse et al.:
Gill, P. M. W., & Adamson, R. D. (1996). A family of attenuated Coulomb operators. Chem. Phys. Lett., 261(1-2), 105–110. http://doi.org/10.1016/0009-2614(96)00931-1
Toulouse, J., Colonna, F., & Savin, A. (2004). Long-range-short-range separation of the electron-electron interaction in density-functional theory. Phys. Rev. A, 70, 62505. http://doi.org/10.1103/PhysRevA.70.062505
See base class for more details.
-
const double
get_c
() const¶ Coefficient of the gaussian.
-
const double
get_alpha
() const¶ Exponential parameter of the gaussian.
-
-
class
GB4RAlphaIntegralLibInt
¶ - #include <ints.h>
Gaussian electron repulsion four-center integrals.
The potential is r^alpha.
Inherits from GB4IntegralLibInt
Public Functions
-
GB4RAlphaIntegralLibInt
(long max_shell_type, double alpha)¶ Initialize a GB4RAlphaIntegralLibInt object.
- Parameters
max_shell_type
: Highest angular momentum index to be expected in the reset method.alpha
: The power of r in the potential.
-
void
laplace_of_potential
(double prefac, double rho, double t, long mmax, double *output)¶ Evaluate the Laplace transform of the r^alpha potential. See Eq. (49) in Ahlrichs’ paper.
See base class for more details.
-
const double
get_alpha
() const¶ The power of r.
Private Members
-
double
alpha
¶ The power of r.
-