Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooRombergIntegrator.cxx
Go to the documentation of this file.
1/// \cond ROOFIT_INTERNAL
2
3/*****************************************************************************
4 * Project: RooFit *
5 * Package: RooFitCore *
6 * @(#)root/roofitcore:$Id$
7 * Authors: *
8 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
9 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
10 * *
11 * Copyright (c) 2000-2005, Regents of the University of California *
12 * and Stanford University. All rights reserved. *
13 * *
14 * Redistribution and use in source and binary forms, *
15 * with or without modification, are permitted according to the terms *
16 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
17 *****************************************************************************/
18
19/**
20\file RooRombergIntegrator.cxx
21\class RooRombergIntegrator
22\ingroup Roofitcore
23
24Adaptive numerical integration algorithm.
25
26It uses Romberg's method with trapezoids or midpoints.
27The integrand is approximated by \f$ 1, 2, 4, 8, \ldots, 2^n \f$ trapezoids, and
28Richardson series acceleration is applied to the series. For refinement step \f$ n \f$, that is
29\f[
30 R(n,m) = \frac{1}{4^m - 1} \left( 4^m R(n, m-1) - R(n-1, m-1) \right)
31\f]
32The integrator will evaluate the first six refinements (i.e. 64 points) in one go,
33apply a five-point series acceleration, and successively add more steps until the
34desired precision is reached.
35**/
36
37#include "Riostream.h"
38
39#include "TClass.h"
41#include "RooArgSet.h"
42#include "RooRealVar.h"
43#include "RooNumber.h"
44#include "RooNumIntConfig.h"
45#include "RooNumIntFactory.h"
46#include "RooMsgService.h"
47
48#include <cassert>
49
50namespace {
51
52constexpr int nPoints = 5;
53
54////////////////////////////////////////////////////////////////////////////////
55/// Calculate the n-th stage of refinement of the Second Euler-Maclaurin
56/// summation rule which has the useful property of not evaluating the
57/// integrand at either of its endpoints but requires more function
58/// evaluations than the trapezoidal rule. This rule can be used with
59/// a suitable change of variables to estimate improper integrals.
60
61double addMidpoints(std::function<double(double)> const &func, double savedResult, int n, double xmin, double xmax)
62{
63 const double range = xmax - xmin;
64
65 if (n == 1) {
66 double xmid = 0.5 * (xmin + xmax);
67 return range * func(xmid);
68 }
69
70 int it = 1;
71 for (int j = 1; j < n - 1; j++) {
72 it *= 3;
73 }
74 double tnm = it;
75 double del = range / (3. * tnm);
76 double ddel = del + del;
77 double x = xmin + 0.5 * del;
78 double sum = 0;
79 for (int j = 1; j <= it; j++) {
80 sum += func(x);
81 x += ddel;
82 sum += func(x);
83 x += del;
84 }
85 return (savedResult + range * sum / tnm) / 3.;
86}
87
88////////////////////////////////////////////////////////////////////////////////
89/// Calculate the n-th stage of refinement of the extended trapezoidal
90/// summation rule. This is the most efficient rule for a well behaved
91/// integrands that can be evaluated over its entire range, including the
92/// endpoints.
93
94double addTrapezoids(std::function<double(double)> const &func, double savedResult, int n, double xmin, double xmax)
95{
96 const double range = xmax - xmin;
97
98 if (n == 1) {
99 // use a single trapezoid to cover the full range
100 return 0.5 * range * (func(xmin) + func(xmax));
101 }
102
103 // break the range down into several trapezoids using 2**(n-2)
104 // equally-spaced interior points
105 const int nInt = 1 << (n - 2);
106 const double del = range / nInt;
107
108 double sum = 0.;
109 // TODO Replace by batch computation
110 for (int j = 0; j < nInt; ++j) {
111 double x = xmin + (0.5 + j) * del;
112 sum += func(x);
113 }
114
115 return 0.5 * (savedResult + range * sum / nInt);
116}
117
118////////////////////////////////////////////////////////////////////////////////
119/// Extrapolate result to final value.
120
121std::pair<double, double> extrapolate(int n, double const *h, double const *s, double *c, double *d)
122{
123 double const *xa = &h[n - nPoints];
124 double const *ya = &s[n - nPoints];
125 int ns = 1;
126
127 double dif = std::abs(xa[1]);
128 for (int i = 1; i <= nPoints; i++) {
129 double dift = std::abs(xa[i]);
130 if (dift < dif) {
131 ns = i;
132 dif = dift;
133 }
134 c[i] = ya[i];
135 d[i] = ya[i];
136 }
137 double extrapError = 0.0;
138 double extrapValue = ya[ns--];
139 for (int m = 1; m < nPoints; m++) {
140 for (int i = 1; i <= nPoints - m; i++) {
141 double ho = xa[i];
142 double hp = xa[i + m];
143 double w = c[i + 1] - d[i];
144 double den = ho - hp;
145 if (den == 0.0) {
146 throw std::runtime_error("RooRombergIntegrator::extrapolate: internal error");
147 }
148 den = w / den;
149 d[i] = hp * den;
150 c[i] = ho * den;
151 }
152 extrapError = 2 * ns < (nPoints - m) ? c[ns + 1] : d[ns--];
154 }
155 return {extrapValue, extrapError};
156}
157
158} // namespace
159
160namespace RooFit {
161namespace Detail {
162
163std::pair<double, int> integrate1d(std::function<double(double)> func, bool doTrapezoid, int maxSteps, int minStepsZero,
164 int fixSteps, double epsAbs, double epsRel, bool doExtrap, double xmin, double xmax,
165 std::span<double> hArr, std::span<double> sArr)
166{
167 assert(int(hArr.size()) == maxSteps + 2);
168 assert(int(sArr.size()) == maxSteps + 2);
169
170 const double range = xmax - xmin;
171
172 // Small working arrays can be on the stack
173 std::array<double, nPoints + 1> cArr = {};
174 std::array<double, nPoints + 1> dArr = {};
175
176 hArr[1] = 1.0;
177 double zeroThresh = epsAbs / range;
178 for (int j = 1; j <= maxSteps; ++j) {
179 // refine our estimate using the appropriate summation rule
180 sArr[j] =
181 doTrapezoid ? addTrapezoids(func, sArr[j - 1], j, xmin, xmax) : addMidpoints(func, sArr[j - 1], j, xmin, xmax);
182
183 if (j >= minStepsZero) {
184 bool allZero(true);
185 for (int jj = 0; jj <= j; jj++) {
186 if (sArr[j] >= zeroThresh) {
187 allZero = false;
188 }
189 }
190 if (allZero) {
191 // std::cout << "Roo1DIntegrator(" << this << "): zero convergence at step " << j << ", value = " << 0 <<
192 // std::endl ;
193 return {0, j};
194 }
195 }
196
197 if (fixSteps > 0) {
198
199 // Fixed step mode, return result after fixed number of steps
200 if (j == fixSteps) {
201 // std::cout << "returning result at fixed step " << j << std::endl ;
202 return {sArr[j], j};
203 }
204
205 } else if (j >= nPoints) {
206
207 double extrapValue = 0.0;
208 double extrapError = 0.0;
209
210 // extrapolate the results of recent refinements and check for a stable result
211 if (doExtrap) {
212 std::tie(extrapValue, extrapError) = extrapolate(j, hArr.data(), sArr.data(), cArr.data(), dArr.data());
213 } else {
214 extrapValue = sArr[j];
215 extrapError = sArr[j] - sArr[j - 1];
216 }
217
218 if (std::abs(extrapError) <= epsRel * std::abs(extrapValue)) {
219 return {extrapValue, j};
220 }
221 if (std::abs(extrapError) <= epsAbs) {
222 return {extrapValue, j};
223 }
224 }
225 // update the step size for the next refinement of the summation
226 hArr[j + 1] = doTrapezoid ? hArr[j] / 4. : hArr[j] / 9.;
227 }
228
229 return {sArr[maxSteps], maxSteps};
230}
231
232} // namespace Detail
233} // namespace RooFit
234
235
236// Register this class with RooNumIntConfig
237
238////////////////////////////////////////////////////////////////////////////////
239/// Register integrator plugins, their parameters and capabilities with RooNumIntFactory.
240
241void RooRombergIntegrator::registerIntegrator(RooNumIntFactory &fact)
242{
243 RooCategory sumRule("sumRule", "Summation Rule");
244 sumRule.defineType("Trapezoid", RooRombergIntegrator::Trapezoid);
245 sumRule.defineType("Midpoint", RooRombergIntegrator::Midpoint);
246 sumRule.setLabel("Trapezoid");
247 RooCategory extrap("extrapolation", "Extrapolation procedure");
248 extrap.defineType("None", 0);
249 extrap.defineType("Wynn-Epsilon", 1);
250 extrap.setLabel("Wynn-Epsilon");
251 RooRealVar maxSteps("maxSteps", "Maximum number of steps", 20);
252 RooRealVar minSteps("minSteps", "Minimum number of steps", 999);
253 RooRealVar fixSteps("fixSteps", "Fixed number of steps", 0);
254 RooRealVar numSeg("numSeg", "Number of segments", 3); // only for the segmented integrators
255
256 std::string name = "RooIntegrator1D";
257
258 auto creator = [](const RooAbsFunc &function, const RooNumIntConfig &config) {
259 return std::make_unique<RooRombergIntegrator>(function, config, 1, /*doSegmentation=*/false);
260 };
261
262 fact.registerPlugin(name, creator, {sumRule, extrap, maxSteps, minSteps, fixSteps},
263 /*canIntegrate1D=*/true,
264 /*canIntegrate2D=*/false,
265 /*canIntegrateND=*/false,
266 /*canIntegrateOpenEnded=*/false);
267
268 RooNumIntConfig::defaultConfig().method1D().setLabel(name);
269
270 auto creator2d = [](const RooAbsFunc &function, const RooNumIntConfig &config) {
271 return std::make_unique<RooRombergIntegrator>(function, config, 2, /*doSegmentation=*/false);
272 };
273 std::string name2d = "RooIntegrator2D";
274 fact.registerPlugin(name2d, creator2d, {},
275 /*canIntegrate1D=*/false,
276 /*canIntegrate2D=*/true,
277 /*canIntegrateND=*/false,
278 /*canIntegrateOpenEnded=*/false,
279 /*depName=*/"RooIntegrator1D");
280 RooNumIntConfig::defaultConfig().method2D().setLabel(name2d);
281
282 auto creatorSeg = [](const RooAbsFunc &function, const RooNumIntConfig &config) {
283 return std::make_unique<RooRombergIntegrator>(function, config, 1, /*doSegmentation=*/true);
284 };
285
286 fact.registerPlugin("RooSegmentedIntegrator1D", creatorSeg, {numSeg},
287 /*canIntegrate1D=*/true,
288 /*canIntegrate2D=*/false,
289 /*canIntegrateND=*/false,
290 /*canIntegrateOpenEnded=*/false,
291 /*depName=*/"RooIntegrator1D");
292
293 auto creatorSeg2d = [](const RooAbsFunc &function, const RooNumIntConfig &config) {
294 return std::make_unique<RooRombergIntegrator>(function, config, 2, /*doSegmentation=*/true);
295 };
296
297 fact.registerPlugin("RooSegmentedIntegrator2D", creatorSeg2d, {},
298 /*canIntegrate1D=*/false,
299 /*canIntegrate2D=*/true,
300 /*canIntegrateND=*/false,
301 /*canIntegrateOpenEnded=*/false,
302 /*depName=*/"RooSegmentedIntegrator1D");
303}
304
305////////////////////////////////////////////////////////////////////////////////
306/// Construct integrator on given function binding, using specified summation
307/// rule, maximum number of steps and conversion tolerance. The integration
308/// limits are taken from the function binding.
309
310RooRombergIntegrator::RooRombergIntegrator(const RooAbsFunc &function, SummationRule rule, int maxSteps, double eps)
311 : RooAbsIntegrator(function), _useIntegrandLimits(true), _rule(rule), _maxSteps(maxSteps), _epsAbs(eps), _epsRel(eps)
312{
313 _valid = initialize();
314}
315
316////////////////////////////////////////////////////////////////////////////////
317/// Construct integrator on given function binding for given range,
318/// using specified summation rule, maximum number of steps and
319/// conversion tolerance. The integration limits are taken from the
320/// function binding.
321
322RooRombergIntegrator::RooRombergIntegrator(const RooAbsFunc &function, double xmin, double xmax, SummationRule rule,
323 int maxSteps, double eps)
325 _useIntegrandLimits(false),
326 _rule(rule),
328 _epsAbs(eps),
329 _epsRel(eps)
330{
331 _xmin.push_back(xmin);
332 _xmax.push_back(xmax);
333 _valid = initialize();
334}
335
336////////////////////////////////////////////////////////////////////////////////
337/// Construct integrator on given function binding, using specified
338/// configuration object. The integration limits are taken from the
339/// function binding
340
341RooRombergIntegrator::RooRombergIntegrator(const RooAbsFunc &function, const RooNumIntConfig &config, int nDim,
342 bool doSegmentation)
343 : RooAbsIntegrator(function, config.printEvalCounter()),
344 _nDim{nDim},
345 _epsAbs(config.epsAbs()),
346 _epsRel(config.epsRel())
347{
348 // Extract parameters from config object
349 const RooArgSet &configSet = config.getConfigSection("RooIntegrator1D");
350 _rule = (SummationRule)configSet.getCatIndex("sumRule", Trapezoid);
351 _maxSteps = (int)configSet.getRealValue("maxSteps", 20);
352 _minStepsZero = (int)configSet.getRealValue("minSteps", 999);
353 _fixSteps = (int)configSet.getRealValue("fixSteps", 0);
354 _doExtrap = (bool)configSet.getCatIndex("extrapolation", 1);
355 if (doSegmentation) {
356 _nSeg = (int)config.getConfigSection("RooSegmentedIntegrator1D").getRealValue("numSeg", 3);
357 _epsAbs /= std::sqrt(_nSeg);
358 _epsRel /= std::sqrt(_nSeg);
359 }
360
361 if (_fixSteps > _maxSteps) {
362 oocoutE(nullptr, Integration) << "RooRombergIntegrator::ctor() ERROR: fixSteps>maxSteps, fixSteps set to maxSteps"
363 << std::endl;
365 }
366
367 _useIntegrandLimits = true;
368 _valid = initialize();
369}
370
371////////////////////////////////////////////////////////////////////////////////
372/// Construct integrator on given function binding, using specified
373/// configuration object and integration range
374
375RooRombergIntegrator::RooRombergIntegrator(const RooAbsFunc &function, double xmin, double xmax,
376 const RooNumIntConfig &config, int nDim)
377 : RooAbsIntegrator(function, config.printEvalCounter()),
378 _useIntegrandLimits(false),
379 _nDim{nDim},
380 _epsAbs(config.epsAbs()),
381 _epsRel(config.epsRel())
382{
383 // Extract parameters from config object
384 const RooArgSet &configSet = config.getConfigSection("RooIntegrator1D");
385 _rule = (SummationRule)configSet.getCatIndex("sumRule", Trapezoid);
386 _maxSteps = (int)configSet.getRealValue("maxSteps", 20);
387 _minStepsZero = (int)configSet.getRealValue("minSteps", 999);
388 _fixSteps = (int)configSet.getRealValue("fixSteps", 0);
389 _doExtrap = (bool)configSet.getCatIndex("extrapolation", 1);
390
391 _xmin.push_back(xmin);
392 _xmax.push_back(xmax);
393 _valid = initialize();
394}
395
396////////////////////////////////////////////////////////////////////////////////
397/// Initialize the integrator
398
399bool RooRombergIntegrator::initialize()
400{
401 // apply defaults if necessary
402 if (_maxSteps <= 0) {
403 _maxSteps = (_rule == Trapezoid) ? 20 : 14;
404 }
405
406 if (_epsRel <= 0)
407 _epsRel = 1e-6;
408 if (_epsAbs <= 0)
409 _epsAbs = 1e-6;
410
411 // check that the integrand is a valid function
412 if (!isValid()) {
413 oocoutE(nullptr, Integration) << "RooRombergIntegrator::initialize: cannot integrate invalid function"
414 << std::endl;
415 return false;
416 }
417
418 // Allocate coordinate buffer size after number of function dimensions
419 _x.resize(_function->getDimension());
420
421 // Allocate workspace for numerical integration engine
422 _wksp.resize(_nDim * 2 * _maxSteps + 4);
423
424 return checkLimits();
425}
426
427////////////////////////////////////////////////////////////////////////////////
428/// Change our integration limits. Return true if the new limits are
429/// ok, or otherwise false. Always returns false and does nothing
430/// if this object was constructed to always use our integrand's limits.
431
432bool RooRombergIntegrator::setLimits(double *xmin, double *xmax)
433{
434 if (_useIntegrandLimits) {
435 oocoutE(nullptr, Integration) << "RooRombergIntegrator::setLimits: cannot override integrand's limits"
436 << std::endl;
437 return false;
438 }
439 _xmin.resize(_nDim);
440 _xmax.resize(_nDim);
441 for (int iDim = 0; iDim < _nDim; ++iDim) {
442 _xmin[iDim] = xmin[iDim];
443 _xmax[iDim] = xmax[iDim];
444 }
445 return checkLimits();
446}
447
448////////////////////////////////////////////////////////////////////////////////
449/// Check that our integration range is finite and otherwise return false.
450/// Update the limits from the integrand if requested.
451
452bool RooRombergIntegrator::checkLimits() const
453{
454 if (_useIntegrandLimits) {
455 assert(nullptr != integrand() && integrand()->isValid());
456 const_cast<std::vector<double> &>(_xmin).resize(_nDim);
457 const_cast<std::vector<double> &>(_xmax).resize(_nDim);
458 for (int iDim = 0; iDim < _nDim; ++iDim) {
459 const_cast<double &>(_xmin[iDim]) = integrand()->getMinLimit(iDim);
460 const_cast<double &>(_xmax[iDim]) = integrand()->getMaxLimit(iDim);
461 }
462 }
463 for (int iDim = 0; iDim < _nDim; ++iDim) {
464 const double xmin = _xmin[iDim];
465 const double xmax = _xmax[iDim];
466 const double range = xmax - xmin;
467 if (range < 0.) {
468 oocoutE(nullptr, Integration) << "RooRombergIntegrator::checkLimits: bad range with min > max (_xmin[" << iDim
469 << "] = " << xmin << " _xmax[" << iDim << "] = " << xmax << ")" << std::endl;
470 return false;
471 }
473 return false;
474 }
475 }
476 return true;
477}
478
479double RooRombergIntegrator::integral(const double *yvec)
480{
481 // Copy yvec to xvec if provided
482 if (yvec) {
483 for (unsigned int i = 0; i < _function->getDimension() - 1; i++) {
484 _x[i + _nDim] = yvec[i];
485 }
486 }
487
488 return integral(_nDim - 1, _nSeg, _wksp);
489}
490
491////////////////////////////////////////////////////////////////////////////////
492/// Calculate numeric integral at given set of function binding parameters.
493
494double RooRombergIntegrator::integral(int iDim, int nSeg, std::span<double> wksp)
495{
496 assert(isValid());
497
498 const double xmin = _xmin[iDim];
499 const double xmax = _xmax[iDim];
500 const double range = xmax - xmin;
501
502 if (range == 0.)
503 return 0.;
504
505 // In case of segmentation, split this integral up in a loop.
506 if (nSeg > 1) {
507 const double segSize = (xmax - xmin) / nSeg;
508 double result = 0.0;
509 for (int iSeg = 0; iSeg < nSeg; iSeg++) {
510 _xmin[iDim] = xmin + iSeg * segSize;
511 _xmax[iDim] = xmin + (iSeg + 1) * segSize;
512 double part = integral(iDim, 1, wksp);
513 result += part;
514 }
515
516 // reset limits
517 _xmin[iDim] = xmin;
518 _xmax[iDim] = xmax;
519
520 return result;
521 }
522
523 // From the working array
524 std::size_t nWorkingArr = _maxSteps + 2;
525 double *hArr = wksp.data();
526 double *sArr = wksp.data() + nWorkingArr;
527
528 double output = 0.0;
529 int steps = 0;
530
531 std::span<double> nextWksp{wksp.data() + 2 * _maxSteps + 4, wksp.data() + wksp.size()};
532
533 auto func = [&](double x) {
534 _x[iDim] = x;
535
536 return iDim == 0 ? integrand(_x.data()) : integral(iDim - 1, _nSeg, nextWksp);
537 };
538
539 std::tie(output, steps) =
540 RooFit::Detail::integrate1d(func, _rule == Trapezoid, _maxSteps, _minStepsZero, _fixSteps, _epsAbs, _epsRel,
542
543 if (steps == _maxSteps) {
544
545 oocoutW(nullptr, Integration) << "RooRombergIntegrator::integral: integral of " << _function->getName()
546 << " over range (" << xmin << "," << xmax << ") did not converge after "
547 << _maxSteps << " steps" << std::endl;
548 for (int j = 1; j <= _maxSteps; ++j) {
549 ooccoutW(nullptr, Integration) << " [" << j << "] h = " << hArr[j] << " , s = " << sArr[j] << std::endl;
550 }
551 }
552
553 return output;
554}
555
556/// \endcond
#define d(i)
Definition RSha256.hxx:102
#define c(i)
Definition RSha256.hxx:101
#define h(i)
Definition RSha256.hxx:106
#define e(i)
Definition RSha256.hxx:103
#define oocoutW(o, a)
#define oocoutE(o, a)
#define ooccoutW(o, a)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t del
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
char name[80]
Definition TGX11.cxx:110
float xmin
float xmax
Abstract interface for evaluating a real-valued function of one real variable and performing numerica...
Definition RooAbsFunc.h:27
Abstract interface for integrators of real-valued functions that implement the RooAbsFunc interface.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Object to represent discrete states.
Definition RooCategory.h:28
Holds the configuration parameters of the various numeric integrators used by RooRealIntegral.
const RooArgSet & getConfigSection(const char *name) const
Retrieve configuration information specific to integrator with given name.
static RooNumIntConfig & defaultConfig()
Return reference to instance of default numeric integrator configuration object.
Factory to instantiate numeric integrators from a given function binding and a given configuration.
static constexpr int isInfinite(double x)
Return true if x is infinite by RooNumber internal specification.
Definition RooNumber.h:27
Variable that can be changed from the outside.
Definition RooRealVar.h:37
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
void function(const Char_t *name_, T fun, const Char_t *docstring=0)
Definition RExports.h:167
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:64
void initialize(typename Architecture_t::Matrix_t &A, EInitialization m)
Definition Functions.h:282
TMarker m
Definition textangle.C:8
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345
static void output()