69 struct gsl_function_struct
71 double (*
function) (
double x,
void * params);
75 #define GSL_FN_EVAL(F,x) (*((F)->function))(x,(F)->params) 92 gsl_integration_workspace;
94 gsl_integration_workspace *
102 double epsabs,
double epsrel,
size_t limit,
104 gsl_integration_workspace * workspace,
105 double *
result,
double *abserr);
110 double epsabs,
double epsrel,
size_t limit,
111 gsl_integration_workspace * workspace,
112 double *
result,
double * abserr) ;
116 double epsabs,
double epsrel,
size_t limit,
117 gsl_integration_workspace * workspace,
118 double *
result,
double *abserr) ;
123 double epsabs,
double epsrel,
size_t limit,
124 gsl_integration_workspace * workspace,
125 double *
result,
double *abserr) ;
130 double epsabs,
double epsrel,
size_t limit,
131 gsl_integration_workspace * workspace,
132 double *
result,
double *abserr) ;
144 RooRealVar maxSeg(
"maxSeg",
"maximum number of segments",100) ;
145 RooCategory method(
"method",
"Integration method for each segment") ;
262 coutE(
Integration) <<
"RooAdaptiveGaussKronrodIntegrator1D::setLimits: cannot override integrand's limits" << endl;
289 if (!infLo && !infHi) {
291 }
else if (infLo && infHi) {
293 }
else if (infLo && !infHi) {
385 #define GSL_SUCCESS 0 388 #define GSL_EBADTOL 13 390 #define GSL_ERROR(a,b) oocoutE((TObject*)0,Integration) << "RooAdaptiveGaussKronrodIntegrator1D::integral() ERROR: " << a << endl ; return b ; 391 #define GSL_DBL_MIN 2.2250738585072014e-308 392 #define GSL_DBL_MAX 1.7976931348623157e+308 393 #define GSL_DBL_EPSILON 2.2204460492503131e-16 396 #define GSL_EMAXITER 3 398 #define GSL_EFAILED 5 399 #define GSL_EDIVERGE 6 402 #define GSL_ERROR_VAL(reason, gsl_errno, value) return value ; 404 #define GSL_MAX(a,b) ((a) > (b) ? (a) : (b)) 420 #define GSL_COERCE_DBL(x) (gsl_coerce_double(x)) 433 double *
result,
double *abserr,
434 double *defabs,
double *resabs);
437 double *
result,
double *abserr,
438 double *resabs,
double *resasc);
441 double *
result,
double *abserr,
442 double *resabs,
double *resasc);
445 double *
result,
double *abserr,
446 double *resabs,
double *resasc);
449 double *
result,
double *abserr,
450 double *resabs,
double *resasc);
453 double *
result,
double *abserr,
454 double *resabs,
double *resasc);
457 double *
result,
double *abserr,
458 double *resabs,
double *resasc);
461 double *cheb12,
double *cheb24);
479 const double wg[],
const double wgk[],
480 double fv1[],
double fv2[],
482 double *
result,
double * abserr,
483 double * resabs,
double * resasc);
488 double epsabs,
double epsrel,
size_t limit,
490 gsl_integration_workspace * workspace,
491 double *
result,
double *abserr);
496 void initialise (gsl_integration_workspace * workspace,
double a,
double b);
499 void initialise (gsl_integration_workspace * workspace,
double a,
double b)
502 workspace->nrmax = 0;
504 workspace->alist[0] =
a;
505 workspace->blist[0] =
b;
506 workspace->rlist[0] = 0.0;
507 workspace->elist[0] = 0.0;
508 workspace->order[0] = 0;
509 workspace->level[0] = 0;
511 workspace->maximum_level = 0;
517 double result,
double error);
521 double result,
double error)
524 workspace->rlist[0] =
result;
525 workspace->elist[0] = error;
530 qpsrt (gsl_integration_workspace * workspace);
533 void qpsrt (gsl_integration_workspace * workspace)
535 const size_t last = workspace->size - 1;
536 const size_t limit = workspace->limit;
538 double * elist = workspace->elist;
539 size_t * order = workspace->order;
545 size_t i_nrmax = workspace->nrmax;
546 size_t i_maxerr = order[i_nrmax] ;
554 workspace->i = i_maxerr ;
558 errmax = elist[i_maxerr] ;
565 while (i_nrmax > 0 && errmax > elist[order[i_nrmax - 1]])
567 order[i_nrmax] = order[i_nrmax - 1] ;
575 if(last < (limit/2 + 2))
581 top = limit - last + 1;
592 while (i < top && errmax < elist[order[i]])
594 order[i-1] = order[i] ;
598 order[i-1] = i_maxerr ;
602 errmin = elist[last] ;
606 while (k > i - 2 && errmin >= elist[order[k]])
608 order[k+1] = order[k] ;
616 i_maxerr = order[i_nrmax] ;
618 workspace->i = i_maxerr ;
619 workspace->nrmax = i_nrmax ;
626 void update (gsl_integration_workspace * workspace,
627 double a1,
double b1,
double area1,
double error1,
628 double a2,
double b2,
double area2,
double error2);
631 retrieve (
const gsl_integration_workspace * workspace,
632 double *
a,
double *
b,
double *
r,
double *
e);
637 void update (gsl_integration_workspace * workspace,
638 double a1,
double b1,
double area1,
double error1,
639 double a2,
double b2,
double area2,
double error2)
641 double * alist = workspace->alist ;
642 double * blist = workspace->blist ;
643 double * rlist = workspace->rlist ;
644 double * elist = workspace->elist ;
645 size_t * level = workspace->level ;
647 const size_t i_max = workspace->i ;
648 const size_t i_new = workspace->size ;
650 const size_t new_level = workspace->level[i_max] + 1;
657 rlist[i_max] = area2;
658 elist[i_max] = error2;
659 level[i_max] = new_level;
663 rlist[i_new] = area1;
664 elist[i_new] = error1;
665 level[i_new] = new_level;
670 rlist[i_max] = area1;
671 elist[i_max] = error1;
672 level[i_max] = new_level;
676 rlist[i_new] = area2;
677 elist[i_new] = error2;
678 level[i_new] = new_level;
683 if (new_level > workspace->maximum_level)
685 workspace->maximum_level = new_level;
692 retrieve (
const gsl_integration_workspace * workspace,
693 double *
a,
double *
b,
double *
r,
double *
e)
695 const size_t i = workspace->i;
696 double * alist = workspace->alist;
697 double * blist = workspace->blist;
698 double * rlist = workspace->rlist;
699 double * elist = workspace->elist;
708 sum_results (
const gsl_integration_workspace * workspace);
713 const double *
const rlist = workspace->rlist ;
714 const size_t n = workspace->size;
717 double result_sum = 0;
719 for (k = 0; k <
n; k++)
721 result_sum += rlist[k];
736 double tmp = (1 + 100 *
e) * (
fabs (a2) + 1000 * u);
738 int status =
fabs (a1) <= tmp &&
fabs (b2) <= tmp;
746 const double a,
const double b,
747 const double epsabs,
const double epsrel,
749 gsl_integration_workspace * workspace,
750 double *
result,
double * abserr,
756 double epsabs,
double epsrel,
size_t limit,
758 gsl_integration_workspace * workspace,
759 double *
result,
double * abserr)
795 status =
qag (f, a, b, epsabs, epsrel, limit,
805 const double a,
const double b,
806 const double epsabs,
const double epsrel,
808 gsl_integration_workspace * workspace,
809 double *
result,
double *abserr,
813 double result0, abserr0, resabs0, resasc0;
815 size_t iteration = 0;
816 int roundoff_type1 = 0, roundoff_type2 = 0, error_type = 0;
827 if (limit > workspace->limit)
834 GSL_ERROR (
"tolerance cannot be acheived with given epsabs and epsrel",
840 q (f, a, b, &result0, &abserr0, &resabs0, &resasc0);
852 if (abserr0 <= round_off && abserr0 > tolerance)
857 GSL_ERROR (
"cannot reach tolerance because of roundoff error " 860 else if ((abserr0 <= tolerance && abserr0 != resasc0) || abserr0 == 0.0)
882 double a1, b1, a2, b2;
883 double a_i, b_i, r_i, e_i;
884 double area1 = 0, area2 = 0, area12 = 0;
885 double error1 = 0, error2 = 0, error12 = 0;
886 double resasc1, resasc2;
887 double resabs1, resabs2;
891 retrieve (workspace, &a_i, &b_i, &r_i, &e_i);
894 b1 = 0.5 * (a_i + b_i);
898 q (f, a1, b1, &area1, &error1, &resabs1, &resasc1);
899 q (f, a2, b2, &area2, &error2, &resabs2, &resasc2);
901 area12 = area1 + area2;
902 error12 = error1 + error2;
904 errsum += (error12 - e_i);
905 area += area12 - r_i;
907 if (resasc1 != error1 && resasc2 != error2)
909 double delta = r_i - area12;
911 if (
fabs (delta) <= 1.0e-5 *
fabs (area12) && error12 >= 0.99 * e_i)
915 if (iteration >= 10 && error12 > e_i)
923 if (errsum > tolerance)
925 if (roundoff_type1 >= 6 || roundoff_type2 >= 20)
939 update (workspace, a1, b1, area1, error1, a2, b2, area2, error2);
941 retrieve (workspace, &a_i, &b_i, &r_i, &e_i);
946 while (iteration < limit && !error_type && errsum > tolerance);
951 if (errsum <= tolerance)
955 else if (error_type == 2)
957 GSL_ERROR (
"roundoff error prevents tolerance from being achieved",
960 else if (error_type == 3)
962 GSL_ERROR (
"bad integrand behavior found in the integration interval",
965 else if (iteration == limit)
975 static double rescale_error (
double err,
const double result_abs,
const double result_asc) ;
978 rescale_error (
double err,
const double result_abs,
const double result_asc)
982 if (result_asc != 0 && err != 0)
984 double scale =
TMath::Power((200 * err / result_asc), 1.5) ;
988 err = result_asc * scale ;
1011 const double xgk[],
const double wg[],
const double wgk[],
1012 double fv1[],
double fv2[],
1014 double *
result,
double *abserr,
1015 double *resabs,
double *resasc)
1018 const double center = 0.5 * (a +
b);
1019 const double half_length = 0.5 * (b -
a);
1020 const double abs_half_length =
fabs (half_length);
1023 double result_gauss = 0;
1024 double result_kronrod = f_center * wgk[n - 1];
1026 double result_abs =
fabs (result_kronrod);
1027 double result_asc = 0;
1028 double mean = 0, err = 0;
1034 result_gauss = f_center * wg[n / 2 - 1];
1037 for (j = 0; j < (n - 1) / 2; j++)
1039 const int jtw = j * 2 + 1;
1040 const double abscissa = half_length * xgk[jtw];
1041 const double fval1 =
GSL_FN_EVAL (f, center - abscissa);
1042 const double fval2 =
GSL_FN_EVAL (f, center + abscissa);
1043 const double fsum = fval1 + fval2;
1046 result_gauss += wg[j] * fsum;
1047 result_kronrod += wgk[jtw] * fsum;
1048 result_abs += wgk[jtw] * (
fabs (fval1) +
fabs (fval2));
1051 for (j = 0; j < n / 2; j++)
1054 const double abscissa = half_length * xgk[jtwm1];
1055 const double fval1 =
GSL_FN_EVAL (f, center - abscissa);
1056 const double fval2 =
GSL_FN_EVAL (f, center + abscissa);
1059 result_kronrod += wgk[jtwm1] * (fval1 + fval2);
1060 result_abs += wgk[jtwm1] * (
fabs (fval1) +
fabs (fval2));
1063 mean = result_kronrod * 0.5;
1065 result_asc = wgk[n - 1] *
fabs (f_center - mean);
1067 for (j = 0; j < n - 1; j++)
1069 result_asc += wgk[j] * (
fabs (fv1[j] - mean) +
fabs (fv2[j] - mean));
1074 err = (result_kronrod - result_gauss) * half_length;
1076 result_kronrod *= half_length;
1077 result_abs *= abs_half_length;
1078 result_asc *= abs_half_length;
1080 *result = result_kronrod;
1081 *resabs = result_abs;
1082 *resasc = result_asc;
1093 0.991455371120812639206854697526329,
1094 0.949107912342758524526189684047851,
1095 0.864864423359769072789712788640926,
1096 0.741531185599394439863864773280788,
1097 0.586087235467691130294144838258730,
1098 0.405845151377397166906606412076961,
1099 0.207784955007898467600689403773245,
1100 0.000000000000000000000000000000000
1108 0.129484966168869693270611432679082,
1109 0.279705391489276667901467771423780,
1110 0.381830050505118944950369775488975,
1111 0.417959183673469387755102040816327
1116 0.022935322010529224963732008058970,
1117 0.063092092629978553290700663189204,
1118 0.104790010322250183839876322541518,
1119 0.140653259715525918745189590510238,
1120 0.169004726639267902826583426598550,
1121 0.190350578064785409913256402421014,
1122 0.204432940075298892414161999234649,
1123 0.209482141084727828012999174891714
1128 double *
result,
double *abserr,
1129 double *resabs,
double *resasc)
1131 double fv1[8], fv2[8];
1133 gsl_integration_qk (8,
xgkA,
wgA,
wgkA, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
1142 0.995657163025808080735527280689003,
1143 0.973906528517171720077964012084452,
1144 0.930157491355708226001207180059508,
1145 0.865063366688984510732096688423493,
1146 0.780817726586416897063717578345042,
1147 0.679409568299024406234327365114874,
1148 0.562757134668604683339000099272694,
1149 0.433395394129247190799265943165784,
1150 0.294392862701460198131126603103866,
1151 0.148874338981631210884826001129720,
1152 0.000000000000000000000000000000000
1160 0.066671344308688137593568809893332,
1161 0.149451349150580593145776339657697,
1162 0.219086362515982043995534934228163,
1163 0.269266719309996355091226921569469,
1164 0.295524224714752870173892994651338
1169 0.011694638867371874278064396062192,
1170 0.032558162307964727478818972459390,
1171 0.054755896574351996031381300244580,
1172 0.075039674810919952767043140916190,
1173 0.093125454583697605535065465083366,
1174 0.109387158802297641899210590325805,
1175 0.123491976262065851077958109831074,
1176 0.134709217311473325928054001771707,
1177 0.142775938577060080797094273138717,
1178 0.147739104901338491374841515972068,
1179 0.149445554002916905664936468389821
1185 double *
result,
double *abserr,
1186 double *resabs,
double *resasc)
1188 double fv1[11], fv2[11];
1190 gsl_integration_qk (11,
xgkB,
wgB,
wgkB, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
1199 0.998002298693397060285172840152271,
1200 0.987992518020485428489565718586613,
1201 0.967739075679139134257347978784337,
1202 0.937273392400705904307758947710209,
1203 0.897264532344081900882509656454496,
1204 0.848206583410427216200648320774217,
1205 0.790418501442465932967649294817947,
1206 0.724417731360170047416186054613938,
1207 0.650996741297416970533735895313275,
1208 0.570972172608538847537226737253911,
1209 0.485081863640239680693655740232351,
1210 0.394151347077563369897207370981045,
1211 0.299180007153168812166780024266389,
1212 0.201194093997434522300628303394596,
1213 0.101142066918717499027074231447392,
1214 0.000000000000000000000000000000000
1222 0.030753241996117268354628393577204,
1223 0.070366047488108124709267416450667,
1224 0.107159220467171935011869546685869,
1225 0.139570677926154314447804794511028,
1226 0.166269205816993933553200860481209,
1227 0.186161000015562211026800561866423,
1228 0.198431485327111576456118326443839,
1229 0.202578241925561272880620199967519
1234 0.005377479872923348987792051430128,
1235 0.015007947329316122538374763075807,
1236 0.025460847326715320186874001019653,
1237 0.035346360791375846222037948478360,
1238 0.044589751324764876608227299373280,
1239 0.053481524690928087265343147239430,
1240 0.062009567800670640285139230960803,
1241 0.069854121318728258709520077099147,
1242 0.076849680757720378894432777482659,
1243 0.083080502823133021038289247286104,
1244 0.088564443056211770647275443693774,
1245 0.093126598170825321225486872747346,
1246 0.096642726983623678505179907627589,
1247 0.099173598721791959332393173484603,
1248 0.100769845523875595044946662617570,
1249 0.101330007014791549017374792767493
1254 double *
result,
double *abserr,
1255 double *resabs,
double *resasc)
1257 double fv1[16], fv2[16];
1259 gsl_integration_qk (16,
xgkC,
wgC,
wgkC, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
1268 0.998859031588277663838315576545863,
1269 0.993128599185094924786122388471320,
1270 0.981507877450250259193342994720217,
1271 0.963971927277913791267666131197277,
1272 0.940822633831754753519982722212443,
1273 0.912234428251325905867752441203298,
1274 0.878276811252281976077442995113078,
1275 0.839116971822218823394529061701521,
1276 0.795041428837551198350638833272788,
1277 0.746331906460150792614305070355642,
1278 0.693237656334751384805490711845932,
1279 0.636053680726515025452836696226286,
1280 0.575140446819710315342946036586425,
1281 0.510867001950827098004364050955251,
1282 0.443593175238725103199992213492640,
1283 0.373706088715419560672548177024927,
1284 0.301627868114913004320555356858592,
1285 0.227785851141645078080496195368575,
1286 0.152605465240922675505220241022678,
1287 0.076526521133497333754640409398838,
1288 0.000000000000000000000000000000000
1296 0.017614007139152118311861962351853,
1297 0.040601429800386941331039952274932,
1298 0.062672048334109063569506535187042,
1299 0.083276741576704748724758143222046,
1300 0.101930119817240435036750135480350,
1301 0.118194531961518417312377377711382,
1302 0.131688638449176626898494499748163,
1303 0.142096109318382051329298325067165,
1304 0.149172986472603746787828737001969,
1305 0.152753387130725850698084331955098
1310 0.003073583718520531501218293246031,
1311 0.008600269855642942198661787950102,
1312 0.014626169256971252983787960308868,
1313 0.020388373461266523598010231432755,
1314 0.025882133604951158834505067096153,
1315 0.031287306777032798958543119323801,
1316 0.036600169758200798030557240707211,
1317 0.041668873327973686263788305936895,
1318 0.046434821867497674720231880926108,
1319 0.050944573923728691932707670050345,
1320 0.055195105348285994744832372419777,
1321 0.059111400880639572374967220648594,
1322 0.062653237554781168025870122174255,
1323 0.065834597133618422111563556969398,
1324 0.068648672928521619345623411885368,
1325 0.071054423553444068305790361723210,
1326 0.073030690332786667495189417658913,
1327 0.074582875400499188986581418362488,
1328 0.075704497684556674659542775376617,
1329 0.076377867672080736705502835038061,
1330 0.076600711917999656445049901530102
1335 double *
result,
double *abserr,
1336 double *resabs,
double *resasc)
1338 double fv1[21], fv2[21];
1340 gsl_integration_qk (21,
xgkD,
wgD,
wgkD, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
1349 0.999262104992609834193457486540341,
1350 0.995556969790498097908784946893902,
1351 0.988035794534077247637331014577406,
1352 0.976663921459517511498315386479594,
1353 0.961614986425842512418130033660167,
1354 0.942974571228974339414011169658471,
1355 0.920747115281701561746346084546331,
1356 0.894991997878275368851042006782805,
1357 0.865847065293275595448996969588340,
1358 0.833442628760834001421021108693570,
1359 0.797873797998500059410410904994307,
1360 0.759259263037357630577282865204361,
1361 0.717766406813084388186654079773298,
1362 0.673566368473468364485120633247622,
1363 0.626810099010317412788122681624518,
1364 0.577662930241222967723689841612654,
1365 0.526325284334719182599623778158010,
1366 0.473002731445714960522182115009192,
1367 0.417885382193037748851814394594572,
1368 0.361172305809387837735821730127641,
1369 0.303089538931107830167478909980339,
1370 0.243866883720988432045190362797452,
1371 0.183718939421048892015969888759528,
1372 0.122864692610710396387359818808037,
1373 0.061544483005685078886546392366797,
1374 0.000000000000000000000000000000000
1382 0.011393798501026287947902964113235,
1383 0.026354986615032137261901815295299,
1384 0.040939156701306312655623487711646,
1385 0.054904695975835191925936891540473,
1386 0.068038333812356917207187185656708,
1387 0.080140700335001018013234959669111,
1388 0.091028261982963649811497220702892,
1389 0.100535949067050644202206890392686,
1390 0.108519624474263653116093957050117,
1391 0.114858259145711648339325545869556,
1392 0.119455763535784772228178126512901,
1393 0.122242442990310041688959518945852,
1394 0.123176053726715451203902873079050
1399 0.001987383892330315926507851882843,
1400 0.005561932135356713758040236901066,
1401 0.009473973386174151607207710523655,
1402 0.013236229195571674813656405846976,
1403 0.016847817709128298231516667536336,
1404 0.020435371145882835456568292235939,
1405 0.024009945606953216220092489164881,
1406 0.027475317587851737802948455517811,
1407 0.030792300167387488891109020215229,
1408 0.034002130274329337836748795229551,
1409 0.037116271483415543560330625367620,
1410 0.040083825504032382074839284467076,
1411 0.042872845020170049476895792439495,
1412 0.045502913049921788909870584752660,
1413 0.047982537138836713906392255756915,
1414 0.050277679080715671963325259433440,
1415 0.052362885806407475864366712137873,
1416 0.054251129888545490144543370459876,
1417 0.055950811220412317308240686382747,
1418 0.057437116361567832853582693939506,
1419 0.058689680022394207961974175856788,
1420 0.059720340324174059979099291932562,
1421 0.060539455376045862945360267517565,
1422 0.061128509717053048305859030416293,
1423 0.061471189871425316661544131965264,
1424 0.061580818067832935078759824240066
1431 double *
result,
double *abserr,
1432 double *resabs,
double *resasc)
1434 double fv1[26], fv2[26];
1436 gsl_integration_qk (26,
xgkE,
wgE,
wgkE, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
1445 0.999484410050490637571325895705811,
1446 0.996893484074649540271630050918695,
1447 0.991630996870404594858628366109486,
1448 0.983668123279747209970032581605663,
1449 0.973116322501126268374693868423707,
1450 0.960021864968307512216871025581798,
1451 0.944374444748559979415831324037439,
1452 0.926200047429274325879324277080474,
1453 0.905573307699907798546522558925958,
1454 0.882560535792052681543116462530226,
1455 0.857205233546061098958658510658944,
1456 0.829565762382768397442898119732502,
1457 0.799727835821839083013668942322683,
1458 0.767777432104826194917977340974503,
1459 0.733790062453226804726171131369528,
1460 0.697850494793315796932292388026640,
1461 0.660061064126626961370053668149271,
1462 0.620526182989242861140477556431189,
1463 0.579345235826361691756024932172540,
1464 0.536624148142019899264169793311073,
1465 0.492480467861778574993693061207709,
1466 0.447033769538089176780609900322854,
1467 0.400401254830394392535476211542661,
1468 0.352704725530878113471037207089374,
1469 0.304073202273625077372677107199257,
1470 0.254636926167889846439805129817805,
1471 0.204525116682309891438957671002025,
1472 0.153869913608583546963794672743256,
1473 0.102806937966737030147096751318001,
1474 0.051471842555317695833025213166723,
1475 0.000000000000000000000000000000000
1483 0.007968192496166605615465883474674,
1484 0.018466468311090959142302131912047,
1485 0.028784707883323369349719179611292,
1486 0.038799192569627049596801936446348,
1487 0.048402672830594052902938140422808,
1488 0.057493156217619066481721689402056,
1489 0.065974229882180495128128515115962,
1490 0.073755974737705206268243850022191,
1491 0.080755895229420215354694938460530,
1492 0.086899787201082979802387530715126,
1493 0.092122522237786128717632707087619,
1494 0.096368737174644259639468626351810,
1495 0.099593420586795267062780282103569,
1496 0.101762389748405504596428952168554,
1497 0.102852652893558840341285636705415
1502 0.001389013698677007624551591226760,
1503 0.003890461127099884051267201844516,
1504 0.006630703915931292173319826369750,
1505 0.009273279659517763428441146892024,
1506 0.011823015253496341742232898853251,
1507 0.014369729507045804812451432443580,
1508 0.016920889189053272627572289420322,
1509 0.019414141193942381173408951050128,
1510 0.021828035821609192297167485738339,
1511 0.024191162078080601365686370725232,
1512 0.026509954882333101610601709335075,
1513 0.028754048765041292843978785354334,
1514 0.030907257562387762472884252943092,
1515 0.032981447057483726031814191016854,
1516 0.034979338028060024137499670731468,
1517 0.036882364651821229223911065617136,
1518 0.038678945624727592950348651532281,
1519 0.040374538951535959111995279752468,
1520 0.041969810215164246147147541285970,
1521 0.043452539701356069316831728117073,
1522 0.044814800133162663192355551616723,
1523 0.046059238271006988116271735559374,
1524 0.047185546569299153945261478181099,
1525 0.048185861757087129140779492298305,
1526 0.049055434555029778887528165367238,
1527 0.049795683427074206357811569379942,
1528 0.050405921402782346840893085653585,
1529 0.050881795898749606492297473049805,
1530 0.051221547849258772170656282604944,
1531 0.051426128537459025933862879215781,
1532 0.051494729429451567558340433647099
1537 double *
result,
double *abserr,
1538 double *resabs,
double *resasc)
1540 double fv1[31], fv2[31];
1542 gsl_integration_qk (31,
xgkF,
wgF,
wgkF, fv1, fv2, f, a, b, result, abserr, resabs, resasc);
1545 gsl_integration_workspace*
1548 gsl_integration_workspace* w ;
1552 GSL_ERROR_VAL (
"workspace length n must be positive integer",
1556 w = (gsl_integration_workspace *)
1557 malloc (
sizeof (gsl_integration_workspace));
1561 GSL_ERROR_VAL (
"failed to allocate space for workspace struct",
1565 w->alist = (
double *)
malloc (n *
sizeof (
double));
1575 w->blist = (
double *)
malloc (n *
sizeof (
double));
1586 w->rlist = (
double *)
malloc (n *
sizeof (
double));
1599 w->elist = (
double *)
malloc (n *
sizeof (
double));
1612 w->order = (
size_t *)
malloc (n *
sizeof (
size_t));
1626 w->level = (
size_t *)
malloc (n *
sizeof (
size_t));
1643 w->maximum_level = 0 ;
1664 reset_nrmax (gsl_integration_workspace * workspace);
1669 workspace->nrmax = 0;
1670 workspace->i = workspace->order[0] ;
1686 int id = workspace->nrmax;
1689 const size_t * level = workspace->level;
1690 const size_t * order = workspace->order;
1692 size_t limit = workspace->limit ;
1693 size_t last = workspace->size - 1 ;
1695 if (last > (1 + limit / 2))
1697 jupbnd = limit + 1 - last;
1704 for (k =
id; k <= jupbnd; k++)
1706 size_t i_max = order[workspace->nrmax];
1708 workspace->i = i_max ;
1710 if (level[i_max] < workspace->maximum_level)
1724 size_t i = workspace->i ;
1725 const size_t * level = workspace->level;
1727 if (level[i] < workspace->maximum_level)
1739 struct extrapolation_table
1764 table->rlist2[0] =
y;
1773 table->rlist2[
n] =
y;
1783 qelg (
struct extrapolation_table *table,
double *
result,
double *abserr);
1786 qelg (
struct extrapolation_table *table,
double *
result,
double *abserr)
1788 double *epstab = table->rlist2;
1789 double *res3la = table->res3la;
1790 const size_t n = table->n - 1;
1792 const double current = epstab[
n];
1797 const size_t newelm = n / 2;
1798 const size_t n_orig =
n;
1802 const size_t nres_orig = table->nres;
1814 epstab[n + 2] = epstab[
n];
1817 for (i = 0; i < newelm; i++)
1819 double res = epstab[n - 2 * i + 2];
1820 double e0 = epstab[n - 2 * i - 2];
1821 double e1 = epstab[n - 2 * i - 1];
1824 double e1abs =
fabs (e1);
1825 double delta2 = e2 - e1;
1826 double err2 =
fabs (delta2);
1828 double delta3 = e1 - e0;
1829 double err3 =
fabs (delta3);
1832 double e3, delta1, err1, tol1, ss;
1834 if (err2 <= tol2 && err3 <= tol3)
1840 absolute = err2 + err3;
1846 e3 = epstab[n - 2 * i];
1847 epstab[n - 2 * i] = e1;
1849 err1 =
fabs (delta1);
1855 if (err1 <= tol1 || err2 <= tol2 || err3 <= tol3)
1861 ss = (1 / delta1 + 1 / delta2) - 1 / delta3;
1867 if (fabs (ss * e1) <= 0.0001)
1877 epstab[n - 2 * i] = res;
1880 const double error = err2 +
fabs (res - e2) + err3;
1882 if (error <= *abserr)
1893 const size_t limexp = 50 - 1;
1895 if (n_final == limexp)
1897 n_final = 2 * (limexp / 2);
1901 if (n_orig % 2 == 1)
1903 for (i = 0; i <= newelm; i++)
1905 epstab[1 + i * 2] = epstab[i * 2 + 3];
1910 for (i = 0; i <= newelm; i++)
1912 epstab[i * 2] = epstab[i * 2 + 2];
1916 if (n_orig != n_final)
1918 for (i = 0; i <= n_final; i++)
1920 epstab[i] = epstab[n_orig - n_final + i];
1924 table->n = n_final + 1;
1928 res3la[nres_orig] = *
result;
1933 *abserr = (
fabs (*result - res3la[2]) +
fabs (*result - res3la[1])
1934 +
fabs (*result - res3la[0]));
1936 res3la[0] = res3la[1];
1937 res3la[1] = res3la[2];
1947 table->nres = nres_orig + 1;
1971 b,
const double epsabs,
const double epsrel,
const size_t limit,
1972 gsl_integration_workspace * workspace,
double *
result,
double *abserr,
1978 double epsabs,
double epsrel,
size_t limit,
1979 gsl_integration_workspace * workspace,
1980 double *
result,
double * abserr)
1982 int status =
qags (f, a, b, epsabs, epsrel, limit,
1996 static double i_transform (
double t,
void *params);
2000 double epsabs,
double epsrel,
size_t limit,
2001 gsl_integration_workspace * workspace,
2002 double *
result,
double *abserr)
2009 f_transform.params =
f;
2011 status =
qags (&f_transform, 0.0, 1.0,
2012 epsabs, epsrel, limit,
2024 double x = (1 - t) / t;
2044 double epsabs,
double epsrel,
size_t limit,
2045 gsl_integration_workspace * workspace,
2046 double *
result,
double *abserr)
2051 struct il_params transform_params ;
2053 transform_params.b =
b ;
2054 transform_params.f =
f ;
2057 f_transform.params = &transform_params;
2059 status =
qags (&f_transform, 0.0, 1.0,
2060 epsabs, epsrel, limit,
2071 struct il_params *p = (
struct il_params *) params;
2074 double x = b - (1 - t) / t;
2093 double epsabs,
double epsrel,
size_t limit,
2094 gsl_integration_workspace * workspace,
2095 double *
result,
double *abserr)
2100 struct iu_params transform_params ;
2102 transform_params.a =
a ;
2103 transform_params.f =
f ;
2106 f_transform.params = &transform_params;
2108 status =
qags (&f_transform, 0.0, 1.0,
2109 epsabs, epsrel, limit,
2120 struct iu_params *p = (
struct iu_params *) params;
2123 double x = a + (1 - t) / t;
2132 const double a,
const double b,
2133 const double epsabs,
const double epsrel,
2135 gsl_integration_workspace * workspace,
2136 double *
result,
double *abserr,
2139 double area, errsum;
2140 double res_ext, err_ext;
2141 double result0, abserr0, resabs0, resasc0;
2145 double error_over_large_intervals = 0;
2146 double reseps = 0, abseps = 0, correc = 0;
2148 int roundoff_type1 = 0, roundoff_type2 = 0, roundoff_type3 = 0;
2149 int error_type = 0, error_type2 = 0;
2151 size_t iteration = 0;
2153 int positive_integrand = 0;
2154 int extrapolate = 0;
2155 int disallow_extrapolation = 0;
2157 struct extrapolation_table table;
2166 if (limit > workspace->limit)
2175 GSL_ERROR (
"tolerance cannot be acheived with given epsabs and epsrel",
2181 q (f, a, b, &result0, &abserr0, &resabs0, &resasc0);
2187 if (abserr0 <= 100 * GSL_DBL_EPSILON * resabs0 && abserr0 > tolerance)
2192 GSL_ERROR (
"cannot reach tolerance because of roundoff error" 2195 else if ((abserr0 <= tolerance && abserr0 != resasc0) || abserr0 == 0.0)
2202 else if (limit == 1)
2227 size_t current_level;
2228 double a1, b1, a2, b2;
2229 double a_i, b_i, r_i, e_i;
2230 double area1 = 0, area2 = 0, area12 = 0;
2231 double error1 = 0, error2 = 0, error12 = 0;
2232 double resasc1, resasc2;
2233 double resabs1, resabs2;
2238 retrieve (workspace, &a_i, &b_i, &r_i, &e_i);
2240 current_level = workspace->level[workspace->i] + 1;
2243 b1 = 0.5 * (a_i + b_i);
2249 q (f, a1, b1, &area1, &error1, &resabs1, &resasc1);
2250 q (f, a2, b2, &area2, &error2, &resabs2, &resasc2);
2252 area12 = area1 + area2;
2253 error12 = error1 + error2;
2263 errsum = errsum + error12 - e_i;
2264 area = area + area12 - r_i;
2268 if (resasc1 != error1 && resasc2 != error2)
2270 double delta = r_i - area12;
2272 if (
fabs (delta) <= 1.0
e-5 *
fabs (area12) && error12 >= 0.99 * e_i)
2283 if (iteration > 10 && error12 > e_i)
2291 if (roundoff_type1 + roundoff_type2 >= 10 || roundoff_type3 >= 20)
2296 if (roundoff_type2 >= 5)
2311 update (workspace, a1, b1, area1, error1, a2, b2, area2, error2);
2313 if (errsum <= tolerance)
2315 goto compute_result;
2323 if (iteration >= limit - 1)
2331 error_over_large_intervals = errsum;
2337 if (disallow_extrapolation)
2342 error_over_large_intervals += -last_e_i;
2344 if (current_level < workspace->maximum_level)
2346 error_over_large_intervals += error12;
2358 workspace->nrmax = 1;
2361 if (!error_type2 && error_over_large_intervals > ertest)
2371 qelg (&table, &reseps, &abseps);
2375 if (ktmin > 5 && err_ext < 0.001 * errsum)
2380 if (abseps < err_ext)
2385 correc = error_over_large_intervals;
2387 if (err_ext <= ertest)
2395 disallow_extrapolation = 1;
2398 if (error_type == 5)
2407 error_over_large_intervals = errsum;
2410 while (iteration < limit);
2416 goto compute_result;
2418 if (error_type || error_type2)
2428 if (res_ext != 0.0 && area != 0.0)
2430 if (err_ext /
fabs (res_ext) > errsum /
fabs (area))
2431 goto compute_result;
2433 else if (err_ext > errsum)
2435 goto compute_result;
2437 else if (area == 0.0)
2448 if (!positive_integrand && max_area < 0.01 * resabs0)
2453 double ratio = res_ext / area;
2455 if (ratio < 0.01 || ratio > 100.0 || errsum >
fabs (area))
2473 if (error_type == 0)
2477 else if (error_type == 1)
2481 else if (error_type == 2)
2483 GSL_ERROR (
"cannot reach tolerance because of roundoff error",
2486 else if (error_type == 3)
2488 GSL_ERROR (
"bad integrand behavior found in the integration interval",
2491 else if (error_type == 4)
2493 GSL_ERROR (
"roundoff error detected in the extrapolation table",
2496 else if (error_type == 5)
2498 GSL_ERROR (
"integral is divergent, or slowly convergent",
static double sum_results(const gsl_integration_workspace *workspace)
static const double xgkC[16]
RooAbsIntegrator is the abstract interface for integrators of real-valued functions that implement th...
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.
RooNumIntConfig holds the configuration parameters of the various numeric integrators used by RooReal...
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.
Double_t _epsAbs
Current coordinate.
RooNumIntFactory is a factory to instantiate numeric integrators from a given function binding and a ...
static const double wgkB[11]
static const double wgkD[21]
const RooArgSet & getConfigSection(const char *name) const
Retrieve configuration information specific to integrator with given name.
static Int_t isInfinite(Double_t x)
Return true if x is infinite by RooNumBer internal specification.
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)
virtual Bool_t checkLimits() const
Check that our integration range is finite and otherwise return kFALSE.
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)
static const double wgkC[16]
static const double wgkA[8]
Double_t integrand(const Double_t x[]) const
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_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 GSL_COERCE_DBL(x)
RooRealVar represents a fundamental (non-derived) real valued object.
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
virtual ~RooAdaptiveGaussKronrodIntegrator1D()
Destructor.
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)
UInt_t getDimension() const
void gsl_integration_qcheb(gsl_function *f, double a, double b, double *cheb12, double *cheb24)
static void update(gsl_integration_workspace *workspace, double a1, double b1, double area1, double error1, double a2, double b2, double area2, double error2)
Double_t * xvec(Double_t &xx)
static void qpsrt(gsl_integration_workspace *workspace)
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)
static void initialise_table(struct extrapolation_table *table)
RooCategory represents a fundamental (non-derived) discrete value object.
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)
const RooAbsFunc * integrand() const
Bool_t _useIntegrandLimits
static const double wgB[5]
struct gsl_function_struct gsl_function
virtual Double_t getMaxLimit(UInt_t dimension) const =0
static double iu_transform(double t, void *params)
virtual RooAbsIntegrator * clone(const RooAbsFunc &function, const RooNumIntConfig &config) const
Virtual constructor.
static int test_positivity(double result, double resabs)
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)
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
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)
static const double xgkE[26]
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...
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.
void gsl_integration_qk31(const gsl_function *f, double a, double b, double *result, double *abserr, double *resabs, double *resasc)
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
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]
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)
RooAdaptiveGaussKronrodIntegrator1D implements the Gauss-Kronrod integration algorithm.
virtual const char * GetName() const
Returns name of object.
static const double wgD[11]
Abstract interface for evaluating a real-valued function of one real variable and performing numerica...
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 _xmax
Lower integration bound.
RooAdaptiveGaussKronrodIntegrator1D()
coverity[UNINIT_CTOR] Default constructor
Bool_t initialize()
Initialize integrator allocate buffers and setup GSL workspace.