Logo ROOT   6.16/01
Reference Guide
TKDE.cxx
Go to the documentation of this file.
1// // @(#)root/hist:$Id$
2// Authors: Bartolomeu Rabacal 07/2010
3/**********************************************************************
4 * *
5 * Copyright (c) 2006 , ROOT MathLib Team *
6 * *
7 * For the licensing terms see $ROOTSYS/LICENSE. *
8 * For the list of contributors see $ROOTSYS/README/CREDITS. *
9 * *
10 **********************************************************************/
11//////////////////////////////////////////////////////////////////////////////
12/** \class TKDE
13 \ingroup Hist
14 Kernel Density Estimation class.
15 The three main references are:
16 1. "Scott DW, Multivariate Density Estimation. Theory, Practice and Visualization. New York: Wiley",
17 2. "Jann Ben - ETH Zurich, Switzerland -, Univariate kernel density estimation document for KDENS:
18 Stata module for univariate kernel density estimation."
19 3. "Hardle W, Muller M, Sperlich S, Werwatz A, Nonparametric and Semiparametric Models. Springer."
20 4. "Cranmer KS, Kernel Estimation in High-Energy
21 Physics. Computer Physics Communications 136:198-207,2001" - e-Print Archive: hep ex/0011057.
22
23 The algorithm is briefly described in (4). A binned version is also implemented to address the
24 performance issue due to its data size dependance.
25 */
26
27
28#include <functional>
29#include <algorithm>
30#include <numeric>
31#include <limits>
32#include <cassert>
33
34#include "Math/Error.h"
35#include "TMath.h"
36#include "Math/Functor.h"
37#include "Math/Integrator.h"
40#include "TGraphErrors.h"
41#include "TF1.h"
42#include "TH1.h"
43#include "TCanvas.h"
44#include "TKDE.h"
45
46
48
49class TKDE::TKernel {
50 TKDE* fKDE;
51 UInt_t fNWeights; // Number of kernel weights (bandwidth as vectorized for binning)
52 std::vector<Double_t> fWeights; // Kernel weights (bandwidth)
53public:
54 TKernel(Double_t weight, TKDE* kde);
55 void ComputeAdaptiveWeights();
57 Double_t GetWeight(Double_t x) const;
59 const std::vector<Double_t> & GetAdaptiveWeights() const;
60};
61
63 enum EIntegralResult{kNorm, kMu, kSigma2, kUnitIntegration};
64 KernelIntegrand(const TKDE* kde, EIntegralResult intRes);
66private:
67 const TKDE* fKDE;
68 EIntegralResult fIntegralResult;
69};
70
71// TKDE::TKDE(UInt_t events, const Double_t* data, Double_t xMin, Double_t xMax, const Option_t* option, Double_t rho) :
72// fData(events, 0.0),
73// fEvents(events, 0.0),
74// fPDF(0),
75// fUpperPDF(0),
76// fLowerPDF(0),
77// fApproximateBias(0),
78// fGraph(0),
79// fNewData(false),
80// fUseMinMaxFromData((xMin >= xMax)),
81// fNBins(events < 10000 ? 100: events / 10),
82// fNEvents(events),
83// fUseBinsNEvents(10000),
84// fMean(0.0),
85// fSigma(0.0),
86// fXMin(xMin),
87// fXMax(xMax),
88// fAdaptiveBandwidthFactor(1.0),
89// fCanonicalBandwidths(std::vector<Double_t>(kTotalKernels, 0.0)),
90// fKernelSigmas2(std::vector<Double_t>(kTotalKernels, -1.0)),
91// fSettedOptions(std::vector<Bool_t>(4, kFALSE))
92// {
93// //Class constructor
94// SetOptions(option, rho);
95// CheckOptions();
96// SetMirror();
97// SetUseBins();
98// SetKernelFunction();
99// SetData(data);
100// SetCanonicalBandwidths();
101// SetKernelSigmas2();
102// SetKernel();
103// }
104
106 //Class destructor
107 if (fPDF) delete fPDF;
108 if (fUpperPDF) delete fUpperPDF;
109 if (fLowerPDF) delete fLowerPDF;
110 if (fGraph) delete fGraph;
113 if (fKernel) delete fKernel;
114}
115
116void TKDE::Instantiate(KernelFunction_Ptr kernfunc, UInt_t events, const Double_t* data, const Double_t* dataWeights, Double_t xMin, Double_t xMax, const Option_t* option, Double_t rho) {
117 // Template's constructor surrogate
118 fData = std::vector<Double_t>(events, 0.0);
119 fEvents = std::vector<Double_t>(events, 0.0);
120 fPDF = 0;
121 fKernel = 0;
122 fKernelFunction = 0;
123 fUpperPDF = 0;
124 fLowerPDF = 0;
126 fGraph = 0;
127 fNewData = false;
128 fUseMirroring = false; fMirrorLeft = false; fMirrorRight = false;
129 fAsymLeft = false; fAsymRight = false;
130 fNBins = events < 10000 ? 100 : events / 10;
131 fNEvents = events;
132 fUseBinsNEvents = 10000;
133 fMean = 0.0;
134 fSigma = 0.0;
135 fXMin = xMin;
136 fXMax = xMax;
138 fSumOfCounts = 0;
140 fRho = rho;
141 fWeightSize = 0;
142 fCanonicalBandwidths = std::vector<Double_t>(kTotalKernels, 0.0);
143 fKernelSigmas2 = std::vector<Double_t>(kTotalKernels, -1.0);
144 fSettedOptions = std::vector<Bool_t>(4, kFALSE);
145 SetOptions(option, rho);
147 SetMirror();
148 SetUseBins();
149 SetData(data, dataWeights);
150 SetKernelFunction(kernfunc);
151}
152
153void TKDE::SetOptions(const Option_t* option, Double_t rho) {
154 //Sets User defined construction options
155 TString opt = option;
156 opt.ToLower();
157 std::string options = opt.Data();
158 size_t numOpt = 4;
159 std::vector<std::string> voption(numOpt, "");
160 for (std::vector<std::string>::iterator it = voption.begin(); it != voption.end() && !options.empty(); ++it) {
161 size_t pos = options.find_last_of(';');
162 if (pos == std::string::npos) {
163 *it = options;
164 break;
165 }
166 *it = options.substr(pos + 1);
167 options = options.substr(0, pos);
168 }
169 for (std::vector<std::string>::iterator it = voption.begin(); it != voption.end(); ++it) {
170 size_t pos = (*it).find(':');
171 if (pos != std::string::npos) {
172 GetOptions((*it).substr(0, pos) , (*it).substr(pos + 1));
173 }
174 }
176 fRho = rho;
177}
178
179void TKDE::SetDrawOptions(const Option_t* option, TString& plotOpt, TString& drawOpt) {
180 // Sets User defined drawing options
181 size_t numOpt = 2;
182 std::string options = TString(option).Data();
183 std::vector<std::string> voption(numOpt, "");
184 for (std::vector<std::string>::iterator it = voption.begin(); it != voption.end() && !options.empty(); ++it) {
185 size_t pos = options.find_last_of(';');
186 if (pos == std::string::npos) {
187 *it = options;
188 break;
189 }
190 *it = options.substr(pos + 1);
191 options = options.substr(0, pos);
192 }
193 Bool_t foundPlotOPt = kFALSE;
194 Bool_t foundDrawOPt = kFALSE;
195 for (std::vector<std::string>::iterator it = voption.begin(); it != voption.end() && !options.empty(); ++it) {
196 size_t pos = (*it).find(':');
197 if (pos == std::string::npos) break;
198 TString optionType = (*it).substr(0, pos);
199 TString optionInstance = (*it).substr(pos + 1);
200 optionType.ToLower();
201 optionInstance.ToLower();
202 if (optionType.Contains("plot")) {
203 foundPlotOPt = kTRUE;
204 if (optionInstance.Contains("estimate") || optionInstance.Contains("errors") || optionInstance.Contains("confidenceinterval"))
205 plotOpt = optionInstance;
206 else
207 this->Warning("SetDrawOptions", "Unknown plotting option: setting to KDE estimate plot.");
208 } else if (optionType.Contains("drawoptions")) {
209 foundDrawOPt = kTRUE;
210 drawOpt = optionInstance;
211 }
212 }
213 if (!foundPlotOPt) {
214 this->Warning("SetDrawOptions", "No plotting option: setting to KDE estimate plot.");
215 plotOpt = "estimate";
216 }
217 if (!foundDrawOPt) {
218 this->Warning("SetDrawOptions", "No drawing options: setting to default ones.");
219 drawOpt = "apl4";
220 }
221}
222
223void TKDE::GetOptions(std::string optionType, std::string option) {
224 // Gets User defined KDE construction options
225 if (optionType.compare("kerneltype") == 0) {
227 if (option.compare("gaussian") == 0) {
229 } else if (option.compare("epanechnikov") == 0) {
231 } else if (option.compare("biweight") == 0) {
233 } else if (option.compare("cosinearch") == 0) {
235 } else if (option.compare("userdefined") == 0) {
237 } else {
238 this->Warning("GetOptions", "Unknown kernel type option: setting to Gaussian");
240 }
241 } else if (optionType.compare("iteration") == 0) {
243 if (option.compare("adaptive") == 0) {
245 } else if (option.compare("fixed") == 0) {
247 } else {
248 this->Warning("GetOptions", "Unknown iteration option: setting to Adaptive");
250 }
251 } else if (optionType.compare("mirror") == 0) {
253 if (option.compare("nomirror") == 0) {
255 } else if (option.compare("mirrorleft") == 0) {
257 } else if (option.compare("mirrorright") == 0) {
259 } else if (option.compare("mirrorboth") == 0) {
261 } else if (option.compare("mirrorasymleft") == 0) {
263 } else if (option.compare("mirrorasymleftright") == 0) {
265 } else if (option.compare("mirrorasymright") == 0) {
267 } else if (option.compare("mirrorleftasymright") == 0) {
269 } else if (option.compare("mirrorasymboth") == 0) {
271 } else {
272 this->Warning("GetOptions", "Unknown mirror option: setting to NoMirror");
274 }
275 } else if (optionType.compare("binning") == 0) {
277 if (option.compare("unbinned") == 0) {
279 } else if (option.compare("relaxedbinning") == 0) {
281 } else if (option.compare("forcedbinning") == 0) {
283 } else {
284 this->Warning("GetOptions", "Unknown binning option: setting to RelaxedBinning");
286 }
287 }
288}
289
291 // Sets missing construction options to default ones
292 if (!fSettedOptions[0]) {
294 }
295 if (!fSettedOptions[1]) {
297 }
298 if (!fSettedOptions[2]) {
300 }
301 if (!fSettedOptions[3]) {
303 }
304}
305
306void TKDE::CheckOptions(Bool_t isUserDefinedKernel) {
307 // Sets User global options
308 if (!(isUserDefinedKernel) && !(fKernelType >= kGaussian && fKernelType < kUserDefined)) {
309 Error("CheckOptions", "Illegal user kernel type input! Use template constructor for user defined kernel.");
310 }
311 if (fIteration != kAdaptive && fIteration != kFixed) {
312 Warning("CheckOptions", "Illegal user iteration type input - use default value !");
314 }
315 if (!(fMirror >= kNoMirror && fMirror <= kMirrorAsymBoth)) {
316 Warning("CheckOptions", "Illegal user mirroring type input - use default value !");
318 }
319 if (!(fBinning >= kUnbinned && fBinning <= kForcedBinning)) {
320 Warning("CheckOptions", "Illegal user binning type input - use default value !");
322 }
323 if (fRho <= 0.0) {
324 Warning("CheckOptions", "Tuning factor rho cannot be non-positive - use default value !");
325 fRho = 1.0;
326 }
327}
328
330 // Sets User option for the choice of kernel estimator
332 delete fKernelFunction;
333 fKernelFunction = 0;
334 }
335 fKernelType = kern;
336 CheckOptions();
338}
339
341 // Sets User option for fixed or adaptive iteration
342 fIteration = iter;
343 CheckOptions();
344 SetKernel();
345}
346
347
349 // Sets User option for mirroring the data
350 fMirror = mir;
351 CheckOptions();
352 SetMirror();
353 if (fUseMirroring) {
355 }
356 SetKernel();
357}
358
359
361 // Sets User option for binning the weights
362 fBinning = bin;
363 CheckOptions();
364 SetUseBins();
365 SetKernel();
366}
367
369 // Sets User option for number of bins
370 if (!nbins) {
371 Error("SetNBins", "Number of bins must be greater than zero.");
372 return;
373 }
374 fNBins = nbins;
378
379 if (fBinning == kUnbinned) {
380 Warning("SetNBins", "Bin type using SetBinning must set for using a binned evaluation");
381 return;
382 }
383
384 SetKernel();
385}
386
388 // Sets User option for the minimum number of events for allowing automatic binning
389 fUseBinsNEvents = nEvents;
390 SetUseBins();
391 SetKernel();
392}
393
395 // Factor which can be used to tune the smoothing.
396 // It is used as multiplicative factor for the fixed and adaptive bandwidth.
397 // A value < 1 will reproduce better the tails but oversmooth the peak
398 // while a factor > 1 will overestimate the tail
399 fRho = rho;
400 CheckOptions();
401 SetKernel();
402}
403
405 // Sets minimum range value and maximum range value
406 if (xMin >= xMax) {
407 Error("SetRange", "Minimum range cannot be bigger or equal than the maximum range! Present range values remain the same.");
408 return;
409 }
410 fXMin = xMin;
411 fXMax = xMax;
412 fUseMinMaxFromData = false;
413 SetKernel();
414}
415
416// private methods
417
419 // Sets User option for using binned weights
420 switch (fBinning) {
421 default:
422 case kRelaxedBinning:
423 if (fNEvents >= fUseBinsNEvents) {
424 fUseBins = kTRUE;
425 } else {
427 }
428 break;
429 case kForcedBinning:
430 fUseBins = kTRUE;
431 break;
432 case kUnbinned:
434 }
435}
436
438 // Sets the mirroring
444}
445
446void TKDE::SetData(const Double_t* data, const Double_t* wgts) {
447 // Sets the data events input sample or bin centres for binned option and computes basic estimators
448 if (!data) {
449 if (fNEvents) fData.reserve(fNEvents);
450 return;
451 }
452 fEvents.assign(data, data + fNEvents);
453 if (wgts) fEventWeights.assign(wgts, wgts + fNEvents);
454
455 if (fUseMinMaxFromData) {
456 fXMin = *std::min_element(fEvents.begin(), fEvents.end());
457 fXMax = *std::max_element(fEvents.begin(), fEvents.end());
458 }
459
460 if (fUseBins) {
461 if (fNBins >= fNEvents) {
462 this->Warning("SetData", "Default number of bins is greater or equal to number of events. Use SetNBins(UInt_t) to set the appropriate number of bins");
463 }
466 } else {
468 fData = fEvents;
469 }
470 // to set fBinCOunt and fSumOfCounts
472
473
475 if (fUseMirroring) {
477 }
478}
479
481 // re-initialize when new data have been filled in TKDE
482 // re-compute kernel quantities and mean and sigma
483 fNewData = false;
484 fEvents = fData;
485 if (fUseMinMaxFromData) {
486 fXMin = *std::min_element(fEvents.begin(), fEvents.end());
487 fXMax = *std::max_element(fEvents.begin(), fEvents.end());
488 }
490 // if (fUseBins) {
491 // } // bin usage is not supported in this case
492 //
494 if (fUseMirroring) {
496 }
497 SetKernel();
498}
499
501 // Mirrors the data
502 std::vector<Double_t> originalEvents = fEvents;
503 std::vector<Double_t> originalWeights = fEventWeights;
504 if (fMirrorLeft) {
505 fEvents.resize(2 * fNEvents, 0.0);
506 transform(fEvents.begin(), fEvents.begin() + fNEvents, fEvents.begin() + fNEvents,
507 std::bind(std::minus<Double_t>(), 2 * fXMin, std::placeholders::_1));
508 }
509 if (fMirrorRight) {
510 fEvents.resize((fMirrorLeft + 2) * fNEvents, 0.0);
511 transform(fEvents.begin(), fEvents.begin() + fNEvents, fEvents.begin() + (fMirrorLeft + 1) * fNEvents,
512 std::bind(std::minus<Double_t>(), 2 * fXMax, std::placeholders::_1));
513 }
514 if (!fEventWeights.empty() && (fMirrorLeft || fMirrorRight)) {
515 // copy weights too
516 fEventWeights.insert(fEventWeights.end(), fEventWeights.begin(), fEventWeights.end() );
517 }
518
519 if(fUseBins) {
524 } else {
525 fData = fEvents;
526 }
528
529 fEvents = originalEvents;
530 fEventWeights = originalWeights;
531}
532
534 // Computes input data's mean
535 fMean = std::accumulate(fEvents.begin(), fEvents.end(), 0.0) / fEvents.size();
536}
537
539 // Computes input data's sigma
540 fSigma = std::sqrt(1. / (fEvents.size() - 1.) * (std::inner_product(fEvents.begin(), fEvents.end(), fEvents.begin(), 0.0) - fEvents.size() * std::pow(fMean, 2.)));
541 fSigmaRob = std::min(fSigma, R / 1.349); // Sigma's robust estimator
542}
543
545 // Sets the kernel density estimator
546 UInt_t n = fData.size();
547 if (n == 0) return;
548 // Optimal bandwidth (Silverman's rule of thumb with assumed Gaussian density)
549 Double_t weight = fCanonicalBandwidths[kGaussian] * fSigmaRob * std::pow(3. / (8. * std::sqrt(M_PI)) * n, -0.2);
551 if (fKernel) delete fKernel;
552 fKernel = new TKernel(weight, this);
553 if (fIteration == kAdaptive) {
554 fKernel->ComputeAdaptiveWeights();
555 }
556}
557
559
560 assert(fKernelFunction == 0); // to avoid memory leaks
561 switch (fKernelType) {
562 case kGaussian :
564 break;
565 case kEpanechnikov :
567 break;
568 case kBiweight :
570 break;
571 case kCosineArch :
573 break;
574 case kUserDefined :
575 fKernelFunction = kernfunc;
577 break;
578 case kTotalKernels :
579 default:
580 /// for user defined kernels
581 fKernelFunction = kernfunc;
583 }
584
585 if (fKernelType == kUserDefined) {
586 if (fKernelFunction) {
590 }
591 else {
592 Error("SetKernelFunction", "User kernel function is not defined !");
593 return;
594 }
595 }
596 assert(fKernelFunction);
599 SetKernel();
600}
601
603 // Sets the canonical bandwidths according to the kernel type
604 fCanonicalBandwidths[kGaussian] = 0.7764; // Checked in Mathematica
605 fCanonicalBandwidths[kEpanechnikov] = 1.7188; // Checked in Mathematica
606 fCanonicalBandwidths[kBiweight] = 2.03617; // Checked in Mathematica
607 fCanonicalBandwidths[kCosineArch] = 1.7663; // Checked in Mathematica
608 fCanonicalBandwidths[kUserDefined] = 1.0; // To be Checked
609}
610
612 // Sets the kernel sigmas2 according to the kernel type
614 fKernelSigmas2[kEpanechnikov] = 1.0 / 5.0;
615 fKernelSigmas2[kBiweight] = 1.0 / 7.0;
616 fKernelSigmas2[kCosineArch] = 1.0 - 8.0 / std::pow(M_PI, 2);
617}
618
620 // Returns the PDF estimate as a function sampled in npx between xMin and xMax
621 // the KDE is not re-normalized to the xMin/xMax range.
622 // The user manages the returned function
623 // For getting a non-sampled TF1, one can create a TF1 directly from the TKDE by doing
624 // TF1 * f1 = new TF1("f1",kde,xMin,xMax,0);
625 return GetKDEFunction(npx,xMin,xMax);
626}
627
628TF1* TKDE::GetUpperFunction(Double_t confidenceLevel, UInt_t npx, Double_t xMin, Double_t xMax) {
629 // Returns the PDF upper estimate (upper confidence interval limit)
630 return GetPDFUpperConfidenceInterval(confidenceLevel,npx,xMin,xMax);
631}
632
633TF1* TKDE::GetLowerFunction(Double_t confidenceLevel, UInt_t npx, Double_t xMin, Double_t xMax) {
634 // Returns the PDF lower estimate (lower confidence interval limit)
635 return GetPDFLowerConfidenceInterval(confidenceLevel,npx,xMin,xMax);
636}
637
639 // Returns the PDF estimate bias
640 return GetKDEApproximateBias(npx,xMin,xMax);
641}
642
644 // Fills data member with User input data event for the unbinned option
645 if (fUseBins) {
646 this->Warning("Fill", "Cannot fill data with data binned option. Data input ignored.");
647 return;
648 }
649 fData.push_back(data);
650 fNEvents++;
651 fNewData = kTRUE;
652}
653
655 // Fills data member with User input data event for the unbinned option
656 if (fUseBins) {
657 this->Warning("Fill", "Cannot fill data with data binned option. Data input ignored.");
658 return;
659 }
660 fData.push_back(data); // should not be here fEvent ??
661 fEventWeights.push_back(weight);
662 fNEvents++;
663 fNewData = kTRUE;
664}
665
667 // The class's unary function: returns the kernel density estimate
668 return (*this)(*x);
669}
670
672 // The class's unary function: returns the kernel density estimate
673 if (fNewData) (const_cast<TKDE*>(this))->InitFromNewData();
674 return (*fKernel)(x);
675}
676
678 // return the mean of the data
679 if (fNewData) (const_cast<TKDE*>(this))->InitFromNewData();
680 return fMean;
681}
682
684 // return the standard deviation of the data
685 if (fNewData) (const_cast<TKDE*>(this))->InitFromNewData();
686 return fSigma;
687}
688
690 // Returns the Root Asymptotic Mean Integrated Squared Error according to Silverman's rule of thumb with assumed Gaussian density
692 return std::sqrt(result);
693}
694
695TKDE::TKernel::TKernel(Double_t weight, TKDE* kde) :
696// Internal class constructor
697fKDE(kde),
698fNWeights(kde->fData.size()),
699fWeights(fNWeights, weight)
700{}
701
702void TKDE::TKernel::ComputeAdaptiveWeights() {
703 // Gets the adaptive weights (bandwidths) for TKernel internal computation
704 std::vector<Double_t> weights = fWeights;
705 Double_t minWeight = weights[0] * 0.05;
706 unsigned int n = fKDE->fData.size();
707 assert( n == weights.size() );
708 bool useDataWeights = (fKDE->fBinCount.size() == n);
709 Double_t f = 0.0;
710 for (unsigned int i = 0; i < n; ++i) {
711// for (; weight != weights.end(); ++weight, ++data, ++dataW) {
712 if (useDataWeights && fKDE->fBinCount[i] <= 0) continue; // skip negative or null weights
713 f = (*fKDE->fKernel)(fKDE->fData[i]);
714 if (f <= 0)
715 fKDE->Warning("ComputeAdativeWeights","function value is zero or negative for x = %f w = %f",
716 fKDE->fData[i],(useDataWeights) ? fKDE->fBinCount[i] : 1.);
717 weights[i] = std::max(weights[i] /= std::sqrt(f), minWeight);
718 fKDE->fAdaptiveBandwidthFactor += std::log(f);
719 //printf("(f = %f w = %f af = %f ),",f,*weight,fKDE->fAdaptiveBandwidthFactor);
720 }
721 Double_t kAPPROX_GEO_MEAN = 0.241970724519143365; // 1 / TMath::Power(2 * TMath::Pi(), .5) * TMath::Exp(-.5). Approximated geometric mean over pointwise data (the KDE function is substituted by the "real Gaussian" pdf) and proportional to sigma. Used directly when the mirroring is enabled, otherwise computed from the data
722 fKDE->fAdaptiveBandwidthFactor = fKDE->fUseMirroring ? kAPPROX_GEO_MEAN / fKDE->fSigmaRob : std::sqrt(std::exp(fKDE->fAdaptiveBandwidthFactor / fKDE->fData.size()));
723 transform(weights.begin(), weights.end(), fWeights.begin(),
724 std::bind(std::multiplies<Double_t>(), std::placeholders::_1, fKDE->fAdaptiveBandwidthFactor));
725 //printf("adaptive bandwidth factor % f weight 0 %f , %f \n",fKDE->fAdaptiveBandwidthFactor, weights[0],fWeights[0] );
726}
727
728Double_t TKDE::TKernel::GetWeight(Double_t x) const {
729 // Returns the bandwidth
730 return fWeights[fKDE->Index(x)];
731}
732
734 // Returns the bins' centres from the data for using with the binned option
735 fData.assign(fNBins, 0.0);
736 Double_t binWidth((xmax - xmin) / fNBins);
737 for (UInt_t i = 0; i < fNBins; ++i) {
738 fData[i] = xmin + (i + 0.5) * binWidth;
739 }
740}
741
743 // Returns the bins' count from the data for using with the binned option
744 // or set the bin count to the weights in case of weighted data
745 if (fUseBins) {
746 fBinCount.resize(fNBins);
747 fSumOfCounts = 0;
748 // case of weighted events
749 if (!fEventWeights.empty() ) {
750 for (UInt_t i = 0; i < fNEvents; ++i) {
751 if (fEvents[i] >= fXMin && fEvents[i] < fXMax) {
754 //printf("sum of counts %f - bin count %d - %f \n",fSumOfCounts, Index(fEvents[i]), fBinCount[Index(fEvents[i])] );
755 }
756 }
757 }
758 // case of unweighted data
759 else {
760 for (UInt_t i = 0; i < fNEvents; ++i) {
761 if (fEvents[i] >= fXMin && fEvents[i] < fXMax) {
762 fBinCount[Index(fEvents[i])] += 1;
763 fSumOfCounts += 1;
764 }
765 }
766 }
767 }
768 else if (!fEventWeights.empty() ) {
770 fSumOfCounts = 0;
771 for (UInt_t i = 0; i < fNEvents; ++i)
773 }
774 else {
776 fBinCount.clear();
777 }
778}
779
780void TKDE::Draw(const Option_t* opt) {
781 // Draws either the KDE functions or its errors
782 // Possible options:
783 // "" (default) - draw just the kde
784 // "same" draw on top of existing pad
785 // "Errors" draw a TGraphErrors with the point and errors
786 // "confidenceinterval" draw KDE + conf interval functions (default is 95%)
787 // "confidenceinterval@0.90" draw KDE + conf interval functions at 90%
788 // Extra options can be passed in opt for drawing the TF1 or the TGraph
789 //
790 //NOTE: The functions GetDrawnFunction(), GetDrawnUpperFunction(), GetDrawnLowerFunction()
791 // and GetGraphWithErrors() return the corresponding drawn objects (which are maneged by the TKDE)
792 // They can be used to changes style, color, etc...
793
794 // TString plotOpt = "";
795 // TString drawOpt = "";
796 // LM : this is too complicates - skip it - not needed for just
797 // three options
798 // SetDrawOptions(opt, plotOpt, drawOpt);
799 TString plotOpt = opt;
800 plotOpt.ToLower();
801 TString drawOpt = plotOpt;
802 if(gPad && !plotOpt.Contains("same")) {
803 gPad->Clear();
804 }
805 if (plotOpt.Contains("errors")) {
806 drawOpt.ReplaceAll("errors","");
807 DrawErrors(drawOpt);
808 }
809 else if (plotOpt.Contains("confidenceinterval") ||
810 plotOpt.Contains("confinterval")) {
811 // parse level option
812 drawOpt.ReplaceAll("confidenceinterval","");
813 drawOpt.ReplaceAll("confinterval","");
814 Double_t level = 0.95;
815 const char * s = strstr(plotOpt.Data(),"interval@");
816 // coverity [secure_coding : FALSE]
817 if (s != 0) sscanf(s,"interval@%lf",&level);
818 if((level <= 0) || (level >= 1)) {
819 Warning("Draw","given confidence level %.3lf is invalid - use default 0.95",level);
820 level = 0.95;
821 }
822 DrawConfidenceInterval(drawOpt,level);
823 }
824 else {
825 if (fPDF) delete fPDF;
827 fPDF->Draw(drawOpt);
828 }
829}
830
831void TKDE::DrawErrors(TString& drawOpt) {
832 // Draws a TGraphErrors for the KDE errors
833 if (fGraph) delete fGraph;
835 fGraph->Draw(drawOpt.Data());
836}
837
839 if (xmin>= xmax) { xmin = fXMin; xmax = fXMax; }
840 // return a TGraphErrors for the KDE errors
841 UInt_t n = npx;
842 Double_t* x = new Double_t[n + 1];
843 Double_t* ex = new Double_t[n + 1];
844 Double_t* y = new Double_t[n + 1];
845 Double_t* ey = new Double_t[n + 1];
846 for (UInt_t i = 0; i <= n; ++i) {
847 x[i] = xmin + i * (xmax - xmin) / n;
848 y[i] = (*this)(x[i]);
849 ex[i] = 0;
850 ey[i] = this->GetError(x[i]);
851 }
852 TGraphErrors* ge = new TGraphErrors(n, &x[0], &y[0], &ex[0], &ey[0]);
853 ge->SetName("kde_graph_error");
854 ge->SetTitle("Errors");
855 delete [] x;
856 delete [] ex;
857 delete [] y;
858 delete [] ey;
859 return ge;
860}
861
862void TKDE::DrawConfidenceInterval(TString& drawOpt,double cl) {
863 // Draws the KDE and its confidence interval
864 GetKDEFunction()->Draw(drawOpt.Data());
866 upper->SetLineColor(kBlue);
867 upper->Draw(("same" + drawOpt).Data());
869 lower->SetLineColor(kRed);
870 lower->Draw(("same" + drawOpt).Data());
871 if (fUpperPDF) delete fUpperPDF;
872 if (fLowerPDF) delete fLowerPDF;
873 fUpperPDF = upper;
874 fLowerPDF = lower;
875}
876
878 // Returns the bandwidth for the non adaptive KDE
879 Double_t result = -1.;
881 this->Warning("GetFixedWeight()", "Fixed iteration option not enabled. Returning %f.", result);
882 } else {
883 result = fKernel->GetFixedWeight();
884 }
885 return result;
886}
887
889 // Returns the bandwidths for the adaptive KDE
891 this->Warning("GetFixedWeight()", "Adaptive iteration option not enabled. Returning a NULL pointer<");
892 return 0;
893 }
894 if (fNewData) (const_cast<TKDE*>(this))->InitFromNewData();
895 return &(fKernel->GetAdaptiveWeights()).front();
896}
897
898Double_t TKDE::TKernel::GetFixedWeight() const {
899 // Returns the bandwidth for the non adaptive KDE
900 return fWeights[0];
901}
902
903const std::vector<Double_t> & TKDE::TKernel::GetAdaptiveWeights() const {
904 // Returns the bandwidth for the non adaptive KDE
905 return fWeights;
906}
907
909 // The internal class's unary function: returns the kernel density estimate
910 Double_t result(0.0);
911 UInt_t n = fKDE->fData.size();
912 // case of bins or weighted data
913 Bool_t useBins = (fKDE->fBinCount.size() == n);
914 Double_t nSum = (useBins) ? fKDE->fSumOfCounts : fKDE->fNEvents;
915 // double dmin = 1.E10;
916 // double xmin,bmin,wmin;
917 for (UInt_t i = 0; i < n; ++i) {
918 Double_t binCount = (useBins) ? fKDE->fBinCount[i] : 1.0;
919 //printf("data point %i %f count %f weight % f result % f\n",i,fKDE->fData[i],binCount,fWeights[i], result);
920 result += binCount / fWeights[i] * (*fKDE->fKernelFunction)((x - fKDE->fData[i]) / fWeights[i]);
921 if (fKDE->fAsymLeft) {
922 result -= binCount / fWeights[i] * (*fKDE->fKernelFunction)((x - (2. * fKDE->fXMin - fKDE->fData[i])) / fWeights[i]);
923 }
924 if (fKDE->fAsymRight) {
925 result -= binCount / fWeights[i] * (*fKDE->fKernelFunction)((x - (2. * fKDE->fXMax - fKDE->fData[i])) / fWeights[i]);
926 }
927 // if ( TMath::IsNaN(result) ) {
928 // printf("event %i count %f weight %f data % f x %f \n",i,binCount,fWeights[i],fKDE->fData[i],x );
929 // }
930 // if ( result <= 0 ) {
931 // printf("event %i count %f weight %f data % f x %f \n",i,binCount,fWeights[i],fKDE->fData[i],x );
932 // }
933 // if (std::abs(x - fKDE->fData[i]) < dmin ) {
934 // xmin = x;
935 // bmin = binCount;
936 // wmin = fWeights[i];
937 // dmin = std::abs(x - fKDE->fData[i]);
938 // }
939 }
940 if ( TMath::IsNaN(result) ) {
941 fKDE->Warning("operator()","Result is NaN for x %f \n",x);
942 //xmin % f , %f, %f \n",result,x,xmin,bmin,wmin );
943 }
944 return result / nSum;
945}
946
948 // Returns the indices (bins) for the binned weights
949 Int_t bin = Int_t((x - fXMin) * fWeightSize);
950 if (bin == (Int_t)fData.size()) return --bin;
952 bin += fData.size() / (fMirrorLeft + fMirrorRight + 1);
953 }
954 if (bin > (Int_t)fData.size()) {
955 return (Int_t)(fData.size()) - 1;
956 } else if (bin <= 0) {
957 return 0;
958 }
959 return bin;
960}
961
963 // Returns the pointwise upper estimated density
964 Double_t f = (*this)(x);
966 Double_t prob = 1. - (1.-*p)/2; // this is 1.-alpha/2
968 return f + z * sigma;
969}
970
972 // Returns the pointwise lower estimated density
973 Double_t f = (*this)(x);
975 Double_t prob = (1.-*p)/2; // this is alpha/2
977 return f + z * sigma;
978}
979
980
982 // Returns the pointwise approximate estimated density bias
985 rd.SetFunction(kern);
986 Double_t df2 = rd.Derivative2(x);
987 Double_t weight = fKernel->GetWeight(x); // Bandwidth
988 return 0.5 * fKernelSigmas2[fKernelType] * std::pow(weight, 2) * df2;
989}
991 // Returns the pointwise sigma of estimated density
992 Double_t kernelL2Norm = ComputeKernelL2Norm();
993 Double_t f = (*this)(x);
994 Double_t weight = fKernel->GetWeight(x); // Bandwidth
995 Double_t result = f * kernelL2Norm / (fNEvents * weight);
996 return std::sqrt(result);
997}
998
1000 // Checks if kernel has unit integral, mu = 0 and positive finite sigma conditions
1001 Double_t valid = kTRUE;
1003 valid = valid && unity == 1.;
1004 if (!valid) {
1005 Error("CheckKernelValidity", "Kernel's integral is %f",unity);
1006 }
1008 valid = valid && mu == 0.;
1009 if (!valid) {
1010 Error("CheckKernelValidity", "Kernel's mu is %f" ,mu);
1011 }
1012 Double_t sigma2 = ComputeKernelSigma2();
1013 valid = valid && sigma2 > 0 && sigma2 != std::numeric_limits<Double_t>::infinity();
1014 if (!valid) {
1015 Error("CheckKernelValidity", "Kernel's sigma2 is %f",sigma2);
1016 }
1017 if (!valid) {
1018 Error("CheckKernelValidity", "Validation conditions: the kernel's integral must be 1, the kernel's mu must be zero and the kernel's sigma2 must be finite positive to be a suitable kernel.");
1019 //exit(EXIT_FAILURE);
1020 }
1021}
1022
1024 // Computes the kernel's L2 norm
1027 ig.SetFunction(kernel);
1028 Double_t result = ig.Integral();
1029 return result;
1030}
1031
1033 // Computes the kernel's sigma squared
1035 KernelIntegrand kernel(this, TKDE::KernelIntegrand::kSigma2);
1036 ig.SetFunction(kernel);
1037 Double_t result = ig.Integral();
1038 return result;
1039}
1040
1042 // Computes the kernel's mu
1044 KernelIntegrand kernel(this, TKDE::KernelIntegrand::kMu);
1045 ig.SetFunction(kernel);
1046 Double_t result = ig.Integral();
1047 return result;
1048}
1049
1051 // Computes the kernel's integral which ought to be unity
1053 KernelIntegrand kernel(this, TKDE::KernelIntegrand::kUnitIntegration);
1054 ig.SetFunction(kernel);
1055 Double_t result = ig.Integral();
1056 return result;
1057}
1058
1060 /// in case of weights use
1061 if (!fEventWeights.empty() ) {
1062 // weighted data
1063 double x1 = fXMin - 0.001*(fXMax-fXMin);
1064 double x2 = fXMax + 0.001*(fXMax-fXMin);
1065 TH1D h1("temphist","", 500, x1, x2);
1066 h1.FillN(fEvents.size(), fEvents.data(), fEventWeights.data() );
1067 assert (h1.GetSumOfWeights() > 0) ;
1068 fMean = h1.GetMean();
1069 fSigma = h1.GetRMS();
1070 // compute robust sigma using midspread
1071 Double_t quantiles[2] = {0.0, 0.0};
1072 Double_t prob[2] = {0.25, 0.75};
1073 h1.GetQuantiles(2, quantiles, prob);
1074 Double_t midspread = quantiles[1] - quantiles[0];
1075 fSigmaRob = std::min(fSigma, midspread / 1.349); // Sigma's robust estimator
1076 //printf("weight case - stat: m = %f, s= %f, sr = %f \n",fMean, fSigma, midspread);
1077 return;
1078 }
1079 else {
1080 // compute statistics using the data
1081 SetMean();
1082 Double_t midspread = ComputeMidspread();
1083 SetSigma(midspread);
1084 //printf("un-weight case - stat: m = %f, s= %f, sr = %f \n",fMean, fSigma, fSigmaRob);
1085 }
1086}
1087
1089 // Computes the inter-quartile range from the data
1090 std::sort(fEvents.begin(), fEvents.end());
1091 Double_t quantiles[2] = {0.0, 0.0};
1092 Double_t prob[2] = {0.25, 0.75};
1093 TMath::Quantiles(fEvents.size(), 2, &fEvents[0], quantiles, prob);
1094 Double_t lowerquartile = quantiles[0];
1095 Double_t upperquartile = quantiles[1];
1096 return upperquartile - lowerquartile;
1097}
1098
1100 // Computes the user's input kernel function canonical bandwidth
1102}
1103
1105 // Computes the user's input kernel function sigma2
1107}
1108
1109TKDE::KernelIntegrand::KernelIntegrand(const TKDE* kde, EIntegralResult intRes) : fKDE(kde), fIntegralResult(intRes) {
1110 // Internal class constructor
1111}
1112
1114 // Internal class unary function
1115 if (fIntegralResult == kNorm) {
1116 return std::pow((*fKDE->fKernelFunction)(x), 2);
1117 } else if (fIntegralResult == kMu) {
1118 return x * (*fKDE->fKernelFunction)(x);
1119 } else if (fIntegralResult == kSigma2) {
1120 return std::pow(x, 2) * (*fKDE->fKernelFunction)(x);
1121 } else if (fIntegralResult == kUnitIntegration) {
1122 return (*fKDE->fKernelFunction)(x);
1123 } else {
1124 return -1;
1125 }
1126}
1127
1129 //Returns the estimated density
1130 TString name = "KDEFunc_"; name+= GetName();
1131 TString title = "KDE "; title+= GetTitle();
1132 if (xMin >= xMax) { xMin = fXMin; xMax = fXMax; }
1133 //Do not register the TF1 to global list
1134 bool previous = TF1::DefaultAddToGlobalList(kFALSE);
1135 TF1 * pdf = new TF1(name.Data(), this, xMin, xMax, 0);
1137 if (npx > 0) pdf->SetNpx(npx);
1138 pdf->SetTitle(title);
1139 return pdf;
1140}
1141
1143 // Returns the upper estimated density
1144 TString name;
1145 name.Form("KDE_UpperCL%f5.3_%s",confidenceLevel,GetName());
1146 if (xMin >= xMax) { xMin = fXMin; xMax = fXMax; }
1147 TF1 * upperPDF = new TF1(name, this, &TKDE::UpperConfidenceInterval, xMin, xMax, 1);
1148 upperPDF->SetParameter(0, confidenceLevel);
1149 if (npx > 0) upperPDF->SetNpx(npx);
1150 TF1 * f = (TF1*)upperPDF->Clone();
1151 delete upperPDF;
1152 return f;
1153}
1154
1156 // Returns the upper estimated density
1157 TString name;
1158 name.Form("KDE_LowerCL%f5.3_%s",confidenceLevel,GetName());
1159 if (xMin >= xMax) { xMin = fXMin; xMax = fXMax; }
1160 TF1 * lowerPDF = new TF1(name, this, &TKDE::LowerConfidenceInterval, xMin, xMax, 1);
1161 lowerPDF->SetParameter(0, confidenceLevel);
1162 if (npx > 0) lowerPDF->SetNpx(npx);
1163 TF1 * f = (TF1*)lowerPDF->Clone();
1164 delete lowerPDF;
1165 return f;
1166}
1167
1169 // Returns the approximate bias
1170 TString name = "KDE_Bias_"; name += GetName();
1171 if (xMin >= xMax) { xMin = fXMin; xMax = fXMax; }
1172 TF1 * approximateBias = new TF1(name, this, &TKDE::ApproximateBias, xMin, xMax, 0);
1173 if (npx > 0) approximateBias->SetNpx(npx);
1174 TF1 * f = (TF1*)approximateBias->Clone();
1175 delete approximateBias;
1176 return f;
1177}
#define f(i)
Definition: RSha256.hxx:104
#define R(a, b, c, d, e, f, g, h, i)
Definition: RSha256.hxx:110
static const double x2[5]
static const double x1[5]
#define M_PI
Definition: Rotated.cxx:105
int Int_t
Definition: RtypesCore.h:41
unsigned int UInt_t
Definition: RtypesCore.h:42
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
const Bool_t kTRUE
Definition: RtypesCore.h:87
const char Option_t
Definition: RtypesCore.h:62
#define ClassImp(name)
Definition: Rtypes.h:363
@ kRed
Definition: Rtypes.h:63
@ kBlue
Definition: Rtypes.h:63
float xmin
Definition: THbookFile.cxx:93
float xmax
Definition: THbookFile.cxx:93
double pow(double, double)
double sqrt(double)
double exp(double)
double log(double)
TRObject operator()(const T1 &t1) const
#define gPad
Definition: TVirtualPad.h:286
Interface (abstract class) for generic functions objects of one-dimension Provides a method to evalua...
Definition: IFunction.h:135
User Class for performing numerical integration of a function in one dimension.
Definition: Integrator.h:94
void SetFunction(Function &f)
method to set the a generic integration function
Definition: Integrator.h:489
double Integral(Function &f, double a, double b)
evaluate the Integral of a function f over the defined interval (a,b)
Definition: Integrator.h:496
User class for calculating the derivatives of a function.
double Derivative2(double x)
Returns the second derivative of the function at point x, computed by Richardson's extrapolation meth...
void SetFunction(const IGenFunction &f)
Set function for derivative calculation (copy the function if option has been enabled in the construc...
Template class to wrap any C++ callable object which takes one argument i.e.
Template class to wrap any member function of a class taking a double and returning a double in a 1D ...
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
1-Dim function class
Definition: TF1.h:211
virtual void SetTitle(const char *title="")
Set function title if title has the form "fffffff;xxxx;yyyy", it is assumed that the function title i...
Definition: TF1.cxx:3542
virtual void SetNpx(Int_t npx=100)
Set the number of points used to draw the function.
Definition: TF1.cxx:3426
virtual void Draw(Option_t *option="")
Draw this function with its current attributes.
Definition: TF1.cxx:1312
static Bool_t DefaultAddToGlobalList(Bool_t on=kTRUE)
Static method to add/avoid to add automatically functions to the global list (gROOT->GetListOfFunctio...
Definition: TF1.cxx:823
virtual void SetParameter(Int_t param, Double_t value)
Definition: TF1.h:618
A TGraphErrors is a TGraph with error bars.
Definition: TGraphErrors.h:26
virtual void SetName(const char *name="")
Set graph name.
Definition: TGraph.cxx:2223
virtual void SetTitle(const char *title="")
Set graph title.
Definition: TGraph.cxx:2232
virtual void Draw(Option_t *chopt="")
Draw this graph with its current attributes.
Definition: TGraph.cxx:747
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:614
virtual Int_t GetQuantiles(Int_t nprobSum, Double_t *q, const Double_t *probSum=0)
Compute Quantiles for this histogram Quantile x_q of a probability distribution Function F is defined...
Definition: TH1.cxx:4342
virtual Double_t GetMean(Int_t axis=1) const
For axis = 1,2 or 3 returns the mean value of the histogram along X,Y or Z axis.
Definition: TH1.cxx:6958
Double_t GetRMS(Int_t axis=1) const
Definition: TH1.h:311
virtual void FillN(Int_t ntimes, const Double_t *x, const Double_t *w, Int_t stride=1)
Fill this histogram with an array x and weights w.
Definition: TH1.cxx:3354
virtual Double_t GetSumOfWeights() const
Return the sum of weights excluding under/overflows.
Definition: TH1.cxx:7297
Kernel Density Estimation class.
Definition: TKDE.h:31
TF1 * GetPDFUpperConfidenceInterval(Double_t confidenceLevel=0.95, UInt_t npx=100, Double_t xMin=1.0, Double_t xMax=0.0)
Definition: TKDE.cxx:1142
TF1 * GetKDEApproximateBias(UInt_t npx=100, Double_t xMin=1.0, Double_t xMax=0.0)
Definition: TKDE.cxx:1168
void SetData(const Double_t *data, const Double_t *weights)
Definition: TKDE.cxx:446
TF1 * fLowerPDF
Definition: TKDE.h:148
virtual ~TKDE()
Definition: TKDE.cxx:105
std::vector< Double_t > fKernelSigmas2
Definition: TKDE.h:178
Double_t ComputeKernelL2Norm() const
Definition: TKDE.cxx:1023
TF1 * fPDF
Definition: TKDE.h:146
TF1 * GetPDFLowerConfidenceInterval(Double_t confidenceLevel=0.95, UInt_t npx=100, Double_t xMin=1.0, Double_t xMax=0.0)
Definition: TKDE.cxx:1155
void SetKernelType(EKernelType kern)
Definition: TKDE.cxx:329
std::vector< Double_t > fCanonicalBandwidths
Definition: TKDE.h:177
void SetKernelFunction(KernelFunction_Ptr kernfunc=0)
Definition: TKDE.cxx:558
UInt_t fNEvents
Definition: TKDE.h:163
void ComputeDataStats()
Definition: TKDE.cxx:1059
Double_t fXMax
Definition: TKDE.h:171
Double_t UpperConfidenceInterval(const Double_t *x, const Double_t *p) const
Definition: TKDE.cxx:962
friend class TKernel
Definition: TKDE.h:137
Double_t ApproximateBias(const Double_t *x, const Double_t *) const
Definition: TKDE.h:208
Double_t ComputeMidspread()
Definition: TKDE.cxx:1088
Bool_t fNewData
Definition: TKDE.h:159
void SetMirror()
Definition: TKDE.cxx:437
Bool_t fUseMirroring
Definition: TKDE.h:157
void DrawConfidenceInterval(TString &drawOpt, double cl=0.95)
Definition: TKDE.cxx:862
void SetMirroredEvents()
Definition: TKDE.cxx:500
EMirror fMirror
Definition: TKDE.h:154
void SetUserCanonicalBandwidth()
Definition: TKDE.cxx:1099
Bool_t fMirrorLeft
Definition: TKDE.h:157
EIteration fIteration
Definition: TKDE.h:153
void CheckKernelValidity()
Definition: TKDE.cxx:999
const Double_t * GetAdaptiveWeights() const
Definition: TKDE.cxx:888
Double_t fAdaptiveBandwidthFactor
Definition: TKDE.h:173
Double_t LowerConfidenceInterval(const Double_t *x, const Double_t *p) const
Definition: TKDE.cxx:971
Double_t fSigmaRob
Definition: TKDE.h:169
std::vector< Double_t > fBinCount
Definition: TKDE.h:180
EBinning fBinning
Definition: TKDE.h:155
EIteration
Definition: TKDE.h:43
@ kAdaptive
Definition: TKDE.h:44
@ kFixed
Definition: TKDE.h:45
Double_t GetRAMISE() const
Definition: TKDE.cxx:689
void SetIteration(EIteration iter)
Definition: TKDE.cxx:340
Double_t ComputeKernelIntegral() const
Definition: TKDE.cxx:1050
Double_t CosineArchKernel(Double_t x) const
Definition: TKDE.h:202
Double_t fXMin
Definition: TKDE.h:170
Double_t operator()(Double_t x) const
Definition: TKDE.cxx:671
void SetUserKernelSigma2()
Definition: TKDE.cxx:1104
Double_t GetBias(Double_t x) const
Definition: TKDE.cxx:981
std::vector< Double_t > fData
Definition: TKDE.h:142
Double_t fSumOfCounts
Definition: TKDE.h:164
UInt_t fUseBinsNEvents
Definition: TKDE.h:165
void SetKernel()
Definition: TKDE.cxx:544
Double_t fSigma
Definition: TKDE.h:168
TGraphErrors * GetGraphWithErrors(UInt_t npx=100, Double_t xMin=1.0, Double_t xMax=0.0)
Definition: TKDE.cxx:838
Double_t GetMean() const
Definition: TKDE.cxx:677
Double_t fRho
Definition: TKDE.h:172
Bool_t fAsymRight
Definition: TKDE.h:157
Double_t fWeightSize
Definition: TKDE.h:175
TKernel * fKernel
Definition: TKDE.h:140
void SetUseBinsNEvents(UInt_t nEvents)
Definition: TKDE.cxx:387
std::vector< Double_t > fEvents
Definition: TKDE.h:143
Double_t GetError(Double_t x) const
Definition: TKDE.cxx:990
TF1 * GetKDEFunction(UInt_t npx=100, Double_t xMin=1.0, Double_t xMax=0.0)
Definition: TKDE.cxx:1128
void SetBinning(EBinning)
Definition: TKDE.cxx:360
void GetOptions(std::string optionType, std::string option)
Definition: TKDE.cxx:223
Bool_t fAsymLeft
Definition: TKDE.h:157
std::vector< Bool_t > fSettedOptions
Definition: TKDE.h:182
Double_t GaussianKernel(Double_t x) const
Definition: TKDE.h:190
void SetRange(Double_t xMin, Double_t xMax)
Definition: TKDE.cxx:404
void SetMean()
Definition: TKDE.cxx:533
virtual void Draw(const Option_t *option="")
Definition: TKDE.cxx:780
Double_t ComputeKernelSigma2() const
Definition: TKDE.cxx:1032
void SetOptions(const Option_t *option, Double_t rho)
Definition: TKDE.cxx:153
void SetUseBins()
Definition: TKDE.cxx:418
Double_t GetFixedWeight() const
Definition: TKDE.cxx:877
void AssureOptions()
Definition: TKDE.cxx:290
TF1 * GetFunction(UInt_t npx=100, Double_t xMin=1.0, Double_t xMax=0.0)
Definition: TKDE.cxx:619
void SetBinCountData()
Definition: TKDE.cxx:742
TF1 * GetUpperFunction(Double_t confidenceLevel=0.95, UInt_t npx=100, Double_t xMin=1.0, Double_t xMax=0.0)
Definition: TKDE.cxx:628
void InitFromNewData()
Definition: TKDE.cxx:480
void Instantiate(KernelFunction_Ptr kernfunc, UInt_t events, const Double_t *data, const Double_t *weight, Double_t xMin, Double_t xMax, const Option_t *option, Double_t rho)
Definition: TKDE.cxx:116
void SetDrawOptions(const Option_t *option, TString &plotOpt, TString &drawOpt)
Definition: TKDE.cxx:179
EMirror
Definition: TKDE.h:48
@ kMirrorLeft
Definition: TKDE.h:50
@ kMirrorAsymLeftRight
Definition: TKDE.h:54
@ kMirrorAsymRight
Definition: TKDE.h:55
@ kMirrorAsymBoth
Definition: TKDE.h:57
@ kNoMirror
Definition: TKDE.h:49
@ kMirrorRight
Definition: TKDE.h:51
@ kMirrorLeftAsymRight
Definition: TKDE.h:56
@ kMirrorBoth
Definition: TKDE.h:52
@ kMirrorAsymLeft
Definition: TKDE.h:53
TGraphErrors * fGraph
Definition: TKDE.h:150
void SetCanonicalBandwidths()
Definition: TKDE.cxx:602
void SetBinCentreData(Double_t xmin, Double_t xmax)
Definition: TKDE.cxx:733
void SetTuneFactor(Double_t rho)
Definition: TKDE.cxx:394
UInt_t Index(Double_t x) const
Definition: TKDE.cxx:947
void Fill(Double_t data)
Definition: TKDE.cxx:643
TF1 * fUpperPDF
Definition: TKDE.h:147
Bool_t fMirrorRight
Definition: TKDE.h:157
UInt_t fNBins
Definition: TKDE.h:162
Double_t ComputeKernelMu() const
Definition: TKDE.cxx:1041
EKernelType
Definition: TKDE.h:34
@ kGaussian
Definition: TKDE.h:35
@ kEpanechnikov
Definition: TKDE.h:36
@ kCosineArch
Definition: TKDE.h:38
@ kBiweight
Definition: TKDE.h:37
@ kTotalKernels
Definition: TKDE.h:40
@ kUserDefined
Definition: TKDE.h:39
Double_t fMean
Definition: TKDE.h:167
void DrawErrors(TString &drawOpt)
Definition: TKDE.cxx:831
void SetNBins(UInt_t nbins)
Definition: TKDE.cxx:368
void CheckOptions(Bool_t isUserDefinedKernel=kFALSE)
Definition: TKDE.cxx:306
Double_t EpanechnikovKernel(Double_t x) const
Definition: TKDE.h:195
Double_t BiweightKernel(Double_t x) const
Definition: TKDE.h:198
void SetKernelSigmas2()
Definition: TKDE.cxx:611
Bool_t fUseMinMaxFromData
Definition: TKDE.h:160
Double_t GetSigma() const
Definition: TKDE.cxx:683
EKernelType fKernelType
Definition: TKDE.h:152
TF1 * fApproximateBias
Definition: TKDE.h:149
void SetSigma(Double_t R)
Definition: TKDE.cxx:538
std::vector< Double_t > fEventWeights
Definition: TKDE.h:144
TF1 * GetApproximateBias(UInt_t npx=100, Double_t xMin=1.0, Double_t xMax=0.0)
Definition: TKDE.cxx:638
TF1 * GetLowerFunction(Double_t confidenceLevel=0.95, UInt_t npx=100, Double_t xMin=1.0, Double_t xMax=0.0)
Definition: TKDE.cxx:633
Bool_t fUseBins
Definition: TKDE.h:158
KernelFunction_Ptr fKernelFunction
Definition: TKDE.h:135
friend struct KernelIntegrand
Definition: TKDE.h:184
EBinning
Definition: TKDE.h:60
@ kUnbinned
Definition: TKDE.h:61
@ kRelaxedBinning
Definition: TKDE.h:62
@ kForcedBinning
Definition: TKDE.h:63
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
virtual TObject * Clone(const char *newname="") const
Make a clone of an object using the Streamer facility.
Definition: TNamed.cxx:74
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:866
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:880
Basic string class.
Definition: TString.h:131
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1100
const char * Data() const
Definition: TString.h:364
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:687
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:619
double normal_quantile(double z, double sigma)
Inverse ( ) of the cumulative distribution function of the lower tail of the normal (Gaussian) distri...
const Double_t sigma
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
const Int_t n
Definition: legend1.C:16
Double_t ey[n]
Definition: legend1.C:17
for(Int_t i=0;i< n;i++)
Definition: legend1.C:18
Double_t ex[n]
Definition: legend1.C:17
TH1F * h1
Definition: legend1.C:5
double inner_product(const LAVector &, const LAVector &)
static constexpr double s
Bool_t IsNaN(Double_t x)
Definition: TMath.h:880
void Quantiles(Int_t n, Int_t nprob, Double_t *x, Double_t *quantiles, Double_t *prob, Bool_t isSorted=kTRUE, Int_t *index=0, Int_t type=7)
Computes sample quantiles, corresponding to the given probabilities.
Definition: TMath.cxx:1190