GBTOlib: library for evaluation of molecular integrals in mixed Gaussian / B-spline basis
111
|
This module contains elementary routines which are used for evaluation of GTO-only integrals. More...
Data Types | |
type | boys_function_obj |
This object contains routines for evaluation of the Boys function in double precision. The object must be first initialized using the method init and its argument values imax and mmax. imax corresponds to the largest power in series expansion of the Boys function that is needed to evaluate the Boys function in the expected range of T values. For T in the range 0,60 the value of imax has been determined to be 140. The routine eval then calculates at once the string of Boys functions for a range of m values. More... | |
Functions/Subroutines | |
pure real(kind=cfp) function | dngto (l, alpha) |
real(kind=cfp) function | cart_gto_int (alp, bet, i, j, k) |
This function calculates \( I = \int\int\int x^{i}y^{j}z^{k}\exp[-(\alpha+\beta)r^2] = {\frac{\Gamma(i/2+1/2)\Gamma(j/2+1/2)\Gamma(k/2+1/2)}{((\alpha+\beta)^{1/2(i+j+k+3)})}} \). In particular, this function can be used to calculate overlap of two primitive GTOs centered on the same nucleus. More... | |
real(kind=cfp) function | norm_cart_gto (alp, i, j, k) |
Calculates the normalization factor for the primitive cartesian GTO function. More... | |
real(kind=cfp) function | contr_cart_gto_norm (i, j, k, alp, ccf) |
This function calculates the normalization factor for a contracted cartesian GTO. The cartesian exponents are i,j,k. The contraction coefficients and exponents of the primitives are in the arrays alp(:) and ccf(:). More... | |
real(kind=cfp) function | olap_ccart_csph (i, j, k, l, m, alp, ccf) |
This function calculates overlap of a contracted cartesian GTO with a contracted spherical GTO (both sitting on the same nucleus). The exponents of the contracted cartesian GTO are: i,j,k. The L,M values of the spherical GTO are: l,m. The exponents and contraction coefficients of the primitives are in the arrays alp(:), ccf(:). More... | |
subroutine | cart_cf_sph_cf (l, ix, iy, iz, alp, ccf, sph_cgto_norm, cart_cf, sph_cf) |
This routine performs the conversion of the orbital coefficients (for one shell) for the cartesian contracted GTOs to the basis of contracted spherical GTOs. The input are the data for the cartesian GTOs (l,ix,iy,iz,alp,ccf) and the corresponding orbital coefficients (cart_cf) and the overall norm of the spherical contracted GTO (sph_cgto_norm). It is assumed that the coefficients correspond to normalized contracted cartesian GTOs. The output are the orbital coefficients sph_cf for the corresponding shell of contracted spherical GTOs with the overall norm given by sph_cgto_norm. More... | |
subroutine | sph_cf_cart_cf (l, ix, iy, iz, alp, ccf, sph_cgto_norm, sph_cf, cart_cf) |
This routine performs the conversion of the orbital coefficients (for one shell) for the spherical contracted GTOs to the basis of contracted cartesian GTOs. The input are the data for the cartesian GTOs (l,ix,iy,iz,alp,ccf) and the orbital coefficients for the spherical GTOs (sph_cf). It is assumed that the coefficients correspond to contracted spherical GTOs normalized with sph_cgto_norm. The output are the orbital coefficients cart_cf for the corresponding shell of normalized contracted cartesian GTOs. More... | |
real(kind=cfp) function | cms_gto_norm (rmat_radius, l, np, alp, ccf, cnorm, norms) |
Computes the normalization factor for a CMS-centered contracted spherical GTO and radius given by rmat_radius. This is computed from the self-overlap of the GTO inside the sphere of radius rmat_radius. More... | |
real(kind=cfp) function | CGTO_amplitude (r, l, number_of_primitives, norm, norms, contractions, exponents) |
Calculates the reduced boundary amplitude of a CGTO shell. More... | |
subroutine | olap_kei_tail (l, np_a, np_b, alp_a, alp_b, ccf_a, ccf_b, cnorm_a, norms_a, cnorm_b, norms_b, rmat_radius, olap_tail, kei_tail, bloch_ab) |
Overlap and kinetic energy tail integrals for the normalized contracted spherical GTOs centered on the CMS. These integrals are given by: More... | |
subroutine | prop_cms_tail (l_c, l_a, l_b, np_a, np_b, alp_a, alp_b, ccf_a, ccf_b, cnorm_a, norms_a, cnorm_b, norms_b, rmat_radius, prop_tail) |
Multipole property tail integrals for a pair of shells of contracted spherical GTOs centered on the CMS. We assume that the multipole is centered on CMS just like the two GTO shells. More... | |
subroutine | nari_tail (xc, yc, zc, l_a, l_b, np_a, np_b, alp_a, alp_b, ccf_a, ccf_b, cnorm_a, norms_a, cnorm_b, norms_b, rmat_radius, na_tail) |
Calculates nuclear attraction tail integrals for a pair of shells of spherical CGTOs centered on the CMS. More... | |
subroutine | eri_tail (tgt_prop, tgt_pair, l_a_tgt, l_b_tgt, l_a, l_b, np_a, np_b, alp_a, alp_b, ccf_a, ccf_b, norms_a, norms_b, rmat_radius, swap_ab_cd, eri_tail_int) |
Calculates 2-electron tail integrals for 1 electron in the continuum for a pair of shells of spherical CGTOs centered on the CMS and a pair of shells of target CGTOs. This is the driver routine for eri_tail_shell. This routine ensures that the parameters are passed to eri_tail_shell in such an order to ensure that the tail integrals on output are ordered in the same way as the integrals calculated by eri. This allows for easy subtraction of the tails. l_a_tgt,l_b_tgt are angular momenta in the pair of shells of target CGTOs, tgt_prop and tgt_pair specify the property integrals for the pair of shells of the target CGTOs. The rest of the input parameters are related to the continuum and the R-matrix radius. More... | |
subroutine | eri_tail_shell (tgt_prop, tgt_pair, l_a_tgt, l_b_tgt, l_a, l_b, np_a, np_b, alp_a, alp_b, ccf_a, ccf_b, norms_a, norms_b, rmat_radius, swap_ab_cd, eri_tail_int) |
Calculates 2-electron tail integrals for 1 electron in the continuum for a pair of shells of spherical CGTOs centered on the CMS and a pair of shells of target CGTOs. More... | |
subroutine | abcd_to_cdab (batch_in, batch_t, na, nb, nc, nd) |
This routine takes on input the batch of integrals [ab|cd] with dimensions given by na,nb,nc,nd and returns in batch_t the transposed batch [cd|ab]. Note that this procedure is equivalent to transposition of a matrix M with dimensions na*nb,nc*nd. More... | |
integer function | check_real_array_size (a, d) |
Takes on input the linear array 'a' and the dimension 'd' which the array should AT LEAST have, i.e. not exactly. If the array is smaller it is reallocated (loosing its contents) to the size 'd'. More... | |
subroutine | index_1el (la, lb, ind_a, ind_b, np, int_index) |
Calculates indices of 1-electron integrals output in the linear arrays from the integral routines given the starting indices of functions in the a and b shells. ind_a corresponds to the function with m=-la, same for ind_b. We assume that indices of the rest of the functions in the shell are sequential. The indices for each integral are reordered in such a way that the index of the a function is always .ge. index of the b function. Since the 1-electron integrals are symmetrical there is no problem. The number np denotes the number of passive indices which don't enter the ordering, i.e. we think of the integral batch as being a 3d array (2*la+1,2*lb+1,np). The third dimension useded to be used in the ordering of the multipole moment integrals where the first two dimensions run over the GTOs A,B and the third dimension over the components of the multipole moments. The requirement to index all multipole integrals separately has been dropped so this option is effectively not being used in the whole library. More... | |
subroutine | index_2el (la, lb, lc, ld, ind_a, ind_b, ind_c, ind_d, int_index, keep_ab_cd_order, swap_ab_cd) |
Calculates indices corresponding to the 2-electron integrals (ab|cd) in the linear array sph_ints (calculated by the routine eri_shell) given the starting indices of the functions in the a,b,c,d shells. ind_a must correspond to the function with m=-la, same for the b,c,d shells. We assume that indices of the rest of the functions in each shell are sequential. The indices of the functions are ordered, exploiting the permutational symmetry of the (ab|cd) integrals, so that: If keep_ab_cd_order .eq. .false. then the input variable swap_ab_cd is ignored and: ind_a .ge. ind_b, ind_c .ge. ind_d, ind_a+ind_b .ge. ind_c+ind_d, ind_a .ge. ind_c. This order of the functions is needed for the integral indexing scheme. If keep_ab_cd_order .eq. .true.: indices of the pairs of shells (ab|, |cd) are swapped or not depending on the value of the logical variable swap_ab_cd. Within the (ab| and |cd) shells the indices are ordered so that a.ge.b, c.ge.d. More... | |
subroutine | index_2el_drv (la, lb, lc, ld, ind_a, ind_b, ind_c, ind_d, ind, keep_ab_cd_order) |
Orders the shells according to the angular momentum in the same way it is done in eri and then computes the indices for each quartet of (ab|cd) integrals. See index_2el for description of the ordering option given by the input variable keep_ab_cd_order. More... | |
subroutine | find_mapping (ind_orig, n, n_map, map) |
Determine the mapping of indices ind_a,ind_b,ind_c,ind_d in the array ind_orig so that ind_ap.ge.ind_bp, ind_cp.ge.ind_dp, ind_ap.ge.ind_cp. More... | |
subroutine | reorder_and_index_2el (la, lb, lc, ld, ind_a, ind_b, ind_c, ind_d, column, ind, ints) |
Reorder the integrals ints in columns a,b,c,d so that ap.ge.bp,cp.ge.dp,ap.ge.cp. The input values are the angular momenta and the starting indices in columns 1,2,3,4 of the ints and ind arrays. On output the indices in ind(1:4,i) are not ordered from the largest to the smallest as the original ones (a.ge.b,c.ge.d,a.ge.c) but that's OK since they're not needed to save the integrals in the MPI mode when this routine is used. More... | |
pure integer function | index_1p_continuum (ordered_pairs, ind1, ind2, ind3, ind4, is_CCTT, last_CT, n_prec, n_TT_pairs) |
This is an index for two-electron integral keeping only 1p in the continuum. It works also for 2p in the continuum in which case is_CCTT must be set to .false. More... | |
subroutine | normalize_cgto (number_of_primitives, l, exponents, contractions, norms, norm) |
integer function | check_cgto_data (number_of_primitives, l, exponents, contractions, norms, number_of_functions) |
subroutine | read_cgto (number_of_primitives, l, exponents, contractions, norms, norm, center, non_zero_at_boundary, number_of_functions, lunit, posit, pos_after_rw) |
subroutine | write_cgto (number_of_primitives, l, exponents, contractions, norms, norm, center, non_zero_at_boundary, number_of_functions, lunit, posit, pos_after_rw) |
subroutine | print_cgto_data (number_of_primitives, l, exponents, contractions, norms, norm, center, non_zero_at_boundary) |
subroutine | eval_cgto (r, n_points, number_of_primitives, l, exponents, contractions, norms, norm, center, eval_CGTO_shell) |
subroutine | compare_print_1el_ints (tag, integrals_1, int_index_1, integrals_2, int_index_2, n_integrals, column) |
Only for debugging. More... | |
subroutine | compare_print_2el_ints (tag, integrals_1, int_index_1, integrals_2, int_index_2, n_integrals, column) |
Only for debugging. More... | |
This module contains elementary routines which are used for evaluation of GTO-only integrals.
subroutine gto_routines_gbl::abcd_to_cdab | ( | real(kind=cfp), dimension(*), intent(in) | batch_in, |
real(kind=cfp), dimension(*), intent(out) | batch_t, | ||
integer, intent(in) | na, | ||
integer, intent(in) | nb, | ||
integer, intent(in) | nc, | ||
integer, intent(in) | nd | ||
) |
This routine takes on input the batch of integrals [ab|cd] with dimensions given by na,nb,nc,nd and returns in batch_t the transposed batch [cd|ab]. Note that this procedure is equivalent to transposition of a matrix M with dimensions na*nb,nc*nd.
subroutine gto_routines_gbl::cart_cf_sph_cf | ( | integer, intent(in) | l, |
integer, dimension(:), intent(in) | ix, | ||
integer, dimension(:), intent(in) | iy, | ||
integer, dimension(:), intent(in) | iz, | ||
real(kind=cfp), dimension(:), intent(in) | alp, | ||
real(kind=cfp), dimension(:), intent(in) | ccf, | ||
real(kind=cfp), intent(in) | sph_cgto_norm, | ||
real(kind=cfp), dimension(:), intent(in) | cart_cf, | ||
real(kind=cfp), dimension(:), intent(out) | sph_cf | ||
) |
This routine performs the conversion of the orbital coefficients (for one shell) for the cartesian contracted GTOs to the basis of contracted spherical GTOs. The input are the data for the cartesian GTOs (l,ix,iy,iz,alp,ccf) and the corresponding orbital coefficients (cart_cf) and the overall norm of the spherical contracted GTO (sph_cgto_norm). It is assumed that the coefficients correspond to normalized contracted cartesian GTOs. The output are the orbital coefficients sph_cf for the corresponding shell of contracted spherical GTOs with the overall norm given by sph_cgto_norm.
real(kind=cfp) function gto_routines_gbl::cart_gto_int | ( | real(kind=cfp), intent(in) | alp, |
real(kind=cfp), intent(in) | bet, | ||
integer, intent(in) | i, | ||
integer, intent(in) | j, | ||
integer, intent(in) | k | ||
) |
This function calculates \( I = \int\int\int x^{i}y^{j}z^{k}\exp[-(\alpha+\beta)r^2] = {\frac{\Gamma(i/2+1/2)\Gamma(j/2+1/2)\Gamma(k/2+1/2)}{((\alpha+\beta)^{1/2(i+j+k+3)})}} \). In particular, this function can be used to calculate overlap of two primitive GTOs centered on the same nucleus.
real(kind=cfp) function gto_routines_gbl::CGTO_amplitude | ( | real(kind=cfp), intent(in) | r, |
integer, intent(in) | l, | ||
integer, intent(in) | number_of_primitives, | ||
real(kind=cfp), intent(in) | norm, | ||
real(kind=cfp), dimension(number_of_primitives), intent(in) | norms, | ||
real(kind=cfp), dimension(number_of_primitives), intent(in) | contractions, | ||
real(kind=cfp), dimension(number_of_primitives), intent(in) | exponents | ||
) |
Calculates the reduced boundary amplitude of a CGTO shell.
integer function gto_routines_gbl::check_cgto_data | ( | integer, intent(in) | number_of_primitives, |
integer, intent(in) | l, | ||
real(kind=cfp), dimension(:), allocatable | exponents, | ||
real(kind=cfp), dimension(:), allocatable | contractions, | ||
real(kind=cfp), dimension(:), allocatable | norms, | ||
integer, intent(in) | number_of_functions | ||
) |
integer function gto_routines_gbl::check_real_array_size | ( | real(kind=cfp), dimension(:), intent(inout), allocatable | a, |
integer, intent(in) | d | ||
) |
Takes on input the linear array 'a' and the dimension 'd' which the array should AT LEAST have, i.e. not exactly. If the array is smaller it is reallocated (loosing its contents) to the size 'd'.
real(kind=cfp) function gto_routines_gbl::cms_gto_norm | ( | real(kind=cfp), intent(in) | rmat_radius, |
integer, intent(in) | l, | ||
integer, intent(in) | np, | ||
real(kind=cfp), dimension(:), intent(in) | alp, | ||
real(kind=cfp), dimension(:), intent(in) | ccf, | ||
real(kind=cfp), intent(in) | cnorm, | ||
real(kind=cfp), dimension(:), intent(in) | norms | ||
) |
Computes the normalization factor for a CMS-centered contracted spherical GTO and radius given by rmat_radius. This is computed from the self-overlap of the GTO inside the sphere of radius rmat_radius.
subroutine gto_routines_gbl::compare_print_1el_ints | ( | character(len=4), intent(in) | tag, |
real(kind=cfp), dimension(:,:), allocatable | integrals_1, | ||
integer, dimension(:,:), allocatable | int_index_1, | ||
real(kind=cfp), dimension(:,:), allocatable | integrals_2, | ||
integer, dimension(:,:), allocatable | int_index_2, | ||
integer, intent(in) | n_integrals, | ||
integer, intent(in) | column | ||
) |
Only for debugging.
subroutine gto_routines_gbl::compare_print_2el_ints | ( | character(len=4), intent(in) | tag, |
real(kind=cfp), dimension(:,:), allocatable | integrals_1, | ||
integer, dimension(:,:), allocatable | int_index_1, | ||
real(kind=cfp), dimension(:,:), allocatable | integrals_2, | ||
integer, dimension(:,:), allocatable | int_index_2, | ||
integer, intent(in) | n_integrals, | ||
integer, intent(in) | column | ||
) |
Only for debugging.
real(kind=cfp) function gto_routines_gbl::contr_cart_gto_norm | ( | integer, intent(in) | i, |
integer, intent(in) | j, | ||
integer, intent(in) | k, | ||
real(kind=cfp), dimension(:), intent(in) | alp, | ||
real(kind=cfp), dimension(:), intent(in) | ccf | ||
) |
This function calculates the normalization factor for a contracted cartesian GTO. The cartesian exponents are i,j,k. The contraction coefficients and exponents of the primitives are in the arrays alp(:) and ccf(:).
pure real(kind=cfp) function gto_routines_gbl::dngto | ( | integer, intent(in) | l, |
real(kind=cfp), intent(in) | alpha | ||
) |
Calculate the normalization factor for a GTO orbital.
\[ N^{GTO}_{\alpha l} = \sqrt{\frac{2(2\alpha)^{l+3/2}}{\Gamma(l+3/2)}}\sqrt{\frac{2l+1}{4\pi}}. \]
Gamma is the gamma function.
[in] | l,alpha | Integer angular momentum l of the GTO and its exponent alpha (real). |
subroutine gto_routines_gbl::eri_tail | ( | real(kind=cfp), dimension(:,:), intent(in) | tgt_prop, |
integer, intent(in) | tgt_pair, | ||
integer, intent(in) | l_a_tgt, | ||
integer, intent(in) | l_b_tgt, | ||
integer, intent(in) | l_a, | ||
integer, intent(in) | l_b, | ||
integer, intent(in) | np_a, | ||
integer, intent(in) | np_b, | ||
real(kind=cfp), dimension(:), intent(in) | alp_a, | ||
real(kind=cfp), dimension(:), intent(in) | alp_b, | ||
real(kind=cfp), dimension(:), intent(in) | ccf_a, | ||
real(kind=cfp), dimension(:), intent(in) | ccf_b, | ||
real(kind=cfp), dimension(:), intent(in) | norms_a, | ||
real(kind=cfp), dimension(:), intent(in) | norms_b, | ||
real(kind=cfp), intent(in) | rmat_radius, | ||
logical, intent(in) | swap_ab_cd, | ||
real(kind=cfp), dimension(:), allocatable | eri_tail_int | ||
) |
Calculates 2-electron tail integrals for 1 electron in the continuum for a pair of shells of spherical CGTOs centered on the CMS and a pair of shells of target CGTOs. This is the driver routine for eri_tail_shell. This routine ensures that the parameters are passed to eri_tail_shell in such an order to ensure that the tail integrals on output are ordered in the same way as the integrals calculated by eri. This allows for easy subtraction of the tails. l_a_tgt,l_b_tgt are angular momenta in the pair of shells of target CGTOs, tgt_prop and tgt_pair specify the property integrals for the pair of shells of the target CGTOs. The rest of the input parameters are related to the continuum and the R-matrix radius.
subroutine gto_routines_gbl::eri_tail_shell | ( | real(kind=cfp), dimension(:,:), intent(in) | tgt_prop, |
integer, intent(in) | tgt_pair, | ||
integer, intent(in) | l_a_tgt, | ||
integer, intent(in) | l_b_tgt, | ||
integer, intent(in) | l_a, | ||
integer, intent(in) | l_b, | ||
integer, intent(in) | np_a, | ||
integer, intent(in) | np_b, | ||
real(kind=cfp), dimension(:), intent(in) | alp_a, | ||
real(kind=cfp), dimension(:), intent(in) | alp_b, | ||
real(kind=cfp), dimension(:), intent(in) | ccf_a, | ||
real(kind=cfp), dimension(:), intent(in) | ccf_b, | ||
real(kind=cfp), dimension(:), intent(in) | norms_a, | ||
real(kind=cfp), dimension(:), intent(in) | norms_b, | ||
real(kind=cfp), intent(in) | rmat_radius, | ||
logical, intent(in) | swap_ab_cd, | ||
real(kind=cfp), dimension(:), allocatable | eri_tail_int | ||
) |
Calculates 2-electron tail integrals for 1 electron in the continuum for a pair of shells of spherical CGTOs centered on the CMS and a pair of shells of target CGTOs.
\[ ERI_{tail} = (-1)^{m_{a}+m_{b}}\frac{4\pi}{\sqrt{(2l_{a}+1)(2l_{b}+1)}}N_{a}N_{b}\sum_{i=1}^{n_{a}}\sum_{j=1}^{n_{j}}c_{i}c_{j}np_{i}np_{j}\sum_{l=\vert l_a-l_b\vert}^{l_a+l_b} \sqrt{\frac{4\pi}{2l+1}}\frac{1}{2}(\alpha_{i}+\beta_{j})^{\frac{l_a+l_b-l+2}{2}}\Gamma[\frac{l_a+l_b-l+2}{2},a^{2}(\alpha_{i}+\beta_{j})] \sum_{m=-l}^{l}\langle l_a,m_a\vert l,m\vert l_b,m_b\rangle_{R}P_{A,B}^{l,m}. \]
This routine calculates the tail integrals in the order (A,B,C,D) where A,B are the continuum shells. The variable swap_ab_cd controls whether the integrals on output should be ordered as (C,D,A,B). The array tgt_prop is the array of the properties calculated by sph_mult_mom at least for L_max = l_a+l_b for the pair of shells (ab) of the target GTOs. However, here we assume that the values in tgt_prop have been calculated for a number of pairs of shells of GTOs and that these are saved in columns of tgt_prop, where the column corresponding to the current combination of target GTOs is given by the index tgt_pair.
subroutine gto_routines_gbl::eval_cgto | ( | real(kind=cfp), dimension(1:3,n_points), intent(in) | r, |
integer, intent(in) | n_points, | ||
integer, intent(in) | number_of_primitives, | ||
integer, intent(in) | l, | ||
real(kind=cfp), dimension(number_of_primitives), intent(in) | exponents, | ||
real(kind=cfp), dimension(number_of_primitives), intent(in) | contractions, | ||
real(kind=cfp), dimension(number_of_primitives), intent(in) | norms, | ||
real(kind=cfp), intent(in) | norm, | ||
real(kind=cfp), dimension(3), intent(in) | center, | ||
real(kind=cfp), dimension(2*l+1,n_points) | eval_CGTO_shell | ||
) |
subroutine gto_routines_gbl::find_mapping | ( | integer, dimension(4), intent(inout) | ind_orig, |
integer, dimension(4), intent(in) | n, | ||
integer, dimension(3), intent(out) | n_map, | ||
integer, dimension(4), intent(out) | map | ||
) |
Determine the mapping of indices ind_a,ind_b,ind_c,ind_d in the array ind_orig so that ind_ap.ge.ind_bp, ind_cp.ge.ind_dp, ind_ap.ge.ind_cp.
subroutine gto_routines_gbl::index_1el | ( | integer, intent(in) | la, |
integer, intent(in) | lb, | ||
integer, intent(in) | ind_a, | ||
integer, intent(in) | ind_b, | ||
integer, intent(in) | np, | ||
integer, dimension(:,:), intent(out) | int_index | ||
) |
Calculates indices of 1-electron integrals output in the linear arrays from the integral routines given the starting indices of functions in the a and b shells. ind_a corresponds to the function with m=-la, same for ind_b. We assume that indices of the rest of the functions in the shell are sequential. The indices for each integral are reordered in such a way that the index of the a function is always .ge. index of the b function. Since the 1-electron integrals are symmetrical there is no problem. The number np denotes the number of passive indices which don't enter the ordering, i.e. we think of the integral batch as being a 3d array (2*la+1,2*lb+1,np). The third dimension useded to be used in the ordering of the multipole moment integrals where the first two dimensions run over the GTOs A,B and the third dimension over the components of the multipole moments. The requirement to index all multipole integrals separately has been dropped so this option is effectively not being used in the whole library.
pure integer function gto_routines_gbl::index_1p_continuum | ( | integer, dimension(:,:), intent(in) | ordered_pairs, |
integer, intent(in) | ind1, | ||
integer, intent(in) | ind2, | ||
integer, intent(in) | ind3, | ||
integer, intent(in) | ind4, | ||
logical, intent(in) | is_CCTT, | ||
integer, intent(in) | last_CT, | ||
integer, intent(in) | n_prec, | ||
integer, intent(in) | n_TT_pairs | ||
) |
This is an index for two-electron integral keeping only 1p in the continuum. It works also for 2p in the continuum in which case is_CCTT must be set to .false.
subroutine gto_routines_gbl::index_2el | ( | integer, intent(in) | la, |
integer, intent(in) | lb, | ||
integer, intent(in) | lc, | ||
integer, intent(in) | ld, | ||
integer, intent(in) | ind_a, | ||
integer, intent(in) | ind_b, | ||
integer, intent(in) | ind_c, | ||
integer, intent(in) | ind_d, | ||
integer, dimension(:,:), allocatable | int_index, | ||
logical, intent(in) | keep_ab_cd_order, | ||
logical, intent(in) | swap_ab_cd | ||
) |
Calculates indices corresponding to the 2-electron integrals (ab|cd) in the linear array sph_ints (calculated by the routine eri_shell) given the starting indices of the functions in the a,b,c,d shells. ind_a must correspond to the function with m=-la, same for the b,c,d shells. We assume that indices of the rest of the functions in each shell are sequential. The indices of the functions are ordered, exploiting the permutational symmetry of the (ab|cd) integrals, so that: If keep_ab_cd_order .eq. .false. then the input variable swap_ab_cd is ignored and: ind_a .ge. ind_b, ind_c .ge. ind_d, ind_a+ind_b .ge. ind_c+ind_d, ind_a .ge. ind_c. This order of the functions is needed for the integral indexing scheme. If keep_ab_cd_order .eq. .true.: indices of the pairs of shells (ab|, |cd) are swapped or not depending on the value of the logical variable swap_ab_cd. Within the (ab| and |cd) shells the indices are ordered so that a.ge.b, c.ge.d.
subroutine gto_routines_gbl::index_2el_drv | ( | integer, intent(in) | la, |
integer, intent(in) | lb, | ||
integer, intent(in) | lc, | ||
integer, intent(in) | ld, | ||
integer, intent(in) | ind_a, | ||
integer, intent(in) | ind_b, | ||
integer, intent(in) | ind_c, | ||
integer, intent(in) | ind_d, | ||
integer, dimension(:,:), allocatable | ind, | ||
logical, intent(in) | keep_ab_cd_order | ||
) |
Orders the shells according to the angular momentum in the same way it is done in eri and then computes the indices for each quartet of (ab|cd) integrals. See index_2el for description of the ordering option given by the input variable keep_ab_cd_order.
subroutine gto_routines_gbl::nari_tail | ( | real(kind=cfp), intent(in) | xc, |
real(kind=cfp), intent(in) | yc, | ||
real(kind=cfp), intent(in) | zc, | ||
integer, intent(in) | l_a, | ||
integer, intent(in) | l_b, | ||
integer, intent(in) | np_a, | ||
integer, intent(in) | np_b, | ||
real(kind=cfp), dimension(:), intent(in) | alp_a, | ||
real(kind=cfp), dimension(:), intent(in) | alp_b, | ||
real(kind=cfp), dimension(:), intent(in) | ccf_a, | ||
real(kind=cfp), dimension(:), intent(in) | ccf_b, | ||
real(kind=cfp), intent(in) | cnorm_a, | ||
real(kind=cfp), dimension(:), intent(in) | norms_a, | ||
real(kind=cfp), intent(in) | cnorm_b, | ||
real(kind=cfp), dimension(:), intent(in) | norms_b, | ||
real(kind=cfp), intent(in) | rmat_radius, | ||
real(kind=cfp), dimension(:), intent(out) | na_tail | ||
) |
Calculates nuclear attraction tail integrals for a pair of shells of spherical CGTOs centered on the CMS.
\[ NAI_{tail} = (-1)^{m_{a}+m_{b}}\frac{16\pi^2}{\sqrt{(2l_{a}+1)(2l_{b}+1)}}N_{a}N_{b}\sum_{i=1}^{n_{a}}\sum_{j=1}^{n_{j}}c_{i}c_{j}np_{i}np_{j}\sum_{l=\vert l_a-l_b\vert}^{l_a+l_b} \frac{R_{c}^{l}}{2l+1}\frac{1}{2}(\alpha_{i}+\beta_{j})^{\frac{l_a+l_b-l+2}{2}}\Gamma[\frac{l_a+l_b-l+2}{2},a^{2}(\alpha_{i}+\beta_{j})] \sum_{m=-l}^{l}\langle l_a,m_a\vert l,m\vert l_b,m_b\rangle_{R}X_{l,m}(\mathbf{R_{c}}). \]
real(kind=cfp) function gto_routines_gbl::norm_cart_gto | ( | real(kind=cfp), intent(in) | alp, |
integer, intent(in) | i, | ||
integer, intent(in) | j, | ||
integer, intent(in) | k | ||
) |
Calculates the normalization factor for the primitive cartesian GTO function.
subroutine gto_routines_gbl::normalize_cgto | ( | integer | number_of_primitives, |
integer | l, | ||
real(kind=cfp), dimension(number_of_primitives), intent(in) | exponents, | ||
real(kind=cfp), dimension(number_of_primitives), intent(in) | contractions, | ||
real(kind=cfp), dimension(number_of_primitives), intent(out) | norms, | ||
real(kind=cfp), intent(out) | norm | ||
) |
real(kind=cfp) function gto_routines_gbl::olap_ccart_csph | ( | integer, intent(in) | i, |
integer, intent(in) | j, | ||
integer, intent(in) | k, | ||
integer, intent(in) | l, | ||
integer, intent(in) | m, | ||
real(kind=cfp), dimension(:), intent(in) | alp, | ||
real(kind=cfp), dimension(:), intent(in) | ccf | ||
) |
This function calculates overlap of a contracted cartesian GTO with a contracted spherical GTO (both sitting on the same nucleus). The exponents of the contracted cartesian GTO are: i,j,k. The L,M values of the spherical GTO are: l,m. The exponents and contraction coefficients of the primitives are in the arrays alp(:), ccf(:).
subroutine gto_routines_gbl::olap_kei_tail | ( | integer, intent(in) | l, |
integer, intent(in) | np_a, | ||
integer, intent(in) | np_b, | ||
real(kind=cfp), dimension(:), intent(in) | alp_a, | ||
real(kind=cfp), dimension(:), intent(in) | alp_b, | ||
real(kind=cfp), dimension(:), intent(in) | ccf_a, | ||
real(kind=cfp), dimension(:), intent(in) | ccf_b, | ||
real(kind=cfp), intent(in) | cnorm_a, | ||
real(kind=cfp), dimension(:), intent(in) | norms_a, | ||
real(kind=cfp), intent(in) | cnorm_b, | ||
real(kind=cfp), dimension(:), intent(in) | norms_b, | ||
real(kind=cfp), intent(in) | rmat_radius, | ||
real(kind=cfp), intent(out) | olap_tail, | ||
real(kind=cfp), intent(out) | kei_tail, | ||
real(kind=cfp), intent(out) | bloch_ab | ||
) |
Overlap and kinetic energy tail integrals for the normalized contracted spherical GTOs centered on the CMS. These integrals are given by:
\[ S_{tail} = \langle S_{LM}(\mathbf{r})\exp[-\alpha r^2] \vert S_{L'M'}(\mathbf{r})\exp[-\beta r^2] \rangle = \delta_{L'L}\delta_{M'M}\frac{4\pi}{2L+1}\frac{1}{2}(\alpha+\beta)^{-L-3/2}\Gamma[3/2+L, a^2(\alpha+\beta)], \]
\[ K_{tail} = \langle S_{LM}(\mathbf{r})\exp[-\alpha r^2] \vert\Delta\vert S_{L'M'}(\mathbf{r})\exp[-\beta r^2] \rangle = \delta_{L'L}\delta_{M'M}\frac{4\pi}{2L+1}\beta(\alpha+\beta)^{-\frac{5}{2}-L} \left(-(\alpha+\beta)(3+2 L)\Gamma\left[\frac{3}{2}+L,a^2 (\alpha+\beta)\right]+2 \beta\Gamma\left[\frac{5}{2}+L,a^2 (\alpha+\beta)\right]\right). \]
The delta terms are actually not taken into account in this routine - hence it is responsibility of the calling routine to subtract the calculated integrals only when appropriate. We assume that the mass of the particle is that of electron. The formula is for a pair of unnormalized primitives. The actual routine performs contraction over the primitives and multiplies by all the norms. The Bloch terms for the (ab) combination of contracted GTOs are also calculated. These should be added to the KE integrals following the tail subtraction.
subroutine gto_routines_gbl::print_cgto_data | ( | integer, intent(in) | number_of_primitives, |
integer, intent(in) | l, | ||
real(kind=cfp), dimension(number_of_primitives), intent(in) | exponents, | ||
real(kind=cfp), dimension(number_of_primitives), intent(in) | contractions, | ||
real(kind=cfp), dimension(number_of_primitives), intent(in) | norms, | ||
real(kind=cfp), intent(in) | norm, | ||
real(kind=cfp), dimension(3), intent(in) | center, | ||
logical, intent(in) | non_zero_at_boundary | ||
) |
subroutine gto_routines_gbl::prop_cms_tail | ( | integer, intent(in) | l_c, |
integer, intent(in) | l_a, | ||
integer, intent(in) | l_b, | ||
integer, intent(in) | np_a, | ||
integer, intent(in) | np_b, | ||
real(kind=cfp), dimension(:), intent(in) | alp_a, | ||
real(kind=cfp), dimension(:), intent(in) | alp_b, | ||
real(kind=cfp), dimension(:), intent(in) | ccf_a, | ||
real(kind=cfp), dimension(:), intent(in) | ccf_b, | ||
real(kind=cfp), intent(in) | cnorm_a, | ||
real(kind=cfp), dimension(:), intent(in) | norms_a, | ||
real(kind=cfp), intent(in) | cnorm_b, | ||
real(kind=cfp), dimension(:), intent(in) | norms_b, | ||
real(kind=cfp), intent(in) | rmat_radius, | ||
real(kind=cfp), dimension(:), intent(out) | prop_tail | ||
) |
Multipole property tail integrals for a pair of shells of contracted spherical GTOs centered on the CMS. We assume that the multipole is centered on CMS just like the two GTO shells.
\[ M_{tail} = \langle N_{a} S_{La,Ma}(\mathbf{r})\exp[-\alpha r^2]\vert S_{Lc,Mc} \vert N_{b} S_{Lb,Mb}(\mathbf{r})\exp[-\beta r^2] \rangle = = \frac{4\pi^{3/2}}{\sqrt{(2La+1)(2Lb+1)(2Lc+1)}}N_{a}N_{b}\frac{1}{2}(\alpha+\beta)^{-(La+Lb+Lc+3)/2}\Gamma[(La+Lb+Lc+3)/2,a^2(\alpha+\beta)], \]
where the first integration is only over $ r \ge a $. The formula is for a pair of normalized primitives. The actual routine performs contraction over the primitives. The calculated integrals are stored in the output 1D array which emulates the 3D array with dimensions (2*la+1,2*lb+1,(l_c+1)**2) where la = max(l_a,l_b), lb = min(l_a,l_b). Hence the tail integral for the M angular numbers ma,mb,m is stored in (ma+la+1,mb+lb+1,l*l+m+l+1), i.e. in the element of the 1D array with index ma+la+1 + (2*la+1)*(mb+lb) + (2*la+1)*(2*lb+1)*(l*l+m+l). Note that the order of the tail integrals as output by this routine is the same as the one in sph_mult_mom so subtraction of the tails is straightforward: sph_mult_mom(:) - prop_tail(:).
subroutine gto_routines_gbl::read_cgto | ( | integer, intent(out) | number_of_primitives, |
integer, intent(out) | l, | ||
real(kind=cfp), dimension(:), allocatable | exponents, | ||
real(kind=cfp), dimension(:), allocatable | contractions, | ||
real(kind=cfp), dimension(:), allocatable | norms, | ||
real(kind=cfp), intent(out) | norm, | ||
real(kind=cfp), dimension(3), intent(out) | center, | ||
logical, intent(out) | non_zero_at_boundary, | ||
integer, intent(out) | number_of_functions, | ||
integer, intent(in) | lunit, | ||
integer, intent(in) | posit, | ||
integer, intent(out) | pos_after_rw | ||
) |
subroutine gto_routines_gbl::reorder_and_index_2el | ( | integer, intent(in) | la, |
integer, intent(in) | lb, | ||
integer, intent(in) | lc, | ||
integer, intent(in) | ld, | ||
integer, intent(in) | ind_a, | ||
integer, intent(in) | ind_b, | ||
integer, intent(in) | ind_c, | ||
integer, intent(in) | ind_d, | ||
integer, intent(in) | column, | ||
integer, dimension(:,:), allocatable | ind, | ||
real(kind=cfp), dimension(:,:), allocatable | ints | ||
) |
Reorder the integrals ints in columns a,b,c,d so that ap.ge.bp,cp.ge.dp,ap.ge.cp. The input values are the angular momenta and the starting indices in columns 1,2,3,4 of the ints and ind arrays. On output the indices in ind(1:4,i) are not ordered from the largest to the smallest as the original ones (a.ge.b,c.ge.d,a.ge.c) but that's OK since they're not needed to save the integrals in the MPI mode when this routine is used.
subroutine gto_routines_gbl::sph_cf_cart_cf | ( | integer, intent(in) | l, |
integer, dimension(:), intent(in) | ix, | ||
integer, dimension(:), intent(in) | iy, | ||
integer, dimension(:), intent(in) | iz, | ||
real(kind=cfp), dimension(:), intent(in) | alp, | ||
real(kind=cfp), dimension(:), intent(in) | ccf, | ||
real(kind=cfp), intent(in) | sph_cgto_norm, | ||
real(kind=cfp), dimension(:), intent(in) | sph_cf, | ||
real(kind=cfp), dimension(:), intent(out) | cart_cf | ||
) |
This routine performs the conversion of the orbital coefficients (for one shell) for the spherical contracted GTOs to the basis of contracted cartesian GTOs. The input are the data for the cartesian GTOs (l,ix,iy,iz,alp,ccf) and the orbital coefficients for the spherical GTOs (sph_cf). It is assumed that the coefficients correspond to contracted spherical GTOs normalized with sph_cgto_norm. The output are the orbital coefficients cart_cf for the corresponding shell of normalized contracted cartesian GTOs.
subroutine gto_routines_gbl::write_cgto | ( | integer, intent(in) | number_of_primitives, |
integer, intent(in) | l, | ||
real(kind=cfp), dimension(number_of_primitives), intent(in) | exponents, | ||
real(kind=cfp), dimension(number_of_primitives), intent(in) | contractions, | ||
real(kind=cfp), dimension(number_of_primitives), intent(in) | norms, | ||
real(kind=cfp), intent(in) | norm, | ||
real(kind=cfp), dimension(3), intent(in) | center, | ||
logical, intent(in) | non_zero_at_boundary, | ||
integer, intent(in) | number_of_functions, | ||
integer, intent(in) | lunit, | ||
integer, intent(in) | posit, | ||
integer, intent(out) | pos_after_rw | ||
) |