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.qag; 9 10 11 import std.conv; 12 import scid.ports.quadpack.qage; 13 14 version(unittest) 15 { 16 import std.math; 17 import scid.core.testing; 18 } 19 20 21 22 23 /// 24 void qag(Real, Func)(Func f, Real a, Real b, Real epsabs, Real epsrel, int key, 25 out Real result, out Real abserr, out int neval, out int ier, 26 int limit, int lenw, out int last, int* iwork, Real* work) 27 { 28 //***begin prologue dqag 29 //***date written 800101 (yymmdd) 30 //***revision date 830518 (yymmdd) 31 //***category no. h2a1a1 32 //***keywords automatic integrator, general-purpose, 33 // integrand examinator, globally adaptive, 34 // gauss-kronrod 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 // definite integral i = integral of f over (a,b), 39 // hopefully satisfying following claim for accuracy 40 // abs(i-result)le.max(epsabs,epsrel*abs(i)). 41 //***description 42 // 43 // computation of a definite integral 44 // standard fortran subroutine 45 // double precision version 46 // 47 // f - double precision 48 // function subprogam defining the integrand 49 // function f(x). the actual name for f needs to be 50 // declared e x t e r n a l in the driver program. 51 // 52 // a - double precision 53 // lower limit of integration 54 // 55 // b - double precision 56 // upper limit of integration 57 // 58 // epsabs - double precision 59 // absolute accoracy requested 60 // epsrel - double precision 61 // relative accuracy requested 62 // if epsabs.le.0 63 // and epsrel.lt.max(50*rel.mach.acc.,0.5d-28), 64 // the routine will end with ier = 6. 65 // 66 // key - integer 67 // key for choice of local integration rule 68 // a gauss-kronrod pair is used with 69 // 7 - 15 points if key.lt.2, 70 // 10 - 21 points if key = 2, 71 // 15 - 31 points if key = 3, 72 // 20 - 41 points if key = 4, 73 // 25 - 51 points if key = 5, 74 // 30 - 61 points if key.gt.5. 75 // 76 // on return 77 // result - double precision 78 // approximation to the integral 79 // 80 // abserr - double precision 81 // estimate of the modulus of the absolute error, 82 // which should equal or exceed abs(i-result) 83 // 84 // neval - integer 85 // number of integrand evaluations 86 // 87 // ier - integer 88 // ier = 0 normal and reliable termination of the 89 // routine. it is assumed that the requested 90 // accuracy has been achieved. 91 // ier.gt.0 abnormal termination of the routine 92 // the estimates for result and error are 93 // less reliable. it is assumed that the 94 // requested accuracy has not been achieved. 95 // error messages 96 // ier = 1 maximum number of subdivisions allowed 97 // has been achieved. one can allow more 98 // subdivisions by increasing the value of 99 // limit (and taking the according dimension 100 // adjustments into account). however, if 101 // this yield no improvement it is advised 102 // to analyze the integrand in order to 103 // determine the integration difficulaties. 104 // if the position of a local difficulty can 105 // be determined (i.e.singularity, 106 // discontinuity within the interval) one 107 // will probably gain from splitting up the 108 // interval at this point and calling the 109 // integrator on the subranges. if possible, 110 // an appropriate special-purpose integrator 111 // should be used which is designed for 112 // handling the type of difficulty involved. 113 // = 2 the occurrence of roundoff error is 114 // detected, which prevents the requested 115 // tolerance from being achieved. 116 // = 3 extremely bad integrand behaviour occurs 117 // at some points of the integration 118 // interval. 119 // = 6 the input is invalid, because 120 // (epsabs.le.0 and 121 // epsrel.lt.max(50*rel.mach.acc.,0.5d-28)) 122 // or limit.lt.1 or lenw.lt.limit*4. 123 // result, abserr, neval, last are set 124 // to zero. 125 // except when lenw is invalid, iwork(1), 126 // work(limit*2+1) and work(limit*3+1) are 127 // set to zero, work(1) is set to a and 128 // work(limit+1) to b. 129 // 130 // dimensioning parameters 131 // limit - integer 132 // dimensioning parameter for iwork 133 // limit determines the maximum number of subintervals 134 // in the partition of the given integration interval 135 // (a,b), limit.ge.1. 136 // if limit.lt.1, the routine will end with ier = 6. 137 // 138 // lenw - integer 139 // dimensioning parameter for work 140 // lenw must be at least limit*4. 141 // if lenw.lt.limit*4, the routine will end with 142 // ier = 6. 143 // 144 // last - integer 145 // on return, last equals the number of subintervals 146 // produced in the subdiviosion process, which 147 // determines the number of significant elements 148 // actually in the work arrays. 149 // 150 // work arrays 151 // iwork - integer 152 // vector of dimension at least limit, the first k 153 // elements of which contain pointers to the error 154 // estimates over the subintervals, such that 155 // work(limit*3+iwork(1)),... , work(limit*3+iwork(k)) 156 // form a decreasing sequence with k = last if 157 // last.le.(limit/2+2), and k = limit+1-last otherwise 158 // 159 // work - double precision 160 // vector of dimension at least lenw 161 // on return 162 // work(1), ..., work(last) contain the left end 163 // points of the subintervals in the partition of 164 // (a,b), 165 // work(limit+1), ..., work(limit+last) contain the 166 // right end points, 167 // work(limit*2+1), ..., work(limit*2+last) contain 168 // the integral approximations over the subintervals, 169 // work(limit*3+1), ..., work(limit*3+last) contain 170 // the error estimates. 171 // 172 //***references (none) 173 //***routines called dqage,xerror 174 //***end prologue dqag 175 int l1,l2,l3; 176 // 177 // check validity of lenw. 178 // 179 //***first executable statement dqag 180 ier = 6; 181 neval = 0; 182 last = 0; 183 result = 0.0; 184 abserr = 0.0; 185 if(limit < 1 || lenw < limit*4) goto l10; 186 // 187 // prepare call for dqage. 188 // 189 l1 = limit; 190 l2 = limit+l1; 191 l3 = limit+l2; 192 // 193 qage!(Real, Func)(f,a,b,epsabs,epsrel,key,limit,result,abserr,neval, 194 ier,work,work+l1,work+l2,work+l3,iwork,last); 195 // 196 // call error handler if necessary. 197 // 198 l10: 199 if(ier != 0) 200 throw new Exception("abnormal return from qag: "~to!string(ier)); 201 return; 202 } 203 204 unittest 205 { 206 alias qag!(float, float delegate(float)) fqag; 207 alias qag!(double, double delegate(double)) dqag; 208 alias qag!(double, double function(double)) dfqag; 209 alias qag!(real, real delegate(real)) rqag; 210 } 211 212 unittest 213 { 214 double f(double x) { return cos(100*sin(x)); } 215 216 enum : double 217 { 218 a = 0.0, 219 b = PI, 220 epsabs = 0.0 221 } 222 double epsrel=1e-10, result, abserr; 223 int key, neval, ier, last; 224 enum 225 { 226 limit = 500, 227 lenw = 4*limit 228 } 229 230 int[limit] iwork; 231 double[lenw] work; 232 233 // The answer is pi*J_0(100) 234 enum double ans = 0.062787400491492695655; 235 236 // Check that the routine works with all keys. 237 for (key=1; key<=6; key++) 238 { 239 qag(&f, a, b, epsabs, epsrel, key, result, abserr, neval, ier, 240 limit, lenw, last, iwork.ptr, work.ptr); 241 assert (isAccurate(result, abserr, ans, epsrel, epsabs)); 242 } 243 } 244 245 unittest 246 { 247 // Integral 1 in the QUADPACK book. 248 double alpha; 249 double f(double x) { return (x^^alpha)*log(1/x); } 250 251 double a=0.0, b=1.0, epsabs=0.0, epsrel=1e-8; 252 double result, abserr; 253 int neval, ier, last; 254 enum { limit=500, lenw=4*limit } 255 int[limit] iwork; 256 double[lenw] work; 257 258 for (alpha=0.0; alpha<=2.6+10*double.epsilon; 259 alpha+=0.2) 260 { 261 double ans = (alpha+1)^^(-2.0); 262 foreach (key; [1,3,6]) 263 { 264 qag(&f, a, b, epsabs, epsrel, key, result, abserr, neval, ier, 265 limit, lenw, last, iwork.ptr, work.ptr); 266 267 assert (isAccurate(result, abserr, ans, epsrel, epsabs)); 268 } 269 } 270 } 271 272 unittest 273 { 274 // Integral 15 in the QUADPACK book. 275 double alpha; 276 double f(double x) { return (x^^2.0) * exp(-(2.0^^(-alpha))*x); } 277 278 double a = 0.0, epsabs = 0.0, epsrel = 1e-8, result, abserr; 279 int key = 6, neval, ier, last; 280 enum { limit = 500, lenw = 4*limit } 281 282 int[limit] iwork; 283 double[lenw] work; 284 285 for (alpha=0.0; alpha<5.001; alpha+=1.0) 286 { 287 double b = 40*(2.0^^alpha); 288 qag (&f, a, b, epsabs, epsrel, key, result, abserr, neval, ier, 289 limit, lenw, last, iwork.ptr, work.ptr); 290 291 double eps = 1e-15; 292 double ans = (1+eps)*(2.0^^(3*alpha+1.0)); 293 assert (isAccurate(result, abserr, ans, epsrel, epsabs)); 294 } 295 }