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