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