becke_numerical_grid

This module contains all quantities needed to build Becke’s grid used in general for DFT integration. Note that it can be used for whatever integration in R^3 as long as the functions to be integrated are mostly concentrated near the atomic regions.

This grid is built as the reunion of a spherical grid around each atom. Each spherical grid contains a certain number of radial and angular points. No pruning is done on the angular part of the grid.

The main keyword for that module is:

  • becke_numerical_grid grid_type_sgn which controls the precision of the grid according the standard SG-n grids. This keyword controls the two providers n_points_integration_angular n_points_radial_grid.

The main providers of that module are:

  • n_points_integration_angular which is the number of angular integration points. WARNING: it obeys to specific rules so it cannot be any integer number. Some of the possible values are [ 50 | 74 | 170 | 194 | 266 | 302 | 590 | 1202 | 2030 | 5810 ] for instance. See angular.f for more details.

  • n_points_radial_grid which is the number of radial angular points. This can be any strictly positive integer. Nevertheless, a minimum of 50 is in general necessary.

  • final_grid_points which are the (x,y,z) coordinates of the grid points.

  • final_weight_at_r_vector which are the weights at each grid point

For a simple example of how to use the grid, see example.irp.f.

The spherical integration uses Lebedev-Laikov grids, which was used from the code distributed through CCL (http://www.ccl.net/). See next section for explanations and citation policies.

This subroutine is part of a set of subroutines that generate
Lebedev grids [1-6] for integration on a sphere. The original
C-code [1] was kindly provided by Dr. Dmitri N. Laikov and
translated into fortran by Dr. Christoph van Wuellen.
This subroutine was translated using a C to fortran77 conversion
tool written by Dr. Christoph van Wuellen.

Users of this code are asked to include reference [1] in their
publications, and in the user- and programmers-manuals
describing their codes.

This code was distributed through CCL (http://www.ccl.net/).

[1] V.I. Lebedev, and D.N. Laikov
    "A quadrature formula for the sphere of the 131st
     algebraic order of accuracy"
    Doklady Mathematics, Vol. 59, No. 3, 1999, pp. 477-481.

[2] V.I. Lebedev
    "A quadrature formula for the sphere of 59th algebraic
     order of accuracy"
    Russian Acad. Sci. Dokl. Math., Vol. 50, 1995, pp. 283-286.

[3] V.I. Lebedev, and A.L. Skorokhodov
    "Quadrature formulas of orders 41, 47, and 53 for the sphere"
    Russian Acad. Sci. Dokl. Math., Vol. 45, 1992, pp. 587-592.

[4] V.I. Lebedev
    "Spherical quadrature formulas exact to orders 25-29"
    Siberian Mathematical Journal, Vol. 18, 1977, pp. 99-107.

[5] V.I. Lebedev
    "Quadratures on a sphere"
    Computational Mathematics and Mathematical Physics, Vol. 16,
    1976, pp. 10-24.

[6] V.I. Lebedev
    "Values of the nodes and weights of ninth to seventeenth
     order Gauss-Markov quadrature formulae invariant under the
     octahedron group with inversion"
    Computational Mathematics and Mathematical Physics, Vol. 15,
    1975, pp. 44-51.

EZFIO parameters

grid_type_sgn

Type of grid used for the Becke’s numerical grid. Can be, by increasing accuracy: [ 0 | 1 | 2 | 3 ]

Default: 2

n_points_final_grid

Total number of grid points

thresh_grid

threshold on the weight of a given grid point

Default: 1.e-20

my_grid_becke

if True, the number of angular and radial grid points are read from EZFIO

Default: False

my_n_pt_r_grid

Number of radial grid points given from input

Default: 300

my_n_pt_a_grid

Number of angular grid points given from input. Warning, this number cannot be any integer. See file list_angular_grid

Default: 1202

n_points_extra_final_grid

Total number of extra_grid points

extra_grid_type_sgn

Type of extra_grid used for the Becke’s numerical extra_grid. Can be, by increasing accuracy: [ 0 | 1 | 2 | 3 ]

Default: 0

thresh_extra_grid

threshold on the weight of a given extra_grid point

Default: 1.e-20

my_extra_grid_becke

if True, the number of angular and radial extra_grid points are read from EZFIO

Default: False

my_n_pt_r_extra_grid

Number of radial extra_grid points given from input

Default: 300

my_n_pt_a_extra_grid

Number of angular extra_grid points given from input. Warning, this number cannot be any integer. See file list_angular_extra_grid

Default: 1202

rad_grid_type

method used to sample the radial space. Possible choices are [KNOWLES | GILL]

Default: KNOWLES

extra_rad_grid_type

method used to sample the radial space. Possible choices are [KNOWLES | GILL]

Default: KNOWLES

Providers

alpha_knowles

File : becke_numerical_grid/integration_radial.irp.f

double precision, allocatable   :: alpha_knowles        (100)

Recommended values for the alpha parameters according to the paper of Knowles (JCP, 104, 1996) as a function of the nuclear charge

Needed by:

  • final_weight_at_r

  • final_weight_at_r_extra

  • grid_points_extra_per_atom

  • grid_points_per_atom

angular_quadrature_points

File : becke_numerical_grid/angular_grid_pts.irp.f

double precision, allocatable   :: angular_quadrature_points    (n_points_integration_angular,3)
double precision, allocatable   :: weights_angular_points       (n_points_integration_angular)

weights and grid points for the integration on the angular variables on the unit sphere centered on (0,0,0) According to the LEBEDEV scheme

Needs:

  • n_points_radial_grid

Needed by:

  • final_weight_at_r

  • grid_points_per_atom

angular_quadrature_points_extra

File : becke_numerical_grid/angular_extra_grid.irp.f

double precision, allocatable   :: angular_quadrature_points_extra      (n_points_extra_integration_angular,3)
double precision, allocatable   :: weights_angular_points_extra (n_points_extra_integration_angular)

weights and grid points_extra for the integration on the angular variables on the unit sphere centered on (0,0,0) According to the LEBEDEV scheme

Needs:

  • n_points_extra_radial_grid

Needed by:

  • final_weight_at_r_extra

  • grid_points_extra_per_atom

dr_radial_extra_integral

File : becke_numerical_grid/extra_grid.irp.f

double precision, allocatable   :: grid_points_extra_radial     (n_points_extra_radial_grid)
double precision        :: dr_radial_extra_integral

points_extra in [0,1] to map the radial integral [0,infty]

Needs:

  • n_points_extra_radial_grid

Needed by:

  • final_weight_at_r_extra

  • grid_points_extra_per_atom

dr_radial_integral

File : becke_numerical_grid/grid_becke.irp.f

double precision, allocatable   :: grid_points_radial   (n_points_radial_grid)
double precision        :: dr_radial_integral

points in [0,1] to map the radial integral [0,infty]

Needs:

  • n_points_radial_grid

Needed by:

  • final_weight_at_r

  • grid_points_per_atom

final_grid_points

File : becke_numerical_grid/grid_becke_vector.irp.f

   double precision, allocatable   :: final_grid_points    (3,n_points_final_grid)
   double precision, allocatable   :: final_weight_at_r_vector     (n_points_final_grid)
   integer, allocatable    :: index_final_points   (3,n_points_final_grid)
   integer, allocatable    :: index_final_points_reverse   (n_points_integration_angular,n_points_radial_grid,nucl_num)


final_grid_points(1:3,j) = (/ x, y, z /) of the jth grid point

final_weight_at_r_vector(i) = Total weight function of the ith grid point which contains the Lebedev, Voronoi and radial weights contributions

index_final_points(1:3,i) = gives the angular, radial and atomic indices associated to the ith grid point

index_final_points_reverse(i,j,k) = index of the grid point having i as angular, j as radial and l as atomic indices

Needs:

  • final_weight_at_r

  • grid_points_per_atom

  • n_points_final_grid

  • n_points_radial_grid

  • nucl_num

  • thresh_grid

Needed by:

  • ao_abs_int_grid

  • ao_overlap_abs_grid

  • ao_prod_abs_r

  • ao_prod_center

  • ao_prod_dist_grid

  • aos_grad_in_r_array

  • aos_in_r_array

  • aos_lapl_in_r_array

  • aos_sr_vc_alpha_lda_w

  • aos_sr_vxc_alpha_lda_w

  • aos_vc_alpha_lda_w

  • aos_vc_alpha_pbe_w

  • aos_vc_alpha_sr_pbe_w

  • aos_vxc_alpha_lda_w

  • aos_vxc_alpha_pbe_w

  • aos_vxc_alpha_sr_pbe_w

  • elec_beta_num_grid_becke

  • energy_c_lda

  • energy_c_sr_lda

  • energy_x_lda

  • energy_x_pbe

  • energy_x_sr_lda

  • energy_x_sr_pbe

  • f_psi_cas_ab

  • f_psi_hf_ab

  • final_grid_points_transp

  • mo_grad_ints

  • mos_in_r_array

  • mos_in_r_array_omp

  • mu_average_prov

  • mu_grad_rho

  • mu_of_r_dft_average

  • mu_rsc_of_r

  • one_e_dm_and_grad_alpha_in_r

final_grid_points_extra

File : becke_numerical_grid/extra_grid_vector.irp.f

   double precision, allocatable   :: final_grid_points_extra      (3,n_points_extra_final_grid)
   double precision, allocatable   :: final_weight_at_r_vector_extra       (n_points_extra_final_grid)
   integer, allocatable    :: index_final_points_extra     (3,n_points_extra_final_grid)
   integer, allocatable    :: index_final_points_extra_reverse     (n_points_extra_integration_angular,n_points_extra_radial_grid,nucl_num)


final_grid_points_extra(1:3,j) = (/ x, y, z /) of the jth grid point

final_weight_at_r_vector_extra(i) = Total weight function of the ith grid point which contains the Lebedev, Voronoi and radial weights contributions

index_final_points_extra(1:3,i) = gives the angular, radial and atomic indices associated to the ith grid point

index_final_points_extra_reverse(i,j,k) = index of the grid point having i as angular, j as radial and l as atomic indices

Needs:

  • final_weight_at_r_extra

  • grid_points_extra_per_atom

  • n_points_extra_final_grid

  • n_points_extra_radial_grid

  • nucl_num

  • thresh_extra_grid

Needed by:

  • aos_in_r_array_extra

final_grid_points_per_atom

File : becke_numerical_grid/grid_becke_per_atom.irp.f

double precision, allocatable   :: final_grid_points_per_atom   (3,n_pts_max_per_atom,nucl_num)
double precision, allocatable   :: final_weight_at_r_vector_per_atom    (n_pts_max_per_atom,nucl_num)
integer, allocatable    :: index_final_points_per_atom  (3,n_pts_max_per_atom,nucl_num)
integer, allocatable    :: index_final_points_per_atom_reverse  (n_points_integration_angular,n_points_radial_grid,nucl_num)

Needs:

  • final_weight_at_r

  • grid_points_per_atom

  • n_points_radial_grid

  • n_pts_per_atom

  • nucl_num

  • thresh_grid

final_grid_points_transp

File : becke_numerical_grid/grid_becke_vector.irp.f

double precision, allocatable   :: final_grid_points_transp     (n_points_final_grid,3)

Transposed final_grid_points

Needs:

  • final_grid_points

  • n_points_final_grid

final_weight_at_r

File : becke_numerical_grid/grid_becke.irp.f

double precision, allocatable   :: final_weight_at_r    (n_points_integration_angular,n_points_radial_grid,nucl_num)

Total weight on each grid point which takes into account all Lebedev, Voronoi and radial weights.

Needs:

  • alpha_knowles

  • angular_quadrature_points

  • grid_atomic_number

  • grid_points_radial

  • m_knowles

  • n_points_radial_grid

  • nucl_num

  • r_gill

  • rad_grid_type

  • weight_at_r

Needed by:

  • final_grid_points

  • final_grid_points_per_atom

  • n_points_final_grid

  • n_pts_per_atom

final_weight_at_r_extra

File : becke_numerical_grid/extra_grid.irp.f

double precision, allocatable   :: final_weight_at_r_extra      (n_points_extra_integration_angular,n_points_extra_radial_grid,nucl_num)

Total weight on each grid point which takes into account all Lebedev, Voronoi and radial weights.

Needs:

  • alpha_knowles

  • angular_quadrature_points_extra

  • extra_rad_grid_type

  • grid_atomic_number

  • grid_points_extra_radial

  • m_knowles

  • n_points_extra_radial_grid

  • nucl_num

  • r_gill

  • weight_at_r_extra

Needed by:

  • final_grid_points_extra

  • n_points_extra_final_grid

final_weight_at_r_vector

File : becke_numerical_grid/grid_becke_vector.irp.f

   double precision, allocatable   :: final_grid_points    (3,n_points_final_grid)
   double precision, allocatable   :: final_weight_at_r_vector     (n_points_final_grid)
   integer, allocatable    :: index_final_points   (3,n_points_final_grid)
   integer, allocatable    :: index_final_points_reverse   (n_points_integration_angular,n_points_radial_grid,nucl_num)


final_grid_points(1:3,j) = (/ x, y, z /) of the jth grid point

final_weight_at_r_vector(i) = Total weight function of the ith grid point which contains the Lebedev, Voronoi and radial weights contributions

index_final_points(1:3,i) = gives the angular, radial and atomic indices associated to the ith grid point

index_final_points_reverse(i,j,k) = index of the grid point having i as angular, j as radial and l as atomic indices

Needs:

  • final_weight_at_r

  • grid_points_per_atom

  • n_points_final_grid

  • n_points_radial_grid

  • nucl_num

  • thresh_grid

Needed by:

  • ao_abs_int_grid

  • ao_overlap_abs_grid

  • ao_prod_abs_r

  • ao_prod_center

  • ao_prod_dist_grid

  • aos_grad_in_r_array

  • aos_in_r_array

  • aos_lapl_in_r_array

  • aos_sr_vc_alpha_lda_w

  • aos_sr_vxc_alpha_lda_w

  • aos_vc_alpha_lda_w

  • aos_vc_alpha_pbe_w

  • aos_vc_alpha_sr_pbe_w

  • aos_vxc_alpha_lda_w

  • aos_vxc_alpha_pbe_w

  • aos_vxc_alpha_sr_pbe_w

  • elec_beta_num_grid_becke

  • energy_c_lda

  • energy_c_sr_lda

  • energy_x_lda

  • energy_x_pbe

  • energy_x_sr_lda

  • energy_x_sr_pbe

  • f_psi_cas_ab

  • f_psi_hf_ab

  • final_grid_points_transp

  • mo_grad_ints

  • mos_in_r_array

  • mos_in_r_array_omp

  • mu_average_prov

  • mu_grad_rho

  • mu_of_r_dft_average

  • mu_rsc_of_r

  • one_e_dm_and_grad_alpha_in_r

final_weight_at_r_vector_extra

File : becke_numerical_grid/extra_grid_vector.irp.f

   double precision, allocatable   :: final_grid_points_extra      (3,n_points_extra_final_grid)
   double precision, allocatable   :: final_weight_at_r_vector_extra       (n_points_extra_final_grid)
   integer, allocatable    :: index_final_points_extra     (3,n_points_extra_final_grid)
   integer, allocatable    :: index_final_points_extra_reverse     (n_points_extra_integration_angular,n_points_extra_radial_grid,nucl_num)


final_grid_points_extra(1:3,j) = (/ x, y, z /) of the jth grid point

final_weight_at_r_vector_extra(i) = Total weight function of the ith grid point which contains the Lebedev, Voronoi and radial weights contributions

index_final_points_extra(1:3,i) = gives the angular, radial and atomic indices associated to the ith grid point

index_final_points_extra_reverse(i,j,k) = index of the grid point having i as angular, j as radial and l as atomic indices

Needs:

  • final_weight_at_r_extra

  • grid_points_extra_per_atom

  • n_points_extra_final_grid

  • n_points_extra_radial_grid

  • nucl_num

  • thresh_extra_grid

Needed by:

  • aos_in_r_array_extra

final_weight_at_r_vector_per_atom

File : becke_numerical_grid/grid_becke_per_atom.irp.f

double precision, allocatable   :: final_grid_points_per_atom   (3,n_pts_max_per_atom,nucl_num)
double precision, allocatable   :: final_weight_at_r_vector_per_atom    (n_pts_max_per_atom,nucl_num)
integer, allocatable    :: index_final_points_per_atom  (3,n_pts_max_per_atom,nucl_num)
integer, allocatable    :: index_final_points_per_atom_reverse  (n_points_integration_angular,n_points_radial_grid,nucl_num)

Needs:

  • final_weight_at_r

  • grid_points_per_atom

  • n_points_radial_grid

  • n_pts_per_atom

  • nucl_num

  • thresh_grid

grid_atomic_number

File : becke_numerical_grid/atomic_number.irp.f

integer, allocatable    :: grid_atomic_number   (nucl_num)

Atomic number used to adjust the grid

Needs:

  • nucl_charge

  • nucl_num

Needed by:

  • final_weight_at_r

  • final_weight_at_r_extra

  • grid_points_extra_per_atom

  • grid_points_per_atom

grid_points_extra_per_atom

File : becke_numerical_grid/extra_grid.irp.f

double precision, allocatable   :: grid_points_extra_per_atom   (3,n_points_extra_integration_angular,n_points_extra_radial_grid,nucl_num)

x,y,z coordinates of grid points_extra used for integration in 3d space

Needs:

  • alpha_knowles

  • angular_quadrature_points_extra

  • extra_rad_grid_type

  • grid_atomic_number

  • grid_points_extra_radial

  • m_knowles

  • n_points_extra_radial_grid

  • nucl_coord

  • nucl_num

  • r_gill

Needed by:

  • final_grid_points_extra

  • weight_at_r_extra

grid_points_extra_radial

File : becke_numerical_grid/extra_grid.irp.f

double precision, allocatable   :: grid_points_extra_radial     (n_points_extra_radial_grid)
double precision        :: dr_radial_extra_integral

points_extra in [0,1] to map the radial integral [0,infty]

Needs:

  • n_points_extra_radial_grid

Needed by:

  • final_weight_at_r_extra

  • grid_points_extra_per_atom

grid_points_per_atom

File : becke_numerical_grid/grid_becke.irp.f

double precision, allocatable   :: grid_points_per_atom (3,n_points_integration_angular,n_points_radial_grid,nucl_num)

x,y,z coordinates of grid points used for integration in 3d space

Needs:

  • alpha_knowles

  • angular_quadrature_points

  • grid_atomic_number

  • grid_points_radial

  • m_knowles

  • n_points_radial_grid

  • nucl_coord

  • nucl_num

  • r_gill

  • rad_grid_type

Needed by:

  • final_grid_points

  • final_grid_points_per_atom

  • weight_at_r

grid_points_radial

File : becke_numerical_grid/grid_becke.irp.f

double precision, allocatable   :: grid_points_radial   (n_points_radial_grid)
double precision        :: dr_radial_integral

points in [0,1] to map the radial integral [0,infty]

Needs:

  • n_points_radial_grid

Needed by:

  • final_weight_at_r

  • grid_points_per_atom

index_final_points

File : becke_numerical_grid/grid_becke_vector.irp.f

   double precision, allocatable   :: final_grid_points    (3,n_points_final_grid)
   double precision, allocatable   :: final_weight_at_r_vector     (n_points_final_grid)
   integer, allocatable    :: index_final_points   (3,n_points_final_grid)
   integer, allocatable    :: index_final_points_reverse   (n_points_integration_angular,n_points_radial_grid,nucl_num)


final_grid_points(1:3,j) = (/ x, y, z /) of the jth grid point

final_weight_at_r_vector(i) = Total weight function of the ith grid point which contains the Lebedev, Voronoi and radial weights contributions

index_final_points(1:3,i) = gives the angular, radial and atomic indices associated to the ith grid point

index_final_points_reverse(i,j,k) = index of the grid point having i as angular, j as radial and l as atomic indices

Needs:

  • final_weight_at_r

  • grid_points_per_atom

  • n_points_final_grid

  • n_points_radial_grid

  • nucl_num

  • thresh_grid

Needed by:

  • ao_abs_int_grid

  • ao_overlap_abs_grid

  • ao_prod_abs_r

  • ao_prod_center

  • ao_prod_dist_grid

  • aos_grad_in_r_array

  • aos_in_r_array

  • aos_lapl_in_r_array

  • aos_sr_vc_alpha_lda_w

  • aos_sr_vxc_alpha_lda_w

  • aos_vc_alpha_lda_w

  • aos_vc_alpha_pbe_w

  • aos_vc_alpha_sr_pbe_w

  • aos_vxc_alpha_lda_w

  • aos_vxc_alpha_pbe_w

  • aos_vxc_alpha_sr_pbe_w

  • elec_beta_num_grid_becke

  • energy_c_lda

  • energy_c_sr_lda

  • energy_x_lda

  • energy_x_pbe

  • energy_x_sr_lda

  • energy_x_sr_pbe

  • f_psi_cas_ab

  • f_psi_hf_ab

  • final_grid_points_transp

  • mo_grad_ints

  • mos_in_r_array

  • mos_in_r_array_omp

  • mu_average_prov

  • mu_grad_rho

  • mu_of_r_dft_average

  • mu_rsc_of_r

  • one_e_dm_and_grad_alpha_in_r

index_final_points_extra

File : becke_numerical_grid/extra_grid_vector.irp.f

   double precision, allocatable   :: final_grid_points_extra      (3,n_points_extra_final_grid)
   double precision, allocatable   :: final_weight_at_r_vector_extra       (n_points_extra_final_grid)
   integer, allocatable    :: index_final_points_extra     (3,n_points_extra_final_grid)
   integer, allocatable    :: index_final_points_extra_reverse     (n_points_extra_integration_angular,n_points_extra_radial_grid,nucl_num)


final_grid_points_extra(1:3,j) = (/ x, y, z /) of the jth grid point

final_weight_at_r_vector_extra(i) = Total weight function of the ith grid point which contains the Lebedev, Voronoi and radial weights contributions

index_final_points_extra(1:3,i) = gives the angular, radial and atomic indices associated to the ith grid point

index_final_points_extra_reverse(i,j,k) = index of the grid point having i as angular, j as radial and l as atomic indices

Needs:

  • final_weight_at_r_extra

  • grid_points_extra_per_atom

  • n_points_extra_final_grid

  • n_points_extra_radial_grid

  • nucl_num

  • thresh_extra_grid

Needed by:

  • aos_in_r_array_extra

index_final_points_extra_reverse

File : becke_numerical_grid/extra_grid_vector.irp.f

   double precision, allocatable   :: final_grid_points_extra      (3,n_points_extra_final_grid)
   double precision, allocatable   :: final_weight_at_r_vector_extra       (n_points_extra_final_grid)
   integer, allocatable    :: index_final_points_extra     (3,n_points_extra_final_grid)
   integer, allocatable    :: index_final_points_extra_reverse     (n_points_extra_integration_angular,n_points_extra_radial_grid,nucl_num)


final_grid_points_extra(1:3,j) = (/ x, y, z /) of the jth grid point

final_weight_at_r_vector_extra(i) = Total weight function of the ith grid point which contains the Lebedev, Voronoi and radial weights contributions

index_final_points_extra(1:3,i) = gives the angular, radial and atomic indices associated to the ith grid point

index_final_points_extra_reverse(i,j,k) = index of the grid point having i as angular, j as radial and l as atomic indices

Needs:

  • final_weight_at_r_extra

  • grid_points_extra_per_atom

  • n_points_extra_final_grid

  • n_points_extra_radial_grid

  • nucl_num

  • thresh_extra_grid

Needed by:

  • aos_in_r_array_extra

index_final_points_per_atom

File : becke_numerical_grid/grid_becke_per_atom.irp.f

double precision, allocatable   :: final_grid_points_per_atom   (3,n_pts_max_per_atom,nucl_num)
double precision, allocatable   :: final_weight_at_r_vector_per_atom    (n_pts_max_per_atom,nucl_num)
integer, allocatable    :: index_final_points_per_atom  (3,n_pts_max_per_atom,nucl_num)
integer, allocatable    :: index_final_points_per_atom_reverse  (n_points_integration_angular,n_points_radial_grid,nucl_num)

Needs:

  • final_weight_at_r

  • grid_points_per_atom

  • n_points_radial_grid

  • n_pts_per_atom

  • nucl_num

  • thresh_grid

index_final_points_per_atom_reverse

File : becke_numerical_grid/grid_becke_per_atom.irp.f

double precision, allocatable   :: final_grid_points_per_atom   (3,n_pts_max_per_atom,nucl_num)
double precision, allocatable   :: final_weight_at_r_vector_per_atom    (n_pts_max_per_atom,nucl_num)
integer, allocatable    :: index_final_points_per_atom  (3,n_pts_max_per_atom,nucl_num)
integer, allocatable    :: index_final_points_per_atom_reverse  (n_points_integration_angular,n_points_radial_grid,nucl_num)

Needs:

  • final_weight_at_r

  • grid_points_per_atom

  • n_points_radial_grid

  • n_pts_per_atom

  • nucl_num

  • thresh_grid

index_final_points_reverse

File : becke_numerical_grid/grid_becke_vector.irp.f

   double precision, allocatable   :: final_grid_points    (3,n_points_final_grid)
   double precision, allocatable   :: final_weight_at_r_vector     (n_points_final_grid)
   integer, allocatable    :: index_final_points   (3,n_points_final_grid)
   integer, allocatable    :: index_final_points_reverse   (n_points_integration_angular,n_points_radial_grid,nucl_num)


final_grid_points(1:3,j) = (/ x, y, z /) of the jth grid point

final_weight_at_r_vector(i) = Total weight function of the ith grid point which contains the Lebedev, Voronoi and radial weights contributions

index_final_points(1:3,i) = gives the angular, radial and atomic indices associated to the ith grid point

index_final_points_reverse(i,j,k) = index of the grid point having i as angular, j as radial and l as atomic indices

Needs:

  • final_weight_at_r

  • grid_points_per_atom

  • n_points_final_grid

  • n_points_radial_grid

  • nucl_num

  • thresh_grid

Needed by:

  • ao_abs_int_grid

  • ao_overlap_abs_grid

  • ao_prod_abs_r

  • ao_prod_center

  • ao_prod_dist_grid

  • aos_grad_in_r_array

  • aos_in_r_array

  • aos_lapl_in_r_array

  • aos_sr_vc_alpha_lda_w

  • aos_sr_vxc_alpha_lda_w

  • aos_vc_alpha_lda_w

  • aos_vc_alpha_pbe_w

  • aos_vc_alpha_sr_pbe_w

  • aos_vxc_alpha_lda_w

  • aos_vxc_alpha_pbe_w

  • aos_vxc_alpha_sr_pbe_w

  • elec_beta_num_grid_becke

  • energy_c_lda

  • energy_c_sr_lda

  • energy_x_lda

  • energy_x_pbe

  • energy_x_sr_lda

  • energy_x_sr_pbe

  • f_psi_cas_ab

  • f_psi_hf_ab

  • final_grid_points_transp

  • mo_grad_ints

  • mos_in_r_array

  • mos_in_r_array_omp

  • mu_average_prov

  • mu_grad_rho

  • mu_of_r_dft_average

  • mu_rsc_of_r

  • one_e_dm_and_grad_alpha_in_r

m_knowles

File : becke_numerical_grid/grid_becke.irp.f

integer :: m_knowles

value of the “m” parameter in the equation (7) of the paper of Knowles (JCP, 104, 1996)

Needed by:

  • final_weight_at_r

  • final_weight_at_r_extra

  • grid_points_extra_per_atom

  • grid_points_per_atom

n_points_extra_final_grid

File : becke_numerical_grid/extra_grid_vector.irp.f

integer :: n_points_extra_final_grid

Number of points_extra which are non zero

Needs:

  • final_weight_at_r_extra

  • n_points_extra_radial_grid

  • nucl_num

  • thresh_extra_grid

Needed by:

  • aos_in_r_array_extra

  • aos_in_r_array_extra_transp

  • final_grid_points_extra

n_points_extra_grid_per_atom

File : becke_numerical_grid/extra_grid.irp.f

integer :: n_points_extra_grid_per_atom

Number of grid points_extra per atom

Needs:

  • n_points_extra_radial_grid

n_points_extra_integration_angular

File : becke_numerical_grid/extra_grid.irp.f

integer :: n_points_extra_radial_grid
integer :: n_points_extra_integration_angular

n_points_extra_radial_grid = number of radial grid points_extra per atom

n_points_extra_integration_angular = number of angular grid points_extra per atom

These numbers are automatically set by setting the grid_type_sgn parameter

Needs:

  • extra_grid_type_sgn

  • my_extra_grid_becke

  • my_n_pt_a_extra_grid

  • my_n_pt_r_extra_grid

Needed by:

  • angular_quadrature_points_extra

  • final_grid_points_extra

  • final_weight_at_r_extra

  • grid_points_extra_per_atom

  • grid_points_extra_radial

  • n_points_extra_final_grid

  • n_points_extra_grid_per_atom

  • weight_at_r_extra

n_points_extra_radial_grid

File : becke_numerical_grid/extra_grid.irp.f

integer :: n_points_extra_radial_grid
integer :: n_points_extra_integration_angular

n_points_extra_radial_grid = number of radial grid points_extra per atom

n_points_extra_integration_angular = number of angular grid points_extra per atom

These numbers are automatically set by setting the grid_type_sgn parameter

Needs:

  • extra_grid_type_sgn

  • my_extra_grid_becke

  • my_n_pt_a_extra_grid

  • my_n_pt_r_extra_grid

Needed by:

  • angular_quadrature_points_extra

  • final_grid_points_extra

  • final_weight_at_r_extra

  • grid_points_extra_per_atom

  • grid_points_extra_radial

  • n_points_extra_final_grid

  • n_points_extra_grid_per_atom

  • weight_at_r_extra

n_points_final_grid

File : becke_numerical_grid/grid_becke_vector.irp.f

integer :: n_points_final_grid

Number of points which are non zero

Needs:

  • final_weight_at_r

  • n_points_radial_grid

  • nucl_num

  • thresh_grid

Needed by:

  • act_mos_in_r_array

  • alpha_dens_kin_in_r

  • ao_abs_int_grid

  • ao_overlap_abs_grid

  • ao_prod_abs_r

  • ao_prod_center

  • ao_prod_dist_grid

  • aos_grad_in_r_array

  • aos_grad_in_r_array_transp

  • aos_grad_in_r_array_transp_3

  • aos_grad_in_r_array_transp_bis

  • aos_in_r_array

  • aos_in_r_array_transp

  • aos_lapl_in_r_array

  • aos_lapl_in_r_array_transp

  • aos_sr_vc_alpha_lda_w

  • aos_sr_vxc_alpha_lda_w

  • aos_vc_alpha_lda_w

  • aos_vc_alpha_pbe_w

  • aos_vc_alpha_sr_pbe_w

  • aos_vxc_alpha_lda_w

  • aos_vxc_alpha_pbe_w

  • aos_vxc_alpha_sr_pbe_w

  • basis_mos_in_r_array

  • core_density

  • core_inact_act_mos_grad_in_r_array

  • core_inact_act_mos_in_r_array

  • core_inact_act_v_kl_contracted

  • core_mos_in_r_array

  • effective_alpha_dm

  • effective_spin_dm

  • elec_beta_num_grid_becke

  • energy_c_lda

  • energy_c_sr_lda

  • energy_x_lda

  • energy_x_pbe

  • energy_x_sr_lda

  • energy_x_sr_pbe

  • f_psi_cas_ab

  • f_psi_cas_ab_old

  • f_psi_hf_ab

  • final_grid_points

  • final_grid_points_transp

  • full_occ_2_rdm_cntrctd

  • full_occ_2_rdm_cntrctd_trans

  • full_occ_v_kl_cntrctd

  • grad_total_cas_on_top_density

  • inact_density

  • inact_mos_in_r_array

  • kinetic_density_generalized

  • mo_grad_ints

  • mos_grad_in_r_array

  • mos_grad_in_r_array_tranp

  • mos_grad_in_r_array_transp_3

  • mos_grad_in_r_array_transp_bis

  • mos_in_r_array

  • mos_in_r_array_omp

  • mos_in_r_array_transp

  • mos_lapl_in_r_array

  • mos_lapl_in_r_array_tranp

  • mu_average_prov

  • mu_grad_rho

  • mu_of_r_dft

  • mu_of_r_dft_average

  • mu_of_r_hf

  • mu_of_r_prov

  • mu_of_r_psi_cas

  • mu_rsc_of_r

  • one_e_act_density_alpha

  • one_e_act_density_beta

  • one_e_cas_total_density

  • one_e_dm_and_grad_alpha_in_r

  • pot_grad_x_alpha_ao_pbe

  • pot_grad_x_alpha_ao_sr_pbe

  • pot_grad_xc_alpha_ao_pbe

  • pot_grad_xc_alpha_ao_sr_pbe

  • pot_scal_x_alpha_ao_pbe

  • pot_scal_x_alpha_ao_sr_pbe

  • pot_scal_xc_alpha_ao_pbe

  • pot_scal_xc_alpha_ao_sr_pbe

  • potential_c_alpha_ao_lda

  • potential_c_alpha_ao_sr_lda

  • potential_x_alpha_ao_lda

  • potential_x_alpha_ao_sr_lda

  • potential_xc_alpha_ao_lda

  • potential_xc_alpha_ao_sr_lda

  • total_cas_on_top_density

  • virt_mos_in_r_array

n_points_grid_per_atom

File : becke_numerical_grid/grid_becke.irp.f

integer :: n_points_grid_per_atom

Number of grid points per atom

Needs:

  • n_points_radial_grid

n_points_integration_angular

File : becke_numerical_grid/grid_becke.irp.f

integer :: n_points_radial_grid
integer :: n_points_integration_angular

n_points_radial_grid = number of radial grid points per atom

n_points_integration_angular = number of angular grid points per atom

These numbers are automatically set by setting the grid_type_sgn parameter

Needs:

  • grid_type_sgn

  • my_grid_becke

  • my_n_pt_a_grid

  • my_n_pt_r_grid

Needed by:

  • angular_quadrature_points

  • final_grid_points

  • final_grid_points_per_atom

  • final_weight_at_r

  • grid_points_per_atom

  • grid_points_radial

  • n_points_final_grid

  • n_points_grid_per_atom

  • n_pts_per_atom

  • weight_at_r

n_points_radial_grid

File : becke_numerical_grid/grid_becke.irp.f

integer :: n_points_radial_grid
integer :: n_points_integration_angular

n_points_radial_grid = number of radial grid points per atom

n_points_integration_angular = number of angular grid points per atom

These numbers are automatically set by setting the grid_type_sgn parameter

Needs:

  • grid_type_sgn

  • my_grid_becke

  • my_n_pt_a_grid

  • my_n_pt_r_grid

Needed by:

  • angular_quadrature_points

  • final_grid_points

  • final_grid_points_per_atom

  • final_weight_at_r

  • grid_points_per_atom

  • grid_points_radial

  • n_points_final_grid

  • n_points_grid_per_atom

  • n_pts_per_atom

  • weight_at_r

n_pts_max_per_atom

File : becke_numerical_grid/grid_becke_per_atom.irp.f

integer, allocatable    :: n_pts_per_atom       (nucl_num)
integer :: n_pts_max_per_atom

Number of points which are non zero

Needs:

  • final_weight_at_r

  • n_points_radial_grid

  • nucl_num

  • thresh_grid

Needed by:

  • final_grid_points_per_atom

n_pts_per_atom

File : becke_numerical_grid/grid_becke_per_atom.irp.f

integer, allocatable    :: n_pts_per_atom       (nucl_num)
integer :: n_pts_max_per_atom

Number of points which are non zero

Needs:

  • final_weight_at_r

  • n_points_radial_grid

  • nucl_num

  • thresh_grid

Needed by:

  • final_grid_points_per_atom

r_gill

File : becke_numerical_grid/grid_becke.irp.f

double precision        :: r_gill

Needed by:

  • final_weight_at_r

  • final_weight_at_r_extra

  • grid_points_extra_per_atom

  • grid_points_per_atom

weight_at_r

File : becke_numerical_grid/grid_becke.irp.f

double precision, allocatable   :: weight_at_r  (n_points_integration_angular,n_points_radial_grid,nucl_num)

Weight function at grid points : w_n(r) according to the equation (22) of Becke original paper (JCP, 88, 1988)

The “n” discrete variable represents the nucleis which in this array is represented by the last dimension and the points are labelled by the other dimensions.

Needs:

  • grid_points_per_atom

  • n_points_radial_grid

  • nucl_coord_transp

  • nucl_dist_inv

  • nucl_num

  • slater_bragg_type_inter_distance_ua

Needed by:

  • final_weight_at_r

weight_at_r_extra

File : becke_numerical_grid/extra_grid.irp.f

double precision, allocatable   :: weight_at_r_extra    (n_points_extra_integration_angular,n_points_extra_radial_grid,nucl_num)

Weight function at grid points_extra : w_n(r) according to the equation (22) of Becke original paper (JCP, 88, 1988)

The “n” discrete variable represents the nucleis which in this array is represented by the last dimension and the points_extra are labelled by the other dimensions.

Needs:

  • grid_points_extra_per_atom

  • n_points_extra_radial_grid

  • nucl_coord_transp

  • nucl_dist_inv

  • nucl_num

  • slater_bragg_type_inter_distance_ua

Needed by:

  • final_weight_at_r_extra

weights_angular_points

File : becke_numerical_grid/angular_grid_pts.irp.f

double precision, allocatable   :: angular_quadrature_points    (n_points_integration_angular,3)
double precision, allocatable   :: weights_angular_points       (n_points_integration_angular)

weights and grid points for the integration on the angular variables on the unit sphere centered on (0,0,0) According to the LEBEDEV scheme

Needs:

  • n_points_radial_grid

Needed by:

  • final_weight_at_r

  • grid_points_per_atom

weights_angular_points_extra

File : becke_numerical_grid/angular_extra_grid.irp.f

double precision, allocatable   :: angular_quadrature_points_extra      (n_points_extra_integration_angular,3)
double precision, allocatable   :: weights_angular_points_extra (n_points_extra_integration_angular)

weights and grid points_extra for the integration on the angular variables on the unit sphere centered on (0,0,0) According to the LEBEDEV scheme

Needs:

  • n_points_extra_radial_grid

Needed by:

  • final_weight_at_r_extra

  • grid_points_extra_per_atom

Subroutines / functions

cell_function_becke:

File : becke_numerical_grid/step_function_becke.irp.f

double precision function cell_function_becke(r, atom_number)

atom_number :: atom on which the cell function of Becke (1988, JCP,88(4)) r(1:3) :: x,y,z coordinantes of the current point

Needs:

  • nucl_coord_transp

  • nucl_dist_inv

  • nucl_num

  • slater_bragg_type_inter_distance_ua

derivative_knowles_function:

File : becke_numerical_grid/integration_radial.irp.f

double precision function derivative_knowles_function(alpha, m, x)

Derivative of the function proposed by Knowles (JCP, 104, 1996) for distributing the radial points

example_becke_numerical_grid:

File : becke_numerical_grid/example.irp.f

subroutine example_becke_numerical_grid

subroutine that illustrates the main features available in becke_numerical_grid

Needs:

  • final_grid_points

  • final_weight_at_r

  • grid_points_per_atom

  • n_points_final_grid

  • n_points_radial_grid

  • nucl_coord

  • nucl_num

f_function_becke:

File : becke_numerical_grid/step_function_becke.irp.f

double precision function f_function_becke(x)
knowles_function:

File : becke_numerical_grid/integration_radial.irp.f

double precision function knowles_function(alpha, m, x)

Function proposed by Knowles (JCP, 104, 1996) for distributing the radial points : the Log “m” function ( equation (7) in the paper )

step_function_becke:

File : becke_numerical_grid/step_function_becke.irp.f

double precision function step_function_becke(x)

Step function of the Becke paper (1988, JCP,88(4))