ROOT  6.06/09
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 #include <iostream>
13 #include <algorithm>
14 #include <string>
15 
16 #include "TMath.h"
17 #include "TMatrixDSymEigen.h"
18 #include "RooNDKeysPdf.h"
19 #include "RooAbsReal.h"
20 #include "RooRealVar.h"
21 #include "RooRandom.h"
22 #include "RooDataSet.h"
23 #include "RooHist.h"
24 #include "RooMsgService.h"
25 
26 //////////////////////////////////////////////////////////////////////////////
27 //
28 // BEGIN_HTML
29 // Generic N-dimensional implementation of a kernel estimation p.d.f. This p.d.f. models the distribution
30 // of an arbitrary input dataset as a superposition of Gaussian kernels, one for each data point,
31 // each contributing 1/N to the total integral of the p.d.f.
32 // <p>
33 // If the 'adaptive mode' is enabled, the width of the Gaussian is adaptively calculated from the
34 // local density of events, i.e. narrow for regions with high event density to preserve details and
35 // wide for regions with log event density to promote smoothness. The details of the general algorithm
36 // are described in the following paper:
37 // <p>
38 // Cranmer KS, Kernel Estimation in High-Energy Physics.
39 // Computer Physics Communications 136:198-207,2001 - e-Print Archive: hep ex/0011057
40 // <p>
41 // For multi-dimensional datasets, the kernels are modeled by multidimensional Gaussians. The kernels are
42 // constructed such that they reflect the correlation coefficients between the observables
43 // in the input dataset.
44 // END_HTML
45 //
46 
47 #include "TError.h"
48 
49 using namespace std;
50 
52 
53 
54 
55 ////////////////////////////////////////////////////////////////////////////////
56 /// Construct N-dimensional kernel estimation p.d.f. in observables 'varList'
57 /// from dataset 'data'. Options can be
58 ///
59 /// 'a' = Use adaptive kernels (width varies with local event density)
60 /// 'm' = Mirror data points over observable boundaries. Improves modeling
61 /// behavior at edges for distributions that are not close to zero
62 /// at edge
63 /// 'd' = Debug flag
64 /// 'v' = Verbose flag
65 ///
66 /// The parameter rho (default = 1) provides an overall scale factor that can
67 /// be applied to the bandwith calculated for each kernel. The nSigma parameter
68 /// determines the size of the box that is used to search for contributing kernels
69 /// around a given point in observable space. The nSigma parameters is used
70 /// in case of non-adaptive bandwidths and for the 1st non-adaptive pass for
71 /// the calculation of adaptive keys p.d.f.s.
72 ///
73 /// The optional weight arguments allows to specify an observable or function
74 /// expression in observables that specifies the weight of each event.
75 
76 RooNDKeysPdf::RooNDKeysPdf(const char *name, const char *title,
77  const RooArgList& varList, RooDataSet& data,
78  TString options, Double_t rho, Double_t nSigma, Bool_t rotate) :
79  RooAbsPdf(name,title),
80  _varList("varList","List of variables",this),
81  _data(data),
82  _options(options),
83  _widthFactor(rho),
84  _nSigma(nSigma),
85  _weights(&_weights0),
86  _rotate(rotate)
87 {
88  // Constructor
89  _varItr = _varList.createIterator() ;
90 
91  TIterator* varItr = varList.createIterator() ;
92  RooAbsArg* var ;
93  for (Int_t i=0; (var = (RooAbsArg*)varItr->Next()); ++i) {
94  if (!dynamic_cast<RooAbsReal*>(var)) {
95  coutE(InputArguments) << "RooNDKeysPdf::ctor(" << GetName() << ") ERROR: variable " << var->GetName()
96  << " is not of type RooAbsReal" << endl ;
97  R__ASSERT(0) ;
98  }
99  _varList.add(*var) ;
100  _varName.push_back(var->GetName());
101  }
102  delete varItr ;
103 
104  createPdf();
105 }
106 
107 ////////////////////////////////////////////////////////////////////////////////
108 /// Constructor
109 
110 RooNDKeysPdf::RooNDKeysPdf(const char *name, const char *title,
111  const RooArgList& varList, RooDataSet& data, const TVectorD& rho,
112  TString options, Double_t nSigma, Bool_t rotate) :
113  RooAbsPdf(name,title),
114  _varList("varList","List of variables",this),
115  _data(data),
116  _options(options),
117  _widthFactor(-1.0),
118  _nSigma(nSigma),
119  _weights(&_weights0),
120  _rotate(rotate)
121 {
123 
124  TIterator* varItr = varList.createIterator() ;
125  RooAbsArg* var ;
126  for (Int_t i=0; (var = (RooAbsArg*)varItr->Next()); ++i) {
127  if (!dynamic_cast<RooAbsReal*>(var)) {
128  coutE(InputArguments) << "RhhNDKeysPdf::ctor(" << GetName() << ") ERROR: variable " << var->GetName()
129  << " is not of type RooAbsReal" << endl ;
130  R__ASSERT(0) ;
131  }
132  _varList.add(*var) ;
133  _varName.push_back(var->GetName());
134  }
135  delete varItr ;
136 
137  // copy rho widths
138  if( _varList.getSize() != rho.GetNrows() ) {
139  coutE(InputArguments) << "ERROR: RhhNDKeysPdf::RhhNDKeysPdf() : The vector-size of rho is different from that of varList."
140  << "Unable to create the PDF." << endl;
141  R__ASSERT ( _varList.getSize()==rho.GetNrows() );
142  }
143 
144  // negative width factor will serve as a switch in initialize()
145  // negative value means that a vector has been provided as input,
146  // and that _rho has already been set ...
147  _rho.resize( rho.GetNrows() );
148  for (Int_t j=0; j<rho.GetNrows(); j++) { _rho[j]=rho[j]; /*cout<<"RooNDKeysPdf c'tor, _rho["<<j<<"]="<<_rho[j]<<endl;*/ }
149 
150  createPdf(); // calls initialize ...
151 }
152 
153 
154 ////////////////////////////////////////////////////////////////////////////////
155 /// Backward compatibility constructor for (1-dim) RooKeysPdf. If you are a new user,
156 /// please use the first constructor form.
157 
158 RooNDKeysPdf::RooNDKeysPdf(const char *name, const char *title,
159  RooAbsReal& x, RooDataSet& data,
160  Mirror mirror, Double_t rho, Double_t nSigma, Bool_t rotate) :
161  RooAbsPdf(name,title),
162  _varList("varList","List of variables",this),
163  _data(data),
164  _options("a"),
165  _widthFactor(rho),
166  _nSigma(nSigma),
167  _weights(&_weights0),
168  _rotate(rotate)
169 {
171 
172  _varList.add(x) ;
173  _varName.push_back(x.GetName());
174 
175  if (mirror!=NoMirror) {
176  if (mirror!=MirrorBoth)
177  coutW(InputArguments) << "RooNDKeysPdf::RooNDKeysPdf() : Warning : asymmetric mirror(s) no longer supported." << endl;
178  _options="m";
179  }
180 
181  createPdf();
182 }
183 
184 
185 
186 ////////////////////////////////////////////////////////////////////////////////
187 /// Backward compatibility constructor for Roo2DKeysPdf. If you are a new user,
188 /// please use the first constructor form.
189 
190 RooNDKeysPdf::RooNDKeysPdf(const char *name, const char *title, RooAbsReal& x, RooAbsReal & y,
191  RooDataSet& data, TString options, Double_t rho, Double_t nSigma, Bool_t rotate) :
192  RooAbsPdf(name,title),
193  _varList("varList","List of variables",this),
194  _data(data),
195  _options(options),
196  _widthFactor(rho),
197  _nSigma(nSigma),
198  _weights(&_weights0),
199  _rotate(rotate)
200 {
202 
203  _varList.add(RooArgSet(x,y)) ;
204  _varName.push_back(x.GetName());
205  _varName.push_back(y.GetName());
206 
207  createPdf();
208 }
209 
210 
211 
212 ////////////////////////////////////////////////////////////////////////////////
213 /// Constructor
214 
215 RooNDKeysPdf::RooNDKeysPdf(const RooNDKeysPdf& other, const char* name) :
216  RooAbsPdf(other,name),
217  _varList("varList",this,other._varList),
218  _data(other._data),
219  _options(other._options),
220  _widthFactor(other._widthFactor),
221  _nSigma(other._nSigma),
222  _weights(&_weights0),
223  _rotate(other._rotate)
224 {
226 
227  _fixedShape = other._fixedShape;
228  _mirror = other._mirror;
229  _debug = other._debug;
230  _verbose = other._verbose;
231  _sqrt2pi = other._sqrt2pi;
232  _nDim = other._nDim;
233  _nEvents = other._nEvents;
234  _nEventsM = other._nEventsM;
235  _nEventsW = other._nEventsW;
236  _d = other._d;
237  _n = other._n;
238  _dataPts = other._dataPts;
239  _dataPtsR = other._dataPtsR;
240  _weights0 = other._weights0;
241  _weights1 = other._weights1;
242  if (_options.Contains("a")) { _weights = &_weights1; }
243  //_sortIdcs = other._sortIdcs;
244  _sortTVIdcs = other._sortTVIdcs;
245  _varName = other._varName;
246  _rho = other._rho;
247  _x = other._x;
248  _x0 = other._x0 ;
249  _x1 = other._x1 ;
250  _x2 = other._x2 ;
251  _xDatLo = other._xDatLo;
252  _xDatHi = other._xDatHi;
253  _xDatLo3s = other._xDatLo3s;
254  _xDatHi3s = other._xDatHi3s;
255  _mean = other._mean;
256  _sigma = other._sigma;
257 
258  // BoxInfo
259  _netFluxZ = other._netFluxZ;
260  _nEventsBW = other._nEventsBW;
261  _nEventsBMSW = other._nEventsBMSW;
262  _xVarLo = other._xVarLo;
263  _xVarHi = other._xVarHi;
264  _xVarLoM3s = other._xVarLoM3s;
265  _xVarLoP3s = other._xVarLoP3s;
266  _xVarHiM3s = other._xVarHiM3s;
267  _xVarHiP3s = other._xVarHiP3s;
268  _bpsIdcs = other._bpsIdcs;
269  _sIdcs = other._sIdcs;
270  _bIdcs = other._bIdcs;
271  _bmsIdcs = other._bmsIdcs;
272 
274  _fullBoxInfo = other._fullBoxInfo ;
275 
276  _idx = other._idx;
277  _minWeight = other._minWeight;
278  _maxWeight = other._maxWeight;
279  _wMap = other._wMap;
280 
281  _covMat = new TMatrixDSym(*other._covMat);
282  _corrMat = new TMatrixDSym(*other._corrMat);
283  _rotMat = new TMatrixD(*other._rotMat);
284  _sigmaR = new TVectorD(*other._sigmaR);
285  _dx = new TVectorD(*other._dx);
286  _sigmaAvgR = other._sigmaAvgR;
287 }
288 
289 
290 
291 ////////////////////////////////////////////////////////////////////////////////
292 
294 {
295  if (_varItr) delete _varItr;
296  if (_covMat) delete _covMat;
297  if (_corrMat) delete _corrMat;
298  if (_rotMat) delete _rotMat;
299  if (_sigmaR) delete _sigmaR;
300  if (_dx) delete _dx;
301 
302  // delete all the boxinfos map
303  while ( !_rangeBoxInfo.empty() ) {
304  map<pair<string,int>,BoxInfo*>::iterator iter = _rangeBoxInfo.begin();
305  BoxInfo* box= (*iter).second;
306  _rangeBoxInfo.erase(iter);
307  delete box;
308  }
309 
310  _dataPts.clear();
311  _dataPtsR.clear();
312  _weights0.clear();
313  _weights1.clear();
314  //_sortIdcs.clear();
315  _sortTVIdcs.clear();
316 }
317 
318 
319 void
320 
321 ////////////////////////////////////////////////////////////////////////////////
322 /// evaluation order of constructor.
323 
325 {
326  if (firstCall) {
327  // set options
328  setOptions();
329  // initialization
330  initialize();
331  }
332 
333 
334  // copy dataset, calculate sigma_i's, determine min and max event weight
335  loadDataSet(firstCall);
336 
337  // mirror dataset around dataset boundaries -- does not depend on event weights
338  if (_mirror) mirrorDataSet();
339 
340  // store indices and weights of events with high enough weights
341  loadWeightSet();
342 
343 
344  // store indices of events in variable boundaries and box shell.
345 //calculateShell(&_fullBoxInfo);
346  // calculate normalization needed in analyticalIntegral()
347 //calculatePreNorm(&_fullBoxInfo);
348 
349  // lookup table for determining which events contribute to a certain coordinate
350  sortDataIndices();
351 
352  // determine static and/or adaptive bandwidth
354 
355 }
356 
357 
358 void
359 
360 ////////////////////////////////////////////////////////////////////////////////
361 /// set the configuration
362 
364 {
365  _options.ToLower();
366 
367  if( _options.Contains("a") ) { _weights = &_weights1; }
368  else { _weights = &_weights0; }
369  if( _options.Contains("m") ) _mirror = true;
370  else _mirror = false;
371  if( _options.Contains("d") ) _debug = true;
372  else _debug = false;
373  if( _options.Contains("v") ) { _debug = true; _verbose = true; }
374  else { _debug = false; _verbose = false; }
375 
376  cxcoutD(InputArguments) << "RooNDKeysPdf::setOptions() options = " << _options
377  << "\n\tbandWidthType = " << _options.Contains("a")
378  << "\n\tmirror = " << _mirror
379  << "\n\tdebug = " << _debug
380  << "\n\tverbose = " << _verbose
381  << endl;
382 
383  if (_nSigma<2.0) {
384  coutW(InputArguments) << "RooNDKeysPdf::setOptions() : Warning : nSigma = " << _nSigma << " < 2.0. "
385  << "Calculated normalization could be too large."
386  << endl;
387  }
388 }
389 
390 
391 void
392 
393 ////////////////////////////////////////////////////////////////////////////////
394 /// initialization
395 
397 {
398  _sqrt2pi = sqrt(2.0*TMath::Pi()) ;
399  _nDim = _varList.getSize();
403 
404  _netFluxZ = kFALSE;
405  _nEventsBW = 0;
406  _nEventsBMSW = 0;
407 
408  if(_nDim==0) {
409  coutE(InputArguments) << "ERROR: RooNDKeysPdf::initialize() : The observable list is empty. "
410  << "Unable to begin generating the PDF." << endl;
411  R__ASSERT (_nDim!=0);
412  }
413 
414  if(_nEvents==0) {
415  coutE(InputArguments) << "ERROR: RooNDKeysPdf::initialize() : The input data set is empty. "
416  << "Unable to begin generating the PDF." << endl;
417  R__ASSERT (_nEvents!=0);
418  }
419 
420  _d = static_cast<Double_t>(_nDim);
421 
422  vector<Double_t> dummy(_nDim,0.);
423  _dataPts.resize(_nEvents,dummy);
424  _weights0.resize(_nEvents,dummy);
425  //_sortIdcs.resize(_nDim);
426  _sortTVIdcs.resize(_nDim);
427 
428  //rdh _rho.resize(_nDim,_widthFactor);
429 
430  if (_widthFactor>0) { _rho.resize(_nDim,_widthFactor); }
431  // else: _rho has been provided as external input
432 
433  _x.resize(_nDim,0.);
434  _x0.resize(_nDim,0.);
435  _x1.resize(_nDim,0.);
436  _x2.resize(_nDim,0.);
437 
438  _mean.resize(_nDim,0.);
439  _sigma.resize(_nDim,0.);
440 
441  _xDatLo.resize(_nDim,0.);
442  _xDatHi.resize(_nDim,0.);
443  _xDatLo3s.resize(_nDim,0.);
444  _xDatHi3s.resize(_nDim,0.);
445 
446  boxInfoInit(&_fullBoxInfo,0,0xFFFF);
447 
448  _minWeight=0;
449  _maxWeight=0;
450  _wMap.clear();
451 
452  _covMat = 0;
453  _corrMat= 0;
454  _rotMat = 0;
455  _sigmaR = 0;
456  _dx = new TVectorD(_nDim); _dx->Zero();
457  _dataPtsR.resize(_nEvents,*_dx);
458 
459  _varItr->Reset() ;
460  RooRealVar* var ;
461  for(Int_t j=0; (var=(RooRealVar*)_varItr->Next()); ++j) {
462  _xDatLo[j] = var->getMin();
463  _xDatHi[j] = var->getMax();
464  }
465 }
466 
467 
468 void
469 
470 ////////////////////////////////////////////////////////////////////////////////
471 /// copy the dataset and calculate some useful variables
472 
474 {
475  // first some initialization
476  _nEventsW=0.;
477 
478  TMatrixD mat(_nDim,_nDim);
479  if (!_covMat) _covMat = new TMatrixDSym(_nDim);
480  if (!_corrMat) _corrMat= new TMatrixDSym(_nDim);
481  if (!_rotMat) _rotMat = new TMatrixD(_nDim,_nDim);
482  if (!_sigmaR) _sigmaR = new TVectorD(_nDim);
483 
484  mat.Zero();
485  _covMat->Zero();
486  _corrMat->Zero();
487  _rotMat->Zero();
488  _sigmaR->Zero();
489 
490  const RooArgSet* values= _data.get();
491  vector<RooRealVar*> dVars(_nDim);
492  for (Int_t j=0; j<_nDim; j++) {
493  dVars[j] = (RooRealVar*)values->find(_varName[j].c_str());
494  _x0[j]=_x1[j]=_x2[j]=0.;
495  }
496 
497  _idx.clear();
498  for (Int_t i=0; i<_nEvents; i++) {
499 
500  _data.get(i); // fills dVars
501  _idx.push_back(i);
502  vector<Double_t>& point = _dataPts[i];
503  TVectorD& pointV = _dataPtsR[i];
504 
505  Double_t myweight = _data.weight(); // default is one?
506  if ( TMath::Abs(myweight)>_maxWeight ) { _maxWeight = TMath::Abs(myweight); }
507  _nEventsW += myweight;
508 
509  for (Int_t j=0; j<_nDim; j++) {
510  for (Int_t k=0; k<_nDim; k++) {
511  mat(j,k) += dVars[j]->getVal() * dVars[k]->getVal() * myweight;
512  }
513  // only need to do once
514  if (firstCall)
515  point[j] = pointV[j] = dVars[j]->getVal();
516 
517  _x0[j] += 1. * myweight;
518  _x1[j] += point[j] * myweight ;
519  _x2[j] += point[j] * point[j] * myweight ;
520  if (_x2[j]!=_x2[j]) exit(3);
521 
522  // only need to do once
523  if (firstCall) {
524  if (point[j]<_xDatLo[j]) { _xDatLo[j]=point[j]; }
525  if (point[j]>_xDatHi[j]) { _xDatHi[j]=point[j]; }
526  }
527  }
528  }
529 
530  _n = TMath::Power(4./(_nEventsW*(_d+2.)), 1./(_d+4.)) ;
531  // = (4/[n(dim(R) + 2)])^1/(dim(R)+4); dim(R) = 2
532  _minWeight = (0.5 - TMath::Erf(_nSigma/sqrt(2.))/2.) * _maxWeight;
533 
534  for (Int_t j=0; j<_nDim; j++) {
535  _mean[j] = _x1[j]/_x0[j];
536  _sigma[j] = sqrt(_x2[j]/_x0[j]-_mean[j]*_mean[j]);
537  }
538 
539  TMatrixDSym covMatRho(_nDim); // covariance matrix times rho parameters
540  for (Int_t j=0; j<_nDim; j++) {
541  for (Int_t k=0; k<_nDim; k++) {
542  (*_covMat)(j,k) = mat(j,k)/_x0[j] - _mean[j]*_mean[k];
543  covMatRho(j,k) = (mat(j,k)/_x0[j] - _mean[j]*_mean[k]) * _rho[j] * _rho[k];
544  }
545  }
546 
547  // find decorrelation matrix and eigenvalues (R)
548  TMatrixDSymEigen evCalculatorRho(covMatRho);
549  *_rotMat = evCalculatorRho.GetEigenVectors();
550  *_rotMat = _rotMat->T(); // transpose
551  *_sigmaR = evCalculatorRho.GetEigenValues();
552 
553 
554  // set rho = 1 because sigmaR now contains rho
555  for (Int_t j=0; j<_nDim; j++) { _rho[j] = 1.; }
556 
557  for (Int_t j=0; j<_nDim; j++) { (*_sigmaR)[j] = sqrt((*_sigmaR)[j]); }
558 
559  for (Int_t j=0; j<_nDim; j++) {
560  for (Int_t k=0; k<_nDim; k++)
561  (*_corrMat)(j,k) = (*_covMat)(j,k)/(_sigma[j]*_sigma[k]) ;
562  }
563 
564  if (_verbose) {
565  //_covMat->Print();
566  _rotMat->Print();
567  _corrMat->Print();
568  _sigmaR->Print();
569  }
570 
571  if (!_rotate) {
572  _rotMat->Print();
573  _sigmaR->Print();
574  TMatrixD haar(_nDim,_nDim);
575  TMatrixD unit(TMatrixD::kUnit,haar);
576  *_rotMat = unit;
577  for (Int_t j=0; j<_nDim; j++) { (*_sigmaR)[j] = _sigma[j]; }
578  _rotMat->Print();
579  _sigmaR->Print();
580  }
581 
582  // use raw sigmas (without rho) for sigmaAvgR
583  TMatrixDSymEigen evCalculator(*_covMat);
584  TVectorD sigmaRraw = evCalculator.GetEigenValues();
585  for (Int_t j=0; j<_nDim; j++) { sigmaRraw[j] = sqrt(sigmaRraw[j]); }
586 
587  _sigmaAvgR=1.;
588  //for (Int_t j=0; j<_nDim; j++) { _sigmaAvgR *= (*_sigmaR)[j]; }
589  for (Int_t j=0; j<_nDim; j++) { _sigmaAvgR *= sigmaRraw[j]; }
591 
592  for (Int_t i=0; i<_nEvents; i++) {
593 
594  TVectorD& pointR = _dataPtsR[i];
595  pointR *= *_rotMat;
596  }
597 
598  coutI(Contents) << "RooNDKeysPdf::loadDataSet(" << this << ")"
599  << "\n Number of events in dataset: " << _nEvents
600  << "\n Weighted number of events in dataset: " << _nEventsW
601  << endl;
602 }
603 
604 
605 void
606 
607 ////////////////////////////////////////////////////////////////////////////////
608 /// determine mirror dataset.
609 /// mirror points are added around the physical boundaries of the dataset
610 /// Two steps:
611 /// 1. For each entry, determine if it should be mirrored (the mirror configuration).
612 /// 2. For each mirror configuration, make the mirror points.
613 
615 {
616  for (Int_t j=0; j<_nDim; j++) {
617  _xDatLo3s[j] = _xDatLo[j] + _nSigma * (_rho[j] * _n * _sigma[j]);
618  _xDatHi3s[j] = _xDatHi[j] - _nSigma * (_rho[j] * _n * _sigma[j]);
619 
620  //cout<<"xDatLo3s["<<j<<"]="<<_xDatLo3s[j]<<endl;
621  //cout<<"xDatHi3s["<<j<<"]="<<_xDatHi3s[j]<<endl;
622  }
623 
624  vector<Double_t> dummy(_nDim,0.);
625 
626  // 1.
627  for (Int_t i=0; i<_nEvents; i++) {
628  vector<Double_t>& x = _dataPts[i];
629 
630  Int_t size = 1;
631  vector<vector<Double_t> > mpoints(size,dummy);
632  vector<vector<Int_t> > mjdcs(size);
633 
634  // create all mirror configurations for event i
635  for (Int_t j=0; j<_nDim; j++) {
636 
637  vector<Int_t>& mjdxK = mjdcs[0];
638  vector<Double_t>& mpointK = mpoints[0];
639 
640  // single mirror *at physical boundaries*
641  if ((x[j]>_xDatLo[j] && x[j]<_xDatLo3s[j]) && x[j]<(_xDatLo[j]+_xDatHi[j])/2.) {
642  mpointK[j] = 2.*_xDatLo[j]-x[j];
643  mjdxK.push_back(j);
644  } else if ((x[j]<_xDatHi[j] && x[j]>_xDatHi3s[j]) && x[j]>(_xDatLo[j]+_xDatHi[j])/2.) {
645  mpointK[j] = 2.*_xDatHi[j]-x[j];
646  mjdxK.push_back(j);
647  }
648  }
649 
650  vector<Int_t>& mjdx0 = mjdcs[0];
651  // no mirror point(s) for this event
652  if (size==1 && mjdx0.size()==0) continue;
653 
654  // 2.
655  // generate all mirror points for event i
656  vector<Int_t>& mjdx = mjdcs[0];
657  vector<Double_t>& mpoint = mpoints[0];
658 
659  // number of mirror points for this mirror configuration
660  Int_t eMir = 1 << mjdx.size();
661  vector<vector<Double_t> > epoints(eMir,x);
662 
663  for (Int_t m=0; m<Int_t(mjdx.size()); m++) {
664  Int_t size1 = 1 << m;
665  Int_t size2 = 1 << (m+1);
666  // copy all previous mirror points
667  for (Int_t l=size1; l<size2; ++l) {
668  epoints[l] = epoints[l-size1];
669  // fill high mirror points
670  vector<Double_t>& epoint = epoints[l];
671  epoint[mjdx[Int_t(mjdx.size()-1)-m]] = mpoint[mjdx[Int_t(mjdx.size()-1)-m]];
672  }
673  }
674 
675  // remove duplicate mirror points
676  // note that: first epoint == x
677  epoints.erase(epoints.begin());
678 
679  // add mirror points of event i to total dataset
680  TVectorD pointR(_nDim);
681 
682  for (Int_t m=0; m<Int_t(epoints.size()); m++) {
683  _idx.push_back(i);
684  _dataPts.push_back(epoints[m]);
685  //_weights0.push_back(_weights0[i]);
686  for (Int_t j=0; j<_nDim; j++) { pointR[j] = (epoints[m])[j]; }
687  pointR *= *_rotMat;
688  _dataPtsR.push_back(pointR);
689  }
690 
691  epoints.clear();
692  mpoints.clear();
693  mjdcs.clear();
694  } // end of event loop
695 
696  _nEventsM = Int_t(_dataPts.size());
697 }
698 
699 
700 void
701 ////////////////////////////////////////////////////////////////////////////////
702 
704 {
705  _wMap.clear();
706 
707  for (Int_t i=0; i<_nEventsM; i++) {
708  _data.get(_idx[i]);
709  Double_t myweight = _data.weight();
710  //if ( TMath::Abs(myweight)>_minWeight ) {
711  _wMap[i] = myweight;
712  //}
713  }
714 
715  coutI(Contents) << "RooNDKeysPdf::loadWeightSet(" << this << ") : Number of weighted events : " << _wMap.size() << endl;
716 }
717 
718 
719 void
720 ////////////////////////////////////////////////////////////////////////////////
721 /// determine points in +/- nSigma shell around the box determined by the variable
722 /// ranges. These points are needed in the normalization, to determine probability
723 /// leakage in and out of the box.
724 
726 {
727  for (Int_t j=0; j<_nDim; j++) {
728  if (bi->xVarLo[j]==_xDatLo[j] && bi->xVarHi[j]==_xDatHi[j]) {
729  bi->netFluxZ = bi->netFluxZ && kTRUE;
730  } else { bi->netFluxZ = kFALSE; }
731 
732  bi->xVarLoM3s[j] = bi->xVarLo[j] - _nSigma * (_rho[j] * _n * _sigma[j]);
733  bi->xVarLoP3s[j] = bi->xVarLo[j] + _nSigma * (_rho[j] * _n * _sigma[j]);
734  bi->xVarHiM3s[j] = bi->xVarHi[j] - _nSigma * (_rho[j] * _n * _sigma[j]);
735  bi->xVarHiP3s[j] = bi->xVarHi[j] + _nSigma * (_rho[j] * _n * _sigma[j]);
736 
737  //cout<<"bi->xVarLoM3s["<<j<<"]="<<bi->xVarLoM3s[j]<<endl;
738  //cout<<"bi->xVarLoP3s["<<j<<"]="<<bi->xVarLoP3s[j]<<endl;
739  //cout<<"bi->xVarHiM3s["<<j<<"]="<<bi->xVarHiM3s[j]<<endl;
740  //cout<<"bi->xVarHiM3s["<<j<<"]="<<bi->xVarHiM3s[j]<<endl;
741  }
742 
743  map<Int_t,Double_t>::iterator wMapItr = _wMap.begin();
744 
745  //for (Int_t i=0; i<_nEventsM; i++) {
746  for (; wMapItr!=_wMap.end(); ++wMapItr) {
747  Int_t i = (*wMapItr).first;
748 
749  const vector<Double_t>& x = _dataPts[i];
750  Bool_t inVarRange(kTRUE);
751  Bool_t inVarRangePlusShell(kTRUE);
752 
753  for (Int_t j=0; j<_nDim; j++) {
754 
755  if (x[j]>bi->xVarLo[j] && x[j]<bi->xVarHi[j]) {
756  inVarRange = inVarRange && kTRUE;
757  } else { inVarRange = inVarRange && kFALSE; }
758 
759  if (x[j]>bi->xVarLoM3s[j] && x[j]<bi->xVarHiP3s[j]) {
760  inVarRangePlusShell = inVarRangePlusShell && kTRUE;
761  } else { inVarRangePlusShell = inVarRangePlusShell && kFALSE; }
762  }
763 
764  // event in range?
765  if (inVarRange) {
766  bi->bIdcs.push_back(i);
767  }
768 
769  // event in shell?
770  if (inVarRangePlusShell) {
771  bi->bpsIdcs[i] = kTRUE;
772  Bool_t inShell(kFALSE);
773  for (Int_t j=0; j<_nDim; j++) {
774  if ((x[j]>bi->xVarLoM3s[j] && x[j]<bi->xVarLoP3s[j]) && x[j]<(bi->xVarLo[j]+bi->xVarHi[j])/2.) {
775  inShell = kTRUE;
776  } else if ((x[j]>bi->xVarHiM3s[j] && x[j]<bi->xVarHiP3s[j]) && x[j]>(bi->xVarLo[j]+bi->xVarHi[j])/2.) {
777  inShell = kTRUE;
778  }
779  }
780  if (inShell) bi->sIdcs.push_back(i); // needed for normalization
781  else {
782  bi->bmsIdcs.push_back(i); // idem
783  }
784  }
785  }
786 
787  coutI(Contents) << "RooNDKeysPdf::calculateShell() : "
788  << "\n Events in shell " << bi->sIdcs.size()
789  << "\n Events in box " << bi->bIdcs.size()
790  << "\n Events in box and shell " << bi->bpsIdcs.size()
791  << endl;
792 }
793 
794 
795 void
796 ////////////////////////////////////////////////////////////////////////////////
797 ///bi->nEventsBMSW=0.;
798 ///bi->nEventsBW=0.;
799 
801 {
802  // box minus shell
803  for (Int_t i=0; i<Int_t(bi->bmsIdcs.size()); i++)
804  bi->nEventsBMSW += _wMap[bi->bmsIdcs[i]];
805 
806  // box
807  for (Int_t i=0; i<Int_t(bi->bIdcs.size()); i++)
808  bi->nEventsBW += _wMap[bi->bIdcs[i]];
809 
810  cxcoutD(Eval) << "RooNDKeysPdf::calculatePreNorm() : "
811  << "\n nEventsBMSW " << bi->nEventsBMSW
812  << "\n nEventsBW " << bi->nEventsBW
813  << endl;
814 }
815 
816 
817 void
818 ////////////////////////////////////////////////////////////////////////////////
819 /// sort entries, as needed for loopRange()
820 
822 {
823  itVec itrVecR;
824  vector<TVectorD>::iterator dpRItr = _dataPtsR.begin();
825  for (Int_t i=0; dpRItr!=_dataPtsR.end(); ++dpRItr, ++i) {
826  if (bi) {
827  if (bi->bpsIdcs.find(i)!=bi->bpsIdcs.end())
828  //if (_wMap.find(i)!=_wMap.end())
829  itrVecR.push_back(itPair(i,dpRItr));
830  } else itrVecR.push_back(itPair(i,dpRItr));
831  }
832 
833  for (Int_t j=0; j<_nDim; j++) {
834  _sortTVIdcs[j].clear();
835  sort(itrVecR.begin(),itrVecR.end(),SorterTV_L2H(j));
836  _sortTVIdcs[j] = itrVecR;
837  }
838 
839  for (Int_t j=0; j<_nDim; j++) {
840  cxcoutD(Eval) << "RooNDKeysPdf::sortDataIndices() : Number of sorted events : " << _sortTVIdcs[j].size() << endl;
841  }
842 }
843 
844 
845 void
846 ////////////////////////////////////////////////////////////////////////////////
847 
849 {
850  cxcoutD(Eval) << "RooNDKeysPdf::calculateBandWidth()" << endl;
851 
852  // non-adaptive bandwidth
853  // (default, and needed to calculate adaptive bandwidth)
854 
855  if(!_options.Contains("a")) {
856  cxcoutD(Eval) << "RooNDKeysPdf::calculateBandWidth() Using static bandwidth." << endl;
857  }
858 
859  for (Int_t i=0; i<_nEvents; i++) {
860  vector<Double_t>& weight = _weights0[i];
861  for (Int_t j=0; j<_nDim; j++) {
862  weight[j] = _rho[j] * _n * (*_sigmaR)[j] ;
863  //cout<<"j: "<<j<<", rho="<<_rho[j]<<", _n: "<<_n<<", sigmaR="<<(*_sigmaR)[j]<<", weight="<<weight[j]<<endl;
864  }
865  }
866 
867 
868  // adaptive width
869  if(_options.Contains("a")) {
870  cxcoutD(Eval) << "RooNDKeysPdf::calculateBandWidth() Using adaptive bandwidth." << endl;
871 
872  double sqrt12=sqrt(12.);
873  double sqrtSigmaAvgR=sqrt(_sigmaAvgR);
874 
875  vector<Double_t> dummy(_nDim,0.);
876  _weights1.resize(_nEvents,dummy);
877 
878  for(Int_t i=0; i<_nEvents; ++i) {
879 
880  vector<Double_t>& x = _dataPts[i];
881  Double_t f = TMath::Power( gauss(x,_weights0)/_nEventsW , -1./(2.*_d) ) ;
882 
883  vector<Double_t>& weight = _weights1[i];
884  for (Int_t j=0; j<_nDim; j++) {
885  Double_t norm = (_rho[j]*_n*(*_sigmaR)[j]) / sqrtSigmaAvgR ;
886  //cout<<"norm["<<j<<"]="<<norm<<endl;
887  weight[j] = norm * f / sqrt12 ; // note additional factor of sqrt(12) compared with HEP-EX/0011057
888  }
889  }
890  _weights = &_weights1;
891  }
892 }
893 
894 
895 Double_t
896 ////////////////////////////////////////////////////////////////////////////////
897 /// loop over all closest point to x, as determined by loopRange()
898 
899 RooNDKeysPdf::gauss(vector<Double_t>& x, vector<vector<Double_t> >& weights) const
900 {
901  if(_nEvents==0) return 0.;
902 
903  Double_t z=0.;
904  map<Int_t,Bool_t> ibMap;
905  ibMap.clear();
906 
907  // determine loop range for event x
908  loopRange(x,ibMap);
909 
910  map<Int_t,Bool_t>::iterator ibMapItr = ibMap.begin();
911 
912  for (; ibMapItr!=ibMap.end(); ++ibMapItr) {
913  Int_t i = (*ibMapItr).first;
914 
915  Double_t g(1.);
916 
917  if(i>=(Int_t)_idx.size()) {continue;} //---> 1.myline
918 
919  const vector<Double_t>& point = _dataPts[i];
920  const vector<Double_t>& weight = weights[_idx[i]];
921 
922  for (Int_t j=0; j<_nDim; j++) {
923  (*_dx)[j] = x[j]-point[j];
924  }
925 
926  if (_nDim>1) {
927  *_dx *= *_rotMat; // rotate to decorrelated frame!
928  }
929 
930  for (Int_t j=0; j<_nDim; j++) {
931  Double_t r = (*_dx)[j]; //x[j] - point[j];
932  Double_t c = 1./(2.*weight[j]*weight[j]);
933 
934  g *= exp( -c*r*r );
935  g *= 1./(_sqrt2pi*weight[j]);
936  }
937  z += (g*_wMap[_idx[i]]);
938  }
939  return z;
940 }
941 
942 
943 void
944 ////////////////////////////////////////////////////////////////////////////////
945 /// determine closest points to x, to loop over in evaluate()
946 
947 RooNDKeysPdf::loopRange(vector<Double_t>& x, map<Int_t,Bool_t>& ibMap) const
948 {
949  TVectorD xRm(_nDim);
950  TVectorD xRp(_nDim);
951 
952  for (Int_t j=0; j<_nDim; j++) { xRm[j] = xRp[j] = x[j]; }
953 
954  xRm *= *_rotMat;
955  xRp *= *_rotMat;
956  for (Int_t j=0; j<_nDim; j++) {
957  xRm[j] -= _nSigma * (_rho[j] * _n * (*_sigmaR)[j]);
958  xRp[j] += _nSigma * (_rho[j] * _n * (*_sigmaR)[j]);
959  //cout<<"xRm["<<j<<"]="<<xRm[j]<<endl;
960  //cout<<"xRp["<<j<<"]="<<xRp[j]<<endl;
961  }
962 
963  vector<TVectorD> xvecRm(1,xRm);
964  vector<TVectorD> xvecRp(1,xRp);
965 
966  map<Int_t,Bool_t> ibMapRT;
967 
968  for (Int_t j=0; j<_nDim; j++) {
969  ibMap.clear();
970  itVec::iterator lo = lower_bound(_sortTVIdcs[j].begin(), _sortTVIdcs[j].end(),
971  itPair(0,xvecRm.begin()), SorterTV_L2H(j));
972  itVec::iterator hi = upper_bound(_sortTVIdcs[j].begin(), _sortTVIdcs[j].end(),
973  itPair(0,xvecRp.begin()), SorterTV_L2H(j));
974  itVec::iterator it=lo;
975  if (j==0) {
976  if (_nDim==1) { for (it=lo; it!=hi; ++it) ibMap[(*it).first] = kTRUE; }
977  else { for (it=lo; it!=hi; ++it) ibMapRT[(*it).first] = kTRUE; }
978  continue;
979  }
980 
981  for (it=lo; it!=hi; ++it)
982  if (ibMapRT.find((*it).first)!=ibMapRT.end()) { ibMap[(*it).first] = kTRUE; }
983 
984  ibMapRT.clear();
985  if (j!=_nDim-1) { ibMapRT = ibMap; }
986  }
987 }
988 
989 
990 void
991 ////////////////////////////////////////////////////////////////////////////////
992 
993 RooNDKeysPdf::boxInfoInit(BoxInfo* bi, const char* rangeName, Int_t /*code*/) const
994 {
995  vector<Bool_t> doInt(_nDim,kTRUE);
996 
997  bi->filled = kFALSE;
998 
999  bi->xVarLo.resize(_nDim,0.);
1000  bi->xVarHi.resize(_nDim,0.);
1001  bi->xVarLoM3s.resize(_nDim,0.);
1002  bi->xVarLoP3s.resize(_nDim,0.);
1003  bi->xVarHiM3s.resize(_nDim,0.);
1004  bi->xVarHiP3s.resize(_nDim,0.);
1005 
1006  bi->netFluxZ = kTRUE;
1007  bi->bpsIdcs.clear();
1008  bi->bIdcs.clear();
1009  bi->sIdcs.clear();
1010  bi->bmsIdcs.clear();
1011 
1012  bi->nEventsBMSW=0.;
1013  bi->nEventsBW=0.;
1014 
1015  _varItr->Reset() ;
1016  RooRealVar* var ;
1017  for(Int_t j=0; (var=(RooRealVar*)_varItr->Next()); ++j) {
1018  if (doInt[j]) {
1019  bi->xVarLo[j] = var->getMin(rangeName);
1020  bi->xVarHi[j] = var->getMax(rangeName);
1021  } else {
1022  bi->xVarLo[j] = var->getVal() ;
1023  bi->xVarHi[j] = var->getVal() ;
1024  }
1025  }
1026 }
1027 
1028 
1029 Double_t
1030 ////////////////////////////////////////////////////////////////////////////////
1031 
1033 {
1034  _varItr->Reset() ;
1035  RooAbsReal* var ;
1036  const RooArgSet* nset = _varList.nset() ;
1037  for(Int_t j=0; (var=(RooAbsReal*)_varItr->Next()); ++j) {
1038  _x[j] = var->getVal(nset);
1039  }
1040 
1041  Double_t val = gauss(_x,*_weights);
1042  //cout<<"returning "<<val<<endl;
1043 
1044  if (val>=1E-20)
1045  return val ;
1046  else
1047  return (1E-20) ;
1048 }
1049 
1050 
1051 Int_t
1052 ////////////////////////////////////////////////////////////////////////////////
1053 
1054 RooNDKeysPdf::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
1055 {
1056  if (rangeName) return 0 ;
1057 
1058  Int_t code=0;
1059  if (matchArgs(allVars,analVars,RooArgSet(_varList))) { code=1; }
1060 
1061  return code;
1062 
1063 }
1064 
1065 
1066 Double_t
1067 ////////////////////////////////////////////////////////////////////////////////
1068 
1069 RooNDKeysPdf::analyticalIntegral(Int_t code, const char* rangeName) const
1070 {
1071  cxcoutD(Eval) << "Calling RooNDKeysPdf::analyticalIntegral(" << GetName() << ") with code " << code
1072  << " and rangeName " << (rangeName?rangeName:"<none>") << endl;
1073 
1074  // determine which observables need to be integrated over ...
1075  Int_t nComb = 1 << (_nDim);
1076  R__ASSERT(code>=1 && code<nComb) ;
1077 
1078  vector<Bool_t> doInt(_nDim,kTRUE);
1079 
1080  // get BoxInfo
1081  BoxInfo* bi(0);
1082 
1083  if (rangeName) {
1084  string rangeNameStr(rangeName) ;
1085  bi = _rangeBoxInfo[make_pair(rangeNameStr,code)] ;
1086  if (!bi) {
1087  bi = new BoxInfo ;
1088  _rangeBoxInfo[make_pair(rangeNameStr,code)] = bi ;
1089  boxInfoInit(bi,rangeName,code);
1090  }
1091  } else bi= &_fullBoxInfo ;
1092 
1093  // have boundaries changed?
1094  Bool_t newBounds(kFALSE);
1095  _varItr->Reset() ;
1096  RooRealVar* var ;
1097  for(Int_t j=0; (var=(RooRealVar*)_varItr->Next()); ++j) {
1098  if ((var->getMin(rangeName)-bi->xVarLo[j]!=0) ||
1099  (var->getMax(rangeName)-bi->xVarHi[j]!=0)) {
1100  newBounds = kTRUE;
1101  }
1102  }
1103 
1104  // reset
1105  if (newBounds) {
1106  cxcoutD(Eval) << "RooNDKeysPdf::analyticalIntegral() : Found new boundaries ... " << (rangeName?rangeName:"<none>") << endl;
1107  boxInfoInit(bi,rangeName,code);
1108  }
1109 
1110  // recalculates netFluxZero and nEventsIR
1111  if (!bi->filled || newBounds) {
1112  // Fill box info with contents
1113  calculateShell(bi);
1114  calculatePreNorm(bi);
1115  bi->filled = kTRUE;
1116  sortDataIndices(bi);
1117  }
1118 
1119  // first guess
1120  Double_t norm=bi->nEventsBW;
1121 
1122  if (_mirror && bi->netFluxZ) {
1123  // KEYS expression is self-normalized
1124  cxcoutD(Eval) << "RooNDKeysPdf::analyticalIntegral() : Using mirrored normalization : " << bi->nEventsBW << endl;
1125  return bi->nEventsBW;
1126  }
1127  // calculate leakage in and out of variable range box
1128  else
1129  {
1130  norm = bi->nEventsBMSW;
1131  if (norm<0.) norm=0.;
1132 
1133  for (Int_t i=0; i<Int_t(bi->sIdcs.size()); ++i) {
1134  Double_t prob=1.;
1135  const vector<Double_t>& x = _dataPts[bi->sIdcs[i]];
1136  const vector<Double_t>& weight = (*_weights)[_idx[bi->sIdcs[i]]];
1137 
1138  vector<Double_t> chi(_nDim,100.);
1139 
1140  for (Int_t j=0; j<_nDim; j++) {
1141  if(!doInt[j]) continue;
1142 
1143  if ((x[j]>bi->xVarLoM3s[j] && x[j]<bi->xVarLoP3s[j]) && x[j]<(bi->xVarLo[j]+bi->xVarHi[j])/2.)
1144  chi[j] = (x[j]-bi->xVarLo[j])/weight[j];
1145  else if ((x[j]>bi->xVarHiM3s[j] && x[j]<bi->xVarHiP3s[j]) && x[j]>(bi->xVarLo[j]+bi->xVarHi[j])/2.)
1146  chi[j] = (bi->xVarHi[j]-x[j])/weight[j];
1147 
1148  if (chi[j]>0) // inVarRange
1149  prob *= (0.5 + TMath::Erf(fabs(chi[j])/sqrt(2.))/2.);
1150  else // outside Var range
1151  prob *= (0.5 - TMath::Erf(fabs(chi[j])/sqrt(2.))/2.);
1152  }
1153 
1154  norm += prob * _wMap[_idx[bi->sIdcs[i]]];
1155  }
1156 
1157  cxcoutD(Eval) << "RooNDKeysPdf::analyticalIntegral() : Final normalization : " << norm << " " << bi->nEventsBW << endl;
1158  return norm;
1159  }
1160 }
1161 
1162 
Double_t _n
Definition: RooNDKeysPdf.h:135
Bool_t _rotate
Definition: RooNDKeysPdf.h:184
for(Int_t i=0;i< n;i++)
Definition: legend1.C:18
const TVectorD & GetEigenValues() const
std::vector< itPair > itVec
Definition: RooNDKeysPdf.h:40
const RooArgSet * nset() const
Definition: RooAbsProxy.h:47
#define coutE(a)
Definition: RooMsgService.h:35
RooListProxy _varList
Definition: RooNDKeysPdf.h:99
std::vector< Double_t > xVarHi
Definition: RooNDKeysPdf.h:89
std::vector< Double_t > _xVarLo
Definition: RooNDKeysPdf.h:162
Double_t _nEventsW
Definition: RooNDKeysPdf.h:133
TMatrixDSym * _corrMat
Definition: RooNDKeysPdf.h:178
Double_t _nSigma
Definition: RooNDKeysPdf.h:122
virtual void Reset()=0
virtual TMatrixTBase< Element > & Zero()
Set matrix elements to zero.
void variables(TString fin="TMVA.root", TString dirName="InputVariables_Id", TString title="TMVA Input Variables", Bool_t isRegression=kFALSE, Bool_t useTMVAStyle=kTRUE)
Definition: variables.cxx:10
std::vector< std::string > _varName
Definition: RooNDKeysPdf.h:150
#define coutI(a)
Definition: RooMsgService.h:32
std::vector< Int_t > bIdcs
Definition: RooNDKeysPdf.h:93
#define cxcoutD(a)
Definition: RooMsgService.h:80
void boxInfoInit(BoxInfo *bi, const char *rangeName, Int_t code) const
Int_t GetNrows() const
Definition: TVectorT.h:81
TVectorD * _dx
Definition: RooNDKeysPdf.h:181
Bool_t _netFluxZ
Definition: RooNDKeysPdf.h:159
#define R__ASSERT(e)
Definition: TError.h:98
std::vector< Int_t > _bIdcs
Definition: RooNDKeysPdf.h:166
Basic string class.
Definition: TString.h:137
virtual Double_t getMin(const char *name=0) const
std::vector< Double_t > _xDatLo
Definition: RooNDKeysPdf.h:156
TMatrixT< Element > & T()
Definition: TMatrixT.h:151
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1088
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
Double_t evaluate() const
do not persist
STL namespace.
RooDataSet & _data
Definition: RooNDKeysPdf.h:119
#define coutW(a)
Definition: RooMsgService.h:34
TMatrixD * _rotMat
Definition: RooNDKeysPdf.h:179
virtual ~RooNDKeysPdf()
void box(Int_t pat, Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Definition: fillpatterns.C:1
std::vector< Double_t > _x1
Definition: RooNDKeysPdf.h:154
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
Double_t _widthFactor
Definition: RooNDKeysPdf.h:121
Iterator abstract base class.
Definition: TIterator.h:32
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported...
Bool_t _verbose
Definition: RooNDKeysPdf.h:127
std::vector< Double_t > _xDatHi
Definition: RooNDKeysPdf.h:156
std::vector< Double_t > _xVarLoM3s
Definition: RooNDKeysPdf.h:163
std::vector< std::vector< Double_t > > _weights1
Definition: RooNDKeysPdf.h:142
double sqrt(double)
std::vector< Double_t > _xVarHiM3s
Definition: RooNDKeysPdf.h:163
Double_t x[n]
Definition: legend1.C:17
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()
std::list< Elem > List
Definition: ModulekNN.h:107
std::map< Int_t, Double_t > _wMap
Definition: RooNDKeysPdf.h:175
std::map< std::string, std::string >::const_iterator iter
Definition: TAlienJob.cxx:54
void calculateBandWidth() const
TIterator * createIterator(Bool_t dir=kIterForward) const
std::map< std::pair< std::string, int >, BoxInfo * > _rangeBoxInfo
Definition: RooNDKeysPdf.h:169
friend class RooArgSet
Definition: RooAbsArg.h:469
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
TVectorD * _sigmaR
Definition: RooNDKeysPdf.h:180
std::vector< Double_t > _xDatLo3s
Definition: RooNDKeysPdf.h:157
TVectorT< Double_t > TVectorD
Definition: TVectorDfwd.h:24
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:143
TMatrixT< Double_t > TMatrixD
Definition: TMatrixDfwd.h:24
void sortDataIndices(BoxInfo *bi=0) const
sort entries, as needed for loopRange()
void Print(Option_t *name="") const
Print the matrix as a table of elements.
void Print(Option_t *option="") const
Print the vector as a list of elements.
Definition: TVectorT.cxx:1360
std::vector< Double_t > xVarHiP3s
Definition: RooNDKeysPdf.h:90
void calculatePreNorm(BoxInfo *bi) const
bi->nEventsBMSW=0.
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
std::vector< Double_t > _x
Definition: RooNDKeysPdf.h:153
ROOT::R::TRInterface & r
Definition: Object.C:4
RooAbsArg * find(const char *name) const
Find object with given name in list.
virtual Double_t weight() const
Return event weight of current event.
std::vector< itVec > _sortTVIdcs
Definition: RooNDKeysPdf.h:147
Double_t Erf(Double_t x)
Computation of the error function erf(x).
Definition: TMath.cxx:187
void setOptions() const
set the configuration
virtual Int_t numEntries() const
Definition: RooAbsData.cxx:291
Double_t _maxWeight
Definition: RooNDKeysPdf.h:174
void createPdf(Bool_t firstCall=kTRUE) const
evaluation order of constructor.
TVectorT< Element > & Zero()
Set vector elements to zero.
Definition: TVectorT.cxx:450
TMarker * m
Definition: textangle.C:8
Double_t _minWeight
Definition: RooNDKeysPdf.h:173
Double_t _sigmaAvgR
Definition: RooNDKeysPdf.h:182
Double_t E()
Definition: TMath.h:54
TLine * l
Definition: textangle.C:4
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
TMatrixDSym * _covMat
Definition: RooNDKeysPdf.h:177
std::vector< Double_t > xVarHiM3s
Definition: RooNDKeysPdf.h:90
RooNDKeysPdf(const char *name, const char *title, const RooArgList &varList, RooDataSet &data, TString options="a", Double_t rho=1, Double_t nSigma=3, Bool_t rotate=kTRUE)
std::vector< Double_t > _x0
Definition: RooNDKeysPdf.h:154
Double_t _nEventsBW
Definition: RooNDKeysPdf.h:160
TString _options
Definition: RooNDKeysPdf.h:120
BoxInfo _fullBoxInfo
Definition: RooNDKeysPdf.h:170
Double_t Pi()
Definition: TMath.h:44
std::vector< Double_t > _mean
Definition: RooNDKeysPdf.h:155
std::vector< Double_t > _xVarHiP3s
Definition: RooNDKeysPdf.h:163
Double_t _nEventsBMSW
Definition: RooNDKeysPdf.h:161
double f(double x)
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
const TMatrixD & GetEigenVectors() const
std::map< Int_t, Bool_t > _bpsIdcs
Definition: RooNDKeysPdf.h:164
std::vector< Int_t > sIdcs
Definition: RooNDKeysPdf.h:92
void mirrorDataSet() const
determine mirror dataset.
TMatrixTSym< Double_t > TMatrixDSym
static RooMathCoreReg dummy
std::vector< Double_t > xVarLoP3s
Definition: RooNDKeysPdf.h:90
Double_t y[n]
Definition: legend1.C:17
std::vector< Double_t > _xVarHi
Definition: RooNDKeysPdf.h:162
std::vector< Int_t > _idx
Definition: RooNDKeysPdf.h:172
std::pair< Int_t, VecTVecDouble::iterator > itPair
Definition: RooNDKeysPdf.h:39
std::vector< Double_t > _xDatHi3s
Definition: RooNDKeysPdf.h:157
std::vector< Int_t > _sIdcs
Definition: RooNDKeysPdf.h:165
#define name(a, b)
Definition: linkTestLib0.cpp:5
std::map< Int_t, Bool_t > bpsIdcs
Definition: RooNDKeysPdf.h:91
virtual Double_t getMax(const char *name=0) const
std::vector< Int_t > bmsIdcs
Definition: RooNDKeysPdf.h:94
std::vector< Double_t > _rho
Definition: RooNDKeysPdf.h:151
void initialize() const
initialization
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral. ...
Bool_t _fixedShape
Definition: RooNDKeysPdf.h:124
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
virtual TObject * Next()=0
void calculateShell(BoxInfo *bi) const
determine points in +/- nSigma shell around the box determined by the variable ranges.
Double_t _sqrt2pi
Definition: RooNDKeysPdf.h:129
std::vector< TVectorD > _dataPtsR
Definition: RooNDKeysPdf.h:140
std::vector< std::vector< Double_t > > _dataPts
Definition: RooNDKeysPdf.h:139
std::vector< std::vector< Double_t > > _weights0
Definition: RooNDKeysPdf.h:141
void loadWeightSet() const
TIterator * _varItr
Definition: RooNDKeysPdf.h:100
float type_of_call hi(const int &, const int &)
void loadDataSet(Bool_t firstCall) const
copy the dataset and calculate some useful variables
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
Int_t getSize() const
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
double exp(double)
const Bool_t kTRUE
Definition: Rtypes.h:91
std::vector< Double_t > _xVarLoP3s
Definition: RooNDKeysPdf.h:163
Bool_t _mirror
Definition: RooNDKeysPdf.h:125
double norm(double *x, double *p)
Definition: unuranDistr.cxx:40
std::vector< Double_t > _sigma
Definition: RooNDKeysPdf.h:155
virtual const RooArgSet * get(Int_t index) const
Return RooArgSet with coordinates of event 'index'.
std::vector< Double_t > xVarLo
Definition: RooNDKeysPdf.h:89
std::vector< Int_t > _bmsIdcs
Definition: RooNDKeysPdf.h:167
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()
std::vector< Double_t > xVarLoM3s
Definition: RooNDKeysPdf.h:90
Double_t _d
Definition: RooNDKeysPdf.h:134
std::vector< Double_t > _x2
Definition: RooNDKeysPdf.h:154
ClassImp(RooNDKeysPdf) RooNDKeysPdf
Construct N-dimensional kernel estimation p.d.f.