2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved.
3 * See cxx source for full Copyright notice
7 /************************************
8 * flow analysis with multi-particle *
11 * author: Ante Bilandzic *
12 * (abilandzic@gmail.com) *
13 ************************************/
15 #ifndef ALIFLOWANALYSISWITHMULTIPARTICLECORRELATIONS_H
16 #define ALIFLOWANALYSISWITHMULTIPARTICLECORRELATIONS_H
21 #include "TProfile2D.h"
24 #include "TDirectoryFile.h"
25 #include "Riostream.h"
28 #include "TGraphErrors.h"
29 #include "TStopwatch.h"
30 #include "AliFlowEventSimple.h"
31 #include "AliFlowTrackSimple.h"
33 class AliFlowAnalysisWithMultiparticleCorrelations{
35 AliFlowAnalysisWithMultiparticleCorrelations();
36 virtual ~AliFlowAnalysisWithMultiparticleCorrelations();
37 // Member functions are grouped as:
38 // 0.) Methods called in the constructor;
39 // 1.) Method Init() and methods called in it (!);
40 // 2.) Method Make() and methods called in it;
41 // 3.) Method Finish() and methods called in it;
42 // 4.) Method GetOutputHistograms() and methods called in it;
43 // 5.) Setters and getters;
46 // 0.) Methods called in the constructor:
47 virtual void InitializeArraysForControlHistograms();
48 virtual void InitializeArraysForQvector();
49 virtual void InitializeArraysForCorrelations();
50 virtual void InitializeArraysForEbECumulants();
51 virtual void InitializeArraysForWeights();
52 virtual void InitializeArraysForQcumulants();
53 virtual void InitializeArraysForDiffCorrelations();
54 virtual void InitializeArraysForNestedLoops();
56 // 1.) Method Init() and methods called in it (!):
58 virtual void CrossCheckSettings();
59 virtual void BookAndNestAllLists();
60 virtual void BookEverythingForBase();
61 virtual void BookEverythingForControlHistograms();
62 virtual void BookEverythingForQvector();
63 virtual void BookEverythingForWeights();
64 virtual void BookEverythingForCorrelations();
65 virtual void BookEverythingForEbECumulants();
66 virtual void BookEverythingForNestedLoops();
67 virtual void BookEverythingForStandardCandles();
68 virtual void BookEverythingForQcumulants();
69 virtual void BookEverythingForDiffCorrelations();
71 // 2.) Method Make() and methods called in it:
72 virtual void Make(AliFlowEventSimple *anEvent);
73 virtual Bool_t CrossCheckInternalFlags(AliFlowEventSimple *anEvent);
74 virtual void CrossCheckPointersUsedInMake();
75 virtual void FillControlHistograms(AliFlowEventSimple *anEvent);
76 virtual void FillQvector(AliFlowEventSimple *anEvent);
77 virtual void CalculateCorrelations(AliFlowEventSimple *anEvent);
78 virtual void CalculateDiffCorrelations(AliFlowEventSimple *anEvent);
79 virtual void CalculateEbECumulants(AliFlowEventSimple *anEvent);
80 virtual void ResetQvector();
81 virtual void CrossCheckWithNestedLoops(AliFlowEventSimple *anEvent);
82 virtual void CrossCheckDiffWithNestedLoops(AliFlowEventSimple *anEvent);
84 // 3.) Method Finish() and methods called in it:
85 virtual void Finish();
86 virtual void CrossCheckPointersUsedInFinish();
87 virtual void CalculateStandardCandles();
88 virtual void CalculateQcumulants();
89 virtual void CalculateReferenceFlow();
91 // 4.) Method GetOutputHistograms() and methods called in it:
92 virtual void GetOutputHistograms(TList *histList);
93 virtual void GetPointersForControlHistograms();
94 virtual void GetPointersForQvector();
95 virtual void GetPointersForCorrelations();
96 virtual void GetPointersForStandardCandles();
97 virtual void GetPointersForQcumulants();
98 virtual void GetPointersForDiffCorrelations();
100 // 5.) Setters and getters:
101 // 5.0.) Base list and internal flags:
102 void SetHistList(TList* const hlist) {this->fHistList = hlist;}
103 TList* GetHistList() const {return this->fHistList;}
104 void SetInternalFlagsPro(TProfile* const ifp) {this->fInternalFlagsPro = ifp;};
105 TProfile* GetInternalFlagsPro() const {return this->fInternalFlagsPro;};
106 void SetMinNoRPs(Int_t const min) {fUseInternalFlags = kTRUE; this->fMinNoRPs = min;};
107 Int_t GetMinNoRPs() const {return this->fMinNoRPs;};
108 void SetMaxNoRPs(Int_t const max) {fUseInternalFlags = kTRUE; this->fMaxNoRPs = max;};
109 Int_t GetMaxNoRPs() const {return this->fMaxNoRPs;};
110 void SetExactNoRPs(Int_t const exact) {fUseInternalFlags = kTRUE; this->fExactNoRPs = exact;};
111 Int_t GetExactNoRPs() const {return this->fExactNoRPs;};
113 // 5.1.) Control histograms:
114 void SetControlHistogramsList(TList* const chl) {this->fControlHistogramsList = chl;};
115 TList* GetControlHistogramsList() const {return this->fControlHistogramsList;}
116 void SetControlHistogramsFlagsPro(TProfile* const chfp) {this->fControlHistogramsFlagsPro = chfp;};
117 TProfile* GetControlHistogramsFlagsPro() const {return this->fControlHistogramsFlagsPro;};
118 void SetFillControlHistograms(Bool_t const fch) {this->fFillControlHistograms = fch;};
119 Bool_t GetFillControlHistograms() const {return this->fFillControlHistograms;};
120 void SetFillKinematicsHist(Bool_t const fkh) {this->fFillKinematicsHist = fkh;};
121 Bool_t GetFillKinematicsHist() const {return this->fFillKinematicsHist;};
122 void SetFillMultDistributionsHist(Bool_t const mdh) {this->fFillMultDistributionsHist = mdh;};
123 Bool_t GetFillMultDistributionsHist() const {return this->fFillMultDistributionsHist;};
124 void SetFillMultCorrelationsHist(Bool_t const mch) {this->fFillMultCorrelationsHist = mch;};
125 Bool_t GetFillMultCorrelationsHist() const {return this->fFillMultCorrelationsHist;};
126 void SetnBins(const char *type, const char *variable, const Int_t nBins); // .cxx
127 void SetMin(const char *type, const char *variable, const Double_t min); // .cxx
128 void SetMax(const char *type, const char *variable, const Double_t max); // .cxx
129 void SetnBinsMult(const char *type, const Int_t nBinsMult); // .cxx
130 void SetMinMult(const char *type, const Double_t minMult); // .cxx
131 void SetMaxMult(const char *type, const Double_t maxMult); // .cxx
134 void SetQvectorList(TList* const qvl) {this->fQvectorList = qvl;};
135 TList* GetQvectorList() const {return this->fQvectorList;}
136 void SetQvectorFlagsPro(TProfile* const qvfp) {this->fQvectorFlagsPro = qvfp;};
137 TProfile* GetQvectorFlagsPro() const {return this->fQvectorFlagsPro;};
138 void SetCalculateQvector(Bool_t const cqv) {this->fCalculateQvector = cqv;};
139 Bool_t GetCalculateQvector() const {return this->fCalculateQvector;};
140 void SetCalculateDiffQvectors(Bool_t const cdqv) {this->fCalculateDiffQvectors = cdqv;};
141 Bool_t GetCalculateDiffQvectors() const {return this->fCalculateDiffQvectors;};
143 // 5.3.) Correlations:
144 void SetCorrelationsList(TList* const cl) {this->fCorrelationsList = cl;};
145 TList* GetCorrelationsList() const {return this->fCorrelationsList;}
146 void SetCorrelationsFlagsPro(TProfile* const cfp) {this->fCorrelationsFlagsPro = cfp;};
147 TProfile* GetCorrelationsFlagsPro() const {return this->fCorrelationsFlagsPro;};
148 void SetCalculateCorrelations(Bool_t const cc) {this->fCalculateCorrelations = cc;};
149 Bool_t GetCalculateCorrelations() const {return this->fCalculateCorrelations;};
150 void SetCalculateIsotropic(Bool_t const ci) {this->fCalculateIsotropic = ci;};
151 Bool_t GetCalculateIsotropic() const {return this->fCalculateIsotropic;};
152 void SetCalculateSame(Bool_t const cs) {this->fCalculateSame = cs;};
153 Bool_t GetCalculateSame() const {return this->fCalculateSame;};
154 void SetSkipZeroHarmonics(Bool_t const szh) {this->fSkipZeroHarmonics = szh;};
155 Bool_t GetSkipZeroHarmonics() const {return this->fSkipZeroHarmonics;};
156 void SetCalculateSameIsotropic(Bool_t const csi) {this->fCalculateSameIsotropic = csi;};
157 Bool_t GetCalculateSameIsotropic() const {return this->fCalculateSameIsotropic;};
158 void SetCalculateAll(Bool_t const ca) {this->fCalculateAll = ca;};
159 Bool_t GetCalculateAll() const {return this->fCalculateAll;};
160 void SetDontGoBeyond(Int_t const dgb) {this->fDontGoBeyond = dgb;};
161 Int_t GetDontGoBeyond() const {return this->fDontGoBeyond;};
162 void SetCalculateOnlyForHarmonicQC(Bool_t const cofhqc) {this->fCalculateOnlyForHarmonicQC = cofhqc;};
163 Bool_t GetCalculateOnlyForHarmonicQC() const {return this->fCalculateOnlyForHarmonicQC;};
164 void SetCalculateOnlyForSC(Bool_t const cofsc) {this->fCalculateOnlyForSC = cofsc;};
165 Bool_t GetCalculateOnlyForSC() const {return this->fCalculateOnlyForSC;};
166 void SetCalculateOnlyCos(Bool_t const coc) {this->fCalculateOnlyCos = coc;};
167 Bool_t GetCalculateOnlyCos() const {return this->fCalculateOnlyCos;};
168 void SetCalculateOnlySin(Bool_t const cos) {this->fCalculateOnlySin = cos;};
169 Bool_t GetCalculateOnlySin() const {return this->fCalculateOnlySin;};
171 // 5.4.) Event-by-event cumulants:
172 void SetEbECumulantsList(TList* const ebecl) {this->fEbECumulantsList = ebecl;};
173 TList* GetEbECumulantsList() const {return this->fEbECumulantsList;}
174 void SetEbECumulantsFlagsPro(TProfile* const ebecfp) {this->fEbECumulantsFlagsPro = ebecfp;};
175 TProfile* GetEbECumulantsFlagsPro() const {return this->fEbECumulantsFlagsPro;};
176 void SetCalculateEbECumulants(Bool_t const cebec) {this->fCalculateEbECumulants = cebec;};
177 Bool_t GetCalculateEbECumulants() const {return this->fCalculateEbECumulants;};
180 void SetWeightsList(TList* const wl) {this->fWeightsList = (TList*)wl->Clone();};
181 TList* GetWeightsList() const {return this->fWeightsList;}
182 void SetWeightsFlagsPro(TProfile* const wfp) {this->fWeightsFlagsPro = wfp;};
183 TProfile* GetWeightsFlagsPro() const {return this->fWeightsFlagsPro;};
184 void SetWeightsHist(TH1D* const hist, const char *type, const char *variable); // .cxx
186 // 5.6.) Nested loops:
187 void SetNestedLoopsList(TList* const nll) {this->fNestedLoopsList = nll;}
188 TList* GetNestedLoopsList() const {return this->fNestedLoopsList;}
189 void SetNestedLoopsFlagsPro(TProfile* const nlfp) {this->fNestedLoopsFlagsPro = nlfp;};
190 TProfile* GetNestedLoopsFlagsPro() const {return this->fNestedLoopsFlagsPro;};
191 void SetCrossCheckWithNestedLoops(Bool_t const ccwnl) {this->fCrossCheckWithNestedLoops = ccwnl;};
192 Bool_t GetCrossCheckWithNestedLoops() const {return this->fCrossCheckWithNestedLoops;};
193 void SetCrossCheckDiffWithNestedLoops(Bool_t const ccdwnl) {this->fCrossCheckDiffWithNestedLoops = ccdwnl;};
194 Bool_t GetCrossCheckDiffWithNestedLoops() const {return this->fCrossCheckDiffWithNestedLoops;};
195 void SetCrossCheckDiffCSCOBN(Int_t const cs, Int_t const co, Int_t const bn)
197 this->fCrossCheckDiffCSCOBN[0] = cs; // cos/sin
198 this->fCrossCheckDiffCSCOBN[1] = co; // correlator order [1p,2p,3p,4p]
199 this->fCrossCheckDiffCSCOBN[2] = bn; // bin number
202 // 5.7.) 'Standard candles':
203 void SetStandardCandlesList(TList* const scl) {this->fStandardCandlesList = scl;}
204 TList* GetStandardCandlesList() const {return this->fStandardCandlesList;}
205 void SetStandardCandlesFlagsPro(TProfile* const scfp) {this->fStandardCandlesFlagsPro = scfp;};
206 TProfile* GetStandardCandlesFlagsPro() const {return this->fStandardCandlesFlagsPro;};
207 void SetCalculateStandardCandles(Bool_t const csc) {this->fCalculateStandardCandles = csc;};
208 Bool_t GetCalculateStandardCandles() const {return this->fCalculateStandardCandles;};
209 void SetPropagateErrorSC(Bool_t const pesc) {this->fPropagateErrorSC = pesc;};
210 Bool_t GetPropagateErrorSC() const {return this->fPropagateErrorSC;};
211 void SetStandardCandlesHist(TH1D* const sch) {this->fStandardCandlesHist = sch;};
212 TH1D* GetStandardCandlesHist() const {return this->fStandardCandlesHist;};
213 void SetProductsSCPro(TProfile2D* const psc) {this->fProductsSCPro = psc;};
214 TProfile2D* GetProductsSCPro() const {return this->fProductsSCPro;};
216 // 5.8.) Q-cumulants:
217 void SetQcumulantsList(TList* const qcl) {this->fQcumulantsList = qcl;};
218 TList* GetQcumulantsList() const {return this->fQcumulantsList;}
219 void SetQcumulantsFlagsPro(TProfile* const qcfp) {this->fQcumulantsFlagsPro = qcfp;};
220 TProfile* GetQcumulantsFlagsPro() const {return this->fQcumulantsFlagsPro;};
221 void SetCalculateQcumulants(Bool_t const cqc) {this->fCalculateQcumulants = cqc;};
222 Bool_t GetCalculateQcumulants() const {return this->fCalculateQcumulants;};
223 void SetHarmonicQC(Int_t const hqc) {this->fHarmonicQC = hqc;};
224 Int_t GetHarmonicQC() const {return this->fHarmonicQC;};
225 void SetPropagateErrorQC(Bool_t const peqc) {this->fPropagateErrorQC = peqc;};
226 Bool_t GetPropagateErrorQC() const {return this->fPropagateErrorQC;};
227 void SetQcumulantsHist(TH1D* const qch) {this->fQcumulantsHist = qch;};
228 TH1D* GetQcumulantsHist() const {return this->fQcumulantsHist;};
229 void SetReferenceFlowHist(TH1D* const rfh) {this->fReferenceFlowHist = rfh;};
230 TH1D* GetReferenceFlowHist() const {return this->fReferenceFlowHist;};
231 void SetProductsQCPro(TProfile2D* const pqc) {this->fProductsQCPro = pqc;};
232 TProfile2D* GetProductsQCPro() const {return this->fProductsQCPro;};
234 // 5.9.) Differential correlations:
235 void SetDiffCorrelationsList(TList* const dcl) {this->fDiffCorrelationsList = dcl;};
236 TList* GetDiffCorrelationsList() const {return this->fDiffCorrelationsList;}
237 void SetDiffCorrelationsFlagsPro(TProfile* const cdfp) {this->fDiffCorrelationsFlagsPro = cdfp;};
238 TProfile* GetDiffCorrelationsFlagsPro() const {return this->fDiffCorrelationsFlagsPro;};
239 void SetCalculateDiffCorrelations(Bool_t const cdc) {this->fCalculateDiffCorrelations = cdc;};
240 Bool_t GetCalculateDiffCorrelations() const {return this->fCalculateDiffCorrelations;};
241 void SetDiffHarmonics(Int_t order, Int_t *harmonics); // see implementation in .cxx file
242 void SetCalculateDiffCos(Bool_t const cdc) {this->fCalculateDiffCos = cdc;};
243 Bool_t GetCalculateDiffCos() const {return this->fCalculateDiffCos;};
244 void SetCalculateDiffSin(Bool_t const cds) {this->fCalculateDiffSin = cds;};
245 Bool_t GetCalculateDiffSin() const {return this->fCalculateDiffSin;};
246 void SetCalculateDiffCorrelationsVsPt(Bool_t const cdcvspt) {this->fCalculateDiffCorrelationsVsPt = cdcvspt;};
247 Bool_t GetCalculateDiffCorrelationsVsPt() const {return this->fCalculateDiffCorrelationsVsPt;};
248 void SetUseDefaultBinning(Bool_t const udb) {this->fUseDefaultBinning = udb;};
249 Bool_t GetUseDefaultBinning() const {return this->fUseDefaultBinning;};
250 void SetnDiffBins(Int_t const ndb) {this->fnDiffBins = ndb;};
251 Int_t GetnDiffBins() const {return this->fnDiffBins;};
252 void SetRangesDiffBins(Double_t* const rdb) {this->fRangesDiffBins = rdb;};
253 Double_t* GetRangesDiffBins() const {return this->fRangesDiffBins;};
256 virtual void WriteHistograms(TString outputFileName);
257 virtual void WriteHistograms(TDirectoryFile *outputFileName);
258 virtual TComplex Q(Int_t n, Int_t p);
259 virtual TComplex p(Int_t n, Int_t p);
260 virtual TComplex q(Int_t n, Int_t p);
261 virtual TComplex One(Int_t n1);
262 virtual TComplex Two(Int_t n1, Int_t n2);
263 virtual TComplex Three(Int_t n1, Int_t n2, Int_t n3);
264 virtual TComplex Four(Int_t n1, Int_t n2, Int_t n3, Int_t n4);
265 virtual TComplex Five(Int_t n1, Int_t n2, Int_t n3, Int_t n4, Int_t n5);
266 virtual TComplex Six(Int_t n1, Int_t n2, Int_t n3, Int_t n4, Int_t n5, Int_t n6);
267 virtual TComplex Seven(Int_t n1, Int_t n2, Int_t n3, Int_t n4, Int_t n5, Int_t n6, Int_t n7);
268 virtual TComplex Eight(Int_t n1, Int_t n2, Int_t n3, Int_t n4, Int_t n5, Int_t n6, Int_t n7, Int_t n8);
269 virtual TComplex OneDiff(Int_t n1);
270 virtual TComplex TwoDiff(Int_t n1, Int_t n2);
271 virtual TComplex ThreeDiff(Int_t n1, Int_t n2, Int_t n3);
272 virtual TComplex FourDiff(Int_t n1, Int_t n2, Int_t n3, Int_t n4);
273 virtual Double_t Weight(const Double_t &value, const char *type, const char *variable); // value, [RP,POI], [phi,pt,eta]
274 virtual Double_t CastStringToCorrelation(const char *string, Bool_t numerator);
275 virtual Double_t Covariance(const char *x, const char *y, TProfile2D *profile2D, Bool_t bUnbiasedEstimator = kFALSE);
276 virtual TComplex Recursion(Int_t n, Int_t* harmonic, Int_t mult = 1, Int_t skip = 0); // Credits: Kristjan Gulbrandsen (gulbrand@nbi.dk)
277 virtual void CalculateProductsOfCorrelations(AliFlowEventSimple *anEvent, TProfile2D *profile2D);
278 static void DumpPointsForDurham(TGraphErrors *ge);
279 static void DumpPointsForDurham(TH1D *h);
280 static void DumpPointsForDurham(TH1F *h);
283 AliFlowAnalysisWithMultiparticleCorrelations(const AliFlowAnalysisWithMultiparticleCorrelations& afawQc);
284 AliFlowAnalysisWithMultiparticleCorrelations& operator=(const AliFlowAnalysisWithMultiparticleCorrelations& afawQc);
285 // Data members are grouped as:
286 // 0.) Base list and internal flags;
287 // 1.) Control histograms;
290 // 4.) Event-by-event cumulants;
293 // 7.) 'Standard candles';
295 // 9.) Differential correlations.
297 // 0.) Base list and internal flags:
298 TList* fHistList; // base list to hold all output object (a.k.a. grandmother of all lists)
299 TProfile *fInternalFlagsPro; // profile to hold all internal flags and settings
300 Bool_t fUseInternalFlags; // use internal flags (automatically set if some internal flag is used)
301 Int_t fMinNoRPs; // minimum number of RPs required for the analysis
302 Int_t fMaxNoRPs; // maximum number of RPs allowed for the analysis
303 Int_t fExactNoRPs; // exact (randomly shuffled) number of RPs selected for the analysis
304 Bool_t fPropagateError; // prevent error propagation if something strange happens during calculations
306 // 1.) Control histograms:
307 TList *fControlHistogramsList; // list to hold all 'control histograms' objects
308 TProfile *fControlHistogramsFlagsPro; // profile to hold all flags for control histograms
309 Bool_t fFillControlHistograms; // fill or not control histograms (by default they are all filled)
310 Bool_t fFillKinematicsHist; // fill or not fKinematicsHist[2][3]
311 Bool_t fFillMultDistributionsHist; // fill or not TH1D *fMultDistributionsHist[3]
312 Bool_t fFillMultCorrelationsHist; // fill or not TH2D *fMultCorrelationsHist[3]
313 TH1D *fKinematicsHist[2][3]; // [RP,POI][phi,pt,eta] distributions
314 TH1D *fMultDistributionsHist[3]; // multiplicity distribution [RP,POI,reference multiplicity]
315 TH2D *fMultCorrelationsHist[3]; // [RP vs. POI, RP vs. refMult, POI vs. refMult]
316 Int_t fnBins[2][3]; // [RP,POI][phi,pt,eta], corresponds to fKinematicsHist[2][3]
317 Double_t fMin[2][3]; // [RP,POI][phi,pt,eta], corresponds to fKinematicsHist[2][3]
318 Double_t fMax[2][3]; // [RP,POI][phi,pt,eta], corresponds to fKinematicsHist[2][3]
319 Int_t fnBinsMult[3]; // [RP,POI,REF], corresponds to fMultDistributionsHist[3]
320 Double_t fMinMult[3]; // [RP,POI,REF], corresponds to fMultDistributionsHist[3]
321 Double_t fMaxMult[3]; // [RP,POI,REF], corresponds to fMultDistributionsHist[3]
324 TList *fQvectorList; // list to hold all Q-vector objects
325 TProfile *fQvectorFlagsPro; // profile to hold all flags for Q-vector
326 Bool_t fCalculateQvector; // to calculate or not to calculate Q-vector components, that's a Boolean...
327 TComplex fQvector[49][9]; // Q-vector components [fMaxHarmonic*fMaxCorrelator+1][fMaxCorrelator+1] = [6*8+1][8+1]
328 Bool_t fCalculateDiffQvectors; // to calculate or not to calculate p- and q-vector components, that's a Boolean...
329 TComplex fpvector[100][49][9]; // p-vector components [bin][fMaxHarmonic*fMaxCorrelator+1][fMaxCorrelator+1] = [6*8+1][8+1] TBI hardwired 100
330 TComplex fqvector[100][49][9]; // q-vector components [bin][fMaxHarmonic*fMaxCorrelator+1][fMaxCorrelator+1] = [6*8+1][8+1] TBI hardwired 100
333 TList *fCorrelationsList; // list to hold all correlations objects
334 TProfile *fCorrelationsFlagsPro; // profile to hold all flags for correlations
335 TProfile *fCorrelationsPro[2][8]; // multi-particle correlations [0=cos,1=sin][1p,2p,...,8p]
336 Bool_t fCalculateCorrelations; // calculate and store correlations
337 Int_t fMaxHarmonic; // 6 (not going beyond v6, if you change this value, change also fQvector[49][9])
338 Int_t fMaxCorrelator; // 8 (not going beyond 8-p correlations, if you change this value, change also fQvector[49][9])
339 Bool_t fCalculateIsotropic; // calculate only isotropic correlations
340 Bool_t fCalculateSame; // calculate only 'same abs harmonics' correlations TBI
341 Bool_t fSkipZeroHarmonics; // skip correlations which have some of the harmonicc equal to zero
342 Bool_t fCalculateSameIsotropic; // calculate all isotropic correlations in 'same abs harmonic' TBI this can be implemented better
343 Bool_t fCalculateAll; // calculate all possible correlations
344 Int_t fDontGoBeyond; // do not go beyond fDontGoBeyond-p correlators
345 Bool_t fCalculateOnlyForHarmonicQC; // calculate only isotropic correlations in |fHarmonicQC|
346 Bool_t fCalculateOnlyForSC; // calculate only correlations needed for 'standard candles'
347 Bool_t fCalculateOnlyCos; // calculate only 'cos' correlations
348 Bool_t fCalculateOnlySin; // calculate only 'sin' correlations
350 // 4.) Event-by-event cumulants:
351 TList *fEbECumulantsList; // list to hold all e-b-e cumulants objects
352 TProfile *fEbECumulantsFlagsPro; // profile to hold all flags for e-b-e cumulants
353 TProfile *fEbECumulantsPro[2][8]; // multi-particle e-b-e cumulants [0=cos,1=sin][1p,2p,...,8p]
354 Bool_t fCalculateEbECumulants; // calculate and store e-b-e cumulants
357 TList *fWeightsList; // list to hold all weights objects
358 TProfile *fWeightsFlagsPro; // profile to hold all flags for weights
359 Bool_t fUseWeights[2][3]; // use weights [RP,POI][phi,pt,eta]
360 TH1D *fWeightsHist[2][3]; // histograms holding weights [RP,POI][phi,pt,eta]
363 TList *fNestedLoopsList; // list to hold all nested loops objects
364 TProfile *fNestedLoopsFlagsPro; // profile to hold all flags for nested loops
365 Bool_t fCrossCheckWithNestedLoops; // cross-check results with nested loops
366 Bool_t fCrossCheckDiffWithNestedLoops; // cross-check differential correlators with nested loops
367 Int_t fCrossCheckDiffCSCOBN[3]; // [0=cos,1=sin][1p,2p,...,4p][binNo]
368 TProfile *fNestedLoopsResultsCosPro; // profile to hold nested loops results (cosine)
369 TProfile *fNestedLoopsResultsSinPro; // profile to hold nested loops results (sinus)
370 TProfile *fNestedLoopsDiffResultsPro; // profile to hold differential nested loops results // TBI
372 // 7.) 'Standard candles':
373 TList *fStandardCandlesList; // list to hold all 'standard candles' objects
374 TProfile *fStandardCandlesFlagsPro; // profile to hold all flags fo 'standard candles'
375 Bool_t fCalculateStandardCandles; // calculate and store 'standard candles'
376 Bool_t fPropagateErrorSC; // propagate and store error for 'standard candles'
377 TH1D *fStandardCandlesHist; // histogram to hold results for 'standard candles'
378 TProfile2D *fProductsSCPro; // 2D profile to hold products of correlations needed for SC error propagation
381 TList *fQcumulantsList; // list to hold all Q-cumulants objects
382 TProfile *fQcumulantsFlagsPro; // profile to hold all flags for Q-cumulants
383 Bool_t fCalculateQcumulants; // calculate and store Q-cumulants
384 Int_t fHarmonicQC; // calculate Q-cumulants in this harmonic (default is 2)
385 Bool_t fPropagateErrorQC; // propagate and store error for Q-cumulants
386 TH1D *fQcumulantsHist; // two- and multi-particle Q-cumulants
387 TH1D *fReferenceFlowHist; // reference flow from two- and multi-particle Q-cumulants
388 TProfile2D *fProductsQCPro; // 2D profile to hold products of correlations needed for QC error propagation
390 // 9.) Differential correlations:
391 TList *fDiffCorrelationsList; // list to hold all correlations objects
392 TProfile *fDiffCorrelationsFlagsPro; // profile to hold all flags for correlations
393 Bool_t fCalculateDiffCorrelations; // calculate and store differential correlations
394 Bool_t fCalculateDiffCos; // calculate and store differential cosine correlations (kTRUE by default)
395 Bool_t fCalculateDiffSin; // calculate and store differential sinus correlations (kFALSE by default)
396 Bool_t fCalculateDiffCorrelationsVsPt; // calculate differential correlations vs pt (default), or vs eta
397 Bool_t fUseDefaultBinning; // use default binning in pt or in eta
398 Int_t fnDiffBins; // number of differential bins in pt or in eta (when non-default binning is used)
399 Double_t *fRangesDiffBins; // ranges for differential bins in pt or in eta (when non-default binning is used)
400 Int_t fDiffHarmonics[4][4]; // harmonics for differential correlations [order][{n1},{n1,n2},...,{n1,n2,n3,n4}]
401 TProfile *fDiffCorrelationsPro[2][4]; // multi-particle correlations [0=cos,1=sin][1p,2p,3p,4p]
402 UInt_t fDiffBinNo; // differential bin number
404 ClassDef(AliFlowAnalysisWithMultiparticleCorrelations,1);
408 //================================================================================================================