ROOT  6.06/09
Reference Guide
RooAdaptiveGaussKronrodIntegrator1D.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitCore *
4  * @(#)root/roofitcore:$Id$
5  * Authors: *
6  * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7  * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8  * *
9  * Copyright (c) 2000-2005, Regents of the University of California *
10  * and Stanford University. All rights reserved. *
11  * *
12  * Redistribution and use in source and binary forms, *
13  * with or without modification, are permitted according to the terms *
14  * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15  *****************************************************************************/
16 
17 //////////////////////////////////////////////////////////////////////////////
18 //
19 // BEGIN_HTML
20 // RooAdaptiveGaussKronrodIntegrator1D implements the Gauss-Kronrod integration algorithm.
21 //
22 // An adaptive Gaussian quadrature method for numerical integration in
23 // which error is estimation based on evaluation at special points
24 // known as "Kronrod points." By suitably picking these points,
25 // abscissas from previous iterations can be reused as part of the new
26 // set of points, whereas usual Gaussian quadrature would require
27 // recomputation of all abscissas at each iteration.
28 //
29 // This class automatically handles (-inf,+inf) integrals by dividing
30 // the integration in three regions (-inf,-1), (-1,1), (1,inf) and
31 // calculating the 1st and 3rd term using a x -> 1/x coordinate
32 // transformation
33 //
34 // This class embeds the adaptive Gauss-Kronrod integrator from the
35 // GNU Scientific Library version 1.5 and applies a chosen rule ( 10-,
36 // 21-, 31-, 41, 51- or 61-point). The integration domain is
37 // subdivided and recursively integrated until the required precision
38 // is reached.
39 //
40 // For integrands with integrable singulaties the Wynn epsilon rule
41 // can be selected to speed up the converges of these integrals
42 // END_HTML
43 //
44 
45 #include "RooFit.h"
46 
47 #include <assert.h>
48 #include <stdlib.h>
49 #include "TClass.h"
50 #include "Riostream.h"
52 #include "RooArgSet.h"
53 #include "RooRealVar.h"
54 #include "RooNumber.h"
55 #include "RooNumIntFactory.h"
56 #include "RooIntegratorBinding.h"
57 #include "TMath.h"
58 #include "RooMsgService.h"
59 
60 using namespace std ;
61 
62 
64 ;
65 
66 
67 // --- From GSL_MATH.h -------------------------------------------
68 struct gsl_function_struct
69 {
70  double (* function) (double x, void * params);
71  void * params;
72 };
73 typedef struct gsl_function_struct gsl_function ;
74 #define GSL_FN_EVAL(F,x) (*((F)->function))(x,(F)->params)
75 
76 //----From GSL_INTEGRATION.h ---------------------------------------
77 typedef struct
78  {
79  size_t limit;
80  size_t size;
81  size_t nrmax;
82  size_t i;
83  size_t maximum_level;
84  double *alist;
85  double *blist;
86  double *rlist;
87  double *elist;
88  size_t *order;
89  size_t *level;
90  }
91 gsl_integration_workspace;
92 
93 gsl_integration_workspace *
94  gsl_integration_workspace_alloc (const size_t n);
95 
96 void
97  gsl_integration_workspace_free (gsl_integration_workspace * w);
98 
100  double a, double b,
101  double epsabs, double epsrel, size_t limit,
102  int key,
103  gsl_integration_workspace * workspace,
104  double *result, double *abserr);
105 
106 int
108  double a, double b,
109  double epsabs, double epsrel, size_t limit,
110  gsl_integration_workspace * workspace,
111  double * result, double * abserr) ;
112 
113 int
115  double epsabs, double epsrel, size_t limit,
116  gsl_integration_workspace * workspace,
117  double *result, double *abserr) ;
118 
119 int
121  double b,
122  double epsabs, double epsrel, size_t limit,
123  gsl_integration_workspace * workspace,
124  double *result, double *abserr) ;
125 
126 int
128  double a,
129  double epsabs, double epsrel, size_t limit,
130  gsl_integration_workspace * workspace,
131  double *result, double *abserr) ;
132 
133 
134 //-------------------------------------------------------------------
135 
136 
137 ////////////////////////////////////////////////////////////////////////////////
138 /// Register this class with RooNumIntConfig as a possible choice of numeric
139 /// integrator for one-dimensional integrals over finite and infinite domains
140 
142 {
143  RooRealVar maxSeg("maxSeg","maximum number of segments",100) ;
144  RooCategory method("method","Integration method for each segment") ;
145  method.defineType("WynnEpsilon",0) ;
146  method.defineType("15Points",1) ;
147  method.defineType("21Points",2) ;
148  method.defineType("31Points",3) ;
149  method.defineType("41Points",4) ;
150  method.defineType("51Points",5) ;
151  method.defineType("61Points",6) ;
152  method.setIndex(2) ;
154 }
155 
156 
157 
158 ////////////////////////////////////////////////////////////////////////////////
159 /// coverity[UNINIT_CTOR]
160 /// Default constructor
161 
163 {
164 }
165 
166 
167 
168 
169 ////////////////////////////////////////////////////////////////////////////////
170 /// Constructor taking a function binding and a configuration object
171 
173  const RooNumIntConfig& config) :
175  _epsAbs(config.epsRel()),
176  _epsRel(config.epsAbs()),
177  _workspace(0)
178 {
179  // Use this form of the constructor to integrate over the function's default range.
180  const RooArgSet& confSet = config.getConfigSection(IsA()->GetName()) ;
181  _maxSeg = (Int_t) confSet.getRealValue("maxSeg",100) ;
182  _methodKey = confSet.getCatIndex("method",2) ;
183 
185  _valid= initialize();
186 }
187 
188 
189 
190 ////////////////////////////////////////////////////////////////////////////////
191 /// Constructor taking a function binding, an integration range and a configuration object
192 
195  const RooNumIntConfig& config) :
197  _epsAbs(config.epsRel()),
198  _epsRel(config.epsAbs()),
199  _workspace(0),
200  _xmin(xmin),
201  _xmax(xmax)
202 {
203  // Use this form of the constructor to integrate over the function's default range.
204  const RooArgSet& confSet = config.getConfigSection(IsA()->GetName()) ;
205  _maxSeg = (Int_t) confSet.getRealValue("maxSeg",100) ;
206  _methodKey = confSet.getCatIndex("method",2) ;
207 
209  _valid= initialize();
210 }
211 
212 
213 
214 ////////////////////////////////////////////////////////////////////////////////
215 /// Virtual constructor
216 
218 {
219  return new RooAdaptiveGaussKronrodIntegrator1D(function,config) ;
220 }
221 
222 
223 
224 ////////////////////////////////////////////////////////////////////////////////
225 /// Initialize integrator allocate buffers and setup GSL workspace
226 
228 {
229  // Allocate coordinate buffer size after number of function dimensions
230  _x = new Double_t[_function->getDimension()] ;
232 
233  return checkLimits();
234 }
235 
236 
237 
238 ////////////////////////////////////////////////////////////////////////////////
239 /// Destructor.
240 
242 {
243  if (_workspace) {
244  gsl_integration_workspace_free ((gsl_integration_workspace*) _workspace) ;
245  }
246  if (_x) {
247  delete[] _x ;
248  }
249 }
250 
251 
252 
253 ////////////////////////////////////////////////////////////////////////////////
254 /// Change our integration limits. Return kTRUE if the new limits are
255 /// ok, or otherwise kFALSE. Always returns kFALSE and does nothing
256 /// if this object was constructed to always use our integrand's limits.
257 
259 {
260  if(_useIntegrandLimits) {
261  coutE(Integration) << "RooAdaptiveGaussKronrodIntegrator1D::setLimits: cannot override integrand's limits" << endl;
262  return kFALSE;
263  }
264 
265  _xmin= *xmin;
266  _xmax= *xmax;
267  return checkLimits();
268 }
269 
270 
271 
272 ////////////////////////////////////////////////////////////////////////////////
273 /// Check that our integration range is finite and otherwise return kFALSE.
274 /// Update the limits from the integrand if requested.
275 
277 {
278  if(_useIntegrandLimits) {
279  assert(0 != integrand() && integrand()->isValid());
280  _xmin= integrand()->getMinLimit(0);
281  _xmax= integrand()->getMaxLimit(0);
282  }
283 
284  // Determine domain type
287 
288  if (!infLo && !infHi) {
289  _domainType = Closed ;
290  } else if (infLo && infHi) {
291  _domainType = Open ;
292  } else if (infLo && !infHi) {
293  _domainType = OpenLo ;
294  } else {
295  _domainType = OpenHi ;
296  }
297 
298 
299  return kTRUE ;
300 }
301 
302 
303 
304 ////////////////////////////////////////////////////////////////////////////////
305 /// Glue function interacing to GSL code
306 
308 {
310  return instance->integrand(instance->xvec(x)) ;
311 }
312 
313 
314 
315 ////////////////////////////////////////////////////////////////////////////////
316 /// Calculate and return integral at at given parameter values
317 
319 {
320  assert(isValid());
321 
322  // Copy yvec to xvec if provided
323  if (yvec) {
324  UInt_t i ; for (i=0 ; i<_function->getDimension()-1 ; i++) {
325  _x[i+1] = yvec[i] ;
326  }
327  }
328 
329  // Setup glue function
330  gsl_function F;
332  F.params = this ;
333 
334  // Return values
335  double result, error;
336 
337  // Call GSL implementation of integeator
338  switch(_domainType) {
339  case Closed:
340  if (_methodKey==0) {
341  gsl_integration_qags (&F, _xmin, _xmax, _epsAbs, _epsRel, _maxSeg, (gsl_integration_workspace*)_workspace,&result, &error);
342  } else {
343  gsl_integration_qag (&F, _xmin, _xmax, _epsAbs, _epsRel, _maxSeg, _methodKey, (gsl_integration_workspace*)_workspace,&result, &error);
344  }
345  break ;
346  case OpenLo:
347  gsl_integration_qagil (&F, _xmax, _epsAbs, _epsRel, _maxSeg, (gsl_integration_workspace*)_workspace,&result, &error);
348  break ;
349  case OpenHi:
350  gsl_integration_qagiu (&F, _xmin, _epsAbs, _epsRel, _maxSeg, (gsl_integration_workspace*)_workspace,&result, &error);
351  break ;
352  case Open:
353  gsl_integration_qagi (&F, _epsAbs, _epsRel, _maxSeg, (gsl_integration_workspace*)_workspace,&result, &error);
354  break ;
355  }
356 
357  return result;
358 }
359 
360 
361 // ----------------------------------------------------------------------------
362 // ---------- Code below imported from GSL version 1.5 ------------------------
363 // ----------------------------------------------------------------------------
364 
365 /*
366  *
367  * Copyright (C) 1996, 1997, 1998, 1999, 2000 Brian Gough
368  *
369  * This program is free software; you can redistribute it and/or modify
370  * it under the terms of the GNU General Public License as published by
371  * the Free Software Foundation; either version 2 of the License, or (at
372  * your option) any later version.
373  *
374  * This program is distributed in the hope that it will be useful, but
375  * WITHOUT ANY WARRANTY; without even the implied warranty of
376  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
377  * General Public License for more details.
378  *
379  * You should have received a copy of the GNU General Public License
380  * along with this program; if not, write to the Free Software
381  * Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
382  */
383 
384 #define GSL_SUCCESS 0
385 #define GSL_EDOM 1 /* input domain error, e.g sqrt(-1) */
386 #define GSL_ENOMEM 8 /* malloc failed */
387 #define GSL_EBADTOL 13 /* user specified an invalid tolerance */
388 #define GSL_ETOL 14 /* failed to reach the specified tolerance */
389 #define GSL_ERROR(a,b) oocoutE((TObject*)0,Integration) << "RooAdaptiveGaussKronrodIntegrator1D::integral() ERROR: " << a << endl ; return b ;
390 #define GSL_DBL_MIN 2.2250738585072014e-308
391 #define GSL_DBL_MAX 1.7976931348623157e+308
392 #define GSL_DBL_EPSILON 2.2204460492503131e-16
393 
394 #define GSL_EINVAL 2
395 #define GSL_EMAXITER 3
396 #define GSL_ESING 4
397 #define GSL_EFAILED 5
398 #define GSL_EDIVERGE 6
399 #define GSL_EROUND 7
400 
401 #define GSL_ERROR_VAL(reason, gsl_errno, value) return value ;
402 
403 #define GSL_MAX(a,b) ((a) > (b) ? (a) : (b))
404 extern inline double
405 GSL_MAX_DBL (double a, double b)
406 {
407  return GSL_MAX (a, b);
408 }
409 
410 double gsl_coerce_double (const double x);
411 
412 double
413 gsl_coerce_double (const double x)
414 {
415  volatile double y;
416  y = x;
417  return y;
418 }
419 #define GSL_COERCE_DBL(x) (gsl_coerce_double(x))
420 
421 // INCLUDED BELOW #include <gsl/gsl_integration.h>
422 
423 
424 /* Workspace for adaptive integrators */
425 // WVE MOVED TO HEAD OF FILE
426 
427 
428 /* Definition of an integration rule */
429 
430 typedef void gsl_integration_rule (const gsl_function * f,
431  double a, double b,
432  double *result, double *abserr,
433  double *defabs, double *resabs);
434 
435 void gsl_integration_qk15 (const gsl_function * f, double a, double b,
436  double *result, double *abserr,
437  double *resabs, double *resasc);
438 
439 void gsl_integration_qk21 (const gsl_function * f, double a, double b,
440  double *result, double *abserr,
441  double *resabs, double *resasc);
442 
443 void gsl_integration_qk31 (const gsl_function * f, double a, double b,
444  double *result, double *abserr,
445  double *resabs, double *resasc);
446 
447 void gsl_integration_qk41 (const gsl_function * f, double a, double b,
448  double *result, double *abserr,
449  double *resabs, double *resasc);
450 
451 void gsl_integration_qk51 (const gsl_function * f, double a, double b,
452  double *result, double *abserr,
453  double *resabs, double *resasc);
454 
455 void gsl_integration_qk61 (const gsl_function * f, double a, double b,
456  double *result, double *abserr,
457  double *resabs, double *resasc);
458 
459 void gsl_integration_qcheb (gsl_function * f, double a, double b,
460  double *cheb12, double *cheb24);
461 
462 /* The low-level integration rules in QUADPACK are identified by small
463  integers (1-6). We'll use symbolic constants to refer to them. */
464 
465 enum
466  {
467  GSL_INTEG_GAUSS15 = 1, /* 15 point Gauss-Kronrod rule */
468  GSL_INTEG_GAUSS21 = 2, /* 21 point Gauss-Kronrod rule */
469  GSL_INTEG_GAUSS31 = 3, /* 31 point Gauss-Kronrod rule */
470  GSL_INTEG_GAUSS41 = 4, /* 41 point Gauss-Kronrod rule */
471  GSL_INTEG_GAUSS51 = 5, /* 51 point Gauss-Kronrod rule */
472  GSL_INTEG_GAUSS61 = 6 /* 61 point Gauss-Kronrod rule */
473  };
474 
475 
476 void
477 gsl_integration_qk (const int n, const double xgk[],
478  const double wg[], const double wgk[],
479  double fv1[], double fv2[],
480  const gsl_function *f, double a, double b,
481  double * result, double * abserr,
482  double * resabs, double * resasc);
483 
484 
485 int gsl_integration_qag (const gsl_function * f,
486  double a, double b,
487  double epsabs, double epsrel, size_t limit,
488  int key,
489  gsl_integration_workspace * workspace,
490  double *result, double *abserr);
491 
492 
493 // INCLUDED BELOW #include "initialise.c"
494 static inline
495 void initialise (gsl_integration_workspace * workspace, double a, double b);
496 
497 static inline
498 void initialise (gsl_integration_workspace * workspace, double a, double b)
499 {
500  workspace->size = 0;
501  workspace->nrmax = 0;
502  workspace->i = 0;
503  workspace->alist[0] = a;
504  workspace->blist[0] = b;
505  workspace->rlist[0] = 0.0;
506  workspace->elist[0] = 0.0;
507  workspace->order[0] = 0;
508  workspace->level[0] = 0;
509 
510  workspace->maximum_level = 0;
511 }
512 
513 // INCLUDED BELOW #include "set_initial.c"
514 static inline
515 void set_initial_result (gsl_integration_workspace * workspace,
516  double result, double error);
517 
518 static inline
519 void set_initial_result (gsl_integration_workspace * workspace,
520  double result, double error)
521 {
522  workspace->size = 1;
523  workspace->rlist[0] = result;
524  workspace->elist[0] = error;
525 }
526 
527 // INCLUDED BELOW #include "qpsrt.c"
528 static inline void
529 qpsrt (gsl_integration_workspace * workspace);
530 
531 static inline
532 void qpsrt (gsl_integration_workspace * workspace)
533 {
534  const size_t last = workspace->size - 1;
535  const size_t limit = workspace->limit;
536 
537  double * elist = workspace->elist;
538  size_t * order = workspace->order;
539 
540  double errmax ;
541  double errmin ;
542  int i, k, top;
543 
544  size_t i_nrmax = workspace->nrmax;
545  size_t i_maxerr = order[i_nrmax] ;
546 
547  /* Check whether the list contains more than two error estimates */
548 
549  if (last < 2)
550  {
551  order[0] = 0 ;
552  order[1] = 1 ;
553  workspace->i = i_maxerr ;
554  return ;
555  }
556 
557  errmax = elist[i_maxerr] ;
558 
559  /* This part of the routine is only executed if, due to a difficult
560  integrand, subdivision increased the error estimate. In the normal
561  case the insert procedure should start after the nrmax-th largest
562  error estimate. */
563 
564  while (i_nrmax > 0 && errmax > elist[order[i_nrmax - 1]])
565  {
566  order[i_nrmax] = order[i_nrmax - 1] ;
567  i_nrmax-- ;
568  }
569 
570  /* Compute the number of elements in the list to be maintained in
571  descending order. This number depends on the number of
572  subdivisions still allowed. */
573 
574  if(last < (limit/2 + 2))
575  {
576  top = last ;
577  }
578  else
579  {
580  top = limit - last + 1;
581  }
582 
583  /* Insert errmax by traversing the list top-down, starting
584  comparison from the element elist(order(i_nrmax+1)). */
585 
586  i = i_nrmax + 1 ;
587 
588  /* The order of the tests in the following line is important to
589  prevent a segmentation fault */
590 
591  while (i < top && errmax < elist[order[i]])
592  {
593  order[i-1] = order[i] ;
594  i++ ;
595  }
596 
597  order[i-1] = i_maxerr ;
598 
599  /* Insert errmin by traversing the list bottom-up */
600 
601  errmin = elist[last] ;
602 
603  k = top - 1 ;
604 
605  while (k > i - 2 && errmin >= elist[order[k]])
606  {
607  order[k+1] = order[k] ;
608  k-- ;
609  }
610 
611  order[k+1] = last ;
612 
613  /* Set i_max and e_max */
614 
615  i_maxerr = order[i_nrmax] ;
616 
617  workspace->i = i_maxerr ;
618  workspace->nrmax = i_nrmax ;
619 }
620 
621 
622 
623 // INCLUDED BELOW #include "util.c"
624 static inline
625 void update (gsl_integration_workspace * workspace,
626  double a1, double b1, double area1, double error1,
627  double a2, double b2, double area2, double error2);
628 
629 static inline void
630 retrieve (const gsl_integration_workspace * workspace,
631  double * a, double * b, double * r, double * e);
632 
633 
634 
635 static inline
636 void update (gsl_integration_workspace * workspace,
637  double a1, double b1, double area1, double error1,
638  double a2, double b2, double area2, double error2)
639 {
640  double * alist = workspace->alist ;
641  double * blist = workspace->blist ;
642  double * rlist = workspace->rlist ;
643  double * elist = workspace->elist ;
644  size_t * level = workspace->level ;
645 
646  const size_t i_max = workspace->i ;
647  const size_t i_new = workspace->size ;
648 
649  const size_t new_level = workspace->level[i_max] + 1;
650 
651  /* append the newly-created intervals to the list */
652 
653  if (error2 > error1)
654  {
655  alist[i_max] = a2; /* blist[maxerr] is already == b2 */
656  rlist[i_max] = area2;
657  elist[i_max] = error2;
658  level[i_max] = new_level;
659 
660  alist[i_new] = a1;
661  blist[i_new] = b1;
662  rlist[i_new] = area1;
663  elist[i_new] = error1;
664  level[i_new] = new_level;
665  }
666  else
667  {
668  blist[i_max] = b1; /* alist[maxerr] is already == a1 */
669  rlist[i_max] = area1;
670  elist[i_max] = error1;
671  level[i_max] = new_level;
672 
673  alist[i_new] = a2;
674  blist[i_new] = b2;
675  rlist[i_new] = area2;
676  elist[i_new] = error2;
677  level[i_new] = new_level;
678  }
679 
680  workspace->size++;
681 
682  if (new_level > workspace->maximum_level)
683  {
684  workspace->maximum_level = new_level;
685  }
686 
687  qpsrt (workspace) ;
688 }
689 
690 static inline void
691 retrieve (const gsl_integration_workspace * workspace,
692  double * a, double * b, double * r, double * e)
693 {
694  const size_t i = workspace->i;
695  double * alist = workspace->alist;
696  double * blist = workspace->blist;
697  double * rlist = workspace->rlist;
698  double * elist = workspace->elist;
699 
700  *a = alist[i] ;
701  *b = blist[i] ;
702  *r = rlist[i] ;
703  *e = elist[i] ;
704 }
705 
706 static inline double
707 sum_results (const gsl_integration_workspace * workspace);
708 
709 static inline double
710 sum_results (const gsl_integration_workspace * workspace)
711 {
712  const double * const rlist = workspace->rlist ;
713  const size_t n = workspace->size;
714 
715  size_t k;
716  double result_sum = 0;
717 
718  for (k = 0; k < n; k++)
719  {
720  result_sum += rlist[k];
721  }
722 
723  return result_sum;
724 }
725 
726 static inline int
727 subinterval_too_small (double a1, double a2, double b2);
728 
729 static inline int
730 subinterval_too_small (double a1, double a2, double b2)
731 {
732  const double e = GSL_DBL_EPSILON;
733  const double u = GSL_DBL_MIN;
734 
735  double tmp = (1 + 100 * e) * (fabs (a2) + 1000 * u);
736 
737  int status = fabs (a1) <= tmp && fabs (b2) <= tmp;
738 
739  return status;
740 }
741 
742 
743 static int
744 qag (const gsl_function *f,
745  const double a, const double b,
746  const double epsabs, const double epsrel,
747  const size_t limit,
748  gsl_integration_workspace * workspace,
749  double * result, double * abserr,
751 
752 int
754  double a, double b,
755  double epsabs, double epsrel, size_t limit,
756  int key,
757  gsl_integration_workspace * workspace,
758  double * result, double * abserr)
759 {
760  int status ;
761  gsl_integration_rule * integration_rule = gsl_integration_qk15 ;
762 
763  if (key < GSL_INTEG_GAUSS15)
764  {
765  key = GSL_INTEG_GAUSS15 ;
766  }
767  else if (key > GSL_INTEG_GAUSS61)
768  {
769  key = GSL_INTEG_GAUSS61 ;
770  }
771 
772  switch (key)
773  {
774  case GSL_INTEG_GAUSS15:
775  integration_rule = gsl_integration_qk15 ;
776  break ;
777  case GSL_INTEG_GAUSS21:
778  integration_rule = gsl_integration_qk21 ;
779  break ;
780  case GSL_INTEG_GAUSS31:
781  integration_rule = gsl_integration_qk31 ;
782  break ;
783  case GSL_INTEG_GAUSS41:
784  integration_rule = gsl_integration_qk41 ;
785  break ;
786  case GSL_INTEG_GAUSS51:
787  integration_rule = gsl_integration_qk51 ;
788  break ;
789  case GSL_INTEG_GAUSS61:
790  integration_rule = gsl_integration_qk61 ;
791  break ;
792  }
793 
794  status = qag (f, a, b, epsabs, epsrel, limit,
795  workspace,
796  result, abserr,
797  integration_rule) ;
798 
799  return status ;
800 }
801 
802 static int
803 qag (const gsl_function * f,
804  const double a, const double b,
805  const double epsabs, const double epsrel,
806  const size_t limit,
807  gsl_integration_workspace * workspace,
808  double *result, double *abserr,
810 {
811  double area, errsum;
812  double result0, abserr0, resabs0, resasc0;
813  double tolerance;
814  size_t iteration = 0;
815  int roundoff_type1 = 0, roundoff_type2 = 0, error_type = 0;
816 
817  double round_off;
818 
819  /* Initialize results */
820 
821  initialise (workspace, a, b);
822 
823  *result = 0;
824  *abserr = 0;
825 
826  if (limit > workspace->limit)
827  {
828  GSL_ERROR ("iteration limit exceeds available workspace", GSL_EINVAL) ;
829  }
830 
831  if (epsabs <= 0 && (epsrel < 50 * GSL_DBL_EPSILON || epsrel < 0.5e-28))
832  {
833  GSL_ERROR ("tolerance cannot be acheived with given epsabs and epsrel",
834  GSL_EBADTOL);
835  }
836 
837  /* perform the first integration */
838 
839  q (f, a, b, &result0, &abserr0, &resabs0, &resasc0);
840 
841  set_initial_result (workspace, result0, abserr0);
842 
843  /* Test on accuracy */
844 
845  tolerance = GSL_MAX_DBL (epsabs, epsrel * fabs (result0));
846 
847  /* need IEEE rounding here to match original quadpack behavior */
848 
849  round_off = GSL_COERCE_DBL (50 * GSL_DBL_EPSILON * resabs0);
850 
851  if (abserr0 <= round_off && abserr0 > tolerance)
852  {
853  *result = result0;
854  *abserr = abserr0;
855 
856  GSL_ERROR ("cannot reach tolerance because of roundoff error "
857  "on first attempt", GSL_EROUND);
858  }
859  else if ((abserr0 <= tolerance && abserr0 != resasc0) || abserr0 == 0.0)
860  {
861  *result = result0;
862  *abserr = abserr0;
863 
864  return GSL_SUCCESS;
865  }
866  else if (limit == 1)
867  {
868  *result = result0;
869  *abserr = abserr0;
870 
871  GSL_ERROR ("a maximum of one iteration was insufficient", GSL_EMAXITER);
872  }
873 
874  area = result0;
875  errsum = abserr0;
876 
877  iteration = 1;
878 
879  do
880  {
881  double a1, b1, a2, b2;
882  double a_i, b_i, r_i, e_i;
883  double area1 = 0, area2 = 0, area12 = 0;
884  double error1 = 0, error2 = 0, error12 = 0;
885  double resasc1, resasc2;
886  double resabs1, resabs2;
887 
888  /* Bisect the subinterval with the largest error estimate */
889 
890  retrieve (workspace, &a_i, &b_i, &r_i, &e_i);
891 
892  a1 = a_i;
893  b1 = 0.5 * (a_i + b_i);
894  a2 = b1;
895  b2 = b_i;
896 
897  q (f, a1, b1, &area1, &error1, &resabs1, &resasc1);
898  q (f, a2, b2, &area2, &error2, &resabs2, &resasc2);
899 
900  area12 = area1 + area2;
901  error12 = error1 + error2;
902 
903  errsum += (error12 - e_i);
904  area += area12 - r_i;
905 
906  if (resasc1 != error1 && resasc2 != error2)
907  {
908  double delta = r_i - area12;
909 
910  if (fabs (delta) <= 1.0e-5 * fabs (area12) && error12 >= 0.99 * e_i)
911  {
912  roundoff_type1++;
913  }
914  if (iteration >= 10 && error12 > e_i)
915  {
916  roundoff_type2++;
917  }
918  }
919 
920  tolerance = GSL_MAX_DBL (epsabs, epsrel * fabs (area));
921 
922  if (errsum > tolerance)
923  {
924  if (roundoff_type1 >= 6 || roundoff_type2 >= 20)
925  {
926  error_type = 2; /* round off error */
927  }
928 
929  /* set error flag in the case of bad integrand behaviour at
930  a point of the integration range */
931 
932  if (subinterval_too_small (a1, a2, b2))
933  {
934  error_type = 3;
935  }
936  }
937 
938  update (workspace, a1, b1, area1, error1, a2, b2, area2, error2);
939 
940  retrieve (workspace, &a_i, &b_i, &r_i, &e_i);
941 
942  iteration++;
943 
944  }
945  while (iteration < limit && !error_type && errsum > tolerance);
946 
947  *result = sum_results (workspace);
948  *abserr = errsum;
949 
950  if (errsum <= tolerance)
951  {
952  return GSL_SUCCESS;
953  }
954  else if (error_type == 2)
955  {
956  GSL_ERROR ("roundoff error prevents tolerance from being achieved",
957  GSL_EROUND);
958  }
959  else if (error_type == 3)
960  {
961  GSL_ERROR ("bad integrand behavior found in the integration interval",
962  GSL_ESING);
963  }
964  else if (iteration == limit)
965  {
966  GSL_ERROR ("maximum number of subdivisions reached", GSL_EMAXITER);
967  }
968 
969  GSL_ERROR ("could not integrate function", GSL_EFAILED);
970 }
971 
972 
973 // INCLUDED BELOW #include "err.c"
974 static double rescale_error (double err, const double result_abs, const double result_asc) ;
975 
976 static double
977 rescale_error (double err, const double result_abs, const double result_asc)
978 {
979  err = fabs(err) ;
980 
981  if (result_asc != 0 && err != 0)
982  {
983  double scale = TMath::Power((200 * err / result_asc), 1.5) ;
984 
985  if (scale < 1)
986  {
987  err = result_asc * scale ;
988  }
989  else
990  {
991  err = result_asc ;
992  }
993  }
994  if (result_abs > GSL_DBL_MIN / (50 * GSL_DBL_EPSILON))
995  {
996  double min_err = 50 * GSL_DBL_EPSILON * result_abs ;
997 
998  if (min_err > err)
999  {
1000  err = min_err ;
1001  }
1002  }
1003 
1004  return err ;
1005 }
1006 
1007 
1008 void
1009 gsl_integration_qk (const int n,
1010  const double xgk[], const double wg[], const double wgk[],
1011  double fv1[], double fv2[],
1012  const gsl_function * f, double a, double b,
1013  double *result, double *abserr,
1014  double *resabs, double *resasc)
1015 {
1016 
1017  const double center = 0.5 * (a + b);
1018  const double half_length = 0.5 * (b - a);
1019  const double abs_half_length = fabs (half_length);
1020  const double f_center = GSL_FN_EVAL (f, center);
1021 
1022  double result_gauss = 0;
1023  double result_kronrod = f_center * wgk[n - 1];
1024 
1025  double result_abs = fabs (result_kronrod);
1026  double result_asc = 0;
1027  double mean = 0, err = 0;
1028 
1029  int j;
1030 
1031  if (n % 2 == 0)
1032  {
1033  result_gauss = f_center * wg[n / 2 - 1];
1034  }
1035 
1036  for (j = 0; j < (n - 1) / 2; j++)
1037  {
1038  const int jtw = j * 2 + 1; /* j=1,2,3 jtw=2,4,6 */
1039  const double abscissa = half_length * xgk[jtw];
1040  const double fval1 = GSL_FN_EVAL (f, center - abscissa);
1041  const double fval2 = GSL_FN_EVAL (f, center + abscissa);
1042  const double fsum = fval1 + fval2;
1043  fv1[jtw] = fval1;
1044  fv2[jtw] = fval2;
1045  result_gauss += wg[j] * fsum;
1046  result_kronrod += wgk[jtw] * fsum;
1047  result_abs += wgk[jtw] * (fabs (fval1) + fabs (fval2));
1048  }
1049 
1050  for (j = 0; j < n / 2; j++)
1051  {
1052  int jtwm1 = j * 2;
1053  const double abscissa = half_length * xgk[jtwm1];
1054  const double fval1 = GSL_FN_EVAL (f, center - abscissa);
1055  const double fval2 = GSL_FN_EVAL (f, center + abscissa);
1056  fv1[jtwm1] = fval1;
1057  fv2[jtwm1] = fval2;
1058  result_kronrod += wgk[jtwm1] * (fval1 + fval2);
1059  result_abs += wgk[jtwm1] * (fabs (fval1) + fabs (fval2));
1060  };
1061 
1062  mean = result_kronrod * 0.5;
1063 
1064  result_asc = wgk[n - 1] * fabs (f_center - mean);
1065 
1066  for (j = 0; j < n - 1; j++)
1067  {
1068  result_asc += wgk[j] * (fabs (fv1[j] - mean) + fabs (fv2[j] - mean));
1069  }
1070 
1071  /* scale by the width of the integration region */
1072 
1073  err = (result_kronrod - result_gauss) * half_length;
1074 
1075  result_kronrod *= half_length;
1076  result_abs *= abs_half_length;
1077  result_asc *= abs_half_length;
1078 
1079  *result = result_kronrod;
1080  *resabs = result_abs;
1081  *resasc = result_asc;
1082  *abserr = rescale_error (err, result_abs, result_asc);
1083 
1084 }
1085 
1086 /* Gauss quadrature weights and kronrod quadrature abscissae and
1087  weights as evaluated with 80 decimal digit arithmetic by
1088  L. W. Fullerton, Bell Labs, Nov. 1981. */
1089 
1090 static const double xgkA[8] = /* abscissae of the 15-point kronrod rule */
1091 {
1092  0.991455371120812639206854697526329,
1093  0.949107912342758524526189684047851,
1094  0.864864423359769072789712788640926,
1095  0.741531185599394439863864773280788,
1096  0.586087235467691130294144838258730,
1097  0.405845151377397166906606412076961,
1098  0.207784955007898467600689403773245,
1099  0.000000000000000000000000000000000
1100 };
1101 
1102 /* xgk[1], xgk[3], ... abscissae of the 7-point gauss rule.
1103  xgk[0], xgk[2], ... abscissae to optimally extend the 7-point gauss rule */
1104 
1105 static const double wgA[4] = /* weights of the 7-point gauss rule */
1106 {
1107  0.129484966168869693270611432679082,
1108  0.279705391489276667901467771423780,
1109  0.381830050505118944950369775488975,
1110  0.417959183673469387755102040816327
1111 };
1112 
1113 static const double wgkA[8] = /* weights of the 15-point kronrod rule */
1114 {
1115  0.022935322010529224963732008058970,
1116  0.063092092629978553290700663189204,
1117  0.104790010322250183839876322541518,
1118  0.140653259715525918745189590510238,
1119  0.169004726639267902826583426598550,
1120  0.190350578064785409913256402421014,
1121  0.204432940075298892414161999234649,
1122  0.209482141084727828012999174891714
1123 };
1124 
1125 void
1126 gsl_integration_qk15 (const gsl_function * f, double a, double b,
1127  double *result, double *abserr,
1128  double *resabs, double *resasc)
1129 {
1130  double fv1[8], fv2[8];
1131  // coverity[UNINIT_CTOR]
1132  gsl_integration_qk (8, xgkA, wgA, wgkA, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
1133 }
1134 
1135 /* Gauss quadrature weights and kronrod quadrature abscissae and
1136  weights as evaluated with 80 decimal digit arithmetic by
1137  L. W. Fullerton, Bell Labs, Nov. 1981. */
1138 
1139 static const double xgkB[11] = /* abscissae of the 21-point kronrod rule */
1140 {
1141  0.995657163025808080735527280689003,
1142  0.973906528517171720077964012084452,
1143  0.930157491355708226001207180059508,
1144  0.865063366688984510732096688423493,
1145  0.780817726586416897063717578345042,
1146  0.679409568299024406234327365114874,
1147  0.562757134668604683339000099272694,
1148  0.433395394129247190799265943165784,
1149  0.294392862701460198131126603103866,
1150  0.148874338981631210884826001129720,
1151  0.000000000000000000000000000000000
1152 };
1153 
1154 /* xgk[1], xgk[3], ... abscissae of the 10-point gauss rule.
1155  xgk[0], xgk[2], ... abscissae to optimally extend the 10-point gauss rule */
1156 
1157 static const double wgB[5] = /* weights of the 10-point gauss rule */
1158 {
1159  0.066671344308688137593568809893332,
1160  0.149451349150580593145776339657697,
1161  0.219086362515982043995534934228163,
1162  0.269266719309996355091226921569469,
1163  0.295524224714752870173892994651338
1164 };
1165 
1166 static const double wgkB[11] = /* weights of the 21-point kronrod rule */
1167 {
1168  0.011694638867371874278064396062192,
1169  0.032558162307964727478818972459390,
1170  0.054755896574351996031381300244580,
1171  0.075039674810919952767043140916190,
1172  0.093125454583697605535065465083366,
1173  0.109387158802297641899210590325805,
1174  0.123491976262065851077958109831074,
1175  0.134709217311473325928054001771707,
1176  0.142775938577060080797094273138717,
1177  0.147739104901338491374841515972068,
1178  0.149445554002916905664936468389821
1179 };
1180 
1181 
1182 void
1183 gsl_integration_qk21 (const gsl_function * f, double a, double b,
1184  double *result, double *abserr,
1185  double *resabs, double *resasc)
1186 {
1187  double fv1[11], fv2[11];
1188  // coverity[UNINIT_CTOR]
1189  gsl_integration_qk (11, xgkB, wgB, wgkB, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
1190 }
1191 
1192 /* Gauss quadrature weights and kronrod quadrature abscissae and
1193  weights as evaluated with 80 decimal digit arithmetic by
1194  L. W. Fullerton, Bell Labs, Nov. 1981. */
1195 
1196 static const double xgkC[16] = /* abscissae of the 31-point kronrod rule */
1197 {
1198  0.998002298693397060285172840152271,
1199  0.987992518020485428489565718586613,
1200  0.967739075679139134257347978784337,
1201  0.937273392400705904307758947710209,
1202  0.897264532344081900882509656454496,
1203  0.848206583410427216200648320774217,
1204  0.790418501442465932967649294817947,
1205  0.724417731360170047416186054613938,
1206  0.650996741297416970533735895313275,
1207  0.570972172608538847537226737253911,
1208  0.485081863640239680693655740232351,
1209  0.394151347077563369897207370981045,
1210  0.299180007153168812166780024266389,
1211  0.201194093997434522300628303394596,
1212  0.101142066918717499027074231447392,
1213  0.000000000000000000000000000000000
1214 };
1215 
1216 /* xgk[1], xgk[3], ... abscissae of the 15-point gauss rule.
1217  xgk[0], xgk[2], ... abscissae to optimally extend the 15-point gauss rule */
1218 
1219 static const double wgC[8] = /* weights of the 15-point gauss rule */
1220 {
1221  0.030753241996117268354628393577204,
1222  0.070366047488108124709267416450667,
1223  0.107159220467171935011869546685869,
1224  0.139570677926154314447804794511028,
1225  0.166269205816993933553200860481209,
1226  0.186161000015562211026800561866423,
1227  0.198431485327111576456118326443839,
1228  0.202578241925561272880620199967519
1229 };
1230 
1231 static const double wgkC[16] = /* weights of the 31-point kronrod rule */
1232 {
1233  0.005377479872923348987792051430128,
1234  0.015007947329316122538374763075807,
1235  0.025460847326715320186874001019653,
1236  0.035346360791375846222037948478360,
1237  0.044589751324764876608227299373280,
1238  0.053481524690928087265343147239430,
1239  0.062009567800670640285139230960803,
1240  0.069854121318728258709520077099147,
1241  0.076849680757720378894432777482659,
1242  0.083080502823133021038289247286104,
1243  0.088564443056211770647275443693774,
1244  0.093126598170825321225486872747346,
1245  0.096642726983623678505179907627589,
1246  0.099173598721791959332393173484603,
1247  0.100769845523875595044946662617570,
1248  0.101330007014791549017374792767493
1249 };
1250 
1251 void
1252 gsl_integration_qk31 (const gsl_function * f, double a, double b,
1253  double *result, double *abserr,
1254  double *resabs, double *resasc)
1255 {
1256  double fv1[16], fv2[16];
1257  // coverity[UNINIT_CTOR]
1258  gsl_integration_qk (16, xgkC, wgC, wgkC, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
1259 }
1260 
1261 /* Gauss quadrature weights and kronrod quadrature abscissae and
1262  weights as evaluated with 80 decimal digit arithmetic by
1263  L. W. Fullerton, Bell Labs, Nov. 1981. */
1264 
1265 static const double xgkD[21] = /* abscissae of the 41-point kronrod rule */
1266 {
1267  0.998859031588277663838315576545863,
1268  0.993128599185094924786122388471320,
1269  0.981507877450250259193342994720217,
1270  0.963971927277913791267666131197277,
1271  0.940822633831754753519982722212443,
1272  0.912234428251325905867752441203298,
1273  0.878276811252281976077442995113078,
1274  0.839116971822218823394529061701521,
1275  0.795041428837551198350638833272788,
1276  0.746331906460150792614305070355642,
1277  0.693237656334751384805490711845932,
1278  0.636053680726515025452836696226286,
1279  0.575140446819710315342946036586425,
1280  0.510867001950827098004364050955251,
1281  0.443593175238725103199992213492640,
1282  0.373706088715419560672548177024927,
1283  0.301627868114913004320555356858592,
1284  0.227785851141645078080496195368575,
1285  0.152605465240922675505220241022678,
1286  0.076526521133497333754640409398838,
1287  0.000000000000000000000000000000000
1288 };
1289 
1290 /* xgk[1], xgk[3], ... abscissae of the 20-point gauss rule.
1291  xgk[0], xgk[2], ... abscissae to optimally extend the 20-point gauss rule */
1292 
1293 static const double wgD[11] = /* weights of the 20-point gauss rule */
1294 {
1295  0.017614007139152118311861962351853,
1296  0.040601429800386941331039952274932,
1297  0.062672048334109063569506535187042,
1298  0.083276741576704748724758143222046,
1299  0.101930119817240435036750135480350,
1300  0.118194531961518417312377377711382,
1301  0.131688638449176626898494499748163,
1302  0.142096109318382051329298325067165,
1303  0.149172986472603746787828737001969,
1304  0.152753387130725850698084331955098
1305 };
1306 
1307 static const double wgkD[21] = /* weights of the 41-point kronrod rule */
1308 {
1309  0.003073583718520531501218293246031,
1310  0.008600269855642942198661787950102,
1311  0.014626169256971252983787960308868,
1312  0.020388373461266523598010231432755,
1313  0.025882133604951158834505067096153,
1314  0.031287306777032798958543119323801,
1315  0.036600169758200798030557240707211,
1316  0.041668873327973686263788305936895,
1317  0.046434821867497674720231880926108,
1318  0.050944573923728691932707670050345,
1319  0.055195105348285994744832372419777,
1320  0.059111400880639572374967220648594,
1321  0.062653237554781168025870122174255,
1322  0.065834597133618422111563556969398,
1323  0.068648672928521619345623411885368,
1324  0.071054423553444068305790361723210,
1325  0.073030690332786667495189417658913,
1326  0.074582875400499188986581418362488,
1327  0.075704497684556674659542775376617,
1328  0.076377867672080736705502835038061,
1329  0.076600711917999656445049901530102
1330 };
1331 
1332 void
1333 gsl_integration_qk41 (const gsl_function * f, double a, double b,
1334  double *result, double *abserr,
1335  double *resabs, double *resasc)
1336 {
1337  double fv1[21], fv2[21];
1338  // coverity[UNINIT]
1339  gsl_integration_qk (21, xgkD, wgD, wgkD, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
1340 }
1341 
1342 /* Gauss quadrature weights and kronrod quadrature abscissae and
1343  weights as evaluated with 80 decimal digit arithmetic by
1344  L. W. Fullerton, Bell Labs, Nov. 1981. */
1345 
1346 static const double xgkE[26] = /* abscissae of the 51-point kronrod rule */
1347 {
1348  0.999262104992609834193457486540341,
1349  0.995556969790498097908784946893902,
1350  0.988035794534077247637331014577406,
1351  0.976663921459517511498315386479594,
1352  0.961614986425842512418130033660167,
1353  0.942974571228974339414011169658471,
1354  0.920747115281701561746346084546331,
1355  0.894991997878275368851042006782805,
1356  0.865847065293275595448996969588340,
1357  0.833442628760834001421021108693570,
1358  0.797873797998500059410410904994307,
1359  0.759259263037357630577282865204361,
1360  0.717766406813084388186654079773298,
1361  0.673566368473468364485120633247622,
1362  0.626810099010317412788122681624518,
1363  0.577662930241222967723689841612654,
1364  0.526325284334719182599623778158010,
1365  0.473002731445714960522182115009192,
1366  0.417885382193037748851814394594572,
1367  0.361172305809387837735821730127641,
1368  0.303089538931107830167478909980339,
1369  0.243866883720988432045190362797452,
1370  0.183718939421048892015969888759528,
1371  0.122864692610710396387359818808037,
1372  0.061544483005685078886546392366797,
1373  0.000000000000000000000000000000000
1374 };
1375 
1376 /* xgk[1], xgk[3], ... abscissae of the 25-point gauss rule.
1377  xgk[0], xgk[2], ... abscissae to optimally extend the 25-point gauss rule */
1378 
1379 static const double wgE[13] = /* weights of the 25-point gauss rule */
1380 {
1381  0.011393798501026287947902964113235,
1382  0.026354986615032137261901815295299,
1383  0.040939156701306312655623487711646,
1384  0.054904695975835191925936891540473,
1385  0.068038333812356917207187185656708,
1386  0.080140700335001018013234959669111,
1387  0.091028261982963649811497220702892,
1388  0.100535949067050644202206890392686,
1389  0.108519624474263653116093957050117,
1390  0.114858259145711648339325545869556,
1391  0.119455763535784772228178126512901,
1392  0.122242442990310041688959518945852,
1393  0.123176053726715451203902873079050
1394 };
1395 
1396 static const double wgkE[26] = /* weights of the 51-point kronrod rule */
1397 {
1398  0.001987383892330315926507851882843,
1399  0.005561932135356713758040236901066,
1400  0.009473973386174151607207710523655,
1401  0.013236229195571674813656405846976,
1402  0.016847817709128298231516667536336,
1403  0.020435371145882835456568292235939,
1404  0.024009945606953216220092489164881,
1405  0.027475317587851737802948455517811,
1406  0.030792300167387488891109020215229,
1407  0.034002130274329337836748795229551,
1408  0.037116271483415543560330625367620,
1409  0.040083825504032382074839284467076,
1410  0.042872845020170049476895792439495,
1411  0.045502913049921788909870584752660,
1412  0.047982537138836713906392255756915,
1413  0.050277679080715671963325259433440,
1414  0.052362885806407475864366712137873,
1415  0.054251129888545490144543370459876,
1416  0.055950811220412317308240686382747,
1417  0.057437116361567832853582693939506,
1418  0.058689680022394207961974175856788,
1419  0.059720340324174059979099291932562,
1420  0.060539455376045862945360267517565,
1421  0.061128509717053048305859030416293,
1422  0.061471189871425316661544131965264,
1423  0.061580818067832935078759824240066
1424 };
1425 
1426 /* wgk[25] was calculated from the values of wgk[0..24] */
1427 
1428 void
1429 gsl_integration_qk51 (const gsl_function * f, double a, double b,
1430  double *result, double *abserr,
1431  double *resabs, double *resasc)
1432 {
1433  double fv1[26], fv2[26];
1434  //coverity[UNINIT]
1435  gsl_integration_qk (26, xgkE, wgE, wgkE, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
1436 }
1437 
1438 /* Gauss quadrature weights and kronrod quadrature abscissae and
1439  weights as evaluated with 80 decimal digit arithmetic by
1440  L. W. Fullerton, Bell Labs, Nov. 1981. */
1441 
1442 static const double xgkF[31] = /* abscissae of the 61-point kronrod rule */
1443 {
1444  0.999484410050490637571325895705811,
1445  0.996893484074649540271630050918695,
1446  0.991630996870404594858628366109486,
1447  0.983668123279747209970032581605663,
1448  0.973116322501126268374693868423707,
1449  0.960021864968307512216871025581798,
1450  0.944374444748559979415831324037439,
1451  0.926200047429274325879324277080474,
1452  0.905573307699907798546522558925958,
1453  0.882560535792052681543116462530226,
1454  0.857205233546061098958658510658944,
1455  0.829565762382768397442898119732502,
1456  0.799727835821839083013668942322683,
1457  0.767777432104826194917977340974503,
1458  0.733790062453226804726171131369528,
1459  0.697850494793315796932292388026640,
1460  0.660061064126626961370053668149271,
1461  0.620526182989242861140477556431189,
1462  0.579345235826361691756024932172540,
1463  0.536624148142019899264169793311073,
1464  0.492480467861778574993693061207709,
1465  0.447033769538089176780609900322854,
1466  0.400401254830394392535476211542661,
1467  0.352704725530878113471037207089374,
1468  0.304073202273625077372677107199257,
1469  0.254636926167889846439805129817805,
1470  0.204525116682309891438957671002025,
1471  0.153869913608583546963794672743256,
1472  0.102806937966737030147096751318001,
1473  0.051471842555317695833025213166723,
1474  0.000000000000000000000000000000000
1475 };
1476 
1477 /* xgk[1], xgk[3], ... abscissae of the 30-point gauss rule.
1478  xgk[0], xgk[2], ... abscissae to optimally extend the 30-point gauss rule */
1479 
1480 static const double wgF[15] = /* weights of the 30-point gauss rule */
1481 {
1482  0.007968192496166605615465883474674,
1483  0.018466468311090959142302131912047,
1484  0.028784707883323369349719179611292,
1485  0.038799192569627049596801936446348,
1486  0.048402672830594052902938140422808,
1487  0.057493156217619066481721689402056,
1488  0.065974229882180495128128515115962,
1489  0.073755974737705206268243850022191,
1490  0.080755895229420215354694938460530,
1491  0.086899787201082979802387530715126,
1492  0.092122522237786128717632707087619,
1493  0.096368737174644259639468626351810,
1494  0.099593420586795267062780282103569,
1495  0.101762389748405504596428952168554,
1496  0.102852652893558840341285636705415
1497 };
1498 
1499 static const double wgkF[31] = /* weights of the 61-point kronrod rule */
1500 {
1501  0.001389013698677007624551591226760,
1502  0.003890461127099884051267201844516,
1503  0.006630703915931292173319826369750,
1504  0.009273279659517763428441146892024,
1505  0.011823015253496341742232898853251,
1506  0.014369729507045804812451432443580,
1507  0.016920889189053272627572289420322,
1508  0.019414141193942381173408951050128,
1509  0.021828035821609192297167485738339,
1510  0.024191162078080601365686370725232,
1511  0.026509954882333101610601709335075,
1512  0.028754048765041292843978785354334,
1513  0.030907257562387762472884252943092,
1514  0.032981447057483726031814191016854,
1515  0.034979338028060024137499670731468,
1516  0.036882364651821229223911065617136,
1517  0.038678945624727592950348651532281,
1518  0.040374538951535959111995279752468,
1519  0.041969810215164246147147541285970,
1520  0.043452539701356069316831728117073,
1521  0.044814800133162663192355551616723,
1522  0.046059238271006988116271735559374,
1523  0.047185546569299153945261478181099,
1524  0.048185861757087129140779492298305,
1525  0.049055434555029778887528165367238,
1526  0.049795683427074206357811569379942,
1527  0.050405921402782346840893085653585,
1528  0.050881795898749606492297473049805,
1529  0.051221547849258772170656282604944,
1530  0.051426128537459025933862879215781,
1531  0.051494729429451567558340433647099
1532 };
1533 
1534 void
1535 gsl_integration_qk61 (const gsl_function * f, double a, double b,
1536  double *result, double *abserr,
1537  double *resabs, double *resasc)
1538 {
1539  double fv1[31], fv2[31];
1540  //coverity[UNINIT]
1541  gsl_integration_qk (31, xgkF, wgF, wgkF, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
1542 }
1543 
1544 gsl_integration_workspace*
1546 {
1547  gsl_integration_workspace* w ;
1548 
1549  if (n == 0)
1550  {
1551  GSL_ERROR_VAL ("workspace length n must be positive integer",
1552  GSL_EDOM, 0);
1553  }
1554 
1555  w = (gsl_integration_workspace *)
1556  malloc (sizeof (gsl_integration_workspace));
1557 
1558  if (w == 0)
1559  {
1560  GSL_ERROR_VAL ("failed to allocate space for workspace struct",
1561  GSL_ENOMEM, 0);
1562  }
1563 
1564  w->alist = (double *) malloc (n * sizeof (double));
1565 
1566  if (w->alist == 0)
1567  {
1568  free (w); /* exception in constructor, avoid memory leak */
1569 
1570  GSL_ERROR_VAL ("failed to allocate space for alist ranges",
1571  GSL_ENOMEM, 0);
1572  }
1573 
1574  w->blist = (double *) malloc (n * sizeof (double));
1575 
1576  if (w->blist == 0)
1577  {
1578  free (w->alist);
1579  free (w); /* exception in constructor, avoid memory leak */
1580 
1581  GSL_ERROR_VAL ("failed to allocate space for blist ranges",
1582  GSL_ENOMEM, 0);
1583  }
1584 
1585  w->rlist = (double *) malloc (n * sizeof (double));
1586 
1587  if (w->rlist == 0)
1588  {
1589  free (w->blist);
1590  free (w->alist);
1591  free (w); /* exception in constructor, avoid memory leak */
1592 
1593  GSL_ERROR_VAL ("failed to allocate space for rlist ranges",
1594  GSL_ENOMEM, 0);
1595  }
1596 
1597 
1598  w->elist = (double *) malloc (n * sizeof (double));
1599 
1600  if (w->elist == 0)
1601  {
1602  free (w->rlist);
1603  free (w->blist);
1604  free (w->alist);
1605  free (w); /* exception in constructor, avoid memory leak */
1606 
1607  GSL_ERROR_VAL ("failed to allocate space for elist ranges",
1608  GSL_ENOMEM, 0);
1609  }
1610 
1611  w->order = (size_t *) malloc (n * sizeof (size_t));
1612 
1613  if (w->order == 0)
1614  {
1615  free (w->elist);
1616  free (w->rlist);
1617  free (w->blist);
1618  free (w->alist);
1619  free (w); /* exception in constructor, avoid memory leak */
1620 
1621  GSL_ERROR_VAL ("failed to allocate space for order ranges",
1622  GSL_ENOMEM, 0);
1623  }
1624 
1625  w->level = (size_t *) malloc (n * sizeof (size_t));
1626 
1627  if (w->level == 0)
1628  {
1629  free (w->order);
1630  free (w->elist);
1631  free (w->rlist);
1632  free (w->blist);
1633  free (w->alist);
1634  free (w); /* exception in constructor, avoid memory leak */
1635 
1636  GSL_ERROR_VAL ("failed to allocate space for order ranges",
1637  GSL_ENOMEM, 0);
1638  }
1639 
1640  w->size = 0 ;
1641  w->limit = n ;
1642  w->maximum_level = 0 ;
1643 
1644  return w ;
1645 }
1646 
1647 void
1648 gsl_integration_workspace_free (gsl_integration_workspace * w)
1649 {
1650  free (w->level) ;
1651  free (w->order) ;
1652  free (w->elist) ;
1653  free (w->rlist) ;
1654  free (w->blist) ;
1655  free (w->alist) ;
1656  free (w) ;
1657 }
1658 
1659 
1660 
1661 // INCLUDED BELOW #include "reset.c"
1662 static inline void
1663 reset_nrmax (gsl_integration_workspace * workspace);
1664 
1665 static inline void
1666 reset_nrmax (gsl_integration_workspace * workspace)
1667 {
1668  workspace->nrmax = 0;
1669  workspace->i = workspace->order[0] ;
1670 }
1671 
1672 
1673 // INCLUDED BELOW #include "qpsrt2.c"
1674 /* The smallest interval has the largest error. Before bisecting
1675  decrease the sum of the errors over the larger intervals
1676  (error_over_large_intervals) and perform extrapolation. */
1677 
1678 static int
1679 increase_nrmax (gsl_integration_workspace * workspace);
1680 
1681 static int
1682 increase_nrmax (gsl_integration_workspace * workspace)
1683 {
1684  int k;
1685  int id = workspace->nrmax;
1686  int jupbnd;
1687 
1688  const size_t * level = workspace->level;
1689  const size_t * order = workspace->order;
1690 
1691  size_t limit = workspace->limit ;
1692  size_t last = workspace->size - 1 ;
1693 
1694  if (last > (1 + limit / 2))
1695  {
1696  jupbnd = limit + 1 - last;
1697  }
1698  else
1699  {
1700  jupbnd = last;
1701  }
1702 
1703  for (k = id; k <= jupbnd; k++)
1704  {
1705  size_t i_max = order[workspace->nrmax];
1706 
1707  workspace->i = i_max ;
1708 
1709  if (level[i_max] < workspace->maximum_level)
1710  {
1711  return 1;
1712  }
1713 
1714  workspace->nrmax++;
1715 
1716  }
1717  return 0;
1718 }
1719 
1720 static int
1721 large_interval (gsl_integration_workspace * workspace)
1722 {
1723  size_t i = workspace->i ;
1724  const size_t * level = workspace->level;
1725 
1726  if (level[i] < workspace->maximum_level)
1727  {
1728  return 1 ;
1729  }
1730  else
1731  {
1732  return 0 ;
1733  }
1734 }
1735 
1736 
1737 // INCLUDED BELOW #include "qelg.c"
1738 struct extrapolation_table
1739  {
1740  size_t n;
1741  double rlist2[52];
1742  size_t nres;
1743  double res3la[3];
1744  };
1745 
1746 static void
1747  initialise_table (struct extrapolation_table *table);
1748 
1749 static void
1750  append_table (struct extrapolation_table *table, double y);
1751 
1752 static void
1753 initialise_table (struct extrapolation_table *table)
1754 {
1755  table->n = 0;
1756  table->nres = 0;
1757 }
1758 #ifdef JUNK
1759 static void
1760 initialise_table (struct extrapolation_table *table, double y)
1761 {
1762  table->n = 0;
1763  table->rlist2[0] = y;
1764  table->nres = 0;
1765 }
1766 #endif
1767 static void
1768 append_table (struct extrapolation_table *table, double y)
1769 {
1770  size_t n;
1771  n = table->n;
1772  table->rlist2[n] = y;
1773  table->n++;
1774 }
1775 
1776 /* static inline void
1777  qelg (size_t * n, double epstab[],
1778  double * result, double * abserr,
1779  double res3la[], size_t * nres); */
1780 
1781 static inline void
1782  qelg (struct extrapolation_table *table, double *result, double *abserr);
1783 
1784 static inline void
1785 qelg (struct extrapolation_table *table, double *result, double *abserr)
1786 {
1787  double *epstab = table->rlist2;
1788  double *res3la = table->res3la;
1789  const size_t n = table->n - 1;
1790 
1791  const double current = epstab[n];
1792 
1793  double absolute = GSL_DBL_MAX;
1794  double relative = 5 * GSL_DBL_EPSILON * fabs (current);
1795 
1796  const size_t newelm = n / 2;
1797  const size_t n_orig = n;
1798  size_t n_final = n;
1799  size_t i;
1800 
1801  const size_t nres_orig = table->nres;
1802 
1803  *result = current;
1804  *abserr = GSL_DBL_MAX;
1805 
1806  if (n < 2)
1807  {
1808  *result = current;
1809  *abserr = GSL_MAX_DBL (absolute, relative);
1810  return;
1811  }
1812 
1813  epstab[n + 2] = epstab[n];
1814  epstab[n] = GSL_DBL_MAX;
1815 
1816  for (i = 0; i < newelm; i++)
1817  {
1818  double res = epstab[n - 2 * i + 2];
1819  double e0 = epstab[n - 2 * i - 2];
1820  double e1 = epstab[n - 2 * i - 1];
1821  double e2 = res;
1822 
1823  double e1abs = fabs (e1);
1824  double delta2 = e2 - e1;
1825  double err2 = fabs (delta2);
1826  double tol2 = GSL_MAX_DBL (fabs (e2), e1abs) * GSL_DBL_EPSILON;
1827  double delta3 = e1 - e0;
1828  double err3 = fabs (delta3);
1829  double tol3 = GSL_MAX_DBL (e1abs, fabs (e0)) * GSL_DBL_EPSILON;
1830 
1831  double e3, delta1, err1, tol1, ss;
1832 
1833  if (err2 <= tol2 && err3 <= tol3)
1834  {
1835  /* If e0, e1 and e2 are equal to within machine accuracy,
1836  convergence is assumed. */
1837 
1838  *result = res;
1839  absolute = err2 + err3;
1840  relative = 5 * GSL_DBL_EPSILON * fabs (res);
1841  *abserr = GSL_MAX_DBL (absolute, relative);
1842  return;
1843  }
1844 
1845  e3 = epstab[n - 2 * i];
1846  epstab[n - 2 * i] = e1;
1847  delta1 = e1 - e3;
1848  err1 = fabs (delta1);
1849  tol1 = GSL_MAX_DBL (e1abs, fabs (e3)) * GSL_DBL_EPSILON;
1850 
1851  /* If two elements are very close to each other, omit a part of
1852  the table by adjusting the value of n */
1853 
1854  if (err1 <= tol1 || err2 <= tol2 || err3 <= tol3)
1855  {
1856  n_final = 2 * i;
1857  break;
1858  }
1859 
1860  ss = (1 / delta1 + 1 / delta2) - 1 / delta3;
1861 
1862  /* Test to detect irregular behaviour in the table, and
1863  eventually omit a part of the table by adjusting the value of
1864  n. */
1865 
1866  if (fabs (ss * e1) <= 0.0001)
1867  {
1868  n_final = 2 * i;
1869  break;
1870  }
1871 
1872  /* Compute a new element and eventually adjust the value of
1873  result. */
1874 
1875  res = e1 + 1 / ss;
1876  epstab[n - 2 * i] = res;
1877 
1878  {
1879  const double error = err2 + fabs (res - e2) + err3;
1880 
1881  if (error <= *abserr)
1882  {
1883  *abserr = error;
1884  *result = res;
1885  }
1886  }
1887  }
1888 
1889  /* Shift the table */
1890 
1891  {
1892  const size_t limexp = 50 - 1;
1893 
1894  if (n_final == limexp)
1895  {
1896  n_final = 2 * (limexp / 2);
1897  }
1898  }
1899 
1900  if (n_orig % 2 == 1)
1901  {
1902  for (i = 0; i <= newelm; i++)
1903  {
1904  epstab[1 + i * 2] = epstab[i * 2 + 3];
1905  }
1906  }
1907  else
1908  {
1909  for (i = 0; i <= newelm; i++)
1910  {
1911  epstab[i * 2] = epstab[i * 2 + 2];
1912  }
1913  }
1914 
1915  if (n_orig != n_final)
1916  {
1917  for (i = 0; i <= n_final; i++)
1918  {
1919  epstab[i] = epstab[n_orig - n_final + i];
1920  }
1921  }
1922 
1923  table->n = n_final + 1;
1924 
1925  if (nres_orig < 3)
1926  {
1927  res3la[nres_orig] = *result;
1928  *abserr = GSL_DBL_MAX;
1929  }
1930  else
1931  { /* Compute error estimate */
1932  *abserr = (fabs (*result - res3la[2]) + fabs (*result - res3la[1])
1933  + fabs (*result - res3la[0]));
1934 
1935  res3la[0] = res3la[1];
1936  res3la[1] = res3la[2];
1937  res3la[2] = *result;
1938  }
1939 
1940  /* In QUADPACK the variable table->nres is incremented at the top of
1941  qelg, so it increases on every call. This leads to the array
1942  res3la being accessed when its elements are still undefined, so I
1943  have moved the update to this point so that its value more
1944  useful. */
1945 
1946  table->nres = nres_orig + 1;
1947 
1948  *abserr = GSL_MAX_DBL (*abserr, 5 * GSL_DBL_EPSILON * fabs (*result));
1949 
1950  return;
1951 }
1952 
1953 
1954 // INCLUDED BELOW #include "positivity.c"
1955 /* Compare the integral of f(x) with the integral of |f(x)|
1956  to determine if f(x) covers both positive and negative values */
1957 
1958 static inline int
1959 test_positivity (double result, double resabs);
1960 
1961 static inline int
1962 test_positivity (double result, double resabs)
1963 {
1964  int status = (fabs (result) >= (1 - 50 * GSL_DBL_EPSILON) * resabs);
1965 
1966  return status;
1967 }
1968 
1969 static int qags (const gsl_function * f, const double a, const double
1970  b, const double epsabs, const double epsrel, const size_t limit,
1971  gsl_integration_workspace * workspace, double *result, double *abserr,
1973 
1974 int
1976  double a, double b,
1977  double epsabs, double epsrel, size_t limit,
1978  gsl_integration_workspace * workspace,
1979  double * result, double * abserr)
1980 {
1981  int status = qags (f, a, b, epsabs, epsrel, limit,
1982  workspace,
1983  result, abserr,
1985  return status ;
1986 }
1987 
1988 /* QAGI: evaluate an integral over an infinite range using the
1989  transformation
1990 
1991  integrate(f(x),-Inf,Inf) = integrate((f((1-t)/t) + f(-(1-t)/t))/t^2,0,1)
1992 
1993  */
1994 
1995 static double i_transform (double t, void *params);
1996 
1997 int
1999  double epsabs, double epsrel, size_t limit,
2000  gsl_integration_workspace * workspace,
2001  double *result, double *abserr)
2002 {
2003  int status;
2004 
2005  gsl_function f_transform;
2006 
2007  f_transform.function = &i_transform;
2008  f_transform.params = f;
2009 
2010  status = qags (&f_transform, 0.0, 1.0,
2011  epsabs, epsrel, limit,
2012  workspace,
2013  result, abserr,
2015 
2016  return status;
2017 }
2018 
2019 static double
2020 i_transform (double t, void *params)
2021 {
2022  gsl_function *f = (gsl_function *) params;
2023  double x = (1 - t) / t;
2024  double y = GSL_FN_EVAL (f, x) + GSL_FN_EVAL (f, -x);
2025  return (y / t) / t;
2026 }
2027 
2028 
2029 /* QAGIL: Evaluate an integral over an infinite range using the
2030  transformation,
2031 
2032  integrate(f(x),-Inf,b) = integrate(f(b-(1-t)/t)/t^2,0,1)
2033 
2034  */
2035 
2036 struct il_params { double b ; gsl_function * f ; } ;
2037 
2038 static double il_transform (double t, void *params);
2039 
2040 int
2042  double b,
2043  double epsabs, double epsrel, size_t limit,
2044  gsl_integration_workspace * workspace,
2045  double *result, double *abserr)
2046 {
2047  int status;
2048 
2049  gsl_function f_transform;
2050  struct il_params transform_params ;
2051 
2052  transform_params.b = b ;
2053  transform_params.f = f ;
2054 
2055  f_transform.function = &il_transform;
2056  f_transform.params = &transform_params;
2057 
2058  status = qags (&f_transform, 0.0, 1.0,
2059  epsabs, epsrel, limit,
2060  workspace,
2061  result, abserr,
2063 
2064  return status;
2065 }
2066 
2067 static double
2068 il_transform (double t, void *params)
2069 {
2070  struct il_params *p = (struct il_params *) params;
2071  double b = p->b;
2072  gsl_function * f = p->f;
2073  double x = b - (1 - t) / t;
2074  double y = GSL_FN_EVAL (f, x);
2075  return (y / t) / t;
2076 }
2077 
2078 /* QAGIU: Evaluate an integral over an infinite range using the
2079  transformation
2080 
2081  integrate(f(x),a,Inf) = integrate(f(a+(1-t)/t)/t^2,0,1)
2082 
2083  */
2084 
2085 struct iu_params { double a ; gsl_function * f ; } ;
2086 
2087 static double iu_transform (double t, void *params);
2088 
2089 int
2091  double a,
2092  double epsabs, double epsrel, size_t limit,
2093  gsl_integration_workspace * workspace,
2094  double *result, double *abserr)
2095 {
2096  int status;
2097 
2098  gsl_function f_transform;
2099  struct iu_params transform_params ;
2100 
2101  transform_params.a = a ;
2102  transform_params.f = f ;
2103 
2104  f_transform.function = &iu_transform;
2105  f_transform.params = &transform_params;
2106 
2107  status = qags (&f_transform, 0.0, 1.0,
2108  epsabs, epsrel, limit,
2109  workspace,
2110  result, abserr,
2112 
2113  return status;
2114 }
2115 
2116 static double
2117 iu_transform (double t, void *params)
2118 {
2119  struct iu_params *p = (struct iu_params *) params;
2120  double a = p->a;
2121  gsl_function * f = p->f;
2122  double x = a + (1 - t) / t;
2123  double y = GSL_FN_EVAL (f, x);
2124  return (y / t) / t;
2125 }
2126 
2127 /* Main integration function */
2128 
2129 static int
2131  const double a, const double b,
2132  const double epsabs, const double epsrel,
2133  const size_t limit,
2134  gsl_integration_workspace * workspace,
2135  double *result, double *abserr,
2137 {
2138  double area, errsum;
2139  double res_ext, err_ext;
2140  double result0, abserr0, resabs0, resasc0;
2141  double tolerance;
2142 
2143  double ertest = 0;
2144  double error_over_large_intervals = 0;
2145  double reseps = 0, abseps = 0, correc = 0;
2146  size_t ktmin = 0;
2147  int roundoff_type1 = 0, roundoff_type2 = 0, roundoff_type3 = 0;
2148  int error_type = 0, error_type2 = 0;
2149 
2150  size_t iteration = 0;
2151 
2152  int positive_integrand = 0;
2153  int extrapolate = 0;
2154  int disallow_extrapolation = 0;
2155 
2156  struct extrapolation_table table;
2157 
2158  /* Initialize results */
2159 
2160  initialise (workspace, a, b);
2161 
2162  *result = 0;
2163  *abserr = 0;
2164 
2165  if (limit > workspace->limit)
2166  {
2167  GSL_ERROR ("iteration limit exceeds available workspace", GSL_EINVAL) ;
2168  }
2169 
2170  /* Test on accuracy */
2171 
2172  if (epsabs <= 0 && (epsrel < 50 * GSL_DBL_EPSILON || epsrel < 0.5e-28))
2173  {
2174  GSL_ERROR ("tolerance cannot be acheived with given epsabs and epsrel",
2175  GSL_EBADTOL);
2176  }
2177 
2178  /* Perform the first integration */
2179 
2180  q (f, a, b, &result0, &abserr0, &resabs0, &resasc0);
2181 
2182  set_initial_result (workspace, result0, abserr0);
2183 
2184  tolerance = GSL_MAX_DBL (epsabs, epsrel * fabs (result0));
2185 
2186  if (abserr0 <= 100 * GSL_DBL_EPSILON * resabs0 && abserr0 > tolerance)
2187  {
2188  *result = result0;
2189  *abserr = abserr0;
2190 
2191  GSL_ERROR ("cannot reach tolerance because of roundoff error"
2192  "on first attempt", GSL_EROUND);
2193  }
2194  else if ((abserr0 <= tolerance && abserr0 != resasc0) || abserr0 == 0.0)
2195  {
2196  *result = result0;
2197  *abserr = abserr0;
2198 
2199  return GSL_SUCCESS;
2200  }
2201  else if (limit == 1)
2202  {
2203  *result = result0;
2204  *abserr = abserr0;
2205 
2206  GSL_ERROR ("a maximum of one iteration was insufficient", GSL_EMAXITER);
2207  }
2208 
2209  /* Initialization */
2210 
2211  initialise_table (&table);
2212  append_table (&table, result0);
2213 
2214  area = result0;
2215  errsum = abserr0;
2216 
2217  res_ext = result0;
2218  err_ext = GSL_DBL_MAX;
2219 
2220  positive_integrand = test_positivity (result0, resabs0);
2221 
2222  iteration = 1;
2223 
2224  do
2225  {
2226  size_t current_level;
2227  double a1, b1, a2, b2;
2228  double a_i, b_i, r_i, e_i;
2229  double area1 = 0, area2 = 0, area12 = 0;
2230  double error1 = 0, error2 = 0, error12 = 0;
2231  double resasc1, resasc2;
2232  double resabs1, resabs2;
2233  double last_e_i;
2234 
2235  /* Bisect the subinterval with the largest error estimate */
2236 
2237  retrieve (workspace, &a_i, &b_i, &r_i, &e_i);
2238 
2239  current_level = workspace->level[workspace->i] + 1;
2240 
2241  a1 = a_i;
2242  b1 = 0.5 * (a_i + b_i);
2243  a2 = b1;
2244  b2 = b_i;
2245 
2246  iteration++;
2247 
2248  q (f, a1, b1, &area1, &error1, &resabs1, &resasc1);
2249  q (f, a2, b2, &area2, &error2, &resabs2, &resasc2);
2250 
2251  area12 = area1 + area2;
2252  error12 = error1 + error2;
2253  last_e_i = e_i;
2254 
2255  /* Improve previous approximations to the integral and test for
2256  accuracy.
2257 
2258  We write these expressions in the same way as the original
2259  QUADPACK code so that the rounding errors are the same, which
2260  makes testing easier. */
2261 
2262  errsum = errsum + error12 - e_i;
2263  area = area + area12 - r_i;
2264 
2265  tolerance = GSL_MAX_DBL (epsabs, epsrel * fabs (area));
2266 
2267  if (resasc1 != error1 && resasc2 != error2)
2268  {
2269  double delta = r_i - area12;
2270 
2271  if (fabs (delta) <= 1.0e-5 * fabs (area12) && error12 >= 0.99 * e_i)
2272  {
2273  if (!extrapolate)
2274  {
2275  roundoff_type1++;
2276  }
2277  else
2278  {
2279  roundoff_type2++;
2280  }
2281  }
2282  if (iteration > 10 && error12 > e_i)
2283  {
2284  roundoff_type3++;
2285  }
2286  }
2287 
2288  /* Test for roundoff and eventually set error flag */
2289 
2290  if (roundoff_type1 + roundoff_type2 >= 10 || roundoff_type3 >= 20)
2291  {
2292  error_type = 2; /* round off error */
2293  }
2294 
2295  if (roundoff_type2 >= 5)
2296  {
2297  error_type2 = 1;
2298  }
2299 
2300  /* set error flag in the case of bad integrand behaviour at
2301  a point of the integration range */
2302 
2303  if (subinterval_too_small (a1, a2, b2))
2304  {
2305  error_type = 4;
2306  }
2307 
2308  /* append the newly-created intervals to the list */
2309 
2310  update (workspace, a1, b1, area1, error1, a2, b2, area2, error2);
2311 
2312  if (errsum <= tolerance)
2313  {
2314  goto compute_result;
2315  }
2316 
2317  if (error_type)
2318  {
2319  break;
2320  }
2321 
2322  if (iteration >= limit - 1)
2323  {
2324  error_type = 1;
2325  break;
2326  }
2327 
2328  if (iteration == 2) /* set up variables on first iteration */
2329  {
2330  error_over_large_intervals = errsum;
2331  ertest = tolerance;
2332  append_table (&table, area);
2333  continue;
2334  }
2335 
2336  if (disallow_extrapolation)
2337  {
2338  continue;
2339  }
2340 
2341  error_over_large_intervals += -last_e_i;
2342 
2343  if (current_level < workspace->maximum_level)
2344  {
2345  error_over_large_intervals += error12;
2346  }
2347 
2348  if (!extrapolate)
2349  {
2350  /* test whether the interval to be bisected next is the
2351  smallest interval. */
2352 
2353  if (large_interval (workspace))
2354  continue;
2355 
2356  extrapolate = 1;
2357  workspace->nrmax = 1;
2358  }
2359 
2360  if (!error_type2 && error_over_large_intervals > ertest)
2361  {
2362  if (increase_nrmax (workspace))
2363  continue;
2364  }
2365 
2366  /* Perform extrapolation */
2367 
2368  append_table (&table, area);
2369 
2370  qelg (&table, &reseps, &abseps);
2371 
2372  ktmin++;
2373 
2374  if (ktmin > 5 && err_ext < 0.001 * errsum)
2375  {
2376  error_type = 5;
2377  }
2378 
2379  if (abseps < err_ext)
2380  {
2381  ktmin = 0;
2382  err_ext = abseps;
2383  res_ext = reseps;
2384  correc = error_over_large_intervals;
2385  ertest = GSL_MAX_DBL (epsabs, epsrel * fabs (reseps));
2386  if (err_ext <= ertest)
2387  break;
2388  }
2389 
2390  /* Prepare bisection of the smallest interval. */
2391 
2392  if (table.n == 1)
2393  {
2394  disallow_extrapolation = 1;
2395  }
2396 
2397  if (error_type == 5)
2398  {
2399  break;
2400  }
2401 
2402  /* work on interval with largest error */
2403 
2404  reset_nrmax (workspace);
2405  extrapolate = 0;
2406  error_over_large_intervals = errsum;
2407 
2408  }
2409  while (iteration < limit);
2410 
2411  *result = res_ext;
2412  *abserr = err_ext;
2413 
2414  if (err_ext == GSL_DBL_MAX)
2415  goto compute_result;
2416 
2417  if (error_type || error_type2)
2418  {
2419  if (error_type2)
2420  {
2421  err_ext += correc;
2422  }
2423 
2424 // if (error_type == 0)
2425 // error_type = 3;
2426 
2427  if (res_ext != 0.0 && area != 0.0)
2428  {
2429  if (err_ext / fabs (res_ext) > errsum / fabs (area))
2430  goto compute_result;
2431  }
2432  else if (err_ext > errsum)
2433  {
2434  goto compute_result;
2435  }
2436  else if (area == 0.0)
2437  {
2438  goto return_error;
2439  }
2440  }
2441 
2442  /* Test on divergence. */
2443 
2444  {
2445  double max_area = GSL_MAX_DBL (fabs (res_ext), fabs (area));
2446 
2447  if (!positive_integrand && max_area < 0.01 * resabs0)
2448  goto return_error;
2449  }
2450 
2451  {
2452  double ratio = res_ext / area;
2453 
2454  if (ratio < 0.01 || ratio > 100.0 || errsum > fabs (area))
2455  error_type = 6;
2456  }
2457 
2458  goto return_error;
2459 
2460 compute_result:
2461 
2462  *result = sum_results (workspace);
2463  *abserr = errsum;
2464 
2465 return_error:
2466 
2467  if (error_type > 2)
2468  error_type--;
2469 
2470 
2471 
2472  if (error_type == 0)
2473  {
2474  return GSL_SUCCESS;
2475  }
2476  else if (error_type == 1)
2477  {
2478  GSL_ERROR ("number of iterations was insufficient", GSL_EMAXITER);
2479  }
2480  else if (error_type == 2)
2481  {
2482  GSL_ERROR ("cannot reach tolerance because of roundoff error",
2483  GSL_EROUND);
2484  }
2485  else if (error_type == 3)
2486  {
2487  GSL_ERROR ("bad integrand behavior found in the integration interval",
2488  GSL_ESING);
2489  }
2490  else if (error_type == 4)
2491  {
2492  GSL_ERROR ("roundoff error detected in the extrapolation table",
2493  GSL_EROUND);
2494  }
2495  else if (error_type == 5)
2496  {
2497  GSL_ERROR ("integral is divergent, or slowly convergent",
2498  GSL_EDIVERGE);
2499  }
2500 
2501  GSL_ERROR ("could not integrate function", GSL_EFAILED);
2502 }
static double sum_results(const gsl_integration_workspace *workspace)
static const double xgkC[16]
#define coutE(a)
Definition: RooMsgService.h:35
float xmin
Definition: THbookFile.cxx:93
static double il_transform(double t, void *params)
double gsl_coerce_double(const double x)
static const double wgF[15]
int gsl_integration_qagiu(gsl_function *f, double a, double epsabs, double epsrel, size_t limit, gsl_integration_workspace *workspace, double *result, double *abserr)
virtual Bool_t setIndex(Int_t index, Bool_t printError=kTRUE)
Set value by specifying the index code of the desired state.
#define assert(cond)
Definition: unittest.h:542
double RooAdaptiveGaussKronrodIntegrator1D_GSL_GlueFunction(double x, void *data)
Glue function interacing to GSL code.
static const double wgkB[11]
TAlienJobStatus * status
Definition: TAlienJob.cxx:51
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
static const double wgkD[21]
const Bool_t kFALSE
Definition: Rtypes.h:92
static Int_t isInfinite(Double_t x)
Return true if x is infinite by RooNumBer internal specification.
Definition: RooNumber.cxx:57
STL namespace.
UInt_t getDimension() const
Definition: RooAbsFunc.h:29
void gsl_integration_qk51(const gsl_function *f, double a, double b, double *result, double *abserr, double *resabs, double *resasc)
int gsl_integration_qag(const gsl_function *f, double a, double b, double epsabs, double epsrel, size_t limit, int key, gsl_integration_workspace *workspace, double *result, double *abserr)
static void qelg(struct extrapolation_table *table, double *result, double *abserr)
static void retrieve(const gsl_integration_workspace *workspace, double *a, double *b, double *r, double *e)
static void append_table(struct extrapolation_table *table, double y)
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
static const double wgkC[16]
static const double wgkA[8]
Double_t x[n]
Definition: legend1.C:17
Vc_ALWAYS_INLINE void free(T *p)
Frees memory that was allocated with Vc::malloc.
Definition: memory.h:94
static const double wgA[4]
void gsl_integration_rule(const gsl_function *f, double a, double b, double *result, double *abserr, double *defabs, double *resabs)
int gsl_integration_qags(const gsl_function *f, double a, double b, double epsabs, double epsrel, size_t limit, gsl_integration_workspace *workspace, double *result, double *abserr)
#define GSL_ERROR(a, b)
#define GSL_FN_EVAL(F, x)
static int subinterval_too_small(double a1, double a2, double b2)
static int increase_nrmax(gsl_integration_workspace *workspace)
#define GSL_ERROR_VAL(reason, gsl_errno, value)
#define F(x, y, z)
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
void function(const char *name_, T fun, const char *docstring=0)
Definition: RExports.h:159
ROOT::R::TRInterface & r
Definition: Object.C:4
return
Definition: TBase64.cxx:62
static void reset_nrmax(gsl_integration_workspace *workspace)
static int large_interval(gsl_integration_workspace *workspace)
void gsl_integration_workspace_free(gsl_integration_workspace *w)
Int_t getCatIndex(const char *name, Int_t defVal=0, Bool_t verbose=kFALSE) const
Get index value of a RooAbsCategory stored in set with given name.
Definition: RooArgSet.cxx:613
void gsl_integration_qcheb(gsl_function *f, double a, double b, double *cheb12, double *cheb24)
TClass * IsA() const
static void update(gsl_integration_workspace *workspace, double a1, double b1, double area1, double error1, double a2, double b2, double area2, double error2)
static void qpsrt(gsl_integration_workspace *workspace)
virtual Bool_t checkLimits() const
Check that our integration range is finite and otherwise return kFALSE.
unsigned int UInt_t
Definition: RtypesCore.h:42
const RooAbsFunc * _function
static const double wgkE[26]
void gsl_integration_qk(const int n, const double xgk[], const double wg[], const double wgk[], double fv1[], double fv2[], const gsl_function *f, double a, double b, double *result, double *abserr, double *resabs, double *resasc)
virtual Double_t getMinLimit(UInt_t dimension) const =0
Bool_t setLimits(Double_t *xmin, Double_t *xmax)
Change our integration limits.
void gsl_integration_qk41(const gsl_function *f, double a, double b, double *result, double *abserr, double *resabs, double *resasc)
void gsl_integration_qk61(const gsl_function *f, double a, double b, double *result, double *abserr, double *resabs, double *resasc)
float xmax
Definition: THbookFile.cxx:93
static void initialise_table(struct extrapolation_table *table)
Double_t integrand(const Double_t x[]) const
double GSL_MAX_DBL(double a, double b)
static void initialise(gsl_integration_workspace *workspace, double a, double b)
static void set_initial_result(gsl_integration_workspace *workspace, double result, double error)
static const double wgB[5]
double f(double x)
struct gsl_function_struct gsl_function
virtual Double_t getMaxLimit(UInt_t dimension) const =0
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:415
double Double_t
Definition: RtypesCore.h:55
static double iu_transform(double t, void *params)
static int test_positivity(double result, double resabs)
Double_t y[n]
Definition: legend1.C:17
static int qag(const gsl_function *f, const double a, const double b, const double epsabs, const double epsrel, const size_t limit, gsl_integration_workspace *workspace, double *result, double *abserr, gsl_integration_rule *q)
virtual RooAbsIntegrator * clone(const RooAbsFunc &function, const RooNumIntConfig &config) const
Virtual constructor.
virtual Double_t integral(const Double_t *yvec=0)
Calculate and return integral at at given parameter values.
static const double wgE[13]
int gsl_integration_qagil(gsl_function *f, double b, double epsabs, double epsrel, size_t limit, gsl_integration_workspace *workspace, double *result, double *abserr)
static const double wgC[8]
static const double xgkD[21]
static int qags(const gsl_function *f, const double a, const double b, const double epsabs, const double epsrel, const size_t limit, gsl_integration_workspace *workspace, double *result, double *abserr, gsl_integration_rule *q)
ClassImp(RooAdaptiveGaussKronrodIntegrator1D)
static const double xgkE[26]
const RooAbsFunc * integrand() const
static const double xgkB[11]
void gsl_integration_qk21(const gsl_function *f, double a, double b, double *result, double *abserr, double *resabs, double *resasc)
static void registerIntegrator(RooNumIntFactory &fact)
Register this class with RooNumIntConfig as a possible choice of numeric integrator for one-dimension...
void gsl_integration_qk31(const gsl_function *f, double a, double b, double *result, double *abserr, double *resabs, double *resasc)
friend double RooAdaptiveGaussKronrodIntegrator1D_GSL_GlueFunction(double x, void *data)
Glue function interacing to GSL code.
static const double xgkF[31]
int gsl_integration_qagi(gsl_function *f, double epsabs, double epsrel, size_t limit, gsl_integration_workspace *workspace, double *result, double *abserr)
Bool_t defineType(const char *label)
Define a state with given name, the lowest available positive integer is assigned as index...
static const double xgkA[8]
double result[121]
static double i_transform(double t, void *params)
void gsl_integration_qk15(const gsl_function *f, double a, double b, double *result, double *abserr, double *resabs, double *resasc)
Bool_t isValid() const
static const double wgD[11]
const Bool_t kTRUE
Definition: Rtypes.h:91
float * q
Definition: THbookFile.cxx:87
Vc_ALWAYS_INLINE_L T *Vc_ALWAYS_INLINE_R malloc(size_t n)
Allocates memory on the Heap with alignment and padding suitable for vectorized access.
Definition: memory.h:67
const Int_t n
Definition: legend1.C:16
const RooArgSet & getConfigSection(const char *name) const
Retrieve configuration information specific to integrator with given name.
gsl_integration_workspace * gsl_integration_workspace_alloc(const size_t n)
static const double wgkF[31]
static double rescale_error(double err, const double result_abs, const double result_asc)
Bool_t storeProtoIntegrator(RooAbsIntegrator *proto, const RooArgSet &defConfig, const char *depName="")
Method accepting registration of a prototype numeric integrator along with a RooArgSet of its default...
Double_t getRealValue(const char *name, Double_t defVal=0, Bool_t verbose=kFALSE) const
Get value of a RooAbsReal stored in set with given name.
Definition: RooArgSet.cxx:527
RooAdaptiveGaussKronrodIntegrator1D()
coverity[UNINIT_CTOR] Default constructor
Bool_t initialize()
Initialize integrator allocate buffers and setup GSL workspace.