Logo ROOT   6.12/07
Reference Guide
RooNDKeysPdf.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitModels *
4  * File: $Id: RooNDKeysPdf.cxx 31258 2009-11-17 22:41:06Z wouter $
5  * Authors: *
6  * Max Baak, CERN, mbaak@cern.ch *
7  * *
8  * Redistribution and use in source and binary forms, *
9  * with or without modification, are permitted according to the terms *
10  * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
11  *****************************************************************************/
12 
13 /** \class RooNDKeysPdf
14  \ingroup Roofit
15 
16 Generic N-dimensional implementation of a kernel estimation p.d.f. This p.d.f. models the distribution
17 of an arbitrary input dataset as a superposition of Gaussian kernels, one for each data point,
18 each contributing 1/N to the total integral of the p.d.f.
19 If the 'adaptive mode' is enabled, the width of the Gaussian is adaptively calculated from the
20 local density of events, i.e. narrow for regions with high event density to preserve details and
21 wide for regions with log event density to promote smoothness. The details of the general algorithm
22 are described in the following paper:
23 Cranmer KS, Kernel Estimation in High-Energy Physics.
24  Computer Physics Communications 136:198-207,2001 - e-Print Archive: hep ex/0011057
25 For multi-dimensional datasets, the kernels are modeled by multidimensional Gaussians. The kernels are
26 constructed such that they reflect the correlation coefficients between the observables
27 in the input dataset.
28 **/
29 
30 #include <iostream>
31 #include <algorithm>
32 #include <string>
33 
34 #include "TMath.h"
35 #include "TMatrixDSymEigen.h"
36 #include "RooNDKeysPdf.h"
37 #include "RooAbsReal.h"
38 #include "RooRealVar.h"
39 #include "RooRandom.h"
40 #include "RooHist.h"
41 #include "RooMsgService.h"
42 #include "RooChangeTracker.h"
43 
44 #include "TError.h"
45 
46 using namespace std;
47 
49 
50 ////////////////////////////////////////////////////////////////////////////////
51 /// Construct N-dimensional kernel estimation p.d.f. in observables 'varList'
52 /// from dataset 'data'. Options can be
53 ///
54 /// - 'a' = Use adaptive kernels (width varies with local event density)
55 /// - 'm' = Mirror data points over observable boundaries. Improves modeling
56 /// behavior at edges for distributions that are not close to zero
57 /// at edge
58 /// - 'd' = Debug flag
59 /// - 'v' = Verbose flag
60 ///
61 /// The parameter rho (default = 1) provides an overall scale factor that can
62 /// be applied to the bandwith calculated for each kernel. The nSigma parameter
63 /// determines the size of the box that is used to search for contributing kernels
64 /// around a given point in observable space. The nSigma parameters is used
65 /// in case of non-adaptive bandwidths and for the 1st non-adaptive pass for
66 /// the calculation of adaptive keys p.d.f.s.
67 ///
68 /// The optional weight arguments allows to specify an observable or function
69 /// expression in observables that specifies the weight of each event.
70 
71 RooNDKeysPdf::RooNDKeysPdf(const char *name, const char *title, const RooArgList &varList, const RooAbsData &data,
72  TString options, Double_t rho, Double_t nSigma, Bool_t rotate, Bool_t sortInput)
73  : RooAbsPdf(name, title), _varList("varList", "List of variables", this),
74  _rhoList("rhoList", "List of rho parameters", this), _dataP(0), _data(data), _options(options), _widthFactor(rho),
75  _nSigma(nSigma), _weights(&_weights0), _rotate(rotate), _sortInput(sortInput), _nAdpt(1), _tracker(0)
76 {
77  // Constructor
80 
81  TIterator* varItr = varList.createIterator() ;
82  RooAbsArg* var ;
83  for (Int_t i=0; (var = (RooAbsArg*)varItr->Next()); ++i) {
84  if (!dynamic_cast<RooAbsReal*>(var)) {
85  coutE(InputArguments) << "RooNDKeysPdf::ctor(" << GetName() << ") ERROR: variable " << var->GetName()
86  << " is not of type RooAbsReal" << endl ;
87  R__ASSERT(0) ;
88  }
89  _varList.add(*var) ;
90  _varName.push_back(var->GetName());
91  }
92  delete varItr ;
93 
94  createPdf();
95 }
96 
97 ////////////////////////////////////////////////////////////////////////////////
98 /// Constructor
99 
100 RooNDKeysPdf::RooNDKeysPdf(const char *name, const char *title, const RooArgList &varList, const TH1 &hist,
101  TString options, Double_t rho, Double_t nSigma, Bool_t rotate, Bool_t sortInput)
102  : RooAbsPdf(name, title), _varList("varList", "List of variables", this),
103  _rhoList("rhoList", "List of rho parameters", this), _dataP(createDatasetFromHist(varList, hist)), _data(*_dataP),
104  _options(options), _widthFactor(rho), _nSigma(nSigma), _weights(&_weights0), _rotate(rotate),
105  _sortInput(sortInput), _nAdpt(1), _tracker(0)
106 {
109 
110  TIterator *varItr = varList.createIterator();
111  RooAbsArg *var;
112  for (Int_t i = 0; (var = (RooAbsArg *)varItr->Next()); ++i) {
113  if (!dynamic_cast<RooAbsReal *>(var)) {
114  coutE(InputArguments) << "RooNDKeysPdf::ctor(" << GetName() << ") ERROR: variable " << var->GetName()
115  << " is not of type RooAbsReal" << endl;
116  assert(0);
117  }
118  _varList.add(*var);
119  _varName.push_back(var->GetName());
120  }
121  delete varItr;
122 
123  createPdf();
124 }
125 
126 ////////////////////////////////////////////////////////////////////////////////
127 /// Constructor
128 
129 RooNDKeysPdf::RooNDKeysPdf(const char *name, const char *title, const RooArgList &varList, const RooAbsData &data,
130  const TVectorD &rho, TString options, Double_t nSigma, Bool_t rotate, Bool_t sortInput)
131  : RooAbsPdf(name, title), _varList("varList", "List of variables", this),
132  _rhoList("rhoList", "List of rho parameters", this), _dataP(0), _data(data), _options(options), _widthFactor(-1.0),
133  _nSigma(nSigma), _weights(&_weights0), _rotate(rotate), _sortInput(sortInput), _nAdpt(1), _tracker(0)
134 {
137 
138  TIterator* varItr = varList.createIterator() ;
139  RooAbsArg* var ;
140  for (Int_t i=0; (var = (RooAbsArg*)varItr->Next()); ++i) {
141  if (!dynamic_cast<RooAbsReal*>(var)) {
142  coutE(InputArguments) << "RooNDKeysPdf::ctor(" << GetName() << ") ERROR: variable " << var->GetName()
143  << " is not of type RooAbsReal" << endl;
144  R__ASSERT(0);
145  }
146  _varList.add(*var) ;
147  _varName.push_back(var->GetName());
148  }
149  delete varItr ;
150 
151  // copy rho widths
152  if( _varList.getSize() != rho.GetNrows() ) {
154  << "ERROR: RooNDKeysPdf::RooNDKeysPdf() : The vector-size of rho is different from that of varList."
155  << "Unable to create the PDF." << endl;
156  R__ASSERT(_varList.getSize() == rho.GetNrows());
157  }
158 
159  // negative width factor will serve as a switch in initialize()
160  // negative value means that a vector has been provided as input,
161  // and that _rho has already been set ...
162  _rho.resize( rho.GetNrows() );
163  for (Int_t j = 0; j < rho.GetNrows(); j++) {
164  _rho[j] = rho[j]; /*cout<<"RooNDKeysPdf ctor, _rho["<<j<<"]="<<_rho[j]<<endl;*/
165  }
166 
167  createPdf(); // calls initialize ...
168 }
169 
170 ////////////////////////////////////////////////////////////////////////////////
171 /// Backward compatibility constructor for (1-dim) RooKeysPdf. If you are a new user,
172 /// please use the first constructor form.
173 
174 RooNDKeysPdf::RooNDKeysPdf(const char *name, const char *title, const RooArgList &varList, const RooAbsData &data,
175  const RooArgList &rhoList, TString options, Double_t nSigma, Bool_t rotate, Bool_t sortInput)
176  : RooAbsPdf(name, title), _varList("varList", "List of variables", this),
177  _rhoList("rhoList", "List of rho parameters", this), _dataP(0), _data(data), _options(options), _widthFactor(-1.0),
178  _nSigma(nSigma), _weights(&_weights0), _rotate(rotate), _sortInput(sortInput), _nAdpt(1)
179 {
182 
183  TIterator *varItr = varList.createIterator();
184  RooAbsArg *var;
185  for (Int_t i = 0; (var = (RooAbsArg *)varItr->Next()); ++i) {
186  if (!dynamic_cast<RooAbsReal *>(var)) {
187  coutE(InputArguments) << "RooNDKeysPdf::ctor(" << GetName() << ") ERROR: variable " << var->GetName()
188  << " is not of type RooAbsReal" << endl;
189  assert(0);
190  }
191  _varList.add(*var);
192  _varName.push_back(var->GetName());
193  }
194  delete varItr;
195 
196  TIterator *rhoItr = rhoList.createIterator();
197  RooAbsArg *rho;
198  _rho.resize(rhoList.getSize(), 1.0);
199 
200  for (Int_t i = 0; (rho = (RooAbsArg *)rhoItr->Next()); ++i) {
201  if (!dynamic_cast<RooAbsReal *>(rho)) {
202  coutE(InputArguments) << "RooNDKeysPdf::ctor(" << GetName() << ") ERROR: parameter " << rho->GetName()
203  << " is not of type RooRealVar" << endl;
204  assert(0);
205  }
206  _rhoList.add(*rho);
207  _rho[i] = (dynamic_cast<RooAbsReal *>(rho))->getVal();
208  }
209  delete rhoItr;
210 
211  // copy rho widths
212  if ((_varList.getSize() != _rhoList.getSize())) {
213  coutE(InputArguments) << "ERROR: RooNDKeysPdf::RooNDKeysPdf() : The size of rhoList is different from varList."
214  << "Unable to create the PDF." << endl;
215  assert(_varList.getSize() == _rhoList.getSize());
216  }
217 
218  // keep track of changes in rho parameters
219  _tracker = new RooChangeTracker("tracker", "track rho parameters", _rhoList, true); // check for value updates
220  (void)_tracker->hasChanged(true); // first evaluation always true for new parameters (?)
221 
222  createPdf();
223 }
224 
225 ////////////////////////////////////////////////////////////////////////////////
226 /// Constructor
227 
228 RooNDKeysPdf::RooNDKeysPdf(const char *name, const char *title, const RooArgList &varList, const TH1 &hist,
229  const RooArgList &rhoList, TString options, Double_t nSigma, Bool_t rotate, Bool_t sortInput)
230  : RooAbsPdf(name, title), _varList("varList", "List of variables", this),
231  _rhoList("rhoList", "List of rho parameters", this), _dataP(createDatasetFromHist(varList, hist)), _data(*_dataP),
232  _options(options), _widthFactor(-1), _nSigma(nSigma), _weights(&_weights0), _rotate(rotate), _sortInput(sortInput),
233  _nAdpt(1)
234 {
237 
238  TIterator *varItr = varList.createIterator();
239  RooAbsArg *var;
240  for (Int_t i = 0; (var = (RooAbsArg *)varItr->Next()); ++i) {
241  if (!dynamic_cast<RooAbsReal *>(var)) {
242  coutE(InputArguments) << "RooNDKeysPdf::ctor(" << GetName() << ") ERROR: variable " << var->GetName()
243  << " is not of type RooAbsReal" << endl;
244  assert(0);
245  }
246  _varList.add(*var);
247  _varName.push_back(var->GetName());
248  }
249  delete varItr;
250 
251  // copy rho widths
252  TIterator *rhoItr = rhoList.createIterator();
253  RooAbsArg *rho;
254  _rho.resize(rhoList.getSize(), 1.0);
255 
256  for (Int_t i = 0; (rho = (RooAbsArg *)rhoItr->Next()); ++i) {
257  if (!dynamic_cast<RooAbsReal *>(rho)) {
258  coutE(InputArguments) << "RooNDKeysPdf::ctor(" << GetName() << ") ERROR: parameter " << rho->GetName()
259  << " is not of type RooRealVar" << endl;
260  assert(0);
261  }
262  _rhoList.add(*rho);
263  _rho[i] = (dynamic_cast<RooAbsReal *>(rho))->getVal();
264  }
265  delete rhoItr;
266 
267  if ((_varList.getSize() != _rhoList.getSize())) {
268  coutE(InputArguments) << "ERROR: RooNDKeysPdf::RooNDKeysPdf() : The size of rhoList is different from varList."
269  << "Unable to create the PDF." << endl;
270  assert(_varList.getSize() == _rhoList.getSize());
271  }
272 
273  // keep track of changes in rho parameters
274  _tracker = new RooChangeTracker("tracker", "track rho parameters", _rhoList, true); // check for value updates
275  (void)_tracker->hasChanged(true); // first evaluation always true for new parameters (?)
276 
277  createPdf();
278 }
279 
280 ////////////////////////////////////////////////////////////////////////////////
281 /// Constructor
282 
283 RooNDKeysPdf::RooNDKeysPdf(const char *name, const char *title, RooAbsReal &x, const RooAbsData &data, Mirror mirror,
284  Double_t rho, Double_t nSigma, Bool_t rotate, Bool_t sortInput)
285  : RooAbsPdf(name, title), _varList("varList", "List of variables", this),
286  _rhoList("rhoList", "List of rho parameters", this), _dataP(0), _data(data), _options("a"), _widthFactor(rho),
287  _nSigma(nSigma), _weights(&_weights0), _rotate(rotate), _sortInput(sortInput), _nAdpt(1), _tracker(0)
288 {
291 
292  _varList.add(x);
293  _varName.push_back(x.GetName());
294 
295  if (mirror != NoMirror) {
296  if (mirror != MirrorBoth)
297  coutW(InputArguments) << "RooNDKeysPdf::RooNDKeysPdf() : Warning : asymmetric mirror(s) no longer supported."
298  << endl;
299  _options = "m";
300  }
301 
302  createPdf();
303 }
304 
305 ////////////////////////////////////////////////////////////////////////////////
306 /// Backward compatibility constructor for Roo2DKeysPdf. If you are a new user,
307 /// please use the first constructor form.
308 
309 RooNDKeysPdf::RooNDKeysPdf(const char *name, const char *title, RooAbsReal &x, RooAbsReal &y, const RooAbsData &data,
310  TString options, Double_t rho, Double_t nSigma, Bool_t rotate, Bool_t sortInput)
311  : RooAbsPdf(name, title), _varList("varList", "List of variables", this),
312  _rhoList("rhoList", "List of rho parameters", this), _dataP(0), _data(data), _options(options), _widthFactor(rho),
313  _nSigma(nSigma), _weights(&_weights0), _rotate(rotate), _sortInput(sortInput), _nAdpt(1), _tracker(0)
314 {
317 
318  _varList.add(RooArgSet(x, y));
319  _varName.push_back(x.GetName());
320  _varName.push_back(y.GetName());
321 
322  createPdf();
323 }
324 
325 ////////////////////////////////////////////////////////////////////////////////
326 /// Constructor
327 
328 RooNDKeysPdf::RooNDKeysPdf(const RooNDKeysPdf &other, const char *name)
329  : RooAbsPdf(other, name), _varList("varList", this, other._varList), _rhoList("rhoList", this, other._rhoList),
330  _dataP(other._dataP != NULL ? new RooDataSet(*other._dataP) : NULL),
331  _data(other._dataP != NULL ? *_dataP : other._data), _options(other._options), _widthFactor(other._widthFactor),
333  _nAdpt(other._nAdpt)
334 {
335  _tracker = (other._tracker != NULL ? new RooChangeTracker(*other._tracker) : NULL);
336  // if (_tracker!=NULL) { _tracker->hasChanged(true); }
337 
340 
341  _fixedShape = other._fixedShape;
342  _mirror = other._mirror;
343  _debug = other._debug;
344  _verbose = other._verbose;
345  _sqrt2pi = other._sqrt2pi;
346  _nDim = other._nDim;
347  _nEvents = other._nEvents;
348  _nEventsM = other._nEventsM;
349  _nEventsW = other._nEventsW;
350  _d = other._d;
351  _n = other._n;
352  _dataPts = other._dataPts;
353  _dataPtsR = other._dataPtsR;
354  _weights0 = other._weights0;
355  _weights1 = other._weights1;
356  if (_options.Contains("a")) {
357  _weights = &_weights1;
358  }
359  //_sortIdcs = other._sortIdcs;
360  _sortTVIdcs = other._sortTVIdcs;
361  _varName = other._varName;
362  _rho = other._rho;
363  _x = other._x;
364  _x0 = other._x0;
365  _x1 = other._x1;
366  _x2 = other._x2;
367  _xDatLo = other._xDatLo;
368  _xDatHi = other._xDatHi;
369  _xDatLo3s = other._xDatLo3s;
370  _xDatHi3s = other._xDatHi3s;
371  _mean = other._mean;
372  _sigma = other._sigma;
373 
374  // BoxInfo
375  _netFluxZ = other._netFluxZ;
376  _nEventsBW = other._nEventsBW;
377  _nEventsBMSW = other._nEventsBMSW;
378  _xVarLo = other._xVarLo;
379  _xVarHi = other._xVarHi;
380  _xVarLoM3s = other._xVarLoM3s;
381  _xVarLoP3s = other._xVarLoP3s;
382  _xVarHiM3s = other._xVarHiM3s;
383  _xVarHiP3s = other._xVarHiP3s;
384  _bpsIdcs = other._bpsIdcs;
385  _ibNoSort = other._ibNoSort;
386  _sIdcs = other._sIdcs;
387  _bIdcs = other._bIdcs;
388  _bmsIdcs = other._bmsIdcs;
389 
391  _fullBoxInfo = other._fullBoxInfo;
392 
393  _idx = other._idx;
394  _minWeight = other._minWeight;
395  _maxWeight = other._maxWeight;
396  _wMap = other._wMap;
397 
398  _covMat = new TMatrixDSym(*other._covMat);
399  _corrMat = new TMatrixDSym(*other._corrMat);
400  _rotMat = new TMatrixD(*other._rotMat);
401  _sigmaR = new TVectorD(*other._sigmaR);
402  _dx = new TVectorD(*other._dx);
403  _sigmaAvgR = other._sigmaAvgR;
404 }
405 
406 ////////////////////////////////////////////////////////////////////////////////
407 
409 {
410  if (_varItr) delete _varItr;
411  if (_rhoItr)
412  delete _rhoItr;
413  if (_covMat) delete _covMat;
414  if (_corrMat) delete _corrMat;
415  if (_rotMat) delete _rotMat;
416  if (_sigmaR) delete _sigmaR;
417  if (_dx) delete _dx;
418  if (_dataP)
419  delete _dataP;
420  if (_tracker)
421  delete _tracker;
422 
423  // delete all the boxinfos map
424  while ( !_rangeBoxInfo.empty() ) {
425  map<pair<string,int>,BoxInfo*>::iterator iter = _rangeBoxInfo.begin();
426  BoxInfo* box= (*iter).second;
427  _rangeBoxInfo.erase(iter);
428  delete box;
429  }
430 
431  _dataPts.clear();
432  _dataPtsR.clear();
433  _weights0.clear();
434  _weights1.clear();
435  //_sortIdcs.clear();
436  _sortTVIdcs.clear();
437 }
438 
439 ////////////////////////////////////////////////////////////////////////////////
440 /// evaluation order of constructor.
441 
442 void RooNDKeysPdf::createPdf(Bool_t firstCall) const
443 {
444  if (firstCall) {
445  // set options
446  setOptions();
447  // initialization
448  initialize();
449  }
450 
451 
452  // copy dataset, calculate sigma_i's, determine min and max event weight
453  loadDataSet(firstCall);
454 
455  // mirror dataset around dataset boundaries -- does not depend on event weights
456  if (_mirror) mirrorDataSet();
457 
458  // store indices and weights of events with high enough weights
459  loadWeightSet();
460 
461  // store indices of events in variable boundaries and box shell.
462 //calculateShell(&_fullBoxInfo);
463  // calculate normalization needed in analyticalIntegral()
464 //calculatePreNorm(&_fullBoxInfo);
465 
466  // lookup table for determining which events contribute to a certain coordinate
467  sortDataIndices();
468 
469  // determine static and/or adaptive bandwidth
471 }
472 
473 ////////////////////////////////////////////////////////////////////////////////
474 /// set the configuration
475 
477 {
478  _options.ToLower();
479 
480  if( _options.Contains("a") ) { _weights = &_weights1; }
481  else { _weights = &_weights0; }
482  if( _options.Contains("m") ) _mirror = true;
483  else _mirror = false;
484  if( _options.Contains("d") ) _debug = true;
485  else _debug = false;
486  if( _options.Contains("v") ) { _debug = true; _verbose = true; }
487  else { _debug = false; _verbose = false; }
488 
489  cxcoutD(InputArguments) << "RooNDKeysPdf::setOptions() options = " << _options
490  << "\n\tbandWidthType = " << _options.Contains("a")
491  << "\n\tmirror = " << _mirror
492  << "\n\tdebug = " << _debug
493  << "\n\tverbose = " << _verbose
494  << endl;
495 
496  if (_nSigma<2.0) {
497  coutW(InputArguments) << "RooNDKeysPdf::setOptions() : Warning : nSigma = " << _nSigma << " < 2.0. "
498  << "Calculated normalization could be too large."
499  << endl;
500  }
501 
502  // number of adaptive width iterations. Default is 1.
503  if (_options.Contains("a")) {
504  if (!sscanf(_options.Data(), "%d%*s", &_nAdpt)) {
505  _nAdpt = 1;
506  }
507  }
508 }
509 
510 ////////////////////////////////////////////////////////////////////////////////
511 /// initialization
512 
514 {
515  _sqrt2pi = sqrt(2.0*TMath::Pi()) ;
516  _nDim = _varList.getSize();
520 
521  _netFluxZ = kFALSE;
522  _nEventsBW = 0;
523  _nEventsBMSW = 0;
524 
525  if(_nDim==0) {
526  coutE(InputArguments) << "ERROR: RooNDKeysPdf::initialize() : The observable list is empty. "
527  << "Unable to begin generating the PDF." << endl;
528  R__ASSERT (_nDim!=0);
529  }
530 
531  if(_nEvents==0) {
532  coutE(InputArguments) << "ERROR: RooNDKeysPdf::initialize() : The input data set is empty. "
533  << "Unable to begin generating the PDF." << endl;
534  R__ASSERT (_nEvents!=0);
535  }
536 
537  _d = static_cast<Double_t>(_nDim);
538 
539  vector<Double_t> dummy(_nDim,0.);
540  _dataPts.resize(_nEvents,dummy);
541  _weights0.resize(_nEvents,dummy);
542  //_sortIdcs.resize(_nDim);
543  _sortTVIdcs.resize(_nDim);
544 
545  //rdh _rho.resize(_nDim,_widthFactor);
546 
547  if (_widthFactor>0) { _rho.resize(_nDim,_widthFactor); }
548  // else: _rho has been provided as external input
549 
550  _x.resize(_nDim,0.);
551  _x0.resize(_nDim,0.);
552  _x1.resize(_nDim,0.);
553  _x2.resize(_nDim,0.);
554 
555  _mean.resize(_nDim,0.);
556  _sigma.resize(_nDim,0.);
557 
558  _xDatLo.resize(_nDim,0.);
559  _xDatHi.resize(_nDim,0.);
560  _xDatLo3s.resize(_nDim,0.);
561  _xDatHi3s.resize(_nDim,0.);
562 
563  boxInfoInit(&_fullBoxInfo,0,0xFFFF);
564 
565  _minWeight=0;
566  _maxWeight=0;
567  _wMap.clear();
568 
569  _covMat = 0;
570  _corrMat= 0;
571  _rotMat = 0;
572  _sigmaR = 0;
573  _dx = new TVectorD(_nDim); _dx->Zero();
574  _dataPtsR.resize(_nEvents,*_dx);
575 
576  _varItr->Reset() ;
577  RooRealVar* var ;
578  for(Int_t j=0; (var=(RooRealVar*)_varItr->Next()); ++j) {
579  _xDatLo[j] = var->getMin();
580  _xDatHi[j] = var->getMax();
581  }
582 }
583 
584 ////////////////////////////////////////////////////////////////////////////////
585 /// copy the dataset and calculate some useful variables
586 
587 void RooNDKeysPdf::loadDataSet(Bool_t firstCall) const
588 {
589  // first some initialization
590  _nEventsW=0.;
591 
592  TMatrixD mat(_nDim,_nDim);
593  if (!_covMat) _covMat = new TMatrixDSym(_nDim);
594  if (!_corrMat) _corrMat= new TMatrixDSym(_nDim);
595  if (!_rotMat) _rotMat = new TMatrixD(_nDim,_nDim);
596  if (!_sigmaR) _sigmaR = new TVectorD(_nDim);
597 
598  mat.Zero();
599  _covMat->Zero();
600  _corrMat->Zero();
601  _rotMat->Zero();
602  _sigmaR->Zero();
603 
604  const RooArgSet* values= _data.get();
605  vector<RooRealVar*> dVars(_nDim);
606  for (Int_t j=0; j<_nDim; j++) {
607  dVars[j] = (RooRealVar*)values->find(_varName[j].c_str());
608  _x0[j]=_x1[j]=_x2[j]=0.;
609  }
610 
611  _idx.clear();
612  for (Int_t i=0; i<_nEvents; i++) {
613 
614  _data.get(i); // fills dVars
615  _idx.push_back(i);
616  vector<Double_t>& point = _dataPts[i];
617  TVectorD& pointV = _dataPtsR[i];
618 
619  Double_t myweight = _data.weight(); // default is one?
620  if ( TMath::Abs(myweight)>_maxWeight ) { _maxWeight = TMath::Abs(myweight); }
621  _nEventsW += myweight;
622 
623  for (Int_t j=0; j<_nDim; j++) {
624  for (Int_t k=0; k<_nDim; k++) {
625  mat(j,k) += dVars[j]->getVal() * dVars[k]->getVal() * myweight;
626  }
627  // only need to do once
628  if (firstCall)
629  point[j] = pointV[j] = dVars[j]->getVal();
630 
631  _x0[j] += 1. * myweight;
632  _x1[j] += point[j] * myweight ;
633  _x2[j] += point[j] * point[j] * myweight ;
634  if (_x2[j]!=_x2[j]) exit(3);
635 
636  // only need to do once
637  if (firstCall) {
638  if (point[j]<_xDatLo[j]) { _xDatLo[j]=point[j]; }
639  if (point[j]>_xDatHi[j]) { _xDatHi[j]=point[j]; }
640  }
641  }
642  }
643 
644  _n = TMath::Power(4./(_nEventsW*(_d+2.)), 1./(_d+4.)) ;
645  // = (4/[n(dim(R) + 2)])^1/(dim(R)+4); dim(R) = 2
646  _minWeight = (0.5 - TMath::Erf(_nSigma/sqrt(2.))/2.) * _maxWeight;
647 
648  for (Int_t j=0; j<_nDim; j++) {
649  _mean[j] = _x1[j]/_x0[j];
650  _sigma[j] = sqrt(_x2[j]/_x0[j]-_mean[j]*_mean[j]);
651  }
652 
653  for (Int_t j=0; j<_nDim; j++) {
654  for (Int_t k=0; k<_nDim; k++) {
655  (*_covMat)(j,k) = mat(j,k)/_x0[j] - _mean[j]*_mean[k];
656  }
657  }
658 
659  for (Int_t j=0; j<_nDim; j++) {
660  for (Int_t k=0; k<_nDim; k++)
661  (*_corrMat)(j,k) = (*_covMat)(j,k)/(_sigma[j]*_sigma[k]) ;
662  }
663 
664  // use raw sigmas (without rho) for sigmaAvgR
665  TMatrixDSymEigen evCalculator(*_covMat);
666  TVectorD sigmaRraw = evCalculator.GetEigenValues();
667  for (Int_t j=0; j<_nDim; j++) { sigmaRraw[j] = sqrt(sigmaRraw[j]); }
668 
669  _sigmaAvgR=1.;
670  for (Int_t j=0; j<_nDim; j++) { _sigmaAvgR *= sigmaRraw[j]; }
672 
673  // find decorrelation matrix and eigenvalues (R)
674  if (_nDim > 1 && _rotate) {
675  // new: rotation matrix now independent of rho evaluation
676  *_rotMat = evCalculator.GetEigenVectors();
677  *_rotMat = _rotMat->T(); // transpose
678  } else {
679  TMatrixD haar(_nDim, _nDim);
680  TMatrixD unit(TMatrixD::kUnit, haar);
681  *_rotMat = unit;
682  }
683 
684  // update sigmas (rho dependent)
685  updateRho();
686 
687  //// rho no longer used after this.
688  //// Now set rho = 1 because sigmaR now contains rho
689  // for (Int_t j=0; j<_nDim; j++) { _rho[j] = 1.; } // reset: important!
690 
691  if (_verbose) {
692  //_covMat->Print();
693  _rotMat->Print();
694  _corrMat->Print();
695  _sigmaR->Print();
696  }
697 
698  if (_nDim > 1 && _rotate) {
699  // apply rotation
700  for (Int_t i = 0; i < _nEvents; i++) {
701  TVectorD &pointR = _dataPtsR[i];
702  pointR *= *_rotMat;
703  }
704  }
705 
706  coutI(Contents) << "RooNDKeysPdf::loadDataSet(" << this << ")"
707  << "\n Number of events in dataset: " << _nEvents
708  << "\n Weighted number of events in dataset: " << _nEventsW << endl;
709 }
710 
711 ////////////////////////////////////////////////////////////////////////////////
712 /// determine mirror dataset.
713 /// mirror points are added around the physical boundaries of the dataset
714 /// Two steps:
715 /// 1. For each entry, determine if it should be mirrored (the mirror configuration).
716 /// 2. For each mirror configuration, make the mirror points.
717 
719 {
720  for (Int_t j=0; j<_nDim; j++) {
721  _xDatLo3s[j] = _xDatLo[j] + _nSigma * (_n * _sigma[j]);
722  _xDatHi3s[j] = _xDatHi[j] - _nSigma * (_n * _sigma[j]);
723 
724  // cout<<"xDatLo3s["<<j<<"]="<<_xDatLo3s[j]<<endl;
725  // cout<<"xDatHi3s["<<j<<"]="<<_xDatHi3s[j]<<endl;
726  }
727 
728  vector<Double_t> dummy(_nDim,0.);
729 
730  // 1.
731  for (Int_t i=0; i<_nEvents; i++) {
732  vector<Double_t>& x = _dataPts[i];
733 
734  Int_t size = 1;
735  vector<vector<Double_t> > mpoints(size,dummy);
736  vector<vector<Int_t> > mjdcs(size);
737 
738  // create all mirror configurations for event i
739  for (Int_t j=0; j<_nDim; j++) {
740 
741  vector<Int_t>& mjdxK = mjdcs[0];
742  vector<Double_t>& mpointK = mpoints[0];
743 
744  // single mirror *at physical boundaries*
745  if ((x[j]>_xDatLo[j] && x[j]<_xDatLo3s[j]) && x[j]<(_xDatLo[j]+_xDatHi[j])/2.) {
746  mpointK[j] = 2.*_xDatLo[j]-x[j];
747  mjdxK.push_back(j);
748  } else if ((x[j]<_xDatHi[j] && x[j]>_xDatHi3s[j]) && x[j]>(_xDatLo[j]+_xDatHi[j])/2.) {
749  mpointK[j] = 2.*_xDatHi[j]-x[j];
750  mjdxK.push_back(j);
751  }
752  }
753 
754  vector<Int_t>& mjdx0 = mjdcs[0];
755  // no mirror point(s) for this event
756  if (size==1 && mjdx0.size()==0) continue;
757 
758  // 2.
759  // generate all mirror points for event i
760  vector<Int_t>& mjdx = mjdcs[0];
761  vector<Double_t>& mpoint = mpoints[0];
762 
763  // number of mirror points for this mirror configuration
764  Int_t eMir = 1 << mjdx.size();
765  vector<vector<Double_t> > epoints(eMir,x);
766 
767  for (Int_t m=0; m<Int_t(mjdx.size()); m++) {
768  Int_t size1 = 1 << m;
769  Int_t size2 = 1 << (m+1);
770  // copy all previous mirror points
771  for (Int_t l=size1; l<size2; ++l) {
772  epoints[l] = epoints[l-size1];
773  // fill high mirror points
774  vector<Double_t>& epoint = epoints[l];
775  epoint[mjdx[Int_t(mjdx.size()-1)-m]] = mpoint[mjdx[Int_t(mjdx.size()-1)-m]];
776  }
777  }
778 
779  // remove duplicate mirror points
780  // note that: first epoint == x
781  epoints.erase(epoints.begin());
782 
783  // add mirror points of event i to total dataset
784  TVectorD pointR(_nDim);
785 
786  for (Int_t m=0; m<Int_t(epoints.size()); m++) {
787  _idx.push_back(i);
788  _dataPts.push_back(epoints[m]);
789  //_weights0.push_back(_weights0[i]);
790  for (Int_t j=0; j<_nDim; j++) { pointR[j] = (epoints[m])[j]; }
791  if (_nDim > 1 && _rotate) {
792  pointR *= *_rotMat;
793  }
794  _dataPtsR.push_back(pointR);
795  }
796 
797  epoints.clear();
798  mpoints.clear();
799  mjdcs.clear();
800  } // end of event loop
801 
802  _nEventsM = Int_t(_dataPts.size());
803 }
804 
805 ////////////////////////////////////////////////////////////////////////////////
806 
808 {
809  _wMap.clear();
810 
811  for (Int_t i=0; i<_nEventsM; i++) {
812  _data.get(_idx[i]);
813  Double_t myweight = _data.weight();
814  //if ( TMath::Abs(myweight)>_minWeight ) {
815  _wMap[i] = myweight;
816  //}
817  }
818 
819  coutI(Contents) << "RooNDKeysPdf::loadWeightSet(" << this << ") : Number of weighted events : " << _wMap.size() << endl;
820 }
821 
822 ////////////////////////////////////////////////////////////////////////////////
823 /// determine points in +/- nSigma shell around the box determined by the variable
824 /// ranges. These points are needed in the normalization, to determine probability
825 /// leakage in and out of the box.
826 
828 {
829  for (Int_t j=0; j<_nDim; j++) {
830  if (bi->xVarLo[j]==_xDatLo[j] && bi->xVarHi[j]==_xDatHi[j]) {
831  bi->netFluxZ = bi->netFluxZ && kTRUE;
832  } else { bi->netFluxZ = kFALSE; }
833 
834  bi->xVarLoM3s[j] = bi->xVarLo[j] - _nSigma * (_n * _sigma[j]);
835  bi->xVarLoP3s[j] = bi->xVarLo[j] + _nSigma * (_n * _sigma[j]);
836  bi->xVarHiM3s[j] = bi->xVarHi[j] - _nSigma * (_n * _sigma[j]);
837  bi->xVarHiP3s[j] = bi->xVarHi[j] + _nSigma * (_n * _sigma[j]);
838 
839  //cout<<"bi->xVarLoM3s["<<j<<"]="<<bi->xVarLoM3s[j]<<endl;
840  //cout<<"bi->xVarLoP3s["<<j<<"]="<<bi->xVarLoP3s[j]<<endl;
841  //cout<<"bi->xVarHiM3s["<<j<<"]="<<bi->xVarHiM3s[j]<<endl;
842  //cout<<"bi->xVarHiM3s["<<j<<"]="<<bi->xVarHiM3s[j]<<endl;
843  }
844 
845  map<Int_t,Double_t>::iterator wMapItr = _wMap.begin();
846 
847  //for (Int_t i=0; i<_nEventsM; i++) {
848  for (; wMapItr!=_wMap.end(); ++wMapItr) {
849  Int_t i = (*wMapItr).first;
850 
851  const vector<Double_t>& x = _dataPts[i];
852  Bool_t inVarRange(kTRUE);
853  Bool_t inVarRangePlusShell(kTRUE);
854 
855  for (Int_t j=0; j<_nDim; j++) {
856 
857  if (x[j]>bi->xVarLo[j] && x[j]<bi->xVarHi[j]) {
858  inVarRange = inVarRange && kTRUE;
859  } else { inVarRange = inVarRange && kFALSE; }
860 
861  if (x[j]>bi->xVarLoM3s[j] && x[j]<bi->xVarHiP3s[j]) {
862  inVarRangePlusShell = inVarRangePlusShell && kTRUE;
863  } else { inVarRangePlusShell = inVarRangePlusShell && kFALSE; }
864  }
865 
866  // event in range?
867  if (inVarRange) {
868  bi->bIdcs.push_back(i);
869  }
870 
871  // event in shell?
872  if (inVarRangePlusShell) {
873  bi->bpsIdcs[i] = kTRUE;
874  Bool_t inShell(kFALSE);
875  for (Int_t j=0; j<_nDim; j++) {
876  if ((x[j]>bi->xVarLoM3s[j] && x[j]<bi->xVarLoP3s[j]) && x[j]<(bi->xVarLo[j]+bi->xVarHi[j])/2.) {
877  inShell = kTRUE;
878  } else if ((x[j]>bi->xVarHiM3s[j] && x[j]<bi->xVarHiP3s[j]) && x[j]>(bi->xVarLo[j]+bi->xVarHi[j])/2.) {
879  inShell = kTRUE;
880  }
881  }
882  if (inShell) bi->sIdcs.push_back(i); // needed for normalization
883  else {
884  bi->bmsIdcs.push_back(i); // idem
885  }
886  }
887  }
888 
889  coutI(Contents) << "RooNDKeysPdf::calculateShell() : "
890  << "\n Events in shell " << bi->sIdcs.size()
891  << "\n Events in box " << bi->bIdcs.size()
892  << "\n Events in box and shell " << bi->bpsIdcs.size()
893  << endl;
894 }
895 
896 ////////////////////////////////////////////////////////////////////////////////
897 ///bi->nEventsBMSW=0.;
898 ///bi->nEventsBW=0.;
899 
901 {
902  // box minus shell
903  for (Int_t i=0; i<Int_t(bi->bmsIdcs.size()); i++)
904  bi->nEventsBMSW += _wMap[bi->bmsIdcs[i]];
905 
906  // box
907  for (Int_t i=0; i<Int_t(bi->bIdcs.size()); i++)
908  bi->nEventsBW += _wMap[bi->bIdcs[i]];
909 
910  cxcoutD(Eval) << "RooNDKeysPdf::calculatePreNorm() : "
911  << "\n nEventsBMSW " << bi->nEventsBMSW
912  << "\n nEventsBW " << bi->nEventsBW
913  << endl;
914 }
915 
916 ////////////////////////////////////////////////////////////////////////////////
917 /// sort entries, as needed for loopRange()
918 
920 {
921  // will loop over all events by default
922  if (!_sortInput) {
923  _ibNoSort.clear();
924  for (unsigned int i = 0; i < _dataPtsR.size(); ++i) {
925  _ibNoSort[i] = kTRUE;
926  }
927  return;
928  }
929 
930  itVec itrVecR;
931  vector<TVectorD>::iterator dpRItr = _dataPtsR.begin();
932  for (Int_t i = 0; dpRItr != _dataPtsR.end(); ++dpRItr, ++i) {
933  if (bi) {
934  if (bi->bpsIdcs.find(i) != bi->bpsIdcs.end())
935  // if (_wMap.find(i)!=_wMap.end())
936  itrVecR.push_back(itPair(i, dpRItr));
937  } else
938  itrVecR.push_back(itPair(i, dpRItr));
939  }
940 
941  for (Int_t j=0; j<_nDim; j++) {
942  _sortTVIdcs[j].clear();
943  sort(itrVecR.begin(),itrVecR.end(),SorterTV_L2H(j));
944  _sortTVIdcs[j] = itrVecR;
945  }
946 
947  for (Int_t j=0; j<_nDim; j++) {
948  cxcoutD(Eval) << "RooNDKeysPdf::sortDataIndices() : Number of sorted events : " << _sortTVIdcs[j].size() << endl;
949  }
950 }
951 
952 ////////////////////////////////////////////////////////////////////////////////
953 
955 {
956  cxcoutD(Eval) << "RooNDKeysPdf::calculateBandWidth()" << endl;
957 
958  // non-adaptive bandwidth
959  // (default, and needed to calculate adaptive bandwidth)
960 
961  if(!_options.Contains("a")) {
962  cxcoutD(Eval) << "RooNDKeysPdf::calculateBandWidth() Using static bandwidth." << endl;
963  }
964 
965  // fixed width approximation
966  for (Int_t i=0; i<_nEvents; i++) {
967  vector<Double_t>& weight = _weights0[i];
968  for (Int_t j = 0; j < _nDim; j++) {
969  weight[j] = _n * (*_sigmaR)[j];
970  // cout<<"j: "<<j<<", _n: "<<_n<<", sigmaR="<<(*_sigmaR)[j]<<", weight="<<weight[j]<<endl;
971  }
972  }
973 
974  // adaptive width
975  if (_options.Contains("a")) {
976  cxcoutD(Eval) << "RooNDKeysPdf::calculateBandWidth() Using adaptive bandwidth." << endl;
977 
978  double sqrt12 = sqrt(12.);
979  double sqrtSigmaAvgR = sqrt(_sigmaAvgR);
980 
981  vector<Double_t> dummy(_nDim, 0.);
982  _weights1.resize(_nEvents, dummy);
983 
984  std::vector<std::vector<Double_t>> *weights_prev(0);
985  std::vector<std::vector<Double_t>> *weights_new(0);
986 
987  // cout << "Number of adaptive iterations: " << _nAdpt << endl;
988 
989  for (Int_t k = 1; k <= _nAdpt; ++k) {
990 
991  // cout << " Cycle: " << k << endl;
992 
993  // if multiple adaptive iterations, need to swap weight sets
994  if (k % 2) {
995  weights_prev = &_weights0;
996  weights_new = &_weights1;
997  } else {
998  weights_prev = &_weights1;
999  weights_new = &_weights0;
1000  }
1001 
1002  for (Int_t i = 0; i < _nEvents; ++i) {
1003  vector<Double_t> &x = _dataPts[i];
1004  Double_t f = TMath::Power(gauss(x, *weights_prev) / _nEventsW, -1. / (2. * _d));
1005 
1006  vector<Double_t> &weight = (*weights_new)[i];
1007  for (Int_t j = 0; j < _nDim; j++) {
1008  Double_t norm = (_n * (*_sigmaR)[j]) / sqrtSigmaAvgR;
1009  weight[j] = norm * f / sqrt12; // note additional factor of sqrt(12) compared with HEP-EX/0011057
1010  }
1011  }
1012  }
1013  // this is the latest updated weights set
1014  _weights = weights_new;
1015  }
1016 }
1017 
1018 ////////////////////////////////////////////////////////////////////////////////
1019 /// loop over all closest point to x, as determined by loopRange()
1020 
1021 Double_t RooNDKeysPdf::gauss(vector<Double_t>& x, vector<vector<Double_t> >& weights) const
1022 {
1023  if(_nEvents==0) return 0.;
1024 
1025  Double_t z=0.;
1026  map<Int_t,Bool_t> ibMap;
1027 
1028  // determine input loop range for event x
1029  if (_sortInput) {
1030  loopRange(x, ibMap);
1031  }
1032 
1033  map<Int_t, Bool_t>::iterator ibMapItr, ibMapEnd;
1034  ibMapItr = (_sortInput ? ibMap.begin() : _ibNoSort.begin());
1035  ibMapEnd = (_sortInput ? ibMap.end() : _ibNoSort.end());
1036 
1037  for (; ibMapItr != ibMapEnd; ++ibMapItr) {
1038  Int_t i = (*ibMapItr).first;
1039 
1040  Double_t g(1.);
1041 
1042  if (i >= (Int_t)_idx.size()) {
1043  continue;
1044  } //---> 1.myline
1045 
1046  const vector<Double_t> &point = _dataPts[i];
1047  const vector<Double_t> &weight = weights[_idx[i]];
1048 
1049  for (Int_t j = 0; j < _nDim; j++) {
1050  (*_dx)[j] = x[j] - point[j];
1051  }
1052 
1053  if (_nDim > 1 && _rotate) {
1054  *_dx *= *_rotMat; // rotate to decorrelated frame!
1055  }
1056 
1057  for (Int_t j = 0; j < _nDim; j++) {
1058  Double_t r = (*_dx)[j]; // x[j] - point[j];
1059  Double_t c = 1. / (2. * weight[j] * weight[j]);
1060 
1061  // cout << "j = " << j << " x[j] = " << point[j] << " w = " << weight[j] << endl;
1062 
1063  g *= exp(-c * r * r);
1064  g *= 1. / (_sqrt2pi * weight[j]);
1065  }
1066  z += (g * _wMap[_idx[i]]);
1067  }
1068  return z;
1069 }
1070 
1071 ////////////////////////////////////////////////////////////////////////////////
1072 /// determine closest points to x, to loop over in evaluate()
1073 
1074 void RooNDKeysPdf::loopRange(vector<Double_t>& x, map<Int_t,Bool_t>& ibMap) const
1075 {
1076  ibMap.clear();
1077  TVectorD xRm(_nDim);
1078  TVectorD xRp(_nDim);
1079 
1080  for (Int_t j = 0; j < _nDim; j++) {
1081  xRm[j] = xRp[j] = x[j];
1082  }
1083 
1084  if (_nDim > 1 && _rotate) {
1085  xRm *= *_rotMat;
1086  xRp *= *_rotMat;
1087  }
1088  for (Int_t j = 0; j < _nDim; j++) {
1089  xRm[j] -= _nSigma * (_n * (*_sigmaR)[j]);
1090  xRp[j] += _nSigma * (_n * (*_sigmaR)[j]);
1091  // cout<<"xRm["<<j<<"]="<<xRm[j]<<endl;
1092  // cout<<"xRp["<<j<<"]="<<xRp[j]<<endl;
1093  }
1094 
1095  vector<TVectorD> xvecRm(1,xRm);
1096  vector<TVectorD> xvecRp(1,xRp);
1097 
1098  map<Int_t,Bool_t> ibMapRT;
1099 
1100  for (Int_t j=0; j<_nDim; j++) {
1101  ibMap.clear();
1102  itVec::iterator lo = lower_bound(_sortTVIdcs[j].begin(), _sortTVIdcs[j].end(),
1103  itPair(0,xvecRm.begin()), SorterTV_L2H(j));
1104  itVec::iterator hi =
1105  upper_bound(_sortTVIdcs[j].begin(), _sortTVIdcs[j].end(), itPair(0, xvecRp.begin()), SorterTV_L2H(j));
1106  itVec::iterator it = lo;
1107 
1108  if (j==0) {
1109  if (_nDim==1) { for (it=lo; it!=hi; ++it) ibMap[(*it).first] = kTRUE; }
1110  else { for (it=lo; it!=hi; ++it) ibMapRT[(*it).first] = kTRUE; }
1111  continue;
1112  }
1113 
1114  for (it=lo; it!=hi; ++it)
1115  if (ibMapRT.find((*it).first)!=ibMapRT.end()) { ibMap[(*it).first] = kTRUE; }
1116 
1117  ibMapRT.clear();
1118  if (j!=_nDim-1) { ibMapRT = ibMap; }
1119  }
1120 }
1121 
1122 ////////////////////////////////////////////////////////////////////////////////
1123 
1124 void RooNDKeysPdf::boxInfoInit(BoxInfo* bi, const char* rangeName, Int_t /*code*/) const
1125 {
1126  vector<Bool_t> doInt(_nDim,kTRUE);
1127 
1128  bi->filled = kFALSE;
1129 
1130  bi->xVarLo.resize(_nDim,0.);
1131  bi->xVarHi.resize(_nDim,0.);
1132  bi->xVarLoM3s.resize(_nDim,0.);
1133  bi->xVarLoP3s.resize(_nDim,0.);
1134  bi->xVarHiM3s.resize(_nDim,0.);
1135  bi->xVarHiP3s.resize(_nDim,0.);
1136 
1137  bi->netFluxZ = kTRUE;
1138  bi->bpsIdcs.clear();
1139  bi->bIdcs.clear();
1140  bi->sIdcs.clear();
1141  bi->bmsIdcs.clear();
1142 
1143  bi->nEventsBMSW=0.;
1144  bi->nEventsBW=0.;
1145 
1146  _varItr->Reset() ;
1147  RooRealVar* var ;
1148  for(Int_t j=0; (var=(RooRealVar*)_varItr->Next()); ++j) {
1149  if (doInt[j]) {
1150  bi->xVarLo[j] = var->getMin(rangeName);
1151  bi->xVarHi[j] = var->getMax(rangeName);
1152  } else {
1153  bi->xVarLo[j] = var->getVal() ;
1154  bi->xVarHi[j] = var->getVal() ;
1155  }
1156  }
1157 }
1158 
1159 ////////////////////////////////////////////////////////////////////////////////
1160 
1162 {
1163  if (_tracker != NULL && _tracker->hasChanged(kTRUE)) {
1164  updateRho(); // update internal rho parameters
1165  // redetermine static and/or adaptive bandwidth
1167  }
1168 
1169  _varItr->Reset();
1170  RooAbsReal *var;
1171  const RooArgSet *nset = _varList.nset();
1172  for (Int_t j = 0; (var = (RooAbsReal *)_varItr->Next()); ++j) {
1173  _x[j] = var->getVal(nset);
1174  }
1175 
1176  Double_t val = gauss(_x,*_weights);
1177  //cout<<"returning "<<val<<endl;
1178 
1179  if (val>=1E-20)
1180  return val ;
1181  else
1182  return (1E-20) ;
1183 }
1184 
1185 ////////////////////////////////////////////////////////////////////////////////
1186 
1187 Int_t RooNDKeysPdf::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
1188 {
1189  if (rangeName) return 0 ;
1190 
1191  Int_t code=0;
1192  if (matchArgs(allVars,analVars,RooArgSet(_varList))) { code=1; }
1193 
1194  return code;
1195 
1196 }
1197 
1198 ////////////////////////////////////////////////////////////////////////////////
1199 
1200 Double_t RooNDKeysPdf::analyticalIntegral(Int_t code, const char* rangeName) const
1201 {
1202  cxcoutD(Eval) << "Calling RooNDKeysPdf::analyticalIntegral(" << GetName() << ") with code " << code
1203  << " and rangeName " << (rangeName?rangeName:"<none>") << endl;
1204 
1205  // determine which observables need to be integrated over ...
1206  Int_t nComb = 1 << (_nDim);
1207  R__ASSERT(code>=1 && code<nComb) ;
1208 
1209  vector<Bool_t> doInt(_nDim,kTRUE);
1210 
1211  // get BoxInfo
1212  BoxInfo* bi(0);
1213 
1214  if (rangeName) {
1215  string rangeNameStr(rangeName) ;
1216  bi = _rangeBoxInfo[make_pair(rangeNameStr,code)] ;
1217  if (!bi) {
1218  bi = new BoxInfo ;
1219  _rangeBoxInfo[make_pair(rangeNameStr,code)] = bi ;
1220  boxInfoInit(bi,rangeName,code);
1221  }
1222  } else bi= &_fullBoxInfo ;
1223 
1224  // have boundaries changed?
1225  Bool_t newBounds(kFALSE);
1226  _varItr->Reset() ;
1227  RooRealVar* var ;
1228  for(Int_t j=0; (var=(RooRealVar*)_varItr->Next()); ++j) {
1229  if ((var->getMin(rangeName)-bi->xVarLo[j]!=0) ||
1230  (var->getMax(rangeName)-bi->xVarHi[j]!=0)) {
1231  newBounds = kTRUE;
1232  }
1233  }
1234 
1235  // reset
1236  if (newBounds) {
1237  cxcoutD(Eval) << "RooNDKeysPdf::analyticalIntegral() : Found new boundaries ... " << (rangeName?rangeName:"<none>") << endl;
1238  boxInfoInit(bi,rangeName,code);
1239  }
1240 
1241  // recalculates netFluxZero and nEventsIR
1242  if (!bi->filled || newBounds) {
1243  // Fill box info with contents
1244  calculateShell(bi);
1245  calculatePreNorm(bi);
1246  bi->filled = kTRUE;
1247  sortDataIndices(bi);
1248  }
1249 
1250  // first guess
1251  Double_t norm=bi->nEventsBW;
1252 
1253  if (_mirror && bi->netFluxZ) {
1254  // KEYS expression is self-normalized
1255  cxcoutD(Eval) << "RooNDKeysPdf::analyticalIntegral() : Using mirrored normalization : " << bi->nEventsBW << endl;
1256  return bi->nEventsBW;
1257  }
1258  // calculate leakage in and out of variable range box
1259  else
1260  {
1261  norm = bi->nEventsBMSW;
1262  if (norm<0.) norm=0.;
1263 
1264  for (Int_t i=0; i<Int_t(bi->sIdcs.size()); ++i) {
1265  Double_t prob=1.;
1266  const vector<Double_t>& x = _dataPts[bi->sIdcs[i]];
1267  const vector<Double_t>& weight = (*_weights)[_idx[bi->sIdcs[i]]];
1268 
1269  vector<Double_t> chi(_nDim,100.);
1270 
1271  for (Int_t j=0; j<_nDim; j++) {
1272  if(!doInt[j]) continue;
1273 
1274  if ((x[j]>bi->xVarLoM3s[j] && x[j]<bi->xVarLoP3s[j]) && x[j]<(bi->xVarLo[j]+bi->xVarHi[j])/2.)
1275  chi[j] = (x[j]-bi->xVarLo[j])/weight[j];
1276  else if ((x[j]>bi->xVarHiM3s[j] && x[j]<bi->xVarHiP3s[j]) && x[j]>(bi->xVarLo[j]+bi->xVarHi[j])/2.)
1277  chi[j] = (bi->xVarHi[j]-x[j])/weight[j];
1278 
1279  if (chi[j]>0) // inVarRange
1280  prob *= (0.5 + TMath::Erf(fabs(chi[j])/sqrt(2.))/2.);
1281  else // outside Var range
1282  prob *= (0.5 - TMath::Erf(fabs(chi[j])/sqrt(2.))/2.);
1283  }
1284 
1285  norm += prob * _wMap[_idx[bi->sIdcs[i]]];
1286  }
1287 
1288  cxcoutD(Eval) << "RooNDKeysPdf::analyticalIntegral() : Final normalization : " << norm << " " << bi->nEventsBW << endl;
1289  return norm;
1290  }
1291 }
1292 
1293 RooDataSet *
1294  ////////////////////////////////////////////////////////////////////////////////
1295 
1296  RooNDKeysPdf::createDatasetFromHist(const RooArgList &varList, const TH1 &hist) const
1297 {
1298  std::vector<RooRealVar *> varVec;
1299  RooArgSet varsAndWeightSet;
1300 
1301  TIterator *varItr = varList.createIterator();
1302  RooAbsArg *var;
1303  for (Int_t i = 0; (var = (RooAbsArg *)varItr->Next()); ++i) {
1304  if (!dynamic_cast<RooRealVar *>(var)) {
1305  coutE(InputArguments) << "RooNDKeysPdf::createDatasetFromHist(" << GetName() << ") WARNING: variable "
1306  << var->GetName() << " is not of type RooRealVar. Skip." << endl;
1307  continue;
1308  }
1309  varsAndWeightSet.add(*var); // used for dataset creation
1310  varVec.push_back(static_cast<RooRealVar *>(var)); // used for setting the variables.
1311  }
1312  delete varItr;
1313 
1314  /// Add weight
1315  RooRealVar weight("weight", "event weight", 0);
1316  varsAndWeightSet.add(weight);
1317 
1318  /// determine histogram dimensionality
1319  unsigned int histndim(0);
1320  std::string classname = hist.ClassName();
1321  if (classname.find("TH1") == 0) {
1322  histndim = 1;
1323  } else if (classname.find("TH2") == 0) {
1324  histndim = 2;
1325  } else if (classname.find("TH3") == 0) {
1326  histndim = 3;
1327  }
1328  assert(histndim == varVec.size());
1329 
1330  if (histndim > 3 || histndim <= 0) {
1331  coutE(InputArguments) << "RooNDKeysPdf::createDatasetFromHist(" << GetName()
1332  << ") ERROR: input histogram dimension not between [1-3]: " << histndim << endl;
1333  assert(0);
1334  }
1335 
1336  /// dataset creation
1337  RooDataSet *dataFromHist = new RooDataSet("datasetFromHist", "datasetFromHist", varsAndWeightSet, weight.GetName());
1338 
1339  /// dataset filling
1340  for (int i = 1; i <= hist.GetXaxis()->GetNbins(); ++i) {
1341  // 1 or more dimension
1342 
1343  Double_t xval = hist.GetXaxis()->GetBinCenter(i);
1344  varVec[0]->setVal(xval);
1345 
1346  if (varVec.size() == 1) {
1347  Double_t fval = hist.GetBinContent(i);
1348  weight.setVal(fval);
1349  dataFromHist->add(varsAndWeightSet, fval);
1350  } else { // 2 or more dimensions
1351 
1352  for (int j = 1; j <= hist.GetYaxis()->GetNbins(); ++j) {
1353  Double_t yval = hist.GetYaxis()->GetBinCenter(j);
1354  varVec[1]->setVal(yval);
1355 
1356  if (varVec.size() == 2) {
1357  Double_t fval = hist.GetBinContent(i, j);
1358  weight.setVal(fval);
1359  dataFromHist->add(varsAndWeightSet, fval);
1360  } else { // 3 dimensions
1361 
1362  for (int k = 1; k <= hist.GetZaxis()->GetNbins(); ++k) {
1363  Double_t zval = hist.GetZaxis()->GetBinCenter(k);
1364  varVec[2]->setVal(zval);
1365 
1366  Double_t fval = hist.GetBinContent(i, j, k);
1367  weight.setVal(fval);
1368  dataFromHist->add(varsAndWeightSet, fval);
1369  }
1370  }
1371  }
1372  }
1373  }
1374 
1375  return dataFromHist;
1376 }
1377 
1378 TMatrixD
1379  ////////////////////////////////////////////////////////////////////////////////
1380  /// Return evaluated weights
1381 
1382  RooNDKeysPdf::getWeights(const int &k) const
1383 {
1384  TMatrixD mref(_nEvents, _nDim + 1);
1385 
1386  cxcoutD(Eval) << "RooNDKeysPdf::getWeights() Return evaluated weights." << endl;
1387 
1388  for (Int_t i = 0; i < _nEvents; ++i) {
1389  vector<Double_t> &x = _dataPts[i];
1390  for (Int_t j = 0; j < _nDim; j++) {
1391  mref(i, j) = x[j];
1392  }
1393 
1394  vector<Double_t> &weight = (*_weights)[i];
1395  mref(i, _nDim) = weight[k];
1396  }
1397 
1398  return mref;
1399 }
1400 
1401 void
1402  ////////////////////////////////////////////////////////////////////////////////
1403 
1405 {
1406  _rhoItr->Reset();
1407  RooAbsReal *rho(0);
1408  for (Int_t j = 0; (rho = (RooAbsReal *)_rhoItr->Next()); ++j) {
1409  _rho[j] = rho->getVal();
1410  }
1411 
1412  if (_nDim > 1 && _rotate) {
1413  TMatrixDSym covMatRho(_nDim); // covariance matrix times rho parameters
1414  for (Int_t j = 0; j < _nDim; j++) {
1415  for (Int_t k = 0; k < _nDim; k++) {
1416  covMatRho(j, k) = (*_covMat)(j, k) * _rho[j] * _rho[k];
1417  }
1418  }
1419  // find decorrelation matrix and eigenvalues (R)
1420  TMatrixDSymEigen evCalculatorRho(covMatRho);
1421  *_sigmaR = evCalculatorRho.GetEigenValues();
1422  for (Int_t j = 0; j < _nDim; j++) {
1423  (*_sigmaR)[j] = sqrt((*_sigmaR)[j]);
1424  }
1425  } else {
1426  for (Int_t j = 0; j < _nDim; j++) {
1427  (*_sigmaR)[j] = (_sigma[j] * _rho[j]);
1428  } // * rho
1429  }
1430 }
Double_t _n
Definition: RooNDKeysPdf.h:154
Bool_t _rotate
Definition: RooNDKeysPdf.h:204
for(Int_t i=0;i< n;i++)
Definition: legend1.C:18
virtual Double_t getMin(const char *name=0) const
Double_t gauss(std::vector< Double_t > &x, std::vector< std::vector< Double_t >> &weights) const
loop over all closest point to x, as determined by loopRange()
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
std::vector< itPair > itVec
Definition: RooNDKeysPdf.h:43
TIterator * createIterator(Bool_t dir=kIterForward) const
#define coutE(a)
Definition: RooMsgService.h:34
void setOptions() const
set the configuration
RooListProxy _varList
Definition: RooNDKeysPdf.h:113
std::vector< Double_t > xVarHi
Definition: RooNDKeysPdf.h:103
std::vector< Double_t > _xVarLo
Definition: RooNDKeysPdf.h:181
void loopRange(std::vector< Double_t > &x, std::map< Int_t, Bool_t > &ibMap) const
determine closest points to x, to loop over in evaluate()
virtual Double_t getMax(const char *name=0) const
virtual Bool_t add(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling add() for each element in the source coll...
Definition: RooArgSet.h:86
Double_t _nEventsW
Definition: RooNDKeysPdf.h:152
auto * m
Definition: textangle.C:8
void updateRho() const
TMatrixDSym * _corrMat
Definition: RooNDKeysPdf.h:198
Double_t _nSigma
Definition: RooNDKeysPdf.h:141
virtual void Reset()=0
void calculateBandWidth() const
void boxInfoInit(BoxInfo *bi, const char *rangeName, Int_t code) const
virtual TMatrixTBase< Element > & Zero()
Set matrix elements to zero.
void loadDataSet(Bool_t firstCall) const
copy the dataset and calculate some useful variables
virtual const RooArgSet * get() const
Definition: RooAbsData.h:79
std::vector< std::string > _varName
Definition: RooNDKeysPdf.h:169
RooNDKeysPdf(const char *name, const char *title, const RooArgList &varList, const RooAbsData &data, TString options="ma", Double_t rho=1, Double_t nSigma=3, Bool_t rotate=kTRUE, Bool_t sortInput=kTRUE)
Construct N-dimensional kernel estimation p.d.f.
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
#define coutI(a)
Definition: RooMsgService.h:31
std::vector< Int_t > bIdcs
Definition: RooNDKeysPdf.h:107
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
#define cxcoutD(a)
Definition: RooMsgService.h:79
TVectorD * _dx
Definition: RooNDKeysPdf.h:201
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition: TH1.cxx:4763
Bool_t _netFluxZ
Definition: RooNDKeysPdf.h:178
Bool_t hasChanged(Bool_t clearState)
Returns true if state has changes since last call with clearState=kTRUE If clearState is true...
void mirrorDataSet() const
determine mirror dataset.
#define R__ASSERT(e)
Definition: TError.h:96
std::vector< Int_t > _bIdcs
Definition: RooNDKeysPdf.h:186
Int_t GetNrows() const
Definition: TVectorT.h:75
std::vector< Double_t > _xDatLo
Definition: RooNDKeysPdf.h:175
TMatrixT< Element > & T()
Definition: TMatrixT.h:150
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
STL namespace.
#define coutW(a)
Definition: RooMsgService.h:33
TMatrixD * _rotMat
Definition: RooNDKeysPdf.h:199
virtual ~RooNDKeysPdf()
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Definition: fillpatterns.C:1
void Print(Option_t *option="") const
Print the vector as a list of elements.
Definition: TVectorT.cxx:1361
std::vector< Double_t > _x1
Definition: RooNDKeysPdf.h:173
Short_t Abs(Short_t d)
Definition: TMathBase.h:108
Double_t _widthFactor
Definition: RooNDKeysPdf.h:140
Iterator abstract base class.
Definition: TIterator.h:30
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:627
TMatrixDSymEigen.
Bool_t _verbose
Definition: RooNDKeysPdf.h:146
std::vector< Double_t > _xDatHi
Definition: RooNDKeysPdf.h:175
std::vector< Double_t > _xVarLoM3s
Definition: RooNDKeysPdf.h:182
std::vector< std::vector< Double_t > > _weights1
Definition: RooNDKeysPdf.h:161
double sqrt(double)
std::vector< Double_t > _xVarHiM3s
Definition: RooNDKeysPdf.h:182
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition: TObject.cxx:128
Double_t x[n]
Definition: legend1.C:17
void initialize() const
initialization
std::map< Int_t, Double_t > _wMap
Definition: RooNDKeysPdf.h:195
std::map< Int_t, Bool_t > _ibNoSort
Definition: RooNDKeysPdf.h:184
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition: TAxis.cxx:464
std::map< std::pair< std::string, int >, BoxInfo * > _rangeBoxInfo
Definition: RooNDKeysPdf.h:189
friend class RooArgSet
Definition: RooAbsArg.h:471
TVectorD * _sigmaR
Definition: RooNDKeysPdf.h:200
std::vector< Double_t > _xDatLo3s
Definition: RooNDKeysPdf.h:176
constexpr Double_t Pi()
Definition: TMath.h:40
RooDataSet * createDatasetFromHist(const RooArgList &varList, const TH1 &hist) const
TVectorT< Double_t > TVectorD
Definition: TVectorDfwd.h:22
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Reimplementation of standard RooArgList::add()
std::vector< std::vector< Double_t > > * _weights
Definition: RooNDKeysPdf.h:162
const RooArgSet * nset() const
Definition: RooAbsProxy.h:46
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
TMatrixT< Double_t > TMatrixD
Definition: TMatrixDfwd.h:22
RooListProxy _rhoList
do not persist
Definition: RooNDKeysPdf.h:116
TMatrixD getWeights(const int &k) const
Return evaluated weights.
virtual void setVal(Double_t value)
Set value of variable to &#39;value&#39;.
Definition: RooRealVar.cxx:205
std::vector< Double_t > xVarHiP3s
Definition: RooNDKeysPdf.h:104
Int_t getSize() const
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
Bool_t _sortInput
Definition: RooNDKeysPdf.h:205
std::vector< Double_t > _x
Definition: RooNDKeysPdf.h:172
ROOT::R::TRInterface & r
Definition: Object.C:4
std::vector< itVec > _sortTVIdcs
Definition: RooNDKeysPdf.h:166
Double_t Erf(Double_t x)
Computation of the error function erf(x).
Definition: TMath.cxx:187
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported...
Double_t evaluate() const
do not persist
Double_t _maxWeight
Definition: RooNDKeysPdf.h:194
TVectorT< Element > & Zero()
Set vector elements to zero.
Definition: TVectorT.cxx:451
Double_t _minWeight
Definition: RooNDKeysPdf.h:193
Double_t _sigmaAvgR
Definition: RooNDKeysPdf.h:202
TMatrixDSym * _covMat
Definition: RooNDKeysPdf.h:197
std::vector< Double_t > xVarHiM3s
Definition: RooNDKeysPdf.h:104
TAxis * GetYaxis()
Definition: TH1.h:316
void calculatePreNorm(BoxInfo *bi) const
bi->nEventsBMSW=0.
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition: RooAbsData.h:37
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:29
std::vector< Double_t > _x0
Definition: RooNDKeysPdf.h:173
const TMatrixD & GetEigenVectors() const
Generic N-dimensional implementation of a kernel estimation p.d.f.
Definition: RooNDKeysPdf.h:48
void createPdf(Bool_t firstCall=kTRUE) const
evaluation order of constructor.
constexpr Double_t E()
Definition: TMath.h:74
virtual void add(const RooArgSet &row, Double_t weight=1.0, Double_t weightError=0)
Add a data point, with its coordinates specified in the &#39;data&#39; argset, to the data set...
Double_t _nEventsBW
Definition: RooNDKeysPdf.h:179
void calculateShell(BoxInfo *bi) const
determine points in +/- nSigma shell around the box determined by the variable ranges.
const Bool_t kFALSE
Definition: RtypesCore.h:88
TString _options
Definition: RooNDKeysPdf.h:139
BoxInfo _fullBoxInfo
Definition: RooNDKeysPdf.h:190
std::vector< Double_t > _mean
Definition: RooNDKeysPdf.h:174
void loadWeightSet() const
std::vector< Double_t > _xVarHiP3s
Definition: RooNDKeysPdf.h:182
Double_t _nEventsBMSW
Definition: RooNDKeysPdf.h:180
const TVectorD & GetEigenValues() const
#define ClassImp(name)
Definition: Rtypes.h:359
RooAbsArg * find(const char *name) const
Find object with given name in list.
double Double_t
Definition: RtypesCore.h:55
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
friend class RooDataSet
Definition: RooAbsArg.h:539
std::map< Int_t, Bool_t > _bpsIdcs
Definition: RooNDKeysPdf.h:183
std::vector< Int_t > sIdcs
Definition: RooNDKeysPdf.h:106
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral. ...
TMatrixTSym< Double_t > TMatrixDSym
static RooMathCoreReg dummy
std::vector< Double_t > xVarLoP3s
Definition: RooNDKeysPdf.h:104
Double_t y[n]
Definition: legend1.C:17
The TH1 histogram class.
Definition: TH1.h:56
std::vector< Double_t > _xVarHi
Definition: RooNDKeysPdf.h:181
std::vector< Int_t > _idx
Definition: RooNDKeysPdf.h:192
std::pair< Int_t, VecTVecDouble::iterator > itPair
Definition: RooNDKeysPdf.h:42
TIterator * _rhoItr
Definition: RooNDKeysPdf.h:117
std::vector< Double_t > _xDatHi3s
Definition: RooNDKeysPdf.h:176
RooChangeTracker is a meta object that tracks value changes in a given set of RooAbsArgs by registeri...
std::vector< Int_t > _sIdcs
Definition: RooNDKeysPdf.h:185
void Print(Option_t *name="") const
Print the matrix as a table of elements.
RooChangeTracker * _tracker
Definition: RooNDKeysPdf.h:208
TAxis * GetZaxis()
Definition: TH1.h:317
you should not use this method at all Int_t Int_t z
Definition: TRolke.cxx:630
std::map< Int_t, Bool_t > bpsIdcs
Definition: RooNDKeysPdf.h:105
std::vector< Int_t > bmsIdcs
Definition: RooNDKeysPdf.h:108
typedef void((*Func_t)())
const RooAbsData & _data
do not persist
Definition: RooNDKeysPdf.h:138
std::vector< Double_t > _rho
Definition: RooNDKeysPdf.h:170
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
RooDataSet * _dataP
Definition: RooNDKeysPdf.h:137
auto * l
Definition: textangle.C:4
Bool_t _fixedShape
Definition: RooNDKeysPdf.h:143
virtual TObject * Next()=0
Double_t _sqrt2pi
Definition: RooNDKeysPdf.h:148
virtual Double_t weight() const =0
std::vector< TVectorD > _dataPtsR
Definition: RooNDKeysPdf.h:159
std::vector< std::vector< Double_t > > _dataPts
Definition: RooNDKeysPdf.h:158
std::vector< std::vector< Double_t > > _weights0
Definition: RooNDKeysPdf.h:160
TIterator * _varItr
Definition: RooNDKeysPdf.h:114
float type_of_call hi(const int &, const int &)
void sortDataIndices(BoxInfo *bi=0) const
sort entries, as needed for loopRange()
Int_t GetNbins() const
Definition: TAxis.h:121
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
double exp(double)
std::vector< Double_t > _xVarLoP3s
Definition: RooNDKeysPdf.h:182
Bool_t _mirror
Definition: RooNDKeysPdf.h:144
const Bool_t kTRUE
Definition: RtypesCore.h:87
std::vector< Double_t > _sigma
Definition: RooNDKeysPdf.h:174
std::vector< Double_t > xVarLo
Definition: RooNDKeysPdf.h:103
std::vector< Int_t > _bmsIdcs
Definition: RooNDKeysPdf.h:187
char name[80]
Definition: TGX11.cxx:109
TAxis * GetXaxis()
Get the behaviour adopted by the object about the statoverflows. See EStatOverflows for more informat...
Definition: TH1.h:315
virtual Int_t numEntries() const
Definition: RooAbsData.cxx:285
std::vector< Double_t > xVarLoM3s
Definition: RooNDKeysPdf.h:104
Double_t _d
Definition: RooNDKeysPdf.h:153
std::vector< Double_t > _x2
Definition: RooNDKeysPdf.h:173
static constexpr double g