]>
Commit | Line | Data |
---|---|---|
ce4a88f5 | 1 | /* |
2 | * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. | |
3 | * See cxx source for full Copyright notice | |
4 | * $Id$ | |
5 | */ | |
6 | ||
7 | /******************************** | |
8 | * integrated flow estimate by * | |
9 | * fitting q-distribution * | |
10 | * * | |
11 | * author: Ante Bilandzic * | |
12 | * (anteb@nikhef.nl) * | |
13 | * * | |
14 | * based on the macro written * | |
15 | * by Sergei Voloshin * | |
16 | *******************************/ | |
17 | ||
18 | #ifndef ALIFLOWANALYSISWITHFITTINGQDISTRIBUTION_H | |
19 | #define ALIFLOWANALYSISWITHFITTINGQDISTRIBUTION_H | |
20 | ||
21 | #include "AliFlowCommonConstants.h" | |
22 | ||
23 | class TObjArray; | |
24 | class TList; | |
25 | class TFile; | |
929098e4 | 26 | class TDirectoryFile; |
ce4a88f5 | 27 | |
929098e4 | 28 | class TH1F; |
29 | class TH1D; | |
ee0860e8 | 30 | class TLegend; |
ce4a88f5 | 31 | class TProfile; |
32 | ||
33 | class AliFlowEventSimple; | |
34 | class AliFlowTrackSimple; | |
35 | class AliFlowCommonHist; | |
36 | class AliFlowCommonHistResults; | |
37 | class AliFlowVector; | |
38 | ||
39 | //================================================================================================================ | |
40 | ||
41 | class AliFlowAnalysisWithFittingQDistribution{ | |
42 | public: | |
43 | AliFlowAnalysisWithFittingQDistribution(); | |
44 | virtual ~AliFlowAnalysisWithFittingQDistribution(); | |
45 | // 0.) methods called in the constructor: | |
46 | virtual void InitializeArrays(); | |
47 | // 1.) method Init() and methods called within Init(): | |
48 | virtual void Init(); | |
49 | virtual void AccessConstants(); | |
50 | virtual void BookCommonHistograms(); | |
51 | virtual void BookAndFillWeightsHistograms(); | |
52 | virtual void BookEverythingForDistributions(); | |
ee0860e8 | 53 | virtual void StoreFittingParameters(); |
54 | virtual void AccessFittingParameters(); | |
ce4a88f5 | 55 | // 2.) method Make() and methods called within Make(): |
56 | virtual void Make(AliFlowEventSimple* anEvent); | |
57 | // 3.) method Finish() and methods called within Finish(): | |
58 | virtual void Finish(Bool_t doFit = kTRUE); | |
97ab16dd | 59 | virtual void DoFit(Bool_t useParticleWeights, Bool_t sigma2Fitted); |
60 | virtual void FillCommonHistResultsIntFlow(Bool_t useParticleWeights, Bool_t sigma2Fitted); | |
ce4a88f5 | 61 | virtual void PrintFinalResultsForIntegratedFlow(); |
62 | // 4.) other methods: | |
63 | virtual void GetOutputHistograms(TList *outputListHistos); | |
64 | virtual void WriteHistograms(TString *outputFileName); | |
65 | virtual void WriteHistograms(TString outputFileName); | |
ad87ae62 | 66 | virtual void WriteHistograms(TDirectoryFile *outputFileName); |
ce4a88f5 | 67 | |
68 | // **** SETTERS and GETTERS **** | |
69 | ||
70 | // 0.) base: | |
71 | TList* GetHistList() const {return this->fHistList;} | |
72 | // 1.) common: | |
73 | void SetCommonHists(AliFlowCommonHist* const ch) {this->fCommonHists = ch;}; | |
74 | AliFlowCommonHist* GetCommonHists() const {return this->fCommonHists;}; | |
75 | void SetCommonHistsResults(AliFlowCommonHistResults* const chr) {this->fCommonHistsResults = chr;}; | |
76 | AliFlowCommonHistResults* GetCommonHistsResults() const {return this->fCommonHistsResults;}; | |
77 | void SetHarmonic(Int_t const harmonic) {this->fHarmonic = harmonic;}; | |
78 | Int_t GetHarmonic() const {return this->fHarmonic;}; | |
79 | void SetAnalysisLabel(const char *aLabel) {this->fAnalysisLabel->Append(*aLabel);}; | |
80 | TString *GetAnalysisLabel() const {return this->fAnalysisLabel;}; | |
ce4a88f5 | 81 | // 2.) weights: |
82 | void SetWeightsList(TList* wlist) {this->fWeightsList = (TList*)wlist->Clone();}; | |
83 | TList* GetWeightsList() const {return this->fWeightsList;} | |
84 | void SetUsePhiWeights(Bool_t const uPhiW) {this->fUsePhiWeights = uPhiW;}; | |
85 | Bool_t GetUsePhiWeights() const {return this->fUsePhiWeights;}; | |
86 | void SetUsePtWeights(Bool_t const uPtW) {this->fUsePtWeights = uPtW;}; | |
87 | Bool_t GetUsePtWeights() const {return this->fUsePtWeights;}; | |
88 | void SetUseEtaWeights(Bool_t const uEtaW) {this->fUseEtaWeights = uEtaW;}; | |
89 | Bool_t GetUseEtaWeights() const {return this->fUseEtaWeights;}; | |
90 | void SetUseParticleWeights(TProfile* const uPW) {this->fUseParticleWeights = uPW;}; | |
91 | TProfile* GetUseParticleWeights() const {return this->fUseParticleWeights;}; | |
92 | void SetPhiWeights(TH1F* const histPhiWeights) {this->fPhiWeights = histPhiWeights;}; | |
93 | TH1F* GetPhiWeights() const {return this->fPhiWeights;}; | |
94 | void SetPtWeights(TH1D* const histPtWeights) {this->fPtWeights = histPtWeights;}; | |
95 | TH1D* GetPtWeights() const {return this->fPtWeights;}; | |
96 | void SetEtaWeights(TH1D* const histEtaWeights) {this->fEtaWeights = histEtaWeights;}; | |
97 | TH1D* GetEtaWeights() const {return this->fEtaWeights;}; | |
98 | // 3.) distributions: | |
99 | void SetSumOfParticleWeights(TH1D* const sopW, Int_t pW) {this->fSumOfParticleWeights[pW] = sopW;}; | |
100 | TH1D* GetSumOfParticleWeights(Int_t pW) const {return this->fSumOfParticleWeights[pW];}; | |
101 | void SetqDistribution(TH1D* const qd, Int_t pW) {this->fqDistribution[pW] = qd;}; | |
102 | TH1D* GetqDistribution(Int_t pW) const {return this->fqDistribution[pW];}; | |
103 | // 4.) final results of fitting: | |
97ab16dd | 104 | void SetIntFlow(TH1D* const intFlow, Int_t pW, Int_t sigmaFitted) {this->fIntFlow[pW][sigmaFitted] = intFlow;}; |
105 | TH1D* GetIntFlow(Int_t pW, Int_t sigmaFitted) const {return this->fIntFlow[pW][sigmaFitted];}; | |
106 | void SetSigma2(TH1D* const sigma2, Int_t pW, Int_t sigmaFitted) {this->fSigma2[pW][sigmaFitted] = sigma2;}; | |
107 | TH1D* GetSigma2(Int_t pW, Int_t sigmaFitted) const {return this->fSigma2[pW][sigmaFitted];}; | |
108 | void SetChi2(TH1D* const chi2, Int_t pW, Int_t sigmaFitted) {this->fChi2[pW][sigmaFitted] = chi2;}; | |
109 | TH1D* GetChi2(Int_t pW, Int_t sigmaFitted) const {return this->fChi2[pW][sigmaFitted];}; | |
ee0860e8 | 110 | // 5.) fitting parameters: |
111 | void SetFittingParameters(TProfile* const fp) {this->fFittingParameters = fp;}; | |
112 | TProfile* GetFittingParameters() const {return this->fFittingParameters;}; | |
113 | void SetTreshold(Double_t const treshold) {this->fTreshold = treshold;}; | |
114 | Double_t GetTreshold() const {return this->fTreshold;}; | |
115 | void SetvStart(Double_t const vStart) {this->fvStart = vStart;}; | |
116 | Double_t GetvStart() const {return this->fvStart;}; | |
117 | void SetvMin(Double_t const vMin) {this->fvMin = vMin;}; | |
118 | Double_t GetvMin() const {return this->fvMin;}; | |
119 | void SetvMax(Double_t const vMax) {this->fvMax = vMax;}; | |
120 | Double_t GetvMax() const {return this->fvMax;}; | |
121 | void SetSigma2Start(Double_t const Sigma2Start) {this->fSigma2Start = Sigma2Start;}; | |
122 | Double_t GetSigma2Start() const {return this->fSigma2Start;}; | |
123 | void SetSigma2Min(Double_t const Sigma2Min) {this->fSigma2Min = Sigma2Min;}; | |
124 | Double_t GetSigma2Min() const {return this->fSigma2Min;}; | |
125 | void SetSigma2Max(Double_t const Sigma2Max) {this->fSigma2Max = Sigma2Max;}; | |
126 | Double_t GetSigma2Max() const {return this->fSigma2Max;}; | |
127 | void SetPlotResults(Bool_t const pr) {this->fPlotResults = pr;}; | |
128 | Double_t GetPlotResults() const {return this->fPlotResults;}; | |
ce4a88f5 | 129 | |
130 | private: | |
131 | AliFlowAnalysisWithFittingQDistribution(const AliFlowAnalysisWithFittingQDistribution &afawfqd); | |
47426051 | 132 | AliFlowAnalysisWithFittingQDistribution& operator=(const AliFlowAnalysisWithFittingQDistribution &afawfqd); |
ce4a88f5 | 133 | // 0.) base: |
134 | TList *fHistList; // base list to hold all output object | |
135 | // 1.) common: | |
136 | AliFlowCommonHist *fCommonHists; // common control histograms | |
137 | AliFlowCommonHistResults *fCommonHistsResults; // final results in common histograms | |
138 | Int_t fnBinsPhi; // number of phi bins | |
139 | Double_t fPhiMin; // minimum phi | |
140 | Double_t fPhiMax; // maximum phi | |
141 | Double_t fPhiBinWidth; // bin width for phi histograms | |
142 | Int_t fnBinsPt; // number of pt bins | |
143 | Double_t fPtMin; // minimum pt | |
144 | Double_t fPtMax; // maximum pt | |
145 | Double_t fPtBinWidth; // bin width for pt histograms | |
146 | Int_t fnBinsEta; // number of eta bins | |
147 | Double_t fEtaMin; // minimum eta | |
148 | Double_t fEtaMax; // maximum eta | |
149 | Double_t fEtaBinWidth; // bin width for eta histograms | |
150 | Int_t fHarmonic; // harmonic | |
151 | TString *fAnalysisLabel; // analysis label (all histograms and output file will have this label) | |
152 | // 2.) particle weights (abbreviated to 'pWeights' or even to 'pW' throughout the code): | |
153 | TList *fWeightsList; // list to hold all histograms with particle weights: fUseParticleWeights, fPhiWeights, fPtWeights and fEtaWeights | |
154 | Bool_t fUsePhiWeights; // use phi weights | |
155 | Bool_t fUsePtWeights; // use pt weights | |
156 | Bool_t fUseEtaWeights; // use eta weights | |
157 | TProfile *fUseParticleWeights; // profile with three bins to hold values of fUsePhiWeights, fUsePtWeights and fUseEtaWeights | |
158 | TH1F *fPhiWeights; // histogram holding phi weights | |
47426051 | 159 | TH1D *fPtWeights; // histogram holding pt weights |
160 | TH1D *fEtaWeights; // histogram holding eta weights | |
ce4a88f5 | 161 | // 3.) distributions: |
162 | TH1D *fSumOfParticleWeights[2]; // [0=particle weights are unit (not used), 1=particle weights are used] | |
163 | TH1D *fqDistribution[2]; // distribution of Q/sqrt{sum of particle weights} [0=particle weights are unit (not used), 1=particle weights are used] | |
164 | // 4.) final results of fitting: | |
97ab16dd | 165 | TH1D *fIntFlow[2][2]; // final result for integrated flow [0=pWeights are unit (not used), 1=pWeights are used][0=sigma^2 not fitted, 1=sigma^2 fitted] |
166 | TH1D *fSigma2[2][2]; // final results for sigma^2 [0=pWeights are unit (not used), 1=pWeights are used][0=sigma^2 not fitted, 1=sigma^2 fitted] | |
167 | TH1D *fChi2[2][2]; // final results for chi^2 from Minuit [0=pWeights are unit (not used), 1=pWeights are used][0=sigma^2 not fitted, 1=sigma^2 fitted] | |
ee0860e8 | 168 | // 5.) fitting parameters: |
169 | TProfile *fFittingParameters; // profile to hold all fitting parameters | |
47426051 | 170 | Double_t fTreshold; // the first bin taken for the fitting is the first bin with nEntries >= fTreshold (analogously for the last bin) |
ee0860e8 | 171 | Double_t fvStart; // fitting of v will start from this point |
172 | Double_t fvMin; // v range, lower boundary | |
173 | Double_t fvMax; // v range, upper boundary | |
174 | Double_t fSigma2Start; // fitting of sigma2 will start from this point | |
47426051 | 175 | Double_t fSigma2Min; // sigma2 range, lower boundary (this should be kept above 0.5 according to theorists...) |
ee0860e8 | 176 | Double_t fSigma2Max; // sigma2 range, upper boundary |
47426051 | 177 | Bool_t fPlotResults; // plot or not q-distribution and resulting fitting function |
ee0860e8 | 178 | // 6.) rest: |
97ab16dd | 179 | TLegend *fLegend; // legend // to be improved (do I need this as data member?) |
ee0860e8 | 180 | |
ce4a88f5 | 181 | ClassDef(AliFlowAnalysisWithFittingQDistribution, 0); |
182 | }; | |
183 | ||
184 | //================================================================================================================ | |
185 | ||
186 | #endif | |
187 | ||
188 | ||
189 | ||
190 | ||
191 |