Logo ROOT   6.18/05
Reference Guide
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 endpoints chosen to give
24a "good" approximation to the true curve. The goodness of the approximation is
25controlled by a precision and a resolution parameter. To view the points where
26a function y(x) is actually evaluated to approximate a smooth curve, use the fact
27that a RooCurve is a TGraph:
28```
29RooPlot *p = y.plotOn(x.frame());
30p->getAttMarker("curve_y")->SetMarkerStyle(20);
31p->setDrawOptions("curve_y","PL");
32p->Draw();
33```
34
35To retrieve a RooCurve from a RooPlot, use RooPlot::getCurve().
36**/
37
38#include "RooFit.h"
39
40#include "RooCurve.h"
41#include "RooHist.h"
42#include "RooAbsReal.h"
43#include "RooArgSet.h"
44#include "RooRealVar.h"
45#include "RooRealIntegral.h"
46#include "RooRealBinding.h"
47#include "RooScaledFunc.h"
48#include "RooMsgService.h"
49
50#include "Riostream.h"
51#include "TClass.h"
52#include "TMath.h"
53#include "TMatrixD.h"
54#include "TVectorD.h"
55#include <iomanip>
56#include <deque>
57
58using namespace std ;
59
61
62
63////////////////////////////////////////////////////////////////////////////////
64/// Default constructor
65
66RooCurve::RooCurve() : _showProgress(kFALSE)
67{
68 initialize();
69}
70
71
72////////////////////////////////////////////////////////////////////////////////
73/// Create a 1-dim curve of the value of the specified real-valued expression
74/// as a function of x. Use the optional precision parameter to control
75/// how precisely the smooth curve is rasterized. Use the optional argument set
76/// to specify how the expression should be normalized. Use the optional scale
77/// factor to rescale the expression after normalization.
78/// If shiftToZero is set, the entire curve is shift down to make the lowest
79/// point in of the curve go through zero.
80
82 Double_t scaleFactor, const RooArgSet *normVars, Double_t prec, Double_t resolution,
83 Bool_t shiftToZero, WingMode wmode, Int_t nEvalError, Int_t doEEVal, Double_t eeVal,
84 Bool_t showProg) : _showProgress(showProg)
85{
86
87 // grab the function's name and title
88 TString name(f.GetName());
89 SetName(name.Data());
90 TString title(f.GetTitle());
91 SetTitle(title.Data());
92 // append " ( [<funit> ][/ <xunit> ])" to our y-axis label if necessary
93 if(0 != strlen(f.getUnit()) || 0 != strlen(x.getUnit())) {
94 title.Append(" ( ");
95 if(0 != strlen(f.getUnit())) {
96 title.Append(f.getUnit());
97 title.Append(" ");
98 }
99 if(0 != strlen(x.getUnit())) {
100 title.Append("/ ");
101 title.Append(x.getUnit());
102 title.Append(" ");
103 }
104 title.Append(")");
105 }
106 setYAxisLabel(title.Data());
107
108 RooAbsFunc *funcPtr = 0;
109 RooAbsFunc *rawPtr = 0;
110 funcPtr= f.bindVars(x,normVars,kTRUE);
111
112 // apply a scale factor if necessary
113 if(scaleFactor != 1) {
114 rawPtr= funcPtr;
115 funcPtr= new RooScaledFunc(*rawPtr,scaleFactor);
116 }
117 assert(0 != funcPtr);
118
119 // calculate the points to add to our curve
120 Double_t prevYMax = getYAxisMax() ;
121 list<Double_t>* hint = f.plotSamplingHint(x,xlo,xhi) ;
122 addPoints(*funcPtr,xlo,xhi,xbins+1,prec,resolution,wmode,nEvalError,doEEVal,eeVal,hint);
123 if (_showProgress) {
124 ccoutP(Plotting) << endl ;
125 }
126 if (hint) {
127 delete hint ;
128 }
129 initialize();
130
131 // cleanup
132 delete funcPtr;
133 if(rawPtr) delete rawPtr;
134 if (shiftToZero) shiftCurveToZero(prevYMax) ;
135
136 // Adjust limits
137 Int_t i ;
138 for (i=0 ; i<GetN() ; i++) {
139 Double_t x2,y2 ;
140 GetPoint(i,x2,y2) ;
142 }
143}
144
145
146
147////////////////////////////////////////////////////////////////////////////////
148/// Create a 1-dim curve of the value of the specified real-valued
149/// expression as a function of x. Use the optional precision
150/// parameter to control how precisely the smooth curve is
151/// rasterized. If shiftToZero is set, the entire curve is shift
152/// down to make the lowest point in of the curve go through zero.
153
154RooCurve::RooCurve(const char *name, const char *title, const RooAbsFunc &func,
155 Double_t xlo, Double_t xhi, UInt_t minPoints, Double_t prec, Double_t resolution,
156 Bool_t shiftToZero, WingMode wmode, Int_t nEvalError, Int_t doEEVal, Double_t eeVal) :
157 _showProgress(kFALSE)
158{
159 SetName(name);
160 SetTitle(title);
161 Double_t prevYMax = getYAxisMax() ;
162 addPoints(func,xlo,xhi,minPoints+1,prec,resolution,wmode,nEvalError,doEEVal,eeVal);
163 initialize();
164 if (shiftToZero) shiftCurveToZero(prevYMax) ;
165
166 // Adjust limits
167 Int_t i ;
168 for (i=0 ; i<GetN() ; i++) {
169 Double_t x,y ;
170 GetPoint(i,x,y) ;
172 }
173}
174
175
176
177////////////////////////////////////////////////////////////////////////////////
178/// Constructor of curve as sum of two other curves
179///
180/// Csum = scale1*c1 + scale2*c2
181///
182
183RooCurve::RooCurve(const char* name, const char* title, const RooCurve& c1, const RooCurve& c2, Double_t scale1, Double_t scale2) :
184 _showProgress(kFALSE)
185{
186 initialize() ;
187 SetName(name) ;
188 SetTitle(title) ;
189
190 // Make deque of points in X
191 deque<Double_t> pointList ;
192 Double_t x,y ;
193
194 // Add X points of C1
195 Int_t i1,n1 = c1.GetN() ;
196 for (i1=0 ; i1<n1 ; i1++) {
197 const_cast<RooCurve&>(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 const_cast<RooCurve&>(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 deque<double>::iterator iter ;
214 for (iter=pointList.begin() ; iter!=pointList.end() ; ++iter) {
215
216 if ((*iter-last)>1e-10) {
217 // Add OR of points to new curve, skipping duplicate points within tolerance
218 addPoint(*iter,scale1*c1.interpolate(*iter)+scale2*c2.interpolate(*iter)) ;
219 }
220 last = *iter ;
221 }
222
223}
224
225
226
227////////////////////////////////////////////////////////////////////////////////
228/// Destructor
229
231{
232}
233
234
235
236////////////////////////////////////////////////////////////////////////////////
237/// Perform initialization that is common to all curves
238
240{
241 // set default line width in pixels
242 SetLineWidth(3);
243 // set default line color
245}
246
247
248
249////////////////////////////////////////////////////////////////////////////////
250/// Find lowest point in curve and move all points in curve so that
251/// lowest point will go exactly through zero
252
254{
255 Int_t i ;
256 Double_t minVal(1e30) ;
257 Double_t maxVal(-1e30) ;
258
259 // First iteration, find current lowest point
260 for (i=1 ; i<GetN()-1 ; i++) {
261 Double_t x,y ;
262 GetPoint(i,x,y) ;
263 if (y<minVal) minVal=y ;
264 if (y>maxVal) maxVal=y ;
265 }
266
267 // Second iteration, lower all points by minVal
268 for (i=1 ; i<GetN()-1 ; i++) {
269 Double_t x,y ;
270 GetPoint(i,x,y) ;
271 SetPoint(i,x,y-minVal) ;
272 }
273
274 // Check if y-axis range needs readjustment
275 if (getYAxisMax()>prevYMax) {
276 Double_t newMax = maxVal - minVal ;
277 setYAxisLimits(getYAxisMin(), newMax<prevYMax ? prevYMax : newMax) ;
278 }
279}
280
281
282
283////////////////////////////////////////////////////////////////////////////////
284/// Add points calculated with the specified function, over the range (xlo,xhi).
285/// Add at least minPoints equally spaced points, and add sufficient points so that
286/// the maximum deviation from the final straight-line segements is prec*(ymax-ymin),
287/// down to a minimum horizontal spacing of resolution*(xhi-xlo).
288
290 Int_t minPoints, Double_t prec, Double_t resolution, WingMode wmode,
291 Int_t numee, Bool_t doEEVal, Double_t eeVal, list<Double_t>* samplingHint)
292{
293 // check the inputs
294 if(!func.isValid()) {
295 coutE(InputArguments) << fName << "::addPoints: input function is not valid" << endl;
296 return;
297 }
298 if(minPoints <= 0 || xhi <= xlo) {
299 coutE(InputArguments) << fName << "::addPoints: bad input (nothing added)" << endl;
300 return;
301 }
302
303 // Perform a coarse scan of the function to estimate its y range.
304 // Save the results so we do not have to re-evaluate at the scan points.
305
306 // Adjust minimum number of points to external sampling hint if used
307 if (samplingHint) {
308 minPoints = samplingHint->size() ;
309 }
310
311 Int_t step;
312 Double_t dx= (xhi-xlo)/(minPoints-1.);
313 Double_t *yval= new Double_t[minPoints];
314
315 // Get list of initial x values. If function provides sampling hint use that,
316 // otherwise use default binning of frame
317 list<Double_t>* xval = samplingHint ;
318 if (!xval) {
319 xval = new list<Double_t> ;
320 for(step= 0; step < minPoints; step++) {
321 xval->push_back(xlo + step*dx) ;
322 }
323 }
324
325
326 Double_t ymax(-1e30), ymin(1e30) ;
327
328 step=0 ;
329 for(list<Double_t>::iterator iter = xval->begin() ; iter!=xval->end() ; ++iter,++step) {
330 Double_t xx = *iter ;
331
332 if (step==minPoints-1) xx-=1e-15 ;
333
334 yval[step]= func(&xx);
335 if (_showProgress) {
336 ccoutP(Plotting) << "." ;
337 cout.flush() ;
338 }
339
341 if (numee>=0) {
342 coutW(Plotting) << "At observable [x]=" << xx << " " ;
344 }
345 if (doEEVal) {
346 yval[step]=eeVal ;
347 }
348 }
350
351
352 if (yval[step]>ymax) ymax=yval[step] ;
353 if (yval[step]<ymin) ymin=yval[step] ;
354 }
355 Double_t yrangeEst=(ymax-ymin) ;
356
357 // store points of the coarse scan and calculate any refinements necessary
358 Double_t minDx= resolution*(xhi-xlo);
359 Double_t x1,x2= xlo;
360
361 if (wmode==Extended) {
362 addPoint(xlo-dx,0) ;
363 addPoint(xlo-dx,yval[0]) ;
364 } else if (wmode==Straight) {
365 addPoint(xlo,0) ;
366 }
367
368 addPoint(xlo,yval[0]);
369
370 list<Double_t>::iterator iter2 = xval->begin() ;
371 x1 = *iter2 ;
372 step=1 ;
373 while(true) {
374 x1= x2;
375 ++iter2 ;
376 if (iter2==xval->end()) {
377 break ;
378 }
379 x2= *iter2 ;
380 if (prec<0) {
381 // If precision is <0, no attempt at recursive interpolation is made
382 addPoint(x2,yval[step]) ;
383 } else {
384 addRange(func,x1,x2,yval[step-1],yval[step],prec*yrangeEst,minDx,numee,doEEVal,eeVal);
385 }
386 step++ ;
387 }
388 addPoint(xhi,yval[minPoints-1]) ;
389
390 if (wmode==Extended) {
391 addPoint(xhi+dx,yval[minPoints-1]) ;
392 addPoint(xhi+dx,0) ;
393 } else if (wmode==Straight) {
394 addPoint(xhi,0) ;
395 }
396
397 // cleanup
398 delete [] yval;
399 if (xval != samplingHint) {
400 delete xval ;
401 }
402
403}
404
405
406////////////////////////////////////////////////////////////////////////////////
407/// Fill the range (x1,x2) with points calculated using func(&x). No point will
408/// be added at x1, and a point will always be added at x2. The density of points
409/// will be calculated so that the maximum deviation from a straight line
410/// approximation is prec*(ymax-ymin) down to the specified minimum horizontal spacing.
411
413 Double_t y1, Double_t y2, Double_t minDy, Double_t minDx,
414 Int_t numee, Bool_t doEEVal, Double_t eeVal)
415{
416 // Explicitly skip empty ranges to eliminate point duplication
417 if (fabs(x2-x1)<1e-20) {
418 return ;
419 }
420
421 // calculate our value at the midpoint of this range
422 Double_t xmid= 0.5*(x1+x2);
423 Double_t ymid= func(&xmid);
424 if (_showProgress) {
425 ccoutP(Plotting) << "." ;
426 cout.flush() ;
427 }
428
430 if (numee>=0) {
431 coutW(Plotting) << "At observable [x]=" << xmid << " " ;
433 }
434 if (doEEVal) {
435 ymid=eeVal ;
436 }
437 }
439
440 // test if the midpoint is sufficiently close to a straight line across this interval
441 Double_t dy= ymid - 0.5*(y1+y2);
442 if((xmid - x1 >= minDx) && fabs(dy)>0 && fabs(dy) >= minDy) {
443 // fill in each subrange
444 addRange(func,x1,xmid,y1,ymid,minDy,minDx,numee,doEEVal,eeVal);
445 addRange(func,xmid,x2,ymid,y2,minDy,minDx,numee,doEEVal,eeVal);
446 }
447 else {
448 // add the endpoint
449 addPoint(x2,y2);
450 }
451}
452
453
454////////////////////////////////////////////////////////////////////////////////
455/// Add a point with the specified coordinates. Update our y-axis limits.
456
458{
459// cout << "RooCurve("<< GetName() << ") adding point at (" << x << "," << y << ")" << endl ;
460 Int_t next= GetN();
461 SetPoint(next, x, y);
463}
464
465
466////////////////////////////////////////////////////////////////////////////////
467/// Return the number of events associated with the plotable object,
468/// it is always 1 for curves
469
471 return 1;
472}
473
474
475////////////////////////////////////////////////////////////////////////////////
476/// Return the number of events associated with the plotable object,
477/// in the given range. It is always 1 for curves
478
480{
481 return 1 ;
482}
483
484
485////////////////////////////////////////////////////////////////////////////////
486/// Get the bin width associated with this plotable object.
487/// It is alwats zero for curves
488
490 return 0 ;
491}
492
493
494
495////////////////////////////////////////////////////////////////////////////////
496
497void RooCurve::printName(ostream& os) const
498//
499{
500 // Print the name of this curve
501 os << GetName() ;
502}
503
504
505////////////////////////////////////////////////////////////////////////////////
506/// Print the title of this curve
507
508void RooCurve::printTitle(ostream& os) const
509{
510 os << GetTitle() ;
511}
512
513
514////////////////////////////////////////////////////////////////////////////////
515/// Print the class name of this curve
516
517void RooCurve::printClassName(ostream& os) const
518{
519 os << IsA()->GetName() ;
520}
521
522
523
524////////////////////////////////////////////////////////////////////////////////
525/// Print the details of this curve
526
527void RooCurve::printMultiline(ostream& os, Int_t /*contents*/, Bool_t /*verbose*/, TString indent) const
528{
529 os << indent << "--- RooCurve ---" << endl ;
530 Int_t n= GetN();
531 os << indent << " Contains " << n << " points" << endl;
532 os << indent << " Graph points:" << endl;
533 for(Int_t i= 0; i < n; i++) {
534 os << indent << setw(3) << i << ") x = " << fX[i] << " , y = " << fY[i] << endl;
535 }
536}
537
538
539
540////////////////////////////////////////////////////////////////////////////////
541/// Calculate the chi^2/NDOF of this curve with respect to the histogram
542/// 'hist' accounting nFitParam floating parameters in case the curve
543/// was the result of a fit
544
545Double_t RooCurve::chiSquare(const RooHist& hist, Int_t nFitParam) const
546{
547 Int_t i,np = hist.GetN() ;
548 Double_t x,y,eyl,eyh,exl,exh ;
549
550 // Find starting and ending bin of histogram based on range of RooCurve
551 Double_t xstart,xstop ;
552
553#if ROOT_VERSION_CODE >= ROOT_VERSION(4,0,1)
554 GetPoint(0,xstart,y) ;
555 GetPoint(GetN()-1,xstop,y) ;
556#else
557 const_cast<RooCurve*>(this)->GetPoint(0,xstart,y) ;
558 const_cast<RooCurve*>(this)->GetPoint(GetN()-1,xstop,y) ;
559#endif
560
561 Int_t nbin(0) ;
562
563 Double_t chisq(0) ;
564 for (i=0 ; i<np ; i++) {
565
566 // Retrieve histogram contents
567 ((RooHist&)hist).GetPoint(i,x,y) ;
568
569 // Check if point is in range of curve
570 if (x<xstart || x>xstop) continue ;
571
572 eyl = hist.GetEYlow()[i] ;
573 eyh = hist.GetEYhigh()[i] ;
574 exl = hist.GetEXlow()[i] ;
575 exh = hist.GetEXhigh()[i] ;
576
577 // Integrate function over this bin
578 Double_t avg = average(x-exl,x+exh) ;
579
580 // Add pull^2 to chisq
581 if (y!=0) {
582 Double_t pull = (y>avg) ? ((y-avg)/eyl) : ((y-avg)/eyh) ;
583 chisq += pull*pull ;
584 nbin++ ;
585 }
586 }
587
588 // Return chisq/nDOF
589 return chisq / (nbin-nFitParam) ;
590}
591
592
593
594////////////////////////////////////////////////////////////////////////////////
595/// Return average curve value in [xFirst,xLast] by integrating curve between points
596/// and dividing by xLast-xFirst
597
599{
600 if (xFirst>=xLast) {
601 coutE(InputArguments) << "RooCurve::average(" << GetName()
602 << ") invalid range (" << xFirst << "," << xLast << ")" << endl ;
603 return 0 ;
604 }
605
606 // Find Y values and begin and end points
607 Double_t yFirst = interpolate(xFirst,1e-10) ;
608 Double_t yLast = interpolate(xLast,1e-10) ;
609
610 // Find first and last mid points
611 Int_t ifirst = findPoint(xFirst,1e10) ;
612 Int_t ilast = findPoint(xLast,1e10) ;
613 Double_t xFirstPt,yFirstPt,xLastPt,yLastPt ;
614 const_cast<RooCurve&>(*this).GetPoint(ifirst,xFirstPt,yFirstPt) ;
615 const_cast<RooCurve&>(*this).GetPoint(ilast,xLastPt,yLastPt) ;
616
617 Double_t tolerance=1e-3*(xLast-xFirst) ;
618
619 // Handle trivial scenario -- no midway points, point only at or outside given range
620 if (ilast-ifirst==1 &&(xFirstPt-xFirst)<-1*tolerance && (xLastPt-xLast)>tolerance) {
621 return 0.5*(yFirst+yLast) ;
622 }
623
624 // If first point closest to xFirst is at xFirst or before xFirst take the next point
625 // as the first midway point
626 if ((xFirstPt-xFirst)<-1*tolerance) {
627 ifirst++ ;
628 const_cast<RooCurve&>(*this).GetPoint(ifirst,xFirstPt,yFirstPt) ;
629 }
630
631 // If last point closest to yLast is at yLast or beyond yLast the the previous point
632 // as the last midway point
633 if ((xLastPt-xLast)>tolerance) {
634 ilast-- ;
635 const_cast<RooCurve&>(*this).GetPoint(ilast,xLastPt,yLastPt) ;
636 }
637
638 Double_t sum(0),x1,y1,x2,y2 ;
639
640 // Trapezoid integration from lower edge to first midpoint
641 sum += (xFirstPt-xFirst)*(yFirst+yFirstPt)/2 ;
642
643 // Trapezoid integration between midpoints
644 Int_t i ;
645 for (i=ifirst ; i<ilast ; i++) {
646 const_cast<RooCurve&>(*this).GetPoint(i,x1,y1) ;
647 const_cast<RooCurve&>(*this).GetPoint(i+1,x2,y2) ;
648 sum += (x2-x1)*(y1+y2)/2 ;
649 }
650
651 // Trapezoid integration from last midpoint to upper edge
652 sum += (xLast-xLastPt)*(yLastPt+yLast)/2 ;
653 return sum/(xLast-xFirst) ;
654}
655
656
657
658////////////////////////////////////////////////////////////////////////////////
659/// Find the nearest point to xvalue. Return -1 if distance
660/// exceeds tolerance
661
663{
664 Double_t delta(std::numeric_limits<double>::max()),x,y ;
665 Int_t i,n = GetN() ;
666 Int_t ibest(-1) ;
667 for (i=0 ; i<n ; i++) {
668 ((RooCurve&)*this).GetPoint(i,x,y) ;
669 if (fabs(xvalue-x)<delta) {
670 delta = fabs(xvalue-x) ;
671 ibest = i ;
672 }
673 }
674
675 return (delta<tolerance)?ibest:-1 ;
676}
677
678
679////////////////////////////////////////////////////////////////////////////////
680/// Return linearly interpolated value of curve at xvalue. If distance
681/// to nearest point is less than tolerance, return nearest point value
682/// instead
683
685{
686 // Find best point
687 int n = GetN() ;
688 int ibest = findPoint(xvalue,1e10) ;
689
690 // Get position of best point
691 Double_t xbest, ybest ;
692 const_cast<RooCurve*>(this)->GetPoint(ibest,xbest,ybest) ;
693
694 // Handle trivial case of being dead on
695 if (fabs(xbest-xvalue)<tolerance) {
696 return ybest ;
697 }
698
699 // Get nearest point on other side w.r.t. xvalue
700 Double_t xother,yother, retVal(0) ;
701 if (xbest<xvalue) {
702 if (ibest==n-1) {
703 // Value beyond end requested -- return value of last point
704 return ybest ;
705 }
706 const_cast<RooCurve*>(this)->GetPoint(ibest+1,xother,yother) ;
707 if (xother==xbest) return ybest ;
708 retVal = ybest + (yother-ybest)*(xvalue-xbest)/(xother-xbest) ;
709
710 } else {
711 if (ibest==0) {
712 // Value before 1st point requested -- return value of 1st point
713 return ybest ;
714 }
715 const_cast<RooCurve*>(this)->GetPoint(ibest-1,xother,yother) ;
716 if (xother==xbest) return ybest ;
717 retVal = yother + (ybest-yother)*(xvalue-xother)/(xbest-xother) ;
718 }
719
720 return retVal ;
721}
722
723
724
725
726////////////////////////////////////////////////////////////////////////////////
727/// Construct filled RooCurve represented error band that captures alpha% of the variations
728/// of the curves passed through argument variations, where the percentage alpha corresponds to
729/// the central interval fraction of a significance Z
730
731RooCurve* RooCurve::makeErrorBand(const vector<RooCurve*>& variations, Double_t Z) const
732{
733 RooCurve* band = new RooCurve ;
734 band->SetName(Form("%s_errorband",GetName())) ;
735 band->SetLineWidth(1) ;
736 band->SetFillColor(kCyan) ;
737 band->SetLineColor(kCyan) ;
738
739 vector<double> bandLo(GetN()) ;
740 vector<double> bandHi(GetN()) ;
741 for (int i=0 ; i<GetN() ; i++) {
742 calcBandInterval(variations,i,Z,bandLo[i],bandHi[i],kFALSE) ;
743 }
744
745 for (int i=0 ; i<GetN() ; i++) {
746 band->addPoint(GetX()[i],bandLo[i]) ;
747 }
748 for (int i=GetN()-1 ; i>=0 ; i--) {
749 band->addPoint(GetX()[i],bandHi[i]) ;
750 }
751
752 return band ;
753}
754
755
756
757
758////////////////////////////////////////////////////////////////////////////////
759/// Construct filled RooCurve represented error band represent the error added in quadrature defined by the curves arguments
760/// plusVar and minusVar corresponding to one-sigma variations of each parameter. The resulting error band, combined used the correlation matrix C
761/// is multiplied with the significance parameter Z to construct the equivalent of a Z sigma error band (in Gaussian approximation)
762
763RooCurve* RooCurve::makeErrorBand(const vector<RooCurve*>& plusVar, const vector<RooCurve*>& minusVar, const TMatrixD& C, Double_t Z) const
764{
765 RooCurve* band = new RooCurve ;
766 band->SetName(Form("%s_errorband",GetName())) ;
767 band->SetLineWidth(1) ;
768 band->SetFillColor(kCyan) ;
769 band->SetLineColor(kCyan) ;
770
771 vector<double> bandLo(GetN()) ;
772 vector<double> bandHi(GetN()) ;
773 for (int i=0 ; i<GetN() ; i++) {
774 calcBandInterval(plusVar,minusVar,i,C,Z,bandLo[i],bandHi[i]) ;
775 }
776
777 for (int i=0 ; i<GetN() ; i++) {
778 band->addPoint(GetX()[i],bandLo[i]) ;
779 }
780 for (int i=GetN()-1 ; i>=0 ; i--) {
781 band->addPoint(GetX()[i],bandHi[i]) ;
782 }
783
784 return band ;
785}
786
787
788
789
790
791////////////////////////////////////////////////////////////////////////////////
792/// Retrieve variation points from curves
793
794void RooCurve::calcBandInterval(const vector<RooCurve*>& plusVar, const vector<RooCurve*>& minusVar,Int_t i, const TMatrixD& C, Double_t /*Z*/, Double_t& lo, Double_t& hi) const
795{
796 vector<double> y_plus(plusVar.size()), y_minus(minusVar.size()) ;
797 Int_t j(0) ;
798 for (vector<RooCurve*>::const_iterator iter=plusVar.begin() ; iter!=plusVar.end() ; ++iter) {
799 y_plus[j++] = (*iter)->interpolate(GetX()[i]) ;
800 }
801 j=0 ;
802 for (vector<RooCurve*>::const_iterator iter=minusVar.begin() ; iter!=minusVar.end() ; ++iter) {
803 y_minus[j++] = (*iter)->interpolate(GetX()[i]) ;
804 }
805 Double_t y_cen = GetY()[i] ;
806 Int_t n = j ;
807
808 // Make vector of variations
809 TVectorD F(plusVar.size()) ;
810 for (j=0 ; j<n ; j++) {
811 F[j] = (y_plus[j]-y_minus[j])/2 ;
812 }
813
814 // Calculate error in linear approximation from variations and correlation coefficient
815 Double_t sum = F*(C*F) ;
816
817 lo= y_cen + sqrt(sum) ;
818 hi= y_cen - sqrt(sum) ;
819}
820
821
822
823////////////////////////////////////////////////////////////////////////////////
824
825void RooCurve::calcBandInterval(const vector<RooCurve*>& variations,Int_t i,Double_t Z, Double_t& lo, Double_t& hi, Bool_t approxGauss) const
826{
827 vector<double> y(variations.size()) ;
828 Int_t j(0) ;
829 for (vector<RooCurve*>::const_iterator iter=variations.begin() ; iter!=variations.end() ; ++iter) {
830 y[j++] = (*iter)->interpolate(GetX()[i]) ;
831}
832
833 if (!approxGauss) {
834 // Construct central 68% interval from variations collected at each point
835 Double_t pvalue = TMath::Erfc(Z/sqrt(2.)) ;
836 Int_t delta = Int_t( y.size()*(pvalue)/2 + 0.5) ;
837 sort(y.begin(),y.end()) ;
838 lo = y[delta] ;
839 hi = y[y.size()-delta] ;
840 } else {
841 // Estimate R.M.S of variations at each point and use that as Gaussian sigma
842 Double_t sum_y(0), sum_ysq(0) ;
843 for (unsigned int k=0 ; k<y.size() ; k++) {
844 sum_y += y[k] ;
845 sum_ysq += y[k]*y[k] ;
846 }
847 sum_y /= y.size() ;
848 sum_ysq /= y.size() ;
849
850 Double_t rms = sqrt(sum_ysq - (sum_y*sum_y)) ;
851 lo = GetY()[i] - Z*rms ;
852 hi = GetY()[i] + Z*rms ;
853 }
854}
855
856
857
858
859////////////////////////////////////////////////////////////////////////////////
860/// Return true if curve is identical to other curve allowing for given
861/// absolute tolerance on each point compared point.
862
864{
865 // Determine X range and Y range
866 Int_t n= min(GetN(),other.GetN());
867 Double_t xmin(1e30), xmax(-1e30), ymin(1e30), ymax(-1e30) ;
868 for(Int_t i= 0; i < n; i++) {
869 if (fX[i]<xmin) xmin=fX[i] ;
870 if (fX[i]>xmax) xmax=fX[i] ;
871 if (fY[i]<ymin) ymin=fY[i] ;
872 if (fY[i]>ymax) ymax=fY[i] ;
873 }
874 Double_t Yrange=ymax-ymin ;
875
876 Bool_t ret(kTRUE) ;
877 for(Int_t i= 2; i < n-2; i++) {
878 Double_t yTest = interpolate(other.fX[i],1e-10) ;
879 Double_t rdy = fabs(yTest-other.fY[i])/Yrange ;
880 if (rdy>tol) {
881
882// cout << "xref = " << other.fX[i] << " yref = " << other.fY[i] << " xtest = " << fX[i] << " ytest = " << fY[i]
883// << " ytestInt[other.fX] = " << interpolate(other.fX[i],1e-10) << endl ;
884
885 cout << "RooCurve::isIdentical[" << i << "] Y tolerance exceeded (" << rdy << ">" << tol
886 << "), X=" << other.fX[i] << "(" << fX[i] << ")" << " Ytest=" << yTest << " Yref=" << other.fY[i] << " range = " << Yrange << endl ;
887 ret=kFALSE ;
888 }
889 }
890
891 return ret ;
892}
893
894
#define f(i)
Definition: RSha256.hxx:104
#define e(i)
Definition: RSha256.hxx:103
static const double x2[5]
static const double x1[5]
#define coutW(a)
Definition: RooMsgService.h:33
#define ccoutP(a)
Definition: RooMsgService.h:39
#define coutE(a)
Definition: RooMsgService.h:34
#define ccoutW(a)
Definition: RooMsgService.h:40
int Int_t
Definition: RtypesCore.h:41
unsigned int UInt_t
Definition: RtypesCore.h:42
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
const Bool_t kTRUE
Definition: RtypesCore.h:87
#define ClassImp(name)
Definition: Rtypes.h:365
@ kCyan
Definition: Rtypes.h:64
@ kBlue
Definition: Rtypes.h:64
static void indent(ostringstream &buf, int indent_level)
char name[80]
Definition: TGX11.cxx:109
float xmin
Definition: THbookFile.cxx:93
float ymin
Definition: THbookFile.cxx:93
float xmax
Definition: THbookFile.cxx:93
float type_of_call hi(const int &, const int &)
float ymax
Definition: THbookFile.cxx:93
double sqrt(double)
char * Form(const char *fmt,...)
Abstract interface for evaluating a real-valued function of one real variable and performing numerica...
Definition: RooAbsFunc.h:23
Bool_t isValid() const
Definition: RooAbsFunc.h:33
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
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:28
A RooCurve is a one-dimensional graphical representation of a real-valued function.
Definition: RooCurve.h:32
void addPoints(const RooAbsFunc &func, Double_t xlo, Double_t xhi, Int_t minPoints, Double_t prec, Double_t resolution, WingMode wmode, Int_t numee=0, Bool_t doEEVal=kFALSE, Double_t eeVal=0., std::list< Double_t > *samplingHint=0)
Add points calculated with the specified function, over the range (xlo,xhi).
Definition: RooCurve.cxx:289
void initialize()
Perform initialization that is common to all curves.
Definition: RooCurve.cxx:239
void addPoint(Double_t x, Double_t y)
Add a point with the specified coordinates. Update our y-axis limits.
Definition: RooCurve.cxx:457
Double_t 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:545
Double_t average(Double_t lo, Double_t hi) const
Return average curve value in [xFirst,xLast] by integrating curve between points and dividing by xLas...
Definition: RooCurve.cxx:598
void addRange(const RooAbsFunc &func, Double_t x1, Double_t x2, Double_t y1, Double_t y2, Double_t minDy, Double_t minDx, Int_t numee=0, Bool_t doEEVal=kFALSE, Double_t eeVal=0.)
Fill the range (x1,x2) with points calculated using func(&x).
Definition: RooCurve.cxx:412
Bool_t _showProgress
Definition: RooCurve.h:91
void shiftCurveToZero(Double_t prevYMax)
Find lowest point in curve and move all points in curve so that lowest point will go exactly through ...
Definition: RooCurve.cxx:253
virtual ~RooCurve()
Destructor.
Definition: RooCurve.cxx:230
RooCurve * makeErrorBand(const std::vector< RooCurve * > &variations, Double_t Z=1) const
Construct filled RooCurve represented error band that captures alpha% of the variations of the curves...
Definition: RooCurve.cxx:731
virtual void printName(std::ostream &os) const
Print name of object.
Definition: RooCurve.cxx:497
@ Extended
Definition: RooCurve.h:35
@ Straight
Definition: RooCurve.h:35
RooCurve()
Default constructor.
Definition: RooCurve.cxx:66
Bool_t isIdentical(const RooCurve &other, Double_t tol=1e-6) const
Return true if curve is identical to other curve allowing for given absolute tolerance on each point ...
Definition: RooCurve.cxx:863
Double_t getFitRangeNEvt() const
Return the number of events associated with the plotable object, it is always 1 for curves.
Definition: RooCurve.cxx:470
virtual void printClassName(std::ostream &os) const
Print the class name of this curve.
Definition: RooCurve.cxx:517
virtual void printMultiline(std::ostream &os, Int_t contents, Bool_t verbose=kFALSE, TString indent="") const
Print the details of this curve.
Definition: RooCurve.cxx:527
Double_t interpolate(Double_t x, Double_t tolerance=1e-10) const
Return linearly interpolated value of curve at xvalue.
Definition: RooCurve.cxx:684
virtual void printTitle(std::ostream &os) const
Print the title of this curve.
Definition: RooCurve.cxx:508
Int_t findPoint(Double_t value, Double_t tolerance=1e-10) const
Find the nearest point to xvalue.
Definition: RooCurve.cxx:662
Double_t getFitRangeBinW() const
Get the bin width associated with this plotable object.
Definition: RooCurve.cxx:489
void calcBandInterval(const std::vector< RooCurve * > &variations, Int_t i, Double_t Z, Double_t &lo, Double_t &hi, Bool_t approxGauss) const
Definition: RooCurve.cxx:825
A RooHist is a graphical representation of binned data based on the TGraphAsymmErrors class.
Definition: RooHist.h:26
static Double_t infinity()
Return internal infinity representation.
Definition: RooNumber.cxx:49
Double_t getYAxisMin() const
Definition: RooPlotable.h:41
void updateYAxisLimits(Double_t y)
Definition: RooPlotable.h:33
Double_t getYAxisMax() const
Definition: RooPlotable.h:42
void setYAxisLimits(Double_t ymin, Double_t ymax)
Definition: RooPlotable.h:37
void setYAxisLabel(const char *label)
Definition: RooPlotable.h:32
Lightweight RooAbsFunction implementation that applies a constant scale factor to another RooAbsFunc.
Definition: RooScaledFunc.h:21
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
Double_t * GetEXhigh() const
Double_t * GetEYlow() const
Double_t * GetEXlow() const
Double_t * GetEYhigh() const
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Set x and y values for point number i.
Definition: TGraph.cxx:2198
Double_t * GetY() const
Definition: TGraph.h:131
virtual void SetName(const char *name="")
Set graph name.
Definition: TGraph.cxx:2221
Int_t GetN() const
Definition: TGraph.h:123
virtual void SetTitle(const char *title="")
Change (i.e.
Definition: TGraph.cxx:2237
Double_t * fY
[fNpoints] array of Y points
Definition: TGraph.h:48
Double_t * GetX() const
Definition: TGraph.h:130
Double_t * fX
[fNpoints] array of X points
Definition: TGraph.h:47
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:1580
TString fName
Definition: TNamed.h:32
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
Basic string class.
Definition: TString.h:131
const char * Data() const
Definition: TString.h:364
TString & Append(const char *cs)
Definition: TString.h:559
return c1
Definition: legend1.C:41
Double_t y[n]
Definition: legend1.C:17
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)
std::string GetName(const std::string &scope_name)
Definition: Cppyy.cxx:146
static double C[]
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
@ InputArguments
Definition: RooGlobalFunc.h:58
Double_t Erfc(Double_t x)
Compute the complementary error function erfc(x).
Definition: TMath.cxx:194
static long int sum(long int i)
Definition: Factory.cxx:2258