KaliVeda
Toolkit for HIC analysis
Loading...
Searching...
No Matches
bayesian_estimator.h
1
3
4#ifndef __bayesian_estimator_H
5#define __bayesian_estimator_H
6
7#include "KVBase.h"
8#include "impact_parameter_distribution.h"
9
10#include <vector>
11#include <numeric>
12#include <TF2.h>
13#include <TH1.h>
14#include <TCanvas.h>
15#include <TGraph.h>
16#include <KVValueRange.h>
17#include "TMath.h"
18#include "TNamed.h"
19
21
22namespace KVImpactParameters {
39 public:
40 double operator()(double X, double mean, double reduced_variance)
41 {
43 return ROOT::Math::gamma_pdf(X, mean / reduced_variance, reduced_variance);
44 }
45
47 };
48
64 class BD_kernel {
65 public:
66 double operator()(double X, double mean, double reduced_variance)
67 {
69 return ROOT::Math::binomial_pdf(TMath::Nint(X), 1. - reduced_variance, mean / (1. - reduced_variance));
70 }
71
73 };
74
90 class NBD_kernel {
91 public:
92 double operator()(double X, double mean, double reduced_variance)
93 {
95 return ROOT::Math::negative_binomial_pdf(TMath::Nint(X), 1. / reduced_variance, mean / (reduced_variance - 1.));
96 }
97
99 };
100
101
216 template
217 <class FittingFunction, class FluctuationKernel>
218 class bayesian_estimator : public KVBase {
219
220 FittingFunction theFitter;
221 FluctuationKernel theKernel;
222
225 std::vector<double> sel_rapp;
238
240
241 impact_parameter_distribution fIPDist;// used to convert centrality to impact parameter
242
243 double mean_X_vs_cb(double* x, double* par)
244 {
249
250 theFitter.fill_params_from_array(par);
251 return theFitter.meanX(x[0]);
252 }
253 double mean_X_vs_b(double* x, double* par)
254 {
259
260 theFitter.fill_params_from_array(par);
261 return theFitter.meanX(fIPDist.GetCentrality().Eval(x[0]));
262 }
263 double P_X_cb(double X, double cb)
264 {
268 return fIntegerVariable ? theKernel(TMath::Nint(X), theFitter.meanX(cb), theFitter.redVar(cb))
269 : theKernel(X, theFitter.meanX(cb), theFitter.redVar(cb));
270 }
271 double P_X_cb_for_TF2_obs_vs_b(double* x, double*)
272 {
276
277 auto X = x[1];
278 auto cb = fIPDist.GetCentrality().Eval(x[0]);
279 return x[0] * theKernel(X, theFitter.meanX(cb), theFitter.redVar(cb));
280 }
281 double P_X_cb_for_integral(double* x, double* par)
282 {
289 return P_X_cb(par[0], x[0]);
290 }
291 double P_X_cb_for_X_integral(double* x, double* par)
292 {
300 return P_X_cb(x[0], par[0]);
301 }
302 double P_X_cb_for_X_integral_with_selection(double* x, double* par)
303 {
311
312 int bin = histo->FindBin(x[0]);
313 if (bin < 1 || bin > histo->GetNbinsX()) return 0;
314 return sel_rapp[bin - 1] * P_X_cb(x[0], par[0]);
315 }
316 double cb_integrated_P_X(double* x, double* p)
317 {
325
326 theFitter.fill_params_from_array(p);
327 p_X_cb_integrator.SetParameter(0, x[0]); // set value of X in integrand
328 return p_X_cb_integrator.Integral(0, 1, 1.e-4);
329 }
330 double P_X_from_fit(double* x, double* par)
331 {
336 p_X_cb_integrator.SetParameter(0, x[0]); // set value of X in integrand
337 return par[0] * p_X_cb_integrator.Integral(0, 1, 1.e-4);
338 }
339 double cb_dist_for_X_selection(double* x, double* p)
340 {
350
352 double num = p_X_X_integrator.Integral(p[0], p[1], 1.e-4);
353 double den = fitted_P_X.Integral(p[0], p[1], 1.e-4);
354 if (den > 0) return num / den;
355 return 0;
356 }
357 double cb_dist_for_arb_X_selection(double* x, double* p)
358 {
367
369 double num = p_X_X_integrator_with_selection.Integral(p[0], p[1], 1.e-4);
371
372 return num;
373 }
374 double b_dist_for_X_selection(double* x, double* p)
375 {
385
387 double num = p_X_X_integrator.Integral(p[0], p[1], 1.e-4);
388 return num * fIPDist.GetDifferentialCrossSection(x[0]);
389 }
390 double b_dist_for_arb_X_selection(double* x, double* p)
391 {
400
402 double num = p_X_X_integrator_with_selection.Integral(p[0], p[1], 1.e-4);
403 return num * fIPDist.GetDifferentialCrossSection(x[0]);
404 }
405 double b_dist_for_arb_X_selection_from_histo(double* x, double* p)
406 {
415 TH1 *hh_centb = fIPDist.GetCentralityFromHisto();
416 Double_t centb = hh_centb->GetBinContent(hh_centb->FindBin(x[0]));
417
419 double num = p_X_X_integrator_with_selection.Integral(p[0], p[1], 1.e-4);
421 }
422
423 public:
424 bayesian_estimator(Bool_t integer_variable = false)
425 : KVBase(),
426 theFitter(),
427 p_X_cb_integrator("p_X_cb_integrator", this, &bayesian_estimator::P_X_cb_for_integral, 0, 1, 1),
428 P_X_fit_function("P_X_fit_function", this, &bayesian_estimator::cb_integrated_P_X, 0, 1, theFitter.npar()),
429 mean_X_vs_cb_function("mean_X_vs_cb", this, &bayesian_estimator::mean_X_vs_cb, 0, 1, theFitter.npar()),
430 mean_X_vs_b_function("mean_X_vs_b", this, &bayesian_estimator::mean_X_vs_b, 0, 20, theFitter.npar()),
431 p_X_X_integrator("p_X_X_integrator", this, &bayesian_estimator::P_X_cb_for_X_integral, 0, 1000, 1),
432 p_X_X_integrator_with_selection("p_X_X_integrator_with_selection", this, &bayesian_estimator::P_X_cb_for_X_integral_with_selection, 0, 1000, 1),
433 fitted_P_X("fitted_P_X", this, &bayesian_estimator::P_X_from_fit, 0, 1000, 1),
434 Cb_dist_for_X_select("Cb_dist_for_X_select", this, &bayesian_estimator::cb_dist_for_X_selection, 0, 1, 2),
435 Cb_dist_for_arb_X_select("Cb_dist_for_arb_X_select", this, &bayesian_estimator::cb_dist_for_arb_X_selection, 0, 1, 2),
436 B_dist_for_X_select("b_dist_for_X_select", this, &bayesian_estimator::b_dist_for_X_selection, 0, 20, 2),
437 B_dist_for_arb_X_select("b_dist_for_arb_X_select", this, &bayesian_estimator::b_dist_for_arb_X_selection, 0, 20, 2),
438 B_dist_for_arb_X_select_from_histo("b_dist_for_arb_X_select_from_histo", this, &bayesian_estimator::b_dist_for_arb_X_selection_from_histo, 0, 20, 2),
439 fIntegerVariable(integer_variable)
440 {
445 theFitter.set_par_names(P_X_fit_function);
446 theFitter.set_par_names(mean_X_vs_cb_function);
447 theFitter.set_par_names(mean_X_vs_b_function);
448 }
449
450 bayesian_estimator(const FittingFunction& previous_fit, Bool_t integer_variable = false)
451 : KVBase(),
452 theFitter(previous_fit),
453 p_X_cb_integrator("p_X_cb_integrator", this, &bayesian_estimator::P_X_cb_for_integral, 0, 1, 1),
454 P_X_fit_function("P_X_fit_function", this, &bayesian_estimator::cb_integrated_P_X, 0, 1, theFitter.npar()),
455 mean_X_vs_cb_function("mean_X_vs_cb", this, &bayesian_estimator::mean_X_vs_cb, 0, 1, theFitter.npar()),
456 mean_X_vs_b_function("mean_X_vs_b", this, &bayesian_estimator::mean_X_vs_b, 0, 20, theFitter.npar()),
457 p_X_X_integrator("p_X_X_integrator", this, &bayesian_estimator::P_X_cb_for_X_integral, 0, 1000, 1),
458 p_X_X_integrator_with_selection("p_X_X_integrator_with_selection", this, &bayesian_estimator::P_X_cb_for_X_integral_with_selection, 0, 1000, 1),
459 fitted_P_X("fitted_P_X", this, &bayesian_estimator::P_X_from_fit, 0, 1000, 1),
460 Cb_dist_for_X_select("Cb_dist_for_X_select", this, &bayesian_estimator::cb_dist_for_X_selection, 0, 1, 2),
461 Cb_dist_for_arb_X_select("Cb_dist_for_arb_X_select", this, &bayesian_estimator::cb_dist_for_arb_X_selection, 0, 1, 2),
462 B_dist_for_X_select("b_dist_for_X_select", this, &bayesian_estimator::b_dist_for_X_selection, 0, 20, 2),
463 B_dist_for_arb_X_select("b_dist_for_arb_X_select", this, &bayesian_estimator::b_dist_for_arb_X_selection, 0, 20, 2),
464 B_dist_for_arb_X_select_from_histo("b_dist_for_arb_X_select_from_histo", this, &bayesian_estimator::b_dist_for_arb_X_selection_from_histo, 0, 20, 2),
465 fIntegerVariable(integer_variable)
466 {
472 theFitter.set_par_names(P_X_fit_function);
473 theFitter.set_par_names(mean_X_vs_cb_function);
474 theFitter.set_par_names(mean_X_vs_b_function);
475 }
476
478
480 {
486
487 histo = h;
488 histo->Scale(1. / h->Integral("width"));
489 }
490
491 void FitHisto(TH1* h = nullptr)
492 {
503
504 if (h) histo = h;
506 theFitter.set_initial_parameters(histo, P_X_fit_function);
508 }
509
510 void SetIPDistParams(double sigmaR, double deltab)
511 {
528
530 }
531
532 void SetIPDistFromHisto(TH1 *ip_histo)
533 {
547
548 fIPDist.FitIPDist(ip_histo);
549 }
550
561
562 void DrawMeanXvsCb(const TString& title = "", Color_t color = -1, Option_t* opt = "")
563 {
572
573 GetMeanXvsCb().Draw();
574 }
576 {
583
584 Double_t par[theFitter.npar()];
585 theFitter.fill_array_from_params(par);
588 }
589 TGraph* GraphMeanXvsb(int npts = 500)
590 {
600
601 Double_t par[theFitter.npar()];
602 theFitter.fill_array_from_params(par);
603 KVValueRange<double> brange(0, GetIPDist().GetB0() + 4 * GetIPDist().GetDeltaB());
605 auto gr = new TGraph;
606 for (int i = 0; i <= npts; ++i) gr->SetPoint(i, brange.ValueIofN(i, npts), mean_X_vs_b_function.Eval(brange.ValueIofN(i, npts)));
607 return gr;
608 }
610 {
615
616 if (!histo) {
617 Warning("update_fit_params", "no histogram set with FitHisto(TH1*)");
618 return;
619 }
620 TF1* fit = (TF1*)histo->FindObject("P_X_fit_function");
621 if (!fit) {
622 Warning("update_fit_params", "no fit function found in histogram");
623 return;
624 }
625 theFitter.fill_params_from_array(fit->GetParameters());
626 }
627 double DrawCbDistForXSelection(double X1, double X2, Option_t* opt = "", Color_t color = kRed, const TString& title = "")
628 {
643 f->SetNpx(500);
644 f->SetLineColor(color);
645 f->SetMarkerColor(color);
646 f->SetLineWidth(2);
647 f->SetTitle(title);
648 return f->GetMaximum();
649 }
650 void DrawCbDistForSelection(TH1* sel, TH1* incl, double& mean_cb, double& sigma_cb, Option_t* opt = "", Color_t color = kRed, const TString& title = "")
651 {
666
667 assert(sel->GetNbinsX() == incl->GetNbinsX());
668
671 sel_rapp.assign((std::vector<double>::size_type)incl->GetNbinsX(), 0.0);
672 int first_bin(0), last_bin(0);
673 for (int i = 1; i <= incl->GetNbinsX(); ++i) {
674 if (incl->GetBinContent(i) > 0) {
675 sel_rapp[i - 1] = sel->GetBinContent(i) / incl->GetBinContent(i);
676 if (sel->GetBinContent(i) > 0) {
677 if (!first_bin) first_bin = i;
678 last_bin = i;
679 }
680 }
681 }
682 double Xmin = incl->GetBinLowEdge(first_bin);
683 double Xmax = incl->GetXaxis()->GetBinUpEdge(last_bin);
684
685 histo = incl;
687
692
694 double cb_mean(0), cb_sqrmean(0), sum_pcb(0);
695 TGraph* f = new TGraph;
696 for (int i = 0; i < 500; ++i) {//graph with 500 points
697 double cb = i / 499.;
698 double p_cb = Cb_dist_for_arb_X_select.Eval(cb);
699
700 cb_mean += p_cb * cb;
701 cb_sqrmean += p_cb * cb * cb;
702 sum_pcb += p_cb;
703
704 f->SetPoint(i, cb, p_cb);
705 }
706
708 cb_sqrmean /= sum_pcb;
709 cb_mean /= sum_pcb;
710 mean_cb = cb_mean;
711 sigma_cb = TMath::Sqrt(cb_sqrmean - cb_mean * cb_mean);
712
714 f->SetLineColor(color);
715 f->SetMarkerColor(color);
716 f->SetLineWidth(2);
717 f->SetTitle(title);
718 if (TString(opt) == "same") f->Draw("l");
719 else f->Draw("al");
720 }
721 double DrawBDistForXSelection(KVValueRange<double> Xrange, Option_t* opt = "", Color_t color = kRed, const TString& title = "")
722 {
739
740 B_dist_for_X_select.SetParameters(Xrange.Min(), Xrange.Max());
741 TGraph* f = new TGraph;
742 double maxS = 0;
743 for (int i = 0; i < 500; ++i) {
744 double b = 2 * i * GetIPDist().GetB0() / 499.;
745 double sig = B_dist_for_X_select.Eval(b);
746 if (sig > maxS) maxS = sig;
747 f->SetPoint(i, b, sig);
748 }
749 f->SetLineColor(color);
750 f->SetMarkerColor(color);
751 f->SetLineWidth(2);
752 f->SetTitle(title);
753 if (TString(opt) == "same") f->Draw("l");
754 else f->Draw("al");
755 return maxS;
756 }
762 {
775
776
777 B_dist_for_X_select.SetParameters(Xrange.Min(), Xrange.Max());
778 TGraph* f = new TGraph;
779 for (int i = 0; i < npts; ++i) {
780 double b = 2 * i * GetIPDist().GetB0() / 499.;
781 double sig = B_dist_for_X_select.Eval(b);
782 f->SetPoint(i, b, sig);
783 }
784 f->SetLineWidth(2);
785 return f;
786 }
787 void DrawBDistForSelection(TH1* sel, TH1* incl, double& mean, double& sigma, Option_t* opt = "", Color_t color = kRed, const TString& title = "")
788 {
808
809
810 assert(sel->GetNbinsX() == incl->GetNbinsX());
811
812 sel_rapp.assign((std::vector<double>::size_type)incl->GetNbinsX(), 0.0);
813 int first_bin(0), last_bin(0);
814 for (int i = 1; i <= incl->GetNbinsX(); ++i) {
815 if (incl->GetBinContent(i) > 0) {
816 sel_rapp[i - 1] = sel->GetBinContent(i) / incl->GetBinContent(i);
817 if (sel->GetBinContent(i) > 0) {
818 if (!first_bin) first_bin = i;
819 last_bin = i;
820 }
821 }
822 }
823 double Xmin = incl->GetBinLowEdge(first_bin);
824 double Xmax = incl->GetXaxis()->GetBinUpEdge(last_bin);
825
826 histo = incl;
828
832
833 double bmean(0), bsqrmean(0), sigtot(0);
834 TGraph* f = new TGraph;
835 for (int i = 0; i < 500; ++i) {
836 double b = 2 * i * GetIPDist().GetB0() / 499.;
837 double sig = B_dist_for_arb_X_select.Eval(b);
838 bmean += sig * b;
839 bsqrmean += sig * b * b;
840 sigtot += sig;
841 f->SetPoint(i, b, sig);
842 }
843 f->SetLineColor(color);
844 f->SetMarkerColor(color);
845 f->SetLineWidth(2);
846 f->SetTitle(title);
847 mean = bmean / sigtot;
848 bsqrmean /= sigtot;
849 sigma = TMath::Sqrt(bsqrmean - mean * mean);
850 if (TString(opt) == "same") f->Draw("l");
851 else f->Draw("al");
852 }
853
854 void DrawBDistForSelectionFromIPHisto(TH1* sel, TH1* incl, double& mean, double& sigma, Option_t* opt = "", Color_t color = kRed, const TString& title = "")
855 {
873
874 assert(sel->GetNbinsX() == incl->GetNbinsX());
875
876 sel_rapp.assign((std::vector<double>::size_type)incl->GetNbinsX(), 0.0);
877 int first_bin(0), last_bin(0);
878 for (int i = 1; i <= incl->GetNbinsX(); ++i) {
879 if (incl->GetBinContent(i) > 0) {
880 sel_rapp[i - 1] = sel->GetBinContent(i) / incl->GetBinContent(i);
881 if (sel->GetBinContent(i) > 0) {
882 if (!first_bin) first_bin = i;
883 last_bin = i;
884 }
885 }
886 }
887 double Xmin = incl->GetBinLowEdge(first_bin);
888 double Xmax = incl->GetXaxis()->GetBinUpEdge(last_bin);
889
890 histo = incl;
892
894 B_dist_for_arb_X_select_from_histo.SetRange(0, 2 * GetIPDist().GetB0()); //if SetIPDistFromHisto() was called, a fit have been applied thus B0 was estimated
896
897 double bmean(0), bsqrmean(0), sigtot(0);
898 TGraph* f = new TGraph;
899 for (int i = 0; i < 500; ++i) {
900 double b = 2 * i * GetIPDist().GetB0() / 499.;
902 bmean += sig * b;
903 bsqrmean += sig * b * b;
904 sigtot += sig;
905 f->SetPoint(i, b, sig);
906 }
907 f->SetLineColor(color);
908 f->SetMarkerColor(color);
909 f->SetLineWidth(2);
910 f->SetTitle(title);
911 mean = bmean / sigtot;
912 bsqrmean /= sigtot;
913 sigma = TMath::Sqrt(bsqrmean - mean * mean);
914 if (TString(opt) == "same") f->Draw("l");
915 else f->Draw("al");
916 }
917
918 void GetMeanAndSigmaBDistForSelection(TH1* sel, TH1* incl, double& mean, double& sigma)
919 {
936
937
938 assert(sel->GetNbinsX() == incl->GetNbinsX());
939
940 sel_rapp.assign((std::vector<double>::size_type)incl->GetNbinsX(), 0.0);
941 int first_bin(0), last_bin(0);
942 for (int i = 1; i <= incl->GetNbinsX(); ++i) {
943 if (incl->GetBinContent(i) > 0) {
944 sel_rapp[i - 1] = sel->GetBinContent(i) / incl->GetBinContent(i);
945 if (sel->GetBinContent(i) > 0) {
946 if (!first_bin) first_bin = i;
947 last_bin = i;
948 }
949 }
950 }
951 double Xmin = incl->GetBinLowEdge(first_bin);
952 double Xmax = incl->GetXaxis()->GetBinUpEdge(last_bin);
953
954 histo = incl;
956
958 mean = B_dist_for_arb_X_select.Mean(0, 20);
959 double var = B_dist_for_arb_X_select.Variance(0, 20);
960 sigma = TMath::Sqrt(var);
961 }
962
963 void GetMeanAndSigmaBDistForXSelection(double X1, double X2, double& mean, double& sigma)
964 {
979
980
982 mean = B_dist_for_X_select.Mean(0, 20);
984 }
985 void DrawFittedP_X(double norm = 1.0, Option_t* opt = "", Color_t color = kRed, const TString& title = "")
986 {
993
994 TF1* f = GetFittedP_X(norm)->DrawCopy(opt);
995 f->SetNpx(500);
996 f->SetLineColor(color);
997 f->SetLineWidth(2);
998 f->SetTitle(title);
999 }
1000 TF1* GetFittedP_X(double norm = 1.0)
1001 {
1006 fitted_P_X.SetParameter(0, norm);
1007 return &fitted_P_X;
1008 }
1009 TGraph* GraphP_XForGivenB(double b, KVValueRange<double> Xrange, int npts = 500)
1010 {
1017 auto graph = new TGraph;
1018 for (int i = 0; i < npts; ++i) graph->SetPoint(i, Xrange.ValueIofN(i, npts), p_X_X_integrator.Eval(Xrange.ValueIofN(i, npts)));
1019 return graph;
1020 }
1021 void Print(Option_t* = "") const
1022 {
1024 theFitter.print_fit_params();
1025 }
1026 void DrawNormalisedMeanXvsb(const TString& title, Color_t color, Option_t* opt)
1027 {
1043
1044
1045 theFitter.backup_params();
1046 theFitter.normalise_shape_function();
1047 Double_t par[5];
1048 theFitter.fill_array_from_params(par);
1049 mean_X_vs_b_function.SetRange(0, GetIPDist().GetB0() + 2 * GetIPDist().GetDeltaB());
1051 TF1* copy = mean_X_vs_b_function.DrawCopy(opt);
1052 if (!title.IsNull()) copy->SetTitle(title);
1053 if (color >= 0) copy->SetLineColor(color);
1054 theFitter.restore_params();
1055 }
1056
1058 {
1061
1062 auto jpd = new TF2("joint_proba_dist", this, &bayesian_estimator::P_X_cb_for_TF2_obs_vs_b,
1063 b_range.Min(), b_range.Max(), X_range.Min(), X_range.Max(), 0);
1064 return jpd;
1065 }
1066
1067 ClassDef(bayesian_estimator, 1) //Estimate impact parameter distribution by fits to data
1068 };
1069
1070}
1071
1072#endif
#define f(i)
bool Bool_t
short Color_t
double Double_t
const char Option_t
#define ClassDef(name, id)
#define X(type, name)
winID h TVirtualViewer3D TVirtualGLPainter p
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 Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t sel
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 b
Base class for KaliVeda framework.
Definition KVBase.h:142
Fluctuation kernel using binomial distribution for use with bayesian_estimator.
double operator()(double X, double mean, double reduced_variance)
Fluctuation kernel using negative binomial distribution for use with bayesian_estimator.
double operator()(double X, double mean, double reduced_variance)
Impact parameter distribution reconstruction from experimental data.
double P_X_cb_for_integral(double *x, double *par)
impact_parameter_distribution fIPDist
void DrawBDistForSelection(TH1 *sel, TH1 *incl, double &mean, double &sigma, Option_t *opt="", Color_t color=kRed, const TString &title="")
void DrawCbDistForSelection(TH1 *sel, TH1 *incl, double &mean_cb, double &sigma_cb, Option_t *opt="", Color_t color=kRed, const TString &title="")
double mean_X_vs_b(double *x, double *par)
void GetMeanAndSigmaBDistForXSelection(double X1, double X2, double &mean, double &sigma)
double DrawBDistForXSelection(KVValueRange< double > Xrange, Option_t *opt="", Color_t color=kRed, const TString &title="")
void GetMeanAndSigmaBDistForSelection(TH1 *sel, TH1 *incl, double &mean, double &sigma)
bayesian_estimator(Bool_t integer_variable=false)
double cb_dist_for_X_selection(double *x, double *p)
void SetIPDistParams(double sigmaR, double deltab)
double b_dist_for_arb_X_selection_from_histo(double *x, double *p)
double P_X_cb_for_TF2_obs_vs_b(double *x, double *)
TGraph * GraphBDistForXSelection(KVValueRange< double > Xrange, int npts=500)
double P_X_from_fit(double *x, double *par)
double b_dist_for_arb_X_selection(double *x, double *p)
double DrawCbDistForXSelection(double X1, double X2, Option_t *opt="", Color_t color=kRed, const TString &title="")
bayesian_estimator(const FittingFunction &previous_fit, Bool_t integer_variable=false)
impact_parameter_distribution & GetIPDist()
double b_dist_for_X_selection(double *x, double *p)
TF2 * GetJointProbabilityDistribution(KVValueRange< double > b_range, KVValueRange< double > X_range)
double cb_integrated_P_X(double *x, double *p)
double P_X_cb_for_X_integral(double *x, double *par)
void DrawNormalisedMeanXvsb(const TString &title, Color_t color, Option_t *opt)
double cb_dist_for_arb_X_selection(double *x, double *p)
double mean_X_vs_cb(double *x, double *par)
void DrawBDistForSelectionFromIPHisto(TH1 *sel, TH1 *incl, double &mean, double &sigma, Option_t *opt="", Color_t color=kRed, const TString &title="")
double P_X_cb_for_X_integral_with_selection(double *x, double *par)
void DrawMeanXvsCb(const TString &title="", Color_t color=-1, Option_t *opt="")
TGraph * GraphP_XForGivenB(double b, KVValueRange< double > Xrange, int npts=500)
void DrawFittedP_X(double norm=1.0, Option_t *opt="", Color_t color=kRed, const TString &title="")
Fluctuation kernel using gamma distribution for use with bayesian_estimator.
double operator()(double X, double mean, double reduced_variance)
Class implementing parametrizable impact parameter distributions.
void SetDeltaB_WithConstantCrossSection(Double_t deltab, Double_t sigmaR=0)
Range of values specified by minimum, maximum.
ValueType Max() const
ValueType ValueIofN(Int_t i, Int_t n) const
ValueType Min() const
virtual void SetLineColor(Color_t lcolor)
virtual Double_t GetBinLowEdge(Int_t bin) const
virtual Double_t GetBinUpEdge(Int_t bin) const
virtual Double_t Mean(Double_t a, Double_t b, const Double_t *params=nullptr, Double_t epsilon=0.000001)
virtual TH1 * GetHistogram() const
virtual void SetRange(Double_t xmin, Double_t xmax)
virtual Double_t Integral(Double_t a, Double_t b, Double_t epsrel=1.e-12)
void SetTitle(const char *title="") override
void Draw(Option_t *option="") override
virtual Double_t Variance(Double_t a, Double_t b, const Double_t *params=nullptr, Double_t epsilon=0.000001)
virtual TF1 * DrawCopy(Option_t *option="") const
virtual void SetParameters(const Double_t *params)
virtual void SetParNames(const char *name0="p0", const char *name1="p1", const char *name2="p2", const char *name3="p3", const char *name4="p4", const char *name5="p5", const char *name6="p6", const char *name7="p7", const char *name8="p8", const char *name9="p9", const char *name10="p10")
virtual Double_t Eval(Double_t x, Double_t y=0, Double_t z=0, Double_t t=0) const
virtual void SetParameter(const TString &name, Double_t value)
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
TAxis * GetXaxis()
TObject * FindObject(const char *name) const override
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
virtual Int_t GetNbinsX() const
virtual Double_t GetBinLowEdge(Int_t bin) const
virtual Double_t Integral(Int_t binx1, Int_t binx2, Option_t *option="") const
virtual Double_t GetBinContent(Int_t bin) const
virtual void Scale(Double_t c1=1, Option_t *option="")
virtual Int_t FindBin(Double_t x, Double_t y=0, Double_t z=0)
virtual void Warning(const char *method, const char *msgfmt,...) const
Bool_t IsNull() const
double binomial_pdf(unsigned int k, double p, unsigned int n)
double negative_binomial_pdf(unsigned int k, double p, double n)
double gamma_pdf(double x, double alpha, double theta, double x0=0)
const Double_t sigma
Double_t x[n]
TGraphErrors * gr
TH1 * h
fit(model, train_loader, val_loader, num_epochs, batch_size, optimizer, criterion, save_best, scheduler)
Int_t Nint(T x)
Double_t Sqrt(Double_t x)