Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooIntegralMorph.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2 * Project: RooFit *
3 * *
4 * Copyright (c) 2000-2005, Regents of the University of California *
5 * and Stanford University. All rights reserved. *
6 * *
7 * Redistribution and use in source and binary forms, *
8 * with or without modification, are permitted according to the terms *
9 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
10 *****************************************************************************/
11
12/** \class RooIntegralMorph
13 \ingroup Roofit
14
15Class RooIntegralMorph is an implementation of the histogram interpolation
16technique described by Alex Read in 'NIM A 425 (1999) 357-369 'Linear interpolation of histograms'
17for continuous functions rather than histograms. The interpolation method, in short,
18works as follows.
19
20 - Given a p.d.f f1(x) with c.d.f F1(x) and p.d.f f2(x) with c.d.f F2(x)
21
22 - One finds takes a value 'y' of both c.d.fs and determines the corresponding x
23 values x(1,2) at which F(1,2)(x)==y.
24
25 - The value of the interpolated p.d.f fbar(x) is then calculated as
26 fbar(alpha*x1+(1-alpha)*x2) = f1(x1)*f2(x2) / ( alpha*f2(x2) + (1-alpha)*f1(x1) ) ;
27
28From a technical point of view class RooIntegralMorph is a p.d.f that takes
29two input p.d.fs f1(x,p) an f2(x,q) and an interpolation parameter to
30make a p.d.f fbar(x,p,q,alpha). The shapes f1 and f2 are always taken
31to be end the end-points of the parameter alpha, regardless of what
32the those numeric values are.
33
34Since the value of fbar(x) cannot be easily calculated for a given value
35of x, class RooIntegralMorph is an implementation of RooAbsCachedPdf and
36calculates the shape of the interpolated p.d.f. fbar(x) for all values
37of x for a given value of alpha,p,q and caches these values in a histogram
38(as implemented by RooAbsCachedPdf). The binning granularity of the cache
39can be controlled by the binning named "cache" on the RooRealVar representing
40the observable x. The fbar sampling algorithm is based on a recursive division
41mechanism with a built-in precision cutoff: First an initial sampling in
4264 equally spaced bins is made. Then the value of fbar is calculated in
43the center of each gap. If the calculated value deviates too much from
44the value obtained by linear interpolation from the edge bins, gap
45is recursively divided. This strategy makes it possible to define a very
46fine cache sampling (e.g. 1000 or 10000) bins without incurring a
47corresponding CPU penalty.
48
49Note on numeric stability of the algorithm. Since the algorithm relies
50on a numeric inversion of cumulative distributions functions, some precision
51may be lost at the 'edges' of the same (i.e. at regions in x where the
52c.d.f. value is close to zero or one). The general sampling strategy is
53to start with 64 equally spaces samples in the range y=(0.01-0.99).
54Then the y ranges are pushed outward by reducing y (or the distance of y to 1.0)
55by a factor of sqrt(10) iteratively up to the point where the corresponding
56x value no longer changes significantly. For p.d.f.s with very flat tails
57such as Gaussians some part of the tail may be lost due to limitations
58in numeric precision in the CDF inversion step.
59
60An effect related to the above limitation in numeric precision should
61be anticipated when floating the alpha parameter in a fit. If a p.d.f
62with such flat tails is fitted, it is likely that the dataset contains
63events in the flat tail region. If the alpha parameter is varied, the
64likelihood contribution from such events may exhibit discontinuities
65in alpha, causing discontinuities in the summed likelihood as well
66that will cause convergence problems in MINUIT. To mitigate this effect
67one can use the setCacheAlpha() method to instruct RooIntegralMorph
68to construct a two-dimensional cache for its output values in both
69x and alpha. If linear interpolation is requested on the resulting
70output histogram, the resulting interpolation of the p.d.f in the
71alpha dimension will smooth out the discontinuities in the tail regions
72result in a continuous likelihood distribution that can be fitted.
73An added advantage of the cacheAlpha option is that if parameters
74p,q of f1,f2 are fixed, the cached values in RooIntegralMorph are
75valid for the entire fit session and do not need to be recalculated
76for each change in alpha, which may result an considerable increase
77in calculation speed.
78
79**/
80
81#include "Riostream.h"
82
83#include "RooIntegralMorph.h"
84#include "RooAbsCategory.h"
85#include "RooBrentRootFinder.h"
86#include "RooAbsFunc.h"
87#include "RooRealVar.h"
88#include "RooDataHist.h"
89#include "TH1.h"
90
91using namespace std;
92
94
95////////////////////////////////////////////////////////////////////////////////
96/// Constructor with observables x, pdf shapes pdf1 and pdf2 which represent
97/// the shapes at the end points of the interpolation parameter alpha
98/// If doCacheAlpha is true, a two-dimensional cache is constructed in
99/// both alpha and x
100
101RooIntegralMorph::RooIntegralMorph(const char *name, const char *title,
102 RooAbsReal& _pdf1,
103 RooAbsReal& _pdf2,
104 RooAbsReal& _x,
105 RooAbsReal& _alpha,
106 bool doCacheAlpha) :
107 RooAbsCachedPdf(name,title,2),
108 pdf1("pdf1","pdf1",this,_pdf1),
109 pdf2("pdf2","pdf2",this,_pdf2),
110 x("x","x",this,_x),
111 alpha("alpha","alpha",this,_alpha),
112 _cacheAlpha(doCacheAlpha)
113{
114}
115
116////////////////////////////////////////////////////////////////////////////////
117/// Copy constructor
118
120 RooAbsCachedPdf(other,name),
121 pdf1("pdf1",this,other.pdf1),
122 pdf2("pdf2",this,other.pdf2),
123 x("x",this,other.x),
124 alpha("alpha",this,other.alpha),
125 _cacheAlpha(other._cacheAlpha)
126{
127}
128
129////////////////////////////////////////////////////////////////////////////////
130/// Observable to be cached for given choice of normalization.
131/// Returns the 'x' observable unless doCacheAlpha is set in which
132/// case a set with both x and alpha
133
135{
136 RooArgSet* obs = new RooArgSet ;
137 if (_cacheAlpha) {
138 obs->add(alpha.arg()) ;
139 }
140 obs->add(x.arg()) ;
142}
143
144////////////////////////////////////////////////////////////////////////////////
145/// Parameters of the cache. Returns parameters of both pdf1 and pdf2
146/// and parameter cache, in case doCacheAlpha is not set.
147
149{
150 auto par1 = pdf1->getParameters(static_cast<RooArgSet*>(nullptr));
151 RooArgSet par2;
152 pdf2->getParameters(nullptr, par2);
153 par1->add(par2,true) ;
154 par1->remove(x.arg(),true,true) ;
155 if (!_cacheAlpha) {
156 par1->add(alpha.arg()) ;
157 }
158 return RooFit::OwningPtr<RooArgSet>{std::move(par1)};
159}
160
161////////////////////////////////////////////////////////////////////////////////
162/// Return base name component for cache components in this case
163/// a string encoding the names of both end point p.d.f.s
164
166{
167 static TString name ;
168
169 name = pdf1.arg().GetName() ;
170 name.Append("_MORPH_") ;
171 name.Append(pdf2.arg().GetName()) ;
172 return name.Data() ;
173}
174
175////////////////////////////////////////////////////////////////////////////////
176/// Fill the cache with the interpolated shape.
177
179{
180 MorphCacheElem& mcache = static_cast<MorphCacheElem&>(cache) ;
181
182 // If cacheAlpha is true employ slice iterator here to fill all slices
183
184 if (!_cacheAlpha) {
185
186 std::unique_ptr<TIterator> dIter{cache.hist()->sliceIterator(const_cast<RooAbsReal&>(x.arg()),RooArgSet())};
187 mcache.calculate(dIter.get());
188
189 } else {
190 std::unique_ptr<TIterator> slIter{cache.hist()->sliceIterator(const_cast<RooAbsReal&>(alpha.arg()),RooArgSet())};
191
192 double alphaSave = alpha ;
193 RooArgSet alphaSet(alpha.arg()) ;
194 coutP(Eval) << "RooIntegralMorph::fillCacheObject(" << GetName() << ") filling multi-dimensional cache" ;
195 while(slIter->Next()) {
196 alphaSet.assign(*cache.hist()->get()) ;
197 std::unique_ptr<TIterator> dIter{cache.hist()->sliceIterator(const_cast<RooAbsReal&>(x.arg()),RooArgSet(alpha.arg()))};
198 mcache.calculate(dIter.get());
199 ccoutP(Eval) << "." << flush;
200 }
201 ccoutP(Eval) << std::endl;
202
203 const_cast<RooIntegralMorph*>(this)->alpha = alphaSave ;
204 }
205}
206
207////////////////////////////////////////////////////////////////////////////////
208/// Create and return a derived MorphCacheElem.
209
211{
212 return new MorphCacheElem(const_cast<RooIntegralMorph&>(*this),nset) ;
213}
214
215////////////////////////////////////////////////////////////////////////////////
216/// Return all RooAbsArg components contained in this cache
217
219{
220 RooArgList ret ;
221 ret.add(PdfCacheElem::containedArgs(action)) ;
222 ret.add(*_self) ;
223 ret.add(*_pdf1) ;
224 ret.add(*_pdf2) ;
225 ret.add(*_x ) ;
226 ret.add(*_alpha) ;
227 ret.add(*_c1) ;
228 ret.add(*_c2) ;
229
230 return ret ;
231}
232
233////////////////////////////////////////////////////////////////////////////////
234/// Construct of cache element, copy relevant input from RooIntegralMorph,
235/// create the cdfs from the input p.d.fs and instantiate the root finders
236/// on the cdfs to perform the inversion.
237
239{
240 // Mark in base class that normalization of cached pdf is invariant under pdf parameters
241 _x = (RooRealVar*)self.x.absArg() ;
242 _nset = std::make_unique<RooArgSet>(*_x);
243
244 _alpha = (RooAbsReal*)self.alpha.absArg() ;
245 _pdf1 = (RooAbsPdf*)(self.pdf1.absArg()) ;
246 _pdf2 = (RooAbsPdf*)(self.pdf2.absArg()) ;
247 _c1 = std::unique_ptr<RooAbsReal>{_pdf1->createCdf(*_x)};
248 _c2 = std::unique_ptr<RooAbsReal>{_pdf2->createCdf(*_x)};
249 _cb1 = std::unique_ptr<RooAbsFunc>{_c1->bindVars(*_x,_nset.get())};
250 _cb2 = std::unique_ptr<RooAbsFunc>{_c2->bindVars(*_x,_nset.get())};
251 _self = &self ;
252
253 _rf1 = std::make_unique<RooBrentRootFinder>(*_cb1);
254 _rf2 = std::make_unique<RooBrentRootFinder>(*_cb2);
255 _ccounter = 0 ;
256
257 _rf1->setTol(1e-12) ;
258 _rf2->setTol(1e-12) ;
259 _ycutoff = 1e-7 ;
260
261 // _yatX = 0 ;
262 // _calcX = 0 ;
263
264 // Must do this here too: fillCache() may not be called if cache contents is retrieved from EOcache
265 pdf()->setUnitNorm(true) ;
266
267 _yatXmax = 0 ;
268 _yatXmin = 0 ;
269}
270
271////////////////////////////////////////////////////////////////////////////////
272/// Destructor
273
275{
276}
277
278////////////////////////////////////////////////////////////////////////////////
279/// Calculate the x value of the output p.d.f at the given cdf value y.
280/// The ok boolean is filled with the success status of the operation.
281
283{
284 if (y<0 || y>1) {
285 oocoutW(_self,Eval) << "RooIntegralMorph::MorphCacheElem::calcX() WARNING: requested root finding for unphysical CDF value " << y << endl ;
286 }
287 double x1,x2 ;
288
289 double xmax = _x->getMax("cache") ;
290 double xmin = _x->getMin("cache") ;
291
292 ok=true ;
293 ok &= _rf1->findRoot(x1,xmin,xmax,y) ;
294 ok &= _rf2->findRoot(x2,xmin,xmax,y) ;
295 if (!ok) return 0 ;
296 _ccounter++ ;
297
298 return _alpha->getVal()*x1 + (1-_alpha->getVal())*x2 ;
299}
300
301////////////////////////////////////////////////////////////////////////////////
302/// Return the bin number enclosing the given x value
303
305{
306 double xmax = _x->getMax("cache") ;
307 double xmin = _x->getMin("cache") ;
308 return (Int_t)(_x->numBins("cache")*(X-xmin)/(xmax-xmin)) ;
309}
310
311////////////////////////////////////////////////////////////////////////////////
312/// Calculate shape of p.d.f for x,alpha values
313/// defined by dIter iterator over cache histogram
314
316{
317 double xsave = _self->x ;
318
319 // if (!_yatX) {
320 // _yatX = new double[_x->numBins("cache")+1] ;
321 // _calcX = new double[_x->numBins("cache")+1] ;
322 // }
323
324 _yatX.resize(_x->numBins("cache")+1);
325 _calcX.resize(_x->numBins("cache")+1);
326
327 _ccounter = 0 ;
328
329 // Get number of bins from PdfCacheElem histogram
330 Int_t nbins = _x->numBins("cache") ;
331
332 // Initialize yatX array to 'un-calculated values (-1)'
333 for (int i=0 ; i<nbins ; i++) {
334 _yatX[i] = -1 ;
335 _calcX[i] = 0 ;
336 }
337
338 // Find low and high point
339 findRange() ;
340
341 // Perform initial scan of 100 points
342 for (int i=0 ; i<10 ; i++) {
343
344 // Take a point in y
345 double offset = _yatX[_yatXmin] ;
346 double delta = (_yatX[_yatXmax] - _yatX[_yatXmin])/10 ;
347 double y = offset + i*delta ;
348
349 // Calculate corresponding X
350 bool ok ;
351 double X = calcX(y,ok) ;
352 if (ok) {
353 Int_t iX = binX(X) ;
354 _yatX[iX] = y ;
355 _calcX[iX] = X ;
356 }
357 }
358
359 // Now take an iteration filling the 'gaps'
360 Int_t igapLow = _yatXmin+1 ;
361 while(true) {
362 // Find next gap
363 Int_t igapHigh = igapLow+1 ;
364 while(igapHigh<(_yatXmax) && _yatX[igapHigh]<0) igapHigh++ ;
365
366 // Fill the gap (iteratively and/or using interpolation)
367 fillGap(igapLow-1,igapHigh) ;
368
369 // Terminate after processing of last gap
370 if (igapHigh>=_yatXmax-1) break ;
371 igapLow = igapHigh+1 ;
372 }
373
374 // Make one more iteration to recalculate Y value at bin centers
375 double xmax = _x->getMax("cache") ;
376 double xmin = _x->getMin("cache") ;
377 double binw = (xmax-xmin)/_x->numBins("cache") ;
378 for (int i=_yatXmin+1 ; i<_yatXmax-1 ; i++) {
379
380 // Calculate additional offset to apply if bin ixlo does not have X value calculated at bin center
381 double xBinC = xmin + (i+0.5)*binw ;
382 double xOffset = xBinC-_calcX[i] ;
383 if (std::abs(xOffset/binw)>1e-3) {
384 double slope = (_yatX[i+1]-_yatX[i-1])/(_calcX[i+1]-_calcX[i-1]) ;
385 double newY = _yatX[i] + slope*xOffset ;
386 //cout << "bin " << i << " needs to be re-centered " << xOffset/binw << " slope = " << slope << " origY = " << _yatX[i] << " newY = " << newY << endl ;
387 _yatX[i] = newY ;
388 }
389 }
390
391 // Zero output histogram below lowest calculable X value
392 for (int i=0; i<_yatXmin ; i++) {
393 dIter->Next() ;
394 //_hist->get(i) ;
395 hist()->set(0) ;
396 }
397
398 double x1 = _x->getMin("cache");
399 double x2 = _x->getMin("cache");
400
401 double xMax = _x->getMax("cache");
402
403 // Transfer calculated values to histogram
404 for (int i=_yatXmin ; i<_yatXmax ; i++) {
405
406 double y = _yatX[i] ;
407
408 // Little optimization here exploiting the fact that th cumulative
409 // distribution functions increase monotonically, so we already know that
410 // the next x-value must be higher than the last one as y is increasing. So
411 // we can use the previous x values as lower bounds.
412 _rf1->findRoot(x1,x1,xMax,y) ;
413 _rf2->findRoot(x2,x2,xMax,y) ;
414
415 _x->setVal(x1);
416 double f1x1 = _pdf1->getVal(_nset.get());
417 _x->setVal(x2);
418 double f2x2 = _pdf2->getVal(_nset.get());
419 double fbarX = f1x1*f2x2 / ( _alpha->getVal()*f2x2 + (1-_alpha->getVal())*f1x1 ) ;
420
421 dIter->Next() ;
422 //_hist->get(i) ;
423 hist()->set(fbarX) ;
424 }
425 // Zero output histogram above highest calculable X value
426 for (int i=_yatXmax+1 ; i<nbins ; i++) {
427 dIter->Next() ;
428 //_hist->get(i) ;
429 hist()->set(0) ;
430 }
431
432 pdf()->setUnitNorm(true) ;
433 _self->x = xsave ;
434
435 oocxcoutD(_self,Eval) << "RooIntegralMorph::MorphCacheElem::calculate(" << _self->GetName() << ") calculation required " << _ccounter << " samplings of cdfs" << endl ;
436}
437
438////////////////////////////////////////////////////////////////////////////////
439/// Fill all empty histogram bins between bins ixlo and ixhi. The value of 'splitPoint'
440/// defines the split point for the recursive division strategy to fill the gaps
441/// If the midpoint value of y is very close to the midpoint in x, use interpolation
442/// to fill the gaps, otherwise the intervals again.
443
444void RooIntegralMorph::MorphCacheElem::fillGap(Int_t ixlo, Int_t ixhi, double splitPoint)
445{
446 // CONVENTION: _yatX[ixlo] is filled, _yatX[ixhi] is filled, elements in between are empty
447 // cout << "fillGap: gap from _yatX[" << ixlo << "]=" << _yatX[ixlo] << " to _yatX[" << ixhi << "]=" << _yatX[ixhi] << ", size = " << ixhi-ixlo << endl ;
448
449 if (_yatX[ixlo]<0) {
450 oocoutE(_self,Eval) << "RooIntegralMorph::MorphCacheElme::fillGap(" << _self->GetName() << "): ERROR in fillgap " << ixlo << " = " << ixhi
451 << " splitPoint= " << splitPoint << " _yatX[ixlo] = " << _yatX[ixlo] << endl ;
452 }
453 if (_yatX[ixhi]<0) {
454 oocoutE(_self,Eval) << "RooIntegralMorph::MorphCacheElme::fillGap(" << _self->GetName() << "): ERROR in fillgap " << ixlo << " = " << ixhi
455 << " splitPoint " << splitPoint << " _yatX[ixhi] = " << _yatX[ixhi] << endl ;
456 }
457
458 // Determine where half-way Y value lands
459 double ymid = _yatX[ixlo]*splitPoint + _yatX[ixhi]*(1-splitPoint) ;
460 bool ok ;
461 double Xmid = calcX(ymid,ok) ;
462 if (!ok) {
463 oocoutW(_self,Eval) << "RooIntegralMorph::MorphCacheElem::fillGap(" << _self->GetName() << ") unable to calculate midpoint in gap ["
464 << ixlo << "," << ixhi << "], resorting to interpolation" << endl ;
465 interpolateGap(ixlo,ixhi) ;
466 }
467
468 Int_t iX = binX(Xmid) ;
469 double cq = (Xmid-_calcX[ixlo])/(_calcX[ixhi]-_calcX[ixlo])-0.5 ;
470
471 // Store midway point
472 _yatX[iX] = ymid ;
473 _calcX[iX] = Xmid ;
474
475
476 // Policy: If centration quality is better than 1% OR better than 1/10 of a bin, fill interval with linear interpolation
477 if (std::abs(cq)<0.01 || std::abs(cq*(ixhi-ixlo))<0.1 || ymid<_ycutoff ) {
478
479 // Fill remaining gaps on either side with linear interpolation
480 if (iX-ixlo>1) {
481 interpolateGap(ixlo,iX) ;
482 }
483 if (ixhi-iX>1) {
484 interpolateGap(iX,ixhi) ;
485 }
486
487 } else {
488
489 if (iX==ixlo) {
490
491 if (splitPoint<0.95) {
492 // Midway value lands on lowest bin, retry split with higher split point
493 double newSplit = splitPoint + 0.5*(1-splitPoint) ;
494 fillGap(ixlo,ixhi,newSplit) ;
495 } else {
496 // Give up and resort to interpolation
497 interpolateGap(ixlo,ixhi) ;
498 }
499
500 } else if (iX==ixhi) {
501
502 // Midway value lands on highest bin, retry split with lower split point
503 if (splitPoint>0.05) {
504 double newSplit = splitPoint/2 ;
505 fillGap(ixlo,ixhi,newSplit) ;
506 } else {
507 // Give up and resort to interpolation
508 interpolateGap(ixlo,ixhi) ;
509 }
510
511 } else {
512
513 // Midway point reasonable, iterate on interval on both sides
514 if (iX-ixlo>1) {
515 fillGap(ixlo,iX) ;
516 }
517 if (ixhi-iX>1) {
518 fillGap(iX,ixhi) ;
519 }
520 }
521 }
522}
523
524////////////////////////////////////////////////////////////////////////////////
525/// Fill empty histogram bins between ixlo and ixhi with values obtained
526/// from linear interpolation of ixlo,ixhi elements.
527
529{
530 //cout << "filling gap with linear interpolation ixlo=" << ixlo << " ixhi=" << ixhi << endl ;
531
532 double xmax = _x->getMax("cache") ;
533 double xmin = _x->getMin("cache") ;
534 double binw = (xmax-xmin)/_x->numBins("cache") ;
535
536 // Calculate deltaY in terms of actual X difference calculate, not based on nominal bin width
537 double deltaY = (_yatX[ixhi]-_yatX[ixlo])/((_calcX[ixhi]-_calcX[ixlo])/binw) ;
538
539 // Calculate additional offset to apply if bin ixlo does not have X value calculated at bin center
540 double xBinC = xmin + (ixlo+0.5)*binw ;
541 double xOffset = xBinC-_calcX[ixlo] ;
542
543 for (int j=ixlo+1 ; j<ixhi ; j++) {
544 _yatX[j] = _yatX[ixlo]+(xOffset+(j-ixlo))*deltaY ;
545 _calcX[j] = xmin + (j+0.5)*binw ;
546 }
547
548}
549
550////////////////////////////////////////////////////////////////////////////////
551/// Determine which range of y values can be mapped to x values
552/// from the numeric inversion of the input c.d.fs.
553/// Start with a y range of [0.1-0.9] and push boundaries
554/// outward with a factor of 1/sqrt(10). Stop iteration if
555/// inverted x values no longer change
556
558{
559 double xmin = _x->getMin("cache") ;
560 double xmax = _x->getMax("cache") ;
561 Int_t nbins = _x->numBins("cache") ;
562
563 double x1,x2 ;
564 bool ok = true ;
565 double ymin=0.1,yminSave(-1) ;
566 double Xsave(-1),Xlast=xmax ;
567
568 // Find lowest Y value that can be measured
569 // Start at 0.1 and iteratively lower limit by sqrt(10)
570 while(true) {
571 ok &= _rf1->findRoot(x1,xmin,xmax,ymin) ;
572 ok &= _rf2->findRoot(x2,xmin,xmax,ymin) ;
573 oocxcoutD(_self,Eval) << "RooIntegralMorph::MorphCacheElem::findRange(" << _self->GetName() << ") findMin: x1 = " << x1 << " x2 = " << x2 << " ok = " << (ok?"T":"F") << endl ;
574
575 // Terminate in case of non-convergence
576 if (!ok) break ;
577
578 // Terminate if value of X no longer moves by >0.1 bin size
579 double X = _alpha->getVal()*x1 + (1-_alpha->getVal())*x2 ;
580 if (std::abs(X-Xlast)/(xmax-xmin)<0.0001) {
581 break ;
582 }
583 Xlast=X ;
584
585 // Store new Y value
586 _yatXmin = (Int_t)(nbins*(X-xmin)/(xmax-xmin)) ;
587 _yatX[_yatXmin] = ymin ;
588 _calcX[_yatXmin] = X ;
589 yminSave = ymin ;
590 Xsave=X ;
591
592 // Reduce ymin by half an order of magnitude
593 ymin /=sqrt(10.) ;
594
595 // Emergency break
596 if (ymin<_ycutoff) break ;
597 }
598 _yatX[_yatXmin] = yminSave ;
599 _calcX[_yatXmin] = Xsave ;
600
601 // Find highest Y value that can be measured
602 // Start at 1 - 0.1 and iteratively lower delta by sqrt(10)
603 ok = true ;
604 double deltaymax=0.1, deltaymaxSave(-1) ;
605 Xlast=xmin ;
606 while(true) {
607 ok &= _rf1->findRoot(x1,xmin,xmax,1-deltaymax) ;
608 ok &= _rf2->findRoot(x2,xmin,xmax,1-deltaymax) ;
609
610 oocxcoutD(_self,Eval) << "RooIntegralMorph::MorphCacheElem::findRange(" << _self->GetName() << ") findMax: x1 = " << x1 << " x2 = " << x2 << " ok = " << (ok?"T":"F") << endl ;
611
612 // Terminate in case of non-convergence
613 if (!ok) break ;
614
615 // Terminate if value of X no longer moves by >0.1 bin size
616 double X = _alpha->getVal()*x1 + (1-_alpha->getVal())*x2 ;
617 if (std::abs(X-Xlast)/(xmax-xmin)<0.0001) {
618 break ;
619 }
620 Xlast=X ;
621
622 // Store new Y value
623 _yatXmax = (Int_t)(nbins*(X-xmin)/(xmax-xmin)) ;
624 _yatX[_yatXmax] = 1-deltaymax ;
625 _calcX[_yatXmax] = X ;
626 deltaymaxSave = deltaymax ;
627 Xsave=X ;
628
629 // Reduce ymin by half an order of magnitude
630 deltaymax /=sqrt(10.) ;
631
632 // Emergency break
633 if (deltaymax<_ycutoff) break ;
634 }
635
636 _yatX[_yatXmax] = 1-deltaymaxSave ;
637 _calcX[_yatXmax] = Xsave ;
638
639
640 // Initialize values out of range to 'out-of-range' (-2)
641 for (int i=0 ; i<_yatXmin ; i++) _yatX[i] = -2 ;
642 for (int i=_yatXmax+1 ; i<nbins; i++) _yatX[i] = -2 ;
643 oocxcoutD(_self,Eval) << "RooIntegralMorph::findRange(" << _self->GetName() << "): ymin = " << _yatX[_yatXmin] << " ymax = " << _yatX[_yatXmax] << endl;
644 oocxcoutD(_self,Eval) << "RooIntegralMorph::findRange(" << _self->GetName() << "): xmin = " << _calcX[_yatXmin] << " xmax = " << _calcX[_yatXmax] << endl;
645}
646
647////////////////////////////////////////////////////////////////////////////////
648/// Dummy
649
651{
652 return 0 ;
653}
654
655////////////////////////////////////////////////////////////////////////////////
656/// Indicate to the RooAbsCachedPdf base class that for the filling of the
657/// cache the traversal of the x should be in the innermost loop, to minimize
658/// recalculation of the one-dimensional internal cache for a fixed value of alpha
659
661{
662 // Put x last to minimize cache faulting
663 orderedObs.removeAll() ;
664
665 orderedObs.add(obs) ;
666 RooAbsArg* obsX = obs.find(x.arg().GetName()) ;
667 if (obsX) {
668 orderedObs.remove(*obsX) ;
669 orderedObs.add(*obsX) ;
670 }
671}
#define e(i)
Definition RSha256.hxx:103
#define coutP(a)
#define oocoutW(o, a)
#define oocxcoutD(o, a)
#define ccoutP(a)
#define oocoutE(o, a)
int Int_t
Definition RtypesCore.h:45
#define ClassImp(name)
Definition Rtypes.h:377
#define X(type, name)
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
char name[80]
Definition TGX11.cxx:110
float xmin
float ymin
float xmax
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:79
RooFit::OwningPtr< RooArgSet > getParameters(const RooAbsData *data, bool stripDisconnected=true) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
Abstract base class for p.d.f.s that need or want to cache their evaluate() output in a RooHistPdf de...
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
virtual bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false)
Remove the specified argument from our list.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
void assign(const RooAbsCollection &other) const
Sets the value, cache and constant attribute of any argument in our set that also appears in the othe...
RooAbsArg * find(const char *name) const
Find object with given name in list.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
RooFit::OwningPtr< RooAbsReal > createCdf(const RooArgSet &iset, const RooArgSet &nset=RooArgSet())
Create a cumulative distribution function of this p.d.f in terms of the observables listed in iset.
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
virtual double offset() const
Definition RooAbsReal.h:378
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooAbsArg * absArg() const
Return pointer to contained argument.
Definition RooArgProxy.h:46
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
TIterator * sliceIterator(RooAbsArg &sliceArg, const RooArgSet &otherArgs)
Create an iterator over all bins in a slice defined by the subset of observables listed in sliceArg.
const RooArgSet * get() const override
Get bin centre of current bin.
Definition RooDataHist.h:76
void setUnitNorm(bool flag)
Definition RooHistPdf.h:79
std::unique_ptr< RooBrentRootFinder > _rf1
std::unique_ptr< RooAbsReal > _c1
void calculate(TIterator *iter)
Calculate shape of p.d.f for x,alpha values defined by dIter iterator over cache histogram.
void interpolateGap(Int_t ixlo, Int_t ixhi)
Fill empty histogram bins between ixlo and ixhi with values obtained from linear interpolation of ixl...
MorphCacheElem(RooIntegralMorph &self, const RooArgSet *nset)
Construct of cache element, copy relevant input from RooIntegralMorph, create the cdfs from the input...
std::unique_ptr< RooAbsFunc > _cb2
void fillGap(Int_t ixlo, Int_t ixhi, double splitPoint=0.5)
Fill all empty histogram bins between bins ixlo and ixhi.
std::unique_ptr< RooArgSet > _nset
void findRange()
Determine which range of y values can be mapped to x values from the numeric inversion of the input c...
std::unique_ptr< RooAbsFunc > _cb1
std::unique_ptr< RooBrentRootFinder > _rf2
RooArgList containedArgs(Action) override
Return all RooAbsArg components contained in this cache.
std::unique_ptr< RooAbsReal > _c2
double calcX(double y, bool &ok)
Calculate the x value of the output p.d.f at the given cdf value y.
Int_t binX(double x)
Return the bin number enclosing the given x value.
Class RooIntegralMorph is an implementation of the histogram interpolation technique described by Ale...
RooIntegralMorph()=default
RooFit::OwningPtr< RooArgSet > actualObservables(const RooArgSet &nset) const override
Observable to be cached for given choice of normalization.
PdfCacheElem * createCache(const RooArgSet *nset) const override
Create and return a derived MorphCacheElem.
const char * inputBaseName() const override
Return base name component for cache components in this case a string encoding the names of both end ...
void preferredObservableScanOrder(const RooArgSet &obs, RooArgSet &orderedObs) const override
Indicate to the RooAbsCachedPdf base class that for the filling of the cache the traversal of the x s...
void fillCacheObject(PdfCacheElem &cache) const override
Fill the cache with the interpolated shape.
double evaluate() const override
Dummy.
RooFit::OwningPtr< RooArgSet > actualParameters(const RooArgSet &nset) const override
Parameters of the cache.
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:37
const T & arg() const
Return reference to object held in proxy.
Iterator abstract base class.
Definition TIterator.h:30
virtual TObject * Next()=0
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
Basic string class.
Definition TString.h:139
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
T * OwningPtr
An alias for raw pointers for indicating that the return type of a RooFit function is an owning point...
Definition Config.h:43