00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034 #include "RooFit.h"
00035
00036 #include "RooCurve.h"
00037 #include "RooHist.h"
00038 #include "RooAbsReal.h"
00039 #include "RooArgSet.h"
00040 #include "RooRealVar.h"
00041 #include "RooRealIntegral.h"
00042 #include "RooRealBinding.h"
00043 #include "RooScaledFunc.h"
00044 #include "RooMsgService.h"
00045
00046 #include "Riostream.h"
00047 #include "TClass.h"
00048 #include "TMath.h"
00049 #include "TMatrixD.h"
00050 #include "TVectorD.h"
00051 #include <iomanip>
00052 #include <math.h>
00053 #include <assert.h>
00054 #include <deque>
00055 #include <algorithm>
00056
00057 using namespace std ;
00058
00059 ClassImp(RooCurve)
00060
00061
00062
00063 RooCurve::RooCurve() : _showProgress(kFALSE)
00064 {
00065
00066 initialize();
00067 }
00068
00069
00070
00071 RooCurve::RooCurve(const RooAbsReal &f, RooAbsRealLValue &x, Double_t xlo, Double_t xhi, Int_t xbins,
00072 Double_t scaleFactor, const RooArgSet *normVars, Double_t prec, Double_t resolution,
00073 Bool_t shiftToZero, WingMode wmode, Int_t nEvalError, Int_t doEEVal, Double_t eeVal,
00074 Bool_t showProg) : _showProgress(showProg)
00075 {
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086 TString name(f.GetName());
00087 SetName(name.Data());
00088 TString title(f.GetTitle());
00089 SetTitle(title.Data());
00090
00091 if(0 != strlen(f.getUnit()) || 0 != strlen(x.getUnit())) {
00092 title.Append(" ( ");
00093 if(0 != strlen(f.getUnit())) {
00094 title.Append(f.getUnit());
00095 title.Append(" ");
00096 }
00097 if(0 != strlen(x.getUnit())) {
00098 title.Append("/ ");
00099 title.Append(x.getUnit());
00100 title.Append(" ");
00101 }
00102 title.Append(")");
00103 }
00104 setYAxisLabel(title.Data());
00105
00106 RooAbsFunc *funcPtr = 0;
00107 RooAbsFunc *rawPtr = 0;
00108 funcPtr= f.bindVars(x,normVars,kTRUE);
00109
00110
00111 if(scaleFactor != 1) {
00112 rawPtr= funcPtr;
00113 funcPtr= new RooScaledFunc(*rawPtr,scaleFactor);
00114 }
00115 assert(0 != funcPtr);
00116
00117
00118 Double_t prevYMax = getYAxisMax() ;
00119 list<Double_t>* hint = f.plotSamplingHint(x,xlo,xhi) ;
00120 addPoints(*funcPtr,xlo,xhi,xbins+1,prec,resolution,wmode,nEvalError,doEEVal,eeVal,hint);
00121 if (_showProgress) {
00122 ccoutP(Plotting) << endl ;
00123 }
00124 if (hint) {
00125 delete hint ;
00126 }
00127 initialize();
00128
00129
00130 delete funcPtr;
00131 if(rawPtr) delete rawPtr;
00132 if (shiftToZero) shiftCurveToZero(prevYMax) ;
00133
00134
00135 Int_t i ;
00136 for (i=0 ; i<GetN() ; i++) {
00137 Double_t x2,y2 ;
00138 GetPoint(i,x2,y2) ;
00139 updateYAxisLimits(y2);
00140 }
00141 }
00142
00143
00144
00145
00146 RooCurve::RooCurve(const char *name, const char *title, const RooAbsFunc &func,
00147 Double_t xlo, Double_t xhi, UInt_t minPoints, Double_t prec, Double_t resolution,
00148 Bool_t shiftToZero, WingMode wmode, Int_t nEvalError, Int_t doEEVal, Double_t eeVal) :
00149 _showProgress(kFALSE)
00150 {
00151
00152
00153
00154
00155
00156
00157 SetName(name);
00158 SetTitle(title);
00159 Double_t prevYMax = getYAxisMax() ;
00160 addPoints(func,xlo,xhi,minPoints+1,prec,resolution,wmode,nEvalError,doEEVal,eeVal);
00161 initialize();
00162 if (shiftToZero) shiftCurveToZero(prevYMax) ;
00163
00164
00165 Int_t i ;
00166 for (i=0 ; i<GetN() ; i++) {
00167 Double_t x,y ;
00168 GetPoint(i,x,y) ;
00169 updateYAxisLimits(y);
00170 }
00171 }
00172
00173
00174
00175
00176 RooCurve::RooCurve(const char* name, const char* title, const RooCurve& c1, const RooCurve& c2, Double_t scale1, Double_t scale2) :
00177 _showProgress(kFALSE)
00178 {
00179
00180
00181
00182
00183
00184 initialize() ;
00185 SetName(name) ;
00186 SetTitle(title) ;
00187
00188
00189 deque<Double_t> pointList ;
00190 Double_t x,y ;
00191
00192
00193 Int_t i1,n1 = c1.GetN() ;
00194 for (i1=0 ; i1<n1 ; i1++) {
00195 const_cast<RooCurve&>(c1).GetPoint(i1,x,y) ;
00196 pointList.push_back(x) ;
00197 }
00198
00199
00200 Int_t i2,n2 = c2.GetN() ;
00201 for (i2=0 ; i2<n2 ; i2++) {
00202 const_cast<RooCurve&>(c2).GetPoint(i2,x,y) ;
00203 pointList.push_back(x) ;
00204 }
00205
00206
00207 sort(pointList.begin(),pointList.end()) ;
00208
00209
00210 deque<double>::iterator iter ;
00211 Double_t last(-RooNumber::infinity()) ;
00212 for (iter=pointList.begin() ; iter!=pointList.end() ; ++iter) {
00213
00214 if ((*iter-last)>1e-10) {
00215
00216 addPoint(*iter,scale1*c1.interpolate(*iter)+scale2*c2.interpolate(*iter)) ;
00217 }
00218 last = *iter ;
00219 }
00220
00221 }
00222
00223
00224
00225
00226 RooCurve::~RooCurve()
00227 {
00228
00229 }
00230
00231
00232
00233
00234 void RooCurve::initialize()
00235 {
00236
00237
00238
00239 SetLineWidth(3);
00240
00241 SetLineColor(kBlue);
00242 }
00243
00244
00245
00246
00247 void RooCurve::shiftCurveToZero(Double_t prevYMax)
00248 {
00249
00250
00251
00252 Int_t i ;
00253 Double_t minVal(1e30) ;
00254 Double_t maxVal(-1e30) ;
00255
00256
00257 for (i=1 ; i<GetN()-1 ; i++) {
00258 Double_t x,y ;
00259 GetPoint(i,x,y) ;
00260 if (y<minVal) minVal=y ;
00261 if (y>maxVal) maxVal=y ;
00262 }
00263
00264
00265 for (i=1 ; i<GetN()-1 ; i++) {
00266 Double_t x,y ;
00267 GetPoint(i,x,y) ;
00268 SetPoint(i,x,y-minVal) ;
00269 }
00270
00271
00272 if (getYAxisMax()>prevYMax) {
00273 Double_t newMax = maxVal - minVal ;
00274 setYAxisLimits(getYAxisMin(), newMax<prevYMax ? prevYMax : newMax) ;
00275 }
00276 }
00277
00278
00279
00280
00281 void RooCurve::addPoints(const RooAbsFunc &func, Double_t xlo, Double_t xhi,
00282 Int_t minPoints, Double_t prec, Double_t resolution, WingMode wmode,
00283 Int_t numee, Bool_t doEEVal, Double_t eeVal, list<Double_t>* samplingHint)
00284 {
00285
00286
00287
00288
00289
00290
00291 if(!func.isValid()) {
00292 coutE(InputArguments) << fName << "::addPoints: input function is not valid" << endl;
00293 return;
00294 }
00295 if(minPoints <= 0 || xhi <= xlo) {
00296 coutE(InputArguments) << fName << "::addPoints: bad input (nothing added)" << endl;
00297 return;
00298 }
00299
00300
00301
00302
00303
00304 if (samplingHint) {
00305 minPoints = samplingHint->size() ;
00306 }
00307
00308 Int_t step;
00309 Double_t dx= (xhi-xlo)/(minPoints-1.);
00310 Double_t *yval= new Double_t[minPoints];
00311
00312
00313
00314 list<Double_t>* xval = samplingHint ;
00315 if (!xval) {
00316 xval = new list<Double_t> ;
00317 for(step= 0; step < minPoints; step++) {
00318 xval->push_back(xlo + step*dx) ;
00319 }
00320 }
00321
00322
00323 Double_t ymax(-1e30), ymin(1e30) ;
00324
00325 step=0 ;
00326 for(list<Double_t>::iterator iter = xval->begin() ; iter!=xval->end() ; ++iter,++step) {
00327 Double_t xx = *iter ;
00328
00329 if (step==minPoints-1) xx-=1e-15 ;
00330
00331 yval[step]= func(&xx);
00332 if (_showProgress) {
00333 ccoutP(Plotting) << "." ;
00334 cout.flush() ;
00335 }
00336
00337 if (RooAbsReal::numEvalErrors()>0) {
00338 if (numee>=0) {
00339 coutW(Plotting) << "At observable [x]=" << xx << " " ;
00340 RooAbsReal::printEvalErrors(ccoutW(Plotting),numee) ;
00341 }
00342 if (doEEVal) {
00343 yval[step]=eeVal ;
00344 }
00345 }
00346 RooAbsReal::clearEvalErrorLog() ;
00347
00348
00349 if (yval[step]>ymax) ymax=yval[step] ;
00350 if (yval[step]<ymin) ymin=yval[step] ;
00351 }
00352 Double_t yrangeEst=(ymax-ymin) ;
00353
00354
00355 Double_t minDx= resolution*(xhi-xlo);
00356 Double_t x1,x2= xlo;
00357
00358 if (wmode==Extended) {
00359 addPoint(xlo-dx,0) ;
00360 addPoint(xlo-dx,yval[0]) ;
00361 } else if (wmode==Straight) {
00362 addPoint(xlo,0) ;
00363 }
00364
00365 addPoint(xlo,yval[0]);
00366
00367 list<Double_t>::iterator iter2 = xval->begin() ;
00368 x1 = *iter2 ;
00369 step=1 ;
00370 while(true) {
00371 x1= x2;
00372 iter2++ ;
00373 if (iter2==xval->end()) {
00374 break ;
00375 }
00376 x2= *iter2 ;
00377 addRange(func,x1,x2,yval[step-1],yval[step],prec*yrangeEst,minDx,numee,doEEVal,eeVal);
00378 step++ ;
00379 }
00380 addPoint(xhi,yval[minPoints-1]) ;
00381
00382 if (wmode==Extended) {
00383 addPoint(xhi+dx,yval[minPoints-1]) ;
00384 addPoint(xhi+dx,0) ;
00385 } else if (wmode==Straight) {
00386 addPoint(xhi,0) ;
00387 }
00388
00389
00390 delete [] yval;
00391 if (xval != samplingHint) {
00392 delete xval ;
00393 }
00394
00395 }
00396
00397
00398
00399 void RooCurve::addRange(const RooAbsFunc& func, Double_t x1, Double_t x2,
00400 Double_t y1, Double_t y2, Double_t minDy, Double_t minDx,
00401 Int_t numee, Bool_t doEEVal, Double_t eeVal)
00402 {
00403
00404
00405
00406
00407
00408
00409 if (fabs(x2-x1)<1e-20) {
00410 return ;
00411 }
00412
00413
00414 Double_t xmid= 0.5*(x1+x2);
00415 Double_t ymid= func(&xmid);
00416 if (_showProgress) {
00417 ccoutP(Plotting) << "." ;
00418 cout.flush() ;
00419 }
00420
00421 if (RooAbsReal::numEvalErrors()>0) {
00422 if (numee>=0) {
00423 coutW(Plotting) << "At observable [x]=" << xmid << " " ;
00424 RooAbsReal::printEvalErrors(ccoutW(Plotting),numee) ;
00425 }
00426 if (doEEVal) {
00427 ymid=eeVal ;
00428 }
00429 }
00430 RooAbsReal::clearEvalErrorLog() ;
00431
00432
00433 Double_t dy= ymid - 0.5*(y1+y2);
00434 if((xmid - x1 >= minDx) && fabs(dy)>0 && fabs(dy) >= minDy) {
00435
00436 addRange(func,x1,xmid,y1,ymid,minDy,minDx,numee,doEEVal,eeVal);
00437 addRange(func,xmid,x2,ymid,y2,minDy,minDx,numee,doEEVal,eeVal);
00438 }
00439 else {
00440
00441 addPoint(x2,y2);
00442 }
00443 }
00444
00445
00446
00447 void RooCurve::addPoint(Double_t x, Double_t y)
00448 {
00449
00450
00451
00452 Int_t next= GetN();
00453 SetPoint(next, x, y);
00454 updateYAxisLimits(y) ;
00455 }
00456
00457
00458
00459 Double_t RooCurve::getFitRangeNEvt() const {
00460
00461
00462 return 1;
00463 }
00464
00465
00466
00467 Double_t RooCurve::getFitRangeNEvt(Double_t, Double_t) const
00468 {
00469
00470
00471 return 1 ;
00472 }
00473
00474
00475
00476 Double_t RooCurve::getFitRangeBinW() const {
00477
00478
00479 return 0 ;
00480 }
00481
00482
00483
00484
00485 void RooCurve::printName(ostream& os) const
00486
00487 {
00488
00489 os << GetName() ;
00490 }
00491
00492
00493
00494 void RooCurve::printTitle(ostream& os) const
00495 {
00496
00497 os << GetTitle() ;
00498 }
00499
00500
00501
00502 void RooCurve::printClassName(ostream& os) const
00503 {
00504
00505 os << IsA()->GetName() ;
00506 }
00507
00508
00509
00510
00511 void RooCurve::printMultiline(ostream& os, Int_t , Bool_t , TString indent) const
00512 {
00513
00514 os << indent << "--- RooCurve ---" << endl ;
00515 Int_t n= GetN();
00516 os << indent << " Contains " << n << " points" << endl;
00517 os << indent << " Graph points:" << endl;
00518 for(Int_t i= 0; i < n; i++) {
00519 os << indent << setw(3) << i << ") x = " << fX[i] << " , y = " << fY[i] << endl;
00520 }
00521 }
00522
00523
00524
00525
00526 Double_t RooCurve::chiSquare(const RooHist& hist, Int_t nFitParam) const
00527 {
00528
00529
00530
00531
00532 Int_t i,np = hist.GetN() ;
00533 Double_t x,y,eyl,eyh,exl,exh ;
00534
00535
00536 Double_t xstart,xstop ;
00537
00538 #if ROOT_VERSION_CODE >= ROOT_VERSION(4,0,1)
00539 GetPoint(0,xstart,y) ;
00540 GetPoint(GetN()-1,xstop,y) ;
00541 #else
00542 const_cast<RooCurve*>(this)->GetPoint(0,xstart,y) ;
00543 const_cast<RooCurve*>(this)->GetPoint(GetN()-1,xstop,y) ;
00544 #endif
00545
00546 Int_t nbin(0) ;
00547
00548 Double_t chisq(0) ;
00549 for (i=0 ; i<np ; i++) {
00550
00551
00552 ((RooHist&)hist).GetPoint(i,x,y) ;
00553
00554
00555 if (x<xstart || x>xstop) continue ;
00556
00557 nbin++ ;
00558 eyl = hist.GetEYlow()[i] ;
00559 eyh = hist.GetEYhigh()[i] ;
00560 exl = hist.GetEXlow()[i] ;
00561 exh = hist.GetEXhigh()[i] ;
00562
00563
00564 Double_t avg = average(x-exl,x+exh) ;
00565
00566
00567 if (y!=0) {
00568 Double_t pull = (y>avg) ? ((y-avg)/eyl) : ((y-avg)/eyh) ;
00569 chisq += pull*pull ;
00570 }
00571 }
00572
00573
00574 return chisq / (nbin-nFitParam) ;
00575 }
00576
00577
00578
00579
00580 Double_t RooCurve::average(Double_t xFirst, Double_t xLast) const
00581 {
00582
00583
00584
00585 if (xFirst>=xLast) {
00586 coutE(InputArguments) << "RooCurve::average(" << GetName()
00587 << ") invalid range (" << xFirst << "," << xLast << ")" << endl ;
00588 return 0 ;
00589 }
00590
00591
00592 Double_t yFirst = interpolate(xFirst,1e-10) ;
00593 Double_t yLast = interpolate(xLast,1e-10) ;
00594
00595
00596 Int_t ifirst = findPoint(xFirst,1e10) ;
00597 Int_t ilast = findPoint(xLast,1e10) ;
00598 Double_t xFirstPt,yFirstPt,xLastPt,yLastPt ;
00599 const_cast<RooCurve&>(*this).GetPoint(ifirst,xFirstPt,yFirstPt) ;
00600 const_cast<RooCurve&>(*this).GetPoint(ilast,xLastPt,yLastPt) ;
00601
00602 Double_t tolerance=1e-3*(xLast-xFirst) ;
00603
00604
00605 if (ilast-ifirst==1 &&(xFirstPt-xFirst)<-1*tolerance && (xLastPt-xLast)>tolerance) {
00606 return 0.5*(yFirst+yLast) ;
00607 }
00608
00609
00610
00611 if ((xFirstPt-xFirst)<-1*tolerance) {
00612 ifirst++ ;
00613 const_cast<RooCurve&>(*this).GetPoint(ifirst,xFirstPt,yFirstPt) ;
00614 }
00615
00616
00617
00618 if ((xLastPt-xLast)>tolerance) {
00619 ilast-- ;
00620 const_cast<RooCurve&>(*this).GetPoint(ilast,xLastPt,yLastPt) ;
00621 }
00622
00623 Double_t sum(0),x1,y1,x2,y2 ;
00624
00625
00626 sum += (xFirstPt-xFirst)*(yFirst+yFirstPt)/2 ;
00627
00628
00629 Int_t i ;
00630 for (i=ifirst ; i<ilast ; i++) {
00631 const_cast<RooCurve&>(*this).GetPoint(i,x1,y1) ;
00632 const_cast<RooCurve&>(*this).GetPoint(i+1,x2,y2) ;
00633 sum += (x2-x1)*(y1+y2)/2 ;
00634 }
00635
00636
00637 sum += (xLast-xLastPt)*(yLastPt+yLast)/2 ;
00638 return sum/(xLast-xFirst) ;
00639 }
00640
00641
00642
00643
00644 Int_t RooCurve::findPoint(Double_t xvalue, Double_t tolerance) const
00645 {
00646
00647
00648
00649 Double_t delta(999.),x,y ;
00650 Int_t i,n = GetN() ;
00651 Int_t ibest(-1) ;
00652 for (i=0 ; i<n ; i++) {
00653 ((RooCurve&)*this).GetPoint(i,x,y) ;
00654 if (fabs(xvalue-x)<delta) {
00655 delta = fabs(xvalue-x) ;
00656 ibest = i ;
00657 }
00658 }
00659
00660 return (delta<tolerance)?ibest:-1 ;
00661 }
00662
00663
00664
00665 Double_t RooCurve::interpolate(Double_t xvalue, Double_t tolerance) const
00666 {
00667
00668
00669
00670
00671
00672 int n = GetN() ;
00673 int ibest = findPoint(xvalue,1e10) ;
00674
00675
00676 Double_t xbest, ybest ;
00677 const_cast<RooCurve*>(this)->GetPoint(ibest,xbest,ybest) ;
00678
00679
00680 if (fabs(xbest-xvalue)<tolerance) {
00681 return ybest ;
00682 }
00683
00684
00685 Double_t xother,yother, retVal(0) ;
00686 if (xbest<xvalue) {
00687 if (ibest==n-1) {
00688
00689 return ybest ;
00690 }
00691 const_cast<RooCurve*>(this)->GetPoint(ibest+1,xother,yother) ;
00692 if (xother==xbest) return ybest ;
00693 retVal = ybest + (yother-ybest)*(xvalue-xbest)/(xother-xbest) ;
00694
00695 } else {
00696 if (ibest==0) {
00697
00698 return ybest ;
00699 }
00700 const_cast<RooCurve*>(this)->GetPoint(ibest-1,xother,yother) ;
00701 if (xother==xbest) return ybest ;
00702 retVal = yother + (ybest-yother)*(xvalue-xother)/(xbest-xother) ;
00703 }
00704
00705 return retVal ;
00706 }
00707
00708
00709
00710
00711
00712 RooCurve* RooCurve::makeErrorBand(const vector<RooCurve*>& variations, Double_t Z) const
00713 {
00714
00715
00716
00717
00718 RooCurve* band = new RooCurve ;
00719 band->SetName(Form("%s_errorband",GetName())) ;
00720 band->SetLineWidth(1) ;
00721 band->SetFillColor(kCyan) ;
00722 band->SetLineColor(kCyan) ;
00723
00724 vector<double> bandLo(GetN()) ;
00725 vector<double> bandHi(GetN()) ;
00726 for (int i=0 ; i<GetN() ; i++) {
00727 calcBandInterval(variations,i,Z,bandLo[i],bandHi[i],kFALSE) ;
00728 }
00729
00730 for (int i=0 ; i<GetN() ; i++) {
00731 band->addPoint(GetX()[i],bandLo[i]) ;
00732 }
00733 for (int i=GetN()-1 ; i>=0 ; i--) {
00734 band->addPoint(GetX()[i],bandHi[i]) ;
00735 }
00736
00737 return band ;
00738 }
00739
00740
00741
00742
00743
00744 RooCurve* RooCurve::makeErrorBand(const vector<RooCurve*>& plusVar, const vector<RooCurve*>& minusVar, const TMatrixD& C, Double_t Z) const
00745 {
00746
00747
00748
00749
00750 RooCurve* band = new RooCurve ;
00751 band->SetName(Form("%s_errorband",GetName())) ;
00752 band->SetLineWidth(1) ;
00753 band->SetFillColor(kCyan) ;
00754 band->SetLineColor(kCyan) ;
00755
00756 vector<double> bandLo(GetN()) ;
00757 vector<double> bandHi(GetN()) ;
00758 for (int i=0 ; i<GetN() ; i++) {
00759 calcBandInterval(plusVar,minusVar,i,C,Z,bandLo[i],bandHi[i]) ;
00760 }
00761
00762 for (int i=0 ; i<GetN() ; i++) {
00763 band->addPoint(GetX()[i],bandLo[i]) ;
00764 }
00765 for (int i=GetN()-1 ; i>=0 ; i--) {
00766 band->addPoint(GetX()[i],bandHi[i]) ;
00767 }
00768
00769 return band ;
00770 }
00771
00772
00773
00774
00775
00776
00777 void RooCurve::calcBandInterval(const vector<RooCurve*>& plusVar, const vector<RooCurve*>& minusVar,Int_t i, const TMatrixD& C, Double_t , Double_t& lo, Double_t& hi) const
00778 {
00779
00780 vector<double> y_plus(plusVar.size()), y_minus(minusVar.size()) ;
00781 Int_t j(0) ;
00782 for (vector<RooCurve*>::const_iterator iter=plusVar.begin() ; iter!=plusVar.end() ; iter++) {
00783 y_plus[j++] = (*iter)->interpolate(GetX()[i]) ;
00784 }
00785 j=0 ;
00786 for (vector<RooCurve*>::const_iterator iter=minusVar.begin() ; iter!=minusVar.end() ; iter++) {
00787 y_minus[j++] = (*iter)->interpolate(GetX()[i]) ;
00788 }
00789 Double_t y_cen = GetY()[i] ;
00790 Int_t n = j ;
00791
00792
00793 TVectorD F(plusVar.size()) ;
00794 for (j=0 ; j<n ; j++) {
00795 F[j] = (y_plus[j]-y_minus[j])/2 ;
00796 }
00797
00798
00799 Double_t sum = F*(C*F) ;
00800
00801 lo= y_cen + sqrt(sum) ;
00802 hi= y_cen - sqrt(sum) ;
00803 }
00804
00805
00806
00807
00808 void RooCurve::calcBandInterval(const vector<RooCurve*>& variations,Int_t i,Double_t Z, Double_t& lo, Double_t& hi, Bool_t approxGauss) const
00809 {
00810 vector<double> y(variations.size()) ;
00811 Int_t j(0) ;
00812 for (vector<RooCurve*>::const_iterator iter=variations.begin() ; iter!=variations.end() ; iter++) {
00813 y[j++] = (*iter)->interpolate(GetX()[i]) ;
00814 }
00815
00816 if (!approxGauss) {
00817
00818 Double_t pvalue = TMath::Erfc(Z/sqrt(2.)) ;
00819 Int_t delta = Int_t( y.size()*(pvalue)/2 + 0.5) ;
00820 sort(y.begin(),y.end()) ;
00821 lo = y[delta] ;
00822 hi = y[y.size()-delta] ;
00823 } else {
00824
00825 Double_t sum_y(0), sum_ysq(0) ;
00826 for (unsigned int k=0 ; k<y.size() ; k++) {
00827 sum_y += y[k] ;
00828 sum_ysq += y[k]*y[k] ;
00829 }
00830 sum_y /= y.size() ;
00831 sum_ysq /= y.size() ;
00832
00833 Double_t rms = sqrt(sum_ysq - (sum_y*sum_y)) ;
00834 lo = GetY()[i] - Z*rms ;
00835 hi = GetY()[i] + Z*rms ;
00836 }
00837 }
00838
00839
00840
00841
00842
00843 Bool_t RooCurve::isIdentical(const RooCurve& other, Double_t tol) const
00844 {
00845
00846
00847
00848
00849 Int_t n= min(GetN(),other.GetN());
00850 Double_t xmin(1e30), xmax(-1e30), ymin(1e30), ymax(-1e30) ;
00851 for(Int_t i= 0; i < n; i++) {
00852 if (fX[i]<xmin) xmin=fX[i] ;
00853 if (fX[i]>xmax) xmax=fX[i] ;
00854 if (fY[i]<ymin) ymin=fY[i] ;
00855 if (fY[i]>ymax) ymax=fY[i] ;
00856 }
00857 Double_t Yrange=ymax-ymin ;
00858
00859 Bool_t ret(kTRUE) ;
00860 for(Int_t i= 2; i < n-2; i++) {
00861 Double_t yTest = interpolate(other.fX[i],1e-10) ;
00862 Double_t rdy = fabs(yTest-other.fY[i])/Yrange ;
00863 if (rdy>tol) {
00864
00865
00866
00867
00868 cout << "RooCurve::isIdentical[" << i << "] Y tolerance exceeded (" << rdy << ">" << tol
00869 << "), X=" << other.fX[i] << "(" << fX[i] << ")" << " Ytest=" << yTest << " Yref=" << other.fY[i] << " range = " << Yrange << endl ;
00870 ret=kFALSE ;
00871 }
00872 }
00873
00874 return ret ;
00875 }
00876
00877