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.qawc; 9 10 11 import std.conv; 12 import scid.ports.quadpack.qawce; 13 14 version(unittest) 15 { 16 import std.math; 17 import scid.core.testing; 18 } 19 20 21 22 23 /// 24 void qawc(Real, Func)(Func f, Real a, Real b, Real c, Real epsabs, Real epsrel, 25 out Real result, out Real abserr, out int neval, out int ier, int limit, 26 int lenw, out int last, int* iwork, Real* work) 27 { 28 //***begin prologue dqawc 29 //***date written 800101 (yymmdd) 30 //***revision date 830518 (yymmdd) 31 //***category no. h2a2a1,j4 32 //***keywords automatic integrator, special-purpose, 33 // cauchy principal value, 34 // clenshaw-curtis, globally adaptive 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 38 // cauchy principal value i = integral of f*w over (a,b) 39 // (w(x) = 1/((x-c), c.ne.a, c.ne.b), hopefully satisfying 40 // following claim for accuracy 41 // abs(i-result).le.max(epsabe,epsrel*abs(i)). 42 //***description 43 // 44 // computation of a cauchy principal value 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 // under limit of integration 58 // 59 // b - double precision 60 // upper limit of integration 61 // 62 // c - parameter in the weight function, c.ne.a, c.ne.b. 63 // if c = a or c = b, the routine will end with 64 // ier = 6 . 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 // on return 75 // result - double precision 76 // approximation to the integral 77 // 78 // abserr - double precision 79 // estimate or the modulus of the absolute error, 80 // which should equal or exceed abs(i-result) 81 // 82 // neval - integer 83 // number of integrand evaluations 84 // 85 // ier - integer 86 // ier = 0 normal and reliable termination of the 87 // routine. it is assumed that the requested 88 // accuracy has been achieved. 89 // ier.gt.0 abnormal termination of the routine 90 // the estimates for integral and error are 91 // less reliable. it is assumed that the 92 // requested accuracy has not been achieved. 93 // error messages 94 // ier = 1 maximum number of subdivisions allowed 95 // has been achieved. one can allow more sub- 96 // divisions by increasing the value of limit 97 // (and taking the according dimension 98 // adjustments into account). however, if 99 // this yields no improvement it is advised 100 // to analyze the integrand in order to 101 // determine the integration difficulties. 102 // if the position of a local difficulty 103 // can be determined (e.g. singularity, 104 // discontinuity within the interval) one 105 // will probably gain from splitting up the 106 // interval at this point and calling 107 // appropriate integrators on the subranges. 108 // = 2 the occurrence of roundoff error is detec- 109 // ted, which prevents the requested 110 // tolerance from being achieved. 111 // = 3 extremely bad integrand behaviour occurs 112 // at some points of the integration 113 // interval. 114 // = 6 the input is invalid, because 115 // c = a or c = b or 116 // (epsabs.le.0 and 117 // epsrel.lt.max(50*rel.mach.acc.,0.5d-28)) 118 // or limit.lt.1 or lenw.lt.limit*4. 119 // result, abserr, neval, last are set to 120 // zero. exept when lenw or limit is invalid, 121 // iwork(1), work(limit*2+1) and 122 // work(limit*3+1) are set to zero, work(1) 123 // is set to a and work(limit+1) to b. 124 // 125 // dimensioning parameters 126 // limit - integer 127 // dimensioning parameter for iwork 128 // limit determines the maximum number of subintervals 129 // in the partition of the given integration interval 130 // (a,b), limit.ge.1. 131 // if limit.lt.1, the routine will end with ier = 6. 132 // 133 // lenw - integer 134 // dimensioning parameter for work 135 // lenw must be at least limit*4. 136 // if lenw.lt.limit*4, the routine will end with 137 // ier = 6. 138 // 139 // last - integer 140 // on return, last equals the number of subintervals 141 // produced in the subdivision process, which 142 // determines the number of significant elements 143 // actually in the work arrays. 144 // 145 // work arrays 146 // iwork - integer 147 // vector of dimension at least limit, the first k 148 // elements of which contain pointers 149 // to the error estimates over the subintervals, 150 // such that work(limit*3+iwork(1)), ... , 151 // work(limit*3+iwork(k)) form a decreasing 152 // sequence, with k = last if last.le.(limit/2+2), 153 // and k = limit+1-last otherwise 154 // 155 // work - double precision 156 // vector of dimension at least lenw 157 // on return 158 // work(1), ..., work(last) contain the left 159 // end points of the subintervals in the 160 // partition of (a,b), 161 // work(limit+1), ..., work(limit+last) contain 162 // the right end points, 163 // work(limit*2+1), ..., work(limit*2+last) contain 164 // the integral approximations over the subintervals, 165 // work(limit*3+1), ..., work(limit*3+last) 166 // contain the error estimates. 167 // 168 //***references (none) 169 //***routines called dqawce,xerror 170 //***end prologue dqawc 171 // 172 int l1,l2,l3; 173 // 174 // check validity of limit and lenw. 175 // 176 //***first executable statement dqawc 177 ier = 6; 178 neval = 0; 179 last = 0; 180 result = 0.0; 181 abserr = 0.0; 182 if(limit < 1 || lenw < limit*4) goto l10; 183 // 184 // prepare call for dqawce. 185 // 186 l1 = limit; 187 l2 = limit+l1; 188 l3 = limit+l2; 189 qawce!(Real, Func)(f,a,b,c,epsabs,epsrel,limit,result,abserr,neval,ier, 190 work,work+l1,work+l2,work+l3,iwork,last); 191 // 192 // call error handler if necessary. 193 // 194 l10: if(ier != 0) 195 throw new Exception("abnormal return from qawc: "~to!string(ier)); 196 return; 197 } 198 199 200 unittest 201 { 202 alias qawc!(float, float delegate(float)) fqawc; 203 alias qawc!(double, double delegate(double)) dqawc; 204 alias qawc!(double, double function(double)) dfqawc; 205 alias qawc!(real, real delegate(real)) rqawc; 206 } 207 208 209 unittest 210 { 211 // Integral 17 in the QUADPACK book. 212 enum : double 213 { 214 a = 0, 215 b = 5, 216 c = 2, 217 epsabs = 0, 218 epsrel = 1e-8, 219 } 220 double result, abserr; 221 int neval, ier, last; 222 enum limit = 500, lenw = limit*4; 223 int[limit] iwork; 224 double[lenw] work; 225 226 foreach (alpha; 0 .. 10) 227 { 228 double f(double x) 229 { 230 return 2.0^^(-alpha) / ((x-1)^^2 + 4.0^^(-alpha)); 231 } 232 233 qawc(&f, a, b, c, epsabs, epsrel, 234 result, abserr, neval, ier, limit, 235 lenw, last, iwork.ptr, work.ptr); 236 237 double exact = 238 ( 239 2.0^^(-alpha) * log(3.0/2) 240 - 2.0^^(-alpha-1) * log((16+4.0^^(-alpha))/(1+4.0^^(-alpha))) 241 - atan(2.0^^(alpha+2)) 242 - atan(2.0^^alpha) 243 ) 244 / (1 + 4.0^^(-alpha)); 245 246 assert (isAccurate(result, abserr, exact, epsrel, epsabs)); 247 } 248 }