Multidip  1.0
Multi-photon matrix elements
multidip_romberg Module Reference

Romberg quadrature for numerical integration. More...

Functions/Subroutines

subroutine nested_cgreen_correct_romberg (Z, a, b, c, N, sa, sb, m, l, k, v)
 Numerically correct asymptotic approximation of the Coulomb-Green integral (Romberg integration) More...
 
complex(wp) function nested_cgreen_eval (Z, c, N, sa, sb, m, l, k, rs, asy)
 Evaluate the Coulomb-Green's integrand. More...
 

Detailed Description

Romberg quadrature for numerical integration.

Author
J Benda
Date
2021

Function/Subroutine Documentation

◆ nested_cgreen_correct_romberg()

subroutine multidip_romberg::nested_cgreen_correct_romberg ( real(wp), intent(in)  Z,
real(wp), intent(in)  a,
real(wp), intent(in)  b,
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,
complex(wp), intent(inout)  v 
)

Numerically correct asymptotic approximation of the Coulomb-Green integral (Romberg integration)

Author
J Benda
Date
2021 - 2023

This function computes the integral of the Coulomb-Green's integrand over Q = (a,+∞)^N \ (b,+∞)^N numerically. Currently, Romberg integration based on the trapezoidal rule is used.

Parameters
ZResidual ion charge
aLower bound
bUpper bound
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
vValue of the asymptotic integral to correct (updated on exit from this subroutine)
Warning
This is currently only implemented for one-dimensional case containing no Green's function at all.
Todo:
Generalize for arbitrary dimension!

Definition at line 58 of file multidip_romberg.f90.

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

◆ nested_cgreen_eval()

complex(wp) function multidip_romberg::nested_cgreen_eval ( real(wp), intent(in)  Z,
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,
real(wp), dimension(:), intent(in)  rs,
logical, intent(in)  asy 
)

Evaluate the Coulomb-Green's integrand.

Author
J Benda
Date
2021 - 2024
Parameters
ZResidual ion charge
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
rsSingle point in the multidimensional position space where the evaluate the integrand
asyUse the asymptotic form of the Coulomb-Green's functions.
Warning
This is currently only implemented for one-dimensional case containing no Green's function at all.
Todo:
Generalize to arbitrary dimension!

Definition at line 142 of file multidip_romberg.f90.

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