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