GBTOlib: library for evaluation of molecular integrals in mixed Gaussian / B-spline basis  111
Data Types | Functions/Subroutines
special_functions_gbl Module Reference

Data Types

interface  cfp_9gmic
 
interface  cfp_9gmit
 
interface  cfp_9lgic
 
interface  cfp_9lgit
 
interface  cfp_9lgmc
 
interface  cfp_asyik
 
interface  cfp_asyjy
 
interface  cfp_besi
 
interface  cfp_besj
 
interface  cfp_binom
 
interface  cfp_csevl
 
interface  cfp_eval_poly_horner
 
interface  cfp_gamic
 
interface  cfp_gamlm
 
interface  cfp_gamma_fun
 
interface  cfp_gamma_slatec
 
interface  cfp_initds
 
interface  cfp_jairy
 
interface  cfp_lgams
 
interface  cfp_lngam
 
interface  cfp_lnrel
 

Functions/Subroutines

elemental integer function, public ipair (i)
 This function is used to index an ordered pair of values (i,j), where i .ge. j, i .ge. 1. The index of the ordered pair is: ipair(i) + j. This process of indexing can be used in a nested way to index quartets of integers, say (i,j,k,l). We compute ij = ipair(i)+j and kl = ipair(k)+l. The index of the quartet (i,j,k,l) is: ipair(ij)+kl. We assume that i.ge.j, k.ge.l, ij.ge.kl. This nesting (triangularization) is used heavilly to index the 2-particle symmetric integrals. More...
 
subroutine, public unpack_pq (pq, n, p, q)
 Assuming \( pq = p + (q-1)*n \) this routine returns p,q given pq,n. More...
 
subroutine, public cfp_nlm (INorm, L)
 
subroutine, public cfp_resh (SH, X, Y, Z, L)
 
subroutine, public cfp_solh (SH, x, y, z, L)
 
subroutine, public cfp_zhar (SH, X, Y, Z, L)
 
real(kind=ep1) function, dimension(1:mmax+1), public boys_function_quad (T, mmax)
 Quadruple precision version of boys_function. See boys_function for details on the method of evaluation. More...
 
elemental integer function dfact (n)
 Computes the double factorial \((n)!!\) of the integer number n. n must be .ge. -1. More...
 
real(kind=cfp) function cfp_eval_poly_horner_single (n, x, a)
 This function evaluates a polynomial of degree \(n\) at point \(x\) using the Horner form for polynomials. More...
 
subroutine, public cfp_eval_poly_horner_many (n, x, n_x, a, res)
 This subroutine evaluates a set of polynomials of degree \(n\) at a number of points \(x\) using the Horner form for polynomials. More...
 
subroutine, public cfp_sph_to_cart_mapping (l, m, c, i_exp, j_exp, k_exp)
 This routine constructs for a given real solid harmonic L,M the list of cartesian functions which build up this real solid harmonic. On output the array c(:) contains the coefficients with which the individual cartesians contribute to the solid harmonic. Exponents of the cartesians are then given in the arrays i_exp,j_exp,k_exp. More...
 
subroutine, public cfp_sph_shell_to_cart_shell (l, c)
 This routine constructs for a given real solid harmonic L the matrix of coefficients for all M of the transformation from the cartesian harmonics (~x^i*y^j*z^k, i+j+k=L). On output the linear array c(:) emulates 2D matrix with 2*L+1 rows and (L+1)*(L+2)/2 columns. The rows of c correspond to the M values: -L,-L+1,...,0,1,...,L.
The coefficients in each row are ordered so that the n-th coefficient in the given row of c corresponds to the coefficient for the cartesian harmonic with the shell canonical index n. More...
 
subroutine, public cfp_sph_shell_to_cart_lshells (l, nz, c_nz, nz_can)
 This routine constructs for all real solid harmonic with L .le. l the matrix of coefficients for all M of the transformation from the cartesian harmonics (~x^i*y^j*z^k, i+j+k=L). On output the linear array c(:) emulates 2D matrix with 2*L+1 rows and (L+1)*(L+2)/2 columns. The rows of c correspond to the M values: -L,-L+1,...,0,1,...,L.
The coefficients in each row are ordered so that the n-th coefficient in the given row of c corresponds to the coefficient for the cartesian harmonic with the shell canonical index n. More...
 

Function/Subroutine Documentation

◆ boys_function_quad()

real(kind=ep1) function, dimension(1:mmax+1), public special_functions_gbl::boys_function_quad ( real(kind=ep1), intent(in)  T,
integer, intent(in)  mmax 
)

Quadruple precision version of boys_function. See boys_function for details on the method of evaluation.

Here is the call graph for this function:

◆ cfp_eval_poly_horner_many()

subroutine, public special_functions_gbl::cfp_eval_poly_horner_many ( integer, intent(in)  n,
real(kind=cfp), dimension(n_x), intent(in)  x,
integer, intent(in)  n_x,
real(kind=cfp), dimension(n_x,n+1), intent(in)  a,
real(kind=cfp), dimension(n_x), intent(out)  res 
)

This subroutine evaluates a set of polynomials of degree \(n\) at a number of points \(x\) using the Horner form for polynomials.

Parameters
[in]nOrder of the polynomials.
[in]xArray of points at which to evaluate.
[in]n_xThe number of points (dimension of x).
[in]aThe array of the coefficients a(1:n_x,1:n+1): different for each point x.
[out]resThe polynomials evaluated at points x.

◆ cfp_eval_poly_horner_single()

real(kind=cfp) function special_functions_gbl::cfp_eval_poly_horner_single ( integer, intent(in)  n,
real(kind=cfp), intent(in)  x,
real(kind=cfp), dimension(1:n+1), intent(in)  a 
)

This function evaluates a polynomial of degree \(n\) at point \(x\) using the Horner form for polynomials.

Parameters
[in]nOrder of the polynomial
[in]xThe point at which to evaluate.
[in]aThe array of the coefficients a(1:n+1): \(a_{0},a_{1},\dots,a_{n}\).

◆ cfp_nlm()

subroutine, public special_functions_gbl::cfp_nlm ( real(kind=cfp), dimension(-l:l,0:l), intent(out)  INorm,
integer, intent(in)  L 
)
Purpose:
 Calculate the full set of inverse normalization factors Ilm for the real unnormalized solid spherical harmonics with l,m .le. L using recursion.
 The unnormalized solid spherical harmonics are defined in Rico, Lopez, et al. Int. J. Q. Chem.,2012,DOI:10.1002/qua.24356.

\[ I_{lm} = \frac{1}{N_{lm}} = \sqrt{\frac{2\pi(1+\delta_{m,0})}{2l+1}\frac{(l+|m|)!}{(l-|m|)!}}. \]

Todo:
{add error checking for case l > size(INorm(,:))}
  Nlm is the normalization factor for the real unnormalized solid spherical harmonics zlm.
Parameters
[out]INorm
  Real array (-L:L,0:L) containing the values of the inverse normalization factors Ilm for all l .le. L and the corresponding values of m.
[in]L
  Integer angular momentum of the last Ilm to be generated.

◆ cfp_resh()

subroutine, public special_functions_gbl::cfp_resh ( real(kind=cfp), dimension(-l:l,0:l), intent(out)  SH,
real(kind=cfp), intent(in)  X,
real(kind=cfp), intent(in)  Y,
real(kind=cfp), intent(in)  Z,
integer, intent(in)  L 
)
Purpose:
 Calculate the full set of real normalized spherical harmonics with l,m .le. L using recursion. See Helgaker, p.210, p.215-218 (Section 6.4) for details.
 The real spherical harmonics are related to the solid ones using the formula:

\[ X_{lm}(\Omega) = S_{lm}(x/r,y/r,z/r)\sqrt{\frac{2l+1}{4\pi}}. \]

We define the normalization factor:

\[ n_{lm} = \sqrt{\frac{2l+1}{4\pi}} \]

and then use Helgaker's method to calculate the 'scaled' solid harmonics:

\[ X_{lm}(\Omega) = {\overline{S}_{lm}}(x/r,y/r,z/r)=n_{lm}S_{lm}(x/r,y/r,z/r). \]

 The real spherical harmonics defined in this way are equivalent to the real spherical harmonics defined using the standard formula:

\[ X_{lm}(\Omega) = \begin{cases} {\sqrt{2}}(-1)^{m}{\cal{R}}\left(Y_{l\vert m\vert}(\Omega)\right), m>0 \\ {\sqrt{2}}(-1)^{m}{\cal{I}}\left(Y_{l\vert m\vert}(\Omega)\right), m<0 \\ Y_{l0}(\Omega), m=0, \end{cases} \]

where the complex spherical harmonics $Y_{lm}(\Omega)$ satisfy the Condon-Shortley phase convention. Note that the convention of Homeier and Steinborn is not to include the factor (-1)**m in the definition of the real spherical harmonics. However, here we follow the std. definition and include the factor (-1)**m in the definition of the real spherical harmonics.
Todo:
{add error checking for case L > size(SH(,:))}
Parameters
[out]SH
  Real array (-L:L,0:L) containing the values of the spherical harmonics for all l .le. L and the corresponding values of m.
  According to my tests the results are in vast majority of cases accurate to full REAL(kind=cfp). In only a few cases the precision in the last 1 or 2 digits is lost.
[in]X,Y,Z
  Real input coordinates for which the real harmonics will be evaluated. 
  The radius vector r={X,Y,Z} need not be normalized to 1 since this is done automatically during the calculation.
[in]L
  Integer angular momentum of the last spherical harmonic to be generated. For L = 0; SH = 1/sqrt(4*pi)

◆ cfp_solh()

subroutine, public special_functions_gbl::cfp_solh ( real(kind=cfp), dimension(-l:l,0:l), intent(out)  SH,
real(kind=cfp), intent(in)  x,
real(kind=cfp), intent(in)  y,
real(kind=cfp), intent(in)  z,
integer, intent(in)  L 
)
Purpose:
 Calculate the full set of real normalized solid spherical harmonics with l,m .le. L using recursion. See Helgaker, p.218 (Section 6.4) for details.
 For l=m=0: Slm=1.
Todo:
{add error checking for case L > size(SH(,:))}
Parameters
[out]SH
  Real array (-L:L,0:L) containing the values of the solid spherical harmonics for all l .le. L and the corresponding values of m.
[in]X,Y,Z
  real input coordinates at which the solid harmonics will be evaluated.
[in]L
  Integer angular momentum of the last spherical harmonic to be generated.

◆ cfp_sph_shell_to_cart_lshells()

subroutine, public special_functions_gbl::cfp_sph_shell_to_cart_lshells ( integer, intent(in)  l,
integer, dimension(:), intent(out)  nz,
real(kind=cfp), dimension(:), intent(out)  c_nz,
integer, dimension(:), intent(out)  nz_can 
)

This routine constructs for all real solid harmonic with L .le. l the matrix of coefficients for all M of the transformation from the cartesian harmonics (~x^i*y^j*z^k, i+j+k=L). On output the linear array c(:) emulates 2D matrix with 2*L+1 rows and (L+1)*(L+2)/2 columns. The rows of c correspond to the M values: -L,-L+1,...,0,1,...,L.
The coefficients in each row are ordered so that the n-th coefficient in the given row of c corresponds to the coefficient for the cartesian harmonic with the shell canonical index n.

Here is the call graph for this function:

◆ cfp_sph_shell_to_cart_shell()

subroutine, public special_functions_gbl::cfp_sph_shell_to_cart_shell ( integer, intent(in)  l,
real(kind=cfp), dimension(:), intent(out)  c 
)

This routine constructs for a given real solid harmonic L the matrix of coefficients for all M of the transformation from the cartesian harmonics (~x^i*y^j*z^k, i+j+k=L). On output the linear array c(:) emulates 2D matrix with 2*L+1 rows and (L+1)*(L+2)/2 columns. The rows of c correspond to the M values: -L,-L+1,...,0,1,...,L.
The coefficients in each row are ordered so that the n-th coefficient in the given row of c corresponds to the coefficient for the cartesian harmonic with the shell canonical index n.

Here is the call graph for this function:

◆ cfp_sph_to_cart_mapping()

subroutine, public special_functions_gbl::cfp_sph_to_cart_mapping ( integer, intent(in)  l,
integer, intent(in)  m,
real(kind=cfp), dimension(:), intent(out), allocatable  c,
integer, dimension(:), intent(out), allocatable  i_exp,
integer, dimension(:), intent(out), allocatable  j_exp,
integer, dimension(:), intent(out), allocatable  k_exp 
)

This routine constructs for a given real solid harmonic L,M the list of cartesian functions which build up this real solid harmonic. On output the array c(:) contains the coefficients with which the individual cartesians contribute to the solid harmonic. Exponents of the cartesians are then given in the arrays i_exp,j_exp,k_exp.

Todo:
the construction of the coefficient list (the summations) should be included along the lines of sph_shell_to_cart_shell
Here is the call graph for this function:

◆ cfp_zhar()

subroutine, public special_functions_gbl::cfp_zhar ( real(kind=cfp), dimension(-l:l,0:l), intent(out)  SH,
real(kind=cfp), intent(in)  X,
real(kind=cfp), intent(in)  Y,
real(kind=cfp), intent(in)  Z,
integer, intent(in)  L 
)
Purpose:
 Calculate the full set of real unnormalized solid spherical harmonics with l,m .le. L using recursion. See Helgaker, p.210, p.218 (Section 6.4) for details.
 The unnormalized solid spherical harmonics are defined in Rico, Lopez, et al. Int. J. Q. Chem.,2012,DOI:10.1002/qua.24356.
 The mehod of calculation is to calculate the real solid spherical harmonics and then renormalize them to obtain the surface 
 spherical harmonics using the formula:

\[ z_{lm}(x,y,z) = S_{lm}(x,y,z)\sqrt{\frac{1+\delta_{m,0}}{2}\frac{(l+|m|)!}{(l-|m|)!}}. \]

Todo:
{rewrite the algorithm to calculate zlm directly using Helgaker's algorithm for the scaled solid harmonics rather than by calculating Slm first and then renormalizing}
Todo:
{add error checking for case L > size(SH(,:))}
Parameters
[out]SH
  Real array (-L:L,0:L) containing the values of the real unnormalized solidspherical harmonics for all l .le. L and the corresponding values of m.
[in]X,Y,Z
  Real input coordinates for which the harmonics will be evaluated.
[in]L
  Integer angular momentum of the last harmonic to be generated.
Here is the call graph for this function:

◆ dfact()

elemental integer function special_functions_gbl::dfact ( integer, intent(in)  n)

Computes the double factorial \((n)!!\) of the integer number n. n must be .ge. -1.

◆ ipair()

elemental integer function, public special_functions_gbl::ipair ( integer, intent(in)  i)

This function is used to index an ordered pair of values (i,j), where i .ge. j, i .ge. 1. The index of the ordered pair is: ipair(i) + j. This process of indexing can be used in a nested way to index quartets of integers, say (i,j,k,l). We compute ij = ipair(i)+j and kl = ipair(k)+l. The index of the quartet (i,j,k,l) is: ipair(ij)+kl. We assume that i.ge.j, k.ge.l, ij.ge.kl. This nesting (triangularization) is used heavilly to index the 2-particle symmetric integrals.

◆ unpack_pq()

subroutine, public special_functions_gbl::unpack_pq ( integer, intent(in)  pq,
integer, intent(in)  n,
integer, intent(out)  p,
integer, intent(out)  q 
)

Assuming \( pq = p + (q-1)*n \) this routine returns p,q given pq,n.