]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FLOW/AliFlowCommon/AliFlowAnalysisWithCumulants.cxx
prevent adding histograms to directory to enable multiple instances of the same class
[u/mrichter/AliRoot.git] / PWG2 / FLOW / AliFlowCommon / AliFlowAnalysisWithCumulants.cxx
1 /*************************************************************************
2 * Copyright(c) 1998-2008, ALICE Experiment at CERN, All rights reserved. *
3 *                                                                        *
4 * Author: The ALICE Off-line Project.                                    *
5 * Contributors are mentioned in the code where appropriate.              *
6 *                                                                        *
7 * Permission to use, copy, modify and distribute this software and its   *
8 * documentation strictly for non-commercial purposes is hereby granted   *
9 * without fee, provided that the above copyright notice appears in all   *
10 * copies and that both the copyright notice and this permission notice   *
11 * appear in the supporting documentation. The authors make no claims     *
12 * about the suitability of this software for any purpose. It is          *
13 * provided "as is" without express or implied warranty.                  * 
14 **************************************************************************/
15
16 /******************************** 
17  * flow analysis with cumulants * 
18  *                              *
19  * author: Ante Bilandzic       * 
20  *          (anteb@nikhef.nl)   *
21  *******************************/ 
22
23 #define AliFlowAnalysisWithCumulants_cxx
24
25 #include "Riostream.h"
26 #include "AliFlowCommonConstants.h"
27 #include "AliFlowCommonHist.h"
28 #include "AliFlowCommonHistResults.h"
29 #include "TChain.h"
30 #include "TFile.h"
31 #include "TList.h" //NEW
32 #include "TParticle.h"
33 #include "TRandom3.h"
34 #include "TProfile.h"
35 #include "TProfile2D.h" 
36 #include "TProfile3D.h"
37 #include "AliFlowEventSimple.h"
38 #include "AliFlowTrackSimple.h"
39 #include "AliFlowAnalysisWithCumulants.h"
40 #include "AliFlowCumuConstants.h"
41 #include "AliCumulantsFunctions.h"
42
43 class TH1;
44 class TGraph;
45 class TPave;
46 class TLatex;
47 class TMarker;
48 class TRandom3;
49 class TObjArray;
50 class TList;
51 class TCanvas;
52 class TSystem;
53 class TROOT;
54 class AliFlowVector;
55 class TVector;
56
57 //================================================================================================================
58
59 ClassImp(AliFlowAnalysisWithCumulants)
60
61 AliFlowAnalysisWithCumulants::AliFlowAnalysisWithCumulants():  
62  fTrack(NULL),
63  fHistList(NULL),
64  fWeightsList(NULL),
65  fR0(0),
66  fPtMax(0),
67  fPtMin(0),
68  fBinWidthPt(0),
69  fgknBinsPt(0),
70  fEtaMax(0),
71  fEtaMin(0),
72  fBinWidthEta(0),
73  fgknBinsEta(0),
74  fAvQx(0),
75  fAvQy(0),
76  fAvQ2x(0),
77  fAvQ2y(0),
78  fAvMultIntFlowGFC(NULL),
79  fQVectorComponentsGFC(NULL),
80  fIntFlowResultsGFC(NULL),
81  fDiffFlowResults2ndOrderGFC(NULL),
82  fDiffFlowResults4thOrderGFC(NULL),
83  fDiffFlowResults6thOrderGFC(NULL),
84  fDiffFlowResults8thOrderGFC(NULL),
85  fCommonHistsResults2nd(NULL),
86  fCommonHistsResults4th(NULL),
87  fCommonHistsResults6th(NULL),
88  fCommonHistsResults8th(NULL),
89  fIntFlowGenFun(NULL),
90  fIntFlowGenFun4(NULL),//(only for other system of Eq.)
91  fIntFlowGenFun6(NULL),//(only for other system of Eq.)
92  fIntFlowGenFun8(NULL),//(only for other system of Eq.)
93  fIntFlowGenFun16(NULL),//(only for other system of Eq.)
94  fAvMultIntFlow4GFC(NULL),//(only for other system of Eq.)
95  fAvMultIntFlow6GFC(NULL),//(only for other system of Eq.)
96  fAvMultIntFlow8GFC(NULL),//(only for other system of Eq.)
97  fAvMultIntFlow16GFC(NULL),//(only for other system of Eq.)
98  fDiffFlowPtRPGenFunRe(NULL),
99  fDiffFlowPtRPGenFunIm(NULL),
100  fPtBinRPNoOfParticles(NULL),
101  fDiffFlowEtaRPGenFunRe(NULL),
102  fDiffFlowEtaRPGenFunIm(NULL),
103  fEtaBinRPNoOfParticles(NULL),
104  fDiffFlowPtPOIGenFunRe(NULL),
105  fDiffFlowPtPOIGenFunIm(NULL),
106  fPtBinPOINoOfParticles(NULL),
107  fDiffFlowEtaPOIGenFunRe(NULL),
108  fDiffFlowEtaPOIGenFunIm(NULL),
109  fEtaBinPOINoOfParticles(NULL),
110  /*
111  fDiffFlowGenFunRe0(NULL),
112  fDiffFlowGenFunRe1(NULL),
113  fDiffFlowGenFunRe2(NULL),
114  fDiffFlowGenFunRe3(NULL),
115  fDiffFlowGenFunRe4(NULL),
116  fDiffFlowGenFunRe5(NULL),
117  fDiffFlowGenFunRe6(NULL),
118  fDiffFlowGenFunRe7(NULL),
119  fDiffFlowGenFunIm0(NULL),
120  fDiffFlowGenFunIm1(NULL),
121  fDiffFlowGenFunIm2(NULL),
122  fDiffFlowGenFunIm3(NULL),
123  fDiffFlowGenFunIm4(NULL),
124  fDiffFlowGenFunIm5(NULL),
125  fDiffFlowGenFunIm6(NULL),
126  fDiffFlowGenFunIm7(NULL),
127  */
128  fCommonHists(NULL),
129  fOtherEquations(kFALSE),
130  fUsePhiWeights(kFALSE),
131  fUsePtWeights(kFALSE),
132  fUseEtaWeights(kFALSE),
133  fAverageOfSquaredWeight(NULL)
134 {
135  //constructor 
136  fHistList = new TList();
137  fWeightsList = new TList();
138  fWeightsList->SetName("Weights");
139  fWeightsList->SetOwner(kTRUE);
140  fR0=AliFlowCumuConstants::GetMaster()->GetR0();
141  //Pt:
142  fPtMax=AliFlowCommonConstants::GetMaster()->GetPtMax(); 
143  fPtMin=AliFlowCommonConstants::GetMaster()->GetPtMin();
144  fgknBinsPt=AliFlowCommonConstants::GetMaster()->GetNbinsPt();
145  if(fgknBinsPt)
146  {
147   fBinWidthPt=(fPtMax-fPtMin)/fgknBinsPt;  
148  } 
149  //Eta: 
150  fEtaMax=AliFlowCommonConstants::GetMaster()->GetEtaMax(); 
151  fEtaMin=AliFlowCommonConstants::GetMaster()->GetEtaMin();
152  fgknBinsEta=AliFlowCommonConstants::GetMaster()->GetNbinsEta();
153  if(fgknBinsEta)
154  {
155   fBinWidthEta=(fEtaMax-fEtaMin)/fgknBinsEta;  
156  } 
157  
158  fOtherEquations=AliFlowCumuConstants::GetMaster()->GetOtherEquations();
159 }
160
161 AliFlowAnalysisWithCumulants::~AliFlowAnalysisWithCumulants()
162 {
163  //desctructor
164  delete fHistList; 
165  delete fWeightsList;
166 }
167
168 //================================================================================================================
169
170 void AliFlowAnalysisWithCumulants::Init()
171 {
172  //various output histograms
173  
174  //save old value and prevent histograms from being added to directory
175  //to avoid name clashes in case multiple analaysis objects are used
176  //in an analysis
177  Bool_t oldHistAddStatus = TH1::AddDirectoryStatus();
178  TH1::AddDirectory(kFALSE);
179  
180  //average multiplicity 
181  fAvMultIntFlowGFC = new TProfile("fAvMultIntFlowGFC","Average Weighted Multiplicity",1,0,1,"s");
182  fAvMultIntFlowGFC->SetXTitle("");
183  fAvMultIntFlowGFC->SetYTitle("");
184  fAvMultIntFlowGFC->SetLabelSize(0.06);
185  fAvMultIntFlowGFC->SetMarkerStyle(25);
186  fAvMultIntFlowGFC->SetLabelOffset(0.01);
187  (fAvMultIntFlowGFC->GetXaxis())->SetBinLabel(1,"Average Weighted Multiplicity");
188  fHistList->Add(fAvMultIntFlowGFC);
189  
190  //averages of Q-vector components (1st bin: <Q_x>, 2nd bin: <Q_y>, 3rd bin: <(Q_x)^2>, 4th bin: <(Q_y)^2>)
191  fQVectorComponentsGFC = new TProfile("fQVectorComponentsGFC","Average of Q-vector components",4,0.,4.);
192  fQVectorComponentsGFC->SetXTitle("");
193  fQVectorComponentsGFC->SetYTitle("");
194  fQVectorComponentsGFC->SetLabelSize(0.06);
195  fQVectorComponentsGFC->SetMarkerStyle(25);
196  (fQVectorComponentsGFC->GetXaxis())->SetBinLabel(1,"<Q_{x}>");
197  (fQVectorComponentsGFC->GetXaxis())->SetBinLabel(2,"<Q_{y}>");
198  (fQVectorComponentsGFC->GetXaxis())->SetBinLabel(3,"<Q_{x}^{2}>");
199  (fQVectorComponentsGFC->GetXaxis())->SetBinLabel(4,"<Q_{y}^{2}>");
200  fHistList->Add(fQVectorComponentsGFC);
201  
202  //final results for integrated flow (v_n{2}, v_n{4},..., v_n{16}) from cumulants (by default n=2) 
203  fIntFlowResultsGFC = new TH1D("fIntFlowResultsGFC","Integrated Flow From Cumulants (Generating Function)",4,0,4);
204  fIntFlowResultsGFC->SetXTitle("");
205  fIntFlowResultsGFC->SetYTitle("");
206  fIntFlowResultsGFC->SetLabelSize(0.06);
207  //fIntFlowResultsGFC->SetTickLength(1);
208  fIntFlowResultsGFC->SetMarkerStyle(25);
209  (fIntFlowResultsGFC->GetXaxis())->SetBinLabel(1,"v_{n}{2}");
210  (fIntFlowResultsGFC->GetXaxis())->SetBinLabel(2,"v_{n}{4}");
211  (fIntFlowResultsGFC->GetXaxis())->SetBinLabel(3,"v_{n}{6}");
212  (fIntFlowResultsGFC->GetXaxis())->SetBinLabel(4,"v_{n}{8}");
213  fHistList->Add(fIntFlowResultsGFC);
214   
215  //final results for differential flow v_p/n{2} (by default p=n=2)
216  fDiffFlowResults2ndOrderGFC = new TH1D("fDiffFlowResults2ndOrderGFC","v'_2/2{2}",fgknBinsPt,fPtMin,fPtMax);
217  fDiffFlowResults2ndOrderGFC->SetXTitle("pt [GeV]");
218  fDiffFlowResults2ndOrderGFC->SetYTitle("");
219  fHistList->Add(fDiffFlowResults2ndOrderGFC);
220
221  //final results for differential flow v_p/n{4} (by default p=n=2) 
222  fDiffFlowResults4thOrderGFC = new TH1D("fDiffFlowResults4thOrderGFC","v'_2/2{4}",fgknBinsPt,fPtMin,fPtMax);
223  fDiffFlowResults4thOrderGFC->SetXTitle("pt [GeV]");
224  fDiffFlowResults4thOrderGFC->SetYTitle("");
225  fHistList->Add(fDiffFlowResults4thOrderGFC);
226  
227  //final results for differential flow v_p/n{6} (by default p=n=2)  
228  fDiffFlowResults6thOrderGFC = new TH1D("fDiffFlowResults6thOrderGFC","v'_2/2{6}",fgknBinsPt,fPtMin,fPtMax);
229  fDiffFlowResults6thOrderGFC->SetXTitle("pt [GeV]");
230  fDiffFlowResults6thOrderGFC->SetYTitle("");
231  fHistList->Add(fDiffFlowResults6thOrderGFC);
232  
233  //final results for differential flow v_p/n{8} (by default p=n=2)
234  fDiffFlowResults8thOrderGFC = new TH1D("fDiffFlowResults8thOrderGFC","v'_2/2{8}",fgknBinsPt,fPtMin,fPtMax);
235  fDiffFlowResults8thOrderGFC->SetXTitle("pt [GeV]");
236  fDiffFlowResults8thOrderGFC->SetYTitle("");
237  fHistList->Add(fDiffFlowResults8thOrderGFC);
238   
239  //avarage of the generating function for integrated flow <G[p][q]>
240  fIntFlowGenFun = new TProfile2D("fIntFlowGenFun","<G[p][q]>",fgkPmax,0.,(Double_t)fgkPmax,fgkQmax,0.,(Double_t)fgkQmax);
241  fIntFlowGenFun->SetXTitle("p");
242  fIntFlowGenFun->SetYTitle("q");
243  fHistList->Add(fIntFlowGenFun);
244
245  if(fOtherEquations)
246  {
247   //avarage of the generating function for integrated flow <G[p][q]> (only for other system of Eq. - up to 4th order)
248   fIntFlowGenFun4 = new TProfile2D("fIntFlowGenFun4","<G4[p4][q4]>",fgkPmax4,0.,(Double_t)fgkPmax4,fgkQmax4,0.,(Double_t)fgkQmax4);
249   fIntFlowGenFun4->SetXTitle("p4");
250   fIntFlowGenFun4->SetYTitle("q4");
251   fHistList->Add(fIntFlowGenFun4);
252
253   //avarage of the generating function for integrated flow <G[p][q]> (only for other system of Eq. - up to 6th order) 
254   fIntFlowGenFun6 = new TProfile2D("fIntFlowGenFun6","<G6[p6][q6]>",fgkPmax6,0.,(Double_t)fgkPmax6,fgkQmax6,0.,(Double_t)fgkQmax6);
255   fIntFlowGenFun6->SetXTitle("p6");
256   fIntFlowGenFun6->SetYTitle("q6");
257   fHistList->Add(fIntFlowGenFun6);
258
259   //avarage of the generating function for integrated flow <G[p][q]> (only for other system of Eq. - up to 8th order)
260   fIntFlowGenFun8 = new TProfile2D("fIntFlowGenFun8","<G8[p8][q8]>",fgkPmax8,0.,(Double_t)fgkPmax8,fgkQmax8,0.,(Double_t)fgkQmax8);
261   fIntFlowGenFun8->SetXTitle("p8");
262   fIntFlowGenFun8->SetYTitle("q8");
263   fHistList->Add(fIntFlowGenFun8);
264
265   //avarage of the generating function for integrated flow <G[p][q]> (only for other system of Eq. - up to 16th order)
266   fIntFlowGenFun16 = new TProfile2D("fIntFlowGenFun16","<G16[p16][q16]>",fgkPmax16,0.,(Double_t)fgkPmax16,fgkQmax16,0.,(Double_t)fgkQmax16);
267   fIntFlowGenFun16->SetXTitle("p16");
268   fIntFlowGenFun16->SetYTitle("q16");
269   fHistList->Add(fIntFlowGenFun16);
270  
271   //average multiplicity (only for other system of Eq. - up to 4th order)
272   fAvMultIntFlow4GFC = new TProfile("fAvMultIntFlow4GFC","Average Multiplicity",1,0,1,"s");
273   fAvMultIntFlow4GFC->SetXTitle("");
274   fAvMultIntFlow4GFC->SetYTitle("");
275   fAvMultIntFlow4GFC->SetLabelSize(0.06);
276   fAvMultIntFlow4GFC->SetMarkerStyle(25);
277   fAvMultIntFlow4GFC->SetLabelOffset(0.01);
278   (fAvMultIntFlow4GFC->GetXaxis())->SetBinLabel(1,"Average Multiplicity");
279   fHistList->Add(fAvMultIntFlow4GFC);
280  
281   //average multiplicity (only for other system of Eq. - up to 6th order)
282   fAvMultIntFlow6GFC = new TProfile("fAvMultIntFlow6GFC","Average Multiplicity",1,0,1,"s");
283   fAvMultIntFlow6GFC->SetXTitle("");
284   fAvMultIntFlow6GFC->SetYTitle("");
285   fAvMultIntFlow6GFC->SetLabelSize(0.06);
286   fAvMultIntFlow6GFC->SetMarkerStyle(25);
287   fAvMultIntFlow6GFC->SetLabelOffset(0.01);
288   (fAvMultIntFlow6GFC->GetXaxis())->SetBinLabel(1,"Average Multiplicity");
289   fHistList->Add(fAvMultIntFlow6GFC);
290  
291   //average multiplicity (only for other system of Eq. - up to 8th order)
292   fAvMultIntFlow8GFC = new TProfile("fAvMultIntFlow8GFC","Average Multiplicity",1,0,1,"s");
293   fAvMultIntFlow8GFC->SetXTitle("");
294   fAvMultIntFlow8GFC->SetYTitle("");
295   fAvMultIntFlow8GFC->SetLabelSize(0.06);
296   fAvMultIntFlow8GFC->SetMarkerStyle(25);
297   fAvMultIntFlow8GFC->SetLabelOffset(0.01);
298   (fAvMultIntFlow8GFC->GetXaxis())->SetBinLabel(1,"Average Multiplicity");
299   fHistList->Add(fAvMultIntFlow8GFC);
300  
301   //average multiplicity (only for other system of Eq. - up to 16th order)
302   fAvMultIntFlow16GFC = new TProfile("fAvMultIntFlow16GFC","Average Multiplicity",1,0,1,"s");
303   fAvMultIntFlow16GFC->SetXTitle("");
304   fAvMultIntFlow16GFC->SetYTitle("");
305   fAvMultIntFlow16GFC->SetLabelSize(0.06);
306   fAvMultIntFlow16GFC->SetMarkerStyle(25);
307   fAvMultIntFlow16GFC->SetLabelOffset(0.01);
308   (fAvMultIntFlow16GFC->GetXaxis())->SetBinLabel(1,"Average Multiplicity");
309   fHistList->Add(fAvMultIntFlow16GFC);
310  }
311  
312  //avarage of the real part of generating function for differential flow in Pt <Re(D[b][p][q])>
313  fDiffFlowPtRPGenFunRe = new TProfile3D("fDiffFlowPtRPGenFunRe","<Re(D[b][p][q])>",fgknBinsPt,(Double_t)(fPtMin/fBinWidthPt),(Double_t)(fPtMax/fBinWidthPt),fgkPmax,0.,(Double_t)fgkPmax,fgkQmax,0.,(Double_t)fgkQmax);
314  fDiffFlowPtRPGenFunRe->SetXTitle("b");
315  fDiffFlowPtRPGenFunRe->SetYTitle("p");
316  fDiffFlowPtRPGenFunRe->SetZTitle("q");
317  fDiffFlowPtRPGenFunRe->SetTitleOffset(1.44,"X");
318  fDiffFlowPtRPGenFunRe->SetTitleOffset(1.44,"Y");
319  fHistList->Add(fDiffFlowPtRPGenFunRe);
320  
321  //avarage of the imaginary part of generating function for differential flow in Pt <Im(D[b][p][q])>
322  fDiffFlowPtRPGenFunIm = new TProfile3D("fDiffFlowPtRPGenFunIm","<Im(D[b][p][q])>",fgknBinsPt,(Double_t)(fPtMin/fBinWidthPt),(Double_t)(fPtMax/fBinWidthPt),fgkPmax,0.,(Double_t)fgkPmax,fgkQmax,0.,(Double_t)fgkQmax);
323  fDiffFlowPtRPGenFunIm->SetXTitle("b");
324  fDiffFlowPtRPGenFunIm->SetYTitle("p");
325  fDiffFlowPtRPGenFunIm->SetZTitle("q");
326  fDiffFlowPtRPGenFunIm->SetTitleOffset(1.44,"X");
327  fDiffFlowPtRPGenFunIm->SetTitleOffset(1.44,"Y");
328  fHistList->Add(fDiffFlowPtRPGenFunIm);
329  
330  //number of particles per pt bin
331  fPtBinRPNoOfParticles = new TProfile("fPtBinRPNoOfParticles","Number of particles per #p_{t} bin",fgknBinsPt,fPtMin,fPtMax);
332  fPtBinRPNoOfParticles->SetXTitle("pt [GeV]");
333  fPtBinRPNoOfParticles->SetYTitle("");
334  fHistList->Add(fPtBinRPNoOfParticles);
335  
336  //avarage of the real part of generating function for differential flow in Eta <Re(D[b][p][q])>
337  fDiffFlowEtaRPGenFunRe = new TProfile3D("fDiffFlowEtaRPGenFunRe","<Re(D[b][p][q])>",fgknBinsEta,(Double_t)(fEtaMin/fBinWidthEta),(Double_t)(fEtaMax/fBinWidthEta),fgkPmax,0.,(Double_t)fgkPmax,fgkQmax,0.,(Double_t)fgkQmax);
338  fDiffFlowEtaRPGenFunRe->SetXTitle("b");
339  fDiffFlowEtaRPGenFunRe->SetYTitle("p");
340  fDiffFlowEtaRPGenFunRe->SetZTitle("q");
341  fDiffFlowEtaRPGenFunRe->SetTitleOffset(1.44,"X");
342  fDiffFlowEtaRPGenFunRe->SetTitleOffset(1.44,"Y");
343  fHistList->Add(fDiffFlowEtaRPGenFunRe);
344  
345  //avarage of the imaginary part of generating function for differential flow in Eta <Im(D[b][p][q])>
346  fDiffFlowEtaRPGenFunIm = new TProfile3D("fDiffFlowEtaRPGenFunIm","<Im(D[b][p][q])>",fgknBinsEta,(Double_t)(fEtaMin/fBinWidthEta),(Double_t)(fEtaMax/fBinWidthEta),fgkPmax,0.,(Double_t)fgkPmax,fgkQmax,0.,(Double_t)fgkQmax);
347  fDiffFlowEtaRPGenFunIm->SetXTitle("b");
348  fDiffFlowEtaRPGenFunIm->SetYTitle("p");
349  fDiffFlowEtaRPGenFunIm->SetZTitle("q");
350  fDiffFlowEtaRPGenFunIm->SetTitleOffset(1.44,"X");
351  fDiffFlowEtaRPGenFunIm->SetTitleOffset(1.44,"Y");
352  fHistList->Add(fDiffFlowEtaRPGenFunIm);
353  
354  //number of particles per eta bin
355  fEtaBinRPNoOfParticles = new TProfile("fEtaBinRPNoOfParticles","Number of particles per #eta bin",fgknBinsEta,fEtaMin,fEtaMax);
356  fEtaBinRPNoOfParticles->SetXTitle("#eta");
357  fEtaBinRPNoOfParticles->SetYTitle("");
358  fHistList->Add(fEtaBinRPNoOfParticles);
359  
360  //avarage of the real part of generating function for differential flow in Pt <Re(D[b][p][q])>
361  fDiffFlowPtPOIGenFunRe = new TProfile3D("fDiffFlowPtPOIGenFunRe","<Re(D[b][p][q])>",fgknBinsPt,(Double_t)(fPtMin/fBinWidthPt),(Double_t)(fPtMax/fBinWidthPt),fgkPmax,0.,(Double_t)fgkPmax,fgkQmax,0.,(Double_t)fgkQmax);
362  fDiffFlowPtPOIGenFunRe->SetXTitle("b");
363  fDiffFlowPtPOIGenFunRe->SetYTitle("p");
364  fDiffFlowPtPOIGenFunRe->SetZTitle("q");
365  fDiffFlowPtPOIGenFunRe->SetTitleOffset(1.44,"X");
366  fDiffFlowPtPOIGenFunRe->SetTitleOffset(1.44,"Y");
367  fHistList->Add(fDiffFlowPtPOIGenFunRe);
368  
369  //avarage of the imaginary part of generating function for differential flow in Pt <Im(D[b][p][q])>
370  fDiffFlowPtPOIGenFunIm = new TProfile3D("fDiffFlowPtPOIGenFunIm","<Im(D[b][p][q])>",fgknBinsPt,(Double_t)(fPtMin/fBinWidthPt),(Double_t)(fPtMax/fBinWidthPt),fgkPmax,0.,(Double_t)fgkPmax,fgkQmax,0.,(Double_t)fgkQmax);
371  fDiffFlowPtPOIGenFunIm->SetXTitle("b");
372  fDiffFlowPtPOIGenFunIm->SetYTitle("p");
373  fDiffFlowPtPOIGenFunIm->SetZTitle("q");
374  fDiffFlowPtPOIGenFunIm->SetTitleOffset(1.44,"X");
375  fDiffFlowPtPOIGenFunIm->SetTitleOffset(1.44,"Y");
376  fHistList->Add(fDiffFlowPtPOIGenFunIm);
377  
378  //number of particles per pt bin
379  fPtBinPOINoOfParticles = new TProfile("fPtBinPOINoOfParticles","Number of particles per #p_{t} bin",fgknBinsPt,fPtMin,fPtMax);
380  fPtBinPOINoOfParticles->SetXTitle("pt [GeV]");
381  fPtBinPOINoOfParticles->SetYTitle("");
382  fHistList->Add(fPtBinPOINoOfParticles);
383  
384  //avarage of the real part of generating function for differential flow in Eta <Re(D[b][p][q])>
385  fDiffFlowEtaPOIGenFunRe = new TProfile3D("fDiffFlowEtaPOIGenFunRe","<Re(D[b][p][q])>",fgknBinsEta,(Double_t)(fEtaMin/fBinWidthEta),(Double_t)(fEtaMax/fBinWidthEta),fgkPmax,0.,(Double_t)fgkPmax,fgkQmax,0.,(Double_t)fgkQmax);
386  fDiffFlowEtaPOIGenFunRe->SetXTitle("b");
387  fDiffFlowEtaPOIGenFunRe->SetYTitle("p");
388  fDiffFlowEtaPOIGenFunRe->SetZTitle("q");
389  fDiffFlowEtaPOIGenFunRe->SetTitleOffset(1.44,"X");
390  fDiffFlowEtaPOIGenFunRe->SetTitleOffset(1.44,"Y");
391  fHistList->Add(fDiffFlowEtaPOIGenFunRe);
392  
393  //avarage of the imaginary part of generating function for differential flow in Eta <Im(D[b][p][q])>
394  fDiffFlowEtaPOIGenFunIm = new TProfile3D("fDiffFlowEtaPOIGenFunIm","<Im(D[b][p][q])>",fgknBinsEta,(Double_t)(fEtaMin/fBinWidthEta),(Double_t)(fEtaMax/fBinWidthEta),fgkPmax,0.,(Double_t)fgkPmax,fgkQmax,0.,(Double_t)fgkQmax);
395  fDiffFlowEtaPOIGenFunIm->SetXTitle("b");
396  fDiffFlowEtaPOIGenFunIm->SetYTitle("p");
397  fDiffFlowEtaPOIGenFunIm->SetZTitle("q");
398  fDiffFlowEtaPOIGenFunIm->SetTitleOffset(1.44,"X");
399  fDiffFlowEtaPOIGenFunIm->SetTitleOffset(1.44,"Y");
400  fHistList->Add(fDiffFlowEtaPOIGenFunIm);
401  
402  //number of particles per eta bin
403  fEtaBinPOINoOfParticles = new TProfile("fEtaBinPOINoOfParticles","Number of particles per #eta bin",fgknBinsEta,fEtaMin,fEtaMax);
404  fEtaBinPOINoOfParticles->SetXTitle("#eta");
405  fEtaBinPOINoOfParticles->SetYTitle("");
406  fHistList->Add(fEtaBinPOINoOfParticles);
407  
408  /*
409  fDiffFlowGenFunRe0 = new TProfile2D("fDiffFlowGenFunRe0","Re(<D[b][0][q]>)",fgknBins,(Double_t)(fPtMin/fBinWidth),(Double_t)(fPtMax/fBinWidth),fgkQmax,0.,(Double_t)fgkQmax);
410  fDiffFlowGenFunRe0->SetXTitle("b");
411  fDiffFlowGenFunRe0->SetYTitle("q");
412  fHistList->Add(fDiffFlowGenFunRe0);
413
414  fDiffFlowGenFunIm0 = new TProfile2D("fDiffFlowGcout<<"HEY M1"<<endl;enFunIm0","Im(<D[b][0][q]>)",fgknBins,(Double_t)(fPtMin/fBinWidth),(Double_t)(fPtMax/fBinWidth),fgkQmax,0.,(Double_t)fgkQmax);
415  fDiffFlowGenFunIm0->SetXTitle("b");
416  fDiffFlowGenFunIm0->SetYTitle("q");
417  fHistList->Add(fDiffFlowGenFunIm0);
418  
419  fDiffFlowGenFunRe1 = new TProfile2D("fDiffFlowGenFunRe1","Re(<D[b][1][q]>)",fgknBins,(Double_t)(fPtMin/fBinWidth),(Double_t)(fPtMax/fBinWidth),fgkQmax,0.,(Double_t)fgkQmax);
420  fDiffFlowGenFunRe1->SetXTitle("b");
421  fDiffFlowGenFunRe1->SetYTitle("q");
422  fHistList->Add(fDiffFlowGenFunRe1);
423  
424  fDiffFlowGenFunIm1 = new TProfile2D("fDiffFlowGenFunIm1","Im(<D[b][1][q]>)",fgknBins,(Double_t)(fPtMin/fBinWidth),(Double_t)(fPtMax/fBinWidth),fgkQmax,0.,(Double_t)fgkQmax);
425  fDiffFlowGenFunIm1->SetXTitle("b");
426  fDiffFlowGenFunIm1->SetYTitle("q");
427  fHistList->Add(fDiffFlowGenFunIm1);
428  
429  fDiffFlowGenFunRe2 = new TProfile2D("fDiffFlowGenFunRe2","Re(<D[b][2][q]>)",fgknBins,(Double_t)(fPtMin/fBinWidth),(Double_t)(fPtMax/fBinWidth),fgkQmax,0.,(Double_t)fgkQmax);
430  fDiffFlowGenFunRe2->SetXTitle("b");
431  fDiffFlowGenFunRe2->SetYTitle("q");
432  fHistList->Add(fDiffFlowGenFunRe2);
433  
434  fDiffFlowGenFunIm2 = new TProfile2D("fDiffFlowGenFunIm2","Im(<D[b][2][q]>)",fgknBins,(Double_t)(fPtMin/fBinWidth),(Double_t)(fPtMax/fBinWidth),fgkQmax,0.,(Double_t)fgkQmax);
435  fDiffFlowGenFunIm2->SetXTitle("b");
436  fDiffFlowGenFunIm2->SetYTitle("q");
437  fHistList->Add(fDiffFlowGenFunIm2);
438  
439  fDiffFlowGenFunRe3 = new TProfile2D("fDiffFlowGenFunRe3","Re(<D[b][3][q]>)",fgknBins,(Double_t)(fPtMin/fBinWidth),(Double_t)(fPtMax/fBinWidth),fgkQmax,0.,(Double_t)fgkQmax);
440  fDiffFlowGenFunRe3->SetXTitle("b");
441  fDiffFlowGenFunRe3->SetYTitle("q");
442  fHistList->Add(fDiffFlowGenFunRe3);
443  
444  fDiffFlowGenFunIm3 = new TProfile2D("fDiffFlowGenFunIm3","Im(<D[b][3][q]>)",fgknBins,(Double_t)(fPtMin/fBinWidth),(Double_t)(fPtMax/fBinWidth),fgkQmax,0.,(Double_t)fgkQmax);
445  fDiffFlowGenFunIm3->SetXTitle("b");
446  fDiffFlowGenFunIm3->SetYTitle("q");
447  fHistList->Add(fDiffFlowGenFunIm3);
448  
449  fDiffFlowGenFunRe4 = new TProfile2D("fDiffFlowGenFunRe4","Re(<D[b][4][q]>)",fgknBins,(Double_t)(fPtMin/fBinWidth),(Double_t)(fPtMax/fBinWidth),fgkQmax,0.,(Double_t)fgkQmax);
450  fDiffFlowGenFunRe4->SetXTitle("b");
451  fDiffFlowGenFunRe4->SetYTitle("q");
452  fHistList->Add(fDiffFlowGenFunRe4);
453  
454  fDiffFlowGenFunIm4 = new TProfile2D("fDiffFlowGenFunIm4","Im(<D[b][4][q]>)",fgknBins,(Double_t)(fPtMin/fBinWidth),(Double_t)(fPtMax/fBinWidth),fgkQmax,0.,(Double_t)fgkQmax);
455  fDiffFlowGenFunIm4->SetXTitle("b");
456  fDiffFlowGenFunIm4->SetYTitle("q");
457  fHistList->Add(fDiffFlowGenFunIm4);
458  
459  fDiffFlowGenFunRe5 = new TProfile2D("fDiffFlowGenFunRe5","Re(<D[b][5][q]>)",fgkQmax,0.,(Double_t)fgkQmax,fgknBins,(Double_t)(fPtMin/fBinWidth),(Double_t)(fPtMax/fBinWidth));
460  fDiffFlowGenFunRe5->SetXTitle("b");
461  fDiffFlowGenFunRe5->SetYTitle("q");
462  fHistList->Add(fDiffFlowGenFunRe5);
463  
464  fDiffFlowGenFunIm5 = new TProfile2D("fDiffFlowGenFunIm5","Im(<D[b][5][q]>)",fgknBins,(Double_t)(fPtMin/fBinWidth),(Double_t)(fPtMax/fBinWidth),fgkQmax,0.,(Double_t)fgkQmax);
465  fDiffFlowGenFunIm5->SetXTitle("b");
466  fDiffFlowGenFunIm5->SetYTitle("q");
467  fHistList->Add(fDiffFlowGenFunIm5);
468  
469  fDiffFlowGenFunRe6 = new TProfile2D("fDiffFlowGenFunRe6","Re(<D[b][6][q]>)",fgknBins,(Double_t)(fPtMin/fBinWidth),(Double_t)(fPtMax/fBinWidth),fgkQmax,0.,(Double_t)fgkQmax);
470  fDiffFlowGenFunRe6->SetXTitle("b");
471  fDiffFlowGenFunRe6->SetYTitle("q");
472  fHistList->Add(fDiffFlowGenFunRe6);
473  
474  fDiffFlowGenFunIm6 = new TProfile2D("fDiffFlowGenFunIm6","Im(<D[b][6][q]>)",fgknBins,(Double_t)(fPtMin/fBinWidth),(Double_t)(fPtMax/fBinWidth),fgkQmax,0.,(Double_t)fgkQmax);
475  fDiffFlowGenFunIm6->SetXTitle("b");
476  fDiffFlowGenFunIm6->SetYTitle("q");
477  fHistList->Add(fDiffFlowGenFunIm6);
478  
479  fDiffFlowGenFunRe7 = new TProfile2D("fDiffFlowGenFunRe7","Re(<D[b][7][q]>)",fgknBins,(Double_t)(fPtMin/fBinWidth),(Double_t)(fPtMax/fBinWidth),fgkQmax,0.,(Double_t)fgkQmax);
480  fDiffFlowGenFunRe7->SetXTitle("b");
481  fDiffFlowGenFunRe7->SetYTitle("q");
482  fHistList->Add(fDiffFlowGenFunRe7);
483  
484  fDiffFlowGenFunIm7 = new TProfile2D("fDiffFlowGenFunIm7","Im(<D[b][7][q]>)",fgknBins,(Double_t)(fPtMin/fBinWidth),(Double_t)(fPtMax/fBinWidth),fgkQmax,0.,(Double_t)fgkQmax);
485  fDiffFlowGenFunIm7->SetXTitle("b");
486  fDiffFlowGenFunIm7->SetYTitle("q");
487  fHistList->Add(fDiffFlowGenFunIm7);
488  */
489  
490  //common control histograms
491  fCommonHists = new AliFlowCommonHist("AliFlowCommonHistGFC");
492  fHistList->Add(fCommonHists);  
493  
494  //common histograms for final results (2nd order)
495  fCommonHistsResults2nd = new AliFlowCommonHistResults("AliFlowCommonHistResults2ndOrderGFC");
496  fHistList->Add(fCommonHistsResults2nd); 
497  
498  //common histograms for final results (4th order)
499  fCommonHistsResults4th = new AliFlowCommonHistResults("AliFlowCommonHistResults4thOrderGFC");
500  fHistList->Add(fCommonHistsResults4th);
501  
502  //common histograms for final results (6th order)
503  fCommonHistsResults6th = new AliFlowCommonHistResults("AliFlowCommonHistResults6thOrderGFC");
504  fHistList->Add(fCommonHistsResults6th); 
505  
506  //common histograms for final results (8th order)
507  fCommonHistsResults8th = new AliFlowCommonHistResults("AliFlowCommonHistResults8thOrderGFC");
508  fHistList->Add(fCommonHistsResults8th);
509  
510  //<w^2>
511  fAverageOfSquaredWeight = new TProfile("fAverageOfSquaredWeight","<w^{2}>",1,0,1);
512  fAverageOfSquaredWeight->SetLabelSize(0.06);
513  fAverageOfSquaredWeight->SetMarkerStyle(25);
514  fAverageOfSquaredWeight->SetLabelOffset(0.01);
515  (fAverageOfSquaredWeight->GetXaxis())->SetBinLabel(1,"<w^{2}>");
516  fHistList->Add(fAverageOfSquaredWeight);
517  
518  // add list fWeightsList with weights to the main list
519  fHistList->Add(fWeightsList);
520
521  TH1::AddDirectory(oldHistAddStatus);
522 }//end of Init()
523
524 //================================================================================================================
525
526 void AliFlowAnalysisWithCumulants::Make(AliFlowEventSimple* anEvent)
527 {
528  //running over data:
529  Int_t nPrim = anEvent->NumberOfTracks(); //total multiplicity
530  
531  Int_t nEventNSelTracksRP = anEvent->GetEventNSelTracksRP(); //selected multiplicity (particles used for int. flow)
532  
533  Int_t n = 2; // int flow harmonic (to be improved)
534  
535  //---------------------------------------------------------------------------------------------------------
536  // weights:
537  Bool_t useWeights = fUsePhiWeights||fUsePtWeights||fUseEtaWeights;
538  
539  TH1F *phiWeights = NULL; // histogram with phi weights
540  TH1D *ptWeights  = NULL; // histogram with pt weights
541  TH1D *etaWeights = NULL; // histogram with eta weights
542
543  Double_t wPhi = 1.; // phi weight
544  Double_t wPt  = 1.; // pt weight
545  Double_t wEta = 1.; // eta weight
546    
547  if(useWeights)
548  {
549   if(!fWeightsList)
550   {
551    cout<<" WARNING: fWeightsList is NULL pointer in GFC::Make() !!!!"<<endl;
552    exit(0);
553   }
554   if(fUsePhiWeights) 
555   {
556    phiWeights = dynamic_cast<TH1F *>(fWeightsList->FindObject("phi_weights"));
557    if(!phiWeights)
558    {
559     cout<<" WARNING: couldn't access the histogram with phi weights in GFC::Make() !!!!"<<endl;
560     exit(0);
561    } 
562   } 
563   if(fUsePtWeights) 
564   { 
565    ptWeights = dynamic_cast<TH1D *>(fWeightsList->FindObject("pt_weights"));
566    if(!ptWeights) 
567    {
568     cout<<" WARNING: couldn't access the histogram with pt weights in GFC::Make() !!!!"<<endl;
569     exit(0);
570    } 
571   } 
572   if(fUseEtaWeights) 
573   {
574    etaWeights = dynamic_cast<TH1D *>(fWeightsList->FindObject("eta_weights"));
575    if(!etaWeights) 
576    {
577     cout<<" WARNING: couldn't access the histogram with eta weights in GFC::Make() !!!!"<<endl;
578     exit(0);
579    }
580   } 
581  } 
582  //---------------------------------------------------------------------------------------------------------
583   
584  if(nEventNSelTracksRP>9) //generating function formalism applied here make sense only for selected multiplicity >= 10 
585  { 
586  //fill the common control histograms
587  fCommonHists->FillControlHistograms(anEvent);   
588   
589  //initializing the generating function G[p][q] for integrated flow 
590  Double_t genfunG[fgkPmax][fgkQmax];
591  
592  for(Int_t p=0;p<fgkPmax;p++)
593  {
594   for(Int_t q=0;q<fgkQmax;q++)
595   {
596    genfunG[p][q]=1.;
597   }   
598  }
599
600  Int_t nSelTracksRP = 0; //cross-checking the selected multiplicity
601   
602  Double_t dPhi = 0.;
603  Double_t dPt  = 0.;
604  Double_t dEta = 0.;
605  Int_t nBinsPhi = 0;
606  
607  for(Int_t i=0;i<nPrim;i++)
608  {
609   fTrack=anEvent->GetTrack(i);
610   if(fTrack && fTrack->InRPSelection())
611   {
612    nSelTracksRP++;
613    // get azimuthal angle, momentum and pseudorapidity of a particle:
614    dPhi = fTrack->Phi();
615    dPt  = fTrack->Pt();
616    dEta = fTrack->Eta();
617    // phi weights:
618    if(fUsePhiWeights) 
619    {
620     nBinsPhi = phiWeights->GetNbinsX();
621     if(nBinsPhi) 
622     {
623      wPhi = phiWeights->GetBinContent(1+(Int_t)(TMath::Floor(dPhi*nBinsPhi/TMath::TwoPi())));
624     }
625    } 
626    // pt weights:
627    if(fUsePtWeights)
628    {          
629     if(fBinWidthPt) 
630     {
631      wPt = ptWeights->GetBinContent(1+(Int_t)(TMath::Floor((dPt-fPtMin)/fBinWidthPt)));
632     }
633    }             
634    // eta weights:
635    if(fUseEtaWeights)
636    {    
637     if(fBinWidthEta)
638     {
639      wEta=etaWeights->GetBinContent(1+(Int_t)(TMath::Floor((dEta-fEtaMin)/fBinWidthEta))); 
640     }
641    }
642    // evaluate the generating function:
643    for(Int_t p=0;p<fgkPmax;p++)
644    {
645     for(Int_t q=0;q<fgkQmax;q++)
646     {
647      genfunG[p][q]*=(1.+wPhi*wPt*wEta*(2.*fR0*sqrt(p+1.)/nEventNSelTracksRP)*cos(fgkFlow*dPhi-2.*q*TMath::Pi()/fgkQmax));
648     }
649    }
650    // calculate <w^2> 
651    fAverageOfSquaredWeight->Fill(0.,pow(wPhi*wPt*wEta,2.),1); 
652   } // end of if(fTrack && fTrack->InRPSelection())
653  } // end of for(Int_t i=0;i<nPrim;i++) 
654  
655  //storing the value of G[p][q] in 2D profile in order to get automatically the avarage <G[p][q]>
656  for(Int_t p=0;p<fgkPmax;p++)
657  {
658   for(Int_t q=0;q<fgkQmax;q++)
659   {
660    fIntFlowGenFun->Fill((Double_t)p,(Double_t)q,genfunG[p][q],1);
661   }
662  } 
663  
664  //storing the selected multiplicity (if fAvMultIntFlow is not filled here then you had wrongly selected the particles used for integrated flow)
665  if(nSelTracksRP==nEventNSelTracksRP)
666  {
667   fAvMultIntFlowGFC->Fill(0.,nSelTracksRP,1);
668  } else 
669    {
670     cout<<"WARNING: nSelTracksRP != nEventNSelTracksRP in GFC::Make(). Something is wrong with RP flagging !!!!"<<endl;
671    }
672  
673  // calculating Q-vector of event (needed for errors)
674  AliFlowVector fQVector;
675  fQVector.Set(0.,0.);
676  fQVector.SetMult(0);
677  fQVector=anEvent->GetQ(1*n,fWeightsList,fUsePhiWeights,fUsePtWeights,fUseEtaWeights); // get the Q vector for this event
678  fQVectorComponentsGFC->Fill(0.,fQVector.X(),1);         // in the 1st bin fill Q_x
679  fQVectorComponentsGFC->Fill(1.,fQVector.Y(),1);         // in the 2nd bin fill Q_y
680  fQVectorComponentsGFC->Fill(2.,pow(fQVector.X(),2.),1); // in the 3rd bin fill (Q_x)^2
681  fQVectorComponentsGFC->Fill(3.,pow(fQVector.Y(),2.),1); // in the 4th bin fill (Q_y)^2
682  
683  //3D profiles for differential flow in pt and eta
684  //second loop over data: evaluating the generating function D[b][p][q] for differential flow 
685  //remark 0: D[b][p][q] is a complex number => real and imaginary part are calculated separately
686  //remark 1: note that bellow G[p][q] is needed, the value of generating function for integrated flow for the CURRENT event
687  //remark 2: results are stored in 3D profiles in order to automatically get <Re(D[b][p][q])> and <Im(D[b][p][q])>
688  for(Int_t i=0;i<nPrim;i++)
689  {
690   fTrack=anEvent->GetTrack(i);
691   if(fTrack)
692   {
693    if(fTrack->InRPSelection() && fTrack->InPOISelection())
694    {
695     // get azimuthal angle, momentum and pseudorapidity of a particle:
696     dPhi = fTrack->Phi();
697     dPt  = fTrack->Pt();
698     dEta = fTrack->Eta();
699     
700     fPtBinPOINoOfParticles->Fill(dPt,dPt,1.);
701     fEtaBinPOINoOfParticles->Fill(dEta,dEta,1.);
702     
703     // phi weights:
704     if(fUsePhiWeights) 
705     {
706      nBinsPhi = phiWeights->GetNbinsX();
707      if(nBinsPhi) 
708      {
709       wPhi = phiWeights->GetBinContent(1+(Int_t)(TMath::Floor(dPhi*nBinsPhi/TMath::TwoPi())));
710      }
711     } 
712     // pt weights:
713     if(fUsePtWeights)
714     {          
715      if(fBinWidthPt) 
716      {
717       wPt = ptWeights->GetBinContent(1+(Int_t)(TMath::Floor((dPt-fPtMin)/fBinWidthPt)));
718      }
719     }             
720     // eta weights:
721     if(fUseEtaWeights)
722     {    
723      if(fBinWidthEta)
724      {
725       wEta=etaWeights->GetBinContent(1+(Int_t)(TMath::Floor((dEta-fEtaMin)/fBinWidthEta))); 
726      }
727     }
728    
729     for(Int_t p=0;p<fgkPmax;p++)
730     {
731      for(Int_t q=0;q<fgkQmax;q++)
732      {
733       //real part (Pt)
734       fDiffFlowPtPOIGenFunRe->Fill(dPt/fBinWidthPt,(Double_t)p,(Double_t)q,genfunG[p][q]*cos(fgkMltpl*fgkFlow*dPhi)/(1.+wPhi*wPt*wEta*(2.*fR0*sqrt(p+1.)/nSelTracksRP)*cos(fgkFlow*dPhi-2.*q*TMath::Pi()/fgkQmax)),1.);
735       //imaginary part (Pt)
736       fDiffFlowPtPOIGenFunIm->Fill(dPt/fBinWidthPt,(Double_t)p,(Double_t)q,genfunG[p][q]*sin(fgkMltpl*fgkFlow*dPhi)/(1.+wPhi*wPt*wEta*(2.*fR0*sqrt(p+1.)/nSelTracksRP)*cos(fgkFlow*dPhi-2.*q*TMath::Pi()/fgkQmax)),1.); 
737       //real part (Eta)
738       fDiffFlowEtaPOIGenFunRe->Fill(dEta/fBinWidthEta,(Double_t)p,(Double_t)q,genfunG[p][q]*cos(fgkMltpl*fgkFlow*dPhi)/(1.+wPhi*wPt*wEta*(2.*fR0*sqrt(p+1.)/nSelTracksRP)*cos(fgkFlow*dPhi-2.*q*TMath::Pi()/fgkQmax)),1.);
739       //imaginary part (Eta)
740       fDiffFlowEtaPOIGenFunIm->Fill(dEta/fBinWidthEta,(Double_t)p,(Double_t)q,genfunG[p][q]*sin(fgkMltpl*fgkFlow*dPhi)/(1.+wPhi*wPt*wEta*(2.*fR0*sqrt(p+1.)/nSelTracksRP)*cos(fgkFlow*dPhi-2.*q*TMath::Pi()/fgkQmax)),1.); 
741      }
742     }
743    } 
744    else if(fTrack->InPOISelection() && !(fTrack->InRPSelection()))
745    {
746     // get azimuthal angle, momentum and pseudorapidity of a particle:
747     dPhi = fTrack->Phi();
748     dPt  = fTrack->Pt();
749     dEta = fTrack->Eta();
750    
751     fPtBinPOINoOfParticles->Fill(dPt,dPt,1.);
752     fEtaBinPOINoOfParticles->Fill(dEta,dEta,1.);
753     for(Int_t p=0;p<fgkPmax;p++)
754     {
755      for(Int_t q=0;q<fgkQmax;q++)
756      {
757       //real part (Pt)
758       fDiffFlowPtPOIGenFunRe->Fill(dPt/fBinWidthPt,(Double_t)p,(Double_t)q,genfunG[p][q]*cos(fgkMltpl*fgkFlow*dPhi),1.);
759       //imaginary part (Pt)
760       fDiffFlowPtPOIGenFunIm->Fill(dPt/fBinWidthPt,(Double_t)p,(Double_t)q,genfunG[p][q]*sin(fgkMltpl*fgkFlow*dPhi),1.); 
761       //real part (Eta)
762       fDiffFlowEtaPOIGenFunRe->Fill(dEta/fBinWidthEta,(Double_t)p,(Double_t)q,genfunG[p][q]*cos(fgkMltpl*fgkFlow*dPhi),1.);
763       //imaginary part (Eta)
764       fDiffFlowEtaPOIGenFunIm->Fill(dEta/fBinWidthEta,(Double_t)p,(Double_t)q,genfunG[p][q]*sin(fgkMltpl*fgkFlow*dPhi),1.);                     
765      } 
766     }
767    }
768    //RP:
769    if(fTrack->InRPSelection())       
770    {
771     // get azimuthal angle, momentum and pseudorapidity of a particle:
772     dPhi = fTrack->Phi();
773     dPt  = fTrack->Pt();
774     dEta = fTrack->Eta();
775    
776     fPtBinRPNoOfParticles->Fill(dPt,dPt,1.);
777     fEtaBinRPNoOfParticles->Fill(dEta,dEta,1.);
778     
779     // phi weights:
780     if(fUsePhiWeights) 
781     {
782      nBinsPhi = phiWeights->GetNbinsX();
783      if(nBinsPhi) 
784      {
785       wPhi = phiWeights->GetBinContent(1+(Int_t)(TMath::Floor(dPhi*nBinsPhi/TMath::TwoPi())));
786      }
787     } 
788     // pt weights:
789     if(fUsePtWeights)
790     {          
791      if(fBinWidthPt) 
792      {
793       wPt = ptWeights->GetBinContent(1+(Int_t)(TMath::Floor((dPt-fPtMin)/fBinWidthPt)));
794      }
795     }             
796     // eta weights:
797     if(fUseEtaWeights)
798     {    
799      if(fBinWidthEta)
800      {
801       wEta=etaWeights->GetBinContent(1+(Int_t)(TMath::Floor((dEta-fEtaMin)/fBinWidthEta))); 
802      }
803     }
804     
805     for(Int_t p=0;p<fgkPmax;p++)
806     {
807      for(Int_t q=0;q<fgkQmax;q++)
808      {
809       //real part (Pt)
810       fDiffFlowPtRPGenFunRe->Fill(dPt/fBinWidthPt,(Double_t)p,(Double_t)q,genfunG[p][q]*cos(fgkMltpl*fgkFlow*dPhi)/(1.+wPhi*wPt*wEta*(2.*fR0*sqrt(p+1.)/nSelTracksRP)*cos(fgkFlow*dPhi-2.*q*TMath::Pi()/fgkQmax)),1.);
811       //imaginary part (Pt)
812       fDiffFlowPtRPGenFunIm->Fill(dPt/fBinWidthPt,(Double_t)p,(Double_t)q,genfunG[p][q]*sin(fgkMltpl*fgkFlow*dPhi)/(1.+wPhi*wPt*wEta*(2.*fR0*sqrt(p+1.)/nSelTracksRP)*cos(fgkFlow*dPhi-2.*q*TMath::Pi()/fgkQmax)),1.); 
813       //real part (Eta)
814       fDiffFlowEtaRPGenFunRe->Fill(dEta/fBinWidthEta,(Double_t)p,(Double_t)q,genfunG[p][q]*cos(fgkMltpl*fgkFlow*dPhi)/(1.+wPhi*wPt*wEta*(2.*fR0*sqrt(p+1.)/nSelTracksRP)*cos(fgkFlow*dPhi-2.*q*TMath::Pi()/fgkQmax)),1.);
815       //imaginary part (Eta)
816       fDiffFlowEtaRPGenFunIm->Fill(dEta/fBinWidthEta,(Double_t)p,(Double_t)q,genfunG[p][q]*sin(fgkMltpl*fgkFlow*dPhi)/(1.+wPhi*wPt*wEta*(2.*fR0*sqrt(p+1.)/nSelTracksRP)*cos(fgkFlow*dPhi-2.*q*TMath::Pi()/fgkQmax)),1.); 
817      }
818     }
819    }//end of if(fTrack->InRPSelection())                   
820   }//end of if(fTrack)  
821  }//ending the second loop over data    
822  
823  
824  
825  /*
826  //sixteen 2D profiles for differential flow          
827  for(Int_t i=0;i<nPrim;i++){
828   fTrack=anEvent->GetTrack(i);
829   if (fTrack && fTrack->InPOISelection()){
830    //for(Int_t p=0;p<fgkPmax;p++){
831     for(Int_t q=0;q<fgkQmax;q++){
832      //real part
833      fDiffFlowGenFunRe0->Fill(fTrack->Pt()/fBinWidth,(Double_t)q,genfunG[0][q]*cos(fgkMltpl*fgkFlow*fTrack->Phi())/(1.+(2.*fR0*sqrt(0.+1.)/nSelTracksRP)*cos(fgkFlow*fTrack->Phi()-2.*q*TMath::Pi()/fgkQmax)),1.);
834      //imaginary part
835      fDiffFlowGenFunIm0->Fill(fTrack->Pt()/fBinWidth,(Double_t)q,genfunG[0][q]*sin(fgkMltpl*fgkFlow*fTrack->Phi())/(1.+(2.*fR0*sqrt(0.+1.)/nSelTracksRP)*cos(fgkFlow*fTrack->Phi()-2.*q*TMath::Pi()/fgkQmax)),1.); 
836      //-----------------------------------------------------------------------
837      //real part
838      fDiffFlowGenFunRe1->Fill(fTrack->Pt()/fBinWidth,(Double_t)q,genfunG[1][q]*cos(fgkMltpl*fgkFlow*fTrack->Phi())/(1.+(2.*fR0*sqrt(1.+1.)/nSelTracksRP)*cos(fgkFlow*fTrack->Phi()-2.*q*TMath::Pi()/fgkQmax)),1.);
839      //imaginary part
840      fDiffFlowGenFunIm1->Fill(fTrack->Pt()/fBinWidth,(Double_t)q,genfunG[1][q]*sin(fgkMltpl*fgkFlow*fTrack->Phi())/(1.+(2.*fR0*sqrt(1.+1.)/nSelTracksRP)*cos(fgkFlow*fTrack->Phi()-2.*q*TMath::Pi()/fgkQmax)),1.);             
841     //-----------------------------------------------------------------------
842      //real part
843      fDiffFlowGenFunRe2->Fill(fTrack->Pt()/fBinWidth,(Double_t)q,genfunG[2][q]*cos(fgkMltpl*fgkFlow*fTrack->Phi())/(1.+(2.*fR0*sqrt(2.+1.)/nSelTracksRP)*cos(fgkFlow*fTrack->Phi()-2.*q*TMath::Pi()/fgkQmax)),1.);
844      //imaginary part
845      fDiffFlowGenFunIm2->Fill(fTrack->Pt()/fBinWidth,(Double_t)q,genfunG[2][q]*sin(fgkMltpl*fgkFlow*fTrack->Phi())/(1.+(2.*fR0*sqrt(2.+1.)/nSelTracksRP)*cos(fgkFlow*fTrack->Phi()-2.*q*TMath::Pi()/fgkQmax)),1.);  
846      //-----------------------------------------------------------------------
847      //real part
848      fDiffFlowGenFunRe3->Fill(fTrack->Pt()/fBinWidth,(Double_t)q,genfunG[3][q]*cos(fgkMltpl*fgkFlow*fTrack->Phi())/(1.+(2.*fR0*sqrt(3.+1.)/nSelTracksRP)*cos(fgkFlow*fTrack->Phi()-2.*q*TMath::Pi()/fgkQmax)),1.);
849      //imaginary part
850      fDiffFlowGenFunIm3->Fill(fTrack->Pt()/fBinWidth,(Double_t)q,genfunG[3][q]*sin(fgkMltpl*fgkFlow*fTrack->Phi())/(1.+(2.*fR0*sqrt(3.+1.)/nSelTracksRP)*cos(fgkFlow*fTrack->Phi()-2.*q*TMath::Pi()/fgkQmax)),1.);  
851      //-----------------------------------------------------------------------
852      //real part
853      fDiffFlowGenFunRe4->Fill(fTrack->Pt()/fBinWidth,(Double_t)q,genfunG[4][q]*cos(fgkMltpl*fgkFlow*fTrack->Phi())/(1.+(2.*fR0*sqrt(4.+1.)/nSelTracksRP)*cos(fgkFlow*fTrack->Phi()-2.*q*TMath::Pi()/fgkQmax)),1.);
854      //imaginary part
855      fDiffFlowGenFunIm4->Fill(fTrack->Pt()/fBinWidth,(Double_t)q,genfunG[4][q]*sin(fgkMltpl*fgkFlow*fTrack->Phi())/(1.+(2.*fR0*sqrt(4.+1.)/nSelTracksRP)*cos(fgkFlow*fTrack->Phi()-2.*q*TMath::Pi()/fgkQmax)),1.);  
856      //-----------------------------------------------------------------------
857      //real part
858      fDiffFlowGenFunRe5->Fill(fTrack->Pt()/fBinWidth,(Double_t)q,genfunG[5][q]*cos(fgkMltpl*fgkFlow*fTrack->Phi())/(1.+(2.*fR0*sqrt(5.+1.)/nSelTracksRP)*cos(fgkFlow*fTrack->Phi()-2.*q*TMath::Pi()/fgkQmax)),1.);
859      //imaginary part
860      fDiffFlowGenFunIm5->Fill(fTrack->Pt()/fBinWidth,(Double_t)q,genfunG[5][q]*sin(fgkMltpl*fgkFlow*fTrack->Phi())/(1.+(2.*fR0*sqrt(5.+1.)/nSelTracksRP)*cos(fgkFlow*fTrack->Phi()-2.*q*TMath::Pi()/fgkQmax)),1.);  
861      //-----------------------------------------------------------------------
862      //real part
863      fDiffFlowGenFunRe6->Fill(fTrack->Pt()/fBinWidth,(Double_t)q,genfunG[6][q]*cos(fgkMltpl*fgkFlow*fTrack->Phi())/(1.+(2.*fR0*sqrt(6.+1.)/nSelTracksRP)*cos(fgkFlow*fTrack->Phi()-2.*q*TMath::Pi()/fgkQmax)),1.);
864      //imaginary part
865      fDiffFlowGenFunIm6->Fill(fTrack->Pt()/fBinWidth,(Double_t)q,genfunG[6][q]*sin(fgkMltpl*fgkFlow*fTrack->Phi())/(1.+(2.*fR0*sqrt(6.+1.)/nSelTracksRP)*cos(fgkFlow*fTrack->Phi()-2.*q*TMath::Pi()/fgkQmax)),1.);
866      //-----------------------------------------------------------------------
867      //real part
868      fDiffFlowGenFunRe7->Fill(fTrack->Pt()/fBinWidth,(Double_t)q,genfunG[7][q]*cos(fgkMltpl*fgkFlow*fTrack->Phi())/(1.+(2.*fR0*sqrt(7.+1.)/nSelTracksRP)*cos(fgkFlow*fTrack->Phi()-2.*q*TMath::Pi()/fgkQmax)),1.);
869      //imaginary part
870      fDiffFlowGenFunIm7->Fill(fTrack->Pt()/fBinWidth,(Double_t)q,genfunG[7][q]*sin(fgkMltpl*fgkFlow*fTrack->Phi())/(1.+(2.*fR0*sqrt(7.+1.)/nSelTracksRP)*cos(fgkFlow*fTrack->Phi()-2.*q*TMath::Pi()/fgkQmax)),1.);    
871    }
872    //}
873   }  
874  }//ending the second loop over data            
875  */
876       
877  }//end of if(nEventNSelTracksRP>9)                   
878  
879
880  
881   
882    
883     
884      
885       
886        
887         
888          
889           
890            
891  //off the record: numerical equations for cumulants solved up to different highest order  
892  if(fOtherEquations)
893  {
894   //running over data
895   Int_t nPrimOE = anEvent->NumberOfTracks();//total multiplicity 
896   
897   Int_t nEventNSelTracksRPOE = anEvent->GetEventNSelTracksRP();
898   
899   Double_t genfunG4[fgkPmax4][fgkQmax4];
900   Double_t genfunG6[fgkPmax6][fgkQmax6];
901   Double_t genfunG8[fgkPmax8][fgkQmax8];
902   Double_t genfunG16[fgkPmax16][fgkQmax16];
903   for(Int_t p=0;p<fgkPmax16;p++)
904   {
905    for(Int_t q=0;q<fgkQmax16;q++)
906    {
907     genfunG16[p][q]=1.;
908     if(p<fgkPmax8 && q<fgkQmax8)
909     {
910      genfunG8[p][q]=1.;
911      if(p<fgkPmax6 && q<fgkQmax6)
912      {
913       genfunG6[p][q]=1.;
914       if(p<fgkPmax4 && q<fgkQmax4)
915       {
916        genfunG4[p][q]=1.;
917       }
918      }
919     } 
920    }
921   } 
922    
923   //multiplicities: 
924   if(nEventNSelTracksRPOE>15) fAvMultIntFlow16GFC->Fill(0.,nEventNSelTracksRPOE,1);
925   if(nEventNSelTracksRPOE>7) fAvMultIntFlow8GFC->Fill(0.,nEventNSelTracksRPOE,1);
926   if(nEventNSelTracksRPOE>5) fAvMultIntFlow6GFC->Fill(0.,nEventNSelTracksRPOE,1);
927   if(nEventNSelTracksRPOE>3) fAvMultIntFlow4GFC->Fill(0.,nEventNSelTracksRPOE,1);  
928   
929   //first loop over data: evaluating the generating function G[p][q] for integrated flow 
930   for(Int_t i=0;i<nPrimOE;i++)
931   {
932    fTrack=anEvent->GetTrack(i);
933    if(fTrack && fTrack->InRPSelection())
934    {
935     for(Int_t p=0;p<fgkPmax16;p++)
936     {
937      for(Int_t q=0;q<fgkQmax16;q++)
938      {
939       if(nEventNSelTracksRPOE>15)
940       {
941        genfunG16[p][q]*=(1.+(2.*fR0*sqrt(p+1.)/nEventNSelTracksRPOE)*cos(fgkFlow*fTrack->Phi()-2.*q*TMath::Pi()/fgkQmax16));
942       }       
943       if(p<fgkPmax8 && q<fgkQmax8)
944       {
945        if(nEventNSelTracksRPOE>7)
946        { 
947         genfunG8[p][q]*=(1.+(2.*fR0*sqrt(p+1.)/nEventNSelTracksRPOE)*cos(fgkFlow*fTrack->Phi()-2.*q*TMath::Pi()/fgkQmax8));
948        }
949        if(p<fgkPmax6 && q<fgkQmax6)
950        {
951         if(nEventNSelTracksRPOE>5) 
952         {
953          genfunG6[p][q]*=(1.+(2.*fR0*sqrt(p+1.)/nEventNSelTracksRPOE)*cos(fgkFlow*fTrack->Phi()-2.*q*TMath::Pi()/fgkQmax6));
954         }
955         if(p<fgkPmax4 && q<fgkQmax4)
956         {
957          if(nEventNSelTracksRPOE>3)
958          {
959           genfunG4[p][q]*=(1.+(2.*fR0*sqrt(p+1.)/nEventNSelTracksRPOE)*cos(fgkFlow*fTrack->Phi()-2.*q*TMath::Pi()/fgkQmax4));
960          } 
961         }
962        }
963       }
964      } 
965     }  
966    }//end of if(fTrack && fTrack->InRPSelection())
967   }//ending the loop over data
968  
969  //storing the value of G[p][q] in 2D profile in order to get automatically the avarage <G[p][q]>
970  for(Int_t p=0;p<fgkPmax16;p++)
971  {
972   for(Int_t q=0;q<fgkQmax16;q++)
973   {
974    if(nEventNSelTracksRPOE>15) fIntFlowGenFun16->Fill((Double_t)p,(Double_t)q,genfunG16[p][q],1);
975    if(p<fgkPmax8 && q<fgkQmax8)
976    {
977     if(nEventNSelTracksRPOE>7) fIntFlowGenFun8->Fill((Double_t)p,(Double_t)q,genfunG8[p][q],1);
978     if(p<fgkPmax6 && q<fgkQmax6)
979     {
980      if(nEventNSelTracksRPOE>5) fIntFlowGenFun6->Fill((Double_t)p,(Double_t)q,genfunG6[p][q],1);
981      if(p<fgkPmax4 && q<fgkQmax4)
982      {
983       if(nEventNSelTracksRPOE>3) fIntFlowGenFun4->Fill((Double_t)p,(Double_t)q,genfunG4[p][q],1);
984      }
985     }
986    } 
987   }
988  }
989 }//end of if(fOtherEquations)                  
990                                                                                                                                                                                                                                                                                                                                                                                                                                                                        
991 }//end of Make()
992
993 //================================================================================================================
994
995 void AliFlowAnalysisWithCumulants::GetOutputHistograms(TList *outputListHistos) 
996 {
997  // get the pointers to all output histograms before calling Finish()
998  if (outputListHistos) 
999  {
1000   //Get the common histograms from the output list
1001   //histograms to store the final results
1002   TH1D *intFlowResults   = dynamic_cast<TH1D*>(outputListHistos->FindObject("fIntFlowResultsGFC"));
1003   TH1D *diffFlowResults2 = dynamic_cast<TH1D*>(outputListHistos->FindObject("fDiffFlowResults2ndOrderGFC"));
1004   TH1D *diffFlowResults4 = dynamic_cast<TH1D*>(outputListHistos->FindObject("fDiffFlowResults4thOrderGFC"));
1005   TH1D *diffFlowResults6 = dynamic_cast<TH1D*>(outputListHistos->FindObject("fDiffFlowResults6thOrderGFC"));
1006   TH1D *diffFlowResults8 = dynamic_cast<TH1D*>(outputListHistos->FindObject("fDiffFlowResults8thOrderGFC"));
1007                     
1008   //common histograms to store the final results  the integrated and differential flow
1009   AliFlowCommonHistResults *commonHistRes2nd = dynamic_cast<AliFlowCommonHistResults*>(outputListHistos->FindObject("AliFlowCommonHistResults2ndOrderGFC"));
1010   AliFlowCommonHistResults *commonHistRes4th = dynamic_cast<AliFlowCommonHistResults*>(outputListHistos->FindObject("AliFlowCommonHistResults4thOrderGFC"));
1011   AliFlowCommonHistResults *commonHistRes6th = dynamic_cast<AliFlowCommonHistResults*>(outputListHistos->FindObject("AliFlowCommonHistResults6thOrderGFC"));
1012   AliFlowCommonHistResults *commonHistRes8th = dynamic_cast<AliFlowCommonHistResults*>(outputListHistos->FindObject("AliFlowCommonHistResults8thOrderGFC"));
1013   
1014   //common control histogram
1015   AliFlowCommonHist *commonHists = dynamic_cast<AliFlowCommonHist*>(outputListHistos->FindObject("AliFlowCommonHistGFC"));
1016   
1017   //profiles with average values of generating functions for int. and diff. flow
1018   TProfile2D *intFlowGenFun    = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fIntFlowGenFun")); 
1019   
1020   TProfile2D *intFlowGenFun4   = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fIntFlowGenFun4"));  //only for other system of Eq.
1021   TProfile2D *intFlowGenFun6   = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fIntFlowGenFun6"));  //only for other system of Eq. 
1022   TProfile2D *intFlowGenFun8   = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fIntFlowGenFun8"));  //only for other system of Eq.
1023   TProfile2D *intFlowGenFun16  = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fIntFlowGenFun16")); //only for other system of Eq.  
1024   
1025   //RP, Pt:
1026   TProfile3D *diffFlowPtRPGenFunRe = dynamic_cast<TProfile3D*>(outputListHistos->FindObject("fDiffFlowPtRPGenFunRe"));
1027   TProfile3D *diffFlowPtRPGenFunIm = dynamic_cast<TProfile3D*>(outputListHistos->FindObject("fDiffFlowPtRPGenFunIm"));
1028   TProfile *ptBinRPNoOfParticles = dynamic_cast<TProfile*>(outputListHistos->FindObject("fPtBinRPNoOfParticles"));
1029   
1030   //RP, Eta:
1031   TProfile3D *diffFlowEtaRPGenFunRe = dynamic_cast<TProfile3D*>(outputListHistos->FindObject("fDiffFlowEtaRPGenFunRe"));
1032   TProfile3D *diffFlowEtaRPGenFunIm = dynamic_cast<TProfile3D*>(outputListHistos->FindObject("fDiffFlowEtaRPGenFunIm"));
1033   TProfile *etaBinRPNoOfParticles = dynamic_cast<TProfile*>(outputListHistos->FindObject("fEtaBinRPNoOfParticles"));
1034   
1035   //POI, Pt:
1036   TProfile3D *diffFlowPtPOIGenFunRe = dynamic_cast<TProfile3D*>(outputListHistos->FindObject("fDiffFlowPtPOIGenFunRe"));
1037   TProfile3D *diffFlowPtPOIGenFunIm = dynamic_cast<TProfile3D*>(outputListHistos->FindObject("fDiffFlowPtPOIGenFunIm"));
1038   TProfile *ptBinPOINoOfParticles = dynamic_cast<TProfile*>(outputListHistos->FindObject("fPtBinPOINoOfParticles"));
1039   
1040   //POI, Eta:
1041   TProfile3D *diffFlowEtaPOIGenFunRe = dynamic_cast<TProfile3D*>(outputListHistos->FindObject("fDiffFlowEtaPOIGenFunRe"));
1042   TProfile3D *diffFlowEtaPOIGenFunIm = dynamic_cast<TProfile3D*>(outputListHistos->FindObject("fDiffFlowEtaPOIGenFunIm"));
1043   TProfile *etaBinPOINoOfParticles = dynamic_cast<TProfile*>(outputListHistos->FindObject("fEtaBinPOINoOfParticles"));
1044   
1045   //average selected multiplicity (for int. flow) 
1046   TProfile *avMult = dynamic_cast<TProfile*>(outputListHistos->FindObject("fAvMultIntFlowGFC"));
1047   
1048   TProfile *avMult4  = dynamic_cast<TProfile*>(outputListHistos->FindObject("fAvMultIntFlow4GFCA"));  //only for other system of Eq.
1049   TProfile *avMult6  = dynamic_cast<TProfile*>(outputListHistos->FindObject("fAvMultIntFlow6GFCA"));  //only for other system of Eq.
1050   TProfile *avMult8  = dynamic_cast<TProfile*>(outputListHistos->FindObject("fAvMultIntFlow8GFCA"));  //only for other system of Eq.
1051   TProfile *avMult16 = dynamic_cast<TProfile*>(outputListHistos->FindObject("fAvMultIntFlow16GFCA")); //only for other system of Eq.
1052   
1053   //average values of Q-vector components (1st bin: <Q_x>, 2nd bin: <Q_y>, 3rd bin: <(Q_x)^2>, 4th bin: <(Q_y)^2>) 
1054   TProfile *qVectorComponents = dynamic_cast<TProfile*>(outputListHistos->FindObject("fQVectorComponentsGFC"));
1055   
1056   //<w^2> 
1057   TProfile *averageOfSquaredWeight = dynamic_cast<TProfile*>(outputListHistos->FindObject("fAverageOfSquaredWeight"));
1058       
1059   /*
1060   TProfile2D *diffFlowPtGenFunRe0 = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fdiffFlowPtGenFunRe0"));
1061   TProfile2D *diffFlowPtGenFunRe1 = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fdiffFlowPtGenFunRe1")); 
1062   TProfile2D *diffFlowPtGenFunRe2 = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fdiffFlowPtGenFunRe2")); 
1063   TProfile2D *diffFlowPtGenFunRe3 = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fdiffFlowPtGenFunRe3")); 
1064   TProfile2D *diffFlowPtGenFunRe4 = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fdiffFlowPtGenFunRe4")); 
1065   TProfile2D *diffFlowPtGenFunRe5 = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fdiffFlowPtGenFunRe5")); 
1066   TProfile2D *diffFlowPtGenFunRe6 = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fdiffFlowPtGenFunRe6")); 
1067   TProfile2D *diffFlowPtGenFunRe7 = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fdiffFlowPtGenFunRe7")); 
1068   TProfile2D *diffFlowPtGenFunIm0 = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fdiffFlowPtGenFunIm0")); 
1069   TProfile2D *diffFlowPtGenFunIm1 = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fdiffFlowPtGenFunIm1")); 
1070   TProfile2D *diffFlowPtGenFunIm2 = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fdiffFlowPtGenFunIm2")); 
1071   TProfile2D *diffFlowPtGenFunIm3 = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fdiffFlowPtGenFunIm3")); 
1072   TProfile2D *diffFlowPtGenFunIm4 = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fdiffFlowPtGenFunIm4")); 
1073   TProfile2D *diffFlowPtGenFunIm5 = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fdiffFlowPtGenFunIm5")); 
1074   TProfile2D *diffFlowPtGenFunIm6 = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fdiffFlowPtGenFunIm6")); 
1075   TProfile2D *diffFlowPtGenFunIm7 = dynamic_cast<TProfile2D*>(outputListHistos->FindObject("fdiffFlowPtGenFunIm7")); 
1076   */
1077
1078   //profile with avarage selected multiplicity for int. flow 
1079   //TProfile *avMult = dynamic_cast<TProfile*>(outputListHistos->FindObject("fAvMultIntFlow"));
1080   
1081   //profile with avarage values of Q-vector components (1st bin: <Q_x>, 2nd bin: <Q_y>, 3rd bin: <(Q_x)^2>, 4th bin: <(Q_y)^2>) 
1082   //TProfile *QVectorComponents = dynamic_cast<TProfile*>(outputListHistos->FindObject("fQVectorComponents"));
1083   
1084   //q-distribution
1085   //TH1D *qDist = dynamic_cast<TH1D*>(outputListHistos->FindObject("fQDist"));
1086   
1087   //AliCumulantsFunctions finalResults(intFlowGenFun,NULL,NULL, intFlowResults,diffFlowResults2,diffFlowResults4,diffFlowResults6,diffFlowResults8,avMult,QVectorComponents,qDist,diffFlowPtGenFunRe0,diffFlowPtGenFunRe1,diffFlowPtGenFunRe2, diffFlowPtGenFunRe3,diffFlowPtGenFunRe4,diffFlowPtGenFunRe5,diffFlowPtGenFunRe6,diffFlowPtGenFunRe7,diffFlowPtGenFunIm0,diffFlowPtGenFunIm1, diffFlowPtGenFunIm2,diffFlowPtGenFunIm3,diffFlowPtGenFunIm4,diffFlowPtGenFunIm5,diffFlowPtGenFunIm6,diffFlowPtGenFunIm7);
1088   
1089   //AliCumulantsFunctions finalResults(intFlowGenFun,diffFlowPtGenFunRe,diffFlowPtGenFunIm, intFlowResults,diffFlowResults2,diffFlowResults4,diffFlowResults6,diffFlowResults8,avMult,QVectorComponents,qDist);
1090          
1091   //finalResults.Calculate();  
1092   
1093   //---------------------------------------------------
1094   
1095   this->SetIntFlowResults(intFlowResults); 
1096   this->SetDiffFlowResults2nd(diffFlowResults2);
1097   this->SetDiffFlowResults4th(diffFlowResults4);
1098   this->SetDiffFlowResults6th(diffFlowResults6);
1099   this->SetDiffFlowResults8th(diffFlowResults8); 
1100   
1101   this->SetCommonHistsResults2nd(commonHistRes2nd); 
1102   this->SetCommonHistsResults4th(commonHistRes4th);
1103   this->SetCommonHistsResults6th(commonHistRes6th);
1104   this->SetCommonHistsResults8th(commonHistRes8th);
1105   
1106   this->SetCommonHists(commonHists);
1107   
1108   this->SetIntFlowGenFun(intFlowGenFun);
1109   
1110   this->SetIntFlowGenFun4(intFlowGenFun4);   //only for other system of Eq.
1111   this->SetIntFlowGenFun6(intFlowGenFun6);   //only for other system of Eq.
1112   this->SetIntFlowGenFun8(intFlowGenFun8);   //only for other system of Eq.
1113   this->SetIntFlowGenFun16(intFlowGenFun16); //only for other system of Eq. 
1114   
1115   this->SetDiffFlowPtRPGenFunRe(diffFlowPtRPGenFunRe);
1116   this->SetDiffFlowPtRPGenFunIm(diffFlowPtRPGenFunIm);
1117   this->SetNumberOfParticlesPerPtBinRP(ptBinRPNoOfParticles);
1118   
1119   this->SetDiffFlowEtaRPGenFunRe(diffFlowEtaRPGenFunRe);
1120   this->SetDiffFlowEtaRPGenFunIm(diffFlowEtaRPGenFunIm);
1121   this->SetNumberOfParticlesPerEtaBinRP(etaBinRPNoOfParticles);     
1122   
1123   this->SetDiffFlowPtPOIGenFunRe(diffFlowPtPOIGenFunRe);
1124   this->SetDiffFlowPtPOIGenFunIm(diffFlowPtPOIGenFunIm);
1125   this->SetNumberOfParticlesPerPtBinPOI(ptBinPOINoOfParticles);
1126   
1127   this->SetDiffFlowEtaPOIGenFunRe(diffFlowEtaPOIGenFunRe);
1128   this->SetDiffFlowEtaPOIGenFunIm(diffFlowEtaPOIGenFunIm);
1129   this->SetNumberOfParticlesPerEtaBinPOI(etaBinPOINoOfParticles);
1130   
1131   this->SetAverageMultiplicity(avMult);
1132   
1133   this->SetAverageMultiplicity4(avMult4);   //only for other system of Eq.
1134   this->SetAverageMultiplicity6(avMult6);   //only for other system of Eq.
1135   this->SetAverageMultiplicity8(avMult8);   //only for other system of Eq.
1136   this->SetAverageMultiplicity16(avMult16); //only for other system of Eq.
1137   
1138   this->SetQVectorComponents(qVectorComponents);
1139   
1140   this->SetAverageOfSquaredWeight(averageOfSquaredWeight);
1141  } 
1142 }
1143
1144 //================================================================================================================
1145
1146 void AliFlowAnalysisWithCumulants::Finish()
1147 {
1148  //calculate the final results
1149  //AliCumulantsFunctions finalResults(fIntFlowGenFun,NULL,NULL, fIntFlowResults,fDiffFlowResults2,fDiffFlowResults4,fDiffFlowResults6,fDiffFlowResults8,fAvMultIntFlow,fQVectorComponents,  fQDist,fDiffFlowGenFunRe0,fDiffFlowGenFunRe1,fDiffFlowGenFunRe2, fDiffFlowGenFunRe3,fDiffFlowGenFunRe4,fDiffFlowGenFunRe5,fDiffFlowGenFunRe6,fDiffFlowGenFunRe7,fDiffFlowGenFunIm0,fDiffFlowGenFunIm1, fDiffFlowGenFunIm2,fDiffFlowGenFunIm3,fDiffFlowGenFunIm4,fDiffFlowGenFunIm5,fDiffFlowGenFunIm6,fDiffFlowGenFunIm7);
1150
1151  AliCumulantsFunctions finalResults(fIntFlowGenFun,fIntFlowGenFun4,fIntFlowGenFun6,fIntFlowGenFun8,fIntFlowGenFun16,fAvMultIntFlow4GFC, fAvMultIntFlow6GFC,fAvMultIntFlow8GFC,fAvMultIntFlow16GFC,fDiffFlowPtRPGenFunRe,fDiffFlowPtRPGenFunIm,fPtBinRPNoOfParticles, fDiffFlowEtaRPGenFunRe,fDiffFlowEtaRPGenFunIm,fEtaBinRPNoOfParticles,fDiffFlowPtPOIGenFunRe,fDiffFlowPtPOIGenFunIm,fPtBinPOINoOfParticles, fDiffFlowEtaPOIGenFunRe,fDiffFlowEtaPOIGenFunIm,fEtaBinPOINoOfParticles,fIntFlowResultsGFC,fDiffFlowResults2ndOrderGFC,fDiffFlowResults4thOrderGFC,fDiffFlowResults6thOrderGFC,fDiffFlowResults8thOrderGFC, fAvMultIntFlowGFC,fQVectorComponentsGFC,fAverageOfSquaredWeight,fCommonHistsResults2nd, fCommonHistsResults4th,fCommonHistsResults6th,fCommonHistsResults8th,fCommonHists);
1152                            
1153  finalResults.Calculate();     
1154 }
1155
1156 //================================================================================================================
1157
1158 void AliFlowAnalysisWithCumulants::WriteHistograms(TString* outputFileName)
1159 {
1160  //store the final results in output .root file
1161  TFile *output = new TFile(outputFileName->Data(),"RECREATE");
1162  //output->WriteObject(fHistList, "cobjGFC","SingleKey");
1163  fHistList->SetName("cobjGFC");
1164  fHistList->SetOwner(kTRUE);
1165  fHistList->Write(fHistList->GetName(), TObject::kSingleKey);
1166  delete output;
1167 }
1168
1169 //================================================================================================================
1170
1171 //================================================================================================================
1172
1173 void AliFlowAnalysisWithCumulants::WriteHistograms(TString outputFileName)
1174 {
1175  //store the final results in output .root file
1176  TFile *output = new TFile(outputFileName.Data(),"RECREATE");
1177  //output->WriteObject(fHistList, "cobjGFC","SingleKey");
1178  fHistList->SetName("cobjGFC");
1179  fHistList->SetOwner(kTRUE);
1180  fHistList->Write(fHistList->GetName(), TObject::kSingleKey);
1181  delete output;
1182 }
1183
1184 //================================================================================================================
1185
1186 //================================================================================================================
1187
1188 void AliFlowAnalysisWithCumulants::WriteHistograms(TDirectoryFile *outputFileName)
1189 {
1190  //store the final results in output .root file
1191  fHistList->SetName("cobjGFC");
1192  fHistList->SetOwner(kTRUE);
1193  outputFileName->Add(fHistList);
1194  outputFileName->Write(outputFileName->GetName(), TObject::kSingleKey);
1195 }
1196
1197 //================================================================================================================
1198
1199