]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/FLOW/Base/AliFlowAnalysisWithMultiparticleCorrelations.h
d977a624bb3758a6beb4b631ffc4f106c044b3aa
[u/mrichter/AliRoot.git] / PWG / FLOW / Base / AliFlowAnalysisWithMultiparticleCorrelations.h
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  * flow analysis with multi-particle *
9  *           correlations            * 
10  *                                   * 
11  * author: Ante Bilandzic            * 
12  *        (abilandzic@gmail.com)     *
13  ************************************/ 
14
15 #ifndef ALIFLOWANALYSISWITHMULTIPARTICLECORRELATIONS_H
16 #define ALIFLOWANALYSISWITHMULTIPARTICLECORRELATIONS_H
17
18 #include "TH1D.h"
19 #include "TH2D.h"
20 #include "TProfile.h"
21 #include "TProfile2D.h"
22 #include "TFile.h"
23 #include "TComplex.h"
24 #include "TDirectoryFile.h"
25 #include "Riostream.h"
26 #include "TRandom3.h"
27 #include "TSystem.h"
28 #include "TGraphErrors.h"
29 #include "TStopwatch.h"
30 #include "AliFlowEventSimple.h"
31 #include "AliFlowTrackSimple.h"
32
33 class AliFlowAnalysisWithMultiparticleCorrelations{
34  public:
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;
44   // 6.) The rest.
45
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(); 
55
56   // 1.) Method Init() and methods called in it (!):
57   virtual void Init();
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();
70    
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);
83
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();
90
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(); 
99
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 min) {fUseInternalFlags = kTRUE; this->fMinNoRPs = min;};
107   Int_t GetMinNoRPs() const {return this->fMinNoRPs;};
108   void SetMaxNoRPs(Int_t max) {fUseInternalFlags = kTRUE; this->fMaxNoRPs = max;};
109   Int_t GetMaxNoRPs() const {return this->fMaxNoRPs;};
110   void SetExactNoRPs(Int_t exact) {fUseInternalFlags = kTRUE; this->fExactNoRPs = exact;};
111   Int_t GetExactNoRPs() const {return this->fExactNoRPs;};
112   void SetAnalysisTag(const char *at) {this->fAnalysisTag = TString(at);};
113   TString GetAnalysisTag() const {return this->fAnalysisTag;};
114   void SetDumpThePoints(Bool_t dtp, Int_t max) {this->fDumpThePoints = dtp; this->fMaxNoEventsPerFile = max;};
115
116   //  5.1.) Control histograms:  
117   void SetControlHistogramsList(TList* const chl) {this->fControlHistogramsList = chl;};
118   TList* GetControlHistogramsList() const {return this->fControlHistogramsList;} 
119   void SetControlHistogramsFlagsPro(TProfile* const chfp) {this->fControlHistogramsFlagsPro = chfp;};
120   TProfile* GetControlHistogramsFlagsPro() const {return this->fControlHistogramsFlagsPro;}; 
121   void SetFillControlHistograms(Bool_t fch) {this->fFillControlHistograms = fch;};
122   Bool_t GetFillControlHistograms() const {return this->fFillControlHistograms;};
123   void SetFillKinematicsHist(Bool_t fkh) {this->fFillKinematicsHist = fkh;};
124   Bool_t GetFillKinematicsHist() const {return this->fFillKinematicsHist;};
125   void SetFillMultDistributionsHist(Bool_t mdh) {this->fFillMultDistributionsHist = mdh;};
126   Bool_t GetFillMultDistributionsHist() const {return this->fFillMultDistributionsHist;};
127   void SetFillMultCorrelationsHist(Bool_t mch) {this->fFillMultCorrelationsHist = mch;};
128   Bool_t GetFillMultCorrelationsHist() const {return this->fFillMultCorrelationsHist;};
129   void SetDontFill(const char *type) 
130   {
131    if(TString(type).EqualTo("RP")){this->fDontFill[0] = kTRUE;} 
132    else if(TString(type).EqualTo("POI")){this->fDontFill[1] = kTRUE;} 
133    else if(TString(type).EqualTo("REF")){this->fDontFill[2] = kTRUE;} 
134    else{Fatal("void SetDontFill(const char *type)","type = %s ???? Allowed: RP, POI and REF.",type);}
135   }; // void SetDontFill(const char *type) 
136   void SetnBins(const char *type, const char *variable, Int_t nBins); // .cxx
137   void SetMin(const char *type, const char *variable, Double_t min); // .cxx
138   void SetMax(const char *type, const char *variable, Double_t max); // .cxx
139   void SetnBinsMult(const char *type, Int_t nBinsMult); // .cxx
140   void SetMinMult(const char *type, Double_t minMult); // .cxx
141   void SetMaxMult(const char *type, Double_t maxMult); // .cxx
142
143   //  5.2.) Q-vectors:
144   void SetQvectorList(TList* const qvl) {this->fQvectorList = qvl;};
145   TList* GetQvectorList() const {return this->fQvectorList;} 
146   void SetQvectorFlagsPro(TProfile* const qvfp) {this->fQvectorFlagsPro = qvfp;};
147   TProfile* GetQvectorFlagsPro() const {return this->fQvectorFlagsPro;}; 
148   void SetCalculateQvector(Bool_t cqv) {this->fCalculateQvector = cqv;};
149   Bool_t GetCalculateQvector() const {return this->fCalculateQvector;};
150   void SetCalculateDiffQvectors(Bool_t cdqv) {this->fCalculateDiffQvectors = cdqv;};
151   Bool_t GetCalculateDiffQvectors() const {return this->fCalculateDiffQvectors;};
152
153   //  5.3.) Correlations:
154   void SetCorrelationsList(TList* const cl) {this->fCorrelationsList = cl;};
155   TList* GetCorrelationsList() const {return this->fCorrelationsList;} 
156   void SetCorrelationsFlagsPro(TProfile* const cfp) {this->fCorrelationsFlagsPro = cfp;};
157   TProfile* GetCorrelationsFlagsPro() const {return this->fCorrelationsFlagsPro;}; 
158   void SetCalculateCorrelations(Bool_t cc) {this->fCalculateCorrelations = cc;};
159   Bool_t GetCalculateCorrelations() const {return this->fCalculateCorrelations;};
160   void SetCalculateIsotropic(Bool_t ci) {this->fCalculateIsotropic = ci;};
161   Bool_t GetCalculateIsotropic() const {return this->fCalculateIsotropic;};
162   void SetCalculateSame(Bool_t cs) {this->fCalculateSame = cs;};
163   Bool_t GetCalculateSame() const {return this->fCalculateSame;};
164   void SetSkipZeroHarmonics(Bool_t szh) {this->fSkipZeroHarmonics = szh;};
165   Bool_t GetSkipZeroHarmonics() const {return this->fSkipZeroHarmonics;};
166   void SetCalculateSameIsotropic(Bool_t csi) {this->fCalculateSameIsotropic = csi;};
167   Bool_t GetCalculateSameIsotropic() const {return this->fCalculateSameIsotropic;};
168   void SetCalculateAll(Bool_t ca) {this->fCalculateAll = ca;};
169   Bool_t GetCalculateAll() const {return this->fCalculateAll;};
170   void SetDontGoBeyond(Int_t dgb) {this->fDontGoBeyond = dgb;};
171   Int_t GetDontGoBeyond() const {return this->fDontGoBeyond;};
172   void SetCalculateOnlyForHarmonicQC(Bool_t cofhqc) {this->fCalculateOnlyForHarmonicQC = cofhqc;};
173   Bool_t GetCalculateOnlyForHarmonicQC() const {return this->fCalculateOnlyForHarmonicQC;};
174   void SetCalculateOnlyForSC(Bool_t cofsc) {this->fCalculateOnlyForSC = cofsc;};
175   Bool_t GetCalculateOnlyForSC() const {return this->fCalculateOnlyForSC;};
176   void SetCalculateOnlyCos(Bool_t coc) {this->fCalculateOnlyCos = coc;};
177   Bool_t GetCalculateOnlyCos() const {return this->fCalculateOnlyCos;};
178   void SetCalculateOnlySin(Bool_t cos) {this->fCalculateOnlySin = cos;};
179   Bool_t GetCalculateOnlySin() const {return this->fCalculateOnlySin;};
180
181   //  5.4.) Event-by-event cumulants:
182   void SetEbECumulantsList(TList* const ebecl) {this->fEbECumulantsList = ebecl;};
183   TList* GetEbECumulantsList() const {return this->fEbECumulantsList;} 
184   void SetEbECumulantsFlagsPro(TProfile* const ebecfp) {this->fEbECumulantsFlagsPro = ebecfp;};
185   TProfile* GetEbECumulantsFlagsPro() const {return this->fEbECumulantsFlagsPro;}; 
186   void SetCalculateEbECumulants(Bool_t cebec) {this->fCalculateEbECumulants = cebec;};
187   Bool_t GetCalculateEbECumulants() const {return this->fCalculateEbECumulants;};
188
189   //  5.5.) Weights: 
190   void SetWeightsList(TList* const wl) {this->fWeightsList = (TList*)wl->Clone();};
191   TList* GetWeightsList() const {return this->fWeightsList;} 
192   void SetWeightsFlagsPro(TProfile* const wfp) {this->fWeightsFlagsPro = wfp;};
193   TProfile* GetWeightsFlagsPro() const {return this->fWeightsFlagsPro;}; 
194   void SetWeightsHist(TH1D* const hist, const char *type, const char *variable); // .cxx
195   
196   //  5.6.) Nested loops:
197   void SetNestedLoopsList(TList* const nll) {this->fNestedLoopsList = nll;} 
198   TList* GetNestedLoopsList() const {return this->fNestedLoopsList;} 
199   void SetNestedLoopsFlagsPro(TProfile* const nlfp) {this->fNestedLoopsFlagsPro = nlfp;};
200   TProfile* GetNestedLoopsFlagsPro() const {return this->fNestedLoopsFlagsPro;}; 
201   void SetCrossCheckWithNestedLoops(Bool_t ccwnl) {this->fCrossCheckWithNestedLoops = ccwnl;};
202   Bool_t GetCrossCheckWithNestedLoops() const {return this->fCrossCheckWithNestedLoops;};
203   void SetCrossCheckDiffWithNestedLoops(Bool_t ccdwnl) {this->fCrossCheckDiffWithNestedLoops = ccdwnl;};
204   Bool_t GetCrossCheckDiffWithNestedLoops() const {return this->fCrossCheckDiffWithNestedLoops;};
205   void SetCrossCheckDiffCSCOBN(Int_t cs, Int_t co, Int_t bn)  
206   {
207    this->fCrossCheckDiffCSCOBN[0] = cs; // cos/sin
208    this->fCrossCheckDiffCSCOBN[1] = co; // correlator order [1p,2p,3p,4p]
209    this->fCrossCheckDiffCSCOBN[2] = bn; // bin number
210   };
211
212   // 5.7.) 'Standard candles':
213   void SetStandardCandlesList(TList* const scl) {this->fStandardCandlesList = scl;} 
214   TList* GetStandardCandlesList() const {return this->fStandardCandlesList;} 
215   void SetStandardCandlesFlagsPro(TProfile* const scfp) {this->fStandardCandlesFlagsPro = scfp;};
216   TProfile* GetStandardCandlesFlagsPro() const {return this->fStandardCandlesFlagsPro;}; 
217   void SetCalculateStandardCandles(Bool_t csc) {this->fCalculateStandardCandles = csc;};
218   Bool_t GetCalculateStandardCandles() const {return this->fCalculateStandardCandles;};
219   void SetPropagateErrorSC(Bool_t pesc) {this->fPropagateErrorSC = pesc;};
220   Bool_t GetPropagateErrorSC() const {return this->fPropagateErrorSC;};
221   void SetStandardCandlesHist(TH1D* const sch) {this->fStandardCandlesHist = sch;};
222   TH1D* GetStandardCandlesHist() const {return this->fStandardCandlesHist;}; 
223   void SetProductsSCPro(TProfile2D* const psc) {this->fProductsSCPro = psc;};
224   TProfile2D* GetProductsSCPro() const {return this->fProductsSCPro;}; 
225
226   //  5.8.) Q-cumulants:
227   void SetQcumulantsList(TList* const qcl) {this->fQcumulantsList = qcl;};
228   TList* GetQcumulantsList() const {return this->fQcumulantsList;} 
229   void SetQcumulantsFlagsPro(TProfile* const qcfp) {this->fQcumulantsFlagsPro = qcfp;};
230   TProfile* GetQcumulantsFlagsPro() const {return this->fQcumulantsFlagsPro;}; 
231   void SetCalculateQcumulants(Bool_t cqc) {this->fCalculateQcumulants = cqc;};
232   Bool_t GetCalculateQcumulants() const {return this->fCalculateQcumulants;};
233   void SetHarmonicQC(Int_t hqc) {this->fHarmonicQC = hqc;};
234   Int_t GetHarmonicQC() const {return this->fHarmonicQC;};
235   void SetPropagateErrorQC(Bool_t peqc) {this->fPropagateErrorQC = peqc;};
236   Bool_t GetPropagateErrorQC() const {return this->fPropagateErrorQC;};
237   void SetQcumulantsHist(TH1D* const qch) {this->fQcumulantsHist = qch;};
238   TH1D* GetQcumulantsHist() const {return this->fQcumulantsHist;}; 
239   void SetReferenceFlowHist(TH1D* const rfh) {this->fReferenceFlowHist = rfh;};
240   TH1D* GetReferenceFlowHist() const {return this->fReferenceFlowHist;}; 
241   void SetProductsQCPro(TProfile2D* const pqc) {this->fProductsQCPro = pqc;};
242   TProfile2D* GetProductsQCPro() const {return this->fProductsQCPro;}; 
243
244   //  5.9.) Differential correlations:
245   void SetDiffCorrelationsList(TList* const dcl) {this->fDiffCorrelationsList = dcl;};
246   TList* GetDiffCorrelationsList() const {return this->fDiffCorrelationsList;} 
247   void SetDiffCorrelationsFlagsPro(TProfile* const cdfp) {this->fDiffCorrelationsFlagsPro = cdfp;};
248   TProfile* GetDiffCorrelationsFlagsPro() const {return this->fDiffCorrelationsFlagsPro;}; 
249   void SetCalculateDiffCorrelations(Bool_t cdc) {this->fCalculateDiffCorrelations = cdc;};
250   Bool_t GetCalculateDiffCorrelations() const {return this->fCalculateDiffCorrelations;};
251   void SetDiffHarmonics(Int_t order, Int_t *harmonics); // see implementation in .cxx file 
252   void SetCalculateDiffCos(Bool_t cdc) {this->fCalculateDiffCos = cdc;};
253   Bool_t GetCalculateDiffCos() const {return this->fCalculateDiffCos;};
254   void SetCalculateDiffSin(Bool_t cds) {this->fCalculateDiffSin = cds;};
255   Bool_t GetCalculateDiffSin() const {return this->fCalculateDiffSin;};
256   void SetCalculateDiffCorrelationsVsPt(Bool_t cdcvspt) {this->fCalculateDiffCorrelationsVsPt = cdcvspt;};
257   Bool_t GetCalculateDiffCorrelationsVsPt() const {return this->fCalculateDiffCorrelationsVsPt;};
258   void SetUseDefaultBinning(Bool_t udb) {this->fUseDefaultBinning = udb;};
259   Bool_t GetUseDefaultBinning() const {return this->fUseDefaultBinning;};
260   void SetnDiffBins(Int_t ndb) {this->fnDiffBins = ndb;};
261   Int_t GetnDiffBins() const {return this->fnDiffBins;};
262   void SetRangesDiffBins(Double_t* const rdb) {this->fRangesDiffBins = rdb;};
263   Double_t* GetRangesDiffBins() const {return this->fRangesDiffBins;};
264
265   // 6.) The rest:
266   virtual void WriteHistograms(TString outputFileName);
267   virtual void WriteHistograms(TDirectoryFile *outputFileName);
268   virtual TComplex Q(Int_t n, Int_t p);
269   virtual TComplex p(Int_t n, Int_t p);
270   virtual TComplex q(Int_t n, Int_t p);
271   virtual TComplex One(Int_t n1);
272   virtual TComplex Two(Int_t n1, Int_t n2);
273   virtual TComplex Three(Int_t n1, Int_t n2, Int_t n3);
274   virtual TComplex Four(Int_t n1, Int_t n2, Int_t n3, Int_t n4);
275   virtual TComplex Five(Int_t n1, Int_t n2, Int_t n3, Int_t n4, Int_t n5);
276   virtual TComplex Six(Int_t n1, Int_t n2, Int_t n3, Int_t n4, Int_t n5, Int_t n6);
277   virtual TComplex Seven(Int_t n1, Int_t n2, Int_t n3, Int_t n4, Int_t n5, Int_t n6, Int_t n7);
278   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);
279   virtual TComplex OneDiff(Int_t n1);
280   virtual TComplex TwoDiff(Int_t n1, Int_t n2);
281   virtual TComplex ThreeDiff(Int_t n1, Int_t n2, Int_t n3);
282   virtual TComplex FourDiff(Int_t n1, Int_t n2, Int_t n3, Int_t n4);
283   virtual Double_t Weight(const Double_t &value, const char *type, const char *variable); // value, [RP,POI], [phi,pt,eta]
284   virtual Double_t CastStringToCorrelation(const char *string, Bool_t numerator);
285   virtual Double_t Covariance(const char *x, const char *y, TProfile2D *profile2D, Bool_t bUnbiasedEstimator = kFALSE);
286   virtual TComplex Recursion(Int_t n, Int_t* harmonic, Int_t mult = 1, Int_t skip = 0); // Credits: Kristjan Gulbrandsen (gulbrand@nbi.dk) 
287   virtual void CalculateProductsOfCorrelations(AliFlowEventSimple *anEvent, TProfile2D *profile2D);
288   static void DumpPointsForDurham(TGraphErrors *ge);
289   static void DumpPointsForDurham(TH1D *h);
290   static void DumpPointsForDurham(TH1F *h);
291   virtual void DumpThePoints(AliFlowEventSimple *anEvent);
292   TH1D* GetHistogramWithWeights(const char *filePath, const char *listName, const char *type, const char *variable);
293
294  private:
295   AliFlowAnalysisWithMultiparticleCorrelations(const AliFlowAnalysisWithMultiparticleCorrelations& afawQc);
296   AliFlowAnalysisWithMultiparticleCorrelations& operator=(const AliFlowAnalysisWithMultiparticleCorrelations& afawQc); 
297   // Data members are grouped as:
298   // 0.) Base list and internal flags;
299   // 1.) Control histograms;  
300   // 2.) Q-vectors;
301   // 3.) Correlations;
302   // 4.) Event-by-event cumulants;
303   // 5.) Weights;
304   // 6.) Nested loops;
305   // 7.) 'Standard candles';
306   // 8.) Q-cumulants;
307   // 9.) Differential correlations. 
308
309   // 0.) Base list and internal flags:
310   TList* fHistList;            // base list to hold all output object (a.k.a. grandmother of all lists)
311   TProfile *fInternalFlagsPro; // profile to hold all internal flags and settings
312   Bool_t fUseInternalFlags;    // use internal flags (automatically set if some internal flag is used)
313   Int_t fMinNoRPs;             // minimum number of RPs required for the analysis 
314   Int_t fMaxNoRPs;             // maximum number of RPs allowed for the analysis 
315   Int_t fExactNoRPs;           // exact (randomly shuffled) number of RPs selected for the analysis 
316   Bool_t fPropagateError;      // prevent error propagation if something strange happens during calculations 
317   TString fAnalysisTag;        // tag internally this analysis
318   Bool_t fDumpThePoints;       // dump the data points into the external file 
319   Int_t fMaxNoEventsPerFile;   // maximum number of events to be dumped in a single file
320
321   // 1.) Control histograms:  
322   TList *fControlHistogramsList;        // list to hold all 'control histograms' objects
323   TProfile *fControlHistogramsFlagsPro; // profile to hold all flags for control histograms
324   Bool_t fFillControlHistograms;        // fill or not control histograms (by default they are all filled)
325   Bool_t fFillKinematicsHist;           // fill or not fKinematicsHist[2][3]
326   Bool_t fFillMultDistributionsHist;    // fill or not TH1D *fMultDistributionsHist[3]    
327   Bool_t fFillMultCorrelationsHist;     // fill or not TH2D *fMultCorrelationsHist[3]
328   Bool_t fDontFill[3];                  // don't fill control histograms [0=RP,1=POI,2=REF]  
329   TH1D *fKinematicsHist[2][3];          // [RP,POI][phi,pt,eta] distributions
330   TH1D *fMultDistributionsHist[3];      // multiplicity distribution [RP,POI,reference multiplicity]
331   TH2I *fMultCorrelationsHist[3];       // [RP vs. POI, RP vs. refMult, POI vs. refMult]  
332   Int_t fnBins[2][3];                   // [RP,POI][phi,pt,eta], corresponds to fKinematicsHist[2][3]
333   Double_t fMin[2][3];                  // [RP,POI][phi,pt,eta], corresponds to fKinematicsHist[2][3]
334   Double_t fMax[2][3];                  // [RP,POI][phi,pt,eta], corresponds to fKinematicsHist[2][3]
335   Int_t fnBinsMult[3];                  // [RP,POI,REF], corresponds to fMultDistributionsHist[3]   
336   Double_t fMinMult[3];                 // [RP,POI,REF], corresponds to fMultDistributionsHist[3]   
337   Double_t fMaxMult[3];                 // [RP,POI,REF], corresponds to fMultDistributionsHist[3]   
338
339   // 2.) Q-vectors:
340   TList *fQvectorList;           // list to hold all Q-vector objects       
341   TProfile *fQvectorFlagsPro;    // profile to hold all flags for Q-vector
342   Bool_t fCalculateQvector;      // to calculate or not to calculate Q-vector components, that's a Boolean...
343   TComplex fQvector[49][9];      // Q-vector components [fMaxHarmonic*fMaxCorrelator+1][fMaxCorrelator+1] = [6*8+1][8+1]  
344   Bool_t fCalculateDiffQvectors; // to calculate or not to calculate p- and q-vector components, that's a Boolean...  
345   TComplex fpvector[100][49][9]; // p-vector components [bin][fMaxHarmonic*fMaxCorrelator+1][fMaxCorrelator+1] = [6*8+1][8+1] TBI hardwired 100
346   TComplex fqvector[100][49][9]; // q-vector components [bin][fMaxHarmonic*fMaxCorrelator+1][fMaxCorrelator+1] = [6*8+1][8+1] TBI hardwired 100
347
348   // 3.) Correlations:
349   TList *fCorrelationsList;           // list to hold all correlations objects
350   TProfile *fCorrelationsFlagsPro;    // profile to hold all flags for correlations
351   TProfile *fCorrelationsPro[2][8];   // multi-particle correlations [0=cos,1=sin][1p,2p,...,8p]
352   Bool_t fCalculateCorrelations;      // calculate and store correlations
353   Int_t fMaxHarmonic;                 // 6 (not going beyond v6, if you change this value, change also fQvector[49][9]) 
354   Int_t fMaxCorrelator;               // 8 (not going beyond 8-p correlations, if you change this value, change also fQvector[49][9]) 
355   Bool_t fCalculateIsotropic;         // calculate only isotropic correlations
356   Bool_t fCalculateSame;              // calculate only 'same abs harmonics' correlations TBI 
357   Bool_t fSkipZeroHarmonics;          // skip correlations which have some of the harmonicc equal to zero
358   Bool_t fCalculateSameIsotropic;     // calculate all isotropic correlations in 'same abs harmonic' TBI this can be implemented better
359   Bool_t fCalculateAll;               // calculate all possible correlations 
360   Int_t  fDontGoBeyond;               // do not go beyond fDontGoBeyond-p correlators
361   Bool_t fCalculateOnlyForHarmonicQC; // calculate only isotropic correlations in |fHarmonicQC|
362   Bool_t fCalculateOnlyForSC;         // calculate only correlations needed for 'standard candles'
363   Bool_t fCalculateOnlyCos;           // calculate only 'cos' correlations
364   Bool_t fCalculateOnlySin;           // calculate only 'sin' correlations
365
366   // 4.) Event-by-event cumulants:
367   TList *fEbECumulantsList;         // list to hold all e-b-e cumulants objects
368   TProfile *fEbECumulantsFlagsPro;  // profile to hold all flags for e-b-e cumulants
369   TProfile *fEbECumulantsPro[2][8]; // multi-particle e-b-e cumulants [0=cos,1=sin][1p,2p,...,8p]
370   Bool_t fCalculateEbECumulants;    // calculate and store e-b-e cumulants
371  
372   // 5.) Weights:
373   TList *fWeightsList;        // list to hold all weights objects
374   TProfile *fWeightsFlagsPro; // profile to hold all flags for weights
375   Bool_t fUseWeights[2][3];   // use weights [RP,POI][phi,pt,eta]
376   TH1D *fWeightsHist[2][3];   // histograms holding weights [RP,POI][phi,pt,eta]
377
378   // 6.) Nested loops:
379   TList *fNestedLoopsList;               // list to hold all nested loops objects
380   TProfile *fNestedLoopsFlagsPro;        // profile to hold all flags for nested loops
381   Bool_t fCrossCheckWithNestedLoops;     // cross-check results with nested loops
382   Bool_t fCrossCheckDiffWithNestedLoops; // cross-check differential correlators with nested loops
383   Int_t fCrossCheckDiffCSCOBN[3];        // [0=cos,1=sin][1p,2p,...,4p][binNo]
384   TProfile *fNestedLoopsResultsCosPro;   // profile to hold nested loops results (cosine)
385   TProfile *fNestedLoopsResultsSinPro;   // profile to hold nested loops results (sinus)
386   TProfile *fNestedLoopsDiffResultsPro;  // profile to hold differential nested loops results // TBI
387
388   // 7.) 'Standard candles':
389   TList *fStandardCandlesList;        // list to hold all 'standard candles' objects
390   TProfile *fStandardCandlesFlagsPro; // profile to hold all flags fo 'standard candles'
391   Bool_t fCalculateStandardCandles;   // calculate and store 'standard candles'
392   Bool_t fPropagateErrorSC;           // propagate and store error for 'standard candles'
393   TH1D *fStandardCandlesHist;         // histogram to hold results for 'standard candles'
394   TProfile2D *fProductsSCPro;         // 2D profile to hold products of correlations needed for SC error propagation
395
396   // 8.) Q-cumulants:
397   TList *fQcumulantsList;        // list to hold all Q-cumulants objects
398   TProfile *fQcumulantsFlagsPro; // profile to hold all flags for Q-cumulants
399   Bool_t fCalculateQcumulants;   // calculate and store Q-cumulants
400   Int_t fHarmonicQC;             // calculate Q-cumulants in this harmonic (default is 2) 
401   Bool_t fPropagateErrorQC;      // propagate and store error for Q-cumulants
402   TH1D *fQcumulantsHist;         // two- and multi-particle Q-cumulants
403   TH1D *fReferenceFlowHist;      // reference flow from two- and multi-particle Q-cumulants
404   TProfile2D *fProductsQCPro;    // 2D profile to hold products of correlations needed for QC error propagation
405
406   // 9.) Differential correlations:
407   TList *fDiffCorrelationsList;          // list to hold all correlations objects
408   TProfile *fDiffCorrelationsFlagsPro;   // profile to hold all flags for correlations
409   Bool_t fCalculateDiffCorrelations;     // calculate and store differential correlations
410   Bool_t fCalculateDiffCos;              // calculate and store differential cosine correlations (kTRUE by default)
411   Bool_t fCalculateDiffSin;              // calculate and store differential sinus correlations (kFALSE by default)
412   Bool_t fCalculateDiffCorrelationsVsPt; // calculate differential correlations vs pt (default), or vs eta
413   Bool_t fUseDefaultBinning;             // use default binning in pt or in eta
414   Int_t fnDiffBins;                      // number of differential bins in pt or in eta (when non-default binning is used)
415   Double_t *fRangesDiffBins;             // ranges for differential bins in pt or in eta (when non-default binning is used)
416   Int_t fDiffHarmonics[4][4];            // harmonics for differential correlations [order][{n1},{n1,n2},...,{n1,n2,n3,n4}] 
417   TProfile *fDiffCorrelationsPro[2][4];  // multi-particle correlations [0=cos,1=sin][1p,2p,3p,4p]
418   UInt_t fDiffBinNo;                     // differential bin number
419
420   ClassDef(AliFlowAnalysisWithMultiparticleCorrelations,1);
421
422 };
423
424 //================================================================================================================
425
426 #endif
427
428
429
430
431