Multidip 1.0
Multi-photon matrix elements
Loading...
Searching...
No Matches
multidip_integ Module Reference

Special integrals needed by MULTIDIP. More...

Functions/Subroutines

complex(wp) function nested_exp_integ (z, a, c, n, m, s, k)
 Multi-dimensional triangular integral of exponentials and powers.
complex(wp) function nested_coul_integ (z, a, c, n, m, s, l, k)
 Multi-dimensional triangular integral of Coulomb-Hankel functions and powers.
complex(wp) function nested_cgreen_integ (z, a, r0, c, n, sa, sb, m, l, k)
 Multi-dimensional integral of Coulomb-Hankel and Coulomb-Green functions and powers.

Detailed Description

Special integrals needed by MULTIDIP.

Author
J Benda
Date
2020 - 2023

Module containing routines for calculation of the multi-dimensional dipole integrals used in outer correction of transition dipole elements in MULTIDIP.

Function/Subroutine Documentation

◆ nested_cgreen_integ()

complex(wp) function multidip_integ::nested_cgreen_integ ( real(wp), intent(in) z,
real(wp), intent(in) a,
real(wp), intent(in) r0,
real(wp), intent(in) c,
integer, intent(in) n,
integer, intent(in) sa,
integer, intent(in) sb,
integer, dimension(:), intent(in) m,
integer, dimension(:), intent(in) l,
complex(wp), dimension(:), intent(in) k )

Multi-dimensional integral of Coulomb-Hankel and Coulomb-Green functions and powers.

Author
J Benda
Date
2020 - 2024

Evaluate the many-dimensional integral

\[ \int\limits_a^{+\infty} \dots \int\limits_a^{+\infty} H_N(r_N) \dots g_3(r_3, r_2) r_2^{m_2} g_2(r_2, r_1) r_1^{m_1} H_1(r_1) \mathrm{d}r_1 \dots \]

using the asymptotic form of Coulomb-Hankel functions. The arrays passed to this function need to be ordered from the right-most integral to the right. This function iterates over all possible orderings of \( (r_1, r_2, \dots, r_N) \) and for each of these it splits the integral into (hyper-)triangular integrals and integrates those using nested_coul_integ.

For small values of "a" it may increase the accuracy if the leading interval (a,r0) (or more specifically the set Q = (a,+∞)^N \ (b,+∞)^N) is integrated numerically. For such use, provide r0 > a. Otherwise set r0 to zero.

When the global paramter coulomb_check is on, the integral is not attempted if any of the Coulomb functions in the integrand is not sufficiently well approximated by the asymptotic form at the R-matrix radius (or at the asymptotic radius if given). In such a case, the function returns the NaN constant.

The one-dimensional case with c = 0, m = 0, l1 = l2 and real momenta is treated in a special way using a closed-form formula from the article: H. F. Arnoldus, T. F. George: Analytical evaluation of elastic Coulomb integrals, J. Math. Phys. 33 (1992) 578–583.

The (adapted) formula reads:

\[ \lim_{c \rightarrow 0+} \int\limits_a^{+\infty} X_l(k_2; r) Y_l(k_1; r) \exp(-cr) dr = \frac{X_l'(k_2; a) Y_l(k_1; a) - X_l(k_2; a) Y_l'(k_1; a)}{k_2^2 - k_1^2} \,, \]

where \( X_l \) and \( Y_l \) stand for arbitrary Coulomb functions and primes denote derivatives with respect to the radius.

Parameters
ZResidual ion charge
aLower bound of the integration
r0Optional radius from which to apply asymptotic integrals
cDamping factor (additional exp(-c*r) added to all r^m functions)
NDimension (number of integration variables)
saSign of the right-most Coulomb-Hankel function
sbSign of the left-most Coulomb-Hankel function
mArray of integer powers of length N
lArray of angular momenta of length N + 1
kArray of linear momenta of length N + 1

Definition at line 355 of file multidip_integ.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ nested_coul_integ()

complex(wp) function multidip_integ::nested_coul_integ ( real(wp), intent(in) z,
real(wp), intent(in) a,
real(wp), intent(in) c,
integer, intent(in) n,
integer, dimension(1:n), intent(in) m,
integer, dimension(1:2*n), intent(in) s,
integer, dimension(1:2*n), intent(in) l,
complex(wp), dimension(1:2*n), intent(in) k )

Multi-dimensional triangular integral of Coulomb-Hankel functions and powers.

Author
J Benda
Date
2020 - 2023

Evaluate the nested many-dimensional triangular integral

\[ \int\limits_a^{+\infty} \dots \int\limits_{r_3}^{+\infty} H_4(r_2) r_2^{m_2} H_3(r_2) \int\limits_{r_2}^{+\infty} H_2(r_1) r_1^{m_1} H_1(r_1) \mathrm{d}r_1 \mathrm{d}r_2 \mathrm{d}r_3 \dots \]

using the asymptotic form of Coulomb-Hankel functions. The arrays passed to this function need to be ordered from the inner-most (right-most) integral outward. For each term, use nested_exp_integ.

Parameters
ZResidual ion charge
aLower bound
cDamping factor (additional exp(-c*r) added to all r^m functions)
NDimension (number of integration variables)
mArray of integer powers of length N
sArray of integer signs (+1 or -1) in exponents of length 2*N
lArray of angular momenta of length 2*N
kArray of linear momenta of length 2*N

Definition at line 227 of file multidip_integ.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ nested_exp_integ()

complex(wp) function multidip_integ::nested_exp_integ ( real(wp), intent(in) z,
real(wp), intent(in) a,
real(wp), intent(in) c,
integer, intent(in) n,
integer, dimension(0:n-1), intent(in) m,
integer, dimension(0:2*n-1), intent(in) s,
complex(wp), dimension(0:2*n-1), intent(in) k )

Multi-dimensional triangular integral of exponentials and powers.

Author
J Benda
Date
2020 - 2023

Evaluate the nested many-dimensional triangular integral

\[ \int\limits_a^{+\infty} \dots \int\limits_{r_3}^{+\infty} \mathrm{e}^{\mathrm{i}s_4 \theta_4(r_2)} r_2^{m_2} \mathrm{e}^{\mathrm{i}s_3 \theta_3(r_2)} \int\limits_{r_2}^{+\infty} \mathrm{e}^{\mathrm{i}s_2 \theta_2(r_1)} r_1^{m_1} \mathrm{e}^{\mathrm{i}s_1 \theta_1(r_1)} \mathrm{d}r_1 \mathrm{d}r_2 \mathrm{d}r_3 \dots \]

where

\[ \theta_1(r_1) = k_1 r_1 + \log(2 k_1 r_1)/k_1 - \pi l_1/2 + \sigma_{l_1}(k_1) \]

is the asymptotic phase of a Coulomb function. The arrays passed to this function need to be ordered from the inner-most (right-most) integral outward. The function result does not contain the overall phase and damp factor, which needs to be added manually.

Parameters
ZResidual ion charge
aLower bound
cDamping factor (additional exp(-c*r) added to all r^m functions)
NDimension (number of integration variables)
mArray of integer powers of length N
sArray of integer signs (+1 or -1) in exponents of length 2*N
kArray of linear momenta of length 2*N

Definition at line 68 of file multidip_integ.f90.

Here is the caller graph for this function: