This module contains routines to calculate the 1-electron integrals in the B-spline basis.
More...
|
subroutine, public | BB_shell_integrals (shell_A, shell_B, starting_index_A, starting_index_B, a, max_property_l, property_center, symmetry_data, olap_column, kei_column, prop_column, nai_column, one_elham_column, int_index, integrals) |
|
subroutine, public | BB_initialize (bspline_grid, max_bspline_l, max_prop_l, dipole_damp_factor, symmetry_data) |
| Calculates: r_points, weights, B_vals, temp_r, bspline_boundary_val, bspline_start_end, prop_on_grid, leg_on_grid, re_sph_harm_center. More...
|
|
subroutine, public | construct_bspline_quadrature_grid (knots, x, w, n, n_rng_knot, r_points, weights, n_total_points) |
| Constructs a quadrature grid for the B-spline basis described by a given knot sequence. The final quadrature grid comprises effectively a series of quadratures over subintervals generated in between each pair of distinct knots. More...
|
|
subroutine | index_1el (sph_shell_A, starting_index_A, sph_shell_B, starting_index_B, n_repeat, int_index) |
|
This module contains routines to calculate the 1-electron integrals in the B-spline basis.
- Todo:
- Add checking that the radial quadratures have high enough order to get exact results. Integrals involving radial B-splines with a non-zero derivative at r=A (i.e. at the origin of the B-spline grid) are not evaluated.
◆ BB_initialize()
subroutine, public bto_integrals_gbl::BB_initialize |
( |
type(bspline_grid_obj), intent(inout) |
bspline_grid, |
|
|
integer, intent(in), optional |
max_bspline_l, |
|
|
integer, intent(in), optional |
max_prop_l, |
|
|
real(kind=cfp), intent(in) |
dipole_damp_factor, |
|
|
type(symmetry_obj), intent(in) |
symmetry_data |
|
) |
| |
Calculates: r_points, weights, B_vals, temp_r, bspline_boundary_val, bspline_start_end, prop_on_grid, leg_on_grid, re_sph_harm_center.
◆ BB_shell_integrals()
subroutine, public bto_integrals_gbl::BB_shell_integrals |
( |
type(bto_shell_data_obj), intent(in), target |
shell_A, |
|
|
type(bto_shell_data_obj), intent(in), target |
shell_B, |
|
|
integer, intent(in) |
starting_index_A, |
|
|
integer, intent(in) |
starting_index_B, |
|
|
real(kind=cfp), intent(in) |
a, |
|
|
integer, intent(in) |
max_property_l, |
|
|
real(kind=cfp), dimension(3), intent(in) |
property_center, |
|
|
type(symmetry_obj), intent(in) |
symmetry_data, |
|
|
integer, intent(in) |
olap_column, |
|
|
integer, intent(in) |
kei_column, |
|
|
integer, intent(in) |
prop_column, |
|
|
integer, intent(in) |
nai_column, |
|
|
integer, intent(in) |
one_elham_column, |
|
|
integer, dimension(:,:), allocatable |
int_index, |
|
|
real(kind=cfp), dimension(:,:), allocatable |
integrals |
|
) |
| |
◆ construct_bspline_quadrature_grid()
subroutine, public bto_integrals_gbl::construct_bspline_quadrature_grid |
( |
real(kind=cfp), dimension(:), intent(in) |
knots, |
|
|
real(kind=cfp), dimension(2*n+1), intent(in) |
x, |
|
|
real(kind=cfp), dimension(2*n+1), intent(in) |
w, |
|
|
integer, intent(in) |
n, |
|
|
integer, intent(in) |
n_rng_knot, |
|
|
real(kind=cfp), dimension(:), allocatable |
r_points, |
|
|
real(kind=cfp), dimension(:), allocatable |
weights, |
|
|
integer, intent(out) |
n_total_points |
|
) |
| |
Constructs a quadrature grid for the B-spline basis described by a given knot sequence. The final quadrature grid comprises effectively a series of quadratures over subintervals generated in between each pair of distinct knots.
◆ index_1el()
subroutine bto_integrals_gbl::index_1el |
( |
integer, intent(in) |
sph_shell_A, |
|
|
integer, intent(in) |
starting_index_A, |
|
|
integer, intent(in) |
sph_shell_B, |
|
|
integer, intent(in) |
starting_index_B, |
|
|
integer, intent(in) |
n_repeat, |
|
|
integer, dimension(:,:), allocatable |
int_index |
|
) |
| |