AliPhysics  2b88e80 (2b88e80)
AliFlowAnalysisWithMultiparticleCorrelations.h
Go to the documentation of this file.
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 "TArrayI.h"
29 #include "TGraphErrors.h"
30 #include "TStopwatch.h"
31 #include "AliFlowEventSimple.h"
32 #include "AliFlowTrackSimple.h"
33 
35  public:
38  // Member functions are grouped as:
39  // 0.) Methods called in the constructor;
40  // 1.) Method Init() and methods called in it (!);
41  // 2.) Method Make() and methods called in it;
42  // 3.) Method Finish() and methods called in it;
43  // 4.) Method GetOutputHistograms() and methods called in it;
44  // 5.) Setters and getters;
45  // 6.) The rest.
46 
47  // 0.) Methods called in the constructor:
49  virtual void InitializeArraysForQvector();
50  virtual void InitializeArraysForCorrelations();
51  virtual void InitializeArraysForEbECumulants();
52  virtual void InitializeArraysForWeights();
53  virtual void InitializeArraysForQcumulants();
55  virtual void InitializeArraysForSymmetryPlanes();
56  virtual void InitializeArraysForNestedLoops();
57  virtual void InitializeArraysForEtaGaps();
58 
59  // 1.) Method Init() and methods called in it (!):
60  virtual void Init();
61  virtual void CrossCheckSettings();
62  virtual void BookAndNestAllLists();
63  virtual void BookEverythingForBase();
65  virtual void BookEverythingForQvector();
66  virtual void BookEverythingForWeights();
67  virtual void BookEverythingForCorrelations();
68  virtual void BookEverythingForEbECumulants();
69  virtual void BookEverythingForNestedLoops();
70  virtual void BookEverythingForStandardCandles();
71  virtual void BookEverythingForQcumulants();
72  virtual void BookEverythingForDiffCorrelations();
73  virtual void BookEverythingForSymmetryPlanes();
74  virtual void BookEverythingForEtaGaps();
75 
76  // 2.) Method Make() and methods called in it:
77  virtual void Make(AliFlowEventSimple *anEvent);
79  virtual void CrossCheckPointersUsedInMake();
80  virtual void DetermineRandomIndices(AliFlowEventSimple *anEvent);
81  virtual void FillControlHistograms(AliFlowEventSimple *anEvent);
82  virtual void FillQvector(AliFlowEventSimple *anEvent);
83  virtual void CalculateCorrelations(AliFlowEventSimple *anEvent);
84  virtual void CalculateDiffCorrelations(AliFlowEventSimple *anEvent);
85  virtual void CalculateEbECumulants(AliFlowEventSimple *anEvent);
86  virtual void CalculateSymmetryPlanes(AliFlowEventSimple *anEvent);
87  virtual void CalculateEtaGaps(AliFlowEventSimple *anEvent);
88  virtual void ResetQvector();
89  virtual void CrossCheckWithNestedLoops(AliFlowEventSimple *anEvent);
91 
92  // 3.) Method Finish() and methods called in it:
93  virtual void Finish();
94  virtual void CrossCheckPointersUsedInFinish();
95  virtual void CalculateStandardCandles();
96  virtual void CalculateQcumulants();
97  virtual void CalculateReferenceFlow();
98 
99  // 4.) Method GetOutputHistograms() and methods called in it:
100  virtual void GetOutputHistograms(TList *histList);
101  virtual void GetPointersForControlHistograms();
102  virtual void GetPointersForQvector();
103  virtual void GetPointersForCorrelations();
104  virtual void GetPointersForStandardCandles();
105  virtual void GetPointersForQcumulants();
106  virtual void GetPointersForDiffCorrelations();
107  virtual void GetPointersForSymmetryPlanes();
108 
109  // 5.) Setters and getters:
110  // 5.0.) Base list and internal flags:
111  void SetHistList(TList* const hlist) {this->fHistList = hlist;}
112  TList* GetHistList() const {return this->fHistList;}
113  void SetInternalFlagsPro(TProfile* const ifp) {this->fInternalFlagsPro = ifp;};
114  TProfile* GetInternalFlagsPro() const {return this->fInternalFlagsPro;};
115  void SetMinNoRPs(Int_t min) {fUseInternalFlags = kTRUE; this->fMinNoRPs = min;};
116  Int_t GetMinNoRPs() const {return this->fMinNoRPs;};
117  void SetMaxNoRPs(Int_t max) {fUseInternalFlags = kTRUE; this->fMaxNoRPs = max;};
118  Int_t GetMaxNoRPs() const {return this->fMaxNoRPs;};
119  void SetExactNoRPs(Int_t exact) {fUseInternalFlags = kTRUE; this->fExactNoRPs = exact;};
120  Int_t GetExactNoRPs() const {return this->fExactNoRPs;};
121  void SetAnalysisTag(const char *at) {this->fAnalysisTag = TString(at);};
122  TString GetAnalysisTag() const {return this->fAnalysisTag;};
123  void SetDumpThePoints(Bool_t dtp, Int_t max) {this->fDumpThePoints = dtp; this->fMaxNoEventsPerFile = max;};
124  void SetSelectRandomlyRPs(Int_t nSelectedRandomlyRPs) {this->fSelectRandomlyRPs = kTRUE; this->fnSelectedRandomlyRPs = nSelectedRandomlyRPs;};
125 
126  // 5.1.) Control histograms:
127  void SetControlHistogramsList(TList* const chl) {this->fControlHistogramsList = chl;};
129  void SetControlHistogramsFlagsPro(TProfile* const chfp) {this->fControlHistogramsFlagsPro = chfp;};
130  TProfile* GetControlHistogramsFlagsPro() const {return this->fControlHistogramsFlagsPro;};
139  void SetDontFill(const char *type)
140  {
141  if(TString(type).EqualTo("RP")){this->fDontFill[0] = kTRUE;}
142  else if(TString(type).EqualTo("POI")){this->fDontFill[1] = kTRUE;}
143  else if(TString(type).EqualTo("REF")){this->fDontFill[2] = kTRUE;}
144  else{Fatal("void SetDontFill(const char *type)","type = %s ???? Allowed: RP, POI and REF.",type);}
145  }; // void SetDontFill(const char *type)
146  void SetnBins(const char *type, const char *variable, Int_t nBins); // .cxx
147  void SetMin(const char *type, const char *variable, Double_t min); // .cxx
148  void SetMax(const char *type, const char *variable, Double_t max); // .cxx
149  void SetnBinsMult(const char *type, Int_t nBinsMult); // .cxx
150  void SetMinMult(const char *type, Double_t minMult); // .cxx
151  void SetMaxMult(const char *type, Double_t maxMult); // .cxx
152  void SetIntervalsToSkip(const char *ppe, Int_t n, Double_t *boundaries); // .cxx
153 
154  // 5.2.) Q-vectors:
155  void SetQvectorList(TList* const qvl) {this->fQvectorList = qvl;};
156  TList* GetQvectorList() const {return this->fQvectorList;}
157  void SetQvectorFlagsPro(TProfile* const qvfp) {this->fQvectorFlagsPro = qvfp;};
158  TProfile* GetQvectorFlagsPro() const {return this->fQvectorFlagsPro;};
159  void SetCalculateQvector(Bool_t cqv) {this->fCalculateQvector = cqv;};
163 
164  // 5.3.) Correlations:
165  void SetCorrelationsList(TList* const cl) {this->fCorrelationsList = cl;};
166  TList* GetCorrelationsList() const {return this->fCorrelationsList;}
167  void SetCorrelationsFlagsPro(TProfile* const cfp) {this->fCorrelationsFlagsPro = cfp;};
168  TProfile* GetCorrelationsFlagsPro() const {return this->fCorrelationsFlagsPro;};
173  void SetCalculateSame(Bool_t cs) {this->fCalculateSame = cs;};
174  Bool_t GetCalculateSame() const {return this->fCalculateSame;};
179  void SetCalculateAll(Bool_t ca) {this->fCalculateAll = ca;};
180  Bool_t GetCalculateAll() const {return this->fCalculateAll;};
181  void SetDontGoBeyond(Int_t dgb) {this->fDontGoBeyond = dgb;};
182  Int_t GetDontGoBeyond() const {return this->fDontGoBeyond;};
185  void SetCalculateOnlyForSC(Bool_t cofsc) {this->fCalculateOnlyForSC = cofsc;};
187  void SetCalculateOnlyCos(Bool_t coc) {this->fCalculateOnlyCos = coc;};
189  void SetCalculateOnlySin(Bool_t cos) {this->fCalculateOnlySin = cos;};
191 
192  // 5.4.) Event-by-event cumulants:
193  void SetEbECumulantsList(TList* const ebecl) {this->fEbECumulantsList = ebecl;};
194  TList* GetEbECumulantsList() const {return this->fEbECumulantsList;}
195  void SetEbECumulantsFlagsPro(TProfile* const ebecfp) {this->fEbECumulantsFlagsPro = ebecfp;};
196  TProfile* GetEbECumulantsFlagsPro() const {return this->fEbECumulantsFlagsPro;};
199 
200  // 5.5.) Weights:
201  void SetWeightsList(TList* const wl) {this->fWeightsList = (TList*)wl->Clone();};
202  TList* GetWeightsList() const {return this->fWeightsList;}
203  void SetWeightsFlagsPro(TProfile* const wfp) {this->fWeightsFlagsPro = wfp;};
204  TProfile* GetWeightsFlagsPro() const {return this->fWeightsFlagsPro;};
205  void SetWeightsHist(TH1D* const hist, const char *type, const char *variable); // .cxx
206 
207  // 5.6.) Nested loops:
208  void SetNestedLoopsList(TList* const nll) {this->fNestedLoopsList = nll;}
209  TList* GetNestedLoopsList() const {return this->fNestedLoopsList;}
210  void SetNestedLoopsFlagsPro(TProfile* const nlfp) {this->fNestedLoopsFlagsPro = nlfp;};
211  TProfile* GetNestedLoopsFlagsPro() const {return this->fNestedLoopsFlagsPro;};
217  {
218  this->fCrossCheckDiffCSCOBN[0] = cs; // cos/sin
219  this->fCrossCheckDiffCSCOBN[1] = co; // correlator order [1p,2p,3p,4p]
220  this->fCrossCheckDiffCSCOBN[2] = bn; // bin number
221  };
222 
223  // 5.7.) 'Standard candles':
224  void SetStandardCandlesList(TList* const scl) {this->fStandardCandlesList = scl;}
226  void SetStandardCandlesFlagsPro(TProfile* const scfp) {this->fStandardCandlesFlagsPro = scfp;};
227  TProfile* GetStandardCandlesFlagsPro() const {return this->fStandardCandlesFlagsPro;};
230  void SetPropagateErrorSC(Bool_t pesc) {this->fPropagateErrorSC = pesc;};
232  void SetStandardCandlesHist(TH1D* const sch) {this->fStandardCandlesHist = sch;};
234  void SetProductsSCPro(TProfile2D* const psc) {this->fProductsSCPro = psc;};
235  TProfile2D* GetProductsSCPro() const {return this->fProductsSCPro;};
236 
237  // 5.8.) Q-cumulants:
238  void SetQcumulantsList(TList* const qcl) {this->fQcumulantsList = qcl;};
239  TList* GetQcumulantsList() const {return this->fQcumulantsList;}
240  void SetQcumulantsFlagsPro(TProfile* const qcfp) {this->fQcumulantsFlagsPro = qcfp;};
241  TProfile* GetQcumulantsFlagsPro() const {return this->fQcumulantsFlagsPro;};
244  void SetHarmonicQC(Int_t hqc) {this->fHarmonicQC = hqc;};
245  Int_t GetHarmonicQC() const {return this->fHarmonicQC;};
246  void SetPropagateErrorQC(Bool_t peqc) {this->fPropagateErrorQC = peqc;};
248  void SetQcumulantsHist(TH1D* const qch) {this->fQcumulantsHist = qch;};
249  TH1D* GetQcumulantsHist() const {return this->fQcumulantsHist;};
250  void SetReferenceFlowHist(TH1D* const rfh) {this->fReferenceFlowHist = rfh;};
251  TH1D* GetReferenceFlowHist() const {return this->fReferenceFlowHist;};
252  void SetProductsQCPro(TProfile2D* const pqc) {this->fProductsQCPro = pqc;};
253  TProfile2D* GetProductsQCPro() const {return this->fProductsQCPro;};
254 
255  // 5.9.) Differential correlations:
256  void SetDiffCorrelationsList(TList* const dcl) {this->fDiffCorrelationsList = dcl;};
258  void SetDiffCorrelationsFlagsPro(TProfile* const cdfp) {this->fDiffCorrelationsFlagsPro = cdfp;};
259  TProfile* GetDiffCorrelationsFlagsPro() const {return this->fDiffCorrelationsFlagsPro;};
262  void SetDiffHarmonics(Int_t order, Int_t *harmonics); // see implementation in .cxx file
263  void SetCalculateDiffCos(Bool_t cdc) {this->fCalculateDiffCos = cdc;};
265  void SetCalculateDiffSin(Bool_t cds) {this->fCalculateDiffSin = cds;};
271  void SetnDiffBins(Int_t ndb) {this->fnDiffBins = ndb;};
272  Int_t GetnDiffBins() const {return this->fnDiffBins;};
273  void SetRangesDiffBins(Double_t* const rdb) {this->fRangesDiffBins = rdb;};
274  Double_t* GetRangesDiffBins() const {return this->fRangesDiffBins;};
275 
276  // 5.10.) Symmetry plane correlations:
277  void SetSymmetryPlanesList(TList* const spl) {this->fSymmetryPlanesList = spl;};
279  void SetSymmetryPlanesFlagsPro(TProfile* const spfp) {this->fSymmetryPlanesFlagsPro = spfp;};
280  TProfile* GetSymmetryPlanesFlagsPro() const {return this->fSymmetryPlanesFlagsPro;};
283 
284  // 5.11.) Eta gaps:
285  void SetEtaGapsList(TList* const egl) {this->fEtaGapsList = egl;};
286  TList* GetEtaGapsList() const {return this->fEtaGapsList;}
287  void SetEtaGapsFlagsPro(TProfile* const egfp) {this->fEtaGapsFlagsPro = egfp;};
288  TProfile* GetEtaGapsFlagsPro() const {return this->fEtaGapsFlagsPro;};
289  void SetCalculateEtaGaps(Bool_t ceg) {this->fCalculateEtaGaps = ceg;};
295 
296  // 6.) The rest:
297  virtual void WriteHistograms(TString outputFileName);
298  virtual void WriteHistograms(TDirectoryFile *outputFileName);
299  virtual TComplex Q(Int_t n, Int_t p);
300  virtual TComplex p(Int_t n, Int_t p);
301  virtual TComplex q(Int_t n, Int_t p);
302  virtual TComplex One(Int_t n1);
303  virtual TComplex Two(Int_t n1, Int_t n2);
304  virtual TComplex Three(Int_t n1, Int_t n2, Int_t n3);
305  virtual TComplex Four(Int_t n1, Int_t n2, Int_t n3, Int_t n4);
306  virtual TComplex Five(Int_t n1, Int_t n2, Int_t n3, Int_t n4, Int_t n5);
307  virtual TComplex Six(Int_t n1, Int_t n2, Int_t n3, Int_t n4, Int_t n5, Int_t n6);
308  virtual TComplex Seven(Int_t n1, Int_t n2, Int_t n3, Int_t n4, Int_t n5, Int_t n6, Int_t n7);
309  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);
310  virtual TComplex OneDiff(Int_t n1);
311  virtual TComplex TwoDiff(Int_t n1, Int_t n2);
312  virtual TComplex ThreeDiff(Int_t n1, Int_t n2, Int_t n3);
313  virtual TComplex FourDiff(Int_t n1, Int_t n2, Int_t n3, Int_t n4);
314  virtual Double_t Weight(const Double_t &value, const char *type, const char *variable); // value, [RP,POI], [phi,pt,eta]
315  virtual Double_t CastStringToCorrelation(const char *string, Bool_t numerator);
316  virtual Double_t Covariance(const char *x, const char *y, TProfile2D *profile2D, Bool_t bUnbiasedEstimator = kFALSE);
317  virtual TComplex Recursion(Int_t n, Int_t* harmonic, Int_t mult = 1, Int_t skip = 0); // Credits: Kristjan Gulbrandsen (gulbrand@nbi.dk)
318  virtual void CalculateProductsOfCorrelations(AliFlowEventSimple *anEvent, TProfile2D *profile2D);
319  static void DumpPointsForDurham(TGraphErrors *ge);
320  static void DumpPointsForDurham(TH1D *h);
321  static void DumpPointsForDurham(TH1F *h);
322  virtual void DumpThePoints(AliFlowEventSimple *anEvent);
323  TH1D* GetHistogramWithWeights(const char *filePath, const char *listName, const char *type, const char *variable, const char *production);
324  virtual Double_t CorrelationPsi2nPsi1n(Int_t n, Int_t k=0);
326 
327  private:
330  // Data members are grouped as:
331  // 0.) Base list and internal flags;
332  // 1.) Control histograms;
333  // 2.) Q-vectors;
334  // 3.) Correlations;
335  // 4.) Event-by-event cumulants;
336  // 5.) Weights;
337  // 6.) Nested loops;
338  // 7.) 'Standard candles';
339  // 8.) Q-cumulants;
340  // 9.) Differential correlations;
341  // 10.) Symmetry plane correlations;
342  // 11.) Eta gaps.
343 
344  // 0.) Base list and internal flags:
345  TList* fHistList; // base list to hold all output object (a.k.a. grandmother of all lists)
346  TProfile *fInternalFlagsPro; // profile to hold all internal flags and settings
347  Bool_t fUseInternalFlags; // use internal flags (automatically set if some internal flag is used)
348  Int_t fMinNoRPs; // minimum number of RPs required for the analysis
349  Int_t fMaxNoRPs; // maximum number of RPs allowed for the analysis
350  Int_t fExactNoRPs; // exact (randomly shuffled) number of RPs selected for the analysis
351  Bool_t fPropagateError; // prevent error propagation if something strange happens during calculations
352  TString fAnalysisTag; // tag internally this analysis
353  Bool_t fDumpThePoints; // dump the data points into the external file
354  Int_t fMaxNoEventsPerFile; // maximum number of events to be dumped in a single file
355  Bool_t fSelectRandomlyRPs; // enable random shuffling to estimate 'fake flow'
356  Int_t fnSelectedRandomlyRPs; // how many RPs will be taken for the analysis after random shuffling?
357  TArrayI *fRandomIndicesRPs; // well, these are random indices...
358 
359  // 1.) Control histograms:
360  TList *fControlHistogramsList; // list to hold all 'control histograms' objects
361  TProfile *fControlHistogramsFlagsPro; // profile to hold all flags for control histograms
362  Bool_t fFillControlHistograms; // fill or not control histograms (by default they are all filled)
363  Bool_t fFillKinematicsHist; // fill or not fKinematicsHist[2][3]
364  Bool_t fFillMultDistributionsHist; // fill or not TH1D *fMultDistributionsHist[3]
365  Bool_t fFillMultCorrelationsHist; // fill or not TH2D *fMultCorrelationsHist[3]
366  Bool_t fDontFill[3]; // don't fill control histograms [0=RP,1=POI,2=REF]
367  TH1D *fKinematicsHist[2][3]; // [RP,POI][phi,pt,eta] distributions
368  TH1D *fMultDistributionsHist[3]; // multiplicity distribution [RP,POI,reference multiplicity]
369  TH2I *fMultCorrelationsHist[3]; // [RP vs. POI, RP vs. refMult, POI vs. refMult]
370  Int_t fnBins[2][3]; // [RP,POI][phi,pt,eta], corresponds to fKinematicsHist[2][3]
371  Double_t fMin[2][3]; // [RP,POI][phi,pt,eta], corresponds to fKinematicsHist[2][3]
372  Double_t fMax[2][3]; // [RP,POI][phi,pt,eta], corresponds to fKinematicsHist[2][3]
373  Int_t fnBinsMult[3]; // [RP,POI,REF], corresponds to fMultDistributionsHist[3]
374  Double_t fMinMult[3]; // [RP,POI,REF], corresponds to fMultDistributionsHist[3]
375  Double_t fMaxMult[3]; // [RP,POI,REF], corresponds to fMultDistributionsHist[3]
376  Bool_t fSkipSomeIntervals; // skip intervals in phi, pt and eta
378  Double_t fSkip[3][10]; // determine intervals in phi, pt and eta to be skipped. TBI hardwired is max 5 intervals. TBI promote this eventually to AFTC class
379 
380  // 2.) Q-vectors:
381  TList *fQvectorList; // list to hold all Q-vector objects
382  TProfile *fQvectorFlagsPro; // profile to hold all flags for Q-vector
383  Bool_t fCalculateQvector; // to calculate or not to calculate Q-vector components, that's a Boolean...
384  TComplex fQvector[49][9]; // Q-vector components [fMaxHarmonic*fMaxCorrelator+1][fMaxCorrelator+1] = [6*8+1][8+1]
385  Bool_t fCalculateDiffQvectors; // to calculate or not to calculate p- and q-vector components, that's a Boolean...
386  TComplex fpvector[100][49][9]; // p-vector components [bin][fMaxHarmonic*fMaxCorrelator+1][fMaxCorrelator+1] = [6*8+1][8+1] TBI hardwired 100
387  TComplex fqvector[100][49][9]; // q-vector components [bin][fMaxHarmonic*fMaxCorrelator+1][fMaxCorrelator+1] = [6*8+1][8+1] TBI hardwired 100
388 
389  // 3.) Correlations:
390  TList *fCorrelationsList; // list to hold all correlations objects
391  TProfile *fCorrelationsFlagsPro; // profile to hold all flags for correlations
392  TProfile *fCorrelationsPro[2][8]; // multi-particle correlations [0=cos,1=sin][1p,2p,...,8p]
393  Bool_t fCalculateCorrelations; // calculate and store correlations
394  Int_t fMaxHarmonic; // 6 (not going beyond v6, if you change this value, change also fQvector[49][9])
395  Int_t fMaxCorrelator; // 8 (not going beyond 8-p correlations, if you change this value, change also fQvector[49][9])
396  Bool_t fCalculateIsotropic; // calculate only isotropic correlations
397  Bool_t fCalculateSame; // calculate only 'same abs harmonics' correlations TBI
398  Bool_t fSkipZeroHarmonics; // skip correlations which have some of the harmonicc equal to zero
399  Bool_t fCalculateSameIsotropic; // calculate all isotropic correlations in 'same abs harmonic' TBI this can be implemented better
400  Bool_t fCalculateAll; // calculate all possible correlations
401  Int_t fDontGoBeyond; // do not go beyond fDontGoBeyond-p correlators
402  Bool_t fCalculateOnlyForHarmonicQC; // calculate only isotropic correlations in |fHarmonicQC|
403  Bool_t fCalculateOnlyForSC; // calculate only correlations needed for 'standard candles'
404  Bool_t fCalculateOnlyCos; // calculate only 'cos' correlations
405  Bool_t fCalculateOnlySin; // calculate only 'sin' correlations
406 
407  // 4.) Event-by-event cumulants:
408  TList *fEbECumulantsList; // list to hold all e-b-e cumulants objects
409  TProfile *fEbECumulantsFlagsPro; // profile to hold all flags for e-b-e cumulants
410  TProfile *fEbECumulantsPro[2][8]; // multi-particle e-b-e cumulants [0=cos,1=sin][1p,2p,...,8p]
411  Bool_t fCalculateEbECumulants; // calculate and store e-b-e cumulants
412 
413  // 5.) Weights:
414  TList *fWeightsList; // list to hold all weights objects
415  TProfile *fWeightsFlagsPro; // profile to hold all flags for weights
416  Bool_t fUseWeights[2][3]; // use weights [RP,POI][phi,pt,eta]
417  TH1D *fWeightsHist[2][3]; // histograms holding weights [RP,POI][phi,pt,eta]
418 
419  // 6.) Nested loops:
420  TList *fNestedLoopsList; // list to hold all nested loops objects
421  TProfile *fNestedLoopsFlagsPro; // profile to hold all flags for nested loops
422  Bool_t fCrossCheckWithNestedLoops; // cross-check results with nested loops
423  Bool_t fCrossCheckDiffWithNestedLoops; // cross-check differential correlators with nested loops
424  Int_t fCrossCheckDiffCSCOBN[3]; // [0=cos,1=sin][1p,2p,...,4p][binNo]
425  TProfile *fNestedLoopsResultsCosPro; // profile to hold nested loops results (cosine)
426  TProfile *fNestedLoopsResultsSinPro; // profile to hold nested loops results (sinus)
427  TProfile *fNestedLoopsDiffResultsPro; // profile to hold differential nested loops results // TBI
428 
429  // 7.) 'Standard candles':
430  TList *fStandardCandlesList; // list to hold all 'standard candles' objects
431  TProfile *fStandardCandlesFlagsPro; // profile to hold all flags fo 'standard candles'
432  Bool_t fCalculateStandardCandles; // calculate and store 'standard candles'
433  Bool_t fPropagateErrorSC; // propagate and store error for 'standard candles'
434  TH1D *fStandardCandlesHist; // histogram to hold results for 'standard candles'
435  TProfile2D *fProductsSCPro; // 2D profile to hold products of correlations needed for SC error propagation
436 
437  // 8.) Q-cumulants:
438  TList *fQcumulantsList; // list to hold all Q-cumulants objects
439  TProfile *fQcumulantsFlagsPro; // profile to hold all flags for Q-cumulants
440  Bool_t fCalculateQcumulants; // calculate and store Q-cumulants
441  Int_t fHarmonicQC; // calculate Q-cumulants in this harmonic (default is 2)
442  Bool_t fPropagateErrorQC; // propagate and store error for Q-cumulants
443  TH1D *fQcumulantsHist; // two- and multi-particle Q-cumulants
444  TH1D *fReferenceFlowHist; // reference flow from two- and multi-particle Q-cumulants
445  TProfile2D *fProductsQCPro; // 2D profile to hold products of correlations needed for QC error propagation
446 
447  // 9.) Differential correlations:
448  TList *fDiffCorrelationsList; // list to hold all correlations objects
449  TProfile *fDiffCorrelationsFlagsPro; // profile to hold all flags for correlations
450  Bool_t fCalculateDiffCorrelations; // calculate and store differential correlations
451  Bool_t fCalculateDiffCos; // calculate and store differential cosine correlations (kTRUE by default)
452  Bool_t fCalculateDiffSin; // calculate and store differential sinus correlations (kFALSE by default)
453  Bool_t fCalculateDiffCorrelationsVsPt; // calculate differential correlations vs pt (default), or vs eta
454  Bool_t fUseDefaultBinning; // use default binning in pt or in eta
455  Int_t fnDiffBins; // number of differential bins in pt or in eta (when non-default binning is used)
456  Double_t *fRangesDiffBins; // ranges for differential bins in pt or in eta (when non-default binning is used)
457  Int_t fDiffHarmonics[4][4]; // harmonics for differential correlations [order][{n1},{n1,n2},...,{n1,n2,n3,n4}]
458  TProfile *fDiffCorrelationsPro[2][4]; // multi-particle correlations [0=cos,1=sin][1p,2p,3p,4p]
459  UInt_t fDiffBinNo; // differential bin number
460 
461  // 10.) Symmetry plane correlations:
462  TList *fSymmetryPlanesList; // list to hold all correlations between symmetry planes
463  TProfile *fSymmetryPlanesFlagsPro; // profile to hold all flags for symmetry plane correlations
464  Bool_t fCalculateSymmetryPlanes; // calculate correlations between symmetry planes
465  //Bool_t ...
466  TProfile *fSymmetryPlanesPro[1][2]; // symmetry planes correlations [[0]:(Psi2n,Psi1n),[1]:...][[0]:n=1,[1]:n=2,...]
467  //Int_t fnHighestCorrelatorSPC; // highest generic correlator to be calculated TBI implement this more differentially
468  //Int_t fnHighestHarmonicSPC; // highest harmonic for evaluation of generic correlators TBI implement this more differentially
469  //Int_t fnHighestOptimizerSPC; // highest optimizer TBI implement this more differentially
470 
471  // 11.) Eta gaps:
472  TList *fEtaGapsList; // list to hold all correlations with eta gaps
473  TProfile *fEtaGapsFlagsPro; // profile to hold all flags for correlations with eta gaps
474  Bool_t fCalculateEtaGaps; // calculate correlations with eta gaps
475  Int_t fLowestHarmonicEtaGaps; // 2-p correlations with eta gaps will be calculated for harmonics [fLowestHarmonicEtaGaps,fHighestHarmonicEtaGaps]
476  Int_t fHighestHarmonicEtaGaps; // 2-p correlations with eta gaps will be calculated for harmonics [fLowestHarmonicEtaGaps,fHighestHarmonicEtaGaps]
477  TProfile *fEtaGapsPro[6]; // [harmonic] different eta gaps are different bins
478 
480 
481 };
482 
483 //================================================================================================================
484 
485 #endif
486 
487 
488 
489 
490 
virtual TComplex Four(Int_t n1, Int_t n2, Int_t n3, Int_t n4)
virtual Double_t CastStringToCorrelation(const char *string, Bool_t numerator)
const Color_t cc[]
Definition: DrawKs.C:1
double Double_t
Definition: External.C:58
virtual TComplex FourDiff(Int_t n1, Int_t n2, Int_t n3, Int_t n4)
void SetMax(const char *type, const char *variable, Double_t max)
void SetMin(const char *type, const char *variable, Double_t min)
virtual TComplex Seven(Int_t n1, Int_t n2, Int_t n3, Int_t n4, Int_t n5, Int_t n6, Int_t n7)
virtual TComplex Six(Int_t n1, Int_t n2, Int_t n3, Int_t n4, Int_t n5, Int_t n6)
TH1D * GetHistogramWithWeights(const char *filePath, const char *listName, const char *type, const char *variable, const char *production)
void SetWeightsHist(TH1D *const hist, const char *type, const char *variable)
virtual void CalculateProductsOfCorrelations(AliFlowEventSimple *anEvent, TProfile2D *profile2D)
virtual TComplex Five(Int_t n1, Int_t n2, Int_t n3, Int_t n4, Int_t n5)
int Int_t
Definition: External.C:63
void SetnBins(const char *type, const char *variable, Int_t nBins)
unsigned int UInt_t
Definition: External.C:33
AliFlowAnalysisWithMultiparticleCorrelations & operator=(const AliFlowAnalysisWithMultiparticleCorrelations &afawQc)
Definition: External.C:212
virtual TComplex Recursion(Int_t n, Int_t *harmonic, Int_t mult=1, Int_t skip=0)
virtual Double_t Weight(const Double_t &value, const char *type, const char *variable)
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)
virtual Double_t Covariance(const char *x, const char *y, TProfile2D *profile2D, Bool_t bUnbiasedEstimator=kFALSE)
void SetIntervalsToSkip(const char *ppe, Int_t n, Double_t *boundaries)
bool Bool_t
Definition: External.C:53
Double_t maxMult