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 #include "RooFit.h"
00028
00029 #include "Roo2DKeysPdf.h"
00030 #include "Roo2DKeysPdf.h"
00031 #include "RooRealVar.h"
00032 #include "TTree.h"
00033 #include "TH2.h"
00034 #include "TFile.h"
00035 #include "TBranch.h"
00036 #include "TMath.h"
00037
00038
00039
00040 ClassImp(Roo2DKeysPdf)
00041
00042
00043
00044 Roo2DKeysPdf::Roo2DKeysPdf(const char *name, const char *title,
00045 RooAbsReal& xx, RooAbsReal & yy, RooDataSet& data, TString options, Double_t widthScaleFactor):
00046 RooAbsPdf(name,title),
00047 x("x", "x dimension",this, xx),
00048 y("y", "y dimension",this, yy)
00049 {
00050 setWidthScaleFactor(widthScaleFactor);
00051 loadDataSet(data, options);
00052 }
00053
00054
00055
00056 Roo2DKeysPdf::Roo2DKeysPdf(const Roo2DKeysPdf & other, const char* name) :
00057 RooAbsPdf(other,name),
00058 x("x", this, other.x),
00059 y("y", this, other.y)
00060 {
00061 if(_verbosedebug) { cout << "Roo2DKeysPdf::Roo2DKeysPdf copy ctor" << endl; }
00062
00063 _xMean = other._xMean;
00064 _xSigma = other._xSigma;
00065 _yMean = other._yMean;
00066 _ySigma = other._ySigma;
00067 _n = other._n;
00068
00069 _BandWidthType = other._BandWidthType;
00070 _MirrorAtBoundary = other._MirrorAtBoundary;
00071 _widthScaleFactor = other._widthScaleFactor;
00072
00073 _2pi = other._2pi;
00074 _sqrt2pi = other._sqrt2pi;
00075 _nEvents = other._nEvents;
00076 _n16 = other._n16;
00077 _debug = other._debug;
00078 _verbosedebug = other._verbosedebug;
00079 _vverbosedebug = other._vverbosedebug;
00080
00081 _lox = other._lox;
00082 _hix = other._hix;
00083 _loy = other._loy;
00084 _hiy = other._hiy;
00085 _xoffset = other._xoffset;
00086 _yoffset = other._yoffset;
00087
00088 _x = new Double_t[_nEvents];
00089 _y = new Double_t[_nEvents];
00090 _hx = new Double_t[_nEvents];
00091 _hy = new Double_t[_nEvents];
00092
00093
00094 for(Int_t iEvt = 0; iEvt< _nEvents; iEvt++)
00095 {
00096 _x[iEvt] = other._x[iEvt];
00097 _y[iEvt] = other._y[iEvt];
00098 _hx[iEvt] = other._hx[iEvt];
00099 _hy[iEvt] = other._hy[iEvt];
00100 }
00101 }
00102
00103
00104
00105 Roo2DKeysPdf::~Roo2DKeysPdf() {
00106 if(_verbosedebug) { cout << "Roo2DKeysPdf::Roo2KeysPdf dtor" << endl; }
00107 delete[] _x;
00108 delete[] _hx;
00109 delete[] _y;
00110 delete[] _hy;
00111 }
00112
00113
00114
00115
00116
00117
00118
00119
00120 Int_t Roo2DKeysPdf::loadDataSet(RooDataSet& data, TString options)
00121 {
00122 if(_verbosedebug) { cout << "Roo2DKeysPdf::loadDataSet" << endl; }
00123
00124 setOptions(options);
00125
00126 if(_verbosedebug) { cout << "Roo2DKeysPdf::loadDataSet(RooDataSet& data, TString options)" << endl; }
00127
00128 _2pi = 2.0*TMath::Pi() ;
00129 _sqrt2pi = sqrt(_2pi);
00130 _nEvents = (Int_t)data.numEntries();
00131 if(_nEvents == 0)
00132 {
00133 cout << "ERROR: Roo2DKeysPdf::loadDataSet The input data set is empty. Unable to begin generating the PDF" << endl;
00134 return 1;
00135 }
00136 _n16 = TMath::Power(_nEvents, -0.166666666);
00137
00138 _lox = x.min();
00139 _hix = x.max();
00140 _loy = y.min();
00141 _hiy = y.max();
00142
00143 _x = new Double_t[_nEvents];
00144 _y = new Double_t[_nEvents];
00145 _hx = new Double_t[_nEvents];
00146 _hy = new Double_t[_nEvents];
00147
00148 Double_t x0 = 0.0;
00149 Double_t x1 = 0.0;
00150 Double_t x_2 = 0.0;
00151 Double_t y0 = 0.0;
00152 Double_t y1 = 0.0;
00153 Double_t y_2 = 0.0;
00154
00155
00156 Int_t bad = 0;
00157 const RooAbsReal & xx = x.arg();
00158 const RooAbsReal & yy = y.arg();
00159 if(! (RooRealVar*)( (RooArgSet*)data.get(0) )->find( xx.GetName() ) )
00160 {
00161 cout << "Roo2DKeysPdf::Roo2DKeysPdf invalid RooAbsReal name: "<<xx.GetName()<<" not in the data set" <<endl;
00162 bad = 1;
00163 }
00164 if(! (RooRealVar*)( (RooArgSet*)data.get(0) )->find( yy.GetName() ) )
00165 {
00166 cout << "Roo2DKeysPdf::Roo2DKeysPdf invalid RooAbsReal name: "<<yy.GetName()<<" not in the data set" << endl;
00167 bad = 1;
00168 }
00169 if(bad)
00170 {
00171 cout << "Roo2DKeysPdf::Roo2DKeysPdf Unable to initilize object; incompatible RooDataSet doesn't contain"<<endl;
00172 cout << " all of the RooAbsReal arguments"<<endl;
00173 return 1;
00174 }
00175
00176
00177 const RooArgSet * values = data.get();
00178 const RooRealVar* X = ((RooRealVar*)(values->find(xx.GetName())) ) ;
00179 const RooRealVar* Y = ((RooRealVar*)(values->find(yy.GetName())) ) ;
00180
00181 for (Int_t j=0;j<_nEvents;++j)
00182 {
00183 data.get(j) ;
00184
00185 _x[j] = X->getVal() ;
00186 _y[j] = Y->getVal() ;
00187
00188 x0+=1; x1+=_x[j]; x_2+=_x[j]*_x[j];
00189 y0+=1; y1+=_y[j]; y_2+=_y[j]*_y[j];
00190 }
00191
00192
00193
00194
00195 if(_nEvents == 0)
00196 {
00197 cout << "Roo2DKeysPdf::Roo2DKeysPdf Empty data set was used; can't generate a PDF"<<endl;
00198 }
00199
00200 _xMean = x1/x0;
00201 _xSigma = sqrt(x_2/_nEvents-_xMean*_xMean);
00202
00203 _yMean = y1/y0;
00204 _ySigma = sqrt(y_2/_nEvents-_yMean*_yMean);
00205
00206 _n = Double_t(1)/(_2pi*_nEvents*_xSigma*_ySigma);
00207
00208
00209 return calculateBandWidth(_BandWidthType);
00210 }
00211
00212
00213
00214 void Roo2DKeysPdf::setOptions(TString options)
00215 {
00216 if(_verbosedebug) { cout << "Roo2DKeysPdf::setOptions" << endl; }
00217
00218 options.ToLower();
00219 if( options.Contains("a") ) _BandWidthType = 0;
00220 else _BandWidthType = 1;
00221 if( options.Contains("n") ) _BandWidthType = 1;
00222 else _BandWidthType = 0;
00223 if( options.Contains("m") ) _MirrorAtBoundary = 1;
00224 else _MirrorAtBoundary = 0;
00225 if( options.Contains("d") ) _debug = 1;
00226 else _debug = 0;
00227 if( options.Contains("v") ) { _debug = 1; _verbosedebug = 1; }
00228 else _verbosedebug = 0;
00229 if( options.Contains("vv") ) { _vverbosedebug = 1; }
00230 else _vverbosedebug = 0;
00231
00232 if( _debug )
00233 {
00234 cout << "Roo2DKeysPdf::setOptions(TString options) options = "<< options << endl;
00235 cout << "\t_BandWidthType = " << _BandWidthType << endl;
00236 cout << "\t_MirrorAtBoundary = " << _MirrorAtBoundary << endl;
00237 cout << "\t_debug = " << _debug << endl;
00238 cout << "\t_verbosedebug = " << _verbosedebug << endl;
00239 cout << "\t_vverbosedebug = " << _vverbosedebug << endl;
00240 }
00241 }
00242
00243
00244
00245 void Roo2DKeysPdf::getOptions(void) const
00246 {
00247 cout << "Roo2DKeysPdf::getOptions(void)" << endl;
00248 cout << "\t_BandWidthType = " << _BandWidthType << endl;
00249 cout << "\t_MirrorAtBoundary = " << _MirrorAtBoundary << endl;
00250 cout << "\t_debug = " << _debug << endl;
00251 cout << "\t_verbosedebug = " << _verbosedebug << endl;
00252 cout << "\t_vverbosedebug = " << _vverbosedebug << endl;
00253 }
00254
00255
00256
00257
00258
00259
00260
00261 Int_t Roo2DKeysPdf::calculateBandWidth(Int_t kernel)
00262 {
00263 if(_verbosedebug) { cout << "Roo2DKeysPdf::calculateBandWidth(Int_t kernel)" << endl; }
00264 if(kernel != -999)
00265 {
00266 _BandWidthType = kernel;
00267 }
00268
00269 Double_t h = 0.0;
00270
00271 Double_t sigSum = _xSigma*_xSigma + _ySigma*_ySigma;
00272 Double_t sqrtSum = sqrt( sigSum );
00273 Double_t sigProd = _ySigma*_xSigma;
00274 if(sigProd != 0.0) h = _n16*sqrt( sigSum/sigProd );
00275 if(sqrtSum == 0)
00276 {
00277 cout << "Roo2DKeysPdf::calculateBandWidth The sqr(variance sum) == 0.0. " << " Your dataset represents a delta function."<<endl;
00278 return 1;
00279 }
00280
00281 Double_t hXSigma = h * _xSigma;
00282 Double_t hYSigma = h * _ySigma;
00283 Double_t xhmin = hXSigma * sqrt(2.)/10;
00284 Double_t yhmin = hYSigma * sqrt(2.)/10;
00285
00286
00287
00288
00289 if(_BandWidthType == 1)
00290 {
00291 cout << "Roo2DKeysPdf::calculateBandWidth Using a normal bandwith (same for a given dimension) based on"<<endl;
00292 cout << " h_j = n^{-1/6}*sigma_j for the j^th dimension and n events * "<<_widthScaleFactor<<endl;
00293 Double_t hxGaussian = _n16 * _xSigma * _widthScaleFactor;
00294 Double_t hyGaussian = _n16 * _ySigma * _widthScaleFactor;
00295 for(Int_t j=0;j<_nEvents;++j)
00296 {
00297 _hx[j] = hxGaussian;
00298 _hy[j] = hyGaussian;
00299 if(_hx[j]<xhmin) _hx[j] = xhmin;
00300 if(_hy[j]<yhmin) _hy[j] = yhmin;
00301 }
00302 }
00303 else
00304 {
00305 cout << "Roo2DKeysPdf::calculateBandWidth Using an adaptive bandwith (in general different for all events) [default]"<<endl;
00306 cout << " scaled by a factor of "<<_widthScaleFactor<<endl;
00307 Double_t xnorm = h * TMath::Power(_xSigma/sqrtSum, 1.5) * _widthScaleFactor;
00308 Double_t ynorm = h * TMath::Power(_ySigma/sqrtSum, 1.5) * _widthScaleFactor;
00309 for(Int_t j=0;j<_nEvents;++j)
00310 {
00311 Double_t f_ti = TMath::Power( g(_x[j], _x, hXSigma, _y[j], _y, hYSigma), -0.25 ) ;
00312 _hx[j] = xnorm * f_ti;
00313 _hy[j] = ynorm * f_ti;
00314 if(_hx[j]<xhmin) _hx[j] = xhmin;
00315 if(_hy[j]<yhmin) _hy[j] = yhmin;
00316 }
00317 }
00318
00319 return 0;
00320 }
00321
00322
00323
00324
00325
00326
00327 Double_t Roo2DKeysPdf::evaluate() const
00328 {
00329
00330
00331
00332 if(_vverbosedebug) { cout << "Roo2DKeysPdf::evaluate()" << endl; }
00333 return evaluateFull(x,y);
00334 }
00335
00336
00337
00338
00339
00340
00341
00342
00343
00344 Double_t Roo2DKeysPdf::evaluateFull(Double_t thisX, Double_t thisY) const
00345 {
00346 if( _vverbosedebug ) { cout << "Roo2DKeysPdf::evaluateFull()" << endl; }
00347
00348 Double_t f=0.0;
00349
00350 Double_t rx2, ry2, zx, zy;
00351 if( _MirrorAtBoundary )
00352 {
00353 for (Int_t j = 0; j < _nEvents; ++j)
00354 {
00355 rx2 = 0.0; ry2 = 0.0; zx = 0.0; zy = 0.0;
00356 if(_hx[j] != 0.0) rx2 = (thisX - _x[j])/_hx[j];
00357 if(_hy[j] != 0.0) ry2 = (thisY - _y[j])/_hy[j];
00358
00359 if(_hx[j] != 0.0) zx = exp(-0.5*rx2*rx2)/_hx[j];
00360 if(_hy[j] != 0.0) zy = exp(-0.5*ry2*ry2)/_hy[j];
00361
00362 zx += highBoundaryCorrection(thisX, _hx[j], x.max(), _x[j])
00363 + lowBoundaryCorrection(thisX, _hx[j], x.min(), _x[j]);
00364 zy += highBoundaryCorrection(thisY, _hy[j], y.max(), _y[j])
00365 + lowBoundaryCorrection(thisY, _hy[j], y.min(), _y[j]);
00366 f += zy * zx;
00367
00368 }
00369 }
00370 else
00371 {
00372 for (Int_t j = 0; j < _nEvents; ++j)
00373 {
00374 rx2 = 0.0; ry2 = 0.0; zx = 0.0; zy = 0.0;
00375 if(_hx[j] != 0.0) rx2 = (thisX - _x[j])/_hx[j];
00376 if(_hy[j] != 0.0) ry2 = (thisY - _y[j])/_hy[j];
00377
00378 if(_hx[j] != 0.0) zx = exp(-0.5*rx2*rx2)/_hx[j];
00379 if(_hy[j] != 0.0) zy = exp(-0.5*ry2*ry2)/_hy[j];
00380 f += zy * zx;
00381
00382 }
00383 }
00384 return f;
00385 }
00386
00387
00388
00389
00390
00391
00392 Double_t Roo2DKeysPdf::highBoundaryCorrection(Double_t thisVar, Double_t thisH, Double_t high, Double_t tVar) const
00393 {
00394 if(_vverbosedebug) { cout << "Roo2DKeysPdf::highBoundaryCorrection" << endl; }
00395
00396 if(thisH == 0.0) return 0.0;
00397 Double_t correction = (thisVar + tVar - 2.0* high )/thisH;
00398 return exp(-0.5*correction*correction)/thisH;
00399 }
00400
00401
00402
00403 Double_t Roo2DKeysPdf::lowBoundaryCorrection(Double_t thisVar, Double_t thisH, Double_t low, Double_t tVar) const
00404 {
00405 if(_vverbosedebug) { cout << "Roo2DKeysPdf::lowBoundaryCorrection" << endl; }
00406
00407 if(thisH == 0.0) return 0.0;
00408 Double_t correction = (thisVar + tVar - 2.0* low )/thisH;
00409 return exp(-0.5*correction*correction)/thisH;
00410 }
00411
00412
00413
00414
00415
00416
00417
00418
00419
00420 Double_t Roo2DKeysPdf::g(Double_t varMean1, Double_t * _var1, Double_t sigma1, Double_t varMean2, Double_t * _var2, Double_t sigma2) const
00421 {
00422 if((_nEvents == 0.0) || (sigma1 == 0.0) || (sigma2 == 0)) return 0.0;
00423
00424 Double_t c1 = -1.0/(2.0*sigma1*sigma1);
00425 Double_t c2 = -1.0/(2.0*sigma2*sigma2);
00426 Double_t d = 4.0*c1*c2 /(_sqrt2pi*_nEvents);
00427 Double_t z = 0.0;
00428
00429 for (Int_t i = 0; i < _nEvents; ++i)
00430 {
00431 Double_t r1 = _var1[i] - varMean1;
00432 Double_t r2 = _var2[i] - varMean2;
00433 z += exp( c1 * r1*r1 ) * exp( c2 * r2*r2 );
00434 }
00435 z = z*d;
00436 return z;
00437 }
00438
00439
00440
00441 Int_t Roo2DKeysPdf::getBandWidthType() const
00442 {
00443 if(_BandWidthType == 1) cout << "The Bandwidth Type selected is Trivial" << endl;
00444 else cout << "The Bandwidth Type selected is Adaptive" << endl;
00445
00446 return _BandWidthType;
00447 }
00448
00449
00450
00451 Double_t Roo2DKeysPdf::getMean(const char * axis) const
00452 {
00453 if(!strcmp(axis,x.GetName()) || !strcmp(axis,"x") || !strcmp(axis,"X")) return _xMean;
00454 else if(!strcmp(axis,y.GetName()) || !strcmp(axis,"y") || !strcmp(axis,"Y")) return _yMean;
00455 else
00456 {
00457 cout << "Roo2DKeysPdf::getMean unknown axis "<<axis<<endl;
00458 }
00459 return 0.0;
00460 }
00461
00462
00463
00464 Double_t Roo2DKeysPdf::getSigma(const char * axis) const
00465 {
00466 if(!strcmp(axis,x.GetName()) || !strcmp(axis,"x") || !strcmp(axis,"X")) return _xSigma;
00467 else if(!strcmp(axis,y.GetName()) || !strcmp(axis,"y") || !strcmp(axis,"Y")) return _ySigma;
00468 else
00469 {
00470 cout << "Roo2DKeysPdf::getSigma unknown axis "<<axis<<endl;
00471 }
00472 return 0.0;
00473 }
00474
00475
00476
00477
00478 void Roo2DKeysPdf::writeToFile(char * outputFile, const char * name) const
00479 {
00480 TString histName = name;
00481 histName += "_hist";
00482 TString nName = name;
00483 nName += "_Ntuple";
00484 writeHistToFile( outputFile, histName);
00485 writeNTupleToFile( outputFile, nName);
00486 }
00487
00488
00489
00490
00491
00492
00493 void Roo2DKeysPdf::writeHistToFile(char * outputFile, const char * histName) const
00494 {
00495 TFile * file = 0;
00496 cout << "Roo2DKeysPdf::writeHistToFile This member function is temporarily disabled" <<endl;
00497
00498 file = new TFile(outputFile, "UPDATE");
00499 if (!file)
00500 {
00501 cout << "Roo2DKeysPdf::writeHistToFile unable to open file "<< outputFile <<endl;
00502 return;
00503 }
00504
00505
00506 const RooAbsReal & xx = x.arg();
00507 const RooAbsReal & yy = y.arg();
00508 RooArgSet values( RooArgList( xx, yy ));
00509 RooRealVar * xArg = ((RooRealVar*)(values.find(xx.GetName())) ) ;
00510 RooRealVar * yArg = ((RooRealVar*)(values.find(yy.GetName())) ) ;
00511
00512 TH2F * hist = (TH2F*)xArg->createHistogram("hist", *yArg);
00513 hist = (TH2F*)this->fillHistogram(hist, RooArgList(*xArg, *yArg) );
00514 hist->SetName(histName);
00515
00516 file->Write();
00517 file->Close();
00518 }
00519
00520
00521
00522
00523
00524
00525 void Roo2DKeysPdf::writeNTupleToFile(char * outputFile, const char * name) const
00526 {
00527 TFile * file = 0;
00528
00529
00530 file = new TFile(outputFile, "UPDATE");
00531 if (!file)
00532 {
00533 cout << "Roo2DKeysPdf::writeNTupleToFile unable to open file "<< outputFile <<endl;
00534 return;
00535 }
00536 RooAbsReal & xArg = (RooAbsReal&)x.arg();
00537 RooAbsReal & yArg = (RooAbsReal&)y.arg();
00538
00539 Double_t theX, theY, hx;
00540 TString label = name;
00541 label += " the source data for 2D Keys PDF";
00542 TTree * _theTree = new TTree(name, label);
00543 if(!_theTree) { cout << "Unable to get a TTree for output" << endl; return; }
00544 _theTree->SetAutoSave(1000000000);
00545
00546
00547 const char * xname = xArg.GetName();
00548 const char * yname = yArg.GetName();
00549 if (!strcmp(xname,"")) xname = "x";
00550 if (!strcmp(yname,"")) yname = "y";
00551
00552 _theTree->Branch(xname, &theX, " x/D");
00553 _theTree->Branch(yname, &theY, " y/D");
00554 _theTree->Branch("hx", &hx, " hx/D");
00555 _theTree->Branch("hy", &hx, " hy/D");
00556
00557 for(Int_t iEvt = 0; iEvt < _nEvents; iEvt++)
00558 {
00559 theX = _x[iEvt];
00560 theY = _y[iEvt];
00561 hx = _hx[iEvt];
00562 hx = _hy[iEvt];
00563 _theTree->Fill();
00564 }
00565 file->Write();
00566 file->Close();
00567 }
00568
00569
00570
00571
00572
00573
00574
00575 void Roo2DKeysPdf::PrintInfo(ostream & out) const
00576 {
00577 out << "Roo2DKeysPDF instance domain information:"<<endl;
00578 out << "\tX_min = " << _lox <<endl;
00579 out << "\tX_max = " << _hix <<endl;
00580 out << "\tY_min = " << _loy <<endl;
00581 out << "\tY_max = " << _hiy <<endl;
00582
00583 out << "Data information:" << endl;
00584 out << "\t<x> = " << _xMean <<endl;
00585 out << "\tsigma(x) = " << _xSigma <<endl;
00586 out << "\t<y> = " << _yMean <<endl;
00587 out << "\tsigma(y) = " << _ySigma <<endl;
00588
00589 out << "END of info for Roo2DKeys pdf instance"<< endl;
00590 }
00591
00592