]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/FLOW/Tasks/AliAnalysisTaskMultiparticleCorrelations.h
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWG / FLOW / Tasks / AliAnalysisTaskMultiparticleCorrelations.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  * analysis task for flow analysis with *
9  *     multi-particle correlations      * 
10  *                                      * 
11  * author: Ante Bilandzic               *
12  *         (abilandzic@gmail.com)       * 
13  ***************************************/
14
15 #ifndef ALIANALYSISTASKMULTIPARTICLECORRELATIONS_H
16 #define ALIANALYSISTASKMULTIPARTICLECORRELATIONS_H
17
18 #include "AliAnalysisTaskSE.h"
19 #include "AliFlowAnalysisWithMultiparticleCorrelations.h"
20 #include "AliFlowEventSimple.h"
21
22 //================================================================================================================
23
24 class AliAnalysisTaskMultiparticleCorrelations : public AliAnalysisTaskSE{
25  public:
26   AliAnalysisTaskMultiparticleCorrelations();
27   AliAnalysisTaskMultiparticleCorrelations(const char *name, Bool_t useParticleWeights=kFALSE);
28   virtual ~AliAnalysisTaskMultiparticleCorrelations(){}; 
29   
30   virtual void UserCreateOutputObjects();
31   virtual void UserExec(Option_t *option);
32   virtual void Terminate(Option_t *);
33   
34   // Internal flags:
35   void SetMinNoRPs(Int_t min) {fUseInternalFlags = kTRUE; this->fMinNoRPs = min;};
36   Int_t GetMinNoRPs() const {return this->fMinNoRPs;};
37   void SetMaxNoRPs(Int_t max) {fUseInternalFlags = kTRUE; this->fMaxNoRPs = max;};
38   Int_t GetMaxNoRPs() const {return this->fMaxNoRPs;};
39   void SetExactNoRPs(Int_t exact) {fUseInternalFlags = kTRUE; this->fExactNoRPs = exact;};
40   Int_t GetExactNoRPs() const {return this->fExactNoRPs;};
41   void SetAnalysisTag(const char *at) {this->fAnalysisTag = TString(at);};
42   TString GetAnalysisTag() const {return this->fAnalysisTag;};
43   void SetDumpThePoints(Bool_t dtp, Int_t max) {this->fDumpThePoints = dtp; this->fMaxNoEventsPerFile = max;};
44
45   // Control histograms:
46   void SetFillControlHistograms(Bool_t fch) {this->fFillControlHistograms = fch;};
47   Bool_t GetFillControlHistograms() const {return this->fFillControlHistograms;};
48   void SetFillKinematicsHist(Bool_t fkh) {this->fFillKinematicsHist = fkh;};
49   Bool_t GetFillKinematicsHist() const {return this->fFillKinematicsHist;};
50   void SetFillMultDistributionsHist(Bool_t mdh) {this->fFillMultDistributionsHist = mdh;};
51   Bool_t GetFillMultDistributionsHist() const {return this->fFillMultDistributionsHist;};
52   void SetFillMultCorrelationsHist(Bool_t mch) {this->fFillMultCorrelationsHist = mch;};
53   Bool_t GetFillMultCorrelationsHist() const {return this->fFillMultCorrelationsHist;};
54   void SetDontFill(const char *type) 
55   {
56    if(TString(type).EqualTo("RP")){this->fDontFill[0] = kTRUE;} 
57    else if(TString(type).EqualTo("POI")){this->fDontFill[1] = kTRUE;} 
58    else if(TString(type).EqualTo("REF")){this->fDontFill[2] = kTRUE;} 
59    else{Fatal("void SetDontFill(const char *type)","type = %s ???? Allowed: RP, POI and REF.",type);}
60   }; // void SetDontFill(const char *type)
61   void SetnBins(const char *type, const char *variable, Int_t nBins); // .cxx
62   void SetMin(const char *type, const char *variable, Double_t min); // .cxx
63   void SetMax(const char *type, const char *variable, Double_t max); // .cxx
64   void SetnBinsMult(const char *type, Int_t nBinsMult); // .cxx
65   void SetMinMult(const char *type, Double_t minMult); // .cxx
66   void SetMaxMult(const char *type, Double_t maxMult); // .cxx
67
68   // Q-vectors:
69   void SetCalculateQvector(Bool_t cqv) {this->fCalculateQvector = cqv;};
70   Bool_t GetCalculateQvector() const {return this->fCalculateQvector;};
71   void SetCalculateDiffQvectors(Bool_t cdqv) {this->fCalculateDiffQvectors = cdqv;};
72   Bool_t GetCalculateDiffQvectors() const {return this->fCalculateDiffQvectors;};
73
74   // Weights:              
75   void SetWeightsHist(TH1D* const hist, const char *type, const char *variable); // .cxx
76   TH1D* GetHistogramWithWeights(const char *filePath, const char *listName, const char *type, const char *variable, const char *production)
77   {
78    fProduction = TString(production);
79    AliFlowAnalysisWithMultiparticleCorrelations *mpc = new AliFlowAnalysisWithMultiparticleCorrelations();
80    return mpc->GetHistogramWithWeights(filePath,listName,type,variable,production);
81   };
82
83   // Correlations:
84   void SetCalculateCorrelations(Bool_t cc) {this->fCalculateCorrelations = cc;};
85   Bool_t GetCalculateCorrelations() const {return this->fCalculateCorrelations;};
86   void SetCalculateIsotropic(Bool_t ci) {this->fCalculateIsotropic = ci;};
87   Bool_t GetCalculateIsotropic() const {return this->fCalculateIsotropic;};
88   void SetCalculateSame(Bool_t csh) {this->fCalculateSame = csh;};
89   Bool_t GetCalculateSame() const {return this->fCalculateSame;};
90   void SetSkipZeroHarmonics(Bool_t szh) {this->fSkipZeroHarmonics = szh;};
91   Bool_t GetSkipZeroHarmonics() const {return this->fSkipZeroHarmonics;};
92   void SetCalculateSameIsotropic(Bool_t csi) {this->fCalculateSameIsotropic = csi;};
93   Bool_t GetCalculateSameIsotropic() const {return this->fCalculateSameIsotropic;};
94   void SetCalculateAll(Bool_t ca) {this->fCalculateAll = ca;};
95   Bool_t GetCalculateAll() const {return this->fCalculateAll;};
96   void SetDontGoBeyond(Int_t dgb) {this->fDontGoBeyond = dgb;};
97   Int_t GetDontGoBeyond() const {return this->fDontGoBeyond;};
98   void SetCalculateOnlyForHarmonicQC(Bool_t cofhqc) {this->fCalculateOnlyForHarmonicQC = cofhqc;};
99   Bool_t GetCalculateOnlyForHarmonicQC() const {return this->fCalculateOnlyForHarmonicQC;};
100   void SetCalculateOnlyForSC(Bool_t cofsc) {this->fCalculateOnlyForSC = cofsc;};
101   Bool_t GetCalculateOnlyForSC() const {return this->fCalculateOnlyForSC;};
102   void SetCalculateOnlyCos(Bool_t coc) {this->fCalculateOnlyCos = coc;};
103   Bool_t GetCalculateOnlyCos() const {return this->fCalculateOnlyCos;};
104   void SetCalculateOnlySin(Bool_t cos) {this->fCalculateOnlySin = cos;};
105   Bool_t GetCalculateOnlySin() const {return this->fCalculateOnlySin;};
106
107   // Event-by-event cumulants:
108   void SetCalculateEbECumulants(Bool_t cebec) {this->fCalculateEbECumulants = cebec;};
109   Bool_t GetCalculateEbECumulants() const {return this->fCalculateEbECumulants;};
110
111   // Nested loops:
112   void SetCrossCheckWithNestedLoops(Bool_t ccwnl) {this->fCrossCheckWithNestedLoops = ccwnl;};
113   Bool_t GetCrossCheckWithNestedLoops() const {return this->fCrossCheckWithNestedLoops;};
114   void SetCrossCheckDiffWithNestedLoops(Bool_t ccdwnl) {this->fCrossCheckDiffWithNestedLoops = ccdwnl;};
115   Bool_t GetCrossCheckDiffWithNestedLoops() const {return this->fCrossCheckDiffWithNestedLoops;};
116   void SetCrossCheckDiffCSCOBN(Int_t cs, Int_t co, Int_t bn)  
117   {
118    this->fCrossCheckDiffCSCOBN[0] = cs; // cos/sin
119    this->fCrossCheckDiffCSCOBN[1] = co; // correlator order [1p,2p,3p,4p]
120    this->fCrossCheckDiffCSCOBN[2] = bn; // bin number
121   };
122
123   // 'Standard candles':
124   void SetCalculateStandardCandles(Bool_t csc) {this->fCalculateStandardCandles = csc;};
125   Bool_t GetCalculateStandardCandles() const {return this->fCalculateStandardCandles;};
126   void SetPropagateErrorSC(Bool_t pesc) {this->fPropagateErrorSC = pesc;};
127   Bool_t GetPropagateErrorSC() const {return this->fPropagateErrorSC;};
128
129   // Q-cumulants:
130   void SetCalculateQcumulants(Bool_t cqc) {this->fCalculateQcumulants = cqc;};
131   Bool_t GetCalculateQcumulants() const {return this->fCalculateQcumulants;};
132   void SetHarmonicQC(Int_t hqc) {this->fHarmonicQC = hqc;};
133   Int_t GetHarmonicQC() const {return this->fHarmonicQC;};
134   void SetPropagateErrorQC(Bool_t peqc) {this->fPropagateErrorQC = peqc;};
135   Bool_t GetPropagateErrorQC() const {return this->fPropagateErrorQC;};
136
137   // Differential correlations:
138   void SetCalculateDiffCorrelations(Bool_t cdc) {this->fCalculateDiffCorrelations = cdc;};
139   Bool_t GetCalculateDiffCorrelations() const {return this->fCalculateDiffCorrelations;};
140   void SetDiffHarmonics(Int_t order, Int_t *harmonics){;} // not implemented 
141   void SetCalculateDiffCos(Bool_t cdc) {this->fCalculateDiffCos = cdc;};
142   Bool_t GetCalculateDiffCos() const {return this->fCalculateDiffCos;};
143   void SetCalculateDiffSin(Bool_t cds) {this->fCalculateDiffSin = cds;};
144   Bool_t GetCalculateDiffSin() const {return this->fCalculateDiffSin;};
145   void SetCalculateDiffCorrelationsVsPt(Bool_t cdcvspt) {this->fCalculateDiffCorrelationsVsPt = cdcvspt;};
146   Bool_t GetCalculateDiffCorrelationsVsPt() const {return this->fCalculateDiffCorrelationsVsPt;};
147   void SetUseDefaultBinning(Bool_t udb) {this->fUseDefaultBinning = udb;};
148   Bool_t GetUseDefaultBinning() const {return this->fUseDefaultBinning;};
149   void SetnDiffBins(Int_t ndb) {this->fnDiffBins = ndb;};
150   Int_t GetnDiffBins() const {return this->fnDiffBins;};
151   void SetRangesDiffBins(Double_t* const rdb) {this->fRangesDiffBins = rdb;};
152   Double_t* GetRangesDiffBins() const {return this->fRangesDiffBins;};
153
154  private:
155   AliAnalysisTaskMultiparticleCorrelations(const AliAnalysisTaskMultiparticleCorrelations& aatqc);
156   AliAnalysisTaskMultiparticleCorrelations& operator=(const AliAnalysisTaskMultiparticleCorrelations& aatqc);
157   
158   AliFlowEventSimple *fEvent; // the input event
159   AliFlowAnalysisWithMultiparticleCorrelations *fMPC; // "multi-particle correlations" object
160   TList *fHistList; // base list to hold all output object (a.k.a. grandmother of all lists)
161
162   // Internal flags:
163   Bool_t fUseInternalFlags;  // use internal flags (automatically set if some internal flag is used)
164   Int_t fMinNoRPs;           // minimum number of RPs required for the analysis 
165   Int_t fMaxNoRPs;           // maximum number of RPs allowed for the analysis 
166   Int_t fExactNoRPs;         // exact (randomly shuffled) number of RPs selected for the analysis 
167   TString fAnalysisTag;      // tag internally this analysis
168   Bool_t fDumpThePoints;     // dump the data points into the external file
169   Int_t fMaxNoEventsPerFile; // maximum number of events to be dumped in a single file
170
171   // Control histograms:
172   Bool_t fFillControlHistograms;     // fill or not control histograms (by default they are filled)
173   Bool_t fFillKinematicsHist;        // fill or not fKinematicsHist[2][3]
174   Bool_t fFillMultDistributionsHist; // fill or not TH1D *fMultDistributionsHist[3]    
175   Bool_t fFillMultCorrelationsHist;  // fill or not TH2D *fMultCorrelationsHist[3] 
176   Bool_t fDontFill[3];               // don't fill control histograms [0=RP,1=POI,2=REF]  
177   Int_t fnBins[2][3];                // [RP,POI][phi,pt,eta], corresponds to fKinematicsHist[2][3]
178   Double_t fMin[2][3];               // [RP,POI][phi,pt,eta], corresponds to fKinematicsHist[2][3]
179   Double_t fMax[2][3];               // [RP,POI][phi,pt,eta], corresponds to fKinematicsHist[2][3]
180   Int_t fnBinsMult[3];               // [RP,POI,REF], corresponds to fMultDistributionsHist[3]   
181   Double_t fMinMult[3];              // [RP,POI,REF], corresponds to fMultDistributionsHist[3]   
182   Double_t fMaxMult[3];              // [RP,POI,REF], corresponds to fMultDistributionsHist[3]  
183
184   // Q-vectors:
185   Bool_t fCalculateQvector;      // to calculate or not to calculate Q-vector components, that's a Boolean...
186   Bool_t fCalculateDiffQvectors; // to calculate or not to calculate p- and q-vector components, that's a Boolean...
187
188   // Weights:
189   Bool_t fUseWeights[2][3]; // use weights [RP,POI][phi,pt,eta]
190   TH1D *fWeightsHist[2][3]; // histograms holding weights [RP,POI][phi,pt,eta]
191   TString fProduction;      // TBI most likely an overkill
192
193   // Correlations:
194   Bool_t fCalculateCorrelations;      // calculate and store correlations, or perhaps not, if the weather is bad...
195   Bool_t fCalculateIsotropic;         // calculate only isotropic correlations
196   Bool_t fCalculateSame;              // calculate only 'same abs harmonics' correlations TBI 
197   Bool_t fSkipZeroHarmonics;          // skip correlations which have some of the harmonicc equal to zero
198   Bool_t fCalculateSameIsotropic;     // calculate all isotropic correlations in 'same abs harmonic' TBI this can be implemented better
199   Bool_t fCalculateAll;               // calculate all possible correlations 
200   Int_t  fDontGoBeyond;               // do not go beyond fDontGoBeyond-p correlators
201   Bool_t fCalculateOnlyForHarmonicQC; // calculate only isotropic correlations in |fHarmonicQC|
202   Bool_t fCalculateOnlyForSC;         // calculate only correlations needed for 'standard candles'
203   Bool_t fCalculateOnlyCos;           // calculate only 'cos' correlations
204   Bool_t fCalculateOnlySin;           // calculate only 'sin' correlations
205
206   // Event-by-event cumulants:
207   Bool_t fCalculateEbECumulants; // calculate and store event-by-event cumulants
208
209   // Nested loops:
210   Bool_t fCrossCheckWithNestedLoops;     // cross-check results with nested loops
211   Bool_t fCrossCheckDiffWithNestedLoops; // cross-check differential correlators with nested loops
212   Int_t fCrossCheckDiffCSCOBN[3];        // [0=cos,1=sin][1p,2p,...,4p][binNo]
213
214   // 'Standard candles':
215   Bool_t fCalculateStandardCandles; // calculate and store 'standard candles'
216   Bool_t fPropagateErrorSC;         // propagate and store error for 'standard candles'
217
218   // Q-cumulants:
219   Bool_t fCalculateQcumulants; // calculate and store Q-cumulants
220   Int_t fHarmonicQC;           // calculate Q-cumulants in this harmonic (default is 2) 
221   Bool_t fPropagateErrorQC;    // propagate and store error for Q-cumulants
222
223   // Differential correlations:
224   Bool_t fCalculateDiffCorrelations;     // calculate and store differential correlations
225   Bool_t fCalculateDiffCos;              // calculate and store differential cosine correlations (kTRUE by default)
226   Bool_t fCalculateDiffSin;              // calculate and store differential sinus correlations (kFALSE by default)
227   Bool_t fCalculateDiffCorrelationsVsPt; // calculate differential correlations vs pt (default), or vs eta
228   Bool_t fUseDefaultBinning;             // use default binning in pt or in eta
229   Int_t fnDiffBins;                      // number of differential bins in pt or in eta (when non-default binning is used)
230   Double_t *fRangesDiffBins;             // ranges for differential bins in pt or in eta (when non-default binning is used)
231
232   ClassDef(AliAnalysisTaskMultiparticleCorrelations,2); 
233
234 };
235
236 //================================================================================================================
237
238 #endif
239
240
241
242
243
244
245
246
247
248
249