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.qawfe;
9 
10 
11 import std.algorithm: max;
12 import std.conv;
13 import std.math;
14 
15 import scid.core.fortran;
16 import scid.ports.quadpack.qagie;
17 import scid.ports.quadpack.qawoe;
18 import scid.ports.quadpack.qelg;
19 
20 
21 
22 
23 ///
24 void qawfe(Real, Func)(Func f, Real a, Real omega, int integr, Real epsabs,
25     int limlst, int limit, int maxp1, out Real result, out Real abserr,
26     out int neval, out int ier, Real* rslst_, Real* erlst_, int* ierlst_,
27     out int lst, Real* alist_, Real* blist_, Real* rlist_, Real* elist_,
28     int* iord_, int* nnlog_, Real* chebmo_)
29 {
30 //***begin prologue  dqawfe
31 //***date written   800101   (yymmdd)
32 //***revision date  830518   (yymmdd)
33 //***category no.  h2a3a1
34 //***keywords  automatic integrator, special-purpose,
35 //             fourier integrals,
36 //             integration between zeros with dqawoe,
37 //             convergence acceleration with dqelg
38 //***author  piessens,robert,appl. math. & progr. div. - k.u.leuven
39 //           dedoncker,elise,appl. math. & progr. div. - k.u.leuven
40 //***purpose  the routine calculates an approximation result to a
41 //            given fourier integal
42 //            i = integral of f(x)*w(x) over (a,infinity)
43 //            where w(x)=cos(omega*x) or w(x)=sin(omega*x),
44 //            hopefully satisfying following claim for accuracy
45 //            abs(i-result).le.epsabs.
46 //***description
47 //
48 //        computation of fourier integrals
49 //        standard fortran subroutine
50 //        double precision version
51 //
52 //        parameters
53 //         on entry
54 //            f      - double precision
55 //                     function subprogram defining the integrand
56 //                     function f(x). the actual name for f needs to
57 //                     be declared e x t e r n a l in the driver program.
58 //
59 //            a      - double precision
60 //                     lower limit of integration
61 //
62 //            omega  - double precision
63 //                     parameter in the weight function
64 //
65 //            integr - integer
66 //                     indicates which weight function is used
67 //                     integr = 1      w(x) = cos(omega*x)
68 //                     integr = 2      w(x) = sin(omega*x)
69 //                     if integr.ne.1.and.integr.ne.2, the routine will
70 //                     end with ier = 6.
71 //
72 //            epsabs - double precision
73 //                     absolute accuracy requested, epsabs.gt.0
74 //                     if epsabs.le.0, the routine will end with ier = 6.
75 //
76 //            limlst - integer
77 //                     limlst gives an upper bound on the number of
78 //                     cycles, limlst.ge.1.
79 //                     if limlst.lt.3, the routine will end with ier = 6.
80 //
81 //            limit  - integer
82 //                     gives an upper bound on the number of subintervals
83 //                     allowed in the partition of each cycle, limit.ge.1
84 //                     each cycle, limit.ge.1.
85 //
86 //            maxp1  - integer
87 //                     gives an upper bound on the number of
88 //                     chebyshev moments which can be stored, i.e.
89 //                     for the intervals of lengths abs(b-a)*2**(-l),
90 //                     l=0,1, ..., maxp1-2, maxp1.ge.1
91 //
92 //         on return
93 //            result - double precision
94 //                     approximation to the integral x
95 //
96 //            abserr - double precision
97 //                     estimate of the modulus of the absolute error,
98 //                     which should equal or exceed abs(i-result)
99 //
100 //            neval  - integer
101 //                     number of integrand evaluations
102 //
103 //            ier    - ier = 0 normal and reliable termination of
104 //                             the routine. it is assumed that the
105 //                             requested accuracy has been achieved.
106 //                     ier.gt.0 abnormal termination of the routine. the
107 //                             estimates for integral and error are less
108 //                             reliable. it is assumed that the requested
109 //                             accuracy has not been achieved.
110 //            error messages
111 //                    if omega.ne.0
112 //                     ier = 1 maximum number of  cycles  allowed
113 //                             has been achieved., i.e. of subintervals
114 //                             (a+(k-1)c,a+kc) where
115 //                             c = (2*int(abs(omega))+1)*pi/abs(omega),
116 //                             for k = 1, 2, ..., lst.
117 //                             one can allow more cycles by increasing
118 //                             the value of limlst (and taking the
119 //                             according dimension adjustments into
120 //                             account).
121 //                             examine the array iwork which contains
122 //                             the error flags on the cycles, in order to
123 //                             look for eventual local integration
124 //                             difficulties. if the position of a local
125 //                             difficulty can be determined (e.g.
126 //                             singularity, discontinuity within the
127 //                             interval) one will probably gain from
128 //                             splitting up the interval at this point
129 //                             and calling appropriate integrators on
130 //                             the subranges.
131 //                         = 4 the extrapolation table constructed for
132 //                             convergence acceleration of the series
133 //                             formed by the integral contributions over
134 //                             the cycles, does not converge to within
135 //                             the requested accuracy. as in the case of
136 //                             ier = 1, it is advised to examine the
137 //                             array iwork which contains the error
138 //                             flags on the cycles.
139 //                         = 6 the input is invalid because
140 //                             (integr.ne.1 and integr.ne.2) or
141 //                              epsabs.le.0 or limlst.lt.3.
142 //                              result, abserr, neval, lst are set
143 //                              to zero.
144 //                         = 7 bad integrand behaviour occurs within one
145 //                             or more of the cycles. location and type
146 //                             of the difficulty involved can be
147 //                             determined from the vector ierlst. here
148 //                             lst is the number of cycles actually
149 //                             needed (see below).
150 //                             ierlst(k) = 1 the maximum number of
151 //                                           subdivisions (= limit) has
152 //                                           been achieved on the k th
153 //                                           cycle.
154 //                                       = 2 occurrence of roundoff error
155 //                                           is detected and prevents the
156 //                                           tolerance imposed on the
157 //                                           k th cycle, from being
158 //                                           achieved.
159 //                                       = 3 extremely bad integrand
160 //                                           behaviour occurs at some
161 //                                           points of the k th cycle.
162 //                                       = 4 the integration procedure
163 //                                           over the k th cycle does
164 //                                           not converge (to within the
165 //                                           required accuracy) due to
166 //                                           roundoff in the
167 //                                           extrapolation procedure
168 //                                           invoked on this cycle. it
169 //                                           is assumed that the result
170 //                                           on this interval is the
171 //                                           best which can be obtained.
172 //                                       = 5 the integral over the k th
173 //                                           cycle is probably divergent
174 //                                           or slowly convergent. it
175 //                                           must be noted that
176 //                                           divergence can occur with
177 //                                           any other value of
178 //                                           ierlst(k).
179 //                    if omega = 0 and integr = 1,
180 //                    the integral is calculated by means of dqagie
181 //                    and ier = ierlst(1) (with meaning as described
182 //                    for ierlst(k), k = 1).
183 //
184 //            rslst  - double precision
185 //                     vector of dimension at least limlst
186 //                     rslst(k) contains the integral contribution
187 //                     over the interval (a+(k-1)c,a+kc) where
188 //                     c = (2*int(abs(omega))+1)*pi/abs(omega),
189 //                     k = 1, 2, ..., lst.
190 //                     note that, if omega = 0, rslst(1) contains
191 //                     the value of the integral over (a,infinity).
192 //
193 //            erlst  - double precision
194 //                     vector of dimension at least limlst
195 //                     erlst(k) contains the error estimate corresponding
196 //                     with rslst(k).
197 //
198 //            ierlst - integer
199 //                     vector of dimension at least limlst
200 //                     ierlst(k) contains the error flag corresponding
201 //                     with rslst(k). for the meaning of the local error
202 //                     flags see description of output parameter ier.
203 //
204 //            lst    - integer
205 //                     number of subintervals needed for the integration
206 //                     if omega = 0 then lst is set to 1.
207 //
208 //            alist, blist, rlist, elist - double precision
209 //                     vector of dimension at least limit,
210 //
211 //            iord, nnlog - integer
212 //                     vector of dimension at least limit, providing
213 //                     space for the quantities needed in the subdivision
214 //                     process of each cycle
215 //
216 //            chebmo - double precision
217 //                     array of dimension at least (maxp1,25), providing
218 //                     space for the chebyshev moments needed within the
219 //                     cycles
220 //
221 //***references  (none)
222 //***routines called  d1mach,dqagie,dqawoe,dqelg
223 //***end prologue  dqawfe
224 //
225       Real abseps,correc,cycle,
226         c1,c2,dl,dla,drl,ep,eps,epsa,
227         errsum,fact,p,p1,reseps,
228         uflow;
229       Real[52] psum_;
230       Real[3] res3la_;
231       int ktmin,l,last,ll,
232           momcom,nev,nres,numrl2;
233 //
234       auto alist = dimension(alist_, limit);
235       auto blist = dimension(blist_, limit);
236       auto chebmo = dimension(chebmo_, maxp1, 25);
237       auto elist = dimension(elist_, limit);
238       auto erlst = dimension(erlst_, limlst);
239       auto ierlst = dimension(ierlst_, limlst);
240       auto iord = dimension(iord_, limit);
241       auto nnlog = dimension(nnlog_, limit);
242       auto psum = dimension(psum_.ptr, 52);
243       auto res3la = dimension(res3la_.ptr, 3);
244       auto rlist = dimension(rlist_, limit);
245       auto rslst = dimension(rslst_, limlst);
246 //
247 //
248 //            the dimension of  psum  is determined by the value of
249 //            limexp in subroutine dqelg (psum must be of dimension
250 //            (limexp+2) at least).
251 //
252 //           list of major variables
253 //           -----------------------
254 //
255 //           c1, c2    - end points of subinterval (of length cycle)
256 //           cycle     - (2*int(abs(omega))+1)*pi/abs(omega)
257 //           psum      - vector of dimension at least (limexp+2)
258 //                       (see routine dqelg)
259 //                       psum contains the part of the epsilon table
260 //                       which is still needed for further computations.
261 //                       each element of psum is a partial sum of the
262 //                       series which should sum to the value of the
263 //                       integral.
264 //           errsum    - sum of error estimates over the subintervals,
265 //                       calculated cumulatively
266 //           epsa      - absolute tolerance requested over current
267 //                       subinterval
268 //           chebmo    - array containing the modified chebyshev
269 //                       moments (see also routine dqc25f)
270 //
271       p = 0.9;
272 //
273 //           test on validity of parameters
274 //           ------------------------------
275 //
276 //***first executable statement  dqawfe
277       result = 0.0;
278       abserr = 0.0;
279       neval = 0;
280       lst = 0;
281       ier = 0;
282       if((integr != 1 && integr != 2) || epsabs <= 0.0 || 
283         limlst < 3) ier = 6;
284       if(ier == 6) goto l999;
285       if(omega != 0.0) goto l10;
286 //
287 //           integration by dqagie if omega is zero
288 //           --------------------------------------
289 //
290       if(integr == 1) qagie!(Real,Func)(f,0.0,1,epsabs,0.0,limit,
291         result,abserr,neval,ier,alist.ptr,blist.ptr,rlist.ptr,elist.ptr,iord.ptr,last);
292       rslst[1] = result;
293       erlst[1] = abserr;
294       ierlst[1] = ier;
295       lst = 1;
296       goto l999;
297 //
298 //           initializations
299 //           ---------------
300 //
301  l10: l = to!int(fabs(omega));
302       dl = 2*l+1;
303       cycle = dl*PI/fabs(omega);
304       ier = 0;
305       ktmin = 0;
306       neval = 0;
307       numrl2 = 0;
308       nres = 0;
309       c1 = a;
310       c2 = cycle+a;
311       p1 = 0.1e1-p;
312       uflow = Real.min_normal;
313       eps = epsabs;
314       if(epsabs > uflow/p1) eps = epsabs*p1;
315       ep = eps;
316       fact = 0.1e1;
317       correc = 0.0;
318       abserr = 0.0;
319       errsum = 0.0;
320 //
321 //           main do-loop
322 //           ------------
323 //
324       for (lst=1; lst<=limlst; lst++) { //do 50 lst = 1,limlst
325 //
326 //           integrate over current subinterval.
327 //
328         dla = lst;
329         epsa = eps*fact;
330         qawoe!(Real,Func)(f,c1,c2,omega,integr,epsa,0.0,limit,lst,maxp1,
331           rslst[lst],erlst[lst],nev,ierlst[lst],last,alist.ptr,blist.ptr,rlist.ptr,
332           elist.ptr,iord.ptr,nnlog.ptr,momcom,chebmo.ptr);
333         neval = neval+nev;
334         fact = fact*p;
335         errsum = errsum+erlst[lst];
336         drl = 0.5e2*fabs(rslst[lst]);
337 //
338 //           test on accuracy with partial sum
339 //
340         if((errsum+drl) <= epsabs && lst >= 6) goto l80;
341         correc = max(correc,erlst[lst]);
342         if(ierlst[lst] != 0) eps = max(ep,correc*p1);
343         if(ierlst[lst] != 0) ier = 7;
344         if(ier == 7 && (errsum+drl) <= correc*0.1e2 && 
345           lst > 5) goto l80;
346         numrl2 = numrl2+1;
347         if(lst > 1) goto l20;
348         psum[1] = rslst[1];
349         goto l40;
350  l20:   psum[numrl2] = psum[ll]+rslst[lst];
351         if(lst == 2) goto l40;
352 //
353 //           test on maximum number of subintervals
354 //
355         if(lst == limlst) ier = 1;
356 //
357 //           perform new extrapolation
358 //
359         qelg!Real(numrl2,psum.ptr,reseps,abseps,res3la.ptr,nres);
360 //
361 //           test whether extrapolated result is influenced by roundoff
362 //
363         ktmin = ktmin+1;
364         if(ktmin >= 15 && abserr <= 0.1e-2*(errsum+drl)) ier = 4;
365         if(abseps > abserr && lst != 3) goto l30;
366         abserr = abseps;
367         result = reseps;
368         ktmin = 0;
369 //
370 //           if ier is not 0, check whether direct result (partial sum)
371 //           or extrapolated result yields the best integral
372 //           approximation
373 //
374         if((abserr+0.1e2*correc) <= epsabs || 
375           (abserr <= epsabs && 0.1e2*correc >= epsabs)) goto l60;
376  l30:   if(ier != 0 && ier != 7) goto l60;
377  l40:   ll = numrl2;
378         c1 = c2;
379         c2 = c2+cycle;
380  l50: ;}
381 //
382 //         set final result and error estimate
383 //         -----------------------------------
384 //
385  l60: abserr = abserr+0.1e2*correc;
386       if(ier == 0) goto l999;
387       if(result != 0.0 && psum[numrl2] != 0.0) goto l70;
388       if(abserr > errsum) goto l80;
389       if(psum[numrl2] == 0.0) goto l999;
390  l70: if(abserr/fabs(result) > (errsum+drl)/fabs(psum[numrl2]))
391         goto l80;
392       if(ier >= 1 && ier != 7) abserr = abserr+drl;
393       goto l999;
394  l80: result = psum[numrl2];
395       abserr = errsum+drl;
396 l999: return;
397 }
398 
399 
400 unittest
401 {
402     alias qawfe!(float, float delegate(float)) fqawfe;
403     alias qawfe!(double, double delegate(double)) dqawfe;
404     alias qawfe!(double, double function(double)) dfqawfe;
405     alias qawfe!(real, real delegate(real)) rqawfe;
406 }