1 /*************************************************************************
2 * Copyright(c) 1998-2008, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 /**********************************
17 * flow analysis with Q-cumulants *
19 * author: Ante Bilandzic *
21 *********************************/
23 #define AliFlowAnalysisWithQCumulants_cxx
25 #include "Riostream.h"
26 #include "AliFlowCommonConstants.h"
27 #include "AliFlowCommonHist.h"
28 #include "AliFlowCommonHistResults.h"
32 #include "TParticle.h"
36 #include "TProfile2D.h"
37 #include "TProfile3D.h"
38 #include "AliFlowEventSimple.h"
39 #include "AliFlowTrackSimple.h"
40 #include "AliFlowAnalysisWithQCumulants.h"
41 #include "AliQCumulantsFunctions.h"
59 //================================================================================================================
61 ClassImp(AliFlowAnalysisWithQCumulants)
63 AliFlowAnalysisWithQCumulants::AliFlowAnalysisWithQCumulants():
66 fAvMultIntFlowQC(NULL),
67 fQvectorComponents(NULL),
68 fIntFlowResultsQC(NULL),
69 fDiffFlowResults2ndOrderQC(NULL),
70 fDiffFlowResults4thOrderQC(NULL),
74 fDirectCorrelations(NULL),
83 f4PerBin1n1n1n1n(NULL),
85 fCommonHistsResults2nd(NULL),
86 fCommonHistsResults4th(NULL),
87 fCommonHistsResults6th(NULL),
88 fCommonHistsResults8th(NULL),
89 f2pDistribution(NULL),
90 f4pDistribution(NULL),
91 f6pDistribution(NULL),
92 f8pDistribution(NULL),
98 fHistList = new TList();
100 fnBinsPt = AliFlowCommonConstants::GetNbinsPt();
101 fPtMin = AliFlowCommonConstants::GetPtMin();
102 fPtMax = AliFlowCommonConstants::GetPtMax();
106 AliFlowAnalysisWithQCumulants::~AliFlowAnalysisWithQCumulants()
112 //================================================================================================================
114 void AliFlowAnalysisWithQCumulants::CreateOutputObjects()
116 //various output histograms
118 //avarage multiplicity
119 fAvMultIntFlowQC = new TProfile("fAvMultIntFlowQC","Average Multiplicity",1,0,1,"s");
120 fAvMultIntFlowQC->SetXTitle("");
121 fAvMultIntFlowQC->SetYTitle("");
122 fAvMultIntFlowQC->SetLabelSize(0.06);
123 fAvMultIntFlowQC->SetMarkerStyle(25);
124 fAvMultIntFlowQC->SetLabelOffset(0.01);
125 (fAvMultIntFlowQC->GetXaxis())->SetBinLabel(1,"Average Multiplicity");
126 fHistList->Add(fAvMultIntFlowQC);
129 fQvectorComponents = new TProfile("fQvectorComponents","Avarage of Q-vector components",44,0.,44.,"s");
130 fQvectorComponents->SetXTitle("");
131 fQvectorComponents->SetYTitle("");
132 //fHistList->Add(fQvectorComponents);
134 //final results for integrated flow from Q-cumulants
135 fIntFlowResultsQC = new TH1D("fIntFlowResultsQC","Integrated Flow from Q-cumulants",4,0,4);
136 //fIntFlowResults->SetXTitle("");
137 //fIntFlowResultsQC->SetYTitle("Integrated Flow");
138 fIntFlowResultsQC->SetLabelSize(0.06);
139 //fIntFlowResultsQC->SetTickLength(1);
140 fIntFlowResultsQC->SetMarkerStyle(25);
141 (fIntFlowResultsQC->GetXaxis())->SetBinLabel(1,"v_{n}{2}");
142 (fIntFlowResultsQC->GetXaxis())->SetBinLabel(2,"v_{n}{4}");
143 (fIntFlowResultsQC->GetXaxis())->SetBinLabel(3,"v_{n}{6}");
144 (fIntFlowResultsQC->GetXaxis())->SetBinLabel(4,"v_{n}{8}");
145 fHistList->Add(fIntFlowResultsQC);
147 //final results for differential flow from 2nd order Q-cumulant
148 fDiffFlowResults2ndOrderQC = new TH1D("fDiffFlowResults2ndOrderQC","Differential Flow from 2nd Order Q-cumulant",fnBinsPt,fPtMin,fPtMax);
149 fDiffFlowResults2ndOrderQC->SetXTitle("p_{t} [GeV]");
150 //fDiffFlowResults2ndOrderQC->SetYTitle("Differential Flow");
151 fHistList->Add(fDiffFlowResults2ndOrderQC);
153 //final results for differential flow from 4th order Q-cumulant
154 fDiffFlowResults4thOrderQC = new TH1D("fDiffFlowResults4thOrderQC","Differential Flow from 4th Order Q-cumulant",fnBinsPt,fPtMin,fPtMax);
155 fDiffFlowResults4thOrderQC->SetXTitle("p_{t} [GeV]");
156 //fDiffFlowResults4thOrderQC->SetYTitle("Differential Flow");
157 fHistList->Add(fDiffFlowResults4thOrderQC);
159 //final results for covariances (1st bin: <2*4>-<2>*<4>, 2nd bin: <2*6>-<2>*<6>, ...)
160 fCovariances = new TH1D("fCovariances","Covariances",6,0,6);
161 //fCovariances->SetXTitle("");
162 //fCovariances->SetYTitle("<covariance>");
163 fCovariances->SetLabelSize(0.04);
164 fCovariances->SetTickLength(1);
165 fCovariances->SetMarkerStyle(25);
166 (fCovariances->GetXaxis())->SetBinLabel(1,"Cov(2,4)");
167 (fCovariances->GetXaxis())->SetBinLabel(2,"Cov(2,6)");
168 (fCovariances->GetXaxis())->SetBinLabel(3,"Cov(2,8)");
169 (fCovariances->GetXaxis())->SetBinLabel(4,"Cov(4,6)");
170 (fCovariances->GetXaxis())->SetBinLabel(5,"Cov(4,8)");
171 (fCovariances->GetXaxis())->SetBinLabel(6,"Cov(6,8)");
172 fHistList->Add(fCovariances);
174 //multi-particle correlations calculated from Q-vectors
175 fQCorrelations = new TProfile("fQCorrelations","multi-particle correlations from Q-vectors",32,0,32,"s");
176 //fQCorrelations->SetXTitle("correlations");
177 //fQCorrelations->SetYTitle("");
178 fQCorrelations->SetTickLength(-0.01,"Y");
179 fQCorrelations->SetMarkerStyle(25);
180 fQCorrelations->SetLabelSize(0.03);
181 fQCorrelations->SetLabelOffset(0.01,"Y");
183 (fQCorrelations->GetXaxis())->SetBinLabel(1,"<<2>>_{n|n}");
184 (fQCorrelations->GetXaxis())->SetBinLabel(2,"<<2>>_{2n|2n}");
185 (fQCorrelations->GetXaxis())->SetBinLabel(3,"<<2>>_{3n|3n}");
186 (fQCorrelations->GetXaxis())->SetBinLabel(4,"<<2>>_{4n|4n}");
188 (fQCorrelations->GetXaxis())->SetBinLabel(6,"<<3>>_{2n|n,n}");
189 (fQCorrelations->GetXaxis())->SetBinLabel(7,"<<3>>_{3n|2n,n}");
190 (fQCorrelations->GetXaxis())->SetBinLabel(8,"<<3>>_{4n|2n,2n}");
191 (fQCorrelations->GetXaxis())->SetBinLabel(9,"<<3>>_{4n|3n,n}");
193 (fQCorrelations->GetXaxis())->SetBinLabel(11,"<<4>>_{n,n|n,n}");
194 (fQCorrelations->GetXaxis())->SetBinLabel(12,"<<4>>_{2n,n|2n,n}");
195 (fQCorrelations->GetXaxis())->SetBinLabel(13,"<<4>>_{2n,2n|2n,2n}");
196 (fQCorrelations->GetXaxis())->SetBinLabel(14,"<<4>>_{3n|n,n,n}");
197 (fQCorrelations->GetXaxis())->SetBinLabel(15,"<<4>>_{3n,n|3n,n}");
198 (fQCorrelations->GetXaxis())->SetBinLabel(16,"<<4>>_{3n,n|2n,2n}");
199 (fQCorrelations->GetXaxis())->SetBinLabel(17,"<<4>>_{4n|2n,n,n}");
201 (fQCorrelations->GetXaxis())->SetBinLabel(19,"<<5>>_{2n|n,n,n,n}");
202 (fQCorrelations->GetXaxis())->SetBinLabel(20,"<<5>>_{2n,2n|2n,n,n}");
203 (fQCorrelations->GetXaxis())->SetBinLabel(21,"<<5>>_{3n,n|2n,n,n}");
204 (fQCorrelations->GetXaxis())->SetBinLabel(22,"<<5>>_{4n|n,n,n,n}");
206 (fQCorrelations->GetXaxis())->SetBinLabel(24,"<<6>>_{n,n,n|n,n,n}");
207 (fQCorrelations->GetXaxis())->SetBinLabel(25,"<<6>>_{2n,n,n|2n,n,n}");
208 (fQCorrelations->GetXaxis())->SetBinLabel(26,"<<6>>_{2n,2n|n,n,n,n}");
209 (fQCorrelations->GetXaxis())->SetBinLabel(27,"<<6>>_{3n,n|n,n,n,n}");
211 (fQCorrelations->GetXaxis())->SetBinLabel(29,"<<7>>_{2n,n,n|n,n,n,n}");
213 (fQCorrelations->GetXaxis())->SetBinLabel(31,"<<8>>_{n,n,n,n|n,n,n,n}");
215 fHistList->Add(fQCorrelations);
218 fQProduct = new TProfile("fQProduct","average of products",6,0,6,"s");
219 fQProduct->SetTickLength(-0.01,"Y");
220 fQProduct->SetMarkerStyle(25);
221 fQProduct->SetLabelSize(0.03);
222 fQProduct->SetLabelOffset(0.01,"Y");
223 (fQProduct->GetXaxis())->SetBinLabel(1,"<<2*4>>");
224 (fQProduct->GetXaxis())->SetBinLabel(2,"<<2*6>>");
225 (fQProduct->GetXaxis())->SetBinLabel(3,"<<2*8>>");
226 (fQProduct->GetXaxis())->SetBinLabel(4,"<<4*6>>");
227 (fQProduct->GetXaxis())->SetBinLabel(5,"<<4*8>>");
228 (fQProduct->GetXaxis())->SetBinLabel(6,"<<6*8>>");
229 fQProduct->SetXTitle("");
230 fQProduct->SetYTitle("");
231 fHistList->Add(fQProduct);
233 //multi-particle correlations calculated with nested loops (0..40 integrated flow; 40..80 differential flow)
234 fDirectCorrelations = new TProfile("fDirectCorrelations","multi-particle correlations with nested loops",80,0,80,"s");
235 fDirectCorrelations->SetXTitle("");
236 fDirectCorrelations->SetYTitle("correlations");
237 fHistList->Add(fDirectCorrelations);
240 fReq1n = new TProfile("fReq1n","Re[q_n]",fnBinsPt,fPtMin,fPtMax,"s");
241 fReq1n->SetXTitle("p_{t} [GeV]");
242 fReq1n->SetYTitle("Re[q_n]");
243 //fHistList->Add(fReq1n);
246 fImq1n = new TProfile("fImq1n","Im[q_n]",fnBinsPt,fPtMin,fPtMax,"s");
247 fImq1n->SetXTitle("p_{t} [GeV]");
248 fImq1n->SetYTitle("Im[q_n]");
249 //fHistList->Add(fImq1n);
252 fReq2n = new TProfile("fReq2n","Re[q_2n]",fnBinsPt,fPtMin,fPtMax,"s");
253 fReq2n->SetXTitle("p_{t} [GeV]");
254 fReq2n->SetYTitle("Im[D]");
255 //fHistList->Add(fReq2n);
258 fImq2n = new TProfile("fImq2n","Im[q_2n]",fnBinsPt,fPtMin,fPtMax,"s");
259 fImq2n->SetXTitle("p_{t} [GeV]");
260 fImq2n->SetYTitle("Im[q_2n]");
261 //fHistList->Add(fImq2n);
264 f2PerBin1n1n = new TProfile("f2PerBin1n1n","<2'>_{n|n}",fnBinsPt,fPtMin,fPtMax,"s");
265 f2PerBin1n1n->SetXTitle("p_{t} [GeV]");
266 //f2PerBin1n1n->SetYTitle("<2'>_{n|n}");
267 fHistList->Add(f2PerBin1n1n);
270 f2PerBin2n2n = new TProfile("f2PerBin2n2n","<2'>_{2n|2n}",fnBinsPt,fPtMin,fPtMax,"s");
271 f2PerBin2n2n->SetXTitle("p_{t} [GeV]");
272 //f2PerBin2n2n->SetYTitle("<2'>_{2n|2n}");
273 fHistList->Add(f2PerBin2n2n);
276 f3PerBin2n1n1n = new TProfile("f3PerBin2n1n1n","<3'>_{2n|n,n}",fnBinsPt,fPtMin,fPtMax,"s");
277 f3PerBin2n1n1n->SetXTitle("p_{t} [GeV]");
278 //f3PerBin2n1n1n->SetYTitle("<3'>_{2n|n,n}");
279 fHistList->Add(f3PerBin2n1n1n);
282 f3PerBin1n1n2n = new TProfile("f3PerBin1n1n2n","<3'>_{n,n|2n}",fnBinsPt,fPtMin,fPtMax,"s");
283 f3PerBin1n1n2n->SetXTitle("p_{t} [GeV]");
284 //f3PerBin1n1n2n->SetYTitle("<3'>_{n,n|2n}");
285 fHistList->Add(f3PerBin1n1n2n);
288 f4PerBin1n1n1n1n = new TProfile("f4PerBin1n1n1n1n","<4'>_{n,n|n,n}",fnBinsPt,fPtMin,fPtMax,"s");
289 f4PerBin1n1n1n1n->SetXTitle("p_{t} [GeV]");
290 //f4PerBin1n1n1n1n->SetYTitle("<4'>_{n,n|n,n}");
291 fHistList->Add(f4PerBin1n1n1n1n);
293 //common control histograms
294 fCommonHists = new AliFlowCommonHist("AliFlowCommonHistQC");
295 fHistList->Add(fCommonHists);
297 //common histograms for final results (2nd order)
298 fCommonHistsResults2nd = new AliFlowCommonHistResults("AliFlowCommonHistResults2ndOrderQC");
299 fHistList->Add(fCommonHistsResults2nd);
301 //common histograms for final results (4th order)
302 fCommonHistsResults4th = new AliFlowCommonHistResults("AliFlowCommonHistResults4thOrderQC");
303 fHistList->Add(fCommonHistsResults4th);
305 //common histograms for final results (6th order)
306 fCommonHistsResults6th = new AliFlowCommonHistResults("AliFlowCommonHistResults6thOrderQC");
307 fHistList->Add(fCommonHistsResults6th);
309 //common histograms for final results (8th order)
310 fCommonHistsResults8th = new AliFlowCommonHistResults("AliFlowCommonHistResults8thOrderQC");
311 fHistList->Add(fCommonHistsResults8th);
313 //weighted <2>_{n|n} distribution
314 f2pDistribution = new TH1D("f2pDistribution","<2>_{n|n} distribution",100000,-0.02,0.1);
315 f2pDistribution->SetXTitle("<2>_{n|n}");
316 f2pDistribution->SetYTitle("Counts");
317 fHistList->Add(f2pDistribution);
319 //weighted <4>_{n,n|n,n} distribution
320 f4pDistribution = new TH1D("f4pDistribution","<4>_{n,n|n,n} distribution",100000,-0.00025,0.002);
321 f4pDistribution->SetXTitle("<4>_{n,n|n,n}");
322 f4pDistribution->SetYTitle("Counts");
323 fHistList->Add(f4pDistribution);
325 //weighted <6>_{n,n,n|n,n,n} distribution
326 f6pDistribution = new TH1D("f6pDistribution","<6>_{n,n,n|n,n,n} distribution",100000,-0.000005,0.000025);
327 f6pDistribution->SetXTitle("<6>_{n,n,n|n,n,n}");
328 f6pDistribution->SetYTitle("Counts");
329 fHistList->Add(f6pDistribution);
331 //weighted <8>_{n,n,n,n|n,n,n,n} distribution
332 f8pDistribution = new TH1D("f8pDistribution","<8>_{n,n,n,n|n,n,n,n} distribution",100000,-0.000000001,0.00000001);
333 f8pDistribution->SetXTitle("<8>_{n,n,n,n|n,n,n,n}");
334 f8pDistribution->SetYTitle("Counts");
335 fHistList->Add(f8pDistribution);
337 }//end of CreateOutputObjects()
339 //================================================================================================================
341 void AliFlowAnalysisWithQCumulants::Make(AliFlowEventSimple* anEvent)
345 //get the total multiplicity of event:
346 //Int_t nPrim = anEvent->NumberOfTracks();//line needed only for nested loops
348 //if(nPrim>8&&nPrim<12)//line needed only for nested loops
349 //{ //line needed only for nested loops
351 //fill the common control histograms:
352 fCommonHists->FillControlHistograms(anEvent);
354 //get the selected multiplicity (i.e. number of particles used for int. flow):
355 //Int_t nEventNSelTracksIntFlow = anEvent->GetEventNSelTracksIntFlow();
357 Int_t n=2; //int flow harmonic (to be improved)
359 //---------------------------------------------------------------------------------------------------------
360 //Q-vectors of an event evaluated in harmonics n, 2n, 3n and 4n:
361 AliFlowVector xQvector1n, xQvector2n, xQvector3n, xQvector4n;
363 xQvector1n.Set(0.,0.);
364 xQvector1n.SetMult(0);
365 xQvector1n=anEvent->GetQ(1*n);
367 xQvector2n.Set(0.,0.);
368 xQvector2n.SetMult(0);
369 xQvector2n=anEvent->GetQ(2*n);
371 xQvector3n.Set(0.,0.);
372 xQvector3n.SetMult(0);
373 xQvector3n=anEvent->GetQ(3*n);
375 xQvector4n.Set(0.,0.);
376 xQvector4n.SetMult(0);
377 xQvector4n=anEvent->GetQ(4*n);
378 //---------------------------------------------------------------------------------------------------------
380 //multiplicity (to be improved, because I already have nEventNSelTracksIntFlow and nPrim)
381 Double_t xMult = xQvector1n.GetMult();
383 fAvMultIntFlowQC->Fill(0.,xMult,1.);
385 //---------------------------------------------------------------------------------------------------------
387 // *******************
388 // **** Q-vectors ****
389 // *******************
391 Double_t reQ2nQ1nstarQ1nstar = pow(xQvector1n.X(),2.)*xQvector2n.X()+2.*xQvector1n.X()*xQvector1n.Y()*xQvector2n.Y()-pow(xQvector1n.Y(),2.)*xQvector2n.X();//Re[Q_{2n} Q_{n}^* Q_{n}^*]
392 //Double_t imQ2nQ1nstarQ1nstar = pow(Qvector1n.X(),2.)*Qvector2n.Y()-2.*Qvector1n.X()*Qvector1n.Y()*Qvector2n.X()-pow(Qvector1n.Y(),2.)*Qvector2n.Y();//Im[Q_{2n} Q_{n}^* Q_{n}^*]
393 Double_t reQ1nQ1nQ2nstar = reQ2nQ1nstarQ1nstar;//Re[Q_{n} Q_{n} Q_{2n}^*] = Re[Q_{2n} Q_{n}^* Q_{n}^*]
394 Double_t reQ3nQ1nQ2nstarQ2nstar = (pow(xQvector2n.X(),2.)-pow(xQvector2n.Y(),2.))*(xQvector3n.X()*xQvector1n.X()-xQvector3n.Y()*xQvector1n.Y())+2.*xQvector2n.X()*xQvector2n.Y()*(xQvector3n.X()*xQvector1n.Y()+xQvector3n.Y()*xQvector1n.X());
395 //Double_t imQ3nQ1nQ2nstarQ2nstar = calculate and implement this (deleteMe)
396 Double_t reQ2nQ2nQ3nstarQ1nstar = reQ3nQ1nQ2nstarQ2nstar;
397 Double_t reQ4nQ2nstarQ2nstar = pow(xQvector2n.X(),2.)*xQvector4n.X()+2.*xQvector2n.X()*xQvector2n.Y()*xQvector4n.Y()-pow(xQvector2n.Y(),2.)*xQvector4n.X();//Re[Q_{4n} Q_{2n}^* Q_{2n}^*]
398 //Double_t imQ4nQ2nstarQ2nstar = calculate and implement this (deleteMe)
399 Double_t reQ2nQ2nQ4nstar = reQ4nQ2nstarQ2nstar;
400 Double_t reQ4nQ3nstarQ1nstar = xQvector4n.X()*(xQvector3n.X()*xQvector1n.X()-xQvector3n.Y()*xQvector1n.Y())+xQvector4n.Y()*(xQvector3n.X()*xQvector1n.Y()+xQvector3n.Y()*xQvector1n.X());//Re[Q_{4n} Q_{3n}^* Q_{n}^*]
401 Double_t reQ3nQ1nQ4nstar = reQ4nQ3nstarQ1nstar;//Re[Q_{3n} Q_{n} Q_{4n}^*] = Re[Q_{4n} Q_{3n}^* Q_{n}^*]
402 //Double_t imQ4nQ3nstarQ1nstar = calculate and implement this (deleteMe)
403 Double_t reQ3nQ2nstarQ1nstar = xQvector3n.X()*xQvector2n.X()*xQvector1n.X()-xQvector3n.X()*xQvector2n.Y()*xQvector1n.Y()+xQvector3n.Y()*xQvector2n.X()*xQvector1n.Y()+xQvector3n.Y()*xQvector2n.Y()*xQvector1n.X();//Re[Q_{3n} Q_{2n}^* Q_{n}^*]
404 Double_t reQ2nQ1nQ3nstar = reQ3nQ2nstarQ1nstar;//Re[Q_{2n} Q_{n} Q_{3n}^*] = Re[Q_{3n} Q_{2n}^* Q_{n}^*]
405 //Double_t imQ3nQ2nstarQ1nstar; //calculate and implement this (deleteMe)
406 Double_t reQ3nQ1nstarQ1nstarQ1nstar = xQvector3n.X()*pow(xQvector1n.X(),3)-3.*xQvector1n.X()*xQvector3n.X()*pow(xQvector1n.Y(),2)+3.*xQvector1n.Y()*xQvector3n.Y()*pow(xQvector1n.X(),2)-xQvector3n.Y()*pow(xQvector1n.Y(),3);//Re[Q_{3n} Q_{n}^* Q_{n}^* Q_{n}^*]
407 //Double_t imQ3nQ1nstarQ1nstarQ1nstar; //calculate and implement this (deleteMe)
408 Double_t xQ2nQ1nQ2nstarQ1nstar = pow(xQvector2n.Mod()*xQvector1n.Mod(),2);//|Q_{2n}|^2 |Q_{n}|^2
409 Double_t reQ4nQ2nstarQ1nstarQ1nstar = (xQvector4n.X()*xQvector2n.X()+xQvector4n.Y()*xQvector2n.Y())*(pow(xQvector1n.X(),2)-pow(xQvector1n.Y(),2))+2.*xQvector1n.X()*xQvector1n.Y()*(xQvector4n.Y()*xQvector2n.X()-xQvector4n.X()*xQvector2n.Y());//Re[Q_{4n} Q_{2n}^* Q_{n}^* Q_{n}^*]
410 //Double_t imQ4nQ2nstarQ1nstarQ1nstar; //calculate and implement this (deleteMe)
411 Double_t reQ2nQ1nQ1nstarQ1nstarQ1nstar = (xQvector2n.X()*xQvector1n.X()-xQvector2n.Y()*xQvector1n.Y())*(pow(xQvector1n.X(),3)-3.*xQvector1n.X()*pow(xQvector1n.Y(),2))+(xQvector2n.X()*xQvector1n.Y()+xQvector1n.X()*xQvector2n.Y())*(3.*xQvector1n.Y()*pow(xQvector1n.X(),2)-pow(xQvector1n.Y(),3));//Re[Q_{2n} Q_{n} Q_{n}^* Q_{n}^* Q_{n}^*]
412 //Double_t imQ2nQ1nQ1nstarQ1nstarQ1nstar; //calculate and implement this (deleteMe)
413 Double_t reQ2nQ2nQ2nstarQ1nstarQ1nstar = pow(xQvector2n.Mod(),2.)*(xQvector2n.X()*(pow(xQvector1n.X(),2.)-pow(xQvector1n.Y(),2.))+2.*xQvector2n.Y()*xQvector1n.X()*xQvector1n.Y());//Re[Q_{2n} Q_{2n} Q_{2n}^* Q_{n}^* Q_{n}^*]
414 //Double_t imQ2nQ2nQ2nstarQ1nstarQ1nstar = pow(Qvector2n.Mod(),2.)*(Qvector2n.Y()*(pow(Qvector1n.X(),2.)-pow(Qvector1n.Y(),2.))-2.*Qvector2n.X()*Qvector1n.X()*Qvector1n.Y());//Im[Q_{2n} Q_{2n} Q_{2n}^* Q_{n}^* Q_{n}^*]
415 Double_t reQ4nQ1nstarQ1nstarQ1nstarQ1nstar = pow(xQvector1n.X(),4.)*xQvector4n.X()-6.*pow(xQvector1n.X(),2.)*xQvector4n.X()*pow(xQvector1n.Y(),2.)+pow(xQvector1n.Y(),4.)*xQvector4n.X()+4.*pow(xQvector1n.X(),3.)*xQvector1n.Y()*xQvector4n.Y()-4.*pow(xQvector1n.Y(),3.)*xQvector1n.X()*xQvector4n.Y();//Re[Q_{4n} Q_{n}^* Q_{n}^* Q_{n}^* Q_{n}^*]
416 //Double_t imQ4nQ1nstarQ1nstarQ1nstarQ1nstar = pow(Qvector1n.X(),4.)*Qvector4n.Y()-6.*pow(Qvector1n.X(),2.)*Qvector4n.Y()*pow(Qvector1n.Y(),2.)+pow(Qvector1n.Y(),4.)*Qvector4n.Y()+4.*pow(Qvector1n.Y(),3.)*Qvector1n.X()*Qvector4n.X()-4.*pow(Qvector1n.X(),3.)*Qvector1n.Y()*Qvector4n.X();//Im[Q_{4n} Q_{n}^* Q_{n}^* Q_{n}^* Q_{n}^*]
417 Double_t reQ3nQ1nQ2nstarQ1nstarQ1nstar = pow(xQvector1n.Mod(),2.)*(xQvector1n.X()*xQvector2n.X()*xQvector3n.X()-xQvector3n.X()*xQvector1n.Y()*xQvector2n.Y()+xQvector2n.X()*xQvector1n.Y()*xQvector3n.Y()+xQvector1n.X()*xQvector2n.Y()*xQvector3n.Y());//Re[Q_{3n} Q_{n} Q_{2n}^* Q_{n}^* Q_{n}^*]
418 //Double_t imQ3nQ1nQ2nstarQ1nstarQ1nstar = pow(xQvector1n.Mod(),2.)*(-xQvector2n.X()*xQvector3n.X()*xQvector1n.Y()-xQvector1n.X()*xQvector3n.X()*xQvector2n.Y()+xQvector1n.X()*xQvector2n.X()*xQvector3n.Y()-xQvector1n.Y()*xQvector2n.Y()*xQvector3n.Y());//Im[Q_{3n} Q_{n} Q_{2n}^* Q_{n}^* Q_{n}^*]
419 Double_t reQ2nQ2nQ1nstarQ1nstarQ1nstarQ1nstar = (pow(xQvector1n.X(),2.)*xQvector2n.X()-2.*xQvector1n.X()*xQvector2n.X()*xQvector1n.Y()-xQvector2n.X()*pow(xQvector1n.Y(),2.)+xQvector2n.Y()*pow(xQvector1n.X(),2.)+2.*xQvector1n.X()*xQvector1n.Y()*xQvector2n.Y()-pow(xQvector1n.Y(),2.)*xQvector2n.Y())*(pow(xQvector1n.X(),2.)*xQvector2n.X()+2.*xQvector1n.X()*xQvector2n.X()*xQvector1n.Y()-xQvector2n.X()*pow(xQvector1n.Y(),2.)-xQvector2n.Y()*pow(xQvector1n.X(),2.)+2.*xQvector1n.X()*xQvector1n.Y()*xQvector2n.Y()+pow(xQvector1n.Y(),2.)*xQvector2n.Y());//Re[Q_{2n} Q_{2n} Q_{n}^* Q_{n}^* Q_{n}^* Q_{n}^*]
420 //Double_t imQ2nQ2nQ1nstarQ1nstarQ1nstarQ1nstar = 2.*(pow(xQvector1n.X(),2.)*xQvector2n.X()-xQvector2n.X()*pow(xQvector1n.Y(),2.)+2.*xQvector1n.X()*xQvector1n.Y()*xQvector2n.Y())*(pow(xQvector1n.X(),2.)*xQvector2n.Y()-2.*xQvector1n.X()*xQvector1n.Y()*xQvector2n.X()-pow(xQvector1n.Y(),2.)*xQvector2n.Y());//Im[Q_{2n} Q_{2n} Q_{n}^* Q_{n}^* Q_{n}^* Q_{n}^*]
421 Double_t reQ3nQ1nQ1nstarQ1nstarQ1nstarQ1nstar = pow(xQvector1n.Mod(),2.)*(pow(xQvector1n.X(),3.)*xQvector3n.X()-3.*xQvector1n.X()*xQvector3n.X()*pow(xQvector1n.Y(),2.)+3.*pow(xQvector1n.X(),2.)*xQvector1n.Y()*xQvector3n.Y()-pow(xQvector1n.Y(),3.)*xQvector3n.Y());//Re[Q_{3n} Q_{n} Q_{n}^* Q_{n}^* Q_{n}^* Q_{n}^*]
422 //Double_t imQ3nQ1nQ1nstarQ1nstarQ1nstarQ1nstar = pow(xQvector1n.Mod(),2.)*(pow(xQvector1n.Y(),3.)*xQvector3n.X()-3.*xQvector1n.Y()*xQvector3n.X()*pow(xQvector1n.X(),2.)-3.*pow(xQvector1n.Y(),2.)*xQvector1n.X()*xQvector3n.Y()+pow(xQvector1n.X(),3.)*xQvector3n.Y());//Im[Q_{3n} Q_{n} Q_{n}^* Q_{n}^* Q_{n}^* Q_{n}^*]
423 Double_t xQ2nQ1nQ1nQ2nstarQ1nstarQ1nstar = pow(xQvector2n.Mod(),2.)*pow(xQvector1n.Mod(),4.);//|Q_{2n}|^2 |Q_{n}|^4
424 Double_t reQ2nQ1nQ1nQ1nstarQ1nstarQ1nstarQ1nstar = pow(xQvector1n.Mod(),4.)*(pow(xQvector1n.X(),2.)*xQvector2n.X()-xQvector2n.X()*pow(xQvector1n.Y(),2.)+2.*xQvector1n.X()*xQvector1n.Y()*xQvector2n.Y());//Re[Q_{2n} Q_{n} Q_{n} Q_{n}^* Q_{n}^* Q_{n}^* Q_{n}^*]
425 //Double_t imQ2nQ1nQ1nQ1nstarQ1nstarQ1nstarQ1nstar = pow(xQvector1n.Mod(),4.)*(pow(xQvector1n.X(),2.)*xQvector2n.Y()-xQvector2n.Y()*pow(xQvector1n.Y(),2.)-2.*xQvector1n.X()*xQvector2n.X()*xQvector1n.Y());//Re[Q_{2n} Q_{n} Q_{n} Q_{n}^* Q_{n}^* Q_{n}^* Q_{n}^*]
426 //---------------------------------------------------------------------------------------------------------
428 //---------------------------------------------------------------------------------------------------------
430 // **************************************
431 // **** multi-particle correlations: ****
432 // **************************************
434 // Remark 1: multi-particle correlations calculated with Q-vectors are stored in fQCorrelations.
435 // Remark 2: binning of fQCorrelations is organized as follows:
437 // 1st bin: <2>_{n|n} = two1n1n
438 // 2nd bin: <2>_{2n|2n} = two2n2n
439 // 3rd bin: <2>_{3n|3n} = two3n3n
440 // 4th bin: <2>_{4n|4n} = two4n4n
441 // 5th bin: -- EMPTY --
442 // 6th bin: <3>_{2n|n,n} = three2n1n1n
443 // 7th bin: <3>_{3n|2n,n} = three3n2n1n
444 // 8th bin: <3>_{4n|2n,2n} = three4n2n2n
445 // 9th bin: <3>_{4n|3n,n} = three4n3n1n
446 //10th bin: -- EMPTY --
447 //11th bin: <4>_{n,n|n,n} = four1n1n1n1n
448 //12th bin: <4>_{2n,n|2n,n} = four2n1n2n1n
449 //13th bin: <4>_{2n,2n|2n,2n} = four2n2n2n2n
450 //14th bin: <4>_{3n|n,n,n} = four3n1n1n1n
451 //15th bin: <4>_{3n,n|3n,n} = four3n1n3n1n
452 //16th bin: <4>_{3n,n|2n,2n} = four3n1n2n2n
453 //17th bin: <4>_{4n|2n,n,n} = four4n2n1n1n
454 //18th bin: -- EMPTY --
455 //19th bin: <5>_{2n|n,n,n,n} = five2n1n1n1n1n
456 //20th bin: <5>_{2n,2n|2n,n,n} = five2n2n2n1n1n
457 //21st bin: <5>_{3n,n|2n,n,n} = five3n1n2n1n1n
458 //22nd bin: <5>_{4n|n,n,n,n} = five4n1n1n1n1n
459 //23rd bin: -- EMPTY --
460 //24th bin: <6>_{n,n,n|n,n,n} = six1n1n1n1n1n1n
461 //25th bin: <6>_{2n,n,n|2n,n,n} = six2n1n1n2n1n1n
462 //26th bin: <6>_{2n,2n|n,n,n,n} = six2n2n1n1n1n1n
463 //27th bin: <6>_{3n,n|n,n,n,n} = six3n1n1n1n1n1n
464 //28th bin: -- EMPTY --
465 //29th bin: <7>_{2n,n,n|n,n,n,n} = seven2n1n1n1n1n1n1n
466 //30th bin: -- EMPTY --
467 //31st bin: <8>_{n,n,n,n|n,n,n,n} = eight1n1n1n1n1n1n1n1n
470 // binning of fQProduct (all correlations are evaluated in harmonic n):
479 Double_t two1n1n=0., two2n2n=0., two3n3n=0., two4n4n=0.;
482 two1n1n = (pow(xQvector1n.Mod(),2.)-xMult)/(xMult*(xMult-1.)); //<2>_{n|n} = <cos(n*(phi1-phi2))>
483 two2n2n = (pow(xQvector2n.Mod(),2.)-xMult)/(xMult*(xMult-1.)); //<2>_{2n|2n} = <cos(2n*(phi1-phi2))>
484 two3n3n = (pow(xQvector3n.Mod(),2.)-xMult)/(xMult*(xMult-1.)); //<2>_{3n|3n} = <cos(3n*(phi1-phi2))>
485 two4n4n = (pow(xQvector4n.Mod(),2.)-xMult)/(xMult*(xMult-1.)); //<2>_{4n|4n} = <cos(4n*(phi1-phi2))>
487 fQCorrelations->Fill(0.,two1n1n,xMult*(xMult-1.));
488 fQCorrelations->Fill(1.,two2n2n,xMult*(xMult-1.));
489 fQCorrelations->Fill(2.,two3n3n,xMult*(xMult-1.));
490 fQCorrelations->Fill(3.,two4n4n,xMult*(xMult-1.));
492 f2pDistribution->Fill(two1n1n,xMult*(xMult-1.));
496 Double_t three2n1n1n=0., three3n2n1n=0., three4n2n2n=0., three4n3n1n=0.;
499 three2n1n1n = (reQ2nQ1nstarQ1nstar-2.*pow(xQvector1n.Mod(),2.)-pow(xQvector2n.Mod(),2.)+2.*xMult)/(xMult*(xMult-1.)*(xMult-2.)); //Re[<3>_{2n|n,n}] = Re[<3>_{n,n|2n}] = <cos(n*(2.*phi1-phi2-phi3))>
500 three3n2n1n = (reQ3nQ2nstarQ1nstar-pow(xQvector3n.Mod(),2.)-pow(xQvector2n.Mod(),2.)-pow(xQvector1n.Mod(),2.)+2.*xMult)/(xMult*(xMult-1.)*(xMult-2.)); //Re[<3>_{3n|2n,n}] = Re[<3>_{2n,n|3n}] = <cos(n*(3.*phi1-2.*phi2-phi3))>
501 three4n2n2n = (reQ4nQ2nstarQ2nstar-2.*pow(xQvector2n.Mod(),2.)-pow(xQvector4n.Mod(),2.)+2.*xMult)/(xMult*(xMult-1.)*(xMult-2.)); //Re[<3>_{4n|2n,2n}] = Re[<3>_{2n,2n|4n}] = <cos(n*(4.*phi1-2.*phi2-2.*phi3))>
502 three4n3n1n = (reQ4nQ3nstarQ1nstar-pow(xQvector4n.Mod(),2.)-pow(xQvector3n.Mod(),2.)-pow(xQvector1n.Mod(),2.)+2.*xMult)/(xMult*(xMult-1.)*(xMult-2.)); //Re[<3>_{4n|3n,n}] = Re[<3>_{3n,n|4n}] = <cos(n*(4.*phi1-3.*phi2-phi3))>
504 fQCorrelations->Fill(5.,three2n1n1n,xMult*(xMult-1.)*(xMult-2.));
505 fQCorrelations->Fill(6.,three3n2n1n,xMult*(xMult-1.)*(xMult-2.));
506 fQCorrelations->Fill(7.,three4n2n2n,xMult*(xMult-1.)*(xMult-2.));
507 fQCorrelations->Fill(8.,three4n3n1n,xMult*(xMult-1.)*(xMult-2.));
511 Double_t four1n1n1n1n=0., four2n2n2n2n=0., four2n1n2n1n=0., four3n1n1n1n=0., four4n2n1n1n=0., four3n1n2n2n=0., four3n1n3n1n=0.;
514 four1n1n1n1n = (2.*xMult*(xMult-3.)+pow(xQvector1n.Mod(),4.)-4.*(xMult-2.)*pow(xQvector1n.Mod(),2.)-2.*reQ2nQ1nstarQ1nstar+pow(xQvector2n.Mod(),2.))/(xMult*(xMult-1)*(xMult-2.)*(xMult-3.));//<4>_{n,n|n,n}
515 four2n2n2n2n = (2.*xMult*(xMult-3.)+pow(xQvector2n.Mod(),4.)-4.*(xMult-2.)*pow(xQvector2n.Mod(),2.)-2.*reQ4nQ2nstarQ2nstar+pow(xQvector4n.Mod(),2.))/(xMult*(xMult-1)*(xMult-2.)*(xMult-3.));//<4>_{2n,2n|2n,2n}
516 four2n1n2n1n = (xQ2nQ1nQ2nstarQ1nstar-2.*reQ3nQ2nstarQ1nstar-2.*reQ2nQ1nstarQ1nstar)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.))-((xMult-5.)*pow(xQvector1n.Mod(),2.)+(xMult-4.)*pow(xQvector2n.Mod(),2.)-pow(xQvector3n.Mod(),2.))/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.))+(xMult-6.)/((xMult-1.)*(xMult-2.)*(xMult-3.));//Re[<4>_{2n,n|2n,n}]
517 four3n1n1n1n = (reQ3nQ1nstarQ1nstarQ1nstar-3.*reQ3nQ2nstarQ1nstar-3.*reQ2nQ1nstarQ1nstar)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.))+(2.*pow(xQvector3n.Mod(),2.)+3.*pow(xQvector2n.Mod(),2.)+6.*pow(xQvector1n.Mod(),2.)-6.*xMult)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.));//Re[<4>_{3n|n,n,n}]
518 four4n2n1n1n = (reQ4nQ2nstarQ1nstarQ1nstar-2.*reQ4nQ3nstarQ1nstar-reQ4nQ2nstarQ2nstar-2.*reQ3nQ2nstarQ1nstar)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.))-(reQ2nQ1nstarQ1nstar-2.*pow(xQvector4n.Mod(),2.)-2.*pow(xQvector3n.Mod(),2.)-3.*pow(xQvector2n.Mod(),2.)-4.*pow(xQvector1n.Mod(),2.))/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.))-(6.)/((xMult-1.)*(xMult-2.)*(xMult-3.));//Re[<4>_{4n|2n,n,n}]
519 four3n1n2n2n = (reQ3nQ1nQ2nstarQ2nstar-reQ4nQ2nstarQ2nstar-reQ3nQ1nQ4nstar-2.*reQ3nQ2nstarQ1nstar)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.))-(2.*reQ1nQ1nQ2nstar-pow(xQvector4n.Mod(),2.)-2.*pow(xQvector3n.Mod(),2.)-4.*pow(xQvector2n.Mod(),2.)-4.*pow(xQvector1n.Mod(),2.))/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.))-(6.)/((xMult-1.)*(xMult-2.)*(xMult-3.));//Re[<4>_{3n,n|2n,2n}]
520 four3n1n3n1n = (pow(xQvector3n.Mod(),2.)*pow(xQvector1n.Mod(),2.)-2.*reQ4nQ3nstarQ1nstar-2.*reQ3nQ2nstarQ1nstar)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.))+(pow(xQvector4n.Mod(),2.)-(xMult-4.)*pow(xQvector3n.Mod(),2.)+pow(xQvector2n.Mod(),2.)-(xMult-4.)*pow(xQvector1n.Mod(),2.))/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.))+(xMult-6.)/((xMult-1.)*(xMult-2.)*(xMult-3.));//<4>_{3n,n|3n,n}
521 //four_3n1n3n1n = Q3nQ1nQ3nstarQ1nstar/(M*(M-1.)*(M-2.)*(M-3.))-(2.*three_3n2n1n+2.*three_4n3n1n)/(M-3.)-(two_4n4n+M*two_3n3n+two_2n2n+M*two_1n1n)/((M-2.)*(M-3.))-M/((M-1.)*(M-2.)*(M-3.));//<4>_{3n,n|3n,n}
523 fQCorrelations->Fill(10.,four1n1n1n1n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.));
524 fQCorrelations->Fill(11.,four2n1n2n1n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.));
525 fQCorrelations->Fill(12.,four2n2n2n2n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.));
526 fQCorrelations->Fill(13.,four3n1n1n1n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.));
527 fQCorrelations->Fill(14.,four3n1n3n1n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.));
528 fQCorrelations->Fill(15.,four3n1n2n2n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.));
529 fQCorrelations->Fill(16.,four4n2n1n1n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.));
531 f4pDistribution->Fill(four1n1n1n1n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.));
533 fQProduct->Fill(0.,two1n1n*four1n1n1n1n,xMult*(xMult-1.)*xMult*(xMult-1.)*(xMult-2.)*(xMult-3.));
537 Double_t five2n1n1n1n1n=0., five2n2n2n1n1n=0., five3n1n2n1n1n=0., five4n1n1n1n1n=0.;
540 five2n1n1n1n1n = (reQ2nQ1nQ1nstarQ1nstarQ1nstar-reQ3nQ1nstarQ1nstarQ1nstar+6.*reQ3nQ2nstarQ1nstar)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.))-(reQ2nQ1nQ3nstar+3.*(xMult-6.)*reQ2nQ1nstarQ1nstar+3.*reQ1nQ1nQ2nstar)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.))-(2.*pow(xQvector3n.Mod(),2.)+3.*pow(xQvector2n.Mod()*xQvector1n.Mod(),2.)-3.*(xMult-4.)*pow(xQvector2n.Mod(),2.))/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.))-3.*(pow(xQvector1n.Mod(),4.)-2.*(2*xMult-5.)*pow(xQvector1n.Mod(),2.)+2.*xMult*(xMult-4.))/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.));//Re[<5>_{2n,n|n,n,n}]
542 five2n2n2n1n1n = (reQ2nQ2nQ2nstarQ1nstarQ1nstar-reQ4nQ2nstarQ1nstarQ1nstar-2.*reQ2nQ2nQ3nstarQ1nstar)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.))+2.*(reQ4nQ2nstarQ2nstar+4.*reQ3nQ2nstarQ1nstar+reQ3nQ1nQ4nstar)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.))+(reQ2nQ2nQ4nstar-2.*(xMult-5.)*reQ2nQ1nstarQ1nstar+2.*reQ1nQ1nQ2nstar)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.))-(2.*pow(xQvector4n.Mod(),2.)+4.*pow(xQvector3n.Mod(),2.)+1.*pow(xQvector2n.Mod(),4.)-2.*(3.*xMult-10.)*pow(xQvector2n.Mod(),2.))/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.))-(4.*pow(xQvector1n.Mod(),2.)*pow(xQvector2n.Mod(),2.)-4.*(xMult-5.)*pow(xQvector1n.Mod(),2.)+4.*xMult*(xMult-6.))/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.));//Re[<5>_{2n,2n|2n,n,n}]
544 //five_2n2n2n1n1n = reQ2nQ2nQ2nstarQ1nstarQ1nstar/(M*(M-1.)*(M-2.)*(M-3.)*(M-4.))-(4.*four_2n1n2n1n+2.*four_3n1n2n2n+1.*four_2n2n2n2n+four_4n2n1n1n)/(M-4.)-(2.*three_4n3n1n+three_4n2n2n+three_4n2n2n+2.*three_3n2n1n)/((M-3.)*(M-4.))-(4.*three_3n2n1n+(2.*M-1.)*three_2n1n1n+2.*three_2n1n1n)/((M-3.)*(M-4.))-(two_4n4n+2.*two_3n3n+4.*(M-1.)*two_2n2n+2.*(2.*M-1.)*two_1n1n)/((M-2.)*(M-3.)*(M-4.))-(2.*M-1.)/((M-1.)*(M-2.)*(M-3.)*(M-4.)); //OK!
546 five4n1n1n1n1n = (reQ4nQ1nstarQ1nstarQ1nstarQ1nstar-6.*reQ4nQ2nstarQ1nstarQ1nstar-4.*reQ3nQ1nstarQ1nstarQ1nstar)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.))+(8.*reQ4nQ3nstarQ1nstar+3.*reQ4nQ2nstarQ2nstar+12.*reQ3nQ2nstarQ1nstar+12.*reQ2nQ1nstarQ1nstar)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.))-(6.*pow(xQvector4n.Mod(),2.)+8.*pow(xQvector3n.Mod(),2.)+12.*pow(xQvector2n.Mod(),2.)+24.*pow(xQvector1n.Mod(),2.)-24.*xMult)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.));//Re[<5>_{4n|n,n,n,n}]
548 //five_4n1n1n1n1n = reQ4nQ1nstarQ1nstarQ1nstarQ1nstar/(M*(M-1.)*(M-2.)*(M-3.)*(M-4.)) - (4.*four_3n1n1n1n+6.*four_4n2n1n1n)/(M-4.) - (6.*three_2n1n1n + 12.*three_3n2n1n + 4.*three_4n3n1n + 3.*three_4n2n2n)/((M-3.)*(M-4.)) - (4.*two_1n1n + 6.*two_2n2n + 4.*two_3n3n + 1.*two_4n4n)/((M-2.)*(M-3.)*(M-4.)) - 1./((M-1.)*(M-2.)*(M-3.)*(M-4.)); //OK!
550 five3n1n2n1n1n = (reQ3nQ1nQ2nstarQ1nstarQ1nstar-reQ4nQ2nstarQ1nstarQ1nstar-reQ3nQ1nstarQ1nstarQ1nstar)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.))-(reQ3nQ1nQ2nstarQ2nstar-3.*reQ4nQ3nstarQ1nstar-reQ4nQ2nstarQ2nstar)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.))-((2.*xMult-13.)*reQ3nQ2nstarQ1nstar-reQ3nQ1nQ4nstar-9.*reQ2nQ1nstarQ1nstar)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.))-(2.*reQ1nQ1nQ2nstar+2.*pow(xQvector4n.Mod(),2.)-2.*(xMult-5.)*pow(xQvector3n.Mod(),2.)+2.*pow(xQvector3n.Mod(),2.)*pow(xQvector1n.Mod(),2.))/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.))+(2.*(xMult-6.)*pow(xQvector2n.Mod(),2.)-2.*pow(xQvector2n.Mod(),2.)*pow(xQvector1n.Mod(),2.)-pow(xQvector1n.Mod(),4.)+2.*(3.*xMult-11.)*pow(xQvector1n.Mod(),2.))/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.))-4.*(xMult-6.)/((xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.));//Re[<5>_{3n,n|2n,n,n}]
552 //five3n1n2n1n1n = reQ3nQ1nQ2nstarQ1nstarQ1nstar/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)) - (four4n2n1n1n+four1n1n1n1n+four3n1n1n1n+2.*four2n1n2n1n+2.*three3n2n1n+2.*four3n1n3n1n+four3n1n2n2n)/(xMult-4.) - (2.*three4n3n1n+three4n2n2n+6.*three3n2n1n+three4n3n1n+2.*three3n2n1n+3.*three2n1n1n+2.*three2n1n1n+4.*two1n1n+2.*two2n2n+2.*two3n3n)/((xMult-3.)*(xMult-4.)) - (5.*two1n1n + 4.*two2n2n + 3.*two3n3n + 1.*two4n4n + 2.)/((xMult-2.)*(xMult-3.)*(xMult-4.)) - 1./((xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)) ;//Re[<5>_{3n,n|2n,n,n}] //OK!
554 //five3n1n2n1n1n = reQ3nQ1nQ2nstarQ1nstarQ1nstar/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)) - (four4n2n1n1n+four1n1n1n1n+four3n1n1n1n+2.*four2n1n2n1n+2.*four3n1n3n1n+four3n1n2n2n)/(xMult-4.) - (2.*three4n3n1n+three4n2n2n+2.*xMult*three3n2n1n+three4n3n1n+2.*three3n2n1n+3.*three2n1n1n+2.*three2n1n1n)/((xMult-3.)*(xMult-4.)) - ((4.*xMult-3.)*two1n1n + 2.*xMult*two2n2n + (2.*xMult-1.)*two3n3n + two4n4n)/((xMult-2.)*(xMult-3.)*(xMult-4.)) - (2.*xMult-1.)/((xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)) ;//Re[<5>_{3n,n|2n,n,n}] //OK!
556 fQCorrelations->Fill(18.,five2n1n1n1n1n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.));
557 fQCorrelations->Fill(19.,five2n2n2n1n1n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.));
558 fQCorrelations->Fill(20.,five3n1n2n1n1n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.));
559 fQCorrelations->Fill(21.,five4n1n1n1n1n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.));
563 Double_t six1n1n1n1n1n1n=0., six2n2n1n1n1n1n=0., six3n1n1n1n1n1n=0., six2n1n1n2n1n1n=0.;
566 six1n1n1n1n1n1n = (pow(xQvector1n.Mod(),6.)+9.*xQ2nQ1nQ2nstarQ1nstar-6.*reQ2nQ1nQ1nstarQ1nstarQ1nstar)/(xMult*(xMult-1)*(xMult-2)*(xMult-3)*(xMult-4)*(xMult-5))+4.*(reQ3nQ1nstarQ1nstarQ1nstar-3.*reQ3nQ2nstarQ1nstar)/(xMult*(xMult-1)*(xMult-2)*(xMult-3)*(xMult-4)*(xMult-5))+2.*(9.*(xMult-4.)*reQ2nQ1nstarQ1nstar+2.*pow(xQvector3n.Mod(),2.))/(xMult*(xMult-1)*(xMult-2)*(xMult-3)*(xMult-4)*(xMult-5))-9.*(pow(xQvector1n.Mod(),4.)+pow(xQvector2n.Mod(),2.))/(xMult*(xMult-1)*(xMult-2)*(xMult-3)*(xMult-5))+(18.*pow(xQvector1n.Mod(),2.))/(xMult*(xMult-1)*(xMult-3)*(xMult-4))-(6.)/((xMult-1)*(xMult-2)*(xMult-3));//<6>_{n,n,n|n,n,n}
568 six2n1n1n2n1n1n = (xQ2nQ1nQ1nQ2nstarQ1nstarQ1nstar-xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(2.*five2n2n2n1n1n+4.*five2n1n1n1n1n+4.*five3n1n2n1n1n+4.*four2n1n2n1n+1.*four1n1n1n1n)-xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(4.*four1n1n1n1n+4.*two1n1n+2.*three2n1n1n+2.*three2n1n1n+4.*four3n1n1n1n+8.*three2n1n1n+2.*four4n2n1n1n+4.*four2n1n2n1n+2.*two2n2n+8.*four2n1n2n1n+4.*four3n1n3n1n+8.*three3n2n1n+4.*four3n1n2n2n+4.*four1n1n1n1n+4.*four2n1n2n1n+1.*four2n2n2n2n)-xMult*(xMult-1.)*(xMult-2.)*(2.*three2n1n1n+8.*two1n1n+4.*two1n1n+2.+4.*two1n1n+4.*three2n1n1n+2.*two2n2n+4.*three2n1n1n+8.*three3n2n1n+8.*two2n2n+4.*three4n3n1n+4.*two3n3n+4.*three3n2n1n+4.*two1n1n+8.*three2n1n1n+4.*two1n1n+4.*three3n2n1n+4.*three2n1n1n+2.*two2n2n+4.*three3n2n1n+2.*three4n2n2n)-xMult*(xMult-1.)*(4.*two1n1n+4.+4.*two1n1n+2.*two2n2n+1.+4.*two1n1n+4.*two2n2n+4.*two3n3n+ 1.+2.*two2n2n+1.*two4n4n)-xMult)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(xMult-5.));//<6>_{2n,n,n|2n,n,n}
570 six2n2n1n1n1n1n = (reQ2nQ2nQ1nstarQ1nstarQ1nstarQ1nstar-xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(five4n1n1n1n1n+8.*five2n1n1n1n1n+6.*five2n2n2n1n1n)-xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(4.*four3n1n1n1n+6.*four4n2n1n1n+12.*three2n1n1n+12.*four1n1n1n1n+24.*four2n1n2n1n+4.*four3n1n2n2n+3.*four2n2n2n2n)-xMult*(xMult-1.)*(xMult-2.)*(6.*three2n1n1n+12.*three3n2n1n+4.*three4n3n1n+3.*three4n2n2n+8.*three2n1n1n+24.*two1n1n+12.*two2n2n+12.*three2n1n1n+8.*three3n2n1n+1.*three4n2n2n)-xMult*(xMult-1.)*(4.*two1n1n+6.*two2n2n+4.*two3n3n+1.*two4n4n+2.*two2n2n+8.*two1n1n+6.)-xMult)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(xMult-5.));//<6>_{2n,2n,n|n,n,n}
572 six3n1n1n1n1n1n = (reQ3nQ1nQ1nstarQ1nstarQ1nstarQ1nstar-xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(five4n1n1n1n1n+4.*five2n1n1n1n1n+6.*five3n1n2n1n1n+4.*four3n1n1n1n)-xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(4.*four3n1n1n1n+6.*four4n2n1n1n+6.*four1n1n1n1n+12.*three2n1n1n+12.*four2n1n2n1n+6.*four3n1n1n1n+12.*three3n2n1n+4.*four3n1n3n1n+3.*four3n1n2n2n)-xMult*(xMult-1.)*(xMult-2.)*(6.*three2n1n1n+12.*three3n2n1n+4.*three4n3n1n+3.*three4n2n2n+4.*two1n1n+12.*two1n1n+6.*three2n1n1n+12.*three2n1n1n+4.*three3n2n1n+12.*two2n2n+4.*three3n2n1n+4.*two3n3n+1.*three4n3n1n+6.*three3n2n1n)-xMult*(xMult-1.)*(4.*two1n1n+6.*two2n2n+4.*two3n3n+1.*two4n4n+1.*two1n1n+4.+6.*two1n1n+4.*two2n2n+1.*two3n3n)-xMult)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(xMult-5.));//<6>_{3n,n|n,n,n,n}
574 fQCorrelations->Fill(23.,six1n1n1n1n1n1n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(xMult-5.));
575 fQCorrelations->Fill(24.,six2n1n1n2n1n1n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(xMult-5.));
576 fQCorrelations->Fill(25.,six2n2n1n1n1n1n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(xMult-5.));
577 fQCorrelations->Fill(26.,six3n1n1n1n1n1n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(xMult-5.));
579 f6pDistribution->Fill(six1n1n1n1n1n1n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(xMult-5.));
581 fQProduct->Fill(1.,two1n1n*six1n1n1n1n1n1n,xMult*(xMult-1.)*xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(xMult-5.));
582 fQProduct->Fill(3.,four1n1n1n1n*six1n1n1n1n1n1n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(xMult-5.));
586 Double_t seven2n1n1n1n1n1n1n=0.;
589 seven2n1n1n1n1n1n1n = (reQ2nQ1nQ1nQ1nstarQ1nstarQ1nstarQ1nstar-xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(xMult-5.)*(2.*six3n1n1n1n1n1n+4.*six1n1n1n1n1n1n+1.*six2n2n1n1n1n1n+6.*six2n1n1n2n1n1n+8.*five2n1n1n1n1n)-xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(1.*five4n1n1n1n1n +8.*five2n1n1n1n1n+8.*four3n1n1n1n+12.*five3n1n2n1n1n+4.*five2n1n1n1n1n+3.*five2n2n2n1n1n+6.*five2n2n2n1n1n+6.*four1n1n1n1n+24.*four1n1n1n1n+12.*five2n1n1n1n1n+12.*five2n1n1n1n1n+12.*three2n1n1n+24.*four2n1n2n1n+4.*five3n1n2n1n1n+4.*five2n1n1n1n1n)-xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(4.*four3n1n1n1n+6.*four4n2n1n1n+12.*four1n1n1n1n+24.*three2n1n1n+24.*four2n1n2n1n+12.*four3n1n1n1n+24.*three3n2n1n+8.*four3n1n3n1n+6.*four3n1n2n2n+6.*three2n1n1n+12.*four1n1n1n1n+12.*four2n1n2n1n+6.*three2n1n1n+12.*four2n1n2n1n+4.*four3n1n2n2n+3.*four2n2n2n2n+4.*four1n1n1n1n+6.*three2n1n1n+24.*two1n1n+24.*four1n1n1n1n+4.*four3n1n1n1n+24.*two1n1n+24.*three2n1n1n+12.*two2n2n+24.*three2n1n1n+12.*four2n1n2n1n+8.*three3n2n1n+8.*four2n1n2n1n+1.*four4n2n1n1n)-xMult*(xMult-1.)*(xMult-2.)*(6.*three2n1n1n+1.*three2n1n1n+8.*two1n1n+12.*three3n2n1n+24.*two1n1n+12.*three2n1n1n+4.*three2n1n1n+8.*two1n1n+4.*three4n3n1n+24.*three2n1n1n+8.*three3n2n1n+12.*two1n1n+12.*two1n1n+3.*three4n2n2n+24.*two2n2n+6.*two2n2n+12.+12.*three3n2n1n+8.*two3n3n+12.*three2n1n1n+24.*two1n1n+4.*three3n2n1n+8.*three3n2n1n+2.*three4n3n1n+12.*two1n1n+8.*three2n1n1n+4.*three2n1n1n+2.*three3n2n1n+6.*two2n2n+8.*two2n2n+1.*three4n2n2n+4.*three3n2n1n+6.*three2n1n1n)-xMult*(xMult-1.)*(4.*two1n1n+2.*two1n1n+6.*two2n2n+8.+1.*two2n2n+4.*two3n3n+12.*two1n1n+4.*two1n1n+1.*two4n4n+8.*two2n2n+6.+2.*two3n3n+4.*two1n1n+1.*two2n2n)-xMult)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(xMult-5.)*(xMult-6.));
591 fQCorrelations->Fill(28.,seven2n1n1n1n1n1n1n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(xMult-5.)*(xMult-6.));
595 Double_t eight1n1n1n1n1n1n1n1n=0.;
598 eight1n1n1n1n1n1n1n1n = (pow(xQvector1n.Mod(),8)-xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(xMult-5.)*(xMult-6.)*(12.*seven2n1n1n1n1n1n1n+16.*six1n1n1n1n1n1n)-xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(xMult-5.)*(8.*six3n1n1n1n1n1n+48.*six1n1n1n1n1n1n+6.*six2n2n1n1n1n1n+96.*five2n1n1n1n1n+72.*four1n1n1n1n+36.*six2n1n1n2n1n1n)-xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(2.*five4n1n1n1n1n+32.*five2n1n1n1n1n+36.*four1n1n1n1n+32.*four3n1n1n1n+48.*five2n1n1n1n1n+48.*five3n1n2n1n1n+144.*five2n1n1n1n1n+288.*four1n1n1n1n+36.*five2n2n2n1n1n+144.*three2n1n1n+96.*two1n1n+144.*four2n1n2n1n)-xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(8.*four3n1n1n1n+48.*four1n1n1n1n+12.*four4n2n1n1n+96.*four2n1n2n1n+96.*three2n1n1n+72.*three2n1n1n+144.*two1n1n+16.*four3n1n3n1n+48.*four3n1n1n1n+144.*four1n1n1n1n+72.*four1n1n1n1n+96.*three3n2n1n+24.*four3n1n2n2n+144.*four2n1n2n1n+288.*two1n1n+288.*three2n1n1n+9.*four2n2n2n2n+72.*two2n2n+24.)-xMult*(xMult-1.)*(xMult-2.)*(12.*three2n1n1n+16.*two1n1n+24.*three3n2n1n+48.*three2n1n1n+96.*two1n1n+8.*three4n3n1n+32.*three3n2n1n+96.*three2n1n1n+144.*two1n1n+6.*three4n2n2n+96.*two2n2n+36.*two2n2n+72.+48.*three3n2n1n+16.*two3n3n+72.*three2n1n1n+144.*two1n1n)-xMult*(xMult-1.)*(8.*two1n1n+12.*two2n2n+16.+8.*two3n3n+48.*two1n1n+1.*two4n4n+16.*two2n2n+18.)-xMult)/(xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(xMult-5.)*(xMult-6.)*(xMult-7.));
600 fQCorrelations->Fill(30.,eight1n1n1n1n1n1n1n1n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(xMult-5.)*(xMult-6.)*(xMult-7.));
602 f8pDistribution->Fill(eight1n1n1n1n1n1n1n1n,xMult*(xMult-1.)*(xMult-2.)*(xMult-3.)*(xMult-4.)*(xMult-5.)*(xMult-6.)*(xMult-7.));
604 //---------------------------------------------------------------------------------------------------------
607 //---------------------------------------------------------------------------------------------------------
610 Double_t xQx = xQvector1n.X();
611 Double_t xQy = xQvector1n.Y();
612 Double_t xQ2x = xQvector2n.X();
613 Double_t xQ2y = xQvector2n.Y();
615 Double_t qx=0.,qy=0.,q2x=0.,q2y=0.,m=0.;//add comments for these variables (deleteMe)
617 for(Int_t i=0;i<xMult;i++) //check if nPrim == M
619 fTrack=anEvent->GetTrack(i);
622 fReq1n->Fill(fTrack->Pt(),cos(n*(fTrack->Phi())),1.);
623 fImq1n->Fill(fTrack->Pt(),sin(n*(fTrack->Phi())),1.);
624 fReq2n->Fill(fTrack->Pt(),cos(2.*n*(fTrack->Phi())),1.);
625 fImq2n->Fill(fTrack->Pt(),sin(2.*n*(fTrack->Phi())),1.);
629 Double_t twoDiff1n1n=0.,twoDiff2n2n=0.,threeDiff2n1n1n=0.,threeDiff1n1n2n=0.,fourDiff1n1n1n1n=0.;
631 for(Int_t bin=1;bin<(fnBinsPt+1);bin++)//loop over pt-bins
633 qx = (fReq1n->GetBinContent(bin))*(fReq1n->GetBinEntries(bin));
634 qy = (fImq1n->GetBinContent(bin))*(fImq1n->GetBinEntries(bin));
635 q2x = (fReq2n->GetBinContent(bin))*(fReq2n->GetBinEntries(bin));
636 q2y = (fImq2n->GetBinContent(bin))*(fImq2n->GetBinEntries(bin));
637 m = fReq1n->GetBinEntries(bin);
641 twoDiff1n1n = (qx*xQx+qy*xQy-m)/(m*(xMult-1.));
642 f2PerBin1n1n->Fill((bin-1)*0.1,twoDiff1n1n,m*(xMult-1.));//<2'>_{n|n}
644 twoDiff2n2n = (q2x*xQ2x+q2y*xQ2y-m)/(m*(xMult-1.));
645 f2PerBin2n2n->Fill((bin-1)*0.1,twoDiff2n2n,m*(xMult-1.));//<2'>_{2n|2n}
650 threeDiff2n1n1n = (q2x*(xQx*xQx-xQy*xQy)+2.*q2y*xQx*xQy-2.*(qx*xQx+qy*xQy)-(q2x*xQ2x+q2y*xQ2y)+2.*m)/(m*(xMult-1.)*(xMult-2.));
651 f3PerBin2n1n1n->Fill((bin-1)*0.1,threeDiff2n1n1n,m*(xMult-1.)*(xMult-2.));//Re[<3'>_{2n|n,n}]
653 threeDiff1n1n2n = (xQ2x*(qx*xQx-qy*xQy)+xQ2y*(qx*xQy+qy*xQx)-2.*(qx*xQx+qy*xQy)-(q2x*xQ2x+q2y*xQ2y)+2.*m)/(m*(xMult-1.)*(xMult-2.));
654 f3PerBin1n1n2n->Fill((bin-1)*0.1,threeDiff1n1n2n,m*(xMult-1.)*(xMult-2.));//Re[<3'>_{n,n|2n}]
659 fourDiff1n1n1n1n = ((xQx*xQx+xQy*xQy)*(qx*xQx+qy*xQy)-(q2x*(xQx*xQx-xQy*xQy)+2.*q2y*xQx*xQy)-(xQ2x*(qx*xQx-qy*xQy)+xQ2y*(qx*xQy+qy*xQx))+(q2x*xQ2x+q2y*xQ2y)-2.*(xMult-3.)*(qx*xQx+qy*xQy)-2.*m*(xQx*xQx+xQy*xQy)+2.*(xQx*qx+xQy*qy)+2.*m*(xMult-3.))/(m*(xMult-1.)*(xMult-2.)*(xMult-3.));
660 f4PerBin1n1n1n1n->Fill((bin-1)*0.1,fourDiff1n1n1n1n,m*(xMult-1.)*(xMult-2.)*(xMult-3.));//Re[<4'>_{n,n|n,n}]
669 //---------------------------------------------------------------------------------------------------------
718 //--------------------------------------------------------------------------------------------------------------------------------
720 // **********************
721 // **** NESTED LOOPS ****
722 // **********************
724 // Remark 1: multi-particle correlations calculated with nested loops are stored in fDirectCorrelations.
725 // Remark 2: binning of fDirectCorrelations: bins 0..40 - correlations needed for integrated flow; bins 40..80 - correlations needed for differential flow (taking as an example bin 0.5 < pt < 0.6)
727 // binning details of fDirectCorrelations (integrated flow):
729 // 1st bin: <2>_{n|n} = two1n1n
730 // 2nd bin: <2>_{2n|2n} = two2n2n
731 // 3rd bin: <2>_{3n|3n} = two3n3n
732 // 4th bin: <2>_{4n|4n} = two4n4n
733 // 5th bin: -- EMPTY --
734 // 6th bin: <3>_{2n|n,n} = three2n1n1n
735 // 7th bin: <3>_{3n|2n,n} = three3n2n1n
736 // 8th bin: <3>_{4n|2n,2n} = three4n2n2n
737 // 9th bin: <3>_{4n|3n,n} = three4n3n1n
738 //10th bin: -- EMPTY --
739 //11th bin: <4>_{n,n|n,n} = four1n1n1n1n
740 //12th bin: <4>_{2n,n|2n,n} = four2n1n2n1n
741 //13th bin: <4>_{2n,2n|2n,2n} = four2n2n2n2n
742 //14th bin: <4>_{3n|n,n,n} = four3n1n1n1n
743 //15th bin: <4>_{3n,n|3n,n} = four3n1n3n1n
744 //16th bin: <4>_{3n,n|2n,2n} = four3n1n2n2n
745 //17th bin: <4>_{4n|2n,n,n} = four4n2n1n1n
746 //18th bin: -- EMPTY --
747 //19th bin: <5>_{2n|n,n,n,n} = five2n1n1n1n1n
748 //20th bin: <5>_{2n,2n|2n,n,n} = five2n2n2n1n1n
749 //21st bin: <5>_{3n,n|2n,n,n} = five3n1n2n1n1n
750 //22nd bin: <5>_{4n|n,n,n,n} = five4n1n1n1n1n
751 //23rd bin: -- EMPTY --
752 //24th bin: <6>_{n,n,n|n,n,n} = six1n1n1n1n1n1n
753 //25th bin: <6>_{2n,n,n|2n,n,n} = six2n1n1n2n1n1n
754 //26th bin: <6>_{2n,2n|n,n,n,n} = six2n2n1n1n1n1n
755 //27th bin: <6>_{3n,n|n,n,n,n} = six3n1n1n1n1n1n
756 //28th bin: -- EMPTY --
757 //29th bin: <7>_{2n,n,n|n,n,n,n} = seven2n1n1n1n1n1n1n
758 //30th bin: -- EMPTY --
759 //31st bin: <8>_{n,n,n,n|n,n,n,n} = eight1n1n1n1n1n1n1n1n
761 Double_t phi1=0., phi2=0., phi3=0., phi4=0., phi5=0., phi6=0., phi7=0., phi8=0.;
763 //<2>_{k*n|k*n} (k=1,2,3 and 4)
764 for(Int_t i1=0;i1<xMult;i1++)
766 fTrack=anEvent->GetTrack(i1);
768 for(Int_t i2=0;i2<xMult;i2++)
771 fTrack=anEvent->GetTrack(i2);
773 fDirectCorrelations->Fill(0.,cos(n*(phi1-phi2)),1); //<2>_{n|n}
774 fDirectCorrelations->Fill(1.,cos(2.*n*(phi1-phi2)),1); //<2>_{2n|2n}
775 fDirectCorrelations->Fill(2.,cos(3.*n*(phi1-phi2)),1); //<2>_{3n|3n}
776 fDirectCorrelations->Fill(3.,cos(4.*n*(phi1-phi2)),1); //<2>_{4n|4n}
780 //<3>_{2n|n,n}, <3>_{3n|2n,n}, <3>_{4n|2n,2n} and <3>_{4n|3n,n}
781 for(Int_t i1=0;i1<xMult;i1++)
783 fTrack=anEvent->GetTrack(i1);
785 for(Int_t i2=0;i2<xMult;i2++)
788 fTrack=anEvent->GetTrack(i2);
790 for(Int_t i3=0;i3<xMult;i3++)
792 if(i3==i1||i3==i2)continue;
793 fTrack=anEvent->GetTrack(i3);
795 fDirectCorrelations->Fill(5.,cos(2*n*phi1-n*(phi2+phi3)),1); //<3>_{2n|n,n}
796 fDirectCorrelations->Fill(6.,cos(3.*n*phi1-2.*n*phi2-n*phi3),1); //<3>_{3n|2n,n}
797 fDirectCorrelations->Fill(7.,cos(4.*n*phi1-2.*n*phi2-2.*n*phi3),1); //<3>_{4n|2n,2n}
798 fDirectCorrelations->Fill(8.,cos(4.*n*phi1-3.*n*phi2-n*phi3),1); //<3>_{4n|3n,n}
803 //<4>_{n,n|n,n}, <4>_{2n,n|2n,n}, <4>_{2n,2n|2n,2n}, <4>_{3n|n,n,n}, <4>_{3n,n|3n,n}, <4>_{3n,n|2n,2n} and <4>_{4n|2n,n,n}
804 for(Int_t i1=0;i1<xMult;i1++)
806 fTrack=anEvent->GetTrack(i1);
808 for(Int_t i2=0;i2<xMult;i2++)
811 fTrack=anEvent->GetTrack(i2);
813 for(Int_t i3=0;i3<xMult;i3++)
815 if(i3==i1||i3==i2)continue;
816 fTrack=anEvent->GetTrack(i3);
818 for(Int_t i4=0;i4<xMult;i4++)
820 if(i4==i1||i4==i2||i4==i3)continue;
821 fTrack=anEvent->GetTrack(i4);
823 fDirectCorrelations->Fill(10.,cos(n*phi1+n*phi2-n*phi3-n*phi4),1); //<4>_{n,n|n,n}
824 fDirectCorrelations->Fill(11.,cos(2.*n*phi1+n*phi2-2.*n*phi3-n*phi4),1); //<4>_{2n,n|2n,n}
825 fDirectCorrelations->Fill(12.,cos(2.*n*phi1+2*n*phi2-2.*n*phi3-2.*n*phi4),1); //<4>_{2n,2n|2n,2n}
826 fDirectCorrelations->Fill(13.,cos(3.*n*phi1-n*phi2-n*phi3-n*phi4),1); //<4>_{3n|n,n,n}
827 fDirectCorrelations->Fill(14.,cos(3.*n*phi1+n*phi2-3.*n*phi3-n*phi4),1); //<4>_{3n,n|3n,n}
828 fDirectCorrelations->Fill(15.,cos(3.*n*phi1+n*phi2-2.*n*phi3-2.*n*phi4),1); //<4>_{3n,n|2n,2n}
829 fDirectCorrelations->Fill(16.,cos(4.*n*phi1-2.*n*phi2-n*phi3-n*phi4),1); //<4>_{4n|2n,n,n}
835 //<5>_{2n,n,n,n,n}, //<5>_{2n,2n|2n,n,n}, <5>_{3n,n|2n,n,n} and <5>_{4n|n,n,n,n}
836 for(Int_t i1=0;i1<xMult;i1++)
838 //cout<<"i1 = "<<i1<<endl;
839 fTrack=anEvent->GetTrack(i1);
841 for(Int_t i2=0;i2<xMult;i2++)
844 fTrack=anEvent->GetTrack(i2);
846 for(Int_t i3=0;i3<xMult;i3++)
848 if(i3==i1||i3==i2)continue;
849 fTrack=anEvent->GetTrack(i3);
851 for(Int_t i4=0;i4<xMult;i4++)
853 if(i4==i1||i4==i2||i4==i3)continue;
854 fTrack=anEvent->GetTrack(i4);
856 for(Int_t i5=0;i5<xMult;i5++)
858 if(i5==i1||i5==i2||i5==i3||i5==i4)continue;
859 fTrack=anEvent->GetTrack(i5);
861 fDirectCorrelations->Fill(18.,cos(2.*n*phi1+n*phi2-n*phi3-n*phi4-n*phi5),1); //<5>_{2n,n|n,n,n}
862 fDirectCorrelations->Fill(19.,cos(2.*n*phi1+2.*n*phi2-2.*n*phi3-n*phi4-n*phi5),1); //<5>_{2n,2n|2n,n,n}
863 fDirectCorrelations->Fill(20.,cos(3.*n*phi1+n*phi2-2.*n*phi3-n*phi4-n*phi5),1); //<5>_{3n,n|2n,n,n}
864 fDirectCorrelations->Fill(21.,cos(4.*n*phi1-n*phi2-n*phi3-n*phi4-n*phi5),1); //<5>_{4n|n,n,n,n}
871 //<6>_{n,n,n,n,n,n}, <6>_{2n,n,n|2n,n,n}, <6>_{2n,2n|n,n,n,n} and <6>_{3n,n|n,n,n,n}
872 for(Int_t i1=0;i1<xMult;i1++)
874 //cout<<"i1 = "<<i1<<endl;
875 fTrack=anEvent->GetTrack(i1);
877 for(Int_t i2=0;i2<xMult;i2++)
880 fTrack=anEvent->GetTrack(i2);
882 for(Int_t i3=0;i3<xMult;i3++)
884 if(i3==i1||i3==i2)continue;
885 fTrack=anEvent->GetTrack(i3);
887 for(Int_t i4=0;i4<xMult;i4++)
889 if(i4==i1||i4==i2||i4==i3)continue;
890 fTrack=anEvent->GetTrack(i4);
892 for(Int_t i5=0;i5<xMult;i5++)
894 if(i5==i1||i5==i2||i5==i3||i5==i4)continue;
895 fTrack=anEvent->GetTrack(i5);
897 for(Int_t i6=0;i6<xMult;i6++)
899 if(i6==i1||i6==i2||i6==i3||i6==i4||i6==i5)continue;
900 fTrack=anEvent->GetTrack(i6);
902 fDirectCorrelations->Fill(23.,cos(n*phi1+n*phi2+n*phi3-n*phi4-n*phi5-n*phi6),1); //<6>_{n,n,n|n,n,n}
903 fDirectCorrelations->Fill(24.,cos(2.*n*phi1+n*phi2+n*phi3-2.*n*phi4-n*phi5-n*phi6),1); //<6>_{2n,n,n|2n,n,n}
904 fDirectCorrelations->Fill(25.,cos(2.*n*phi1+2.*n*phi2-n*phi3-n*phi4-n*phi5-n*phi6),1); //<6>_{2n,2n|n,n,n,n}
905 fDirectCorrelations->Fill(26.,cos(3.*n*phi1+n*phi2-n*phi3-n*phi4-n*phi5-n*phi6),1); //<6>_{3n,n|n,n,n,n}
913 //<7>_{2n,n,n|n,n,n,n}
914 for(Int_t i1=0;i1<xMult;i1++)
916 //cout<<"i1 = "<<i1<<endl;
917 fTrack=anEvent->GetTrack(i1);
919 for(Int_t i2=0;i2<xMult;i2++)
922 fTrack=anEvent->GetTrack(i2);
924 for(Int_t i3=0;i3<xMult;i3++)
926 if(i3==i1||i3==i2)continue;
927 fTrack=anEvent->GetTrack(i3);
929 for(Int_t i4=0;i4<xMult;i4++)
931 if(i4==i1||i4==i2||i4==i3)continue;
932 fTrack=anEvent->GetTrack(i4);
934 for(Int_t i5=0;i5<xMult;i5++)
936 if(i5==i1||i5==i2||i5==i3||i5==i4)continue;
937 fTrack=anEvent->GetTrack(i5);
939 for(Int_t i6=0;i6<xMult;i6++)
941 if(i6==i1||i6==i2||i6==i3||i6==i4||i6==i5)continue;
942 fTrack=anEvent->GetTrack(i6);
944 for(Int_t i7=0;i7<xMult;i7++)
946 if(i7==i1||i7==i2||i7==i3||i7==i4||i7==i5||i7==i6)continue;
947 fTrack=anEvent->GetTrack(i7);
949 fDirectCorrelations->Fill(28.,cos(2.*n*phi1+n*phi2+n*phi3-n*phi4-n*phi5-n*phi6-n*phi7),1);//<7>_{2n,n,n|n,n,n,n}
958 //<8>_{n,n,n,n|n,n,n,n}
959 for(Int_t i1=0;i1<xMult;i1++)
961 cout<<"i1 = "<<i1<<endl;
962 fTrack=anEvent->GetTrack(i1);
964 for(Int_t i2=0;i2<xMult;i2++)
967 fTrack=anEvent->GetTrack(i2);
969 for(Int_t i3=0;i3<xMult;i3++)
971 if(i3==i1||i3==i2)continue;
972 fTrack=anEvent->GetTrack(i3);
974 for(Int_t i4=0;i4<xMult;i4++)
976 if(i4==i1||i4==i2||i4==i3)continue;
977 fTrack=anEvent->GetTrack(i4);
979 for(Int_t i5=0;i5<xMult;i5++)
981 if(i5==i1||i5==i2||i5==i3||i5==i4)continue;
982 fTrack=anEvent->GetTrack(i5);
984 for(Int_t i6=0;i6<xMult;i6++)
986 if(i6==i1||i6==i2||i6==i3||i6==i4||i6==i5)continue;
987 fTrack=anEvent->GetTrack(i6);
989 for(Int_t i7=0;i7<xMult;i7++)
991 if(i7==i1||i7==i2||i7==i3||i7==i4||i7==i5||i7==i6)continue;
992 fTrack=anEvent->GetTrack(i7);
994 for(Int_t i8=0;i8<xMult;i8++)
996 if(i8==i1||i8==i2||i8==i3||i8==i4||i8==i5||i8==i6||i8==i7)continue;
997 fTrack=anEvent->GetTrack(i8);
999 fDirectCorrelations->Fill(30.,cos(n*phi1+n*phi2+n*phi3+n*phi4-n*phi5-n*phi6-n*phi7-n*phi8),1);//<8>_{n,n,n,n|n,n,n,n}
1009 // binning details of fDirectCorrelations (differential flow):
1011 //41st bin: <2'>_{n|n}
1012 //42nd bin: <2'>_{2n|2n}
1013 //46th bin: <3'>_{2n|n,n}
1014 //47th bin: <3'>_{n,n|2n}
1015 //51st bin: <4'>_{n,n|n,n}
1018 for(Int_t i1=0;i1<xMult;i1++)
1020 fTrack=anEvent->GetTrack(i1);
1021 if(fTrack->Pt()>=0.5&&fTrack->Pt()<0.6)
1024 for(Int_t i2=0;i2<xMult;i2++)
1027 fTrack=anEvent->GetTrack(i2);
1029 fDirectCorrelations->Fill(40.,cos(1.*n*(phi1-phi2)),1); //<2'>_{n,n}
1030 fDirectCorrelations->Fill(41.,cos(2.*n*(phi1-phi2)),1); //<2'>_{2n,2n}
1036 for(Int_t i1=0;i1<xMult;i1++)
1038 fTrack=anEvent->GetTrack(i1);
1039 if(fTrack->Pt()>=0.5&&fTrack->Pt()<0.6)
1042 for(Int_t i2=0;i2<xMult;i2++)
1045 fTrack=anEvent->GetTrack(i2);
1047 for(Int_t i3=0;i3<xMult;i3++)
1049 if(i3==i1||i3==i2)continue;
1050 fTrack=anEvent->GetTrack(i3);
1052 fDirectCorrelations->Fill(45.,cos(n*(2.*phi1-phi2-phi3)),1); //<3'>_{2n|n,n}
1053 fDirectCorrelations->Fill(46.,cos(n*(phi1+phi2-2.*phi3)),1); //<3'>_{n,n|2n}
1060 for(Int_t i1=0;i1<xMult;i1++)
1062 fTrack=anEvent->GetTrack(i1);
1063 if(fTrack->Pt()>=0.5&&fTrack->Pt()<0.6)
1066 for(Int_t i2=0;i2<xMult;i2++)
1069 fTrack=anEvent->GetTrack(i2);
1071 for(Int_t i3=0;i3<xMult;i3++)
1073 if(i3==i1||i3==i2)continue;
1074 fTrack=anEvent->GetTrack(i3);
1076 for(Int_t i4=0;i4<xMult;i4++)
1078 if(i4==i1||i4==i2||i4==i3)continue;
1079 fTrack=anEvent->GetTrack(i4);
1081 fDirectCorrelations->Fill(50.,cos(n*(phi1+phi2-phi3-phi4)),1); //<4'>_{n,n|n,n}
1087 //--------------------------------------------------------------------------------------------------------------------------------
1097 //}//line needed only for nested loops - end of if(nPrim>8&&nPrim<14)
1101 //================================================================================================================
1103 void AliFlowAnalysisWithQCumulants::Finish()
1105 //calculate the final results
1106 AliQCumulantsFunctions finalResults(fIntFlowResultsQC,fDiffFlowResults2ndOrderQC,fDiffFlowResults4thOrderQC,fCovariances,fAvMultIntFlowQC,fQvectorComponents,fQCorrelations, fQProduct,fDirectCorrelations, f2PerBin1n1n,f2PerBin2n2n,f3PerBin2n1n1n,f3PerBin1n1n2n,f4PerBin1n1n1n1n,fCommonHistsResults2nd, fCommonHistsResults4th,fCommonHistsResults6th,fCommonHistsResults8th);
1108 finalResults.Calculate();
1111 //================================================================================================================
1113 void AliFlowAnalysisWithQCumulants::WriteHistograms(TString* outputFileName)
1115 //store the final results in output .root file
1116 TFile *output = new TFile(outputFileName->Data(),"RECREATE");
1117 output->mkdir("cobjQC","cobjQC");
1118 output->cd("cobjQC");
1123 //================================================================================================================