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.qaws; 9 10 11 import std.conv; 12 import scid.ports.quadpack.qawse; 13 14 version(unittest) 15 { 16 import std.math; 17 import scid.core.testing; 18 } 19 20 21 22 23 /// 24 void qaws(Real, Func)(Func f, Real a, Real b, Real alfa, Real beta, 25 int integr, Real epsabs, Real epsrel, out Real result, out Real abserr, 26 out int neval, out int ier, int limit, int lenw, out int last, 27 int* iwork, Real* work) 28 { 29 //***begin prologue dqaws 30 //***date written 800101 (yymmdd) 31 //***revision date 830518 (yymmdd) 32 //***category no. h2a2a1 33 //***keywords automatic integrator, special-purpose, 34 // algebraico-logarithmic end-point singularities, 35 // clenshaw-curtis, 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 // definite integral i = integral of f*w over (a,b), 40 // (where w shows a singular behaviour at the end points 41 // see parameter integr). 42 // hopefully satisfying following claim for accuracy 43 // abs(i-result).le.max(epsabs,epsrel*abs(i)). 44 //***description 45 // 46 // integration of functions having algebraico-logarithmic 47 // end point singularities 48 // standard fortran subroutine 49 // double precision version 50 // 51 // parameters 52 // on entry 53 // f - double precision 54 // function subprogram defining the integrand 55 // function f(x). the actual name for f needs to be 56 // declared e x t e r n a l in the driver program. 57 // 58 // a - double precision 59 // lower limit of integration 60 // 61 // b - double precision 62 // upper limit of integration, b.gt.a 63 // if b.le.a, the routine will end with ier = 6. 64 // 65 // alfa - double precision 66 // parameter in the integrand function, alfa.gt.(-1) 67 // if alfa.le.(-1), the routine will end with 68 // ier = 6. 69 // 70 // beta - double precision 71 // parameter in the integrand function, beta.gt.(-1) 72 // if beta.le.(-1), the routine will end with 73 // ier = 6. 74 // 75 // integr - integer 76 // indicates which weight function is to be used 77 // = 1 (x-a)**alfa*(b-x)**beta 78 // = 2 (x-a)**alfa*(b-x)**beta*log(x-a) 79 // = 3 (x-a)**alfa*(b-x)**beta*log(b-x) 80 // = 4 (x-a)**alfa*(b-x)**beta*log(x-a)*log(b-x) 81 // if integr.lt.1 or integr.gt.4, the routine 82 // will end with ier = 6. 83 // 84 // epsabs - double precision 85 // absolute accuracy requested 86 // epsrel - double precision 87 // relative accuracy requested 88 // if epsabs.le.0 89 // and epsrel.lt.max(50*rel.mach.acc.,0.5d-28), 90 // the routine will end with ier = 6. 91 // 92 // on return 93 // result - double precision 94 // approximation to the integral 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 - integer 104 // ier = 0 normal and reliable termination of the 105 // routine. it is assumed that the requested 106 // accuracy has been achieved. 107 // ier.gt.0 abnormal termination of the routine 108 // the estimates for the integral and error 109 // are less reliable. it is assumed that the 110 // requested accuracy has not been achieved. 111 // error messages 112 // ier = 1 maximum number of subdivisions allowed 113 // has been achieved. one can allow more 114 // subdivisions by increasing the value of 115 // limit (and taking the according dimension 116 // adjustments into account). however, if 117 // this yields no improvement it is advised 118 // to analyze the integrand, in order to 119 // determine the integration difficulties 120 // which prevent the requested tolerance from 121 // being achieved. in case of a jump 122 // discontinuity or a local singularity 123 // of algebraico-logarithmic type at one or 124 // more interior points of the integration 125 // range, one should proceed by splitting up 126 // the interval at these points and calling 127 // the integrator on the subranges. 128 // = 2 the occurrence of roundoff error is 129 // detected, which prevents the requested 130 // tolerance from being achieved. 131 // = 3 extremely bad integrand behaviour occurs 132 // at some points of the integration 133 // interval. 134 // = 6 the input is invalid, because 135 // b.le.a or alfa.le.(-1) or beta.le.(-1) or 136 // or integr.lt.1 or integr.gt.4 or 137 // (epsabs.le.0 and 138 // epsrel.lt.max(50*rel.mach.acc.,0.5d-28)) 139 // or limit.lt.2 or lenw.lt.limit*4. 140 // result, abserr, neval, last are set to 141 // zero. except when lenw or limit is invalid 142 // iwork(1), work(limit*2+1) and 143 // work(limit*3+1) are set to zero, work(1) 144 // is set to a and work(limit+1) to b. 145 // 146 // dimensioning parameters 147 // limit - integer 148 // dimensioning parameter for iwork 149 // limit determines the maximum number of 150 // subintervals in the partition of the given 151 // integration interval (a,b), limit.ge.2. 152 // if limit.lt.2, the routine will end with ier = 6. 153 // 154 // lenw - integer 155 // dimensioning parameter for work 156 // lenw must be at least limit*4. 157 // if lenw.lt.limit*4, the routine will end 158 // with ier = 6. 159 // 160 // last - integer 161 // on return, last equals the number of 162 // subintervals produced in the subdivision process, 163 // which determines the significant number of 164 // elements actually in the work arrays. 165 // 166 // work arrays 167 // iwork - integer 168 // vector of dimension limit, the first k 169 // elements of which contain pointers 170 // to the error estimates over the subintervals, 171 // such that work(limit*3+iwork(1)), ..., 172 // work(limit*3+iwork(k)) form a decreasing 173 // sequence with k = last if last.le.(limit/2+2), 174 // and k = limit+1-last otherwise 175 // 176 // work - double precision 177 // vector of dimension lenw 178 // on return 179 // work(1), ..., work(last) contain the left 180 // end points of the subintervals in the 181 // partition of (a,b), 182 // work(limit+1), ..., work(limit+last) contain 183 // the right end points, 184 // work(limit*2+1), ..., work(limit*2+last) 185 // contain the integral approximations over 186 // the subintervals, 187 // work(limit*3+1), ..., work(limit*3+last) 188 // contain the error estimates. 189 // 190 //***references (none) 191 //***routines called dqawse,xerror 192 //***end prologue dqaws 193 // 194 int l1,l2,l3; 195 // 196 // check validity of limit and lenw. 197 // 198 //***first executable statement dqaws 199 ier = 6; 200 neval = 0; 201 last = 0; 202 result = 0.0; 203 abserr = 0.0; 204 if(limit < 2 || lenw < limit*4) goto l10; 205 // 206 // prepare call for dqawse. 207 // 208 l1 = limit; 209 l2 = limit+l1; 210 l3 = limit+l2; 211 // 212 qawse!(Real,Func)(f,a,b,alfa,beta,integr,epsabs,epsrel,limit,result, 213 abserr,neval,ier,work,work+l1,work+l2,work+l3,iwork,last); 214 // 215 // call error handler if necessary. 216 // 217 l10: if(ier != 0) 218 throw new Exception("abnormal return from qaws: "~to!string(ier)); 219 return; 220 } 221 222 unittest 223 { 224 alias qaws!(float, float delegate(float)) fqaws; 225 alias qaws!(double, double delegate(double)) dqaws; 226 alias qaws!(double, double function(double)) dfqaws; 227 alias qaws!(real, real delegate(real)) rqaws; 228 } 229 230 unittest 231 { 232 // From QUADPACK book 233 double f(double x) 234 { 235 if (x <= 0.0) return 0.0; 236 double logx2 = log(x)^^2; 237 return (1+logx2)^^(-2.0); 238 } 239 240 double a = 0.0, b = 1.0, alfa = 0.0, beta = 0.0; 241 int integr = 2, neval, ier, last; 242 double epsabs = 0.0, epsrel = 1e-8; 243 double result, abserr; 244 enum 245 { 246 limit = 500, 247 lenw = 4*limit 248 } 249 250 int[limit] iwork; 251 double[lenw] work; 252 253 qaws(&f, a, b, alfa, beta, integr, epsabs, epsrel, result, abserr, 254 neval, ier, limit, lenw, last, iwork.ptr, work.ptr); 255 256 double ans = -0.18927518788; 257 assert (isAccurate(result, abserr, ans, epsrel, epsabs)); 258 }