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.qawf;
9 
10 
11 import std.conv;
12 import scid.ports.quadpack.qawfe;
13 
14 version(unittest)
15 {
16     import std.math;
17     import scid.core.testing;
18 }
19 
20 
21 
22 ///
23 void qawf(Real, Func)(Func f, Real a, Real omega, int integr, Real epsabs,
24     out Real result, out Real abserr, out int neval, out int ier,
25     int limlst, out int lst, int leniw, int maxp1, int lenw,
26     int* iwork, Real* work)
27 {
28 //***begin prologue  dqawf
29 //***date written   800101   (yymmdd)
30 //***revision date  830518   (yymmdd)
31 //***category no.  h2a3a1
32 //***keywords  automatic integrator, special-purpose,fourier
33 //             integral, integration between zeros with dqawoe,
34 //             convergence acceleration with dqelg
35 //***author  piessens,robert ,appl. math. & progr. div. - k.u.leuven
36 //           de doncker,elise,appl. math & progr. div. - k.u.leuven
37 //***purpose  the routine calculates an approximation result to a given
38 //            fourier integral i=integral of f(x)*w(x) over (a,infinity)
39 //            where w(x) = cos(omega*x) or w(x) = sin(omega*x).
40 //            hopefully satisfying following claim for accuracy
41 //            abs(i-result).le.epsabs.
42 //***description
43 //
44 //        computation of fourier integrals
45 //        standard fortran subroutine
46 //        double precision version
47 //
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 //            a      - double precision
57 //                     lower limit of integration
58 //
59 //            omega  - double precision
60 //                     parameter in the integrand weight function
61 //
62 //            integr - integer
63 //                     indicates which of the weight functions is used
64 //                     integr = 1      w(x) = cos(omega*x)
65 //                     integr = 2      w(x) = sin(omega*x)
66 //                     if integr.ne.1.and.integr.ne.2, the routine
67 //                     will end with ier = 6.
68 //
69 //            epsabs - double precision
70 //                     absolute accuracy requested, epsabs.gt.0.
71 //                     if epsabs.le.0, the routine will end with ier = 6.
72 //
73 //         on return
74 //            result - double precision
75 //                     approximation to the integral
76 //
77 //            abserr - double precision
78 //                     estimate of the modulus of the absolute error,
79 //                     which should equal or exceed abs(i-result)
80 //
81 //            neval  - integer
82 //                     number of integrand evaluations
83 //
84 //            ier    - integer
85 //                     ier = 0 normal and reliable termination of the
86 //                             routine. it is assumed that the requested
87 //                             accuracy has been achieved.
88 //                     ier.gt.0 abnormal termination of the routine.
89 //                             the estimates for integral and error are
90 //                             less reliable. it is assumed that the
91 //                             requested accuracy has not been achieved.
92 //            error messages
93 //                    if omega.ne.0
94 //                     ier = 1 maximum number of cycles allowed
95 //                             has been achieved, i.e. of subintervals
96 //                             (a+(k-1)c,a+kc) where
97 //                             c = (2*int(abs(omega))+1)*pi/abs(omega),
98 //                             for k = 1, 2, ..., lst.
99 //                             one can allow more cycles by increasing
100 //                             the value of limlst (and taking the
101 //                             according dimension adjustments into
102 //                             account). examine the array iwork which
103 //                             contains the error flags on the cycles, in
104 //                             order to look for eventual local
105 //                             integration difficulties.
106 //                             if the position of a local difficulty
107 //                             can be determined (e.g. singularity,
108 //                             discontinuity within the interval) one
109 //                             will probably gain from splitting up the
110 //                             interval at this point and calling
111 //                             appropriate integrators on the subranges.
112 //                         = 4 the extrapolation table constructed for
113 //                             convergence accelaration of the series
114 //                             formed by the integral contributions over
115 //                             the cycles, does not converge to within
116 //                             the requested accuracy.
117 //                             as in the case of ier = 1, it is advised
118 //                             to examine the array iwork which contains
119 //                             the error flags on the cycles.
120 //                         = 6 the input is invalid because
121 //                             (integr.ne.1 and integr.ne.2) or
122 //                              epsabs.le.0 or limlst.lt.1 or
123 //                              leniw.lt.(limlst+2) or maxp1.lt.1 or
124 //                              lenw.lt.(leniw*2+maxp1*25).
125 //                              result, abserr, neval, lst are set to
126 //                              zero.
127 //                         = 7 bad integrand behaviour occurs within
128 //                             one or more of the cycles. location and
129 //                             type of the difficulty involved can be
130 //                             determined from the first lst elements of
131 //                             vector iwork.  here lst is the number of
132 //                             cycles actually needed (see below).
133 //                             iwork(k) = 1 the maximum number of
134 //                                          subdivisions (=(leniw-limlst)
135 //                                          /2) has been achieved on the
136 //                                          k th cycle.
137 //                                      = 2 occurrence of roundoff error
138 //                                          is detected and prevents the
139 //                                          tolerance imposed on the k th
140 //                                          cycle, from being achieved
141 //                                          on this cycle.
142 //                                      = 3 extremely bad integrand
143 //                                          behaviour occurs at some
144 //                                          points of the k th cycle.
145 //                                      = 4 the integration procedure
146 //                                          over the k th cycle does
147 //                                          not converge (to within the
148 //                                          required accuracy) due to
149 //                                          roundoff in the extrapolation
150 //                                          procedure invoked on this
151 //                                          cycle. it is assumed that the
152 //                                          result on this interval is
153 //                                          the best which can be
154 //                                          obtained.
155 //                                      = 5 the integral over the k th
156 //                                          cycle is probably divergent
157 //                                          or slowly convergent. it must
158 //                                          be noted that divergence can
159 //                                          occur with any other value of
160 //                                          iwork(k).
161 //                    if omega = 0 and integr = 1,
162 //                    the integral is calculated by means of dqagie,
163 //                    and ier = iwork(1) (with meaning as described
164 //                    for iwork(k),k = 1).
165 //
166 //         dimensioning parameters
167 //            limlst - integer
168 //                     limlst gives an upper bound on the number of
169 //                     cycles, limlst.ge.3.
170 //                     if limlst.lt.3, the routine will end with ier = 6.
171 //
172 //            lst    - integer
173 //                     on return, lst indicates the number of cycles
174 //                     actually needed for the integration.
175 //                     if omega = 0, then lst is set to 1.
176 //
177 //            leniw  - integer
178 //                     dimensioning parameter for iwork. on entry,
179 //                     (leniw-limlst)/2 equals the maximum number of
180 //                     subintervals allowed in the partition of each
181 //                     cycle, leniw.ge.(limlst+2).
182 //                     if leniw.lt.(limlst+2), the routine will end with
183 //                     ier = 6.
184 //
185 //            maxp1  - integer
186 //                     maxp1 gives an upper bound on the number of
187 //                     chebyshev moments which can be stored, i.e. for
188 //                     the intervals of lengths abs(b-a)*2**(-l),
189 //                     l = 0,1, ..., maxp1-2, maxp1.ge.1.
190 //                     if maxp1.lt.1, the routine will end with ier = 6.
191 //            lenw   - integer
192 //                     dimensioning parameter for work
193 //                     lenw must be at least leniw*2+maxp1*25.
194 //                     if lenw.lt.(leniw*2+maxp1*25), the routine will
195 //                     end with ier = 6.
196 //
197 //         work arrays
198 //            iwork  - integer
199 //                     vector of dimension at least leniw
200 //                     on return, iwork(k) for k = 1, 2, ..., lst
201 //                     contain the error flags on the cycles.
202 //
203 //            work   - double precision
204 //                     vector of dimension at least
205 //                     on return,
206 //                     work(1), ..., work(lst) contain the integral
207 //                      approximations over the cycles,
208 //                     work(limlst+1), ..., work(limlst+lst) contain
209 //                      the error extimates over the cycles.
210 //                     further elements of work have no specific
211 //                     meaning for the user.
212 //
213 //***references  (none)
214 //***routines called  dqawfe,xerror
215 //***end prologue  dqawf
216 //
217        int last, limit, ll2,l1,l2,l3,l4,l5,l6;
218 //
219 //         check validity of limlst, leniw, maxp1 and lenw.
220 //
221 //***first executable statement  dqawf
222       ier = 6;
223       neval = 0;
224       last = 0;
225       result = 0.0;
226       abserr = 0.0;
227       if(limlst < 3 || leniw < (limlst+2) || maxp1 < 1 || lenw <
228          (leniw*2+maxp1*25)) goto l10;
229 //
230 //         prepare call for dqawfe
231 //
232       limit = (leniw-limlst)/2;
233       l1 = limlst;
234       l2 = limlst+l1;
235       l3 = limit+l2;
236       l4 = limit+l3;
237       l5 = limit+l4;
238       l6 = limit+l5;
239       ll2 = limit+l1;
240       qawfe!(Real,Func)(f,a,omega,integr,epsabs,limlst,limit,maxp1,result,
241         abserr,neval,ier,work,work+l1,iwork,lst,work+l2,
242         work+l3,work+l4,work+l5,iwork+l1,iwork+ll2,work+l6);
243 //
244 //         call error handler if necessary
245 //
246 l10:  if(ier != 0)
247         throw new Exception("abnormal return from qawf: "~to!string(ier));
248       return;
249 }
250 
251 
252 unittest
253 {
254     alias qawf!(float, float delegate(float)) fqawf;
255     alias qawf!(double, double delegate(double)) dqawf;
256     alias qawf!(double, double function(double)) dfqawf;
257     alias qawf!(real, real delegate(real)) rqawf;
258 }
259 
260 
261 unittest
262 {
263     double f(double x) { return x > 0.0 ? 1/sqrt(x) : 0.0; }
264 
265     enum : double
266     {
267         a = 0.0,
268         omega = PI_2,
269         epsabs = 1e-8,
270     }
271     double result, abserr;
272     int neval, ier, lst;
273     enum
274     {
275         integr = 1, // cos(pi x/2)
276         limlst = 50,
277         leniw = 1050,
278         maxp1 = 21,
279         lenw = leniw*2 + maxp1*25
280     }
281 
282     int[leniw] iwork;
283     double[lenw] work;
284 
285     qawf(&f, a, omega, integr, epsabs, result, abserr, neval, ier, limlst,
286         lst, leniw, maxp1, lenw, iwork.ptr, work.ptr);
287 }
288 
289 
290 unittest
291 {
292     // This is integral 14 in the QUADPACK book.
293     real alpha;
294     real f(real x) { return x <= 0.0 ? 0.0 : exp(-(2.0L^^(-alpha))*x)/sqrt(x); }
295     enum real omega = 1.0;
296     enum integr = 1; // cos(x)
297     
298     enum : real
299     {
300         a = 0.0,
301         epsabs = 1e-8
302     }
303     real result, abserr;
304     int neval, ier, lst;
305     enum
306     {
307         limlst = 50,
308         leniw = 1050,
309         maxp1 = 21,
310         lenw = leniw*2 + maxp1*25
311     }
312 
313     int[leniw] iwork;
314     real[lenw] work;
315 
316     alpha = 0.0;
317     qawf(&f, a, omega, integr, epsabs, result, abserr, neval, ier, limlst,
318         lst, leniw, maxp1, lenw, iwork.ptr, work.ptr);
319     real ans = sqrt(PI) * ((1+(4.0L^^(-alpha)))^^(-0.25L))
320         * cos(atan(2.0L^^alpha)/2);
321     assert (isAccurate(result, abserr, ans, 0.0L, epsabs));
322 }