4 #ifndef __bayesian_estimator_H
5 #define __bayesian_estimator_H
39 double operator()(
double X,
double mean,
double reduced_variance)
64 double operator()(
double X,
double mean,
double reduced_variance)
89 double operator()(
double X,
double mean,
double reduced_variance)
213 <
class FittingFunction,
class FluctuationKernel>
349 if (den > 0)
return num / den;
573 Warning(
"update_fit_params",
"no histogram set with FitHisto(TH1*)");
578 Warning(
"update_fit_params",
"no fit function found in histogram");
600 f->SetLineColor(color);
601 f->SetMarkerColor(color);
604 return f->GetMaximum();
628 int first_bin(0), last_bin(0);
629 for (
int i = 1; i <= incl->
GetNbinsX(); ++i) {
633 if (!first_bin) first_bin = i;
650 double cb_mean(0), cb_sqrmean(0), sum_pcb(0);
652 for (
int i = 0; i < 500; ++i) {
653 double cb = i / 499.;
656 cb_mean += p_cb * cb;
657 cb_sqrmean += p_cb * cb * cb;
660 f->SetPoint(i, cb, p_cb);
664 cb_sqrmean /= sum_pcb;
667 sigma_cb =
TMath::Sqrt(cb_sqrmean - cb_mean * cb_mean);
670 Info(
"DrawCbDistForSelection",
"Xmin=%lf, Xmax=%lf, cb_mean=%lf, sigma_cb=%lf", Xmin, Xmax, cb_mean, sigma_cb);
673 f->SetLineColor(color);
674 f->SetMarkerColor(color);
677 if (
TString(opt) ==
"same")
f->Draw(
"l");
702 for (
int i = 0; i < 500; ++i) {
705 if (sig > maxS) maxS = sig;
706 f->SetPoint(i,
b, sig);
708 f->SetLineColor(color);
709 f->SetMarkerColor(color);
712 if (
TString(opt) ==
"same")
f->Draw(
"l");
738 for (
int i = 0; i < npts; ++i) {
741 f->SetPoint(i,
b, sig);
772 int first_bin(0), last_bin(0);
773 for (
int i = 1; i <= incl->
GetNbinsX(); ++i) {
777 if (!first_bin) first_bin = i;
792 double bmean(0), bsqrmean(0), sigtot(0);
794 for (
int i = 0; i < 500; ++i) {
798 bsqrmean += sig *
b *
b;
800 f->SetPoint(i,
b, sig);
802 f->SetLineColor(color);
803 f->SetMarkerColor(color);
806 mean = bmean / sigtot;
809 if (
TString(opt) ==
"same")
f->Draw(
"l");
836 int first_bin(0), last_bin(0);
837 for (
int i = 1; i <= incl->
GetNbinsX(); ++i) {
841 if (!first_bin) first_bin = i;
891 f->SetLineColor(color);
958 b_range.
Min(), b_range.
Max(), X_range.
Min(), X_range.
Max(), 0);
#define ClassDef(name, id)
Base class for KaliVeda framework.
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="")
void RenormaliseHisto(TH1 *h)
double mean_X_vs_b(double *x, double *par)
TF1 p_X_X_integrator_with_selection
virtual ~bayesian_estimator()
void GetMeanAndSigmaBDistForXSelection(double X1, double X2, double &mean, double &sigma)
double P_X_cb(double X, double cb)
double DrawBDistForXSelection(KVValueRange< double > Xrange, Option_t *opt="", Color_t color=kRed, const TString &title="")
TF1 Cb_dist_for_arb_X_select
void GetMeanAndSigmaBDistForSelection(TH1 *sel, TH1 *incl, double &mean, double &sigma)
TF1 & GetB_dist_for_X_select()
bayesian_estimator(Bool_t integer_variable=false)
TGraph * GraphMeanXvsb(int npts=500)
double cb_dist_for_X_selection(double *x, double *p)
void SetIPDistParams(double sigmaR, double deltab)
FluctuationKernel theKernel
TF1 * GetFittedP_X(double norm=1.0)
double P_X_cb_for_TF2_obs_vs_b(double *x, double *)
double P_X_from_fit(double *x, double *par)
FittingFunction theFitter
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="")
void FitHisto(TH1 *h=nullptr)
bayesian_estimator(const FittingFunction &previous_fit, Bool_t integer_variable=false)
double b_dist_for_X_selection(double *x, double *p)
impact_parameter_distribution & GetIPDist()
double cb_integrated_P_X(double *x, double *p)
TF1 mean_X_vs_cb_function
std::vector< double > sel_rapp
double P_X_cb_for_X_integral(double *x, double *par)
void DrawNormalisedMeanXvsb(const TString &title, Color_t color, Option_t *opt)
void Print(Option_t *="") const
double cb_dist_for_arb_X_selection(double *x, double *p)
double mean_X_vs_cb(double *x, double *par)
TGraph * GraphBDistForXSelection(KVValueRange< double > Xrange, int npts=500)
double P_X_cb_for_X_integral_with_selection(double *x, double *par)
TGraph * GraphP_XForGivenB(double b, KVValueRange< double > Xrange, int npts=500)
TF2 * GetJointProbabilityDistribution(KVValueRange< double > b_range, KVValueRange< double > X_range)
void DrawMeanXvsCb(const TString &title="", Color_t color=-1, Option_t *opt="")
TF1 B_dist_for_arb_X_select
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.
Double_t GetDifferentialCrossSection(double b) const
const TF1 & GetCentrality()
void SetDeltaB_WithConstantCrossSection(Double_t deltab, Double_t sigmaR=0)
Range of values specified by minimum, maximum.
ValueType ValueIofN(Int_t i, Int_t n) const
virtual void SetLineColor(Color_t lcolor)
virtual Double_t GetBinLowEdge(Int_t bin) const
virtual Double_t GetBinUpEdge(Int_t bin) const
virtual TH1 * GetHistogram() const
virtual void SetTitle(const char *title="")
virtual void SetRange(Double_t xmin, Double_t xmax)
virtual Double_t Mean(Double_t a, Double_t b, const Double_t *params=0, Double_t epsilon=0.000001)
virtual Double_t Variance(Double_t a, Double_t b, const Double_t *params=0, Double_t epsilon=0.000001)
virtual Double_t Integral(Double_t a, Double_t b, Double_t epsrel=1.e-12)
virtual void Draw(Option_t *option="")
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)
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 TObject * FindObject(const char *name) 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
virtual void Info(const char *method, const char *msgfmt,...) 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)
def fit(model, train_loader, val_loader, num_epochs, batch_size, optimizer, criterion, save_best, scheduler)
Double_t Sqrt(Double_t x)