|
GBTOlib: library for evaluation of molecular integrals in mixed Gaussian / B-spline basis 111
|


Public Member Functions | |
| procedure | init (this, n, geometry) |
| Allocates space for a given number and type of shells of basis functions. | |
| procedure | final (this) |
| Finalizes the basis set. | |
| procedure | add_cms_gtos (this, min_l, max_l, max_num, no_exps, exponents, n_gto, non_zero_at_boundary) |
| Adds any function which is a primitive gto (single exponent) located at the center of mass, uses add_shell. | |
| procedure | add_shell (this, shell_data) |
| Adds data for one shell into the basis set. All target shells must be added before the continuum ones. This is important for the integral generation algorithms where this situation simplifies a lot what would otherwise be a complicated book-keeping. Internally, the norms of the contracted GTOs are multiplied in with the primitive norms. This strategy is used in the whole basis set to simplify things. | |
| procedure | function_pair_ordered_index (this, i, j) |
| This is an index for a quartet of shells for two-electron integrals. It works for both 1p and 2p in the continuum. It is assumed that on input the quartet of indices i,j,k,l are ordered so that: i.ge.j, k.ge.l – Testing. | |
| procedure | get_function_pair_type (this, i, j) |
| Returns the integer value characterizing the given type of the pair of functions: TT = 1, CT = 2, CC = 3. | |
| procedure | one_electron_integrals (this, integral_storage, integral_options) |
| Calculates and stores 1-electron integrals for all pairs of shells in the basis. | |
| procedure | two_electron_integrals (this, integral_storage, integral_options) |
| Calculates and stores 2-electron integrals for all pairs of shells in the basis. | |
| procedure | integral_index (this, integral_type, bf_indices, two_p_continuum) |
| Calculates index for the given type of integral. The integral type is given by one of the strings definied in the const module and by the number of basis functions on input. The two_p_continuum input variable distinguishes the two indexing methods used. | |
| procedure | assemble_overlap_matrix (this, ao_integrals, overlap_matrix) |
| Takes the p2d_array of integrals on input containing the overlap integrals and returns the 2D array corresponding to the AO overlap matrix. | |
| procedure | get_basis_name (this) |
| Returns the name of the basis set. | |
| procedure | get_shell_name (this, i) |
| Returns the name of the i-th shell in the basis set. | |
| procedure | get_shell_data (this, i, shell_data) |
| Returns the shell data for the i-th shell in the basis set. | |
| procedure | get_all_CGTO_shells get_all_cgto_shells |
| Returns an array containing data for all shells of CGTOs in the basis. | |
| procedure | get_all_BTO_shells get_all_bto_shells |
| Returns an array containing data for all shells of BTOs in the basis. | |
| procedure | get_bspline_grid (this, bspline_grid) |
| Returns the object bspline_grid_obj if the basis contains shells of BTO functions. | |
| procedure | get_number_of_functions_in_shell (this, i) |
| Returns the number of functions in a shell with the given index. | |
| procedure | is_initialized (this) |
| Returns the value of initialized. | |
| procedure | get_symmetry_flags (this, irr, list, from, to) |
| For a given IRR it returns a list of logical values of size equal to the number of functions in the basis. i-th element of the output array is set to .true. if the i-th function in the basis has symmetry irr andcoms from a shell between indices from and to. | |
| procedure | get_continuum_flags (this, irr, list) |
| For a given IRR it returns a list of logical values of size equal to the number of functions in the basis. i-th element of the output array is set to .true. if the i-th function in the basis is a continuum function with the given symmetry. | |
| procedure | get_continuum_l_range (this, min_l, max_l) |
| Returns the values of the smallest and largest L of the continuum functions in the basis. If the basis doesn't contain continuum then it terminates with an error. | |
| procedure | print (this) |
| Prints the basis set data to stdout. | |
| procedure | calculate_amplitudes (this, a, normalize_to_a, amplitudes, continuum_channels) |
| Calculates amplitudes for all basis functions in the basis and builds the channel numbers. The continuum is first normalized to the requested radius. The additional normalization factor for CGTOs is multiplied in with the primitive norms. | |
| procedure | contains_cgtos (this) |
| Returns .true. if the basis contains the CGTOs. | |
| procedure | contains_btos (this) |
| Returns .true. if the basis contains the BTOs. | |
| procedure | normalize_continuum (this, a) |
| Normalizes the atomic continuum functions to the R-matrix sphere with a given radius. The normalization is only done if a > 0.0_cfp. | |
| procedure | eval_basis (this, r, n_points, ao_basis_at_r) |
| Evaluates the AO basis for the given set of points in space. | |
| Public Member Functions inherited from basis_data_generic_gbl::basis_data_generic_obj | |
| procedure(init_bdg_intf), deferred | init (this, n, geometry) |
| Allocates space for a given number and type of shells of basis functions. | |
| procedure(final_bdg_intf), deferred | final (this) |
| Finalizes the basis set. | |
| procedure(add_shell_bdg_intf), deferred | add_shell (this, shell_data) |
| Adds data for one shell into the basis set. | |
| procedure(one_electron_integrals_bdg_intf), deferred | one_electron_integrals (this, integral_storage, integral_options) |
| Calculates and stores 1-electron integrals for all pairs of shells in the basis. | |
| procedure(two_electron_integrals_bdg_intf), deferred | two_electron_integrals (this, integral_storage, integral_options) |
| Calculates and stores 2-electron integrals for all quartets of shells in the basis with possible exclusion of certain classes of integrals as specified on input in integral_options. | |
| procedure(integral_index_bdg_intf), deferred | integral_index (this, integral_type, bf_indices, two_p_continuum) |
| Calculates indices for 1-electron or 2-electron integrals given a list of pairs or quartets of basis functions and specifying the type of the integral to index. The two_p_continuum input variable is used only for AO integrals. | |
| procedure, non_overridable | write (this, path) |
| Writes the basis set data into the given file. | |
| procedure, non_overridable | read (this, path) |
| Reads the basis set data from the given file. | |
| procedure(print_bdg), deferred | print (this) |
| Prints the basis set data to stdout. | |
| procedure(get_basis_name_bdg_intf), deferred | get_basis_name (this) |
| Returns the character string identifying the basis set. This is used for input/output involving disk. | |
| procedure(get_shell_name_bdg_intf), deferred | get_shell_name (this, i) |
| Returns the character string identifying the i-th shell in the basis set. This is used for input/output involving disk. | |
| procedure(get_shell_data_bdg_intf), deferred | get_shell_data (this, i, shell_data) |
| Returns the shell data for the i-th shell in the basis set. This is used for input/output involving disk. | |
| procedure(is_initialized_bdg_intf), deferred | is_initialized (this) |
| Returns .true. if the basis set has been initialized. | |
| procedure(get_continuum_flags_bdg_intf), deferred | get_continuum_flags (this, irr, list) |
| For a given IRR it returns a list of logical values of size equal to the number of functions with that IRR in the basis. i-th element of the output array is set to .true. if the i-th function in the basis is a continuum function. | |
| procedure(amplitudes_intf), deferred | calculate_amplitudes (this, a, normalize_to_a, amplitudes, continuum_channels) |
| Calculates amplitudes for all functions in the basis up to the channel with the given maximum l-value of and radial distance from CMS. | |
Public Attributes | |
| integer, dimension(:,:), allocatable | indices_to_shells |
| Array of dimensions (2,number_of_basis_functions) containing mapping of basis function indices to shell indices (1) and indices within the shells (2). This table is useful for looking up which shells need to be involved in calculation of a given integral and for picking out the requested integral from the integral batch generated over shells of functions. | |
| integer, dimension(:,:), allocatable | shell_descriptor |
| For each shell (2nd dimension) the first dimension contains the following info: row 1: shell type (CGTO: 1, BTO: 2) row 2: index of the shell within its own type row 3: is a continuum shell (1) or not (0) row 4: starting index for the functions in the shell row 5: the number of functions in the shell row 6: the angular momentum of the shell The shell index corresponds to the order in which it was added to the basis. | |
| integer | n_target_fns = -1 |
| Number of target functions, index of the last CT pair of functions and other values needed to index optimally 2-electron 1p continuum integrals. | |
| integer | n_cont_fns = -1 |
| integer | last_ct_fn = -1 |
| integer | n_prec_ints = -1 |
| integer | n_tt_pairs = -1 |
| integer, dimension(:,:), allocatable | ordered_shell_pairs |
| Indices mapping pairs of shells to a single index. This mapping is required to index properly 2p integrals with only 1p in the continuum in case of multiple MPI tasks being used. | |
| integer | n_target_sh = -1 |
| Number of target shells, index of the last CT pair of shells and other values needed to index optimally 2-electron 1p continuum integrals in case of multiple MPI tasks are used. | |
| integer | n_cont_sh = -1 |
| integer | last_ct_sh = -1 |
| integer | n_prec_sh = -1 |
| integer | n_tt_sh_pairs = -1 |
| integer, dimension(:,:), allocatable | shell_pair_indices |
| shell_pair_indices: Mapping between indices of pairs of shells with the indices of the shells within their own type of shells of basis functions. shell_pair_type: Mapping between pairs of shells and their type: GG (1), BG (2), GB (3), BB (4). | |
| integer, dimension(:), allocatable | shell_pair_type |
| integer | n_cgto_functions = -1 |
| Total number of CGTO functions in the basis. | |
| real(kind=cfp), dimension(:,:), allocatable | ao_basis_at_r |
| Auxiliary arrays used by eval_orbital. | |
| real(kind=cfp), dimension(:,:), allocatable | r |
| Public Attributes inherited from basis_data_generic_gbl::basis_data_generic_obj | |
| type(symmetry_obj) | symmetry_data |
| The symmetry specification is fully in control of the user. | |
| integer | number_of_shells = 0 |
| Total number of shells in the basis. | |
| integer | number_of_functions = 0 |
| The number of functions generated by all shells in the basis. | |
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::add_cms_gtos | ( | class(atomic_orbital_basis_obj) | this, |
| integer, intent(in) | min_l, | ||
| integer, intent(in) | max_l, | ||
| integer, intent(in) | max_num, | ||
| integer, dimension(min_l:max_l), intent(in) | no_exps, | ||
| real(kind=cfp), dimension(1:max_num,min_l:max_l), intent(in) | exponents, | ||
| integer, intent(out) | n_gto, | ||
| logical, intent(in) | non_zero_at_boundary ) |
Adds any function which is a primitive gto (single exponent) located at the center of mass, uses add_shell.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::add_shell | ( | class(atomic_orbital_basis_obj) | this, |
| class(shell_data_obj), intent(inout) | shell_data ) |
Adds data for one shell into the basis set. All target shells must be added before the continuum ones. This is important for the integral generation algorithms where this situation simplifies a lot what would otherwise be a complicated book-keeping. Internally, the norms of the contracted GTOs are multiplied in with the primitive norms. This strategy is used in the whole basis set to simplify things.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::assemble_overlap_matrix | ( | class(atomic_orbital_basis_obj) | this, |
| type(p2d_array_obj) | ao_integrals, | ||
| real(kind=cfp), dimension(:,:), allocatable | overlap_matrix ) |
Takes the p2d_array of integrals on input containing the overlap integrals and returns the 2D array corresponding to the AO overlap matrix.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::calculate_amplitudes | ( | class(atomic_orbital_basis_obj) | this, |
| real(kind=cfp), intent(in) | a, | ||
| logical, intent(in) | normalize_to_a, | ||
| real(kind=cfp), dimension(:,:), allocatable | amplitudes, | ||
| integer, dimension(:,:), allocatable | continuum_channels ) |
Calculates amplitudes for all basis functions in the basis and builds the channel numbers. The continuum is first normalized to the requested radius. The additional normalization factor for CGTOs is multiplied in with the primitive norms.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::contains_btos | ( | class(atomic_orbital_basis_obj) | this | ) |
Returns .true. if the basis contains the BTOs.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::contains_cgtos | ( | class(atomic_orbital_basis_obj) | this | ) |
Returns .true. if the basis contains the CGTOs.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::eval_basis | ( | class(atomic_orbital_basis_obj) | this, |
| real(kind=cfp), dimension(3,n_points), intent(in) | r, | ||
| integer, intent(in) | n_points, | ||
| real(kind=cfp), dimension(:,:), allocatable | ao_basis_at_r ) |
Evaluates the AO basis for the given set of points in space.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::final | ( | class(atomic_orbital_basis_obj) | this | ) |
Finalizes the basis set.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::function_pair_ordered_index | ( | class(atomic_orbital_basis_obj), intent(in) | this, |
| integer, intent(in) | i, | ||
| integer, intent(in) | j ) |
This is an index for a quartet of shells for two-electron integrals. It works for both 1p and 2p in the continuum. It is assumed that on input the quartet of indices i,j,k,l are ordered so that: i.ge.j, k.ge.l – Testing.
Returns the ordered index of the given pair of basis functions: the function-pairs are ordered as TT, CT, CC starting from TT.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::get_all_BTO_shells |
Returns an array containing data for all shells of BTOs in the basis.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::get_all_CGTO_shells |
Returns an array containing data for all shells of CGTOs in the basis.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::get_basis_name | ( | class(atomic_orbital_basis_obj) | this | ) |
Returns the name of the basis set.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::get_bspline_grid | ( | class(atomic_orbital_basis_obj) | this, |
| type(bspline_grid_obj) | bspline_grid ) |
Returns the object bspline_grid_obj if the basis contains shells of BTO functions.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::get_continuum_flags | ( | class(atomic_orbital_basis_obj) | this, |
| integer, intent(in) | irr, | ||
| logical, dimension(:), allocatable | list ) |
For a given IRR it returns a list of logical values of size equal to the number of functions in the basis. i-th element of the output array is set to .true. if the i-th function in the basis is a continuum function with the given symmetry.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::get_continuum_l_range | ( | class(atomic_orbital_basis_obj) | this, |
| integer, intent(out) | min_l, | ||
| integer, intent(out) | max_l ) |
Returns the values of the smallest and largest L of the continuum functions in the basis. If the basis doesn't contain continuum then it terminates with an error.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::get_function_pair_type | ( | class(atomic_orbital_basis_obj), intent(in) | this, |
| integer, intent(in) | i, | ||
| integer, intent(in) | j ) |
Returns the integer value characterizing the given type of the pair of functions: TT = 1, CT = 2, CC = 3.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::get_number_of_functions_in_shell | ( | class(atomic_orbital_basis_obj) | this, |
| integer, intent(in) | i ) |
Returns the number of functions in a shell with the given index.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::get_shell_data | ( | class(atomic_orbital_basis_obj) | this, |
| integer, intent(in) | i, | ||
| class(shell_data_obj), intent(out) | shell_data ) |
Returns the shell data for the i-th shell in the basis set.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::get_shell_name | ( | class(atomic_orbital_basis_obj) | this, |
| integer, intent(in) | i ) |
Returns the name of the i-th shell in the basis set.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::get_symmetry_flags | ( | class(atomic_orbital_basis_obj) | this, |
| integer, intent(in) | irr, | ||
| logical, dimension(:), allocatable | list, | ||
| integer, intent(in) | from, | ||
| integer, intent(in) | to ) |
For a given IRR it returns a list of logical values of size equal to the number of functions in the basis. i-th element of the output array is set to .true. if the i-th function in the basis has symmetry irr andcoms from a shell between indices from and to.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::init | ( | class(atomic_orbital_basis_obj) | this, |
| integer, intent(in) | n, | ||
| class(geometry_obj), intent(in) | geometry ) |
Allocates space for a given number and type of shells of basis functions.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::integral_index | ( | class(atomic_orbital_basis_obj) | this, |
| character(len=*), intent(in) | integral_type, | ||
| integer, dimension(:,:), intent(in) | bf_indices, | ||
| logical, intent(in) | two_p_continuum ) |
Calculates index for the given type of integral. The integral type is given by one of the strings definied in the const module and by the number of basis functions on input. The two_p_continuum input variable distinguishes the two indexing methods used.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::is_initialized | ( | class(atomic_orbital_basis_obj) | this | ) |
Returns the value of initialized.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::normalize_continuum | ( | class(atomic_orbital_basis_obj) | this, |
| real(kind=cfp), intent(in) | a ) |
Normalizes the atomic continuum functions to the R-matrix sphere with a given radius. The normalization is only done if a > 0.0_cfp.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::one_electron_integrals | ( | class(atomic_orbital_basis_obj) | this, |
| class(integral_storage_obj), intent(inout) | integral_storage, | ||
| class(integral_options_obj), intent(in) | integral_options ) |
Calculates and stores 1-electron integrals for all pairs of shells in the basis.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::print | ( | class(atomic_orbital_basis_obj) | this | ) |
Prints the basis set data to stdout.
| procedure atomic_basis_gbl::atomic_orbital_basis_obj::two_electron_integrals | ( | class(atomic_orbital_basis_obj) | this, |
| class(integral_storage_obj), intent(inout) | integral_storage, | ||
| class(integral_options_obj), intent(in) | integral_options ) |
Calculates and stores 2-electron integrals for all pairs of shells in the basis.
| real(kind=cfp), dimension(:,:), allocatable atomic_basis_gbl::atomic_orbital_basis_obj::ao_basis_at_r |
Auxiliary arrays used by eval_orbital.
| integer, dimension(:,:), allocatable atomic_basis_gbl::atomic_orbital_basis_obj::indices_to_shells |
Array of dimensions (2,number_of_basis_functions) containing mapping of basis function indices to shell indices (1) and indices within the shells (2). This table is useful for looking up which shells need to be involved in calculation of a given integral and for picking out the requested integral from the integral batch generated over shells of functions.
| integer atomic_basis_gbl::atomic_orbital_basis_obj::last_ct_fn = -1 |
| integer atomic_basis_gbl::atomic_orbital_basis_obj::last_ct_sh = -1 |
| integer atomic_basis_gbl::atomic_orbital_basis_obj::n_cgto_functions = -1 |
Total number of CGTO functions in the basis.
| integer atomic_basis_gbl::atomic_orbital_basis_obj::n_cont_fns = -1 |
| integer atomic_basis_gbl::atomic_orbital_basis_obj::n_cont_sh = -1 |
| integer atomic_basis_gbl::atomic_orbital_basis_obj::n_prec_ints = -1 |
| integer atomic_basis_gbl::atomic_orbital_basis_obj::n_prec_sh = -1 |
| integer atomic_basis_gbl::atomic_orbital_basis_obj::n_target_fns = -1 |
Number of target functions, index of the last CT pair of functions and other values needed to index optimally 2-electron 1p continuum integrals.
| integer atomic_basis_gbl::atomic_orbital_basis_obj::n_target_sh = -1 |
Number of target shells, index of the last CT pair of shells and other values needed to index optimally 2-electron 1p continuum integrals in case of multiple MPI tasks are used.
| integer atomic_basis_gbl::atomic_orbital_basis_obj::n_tt_pairs = -1 |
| integer atomic_basis_gbl::atomic_orbital_basis_obj::n_tt_sh_pairs = -1 |
| integer, dimension(:,:), allocatable atomic_basis_gbl::atomic_orbital_basis_obj::ordered_shell_pairs |
Indices mapping pairs of shells to a single index. This mapping is required to index properly 2p integrals with only 1p in the continuum in case of multiple MPI tasks being used.
| real(kind=cfp), dimension(:,:), allocatable atomic_basis_gbl::atomic_orbital_basis_obj::r |
| integer, dimension(:,:), allocatable atomic_basis_gbl::atomic_orbital_basis_obj::shell_descriptor |
For each shell (2nd dimension) the first dimension contains the following info: row 1: shell type (CGTO: 1, BTO: 2) row 2: index of the shell within its own type row 3: is a continuum shell (1) or not (0) row 4: starting index for the functions in the shell row 5: the number of functions in the shell row 6: the angular momentum of the shell The shell index corresponds to the order in which it was added to the basis.
| integer, dimension(:,:), allocatable atomic_basis_gbl::atomic_orbital_basis_obj::shell_pair_indices |
shell_pair_indices: Mapping between indices of pairs of shells with the indices of the shells within their own type of shells of basis functions. shell_pair_type: Mapping between pairs of shells and their type: GG (1), BG (2), GB (3), BB (4).
| integer, dimension(:), allocatable atomic_basis_gbl::atomic_orbital_basis_obj::shell_pair_type |