|
ABACUS develop
Atomic-orbital Based Ab-initio Computation at UStc
|
A class to compute two-center integrals. More...
#include <two_center_integrator.h>
Public Member Functions | |
| TwoCenterIntegrator () | |
| TwoCenterIntegrator (const TwoCenterIntegrator &)=delete | |
| TwoCenterIntegrator & | operator= (const TwoCenterIntegrator &)=delete |
| ~TwoCenterIntegrator () | |
| void | tabulate (const RadialCollection &bra, const RadialCollection &ket, const char op, const int nr, const double cutoff) |
| Tabulates the radial part of a two-center integral. | |
| void | calculate (const int itype1, const int l1, const int izeta1, const int m1, const int itype2, const int l2, const int izeta2, const int m2, const ModuleBase::Vector3< double > &vR, double *out=nullptr, double *grad_out=nullptr) const |
| Compute the two-center integrals. | |
| void | snap (const int itype1, const int l1, const int izeta1, const int m1, const int itype2, const ModuleBase::Vector3< double > &vR, const bool deriv, std::vector< std::vector< double > > &out) const |
| Compute a batch of two-center integrals. | |
| size_t | table_memory () const |
| Returns the amount of heap memory used by table_ (in bytes). | |
Private Member Functions | |
| int | ylm_index (const int l, const int m) const |
| Returns the index of (l,m) in the array of spherical harmonics. | |
Private Attributes | |
| bool | is_tabulated_ |
| char | op_ |
| TwoCenterTable | table_ |
A class to compute two-center integrals.
This class computes two-center integrals
/
I(R) = | dr phi1(r) (op) phi2(r - R)
/
as well as their gradients, where op is 1 (overlap) or minus Laplacian (kinetic), and phi1, phi2 are "atomic-orbital-like" functions of the form
phi(r) = chi(|r|) * Ylm(r/|r|)
where chi is some numerical radial function and Ylm is some real spherical harmonics.
This class is designed to efficiently compute the two-center integrals between two "collections" of the above functions with various R, e.g., the overlap integrals between all numerical atomic orbitals and all Kleinman-Bylander nonlocal projectors, the overlap & kinetic integrals between all numerical atomic orbitals, etc. This is done by tabulating the radial part of the integrals on an r-space grid and the real Gaunt coefficients in advance.
See the developer's document for more details.
| TwoCenterIntegrator::TwoCenterIntegrator | ( | ) |
|
delete |
|
inline |
| void TwoCenterIntegrator::calculate | ( | const int | itype1, |
| const int | l1, | ||
| const int | izeta1, | ||
| const int | m1, | ||
| const int | itype2, | ||
| const int | l2, | ||
| const int | izeta2, | ||
| const int | m2, | ||
| const ModuleBase::Vector3< double > & | vR, | ||
| double * | out = nullptr, |
||
| double * | grad_out = nullptr |
||
| ) | const |
Compute the two-center integrals.
This function calculates the two-center integral
/
I(R) = | dr phi1(r) (op_) phi2(r - R)
/
or its gradient by using the tabulated radial part and real Gaunt coefficients.
| [in] | itype1 | Element index of orbital 1. |
| [in] | l1 | Angular momentum of orbital 1. |
| [in] | izeta1 | Zeta number of orbital 1. |
| [in] | m1 | Magnetic quantum number of orbital 1. |
| [in] | itype2 | Element index of orbital 2. |
| [in] | l2 | Angular momentum of orbital 2. |
| [in] | izeta2 | Zeta number of orbital 2. |
| [in] | m2 | Magnetic quantum number of orbital 2. |
| [in] | vR | R2 - R1. |
| [out] | out | Two-center integral. The integral will not be computed if out is nullptr. |
| [out] | grad_out | Gradient of the integral. grad_out[0], grad_out[1] and grad_out[2] are the x, y, z components of the gradient. The gradient will not be computed if grad_out is nullptr. |
|
delete |
| void TwoCenterIntegrator::snap | ( | const int | itype1, |
| const int | l1, | ||
| const int | izeta1, | ||
| const int | m1, | ||
| const int | itype2, | ||
| const ModuleBase::Vector3< double > & | vR, | ||
| const bool | deriv, | ||
| std::vector< std::vector< double > > & | out | ||
| ) | const |
Compute a batch of two-center integrals.
This function calculates the two-center integrals (and optionally their gradients) between one orbital and all orbitals of a certain type from the other collection.
|
inline |
Returns the amount of heap memory used by table_ (in bytes).
| void TwoCenterIntegrator::tabulate | ( | const RadialCollection & | bra, |
| const RadialCollection & | ket, | ||
| const char | op, | ||
| const int | nr, | ||
| const double | cutoff | ||
| ) |
Tabulates the radial part of a two-center integral.
| [in] | bra | The radial functions of the first collection. |
| [in] | ket | The radial functions of the second collection. |
| [in] | op | Operator, could be 'S' or 'T'. |
| [in] | nr | Number of r-space grid points. |
| [in] | cutoff | r-space cutoff radius. |
|
private |
Returns the index of (l,m) in the array of spherical harmonics.
Spherical harmonics in ABACUS are stored in the following order:
index 0 1 2 3 4 5 6 7 8 9 10 11 12 ... l 0 1 1 1 2 2 2 2 2 3 3 3 3 ... m 0 0 1 -1 0 1 -1 2 -2 0 1 -1 2 ...
|
private |
|
private |
|
private |