Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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 "Roo2DKeysPdf.h"
26#include "RooRealVar.h"
27#include "TTree.h"
28#include "TH2.h"
29#include "TFile.h"
30#include "TMath.h"
31
32//#include <math.h>
33
34using std::cout, std::endl, std::ostream;
35
37
38
39////////////////////////////////////////////////////////////////////////////////
40/// Constructor.
41/// \param[in] name
42/// \param[in] title
43/// \param[in] xx
44/// \param[in] yy
45/// \param[in] data
46/// \param[in] options
47/// \param[in] widthScaleFactor
48
49Roo2DKeysPdf::Roo2DKeysPdf(const char *name, const char *title,
50 RooAbsReal& xx, RooAbsReal & yy, RooDataSet& data, TString options, double widthScaleFactor):
51 RooAbsPdf(name,title),
52 x("x", "x dimension",this, xx),
53 y("y", "y dimension",this, yy)
54{
55 setWidthScaleFactor(widthScaleFactor);
56 loadDataSet(data, options);
57}
58
59
60////////////////////////////////////////////////////////////////////////////////
61/// Copy constructor.
62/// \param[in] other
63/// \param[in] name
64
65Roo2DKeysPdf::Roo2DKeysPdf(const Roo2DKeysPdf & other, const char* name) :
66 RooAbsPdf(other,name),
67 x("x", this, other.x),
68 y("y", this, other.y)
69{
70 if(_verbosedebug) { cout << "Roo2DKeysPdf::Roo2DKeysPdf copy ctor" << endl; }
71
72 _xMean = other._xMean;
73 _xSigma = other._xSigma;
74 _yMean = other._yMean;
75 _ySigma = other._ySigma;
76 _n = other._n;
77
81
82 _2pi = other._2pi;
83 _sqrt2pi = other._sqrt2pi;
84 _nEvents = other._nEvents;
85 _n16 = other._n16;
86 _debug = other._debug;
89
90 _lox = other._lox;
91 _hix = other._hix;
92 _loy = other._loy;
93 _hiy = other._hiy;
94 _xoffset = other._xoffset;
95 _yoffset = other._yoffset;
96
97 _x = new double[_nEvents];
98 _y = new double[_nEvents];
99 _hx = new double[_nEvents];
100 _hy = new double[_nEvents];
101
102 //copy the data and bandwidths
103 for(Int_t iEvt = 0; iEvt< _nEvents; iEvt++)
104 {
105 _x[iEvt] = other._x[iEvt];
106 _y[iEvt] = other._y[iEvt];
107 _hx[iEvt] = other._hx[iEvt];
108 _hy[iEvt] = other._hy[iEvt];
109 }
110}
111
112
113////////////////////////////////////////////////////////////////////////////////
114/// Destructor.
115
117 if(_verbosedebug) { cout << "Roo2DKeysPdf::Roo2KeysPdf dtor" << endl; }
118 delete[] _x;
119 delete[] _hx;
120 delete[] _y;
121 delete[] _hy;
122}
123
124
125////////////////////////////////////////////////////////////////////////////////
126/// Loads a new data set into the class instance.
127/// Returns 1 in case of error, 0 otherwise.
128/// \param[in] data
129/// \param[in] options
130
132{
133 if(_verbosedebug) { cout << "Roo2DKeysPdf::loadDataSet" << endl; }
134
135 setOptions(options);
136
137 if(_verbosedebug) { cout << "Roo2DKeysPdf::loadDataSet(RooDataSet& data, TString options)" << endl; }
138
139 _2pi = 2.0*TMath::Pi() ; //use pi from math.h
140 _sqrt2pi = sqrt(_2pi);
141 _nEvents = (Int_t)data.numEntries();
142 if(_nEvents == 0)
143 {
144 cout << "ERROR: Roo2DKeysPdf::loadDataSet The input data set is empty. Unable to begin generating the PDF" << endl;
145 return 1;
146 }
147 _n16 = std::pow(_nEvents, -1./6.); // = (4/[n(dim(R) + 2)])^1/(dim(R)+4); dim(R) = 2
148
149 _lox = x.min();
150 _hix = x.max();
151 _loy = y.min();
152 _hiy = y.max();
153
154 _x = new double[_nEvents];
155 _y = new double[_nEvents];
156 _hx = new double[_nEvents];
157 _hy = new double[_nEvents];
158
159 double x0 = 0.0;
160 double x1 = 0.0;
161 double x_2 = 0.0;
162 double y0 = 0.0;
163 double y1 = 0.0;
164 double y_2 = 0.0;
165
166 //check that the data contain the variable we are interested in
167 Int_t bad = 0;
168 const RooAbsReal & xx = x.arg();
169 const RooAbsReal & yy = y.arg();
170 if(! static_cast<RooRealVar*>((const_cast<RooArgSet *>(data.get(0)))->find( xx.GetName() )) )
171 {
172 cout << "Roo2DKeysPdf::Roo2DKeysPdf invalid RooAbsReal name: "<<xx.GetName()<<" not in the data set" <<endl;
173 bad = 1;
174 }
175 if(! static_cast<RooRealVar*>((const_cast<RooArgSet *>(data.get(0)))->find( yy.GetName() )) )
176 {
177 cout << "Roo2DKeysPdf::Roo2DKeysPdf invalid RooAbsReal name: "<<yy.GetName()<<" not in the data set" << endl;
178 bad = 1;
179 }
180 if(bad)
181 {
182 cout << "Roo2DKeysPdf::Roo2DKeysPdf Unable to initialize object; incompatible RooDataSet doesn't contain"<<endl;
183 cout << " all of the RooAbsReal arguments"<<endl;
184 return 1;
185 }
186
187 //copy the data into local arrays
188 const RooArgSet * values = data.get();
189 auto X = static_cast<RooRealVar const*>(values->find(xx.GetName()));
190 auto Y = static_cast<RooRealVar const*>(values->find(yy.GetName()));
191
192 for (Int_t j=0;j<_nEvents;++j)
193 {
194 data.get(j) ;
195
196 _x[j] = X->getVal() ;
197 _y[j] = Y->getVal() ;
198
199 x0+=1; x1+=_x[j]; x_2+=_x[j]*_x[j];
200 y0+=1; y1+=_y[j]; y_2+=_y[j]*_y[j];
201 }
202
203 //==========================================//
204 //calculate the mean and sigma for the data //
205 //==========================================//
206 if(_nEvents == 0)
207 {
208 cout << "Roo2DKeysPdf::Roo2DKeysPdf Empty data set was used; can't generate a PDF"<<endl;
209 }
210
211 _xMean = x1/x0;
212 _xSigma = sqrt(x_2/_nEvents-_xMean*_xMean);
213
214 _yMean = y1/y0;
215 _ySigma = sqrt(y_2/_nEvents-_yMean*_yMean);
216
218
219 //calculate the PDF
221}
222
223
224////////////////////////////////////////////////////////////////////////////////
225
227{
228 if(_verbosedebug) { cout << "Roo2DKeysPdf::setOptions" << endl; }
229
230 options.ToLower();
231 if( options.Contains("a") ) _BandWidthType = 0;
232 else _BandWidthType = 1;
233 if( options.Contains("n") ) _BandWidthType = 1;
234 else _BandWidthType = 0;
235 if( options.Contains("m") ) _MirrorAtBoundary = 1;
236 else _MirrorAtBoundary = 0;
237 if( options.Contains("d") ) _debug = 1;
238 else _debug = 0;
239 if( options.Contains("v") ) { _debug = 1; _verbosedebug = 1; }
240 else _verbosedebug = 0;
241 if( options.Contains("vv") ) { _vverbosedebug = 1; }
242 else _vverbosedebug = 0;
243
244 if( _debug )
245 {
246 cout << "Roo2DKeysPdf::setOptions(TString options) options = "<< options << endl;
247 cout << "\t_BandWidthType = " << _BandWidthType << endl;
248 cout << "\t_MirrorAtBoundary = " << _MirrorAtBoundary << endl;
249 cout << "\t_debug = " << _debug << endl;
250 cout << "\t_verbosedebug = " << _verbosedebug << endl;
251 cout << "\t_vverbosedebug = " << _vverbosedebug << endl;
252 }
253}
254
255
256////////////////////////////////////////////////////////////////////////////////
257
259{
260 cout << "Roo2DKeysPdf::getOptions(void)" << endl;
261 cout << "\t_BandWidthType = " << _BandWidthType << endl;
262 cout << "\t_MirrorAtBoundary = " << _MirrorAtBoundary << endl;
263 cout << "\t_debug = " << _debug << endl;
264 cout << "\t_verbosedebug = " << _verbosedebug << endl;
265 cout << "\t_vverbosedebug = " << _vverbosedebug << endl;
266}
267
268
269////////////////////////////////////////////////////////////////////////////////
270/// Calculates the kernel bandwidth for x & y and the probability look up table _p[i][j]
271/// \param[in] kernel
272
274{
275 if(_verbosedebug) { cout << "Roo2DKeysPdf::calculateBandWidth(Int_t kernel)" << endl; }
276 if(kernel != -999)
277 {
278 _BandWidthType = kernel;
279 }
280
281 double h = 0.0;
282
283 double sigSum = _xSigma*_xSigma + _ySigma*_ySigma;
284 double sqrtSum = sqrt( sigSum );
285 double sigProd = _ySigma*_xSigma;
286 if(sigProd != 0.0) h = _n16*sqrt( sigSum/sigProd );
287 if(sqrtSum == 0)
288 {
289 cout << "Roo2DKeysPdf::calculateBandWidth The sqr(variance sum) == 0.0. " << " Your dataset represents a delta function."<<endl;
290 return 1;
291 }
292
293 double hXSigma = h * _xSigma;
294 double hYSigma = h * _ySigma;
295 double xhmin = hXSigma * sqrt(2.)/10; //smallest anticipated bandwidth
296 double yhmin = hYSigma * sqrt(2.)/10;
297
298 //////////////////////////////////////
299 //calculate bandwidths from the data//
300 //////////////////////////////////////
301 if(_BandWidthType == 1) //calculate a trivial bandwidth
302 {
303 cout << "Roo2DKeysPdf::calculateBandWidth Using a normal bandwidth (same for a given dimension) based on"<<endl;
304 cout << " h_j = n^{-1/6}*sigma_j for the j^th dimension and n events * "<<_widthScaleFactor<<endl;
305 double hxGaussian = _n16 * _xSigma * _widthScaleFactor;
306 double hyGaussian = _n16 * _ySigma * _widthScaleFactor;
307 for(Int_t j=0;j<_nEvents;++j)
308 {
309 _hx[j] = hxGaussian;
310 _hy[j] = hyGaussian;
311 if(_hx[j]<xhmin) _hx[j] = xhmin;
312 if(_hy[j]<yhmin) _hy[j] = yhmin;
313 }
314 }
315 else //use an adaptive bandwidth to reduce the dependence on global data distribution
316 {
317 cout << "Roo2DKeysPdf::calculateBandWidth Using an adaptive bandwidth (in general different for all events) [default]"<<endl;
318 cout << " scaled by a factor of "<<_widthScaleFactor<<endl;
319 double xnorm = h * std::pow(_xSigma/sqrtSum, 1.5) * _widthScaleFactor;
320 double ynorm = h * std::pow(_ySigma/sqrtSum, 1.5) * _widthScaleFactor;
321 for(Int_t j=0;j<_nEvents;++j)
322 {
323 double f_ti = std::pow( g(_x[j], _x, hXSigma, _y[j], _y, hYSigma), -0.25 ) ;
324 _hx[j] = xnorm * f_ti;
325 _hy[j] = ynorm * f_ti;
326 if(_hx[j]<xhmin) _hx[j] = xhmin;
327 if(_hy[j]<yhmin) _hy[j] = yhmin;
328 }
329 }
330
331 return 0;
332}
333
334
335////////////////////////////////////////////////////////////////////////////////
336/// Evaluates the kernel estimation for x,y, interpolating between the points if necessary
337///
338/// Uses the caching intrinsic in RFC to bypass the grid and remove
339/// the grid and extrapolation approximation in the kernel estimation method
340/// implementation.
341
343{
344 if(_vverbosedebug) { cout << "Roo2DKeysPdf::evaluate()" << endl; }
345 return evaluateFull(x,y);
346}
347
348
349////////////////////////////////////////////////////////////////////////////////
350/// Evaluates the sum of the product of the 2D kernels
351/// for use in calculating the fixed kernel estimate, f,
352/// given the bandwidths _hx[j] and _hy[j].
353///
354/// _n is calculated once in the constructor.
355/// \param[in] thisX
356/// \param[in] thisY
357
358double Roo2DKeysPdf::evaluateFull(double thisX, double thisY) const
359{
360 if( _vverbosedebug ) { cout << "Roo2DKeysPdf::evaluateFull()" << endl; }
361
362 double f=0.0;
363
364 double rx2;
365 double ry2;
366 double zx;
367 double 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
414double Roo2DKeysPdf::highBoundaryCorrection(double thisVar, double thisH, double high, double tVar) const
415{
416 if(_vverbosedebug) { cout << "Roo2DKeysPdf::highBoundaryCorrection" << endl; }
417
418 if(thisH == 0.0) return 0.0;
419 double correction = (thisVar + tVar - 2.0* high )/thisH;
420 return exp(-0.5*correction*correction)/thisH;
421}
422
423
424////////////////////////////////////////////////////////////////////////////////
425
426double Roo2DKeysPdf::lowBoundaryCorrection(double thisVar, double thisH, double low, double tVar) const
427{
428 if(_vverbosedebug) { cout << "Roo2DKeysPdf::lowBoundaryCorrection" << endl; }
429
430 if(thisH == 0.0) return 0.0;
431 double 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 Roo2DKeysPdf::g(double varMean1, double * _var1, double sigma1, double varMean2, double * _var2, double sigma2) const
447{
448 if((_nEvents == 0.0) || (sigma1 == 0.0) || (sigma2 == 0)) return 0.0;
449
450 double c1 = -1.0/(2.0*sigma1*sigma1);
451 double c2 = -1.0/(2.0*sigma2*sigma2);
452 double d = 4.0*c1*c2 /(_sqrt2pi*_nEvents);
453 double z = 0.0;
454
455 for (Int_t i = 0; i < _nEvents; ++i)
456 {
457 double r1 = _var1[i] - varMean1;
458 double 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 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 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 cout << "Roo2DKeysPdf::writeHistToFile This member function is temporarily disabled" <<endl;
528 //make sure that any existing file is not over written
529 std::unique_ptr<TFile> file{TFile::Open(outputFile, "UPDATE")};
530 if (!file)
531 {
532 cout << "Roo2DKeysPdf::writeHistToFile unable to open file "<< outputFile <<endl;
533 return;
534 }
535
536
537 const RooAbsReal & xx = x.arg();
538 const RooAbsReal & yy = y.arg();
539 RooArgSet values( RooArgList( xx, yy ));
540 RooRealVar * xArg = (static_cast<RooRealVar*>(values.find(xx.GetName())) ) ;
541 RooRealVar * yArg = (static_cast<RooRealVar*>(values.find(yy.GetName())) ) ;
542
543 TH2F * hist = (TH2F*)xArg->createHistogram("hist", *yArg);
544 hist = static_cast<TH2F*>(this->fillHistogram(hist, RooArgList(*xArg, *yArg) ));
545 hist->SetName(histName);
546
547 file->Write();
548 file->Close();
549}
550
551
552////////////////////////////////////////////////////////////////////////////////
553/// Saves the data and calculated bandwidths to a file,
554/// as a record of what produced the PDF and to give a reduced
555/// data set in order to facilitate re-calculation in the future.
556/// \param[in] outputFile Name of the file where to store the data
557/// \param[in] name Name of the tree which will contain the data
558
559void Roo2DKeysPdf::writeNTupleToFile(char * outputFile, const char * name) const
560{
561 TFile * file = nullptr;
562
563 //make sure that any existing file is not over written
564 file = new TFile(outputFile, "UPDATE");
565 if (!file)
566 {
567 cout << "Roo2DKeysPdf::writeNTupleToFile unable to open file "<< outputFile <<endl;
568 return;
569 }
570 RooAbsReal & xArg = (RooAbsReal&)x.arg();
571 RooAbsReal & yArg = (RooAbsReal&)y.arg();
572
573 double theX;
574 double theY;
575 double hx;
576 TString label = name;
577 label += " the source data for 2D Keys PDF";
578 TTree * _theTree = new TTree(name, label);
579 if(!_theTree) { cout << "Unable to get a TTree for output" << endl; return; }
580 _theTree->SetAutoSave(1000000000); // autosave when 1 Gbyte written
581
582 //name the TBranches the same as the RooAbsReal's
583 const char * xname = xArg.GetName();
584 const char * yname = yArg.GetName();
585 if (!strcmp(xname,"")) xname = "x";
586 if (!strcmp(yname,"")) yname = "y";
587
588 _theTree->Branch(xname, &theX, " x/D");
589 _theTree->Branch(yname, &theY, " y/D");
590 _theTree->Branch("hx", &hx, " hx/D");
591 _theTree->Branch("hy", &hx, " hy/D");
592
593 for(Int_t iEvt = 0; iEvt < _nEvents; iEvt++)
594 {
595 theX = _x[iEvt];
596 theY = _y[iEvt];
597 hx = _hx[iEvt];
598 hx = _hy[iEvt];
599 _theTree->Fill();
600 }
601 file->Write();
602 file->Close();
603}
604
605
606////////////////////////////////////////////////////////////////////////////////
607/// Prints out _p[_nPoints][_nPoints] indicating the domain limits.
608/// \param[out] out Output stream where to print
609
610void Roo2DKeysPdf::PrintInfo(ostream & out) const
611{
612 out << "Roo2DKeysPDF instance domain information:"<<endl;
613 out << "\tX_min = " << _lox <<endl;
614 out << "\tX_max = " << _hix <<endl;
615 out << "\tY_min = " << _loy <<endl;
616 out << "\tY_max = " << _hiy <<endl;
617
618 out << "Data information:" << endl;
619 out << "\t<x> = " << _xMean <<endl;
620 out << "\tsigma(x) = " << _xSigma <<endl;
621 out << "\t<y> = " << _yMean <<endl;
622 out << "\tsigma(y) = " << _ySigma <<endl;
623
624 out << "END of info for Roo2DKeys pdf instance"<< endl;
625}
#define d(i)
Definition RSha256.hxx:102
#define f(i)
Definition RSha256.hxx:104
#define g(i)
Definition RSha256.hxx:105
#define h(i)
Definition RSha256.hxx:106
int Int_t
Definition RtypesCore.h:45
#define ClassImp(name)
Definition Rtypes.h:382
#define X(type, name)
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y1
char name[80]
Definition TGX11.cxx:110
Two-dimensional kernel estimation PDF.
RooRealProxy y
double getMean(const char *axis) const
Roo2DKeysPdf(const char *name, const char *title, RooAbsReal &xx, RooAbsReal &yy, RooDataSet &data, TString options="a", double widthScaleFactor=1.0)
Constructor.
double _widthScaleFactor
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
void writeToFile(char *outputFile, const char *name) const
Int_t _verbosedebug
double * _hy
double lowBoundaryCorrection(double thisVar, double thisH, double low, double tVar) const
double evaluate() const override
Evaluates the kernel estimation for x,y, interpolating between the points if necessary.
Int_t getBandWidthType() const
Int_t _MirrorAtBoundary
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 g(double var1, double *_var1, double sigma1, double var2, double *_var2, double sigma2) const
Calculates f(t_i) for the bandwidths.
Int_t _vverbosedebug
double * _x
double * _y
double * _hx
double getSigma(const char *axis) const
~Roo2DKeysPdf() override
Destructor.
void setWidthScaleFactor(double widthScaleFactor)
double highBoundaryCorrection(double thisVar, double thisH, double high, double tVar) const
Apply the mirror at boundary correction to a dimension given the space position to evaluate at (thisV...
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].
void setOptions(TString options)
void getOptions(void) const
double evaluateFull(double thisX, double thisY) const
Evaluates the sum of the product of the 2D kernels for use in calculating the fixed kernel estimate,...
RooRealProxy x
RooAbsArg * find(const char *name) const
Find object with given name in list.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
TH1 * createHistogram(const char *name, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}) const
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
TH1 * fillHistogram(TH1 *hist, const RooArgList &plotVars, double scaleFactor=1, const RooArgSet *projectedVars=nullptr, bool scaling=true, const RooArgSet *condObs=nullptr, bool setError=true) const
Fill the ROOT histogram 'hist' with values sampled from this function at the bin centers.
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Container class to hold unbinned data.
Definition RooDataSet.h:34
Variable that can be changed from the outside.
Definition RooRealVar.h:37
double max(const char *rname=nullptr) 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.
double min(const char *rname=nullptr) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
Definition TFile.h:53
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
Definition TFile.cxx:4086
Int_t Write(const char *name=nullptr, Int_t opt=0, Int_t bufsiz=0) override
Write memory objects to this file.
Definition TFile.cxx:2433
void Close(Option_t *option="") override
Close a file.
Definition TFile.cxx:947
void SetName(const char *name) override
Change the name of this histogram.
Definition TH1.cxx:8993
2-D histogram with a float per channel (see TH1 documentation)
Definition TH2.h:308
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
Basic string class.
Definition TString.h:139
void ToLower()
Change string to lower-case.
Definition TString.cxx:1182
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:632
A TTree represents a columnar dataset.
Definition TTree.h:79
virtual Int_t Fill()
Fill all branches.
Definition TTree.cxx:4593
virtual void SetAutoSave(Long64_t autos=-300000000)
In case of a program crash, it will be possible to recover the data in the tree up to the last AutoSa...
Definition TTree.cxx:8332
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:353
Double_t y[n]
Definition legend1.C:17
return c1
Definition legend1.C:41
Double_t x[n]
Definition legend1.C:17
return c2
Definition legend2.C:14
constexpr Double_t Pi()
Definition TMath.h:37