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