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(int( _varList.size()) != 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;
138 R__ASSERT(int(_varList.size()) == rho.GetNrows());
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.size(), 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.size() != _rhoList.size())) {
187 coutE(InputArguments) << "ERROR: RooNDKeysPdf::RooNDKeysPdf() : The size of rhoList is different from varList."
188 << "Unable to create the PDF." << endl;
189 assert(_varList.size() == _rhoList.size());
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.size(), 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.size() != _rhoList.size())) {
234 coutE(InputArguments) << "ERROR: RooNDKeysPdf::RooNDKeysPdf() : The size of rhoList is different from varList."
235 << "Unable to create the PDF." << endl;
236 assert(_varList.size() == _rhoList.size());
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 }
263 _options = "m";
264 }
265
266 createPdf(true, data);
267}
268
269////////////////////////////////////////////////////////////////////////////////
270/// Backward compatibility constructor for Roo2DKeysPdf. If you are a new user,
271/// please use the first constructor form.
272
273RooNDKeysPdf::RooNDKeysPdf(const char *name, const char *title, RooAbsReal &x, RooAbsReal &y, const RooDataSet &data,
274 TString options, double rho, double nSigma, bool rotate, bool sortInput)
275 : RooAbsPdf(name, title), _varList("varList", "List of variables", this),
276 _rhoList("rhoList", "List of rho parameters", this), _options(options), _widthFactor(rho),
277 _nSigma(nSigma), _rotate(rotate), _sortInput(sortInput), _nAdpt(1)
278{
280 _varName.push_back(x.GetName());
281 _varName.push_back(y.GetName());
282
283 createPdf(true, data);
284}
285
286////////////////////////////////////////////////////////////////////////////////
287/// Constructor
288
290 : RooAbsPdf(other, name),
291 _varList("varList", this, other._varList),
292 _rhoList("rhoList", this, other._rhoList),
293 _options(other._options),
294 _widthFactor(other._widthFactor),
295 _nSigma(other._nSigma),
296 _fixedShape(other._fixedShape),
297 _mirror(other._mirror),
298 _debug(other._debug),
299 _verbose(other._verbose),
300 _nDim(other._nDim),
301 _nEvents(other._nEvents),
302 _nEventsM(other._nEventsM),
303 _nEventsW(other._nEventsW),
304 _d(other._d),
305 _n(other._n),
306 _dataPts(other._dataPts),
307 _dataPtsR(other._dataPtsR),
308 _weights0(other._weights0),
309 _weights1(other._weights1),
310 _weights(_options.Contains("a") ? &_weights1 : &_weights0),
311 _sortTVIdcs(other._sortTVIdcs),
312 _varName(other._varName),
313 _rho(other._rho),
314 _x(other._x),
315 _x0(other._x0),
316 _x1(other._x1),
317 _x2(other._x2),
318 _mean(other._mean),
319 _sigma(other._sigma),
320 _xDatLo(other._xDatLo),
321 _xDatHi(other._xDatHi),
322 _xDatLo3s(other._xDatLo3s),
323 _xDatHi3s(other._xDatHi3s),
324 _netFluxZ(other._netFluxZ),
325 _nEventsBW(other._nEventsBW),
326 _nEventsBMSW(other._nEventsBMSW),
327 _xVarLo(other._xVarLo),
328 _xVarHi(other._xVarHi),
329 _xVarLoM3s(other._xVarLoM3s),
330 _xVarLoP3s(other._xVarLoP3s),
331 _xVarHiM3s(other._xVarHiM3s),
332 _xVarHiP3s(other._xVarHiP3s),
333 _bpsIdcs(other._bpsIdcs),
334 _ibNoSort(other._ibNoSort),
335 _sIdcs(other._sIdcs),
336 _bIdcs(other._bIdcs),
337 _bmsIdcs(other._bmsIdcs),
338 _rangeBoxInfo(other._rangeBoxInfo),
339 _fullBoxInfo(other._fullBoxInfo),
340 _idx(other._idx),
341 _minWeight(other._minWeight),
342 _maxWeight(other._maxWeight),
343 _wMap(other._wMap),
344 _covMat(new TMatrixDSym(*other._covMat)),
345 _corrMat(new TMatrixDSym(*other._corrMat)),
346 _rotMat(new TMatrixD(*other._rotMat)),
347 _sigmaR(new TVectorD(*other._sigmaR)),
348 _dx(new TVectorD(*other._dx)),
349 _sigmaAvgR(other._sigmaAvgR),
350 _rotate(other._rotate),
351 _sortInput(other._sortInput),
352 _nAdpt(other._nAdpt),
353 _tracker(other._tracker != nullptr ? new RooChangeTracker(*other._tracker) : nullptr)
354{
355}
356
357////////////////////////////////////////////////////////////////////////////////
358
360{
361 if (_covMat) delete _covMat;
362 if (_corrMat) delete _corrMat;
363 if (_rotMat) delete _rotMat;
364 if (_sigmaR) delete _sigmaR;
365 if (_dx) delete _dx;
366 if (_tracker)
367 delete _tracker;
368
369 // delete all the boxinfos map
370 while ( !_rangeBoxInfo.empty() ) {
371 map<pair<string,int>,BoxInfo*>::iterator iter = _rangeBoxInfo.begin();
372 BoxInfo* box= (*iter).second;
373 _rangeBoxInfo.erase(iter);
374 delete box;
375 }
376}
377
378////////////////////////////////////////////////////////////////////////////////
379/// evaluation order of constructor.
380
381void RooNDKeysPdf::createPdf(bool firstCall, RooDataSet const& data)
382{
383 if (firstCall) {
384 // set options
385 setOptions();
386 // initialization
388 }
389
390
391 // copy dataset, calculate sigma_i's, determine min and max event weight
392 loadDataSet(firstCall, data);
393
394 // mirror dataset around dataset boundaries -- does not depend on event weights
395 if (_mirror) mirrorDataSet();
396
397 // store indices and weights of events with high enough weights
399
400 // store indices of events in variable boundaries and box shell.
401//calculateShell(&_fullBoxInfo);
402 // calculate normalization needed in analyticalIntegral()
403//calculatePreNorm(&_fullBoxInfo);
404
405 // lookup table for determining which events contribute to a certain coordinate
406 const_cast<RooNDKeysPdf*>(this)->sortDataIndices();
407
408 // determine static and/or adaptive bandwidth
410}
411
412////////////////////////////////////////////////////////////////////////////////
413/// set the configuration
414
416{
418
419 if( _options.Contains("a") ) { _weights = &_weights1; }
420 else { _weights = &_weights0; }
421 if( _options.Contains("m") ) _mirror = true;
422 else _mirror = false;
423 if( _options.Contains("d") ) _debug = true;
424 else _debug = false;
425 if( _options.Contains("v") ) { _debug = true; _verbose = true; }
426 else { _debug = false; _verbose = false; }
427
428 cxcoutD(InputArguments) << "RooNDKeysPdf::setOptions() options = " << _options
429 << "\n\tbandWidthType = " << _options.Contains("a")
430 << "\n\tmirror = " << _mirror
431 << "\n\tdebug = " << _debug
432 << "\n\tverbose = " << _verbose
433 << endl;
434
435 if (_nSigma<2.0) {
436 coutW(InputArguments) << "RooNDKeysPdf::setOptions() : Warning : nSigma = " << _nSigma << " < 2.0. "
437 << "Calculated normalization could be too large."
438 << endl;
439 }
440
441 // number of adaptive width iterations. Default is 1.
442 if (_options.Contains("a")) {
443 if (!sscanf(_options.Data(), "%d%*s", &_nAdpt)) {
444 _nAdpt = 1;
445 }
446 }
447}
448
449////////////////////////////////////////////////////////////////////////////////
450/// initialization
451
453{
454 _nDim = _varList.size();
455 _nEvents = (Int_t)data.numEntries();
457 _fixedShape= false;
458
459 _netFluxZ = false;
460 _nEventsBW = 0;
461 _nEventsBMSW = 0;
462
463 if(_nDim==0) {
464 coutE(InputArguments) << "ERROR: RooNDKeysPdf::initialize() : The observable list is empty. "
465 << "Unable to begin generating the PDF." << endl;
466 R__ASSERT (_nDim!=0);
467 }
468
469 if(_nEvents==0) {
470 coutE(InputArguments) << "ERROR: RooNDKeysPdf::initialize() : The input data set is empty. "
471 << "Unable to begin generating the PDF." << endl;
472 R__ASSERT (_nEvents!=0);
473 }
474
475 _d = static_cast<double>(_nDim);
476
477 std::vector<double> dummy(_nDim,0.);
478 _dataPts.resize(_nEvents,dummy);
479 _weights0.resize(_nEvents,dummy);
480
481 //rdh _rho.resize(_nDim,_widthFactor);
482
483 if (_widthFactor>0) { _rho.resize(_nDim,_widthFactor); }
484 // else: _rho has been provided as external input
485
486 _x.resize(_nDim,0.);
487 _x0.resize(_nDim,0.);
488 _x1.resize(_nDim,0.);
489 _x2.resize(_nDim,0.);
490
491 _mean.resize(_nDim,0.);
492 _sigma.resize(_nDim,0.);
493
494 _xDatLo.resize(_nDim,0.);
495 _xDatHi.resize(_nDim,0.);
496 _xDatLo3s.resize(_nDim,0.);
497 _xDatHi3s.resize(_nDim,0.);
498
499 boxInfoInit(&_fullBoxInfo,nullptr,0xFFFF);
500
501 _minWeight=0;
502 _maxWeight=0;
503 _wMap.clear();
504
505 _covMat = nullptr;
506 _corrMat= nullptr;
507 _rotMat = nullptr;
508 _sigmaR = nullptr;
509 _dx = new TVectorD(_nDim); _dx->Zero();
510 _dataPtsR.resize(_nEvents,*_dx);
511
512 for(unsigned int j=0; j < _varList.size(); ++j) {
513 auto& var = static_cast<const RooRealVar&>(_varList[j]);
514 _xDatLo[j] = var.getMin();
515 _xDatHi[j] = var.getMax();
516 }
517}
518
519////////////////////////////////////////////////////////////////////////////////
520/// copy the dataset and calculate some useful variables
521
522void RooNDKeysPdf::loadDataSet(bool firstCall, RooDataSet const& data)
523{
524 // first some initialization
525 _nEventsW=0.;
526
527 TMatrixD mat(_nDim,_nDim);
528 if (!_covMat) _covMat = new TMatrixDSym(_nDim);
530 if (!_rotMat) _rotMat = new TMatrixD(_nDim,_nDim);
531 if (!_sigmaR) _sigmaR = new TVectorD(_nDim);
532
533 mat.Zero();
534 _covMat->Zero();
535 _corrMat->Zero();
536 _rotMat->Zero();
537 _sigmaR->Zero();
538
539 const RooArgSet* values= data.get();
540 vector<RooRealVar*> dVars(_nDim);
541 for (Int_t j=0; j<_nDim; j++) {
542 dVars[j] = static_cast<RooRealVar*>(values->find(_varName[j].c_str()));
543 _x0[j]=_x1[j]=_x2[j]=0.;
544 }
545
546 _idx.clear();
547 for (Int_t i=0; i<_nEvents; i++) {
548
549 data.get(i); // fills dVars
550 _idx.push_back(i);
551 std::vector<double>& point = _dataPts[i];
552 TVectorD& pointV = _dataPtsR[i];
553
554 double myweight = data.weight(); // default is one?
555 if ( TMath::Abs(myweight)>_maxWeight ) { _maxWeight = TMath::Abs(myweight); }
556 _nEventsW += myweight;
557
558 for (Int_t j=0; j<_nDim; j++) {
559 for (Int_t k=0; k<_nDim; k++) {
560 mat(j,k) += dVars[j]->getVal() * dVars[k]->getVal() * myweight;
561 }
562 // only need to do once
563 if (firstCall)
564 point[j] = pointV[j] = dVars[j]->getVal();
565
566 _x0[j] += 1. * myweight;
567 _x1[j] += point[j] * myweight ;
568 _x2[j] += point[j] * point[j] * myweight ;
569 if (_x2[j]!=_x2[j]) exit(3);
570
571 // only need to do once
572 if (firstCall) {
573 if (point[j]<_xDatLo[j]) { _xDatLo[j]=point[j]; }
574 if (point[j]>_xDatHi[j]) { _xDatHi[j]=point[j]; }
575 }
576 }
577 }
578
579 _n = TMath::Power(4./(_nEventsW*(_d+2.)), 1./(_d+4.)) ;
580 // = (4/[n(dim(R) + 2)])^1/(dim(R)+4); dim(R) = 2
581 _minWeight = (0.5 - TMath::Erf(_nSigma/sqrt(2.))/2.) * _maxWeight;
582
583 for (Int_t j=0; j<_nDim; j++) {
584 _mean[j] = _x1[j]/_x0[j];
585 _sigma[j] = sqrt(_x2[j]/_x0[j]-_mean[j]*_mean[j]);
586 }
587
588 for (Int_t j=0; j<_nDim; j++) {
589 for (Int_t k=0; k<_nDim; k++) {
590 (*_covMat)(j,k) = mat(j,k)/_x0[j] - _mean[j]*_mean[k];
591 }
592 }
593
594 for (Int_t j=0; j<_nDim; j++) {
595 for (Int_t k=0; k<_nDim; k++)
596 (*_corrMat)(j,k) = (*_covMat)(j,k)/(_sigma[j]*_sigma[k]) ;
597 }
598
599 // use raw sigmas (without rho) for sigmaAvgR
600 TMatrixDSymEigen evCalculator(*_covMat);
601 TVectorD sigmaRraw = evCalculator.GetEigenValues();
602 for (Int_t j=0; j<_nDim; j++) { sigmaRraw[j] = sqrt(sigmaRraw[j]); }
603
604 _sigmaAvgR=1.;
605 for (Int_t j=0; j<_nDim; j++) { _sigmaAvgR *= sigmaRraw[j]; }
607
608 // find decorrelation matrix and eigenvalues (R)
609 if (_nDim > 1 && _rotate) {
610 // new: rotation matrix now independent of rho evaluation
611 *_rotMat = evCalculator.GetEigenVectors();
612 *_rotMat = _rotMat->T(); // transpose
613 } else {
614 TMatrixD haar(_nDim, _nDim);
615 TMatrixD unit(TMatrixD::kUnit, haar);
616 *_rotMat = unit;
617 }
618
619 // update sigmas (rho dependent)
620 updateRho();
621
622 //// rho no longer used after this.
623 //// Now set rho = 1 because sigmaR now contains rho
624 // for (Int_t j=0; j<_nDim; j++) { _rho[j] = 1.; } // reset: important!
625
626 if (_verbose) {
627 //_covMat->Print();
628 _rotMat->Print();
629 _corrMat->Print();
630 _sigmaR->Print();
631 }
632
633 if (_nDim > 1 && _rotate) {
634 // apply rotation
635 for (Int_t i = 0; i < _nEvents; i++) {
636 TVectorD &pointR = _dataPtsR[i];
637 pointR *= *_rotMat;
638 }
639 }
640
641 coutI(Contents) << "RooNDKeysPdf::loadDataSet(" << this << ")"
642 << "\n Number of events in dataset: " << _nEvents
643 << "\n Weighted number of events in dataset: " << _nEventsW << endl;
644}
645
646////////////////////////////////////////////////////////////////////////////////
647/// determine mirror dataset.
648/// mirror points are added around the physical boundaries of the dataset
649/// Two steps:
650/// 1. For each entry, determine if it should be mirrored (the mirror configuration).
651/// 2. For each mirror configuration, make the mirror points.
652
654{
655 for (Int_t j=0; j<_nDim; j++) {
656 _xDatLo3s[j] = _xDatLo[j] + _nSigma * (_n * _sigma[j]);
657 _xDatHi3s[j] = _xDatHi[j] - _nSigma * (_n * _sigma[j]);
658
659 // cout<<"xDatLo3s["<<j<<"]="<<_xDatLo3s[j]<<endl;
660 // cout<<"xDatHi3s["<<j<<"]="<<_xDatHi3s[j]<<endl;
661 }
662
663 vector<double> dummy(_nDim,0.);
664
665 // 1.
666 for (Int_t i=0; i<_nEvents; i++) {
667 vector<double>& x = _dataPts[i];
668
669 Int_t size = 1;
670 vector<vector<double> > mpoints(size,dummy);
671 vector<vector<Int_t> > mjdcs(size);
672
673 // create all mirror configurations for event i
674 for (Int_t j=0; j<_nDim; j++) {
675
676 vector<Int_t>& mjdxK = mjdcs[0];
677 vector<double>& mpointK = mpoints[0];
678
679 // single mirror *at physical boundaries*
680 if ((x[j]>_xDatLo[j] && x[j]<_xDatLo3s[j]) && x[j]<(_xDatLo[j]+_xDatHi[j])/2.) {
681 mpointK[j] = 2.*_xDatLo[j]-x[j];
682 mjdxK.push_back(j);
683 } else if ((x[j]<_xDatHi[j] && x[j]>_xDatHi3s[j]) && x[j]>(_xDatLo[j]+_xDatHi[j])/2.) {
684 mpointK[j] = 2.*_xDatHi[j]-x[j];
685 mjdxK.push_back(j);
686 }
687 }
688
689 vector<Int_t>& mjdx0 = mjdcs[0];
690 // no mirror point(s) for this event
691 if (size==1 && mjdx0.empty()) continue;
692
693 // 2.
694 // generate all mirror points for event i
695 vector<Int_t>& mjdx = mjdcs[0];
696 vector<double>& mpoint = mpoints[0];
697
698 // number of mirror points for this mirror configuration
699 Int_t eMir = 1 << mjdx.size();
700 vector<vector<double> > epoints(eMir,x);
701
702 for (Int_t m=0; m<Int_t(mjdx.size()); m++) {
703 Int_t size1 = 1 << m;
704 Int_t size2 = 1 << (m+1);
705 // copy all previous mirror points
706 for (Int_t l=size1; l<size2; ++l) {
707 epoints[l] = epoints[l-size1];
708 // fill high mirror points
709 vector<double>& epoint = epoints[l];
710 epoint[mjdx[Int_t(mjdx.size()-1)-m]] = mpoint[mjdx[Int_t(mjdx.size()-1)-m]];
711 }
712 }
713
714 // remove duplicate mirror points
715 // note that: first epoint == x
716 epoints.erase(epoints.begin());
717
718 // add mirror points of event i to total dataset
719 TVectorD pointR(_nDim);
720
721 for (Int_t m=0; m<Int_t(epoints.size()); m++) {
722 _idx.push_back(i);
723 _dataPts.push_back(epoints[m]);
724 //_weights0.push_back(_weights0[i]);
725 for (Int_t j=0; j<_nDim; j++) { pointR[j] = (epoints[m])[j]; }
726 if (_nDim > 1 && _rotate) {
727 pointR *= *_rotMat;
728 }
729 _dataPtsR.push_back(pointR);
730 }
731
732 epoints.clear();
733 mpoints.clear();
734 mjdcs.clear();
735 } // end of event loop
736
737 _nEventsM = Int_t(_dataPts.size());
738}
739
740////////////////////////////////////////////////////////////////////////////////
741
743{
744 _wMap.clear();
745
746 for (Int_t i=0; i<_nEventsM; i++) {
747 data.get(_idx[i]);
748 double myweight = data.weight();
749 //if ( TMath::Abs(myweight)>_minWeight ) {
750 _wMap[i] = myweight;
751 //}
752 }
753
754 coutI(Contents) << "RooNDKeysPdf::loadWeightSet(" << this << ") : Number of weighted events : " << _wMap.size() << endl;
755}
756
757////////////////////////////////////////////////////////////////////////////////
758/// determine points in +/- nSigma shell around the box determined by the variable
759/// ranges. These points are needed in the normalization, to determine probability
760/// leakage in and out of the box.
761
763{
764 for (Int_t j=0; j<_nDim; j++) {
765 if (bi->xVarLo[j]==_xDatLo[j] && bi->xVarHi[j]==_xDatHi[j]) {
766 bi->netFluxZ = bi->netFluxZ && true;
767 } else { bi->netFluxZ = false; }
768
769 bi->xVarLoM3s[j] = bi->xVarLo[j] - _nSigma * (_n * _sigma[j]);
770 bi->xVarLoP3s[j] = bi->xVarLo[j] + _nSigma * (_n * _sigma[j]);
771 bi->xVarHiM3s[j] = bi->xVarHi[j] - _nSigma * (_n * _sigma[j]);
772 bi->xVarHiP3s[j] = bi->xVarHi[j] + _nSigma * (_n * _sigma[j]);
773
774 //cout<<"bi->xVarLoM3s["<<j<<"]="<<bi->xVarLoM3s[j]<<endl;
775 //cout<<"bi->xVarLoP3s["<<j<<"]="<<bi->xVarLoP3s[j]<<endl;
776 //cout<<"bi->xVarHiM3s["<<j<<"]="<<bi->xVarHiM3s[j]<<endl;
777 //cout<<"bi->xVarHiM3s["<<j<<"]="<<bi->xVarHiM3s[j]<<endl;
778 }
779
780 //for (Int_t i=0; i<_nEventsM; i++) {
781 for (const auto& wMapItr : _wMap) {
782 Int_t i = wMapItr.first;
783
784 const vector<double>& x = _dataPts[i];
785 bool inVarRange(true);
786 bool inVarRangePlusShell(true);
787
788 for (Int_t j=0; j<_nDim; j++) {
789
790 if (x[j]>bi->xVarLo[j] && x[j]<bi->xVarHi[j]) {
791 inVarRange = inVarRange && true;
792 } else { inVarRange = inVarRange && false; }
793
794 if (x[j]>bi->xVarLoM3s[j] && x[j]<bi->xVarHiP3s[j]) {
795 inVarRangePlusShell = inVarRangePlusShell && true;
796 } else { inVarRangePlusShell = inVarRangePlusShell && false; }
797 }
798
799 // event in range?
800 if (inVarRange) {
801 bi->bIdcs.push_back(i);
802 }
803
804 // event in shell?
805 if (inVarRangePlusShell) {
806 bi->bpsIdcs[i] = true;
807 bool inShell(false);
808 for (Int_t j=0; j<_nDim; j++) {
809 if ((x[j]>bi->xVarLoM3s[j] && x[j]<bi->xVarLoP3s[j]) && x[j]<(bi->xVarLo[j]+bi->xVarHi[j])/2.) {
810 inShell = true;
811 } else if ((x[j]>bi->xVarHiM3s[j] && x[j]<bi->xVarHiP3s[j]) && x[j]>(bi->xVarLo[j]+bi->xVarHi[j])/2.) {
812 inShell = true;
813 }
814 }
815 if (inShell) bi->sIdcs.push_back(i); // needed for normalization
816 else {
817 bi->bmsIdcs.push_back(i); // idem
818 }
819 }
820 }
821
822 coutI(Contents) << "RooNDKeysPdf::calculateShell() : "
823 << "\n Events in shell " << bi->sIdcs.size()
824 << "\n Events in box " << bi->bIdcs.size()
825 << "\n Events in box and shell " << bi->bpsIdcs.size()
826 << endl;
827}
828
829////////////////////////////////////////////////////////////////////////////////
830///bi->nEventsBMSW=0.;
831///bi->nEventsBW=0.;
832
834{
835 // box minus shell
836 for (Int_t i=0; i<Int_t(bi->bmsIdcs.size()); i++)
837 bi->nEventsBMSW += _wMap.at(bi->bmsIdcs[i]);
838
839 // box
840 for (Int_t i=0; i<Int_t(bi->bIdcs.size()); i++)
841 bi->nEventsBW += _wMap.at(bi->bIdcs[i]);
842
843 cxcoutD(Eval) << "RooNDKeysPdf::calculatePreNorm() : "
844 << "\n nEventsBMSW " << bi->nEventsBMSW
845 << "\n nEventsBW " << bi->nEventsBW
846 << endl;
847}
848
849////////////////////////////////////////////////////////////////////////////////
850/// sort entries, as needed for loopRange()
851
853{
854 // will loop over all events by default
855 if (!_sortInput) {
856 _ibNoSort.clear();
857 for (unsigned int i = 0; i < _dataPtsR.size(); ++i) {
858 _ibNoSort[i] = true;
859 }
860 return;
861 }
862
863 itVec itrVecR;
864 vector<TVectorD>::iterator dpRItr = _dataPtsR.begin();
865 for (Int_t i = 0; dpRItr != _dataPtsR.end(); ++dpRItr, ++i) {
866 if (bi) {
867 if (bi->bpsIdcs.find(i) != bi->bpsIdcs.end()) {
868 // if (_wMap.find(i)!=_wMap.end())
869 itrVecR.push_back(itPair(i, dpRItr));
870 }
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////////////////////////////////////////////////////////////////////////////////
1150
1151Int_t RooNDKeysPdf::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
1152{
1153 if (rangeName) return 0 ;
1154
1155 Int_t code=0;
1156 if (matchArgs(allVars,analVars,RooArgSet(_varList))) { code=1; }
1157
1158 return code;
1159
1160}
1161
1162////////////////////////////////////////////////////////////////////////////////
1163
1164double RooNDKeysPdf::analyticalIntegral(Int_t code, const char* rangeName) const
1165{
1167
1168 cxcoutD(Eval) << "Calling RooNDKeysPdf::analyticalIntegral(" << GetName() << ") with code " << code
1169 << " and rangeName " << (rangeName?rangeName:"<none>") << endl;
1170
1171 // determine which observables need to be integrated over ...
1172 Int_t nComb = 1 << (_nDim);
1173 R__ASSERT(code>=1 && code<nComb) ;
1174
1175 vector<bool> doInt(_nDim,true);
1176
1177 // get BoxInfo
1178 BoxInfo* bi(nullptr);
1179
1180 if (rangeName) {
1181 string rangeNameStr(rangeName) ;
1182 bi = _rangeBoxInfo[make_pair(rangeNameStr,code)] ;
1183 if (!bi) {
1184 bi = new BoxInfo ;
1185 _rangeBoxInfo[make_pair(rangeNameStr,code)] = bi ;
1186 boxInfoInit(bi,rangeName,code);
1187 }
1188 } else bi= &_fullBoxInfo ;
1189
1190 // have boundaries changed?
1191 bool newBounds(false);
1192 for (unsigned int j=0; j < _varList.size(); ++j) {
1193 auto var = static_cast<const RooAbsRealLValue*>(_varList.at(j));
1194 if ((var->getMin(rangeName)-bi->xVarLo[j]!=0) ||
1195 (var->getMax(rangeName)-bi->xVarHi[j]!=0)) {
1196 newBounds = true;
1197 }
1198 }
1199
1200 // reset
1201 if (newBounds) {
1202 cxcoutD(Eval) << "RooNDKeysPdf::analyticalIntegral() : Found new boundaries ... " << (rangeName?rangeName:"<none>") << endl;
1203 boxInfoInit(bi,rangeName,code);
1204 }
1205
1206 // recalculates netFluxZero and nEventsIR
1207 if (!bi->filled || newBounds) {
1208 // Fill box info with contents
1209 calculateShell(bi);
1210 calculatePreNorm(bi);
1211 bi->filled = true;
1212 const_cast<RooNDKeysPdf*>(this)->sortDataIndices(bi);
1213 }
1214
1215 // first guess
1216 double norm=bi->nEventsBW;
1217
1218 if (_mirror && bi->netFluxZ) {
1219 // KEYS expression is self-normalized
1220 cxcoutD(Eval) << "RooNDKeysPdf::analyticalIntegral() : Using mirrored normalization : " << bi->nEventsBW << endl;
1221 return bi->nEventsBW;
1222 }
1223 // calculate leakage in and out of variable range box
1224 else
1225 {
1226 norm = bi->nEventsBMSW;
1227 if (norm<0.) norm=0.;
1228
1229 for (Int_t i=0; i<Int_t(bi->sIdcs.size()); ++i) {
1230 double prob=1.;
1231 const vector<double>& x = _dataPts[bi->sIdcs[i]];
1232 const vector<double>& weight = (*_weights)[_idx[bi->sIdcs[i]]];
1233
1234 vector<double> chi(_nDim,100.);
1235
1236 for (Int_t j=0; j<_nDim; j++) {
1237 if(!doInt[j]) continue;
1238
1239 if ((x[j] > bi->xVarLoM3s[j] && x[j] < bi->xVarLoP3s[j]) && x[j] < (bi->xVarLo[j] + bi->xVarHi[j]) / 2.) {
1240 chi[j] = (x[j] - bi->xVarLo[j]) / weight[j];
1241 } else if ((x[j] > bi->xVarHiM3s[j] && x[j] < bi->xVarHiP3s[j]) && x[j] > (bi->xVarLo[j] + bi->xVarHi[j]) / 2.) {
1242 chi[j] = (bi->xVarHi[j] - x[j]) / weight[j];
1243 }
1244
1245 if (chi[j] > 0) { // inVarRange
1246 prob *= (0.5 + TMath::Erf(std::abs(chi[j]) / sqrt(2.)) / 2.);
1247 } else { // outside Var range
1248 prob *= (0.5 - TMath::Erf(std::abs(chi[j]) / sqrt(2.)) / 2.);
1249 }
1250 }
1251
1252 norm += prob * _wMap.at(_idx[bi->sIdcs[i]]);
1253 }
1254
1255 cxcoutD(Eval) << "RooNDKeysPdf::analyticalIntegral() : Final normalization : " << norm << " " << bi->nEventsBW << endl;
1256 return norm;
1257 }
1258}
1259
1260
1261////////////////////////////////////////////////////////////////////////////////
1262
1264{
1265 RooArgSet varSet{varList};
1266 std::vector<RooRealVar *> varVec;
1267
1268 for (const auto var : varList) {
1269 if (!dynamic_cast<RooRealVar *>(var)) {
1270 coutE(InputArguments) << "RooNDKeysPdf::createDatasetFromHist(" << GetName() << ") WARNING: variable "
1271 << var->GetName() << " is not of type RooRealVar. Skip." << endl;
1272 continue;
1273 }
1274 varVec.push_back(static_cast<RooRealVar *>(var)); // used for setting the variables.
1275 }
1276
1277 /// determine histogram dimensionality
1278 unsigned int histndim(0);
1279 std::string classname = hist.ClassName();
1280 if (classname.find("TH1") == 0) {
1281 histndim = 1;
1282 } else if (classname.find("TH2") == 0) {
1283 histndim = 2;
1284 } else if (classname.find("TH3") == 0) {
1285 histndim = 3;
1286 }
1287 assert(histndim == varVec.size());
1288
1289 if (histndim > 3 || histndim <= 0) {
1290 coutE(InputArguments) << "RooNDKeysPdf::createDatasetFromHist(" << GetName()
1291 << ") ERROR: input histogram dimension not between [1-3]: " << histndim << endl;
1292 assert(0);
1293 }
1294
1295 /// dataset creation
1296 RooDataSet *dataFromHist = new RooDataSet("datasetFromHist", "datasetFromHist", varSet, RooFit::WeightVar());
1297
1298 /// dataset filling
1299 for (int i = 1; i <= hist.GetXaxis()->GetNbins(); ++i) {
1300 // 1 or more dimension
1301
1302 double xval = hist.GetXaxis()->GetBinCenter(i);
1303 varVec[0]->setVal(xval);
1304
1305 if (varVec.size() == 1) {
1306 double fval = hist.GetBinContent(i);
1307 dataFromHist->add(varSet, fval);
1308 } else { // 2 or more dimensions
1309
1310 for (int j = 1; j <= hist.GetYaxis()->GetNbins(); ++j) {
1311 double yval = hist.GetYaxis()->GetBinCenter(j);
1312 varVec[1]->setVal(yval);
1313
1314 if (varVec.size() == 2) {
1315 double fval = hist.GetBinContent(i, j);
1316 dataFromHist->add(varSet, fval);
1317 } else { // 3 dimensions
1318
1319 for (int k = 1; k <= hist.GetZaxis()->GetNbins(); ++k) {
1320 double zval = hist.GetZaxis()->GetBinCenter(k);
1321 varVec[2]->setVal(zval);
1322
1323 double fval = hist.GetBinContent(i, j, k);
1324 dataFromHist->add(varSet, fval);
1325 }
1326 }
1327 }
1328 }
1329 }
1330
1331 return dataFromHist;
1332}
1333
1334
1335////////////////////////////////////////////////////////////////////////////////
1336/// Return evaluated weights
1338{
1340
1341 TMatrixD mref(_nEvents, _nDim + 1);
1342
1343 cxcoutD(Eval) << "RooNDKeysPdf::getWeights() Return evaluated weights." << endl;
1344
1345 for (Int_t i = 0; i < _nEvents; ++i) {
1346 const vector<double>& x = _dataPts[i];
1347 for (Int_t j = 0; j < _nDim; j++) {
1348 mref(i, j) = x[j];
1349 }
1350
1351 const vector<double>& weight = (*_weights)[i];
1352 mref(i, _nDim) = weight[k];
1353 }
1354
1355 return mref;
1356}
1357
1358 ////////////////////////////////////////////////////////////////////////////////
1360{
1361 for (unsigned int j = 0; j < _rhoList.size(); ++j) {
1362 auto rho = static_cast<const RooAbsReal*>(_rhoList.at(j));
1363 _rho[j] = rho->getVal();
1364 }
1365
1366 if (_nDim > 1 && _rotate) {
1367 TMatrixDSym covMatRho(_nDim); // covariance matrix times rho parameters
1368 for (Int_t j = 0; j < _nDim; j++) {
1369 for (Int_t k = 0; k < _nDim; k++) {
1370 covMatRho(j, k) = (*_covMat)(j, k) * _rho[j] * _rho[k];
1371 }
1372 }
1373 // find decorrelation matrix and eigenvalues (R)
1374 TMatrixDSymEigen evCalculatorRho(covMatRho);
1375 *_sigmaR = evCalculatorRho.GetEigenValues();
1376 for (Int_t j = 0; j < _nDim; j++) {
1377 (*_sigmaR)[j] = sqrt((*_sigmaR)[j]);
1378 }
1379 } else {
1380 for (Int_t j = 0; j < _nDim; j++) {
1381 (*_sigmaR)[j] = (_sigma[j] * _rho[j]);
1382 } // * rho
1383 }
1384}
#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
TIterator Use end() or range-based loops.")
Definition RooAbsArg.h:138
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
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
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...
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
TVectorD * _dx
std::map< Int_t, bool > _ibNoSort
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
void loadDataSet(bool firstCall, RooDataSet const &data)
copy the dataset and calculate some useful variables
std::vector< TVectorD > _dataPtsR
void updateRho() const
std::vector< double > _mean
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 > _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
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< std::vector< double > > _weights1
std::vector< double > _xDatHi3s
TString _options
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
TMatrixDSym * _covMat
double _minWeight
RooChangeTracker * _tracker
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:325
TAxis * GetXaxis()
Definition TH1.h:323
TAxis * GetYaxis()
Definition TH1.h:324
virtual Double_t GetBinContent(Int_t bin) const
Return content of bin number bin.
Definition TH1.cxx:5010
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:1171
const char * Data() const
Definition TString.h:377
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:633
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