00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045 #include "RooFit.h"
00046
00047 #include <assert.h>
00048 #include <stdlib.h>
00049 #include "TClass.h"
00050 #include "Riostream.h"
00051 #include "RooAdaptiveGaussKronrodIntegrator1D.h"
00052 #include "RooArgSet.h"
00053 #include "RooRealVar.h"
00054 #include "RooNumber.h"
00055 #include "RooNumIntFactory.h"
00056 #include "RooIntegratorBinding.h"
00057 #include "TMath.h"
00058 #include "RooMsgService.h"
00059
00060 using namespace std ;
00061
00062
00063 ClassImp(RooAdaptiveGaussKronrodIntegrator1D)
00064 ;
00065
00066
00067
00068 struct gsl_function_struct
00069 {
00070 double (* function) (double x, void * params);
00071 void * params;
00072 };
00073 typedef struct gsl_function_struct gsl_function ;
00074 #define GSL_FN_EVAL(F,x) (*((F)->function))(x,(F)->params)
00075
00076
00077 typedef struct
00078 {
00079 size_t limit;
00080 size_t size;
00081 size_t nrmax;
00082 size_t i;
00083 size_t maximum_level;
00084 double *alist;
00085 double *blist;
00086 double *rlist;
00087 double *elist;
00088 size_t *order;
00089 size_t *level;
00090 }
00091 gsl_integration_workspace;
00092
00093 gsl_integration_workspace *
00094 gsl_integration_workspace_alloc (const size_t n);
00095
00096 void
00097 gsl_integration_workspace_free (gsl_integration_workspace * w);
00098
00099 int gsl_integration_qag (const gsl_function * f,
00100 double a, double b,
00101 double epsabs, double epsrel, size_t limit,
00102 int key,
00103 gsl_integration_workspace * workspace,
00104 double *result, double *abserr);
00105
00106 int
00107 gsl_integration_qags (const gsl_function *f,
00108 double a, double b,
00109 double epsabs, double epsrel, size_t limit,
00110 gsl_integration_workspace * workspace,
00111 double * result, double * abserr) ;
00112
00113 int
00114 gsl_integration_qagi (gsl_function * f,
00115 double epsabs, double epsrel, size_t limit,
00116 gsl_integration_workspace * workspace,
00117 double *result, double *abserr) ;
00118
00119 int
00120 gsl_integration_qagil (gsl_function * f,
00121 double b,
00122 double epsabs, double epsrel, size_t limit,
00123 gsl_integration_workspace * workspace,
00124 double *result, double *abserr) ;
00125
00126 int
00127 gsl_integration_qagiu (gsl_function * f,
00128 double a,
00129 double epsabs, double epsrel, size_t limit,
00130 gsl_integration_workspace * workspace,
00131 double *result, double *abserr) ;
00132
00133
00134
00135
00136
00137
00138 void RooAdaptiveGaussKronrodIntegrator1D::registerIntegrator(RooNumIntFactory& fact)
00139 {
00140
00141
00142
00143 RooRealVar maxSeg("maxSeg","maximum number of segments",100) ;
00144 RooCategory method("method","Integration method for each segment") ;
00145 method.defineType("WynnEpsilon",0) ;
00146 method.defineType("15Points",1) ;
00147 method.defineType("21Points",2) ;
00148 method.defineType("31Points",3) ;
00149 method.defineType("41Points",4) ;
00150 method.defineType("51Points",5) ;
00151 method.defineType("61Points",6) ;
00152 method.setIndex(2) ;
00153 fact.storeProtoIntegrator(new RooAdaptiveGaussKronrodIntegrator1D(),RooArgSet(maxSeg,method)) ;
00154 }
00155
00156
00157
00158
00159 RooAdaptiveGaussKronrodIntegrator1D::RooAdaptiveGaussKronrodIntegrator1D() : _x(0), _workspace(0)
00160 {
00161
00162 }
00163
00164
00165
00166
00167
00168 RooAdaptiveGaussKronrodIntegrator1D::RooAdaptiveGaussKronrodIntegrator1D(const RooAbsFunc& function,
00169 const RooNumIntConfig& config) :
00170 RooAbsIntegrator(function),
00171 _epsAbs(config.epsRel()),
00172 _epsRel(config.epsAbs()),
00173 _workspace(0)
00174 {
00175
00176
00177
00178 const RooArgSet& confSet = config.getConfigSection(IsA()->GetName()) ;
00179 _maxSeg = (Int_t) confSet.getRealValue("maxSeg",100) ;
00180 _methodKey = confSet.getCatIndex("method",2) ;
00181
00182 _useIntegrandLimits= kTRUE;
00183 _valid= initialize();
00184 }
00185
00186
00187
00188
00189 RooAdaptiveGaussKronrodIntegrator1D::RooAdaptiveGaussKronrodIntegrator1D(const RooAbsFunc& function,
00190 Double_t xmin, Double_t xmax,
00191 const RooNumIntConfig& config) :
00192 RooAbsIntegrator(function),
00193 _epsAbs(config.epsRel()),
00194 _epsRel(config.epsAbs()),
00195 _workspace(0),
00196 _xmin(xmin),
00197 _xmax(xmax)
00198 {
00199
00200
00201
00202 const RooArgSet& confSet = config.getConfigSection(IsA()->GetName()) ;
00203 _maxSeg = (Int_t) confSet.getRealValue("maxSeg",100) ;
00204 _methodKey = confSet.getCatIndex("method",2) ;
00205
00206 _useIntegrandLimits= kFALSE;
00207 _valid= initialize();
00208 }
00209
00210
00211
00212
00213 RooAbsIntegrator* RooAdaptiveGaussKronrodIntegrator1D::clone(const RooAbsFunc& function, const RooNumIntConfig& config) const
00214 {
00215
00216 return new RooAdaptiveGaussKronrodIntegrator1D(function,config) ;
00217 }
00218
00219
00220
00221
00222 Bool_t RooAdaptiveGaussKronrodIntegrator1D::initialize()
00223 {
00224
00225
00226
00227 _x = new Double_t[_function->getDimension()] ;
00228 _workspace = gsl_integration_workspace_alloc (_maxSeg) ;
00229
00230 return checkLimits();
00231 }
00232
00233
00234
00235
00236 RooAdaptiveGaussKronrodIntegrator1D::~RooAdaptiveGaussKronrodIntegrator1D()
00237 {
00238
00239
00240 if (_workspace) {
00241 gsl_integration_workspace_free ((gsl_integration_workspace*) _workspace) ;
00242 }
00243 if (_x) {
00244 delete[] _x ;
00245 }
00246 }
00247
00248
00249
00250
00251 Bool_t RooAdaptiveGaussKronrodIntegrator1D::setLimits(Double_t* xmin, Double_t* xmax)
00252 {
00253
00254
00255
00256
00257 if(_useIntegrandLimits) {
00258 coutE(Integration) << "RooAdaptiveGaussKronrodIntegrator1D::setLimits: cannot override integrand's limits" << endl;
00259 return kFALSE;
00260 }
00261
00262 _xmin= *xmin;
00263 _xmax= *xmax;
00264 return checkLimits();
00265 }
00266
00267
00268
00269
00270 Bool_t RooAdaptiveGaussKronrodIntegrator1D::checkLimits() const
00271 {
00272
00273
00274
00275 if(_useIntegrandLimits) {
00276 assert(0 != integrand() && integrand()->isValid());
00277 _xmin= integrand()->getMinLimit(0);
00278 _xmax= integrand()->getMaxLimit(0);
00279 }
00280
00281
00282 Bool_t infLo= RooNumber::isInfinite(_xmin);
00283 Bool_t infHi= RooNumber::isInfinite(_xmax);
00284
00285 if (!infLo && !infHi) {
00286 _domainType = Closed ;
00287 } else if (infLo && infHi) {
00288 _domainType = Open ;
00289 } else if (infLo && !infHi) {
00290 _domainType = OpenLo ;
00291 } else {
00292 _domainType = OpenHi ;
00293 }
00294
00295
00296 return kTRUE ;
00297 }
00298
00299
00300
00301
00302 double RooAdaptiveGaussKronrodIntegrator1D_GSL_GlueFunction(double x, void *data)
00303 {
00304
00305 RooAdaptiveGaussKronrodIntegrator1D* instance = (RooAdaptiveGaussKronrodIntegrator1D*) data ;
00306 return instance->integrand(instance->xvec(x)) ;
00307 }
00308
00309
00310
00311
00312 Double_t RooAdaptiveGaussKronrodIntegrator1D::integral(const Double_t *yvec)
00313 {
00314
00315
00316 assert(isValid());
00317
00318
00319 if (yvec) {
00320 UInt_t i ; for (i=0 ; i<_function->getDimension()-1 ; i++) {
00321 _x[i+1] = yvec[i] ;
00322 }
00323 }
00324
00325
00326 gsl_function F;
00327 F.function = &RooAdaptiveGaussKronrodIntegrator1D_GSL_GlueFunction ;
00328 F.params = this ;
00329
00330
00331 double result, error;
00332
00333
00334 switch(_domainType) {
00335 case Closed:
00336 if (_methodKey==0) {
00337 gsl_integration_qags (&F, _xmin, _xmax, _epsAbs, _epsRel, _maxSeg, (gsl_integration_workspace*)_workspace,&result, &error);
00338 } else {
00339 gsl_integration_qag (&F, _xmin, _xmax, _epsAbs, _epsRel, _maxSeg, _methodKey, (gsl_integration_workspace*)_workspace,&result, &error);
00340 }
00341 break ;
00342 case OpenLo:
00343 gsl_integration_qagil (&F, _xmax, _epsAbs, _epsRel, _maxSeg, (gsl_integration_workspace*)_workspace,&result, &error);
00344 break ;
00345 case OpenHi:
00346 gsl_integration_qagiu (&F, _xmin, _epsAbs, _epsRel, _maxSeg, (gsl_integration_workspace*)_workspace,&result, &error);
00347 break ;
00348 case Open:
00349 gsl_integration_qagi (&F, _epsAbs, _epsRel, _maxSeg, (gsl_integration_workspace*)_workspace,&result, &error);
00350 break ;
00351 }
00352
00353 return result;
00354 }
00355
00356
00357
00358
00359
00360
00361
00362
00363
00364
00365
00366
00367
00368
00369
00370
00371
00372
00373
00374
00375
00376
00377
00378
00379
00380 #define GSL_SUCCESS 0
00381 #define GSL_EDOM 1
00382 #define GSL_ENOMEM 8
00383 #define GSL_EBADTOL 13
00384 #define GSL_ETOL 14
00385 #define GSL_ERROR(a,b) oocoutE((TObject*)0,Integration) << "RooAdaptiveGaussKronrodIntegrator1D::integral() ERROR: " << a << endl ; return b ;
00386 #define GSL_DBL_MIN 2.2250738585072014e-308
00387 #define GSL_DBL_MAX 1.7976931348623157e+308
00388 #define GSL_DBL_EPSILON 2.2204460492503131e-16
00389
00390 #define GSL_EINVAL 2
00391 #define GSL_EMAXITER 3
00392 #define GSL_ESING 4
00393 #define GSL_EFAILED 5
00394 #define GSL_EDIVERGE 6
00395 #define GSL_EROUND 7
00396
00397 #define GSL_ERROR_VAL(reason, gsl_errno, value) return value ;
00398
00399 #define GSL_MAX(a,b) ((a) > (b) ? (a) : (b))
00400 extern inline double
00401 GSL_MAX_DBL (double a, double b)
00402 {
00403 return GSL_MAX (a, b);
00404 }
00405
00406 double gsl_coerce_double (const double x);
00407
00408 double
00409 gsl_coerce_double (const double x)
00410 {
00411 volatile double y;
00412 y = x;
00413 return y;
00414 }
00415 #define GSL_COERCE_DBL(x) (gsl_coerce_double(x))
00416
00417
00418
00419
00420
00421
00422
00423
00424
00425
00426 typedef void gsl_integration_rule (const gsl_function * f,
00427 double a, double b,
00428 double *result, double *abserr,
00429 double *defabs, double *resabs);
00430
00431 void gsl_integration_qk15 (const gsl_function * f, double a, double b,
00432 double *result, double *abserr,
00433 double *resabs, double *resasc);
00434
00435 void gsl_integration_qk21 (const gsl_function * f, double a, double b,
00436 double *result, double *abserr,
00437 double *resabs, double *resasc);
00438
00439 void gsl_integration_qk31 (const gsl_function * f, double a, double b,
00440 double *result, double *abserr,
00441 double *resabs, double *resasc);
00442
00443 void gsl_integration_qk41 (const gsl_function * f, double a, double b,
00444 double *result, double *abserr,
00445 double *resabs, double *resasc);
00446
00447 void gsl_integration_qk51 (const gsl_function * f, double a, double b,
00448 double *result, double *abserr,
00449 double *resabs, double *resasc);
00450
00451 void gsl_integration_qk61 (const gsl_function * f, double a, double b,
00452 double *result, double *abserr,
00453 double *resabs, double *resasc);
00454
00455 void gsl_integration_qcheb (gsl_function * f, double a, double b,
00456 double *cheb12, double *cheb24);
00457
00458
00459
00460
00461 enum
00462 {
00463 GSL_INTEG_GAUSS15 = 1,
00464 GSL_INTEG_GAUSS21 = 2,
00465 GSL_INTEG_GAUSS31 = 3,
00466 GSL_INTEG_GAUSS41 = 4,
00467 GSL_INTEG_GAUSS51 = 5,
00468 GSL_INTEG_GAUSS61 = 6
00469 };
00470
00471
00472 void
00473 gsl_integration_qk (const int n, const double xgk[],
00474 const double wg[], const double wgk[],
00475 double fv1[], double fv2[],
00476 const gsl_function *f, double a, double b,
00477 double * result, double * abserr,
00478 double * resabs, double * resasc);
00479
00480
00481 int gsl_integration_qag (const gsl_function * f,
00482 double a, double b,
00483 double epsabs, double epsrel, size_t limit,
00484 int key,
00485 gsl_integration_workspace * workspace,
00486 double *result, double *abserr);
00487
00488
00489
00490 static inline
00491 void initialise (gsl_integration_workspace * workspace, double a, double b);
00492
00493 static inline
00494 void initialise (gsl_integration_workspace * workspace, double a, double b)
00495 {
00496 workspace->size = 0;
00497 workspace->nrmax = 0;
00498 workspace->i = 0;
00499 workspace->alist[0] = a;
00500 workspace->blist[0] = b;
00501 workspace->rlist[0] = 0.0;
00502 workspace->elist[0] = 0.0;
00503 workspace->order[0] = 0;
00504 workspace->level[0] = 0;
00505
00506 workspace->maximum_level = 0;
00507 }
00508
00509
00510 static inline
00511 void set_initial_result (gsl_integration_workspace * workspace,
00512 double result, double error);
00513
00514 static inline
00515 void set_initial_result (gsl_integration_workspace * workspace,
00516 double result, double error)
00517 {
00518 workspace->size = 1;
00519 workspace->rlist[0] = result;
00520 workspace->elist[0] = error;
00521 }
00522
00523
00524 static inline void
00525 qpsrt (gsl_integration_workspace * workspace);
00526
00527 static inline
00528 void qpsrt (gsl_integration_workspace * workspace)
00529 {
00530 const size_t last = workspace->size - 1;
00531 const size_t limit = workspace->limit;
00532
00533 double * elist = workspace->elist;
00534 size_t * order = workspace->order;
00535
00536 double errmax ;
00537 double errmin ;
00538 int i, k, top;
00539
00540 size_t i_nrmax = workspace->nrmax;
00541 size_t i_maxerr = order[i_nrmax] ;
00542
00543
00544
00545 if (last < 2)
00546 {
00547 order[0] = 0 ;
00548 order[1] = 1 ;
00549 workspace->i = i_maxerr ;
00550 return ;
00551 }
00552
00553 errmax = elist[i_maxerr] ;
00554
00555
00556
00557
00558
00559
00560 while (i_nrmax > 0 && errmax > elist[order[i_nrmax - 1]])
00561 {
00562 order[i_nrmax] = order[i_nrmax - 1] ;
00563 i_nrmax-- ;
00564 }
00565
00566
00567
00568
00569
00570 if(last < (limit/2 + 2))
00571 {
00572 top = last ;
00573 }
00574 else
00575 {
00576 top = limit - last + 1;
00577 }
00578
00579
00580
00581
00582 i = i_nrmax + 1 ;
00583
00584
00585
00586
00587 while (i < top && errmax < elist[order[i]])
00588 {
00589 order[i-1] = order[i] ;
00590 i++ ;
00591 }
00592
00593 order[i-1] = i_maxerr ;
00594
00595
00596
00597 errmin = elist[last] ;
00598
00599 k = top - 1 ;
00600
00601 while (k > i - 2 && errmin >= elist[order[k]])
00602 {
00603 order[k+1] = order[k] ;
00604 k-- ;
00605 }
00606
00607 order[k+1] = last ;
00608
00609
00610
00611 i_maxerr = order[i_nrmax] ;
00612
00613 workspace->i = i_maxerr ;
00614 workspace->nrmax = i_nrmax ;
00615 }
00616
00617
00618
00619
00620 static inline
00621 void update (gsl_integration_workspace * workspace,
00622 double a1, double b1, double area1, double error1,
00623 double a2, double b2, double area2, double error2);
00624
00625 static inline void
00626 retrieve (const gsl_integration_workspace * workspace,
00627 double * a, double * b, double * r, double * e);
00628
00629
00630
00631 static inline
00632 void update (gsl_integration_workspace * workspace,
00633 double a1, double b1, double area1, double error1,
00634 double a2, double b2, double area2, double error2)
00635 {
00636 double * alist = workspace->alist ;
00637 double * blist = workspace->blist ;
00638 double * rlist = workspace->rlist ;
00639 double * elist = workspace->elist ;
00640 size_t * level = workspace->level ;
00641
00642 const size_t i_max = workspace->i ;
00643 const size_t i_new = workspace->size ;
00644
00645 const size_t new_level = workspace->level[i_max] + 1;
00646
00647
00648
00649 if (error2 > error1)
00650 {
00651 alist[i_max] = a2;
00652 rlist[i_max] = area2;
00653 elist[i_max] = error2;
00654 level[i_max] = new_level;
00655
00656 alist[i_new] = a1;
00657 blist[i_new] = b1;
00658 rlist[i_new] = area1;
00659 elist[i_new] = error1;
00660 level[i_new] = new_level;
00661 }
00662 else
00663 {
00664 blist[i_max] = b1;
00665 rlist[i_max] = area1;
00666 elist[i_max] = error1;
00667 level[i_max] = new_level;
00668
00669 alist[i_new] = a2;
00670 blist[i_new] = b2;
00671 rlist[i_new] = area2;
00672 elist[i_new] = error2;
00673 level[i_new] = new_level;
00674 }
00675
00676 workspace->size++;
00677
00678 if (new_level > workspace->maximum_level)
00679 {
00680 workspace->maximum_level = new_level;
00681 }
00682
00683 qpsrt (workspace) ;
00684 }
00685
00686 static inline void
00687 retrieve (const gsl_integration_workspace * workspace,
00688 double * a, double * b, double * r, double * e)
00689 {
00690 const size_t i = workspace->i;
00691 double * alist = workspace->alist;
00692 double * blist = workspace->blist;
00693 double * rlist = workspace->rlist;
00694 double * elist = workspace->elist;
00695
00696 *a = alist[i] ;
00697 *b = blist[i] ;
00698 *r = rlist[i] ;
00699 *e = elist[i] ;
00700 }
00701
00702 static inline double
00703 sum_results (const gsl_integration_workspace * workspace);
00704
00705 static inline double
00706 sum_results (const gsl_integration_workspace * workspace)
00707 {
00708 const double * const rlist = workspace->rlist ;
00709 const size_t n = workspace->size;
00710
00711 size_t k;
00712 double result_sum = 0;
00713
00714 for (k = 0; k < n; k++)
00715 {
00716 result_sum += rlist[k];
00717 }
00718
00719 return result_sum;
00720 }
00721
00722 static inline int
00723 subinterval_too_small (double a1, double a2, double b2);
00724
00725 static inline int
00726 subinterval_too_small (double a1, double a2, double b2)
00727 {
00728 const double e = GSL_DBL_EPSILON;
00729 const double u = GSL_DBL_MIN;
00730
00731 double tmp = (1 + 100 * e) * (fabs (a2) + 1000 * u);
00732
00733 int status = fabs (a1) <= tmp && fabs (b2) <= tmp;
00734
00735 return status;
00736 }
00737
00738
00739 static int
00740 qag (const gsl_function *f,
00741 const double a, const double b,
00742 const double epsabs, const double epsrel,
00743 const size_t limit,
00744 gsl_integration_workspace * workspace,
00745 double * result, double * abserr,
00746 gsl_integration_rule * q) ;
00747
00748 int
00749 gsl_integration_qag (const gsl_function *f,
00750 double a, double b,
00751 double epsabs, double epsrel, size_t limit,
00752 int key,
00753 gsl_integration_workspace * workspace,
00754 double * result, double * abserr)
00755 {
00756 int status ;
00757 gsl_integration_rule * integration_rule = gsl_integration_qk15 ;
00758
00759 if (key < GSL_INTEG_GAUSS15)
00760 {
00761 key = GSL_INTEG_GAUSS15 ;
00762 }
00763 else if (key > GSL_INTEG_GAUSS61)
00764 {
00765 key = GSL_INTEG_GAUSS61 ;
00766 }
00767
00768 switch (key)
00769 {
00770 case GSL_INTEG_GAUSS15:
00771 integration_rule = gsl_integration_qk15 ;
00772 break ;
00773 case GSL_INTEG_GAUSS21:
00774 integration_rule = gsl_integration_qk21 ;
00775 break ;
00776 case GSL_INTEG_GAUSS31:
00777 integration_rule = gsl_integration_qk31 ;
00778 break ;
00779 case GSL_INTEG_GAUSS41:
00780 integration_rule = gsl_integration_qk41 ;
00781 break ;
00782 case GSL_INTEG_GAUSS51:
00783 integration_rule = gsl_integration_qk51 ;
00784 break ;
00785 case GSL_INTEG_GAUSS61:
00786 integration_rule = gsl_integration_qk61 ;
00787 break ;
00788 }
00789
00790 status = qag (f, a, b, epsabs, epsrel, limit,
00791 workspace,
00792 result, abserr,
00793 integration_rule) ;
00794
00795 return status ;
00796 }
00797
00798 static int
00799 qag (const gsl_function * f,
00800 const double a, const double b,
00801 const double epsabs, const double epsrel,
00802 const size_t limit,
00803 gsl_integration_workspace * workspace,
00804 double *result, double *abserr,
00805 gsl_integration_rule * q)
00806 {
00807 double area, errsum;
00808 double result0, abserr0, resabs0, resasc0;
00809 double tolerance;
00810 size_t iteration = 0;
00811 int roundoff_type1 = 0, roundoff_type2 = 0, error_type = 0;
00812
00813 double round_off;
00814
00815
00816
00817 initialise (workspace, a, b);
00818
00819 *result = 0;
00820 *abserr = 0;
00821
00822 if (limit > workspace->limit)
00823 {
00824 GSL_ERROR ("iteration limit exceeds available workspace", GSL_EINVAL) ;
00825 }
00826
00827 if (epsabs <= 0 && (epsrel < 50 * GSL_DBL_EPSILON || epsrel < 0.5e-28))
00828 {
00829 GSL_ERROR ("tolerance cannot be acheived with given epsabs and epsrel",
00830 GSL_EBADTOL);
00831 }
00832
00833
00834
00835 q (f, a, b, &result0, &abserr0, &resabs0, &resasc0);
00836
00837 set_initial_result (workspace, result0, abserr0);
00838
00839
00840
00841 tolerance = GSL_MAX_DBL (epsabs, epsrel * fabs (result0));
00842
00843
00844
00845 round_off = GSL_COERCE_DBL (50 * GSL_DBL_EPSILON * resabs0);
00846
00847 if (abserr0 <= round_off && abserr0 > tolerance)
00848 {
00849 *result = result0;
00850 *abserr = abserr0;
00851
00852 GSL_ERROR ("cannot reach tolerance because of roundoff error "
00853 "on first attempt", GSL_EROUND);
00854 }
00855 else if ((abserr0 <= tolerance && abserr0 != resasc0) || abserr0 == 0.0)
00856 {
00857 *result = result0;
00858 *abserr = abserr0;
00859
00860 return GSL_SUCCESS;
00861 }
00862 else if (limit == 1)
00863 {
00864 *result = result0;
00865 *abserr = abserr0;
00866
00867 GSL_ERROR ("a maximum of one iteration was insufficient", GSL_EMAXITER);
00868 }
00869
00870 area = result0;
00871 errsum = abserr0;
00872
00873 iteration = 1;
00874
00875 do
00876 {
00877 double a1, b1, a2, b2;
00878 double a_i, b_i, r_i, e_i;
00879 double area1 = 0, area2 = 0, area12 = 0;
00880 double error1 = 0, error2 = 0, error12 = 0;
00881 double resasc1, resasc2;
00882 double resabs1, resabs2;
00883
00884
00885
00886 retrieve (workspace, &a_i, &b_i, &r_i, &e_i);
00887
00888 a1 = a_i;
00889 b1 = 0.5 * (a_i + b_i);
00890 a2 = b1;
00891 b2 = b_i;
00892
00893 q (f, a1, b1, &area1, &error1, &resabs1, &resasc1);
00894 q (f, a2, b2, &area2, &error2, &resabs2, &resasc2);
00895
00896 area12 = area1 + area2;
00897 error12 = error1 + error2;
00898
00899 errsum += (error12 - e_i);
00900 area += area12 - r_i;
00901
00902 if (resasc1 != error1 && resasc2 != error2)
00903 {
00904 double delta = r_i - area12;
00905
00906 if (fabs (delta) <= 1.0e-5 * fabs (area12) && error12 >= 0.99 * e_i)
00907 {
00908 roundoff_type1++;
00909 }
00910 if (iteration >= 10 && error12 > e_i)
00911 {
00912 roundoff_type2++;
00913 }
00914 }
00915
00916 tolerance = GSL_MAX_DBL (epsabs, epsrel * fabs (area));
00917
00918 if (errsum > tolerance)
00919 {
00920 if (roundoff_type1 >= 6 || roundoff_type2 >= 20)
00921 {
00922 error_type = 2;
00923 }
00924
00925
00926
00927
00928 if (subinterval_too_small (a1, a2, b2))
00929 {
00930 error_type = 3;
00931 }
00932 }
00933
00934 update (workspace, a1, b1, area1, error1, a2, b2, area2, error2);
00935
00936 retrieve (workspace, &a_i, &b_i, &r_i, &e_i);
00937
00938 iteration++;
00939
00940 }
00941 while (iteration < limit && !error_type && errsum > tolerance);
00942
00943 *result = sum_results (workspace);
00944 *abserr = errsum;
00945
00946 if (errsum <= tolerance)
00947 {
00948 return GSL_SUCCESS;
00949 }
00950 else if (error_type == 2)
00951 {
00952 GSL_ERROR ("roundoff error prevents tolerance from being achieved",
00953 GSL_EROUND);
00954 }
00955 else if (error_type == 3)
00956 {
00957 GSL_ERROR ("bad integrand behavior found in the integration interval",
00958 GSL_ESING);
00959 }
00960 else if (iteration == limit)
00961 {
00962 GSL_ERROR ("maximum number of subdivisions reached", GSL_EMAXITER);
00963 }
00964
00965 GSL_ERROR ("could not integrate function", GSL_EFAILED);
00966 }
00967
00968
00969
00970 static double rescale_error (double err, const double result_abs, const double result_asc) ;
00971
00972 static double
00973 rescale_error (double err, const double result_abs, const double result_asc)
00974 {
00975 err = fabs(err) ;
00976
00977 if (result_asc != 0 && err != 0)
00978 {
00979 double scale = TMath::Power((200 * err / result_asc), 1.5) ;
00980
00981 if (scale < 1)
00982 {
00983 err = result_asc * scale ;
00984 }
00985 else
00986 {
00987 err = result_asc ;
00988 }
00989 }
00990 if (result_abs > GSL_DBL_MIN / (50 * GSL_DBL_EPSILON))
00991 {
00992 double min_err = 50 * GSL_DBL_EPSILON * result_abs ;
00993
00994 if (min_err > err)
00995 {
00996 err = min_err ;
00997 }
00998 }
00999
01000 return err ;
01001 }
01002
01003
01004 void
01005 gsl_integration_qk (const int n,
01006 const double xgk[], const double wg[], const double wgk[],
01007 double fv1[], double fv2[],
01008 const gsl_function * f, double a, double b,
01009 double *result, double *abserr,
01010 double *resabs, double *resasc)
01011 {
01012
01013 const double center = 0.5 * (a + b);
01014 const double half_length = 0.5 * (b - a);
01015 const double abs_half_length = fabs (half_length);
01016 const double f_center = GSL_FN_EVAL (f, center);
01017
01018 double result_gauss = 0;
01019 double result_kronrod = f_center * wgk[n - 1];
01020
01021 double result_abs = fabs (result_kronrod);
01022 double result_asc = 0;
01023 double mean = 0, err = 0;
01024
01025 int j;
01026
01027 if (n % 2 == 0)
01028 {
01029 result_gauss = f_center * wg[n / 2 - 1];
01030 }
01031
01032 for (j = 0; j < (n - 1) / 2; j++)
01033 {
01034 const int jtw = j * 2 + 1;
01035 const double abscissa = half_length * xgk[jtw];
01036 const double fval1 = GSL_FN_EVAL (f, center - abscissa);
01037 const double fval2 = GSL_FN_EVAL (f, center + abscissa);
01038 const double fsum = fval1 + fval2;
01039 fv1[jtw] = fval1;
01040 fv2[jtw] = fval2;
01041 result_gauss += wg[j] * fsum;
01042 result_kronrod += wgk[jtw] * fsum;
01043 result_abs += wgk[jtw] * (fabs (fval1) + fabs (fval2));
01044 }
01045
01046 for (j = 0; j < n / 2; j++)
01047 {
01048 int jtwm1 = j * 2;
01049 const double abscissa = half_length * xgk[jtwm1];
01050 const double fval1 = GSL_FN_EVAL (f, center - abscissa);
01051 const double fval2 = GSL_FN_EVAL (f, center + abscissa);
01052 fv1[jtwm1] = fval1;
01053 fv2[jtwm1] = fval2;
01054 result_kronrod += wgk[jtwm1] * (fval1 + fval2);
01055 result_abs += wgk[jtwm1] * (fabs (fval1) + fabs (fval2));
01056 };
01057
01058 mean = result_kronrod * 0.5;
01059
01060 result_asc = wgk[n - 1] * fabs (f_center - mean);
01061
01062 for (j = 0; j < n - 1; j++)
01063 {
01064 result_asc += wgk[j] * (fabs (fv1[j] - mean) + fabs (fv2[j] - mean));
01065 }
01066
01067
01068
01069 err = (result_kronrod - result_gauss) * half_length;
01070
01071 result_kronrod *= half_length;
01072 result_abs *= abs_half_length;
01073 result_asc *= abs_half_length;
01074
01075 *result = result_kronrod;
01076 *resabs = result_abs;
01077 *resasc = result_asc;
01078 *abserr = rescale_error (err, result_abs, result_asc);
01079
01080 }
01081
01082
01083
01084
01085
01086 static const double xgkA[8] =
01087 {
01088 0.991455371120812639206854697526329,
01089 0.949107912342758524526189684047851,
01090 0.864864423359769072789712788640926,
01091 0.741531185599394439863864773280788,
01092 0.586087235467691130294144838258730,
01093 0.405845151377397166906606412076961,
01094 0.207784955007898467600689403773245,
01095 0.000000000000000000000000000000000
01096 };
01097
01098
01099
01100
01101 static const double wgA[4] =
01102 {
01103 0.129484966168869693270611432679082,
01104 0.279705391489276667901467771423780,
01105 0.381830050505118944950369775488975,
01106 0.417959183673469387755102040816327
01107 };
01108
01109 static const double wgkA[8] =
01110 {
01111 0.022935322010529224963732008058970,
01112 0.063092092629978553290700663189204,
01113 0.104790010322250183839876322541518,
01114 0.140653259715525918745189590510238,
01115 0.169004726639267902826583426598550,
01116 0.190350578064785409913256402421014,
01117 0.204432940075298892414161999234649,
01118 0.209482141084727828012999174891714
01119 };
01120
01121 void
01122 gsl_integration_qk15 (const gsl_function * f, double a, double b,
01123 double *result, double *abserr,
01124 double *resabs, double *resasc)
01125 {
01126 double fv1[8], fv2[8];
01127 gsl_integration_qk (8, xgkA, wgA, wgkA, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
01128 }
01129
01130
01131
01132
01133
01134 static const double xgkB[11] =
01135 {
01136 0.995657163025808080735527280689003,
01137 0.973906528517171720077964012084452,
01138 0.930157491355708226001207180059508,
01139 0.865063366688984510732096688423493,
01140 0.780817726586416897063717578345042,
01141 0.679409568299024406234327365114874,
01142 0.562757134668604683339000099272694,
01143 0.433395394129247190799265943165784,
01144 0.294392862701460198131126603103866,
01145 0.148874338981631210884826001129720,
01146 0.000000000000000000000000000000000
01147 };
01148
01149
01150
01151
01152 static const double wgB[5] =
01153 {
01154 0.066671344308688137593568809893332,
01155 0.149451349150580593145776339657697,
01156 0.219086362515982043995534934228163,
01157 0.269266719309996355091226921569469,
01158 0.295524224714752870173892994651338
01159 };
01160
01161 static const double wgkB[11] =
01162 {
01163 0.011694638867371874278064396062192,
01164 0.032558162307964727478818972459390,
01165 0.054755896574351996031381300244580,
01166 0.075039674810919952767043140916190,
01167 0.093125454583697605535065465083366,
01168 0.109387158802297641899210590325805,
01169 0.123491976262065851077958109831074,
01170 0.134709217311473325928054001771707,
01171 0.142775938577060080797094273138717,
01172 0.147739104901338491374841515972068,
01173 0.149445554002916905664936468389821
01174 };
01175
01176
01177 void
01178 gsl_integration_qk21 (const gsl_function * f, double a, double b,
01179 double *result, double *abserr,
01180 double *resabs, double *resasc)
01181 {
01182 double fv1[11], fv2[11];
01183 gsl_integration_qk (11, xgkB, wgB, wgkB, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
01184 }
01185
01186
01187
01188
01189
01190 static const double xgkC[16] =
01191 {
01192 0.998002298693397060285172840152271,
01193 0.987992518020485428489565718586613,
01194 0.967739075679139134257347978784337,
01195 0.937273392400705904307758947710209,
01196 0.897264532344081900882509656454496,
01197 0.848206583410427216200648320774217,
01198 0.790418501442465932967649294817947,
01199 0.724417731360170047416186054613938,
01200 0.650996741297416970533735895313275,
01201 0.570972172608538847537226737253911,
01202 0.485081863640239680693655740232351,
01203 0.394151347077563369897207370981045,
01204 0.299180007153168812166780024266389,
01205 0.201194093997434522300628303394596,
01206 0.101142066918717499027074231447392,
01207 0.000000000000000000000000000000000
01208 };
01209
01210
01211
01212
01213 static const double wgC[8] =
01214 {
01215 0.030753241996117268354628393577204,
01216 0.070366047488108124709267416450667,
01217 0.107159220467171935011869546685869,
01218 0.139570677926154314447804794511028,
01219 0.166269205816993933553200860481209,
01220 0.186161000015562211026800561866423,
01221 0.198431485327111576456118326443839,
01222 0.202578241925561272880620199967519
01223 };
01224
01225 static const double wgkC[16] =
01226 {
01227 0.005377479872923348987792051430128,
01228 0.015007947329316122538374763075807,
01229 0.025460847326715320186874001019653,
01230 0.035346360791375846222037948478360,
01231 0.044589751324764876608227299373280,
01232 0.053481524690928087265343147239430,
01233 0.062009567800670640285139230960803,
01234 0.069854121318728258709520077099147,
01235 0.076849680757720378894432777482659,
01236 0.083080502823133021038289247286104,
01237 0.088564443056211770647275443693774,
01238 0.093126598170825321225486872747346,
01239 0.096642726983623678505179907627589,
01240 0.099173598721791959332393173484603,
01241 0.100769845523875595044946662617570,
01242 0.101330007014791549017374792767493
01243 };
01244
01245 void
01246 gsl_integration_qk31 (const gsl_function * f, double a, double b,
01247 double *result, double *abserr,
01248 double *resabs, double *resasc)
01249 {
01250 double fv1[16], fv2[16];
01251 gsl_integration_qk (16, xgkC, wgC, wgkC, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
01252 }
01253
01254
01255
01256
01257
01258 static const double xgkD[21] =
01259 {
01260 0.998859031588277663838315576545863,
01261 0.993128599185094924786122388471320,
01262 0.981507877450250259193342994720217,
01263 0.963971927277913791267666131197277,
01264 0.940822633831754753519982722212443,
01265 0.912234428251325905867752441203298,
01266 0.878276811252281976077442995113078,
01267 0.839116971822218823394529061701521,
01268 0.795041428837551198350638833272788,
01269 0.746331906460150792614305070355642,
01270 0.693237656334751384805490711845932,
01271 0.636053680726515025452836696226286,
01272 0.575140446819710315342946036586425,
01273 0.510867001950827098004364050955251,
01274 0.443593175238725103199992213492640,
01275 0.373706088715419560672548177024927,
01276 0.301627868114913004320555356858592,
01277 0.227785851141645078080496195368575,
01278 0.152605465240922675505220241022678,
01279 0.076526521133497333754640409398838,
01280 0.000000000000000000000000000000000
01281 };
01282
01283
01284
01285
01286 static const double wgD[11] =
01287 {
01288 0.017614007139152118311861962351853,
01289 0.040601429800386941331039952274932,
01290 0.062672048334109063569506535187042,
01291 0.083276741576704748724758143222046,
01292 0.101930119817240435036750135480350,
01293 0.118194531961518417312377377711382,
01294 0.131688638449176626898494499748163,
01295 0.142096109318382051329298325067165,
01296 0.149172986472603746787828737001969,
01297 0.152753387130725850698084331955098
01298 };
01299
01300 static const double wgkD[21] =
01301 {
01302 0.003073583718520531501218293246031,
01303 0.008600269855642942198661787950102,
01304 0.014626169256971252983787960308868,
01305 0.020388373461266523598010231432755,
01306 0.025882133604951158834505067096153,
01307 0.031287306777032798958543119323801,
01308 0.036600169758200798030557240707211,
01309 0.041668873327973686263788305936895,
01310 0.046434821867497674720231880926108,
01311 0.050944573923728691932707670050345,
01312 0.055195105348285994744832372419777,
01313 0.059111400880639572374967220648594,
01314 0.062653237554781168025870122174255,
01315 0.065834597133618422111563556969398,
01316 0.068648672928521619345623411885368,
01317 0.071054423553444068305790361723210,
01318 0.073030690332786667495189417658913,
01319 0.074582875400499188986581418362488,
01320 0.075704497684556674659542775376617,
01321 0.076377867672080736705502835038061,
01322 0.076600711917999656445049901530102
01323 };
01324
01325 void
01326 gsl_integration_qk41 (const gsl_function * f, double a, double b,
01327 double *result, double *abserr,
01328 double *resabs, double *resasc)
01329 {
01330 double fv1[21], fv2[21];
01331 gsl_integration_qk (21, xgkD, wgD, wgkD, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
01332 }
01333
01334
01335
01336
01337
01338 static const double xgkE[26] =
01339 {
01340 0.999262104992609834193457486540341,
01341 0.995556969790498097908784946893902,
01342 0.988035794534077247637331014577406,
01343 0.976663921459517511498315386479594,
01344 0.961614986425842512418130033660167,
01345 0.942974571228974339414011169658471,
01346 0.920747115281701561746346084546331,
01347 0.894991997878275368851042006782805,
01348 0.865847065293275595448996969588340,
01349 0.833442628760834001421021108693570,
01350 0.797873797998500059410410904994307,
01351 0.759259263037357630577282865204361,
01352 0.717766406813084388186654079773298,
01353 0.673566368473468364485120633247622,
01354 0.626810099010317412788122681624518,
01355 0.577662930241222967723689841612654,
01356 0.526325284334719182599623778158010,
01357 0.473002731445714960522182115009192,
01358 0.417885382193037748851814394594572,
01359 0.361172305809387837735821730127641,
01360 0.303089538931107830167478909980339,
01361 0.243866883720988432045190362797452,
01362 0.183718939421048892015969888759528,
01363 0.122864692610710396387359818808037,
01364 0.061544483005685078886546392366797,
01365 0.000000000000000000000000000000000
01366 };
01367
01368
01369
01370
01371 static const double wgE[13] =
01372 {
01373 0.011393798501026287947902964113235,
01374 0.026354986615032137261901815295299,
01375 0.040939156701306312655623487711646,
01376 0.054904695975835191925936891540473,
01377 0.068038333812356917207187185656708,
01378 0.080140700335001018013234959669111,
01379 0.091028261982963649811497220702892,
01380 0.100535949067050644202206890392686,
01381 0.108519624474263653116093957050117,
01382 0.114858259145711648339325545869556,
01383 0.119455763535784772228178126512901,
01384 0.122242442990310041688959518945852,
01385 0.123176053726715451203902873079050
01386 };
01387
01388 static const double wgkE[26] =
01389 {
01390 0.001987383892330315926507851882843,
01391 0.005561932135356713758040236901066,
01392 0.009473973386174151607207710523655,
01393 0.013236229195571674813656405846976,
01394 0.016847817709128298231516667536336,
01395 0.020435371145882835456568292235939,
01396 0.024009945606953216220092489164881,
01397 0.027475317587851737802948455517811,
01398 0.030792300167387488891109020215229,
01399 0.034002130274329337836748795229551,
01400 0.037116271483415543560330625367620,
01401 0.040083825504032382074839284467076,
01402 0.042872845020170049476895792439495,
01403 0.045502913049921788909870584752660,
01404 0.047982537138836713906392255756915,
01405 0.050277679080715671963325259433440,
01406 0.052362885806407475864366712137873,
01407 0.054251129888545490144543370459876,
01408 0.055950811220412317308240686382747,
01409 0.057437116361567832853582693939506,
01410 0.058689680022394207961974175856788,
01411 0.059720340324174059979099291932562,
01412 0.060539455376045862945360267517565,
01413 0.061128509717053048305859030416293,
01414 0.061471189871425316661544131965264,
01415 0.061580818067832935078759824240066
01416 };
01417
01418
01419
01420 void
01421 gsl_integration_qk51 (const gsl_function * f, double a, double b,
01422 double *result, double *abserr,
01423 double *resabs, double *resasc)
01424 {
01425 double fv1[26], fv2[26];
01426 gsl_integration_qk (26, xgkE, wgE, wgkE, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
01427 }
01428
01429
01430
01431
01432
01433 static const double xgkF[31] =
01434 {
01435 0.999484410050490637571325895705811,
01436 0.996893484074649540271630050918695,
01437 0.991630996870404594858628366109486,
01438 0.983668123279747209970032581605663,
01439 0.973116322501126268374693868423707,
01440 0.960021864968307512216871025581798,
01441 0.944374444748559979415831324037439,
01442 0.926200047429274325879324277080474,
01443 0.905573307699907798546522558925958,
01444 0.882560535792052681543116462530226,
01445 0.857205233546061098958658510658944,
01446 0.829565762382768397442898119732502,
01447 0.799727835821839083013668942322683,
01448 0.767777432104826194917977340974503,
01449 0.733790062453226804726171131369528,
01450 0.697850494793315796932292388026640,
01451 0.660061064126626961370053668149271,
01452 0.620526182989242861140477556431189,
01453 0.579345235826361691756024932172540,
01454 0.536624148142019899264169793311073,
01455 0.492480467861778574993693061207709,
01456 0.447033769538089176780609900322854,
01457 0.400401254830394392535476211542661,
01458 0.352704725530878113471037207089374,
01459 0.304073202273625077372677107199257,
01460 0.254636926167889846439805129817805,
01461 0.204525116682309891438957671002025,
01462 0.153869913608583546963794672743256,
01463 0.102806937966737030147096751318001,
01464 0.051471842555317695833025213166723,
01465 0.000000000000000000000000000000000
01466 };
01467
01468
01469
01470
01471 static const double wgF[15] =
01472 {
01473 0.007968192496166605615465883474674,
01474 0.018466468311090959142302131912047,
01475 0.028784707883323369349719179611292,
01476 0.038799192569627049596801936446348,
01477 0.048402672830594052902938140422808,
01478 0.057493156217619066481721689402056,
01479 0.065974229882180495128128515115962,
01480 0.073755974737705206268243850022191,
01481 0.080755895229420215354694938460530,
01482 0.086899787201082979802387530715126,
01483 0.092122522237786128717632707087619,
01484 0.096368737174644259639468626351810,
01485 0.099593420586795267062780282103569,
01486 0.101762389748405504596428952168554,
01487 0.102852652893558840341285636705415
01488 };
01489
01490 static const double wgkF[31] =
01491 {
01492 0.001389013698677007624551591226760,
01493 0.003890461127099884051267201844516,
01494 0.006630703915931292173319826369750,
01495 0.009273279659517763428441146892024,
01496 0.011823015253496341742232898853251,
01497 0.014369729507045804812451432443580,
01498 0.016920889189053272627572289420322,
01499 0.019414141193942381173408951050128,
01500 0.021828035821609192297167485738339,
01501 0.024191162078080601365686370725232,
01502 0.026509954882333101610601709335075,
01503 0.028754048765041292843978785354334,
01504 0.030907257562387762472884252943092,
01505 0.032981447057483726031814191016854,
01506 0.034979338028060024137499670731468,
01507 0.036882364651821229223911065617136,
01508 0.038678945624727592950348651532281,
01509 0.040374538951535959111995279752468,
01510 0.041969810215164246147147541285970,
01511 0.043452539701356069316831728117073,
01512 0.044814800133162663192355551616723,
01513 0.046059238271006988116271735559374,
01514 0.047185546569299153945261478181099,
01515 0.048185861757087129140779492298305,
01516 0.049055434555029778887528165367238,
01517 0.049795683427074206357811569379942,
01518 0.050405921402782346840893085653585,
01519 0.050881795898749606492297473049805,
01520 0.051221547849258772170656282604944,
01521 0.051426128537459025933862879215781,
01522 0.051494729429451567558340433647099
01523 };
01524
01525 void
01526 gsl_integration_qk61 (const gsl_function * f, double a, double b,
01527 double *result, double *abserr,
01528 double *resabs, double *resasc)
01529 {
01530 double fv1[31], fv2[31];
01531 gsl_integration_qk (31, xgkF, wgF, wgkF, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
01532 }
01533
01534 gsl_integration_workspace*
01535 gsl_integration_workspace_alloc (const size_t n)
01536 {
01537 gsl_integration_workspace* w ;
01538
01539 if (n == 0)
01540 {
01541 GSL_ERROR_VAL ("workspace length n must be positive integer",
01542 GSL_EDOM, 0);
01543 }
01544
01545 w = (gsl_integration_workspace *)
01546 malloc (sizeof (gsl_integration_workspace));
01547
01548 if (w == 0)
01549 {
01550 GSL_ERROR_VAL ("failed to allocate space for workspace struct",
01551 GSL_ENOMEM, 0);
01552 }
01553
01554 w->alist = (double *) malloc (n * sizeof (double));
01555
01556 if (w->alist == 0)
01557 {
01558 free (w);
01559
01560 GSL_ERROR_VAL ("failed to allocate space for alist ranges",
01561 GSL_ENOMEM, 0);
01562 }
01563
01564 w->blist = (double *) malloc (n * sizeof (double));
01565
01566 if (w->blist == 0)
01567 {
01568 free (w->alist);
01569 free (w);
01570
01571 GSL_ERROR_VAL ("failed to allocate space for blist ranges",
01572 GSL_ENOMEM, 0);
01573 }
01574
01575 w->rlist = (double *) malloc (n * sizeof (double));
01576
01577 if (w->rlist == 0)
01578 {
01579 free (w->blist);
01580 free (w->alist);
01581 free (w);
01582
01583 GSL_ERROR_VAL ("failed to allocate space for rlist ranges",
01584 GSL_ENOMEM, 0);
01585 }
01586
01587
01588 w->elist = (double *) malloc (n * sizeof (double));
01589
01590 if (w->elist == 0)
01591 {
01592 free (w->rlist);
01593 free (w->blist);
01594 free (w->alist);
01595 free (w);
01596
01597 GSL_ERROR_VAL ("failed to allocate space for elist ranges",
01598 GSL_ENOMEM, 0);
01599 }
01600
01601 w->order = (size_t *) malloc (n * sizeof (size_t));
01602
01603 if (w->order == 0)
01604 {
01605 free (w->elist);
01606 free (w->rlist);
01607 free (w->blist);
01608 free (w->alist);
01609 free (w);
01610
01611 GSL_ERROR_VAL ("failed to allocate space for order ranges",
01612 GSL_ENOMEM, 0);
01613 }
01614
01615 w->level = (size_t *) malloc (n * sizeof (size_t));
01616
01617 if (w->level == 0)
01618 {
01619 free (w->order);
01620 free (w->elist);
01621 free (w->rlist);
01622 free (w->blist);
01623 free (w->alist);
01624 free (w);
01625
01626 GSL_ERROR_VAL ("failed to allocate space for order ranges",
01627 GSL_ENOMEM, 0);
01628 }
01629
01630 w->size = 0 ;
01631 w->limit = n ;
01632 w->maximum_level = 0 ;
01633
01634 return w ;
01635 }
01636
01637 void
01638 gsl_integration_workspace_free (gsl_integration_workspace * w)
01639 {
01640 free (w->level) ;
01641 free (w->order) ;
01642 free (w->elist) ;
01643 free (w->rlist) ;
01644 free (w->blist) ;
01645 free (w->alist) ;
01646 free (w) ;
01647 }
01648
01649
01650
01651
01652 static inline void
01653 reset_nrmax (gsl_integration_workspace * workspace);
01654
01655 static inline void
01656 reset_nrmax (gsl_integration_workspace * workspace)
01657 {
01658 workspace->nrmax = 0;
01659 workspace->i = workspace->order[0] ;
01660 }
01661
01662
01663
01664
01665
01666
01667
01668 static int
01669 increase_nrmax (gsl_integration_workspace * workspace);
01670
01671 static int
01672 increase_nrmax (gsl_integration_workspace * workspace)
01673 {
01674 int k;
01675 int id = workspace->nrmax;
01676 int jupbnd;
01677
01678 const size_t * level = workspace->level;
01679 const size_t * order = workspace->order;
01680
01681 size_t limit = workspace->limit ;
01682 size_t last = workspace->size - 1 ;
01683
01684 if (last > (1 + limit / 2))
01685 {
01686 jupbnd = limit + 1 - last;
01687 }
01688 else
01689 {
01690 jupbnd = last;
01691 }
01692
01693 for (k = id; k <= jupbnd; k++)
01694 {
01695 size_t i_max = order[workspace->nrmax];
01696
01697 workspace->i = i_max ;
01698
01699 if (level[i_max] < workspace->maximum_level)
01700 {
01701 return 1;
01702 }
01703
01704 workspace->nrmax++;
01705
01706 }
01707 return 0;
01708 }
01709
01710 static int
01711 large_interval (gsl_integration_workspace * workspace)
01712 {
01713 size_t i = workspace->i ;
01714 const size_t * level = workspace->level;
01715
01716 if (level[i] < workspace->maximum_level)
01717 {
01718 return 1 ;
01719 }
01720 else
01721 {
01722 return 0 ;
01723 }
01724 }
01725
01726
01727
01728 struct extrapolation_table
01729 {
01730 size_t n;
01731 double rlist2[52];
01732 size_t nres;
01733 double res3la[3];
01734 };
01735
01736 static void
01737 initialise_table (struct extrapolation_table *table);
01738
01739 static void
01740 append_table (struct extrapolation_table *table, double y);
01741
01742 static void
01743 initialise_table (struct extrapolation_table *table)
01744 {
01745 table->n = 0;
01746 table->nres = 0;
01747 }
01748 #ifdef JUNK
01749 static void
01750 initialise_table (struct extrapolation_table *table, double y)
01751 {
01752 table->n = 0;
01753 table->rlist2[0] = y;
01754 table->nres = 0;
01755 }
01756 #endif
01757 static void
01758 append_table (struct extrapolation_table *table, double y)
01759 {
01760 size_t n;
01761 n = table->n;
01762 table->rlist2[n] = y;
01763 table->n++;
01764 }
01765
01766
01767
01768
01769
01770
01771 static inline void
01772 qelg (struct extrapolation_table *table, double *result, double *abserr);
01773
01774 static inline void
01775 qelg (struct extrapolation_table *table, double *result, double *abserr)
01776 {
01777 double *epstab = table->rlist2;
01778 double *res3la = table->res3la;
01779 const size_t n = table->n - 1;
01780
01781 const double current = epstab[n];
01782
01783 double absolute = GSL_DBL_MAX;
01784 double relative = 5 * GSL_DBL_EPSILON * fabs (current);
01785
01786 const size_t newelm = n / 2;
01787 const size_t n_orig = n;
01788 size_t n_final = n;
01789 size_t i;
01790
01791 const size_t nres_orig = table->nres;
01792
01793 *result = current;
01794 *abserr = GSL_DBL_MAX;
01795
01796 if (n < 2)
01797 {
01798 *result = current;
01799 *abserr = GSL_MAX_DBL (absolute, relative);
01800 return;
01801 }
01802
01803 epstab[n + 2] = epstab[n];
01804 epstab[n] = GSL_DBL_MAX;
01805
01806 for (i = 0; i < newelm; i++)
01807 {
01808 double res = epstab[n - 2 * i + 2];
01809 double e0 = epstab[n - 2 * i - 2];
01810 double e1 = epstab[n - 2 * i - 1];
01811 double e2 = res;
01812
01813 double e1abs = fabs (e1);
01814 double delta2 = e2 - e1;
01815 double err2 = fabs (delta2);
01816 double tol2 = GSL_MAX_DBL (fabs (e2), e1abs) * GSL_DBL_EPSILON;
01817 double delta3 = e1 - e0;
01818 double err3 = fabs (delta3);
01819 double tol3 = GSL_MAX_DBL (e1abs, fabs (e0)) * GSL_DBL_EPSILON;
01820
01821 double e3, delta1, err1, tol1, ss;
01822
01823 if (err2 <= tol2 && err3 <= tol3)
01824 {
01825
01826
01827
01828 *result = res;
01829 absolute = err2 + err3;
01830 relative = 5 * GSL_DBL_EPSILON * fabs (res);
01831 *abserr = GSL_MAX_DBL (absolute, relative);
01832 return;
01833 }
01834
01835 e3 = epstab[n - 2 * i];
01836 epstab[n - 2 * i] = e1;
01837 delta1 = e1 - e3;
01838 err1 = fabs (delta1);
01839 tol1 = GSL_MAX_DBL (e1abs, fabs (e3)) * GSL_DBL_EPSILON;
01840
01841
01842
01843
01844 if (err1 <= tol1 || err2 <= tol2 || err3 <= tol3)
01845 {
01846 n_final = 2 * i;
01847 break;
01848 }
01849
01850 ss = (1 / delta1 + 1 / delta2) - 1 / delta3;
01851
01852
01853
01854
01855
01856 if (fabs (ss * e1) <= 0.0001)
01857 {
01858 n_final = 2 * i;
01859 break;
01860 }
01861
01862
01863
01864
01865 res = e1 + 1 / ss;
01866 epstab[n - 2 * i] = res;
01867
01868 {
01869 const double error = err2 + fabs (res - e2) + err3;
01870
01871 if (error <= *abserr)
01872 {
01873 *abserr = error;
01874 *result = res;
01875 }
01876 }
01877 }
01878
01879
01880
01881 {
01882 const size_t limexp = 50 - 1;
01883
01884 if (n_final == limexp)
01885 {
01886 n_final = 2 * (limexp / 2);
01887 }
01888 }
01889
01890 if (n_orig % 2 == 1)
01891 {
01892 for (i = 0; i <= newelm; i++)
01893 {
01894 epstab[1 + i * 2] = epstab[i * 2 + 3];
01895 }
01896 }
01897 else
01898 {
01899 for (i = 0; i <= newelm; i++)
01900 {
01901 epstab[i * 2] = epstab[i * 2 + 2];
01902 }
01903 }
01904
01905 if (n_orig != n_final)
01906 {
01907 for (i = 0; i <= n_final; i++)
01908 {
01909 epstab[i] = epstab[n_orig - n_final + i];
01910 }
01911 }
01912
01913 table->n = n_final + 1;
01914
01915 if (nres_orig < 3)
01916 {
01917 res3la[nres_orig] = *result;
01918 *abserr = GSL_DBL_MAX;
01919 }
01920 else
01921 {
01922 *abserr = (fabs (*result - res3la[2]) + fabs (*result - res3la[1])
01923 + fabs (*result - res3la[0]));
01924
01925 res3la[0] = res3la[1];
01926 res3la[1] = res3la[2];
01927 res3la[2] = *result;
01928 }
01929
01930
01931
01932
01933
01934
01935
01936 table->nres = nres_orig + 1;
01937
01938 *abserr = GSL_MAX_DBL (*abserr, 5 * GSL_DBL_EPSILON * fabs (*result));
01939
01940 return;
01941 }
01942
01943
01944
01945
01946
01947
01948 static inline int
01949 test_positivity (double result, double resabs);
01950
01951 static inline int
01952 test_positivity (double result, double resabs)
01953 {
01954 int status = (fabs (result) >= (1 - 50 * GSL_DBL_EPSILON) * resabs);
01955
01956 return status;
01957 }
01958
01959 static int qags (const gsl_function * f, const double a, const double
01960 b, const double epsabs, const double epsrel, const size_t limit,
01961 gsl_integration_workspace * workspace, double *result, double *abserr,
01962 gsl_integration_rule * q);
01963
01964 int
01965 gsl_integration_qags (const gsl_function *f,
01966 double a, double b,
01967 double epsabs, double epsrel, size_t limit,
01968 gsl_integration_workspace * workspace,
01969 double * result, double * abserr)
01970 {
01971 int status = qags (f, a, b, epsabs, epsrel, limit,
01972 workspace,
01973 result, abserr,
01974 &gsl_integration_qk21) ;
01975 return status ;
01976 }
01977
01978
01979
01980
01981
01982
01983
01984
01985 static double i_transform (double t, void *params);
01986
01987 int
01988 gsl_integration_qagi (gsl_function * f,
01989 double epsabs, double epsrel, size_t limit,
01990 gsl_integration_workspace * workspace,
01991 double *result, double *abserr)
01992 {
01993 int status;
01994
01995 gsl_function f_transform;
01996
01997 f_transform.function = &i_transform;
01998 f_transform.params = f;
01999
02000 status = qags (&f_transform, 0.0, 1.0,
02001 epsabs, epsrel, limit,
02002 workspace,
02003 result, abserr,
02004 &gsl_integration_qk15);
02005
02006 return status;
02007 }
02008
02009 static double
02010 i_transform (double t, void *params)
02011 {
02012 gsl_function *f = (gsl_function *) params;
02013 double x = (1 - t) / t;
02014 double y = GSL_FN_EVAL (f, x) + GSL_FN_EVAL (f, -x);
02015 return (y / t) / t;
02016 }
02017
02018
02019
02020
02021
02022
02023
02024
02025
02026 struct il_params { double b ; gsl_function * f ; } ;
02027
02028 static double il_transform (double t, void *params);
02029
02030 int
02031 gsl_integration_qagil (gsl_function * f,
02032 double b,
02033 double epsabs, double epsrel, size_t limit,
02034 gsl_integration_workspace * workspace,
02035 double *result, double *abserr)
02036 {
02037 int status;
02038
02039 gsl_function f_transform;
02040 struct il_params transform_params ;
02041
02042 transform_params.b = b ;
02043 transform_params.f = f ;
02044
02045 f_transform.function = &il_transform;
02046 f_transform.params = &transform_params;
02047
02048 status = qags (&f_transform, 0.0, 1.0,
02049 epsabs, epsrel, limit,
02050 workspace,
02051 result, abserr,
02052 &gsl_integration_qk15);
02053
02054 return status;
02055 }
02056
02057 static double
02058 il_transform (double t, void *params)
02059 {
02060 struct il_params *p = (struct il_params *) params;
02061 double b = p->b;
02062 gsl_function * f = p->f;
02063 double x = b - (1 - t) / t;
02064 double y = GSL_FN_EVAL (f, x);
02065 return (y / t) / t;
02066 }
02067
02068
02069
02070
02071
02072
02073
02074
02075 struct iu_params { double a ; gsl_function * f ; } ;
02076
02077 static double iu_transform (double t, void *params);
02078
02079 int
02080 gsl_integration_qagiu (gsl_function * f,
02081 double a,
02082 double epsabs, double epsrel, size_t limit,
02083 gsl_integration_workspace * workspace,
02084 double *result, double *abserr)
02085 {
02086 int status;
02087
02088 gsl_function f_transform;
02089 struct iu_params transform_params ;
02090
02091 transform_params.a = a ;
02092 transform_params.f = f ;
02093
02094 f_transform.function = &iu_transform;
02095 f_transform.params = &transform_params;
02096
02097 status = qags (&f_transform, 0.0, 1.0,
02098 epsabs, epsrel, limit,
02099 workspace,
02100 result, abserr,
02101 &gsl_integration_qk15);
02102
02103 return status;
02104 }
02105
02106 static double
02107 iu_transform (double t, void *params)
02108 {
02109 struct iu_params *p = (struct iu_params *) params;
02110 double a = p->a;
02111 gsl_function * f = p->f;
02112 double x = a + (1 - t) / t;
02113 double y = GSL_FN_EVAL (f, x);
02114 return (y / t) / t;
02115 }
02116
02117
02118
02119 static int
02120 qags (const gsl_function * f,
02121 const double a, const double b,
02122 const double epsabs, const double epsrel,
02123 const size_t limit,
02124 gsl_integration_workspace * workspace,
02125 double *result, double *abserr,
02126 gsl_integration_rule * q)
02127 {
02128 double area, errsum;
02129 double res_ext, err_ext;
02130 double result0, abserr0, resabs0, resasc0;
02131 double tolerance;
02132
02133 double ertest = 0;
02134 double error_over_large_intervals = 0;
02135 double reseps = 0, abseps = 0, correc = 0;
02136 size_t ktmin = 0;
02137 int roundoff_type1 = 0, roundoff_type2 = 0, roundoff_type3 = 0;
02138 int error_type = 0, error_type2 = 0;
02139
02140 size_t iteration = 0;
02141
02142 int positive_integrand = 0;
02143 int extrapolate = 0;
02144 int disallow_extrapolation = 0;
02145
02146 struct extrapolation_table table;
02147
02148
02149
02150 initialise (workspace, a, b);
02151
02152 *result = 0;
02153 *abserr = 0;
02154
02155 if (limit > workspace->limit)
02156 {
02157 GSL_ERROR ("iteration limit exceeds available workspace", GSL_EINVAL) ;
02158 }
02159
02160
02161
02162 if (epsabs <= 0 && (epsrel < 50 * GSL_DBL_EPSILON || epsrel < 0.5e-28))
02163 {
02164 GSL_ERROR ("tolerance cannot be acheived with given epsabs and epsrel",
02165 GSL_EBADTOL);
02166 }
02167
02168
02169
02170 q (f, a, b, &result0, &abserr0, &resabs0, &resasc0);
02171
02172 set_initial_result (workspace, result0, abserr0);
02173
02174 tolerance = GSL_MAX_DBL (epsabs, epsrel * fabs (result0));
02175
02176 if (abserr0 <= 100 * GSL_DBL_EPSILON * resabs0 && abserr0 > tolerance)
02177 {
02178 *result = result0;
02179 *abserr = abserr0;
02180
02181 GSL_ERROR ("cannot reach tolerance because of roundoff error"
02182 "on first attempt", GSL_EROUND);
02183 }
02184 else if ((abserr0 <= tolerance && abserr0 != resasc0) || abserr0 == 0.0)
02185 {
02186 *result = result0;
02187 *abserr = abserr0;
02188
02189 return GSL_SUCCESS;
02190 }
02191 else if (limit == 1)
02192 {
02193 *result = result0;
02194 *abserr = abserr0;
02195
02196 GSL_ERROR ("a maximum of one iteration was insufficient", GSL_EMAXITER);
02197 }
02198
02199
02200
02201 initialise_table (&table);
02202 append_table (&table, result0);
02203
02204 area = result0;
02205 errsum = abserr0;
02206
02207 res_ext = result0;
02208 err_ext = GSL_DBL_MAX;
02209
02210 positive_integrand = test_positivity (result0, resabs0);
02211
02212 iteration = 1;
02213
02214 do
02215 {
02216 size_t current_level;
02217 double a1, b1, a2, b2;
02218 double a_i, b_i, r_i, e_i;
02219 double area1 = 0, area2 = 0, area12 = 0;
02220 double error1 = 0, error2 = 0, error12 = 0;
02221 double resasc1, resasc2;
02222 double resabs1, resabs2;
02223 double last_e_i;
02224
02225
02226
02227 retrieve (workspace, &a_i, &b_i, &r_i, &e_i);
02228
02229 current_level = workspace->level[workspace->i] + 1;
02230
02231 a1 = a_i;
02232 b1 = 0.5 * (a_i + b_i);
02233 a2 = b1;
02234 b2 = b_i;
02235
02236 iteration++;
02237
02238 q (f, a1, b1, &area1, &error1, &resabs1, &resasc1);
02239 q (f, a2, b2, &area2, &error2, &resabs2, &resasc2);
02240
02241 area12 = area1 + area2;
02242 error12 = error1 + error2;
02243 last_e_i = e_i;
02244
02245
02246
02247
02248
02249
02250
02251
02252 errsum = errsum + error12 - e_i;
02253 area = area + area12 - r_i;
02254
02255 tolerance = GSL_MAX_DBL (epsabs, epsrel * fabs (area));
02256
02257 if (resasc1 != error1 && resasc2 != error2)
02258 {
02259 double delta = r_i - area12;
02260
02261 if (fabs (delta) <= 1.0e-5 * fabs (area12) && error12 >= 0.99 * e_i)
02262 {
02263 if (!extrapolate)
02264 {
02265 roundoff_type1++;
02266 }
02267 else
02268 {
02269 roundoff_type2++;
02270 }
02271 }
02272 if (iteration > 10 && error12 > e_i)
02273 {
02274 roundoff_type3++;
02275 }
02276 }
02277
02278
02279
02280 if (roundoff_type1 + roundoff_type2 >= 10 || roundoff_type3 >= 20)
02281 {
02282 error_type = 2;
02283 }
02284
02285 if (roundoff_type2 >= 5)
02286 {
02287 error_type2 = 1;
02288 }
02289
02290
02291
02292
02293 if (subinterval_too_small (a1, a2, b2))
02294 {
02295 error_type = 4;
02296 }
02297
02298
02299
02300 update (workspace, a1, b1, area1, error1, a2, b2, area2, error2);
02301
02302 if (errsum <= tolerance)
02303 {
02304 goto compute_result;
02305 }
02306
02307 if (error_type)
02308 {
02309 break;
02310 }
02311
02312 if (iteration >= limit - 1)
02313 {
02314 error_type = 1;
02315 break;
02316 }
02317
02318 if (iteration == 2)
02319 {
02320 error_over_large_intervals = errsum;
02321 ertest = tolerance;
02322 append_table (&table, area);
02323 continue;
02324 }
02325
02326 if (disallow_extrapolation)
02327 {
02328 continue;
02329 }
02330
02331 error_over_large_intervals += -last_e_i;
02332
02333 if (current_level < workspace->maximum_level)
02334 {
02335 error_over_large_intervals += error12;
02336 }
02337
02338 if (!extrapolate)
02339 {
02340
02341
02342
02343 if (large_interval (workspace))
02344 continue;
02345
02346 extrapolate = 1;
02347 workspace->nrmax = 1;
02348 }
02349
02350 if (!error_type2 && error_over_large_intervals > ertest)
02351 {
02352 if (increase_nrmax (workspace))
02353 continue;
02354 }
02355
02356
02357
02358 append_table (&table, area);
02359
02360 qelg (&table, &reseps, &abseps);
02361
02362 ktmin++;
02363
02364 if (ktmin > 5 && err_ext < 0.001 * errsum)
02365 {
02366 error_type = 5;
02367 }
02368
02369 if (abseps < err_ext)
02370 {
02371 ktmin = 0;
02372 err_ext = abseps;
02373 res_ext = reseps;
02374 correc = error_over_large_intervals;
02375 ertest = GSL_MAX_DBL (epsabs, epsrel * fabs (reseps));
02376 if (err_ext <= ertest)
02377 break;
02378 }
02379
02380
02381
02382 if (table.n == 1)
02383 {
02384 disallow_extrapolation = 1;
02385 }
02386
02387 if (error_type == 5)
02388 {
02389 break;
02390 }
02391
02392
02393
02394 reset_nrmax (workspace);
02395 extrapolate = 0;
02396 error_over_large_intervals = errsum;
02397
02398 }
02399 while (iteration < limit);
02400
02401 *result = res_ext;
02402 *abserr = err_ext;
02403
02404 if (err_ext == GSL_DBL_MAX)
02405 goto compute_result;
02406
02407 if (error_type || error_type2)
02408 {
02409 if (error_type2)
02410 {
02411 err_ext += correc;
02412 }
02413
02414
02415
02416
02417 if (res_ext != 0.0 && area != 0.0)
02418 {
02419 if (err_ext / fabs (res_ext) > errsum / fabs (area))
02420 goto compute_result;
02421 }
02422 else if (err_ext > errsum)
02423 {
02424 goto compute_result;
02425 }
02426 else if (area == 0.0)
02427 {
02428 goto return_error;
02429 }
02430 }
02431
02432
02433
02434 {
02435 double max_area = GSL_MAX_DBL (fabs (res_ext), fabs (area));
02436
02437 if (!positive_integrand && max_area < 0.01 * resabs0)
02438 goto return_error;
02439 }
02440
02441 {
02442 double ratio = res_ext / area;
02443
02444 if (ratio < 0.01 || ratio > 100.0 || errsum > fabs (area))
02445 error_type = 6;
02446 }
02447
02448 goto return_error;
02449
02450 compute_result:
02451
02452 *result = sum_results (workspace);
02453 *abserr = errsum;
02454
02455 return_error:
02456
02457 if (error_type > 2)
02458 error_type--;
02459
02460
02461
02462 if (error_type == 0)
02463 {
02464 return GSL_SUCCESS;
02465 }
02466 else if (error_type == 1)
02467 {
02468 GSL_ERROR ("number of iterations was insufficient", GSL_EMAXITER);
02469 }
02470 else if (error_type == 2)
02471 {
02472 GSL_ERROR ("cannot reach tolerance because of roundoff error",
02473 GSL_EROUND);
02474 }
02475 else if (error_type == 3)
02476 {
02477 GSL_ERROR ("bad integrand behavior found in the integration interval",
02478 GSL_ESING);
02479 }
02480 else if (error_type == 4)
02481 {
02482 GSL_ERROR ("roundoff error detected in the extrapolation table",
02483 GSL_EROUND);
02484 }
02485 else if (error_type == 5)
02486 {
02487 GSL_ERROR ("integral is divergent, or slowly convergent",
02488 GSL_EDIVERGE);
02489 }
02490
02491 GSL_ERROR ("could not integrate function", GSL_EFAILED);
02492 }