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