Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooCurve.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * @(#)root/roofitcore:$Id$
5 * Authors: *
6 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * and Stanford University. All rights reserved. *
11 * *
12 * Redistribution and use in source and binary forms, *
13 * with or without modification, are permitted according to the terms *
14 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15 *****************************************************************************/
16
17/**
18\file RooCurve.cxx
19\class RooCurve
20\ingroup Roofitcore
21
22A RooCurve is a one-dimensional graphical representation of a real-valued function.
23A curve is approximated by straight line segments with end points chosen to give
24a "good" approximation to the true curve. The goodness of the approximation is
25controlled by a precision and a resolution parameter.
26
27A RooCurve derives from TGraph, so it can either be drawn as a line (default) or
28as points:
29```
30RooPlot *p = y.plotOn(x.frame());
31p->getAttMarker("curve_y")->SetMarkerStyle(20);
32p->setDrawOptions("curve_y","PL");
33p->Draw();
34```
35
36To retrieve a RooCurve from a RooPlot, use RooPlot::getCurve().
37**/
38
39#include "RooCurve.h"
40#include "RooHist.h"
41#include "RooAbsReal.h"
42#include "RooArgSet.h"
43#include "RooRealVar.h"
44#include "RooRealIntegral.h"
45#include "RooRealBinding.h"
46#include "RooMsgService.h"
47#include "RooProduct.h"
48#include "RooConstVar.h"
49
50#include "Riostream.h"
51#include "TMath.h"
52#include "TAxis.h"
53#include "TMatrixD.h"
54#include "TVectorD.h"
55#include "Math/Util.h"
56#include <iomanip>
57#include <deque>
58#include <algorithm>
59
60using namespace std ;
61
63
64
65////////////////////////////////////////////////////////////////////////////////
66/// Default constructor.
67
69{
70 initialize();
71}
72
73
74////////////////////////////////////////////////////////////////////////////////
75/// Create a 1-dim curve of the value of the specified real-valued expression
76/// as a function of x. Use the optional precision parameter to control
77/// how precisely the smooth curve is rasterized. Use the optional argument set
78/// to specify how the expression should be normalized. Use the optional scale
79/// factor to rescale the expression after normalization.
80/// If shiftToZero is set, the entire curve is shifted down to make the lowest
81/// point of the curve go through zero.
82RooCurve::RooCurve(const RooAbsReal &f, RooAbsRealLValue &x, double xlo, double xhi, Int_t xbins,
83 double scaleFactor, const RooArgSet *normVars, double prec, double resolution,
84 bool shiftToZero, WingMode wmode, Int_t nEvalError, Int_t doEEVal, double eeVal,
85 bool showProg) :
86 TGraph(),
88 _showProgress(showProg)
89{
90
91 // grab the function's name and title
92 TString name(f.GetName());
93 SetName(name.Data());
94 TString title(f.GetTitle());
95 SetTitle(title.Data());
96 // append " ( [<funit> ][/ <xunit> ])" to our y-axis label if necessary
97 if(0 != strlen(f.getUnit()) || 0 != strlen(x.getUnit())) {
98 title.Append(" ( ");
99 if(0 != strlen(f.getUnit())) {
100 title.Append(f.getUnit());
101 title.Append(" ");
102 }
103 if(0 != strlen(x.getUnit())) {
104 title.Append("/ ");
105 title.Append(x.getUnit());
106 title.Append(" ");
107 }
108 title.Append(")");
109 }
110 setYAxisLabel(title.Data());
111
112 RooProduct scaledFunc{"scaled_func", "scaled_func", {f, RooFit::RooConst(scaleFactor)}};
113 std::unique_ptr<RooAbsFunc> funcPtr{scaledFunc.bindVars(x, normVars, true)};
114
115 // calculate the points to add to our curve
116 if(xbins > 0){
117 // regular mode - use the sampling hint to decide where to evaluate the pdf
118 std::unique_ptr<std::list<double>> hint{f.plotSamplingHint(x,xlo,xhi)};
119 addPoints(*funcPtr,xlo,xhi,xbins+1,prec,resolution,wmode,nEvalError,doEEVal,eeVal,hint.get());
120 if (_showProgress) {
121 ccoutP(Plotting) << endl ;
122 }
123 } else {
124 // if number of bins is set to <= 0, skip any interpolation and just evaluate the pdf at the bin centers
125 // this is useful when plotting a pdf like a histogram
126 int nBinsX = x.numBins();
127 for(int i=0; i<nBinsX; ++i){
128 double xval = x.getBinning().binCenter(i);
129 addPoint(xval,(*funcPtr)(&xval)) ;
130 }
131 }
132 initialize();
133
134 if (shiftToZero) shiftCurveToZero() ;
135
136 // Adjust limits
137 for (int i=0 ; i<GetN() ; i++) {
139 }
140 this->Sort();
141}
142
143
144
145////////////////////////////////////////////////////////////////////////////////
146/// Create a 1-dim curve of the value of the specified real-valued
147/// expression as a function of x. Use the optional precision
148/// parameter to control how precisely the smooth curve is
149/// rasterized. If shiftToZero is set, the entire curve is shifted
150/// down to make the lowest point in of the curve go through zero.
151
152RooCurve::RooCurve(const char *name, const char *title, const RooAbsFunc &func,
153 double xlo, double xhi, UInt_t minPoints, double prec, double resolution,
154 bool shiftToZero, WingMode wmode, Int_t nEvalError, Int_t doEEVal, double eeVal)
155{
156 SetName(name);
157 SetTitle(title);
158 addPoints(func,xlo,xhi,minPoints+1,prec,resolution,wmode,nEvalError,doEEVal,eeVal);
159 initialize();
160 if (shiftToZero) shiftCurveToZero() ;
161
162 // Adjust limits
163 for (int i=0 ; i<GetN() ; i++) {
165 }
166 this->Sort();
167}
168
169
170
171////////////////////////////////////////////////////////////////////////////////
172/// Constructor of a curve as sum of two other curves.
173/// \f[
174/// C_\mathrm{sum} = \mathrm{scale1}*c1 + \mathrm{scale2}*c2
175/// \f]
176///
177/// \param[in] name Name of the curve (to retrieve it from a plot)
178/// \param[in] title Title (for plotting).
179/// \param[in] c1 First curve.
180/// \param[in] c2 Second curve.
181/// \param[in] scale1 Scale y values for c1 by this factor.
182/// \param[in] scale2 Scale y values for c2 by this factor.
183
184RooCurve::RooCurve(const char* name, const char* title, const RooCurve& c1, const RooCurve& c2, double scale1, double scale2)
185{
186 initialize() ;
187 SetName(name) ;
188 SetTitle(title) ;
189
190 // Make deque of points in X
191 deque<double> pointList ;
192 double x,y ;
193
194 // Add X points of C1
195 Int_t i1,n1 = c1.GetN() ;
196 for (i1=0 ; i1<n1 ; i1++) {
197 c1.GetPoint(i1,x,y) ;
198 pointList.push_back(x) ;
199 }
200
201 // Add X points of C2
202 Int_t i2,n2 = c2.GetN() ;
203 for (i2=0 ; i2<n2 ; i2++) {
204 c2.GetPoint(i2,x,y) ;
205 pointList.push_back(x) ;
206 }
207
208 // Sort X points
209 sort(pointList.begin(),pointList.end()) ;
210
211 // Loop over X points
212 double last(-RooNumber::infinity()) ;
213 for (auto point : pointList) {
214
215 if ((point-last)>1e-10) {
216 // Add OR of points to new curve, skipping duplicate points within tolerance
217 addPoint(point,scale1*c1.interpolate(point)+scale2*c2.interpolate(point)) ;
218 }
219 last = point ;
220 }
221
222 this->Sort();
223}
224
225
226RooCurve::~RooCurve() = default;
227
228
229////////////////////////////////////////////////////////////////////////////////
230/// Perform initialization that is common to all curves
231
233{
234 // set default line width in pixels
235 SetLineWidth(3);
236 // set default line color
238}
239
240
241
242////////////////////////////////////////////////////////////////////////////////
243/// Find lowest point in curve and move all points in curve so that
244/// lowest point will go exactly through zero
245
247{
248 double minVal = std::numeric_limits<double>::infinity();
249 double maxVal = -std::numeric_limits<double>::infinity();
250
251 // First iteration, find current lowest point
252 for (int i = 1; i < GetN() - 1; i++) {
253 double x;
254 double y;
255 GetPoint(i, x, y);
256 minVal = std::min(y, minVal);
257 maxVal = std::max(y, maxVal);
258 }
259
260 // Second iteration, lower all points by minVal
261 for (int i = 1; i < GetN() - 1; i++) {
262 double x;
263 double y;
264 GetPoint(i, x, y);
265 SetPoint(i, x, y - minVal);
266 }
267
268 setYAxisLimits(0, maxVal - minVal);
269}
270
271
272
273////////////////////////////////////////////////////////////////////////////////
274/// Add points calculated with the specified function, over the range (xlo,xhi).
275/// Add at least minPoints equally spaced points, and add sufficient points so that
276/// the maximum deviation from the final straight-line segments is prec*(ymax-ymin),
277/// down to a minimum horizontal spacing of resolution*(xhi-xlo).
278
279void RooCurve::addPoints(const RooAbsFunc &func, double xlo, double xhi,
280 Int_t minPoints, double prec, double resolution, WingMode wmode,
281 Int_t numee, bool doEEVal, double eeVal, list<double>* samplingHint)
282{
283 // check the inputs
284 if(!func.isValid()) {
285 coutE(InputArguments) << fName << "::addPoints: input function is not valid" << endl;
286 return;
287 }
288 if(minPoints <= 0 || xhi <= xlo) {
289 coutE(InputArguments) << fName << "::addPoints: bad input (nothing added)" << endl;
290 return;
291 }
292
293 // Perform a coarse scan of the function to estimate its y range.
294 // Save the results so we do not have to re-evaluate at the scan points.
295
296 // Adjust minimum number of points to external sampling hint if used
297 if (samplingHint) {
298 minPoints = samplingHint->size() ;
299 }
300
301 double dx= (xhi-xlo)/(minPoints-1.);
302 const double epsilon = (xhi - xlo) * relativeXEpsilon();;
303 std::vector<double> yval(minPoints);
304
305 // Get list of initial x values. If function provides sampling hint use that,
306 // otherwise use default binning of frame
307 std::vector<double> xval;
308 if (!samplingHint) {
309 for(int step= 0; step < minPoints; step++) {
310 xval.push_back(xlo + step*dx) ;
311 }
312 } else {
313 std::copy(samplingHint->begin(), samplingHint->end(), std::back_inserter(xval));
314 }
315
316 for (unsigned int step=0; step < xval.size(); ++step) {
317 double xx = xval[step];
318 if (step == static_cast<unsigned int>(minPoints-1))
319 xx -= 1e-9 * dx;
320
321 yval[step]= func(&xx);
322 if (_showProgress) {
323 ccoutP(Plotting) << "." ;
324 cout.flush() ;
325 }
326
328 if (numee>=0) {
329 coutW(Plotting) << "At observable [x]=" << xx << " " ;
330 RooAbsReal::printEvalErrors(ccoutW(Plotting),numee) ;
331 }
332 if (doEEVal) {
333 yval[step]=eeVal ;
334 }
335 }
337 }
338
339 const double ymax = *std::max_element(yval.begin(), yval.end());
340 const double ymin = *std::min_element(yval.begin(), yval.end());
341 double yrangeEst=(ymax-ymin) ;
342
343 // store points of the coarse scan and calculate any refinements necessary
344 double minDx= resolution*(xhi-xlo);
345 double x1,x2= xlo;
346
347 if (wmode==Extended) {
348 // Add two points to make curve jump from 0 to yval at the left end of the plotting range.
349 // This ensures that filled polygons are drawn properly. The first point needs to be to the
350 // left of the second, so it's shifted by 1/1000 more than the second.
351 addPoint(xlo-dx*1.001, 0);
352 addPoint(xlo-dx,yval[0]) ;
353 } else if (wmode==Straight) {
354 addPoint(xlo-dx*0.001,0) ;
355 }
356
357 addPoint(xlo,yval[0]);
358
359 auto iter2 = xval.begin() ;
360 x1 = *iter2 ;
361 int step=1 ;
362 while(true) {
363 x1= x2;
364 ++iter2 ;
365 if (iter2==xval.end()) {
366 break ;
367 }
368 x2= *iter2 ;
369 if (prec<0) {
370 // If precision is <0, no attempt at recursive interpolation is made
371 addPoint(x2,yval[step]) ;
372 } else {
373 addRange(func,x1,x2,yval[step-1],yval[step],prec*yrangeEst,minDx,numee,doEEVal,eeVal,epsilon);
374 }
375 step++ ;
376 }
377 addPoint(xhi,yval[minPoints-1]) ;
378
379 if (wmode==Extended) {
380 // Add two points to close polygon. The order matters. Since they are sorted in x later, the second
381 // point is shifted by 1/1000 more than the second-to-last point.
382 addPoint(xhi+dx,yval[minPoints-1]) ;
383 addPoint(xhi+dx*1.001, 0);
384 } else if (wmode==Straight) {
385 addPoint(xhi+dx*0.001,0) ;
386 }
387}
388
389
390////////////////////////////////////////////////////////////////////////////////
391/// Fill the range (x1,x2) with points calculated using func(&x). No point will
392/// be added at x1, and a point will always be added at x2. The density of points
393/// will be calculated so that the maximum deviation from a straight line
394/// approximation is prec*(ymax-ymin) down to the specified minimum horizontal spacing.
395
396void RooCurve::addRange(const RooAbsFunc& func, double x1, double x2,
397 double y1, double y2, double minDy, double minDx,
398 int numee, bool doEEVal, double eeVal, double epsilon)
399{
400 // Explicitly skip empty ranges to eliminate point duplication
401 if (std::abs(x2-x1) <= epsilon) {
402 return ;
403 }
404
405 // calculate our value at the midpoint of this range
406 double xmid= 0.5*(x1+x2);
407 double ymid= func(&xmid);
408 if (_showProgress) {
409 ccoutP(Plotting) << "." ;
410 cout.flush() ;
411 }
412
414 if (numee>=0) {
415 coutW(Plotting) << "At observable [x]=" << xmid << " " ;
416 RooAbsReal::printEvalErrors(ccoutW(Plotting),numee) ;
417 }
418 if (doEEVal) {
419 ymid=eeVal ;
420 }
421 }
423
424 // test if the midpoint is sufficiently close to a straight line across this interval
425 double dy= ymid - 0.5*(y1+y2);
426 if((xmid - x1 >= minDx) && std::abs(dy)>0 && std::abs(dy) >= minDy) {
427 // fill in each subrange
428 addRange(func,x1,xmid,y1,ymid,minDy,minDx,numee,doEEVal,eeVal,epsilon);
429 addRange(func,xmid,x2,ymid,y2,minDy,minDx,numee,doEEVal,eeVal,epsilon);
430 }
431 else {
432 // add the endpoint
433 addPoint(x2,y2);
434 }
435}
436
437
438////////////////////////////////////////////////////////////////////////////////
439/// Add a point with the specified coordinates. Update our y-axis limits.
440
441void RooCurve::addPoint(double x, double y)
442{
443// cout << "RooCurve("<< GetName() << ") adding point at (" << x << "," << y << ")" << endl ;
444 Int_t next= GetN();
445 SetPoint(next, x, y);
447}
448
449
450////////////////////////////////////////////////////////////////////////////////
451/// Return the number of events associated with the plotable object,
452/// it is always 1 for curves
453
455 return 1;
456}
457
458
459////////////////////////////////////////////////////////////////////////////////
460/// Return the number of events associated with the plotable object,
461/// in the given range. It is always 1 for curves
462
463double RooCurve::getFitRangeNEvt(double, double) const
464{
465 return 1 ;
466}
467
468
469////////////////////////////////////////////////////////////////////////////////
470/// Get the bin width associated with this plotable object.
471/// It is alwats zero for curves
472
474 return 0 ;
475}
476
477
478
479////////////////////////////////////////////////////////////////////////////////
480
481void RooCurve::printName(ostream& os) const
482//
483{
484 // Print the name of this curve
485 os << GetName() ;
486}
487
488
489////////////////////////////////////////////////////////////////////////////////
490/// Print the title of this curve
491
492void RooCurve::printTitle(ostream& os) const
493{
494 os << GetTitle() ;
495}
496
497
498////////////////////////////////////////////////////////////////////////////////
499/// Print the class name of this curve
500
501void RooCurve::printClassName(ostream& os) const
502{
503 os << ClassName() ;
504}
505
506
507
508////////////////////////////////////////////////////////////////////////////////
509/// Print the details of this curve
510
511void RooCurve::printMultiline(ostream& os, Int_t /*contents*/, bool /*verbose*/, TString indent) const
512{
513 os << indent << "--- RooCurve ---" << endl ;
514 Int_t n= GetN();
515 os << indent << " Contains " << n << " points" << endl;
516 os << indent << " Graph points:" << endl;
517 for(Int_t i= 0; i < n; i++) {
518 os << indent << setw(3) << i << ") x = " << fX[i] << " , y = " << fY[i] << endl;
519 }
520}
521
522
523
524////////////////////////////////////////////////////////////////////////////////
525/// Calculate the chi^2/NDOF of this curve with respect to the histogram
526/// 'hist' accounting nFitParam floating parameters in case the curve
527/// was the result of a fit
528
529double RooCurve::chiSquare(const RooHist& hist, Int_t nFitParam) const
530{
531 Int_t i,np = hist.GetN() ;
532 double x,y,eyl,eyh,exl,exh ;
533
534 // Find starting and ending bin of histogram based on range of RooCurve
535 double xstart,xstop ;
536
537 GetPoint(0,xstart,y) ;
538 GetPoint(GetN()-1,xstop,y) ;
539
540 Int_t nbin(0) ;
541
543 for (i=0 ; i<np ; i++) {
544
545 // Retrieve histogram contents
546 hist.GetPoint(i,x,y) ;
547
548 // Check if point is in range of curve
549 if (x<xstart || x>xstop) continue ;
550
551 eyl = hist.GetEYlow()[i] ;
552 eyh = hist.GetEYhigh()[i] ;
553 exl = hist.GetEXlow()[i] ;
554 exh = hist.GetEXhigh()[i] ;
555
556 // Integrate function over this bin
557 double avg = average(x-exl,x+exh) ;
558
559 // Add pull^2 to chisq
560 if (y!=0) {
561 double pull = (y>avg) ? ((y-avg)/eyl) : ((y-avg)/eyh) ;
562 chisq += pull*pull ;
563 nbin++ ;
564 }
565 }
566
567 // Return chisq/nDOF
568 return chisq.Sum() / (nbin-nFitParam) ;
569}
570
571
572
573////////////////////////////////////////////////////////////////////////////////
574/// Return average curve value in [xFirst,xLast] by integrating curve between points
575/// and dividing by xLast-xFirst
576
577double RooCurve::average(double xFirst, double xLast) const
578{
579 if (xFirst>=xLast) {
580 coutE(InputArguments) << "RooCurve::average(" << GetName()
581 << ") invalid range (" << xFirst << "," << xLast << ")" << endl ;
582 return 0 ;
583 }
584
585 // Find Y values and begin and end points
586 double yFirst = interpolate(xFirst,1e-10) ;
587 double yLast = interpolate(xLast,1e-10) ;
588
589 // Find first and last mid points
590 Int_t ifirst = findPoint(xFirst,1e10) ;
591 Int_t ilast = findPoint(xLast,1e10) ;
592 double xFirstPt,yFirstPt,xLastPt,yLastPt ;
593 GetPoint(ifirst,xFirstPt,yFirstPt) ;
594 GetPoint(ilast,xLastPt,yLastPt) ;
595
596 double tolerance=1e-3*(xLast-xFirst) ;
597
598 // Handle trivial scenario -- no midway points, point only at or outside given range
599 if (ilast-ifirst==1 &&(xFirstPt-xFirst)<-1*tolerance && (xLastPt-xLast)>tolerance) {
600 return 0.5*(yFirst+yLast) ;
601 }
602
603 // If first point closest to xFirst is at xFirst or before xFirst take the next point
604 // as the first midway point
605 if ((xFirstPt-xFirst)<-1*tolerance) {
606 ifirst++ ;
607 const_cast<RooCurve&>(*this).GetPoint(ifirst,xFirstPt,yFirstPt) ;
608 }
609
610 // If last point closest to yLast is at yLast or beyond yLast the previous point
611 // as the last midway point
612 if ((xLastPt-xLast)>tolerance) {
613 ilast-- ;
614 const_cast<RooCurve&>(*this).GetPoint(ilast,xLastPt,yLastPt) ;
615 }
616
617 double sum(0),x1,y1,x2,y2 ;
618
619 // Trapezoid integration from lower edge to first midpoint
620 sum += (xFirstPt-xFirst)*(yFirst+yFirstPt)/2 ;
621
622 // Trapezoid integration between midpoints
623 Int_t i ;
624 for (i=ifirst ; i<ilast ; i++) {
625 const_cast<RooCurve&>(*this).GetPoint(i,x1,y1) ;
626 const_cast<RooCurve&>(*this).GetPoint(i+1,x2,y2) ;
627 sum += (x2-x1)*(y1+y2)/2 ;
628 }
629
630 // Trapezoid integration from last midpoint to upper edge
631 sum += (xLast-xLastPt)*(yLastPt+yLast)/2 ;
632 return sum/(xLast-xFirst) ;
633}
634
635
636
637////////////////////////////////////////////////////////////////////////////////
638/// Find the nearest point to xvalue. Return -1 if distance
639/// exceeds tolerance
640
641Int_t RooCurve::findPoint(double xvalue, double tolerance) const
642{
643 double delta(std::numeric_limits<double>::max()),x,y ;
644 Int_t i,n = GetN() ;
645 Int_t ibest(-1) ;
646 for (i=0 ; i<n ; i++) {
647 GetPoint(i,x,y);
648 if (std::abs(xvalue-x)<delta) {
649 delta = std::abs(xvalue-x) ;
650 ibest = i ;
651 }
652 }
653
654 return (delta<tolerance)?ibest:-1 ;
655}
656
657
658////////////////////////////////////////////////////////////////////////////////
659/// Return linearly interpolated value of curve at xvalue. If distance
660/// to nearest point is less than tolerance, return nearest point value
661/// instead
662
663double RooCurve::interpolate(double xvalue, double tolerance) const
664{
665 // Find best point
666 int n = GetN() ;
667 int ibest = findPoint(xvalue,1e10) ;
668
669 // Get position of best point
670 double xbest, ybest ;
671 const_cast<RooCurve*>(this)->GetPoint(ibest,xbest,ybest) ;
672
673 // Handle trivial case of being dead on
674 if (std::abs(xbest-xvalue)<tolerance) {
675 return ybest ;
676 }
677
678 // Get nearest point on other side w.r.t. xvalue
679 double xother,yother, retVal(0) ;
680 if (xbest<xvalue) {
681 if (ibest==n-1) {
682 // Value beyond end requested -- return value of last point
683 return ybest ;
684 }
685 const_cast<RooCurve*>(this)->GetPoint(ibest+1,xother,yother) ;
686 if (xother==xbest) return ybest ;
687 retVal = ybest + (yother-ybest)*(xvalue-xbest)/(xother-xbest) ;
688
689 } else {
690 if (ibest==0) {
691 // Value before 1st point requested -- return value of 1st point
692 return ybest ;
693 }
694 const_cast<RooCurve*>(this)->GetPoint(ibest-1,xother,yother) ;
695 if (xother==xbest) return ybest ;
696 retVal = yother + (ybest-yother)*(xvalue-xother)/(xbest-xother) ;
697 }
698
699 return retVal ;
700}
701
702
703
704
705////////////////////////////////////////////////////////////////////////////////
706/// Construct filled RooCurve represented error band that captures alpha% of the variations
707/// of the curves passed through argument variations, where the percentage alpha corresponds to
708/// the central interval fraction of a significance Z
709
710RooCurve* RooCurve::makeErrorBand(const vector<RooCurve*>& variations, double Z) const
711{
712 RooCurve* band = new RooCurve ;
713 band->SetName((std::string(GetName()) + "_errorband").c_str());
714 band->SetLineWidth(1) ;
715 band->SetFillColor(kCyan) ;
716 band->SetLineColor(kCyan) ;
717
718 vector<double> bandLo(GetN()) ;
719 vector<double> bandHi(GetN()) ;
720 for (int i=0 ; i<GetN() ; i++) {
721 calcBandInterval(variations,i,Z,bandLo[i],bandHi[i],false) ;
722 }
723
724 for (int i=0 ; i<GetN() ; i++) {
725 band->addPoint(GetX()[i],bandLo[i]) ;
726 }
727 for (int i=GetN()-1 ; i>=0 ; i--) {
728 band->addPoint(GetX()[i],bandHi[i]) ;
729 }
730 // if the axis of the old graph is alphanumeric, copy the labels to the new one as well
731 if(this->GetXaxis() && this->GetXaxis()->IsAlphanumeric()){
732 band->GetXaxis()->Set(this->GetXaxis()->GetNbins(),this->GetXaxis()->GetXmin(),this->GetXaxis()->GetXmax());
733 for(int i=0; i<this->GetXaxis()->GetNbins(); ++i){
734 band->GetXaxis()->SetBinLabel(i+1,this->GetXaxis()->GetBinLabel(i+1));
735 }
736 }
737
738 return band ;
739}
740
741
742
743
744////////////////////////////////////////////////////////////////////////////////
745/// Construct filled RooCurve represented error band represent the error added in quadrature defined by the curves arguments
746/// plusVar and minusVar corresponding to one-sigma variations of each parameter. The resulting error band, combined used the correlation matrix C
747/// is multiplied with the significance parameter Z to construct the equivalent of a Z sigma error band (in Gaussian approximation)
748
749RooCurve* RooCurve::makeErrorBand(const vector<RooCurve*>& plusVar, const vector<RooCurve*>& minusVar, const TMatrixD& C, double Z) const
750{
751
752 RooCurve* band = new RooCurve ;
753 band->SetName((std::string(GetName()) + "_errorband").c_str());
754 band->SetLineWidth(1) ;
755 band->SetFillColor(kCyan) ;
756 band->SetLineColor(kCyan) ;
757
758 vector<double> bandLo(GetN()) ;
759 vector<double> bandHi(GetN()) ;
760 for (int i=0 ; i<GetN() ; i++) {
761 calcBandInterval(plusVar,minusVar,i,C,Z,bandLo[i],bandHi[i]) ;
762 }
763
764 for (int i=0 ; i<GetN() ; i++) {
765 band->addPoint(GetX()[i],bandLo[i]) ;
766 }
767 for (int i=GetN()-1 ; i>=0 ; i--) {
768 band->addPoint(GetX()[i],bandHi[i]) ;
769 }
770
771 // if the axis of the old graph is alphanumeric, copy the labels to the new one as well
772 if(this->GetXaxis() && this->GetXaxis()->IsAlphanumeric()){
773 band->GetXaxis()->Set(this->GetXaxis()->GetNbins(),this->GetXaxis()->GetXmin(),this->GetXaxis()->GetXmax());
774 for(int i=0; i<this->GetXaxis()->GetNbins(); ++i){
775 band->GetXaxis()->SetBinLabel(i+1,this->GetXaxis()->GetBinLabel(i+1));
776 }
777 }
778
779 return band ;
780}
781
782
783
784
785
786////////////////////////////////////////////////////////////////////////////////
787/// Retrieve variation points from curves
788
789void RooCurve::calcBandInterval(const vector<RooCurve*>& plusVar, const vector<RooCurve*>& minusVar,Int_t i, const TMatrixD& C, double /*Z*/, double& lo, double& hi) const
790{
791 vector<double> y_plus(plusVar.size()), y_minus(minusVar.size()) ;
792 Int_t j(0) ;
793 for (vector<RooCurve*>::const_iterator iter=plusVar.begin() ; iter!=plusVar.end() ; ++iter) {
794 y_plus[j++] = (*iter)->interpolate(GetX()[i]) ;
795 }
796 j=0 ;
797 for (vector<RooCurve*>::const_iterator iter=minusVar.begin() ; iter!=minusVar.end() ; ++iter) {
798 y_minus[j++] = (*iter)->interpolate(GetX()[i]) ;
799 }
800 double y_cen = GetY()[i] ;
801 Int_t n = j ;
802
803 // Make vector of variations
804 TVectorD F(plusVar.size()) ;
805 for (j=0 ; j<n ; j++) {
806 F[j] = (y_plus[j]-y_minus[j])/2 ;
807 }
808
809 // Calculate error in linear approximation from variations and correlation coefficient
810 double sum = F*(C*F) ;
811
812 lo= y_cen + sqrt(sum) ;
813 hi= y_cen - sqrt(sum) ;
814}
815
816
817
818////////////////////////////////////////////////////////////////////////////////
819
820void RooCurve::calcBandInterval(const vector<RooCurve*>& variations,Int_t i,double Z, double& lo, double& hi, bool approxGauss) const
821{
822 vector<double> y(variations.size()) ;
823 Int_t j(0) ;
824 for (vector<RooCurve*>::const_iterator iter=variations.begin() ; iter!=variations.end() ; ++iter) {
825 y[j++] = (*iter)->interpolate(GetX()[i]) ;
826}
827
828 if (!approxGauss) {
829 // Construct central 68% interval from variations collected at each point
830 double pvalue = TMath::Erfc(Z/sqrt(2.)) ;
831 Int_t delta = Int_t( y.size()*(pvalue)/2 + 0.5) ;
832 sort(y.begin(),y.end()) ;
833 lo = y[delta] ;
834 hi = y[y.size()-delta] ;
835 } else {
836 // Estimate R.M.S of variations at each point and use that as Gaussian sigma
837 double sum_y(0), sum_ysq(0) ;
838 for (unsigned int k=0 ; k<y.size() ; k++) {
839 sum_y += y[k] ;
840 sum_ysq += y[k]*y[k] ;
841 }
842 sum_y /= y.size() ;
843 sum_ysq /= y.size() ;
844
845 double rms = sqrt(sum_ysq - (sum_y*sum_y)) ;
846 lo = GetY()[i] - Z*rms ;
847 hi = GetY()[i] + Z*rms ;
848 }
849}
850
851
852
853
854////////////////////////////////////////////////////////////////////////////////
855/// Return true if curve is identical to other curve allowing for given
856/// absolute tolerance on each point compared point.
857
858bool RooCurve::isIdentical(const RooCurve& other, double tol, bool verbose) const
859{
860 // Determine X range and Y range
861 Int_t n= min(GetN(),other.GetN());
862 double xmin(1e30), xmax(-1e30), ymin(1e30), ymax(-1e30) ;
863 for(Int_t i= 0; i < n; i++) {
864 if (fX[i]<xmin) xmin=fX[i] ;
865 if (fX[i]>xmax) xmax=fX[i] ;
866 if (fY[i]<ymin) ymin=fY[i] ;
867 if (fY[i]>ymax) ymax=fY[i] ;
868 }
869 const double Yrange=ymax-ymin ;
870
871 bool ret(true) ;
872 for(Int_t i= 2; i < n-2; i++) {
873 double yTest = interpolate(other.fX[i],1e-10) ;
874 double rdy = std::abs(yTest-other.fY[i])/Yrange ;
875 if (rdy>tol) {
876 ret = false;
877 if(!verbose) continue;
878 cout << "RooCurve::isIdentical[" << std::setw(3) << i << "] Y tolerance exceeded (" << std::setprecision(5) << std::setw(10) << rdy << ">" << tol << "),";
879 cout << " x,y=(" << std::right << std::setw(10) << fX[i] << "," << std::setw(10) << fY[i] << ")\tref: y="
880 << std::setw(10) << other.interpolate(fX[i], 1.E-15) << ". [Nearest point from ref: ";
881 auto j = other.findPoint(fX[i], 1.E10);
882 std::cout << "j=" << j << "\tx,y=(" << std::setw(10) << other.fX[j] << "," << std::setw(10) << other.fY[j] << ") ]" << "\trange=" << Yrange << std::endl;
883 }
884 }
885
886 return ret ;
887}
888
889
890
891////////////////////////////////////////////////////////////////////////////////
892/// Returns sampling hints for a histogram with given boundaries. This helper
893/// function is meant to be used by binned RooAbsReals to produce sampling
894/// hints that are working well with RooFits plotting.
895
896std::list<double> *
897RooCurve::plotSamplingHintForBinBoundaries(std::span<const double> boundaries, double xlo, double xhi)
898{
899 auto hint = new std::list<double>;
900
901 // Make sure the difference between two points around a bin boundary is
902 // larger than the relative epsilon for which the RooCurve considers two
903 // points as the same. Otherwise, the points right of the bin boundary would
904 // be skipped.
905 const double delta = (xhi - xlo) * RooCurve::relativeXEpsilon();
906
907 // Sample points right next to the plot limits
908 hint->push_back(xlo + delta);
909 hint->push_back(xhi - delta);
910
911 // Sample points very close to the left and right of the bin boundaries that
912 // are strictly in between the plot limits.
913 for (const double x : boundaries) {
914 if (x - xlo > delta && xhi - x > delta) {
915 hint->push_back(x - delta);
916 hint->push_back(x + delta);
917 }
918 }
919
920 hint->sort();
921
922 return hint;
923}
#define f(i)
Definition RSha256.hxx:104
#define e(i)
Definition RSha256.hxx:103
bool _showProgress
Show progress indication during evaluation if true.
#define coutW(a)
#define ccoutP(a)
#define coutE(a)
#define ccoutW(a)
int Int_t
Definition RtypesCore.h:45
#define ClassImp(name)
Definition Rtypes.h:377
@ kCyan
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
static void indent(ostringstream &buf, int indent_level)
Option_t Option_t SetLineWidth
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t np
Option_t Option_t SetLineColor
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char y1
char name[80]
Definition TGX11.cxx:110
float xmin
#define hi
float ymin
float xmax
float ymax
The Kahan summation is a compensated summation algorithm, which significantly reduces numerical error...
Definition Util.h:122
T Sum() const
Definition Util.h:240
Abstract interface for evaluating a real-valued function of one real variable and performing numerica...
Definition RooAbsFunc.h:27
bool isValid() const
Definition RooAbsFunc.h:37
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
RooFit::OwningPtr< RooAbsFunc > bindVars(const RooArgSet &vars, const RooArgSet *nset=nullptr, bool clipInvalid=false) const
Create an interface adaptor f(vars) that binds us to the specified variables (in arbitrary order).
static Int_t numEvalErrors()
Return the number of logged evaluation errors since the last clearing.
static void printEvalErrors(std::ostream &os=std::cout, Int_t maxPerNode=10000000)
Print all outstanding logged evaluation error on the given ostream.
static void clearEvalErrorLog()
Clear the stack of evaluation error messages.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
A RooCurve is a one-dimensional graphical representation of a real-valued function.
Definition RooCurve.h:36
void addPoints(const RooAbsFunc &func, double xlo, double xhi, Int_t minPoints, double prec, double resolution, WingMode wmode, Int_t numee=0, bool doEEVal=false, double eeVal=0.0, std::list< double > *samplingHint=nullptr)
Add points calculated with the specified function, over the range (xlo,xhi).
Definition RooCurve.cxx:279
void printTitle(std::ostream &os) const override
Print the title of this curve.
Definition RooCurve.cxx:492
void initialize()
Perform initialization that is common to all curves.
Definition RooCurve.cxx:232
~RooCurve() override
double getFitRangeBinW() const override
Get the bin width associated with this plotable object.
Definition RooCurve.cxx:473
static constexpr double relativeXEpsilon()
The distance between two points x1 and x2 relative to the full plot range below which two points are ...
Definition RooCurve.h:104
void addRange(const RooAbsFunc &func, double x1, double x2, double y1, double y2, double minDy, double minDx, int numee, bool doEEVal, double eeVal, double epsilon)
Fill the range (x1,x2) with points calculated using func(&x).
Definition RooCurve.cxx:396
void printName(std::ostream &os) const override
Print name of object.
Definition RooCurve.cxx:481
void printMultiline(std::ostream &os, Int_t contents, bool verbose=false, TString indent="") const override
Print the details of this curve.
Definition RooCurve.cxx:511
double interpolate(double x, double tolerance=1e-10) const
Return linearly interpolated value of curve at xvalue.
Definition RooCurve.cxx:663
static std::list< double > * plotSamplingHintForBinBoundaries(std::span< const double > boundaries, double xlo, double xhi)
Returns sampling hints for a histogram with given boundaries.
Definition RooCurve.cxx:897
void printClassName(std::ostream &os) const override
Print the class name of this curve.
Definition RooCurve.cxx:501
@ Extended
Definition RooCurve.h:39
@ Straight
Definition RooCurve.h:39
RooCurve()
Default constructor.
Definition RooCurve.cxx:68
bool _showProgress
! Show progress indication when adding points
Definition RooCurve.h:97
RooCurve * makeErrorBand(const std::vector< RooCurve * > &variations, double Z=1) const
Construct filled RooCurve represented error band that captures alpha% of the variations of the curves...
Definition RooCurve.cxx:710
double chiSquare(const RooHist &hist, int nFitParam) const
Calculate the chi^2/NDOF of this curve with respect to the histogram 'hist' accounting nFitParam floa...
Definition RooCurve.cxx:529
void calcBandInterval(const std::vector< RooCurve * > &variations, Int_t i, double Z, double &lo, double &hi, bool approxGauss) const
Definition RooCurve.cxx:820
double getFitRangeNEvt() const override
Return the number of events associated with the plotable object, it is always 1 for curves.
Definition RooCurve.cxx:454
void shiftCurveToZero()
Find lowest point in curve and move all points in curve so that lowest point will go exactly through ...
Definition RooCurve.cxx:246
void addPoint(double x, double y)
Add a point with the specified coordinates. Update our y-axis limits.
Definition RooCurve.cxx:441
bool isIdentical(const RooCurve &other, double tol=1e-6, bool verbose=true) const
Return true if curve is identical to other curve allowing for given absolute tolerance on each point ...
Definition RooCurve.cxx:858
Int_t findPoint(double value, double tolerance=1e-10) const
Find the nearest point to xvalue.
Definition RooCurve.cxx:641
double average(double lo, double hi) const
Return average curve value in [xFirst,xLast] by integrating curve between points and dividing by xLas...
Definition RooCurve.cxx:577
A RooHist is a graphical representation of binned data based on the TGraphAsymmErrors class.
Definition RooHist.h:29
static constexpr double infinity()
Return internal infinity representation.
Definition RooNumber.h:25
Class RooPotable is a base class for objects that can be inserted into RooPlots and take advantage of...
Definition RooPlotable.h:26
void updateYAxisLimits(double y)
Definition RooPlotable.h:30
void setYAxisLimits(double ymin, double ymax)
Definition RooPlotable.h:34
void setYAxisLabel(const char *label)
Definition RooPlotable.h:29
A RooProduct represents the product of a given set of RooAbsReal objects.
Definition RooProduct.h:29
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:37
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:43
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
virtual void SetBinLabel(Int_t bin, const char *label)
Set label for bin.
Definition TAxis.cxx:886
Bool_t IsAlphanumeric() const
Definition TAxis.h:88
Double_t GetXmax() const
Definition TAxis.h:140
const char * GetBinLabel(Int_t bin) const
Return label for bin.
Definition TAxis.cxx:440
virtual void Set(Int_t nbins, Double_t xmin, Double_t xmax)
Initialize axis with fix bins.
Definition TAxis.cxx:794
Double_t GetXmin() const
Definition TAxis.h:139
Int_t GetNbins() const
Definition TAxis.h:125
Double_t * GetEXlow() const override
Double_t * GetEYhigh() const override
Double_t * GetEXhigh() const override
Double_t * GetEYlow() const override
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Set x and y values for point number i.
Definition TGraph.cxx:2315
Double_t * GetY() const
Definition TGraph.h:138
Int_t GetN() const
Definition TGraph.h:130
Double_t * fY
[fNpoints] array of Y points
Definition TGraph.h:48
virtual void Sort(Bool_t(*greater)(const TGraph *, Int_t, Int_t)=&TGraph::CompareX, Bool_t ascending=kTRUE, Int_t low=0, Int_t high=-1111)
Sorts the points of this TGraph using in-place quicksort (see e.g.
Definition TGraph.cxx:2459
Double_t * GetX() const
Definition TGraph.h:137
void SetName(const char *name="") override
Set graph name.
Definition TGraph.cxx:2354
TAxis * GetXaxis() const
Get x axis of the graph.
Definition TGraph.cxx:1540
Double_t * fX
[fNpoints] array of X points
Definition TGraph.h:47
void SetTitle(const char *title="") override
Change (i.e.
Definition TGraph.cxx:2370
virtual Int_t GetPoint(Int_t i, Double_t &x, Double_t &y) const
Get x and y values for point number i.
Definition TGraph.cxx:1507
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
TString fName
Definition TNamed.h:32
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:207
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:380
TString & Append(const char *cs)
Definition TString.h:576
void addPoints()
Definition event_demo.C:69
RooConstVar & RooConst(double val)
Double_t y[n]
Definition legend1.C:17
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
return c2
Definition legend2.C:14
#define F(x, y, z)
Double_t Erfc(Double_t x)
Computes the complementary error function erfc(x).
Definition TMath.cxx:199
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345
double epsilon
Definition triangle.c:618