Multidip 1.0
Multi-photon matrix elements
Loading...
Searching...
No Matches
multidip_romberg.f90
Go to the documentation of this file.
1! Copyright 2020
2!
3! For a comprehensive list of the developers that contributed to these codes
4! see the UK-AMOR website.
5!
6! This file is part of UKRmol-out (UKRmol+ suite).
7!
8! UKRmol-out is free software: you can redistribute it and/or modify
9! it under the terms of the GNU General Public License as published by
10! the Free Software Foundation, either version 3 of the License, or
11! (at your option) any later version.
12!
13! UKRmol-out is distributed in the hope that it will be useful,
14! but WITHOUT ANY WARRANTY; without even the implied warranty of
15! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16! GNU General Public License for more details.
17!
18! You should have received a copy of the GNU General Public License
19! along with UKRmol-out (in source/COPYING). Alternatively, you can also visit
20! <https://www.gnu.org/licenses/>.
21!
22!> \brief Romberg quadrature for numerical integration
23!> \author J Benda
24!> \date 2021
25!>
27
28 use precisn_gbl, only: wp
29
30 implicit none
31
32contains
33
34 !> \brief Numerically correct asymptotic approximation of the Coulomb-Green integral (Romberg integration)
35 !> \author J Benda
36 !> \date 2021 - 2023
37 !>
38 !> This function computes the integral of the Coulomb-Green's integrand over Q = (a,+∞)^N \ (b,+∞)^N numerically.
39 !> Currently, Romberg integration based on the trapezoidal rule is used.
40 !>
41 !> \param Z Residual ion charge
42 !> \param a Lower bound
43 !> \param b Upper bound
44 !> \param c Damping factor (additional exp(-c*r) added to all r^m functions)
45 !> \param N Dimension (number of integration variables)
46 !> \param sa Sign of the right-most Coulomb-Hankel function
47 !> \param sb Sign of the left-most Coulomb-Hankel function
48 !> \param m Array of integer powers of length N
49 !> \param l Array of angular momenta of length N + 1
50 !> \param k Array of linear momenta of length N + 1
51 !> \param v Value of the asymptotic integral to correct (updated on exit from this subroutine)
52 !>
53 !> \warning This is currently only implemented for one-dimensional case containing no
54 !> Green's function at all.
55 !>
56 !> \todo Generalize for arbitrary dimension!
57 !>
58 subroutine nested_cgreen_correct_romberg (Z, a, b, c, N, sa, sb, m, l, k, v)
59
62
63 integer, intent(in) :: N, sa, sb, m(:), l(:)
64 real(wp), intent(in) :: Z, a, b, c
65 complex(wp), intent(in) :: k(:)
66 complex(wp), intent(inout) :: v
67
68 integer :: level, i
69 real(wp) :: h, r
70
71 complex(wp) :: romb(max_romb_level) ! sequence of Romberg estimations, with the highest order placed in the first element
72 complex(wp) :: integ, y, ya, yb, error, prev
73
74 h = b - a
75 ya = nested_cgreen_eval(z, c, n, sa, sb, m, l, k, [ a ], .false.)
76 yb = nested_cgreen_eval(z, c, n, sa, sb, m, l, k, [ b ], .false.)
77 romb(1) = h * (ya + yb) / 2
78
79 do level = 2, max_romb_level
80
81 ! remember the previous best Romberg estimate
82 prev = romb(1)
83
84 ! calculate next trapezoidal estimate from the previous one
85 integ = 0
86 error = 0
87 h = h / 2
88 do i = 1, 2**(level - 2)
89 r = a + (2*i - 1)*h
90 y = nested_cgreen_eval(z, c, n, sa, sb, m, l, k, [ r ], .false.)
91 call kahan_add(integ, y, error)
92 end do
93 romb(level) = romb(level - 1)/2 + h*integ
94
95 ! perform the Romberg extrapolation step
96 do i = level - 1, 1, -1
97 romb(i) = (4**(level - 1) * romb(i + 1) - romb(i)) / (4**(level - 1) - 1)
98 end do
99
100 ! abort on NaN
101 if (.not. romb(1) == romb(1)) then
102 print '(a)', 'WARNING: NaN encountered during Romberg quadrature'
103 exit
104 end if
105
106 ! compare the current best estimate to the previous one
107 if (abs(romb(1) - prev) <= epsrel * (abs(romb(1)) + abs(prev))) then
108 exit
109 end if
110
111 end do
112
113 if (level > max_romb_level) then
114 print '(A,I0,A)', 'WARNING: Romberg quadrature did not converge in ', max_romb_level, ' subdivisions'
115 end if
116
117 v = v + romb(1)
118
119 end subroutine nested_cgreen_correct_romberg
120
121
122 !> \brief Evaluate the Coulomb-Green's integrand
123 !> \author J Benda
124 !> \date 2021 - 2024
125 !>
126 !> \param Z Residual ion charge
127 !> \param c Damping factor (additional exp(-c*r) added to all r^m functions)
128 !> \param N Dimension (number of integration variables)
129 !> \param sa Sign of the right-most Coulomb-Hankel function
130 !> \param sb Sign of the left-most Coulomb-Hankel function
131 !> \param m Array of integer powers of length N
132 !> \param l Array of angular momenta of length N + 1
133 !> \param k Array of linear momenta of length N + 1
134 !> \param rs Single point in the multidimensional position space where the evaluate the integrand
135 !> \param asy Use the asymptotic form of the Coulomb-Green's functions.
136 !>
137 !> \warning This is currently only implemented for one-dimensional case containing no
138 !> Green's function at all.
139 !>
140 !> \todo Generalize to arbitrary dimension!
141 !>
142 complex(wp) function nested_cgreen_eval (Z, c, N, sa, sb, m, l, k, rs, asy) result (val)
143
144 use mpi_gbl, only: mpi_xermsg
145 use multidip_special, only: coulh, coulh_asy
146
147 integer, intent(in) :: n, sa, sb, m(:), l(:)
148 real(wp), intent(in) :: z, c, rs(:)
149 complex(wp), intent(in) :: k(:)
150 logical, intent(in) :: asy
151
152 real(wp) :: ek(n+1)
153 complex(wp) :: h1, h2
154
155 ek = real(k(1:n+1)**2, wp)/2
156
157 ! WARNING: Only one-dimensional at the moment (no Green's functions)
158 if (n > 1) then
159 call mpi_xermsg('multidip_romberg', 'nested_cgreen_eval', 'nested_cgreen_eval not implemented for higher orders', 1, 1)
160 end if
161
162 if (asy) then
163 h1 = coulh_asy(z, sa, l(1), ek(1), rs(1))
164 h2 = coulh_asy(z, sb, l(2), ek(2), rs(1))
165 else
166 h1 = coulh(z, sa, l(1), ek(1), rs(1))
167 h2 = coulh(z, sb, l(2), ek(2), rs(1))
168 end if
169
170 val = h2 * rs(1)**m(1) * exp(-c*rs(1)) * h1
171
172 end function nested_cgreen_eval
173
174end module multidip_romberg
Hard-coded parameters of MULTIDIP.
integer, parameter max_romb_level
Maximal nesting level for Romberg integration.
real(wp) epsrel
Romberg integration relative tolerance for convergence.
Romberg quadrature for numerical integration.
complex(wp) function nested_cgreen_eval(z, c, n, sa, sb, m, l, k, rs, asy)
Evaluate the Coulomb-Green's integrand.
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)
Special functions and objects used by MULTIDIP.
subroutine kahan_add(x, dx, err)
Compensated summation.