Logo ROOT   6.08/07
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  _yatX.resize(_x->numBins("cache")+1);
357  _calcX.resize(_x->numBins("cache")+1);
358 
359  RooArgSet nsetTmp(*_x) ;
360  _ccounter = 0 ;
361 
362  // Get number of bins from PdfCacheElem histogram
363  Int_t nbins = _x->numBins("cache") ;
364 
365  // Initialize yatX array to 'uncalculated values (-1)'
366  for (int i=0 ; i<nbins ; i++) {
367  _yatX[i] = -1 ;
368  _calcX[i] = 0 ;
369  }
370 
371  // Find low and high point
372  findRange() ;
373 
374  // Perform initial scan of 100 points
375  for (int i=0 ; i<10 ; i++) {
376 
377  // Take a point in y
379  Double_t delta = (_yatX[_yatXmax] - _yatX[_yatXmin])/10 ;
380  Double_t y = offset + i*delta ;
381 
382  // Calculate corresponding X
383  Bool_t ok ;
384  Double_t X = calcX(y,ok) ;
385  if (ok) {
386  Int_t iX = binX(X) ;
387  _yatX[iX] = y ;
388  _calcX[iX] = X ;
389  }
390  }
391 
392  // Now take an iteration filling the 'gaps'
393  Int_t igapLow = _yatXmin+1 ;
394  while(true) {
395  // Find next gap
396  Int_t igapHigh = igapLow+1 ;
397  while(_yatX[igapHigh]<0 && igapHigh<(_yatXmax)) igapHigh++ ;
398 
399  // Fill the gap (iteratively and/or using interpolation)
400  fillGap(igapLow-1,igapHigh) ;
401 
402  // Terminate after processing of last gap
403  if (igapHigh>=_yatXmax-1) break ;
404  igapLow = igapHigh+1 ;
405  }
406 
407  // Make one more iteration to recalculate Y value at bin centers
408  Double_t xmax = _x->getMax("cache") ;
409  Double_t xmin = _x->getMin("cache") ;
410  Double_t binw = (xmax-xmin)/_x->numBins("cache") ;
411  for (int i=_yatXmin+1 ; i<_yatXmax-1 ; i++) {
412 
413  // Calculate additional offset to apply if bin ixlo does not have X value calculated at bin center
414  Double_t xBinC = xmin + (i+0.5)*binw ;
415  Double_t xOffset = xBinC-_calcX[i] ;
416  if (fabs(xOffset/binw)>1e-3) {
417  Double_t slope = (_yatX[i+1]-_yatX[i-1])/(_calcX[i+1]-_calcX[i-1]) ;
418  Double_t newY = _yatX[i] + slope*xOffset ;
419  //cout << "bin " << i << " needs to be recentered " << xOffset/binw << " slope = " << slope << " origY = " << _yatX[i] << " newY = " << newY << endl ;
420  _yatX[i] = newY ;
421  }
422  }
423 
424  // Zero output histogram below lowest calculable X value
425  for (int i=0; i<_yatXmin ; i++) {
426  dIter->Next() ;
427  //_hist->get(i) ;
428  hist()->set(0) ;
429  }
430  // Transfer calculated values to histogram
431  for (int i=_yatXmin ; i<_yatXmax ; i++) {
432 
433  Double_t y = _yatX[i] ;
434 
435  Double_t x1,x2 ;
436 
437  Double_t xMin = _x->getMin("cache") ;
438  Double_t xMax = _x->getMax("cache") ;
439  _rf1->findRoot(x1,xMin,xMax,y) ;
440  _rf2->findRoot(x2,xMin,xMax,y) ;
441 
442  _x->setVal(x1) ; Double_t f1x1 = _pdf1->getVal(&nsetTmp) ;
443  _x->setVal(x2) ; Double_t f2x2 = _pdf2->getVal(&nsetTmp) ;
444  Double_t fbarX = f1x1*f2x2 / ( _alpha->getVal()*f2x2 + (1-_alpha->getVal())*f1x1 ) ;
445 
446  dIter->Next() ;
447  //_hist->get(i) ;
448  hist()->set(fbarX) ;
449  }
450  // Zero output histogram above highest calculable X value
451  for (int i=_yatXmax+1 ; i<nbins ; i++) {
452  dIter->Next() ;
453  //_hist->get(i) ;
454  hist()->set(0) ;
455  }
456 
457  pdf()->setUnitNorm(kTRUE) ;
458  _self->x = xsave ;
459 
460  oocxcoutD(_self,Eval) << "RooIntegralMorph::MorphCacheElem::calculate(" << _self->GetName() << ") calculation required " << _ccounter << " samplings of cdfs" << endl ;
461 }
462 
463 
464 ////////////////////////////////////////////////////////////////////////////////
465 /// Fill all empty histogram bins between bins ixlo and ixhi. The value of 'splitPoint'
466 /// defines the split point for the recursive division strategy to fill the gaps
467 /// If the midpoint value of y is very close to the midpoint in x, use interpolation
468 /// to fill the gaps, otherwise the intervals again.
469 
471 {
472  // CONVENTION: _yatX[ixlo] is filled, _yatX[ixhi] is filled, elements in between are empty
473  // cout << "fillGap: gap from _yatX[" << ixlo << "]=" << _yatX[ixlo] << " to _yatX[" << ixhi << "]=" << _yatX[ixhi] << ", size = " << ixhi-ixlo << endl ;
474 
475  if (_yatX[ixlo]<0) {
476  oocoutE(_self,Eval) << "RooIntegralMorph::MorphCacheElme::fillGap(" << _self->GetName() << "): ERROR in fillgap " << ixlo << " = " << ixhi
477  << " splitPoint= " << splitPoint << " _yatX[ixlo] = " << _yatX[ixlo] << endl ;
478  }
479  if (_yatX[ixhi]<0) {
480  oocoutE(_self,Eval) << "RooIntegralMorph::MorphCacheElme::fillGap(" << _self->GetName() << "): ERROR in fillgap " << ixlo << " = " << ixhi
481  << " splitPoint " << splitPoint << " _yatX[ixhi] = " << _yatX[ixhi] << endl ;
482  }
483 
484  // Determine where half-way Y value lands
485  Double_t ymid = _yatX[ixlo]*splitPoint + _yatX[ixhi]*(1-splitPoint) ;
486  Bool_t ok ;
487  Double_t Xmid = calcX(ymid,ok) ;
488  if (!ok) {
489  oocoutW(_self,Eval) << "RooIntegralMorph::MorphCacheElem::fillGap(" << _self->GetName() << ") unable to calculate midpoint in gap ["
490  << ixlo << "," << ixhi << "], resorting to interpolation" << endl ;
491  interpolateGap(ixlo,ixhi) ;
492  }
493 
494  Int_t iX = binX(Xmid) ;
495  Double_t cq = (Xmid-_calcX[ixlo])/(_calcX[ixhi]-_calcX[ixlo])-0.5 ;
496 
497  // Store midway point
498  _yatX[iX] = ymid ;
499  _calcX[iX] = Xmid ;
500 
501 
502  // Policy: If centration quality is better than 1% OR better than 1/10 of a bin, fill interval with linear interpolation
503  if (fabs(cq)<0.01 || fabs(cq*(ixhi-ixlo))<0.1 || ymid<_ycutoff ) {
504 
505  // Fill remaining gaps on either side with linear interpolation
506  if (iX-ixlo>1) {
507  interpolateGap(ixlo,iX) ;
508  }
509  if (ixhi-iX>1) {
510  interpolateGap(iX,ixhi) ;
511  }
512 
513  } else {
514 
515  if (iX==ixlo) {
516 
517  if (splitPoint<0.95) {
518  // Midway value lands on lowest bin, retry split with higher split point
519  Double_t newSplit = splitPoint + 0.5*(1-splitPoint) ;
520  fillGap(ixlo,ixhi,newSplit) ;
521  } else {
522  // Give up and resort to interpolation
523  interpolateGap(ixlo,ixhi) ;
524  }
525 
526  } else if (iX==ixhi) {
527 
528  // Midway value lands on highest bin, retry split with lower split point
529  if (splitPoint>0.05) {
530  Double_t newSplit = splitPoint/2 ;
531  fillGap(ixlo,ixhi,newSplit) ;
532  } else {
533  // Give up and resort to interpolation
534  interpolateGap(ixlo,ixhi) ;
535  }
536 
537  } else {
538 
539  // Midway point reasonable, iterate on interval on both sides
540  if (iX-ixlo>1) {
541  fillGap(ixlo,iX) ;
542  }
543  if (ixhi-iX>1) {
544  fillGap(iX,ixhi) ;
545  }
546  }
547  }
548 }
549 
550 
551 ////////////////////////////////////////////////////////////////////////////////
552 /// Fill empty histogram bins between ixlo and ixhi with values obtained
553 /// from linear interpolation of ixlo,ixhi elements.
554 
556 {
557  //cout << "filling gap with linear interpolation ixlo=" << ixlo << " ixhi=" << ixhi << endl ;
558 
559  Double_t xmax = _x->getMax("cache") ;
560  Double_t xmin = _x->getMin("cache") ;
561  Double_t binw = (xmax-xmin)/_x->numBins("cache") ;
562 
563  // Calculate deltaY in terms of actuall X difference calculate, not based on nominal bin width
564  Double_t deltaY = (_yatX[ixhi]-_yatX[ixlo])/((_calcX[ixhi]-_calcX[ixlo])/binw) ;
565 
566  // Calculate additional offset to apply if bin ixlo does not have X value calculated at bin center
567  Double_t xBinC = xmin + (ixlo+0.5)*binw ;
568  Double_t xOffset = xBinC-_calcX[ixlo] ;
569 
570  for (int j=ixlo+1 ; j<ixhi ; j++) {
571  _yatX[j] = _yatX[ixlo]+(xOffset+(j-ixlo))*deltaY ;
572  _calcX[j] = xmin + (j+0.5)*binw ;
573  }
574 
575 }
576 
577 
578 
579 ////////////////////////////////////////////////////////////////////////////////
580 /// Determine which range of y values can be mapped to x values
581 /// from the numeric inversion of the input c.d.fs.
582 /// Start with a y rannge of [0.1-0.9] and push boundaries
583 /// outward with a factor of 1/sqrt(10). Stop iteration if
584 /// inverted x values no longer change
585 
587 {
588  Double_t xmin = _x->getMin("cache") ;
589  Double_t xmax = _x->getMax("cache") ;
590  Int_t nbins = _x->numBins("cache") ;
591 
592  Double_t x1,x2 ;
593  Bool_t ok = kTRUE ;
594  Double_t ymin=0.1,yminSave(-1) ;
595  Double_t Xsave(-1),Xlast=xmax ;
596 
597  // Find lowest Y value that can be measured
598  // Start at 0.1 and iteratively lower limit by sqrt(10)
599  while(true) {
600  ok &= _rf1->findRoot(x1,xmin,xmax,ymin) ;
601  ok &= _rf2->findRoot(x2,xmin,xmax,ymin) ;
602  oocxcoutD(_self,Eval) << "RooIntegralMorph::MorphCacheElem::findRange(" << _self->GetName() << ") findMin: x1 = " << x1 << " x2 = " << x2 << " ok = " << (ok?"T":"F") << endl ;
603 
604  // Terminate in case of non-convergence
605  if (!ok) break ;
606 
607  // Terminate if value of X no longer moves by >0.1 bin size
608  Double_t X = _alpha->getVal()*x1 + (1-_alpha->getVal())*x2 ;
609  if (fabs(X-Xlast)/(xmax-xmin)<0.0001) {
610  break ;
611  }
612  Xlast=X ;
613 
614  // Store new Y value
615  _yatXmin = (Int_t)(nbins*(X-xmin)/(xmax-xmin)) ;
616  _yatX[_yatXmin] = ymin ;
617  _calcX[_yatXmin] = X ;
618  yminSave = ymin ;
619  Xsave=X ;
620 
621  // Reduce ymin by half an order of magnitude
622  ymin /=sqrt(10.) ;
623 
624  // Emergency break
625  if (ymin<_ycutoff) break ;
626  }
627  _yatX[_yatXmin] = yminSave ;
628  _calcX[_yatXmin] = Xsave ;
629 
630  // Find highst Y value that can be measured
631  // Start at 1 - 0.1 and iteratively lower delta by sqrt(10)
632  ok = kTRUE ;
633  Double_t deltaymax=0.1, deltaymaxSave(-1) ;
634  Xlast=xmin ;
635  while(true) {
636  ok &= _rf1->findRoot(x1,xmin,xmax,1-deltaymax) ;
637  ok &= _rf2->findRoot(x2,xmin,xmax,1-deltaymax) ;
638 
639  oocxcoutD(_self,Eval) << "RooIntegralMorph::MorphCacheElem::findRange(" << _self->GetName() << ") findMax: x1 = " << x1 << " x2 = " << x2 << " ok = " << (ok?"T":"F") << endl ;
640 
641  // Terminate in case of non-convergence
642  if (!ok) break ;
643 
644  // Terminate if value of X no longer moves by >0.1 bin size
645  Double_t X = _alpha->getVal()*x1 + (1-_alpha->getVal())*x2 ;
646  if (fabs(X-Xlast)/(xmax-xmin)<0.0001) {
647  break ;
648  }
649  Xlast=X ;
650 
651  // Store new Y value
652  _yatXmax = (Int_t)(nbins*(X-xmin)/(xmax-xmin)) ;
653  _yatX[_yatXmax] = 1-deltaymax ;
654  _calcX[_yatXmax] = X ;
655  deltaymaxSave = deltaymax ;
656  Xsave=X ;
657 
658  // Reduce ymin by half an order of magnitude
659  deltaymax /=sqrt(10.) ;
660 
661  // Emergency break
662  if (deltaymax<_ycutoff) break ;
663  }
664 
665  _yatX[_yatXmax] = 1-deltaymaxSave ;
666  _calcX[_yatXmax] = Xsave ;
667 
668 
669  // Initialize values out of range to 'out-of-range' (-2)
670  for (int i=0 ; i<_yatXmin ; i++) _yatX[i] = -2 ;
671  for (int i=_yatXmax+1 ; i<nbins; i++) _yatX[i] = -2 ;
672  oocxcoutD(_self,Eval) << "RooIntegralMorph::findRange(" << _self->GetName() << "): ymin = " << _yatX[_yatXmin] << " ymax = " << _yatX[_yatXmax] << endl;
673  oocxcoutD(_self,Eval) << "RooIntegralMorph::findRange(" << _self->GetName() << "): xmin = " << _calcX[_yatXmin] << " xmax = " << _calcX[_yatXmax] << endl;
674 }
675 
676 
677 
678 ////////////////////////////////////////////////////////////////////////////////
679 /// Dummy
680 
682 {
683  return 0 ;
684 }
685 
686 
687 
688 ////////////////////////////////////////////////////////////////////////////////
689 /// Indicate to the RooAbsCachedPdf base class that for the filling of the
690 /// cache the traversal of the x should be in the innermost loop, to minimize
691 /// recalculation of the one-dimensional internal cache for a fixed value of alpha
692 
694 {
695  // Put x last to minimize cache faulting
696  orderedObs.removeAll() ;
697 
698  orderedObs.add(obs) ;
699  RooAbsArg* obsX = obs.find(x.arg().GetName()) ;
700  if (obsX) {
701  orderedObs.remove(*obsX) ;
702  orderedObs.add(*obsX) ;
703  }
704 }
705 
706 
707 
virtual Double_t getMin(const char *name=0) const
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
float xmin
Definition: THbookFile.cxx:93
virtual Double_t getMax(const char *name=0) const
virtual Bool_t add(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling add() for each element in the source coll...
Definition: RooArgSet.h:86
virtual const RooArgSet * get() const
Definition: RooDataHist.h:77
std::vector< Double_t > _calcX
std::vector< Double_t > _yatX
virtual RooArgSet * actualObservables(const RooArgSet &nset) const
Observable to be cached for given choice of normalization.
float ymin
Definition: THbookFile.cxx:93
virtual Int_t numBins(const char *rangeName=0) const
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
void set(Double_t weight, Double_t wgtErr=-1)
Increment the weight of the bin enclosing the coordinates given by &#39;row&#39; by the specified amount...
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.
Iterator abstract base class.
Definition: TIterator.h:32
#define coutP(a)
Definition: RooMsgService.h:33
virtual RooArgSet * actualParameters(const RooArgSet &nset) const
Parameters of the cache.
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
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
friend class RooArgSet
Definition: RooAbsArg.h:469
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
virtual Double_t offset() const
Definition: RooAbsReal.h:309
virtual void setVal(Double_t value)
Set value of variable to &#39;value&#39;.
Definition: RooRealVar.cxx:205
virtual Bool_t findRoot(Double_t &result, Double_t xlo, Double_t xhi, Double_t value=0) const
Do the root finding using the Brent-Decker method.
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.
virtual const char * inputBaseName() const
Return base name component for cache components in this case a string encoding the names of both end ...
void findRange()
Determine which range of y values can be mapped to x values from the numeric inversion of the input c...
Implement the abstract 1-dimensional root finding interface using the Brent-Decker method...
RooAbsArg * absArg() const
Definition: RooArgProxy.h:37
void setTol(Double_t tol)
float xmax
Definition: THbookFile.cxx:93
#define ccoutP(a)
Definition: RooMsgService.h:40
virtual void fillCacheObject(PdfCacheElem &cache) const
Fill the cache with the interpolated shape.
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]
RooAbsCachedPdf is the abstract base class for p.d.f.s that need or want to cache their evaluate() ou...
#define ClassImp(name)
Definition: Rtypes.h:279
RooAbsArg * find(const char *name) const
Find object with given name in list.
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
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)...
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&#39;t match any of...
Definition: RooAbsArg.cxx:560
RooRealProxy pdf1
Double_t y[n]
Definition: legend1.C:17
RooRealProxy alpha
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...
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
Definition: TRolke.cxx:630
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 oocoutW(o, a)
Definition: RooMsgService.h:47
virtual RooArgList containedArgs(Action)
Return all RooAbsArg components contained in this cache.
Double_t evaluate() const
Dummy.
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
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
const RooAbsReal & arg() const
Definition: RooRealProxy.h:43
friend class MorphCacheElem
const Bool_t kTRUE
Definition: Rtypes.h:91
virtual PdfCacheElem * createCache(const RooArgSet *nset) const
Create and return a derived MorphCacheElem.
void calculate(TIterator *iter)
Calculate shape of p.d.f for x,alpha values defined by dIter iterator over cache histogram.
char name[80]
Definition: TGX11.cxx:109