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