108 pdf1(
"pdf1",
"pdf1",this,_pdf1),
109 pdf2(
"pdf2",
"pdf2",this,_pdf2),
111 alpha(
"alpha",
"alpha",this,_alpha),
112 _cacheAlpha(doCacheAlpha)
121 pdf1(
"pdf1",this,other.pdf1),
122 pdf2(
"pdf2",this,other.pdf2),
124 alpha(
"alpha",this,other.alpha),
125 _cacheAlpha(other._cacheAlpha)
153 par1->add(par2,
true) ;
154 par1->remove(
x.
arg(),
true,
true) ;
170 name.Append(
"_MORPH_") ;
192 double alphaSave =
alpha ;
194 coutP(Eval) <<
"RooIntegralMorph::fillCacheObject(" <<
GetName() <<
") filling multi-dimensional cache" ;
195 while(slIter->Next()) {
199 ccoutP(Eval) <<
"." << flush;
201 ccoutP(Eval) << std::endl;
221 ret.
add(PdfCacheElem::containedArgs(action)) ;
242 _nset = std::make_unique<RooArgSet>(*
_x);
249 _cb1 = std::unique_ptr<RooAbsFunc>{
_c1->bindVars(*
_x,
_nset.get())};
250 _cb2 = std::unique_ptr<RooAbsFunc>{
_c2->bindVars(*
_x,
_nset.get())};
253 _rf1 = std::make_unique<RooBrentRootFinder>(*
_cb1);
254 _rf2 = std::make_unique<RooBrentRootFinder>(*
_cb2);
257 _rf1->setTol(1
e-12) ;
258 _rf2->setTol(1
e-12) ;
285 oocoutW(_self,Eval) <<
"RooIntegralMorph::MorphCacheElem::calcX() WARNING: requested root finding for unphysical CDF value " <<
y << endl ;
289 double xmax = _x->getMax(
"cache") ;
290 double xmin = _x->getMin(
"cache") ;
298 return _alpha->getVal()*
x1 + (1-_alpha->getVal())*
x2 ;
306 double xmax = _x->getMax(
"cache") ;
307 double xmin = _x->getMin(
"cache") ;
317 double xsave = _self->x ;
324 _yatX.resize(_x->numBins(
"cache")+1);
325 _calcX.resize(_x->numBins(
"cache")+1);
330 Int_t nbins = _x->numBins(
"cache") ;
333 for (
int i=0 ; i<nbins ; i++) {
342 for (
int i=0 ; i<10 ; i++) {
345 double offset = _yatX[_yatXmin] ;
346 double delta = (_yatX[_yatXmax] - _yatX[_yatXmin])/10 ;
351 double X = calcX(
y,ok) ;
360 Int_t igapLow = _yatXmin+1 ;
363 Int_t igapHigh = igapLow+1 ;
364 while(igapHigh<(_yatXmax) && _yatX[igapHigh]<0) igapHigh++ ;
367 fillGap(igapLow-1,igapHigh) ;
370 if (igapHigh>=_yatXmax-1) break ;
371 igapLow = igapHigh+1 ;
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++) {
381 double xBinC =
xmin + (i+0.5)*binw ;
382 double xOffset = xBinC-_calcX[i] ;
383 if (std::abs(xOffset/binw)>1
e-3) {
384 double slope = (_yatX[i+1]-_yatX[i-1])/(_calcX[i+1]-_calcX[i-1]) ;
385 double newY = _yatX[i] + slope*xOffset ;
392 for (
int i=0; i<_yatXmin ; i++) {
398 double x1 = _x->getMin(
"cache");
399 double x2 = _x->getMin(
"cache");
401 double xMax = _x->getMax(
"cache");
404 for (
int i=_yatXmin ; i<_yatXmax ; i++) {
406 double y = _yatX[i] ;
412 _rf1->findRoot(
x1,
x1,xMax,
y) ;
413 _rf2->findRoot(
x2,
x2,xMax,
y) ;
416 double f1x1 = _pdf1->getVal(_nset.get());
418 double f2x2 = _pdf2->getVal(_nset.get());
419 double fbarX = f1x1*f2x2 / ( _alpha->getVal()*f2x2 + (1-_alpha->getVal())*f1x1 ) ;
426 for (
int i=_yatXmax+1 ; i<nbins ; i++) {
432 pdf()->setUnitNorm(
true) ;
435 oocxcoutD(_self,Eval) <<
"RooIntegralMorph::MorphCacheElem::calculate(" << _self->GetName() <<
") calculation required " << _ccounter <<
" samplings of cdfs" << endl ;
450 oocoutE(_self,Eval) <<
"RooIntegralMorph::MorphCacheElme::fillGap(" << _self->GetName() <<
"): ERROR in fillgap " << ixlo <<
" = " << ixhi
451 <<
" splitPoint= " << splitPoint <<
" _yatX[ixlo] = " << _yatX[ixlo] << endl ;
454 oocoutE(_self,Eval) <<
"RooIntegralMorph::MorphCacheElme::fillGap(" << _self->GetName() <<
"): ERROR in fillgap " << ixlo <<
" = " << ixhi
455 <<
" splitPoint " << splitPoint <<
" _yatX[ixhi] = " << _yatX[ixhi] << endl ;
459 double ymid = _yatX[ixlo]*splitPoint + _yatX[ixhi]*(1-splitPoint) ;
461 double Xmid = calcX(ymid,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) ;
468 Int_t iX = binX(Xmid) ;
469 double cq = (Xmid-_calcX[ixlo])/(_calcX[ixhi]-_calcX[ixlo])-0.5 ;
477 if (std::abs(cq)<0.01 || std::abs(cq*(ixhi-ixlo))<0.1 || ymid<_ycutoff ) {
481 interpolateGap(ixlo,iX) ;
484 interpolateGap(iX,ixhi) ;
491 if (splitPoint<0.95) {
493 double newSplit = splitPoint + 0.5*(1-splitPoint) ;
494 fillGap(ixlo,ixhi,newSplit) ;
497 interpolateGap(ixlo,ixhi) ;
500 }
else if (iX==ixhi) {
503 if (splitPoint>0.05) {
504 double newSplit = splitPoint/2 ;
505 fillGap(ixlo,ixhi,newSplit) ;
508 interpolateGap(ixlo,ixhi) ;
532 double xmax = _x->getMax(
"cache") ;
533 double xmin = _x->getMin(
"cache") ;
534 double binw = (
xmax-
xmin)/_x->numBins(
"cache") ;
537 double deltaY = (_yatX[ixhi]-_yatX[ixlo])/((_calcX[ixhi]-_calcX[ixlo])/binw) ;
540 double xBinC =
xmin + (ixlo+0.5)*binw ;
541 double xOffset = xBinC-_calcX[ixlo] ;
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 ;
559 double xmin = _x->getMin(
"cache") ;
560 double xmax = _x->getMax(
"cache") ;
561 Int_t nbins = _x->numBins(
"cache") ;
565 double ymin=0.1,yminSave(-1) ;
566 double Xsave(-1),Xlast=
xmax ;
573 oocxcoutD(_self,Eval) <<
"RooIntegralMorph::MorphCacheElem::findRange(" << _self->GetName() <<
") findMin: x1 = " <<
x1 <<
" x2 = " <<
x2 <<
" ok = " << (ok?
"T":
"F") << endl ;
579 double X = _alpha->getVal()*
x1 + (1-_alpha->getVal())*
x2 ;
580 if (std::abs(
X-Xlast)/(
xmax-
xmin)<0.0001) {
587 _yatX[_yatXmin] =
ymin ;
588 _calcX[_yatXmin] =
X ;
596 if (
ymin<_ycutoff) break ;
598 _yatX[_yatXmin] = yminSave ;
599 _calcX[_yatXmin] = Xsave ;
604 double deltaymax=0.1, deltaymaxSave(-1) ;
607 ok &= _rf1->findRoot(
x1,
xmin,
xmax,1-deltaymax) ;
608 ok &= _rf2->findRoot(
x2,
xmin,
xmax,1-deltaymax) ;
610 oocxcoutD(_self,Eval) <<
"RooIntegralMorph::MorphCacheElem::findRange(" << _self->GetName() <<
") findMax: x1 = " <<
x1 <<
" x2 = " <<
x2 <<
" ok = " << (ok?
"T":
"F") << endl ;
616 double X = _alpha->getVal()*
x1 + (1-_alpha->getVal())*
x2 ;
617 if (std::abs(
X-Xlast)/(
xmax-
xmin)<0.0001) {
624 _yatX[_yatXmax] = 1-deltaymax ;
625 _calcX[_yatXmax] =
X ;
626 deltaymaxSave = deltaymax ;
630 deltaymax /=sqrt(10.) ;
633 if (deltaymax<_ycutoff) break ;
636 _yatX[_yatXmax] = 1-deltaymaxSave ;
637 _calcX[_yatXmax] = Xsave ;
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;
665 orderedObs.
add(obs) ;
668 orderedObs.
remove(*obsX) ;
669 orderedObs.
add(*obsX) ;
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Common abstract base class for objects that represent a value and a "shape" in RooFit.
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.
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...
virtual double offset() const
RooArgList is a container object that can hold multiple RooAbsArg objects.
RooAbsArg * absArg() const
Return pointer to contained argument.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
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.
void setUnitNorm(bool flag)
~MorphCacheElem() override
Destructor.
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.
const T & arg() const
Return reference to object held in proxy.
Iterator abstract base class.
virtual TObject * Next()=0
const char * GetName() const override
Returns name of object.
T * OwningPtr
An alias for raw pointers for indicating that the return type of a RooFit function is an owning point...