GBTOlib: library for evaluation of molecular integrals in mixed Gaussian / B-spline basis 111
Loading...
Searching...
No Matches
grid_gbl::grid_r1_r2_obj Type Reference
Collaboration diagram for grid_gbl::grid_r1_r2_obj:

Public Member Functions

procedure init (this, bspline_grid, first_bspline_index, max_bspline_l, max_prop_l, nuclei, delta_r1, n1, n2, only_on_bto_grid, rmat_radius, dipole_damp_factor, max_l_legendre)
 Initializes the grid given the B-spline grid, the elementary quadrature rules, etc.
procedure construct_r1_r2_grids (this, inp_bspline_grid, first_bspline_index, max_bspline_l, max_prop_l, dipole_damp_factor, max_l_legendre, nuclei, delta_r1, x1, w1, n1, x2, w2, n2, only_on_bto_grid, rmat_radius)
 Constructs the quadrature grids appropriate for the Legendre expansion and evaluates other auxiliary arrays.
procedure eval_lambda_l_BB eval_lambda_l_bb
 Calculates the Y function for each unique pair of radial B-splines.
procedure construct_lebedev_grid (this, n)
 Construct the Lebedev grid of a given (or nearest larger) order.
procedure construct_GK_angular_grid construct_gk_angular_grid
 Constructs angular grid quadrature grid based on the Gauss-Kronrod product rule.
procedure eval_Xlm_on_lebedev_grid eval_xlm_on_lebedev_grid
 Construct the Lebedev grid and evaluate the real spherical harmonics on the angular Lebedev grid.
procedure final (this)
 Deallocates everything.

Public Attributes

real(kind=cfp), dimension(:), allocatable r1
 Quadrature points and weights for r1 and r2 radial grids for the Legendre expansion. The r1 grid spans the radial distance from the start of the B-spline grid up to the larger of: the end of the B-spline grid and the (adjusted) R-matrix radius. The value of the adjusted R-matrix radius is obtained as the end point of the last B-spline that starts in the inner region and extends to the outer region. The r2 grid is restricted to the (adjusted) inner region only: it is not needed in the outer region.
real(kind=cfp), dimension(:), allocatable w1
real(kind=cfp), dimension(:), allocatable r2
real(kind=cfp), dimension(:), allocatable w2
integer n1_total_points
integer n2_total_points
integer max_bspline_l
integer max_prop_l
integer first_bspline_index
integer n_unique_pairs
integer max_l_legendre
type(bspline_grid_objbspline_grid
 Grid of radial B-splines which will be used to calculate the auxiliary functions.
integer last_bspline_inner = 0
 Index of the last B-spline which lies in the inner region (r <= R-matrix sphere).
real(kind=cfp), dimension(:,:), allocatable b_vals_r1
 Normalized B-splines (and pairs of B-splines) evaluated on the r1 and r2 grids. We don't need BB pairs on the r2 grid.
real(kind=cfp), dimension(:,:), allocatable b_vals_r2
real(kind=cfp), dimension(:,:), allocatable bb_vals_r1
real(kind=cfp), dimension(:,:), allocatable bto_radial_olap
 The evaluation of the 1-electron mixed integrals involves quadrature over r1 grid involving a product of GTO-related part and a BTO-related part. The part related to the BTOs is stored in the arrays below. The radial part of the BTO is defined as: B(r)/r. The values stored include multiplication by the radial part of the Jacobian (r**2) and the quadrature weights for the r1 grid.
real(kind=cfp), dimension(:,:,:), allocatable bto_radial_kei
real(kind=cfp), dimension(:,:,:), allocatable bto_radial_prop
real(kind=cfp), dimension(:,:), allocatable bto_end_points
integer, dimension(:,:), allocatable bspline_start_end_r1
 Starting and ending points of the B-splines (and pairs of B-splines) on the r1 and r2 grids. The value in bspline_start_end_*(3,ind) contains the index of the last end-point of the B-spline with index ind such that it is never larger than the last point in the inner region.
integer, dimension(:,:), allocatable bspline_start_end_r2
integer, dimension(:,:), allocatable bb_start_end_r1
integer, dimension(:,:), allocatable pairs_to_unique_pairs
 For a given pair of B-splines it contains the index of the corresponding unique pair of B-splines.
integer, dimension(:,:), allocatable unique_pairs_to_pairs
 For a given index of the unique pair of B-splines it contains the indices of the two corresponding B-splines.
real(kind=cfp), dimension(:,:,:), allocatable lambda_bb_r1_l_ij
 Y_l function on the r1 grid for all pairs of B-splines.
real(kind=cfp), dimension(:), allocatable xlm_lebedev
 Lebedev grid and the real spherical harmonics evaluated at the Lebedev angular points.
real(kind=cfp), dimension(:), allocatable leb_r1
real(kind=cfp), dimension(:), allocatable leb_r2
real(kind=cfp), dimension(:), allocatable leb_r3
real(kind=cfp), dimension(:), allocatable leb_w
real(kind=cfp), dimension(:), allocatable bto_norm
 Normalization factors for the radial B-splines.
integer lebedev_order = -1
 Number of points of the chosen Lebedev rule.

Member Function/Subroutine Documentation

◆ construct_GK_angular_grid()

procedure grid_gbl::grid_r1_r2_obj::construct_GK_angular_grid

Constructs angular grid quadrature grid based on the Gauss-Kronrod product rule.

◆ construct_lebedev_grid()

procedure grid_gbl::grid_r1_r2_obj::construct_lebedev_grid ( class(grid_r1_r2_obj) this,
integer, intent(in) n )

Construct the Lebedev grid of a given (or nearest larger) order.

◆ construct_r1_r2_grids()

procedure grid_gbl::grid_r1_r2_obj::construct_r1_r2_grids ( class(grid_r1_r2_obj) this,
type(bspline_grid_obj) inp_bspline_grid,
integer, intent(in) first_bspline_index,
integer, intent(in) max_bspline_l,
integer, intent(in) max_prop_l,
real(kind=cfp), intent(in) dipole_damp_factor,
integer, intent(in) max_l_legendre,
type(nucleus_type), dimension(:) nuclei,
real(kind=cfp), intent(in) delta_r1,
real(kind=cfp), dimension(n1), intent(in) x1,
real(kind=cfp), dimension(n1), intent(in) w1,
integer, intent(in) n1,
real(kind=cfp), dimension(n2), intent(in) x2,
real(kind=cfp), dimension(n2), intent(in) w2,
integer, intent(in) n2,
logical, intent(in) only_on_bto_grid,
real(kind=cfp), intent(in) rmat_radius )

Constructs the quadrature grids appropriate for the Legendre expansion and evaluates other auxiliary arrays.

◆ eval_lambda_l_BB()

procedure grid_gbl::grid_r1_r2_obj::eval_lambda_l_BB

Calculates the Y function for each unique pair of radial B-splines.

◆ eval_Xlm_on_lebedev_grid()

procedure grid_gbl::grid_r1_r2_obj::eval_Xlm_on_lebedev_grid

Construct the Lebedev grid and evaluate the real spherical harmonics on the angular Lebedev grid.

◆ final()

procedure grid_gbl::grid_r1_r2_obj::final ( class(grid_r1_r2_obj) this)

Deallocates everything.

◆ init()

procedure grid_gbl::grid_r1_r2_obj::init ( class(grid_r1_r2_obj) this,
type(bspline_grid_obj) bspline_grid,
integer, intent(in) first_bspline_index,
integer, intent(in) max_bspline_l,
integer, intent(in) max_prop_l,
type(nucleus_type), dimension(:) nuclei,
real(kind=cfp), intent(in) delta_r1,
integer, intent(in) n1,
integer, intent(in) n2,
logical, intent(in) only_on_bto_grid,
real(kind=cfp), intent(in) rmat_radius,
real(kind=cfp), intent(in) dipole_damp_factor,
integer, intent(in) max_l_legendre )

Initializes the grid given the B-spline grid, the elementary quadrature rules, etc.

Member Data Documentation

◆ b_vals_r1

real(kind=cfp), dimension(:,:), allocatable grid_gbl::grid_r1_r2_obj::b_vals_r1

Normalized B-splines (and pairs of B-splines) evaluated on the r1 and r2 grids. We don't need BB pairs on the r2 grid.

◆ b_vals_r2

real(kind=cfp), dimension(:,:), allocatable grid_gbl::grid_r1_r2_obj::b_vals_r2

◆ bb_start_end_r1

integer, dimension(:,:), allocatable grid_gbl::grid_r1_r2_obj::bb_start_end_r1

◆ bb_vals_r1

real(kind=cfp), dimension(:,:), allocatable grid_gbl::grid_r1_r2_obj::bb_vals_r1

◆ bspline_grid

type(bspline_grid_obj) grid_gbl::grid_r1_r2_obj::bspline_grid

Grid of radial B-splines which will be used to calculate the auxiliary functions.

◆ bspline_start_end_r1

integer, dimension(:,:), allocatable grid_gbl::grid_r1_r2_obj::bspline_start_end_r1

Starting and ending points of the B-splines (and pairs of B-splines) on the r1 and r2 grids. The value in bspline_start_end_*(3,ind) contains the index of the last end-point of the B-spline with index ind such that it is never larger than the last point in the inner region.

◆ bspline_start_end_r2

integer, dimension(:,:), allocatable grid_gbl::grid_r1_r2_obj::bspline_start_end_r2

◆ bto_end_points

real(kind=cfp), dimension(:,:), allocatable grid_gbl::grid_r1_r2_obj::bto_end_points

◆ bto_norm

real(kind=cfp), dimension(:), allocatable grid_gbl::grid_r1_r2_obj::bto_norm

Normalization factors for the radial B-splines.

◆ bto_radial_kei

real(kind=cfp), dimension(:,:,:), allocatable grid_gbl::grid_r1_r2_obj::bto_radial_kei

◆ bto_radial_olap

real(kind=cfp), dimension(:,:), allocatable grid_gbl::grid_r1_r2_obj::bto_radial_olap

The evaluation of the 1-electron mixed integrals involves quadrature over r1 grid involving a product of GTO-related part and a BTO-related part. The part related to the BTOs is stored in the arrays below. The radial part of the BTO is defined as: B(r)/r. The values stored include multiplication by the radial part of the Jacobian (r**2) and the quadrature weights for the r1 grid.

◆ bto_radial_prop

real(kind=cfp), dimension(:,:,:), allocatable grid_gbl::grid_r1_r2_obj::bto_radial_prop

◆ first_bspline_index

integer grid_gbl::grid_r1_r2_obj::first_bspline_index

◆ lambda_bb_r1_l_ij

real(kind=cfp), dimension(:,:,:), allocatable grid_gbl::grid_r1_r2_obj::lambda_bb_r1_l_ij

Y_l function on the r1 grid for all pairs of B-splines.

◆ last_bspline_inner

integer grid_gbl::grid_r1_r2_obj::last_bspline_inner = 0

Index of the last B-spline which lies in the inner region (r <= R-matrix sphere).

◆ leb_r1

real(kind=cfp), dimension(:), allocatable grid_gbl::grid_r1_r2_obj::leb_r1

◆ leb_r2

real(kind=cfp), dimension(:), allocatable grid_gbl::grid_r1_r2_obj::leb_r2

◆ leb_r3

real(kind=cfp), dimension(:), allocatable grid_gbl::grid_r1_r2_obj::leb_r3

◆ leb_w

real(kind=cfp), dimension(:), allocatable grid_gbl::grid_r1_r2_obj::leb_w

◆ lebedev_order

integer grid_gbl::grid_r1_r2_obj::lebedev_order = -1

Number of points of the chosen Lebedev rule.

◆ max_bspline_l

integer grid_gbl::grid_r1_r2_obj::max_bspline_l

◆ max_l_legendre

integer grid_gbl::grid_r1_r2_obj::max_l_legendre

◆ max_prop_l

integer grid_gbl::grid_r1_r2_obj::max_prop_l

◆ n1_total_points

integer grid_gbl::grid_r1_r2_obj::n1_total_points

◆ n2_total_points

integer grid_gbl::grid_r1_r2_obj::n2_total_points

◆ n_unique_pairs

integer grid_gbl::grid_r1_r2_obj::n_unique_pairs

◆ pairs_to_unique_pairs

integer, dimension(:,:), allocatable grid_gbl::grid_r1_r2_obj::pairs_to_unique_pairs

For a given pair of B-splines it contains the index of the corresponding unique pair of B-splines.

◆ r1

real(kind=cfp), dimension(:), allocatable grid_gbl::grid_r1_r2_obj::r1

Quadrature points and weights for r1 and r2 radial grids for the Legendre expansion. The r1 grid spans the radial distance from the start of the B-spline grid up to the larger of: the end of the B-spline grid and the (adjusted) R-matrix radius. The value of the adjusted R-matrix radius is obtained as the end point of the last B-spline that starts in the inner region and extends to the outer region. The r2 grid is restricted to the (adjusted) inner region only: it is not needed in the outer region.

◆ r2

real(kind=cfp), dimension(:), allocatable grid_gbl::grid_r1_r2_obj::r2

◆ unique_pairs_to_pairs

integer, dimension(:,:), allocatable grid_gbl::grid_r1_r2_obj::unique_pairs_to_pairs

For a given index of the unique pair of B-splines it contains the indices of the two corresponding B-splines.

◆ w1

real(kind=cfp), dimension(:), allocatable grid_gbl::grid_r1_r2_obj::w1

◆ w2

real(kind=cfp), dimension(:), allocatable grid_gbl::grid_r1_r2_obj::w2

◆ xlm_lebedev

real(kind=cfp), dimension(:), allocatable grid_gbl::grid_r1_r2_obj::xlm_lebedev

Lebedev grid and the real spherical harmonics evaluated at the Lebedev angular points.


The documentation for this type was generated from the following file: