Centrality Framework
Loading...
Searching...
No Matches
Fitter.hpp
1
8#ifndef GlauberFitter_H
9#define GlauberFitter_H 1
10
11#include "TH1F.h"
12#include "TNamed.h"
13#include "TString.h"
14#include "TTree.h"
15
16#include <vector>
17// #include "TMinuit.h"
18
19namespace Glauber {
20class Fitter {
21
22 public:
24 Fitter() = default;
25 ;
26 explicit Fitter(std::unique_ptr<TTree> tree);
28 virtual ~Fitter() = default;
29 ;
30
31 void Init(int nEntries);
32 void SetGlauberFitHisto(double f, double mu, double k, Int_t n = 10000, Bool_t Norm2Data = true);
33 void NormalizeGlauberFit();
34 // void DrawHistos(Bool_t isSim = true, Bool_t isData = true, Bool_t isGlauber = false, Bool_t isNBD = false) {};
35
36 double FitGlauber(double* par, double f0, Int_t k0, Int_t k1, Int_t nEvents);
37 void FindMuGoldenSection(double* mu,
38 double* chi2,
39 double mu_min,
40 double mu_max,
41 double f,
42 double k,
43 Int_t nEvents = 10000,
44 Int_t nIter = 5);
45
46 double GetChi2() const;
47
48 double NBD(double n, double mu, double k) const;
49 void SetNBDhist(double mu, double k);
50 double Nancestors(double f) const;
51 double NancestorsMax(double f) const;
52
53 std::unique_ptr<TH1F> GetModelHisto(const double range[2], const TString& name, const double par[3], Int_t nEvents);
54
55 //
56 // Setters
57 //
58 void SetInputHisto(const TH1F& h) { fDataHisto = h; }
59 void SetFitMinBin(Int_t min) { fFitMinBin = min; }
60 void SetFitMaxBin(Int_t min) { fFitMaxBin = min; }
61 void SetNormMinBin(Int_t min) { fNormMinBin = min; }
62 void SetBinSize(Int_t size) { fBinSize = size; }
63 void SetOutDirName(const TString& name) { fOutDirName = name; }
64 void SetMode(const TString& mode) { fMode = mode; }
65
66 //
67 // Getters
68 //
69 TH1F GetGlauberFitHisto() const { return fGlauberFitHisto; }
70 TH1F GetDataHisto() const { return fDataHisto; }
71 TH1F GetNBDHisto() const { return fNbdHisto; }
72 TH1F GetNpartHisto() const { return fNpartHisto; }
73 TH1F GetNcollHisto() const { return fNcollHisto; }
74 TH1F GetBestFiHisto() const { return fBestFitHisto; }
75
76 private:
78 TH1F fNpartHisto;
79 TH1F fNcollHisto;
80 TH1F fDataHisto;
81 TH1F fNbdHisto;
82 TH1F fGlauberFitHisto;
83 TH1F fBestFitHisto;
84
85 /* MC data */
86 std::unique_ptr<TTree> fSimTree{nullptr};
87
88 float fNpart{-1.};
89 float fNcoll{-1.};
90
91 double fMaxValue{-1.};
92
93 Int_t fNbins{-1};
94 Int_t fBinSize{-1};
95
96 Int_t fFitMinBin{-1};
97 Int_t fFitMaxBin{-1};
98
99 Int_t fNormMinBin{-1};
100
101 TString fMode{"Default"};
102
103 TString fOutDirName{""};
104 ClassDef(Fitter, 2);
105};
106}// namespace Glauber
107
108#endif
Class to fit histo with Glauber based function.
Definition Fitter.hpp:20
virtual ~Fitter()=default
void FindMuGoldenSection(double *mu, double *chi2, double mu_min, double mu_max, double f, double k, Int_t nEvents=10000, Int_t nIter=5)
Definition Fitter.cpp:148
std::unique_ptr< TH1F > GetModelHisto(const double range[2], const TString &name, const double par[3], Int_t nEvents)
Definition Fitter.cpp:339
double GetChi2() const
Definition Fitter.cpp:267
void SetNBDhist(double mu, double k)
Definition Fitter.cpp:290
double FitGlauber(double *par, double f0, Int_t k0, Int_t k1, Int_t nEvents)
Definition Fitter.cpp:206
double NBD(double n, double mu, double k) const
Definition Fitter.cpp:311
Fitter()=default