Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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 <memory>
31#include <numeric>
32#include <limits>
33#include <cassert>
34
35#include "Math/Error.h"
36#include "TMath.h"
37#include "Math/Functor.h"
38#include "Math/Integrator.h"
41#include "TGraphErrors.h"
42#include "TF1.h"
43#include "TH1.h"
44#include "TVirtualPad.h"
45#include "TKDE.h"
46
48
49
58
59///////////////////////////////////////////////////////////////////////
60/// default constructor used by I/O
62 fKernelFunction(nullptr),
63 fKernel(nullptr),
64 fPDF(nullptr),
65 fUpperPDF(nullptr),
66 fLowerPDF(nullptr),
67 fApproximateBias(nullptr),
68 fGraph(nullptr),
69 fUseMirroring(false), fMirrorLeft(false), fMirrorRight(false), fAsymLeft(false), fAsymRight(false),
70 fUseBins(false), fNewData(false), fUseMinMaxFromData(false),
71 fNBins(0), fNEvents(0), fSumOfCounts(0), fUseBinsNEvents(0),
72 fMean(0.),fSigma(0.), fSigmaRob(0.), fXMin(0.), fXMax(0.),
73 fRho(0.), fAdaptiveBandwidthFactor(0.), fWeightSize(0)
74{
75}
76
78 //Class destructor
79 if (fPDF) delete fPDF;
80 if (fUpperPDF) delete fUpperPDF;
81 if (fLowerPDF) delete fLowerPDF;
82 if (fGraph) delete fGraph;
84 // delete fKernelFunction if it is not user defined
86 delete fKernelFunction;
87}
88
89////////////////////////////////////////////////////////////////////
90//// Internal function to instantiate a TKDE object
92
93 fPDF = nullptr;
94 fKernelFunction = nullptr;
95 fUpperPDF = nullptr;
96 fLowerPDF = nullptr;
97 fApproximateBias = nullptr;
98 fGraph = nullptr;
99 fNewData = false;
100 fUseMirroring = false; fMirrorLeft = false; fMirrorRight = false;
101 fAsymLeft = false; fAsymRight = false;
102 fNBins = events < 10000 ? 1000 : std::min(10000, int(events / 100)*10);
103 fNEvents = events;
104 fUseBinsNEvents = 10000;
105 fMean = 0.0;
106 fSigma = 0.0;
107 fXMin = xMin;
108 fXMax = xMax;
110 fSumOfCounts = 0;
112 fRho = rho;
113 fWeightSize = 0;
114 fCanonicalBandwidths = std::vector<Double_t>(kTotalKernels, 0.0);
115 fKernelSigmas2 = std::vector<Double_t>(kTotalKernels, -1.0);
116 fSettedOptions = std::vector<Bool_t>(4, kFALSE);
117 SetOptions(option, rho);
119 SetMirror();
120 SetUseBins();
123
124}
125
127 //Sets User defined construction options
128 TString opt = option;
129 opt.ToLower();
130 std::string options = opt.Data();
131 size_t numOpt = 4;
132 std::vector<std::string> voption(numOpt, "");
133 for (std::vector<std::string>::iterator it = voption.begin(); it != voption.end() && !options.empty(); ++it) {
134 size_t pos = options.find_last_of(';');
135 if (pos == std::string::npos) {
136 *it = options;
137 break;
138 }
139 *it = options.substr(pos + 1);
140 options = options.substr(0, pos);
141 }
142 for (std::vector<std::string>::iterator it = voption.begin(); it != voption.end(); ++it) {
143 size_t pos = (*it).find(':');
144 if (pos != std::string::npos) {
145 GetOptions((*it).substr(0, pos) , (*it).substr(pos + 1));
146 }
147 }
149 fRho = rho;
150}
151
153 // Sets User defined drawing options
154 size_t numOpt = 2;
155 std::string options = TString(option).Data();
156 std::vector<std::string> voption(numOpt, "");
157 for (std::vector<std::string>::iterator it = voption.begin(); it != voption.end() && !options.empty(); ++it) {
158 size_t pos = options.find_last_of(';');
159 if (pos == std::string::npos) {
160 *it = options;
161 break;
162 }
163 *it = options.substr(pos + 1);
164 options = options.substr(0, pos);
165 }
168 for (std::vector<std::string>::iterator it = voption.begin(); it != voption.end() && !options.empty(); ++it) {
169 size_t pos = (*it).find(':');
170 if (pos == std::string::npos) break;
171 TString optionType = (*it).substr(0, pos);
172 TString optionInstance = (*it).substr(pos + 1);
173 optionType.ToLower();
174 optionInstance.ToLower();
175 if (optionType.Contains("plot")) {
177 if (optionInstance.Contains("estimate") || optionInstance.Contains("errors") || optionInstance.Contains("confidenceinterval"))
179 else {
180 this->Warning("SetDrawOptions", "Unknown plotting option %s: setting to KDE estimate plot.",optionInstance.Data());
181 this->Info("SetDrawOptions", "Possible plotting options are: Estimate,Errors,ConfidenceInterval");
182 }
183 } else if (optionType.Contains("drawoptions")) {
186 }
187 }
188 if (!foundPlotOPt) {
189 this->Warning("SetDrawOptions", "No plotting option: setting to KDE estimate plot.");
190 plotOpt = "estimate";
191 }
192 if (!foundDrawOPt) {
193 this->Warning("SetDrawOptions", "No drawing options: setting to default ones.");
194 drawOpt = "apl4";
195 }
196}
197
198void TKDE::GetOptions(std::string optionType, std::string option) {
199 // Gets User defined KDE construction options
200 if (optionType == "kerneltype") {
202 if (option == "gaussian") {
204 } else if (option == "epanechnikov") {
206 } else if (option == "biweight") {
208 } else if (option == "cosinearch") {
210 } else if (option == "userdefined") {
212 } else {
213 this->Warning("GetOptions", "Unknown kernel type option %s: setting to Gaussian",option.c_str());
214 this->Info("GetOptions", "Possible kernel type options are: Gaussian, Epanechnikov, Biweight, Cosinearch, Userdefined");
216 }
217 } else if (optionType == "iteration") {
219 if (option == "adaptive") {
221 } else if (option == "fixed") {
223 } else {
224 this->Warning("GetOptions", "Unknown iteration option %s: setting to Adaptive",option.c_str());
225 this->Info("GetOptions","Possible iteration type options are: Adaptive, Fixed");
227 }
228 } else if (optionType == "mirror") {
230 if (option == "nomirror") {
232 } else if (option == "mirrorleft") {
234 } else if (option == "mirrorright") {
236 } else if (option == "mirrorboth") {
238 } else if (option == "mirrorasymleft") {
240 } else if (option == "mirrorrightasymleft") {
242 } else if (option == "mirrorasymright") {
244 } else if (option == "mirrorleftasymright") {
246 } else if (option == "mirrorasymboth") {
248 } else {
249 this->Warning("GetOptions", "Unknown mirror option %s: setting to NoMirror",option.c_str());
250 this->Info("GetOptions", "Possible mirror type options are: NoMirror, MirrorLeft, MirrorRight, MirrorAsymLeft,"
251 "MirrorAsymRight, MirrorRightAsymLeft, MirrorLeftAsymRight, MirrorAsymBoth");
253 }
254 } else if (optionType == "binning") {
256 if (option == "unbinned") {
258 } else if (option == "relaxedbinning") {
260 } else if (option == "forcedbinning") {
262 } else {
263 this->Warning("GetOptions", "Unknown binning option %s: setting to RelaxedBinning", option.c_str());
264 this->Info("GetOptions", "Possible binning type options are: Unbinned, ForcedBinning, RelaxedBinning");
266 }
267 }
268}
269
271 // Sets missing construction options to default ones
272 if (!fSettedOptions[0]) {
274 }
275 if (!fSettedOptions[1]) {
277 }
278 if (!fSettedOptions[2]) {
280 }
281 if (!fSettedOptions[3]) {
283 }
284}
285
287 // Sets User global options
289 Error("CheckOptions", "Illegal user kernel type input! Use template constructor for user defined kernel.");
290 }
291 if (fIteration != kAdaptive && fIteration != kFixed) {
292 Warning("CheckOptions", "Illegal user iteration type input - use default value !");
294 }
295 if (!(fMirror >= kNoMirror && fMirror <= kMirrorAsymBoth)) {
296 Warning("CheckOptions", "Illegal user mirroring type input - use default value !");
298 }
299 if (!(fBinning >= kUnbinned && fBinning <= kForcedBinning)) {
300 Warning("CheckOptions", "Illegal user binning type input - use default value !");
302 }
303 if (fRho <= 0.0) {
304 Warning("CheckOptions", "Tuning factor rho cannot be non-positive - use default value !");
305 fRho = 1.0;
306 }
307}
308
310 // Sets User option for the choice of kernel estimator
312 delete fKernelFunction;
313 fKernelFunction = nullptr;
314 }
316 CheckOptions();
317 // resetting fKernel automatically calls ReInit when needed
318 fKernel.reset();
319}
320
322 // Sets User option for fixed or adaptive iteration
323 fIteration = iter;
324 CheckOptions();
325 fKernel.reset();
326}
327
329 // Sets User option for mirroring the data
330 fMirror = mir;
331 CheckOptions();
332 SetMirror();
333 if (fUseMirroring) {
335 }
336 fKernel.reset();
337}
338
340 // Sets User option for binning the weights
341 fBinning = bin;
342 CheckOptions();
343 SetUseBins();
344}
345
347 // Sets User option for number of bins
348 if (!nbins) {
349 Error("SetNBins", "Number of bins must be greater than zero.");
350 return;
351 }
352
353 fNBins = nbins;
354
355 SetUseBins();
356 if (!fUseBins) {
357 if (fBinning == kUnbinned)
358 Warning("SetNBins", "Bin type using SetBinning must be set for using a binned evaluation");
359 else
360 Warning("SetNBins", "Bin type using SetBinning or with SetUseBinsNEvents must be set for using a binned evaluation");
361 }
362}
363
365 // Sets User option for the minimum number of events for allowing automatic binning
366 fUseBinsNEvents = nEvents;
367 SetUseBins();
368}
369
371 // Factor which can be used to tune the smoothing.
372 // It is used as multiplicative factor for the fixed and adaptive bandwidth.
373 // A value < 1 will reproduce better the tails but oversmooth the peak
374 // while a factor > 1 will overestimate the tail
375 fRho = rho;
376 CheckOptions();
377 fKernel.reset();
378}
379
381 // Sets minimum range value and maximum range value
382 if (xMin >= xMax) {
383 Error("SetRange", "Minimum range cannot be bigger or equal than the maximum range! Present range values remain the same.");
384 return;
385 }
386 fXMin = xMin;
387 fXMax = xMax;
388 fUseMinMaxFromData = false;
389 fKernel.reset();
390}
391
392// private methods
393
395 // Sets User option for using binned weights
396 switch (fBinning) {
397 default:
398 case kRelaxedBinning:
399 if (fNEvents >= fUseBinsNEvents) {
400 fUseBins = kTRUE;
401 } else {
403 }
404 break;
405 case kForcedBinning:
406 fUseBins = kTRUE;
407 break;
408 case kUnbinned:
410 }
411
412 // during initialization we don't need to recompute the bins
413 // it is done within TKDE::SetData
414 // in this case fEvents is empty
415 if (fEvents.empty()) return;
416
417 // in case we need to recompute bins structure
418 // needed when called from the TKDE setters method
419 if (fUseBins) {
423 }
424 else {
425 // unbinned case
426 fWeightSize = 0.;
427 fBinCount.clear();
428 fData = fEvents;
429 }
430 fKernel.reset();
431}
432
441
443 // Sets the data events input sample or bin centres for binned option and computes basic estimators
444 if (!data) {
445 if (fNEvents) fData.reserve(fNEvents);
446 return;
447 }
448 fEvents.assign(data, data + fNEvents);
449 if (wgts) fEventWeights.assign(wgts, wgts + fNEvents);
450
451 if (fUseMinMaxFromData) {
452 fXMin = *std::min_element(fEvents.begin(), fEvents.end());
453 fXMax = *std::max_element(fEvents.begin(), fEvents.end());
454 }
455
456 if (fUseBins) {
457 if (fNBins >= fNEvents) {
458 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");
459 }
462 } else {
463 fWeightSize = 0; //fNEvents / (fXMax - fXMin);
464 fData = fEvents;
465 }
466
468 if (fUseMirroring) {
469 // set mirror events called automatically SetBinCount
471 }
472 else {
473 // to set fBinCount and fSumOfCounts
475 }
476}
477
479 // re-initialize the KDE in case of new data or in case of reading from a file
480
481 // in case of new data re-compute the statistics (works only for unbinned case)
482 if (fNewData) {
484 return;
485 }
486
487 if (fEvents.empty()) {
488 Error("ReInit","TKDE does not contain any data !");
489 return;
490 }
491
492 // when of reading from a file fKernelFunction is a nullptr
493 // we need to recreate Kernel class (with adaptive weights if needed) and
494 // recreate kernel function pointer
495 if (!fKernelFunction) SetKernelFunction(nullptr);
496
497 SetKernel();
498}
499
501 // re-initialize when new data have been filled in TKDE
502 // re-compute kernel quantities and statistics
503 if (fUseBins) {
504 Error("InitFromNewData","Re-felling is not supported with binning");
505 return;
506 }
507
508 fNewData = false;
509 // in case of unbinned data
510 if (!fUseBins) fEvents = fData;
511 if (fUseMinMaxFromData) {
512 fXMin = *std::min_element(fEvents.begin(), fEvents.end());
513 fXMax = *std::max_element(fEvents.begin(), fEvents.end());
514 }
516 //
518 if (fUseMirroring) {
520 }
521 else {
522 // to set fBinCount and fSumOfCounts
524 }
525 // in case of I/O reset the kernel
526 if (!fKernelFunction) SetKernelFunction(nullptr);
527
528 SetKernel();
529
530}
531
532///////////////////////////////////////////////////////////////////////////////////////////////////////
533/// Intgernal function to mirror the data
536 if (fUseBins) {
537 // binned case
538 // first fill the bin data in [fXmin,fXmax]
539 // bin center should have been already filled before
542
543 // now mirror the bins
544 assert(fNBins == fData.size());
545 if (fMirrorLeft) {
546 fData.insert(fData.begin(), fNBins, 0.);
547 fBinCount.insert(fBinCount.begin(), fNBins, 0.);
548 for (UInt_t i = 0; i < fNBins; ++i) {
549 fData[i] = fData[i + fNBins] - (fXMax - fXMin);
550 fBinCount[fNBins - i - 1] = fBinCount[i + fNBins];
551 }
552 }
553 if (fMirrorRight) {
554 fData.insert(fData.end(), fNBins, 0.);
555 fBinCount.insert(fBinCount.end(), fNBins, 0.);
556 int i0 = (fMirrorLeft) ? fNBins : 0;
557 int iLast = i0 + 2 * fNBins - 1;
558 for (UInt_t i = 0; i < fNBins; ++i) {
559 fData[i0 + fNBins + i] = fData[i0 + i] + (fXMax - fXMin);
560 fBinCount[iLast - i] = fBinCount[i0 + i];
561 }
562 }
563 }
564 else {
565 // unbinned case (transform events)
566 std::vector<Double_t> originalEvents = fEvents;
567 std::vector<Double_t> originalWeights = fEventWeights;
568 if (fMirrorLeft) {
569 fEvents.resize(2 * fNEvents, 0.0);
570 transform(fEvents.begin(), fEvents.begin() + fNEvents, fEvents.begin() + fNEvents,
571 std::bind(std::minus<Double_t>(), 2 * fXMin, std::placeholders::_1));
572 }
573 if (fMirrorRight) {
574 fEvents.resize((fMirrorLeft + 2) * fNEvents, 0.0);
575 transform(fEvents.begin(), fEvents.begin() + fNEvents, fEvents.begin() + (fMirrorLeft + 1) * fNEvents,
576 std::bind(std::minus<Double_t>(), 2 * fXMax, std::placeholders::_1));
577 }
578 if (!fEventWeights.empty() && (fMirrorLeft || fMirrorRight)) {
579 // copy weights too
580 fEventWeights.insert(fEventWeights.end(), fEventWeights.begin(), fEventWeights.begin() + fNEvents);
582 fEventWeights.insert(fEventWeights.end(), fEventWeights.begin(), fEventWeights.begin() + fNEvents);
583 }
584
585 fData = fEvents;
589 }
590
591}
592
594 // Computes input data's mean
595 fMean = std::accumulate(fEvents.begin(), fEvents.end(), 0.0) / fEvents.size();
596}
597
599 // Computes input data's sigma
600 fSigma = std::sqrt(1. / (fEvents.size() - 1.) * (std::inner_product(fEvents.begin(), fEvents.end(), fEvents.begin(), 0.0) - fEvents.size() * std::pow(fMean, 2.)));
601 fSigmaRob = std::min(fSigma, R / 1.349); // Sigma's robust estimator
602}
603
605 // Sets the kernel density estimator
606 // n here should be fNEvents or the total size including the mirrored events ???
607 // it should not be fData.size() that in binned case is number of bins
609 if (n == 0) return;
610 // Optimal bandwidth (Silverman's rule of thumb with assumed Gaussian density)
611 Double_t weight = fCanonicalBandwidths[kGaussian] * fSigmaRob * std::pow(3. / (8. * std::sqrt(M_PI)) * n, -0.2);
613
614 fKernel = std::make_unique<TKernel>(weight, this);
615
616 if (fIteration == kAdaptive) {
617 fKernel->ComputeAdaptiveWeights();
618 }
619 if (gDebug) {
620 if (fIteration != kAdaptive)
621 Info("SetKernel",
622 "Using a fix kernel - bandwidth = %f - using n = %d, rho = %f , sigmaRob = %f , mean = %f , sigma = %f "
623 ", canonicalBandwidth= %f",
625 else
626 Info("SetKernel",
627 "Using an adaptive kernel - weight = %f - using n = %d, rho = %f , sigmaRob = %f , mean = %f , sigma = %f "
628 ", canonicalBandwidth= %f",
630 }
631}
632
634 // sets the kernel function pointer. It creates and manages the pointer for internal classes
635 if (fKernelFunction) {
636 // kernel function pointer must be set to null before calling SetKernelFunction to avoid memory leaks
637 Fatal("SetKernelFunction", "Kernel function pointer is not null");
638 }
639 switch (fKernelType) {
640 case kGaussian :
642 break;
643 case kEpanechnikov :
645 break;
646 case kBiweight :
648 break;
649 case kCosineArch :
651 break;
652 case kUserDefined :
655 break;
656 case kTotalKernels :
657 default:
658 /// for user defined kernels
661 }
662
663 if (fKernelType == kUserDefined) {
664 if (fKernelFunction) {
668 }
669 else {
670 Error("SetKernelFunction", "User kernel function is not defined !");
671 return;
672 }
673 }
677 SetKernel();
678}
679
681 // Sets the canonical bandwidths according to the kernel type
682 fCanonicalBandwidths[kGaussian] = 0.7764; // Checked in Mathematica
683 fCanonicalBandwidths[kEpanechnikov] = 1.7188; // Checked in Mathematica
684 fCanonicalBandwidths[kBiweight] = 2.03617; // Checked in Mathematica
685 fCanonicalBandwidths[kCosineArch] = 1.7663; // Checked in Mathematica
686 fCanonicalBandwidths[kUserDefined] = 1.0; // To be Checked
687}
688
690 // Sets the kernel sigmas2 according to the kernel type
692 fKernelSigmas2[kEpanechnikov] = 1.0 / 5.0;
693 fKernelSigmas2[kBiweight] = 1.0 / 7.0;
694 fKernelSigmas2[kCosineArch] = 1.0 - 8.0 / std::pow(M_PI, 2);
695}
696
698 // Returns the PDF estimate as a function sampled in npx between xMin and xMax
699 // the KDE is not re-normalized to the xMin/xMax range.
700 // The user manages the returned function
701 // For getting a non-sampled TF1, one can create a TF1 directly from the TKDE by doing
702 // TF1 * f1 = new TF1("f1",kde,xMin,xMax,0);
703 return GetKDEFunction(npx,xMin,xMax);
704}
705
707 // Returns the PDF upper estimate (upper confidence interval limit)
709}
710
712 // Returns the PDF lower estimate (lower confidence interval limit)
714}
715
717 // Returns the PDF estimate bias
719}
720
722 // Fills data member with User input data event for the unbinned option
723 if (fUseBins) {
724 this->Warning("Fill", "Cannot fill data with data binned option. Data input ignored.");
725 return;
726 }
727 fData.push_back(data);
728 fNEvents++;
729 fNewData = kTRUE;
730}
731
733 // Fills data member with User input data event for the unbinned option
734 if (fUseBins) {
735 this->Warning("Fill", "Cannot fill data with data binned option. Data input ignored.");
736 return;
737 }
738 fData.push_back(data); // should not be here fEvent ??
739 fEventWeights.push_back(weight);
740 fNEvents++;
741 fNewData = kTRUE;
742}
743
745 // The class's unary function: returns the kernel density estimate
746 return (*this)(*x);
747}
748
750 // The class's unary function: returns the kernel density estimate
751 if (!fKernel) {
752 (const_cast<TKDE*>(this))->ReInit();
753 // in case of failed re-initialization
754 if (!fKernel) return TMath::QuietNaN();
755 }
756 return (*fKernel)(x);
757}
758
760 // return the mean of the data
761 if (fNewData) (const_cast<TKDE*>(this))->InitFromNewData();
762 return fMean;
763}
764
766 // return the standard deviation of the data
767 if (fNewData) (const_cast<TKDE*>(this))->InitFromNewData();
768 return fSigma;
769}
770
772 // Returns the Root Asymptotic Mean Integrated Squared Error according to Silverman's rule of thumb with assumed Gaussian density
773 Double_t result = 5. / 4. * fKernelSigmas2[fKernelType] * std::pow(fCanonicalBandwidths[fKernelType], 4) * std::pow(3. / (8. * std::sqrt(M_PI)) , -0.2) * fSigmaRob * std::pow(fNEvents, -0.8);
774 return std::sqrt(result);
775}
776
778// Internal class constructor
779fKDE(kde),
780fNWeights(kde->fData.size()),
781fWeights(1, weight)
782{}
783
785 // Gets the adaptive weights (bandwidths) for TKernel internal computation
786 unsigned int n = fKDE->fData.size();
787 Double_t minWeight = fWeights[0] * 0.05;
788 // we will store computed adaptive weights in weights
789 std::vector<Double_t> weights(n, fWeights[0]);
790 bool useDataWeights = (fKDE->fBinCount.size() == n);
791 Double_t f = 0.0;
792 for (unsigned int i = 0; i < n; ++i) {
793 // for negative or null bin contents use the fixed weight value (fWeights[0])
794 if (useDataWeights && fKDE->fBinCount[i] <= 0) {
795 weights[i] = fWeights[0];
796 continue; // skip negative or null weights
797 }
798 f = (*fKDE->fKernel)(fKDE->fData[i]);
799 if (f <= 0) {
800 // this can happen when data are outside range and fAsymLeft or fAsymRight is on
801 fKDE->Warning("ComputeAdativeWeights","function value is zero or negative for x = %f w = %f - set their bandwidth to zero",
802 fKDE->fData[i],(useDataWeights) ? fKDE->fBinCount[i] : 1.);
803 // set bandwidth for these points to zero.
804 weights[i] = 0;
805 continue;
806 }
807 // use weight if it is larger
808 weights[i] = std::max(weights[i] /= std::sqrt(f), minWeight);
809 fKDE->fAdaptiveBandwidthFactor += std::log(f);
810 // printf("(f = %f w = %f af = %f ),",f,*weight,fKDE->fAdaptiveBandwidthFactor);
811 }
812 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
813 // not sure for this special case for mirror. This results in a much smaller bandwidth for mirror case
814 fKDE->fAdaptiveBandwidthFactor = fKDE->fUseMirroring ? kAPPROX_GEO_MEAN / fKDE->fSigmaRob : std::sqrt(std::exp(fKDE->fAdaptiveBandwidthFactor / fKDE->fData.size()));
815 // set adaptive weights in fWeights matrix
816 fWeights.resize(n);
817 transform(weights.begin(), weights.end(), fWeights.begin(),
818 std::bind(std::multiplies<Double_t>(), std::placeholders::_1, fKDE->fAdaptiveBandwidthFactor));
819 //printf("adaptive bandwidth factor % f weight 0 %f , %f \n",fKDE->fAdaptiveBandwidthFactor, weights[0],fWeights[0] );
820}
821
823 // Returns the bandwidth
824 return fWeights[fKDE->Index(x)];
825}
826
828 // Returns the bins' centres from the data for using with the binned option
829 fData.assign(fNBins, 0.0);
830 Double_t binWidth((xmax - xmin) / fNBins);
831 for (UInt_t i = 0; i < fNBins; ++i) {
832 fData[i] = xmin + (i + 0.5) * binWidth;
833 }
834}
835
837 // Returns the bins' count from the data for using with the binned option
838 // or set the bin count to the weights in case of weighted data
839 if (fUseBins) {
840 fBinCount.assign(fNBins, 0.0);
841 fSumOfCounts = 0;
842 // note this function should be called before the data have been mirrored
844 assert(fEvents.size() == nevents);
845 // case of weighted events
846 if (!fEventWeights.empty() ) {
847 assert(nevents == fEventWeights.size());
848 for (UInt_t i = 0; i < nevents; ++i) {
849 if (fEvents[i] >= fXMin && fEvents[i] < fXMax) {
852 //printf("sum of counts %f - bin count %d - %f \n",fSumOfCounts, Index(fEvents[i]), fBinCount[Index(fEvents[i])] );
853 }
854 }
855 }
856 // case of unweighted data
857 else {
858 for (UInt_t i = 0; i < nevents; ++i) {
859 if (fEvents[i] >= fXMin && fEvents[i] < fXMax) {
860 fBinCount[Index(fEvents[i])] += 1;
861 fSumOfCounts += 1;
862 }
863 }
864 }
865 }
866 else if (!fEventWeights.empty() ) {
867 //unbinned weighted data
869 fSumOfCounts = 0;
870 for (UInt_t i = 0; i < fNEvents; ++i) {
871 if (fEvents[i] >= fXMin && fEvents[i] < fXMax) {
873 }
874 }
875 }
876 else {
877 // unbinned unweighted data
878 //fSumOfCounts = fNEvents;
879 fSumOfCounts = 0;
880 for (UInt_t i = 0; i < fNEvents; ++i) {
881 if (fEvents[i] >= fXMin && fEvents[i] < fXMax) {
882 fSumOfCounts += 1;
883 }
884 }
885 fBinCount.clear();
886 }
887 if (fSumOfCounts == 0) {
888 Warning("SetBinCountData()", "Empty sum of counts, might lead to nan/inf errors when normalizing.");
889 }
890}
891
892////////////////////////////////////////////////////////////////////////////////
893/// Draws either the KDE functions or its errors
894// @param opt : Drawing options:
895// - "" (default) - draw just the kde
896// - "same" draw on top of existing pad
897// - "Errors" draw a TGraphErrors with the point and errors
898// -"confidenceinterval" draw KDE + conf interval functions (default is 95%)
899// -"confidenceinterval@0.90" draw KDE + conf interval functions at 90%
900// - Extra options can be passed in opt for the drawing of the corresponding TF1 or TGraph
901//
902// NOTE: The functions GetDrawnFunction(), GetDrawnUpperFunction(), GetDrawnLowerFunction()
903// and GetGraphWithErrors() return the corresponding drawn objects (which are managed by the TKDE)
904// They can be used to changes style, color, etc..
905////
906void TKDE::Draw(const Option_t* opt) {
907
908 TString plotOpt = opt;
909 plotOpt.ToLower();
911 if(gPad && !plotOpt.Contains("same")) {
912 gPad->Clear();
913 }
914 if (plotOpt.Contains("errors")) {
915 drawOpt.ReplaceAll("errors","");
917 }
918 else if (plotOpt.Contains("confidenceinterval") ||
919 plotOpt.Contains("confinterval")) {
920 // parse level option
921 drawOpt.ReplaceAll("confidenceinterval","");
922 drawOpt.ReplaceAll("confinterval","");
923 Double_t level = 0.95;
924 const char * s = strstr(plotOpt.Data(),"interval@");
925 // coverity [secure_coding : FALSE]
926 if (s != nullptr) sscanf(s,"interval@%lf",&level);
927 if((level <= 0) || (level >= 1)) {
928 Warning("Draw","given confidence level %.3lf is invalid - use default 0.95",level);
929 level = 0.95;
930 }
932 }
933 else {
934 if (fPDF) delete fPDF;
936 fPDF->Draw(drawOpt);
937 }
938}
939
940///////////////////////////////////////////////////////////////////////
941/// Draws a TGraphErrors with KDE values and errors
943{
944 if (fGraph) delete fGraph;
946 fGraph->Draw(drawOpt.Data());
947}
948///////////////////////////////////////////////////////////////////////
949/// return a TGraphErrors with the KDE values and errors
950/// The return object is managed by the user
952 if (xmin>= xmax) { xmin = fXMin; xmax = fXMax; }
953 UInt_t n = npx;
954 Double_t* x = new Double_t[n + 1];
955 Double_t* ex = new Double_t[n + 1];
956 Double_t* y = new Double_t[n + 1];
957 Double_t* ey = new Double_t[n + 1];
958 for (UInt_t i = 0; i <= n; ++i) {
959 x[i] = xmin + i * (xmax - xmin) / n;
960 y[i] = (*this)(x[i]);
961 ex[i] = 0;
962 ey[i] = this->GetError(x[i]);
963 }
964 TGraphErrors* ge = new TGraphErrors(n, &x[0], &y[0], &ex[0], &ey[0]);
965 ge->SetName("kde_graph_error");
966 ge->SetTitle("Errors");
967 delete [] x;
968 delete [] ex;
969 delete [] y;
970 delete [] ey;
971 return ge;
972}
973
974//////////////////////////////////////////////////////////////
975/// // Draws the KDE and its confidence interval
977 GetKDEFunction()->Draw(drawOpt.Data());
979 upper->SetLineColor(kBlue);
980 upper->Draw(("same" + drawOpt).Data());
982 lower->SetLineColor(kRed);
983 lower->Draw(("same" + drawOpt).Data());
984 if (fUpperPDF) delete fUpperPDF;
985 if (fLowerPDF) delete fLowerPDF;
988}
989
991 // Returns the bandwidth for the non adaptive KDE
992 Double_t result = -1.;
994 this->Warning("GetFixedWeight()", "Fixed iteration option not enabled. Returning %f.", result);
995 } else {
996 result = fKernel->GetFixedWeight();
997 }
998 return result;
999}
1000
1002 // Returns the bandwidths for the adaptive KDE
1003 if (fIteration != TKDE::kAdaptive) {
1004 this->Warning("GetFixedWeight()", "Adaptive iteration option not enabled. Returning a NULL pointer<");
1005 return nullptr;
1006 }
1007 if (fNewData) (const_cast<TKDE*>(this))->InitFromNewData();
1008 return &(fKernel->GetAdaptiveWeights()).front();
1009}
1010
1012 // Returns the bandwidth for the non adaptive KDE
1013 return fWeights[0];
1014}
1015
1016const std::vector<Double_t> & TKDE::TKernel::GetAdaptiveWeights() const {
1017 // Returns the bandwidth for the adaptive KDE
1018 return fWeights;
1019}
1020
1022 // The internal class's unary function: returns the kernel density estimate
1023 Double_t result(0.0);
1024 UInt_t n = fKDE->fData.size();
1025 // case of bins or weighted data
1026 Bool_t useCount = (fKDE->fBinCount.size() == n);
1027 // also in case of unbinned unweighted data fSumOfCounts is sum of events in range
1028 // events outside range should be used to normalize the TKDE ??
1029 Double_t nSum = fKDE->fSumOfCounts; //(useCount) ? fKDE->fSumOfCounts : fKDE->fNEvents;
1030 //if (!useCount) nSum = fKDE->fNEvents;
1031 // in case of non-adaptive fWeights is a vector of size 1
1032 Bool_t hasAdaptiveWeights = (fWeights.size() == n);
1033 Double_t invWeight = (!hasAdaptiveWeights) ? 1. / fWeights[0] : 0;
1034 for (UInt_t i = 0; i < n; ++i) {
1035 Double_t binCount = (useCount) ? fKDE->fBinCount[i] : 1.0;
1036 // uncommenting following line slows down so keep computation for
1037 // zero bincounts
1038 //if (binCount <= 0) continue;
1039 if (hasAdaptiveWeights) {
1040 // skip data points that have 0 bandwidth (this can happen, see TKernel::ComputeAdaptiveWeight)
1041 if (fWeights[i] == 0) continue;
1042 invWeight = 1. / fWeights[i];
1043 }
1044 result += binCount * invWeight * (*fKDE->fKernelFunction)((x - fKDE->fData[i]) * invWeight );
1045 if (fKDE->fAsymLeft) {
1046 result += binCount * invWeight * (*fKDE->fKernelFunction)((x - (2. * fKDE->fXMin - fKDE->fData[i])) * invWeight);
1047 }
1048 if (fKDE->fAsymRight) {
1049 result += binCount * invWeight * (*fKDE->fKernelFunction)((x - (2. * fKDE->fXMax - fKDE->fData[i])) * invWeight);
1050 }
1051 // printf("data point %i %f %f count %f weight % f result % f\n",i,fKDE->fData[i],fKDE->fEvents[i],binCount,fWeights[i], result);
1052 }
1053 if ( TMath::IsNaN(result) ) {
1054 fKDE->Warning("operator()","Result is NaN for x %f \n",x);
1055 }
1056 return result / nSum;
1057}
1058
1059////////////////////////////////////////////////////
1060/// compute the bin index given a data point x
1062 // Returns the indices (bins) for the data point
1063 // fWeightSize is the inverse of the bin width
1064 Int_t bin = Int_t((x - fXMin) * fWeightSize);
1065 if (bin == (Int_t)fData.size()) return --bin;
1066 if (bin > (Int_t)fData.size()) {
1067 return (Int_t)(fData.size()) - 1;
1068 } else if (bin <= 0) {
1069 return 0;
1070 }
1071 return bin;
1072}
1073
1075 // Returns the pointwise upper estimated density
1076 Double_t f = (*this)(x);
1078 Double_t prob = 1. - (1.-*p)/2; // this is 1.-alpha/2
1080 return f + z * sigma;
1081}
1082
1084 // Returns the pointwise lower estimated density
1085 Double_t f = (*this)(x);
1087 Double_t prob = (1.-*p)/2; // this is alpha/2
1089 return f + z * sigma;
1090}
1091
1092
1094 // Returns the pointwise approximate estimated density bias
1097 rd.SetFunction(kern);
1098 Double_t df2 = rd.Derivative2(x);
1099 Double_t weight = fKernel->GetWeight(x); // Bandwidth
1100 return 0.5 * fKernelSigmas2[fKernelType] * std::pow(weight, 2) * df2;
1101}
1103 // Returns the pointwise sigma of estimated density
1105 Double_t f = (*this)(x);
1106 Double_t weight = fKernel->GetWeight(x); // Bandwidth
1107 Double_t result = f * kernelL2Norm / (fNEvents * weight);
1108 return std::sqrt(result);
1109}
1110
1112 // Checks if kernel has unit integral, mu = 0 and positive finite sigma conditions
1115 valid = valid && unity == 1.;
1116 if (!valid) {
1117 Error("CheckKernelValidity", "Kernel's integral is %f",unity);
1118 }
1120 valid = valid && mu == 0.;
1121 if (!valid) {
1122 Error("CheckKernelValidity", "Kernel's mu is %f" ,mu);
1123 }
1125 valid = valid && sigma2 > 0 && sigma2 != std::numeric_limits<Double_t>::infinity();
1126 if (!valid) {
1127 Error("CheckKernelValidity", "Kernel's sigma2 is %f",sigma2);
1128 }
1129 if (!valid) {
1130 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.");
1131 //exit(EXIT_FAILURE);
1132 }
1133}
1134
1136 // Computes the kernel's L2 norm
1139 ig.SetFunction(kernel);
1140 Double_t result = ig.Integral();
1141 return result;
1142}
1143
1145 // Computes the kernel's sigma squared
1148 ig.SetFunction(kernel);
1149 Double_t result = ig.Integral();
1150 return result;
1151}
1152
1154 // Computes the kernel's mu
1157 ig.SetFunction(kernel);
1158 Double_t result = ig.Integral();
1159 return result;
1160}
1161
1163 // Computes the kernel's integral which ought to be unity
1166 ig.SetFunction(kernel);
1167 Double_t result = ig.Integral();
1168 return result;
1169}
1170
1171//////////////////////////////////////////////////////////////////////////////////////////////////////////////
1172/// Internal function to compute statistics (mean,stddev) using always all the provided data (i.e. no binning)
1174 /// in case of weights use
1175 if (!fEventWeights.empty() ) {
1176 // weighted data
1177 double x1 = fXMin - 0.001*(fXMax-fXMin);
1178 double x2 = fXMax + 0.001*(fXMax-fXMin);
1179 TH1D h1("temphist","", 500, x1, x2);
1180 h1.FillN(fEvents.size(), fEvents.data(), fEventWeights.data() );
1181 assert (h1.GetSumOfWeights() > 0) ;
1182 fMean = h1.GetMean();
1183 fSigma = h1.GetRMS();
1184 // compute robust sigma using midspread
1185 Double_t quantiles[2] = {0.0, 0.0};
1186 Double_t prob[2] = {0.25, 0.75};
1189 fSigmaRob = std::min(fSigma, midspread / 1.349); // Sigma's robust estimator
1190 //printf("weight case - stat: m = %f, s= %f, sr = %f \n",fMean, fSigma, midspread);
1191 return;
1192 }
1193 else {
1194 // compute statistics using the data
1195 SetMean();
1198 //printf("un-weight case - stat: m = %f, s= %f, sr = %f \n",fMean, fSigma, fSigmaRob);
1199 }
1200}
1201
1203 // Computes the inter-quartile range from the data
1204 std::sort(fEvents.begin(), fEvents.end());
1205 Double_t quantiles[2] = {0.0, 0.0};
1206 Double_t prob[2] = {0.25, 0.75};
1207 TMath::Quantiles(fEvents.size(), 2, &fEvents[0], quantiles, prob);
1211}
1212
1214 // Computes the user's input kernel function canonical bandwidth
1215 fCanonicalBandwidths[kUserDefined] = std::pow(ComputeKernelL2Norm() / std::pow(ComputeKernelSigma2(), 2), 1. / 5.);
1216}
1217
1219 // Computes the user's input kernel function sigma2
1221}
1222
1224 // Internal class constructor
1225}
1226
1228 // Internal class unary function
1229 if (fIntegralResult == kNorm) {
1230 return std::pow((*fKDE->fKernelFunction)(x), 2);
1231 } else if (fIntegralResult == kMu) {
1232 return x * (*fKDE->fKernelFunction)(x);
1233 } else if (fIntegralResult == kSigma2) {
1234 return std::pow(x, 2) * (*fKDE->fKernelFunction)(x);
1235 } else if (fIntegralResult == kUnitIntegration) {
1236 return (*fKDE->fKernelFunction)(x);
1237 } else {
1238 return -1;
1239 }
1240}
1241
1243 //Returns the estimated density
1244 TString name = "KDEFunc_"; name+= GetName();
1245 TString title = "KDE "; title+= GetTitle();
1246 if (xMin >= xMax) { xMin = fXMin; xMax = fXMax; }
1247 //Do not register the TF1 to global list
1248 TF1 * pdf = new TF1(name.Data(), this, xMin, xMax, 0, 1, TF1::EAddToList::kNo);
1249 if (npx > 0) pdf->SetNpx(npx);
1250 pdf->SetTitle(title);
1251 return pdf;
1252}
1253
1255 // Returns the upper estimated density
1256 TString name;
1257 name.Form("KDE_UpperCL%f5.3_%s",confidenceLevel,GetName());
1258 if (xMin >= xMax) { xMin = fXMin; xMax = fXMax; }
1260 upperPDF->SetParameter(0, confidenceLevel);
1261 if (npx > 0) upperPDF->SetNpx(npx);
1262 TF1 * f = (TF1*)upperPDF->Clone();
1263 delete upperPDF;
1264 return f;
1265}
1266
1268 // Returns the upper estimated density
1269 TString name;
1270 name.Form("KDE_LowerCL%f5.3_%s",confidenceLevel,GetName());
1271 if (xMin >= xMax) { xMin = fXMin; xMax = fXMax; }
1273 lowerPDF->SetParameter(0, confidenceLevel);
1274 if (npx > 0) lowerPDF->SetNpx(npx);
1275 TF1 * f = (TF1*)lowerPDF->Clone();
1276 delete lowerPDF;
1277 return f;
1278}
1279
1281 // Returns the approximate bias
1282 TString name = "KDE_Bias_"; name += GetName();
1283 if (xMin >= xMax) { xMin = fXMin; xMax = fXMax; }
1285 if (npx > 0) approximateBias->SetNpx(npx);
1286 TF1 * f = (TF1*)approximateBias->Clone();
1287 delete approximateBias;
1288 return f;
1289}
#define f(i)
Definition RSha256.hxx:104
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define M_PI
Definition Rotated.cxx:105
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
#define ClassImp(name)
Definition Rtypes.h:376
@ kRed
Definition Rtypes.h:67
@ kBlue
Definition Rtypes.h:67
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t option
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 Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
char name[80]
Definition TGX11.cxx:110
float xmin
float xmax
Int_t gDebug
Global variable setting the debug level. Set to 0 to disable, increase it in steps of 1 to increase t...
Definition TROOT.cxx:627
#define gPad
Interface (abstract class) for generic functions objects of one-dimension Provides a method to evalua...
Definition IFunction.h:170
User Class for performing numerical integration of a function in one dimension.
Definition Integrator.h:98
User class for calculating the derivatives of a function.
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 ...
const_iterator begin() const
const_iterator end() const
1-Dim function class
Definition TF1.h:182
void SetTitle(const char *title="") override
Set function title if title has the form "fffffff;xxxx;yyyy", it is assumed that the function title i...
Definition TF1.cxx:3590
virtual void SetNpx(Int_t npx=100)
Set the number of points used to draw the function.
Definition TF1.cxx:3465
void Draw(Option_t *option="") override
Draw this function with its current attributes.
Definition TF1.cxx:1341
A TGraphErrors is a TGraph with error bars.
void Draw(Option_t *chopt="") override
Draw this graph with its current attributes.
Definition TGraph.cxx:833
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:927
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:7559
virtual Double_t GetSumOfWeights() const
Return the sum of weights across all bins excluding under/overflows.
Definition TH1.h:560
Double_t GetRMS(Int_t axis=1) const
This function returns the Standard Deviation (Sigma) of the distribution not the Root Mean Square (RM...
Definition TH1.h:567
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:3419
virtual Int_t GetQuantiles(Int_t n, Double_t *xp, const Double_t *p=nullptr)
Compute Quantiles for this histogram.
Definition TH1.cxx:4593
void ComputeAdaptiveWeights()
Definition TKDE.cxx:784
Double_t GetFixedWeight() const
Definition TKDE.cxx:1011
TKernel(Double_t weight, TKDE *kde)
Definition TKDE.cxx:777
Double_t GetWeight(Double_t x) const
Definition TKDE.cxx:822
Double_t operator()(Double_t x) const
Definition TKDE.cxx:1021
const std::vector< Double_t > & GetAdaptiveWeights() const
Definition TKDE.cxx:1016
Kernel Density Estimation class.
Definition TKDE.h:37
~TKDE() override
Definition TKDE.cxx:77
TF1 * GetPDFUpperConfidenceInterval(Double_t confidenceLevel=0.95, UInt_t npx=100, Double_t xMin=1.0, Double_t xMax=0.0)
Definition TKDE.cxx:1254
TF1 * GetKDEApproximateBias(UInt_t npx=100, Double_t xMin=1.0, Double_t xMax=0.0)
Definition TKDE.cxx:1280
void SetData(const Double_t *data, const Double_t *weights)
Definition TKDE.cxx:442
TF1 * fLowerPDF
Output Kernel Density Estimation upper confidence interval PDF function.
Definition TKDE.h:195
std::vector< Double_t > fKernelSigmas2
Definition TKDE.h:226
Double_t ComputeKernelL2Norm() const
Definition TKDE.cxx:1135
TF1 * fPDF
Definition TKDE.h:193
TF1 * GetPDFLowerConfidenceInterval(Double_t confidenceLevel=0.95, UInt_t npx=100, Double_t xMin=1.0, Double_t xMax=0.0)
Definition TKDE.cxx:1267
void SetKernelType(EKernelType kern)
Definition TKDE.cxx:309
std::vector< Double_t > fCanonicalBandwidths
Definition TKDE.h:225
UInt_t fNEvents
Data's number of events.
Definition TKDE.h:211
void ComputeDataStats()
Internal function to compute statistics (mean,stddev) using always all the provided data (i....
Definition TKDE.cxx:1173
Double_t fXMax
Data maximum value.
Definition TKDE.h:219
Double_t UpperConfidenceInterval(const Double_t *x, const Double_t *p) const
Valid if the bandwidth is small compared to nEvents**1/5.
Definition TKDE.cxx:1074
void ReInit()
Definition TKDE.cxx:478
Double_t ApproximateBias(const Double_t *x, const Double_t *) const
Definition TKDE.h:259
Double_t ComputeMidspread()
Definition TKDE.cxx:1202
Bool_t fNewData
Flag to control when new data are given.
Definition TKDE.h:207
void SetMirror()
Definition TKDE.cxx:433
Bool_t fUseMirroring
Definition TKDE.h:205
void DrawConfidenceInterval(TString &drawOpt, double cl=0.95)
// Draws the KDE and its confidence interval
Definition TKDE.cxx:976
void SetMirroredEvents()
Intgernal function to mirror the data.
Definition TKDE.cxx:534
EMirror fMirror
Definition TKDE.h:201
void SetUserCanonicalBandwidth()
Definition TKDE.cxx:1213
Bool_t fMirrorLeft
Definition TKDE.h:205
EIteration fIteration
Definition TKDE.h:200
void CheckKernelValidity()
Definition TKDE.cxx:1111
const Double_t * GetAdaptiveWeights() const
Definition TKDE.cxx:1001
Double_t fAdaptiveBandwidthFactor
Geometric mean of the kernel density estimation from the data for adaptive iteration.
Definition TKDE.h:221
void SetKernelFunction(KernelFunction_Ptr kernfunc=nullptr)
Definition TKDE.cxx:633
Double_t LowerConfidenceInterval(const Double_t *x, const Double_t *p) const
Valid if the bandwidth is small compared to nEvents**1/5.
Definition TKDE.cxx:1083
Double_t fSigmaRob
Data std deviation (robust estimation)
Definition TKDE.h:217
std::vector< Double_t > fBinCount
Number of events per bin for binned data option.
Definition TKDE.h:228
void Draw(const Option_t *option="") override
Draws either the KDE functions or its errors.
Definition TKDE.cxx:906
EBinning fBinning
Definition TKDE.h:202
EIteration
Iteration types. They can be set using SetIteration()
Definition TKDE.h:52
@ kAdaptive
Definition TKDE.h:53
@ kFixed
Definition TKDE.h:54
Double_t GetRAMISE() const
Definition TKDE.cxx:771
void SetIteration(EIteration iter)
Definition TKDE.cxx:321
Double_t ComputeKernelIntegral() const
Definition TKDE.cxx:1162
Double_t CosineArchKernel(Double_t x) const
Returns the kernel evaluation at x.
Definition TKDE.h:254
Double_t fXMin
Data minimum value.
Definition TKDE.h:218
Double_t operator()(Double_t x) const
Definition TKDE.cxx:749
void SetUserKernelSigma2()
Definition TKDE.cxx:1218
Double_t GetBias(Double_t x) const
Definition TKDE.cxx:1093
std::vector< Double_t > fData
Data events.
Definition TKDE.h:189
Double_t fSumOfCounts
Data sum of weights.
Definition TKDE.h:212
UInt_t fUseBinsNEvents
If the algorithm is allowed to use automatic (relaxed) binning this is the minimum number of events t...
Definition TKDE.h:213
void SetKernel()
Definition TKDE.cxx:604
Double_t fSigma
Data std deviation.
Definition TKDE.h:216
TGraphErrors * GetGraphWithErrors(UInt_t npx=100, Double_t xMin=1.0, Double_t xMax=0.0)
return a TGraphErrors with the KDE values and errors The return object is managed by the user
Definition TKDE.cxx:951
Double_t GetMean() const
Definition TKDE.cxx:759
Double_t fRho
Adjustment factor for sigma.
Definition TKDE.h:220
Bool_t fAsymRight
Definition TKDE.h:205
Double_t fWeightSize
Caches the weight size.
Definition TKDE.h:223
void SetUseBinsNEvents(UInt_t nEvents)
Definition TKDE.cxx:364
std::vector< Double_t > fEvents
Original data storage.
Definition TKDE.h:190
Double_t GetError(Double_t x) const
Definition TKDE.cxx:1102
TF1 * GetKDEFunction(UInt_t npx=100, Double_t xMin=1.0, Double_t xMax=0.0)
Definition TKDE.cxx:1242
void SetBinning(EBinning)
Definition TKDE.cxx:339
void GetOptions(std::string optionType, std::string option)
Definition TKDE.cxx:198
Bool_t fAsymLeft
Definition TKDE.h:205
std::vector< Bool_t > fSettedOptions
User input options flag.
Definition TKDE.h:230
Double_t GaussianKernel(Double_t x) const
Returns the kernel evaluation at x.
Definition TKDE.h:239
void SetRange(Double_t xMin, Double_t xMax)
By default computed from the data.
Definition TKDE.cxx:380
void SetMean()
Definition TKDE.cxx:593
Double_t ComputeKernelSigma2() const
Definition TKDE.cxx:1144
void SetOptions(const Option_t *option, Double_t rho)
Definition TKDE.cxx:126
void SetUseBins()
Definition TKDE.cxx:394
Double_t GetFixedWeight() const
Definition TKDE.cxx:990
void AssureOptions()
Definition TKDE.cxx:270
TF1 * GetFunction(UInt_t npx=100, Double_t xMin=1.0, Double_t xMax=0.0)
Definition TKDE.cxx:697
void SetBinCountData()
Definition TKDE.cxx:836
TF1 * GetUpperFunction(Double_t confidenceLevel=0.95, UInt_t npx=100, Double_t xMin=1.0, Double_t xMax=0.0)
Definition TKDE.cxx:706
void InitFromNewData()
Definition TKDE.cxx:500
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:91
void SetDrawOptions(const Option_t *option, TString &plotOpt, TString &drawOpt)
Definition TKDE.cxx:152
EMirror
Data "mirroring" option to address the probability "spill out" boundary effect They can be set using ...
Definition TKDE.h:59
@ kMirrorRightAsymLeft
Definition TKDE.h:65
@ kMirrorLeft
Definition TKDE.h:61
@ kMirrorAsymRight
Definition TKDE.h:66
@ kMirrorAsymBoth
Definition TKDE.h:68
@ kNoMirror
Definition TKDE.h:60
@ kMirrorRight
Definition TKDE.h:62
@ kMirrorLeftAsymRight
Definition TKDE.h:67
@ kMirrorBoth
Definition TKDE.h:63
@ kMirrorAsymLeft
Definition TKDE.h:64
TGraphErrors * fGraph
Output Kernel Density Estimation approximate bias.
Definition TKDE.h:197
void SetCanonicalBandwidths()
Definition TKDE.cxx:680
TKDE()
default constructor used only by I/O
Definition TKDE.cxx:61
void SetBinCentreData(Double_t xmin, Double_t xmax)
Definition TKDE.cxx:827
void SetTuneFactor(Double_t rho)
Definition TKDE.cxx:370
UInt_t Index(Double_t x) const
compute the bin index given a data point x
Definition TKDE.cxx:1061
void Fill(Double_t data)
Definition TKDE.cxx:721
TF1 * fUpperPDF
Output Kernel Density Estimation PDF function.
Definition TKDE.h:194
Bool_t fMirrorRight
Definition TKDE.h:205
UInt_t fNBins
Number of bins for binned data option.
Definition TKDE.h:210
Double_t ComputeKernelMu() const
Definition TKDE.cxx:1153
EKernelType
Types of Kernel functions They can be set using the function SetKernelType() or as a string in the co...
Definition TKDE.h:42
@ kGaussian
Definition TKDE.h:43
@ kEpanechnikov
Definition TKDE.h:44
@ kCosineArch
Definition TKDE.h:46
@ kBiweight
Definition TKDE.h:45
@ kTotalKernels
Internal use only for member initialization.
Definition TKDE.h:48
@ kUserDefined
Internal use only for the class's template constructor.
Definition TKDE.h:47
Double_t fMean
Data mean.
Definition TKDE.h:215
void DrawErrors(TString &drawOpt)
Draws a TGraphErrors with KDE values and errors.
Definition TKDE.cxx:942
void SetNBins(UInt_t nbins)
Definition TKDE.cxx:346
void CheckOptions(Bool_t isUserDefinedKernel=kFALSE)
Definition TKDE.cxx:286
Double_t EpanechnikovKernel(Double_t x) const
Definition TKDE.h:244
Double_t BiweightKernel(Double_t x) const
Returns the kernel evaluation at x.
Definition TKDE.h:249
void SetKernelSigmas2()
Definition TKDE.cxx:689
std::unique_ptr< TKernel > fKernel
! internal kernel class. Transient because it is recreated after reading from a file
Definition TKDE.h:187
Bool_t fUseMinMaxFromData
Flag top control if min and max must be used from data.
Definition TKDE.h:208
Double_t GetSigma() const
Definition TKDE.cxx:765
EKernelType fKernelType
Graph with the errors.
Definition TKDE.h:199
TF1 * fApproximateBias
Output Kernel Density Estimation lower confidence interval PDF function.
Definition TKDE.h:196
void SetSigma(Double_t R)
Definition TKDE.cxx:598
std::vector< Double_t > fEventWeights
Original data weights.
Definition TKDE.h:191
TF1 * GetApproximateBias(UInt_t npx=100, Double_t xMin=1.0, Double_t xMax=0.0)
Definition TKDE.cxx:716
TF1 * GetLowerFunction(Double_t confidenceLevel=0.95, UInt_t npx=100, Double_t xMin=1.0, Double_t xMax=0.0)
Definition TKDE.cxx:711
Bool_t fUseBins
Definition TKDE.h:206
KernelFunction_Ptr fKernelFunction
! pointer to kernel function
Definition TKDE.h:185
friend struct KernelIntegrand
Definition TKDE.h:233
EBinning
Data binning option.
Definition TKDE.h:73
@ kUnbinned
Definition TKDE.h:74
@ kRelaxedBinning
The algorithm is allowed to use binning if the data is large enough.
Definition TKDE.h:75
@ kForcedBinning
Definition TKDE.h:76
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:50
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:1058
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1072
virtual void Fatal(const char *method, const char *msgfmt,...) const
Issue fatal error message.
Definition TObject.cxx:1100
virtual void Info(const char *method, const char *msgfmt,...) const
Issue info message.
Definition TObject.cxx:1046
Basic string class.
Definition TString.h:138
void ToLower()
Change string to lower-case.
Definition TString.cxx:1190
const char * Data() const
Definition TString.h:384
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:640
@ kGAUSS
simple Gauss integration method with fixed rule
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
Double_t ey[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
Double_t ex[n]
Definition legend1.C:17
TH1F * h1
Definition legend1.C:5
Bool_t IsNaN(Double_t x)
Definition TMath.h:898
Double_t QuietNaN()
Returns a quiet NaN as defined by IEEE 754.
Definition TMath.h:908
void Quantiles(Int_t n, Int_t nprob, Double_t *x, Double_t *quantiles, Double_t *prob, Bool_t isSorted=kTRUE, Int_t *index=nullptr, Int_t type=7)
Computes sample quantiles, corresponding to the given probabilities.
Definition TMath.cxx:1207
const TKDE * fKDE
Definition TKDE.cxx:55
KernelIntegrand(const TKDE *kde, EIntegralResult intRes)
Definition TKDE.cxx:1223
Double_t operator()(Double_t x) const
Definition TKDE.cxx:1227
EIntegralResult fIntegralResult
Definition TKDE.cxx:56