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