1 // This code has been mechanically translated from the original FORTRAN
2 // code at http://netlib.org/quadpack.
3 
4 /** Authors:    Lars Tandle Kyllingstad
5     Copyright:  Copyright (c) 2009, Lars T. Kyllingstad. All rights reserved.
6     License:    Boost License 1.0
7 */
8 module scid.ports.quadpack.qagi;
9 
10 
11 import std.conv;
12 import scid.ports.quadpack.qagie;
13 
14 version(unittest)
15 {
16     import std.math;
17     import scid.core.testing;
18 }
19 
20 
21 
22 
23 ///
24 void qagi(Real, Func)(Func f, Real bound, int inf, Real epsabs,
25     Real epsrel, out Real result, out Real abserr, out int neval,
26     out int ier, int limit, int lenw, out int last,
27     int* iwork, Real* work)
28 {
29 //***begin prologue  dqagi
30 //***date written   800101   (yymmdd)
31 //***revision date  830518   (yymmdd)
32 //***category no.  h2a3a1,h2a4a1
33 //***keywords  automatic integrator, infinite intervals,
34 //             general-purpose, transformation, extrapolation,
35 //             globally adaptive
36 //***author  piessens,robert,appl. math. & progr. div. - k.u.leuven
37 //           de doncker,elise,appl. math. & progr. div. -k.u.leuven
38 //***purpose  the routine calculates an approximation result to a given
39 //            integral   i = integral of f over (bound,+infinity)
40 //            or i = integral of f over (-infinity,bound)
41 //            or i = integral of f over (-infinity,+infinity)
42 //            hopefully satisfying following claim for accuracy
43 //            abs(i-result).le.max(epsabs,epsrel*abs(i)).
44 //***description
45 //
46 //        integration over infinite intervals
47 //        standard fortran subroutine
48 //
49 //        parameters
50 //         on entry
51 //            f      - double precision
52 //                     function subprogram defining the integrand
53 //                     function f(x). the actual name for f needs to be
54 //                     declared e x t e r n a l in the driver program.
55 //
56 //            bound  - double precision
57 //                     finite bound of integration range
58 //                     (has no meaning if interval is doubly-infinite)
59 //
60 //            inf    - integer
61 //                     indicating the kind of integration range involved
62 //                     inf = 1 corresponds to  (bound,+infinity),
63 //                     inf = -1            to  (-infinity,bound),
64 //                     inf = 2             to (-infinity,+infinity).
65 //
66 //            epsabs - double precision
67 //                     absolute accuracy requested
68 //            epsrel - double precision
69 //                     relative accuracy requested
70 //                     if  epsabs.le.0
71 //                     and epsrel.lt.max(50*rel.mach.acc.,0.5d-28),
72 //                     the routine will end with ier = 6.
73 //
74 //
75 //         on return
76 //            result - double precision
77 //                     approximation to the integral
78 //
79 //            abserr - double precision
80 //                     estimate of the modulus of the absolute error,
81 //                     which should equal or exceed abs(i-result)
82 //
83 //            neval  - integer
84 //                     number of integrand evaluations
85 //
86 //            ier    - integer
87 //                     ier = 0 normal and reliable termination of the
88 //                             routine. it is assumed that the requested
89 //                             accuracy has been achieved.
90 //                   - ier.gt.0 abnormal termination of the routine. the
91 //                             estimates for result and error are less
92 //                             reliable. it is assumed that the requested
93 //                             accuracy has not been achieved.
94 //            error messages
95 //                     ier = 1 maximum number of subdivisions allowed
96 //                             has been achieved. one can allow more
97 //                             subdivisions by increasing the value of
98 //                             limit (and taking the according dimension
99 //                             adjustments into account). however, if
100 //                             this yields no improvement it is advised
101 //                             to analyze the integrand in order to
102 //                             determine the integration difficulties. if
103 //                             the position of a local difficulty can be
104 //                             determined (e.g. singularity,
105 //                             discontinuity within the interval) one
106 //                             will probably gain from splitting up the
107 //                             interval at this point and calling the
108 //                             integrator on the subranges. if possible,
109 //                             an appropriate special-purpose integrator
110 //                             should be used, which is designed for
111 //                             handling the type of difficulty involved.
112 //                         = 2 the occurrence of roundoff error is
113 //                             detected, which prevents the requested
114 //                             tolerance from being achieved.
115 //                             the error may be under-estimated.
116 //                         = 3 extremely bad integrand behaviour occurs
117 //                             at some points of the integration
118 //                             interval.
119 //                         = 4 the algorithm does not converge.
120 //                             roundoff error is detected in the
121 //                             extrapolation table.
122 //                             it is assumed that the requested tolerance
123 //                             cannot be achieved, and that the returned
124 //                             result is the best which can be obtained.
125 //                         = 5 the integral is probably divergent, or
126 //                             slowly convergent. it must be noted that
127 //                             divergence can occur with any other value
128 //                             of ier.
129 //                         = 6 the input is invalid, because
130 //                             (epsabs.le.0 and
131 //                              epsrel.lt.max(50*rel.mach.acc.,0.5d-28))
132 //                              or limit.lt.1 or leniw.lt.limit*4.
133 //                             result, abserr, neval, last are set to
134 //                             zero. exept when limit or leniw is
135 //                             invalid, iwork(1), work(limit*2+1) and
136 //                             work(limit*3+1) are set to zero, work(1)
137 //                             is set to a and work(limit+1) to b.
138 //
139 //         dimensioning parameters
140 //            limit - integer
141 //                    dimensioning parameter for iwork
142 //                    limit determines the maximum number of subintervals
143 //                    in the partition of the given integration interval
144 //                    (a,b), limit.ge.1.
145 //                    if limit.lt.1, the routine will end with ier = 6.
146 //
147 //            lenw  - integer
148 //                    dimensioning parameter for work
149 //                    lenw must be at least limit*4.
150 //                    if lenw.lt.limit*4, the routine will end
151 //                    with ier = 6.
152 //
153 //            last  - integer
154 //                    on return, last equals the number of subintervals
155 //                    produced in the subdivision process, which
156 //                    determines the number of significant elements
157 //                    actually in the work arrays.
158 //
159 //         work arrays
160 //            iwork - integer
161 //                    vector of dimension at least limit, the first
162 //                    k elements of which contain pointers
163 //                    to the error estimates over the subintervals,
164 //                    such that work(limit*3+iwork(1)),... ,
165 //                    work(limit*3+iwork(k)) form a decreasing
166 //                    sequence, with k = last if last.le.(limit/2+2), and
167 //                    k = limit+1-last otherwise
168 //
169 //            work  - double precision
170 //                    vector of dimension at least lenw
171 //                    on return
172 //                    work(1), ..., work(last) contain the left
173 //                     end points of the subintervals in the
174 //                     partition of (a,b),
175 //                    work(limit+1), ..., work(limit+last) contain
176 //                     the right end points,
177 //                    work(limit*2+1), ...,work(limit*2+last) contain the
178 //                     integral approximations over the subintervals,
179 //                    work(limit*3+1), ..., work(limit*3)
180 //                     contain the error estimates.
181 //***references  (none)
182 //***routines called  dqagie,xerror
183 //***end prologue  dqagi
184 //
185       int lvl=1,l1=1,l2=1,l3=1;
186 //
187 //         check validity of limit and lenw.
188 //
189 //***first executable statement  dqagi
190       ier = 6;
191       neval = 0;
192       last = 0;
193       result = 0.0;
194       abserr = 0.0;
195       if(limit < 1 || lenw < limit*4) goto l10;
196 //
197 //         prepare call for dqagie.
198 //
199       l1 = limit;
200       l2 = limit+l1;
201       l3 = limit+l2;
202 //
203       qagie!(Real,Func)(f,bound,inf,epsabs,epsrel,limit,result,abserr,
204        neval,ier,work,work+l1,work+l2,work+l3,iwork,last);
205 //
206 //         call error handler if necessary.
207 //
208        lvl = 0;
209 l10:  if(ier == 6) lvl = 1;
210       if(ier != 0)
211         throw new Exception("abnormal return from  qagi: "~to!string(ier));
212       return;
213 }
214 
215 
216 unittest
217 {
218     alias qagi!(float, float delegate(float)) fqagi;
219     alias qagi!(double, double delegate(double)) dqagi;
220     alias qagi!(double, double function(double)) dfqagi;
221     alias qagi!(real, real delegate(real)) rqagi;
222 }
223 
224 unittest
225 {
226     real f(real x) { return exp(-x*x); }
227 
228     enum : real
229     {
230         bound = 0.0,
231         epsabs = 0.0,
232         epsrel = 1e-10
233     }
234     real result, abserr;
235     int neval, ier;
236     enum
237     {
238         limit = 500,
239         lenw = 4*limit
240     }
241     int last;
242     
243     int[limit] iwork;
244     real[lenw] work;
245 
246     enum real sqrtPi = sqrt(PI);
247 
248     int inf = 1; // (bound,+inf)
249     qagi(&f, bound, inf, epsabs, epsrel, result, abserr, neval, ier,
250         limit, lenw, last, iwork.ptr, work.ptr);
251     assert (isAccurate(result, abserr, sqrtPi/2, epsrel, epsabs));
252 
253     inf = -1; // (-inf,bound)
254     qagi(&f, bound, inf, epsabs, epsrel, result, abserr, neval, ier,
255         limit, lenw, last, iwork.ptr, work.ptr);
256     assert (isAccurate(result, abserr, sqrtPi/2, epsrel, epsabs));
257     
258     inf = 2; // (-inf,+inf)
259     qagi(&f, bound, inf, epsabs, epsrel, result, abserr, neval, ier,
260         limit, lenw, last, iwork.ptr, work.ptr);
261 
262     assert (isAccurate(result, abserr, sqrtPi, epsrel, epsabs));
263 }
264 
265 unittest
266 {
267     // Integral 15 in the QUADPACK book.
268     double alpha;
269     double f(double x) { return (x^^2.0) * exp(-(2.0^^(-alpha))*x); }
270     
271     double bound = 0.0, epsabs = 0.0, epsrel = 1e-8, result, abserr;
272     int inf = 1, neval, ier, last;
273     enum { limit = 500, lenw = 4*limit }
274 
275     int[limit] iwork;
276     double[lenw] work;
277 
278     for (alpha=0.0; alpha<5.001; alpha+=1.0)
279     {
280         qagi (&f, bound, inf, epsabs, epsrel, result, abserr, neval, ier,
281             limit, lenw, last, iwork.ptr, work.ptr);
282 
283         double ans = 2.0^^(3*alpha+1.0);
284         assert (isAccurate(result, abserr, ans, epsrel, epsabs));
285     }
286 }
287 
288 
289 unittest
290 {
291     // Integral 16 in the QUADPACK book.
292     double alpha;
293     double f(double x) { return (x^^(alpha-1))/((1+10*x)^^2); }
294     
295     double bound = 0.0, epsabs = 0.0, epsrel = 1e-8, result, abserr;
296     int inf = 1, neval, ier, last;
297     enum { limit = 500, lenw = 4*limit }
298 
299     int[limit] iwork;
300     double[lenw] work;
301 
302     for (int i=-9; i<=9; i++)
303     {
304         alpha = 1.0 + i*0.1;
305         qagi (&f, bound, inf, epsabs, epsrel, result, abserr, neval, ier,
306             limit, lenw, last, iwork.ptr, work.ptr);
307 
308         double ans;
309         if (alpha == 1.0)  ans = 0.1;
310         else ans = (10.0^^(-alpha))*(1-alpha)*PI/sin(PI*alpha);
311 
312         assert (isAccurate(result, abserr, ans, epsrel, epsabs));
313     }
314 }