3.3.9. horton/gbasis/fns.h
– Evaluate functions on grids and derive Fock matrices from potentials on grids.¶
Evaluate functions on grids and derive Fock matrices from potentials on grids.
Functions on grid points are calculated in three steps, which are controlled by methods in the GBasis class in gbasis.h. Also code is shared with the Fock build from potentials on grids.
GB1GridFn functions require only a single loop to compute properties of basis functions. See following methods in gbasis.h:
GOBasis::compute_grid1_exp
GOBasis::compute_grid1_dm
GOBasis::compute_grid1_fock
GOBasis::compute_grid_point1 (called by the previous three)
GB2GridFn functions requires double loop to compute properties of basis pairs of basis functions. See following methods in gbasis.h:
GOBasis::compute_grid2_dm
GOBasis::compute_grid_point2 (called by compute_grid2_dm)
The methods reset
, add
, cart_to_pure
and compute_point_from_*
or compute_fock_from_pot
, are called from functions in gbasis.cpp (GOBasis). For a given grid point…
1) The method reset
is called with basic information of the contraction and the position of the grid point.
2) The method add
is called to evaluate properties (function value and/or derivatives) of primitives on the grid point. Code in gbasis.cpp takes care of contractions and calls add
with the proper prefactor coeff
. The add
method adds results for one primitive in work_cart
.
3) cart_to_pure
to transform the results for a contraction to pure functions when needed.
4) After looping over all the (pairs of) contractions (steps 1 to 3), one of the following two happens, still just for one grid point
4a) The compute_point_from_* method is called to convert the properties of the basis functions into the (density) function(s) of interest, making use of either expansion coefficients of the orbitals or first-order density matrix coefficients.
4b) The compute_fock_from_pot method is called to convert a potential on a grid into a Fock operator.
The above work flow is repeated for every grid point.
The Cartesian polynomials in the Gaussian primitives (and/or their derivatives) are computed only once for a given contraction when calling the reset
method. This is done by calling fill_cartesian_polynomials
, which computes all mononomials in alphabetical order up to a given order. It returns an offset, which is the position of the first mononomial of the highest order. This offset is used to quickly find the desired mononomial for evaluating each primitive. The offset is stored as a class attribute and used on the add
method:
- There is only only offset attribute in many cases, i.e. when only mononomials of one order are needed.
- When mononomials of different orders are needed, their offsets are computed in the
reset
method: offset_l1, offset_h1, offset_l2 and offset_h2. The suffixes h and l refer to higher and lower. - Keep in mind that the derivative (toward x) of a primitive (with x^k) includes two terms (one with a mononomial x^{k-1} and one with x^(k+1}). So in the case of GGA and even more so for MGGA, many mononomials are needed and several relevant offsets must be stored.
Given the position of a mononomial, related mononomials (increasing or decreasing one or two powers) are “easily” found. Given a mononomial:
- mono = x**nx * y**ny * z**nz: offset + i
- mono * x: offset_h1 + i
- mono / x: offset_l1 + i
- mono * y: offset_h1 + i + 1 + ny + nz
- mono / y: offset_l1 + i - ny - nz
- mono * z: offset_h1 + i + 2 + ny + nz
- mono / z: offset_l1 + i - 1 - ny - nz
- mono * (x*x): offset_h2 + i
- mono / (x*x): offset_l2 + i
- mono * (x*y): offset_h2 + i + 1 + ny + nz
- mono / (x*y): offset_l2 + i - ny - nz
- mono * (x*z): offset_h2 + i + 2 + ny + nz
- mono / (x*z): offset_l2 + i - 1 - ny - nz
- mono * (y*y): offset_h2 + i + 3 + 2*ny + 2*nz
- mono / (y*y): offset_l2 + i + 1 - 2*ny - 2*nz
- mono * (y*z): offset_h2 + i + 4 + 2*ny + 2*nz
- mono / (y*z): offset_l2 + i - 2*ny - 2*nz
- mono * (z*z): offset_h2 + i + 5 + 2*ny + 2*nz
- mono / (z*z): offset_l2 + i - 1 - 2*ny - 2*nz
These rules are used on the add
methods.
- class
- #include <fns.h>
Base class for grid calculators that require only a single loop over all basis functions.
Inherits from GBCalculator
Subclassed by GB1DMGridFn, GB1ExpGridFn
Public Functions
-
GB1GridFn::
GB1GridFn
(long max_shell_type, long dim_work, long dim_output)¶ Construct a GB1GridFn object.
- Parameters
max_shell_type
-The maximum shell type in the basis set.
dim_work
-A multiplier for the size of the work array, e.g. when multiple results need to be stored, such as the orbitals and its gradient.
dim_output
-The number of results for each grid point, e.g. 3 for a density gradient.
-
void
GB1GridFn::
reset
(long shell_type0, const double *r0, const double *point)¶ Reset calculator for a new contraction.
- Parameters
shell_type0
-Shell type of contraction.
r0
-Center of the contraction. (size=3)
point
-Cartesian coordinates of the grid point. (size=3)
-
void
GB1GridFn::
cart_to_pure
()¶ Convert results in work array from Cartesian to pure functions where needed.
-
const long
GB1GridFn::
get_shell_type0
() const¶ Shell type of the contraction.
-
long
GB1GridFn::
get_dim_work
()¶ Multiplier for the size of the work array.
-
long
GB1GridFn::
get_dim_output
()¶ Number of results per grid point.
-
virtual void
GB1GridFn::
add
(double coeff, double alpha0, const double *scales0)¶
= 0 Add contributions to work array for current grid point and given primitive shell.
- Parameters
coeff
-Contraction coefficient of current primitive shells.
alpha0
-Exponent of the primitive shell.
scales0
-Normalization prefactors for primitives in the shell.
Protected Attributes
-
const long
GB1GridFn::
dim_work
¶ Multiplier for the size of the work array.
-
const long
GB1GridFn::
dim_output
¶ Number of results per grid point.
-
long
GB1GridFn::
shell_type0
¶ Shell type of the current contraction.
-
const double *
GB1GridFn::
r0
¶ Center of the current contraction.
-
const double *
GB1GridFn::
point
¶ Grid point at which the fn is evaluated.
-
- class
- #include <fns.h>
Base class for GB1 grid calculators that use the expansion coefficients of the orbitals.
Inherits from GB1GridFn
Subclassed by GB1ExpGridOrbGradientFn, GB1ExpGridOrbitalFn
Public Functions
-
GB1ExpGridFn::
GB1ExpGridFn
(long max_shell_type, long nfn, long dim_work, long dim_output)¶ Construct a GB1ExpGridFn object.
- Parameters
max_shell_type
-The maximum shell type in the basis set.
nfn
-The number of orbitals (occupied and virtual).
dim_work
-A multiplier for the size of the work array, e.g. when multiple results need to be stored, such as an orbital and its gradient.
dim_output
-The number of results, i.e. elements in output and pot arguments, for each grid point, e.g. 3 for a density gradient.
-
virtual void
GB1ExpGridFn::
compute_point_from_exp
(double *work_basis, double *coeffs, long nbasis, double *output)¶
= 0 Compute (final) results for a given grid point.
- Parameters
work_basis
-Properties of basis functions computed for the current grid point. (Work done by add method.) (size=nbasis*dim_work)
coeffs
-The orbital expansion coefficients. (size=nbasis*nfn)
nbasis
-The number of basis functions.
output
-The output array for the current grid point. (size=dim_output)
Protected Attributes
-
long
GB1ExpGridFn::
nfn
¶ Number of orbitals (occupied and virtual).
-
- class
- #include <fns.h>
Evaluates a selection of orbitals on a grid.
Content of work_basis (at one grid point): [0] Basis function value. Content of the argument ‘output’ (at one grid point): [0-norb] Values of the orbitals.
Inherits from GB1ExpGridFn
Public Functions
-
GB1ExpGridOrbitalFn::
GB1ExpGridOrbitalFn
(long max_shell_type, long nfn, long *iorbs, long norb)¶ Construct a GB1ExpGridOrbitalFn object.
- Parameters
max_shell_type
-The maximum shell type in the basis set.
nfn
-The number of orbitals (occupied and virtual).
iorbs
-An array with orbitals to be computed.
norb
-The number of elements in iorbs.
-
void
GB1ExpGridOrbitalFn::
reset
(long _shell_type0, const double *_r0, const double *_point)¶ Reset calculator for a new contraction. (See base class for details.)
-
void
GB1ExpGridOrbitalFn::
add
(double coeff, double alpha0, const double *scales0)¶ Add contributions to work array for current grid point and given primitive shell. (See base class for more details.)
-
void
GB1ExpGridOrbitalFn::
compute_point_from_exp
(double *work_basis, double *coeffs, long nbasis, double *output)¶ Compute (final) results for a given grid point. (See base class for details.)
Protected Attributes
-
double
GB1ExpGridOrbitalFn::
poly_work
[MAX_NCART_CUMUL]¶ Work array with Cartesian polynomials.
-
long
GB1ExpGridOrbitalFn::
offset
¶ Offset for the polynomials for the density.
-
long *
GB1ExpGridOrbitalFn::
iorbs
¶ Array of indices of orbitals to be evaluated on grid.
-
long
GB1ExpGridOrbitalFn::
norb
¶ The number of elements in iorbs.
-
- class
- #include <fns.h>
Evaluates the gradient of a selection of orbitals on a grid.
Content of work_basis (at one grid point): [0] Basis function derivative toward x. [1] Basis function derivative toward y. [2] Basis function derivative toward z. Content of the argument ‘output’ (at one grid point): [i*norb + 0] derivative of orbital i toward x. [i*norb + 1] derivative of orbital i toward y. [i*norb + 2] derivative of orbital i toward z.
Inherits from GB1ExpGridFn
Public Functions
-
GB1ExpGridOrbGradientFn::
GB1ExpGridOrbGradientFn
(long max_shell_type, long nfn, long *iorbs, long norb)¶ Construct a GB1ExpGridOrbGradientFn object.
- Parameters
max_shell_type
-The maximum shell type in the basis set.
nfn
-The number of orbitals (occupied and virtual).
iorbs
-An array with orbitals to be computed.
norb
-The number of elements in iorbs.
-
void
GB1ExpGridOrbGradientFn::
reset
(long _shell_type0, const double *_r0, const double *_point)¶ Reset calculator for a new contraction. (See base class for details.)
-
void
GB1ExpGridOrbGradientFn::
add
(double coeff, double alpha0, const double *scales0)¶ Add contributions to work array for current grid point and given primitive shell. (See base class for more details.)
-
void
GB1ExpGridOrbGradientFn::
compute_point_from_exp
(double *work_basis, double *coeffs, long nbasis, double *output)¶ Compute the final result on one grid point. (See base class for details.)
Protected Attributes
-
long *
GB1ExpGridOrbGradientFn::
iorbs
¶ Array of indices of orbitals to be evaluated on grid.
-
long
GB1ExpGridOrbGradientFn::
norb
¶ The number of elements in iorbs.
-
double
GB1ExpGridOrbGradientFn::
poly_work
[MAX_NCART_CUMUL_D]¶ Work array with Cartesian polynomials.
-
long
GB1ExpGridOrbGradientFn::
offset
¶ Offset for the polynomials for the density.
-
long
GB1ExpGridOrbGradientFn::
offset_l1
¶ Lower offset for the polynomials for the gradient.
-
long
GB1ExpGridOrbGradientFn::
offset_h1
¶ Higher offset for the polynomials for the gradient.
-
- class
- #include <fns.h>
Base class for GB1 grid calculators that use the first-order density matrix coefficients.
Inherits from GB1GridFn
Subclassed by GB1DMGridDensityFn, GB1DMGridGradientFn, GB1DMGridHessianFn, GB1DMGridKineticFn, GB1DMGridMGGAFn
Public Functions
-
GB1DMGridFn::
GB1DMGridFn
(long max_shell_type, long dim_work, long dim_output)¶ Construct a GB1GridFn object. (See base class for details.)
-
virtual void
GB1DMGridFn::
compute_point_from_dm
(double *work_basis, double *dm, long nbasis, double *output, double epsilon, double *dmmaxrow)¶
= 0 Compute the final result on one grid point.
- Parameters
work_basis
-Properties of basis functions computed for the current grid point. (Work done by add method.) (size=nbasis*dim_work)
dm
-The coefficients of the first-order density matrix. (size=nbasis*nbasis)
nbasis
-The number of basis functions.
output
-The output array for the current grid point. (size=dim_output)
epsilon
-A cutoff value used to discard small contributions.
dmmaxrow
-The maximum value of the density matrix on each row. (size=nbasis)
-
virtual void
GB1DMGridFn::
compute_fock_from_pot
(double *pot, double *work_basis, long nbasis, double *fock)¶
= 0 Add contribution to Fock matrix from one grid point.
The chain rule is used to transform grid potential data (in one point, see reset method) into a Fock matrix contribution.
- Parameters
pot
-The value of the potential at the grid point. This may be multiple values, e.g in case of GGA, this contains four elements: the functional derivative of the energy w.r.t. the density and the components of the gradient. (size=dim_output)
work_basis
-Properties of the orbital basis in the current grid point, typically the value of the basis function and optionally first or second derivatives toward x, y and z, all evaluated in
point
(see reset method). (size=nbasis*dim_work)nbasis
-The number of basis functions.
fock
-The Fock matrix to which the result will be added. (size=nbasis*nbasis)
-
- class
- #include <fns.h>
Compute just the electron density on a grid.
Content of work_basis (at one grid point): [0] Basis function value. Content of the argument ‘output’ and the energy derivative in ‘pot’ (at one grid point): [0] The electron density
Inherits from GB1DMGridFn
Public Functions
-
GB1DMGridDensityFn::
GB1DMGridDensityFn
(long max_shell_type)¶ Construct a GB1DMGridDensityFn object.
- Parameters
max_shell_type
-The maximum shell type in the basis set.
-
void
GB1DMGridDensityFn::
reset
(long _shell_type0, const double *_r0, const double *_point)¶ Reset calculator for a new contraction. (See base class for details.)
-
void
GB1DMGridDensityFn::
add
(double coeff, double alpha0, const double *scales0)¶ Add contributions to work array for current grid point and given primitive shell. (See base class for more details.)
-
void
GB1DMGridDensityFn::
compute_point_from_dm
(double *work_basis, double *dm, long nbasis, double *output, double epsilon, double *dmmaxrow)¶ Compute the final result on one grid point. (See base class for details.)
-
void
GB1DMGridDensityFn::
compute_fock_from_pot
(double *pot, double *work_basis, long nbasis, double *fock)¶ Add contribution to Fock matrix for one grid point. (See base class for details.)
-
- class
- #include <fns.h>
Compute gradient of the electron density on a grid.
Content of work_basis (at one grid point): [0] Basis function value. [1] Basis function derivative toward x. [2] Basis function derivative toward y. [3] Basis function derivative toward z. Content of the argument ‘output’ and the energy derivative in ‘pot’ (at one grid point): [0] Density derivative toward x. [1] Density derivative toward y. [2] Density derivative toward z.
Inherits from GB1DMGridFn
Subclassed by GB1DMGridGGAFn
Public Functions
-
GB1DMGridGradientFn::
GB1DMGridGradientFn
(long max_shell_type)¶ Construct a GB1DMGridGradientFn object.
- Parameters
max_shell_type
-The maximum shell type in the basis set.
-
GB1DMGridGradientFn::
GB1DMGridGradientFn
(long max_shell_type, long dim_output)¶ Construct a GB1DMGridGradientFn object.
- Parameters
max_shell_type
-The maximum shell type in the basis set.
dim_output
-The number of results for each grid point. This may be different from 3 for derived classes.
-
void
GB1DMGridGradientFn::
reset
(long _shell_type0, const double *_r0, const double *_point)¶ Reset calculator for a new contraction. (See base class for details.)
-
void
GB1DMGridGradientFn::
add
(double coeff, double alpha0, const double *scales0)¶ Add contributions to work array for current grid point and given primitive shell. (See base class for more details.)
-
void
GB1DMGridGradientFn::
compute_point_from_dm
(double *work_basis, double *dm, long nbasis, double *output, double epsilon, double *dmmaxrow)¶ Compute the final result on one grid point. (See base class for details.)
-
void
GB1DMGridGradientFn::
compute_fock_from_pot
(double *pot, double *work_basis, long nbasis, double *fock)¶ Add contribution to Fock matrix for one grid point. (See base class for details.)
Protected Attributes
-
double
GB1DMGridGradientFn::
poly_work
[MAX_NCART_CUMUL_D]¶ Work array with Cartesian polynomials.
-
long
GB1DMGridGradientFn::
offset
¶ Offset for the polynomials for the density.
-
long
GB1DMGridGradientFn::
offset_l1
¶ Lower offset for the polynomials for the gradient.
-
long
GB1DMGridGradientFn::
offset_h1
¶ Higher offset for the polynomials for the gradient.
-
- class
- #include <fns.h>
Compute density and gradient on a grid.
Content of work_basis (at one grid point): [0] Basis function value. [1] Basis function derivative toward x. [2] Basis function derivative toward y. [3] Basis function derivative toward z. Content of the argument ‘output’ and the energy derivative in ‘pot’ (at one grid point): [0] Density. [1] Density derivative toward x. [2] Density derivative toward y. [3] Density derivative toward z.
Inherits from GB1DMGridGradientFn
Public Functions
-
GB1DMGridGGAFn::
GB1DMGridGGAFn
(long max_shell_type)¶ Construct a GB1DMGridGGAFn object.
- Parameters
max_shell_type
-The maximum shell type in the basis set.
-
void
GB1DMGridGGAFn::
compute_point_from_dm
(double *work_basis, double *dm, long nbasis, double *output, double epsilon, double *dmmaxrow)¶ Compute the final result on one grid point. (See base class for details.)
-
void
GB1DMGridGGAFn::
compute_fock_from_pot
(double *pot, double *work_basis, long nbasis, double *fock)¶ Add contribution to Fock matrix for one grid point. (See base class for details.)
-
- class
- #include <fns.h>
Compute kinetic energy density on a grid.
Content of work_basis (at one grid point): [0] Basis function derivative toward x. [1] Basis function derivative toward y. [2] Basis function derivative toward z. Content of the argument ‘output’ and the energy derivative in ‘pot’ (at one grid point): [0] Kinetic energy density.
Inherits from GB1DMGridFn
Public Functions
-
GB1DMGridKineticFn::
GB1DMGridKineticFn
(long max_shell_type)¶ Construct a GB1DMGridKineticFn object.
- Parameters
max_shell_type
-The maximum shell type in the basis set.
-
void
GB1DMGridKineticFn::
reset
(long _shell_type0, const double *_r0, const double *_point)¶ Reset calculator for a new contraction. (See base class for details.)
-
void
GB1DMGridKineticFn::
add
(double coeff, double alpha0, const double *scales0)¶ Add contributions to work array for current grid point and given primitive shell. (See base class for more details.)
-
void
GB1DMGridKineticFn::
compute_point_from_dm
(double *work_basis, double *dm, long nbasis, double *output, double epsilon, double *dmmaxrow)¶ Compute the final result on one grid point. (See base class for details.)
-
void
GB1DMGridKineticFn::
compute_fock_from_pot
(double *pot, double *work_basis, long nbasis, double *fock)¶ Add contribution to Fock matrix for one grid point. (See base class for details.)
Private Members
-
double
GB1DMGridKineticFn::
poly_work
[MAX_NCART_CUMUL_D]¶ Work array with Cartesian polynomials.
-
long
GB1DMGridKineticFn::
offset
¶ Offset for the polynomials for the density.
-
long
GB1DMGridKineticFn::
offset_l1
¶ Lower offset for the polynomials for the gradient.
-
long
GB1DMGridKineticFn::
offset_h1
¶ Higher offset for the polynomials for the gradient.
-
- class
- #include <fns.h>
Compute density Hessian on a grid: xx, xy, xz, yy, yz, zz.
Content of work_basis (at one grid point): [0] Basis function value. [1] Basis function derivative toward x. [2] Basis function derivative toward y. [3] Basis function derivative toward z. [4] Basis function derivative toward xx. [5] Basis function derivative toward xy. [6] Basis function derivative toward xz. [7] Basis function derivative toward yy. [8] Basis function derivative toward yz. [9] Basis function derivative toward zz. Content of the argument ‘output’ and the energy derivative in ‘pot’ (at one grid point): [0] Density derivative toward xx. [1] Density derivative toward xy. [2] Density derivative toward xz. [3] Density derivative toward yy. [4] Density derivative toward yz. [5] Density derivative toward zz.
Inherits from GB1DMGridFn
Public Functions
-
GB1DMGridHessianFn::
GB1DMGridHessianFn
(long max_shell_type)¶ Construct a GB1DMGridHessianFn object.
- Parameters
max_shell_type
-The maximum shell type in the basis set.
-
void
GB1DMGridHessianFn::
reset
(long _shell_type0, const double *_r0, const double *_point)¶ Reset calculator for a new contraction. (See base class for details.)
-
void
GB1DMGridHessianFn::
add
(double coeff, double alpha0, const double *scales0)¶ Add contributions to work array for current grid point and given primitive shell. (See base class for more details.)
-
void
GB1DMGridHessianFn::
compute_point_from_dm
(double *work_basis, double *dm, long nbasis, double *output, double epsilon, double *dmmaxrow)¶ Compute the final result on one grid point. (See base class for details.)
-
void
GB1DMGridHessianFn::
compute_fock_from_pot
(double *pot, double *work_basis, long nbasis, double *fock)¶ Add contribution to Fock matrix for one grid point. (See base class for details.)
Private Members
-
double
GB1DMGridHessianFn::
poly_work
[MAX_NCART_CUMUL_DD]¶ Work array with Cartesian polynomials.
-
long
GB1DMGridHessianFn::
offset
¶ Offset for the polynomials for the density.
-
long
GB1DMGridHessianFn::
offset_l1
¶ Lower offset for the polynomials for the gradient.
-
long
GB1DMGridHessianFn::
offset_h1
¶ Higher offset for the polynomials for the gradient.
-
long
GB1DMGridHessianFn::
offset_l2
¶ Lower offset for the polynomials for the hessian.
-
long
GB1DMGridHessianFn::
offset_h2
¶ Higher offset for the polynomials for the hessian.
-
- class
- #include <fns.h>
Compute MGGA properties on a grid: density, gradient, laplacian and kinetic energy density.
Content of work_basis (at one grid point): [0] Basis function value. [1] Basis function derivative toward x. [2] Basis function derivative toward y. [3] Basis function derivative toward z. [4] Basis function Laplacian. Content of the argument ‘output’ and the energy derivative in ‘pot’ (at one grid point): [0] Density. [1] Density derivative toward x. [2] Density derivative toward y. [3] Density derivative toward z. [4] Laplacian of the density. [5] Kinetic energy density.
Inherits from GB1DMGridFn
Public Functions
-
GB1DMGridMGGAFn::
GB1DMGridMGGAFn
(long max_shell_type)¶ Construct a GB1DMGridMGGAFn object.
- Parameters
max_shell_type
-The maximum shell type in the basis set.
-
void
GB1DMGridMGGAFn::
reset
(long _shell_type0, const double *_r0, const double *_point)¶ Reset calculator for a new contraction. (See base class for details.)
-
void
GB1DMGridMGGAFn::
add
(double coeff, double alpha0, const double *scales0)¶ Add contributions to work array for current grid point and given primitive shell. (See base class for more details.)
-
void
GB1DMGridMGGAFn::
compute_point_from_dm
(double *work_basis, double *dm, long nbasis, double *output, double epsilon, double *dmmaxrow)¶ Compute the final result on one grid point. (See base class for details.)
-
void
GB1DMGridMGGAFn::
compute_fock_from_pot
(double *pot, double *work_basis, long nbasis, double *fock)¶ Add contribution to Fock matrix for one grid point. (See base class for details.)
Private Members
-
double
GB1DMGridMGGAFn::
poly_work
[MAX_NCART_CUMUL_DD]¶ Work array with Cartesian polynomials.
-
long
GB1DMGridMGGAFn::
offset
¶ Offset for the polynomials for the density.
-
long
GB1DMGridMGGAFn::
offset_l1
¶ Lower offset for the polynomials for the gradient.
-
long
GB1DMGridMGGAFn::
offset_h1
¶ Higher offset for the polynomials for the gradient.
-
long
GB1DMGridMGGAFn::
offset_l2
¶ Lower offset for the polynomials for the hessian.
-
long
GB1DMGridMGGAFn::
offset_h2
¶ Higher offset for the polynomials for the hessian.
-
- class
- #include <fns.h>
Base class for grid calculators that require a double loop over all basis functions.
Inherits from GBCalculator
Subclassed by GB2DMGridHartreeFn
Public Functions
-
GB2DMGridFn::
GB2DMGridFn
(long max_shell_type)¶ Construct a GB2DMGridFn object. (See base class for details.)
-
void
GB2DMGridFn::
reset
(long shell_type0, long shell_type1, const double *r0, const double *r1, const double *point)¶ Reset calculator for a new pair of contractions.
- Parameters
shell_type0
-Shell type of first contraction.
shell_type1
-Shell type of second contraction.
r0
-Center of the first contraction. (size=3)
r1
-Center of the second contraction. (size=3)
point
-Cartesian coordinates of the grid point. (size=3)
-
void
GB2DMGridFn::
cart_to_pure
()¶ Convert results in work array from Cartesian to pure functions where needed.
-
const long
GB2DMGridFn::
get_shell_type0
() const¶ Shell type of the first contraction.
-
const long
GB2DMGridFn::
get_shell_type1
() const¶ Shell type of the second contraction.
-
virtual void
GB2DMGridFn::
add
(double coeff, double alpha0, double alpha1, const double *scales0, const double *scales1)¶
= 0 Add contributions to work array for current grid point and given pair of primitive shells.
- Parameters
coeff
-Product of contraction coefficients of current primitive shells.
alpha0
-Exponent of the first shell.
alpha1
-Exponent of the second shell.
scales0
-Normalization prefactors for primitives in first shell.
scales1
-Normalization prefactors for primitives in second shell.
Protected Attributes
-
long
GB2DMGridFn::
shell_type0
¶ Shell type of the first contraction.
-
long
GB2DMGridFn::
shell_type1
¶ Shell type of the second contraction.
-
const double *
GB2DMGridFn::
r0
¶ Center of first basis contraction.
-
const double *
GB2DMGridFn::
r1
¶ Center of second basis contraction.
-
const double *
GB2DMGridFn::
point
¶ Grid point at which the fn is evaluated.
-
- class
- #include <fns.h>
Calculator for Hartree potential on a set of grid points.
This requires a double loop over all basis functions because the contributions from the products of basis functions cannot completely be factorized.
Content of work_basis (at one grid point): [0] The Hartree potential due to the product of the two basis functions.
Inherits from GB2DMGridFn
Public Functions
-
GB2DMGridHartreeFn::
GB2DMGridHartreeFn
(long max_shell_type)¶ Construct a GB2DMGridHartreeFn object. (See base class for details.)
-
GB2DMGridHartreeFn::
~GB2DMGridHartreeFn
()¶
-
void
GB2DMGridHartreeFn::
add
(double coeff, double alpha0, double alpha1, const double *scales0, const double *scales1)¶ Add contributions to work array for current grid point and given pair of primitive shells. (See base class for details.)
Private Members
-
double *
GB2DMGridHartreeFn::
work_g0
¶ Work array for results from NAI helper code (x).
-
double *
GB2DMGridHartreeFn::
work_g1
¶ Work array for results from NAI helper code (y).
-
double *
GB2DMGridHartreeFn::
work_g2
¶ Work array for results from NAI helper code (z).
-
double *
GB2DMGridHartreeFn::
work_boys
¶ Work array for Boys function values for different m values.
-