Logo ROOT  
Reference Guide
Roo2DKeysPdf.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * @(#)root/roofit:$Id$
5 * Authors: *
6 * AB, Adrian Bevan, Liverpool University, bevan@slac.stanford.edu *
7 * *
8 * Copyright (c) 2000-2005, Regents of the University of California, *
9 * Liverpool University, *
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/** \class Roo2DKeysPdf
18 \ingroup Roofit
19
20Two-dimensional kernel estimation PDF.
21
22<b>This function has been superseded by the more general RooNDKeysPdf.</b>
23*/
24
25#include "RooFit.h"
26
27#include "Roo2DKeysPdf.h"
28#include "RooRealVar.h"
29#include "TTree.h"
30#include "TH2.h"
31#include "TFile.h"
32#include "TBranch.h"
33#include "TMath.h"
34
35//#include <math.h>
36
37using namespace std;
38
40
41
42////////////////////////////////////////////////////////////////////////////////
43/// Constructor.
44/// \param[in] name
45/// \param[in] title
46/// \param[in] xx
47/// \param[in] yy
48/// \param[in] data
49/// \param[in] options
50/// \param[in] widthScaleFactor
51
52Roo2DKeysPdf::Roo2DKeysPdf(const char *name, const char *title,
53 RooAbsReal& xx, RooAbsReal & yy, RooDataSet& data, TString options, Double_t widthScaleFactor):
54 RooAbsPdf(name,title),
55 x("x", "x dimension",this, xx),
56 y("y", "y dimension",this, yy)
57{
58 setWidthScaleFactor(widthScaleFactor);
59 loadDataSet(data, options);
60}
61
62
63////////////////////////////////////////////////////////////////////////////////
64/// Copy constructor.
65/// \param[in] other
66/// \param[in] name
67
68Roo2DKeysPdf::Roo2DKeysPdf(const Roo2DKeysPdf & other, const char* name) :
69 RooAbsPdf(other,name),
70 x("x", this, other.x),
71 y("y", this, other.y)
72{
73 if(_verbosedebug) { cout << "Roo2DKeysPdf::Roo2DKeysPdf copy ctor" << endl; }
74
75 _xMean = other._xMean;
76 _xSigma = other._xSigma;
77 _yMean = other._yMean;
78 _ySigma = other._ySigma;
79 _n = other._n;
80
84
85 _2pi = other._2pi;
86 _sqrt2pi = other._sqrt2pi;
87 _nEvents = other._nEvents;
88 _n16 = other._n16;
89 _debug = other._debug;
92
93 _lox = other._lox;
94 _hix = other._hix;
95 _loy = other._loy;
96 _hiy = other._hiy;
97 _xoffset = other._xoffset;
98 _yoffset = other._yoffset;
99
100 _x = new Double_t[_nEvents];
101 _y = new Double_t[_nEvents];
102 _hx = new Double_t[_nEvents];
103 _hy = new Double_t[_nEvents];
104
105 //copy the data and bandwidths
106 for(Int_t iEvt = 0; iEvt< _nEvents; iEvt++)
107 {
108 _x[iEvt] = other._x[iEvt];
109 _y[iEvt] = other._y[iEvt];
110 _hx[iEvt] = other._hx[iEvt];
111 _hy[iEvt] = other._hy[iEvt];
112 }
113}
114
115
116////////////////////////////////////////////////////////////////////////////////
117/// Destructor.
118
120 if(_verbosedebug) { cout << "Roo2DKeysPdf::Roo2KeysPdf dtor" << endl; }
121 delete[] _x;
122 delete[] _hx;
123 delete[] _y;
124 delete[] _hy;
125}
126
127
128////////////////////////////////////////////////////////////////////////////////
129/// Loads a new data set into the class instance.
130/// Returns 1 in case of error, 0 otherwise.
131/// \param[in] data
132/// \param[in] options
133
135{
136 if(_verbosedebug) { cout << "Roo2DKeysPdf::loadDataSet" << endl; }
137
138 setOptions(options);
139
140 if(_verbosedebug) { cout << "Roo2DKeysPdf::loadDataSet(RooDataSet& data, TString options)" << endl; }
141
142 _2pi = 2.0*TMath::Pi() ; //use pi from math.h
143 _sqrt2pi = sqrt(_2pi);
144 _nEvents = (Int_t)data.numEntries();
145 if(_nEvents == 0)
146 {
147 cout << "ERROR: Roo2DKeysPdf::loadDataSet The input data set is empty. Unable to begin generating the PDF" << endl;
148 return 1;
149 }
150 _n16 = TMath::Power(_nEvents, -0.166666666); // = (4/[n(dim(R) + 2)])^1/(dim(R)+4); dim(R) = 2
151
152 _lox = x.min();
153 _hix = x.max();
154 _loy = y.min();
155 _hiy = y.max();
156
157 _x = new Double_t[_nEvents];
158 _y = new Double_t[_nEvents];
159 _hx = new Double_t[_nEvents];
160 _hy = new Double_t[_nEvents];
161
162 Double_t x0 = 0.0;
163 Double_t x1 = 0.0;
164 Double_t x_2 = 0.0;
165 Double_t y0 = 0.0;
166 Double_t y1 = 0.0;
167 Double_t y_2 = 0.0;
168
169 //check that the data contain the variable we are interested in
170 Int_t bad = 0;
171 const RooAbsReal & xx = x.arg();
172 const RooAbsReal & yy = y.arg();
173 if(! (RooRealVar*)( (RooArgSet*)data.get(0) )->find( xx.GetName() ) )
174 {
175 cout << "Roo2DKeysPdf::Roo2DKeysPdf invalid RooAbsReal name: "<<xx.GetName()<<" not in the data set" <<endl;
176 bad = 1;
177 }
178 if(! (RooRealVar*)( (RooArgSet*)data.get(0) )->find( yy.GetName() ) )
179 {
180 cout << "Roo2DKeysPdf::Roo2DKeysPdf invalid RooAbsReal name: "<<yy.GetName()<<" not in the data set" << endl;
181 bad = 1;
182 }
183 if(bad)
184 {
185 cout << "Roo2DKeysPdf::Roo2DKeysPdf Unable to initilize object; incompatible RooDataSet doesn't contain"<<endl;
186 cout << " all of the RooAbsReal arguments"<<endl;
187 return 1;
188 }
189
190 //copy the data into local arrays
191 const RooArgSet * values = data.get();
192 const RooRealVar* X = ((RooRealVar*)(values->find(xx.GetName())) ) ;
193 const RooRealVar* Y = ((RooRealVar*)(values->find(yy.GetName())) ) ;
194
195 for (Int_t j=0;j<_nEvents;++j)
196 {
197 data.get(j) ;
198
199 _x[j] = X->getVal() ;
200 _y[j] = Y->getVal() ;
201
202 x0+=1; x1+=_x[j]; x_2+=_x[j]*_x[j];
203 y0+=1; y1+=_y[j]; y_2+=_y[j]*_y[j];
204 }
205
206 //==========================================//
207 //calculate the mean and sigma for the data //
208 //==========================================//
209 if(_nEvents == 0)
210 {
211 cout << "Roo2DKeysPdf::Roo2DKeysPdf Empty data set was used; can't generate a PDF"<<endl;
212 }
213
214 _xMean = x1/x0;
216
217 _yMean = y1/y0;
219
221
222 //calculate the PDF
224}
225
226
227////////////////////////////////////////////////////////////////////////////////
228
230{
231 if(_verbosedebug) { cout << "Roo2DKeysPdf::setOptions" << endl; }
232
233 options.ToLower();
234 if( options.Contains("a") ) _BandWidthType = 0;
235 else _BandWidthType = 1;
236 if( options.Contains("n") ) _BandWidthType = 1;
237 else _BandWidthType = 0;
238 if( options.Contains("m") ) _MirrorAtBoundary = 1;
239 else _MirrorAtBoundary = 0;
240 if( options.Contains("d") ) _debug = 1;
241 else _debug = 0;
242 if( options.Contains("v") ) { _debug = 1; _verbosedebug = 1; }
243 else _verbosedebug = 0;
244 if( options.Contains("vv") ) { _vverbosedebug = 1; }
245 else _vverbosedebug = 0;
246
247 if( _debug )
248 {
249 cout << "Roo2DKeysPdf::setOptions(TString options) options = "<< options << endl;
250 cout << "\t_BandWidthType = " << _BandWidthType << endl;
251 cout << "\t_MirrorAtBoundary = " << _MirrorAtBoundary << endl;
252 cout << "\t_debug = " << _debug << endl;
253 cout << "\t_verbosedebug = " << _verbosedebug << endl;
254 cout << "\t_vverbosedebug = " << _vverbosedebug << endl;
255 }
256}
257
258
259////////////////////////////////////////////////////////////////////////////////
260
262{
263 cout << "Roo2DKeysPdf::getOptions(void)" << endl;
264 cout << "\t_BandWidthType = " << _BandWidthType << endl;
265 cout << "\t_MirrorAtBoundary = " << _MirrorAtBoundary << endl;
266 cout << "\t_debug = " << _debug << endl;
267 cout << "\t_verbosedebug = " << _verbosedebug << endl;
268 cout << "\t_vverbosedebug = " << _vverbosedebug << endl;
269}
270
271
272////////////////////////////////////////////////////////////////////////////////
273/// Calculates the kernel bandwidth for x & y and the probability look up table _p[i][j]
274/// \param[in] kernel
275
277{
278 if(_verbosedebug) { cout << "Roo2DKeysPdf::calculateBandWidth(Int_t kernel)" << endl; }
279 if(kernel != -999)
280 {
281 _BandWidthType = kernel;
282 }
283
284 Double_t h = 0.0;
285
287 Double_t sqrtSum = sqrt( sigSum );
288 Double_t sigProd = _ySigma*_xSigma;
289 if(sigProd != 0.0) h = _n16*sqrt( sigSum/sigProd );
290 if(sqrtSum == 0)
291 {
292 cout << "Roo2DKeysPdf::calculateBandWidth The sqr(variance sum) == 0.0. " << " Your dataset represents a delta function."<<endl;
293 return 1;
294 }
295
296 Double_t hXSigma = h * _xSigma;
297 Double_t hYSigma = h * _ySigma;
298 Double_t xhmin = hXSigma * sqrt(2.)/10; //smallest anticipated bandwidth
299 Double_t yhmin = hYSigma * sqrt(2.)/10;
300
301 //////////////////////////////////////
302 //calculate bandwidths from the data//
303 //////////////////////////////////////
304 if(_BandWidthType == 1) //calculate a trivial bandwidth
305 {
306 cout << "Roo2DKeysPdf::calculateBandWidth Using a normal bandwidth (same for a given dimension) based on"<<endl;
307 cout << " h_j = n^{-1/6}*sigma_j for the j^th dimension and n events * "<<_widthScaleFactor<<endl;
308 Double_t hxGaussian = _n16 * _xSigma * _widthScaleFactor;
309 Double_t hyGaussian = _n16 * _ySigma * _widthScaleFactor;
310 for(Int_t j=0;j<_nEvents;++j)
311 {
312 _hx[j] = hxGaussian;
313 _hy[j] = hyGaussian;
314 if(_hx[j]<xhmin) _hx[j] = xhmin;
315 if(_hy[j]<yhmin) _hy[j] = yhmin;
316 }
317 }
318 else //use an adaptive bandwidth to reduce the dependence on global data distribution
319 {
320 cout << "Roo2DKeysPdf::calculateBandWidth Using an adaptive bandwidth (in general different for all events) [default]"<<endl;
321 cout << " scaled by a factor of "<<_widthScaleFactor<<endl;
322 Double_t xnorm = h * TMath::Power(_xSigma/sqrtSum, 1.5) * _widthScaleFactor;
323 Double_t ynorm = h * TMath::Power(_ySigma/sqrtSum, 1.5) * _widthScaleFactor;
324 for(Int_t j=0;j<_nEvents;++j)
325 {
326 Double_t f_ti = TMath::Power( g(_x[j], _x, hXSigma, _y[j], _y, hYSigma), -0.25 ) ;
327 _hx[j] = xnorm * f_ti;
328 _hy[j] = ynorm * f_ti;
329 if(_hx[j]<xhmin) _hx[j] = xhmin;
330 if(_hy[j]<yhmin) _hy[j] = yhmin;
331 }
332 }
333
334 return 0;
335}
336
337
338////////////////////////////////////////////////////////////////////////////////
339/// Evaluates the kernel estimation for x,y, interpolating between the points if necessary
340///
341/// Uses the caching intrinsic in RFC to bypass the grid and remove
342/// the grid and extrapolation approximation in the kernel estimation method
343/// implementation.
344
346{
347 if(_vverbosedebug) { cout << "Roo2DKeysPdf::evaluate()" << endl; }
348 return evaluateFull(x,y);
349}
350
351
352////////////////////////////////////////////////////////////////////////////////
353/// Evaluates the sum of the product of the 2D kernels
354/// for use in calculating the fixed kernel estimate, f,
355/// given the bandwidths _hx[j] and _hy[j].
356///
357/// _n is calculated once in the constructor.
358/// \param[in] thisX
359/// \param[in] thisY
360
362{
363 if( _vverbosedebug ) { cout << "Roo2DKeysPdf::evaluateFull()" << endl; }
364
365 Double_t f=0.0;
366
367 Double_t rx2, ry2, zx, zy;
369 {
370 for (Int_t j = 0; j < _nEvents; ++j)
371 {
372 rx2 = 0.0; ry2 = 0.0; zx = 0.0; zy = 0.0;
373 if(_hx[j] != 0.0) rx2 = (thisX - _x[j])/_hx[j];
374 if(_hy[j] != 0.0) ry2 = (thisY - _y[j])/_hy[j];
375
376 if(_hx[j] != 0.0) zx = exp(-0.5*rx2*rx2)/_hx[j];
377 if(_hy[j] != 0.0) zy = exp(-0.5*ry2*ry2)/_hy[j];
378
379 zx += highBoundaryCorrection(thisX, _hx[j], x.max(), _x[j])
380 + lowBoundaryCorrection(thisX, _hx[j], x.min(), _x[j]);
381 zy += highBoundaryCorrection(thisY, _hy[j], y.max(), _y[j])
382 + lowBoundaryCorrection(thisY, _hy[j], y.min(), _y[j]);
383 f += zy * zx;
384 // f += _n * zy * zx; // ooops this is a normalisation factor :(
385 }
386 }
387 else
388 {
389 for (Int_t j = 0; j < _nEvents; ++j)
390 {
391 rx2 = 0.0; ry2 = 0.0; zx = 0.0; zy = 0.0;
392 if(_hx[j] != 0.0) rx2 = (thisX - _x[j])/_hx[j];
393 if(_hy[j] != 0.0) ry2 = (thisY - _y[j])/_hy[j];
394
395 if(_hx[j] != 0.0) zx = exp(-0.5*rx2*rx2)/_hx[j];
396 if(_hy[j] != 0.0) zy = exp(-0.5*ry2*ry2)/_hy[j];
397 f += zy * zx;
398 // f += _n * zy * zx; // ooops this is a normalisation factor :(
399 }
400 }
401 return f;
402}
403
404
405////////////////////////////////////////////////////////////////////////////////
406/// Apply the mirror at boundary correction to a dimension given the space position to evaluate
407/// at (thisVar), the bandwidth at this position (thisH), the boundary (high/low) and the
408/// value of the data kernel that this correction is being applied to tVar (i.e. the _x[ix] etc.).
409/// \param[in] thisVar
410/// \param[in] thisH
411/// \param[in] high
412/// \param[in] tVar
413
415{
416 if(_vverbosedebug) { cout << "Roo2DKeysPdf::highBoundaryCorrection" << endl; }
417
418 if(thisH == 0.0) return 0.0;
419 Double_t correction = (thisVar + tVar - 2.0* high )/thisH;
420 return exp(-0.5*correction*correction)/thisH;
421}
422
423
424////////////////////////////////////////////////////////////////////////////////
425
427{
428 if(_vverbosedebug) { cout << "Roo2DKeysPdf::lowBoundaryCorrection" << endl; }
429
430 if(thisH == 0.0) return 0.0;
431 Double_t correction = (thisVar + tVar - 2.0* low )/thisH;
432 return exp(-0.5*correction*correction)/thisH;
433}
434
435
436////////////////////////////////////////////////////////////////////////////////
437/// Calculates f(t_i) for the bandwidths.
438/// \f$ g = 1/(N_{evt} * \sigma_j * \sqrt{2\pi})*\sum_{all evts}{\prod d K[ \exp({-(xd - ti)/\sigma_{j}d^2}) ]}\f$
439/// \param[in] varMean1
440/// \param[in] _var1
441/// \param[in] sigma1
442/// \param[in] varMean2
443/// \param[in] _var2
444/// \param[in] sigma2
445
446Double_t Roo2DKeysPdf::g(Double_t varMean1, Double_t * _var1, Double_t sigma1, Double_t varMean2, Double_t * _var2, Double_t sigma2) const
447{
448 if((_nEvents == 0.0) || (sigma1 == 0.0) || (sigma2 == 0)) return 0.0;
449
450 Double_t c1 = -1.0/(2.0*sigma1*sigma1);
451 Double_t c2 = -1.0/(2.0*sigma2*sigma2);
452 Double_t d = 4.0*c1*c2 /(_sqrt2pi*_nEvents);
453 Double_t z = 0.0;
454
455 for (Int_t i = 0; i < _nEvents; ++i)
456 {
457 Double_t r1 = _var1[i] - varMean1;
458 Double_t r2 = _var2[i] - varMean2;
459 z += exp( c1 * r1*r1 ) * exp( c2 * r2*r2 );
460 }
461 z = z*d;
462 return z;
463}
464
465
466////////////////////////////////////////////////////////////////////////////////
467
469{
470 if(_BandWidthType == 1) cout << "The Bandwidth Type selected is Trivial" << endl;
471 else cout << "The Bandwidth Type selected is Adaptive" << endl;
472
473 return _BandWidthType;
474}
475
476
477////////////////////////////////////////////////////////////////////////////////
478
479Double_t Roo2DKeysPdf::getMean(const char * axis) const
480{
481 if(!strcmp(axis,x.GetName()) || !strcmp(axis,"x") || !strcmp(axis,"X")) return _xMean;
482 else if(!strcmp(axis,y.GetName()) || !strcmp(axis,"y") || !strcmp(axis,"Y")) return _yMean;
483 else
484 {
485 cout << "Roo2DKeysPdf::getMean unknown axis "<<axis<<endl;
486 }
487 return 0.0;
488}
489
490
491////////////////////////////////////////////////////////////////////////////////
492
493Double_t Roo2DKeysPdf::getSigma(const char * axis) const
494{
495 if(!strcmp(axis,x.GetName()) || !strcmp(axis,"x") || !strcmp(axis,"X")) return _xSigma;
496 else if(!strcmp(axis,y.GetName()) || !strcmp(axis,"y") || !strcmp(axis,"Y")) return _ySigma;
497 else
498 {
499 cout << "Roo2DKeysPdf::getSigma unknown axis "<<axis<<endl;
500 }
501 return 0.0;
502}
503
504
505
506////////////////////////////////////////////////////////////////////////////////
507
508void Roo2DKeysPdf::writeToFile(char * outputFile, const char * name) const
509{
510 TString histName = name;
511 histName += "_hist";
512 TString nName = name;
513 nName += "_Ntuple";
514 writeHistToFile( outputFile, histName);
515 writeNTupleToFile( outputFile, nName);
516}
517
518
519////////////////////////////////////////////////////////////////////////////////
520/// Plots the PDF as a histogram and saves it to a file, so that it can be loaded in
521/// as a Roo2DHist PDF in the future to save on calculation time.
522/// \param[in] outputFile Name of the file where to store the PDF
523/// \param[in] histName PDF histogram name
524
525void Roo2DKeysPdf::writeHistToFile(char * outputFile, const char * histName) const
526{
527 TFile * file = 0;
528 cout << "Roo2DKeysPdf::writeHistToFile This member function is temporarily disabled" <<endl;
529 //make sure that any existing file is not over written
530 file = new TFile(outputFile, "UPDATE");
531 if (!file)
532 {
533 cout << "Roo2DKeysPdf::writeHistToFile unable to open file "<< outputFile <<endl;
534 return;
535 }
536
537
538 const RooAbsReal & xx = x.arg();
539 const RooAbsReal & yy = y.arg();
540 RooArgSet values( RooArgList( xx, yy ));
541 RooRealVar * xArg = ((RooRealVar*)(values.find(xx.GetName())) ) ;
542 RooRealVar * yArg = ((RooRealVar*)(values.find(yy.GetName())) ) ;
543
544 TH2F * hist = (TH2F*)xArg->createHistogram("hist", *yArg);
545 hist = (TH2F*)this->fillHistogram(hist, RooArgList(*xArg, *yArg) );
546 hist->SetName(histName);
547
548 file->Write();
549 file->Close();
550}
551
552
553////////////////////////////////////////////////////////////////////////////////
554/// Saves the data and calculated bandwidths to a file,
555/// as a record of what produced the PDF and to give a reduced
556/// data set in order to facilitate re-calculation in the future.
557/// \param[in] outputFile Name of the file where to store the data
558/// \param[in] name Name of the tree which will contain the data
559
560void Roo2DKeysPdf::writeNTupleToFile(char * outputFile, const char * name) const
561{
562 TFile * file = 0;
563
564 //make sure that any existing file is not over written
565 file = new TFile(outputFile, "UPDATE");
566 if (!file)
567 {
568 cout << "Roo2DKeysPdf::writeNTupleToFile unable to open file "<< outputFile <<endl;
569 return;
570 }
571 RooAbsReal & xArg = (RooAbsReal&)x.arg();
572 RooAbsReal & yArg = (RooAbsReal&)y.arg();
573
574 Double_t theX, theY, hx/*, hy*/;
575 TString label = name;
576 label += " the source data for 2D Keys PDF";
577 TTree * _theTree = new TTree(name, label);
578 if(!_theTree) { cout << "Unable to get a TTree for output" << endl; return; }
579 _theTree->SetAutoSave(1000000000); // autosave when 1 Gbyte written
580
581 //name the TBranches the same as the RooAbsReal's
582 const char * xname = xArg.GetName();
583 const char * yname = yArg.GetName();
584 if (!strcmp(xname,"")) xname = "x";
585 if (!strcmp(yname,"")) yname = "y";
586
587 _theTree->Branch(xname, &theX, " x/D");
588 _theTree->Branch(yname, &theY, " y/D");
589 _theTree->Branch("hx", &hx, " hx/D");
590 _theTree->Branch("hy", &hx, " hy/D");
591
592 for(Int_t iEvt = 0; iEvt < _nEvents; iEvt++)
593 {
594 theX = _x[iEvt];
595 theY = _y[iEvt];
596 hx = _hx[iEvt];
597 hx = _hy[iEvt];
598 _theTree->Fill();
599 }
600 file->Write();
601 file->Close();
602}
603
604
605////////////////////////////////////////////////////////////////////////////////
606/// Prints out _p[_nPoints][_nPoints] indicating the domain limits.
607/// \param[out] out Output stream where to print
608
609void Roo2DKeysPdf::PrintInfo(ostream & out) const
610{
611 out << "Roo2DKeysPDF instance domain information:"<<endl;
612 out << "\tX_min = " << _lox <<endl;
613 out << "\tX_max = " << _hix <<endl;
614 out << "\tY_min = " << _loy <<endl;
615 out << "\tY_max = " << _hiy <<endl;
616
617 out << "Data information:" << endl;
618 out << "\t<x> = " << _xMean <<endl;
619 out << "\tsigma(x) = " << _xSigma <<endl;
620 out << "\t<y> = " << _yMean <<endl;
621 out << "\tsigma(y) = " << _ySigma <<endl;
622
623 out << "END of info for Roo2DKeys pdf instance"<< endl;
624}
#define d(i)
Definition: RSha256.hxx:102
#define f(i)
Definition: RSha256.hxx:104
#define h(i)
Definition: RSha256.hxx:106
static const double x1[5]
int Int_t
Definition: RtypesCore.h:43
double Double_t
Definition: RtypesCore.h:57
#define ClassImp(name)
Definition: Rtypes.h:361
char name[80]
Definition: TGX11.cxx:109
double sqrt(double)
double exp(double)
Two-dimensional kernel estimation PDF.
Definition: Roo2DKeysPdf.h:26
RooRealProxy y
Definition: Roo2DKeysPdf.h:79
Double_t _lox
Definition: Roo2DKeysPdf.h:108
Double_t _loy
Definition: Roo2DKeysPdf.h:109
void writeNTupleToFile(char *outputFile, const char *name) const
Saves the data and calculated bandwidths to a file, as a record of what produced the PDF and to give ...
Int_t _BandWidthType
Definition: Roo2DKeysPdf.h:115
Double_t getMean(const char *axis) const
Double_t g(Double_t var1, Double_t *_var1, Double_t sigma1, Double_t var2, Double_t *_var2, Double_t sigma2) const
Calculates f(t_i) for the bandwidths.
Double_t _hix
Definition: Roo2DKeysPdf.h:108
void writeToFile(char *outputFile, const char *name) const
Int_t _verbosedebug
Definition: Roo2DKeysPdf.h:118
Double_t _ySigma
Definition: Roo2DKeysPdf.h:103
Double_t * _x
Definition: Roo2DKeysPdf.h:95
Double_t _n
Definition: Roo2DKeysPdf.h:104
Roo2DKeysPdf(const char *name, const char *title, RooAbsReal &xx, RooAbsReal &yy, RooDataSet &data, TString options="a", Double_t widthScaleFactor=1.0)
Constructor.
Double_t _n16
Definition: Roo2DKeysPdf.h:105
Double_t getSigma(const char *axis) const
Double_t _hiy
Definition: Roo2DKeysPdf.h:109
Double_t _widthScaleFactor
Definition: Roo2DKeysPdf.h:112
Double_t * _hy
Definition: Roo2DKeysPdf.h:98
Double_t _yoffset
Definition: Roo2DKeysPdf.h:111
Int_t getBandWidthType() const
Double_t _xSigma
Definition: Roo2DKeysPdf.h:101
Int_t _MirrorAtBoundary
Definition: Roo2DKeysPdf.h:116
void PrintInfo(std::ostream &) const
Prints out _p[_nPoints][_nPoints] indicating the domain limits.
Int_t loadDataSet(RooDataSet &data, TString options)
Loads a new data set into the class instance.
Double_t * _y
Definition: Roo2DKeysPdf.h:97
Int_t _vverbosedebug
Definition: Roo2DKeysPdf.h:119
Double_t * _hx
Definition: Roo2DKeysPdf.h:96
Double_t evaluate() const
Evaluates the kernel estimation for x,y, interpolating between the points if necessary.
virtual ~Roo2DKeysPdf()
Destructor.
Double_t highBoundaryCorrection(Double_t thisVar, Double_t thisH, Double_t high, Double_t tVar) const
Apply the mirror at boundary correction to a dimension given the space position to evaluate at (thisV...
Double_t _xoffset
Definition: Roo2DKeysPdf.h:110
Double_t _xMean
Definition: Roo2DKeysPdf.h:100
Double_t lowBoundaryCorrection(Double_t thisVar, Double_t thisH, Double_t low, Double_t tVar) const
Double_t _yMean
Definition: Roo2DKeysPdf.h:102
void writeHistToFile(char *outputFile, const char *histName) const
Plots the PDF as a histogram and saves it to a file, so that it can be loaded in as a Roo2DHist PDF i...
Int_t calculateBandWidth(Int_t kernel=-999)
Calculates the kernel bandwidth for x & y and the probability look up table _p[i][j].
Double_t _sqrt2pi
Definition: Roo2DKeysPdf.h:106
void setOptions(TString options)
Double_t evaluateFull(Double_t thisX, Double_t thisY) const
Evaluates the sum of the product of the 2D kernels for use in calculating the fixed kernel estimate,...
void getOptions(void) const
void setWidthScaleFactor(Double_t widthScaleFactor)
Definition: Roo2DKeysPdf.h:124
Double_t _2pi
Definition: Roo2DKeysPdf.h:107
RooRealProxy x
Definition: Roo2DKeysPdf.h:78
RooAbsArg * find(const char *name) const
Find object with given name in list.
virtual Int_t numEntries() const
Definition: RooAbsData.cxx:306
TH1 * createHistogram(const char *name, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:60
TH1 * fillHistogram(TH1 *hist, const RooArgList &plotVars, Double_t scaleFactor=1, const RooArgSet *projectedVars=0, Bool_t scaling=kTRUE, const RooArgSet *condObs=0, Bool_t setError=kTRUE) const
Fill the ROOT histogram 'hist' with values sampled from this function at the bin centers.
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:90
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgList.h:21
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:33
virtual const RooArgSet * get(Int_t index) const override
Return RooArgSet with coordinates of event 'index'.
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
double min(const char *rname=0) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
double max(const char *rname=0) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
const T & arg() const
Return reference to object held in proxy.
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition: TFile.h:53
virtual void SetName(const char *name)
Change the name of this histogram.
Definition: TH1.cxx:8416
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:251
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
Basic string class.
Definition: TString.h:131
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1125
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:619
A TTree represents a columnar dataset.
Definition: TTree.h:78
virtual Int_t Fill()
Fill all branches.
Definition: TTree.cxx:4524
virtual void SetAutoSave(Long64_t autos=-300000000)
This function may be called at the start of a program to change the default value for fAutoSave (and ...
Definition: TTree.cxx:8194
TBranch * Branch(const char *name, T *obj, Int_t bufsize=32000, Int_t splitlevel=99)
Add a new branch, and infer the data type from the type of obj being passed.
Definition: TTree.h:348
return c1
Definition: legend1.C:41
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
return c2
Definition: legend2.C:14
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:725
constexpr Double_t Pi()
Definition: TMath.h:38
Definition: file.py:1