]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/FLOW/AliFlowCommon/AliFlowAnalysisWithMixedHarmonics.cxx
changed behavior cut on integer values
[u/mrichter/AliRoot.git] / PWG2 / FLOW / AliFlowCommon / AliFlowAnalysisWithMixedHarmonics.cxx
CommitLineData
489d5531 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. *
d66d46f7 14**************************************************************************/
15
489d5531 16/* $Id$ */
17
d66d46f7 18/**********************************************************
19 * In this class azimuthal correlators in mixed harmonics *
20 * are implemented in terms of Q-vectors. This approach *
21 * doesn't require evaluation of nested loops. This class *
22 * can be used to: *
23 * *
24 * a) Extract subdominant harmonics (like v1 and v4); *
5c078806 25 * b) Study flow of two-particle resonances; *
26 * c) Study strong parity violation. *
d66d46f7 27 * *
28 * Author: Ante Bilandzic (abilandzic@gmail.com) *
489d5531 29 *********************************************************/
30
31#define AliFlowAnalysisWithMixedHarmonics_cxx
32
33#include "Riostream.h"
34#include "AliFlowCommonConstants.h"
35#include "AliFlowCommonHist.h"
36#include "AliFlowCommonHistResults.h"
37
38#include "TMath.h"
39#include "TFile.h"
40#include "TList.h"
d66d46f7 41#include "TProfile.h"
489d5531 42#include "TProfile2D.h"
43
44#include "AliFlowEventSimple.h"
45#include "AliFlowTrackSimple.h"
46#include "AliFlowAnalysisWithMixedHarmonics.h"
47
48class TH1;
49class TList;
50ClassImp(AliFlowAnalysisWithMixedHarmonics)
51
d66d46f7 52//================================================================================================================
53
489d5531 54AliFlowAnalysisWithMixedHarmonics::AliFlowAnalysisWithMixedHarmonics():
d66d46f7 55fHistList(NULL),
489d5531 56fHistListName(NULL),
e999459d 57fHarmonic(1),
d66d46f7 58fAnalysisLabel(NULL),
59fAnalysisSettings(NULL),
d66d46f7 60fNoOfMultipicityBins(10),
61fMultipicityBinWidth(2),
62fMinMultiplicity(1),
e999459d 63fOppositeChargesPOI(kFALSE),
64fEvaluateDifferential3pCorrelator(kFALSE),
489d5531 65fCorrectForDetectorEffects(kTRUE),
5c078806 66fPrintOnTheScreen(kTRUE),
d66d46f7 67fCommonHists(NULL),
489d5531 68fnBinsPhi(0),
69fPhiMin(0),
70fPhiMax(0),
71fPhiBinWidth(0),
72fnBinsPt(0),
73fPtMin(0),
74fPtMax(0),
75fPtBinWidth(0),
76fnBinsEta(0),
77fEtaMin(0),
78fEtaMax(0),
d66d46f7 79fEtaBinWidth(0),
489d5531 80fWeightsList(NULL),
81fUsePhiWeights(kFALSE),
82fUsePtWeights(kFALSE),
83fUseEtaWeights(kFALSE),
84fUseParticleWeights(NULL),
85fPhiWeights(NULL),
86fPtWeights(NULL),
d66d46f7 87fEtaWeights(NULL),
489d5531 88fReQnk(NULL),
89fImQnk(NULL),
d66d46f7 90fSpk(NULL),
ecac11c2 91fProfileList(NULL),
d66d46f7 92f3pCorrelatorPro(NULL),
93fNonIsotropicTermsPro(NULL),
94f3pCorrelatorVsMPro(NULL),
95fNonIsotropicTermsVsMPro(NULL),
96fResultsList(NULL),
97f3pCorrelatorHist(NULL),
98fDetectorBiasHist(NULL),
e558a733 99f3pCorrelatorVsMHist(NULL),
5c078806 100fDetectorBiasVsMHist(NULL)
489d5531 101{
102 // Constructor.
5c078806 103
104 // Base list to hold all output objects:
d66d46f7 105 fHistList = new TList();
489d5531 106 fHistListName = new TString("cobjMH");
107 fHistList->SetName(fHistListName->Data());
d66d46f7 108 fHistList->SetOwner(kTRUE);
109
489d5531 110 // List to hold histograms with phi, pt and eta weights:
d66d46f7 111 fWeightsList = new TList();
112
ecac11c2 113 // List to hold all all-event profiles:
114 fProfileList = new TList();
115
489d5531 116 // List to hold objects with final results:
d66d46f7 117 fResultsList = new TList();
5c078806 118
119 // Initialize all arrays:
120 this->InitializeArrays();
489d5531 121
122} // AliFlowAnalysisWithMixedHarmonics::AliFlowAnalysisWithMixedHarmonics()
123
124//================================================================================================================
125
126AliFlowAnalysisWithMixedHarmonics::~AliFlowAnalysisWithMixedHarmonics()
127{
128 // Destructor.
129
130 delete fHistList;
131
132} // end of AliFlowAnalysisWithMixedHarmonics::~AliFlowAnalysisWithMixedHarmonics()
133
134//================================================================================================================
135
136void AliFlowAnalysisWithMixedHarmonics::Init()
d66d46f7 137{
138 // Initialize and book all objects.
139
140 // a) Cross check if the user settings make sense before starting;
141 // b) Access all common constants;
489d5531 142 // c) Book and nest all lists in the base list fHistList;
d66d46f7 143 // d) Book common control histograms;
144 // e) Book all event-by-event quantities;
145 // f) Book all all-event quantities;
5c078806 146 // g) Book and fill histograms to hold phi, pt and eta weights;
e999459d 147 // h) Store harmonic n used in cos[n*(phi1+phi2-2phi3)] and cos[n*(psi1+psi2-2phi3)].
5c078806 148
489d5531 149 //save old value and prevent histograms from being added to directory
150 //to avoid name clashes in case multiple analaysis objects are used
151 //in an analysis
152 Bool_t oldHistAddStatus = TH1::AddDirectoryStatus();
153 TH1::AddDirectory(kFALSE);
154
5c078806 155 TH1::SetDefaultSumw2();
156
489d5531 157 this->CrossCheckSettings();
d66d46f7 158 this->AccessConstants();
159 this->BookAndNestAllLists();
489d5531 160 this->BookProfileHoldingSettings();
d66d46f7 161 this->BookCommonHistograms();
162 this->BookAllEventByEventQuantities();
163 this->BookAllAllEventQuantities();
489d5531 164 this->BookAndFillWeightsHistograms();
e999459d 165 this->StoreHarmonic();
489d5531 166
167 TH1::AddDirectory(oldHistAddStatus);
5c078806 168
489d5531 169} // end of void AliFlowAnalysisWithMixedHarmonics::Init()
170
171//================================================================================================================
172
173void AliFlowAnalysisWithMixedHarmonics::Make(AliFlowEventSimple* anEvent)
174{
175 // Running over data only in this method.
d66d46f7 176
5c078806 177 // a) Check all pointers used in this method;
178 // b) Define local variables;
179 // c) Fill common control histograms;
180 // d) Loop over data and calculate e-b-e quantities Q_{n,k} and S_{p,k};
e999459d 181 // e) Calculate 3-p azimuthal correlator cos[n(phi1+phi2-2*phi3)] and non-isotropic terms in terms of Q_{n,k} and S_{p,k};
182 // f) Calculate differential 3-p azimuthal correlator cos[n(psi1+psi2-2*phi3)] in terms of Q_{2n} and p_{n}:
5c078806 183 // g) Reset all event-by-event quantities.
d66d46f7 184
5c078806 185 // a) Check all pointers used in this method:
d66d46f7 186 this->CheckPointersUsedInMake();
187
5c078806 188 // b) Define local variables:
489d5531 189 Double_t dPhi = 0.; // azimuthal angle in the laboratory frame
190 Double_t dPt = 0.; // transverse momentum
191 Double_t dEta = 0.; // pseudorapidity
192 Double_t wPhi = 1.; // phi weight
193 Double_t wPt = 1.; // pt weight
d66d46f7 194 Double_t wEta = 1.; // eta weight
195 AliFlowTrackSimple *aftsTrack = NULL; // simple track
196
5c078806 197 // c) Fill common control histograms:
d66d46f7 198 fCommonHists->FillControlHistograms(anEvent);
199
5c078806 200 // d) Loop over data and calculate e-b-e quantities:
e999459d 201 Int_t nPrim = anEvent->NumberOfTracks(); // nPrim = total number of primary tracks, i.e. nPrim = nRP + nPOI
489d5531 202 // nRP = # of particles used to determine the reaction plane ("Reference Particles");
203 // nPOI = # of particles of interest for a detailed flow analysis ("Particles of Interest");
d66d46f7 204 // Start loop over data:
489d5531 205 for(Int_t i=0;i<nPrim;i++)
206 {
207 aftsTrack=anEvent->GetTrack(i);
208 if(aftsTrack)
209 {
210 if(!(aftsTrack->InRPSelection() || aftsTrack->InPOISelection())) continue; // consider only tracks which are either RPs or POIs
e999459d 211 Int_t n = fHarmonic;
489d5531 212 if(aftsTrack->InRPSelection()) // checking RP condition:
213 {
214 dPhi = aftsTrack->Phi();
215 dPt = aftsTrack->Pt();
216 dEta = aftsTrack->Eta();
217 if(fUsePhiWeights && fPhiWeights && fnBinsPhi) // determine phi-weight for this particle:
218 {
219 wPhi = fPhiWeights->GetBinContent(1+(Int_t)(TMath::Floor(dPhi*fnBinsPhi/TMath::TwoPi())));
220 }
221 if(fUsePtWeights && fPtWeights && fnBinsPt) // determine pt-weight for this particle:
222 {
223 wPt = fPtWeights->GetBinContent(1+(Int_t)(TMath::Floor((dPt-fPtMin)/fPtBinWidth)));
224 }
225 if(fUseEtaWeights && fEtaWeights && fEtaBinWidth) // determine eta-weight for this particle:
226 {
227 wEta = fEtaWeights->GetBinContent(1+(Int_t)(TMath::Floor((dEta-fEtaMin)/fEtaBinWidth)));
228 }
e999459d 229 // Calculate Re[Q_{m,k}] and Im[Q_{m,k}], (m <= n and k = 0,1,2,3) for this event:
230 for(Int_t m=0;m<2;m++) // needs only Q-vectors in harmonincs n and 2n
489d5531 231 {
e999459d 232 for(Int_t k=0;k<4;k++) // to be improved (what is the maximum k that I need?)
489d5531 233 {
234 (*fReQnk)(m,k)+=pow(wPhi*wPt*wEta,k)*TMath::Cos((m+1)*n*dPhi);
235 (*fImQnk)(m,k)+=pow(wPhi*wPt*wEta,k)*TMath::Sin((m+1)*n*dPhi);
236 }
237 }
238 // Calculate partially S_{p,k} for this event (final calculation of S_{p,k} follows after the loop over data bellow):
239 for(Int_t p=0;p<4;p++) // to be improved (what is maximum p that I need?)
240 {
241 for(Int_t k=0;k<4;k++) // to be improved (what is maximum k that I need?)
242 {
243 (*fSpk)(p,k)+=pow(wPhi*wPt*wEta,k);
244 }
245 }
5c078806 246 } // end of if(aftsTrack->InRPSelection())
247 // POIs:
e999459d 248 if(fEvaluateDifferential3pCorrelator)
5c078806 249 {
e999459d 250 if(aftsTrack->InPOISelection()) // 1st POI
5c078806 251 {
e999459d 252 Double_t dPsi1 = aftsTrack->Phi();
253 Double_t dPt1 = aftsTrack->Pt();
254 Int_t iCharge1 = aftsTrack->Charge();
255 Bool_t b1stPOIisAlsoRP = kFALSE;
256 if(aftsTrack->InRPSelection()){b1stPOIisAlsoRP = kTRUE;}
257 for(Int_t j=0;j<nPrim;j++)
5c078806 258 {
e999459d 259 if(j==i){continue;}
260 aftsTrack=anEvent->GetTrack(j);
261 if(aftsTrack->InPOISelection()) // 2nd POI
262 {
263 Double_t dPsi2 = aftsTrack->Phi();
264 Double_t dPt2 = aftsTrack->Pt();
265 Int_t iCharge2 = aftsTrack->Charge();
266 if(fOppositeChargesPOI && iCharge1 == iCharge2){continue;}
267 Bool_t b2ndPOIisAlsoRP = kFALSE;
268 if(aftsTrack->InRPSelection()){b2ndPOIisAlsoRP = kTRUE;}
269 // Fill:
270 fRePEBE[0]->Fill((dPt1+dPt2)/2.,TMath::Cos(n*(dPsi1+dPsi2)),1.);
271 fImPEBE[0]->Fill((dPt1+dPt2)/2.,TMath::Sin(n*(dPsi1+dPsi2)),1.);
272 fRePEBE[1]->Fill(TMath::Abs(dPt1-dPt2),TMath::Cos(n*(dPsi1+dPsi2)),1.);
273 fImPEBE[1]->Fill(TMath::Abs(dPt1-dPt2),TMath::Sin(n*(dPsi1+dPsi2)),1.);
274 if(b1stPOIisAlsoRP)
275 {
276 fOverlapEBE[0][0]->Fill((dPt1+dPt2)/2.,TMath::Cos(n*(dPsi1-dPsi2)),1.);
277 fOverlapEBE[0][1]->Fill(TMath::Abs(dPt1-dPt2),TMath::Cos(n*(dPsi1-dPsi2)),1.);
278 }
279 if(b2ndPOIisAlsoRP)
280 {
281 fOverlapEBE[1][0]->Fill((dPt1+dPt2)/2.,TMath::Cos(n*(dPsi1-dPsi2)),1.);
282 fOverlapEBE[1][1]->Fill(TMath::Abs(dPt1-dPt2),TMath::Cos(n*(dPsi1-dPsi2)),1.);
283 }
284 } // end of if(aftsTrack->InPOISelection()) // 2nd POI
285 } // end of for(Int_t j=i+1;j<nPrim;j++)
286 } // end of if(aftsTrack->InPOISelection()) // 1st POI
287 } // end of if(fEvaluateDifferential3pCorrelator)
489d5531 288 } else // to if(aftsTrack)
289 {
290 cout<<endl;
e999459d 291 cout<<" WARNING (MH): No particle! (i.e. aftsTrack is a NULL pointer in Make().)"<<endl;
489d5531 292 cout<<endl;
293 }
294 } // end of for(Int_t i=0;i<nPrim;i++)
295
296 // Calculate the final expressions for S_{p,k}:
297 for(Int_t p=0;p<4;p++) // to be improved (what is maximum p that I need?)
298 {
e999459d 299 for(Int_t k=0;k<4;k++) // to be improved (what is maximum k that I need?)
489d5531 300 {
301 (*fSpk)(p,k)=pow((*fSpk)(p,k),p+1);
302 }
d66d46f7 303 }
304
e999459d 305 // e) Calculate 3-p correlator cos[n(phi1+phi2-2*phi3)] in terms of Q_{n,k} and S_{p,k}:
d66d46f7 306 if(anEvent->GetEventNSelTracksRP() >= 3)
307 {
308 this->Calculate3pCorrelator();
309 }
310 if(anEvent->GetEventNSelTracksRP() >= 0) // to be improved (is this correct if statement?)
311 {
312 this->CalculateNonIsotropicTerms();
313 }
314
e999459d 315 // f) Calculate differential 3-p azimuthal correlator cos[n(psi1+psi2-2*phi3)] in terms of Q_{2n} and p_{n}:
316 if(fEvaluateDifferential3pCorrelator && anEvent->GetEventNSelTracksRP() >= 1)
317 {
318 this->CalculateDifferential3pCorrelator(); // to be improved - add relevant if statements for the min # POIs as well
319 }
5c078806 320
321 // g) Reset all event-by-event quantities:
d66d46f7 322 this->ResetEventByEventQuantities();
5c078806 323
489d5531 324} // end of AliFlowAnalysisWithMixedHarmonics::Make(AliFlowEventSimple* anEvent)
5c078806 325
489d5531 326//================================================================================================================
327
328void AliFlowAnalysisWithMixedHarmonics::Finish()
329{
d66d46f7 330 // Calculate the final results.
331
332 // a) Check all pointers used in this method;
333 // b) Access settings for analysis with mixed harmonics;
5c078806 334 // c) Correct for detector effects;
335 // d) Print on the screen the final results.
d66d46f7 336
337 this->CheckPointersUsedInFinish();
338 this->AccessSettings();
e558a733 339 if(fCorrectForDetectorEffects)
340 {
341 this->CorrectForDetectorEffects();
342 this->CorrectForDetectorEffectsVsM();
343 }
e999459d 344 if(fPrintOnTheScreen) {this->PrintOnTheScreen();}
489d5531 345
346} // end of AliFlowAnalysisWithMixedHarmonics::Finish()
347
348//================================================================================================================
349
350void AliFlowAnalysisWithMixedHarmonics::GetOutputHistograms(TList *outputListHistos)
351{
d66d46f7 352 // Get pointers to all objects saved in the output file.
353
ecac11c2 354 // a) Get pointers for common control histograms.
489d5531 355 if(outputListHistos)
356 {
d66d46f7 357 this->SetHistList(outputListHistos);
358 if(!fHistList)
359 {
360 cout<<endl;
e999459d 361 cout<<" WARNING (MH): fHistList is NULL in GetOutputHistograms() !!!!"<<endl;
d66d46f7 362 cout<<endl;
363 exit(0);
364 }
ecac11c2 365 this->GetPointersForBaseHistograms();
d66d46f7 366 this->GetPointersForCommonHistograms();
ecac11c2 367 this->GetPointersForAllEventProfiles();
489d5531 368 this->GetPointersForResultsHistograms();
d66d46f7 369 } else
370 {
371 cout<<endl;
e999459d 372 cout<<" WARNING (MH): outputListHistos is NULL in GetOutputHistograms() !!!!"<<endl;
d66d46f7 373 cout<<endl;
374 exit(0);
375 }
489d5531 376
377} // end of void AliFlowAnalysisWithMixedHarmonics::GetOutputHistograms(TList *outputListHistos)
378
d66d46f7 379//================================================================================================================
380
ecac11c2 381void AliFlowAnalysisWithMixedHarmonics::GetPointersForBaseHistograms()
382{
383 // Get pointers to base histograms.
384
385 TString analysisSettingsName = "fAnalysisSettings";
386 TProfile *analysisSettings = dynamic_cast<TProfile*>(fHistList->FindObject(analysisSettingsName.Data()));
387 if(analysisSettings)
388 {
389 this->SetAnalysisSettings(analysisSettings);
390 } else
391 {
392 cout<<endl;
e999459d 393 cout<<" WARNING (MH): analysisSettings is NULL in GetPointersForBaseHistograms() !!!!"<<endl;
ecac11c2 394 cout<<endl;
395 exit(0);
396 }
397
398} // end of void AliFlowAnalysisWithMixedHarmonics::GetPointersForBaseHistograms()
399
400//================================================================================================================
401
489d5531 402void AliFlowAnalysisWithMixedHarmonics::GetPointersForCommonHistograms()
403{
404 // Get pointers to common control histograms.
405
406 TString commonHistsName = "AliFlowCommonHistMH";
407 AliFlowCommonHist *commonHist = dynamic_cast<AliFlowCommonHist*>(fHistList->FindObject(commonHistsName.Data()));
d66d46f7 408 if(commonHist)
409 {
410 this->SetCommonHists(commonHist);
411 } else
412 {
413 cout<<endl;
e999459d 414 cout<<" WARNING (MH): commonHist is NULL in GetPointersForCommonHistograms() !!!!"<<endl;
d66d46f7 415 cout<<endl;
416 exit(0);
417 }
418
419} // end of void AliFlowAnalysisWithMixedHarmonics::GetPointersForCommonHistograms()
420
421//================================================================================================================
422
ecac11c2 423void AliFlowAnalysisWithMixedHarmonics::GetPointersForAllEventProfiles()
424{
5c078806 425 // Get pointers to profiles holding final results.
ecac11c2 426
427 TList *profileList = NULL;
428 profileList = dynamic_cast<TList*>(fHistList->FindObject("Profiles"));
429 if(!profileList)
430 {
e999459d 431 cout<<endl;
432 cout<<" WARNING (MH): profileList is NULL in GetPointersForAllEventProfiles() !!!!"<<endl;
433 cout<<endl;
ecac11c2 434 exit(0);
435 }
436
437 TString s3pCorrelatorProName = "f3pCorrelatorPro";
5c078806 438 TProfile *p3pCorrelatorPro = dynamic_cast<TProfile*>(profileList->FindObject(s3pCorrelatorProName.Data()));
439 if(p3pCorrelatorPro)
ecac11c2 440 {
5c078806 441 this->Set3pCorrelatorPro(p3pCorrelatorPro);
ecac11c2 442 }
443 TString s3pCorrelatorVsMProName = "f3pCorrelatorVsMPro";
5c078806 444 TProfile *p3pCorrelatorVsMPro = dynamic_cast<TProfile*>(profileList->FindObject(s3pCorrelatorVsMProName.Data()));
445 if(p3pCorrelatorVsMPro)
ecac11c2 446 {
5c078806 447 this->Set3pCorrelatorVsMPro(p3pCorrelatorVsMPro);
ecac11c2 448 }
449 TString nonIsotropicTermsProName = "fNonIsotropicTermsPro";
450 TProfile *nonIsotropicTermsPro = dynamic_cast<TProfile*>(profileList->FindObject(nonIsotropicTermsProName.Data()));
451 if(nonIsotropicTermsPro)
452 {
453 this->SetNonIsotropicTermsPro(nonIsotropicTermsPro);
454 }
455 TString nonIsotropicTermsVsMProName = "fNonIsotropicTermsVsMPro";
456 TProfile2D *nonIsotropicTermsVsMPro = dynamic_cast<TProfile2D*>(profileList->FindObject(nonIsotropicTermsVsMProName.Data()));
457 if(nonIsotropicTermsVsMPro)
458 {
459 this->SetNonIsotropicTermsVsMPro(nonIsotropicTermsVsMPro);
5c078806 460 }
461 TString psdFlag[2] = {"PtSum","PtDiff"};
462 for(Int_t sd=0;sd<2;sd++)
463 {
464 TProfile *p3pCorrelatorVsPtSumDiffPro = dynamic_cast<TProfile*>(profileList->FindObject(Form("f3pCorrelatorVs%sPro",psdFlag[sd].Data())));
465 if(p3pCorrelatorVsPtSumDiffPro)
466 {
467 this->Set3pCorrelatorVsPtSumDiffPro(p3pCorrelatorVsPtSumDiffPro,sd);
468 }
469 }
470
ecac11c2 471} // end of void AliFlowAnalysisWithMixedHarmonics::GetPointersForAllEventProfiles()
472
473//================================================================================================================
474
489d5531 475void AliFlowAnalysisWithMixedHarmonics::GetPointersForResultsHistograms()
476{
477 // Get pointers to histograms holding final results.
478
479 TList *resultsList = NULL;
480 resultsList = dynamic_cast<TList*>(fHistList->FindObject("Results"));
481 if(!resultsList)
482 {
e999459d 483 cout<<endl;
484 cout<<" WARNING (MH): resultsList is NULL in GetPointersForResultsHistograms() !!!!"<<endl;
485 cout<<endl;
489d5531 486 exit(0);
d66d46f7 487 }
489d5531 488 TString s3pCorrelatorHistName = "f3pCorrelatorHist";
489 TH1D *h3pCorrelatorHist = dynamic_cast<TH1D*>(resultsList->FindObject(s3pCorrelatorHistName.Data()));
490 if(h3pCorrelatorHist)
491 {
492 this->Set3pCorrelatorHist(h3pCorrelatorHist);
d66d46f7 493 }
e558a733 494 TString s3pCorrelatorVsMHistName = "f3pCorrelatorVsMHist";
495 TH1D *h3pCorrelatorVsMHist = dynamic_cast<TH1D*>(resultsList->FindObject(s3pCorrelatorVsMHistName.Data()));
496 if(h3pCorrelatorVsMHist)
497 {
498 this->Set3pCorrelatorVsMHist(h3pCorrelatorVsMHist);
499 }
489d5531 500 TString detectorBiasHistName = "fDetectorBiasHist";
501 TH1D *detectorBiasHist = dynamic_cast<TH1D*>(resultsList->FindObject(detectorBiasHistName.Data()));
502 if(detectorBiasHist)
503 {
504 this->SetDetectorBiasHist(detectorBiasHist);
d66d46f7 505 }
ecac11c2 506 TString detectorBiasVsMHistName = "fDetectorBiasVsMHist";
507 TH1D *detectorBiasVsMHist = dynamic_cast<TH1D*>(resultsList->FindObject(detectorBiasVsMHistName.Data()));
508 if(detectorBiasVsMHist)
509 {
510 this->SetDetectorBiasVsMHist(detectorBiasVsMHist);
511 }
ecac11c2 512
d66d46f7 513} // end of void AliFlowAnalysisWithMixedHarmonics::GetPointersForResultsHistograms()
489d5531 514
d66d46f7 515//================================================================================================================
516
489d5531 517void AliFlowAnalysisWithMixedHarmonics::WriteHistograms(TString outputFileName)
518{
519 // Store the final results in output .root file.
520 TFile *output = new TFile(outputFileName.Data(),"RECREATE");
521 fHistList->Write(fHistList->GetName(),TObject::kSingleKey);
522 delete output;
523}
524
525//================================================================================================================
d66d46f7 526
527void AliFlowAnalysisWithMixedHarmonics::WriteHistograms(TDirectoryFile *outputFileName)
528{
529 // Store the final results in output .root file.
530 fHistList->SetName("cobjMH");
531 fHistList->SetOwner(kTRUE);
532 outputFileName->Add(fHistList);
533 outputFileName->Write(outputFileName->GetName(),TObject::kSingleKey);
534}
535
536//================================================================================================================
537
e999459d 538void AliFlowAnalysisWithMixedHarmonics::StoreHarmonic()
539{
540 // Store harmonic n used in cos[n*(phi1+phi2-2phi3)] and cos[n*(psi1+psi2-2phi3)].
541
542 (fCommonHists->GetHarmonic())->Fill(0.5,fHarmonic);
543
544} // end of void AliFlowAnalysisWithMixedHarmonics::StoreHarmonic()
545
546//================================================================================================================
547
5c078806 548void AliFlowAnalysisWithMixedHarmonics::InitializeArrays()
549{
550 // Initialize arrays.
551
552 for(Int_t sd=0;sd<2;sd++)
553 {
554 fRePEBE[sd] = NULL;
555 fImPEBE[sd] = NULL;
556 f3pCorrelatorVsPtSumDiffPro[sd] = NULL;
557 }
e999459d 558 for(Int_t fs=0;fs<2;fs++) // 1st/2nd POI which is also RP
559 {
560 for(Int_t sd=0;sd<2;sd++)
561 {
562 fOverlapEBE[fs][sd] = NULL;
563 }
564 }
5c078806 565
566} // end of AliFlowAnalysisWithMixedHarmonics::InitializeArrays()
567
568//================================================================================================================
569
489d5531 570void AliFlowAnalysisWithMixedHarmonics::BookAndNestAllLists()
d66d46f7 571{
572 // Book and nest all list in base list fHistList.
573
574 // Weights:
489d5531 575 fWeightsList->SetName("Weights");
576 fWeightsList->SetOwner(kTRUE);
d66d46f7 577 fHistList->Add(fWeightsList);
ecac11c2 578 // Profiles:
579 fProfileList->SetName("Profiles");
580 fProfileList->SetOwner(kTRUE);
581 fHistList->Add(fProfileList);
d66d46f7 582 // Results:
489d5531 583 fResultsList->SetName("Results");
584 fResultsList->SetOwner(kTRUE);
d66d46f7 585 fHistList->Add(fResultsList);
586
587} // end of void AliFlowAnalysisWithMixedHarmonics::BookAndNestAllLists()
588
589//================================================================================================================
590
489d5531 591void AliFlowAnalysisWithMixedHarmonics::BookProfileHoldingSettings()
d66d46f7 592{
593 // Book profile to hold all analysis settings.
594
489d5531 595 TString analysisSettingsName = "fAnalysisSettings";
e999459d 596 fAnalysisSettings = new TProfile(analysisSettingsName.Data(),"Settings for analysis with mixed harmonics",8,0,8);
597 fAnalysisSettings->SetStats(kFALSE);
598 fAnalysisSettings->GetXaxis()->SetLabelSize(0.03);
599 fAnalysisSettings->GetXaxis()->SetBinLabel(1,"Corr. for det. effects?");
600 fAnalysisSettings->Fill(0.5,(Int_t)fCorrectForDetectorEffects);
601 fAnalysisSettings->GetXaxis()->SetBinLabel(2,"# of mult. bins");
602 fAnalysisSettings->Fill(1.5,fNoOfMultipicityBins);
603 fAnalysisSettings->GetXaxis()->SetBinLabel(3,"Width of mult. bins");
604 fAnalysisSettings->Fill(2.5,fMultipicityBinWidth);
605 fAnalysisSettings->GetXaxis()->SetBinLabel(4,"Minimal mult.");
606 fAnalysisSettings->Fill(3.5,fMinMultiplicity);
607 fAnalysisSettings->GetXaxis()->SetBinLabel(5,"Print on the screen?");
608 fAnalysisSettings->Fill(4.5,(Int_t)fPrintOnTheScreen);
609 fAnalysisSettings->GetXaxis()->SetBinLabel(6,"fHarmonic");
610 fAnalysisSettings->Fill(5.5,(Int_t)fHarmonic);
611 fAnalysisSettings->GetXaxis()->SetBinLabel(7,"fOppositeChargesPOI");
612 fAnalysisSettings->Fill(6.5,(Int_t)fOppositeChargesPOI);
613 fAnalysisSettings->GetXaxis()->SetBinLabel(8,"fEvaluateDifferential3pCorrelator");
614 fAnalysisSettings->Fill(7.5,(Int_t)fOppositeChargesPOI);
d66d46f7 615 fHistList->Add(fAnalysisSettings);
616
617} // end of void AliFlowAnalysisWithMixedHarmonics::BookProfileHoldingSettings()
489d5531 618
d66d46f7 619//================================================================================================================
620
489d5531 621void AliFlowAnalysisWithMixedHarmonics::BookCommonHistograms()
622{
d66d46f7 623 // Book common control histograms and common histograms for final results.
489d5531 624
625 TString commonHistsName = "AliFlowCommonHistMH";
626 fCommonHists = new AliFlowCommonHist(commonHistsName.Data());
627 fHistList->Add(fCommonHists);
628
629} // end of void AliFlowAnalysisWithMixedHarmonics::BookCommonHistograms()
d66d46f7 630
631//================================================================================================================
632
489d5531 633void AliFlowAnalysisWithMixedHarmonics::BookAllEventByEventQuantities()
634{
d66d46f7 635 // Book all event-by-event quantitites.
636
e999459d 637 Int_t n = fHarmonic;
638
d66d46f7 639 // Q_{n,k} and S{p,k}:
e999459d 640 fReQnk = new TMatrixD(2*n,9); // to be improved(bound on k)
641 fImQnk = new TMatrixD(2*n,9); // to be improved(bound on k)
d66d46f7 642 fSpk = new TMatrixD(4,4); // to be improved(bound on p and k)
5c078806 643
644 // p_n vs [(p1+p2)/2,|p1-p2|]
e999459d 645 if(!fEvaluateDifferential3pCorrelator){return;}
5c078806 646 TString psdFlag[2] = {"PtSum","PtDiff"};
e999459d 647 TString fsFlag[2] = {"1st","2nd"};
5c078806 648 for(Int_t sd=0;sd<2;sd++)
649 {
5c078806 650 fRePEBE[sd] = new TProfile(Form("fRePEBE%s",psdFlag[sd].Data()),"",fnBinsPt,0.,fPtMax);
651 fImPEBE[sd] = new TProfile(Form("fImPEBE%s",psdFlag[sd].Data()),"",fnBinsPt,0.,fPtMax);
652 }
e999459d 653 for(Int_t fs=0;fs<2;fs++)
654 {
655 for(Int_t sd=0;sd<2;sd++)
656 {
657 fOverlapEBE[fs][sd] = new TProfile(Form("%s POI, %s",fsFlag[sd].Data(),psdFlag[sd].Data()),"",fnBinsPt,0.,fPtMax);
658 }
659 }
660
489d5531 661} // end fo void AliFlowAnalysisWithMixedHarmonics::BookAllEventByEventQuantities()
5c078806 662
d66d46f7 663//================================================================================================================
664
665void AliFlowAnalysisWithMixedHarmonics::BookAllAllEventQuantities()
666{
667 // Book all all-event quantitites.
668
669 // a) Book quantites without multiplicity binning;
5c078806 670 // b) Book quantites with multiplicity binning;
671 // c) Book quantities with binning in (p1+p2)/2, |p1-p2|, (eta1+eta2)/2 or |eta1-eta2|.
672
d66d46f7 673 // a) Book quantites without multiplicity binning:
e999459d 674 // 3-p correlator <<cos[n*(phi1+phi2-2phi3)]>> for all events:
489d5531 675 TString s3pCorrelatorProName = "f3pCorrelatorPro";
5c078806 676 f3pCorrelatorPro = new TProfile(s3pCorrelatorProName.Data(),"",1,0,1);
e999459d 677 f3pCorrelatorPro->SetStats(kFALSE);
5c078806 678 f3pCorrelatorPro->GetXaxis()->SetLabelOffset(0.01);
679 f3pCorrelatorPro->GetXaxis()->SetLabelSize(0.05);
e999459d 680 if(fHarmonic == 1)
5c078806 681 {
e999459d 682 f3pCorrelatorPro->GetXaxis()->SetBinLabel(1,"#LT#LTcos(#phi_{1}+#phi_{2}-2#phi_{3})#GT#GT");
5c078806 683 } else
684 {
e999459d 685 f3pCorrelatorPro->GetXaxis()->SetBinLabel(1,Form("#LT#LTcos[%i(#phi_{1}+#phi_{2}-2#phi_{3})]#GT#GT",fHarmonic));
5c078806 686 }
ecac11c2 687 fProfileList->Add(f3pCorrelatorPro);
e999459d 688 // Non-isotropic terms in the decomposition of <<cos[n(phi1+phi2-2phi3)]>>:
489d5531 689 TString nonIsotropicTermsProName = "fNonIsotropicTermsPro";
e999459d 690 fNonIsotropicTermsPro = new TProfile(nonIsotropicTermsProName.Data(),"",8,0,8);
691 fNonIsotropicTermsPro->SetStats(kFALSE);
692 if(fHarmonic == 1)
5c078806 693 {
e999459d 694 fNonIsotropicTermsPro->SetTitle("Non-isotropic terms in decomposition of #LT#LTcos(#phi_{1}+#phi_{2}-2#phi_{3})#GT#GT");
5c078806 695 fNonIsotropicTermsPro->GetXaxis()->SetBinLabel(1,"cos(#phi_{1})");
696 fNonIsotropicTermsPro->GetXaxis()->SetBinLabel(2,"sin(#phi_{1})");
697 fNonIsotropicTermsPro->GetXaxis()->SetBinLabel(3,"cos(2#phi_{1})");
698 fNonIsotropicTermsPro->GetXaxis()->SetBinLabel(4,"sin(2#phi_{1})");
699 fNonIsotropicTermsPro->GetXaxis()->SetBinLabel(5,"cos(#phi_{1}+#phi_{2})");
700 fNonIsotropicTermsPro->GetXaxis()->SetBinLabel(6,"sin(#phi_{1}+#phi_{2})");
701 fNonIsotropicTermsPro->GetXaxis()->SetBinLabel(7,"cos(2#phi_{1}-#phi_{2})");
702 fNonIsotropicTermsPro->GetXaxis()->SetBinLabel(8,"sin(2#phi_{1}-#phi_{2})");
e999459d 703 // fNonIsotropicTermsPro->GetXaxis()->SetBinLabel(9,"cos(#phi_{1}-#phi_{2}-#phi_{3})"); // not needed
704 // fNonIsotropicTermsPro->GetXaxis()->SetBinLabel(10,"sin(#phi_{1}-#phi_{2}-#phi_{3})"); // not needed
5c078806 705 } else
706 {
e999459d 707 fNonIsotropicTermsPro->SetTitle(Form("Non-isotropic terms in decomposition of #LT#LTcos[%i(#phi_{1}+#phi_{2}-2#phi_{3})]#GT#GT",fHarmonic));
708 fNonIsotropicTermsPro->GetXaxis()->SetBinLabel(1,Form("cos(%d#phi_{1})",fHarmonic));
709 fNonIsotropicTermsPro->GetXaxis()->SetBinLabel(2,Form("sin(%d#phi_{1})",fHarmonic));
710 fNonIsotropicTermsPro->GetXaxis()->SetBinLabel(3,Form("cos(%d#phi_{1})",2*fHarmonic));
711 fNonIsotropicTermsPro->GetXaxis()->SetBinLabel(4,Form("sin(%d#phi_{1})",2*fHarmonic));
712 fNonIsotropicTermsPro->GetXaxis()->SetBinLabel(5,Form("cos[%d(#phi_{1}+#phi_{2})]",fHarmonic));
713 fNonIsotropicTermsPro->GetXaxis()->SetBinLabel(6,Form("sin[%d(#phi_{1}+#phi_{2})]",fHarmonic));
714 fNonIsotropicTermsPro->GetXaxis()->SetBinLabel(7,Form("cos[%d(2#phi_{1}-#phi_{2})]",fHarmonic));
715 fNonIsotropicTermsPro->GetXaxis()->SetBinLabel(8,Form("sin[%d(2#phi_{1}-#phi_{2})]",fHarmonic));
716 // fNonIsotropicTermsPro->GetXaxis()->SetBinLabel(9,Form("cos(%d(#phi_{1}-#phi_{2}-#phi_{3}))",fHarmonic)); // not needed
717 // fNonIsotropicTermsPro->GetXaxis()->SetBinLabel(10,Form("sin(%d(#phi_{1}-#phi_{2}-#phi_{3}))",fHarmonic)); // not needed
5c078806 718 }
ecac11c2 719 fProfileList->Add(fNonIsotropicTermsPro);
e999459d 720 // 3-p correlator <<cos[n(phi1+phi2-2phi3)]>> corrected for detector effects:
5c078806 721 TString s3pCorrelatorHistName = "f3pCorrelatorHist";
722 f3pCorrelatorHist = new TH1D(s3pCorrelatorHistName.Data(),"",1,0,1);
e999459d 723 f3pCorrelatorHist->SetStats(kFALSE);
5c078806 724 f3pCorrelatorHist->GetXaxis()->SetLabelOffset(0.01);
725 f3pCorrelatorHist->GetXaxis()->SetLabelSize(0.05);
e999459d 726 if(fHarmonic == 1)
5c078806 727 {
e999459d 728 f3pCorrelatorHist->GetXaxis()->SetBinLabel(1,"#LT#LTcos(#phi_{1}+#phi_{2}-2#phi_{3})#GT#GT");
5c078806 729 } else
730 {
e999459d 731 f3pCorrelatorHist->GetXaxis()->SetBinLabel(1,Form("#LT#LTcos[%i(#phi_{1}+#phi_{2}-2#phi_{3})]#GT#GT",fHarmonic));
5c078806 732 }
733 fResultsList->Add(f3pCorrelatorHist);
e999459d 734 // Quantified bias coming from detector inefficiencies to 3-p correlator <<cos[n(phi1+phi2-2phi3)]>>:
489d5531 735 TString detectorBiasHistName = "fDetectorBiasHist";
d66d46f7 736 fDetectorBiasHist = new TH1D(detectorBiasHistName.Data(),"Bias coming from detector inefficiences",1,0,1);
e999459d 737 fDetectorBiasHist->SetStats(kFALSE);
738 if(fHarmonic == 1)
739 {
740 fDetectorBiasHist->GetXaxis()->SetBinLabel(1,"#frac{corrected}{measured} #LT#LTcos(#phi_{1}+#phi_{2}-2#phi_{3})#GT#GT");
741 } else
742 {
743 fDetectorBiasHist->GetXaxis()->SetBinLabel(1,Form("#frac{corrected}{measured} #LT#LTcos[%i(#phi_{1}+#phi_{2}-2#phi_{3})]#GT#GT",fHarmonic));
744 }
d66d46f7 745 fResultsList->Add(fDetectorBiasHist);
746
747 // b) Book quantites with multiplicity binning.
e999459d 748 // 3-p correlator <<cos[n(phi1+phi2-2phi3)]>> versus multiplicity:
489d5531 749 TString s3pCorrelatorVsMProName = "f3pCorrelatorVsMPro";
5c078806 750 f3pCorrelatorVsMPro = new TProfile(s3pCorrelatorVsMProName.Data(),"",fNoOfMultipicityBins+2,0,fNoOfMultipicityBins+2);
e999459d 751 f3pCorrelatorVsMPro->SetStats(kFALSE);
752 if(fHarmonic == 1)
d66d46f7 753 {
e999459d 754 f3pCorrelatorVsMPro->SetTitle("#LT#LTcos(#phi_{1}+#phi_{2}-2#phi_{3})#GT#GT #font[72]{vs} M");
5c078806 755 } else
756 {
e999459d 757 f3pCorrelatorVsMPro->SetTitle(Form("#LT#LTcos[%d(#phi_{1}+#phi_{2}-2#phi_{3})]#GT#GT #font[72]{vs} M",fHarmonic));
5c078806 758 }
759 f3pCorrelatorVsMPro->GetXaxis()->SetBinLabel(1,Form("M < %d",(Int_t)fMinMultiplicity));
760 for(Int_t b=2;b<=fNoOfMultipicityBins+1;b++)
d66d46f7 761 {
5c078806 762 f3pCorrelatorVsMPro->GetXaxis()->SetBinLabel(b,Form("%d #leq M < %d",(Int_t)(fMinMultiplicity+(b-2)*fMultipicityBinWidth),(Int_t)(fMinMultiplicity+(b-1)*fMultipicityBinWidth)));
d66d46f7 763 }
5c078806 764 f3pCorrelatorVsMPro->GetXaxis()->SetBinLabel(fNoOfMultipicityBins+2,Form(" M #geq %d",(Int_t)(fMinMultiplicity+fNoOfMultipicityBins*fMultipicityBinWidth)));
e558a733 765 fProfileList->Add(f3pCorrelatorVsMPro);
e999459d 766 // 3-p correlator <<cos[n(phi1+phi2-2phi3)]>> versus multiplicity corrected for detector effects:
e558a733 767 TString s3pCorrelatorVsMHistName = "f3pCorrelatorVsMHist";
768 f3pCorrelatorVsMHist = new TH1D(s3pCorrelatorVsMHistName.Data(),"",fNoOfMultipicityBins+2,0,fNoOfMultipicityBins+2);
e999459d 769 f3pCorrelatorVsMHist->SetStats(kFALSE);
770 if(fHarmonic == 1)
e558a733 771 {
e999459d 772 f3pCorrelatorVsMHist->SetTitle("cos(#phi_{1}+#phi_{2}-2#phi_{3}) #font[72]{vs} M");
e558a733 773 } else
774 {
e999459d 775 f3pCorrelatorVsMHist->SetTitle(Form("cos[%d(#phi_{1}+#phi_{2}-2#phi_{3})] #font[72]{vs} M",fHarmonic));
e558a733 776 }
777 f3pCorrelatorVsMHist->GetXaxis()->SetBinLabel(1,Form("M < %d",(Int_t)fMinMultiplicity));
778 for(Int_t b=2;b<=fNoOfMultipicityBins+1;b++)
779 {
780 f3pCorrelatorVsMHist->GetXaxis()->SetBinLabel(b,Form("%d #leq M < %d",(Int_t)(fMinMultiplicity+(b-2)*fMultipicityBinWidth),(Int_t)(fMinMultiplicity+(b-1)*fMultipicityBinWidth)));
781 }
782 f3pCorrelatorVsMHist->GetXaxis()->SetBinLabel(fNoOfMultipicityBins+2,Form(" M #geq %d",(Int_t)(fMinMultiplicity+fNoOfMultipicityBins*fMultipicityBinWidth)));
783 fResultsList->Add(f3pCorrelatorVsMHist);
e999459d 784 // Non-isotropic terms in the decomposition of <<cos[n(phi1+phi2-2phi3)]>> vs multiplicity:
489d5531 785 TString nonIsotropicTermsVsMProName = "fNonIsotropicTermsVsMPro";
e999459d 786 fNonIsotropicTermsVsMPro = new TProfile2D(nonIsotropicTermsVsMProName.Data(),"",8,0,8,fNoOfMultipicityBins+2,0,fNoOfMultipicityBins+2);
787 fNonIsotropicTermsVsMPro->SetStats(kFALSE);
788 if(fHarmonic == 1)
d66d46f7 789 {
e999459d 790 fNonIsotropicTermsVsMPro->SetTitle("Non-isotropic terms in decomposition of #LT#LTcos(#phi_{1}+#phi_{2}-2#phi_{3})#GT#GT #font[72]{vs} M");
5c078806 791 fNonIsotropicTermsVsMPro->GetXaxis()->SetBinLabel(1,"cos(#phi_{1})");
792 fNonIsotropicTermsVsMPro->GetXaxis()->SetBinLabel(2,"sin(#phi_{1})");
793 fNonIsotropicTermsVsMPro->GetXaxis()->SetBinLabel(3,"cos(2#phi_{1})");
794 fNonIsotropicTermsVsMPro->GetXaxis()->SetBinLabel(4,"sin(2#phi_{1})");
795 fNonIsotropicTermsVsMPro->GetXaxis()->SetBinLabel(5,"cos(#phi_{1}+#phi_{2})");
796 fNonIsotropicTermsVsMPro->GetXaxis()->SetBinLabel(6,"sin(#phi_{1}+#phi_{2})");
797 fNonIsotropicTermsVsMPro->GetXaxis()->SetBinLabel(7,"cos(2#phi_{1}-#phi_{2})");
798 fNonIsotropicTermsVsMPro->GetXaxis()->SetBinLabel(8,"sin(2#phi_{1}-#phi_{2})");
e999459d 799 // fNonIsotropicTermsVsMPro->GetXaxis()->SetBinLabel(9,"cos(#phi_{1}-#phi_{2}-#phi_{3})"); // not needed
800 // fNonIsotropicTermsVsMPro->GetXaxis()->SetBinLabel(10,"sin(#phi_{1}-#phi_{2}-#phi_{3})"); // not needed
5c078806 801 } else
802 {
e999459d 803 fNonIsotropicTermsVsMPro->SetTitle(Form("Non-isotropic terms in decomposition of #LT#LTcos[%d(#phi_{1}+#phi_{2}-2#phi_{3})]#GT#GT",fHarmonic));
804 fNonIsotropicTermsVsMPro->GetXaxis()->SetBinLabel(1,Form("cos(%d#phi_{1})",fHarmonic));
805 fNonIsotropicTermsVsMPro->GetXaxis()->SetBinLabel(2,Form("sin(%d#phi_{1})",fHarmonic));
806 fNonIsotropicTermsVsMPro->GetXaxis()->SetBinLabel(3,Form("cos(%d#phi_{1})",2*fHarmonic));
807 fNonIsotropicTermsVsMPro->GetXaxis()->SetBinLabel(4,Form("sin(%d#phi_{1})",2*fHarmonic));
808 fNonIsotropicTermsVsMPro->GetXaxis()->SetBinLabel(5,Form("cos[%d(#phi_{1}+#phi_{2})]",fHarmonic));
809 fNonIsotropicTermsVsMPro->GetXaxis()->SetBinLabel(6,Form("sin[%d(#phi_{1}+#phi_{2})]",fHarmonic));
810 fNonIsotropicTermsVsMPro->GetXaxis()->SetBinLabel(7,Form("cos[%d(2#phi_{1}-#phi_{2})]",fHarmonic));
811 fNonIsotropicTermsVsMPro->GetXaxis()->SetBinLabel(8,Form("sin[%d(2#phi_{1}-#phi_{2})]",fHarmonic));
812 // fNonIsotropicTermsVsMPro->GetXaxis()->SetBinLabel(9,Form("cos(%d(#phi_{1}-#phi_{2}-#phi_{3}))",fHarmonic)); // not needed
813 // fNonIsotropicTermsVsMPro->GetXaxis()->SetBinLabel(10,Form("sin(%d(#phi_{1}-#phi_{2}-#phi_{3}))",fHarmonic)); // not needed
5c078806 814 }
815 fNonIsotropicTermsVsMPro->GetYaxis()->SetBinLabel(1,Form("M < %d",(Int_t)fMinMultiplicity));
816 for(Int_t b=2;b<=fNoOfMultipicityBins+1;b++)
d66d46f7 817 {
e558a733 818 fNonIsotropicTermsVsMPro->GetYaxis()->SetBinLabel(b,Form("%d #leq M < %d",(Int_t)(fMinMultiplicity+(b-2)*fMultipicityBinWidth),(Int_t)(fMinMultiplicity+(b-1)*fMultipicityBinWidth)));
d66d46f7 819 }
5c078806 820 fNonIsotropicTermsVsMPro->GetYaxis()->SetBinLabel(fNoOfMultipicityBins+2,Form(" M #geq %d",(Int_t)(fMinMultiplicity+fNoOfMultipicityBins*fMultipicityBinWidth)));
821 fProfileList->Add(fNonIsotropicTermsVsMPro);
e999459d 822 // Quantified bias coming from detector inefficiencies to 3-p correlator <<cos[n(phi1+phi2-2phi3)]>> vs multiplicity
489d5531 823 TString detectorBiasVsMHistName = "fDetectorBiasVsMHist";
e558a733 824 fDetectorBiasVsMHist = new TH1D(detectorBiasVsMHistName.Data(),"",fNoOfMultipicityBins+2,0,fNoOfMultipicityBins+2);
e999459d 825 fDetectorBiasVsMHist->SetStats(kFALSE);
826 if(fHarmonic == 1)
e558a733 827 {
e999459d 828 fDetectorBiasVsMHist->SetTitle("#frac{corrected}{measured} #LT#LTcos(#phi_{1}+#phi_{2}-2#phi_{3})#GT#GT #font[72]{vs} M");
e558a733 829 } else
830 {
e999459d 831 fDetectorBiasVsMHist->SetTitle(Form("#frac{corrected}{measured} #LT#LTcos[%d(#phi_{1}+#phi_{2}-2#phi_{3})]#GT#GT #font[72]{vs} M",fHarmonic));
e558a733 832 }
833 fDetectorBiasVsMHist->GetXaxis()->SetBinLabel(1,Form("M < %d",(Int_t)fMinMultiplicity));
834 for(Int_t b=2;b<=fNoOfMultipicityBins+1;b++)
d66d46f7 835 {
e558a733 836 fDetectorBiasVsMHist->GetXaxis()->SetBinLabel(b,Form("%d #leq M < %d",(Int_t)(fMinMultiplicity+(b-2)*fMultipicityBinWidth),(Int_t)(fMinMultiplicity+(b-1)*fMultipicityBinWidth)));
d66d46f7 837 }
e558a733 838 fDetectorBiasVsMHist->GetXaxis()->SetBinLabel(fNoOfMultipicityBins+2,Form(" M #geq %d",(Int_t)(fMinMultiplicity+fNoOfMultipicityBins*fMultipicityBinWidth)));
839 fResultsList->Add(fDetectorBiasVsMHist);
d66d46f7 840
5c078806 841 // c) Book quantities with binning in (p1+p2)/2, |p1-p2|, (eta1+eta2)/2 or |eta1-eta2|:
e999459d 842 if(!fEvaluateDifferential3pCorrelator){return;}
5c078806 843 TString psdFlag[2] = {"PtSum","PtDiff"};
e999459d 844 TString psdTitleFlag[2] = {"(p_{T,1}+ p_{T,2})/2","#left|p_{T,1}- p_{T,2}#right|"};
5c078806 845 //TString s3pCorrelatorVsPtSumDiffProName = "f3pCorrelatorVsPtSumDiffPro";
846 for(Int_t sd=0;sd<2;sd++)
847 {
5c078806 848 f3pCorrelatorVsPtSumDiffPro[sd] = new TProfile(Form("f3pCorrelatorVs%sPro",psdFlag[sd].Data()),"",fnBinsPt,0.,fPtMax);
e999459d 849 f3pCorrelatorVsPtSumDiffPro[sd]->SetStats(kFALSE);
5c078806 850 //f3pCorrelatorVsPtSumDiffPro[sd]->SetLabelSize(0.05);
851 //f3pCorrelatorVsPtSumDiffPro[sd]->SetMarkerStyle(25);
e999459d 852 if(fHarmonic == 1)
853 {
854 f3pCorrelatorVsPtSumDiffPro[sd]->SetTitle(Form("#LT#LTcos(#psi_{1}+#psi_{2}-2#phi_{3})#GT#GT #font[72]{vs} %s",psdTitleFlag[sd].Data()));
855 } else
856 {
857 f3pCorrelatorVsPtSumDiffPro[sd]->SetTitle(Form("#LT#LTcos[%d(#psi_{1}+#psi_{2}-2#phi_{3})]#GT#GT #font[72]{vs} %s",fHarmonic,psdTitleFlag[sd].Data()));
858 }
5c078806 859 f3pCorrelatorVsPtSumDiffPro[sd]->GetXaxis()->SetTitle(psdTitleFlag[sd].Data());
860 fProfileList->Add(f3pCorrelatorVsPtSumDiffPro[sd]);
861 }
862
d66d46f7 863} // end of void AliFlowAnalysisWithMixedHarmonics::BookAllAllEventQuantities()
489d5531 864
d66d46f7 865//================================================================================================================
866
489d5531 867void AliFlowAnalysisWithMixedHarmonics::AccessConstants()
868{
869 // Access needed common constants from AliFlowCommonConstants.
870
871 fnBinsPhi = AliFlowCommonConstants::GetMaster()->GetNbinsPhi();
872 fPhiMin = AliFlowCommonConstants::GetMaster()->GetPhiMin();
873 fPhiMax = AliFlowCommonConstants::GetMaster()->GetPhiMax();
874 if(fnBinsPhi) fPhiBinWidth = (fPhiMax-fPhiMin)/fnBinsPhi;
875 fnBinsPt = AliFlowCommonConstants::GetMaster()->GetNbinsPt();
876 fPtMin = AliFlowCommonConstants::GetMaster()->GetPtMin();
877 fPtMax = AliFlowCommonConstants::GetMaster()->GetPtMax();
878 if(fnBinsPt) fPtBinWidth = (fPtMax-fPtMin)/fnBinsPt;
879 fnBinsEta = AliFlowCommonConstants::GetMaster()->GetNbinsEta();
880 fEtaMin = AliFlowCommonConstants::GetMaster()->GetEtaMin();
881 fEtaMax = AliFlowCommonConstants::GetMaster()->GetEtaMax();
882 if(fnBinsEta) fEtaBinWidth = (fEtaMax-fEtaMin)/fnBinsEta;
883
d66d46f7 884} // end of void AliFlowAnalysisWithMixedHarmonics::AccessConstants()
885
489d5531 886//================================================================================================================
887
888void AliFlowAnalysisWithMixedHarmonics::CrossCheckSettings()
889{
890 // Cross-check if the user settings make sense.
891
e999459d 892 // ...
893
489d5531 894} // end of void AliFlowAnalysisWithMixedHarmonics::CrossCheckSettings()
895
896//================================================================================================================
897
898void AliFlowAnalysisWithMixedHarmonics::BookAndFillWeightsHistograms()
899{
900 // Book and fill (by accessing file "weights.root") histograms which hold phi, pt and eta weights.
901
902 if(!fWeightsList)
903 {
e999459d 904 cout<<endl;
905 cout<<" WARNING (MH): fWeightsList is NULL in BookAndFillWeightsHistograms() !!!!"<<endl;
906 cout<<endl;
489d5531 907 exit(0);
908 }
909 // Profile to hold flags for weights:
910 TString fUseParticleWeightsName = "fUseParticleWeightsMH";
911 fUseParticleWeights = new TProfile(fUseParticleWeightsName.Data(),"0 = particle weight not used, 1 = particle weight used ",3,0,3);
e999459d 912 fUseParticleWeights->SetStats(kFALSE);
489d5531 913 fUseParticleWeights->SetLabelSize(0.06);
914 (fUseParticleWeights->GetXaxis())->SetBinLabel(1,"w_{#phi}");
915 (fUseParticleWeights->GetXaxis())->SetBinLabel(2,"w_{p_{T}}");
916 (fUseParticleWeights->GetXaxis())->SetBinLabel(3,"w_{#eta}");
917 fUseParticleWeights->Fill(0.5,(Int_t)fUsePhiWeights);
918 fUseParticleWeights->Fill(1.5,(Int_t)fUsePtWeights);
919 fUseParticleWeights->Fill(2.5,(Int_t)fUseEtaWeights);
920 fWeightsList->Add(fUseParticleWeights);
921 // Phi-weights:
922 if(fUsePhiWeights)
923 {
924 if(fWeightsList->FindObject("phi_weights"))
925 {
926 fPhiWeights = dynamic_cast<TH1F*>(fWeightsList->FindObject("phi_weights"));
927 if(TMath::Abs(fPhiWeights->GetBinWidth(1)-fPhiBinWidth)>pow(10.,-6.))
928 {
d66d46f7 929 cout<<endl;
e999459d 930 cout<<" WARNING (MH): Inconsistent binning in histograms for phi-weights throughout the code."<<endl;
d66d46f7 931 cout<<endl;
489d5531 932 exit(0);
933 }
934 } else
935 {
e999459d 936 cout<<endl;
937 cout<<" WARNING (MH): fWeightsList->FindObject(\"phi_weights\") is NULL in BookAndFillWeightsHistograms() !!!!"<<endl;
938 cout<<endl;
489d5531 939 exit(0);
940 }
941 } // end of if(fUsePhiWeights)
942 // Pt-weights:
943 if(fUsePtWeights)
944 {
945 if(fWeightsList->FindObject("pt_weights"))
946 {
947 fPtWeights = dynamic_cast<TH1D*>(fWeightsList->FindObject("pt_weights"));
948 if(TMath::Abs(fPtWeights->GetBinWidth(1)-fPtBinWidth)>pow(10.,-6.))
949 {
d66d46f7 950 cout<<endl;
e999459d 951 cout<<" WARNING (MH): Inconsistent binning in histograms for pt-weights throughout the code."<<endl;
d66d46f7 952 cout<<endl;
489d5531 953 exit(0);
954 }
955 } else
956 {
e999459d 957 cout<<endl;
958 cout<<" WARNING (MH): fWeightsList->FindObject(\"pt_weights\") is NULL in BookAndFillWeightsHistograms() !!!!"<<endl;
959 cout<<endl;
489d5531 960 exit(0);
961 }
962 } // end of if(fUsePtWeights)
963 // Eta-weights:
964 if(fUseEtaWeights)
965 {
966 if(fWeightsList->FindObject("eta_weights"))
967 {
968 fEtaWeights = dynamic_cast<TH1D*>(fWeightsList->FindObject("eta_weights"));
969 if(TMath::Abs(fEtaWeights->GetBinWidth(1)-fEtaBinWidth)>pow(10.,-6.))
970 {
d66d46f7 971 cout<<endl;
e999459d 972 cout<<" WARNING (MH): Inconsistent binning in histograms for eta-weights throughout the code."<<endl;
d66d46f7 973 cout<<endl;
489d5531 974 exit(0);
975 }
976 } else
977 {
e999459d 978 cout<<endl;
979 cout<<" WARNING (MH): fUseEtaWeights && fWeightsList->FindObject(\"eta_weights\") is NULL in BookAndFillWeightsHistograms() !!!!"<<endl;
980 cout<<endl;
489d5531 981 exit(0);
982 }
983 } // end of if(fUseEtaWeights)
984
d66d46f7 985} // end of AliFlowAnalysisWithMixedHarmonics::BookAndFillWeightsHistograms()
986
987//================================================================================================================
489d5531 988
989void AliFlowAnalysisWithMixedHarmonics::CheckPointersUsedInMake()
990{
991 // Check pointers used in method Make().
d66d46f7 992
e558a733 993 if(!fReQnk || !fImQnk || !fSpk )
d66d46f7 994 {
995 cout<<endl;
e999459d 996 cout<<" WARNING (MH): fReQnk || fImQnk || fSpk is NULL in CheckPointersUsedInMake() !!!!"<<endl;
d66d46f7 997 cout<<endl;
998 exit(0);
999 }
1000 if(!f3pCorrelatorPro)
1001 {
1002 cout<<endl;
e999459d 1003 cout<<" WARNING (MH): f3pCorrelatorPro is NULL in CheckPointersUsedInMake() !!!!"<<endl;
d66d46f7 1004 cout<<endl;
1005 exit(0);
5c078806 1006 }
1007 if(!fNonIsotropicTermsPro)
1008 {
1009 cout<<endl;
e999459d 1010 cout<<" WARNING (MH): fNonIsotropicTermsPro is NULL in CheckPointersUsedInMake() !!!!"<<endl;
5c078806 1011 cout<<endl;
1012 exit(0);
1013 }
1014 if(!f3pCorrelatorVsMPro)
1015 {
1016 cout<<endl;
e999459d 1017 cout<<" WARNING (MH): f3pCorrelatorVsMPro is NULL in CheckPointersUsedInMake() !!!!"<<endl;
5c078806 1018 cout<<endl;
1019 exit(0);
1020 }
1021 if(!fNonIsotropicTermsVsMPro)
1022 {
1023 cout<<endl;
e999459d 1024 cout<<" WARNING (MH): fNonIsotropicTermsVsMPro is NULL in CheckPointersUsedInMake() !!!!"<<endl;
5c078806 1025 cout<<endl;
1026 exit(0);
1027 }
e999459d 1028 if(!fEvaluateDifferential3pCorrelator){return;}
5c078806 1029 for(Int_t sd=0;sd<2;sd++)
1030 {
1031 if(!(f3pCorrelatorVsPtSumDiffPro[sd]))
1032 {
1033 cout<<endl;
e999459d 1034 cout<<" WARNING (MH): "<<Form("f3pCorrelatorVsPtSumDiffPro[%d]",sd)<<" is NULL in CheckPointersUsedInMake() !!!!"<<endl;
5c078806 1035 cout<<endl;
1036 exit(0);
1037 }
d66d46f7 1038 }
e999459d 1039 for(Int_t sd=0;sd<2;sd++)
1040 {
1041 if(!fRePEBE[sd]||!fImPEBE[sd])
1042 {
1043 cout<<endl;
1044 cout<<" WARNING (MH): "<<Form("!fRePEBE[%d]||!fImPEBE[%d]",sd)<<" is NULL in CheckPointersUsedInMake() !!!!"<<endl;
1045 cout<<endl;
1046 exit(0);
1047 }
1048 for(Int_t fs=0;fs<2;fs++)
1049 {
1050 if(!fOverlapEBE[fs][sd]||!fOverlapEBE[fs][sd])
1051 {
1052 cout<<endl;
1053 cout<<" WARNING (MH): "<<Form("!fOverlapEBE[%d][%d]||!fOverlapEBE[%d][%d]",fs,sd)<<" is NULL in CheckPointersUsedInMake() !!!!"<<endl;
1054 cout<<endl;
1055 exit(0);
1056 }
1057 } // end of for(Int_t fs=0;fs<2;fs++)
1058 } // end of for(Int_t sd=0;sd<2;sd++)
489d5531 1059
d66d46f7 1060} // end of AliFlowAnalysisWithMixedHarmonics::CheckPointersUsedInMake()
1061
1062//================================================================================================================
489d5531 1063
1064void AliFlowAnalysisWithMixedHarmonics::CheckPointersUsedInFinish()
1065{
d66d46f7 1066 // Check pointers used in method Finish().
ecac11c2 1067
1068 if(!fAnalysisSettings)
1069 {
1070 cout<<endl;
e999459d 1071 cout<<" WARNING (MH): fAnalysisSettings is NULL in CheckPointersUsedInFinish() !!!!"<<endl;
ecac11c2 1072 cout<<endl;
1073 exit(0);
1074 }
1075 if(!f3pCorrelatorPro)
1076 {
1077 cout<<endl;
e999459d 1078 cout<<" WARNING (MH): f3pCorrelatorPro is NULL in CheckPointersUsedInFinish() !!!!"<<endl;
ecac11c2 1079 cout<<endl;
1080 exit(0);
1081 }
1082 if(!fNonIsotropicTermsPro)
1083 {
1084 cout<<endl;
e999459d 1085 cout<<" WARNING (MH): fNonIsotropicTermsPro is NULL in CheckPointersUsedInFinish() !!!!"<<endl;
ecac11c2 1086 cout<<endl;
1087 exit(0);
1088 }
1089 if(!f3pCorrelatorVsMPro)
1090 {
1091 cout<<endl;
e999459d 1092 cout<<" WARNING (MH): f3pCorrelatorVsMPro is NULL in CheckPointersUsedInFinish() !!!!"<<endl;
ecac11c2 1093 cout<<endl;
1094 exit(0);
1095 }
e558a733 1096 if(!f3pCorrelatorVsMHist)
1097 {
1098 cout<<endl;
e999459d 1099 cout<<" WARNING (MH): f3pCorrelatorVsMHist is NULL in CheckPointersUsedInFinish() !!!!"<<endl;
e558a733 1100 cout<<endl;
1101 exit(0);
1102 }
ecac11c2 1103 if(!fNonIsotropicTermsVsMPro)
1104 {
1105 cout<<endl;
e999459d 1106 cout<<" WARNING (MH): fNonIsotropicTermsVsMPro is NULL in CheckPointersUsedInFinish() !!!!"<<endl;
ecac11c2 1107 cout<<endl;
1108 exit(0);
1109 }
1110 if(!f3pCorrelatorHist)
1111 {
1112 cout<<endl;
e999459d 1113 cout<<" WARNING (MH): f3pCorrelatorHist is NULL in CheckPointersUsedInFinish() !!!!"<<endl;
ecac11c2 1114 cout<<endl;
1115 exit(0);
5c078806 1116 }
ecac11c2 1117 if(!fDetectorBiasHist)
1118 {
1119 cout<<endl;
e999459d 1120 cout<<" WARNING (MH): fDetectorBiasHist is NULL in CheckPointersUsedInFinish() !!!!"<<endl;
ecac11c2 1121 cout<<endl;
1122 exit(0);
1123 }
5c078806 1124 /* to be improved - enabled eventually
ecac11c2 1125 if(!fDetectorBiasVsMHist)
1126 {
1127 cout<<endl;
e999459d 1128 cout<<" WARNING (MH): !fDetectorBiasVsMHist is NULL in CheckPointersUsedInFinish() !!!!"<<endl;
ecac11c2 1129 cout<<endl;
1130 exit(0);
5c078806 1131 }
1132 */
e999459d 1133 if(!fEvaluateDifferential3pCorrelator){return;}
5c078806 1134 for(Int_t sd=0;sd<2;sd++)
1135 {
1136 if(!(f3pCorrelatorVsPtSumDiffPro[sd]))
1137 {
1138 cout<<endl;
e999459d 1139 cout<<" WARNING (MH): "<<Form("f3pCorrelatorVsPtSumDiffPro[%d]",sd)<<" is NULL in CheckPointersUsedInFinish() !!!!"<<endl;
5c078806 1140 cout<<endl;
1141 exit(0);
1142 }
1143 }
ecac11c2 1144
489d5531 1145} // end of AliFlowAnalysisWithMixedHarmonics::CheckPointersUsedInFinish()
d66d46f7 1146
1147//================================================================================================================
489d5531 1148
5c078806 1149void AliFlowAnalysisWithMixedHarmonics::PrintOnTheScreen()
1150{
1151 // Print the final results on the screen.
1152
1153 cout<<endl;
1154 cout<<"*******************************************************"<<endl;
1155 cout<<"*******************************************************"<<endl;
1156 cout<<" Mixed Harmonics "<<endl;
1157 cout<<endl;
e999459d 1158 if(fHarmonic!=1)
5c078806 1159 {
e999459d 1160 cout<<" cos["<<fHarmonic<<"(phi1+phi2-2phi3)] = "<<f3pCorrelatorHist->GetBinContent(1)<<
5c078806 1161 " +/- "<<f3pCorrelatorHist->GetBinError(1)<<endl;
1162 } else
1163 {
e999459d 1164 cout<<" cos(phi1+phi2-2phi3) = "<<f3pCorrelatorHist->GetBinContent(1)<<" +/- "<<f3pCorrelatorHist->GetBinError(1)<<endl;
5c078806 1165 }
1166 cout<<" Detector Bias = "<<fDetectorBiasHist->GetBinContent(1)<<endl;
1167 cout<<endl;
e999459d 1168 cout<<" nEvts = "<<(Int_t)fCommonHists->GetHistMultRP()->GetEntries()<<", <M> = "<<(Double_t)fCommonHists->GetHistMultRP()->GetMean()<<endl;
5c078806 1169 cout<<"*******************************************************"<<endl;
1170 cout<<"*******************************************************"<<endl;
1171
1172} // end of void AliFlowAnalysisWithMixedHarmonics::PrintOnTheScreen()
1173
1174//================================================================================================================
1175
489d5531 1176void AliFlowAnalysisWithMixedHarmonics::AccessSettings()
1177{
d66d46f7 1178 // Access the settings for analysis with mixed harmonics.
1179
e999459d 1180 fCorrectForDetectorEffects = (Bool_t)fAnalysisSettings->GetBinContent(1);
1181 fNoOfMultipicityBins = (Int_t)fAnalysisSettings->GetBinContent(2);
1182 fMultipicityBinWidth = (Double_t)fAnalysisSettings->GetBinContent(3);
1183 fMinMultiplicity = (Double_t)fAnalysisSettings->GetBinContent(4);
1184 fPrintOnTheScreen = (Bool_t)fAnalysisSettings->GetBinContent(5);
1185 fHarmonic = (Int_t)fAnalysisSettings->GetBinContent(6);
1186 fOppositeChargesPOI = (Bool_t)fAnalysisSettings->GetBinContent(7);
1187 fEvaluateDifferential3pCorrelator = (Bool_t)fAnalysisSettings->GetBinContent(8);
489d5531 1188
d66d46f7 1189} // end of AliFlowAnalysisWithMixedHarmonics::AccessSettings()
1190
1191//================================================================================================================
489d5531 1192
1193void AliFlowAnalysisWithMixedHarmonics::CorrectForDetectorEffects()
1194{
e999459d 1195 // Correct measured 3-p correlator cos[n(phi1+phi2-2phi3)] for detector effects.
5c078806 1196
d66d46f7 1197 Double_t measured3pCorrelator = f3pCorrelatorPro->GetBinContent(1); // biased by detector effects
1198 Double_t corrected3pCorrelator = 0.; // corrected for detector effects
1199 Double_t nonIsotropicTerms[10] = {0.}; // there are 10 distinct non-isotropic terms
1200 for(Int_t nit=0;nit<10;nit++)
1201 {
1202 nonIsotropicTerms[nit] = fNonIsotropicTermsPro->GetBinContent(nit+1);
1203 }
1204 // Calculate corrected 3-p correlator:
1205 corrected3pCorrelator = measured3pCorrelator
1206 - nonIsotropicTerms[2]*nonIsotropicTerms[4]
1207 - nonIsotropicTerms[3]*nonIsotropicTerms[5]
1208 - 2.*nonIsotropicTerms[0]*nonIsotropicTerms[6]
1209 - 2.*nonIsotropicTerms[1]*nonIsotropicTerms[7]
1210 + 2.*nonIsotropicTerms[2]*(pow(nonIsotropicTerms[0],2.)-pow(nonIsotropicTerms[1],2.))
5c078806 1211 + 4.*nonIsotropicTerms[3]*nonIsotropicTerms[0]*nonIsotropicTerms[1];
1212 // Store corrected correlator:
1213 f3pCorrelatorHist->SetBinContent(1,corrected3pCorrelator);
1214 f3pCorrelatorHist->SetBinError(1,f3pCorrelatorPro->GetBinError(1)); // to be improved (propagate error for non-isotropic terms)
1215 // Quantify bias from detector inefficiences to 3-p correlator. Remark: Bias is quantified as a
1216 // ratio between corrected and measured 3-p correlator:
1217 // bias = corrected/measured
1218 // This bias is stored in histogram fDetectorBias.
1219 Double_t bias = 0.;
1220 if(measured3pCorrelator)
1221 {
1222 bias = corrected3pCorrelator/measured3pCorrelator;
1223 fDetectorBiasHist->SetBinContent(1,bias);
e558a733 1224 }
489d5531 1225
d66d46f7 1226} // end of AliFlowAnalysisWithMixedHarmonics::CorrectForDetectorEffects()
1227
1228//================================================================================================================
489d5531 1229
e558a733 1230void AliFlowAnalysisWithMixedHarmonics::CorrectForDetectorEffectsVsM()
1231{
e999459d 1232 // Correct measured 3-p correlator cos[n(phi1+phi2-2phi3)] vs M for detector effects.
e558a733 1233
1234 for(Int_t b=1;b<=fNoOfMultipicityBins+2;b++)
1235 {
1236 Double_t measured3pCorrelator = f3pCorrelatorVsMPro->GetBinContent(b); // biased by detector effects
1237 Double_t corrected3pCorrelator = 0.; // corrected for detector effects
1238 Double_t nonIsotropicTerms[10] = {0.}; // there are 10 distinct non-isotropic terms
1239 for(Int_t nit=0;nit<10;nit++)
1240 {
1241 nonIsotropicTerms[nit] = fNonIsotropicTermsVsMPro->GetBinContent(fNonIsotropicTermsVsMPro->GetBin(nit+1,b));
1242 }
1243 // Calculate corrected 3-p correlator:
1244 corrected3pCorrelator = measured3pCorrelator
1245 - nonIsotropicTerms[2]*nonIsotropicTerms[4]
1246 - nonIsotropicTerms[3]*nonIsotropicTerms[5]
1247 - 2.*nonIsotropicTerms[0]*nonIsotropicTerms[6]
1248 - 2.*nonIsotropicTerms[1]*nonIsotropicTerms[7]
1249 + 2.*nonIsotropicTerms[2]*(pow(nonIsotropicTerms[0],2.)-pow(nonIsotropicTerms[1],2.))
1250 + 4.*nonIsotropicTerms[3]*nonIsotropicTerms[0]*nonIsotropicTerms[1];
1251 // Store corrected correlator:
1252 f3pCorrelatorVsMHist->SetBinContent(b,corrected3pCorrelator);
1253 f3pCorrelatorVsMHist->SetBinError(b,f3pCorrelatorVsMPro->GetBinError(b)); // to be improved (propagate error for non-isotropic terms)
1254 // Quantify bias from detector inefficiences to 3-p correlator. Remark: Bias is quantified as a
1255 // ratio between corrected and measured 3-p correlator:
1256 // bias = corrected/measured
1257 // This bias is stored in histogram fDetectorBias.
1258 Double_t bias = 0.;
1259 if(measured3pCorrelator)
1260 {
1261 bias = corrected3pCorrelator/measured3pCorrelator;
1262 fDetectorBiasVsMHist->SetBinContent(b,bias);
1263 }
1264 } // end of for(Int_t b=1;b<=fNoOfMultipicityBins;b++)
1265
1266} // end of AliFlowAnalysisWithMixedHarmonics::CorrectForDetectorEffectsVsM()
1267
1268//================================================================================================================
1269
489d5531 1270void AliFlowAnalysisWithMixedHarmonics::ResetEventByEventQuantities()
1271{
1272 // Reset all event by event quantities.
1273
1274 fReQnk->Zero();
1275 fImQnk->Zero();
d66d46f7 1276 fSpk->Zero();
e558a733 1277
e999459d 1278 if(!fEvaluateDifferential3pCorrelator){return;}
5c078806 1279 for(Int_t sd=0;sd<2;sd++)
1280 {
1281 fRePEBE[sd]->Reset();
1282 fImPEBE[sd]->Reset();
1283 }
e999459d 1284 for(Int_t fs=0;fs<2;fs++)
1285 {
1286 for(Int_t sd=0;sd<2;sd++)
1287 {
1288 fOverlapEBE[fs][sd]->Reset();
1289 }
1290 }
5c078806 1291
d66d46f7 1292} // end of void AliFlowAnalysisWithMixedHarmonics::ResetEventByEventQuantities()
1293
1294//================================================================================================================
1295
489d5531 1296void AliFlowAnalysisWithMixedHarmonics::Calculate3pCorrelator()
1297{
e999459d 1298 // Calculate 3-p azimuthal correlator cos[n(phi1+phi2-2phi3)] in terms of Q_{n,k} and S_{p,k}.
d66d46f7 1299
1300 // a) Calculate 3-p correlator without using particle weights;
1301 // b) Calculate 3-p correlator with using particle weights.
1302
489d5531 1303 // a) Calculate 3-p correlator without using particle weights:
d66d46f7 1304 if(!(fUsePhiWeights || fUsePtWeights || fUseEtaWeights))
1305 {
489d5531 1306 // Multiplicity (number of RPs):
1307 Double_t dMult = (*fSpk)(0,0);
1308 // Real and imaginary parts of non-weighted Q-vectors (Q_{n,0}) evaluated in harmonics n and 2n:
1309 Double_t dReQ1n = (*fReQnk)(0,0);
1310 Double_t dReQ2n = (*fReQnk)(1,0);
1311 Double_t dImQ1n = (*fImQnk)(0,0);
d66d46f7 1312 Double_t dImQ2n = (*fImQnk)(1,0);
e999459d 1313 // 3-particle azimuthal correlator <cos(n*(phi1+phi2-2phi3))>:
1314 Double_t three1n1n2n = (pow(dReQ1n,2.)*dReQ2n + 2.*dReQ1n*dImQ1n*dImQ2n - pow(dImQ1n,2.)*dReQ2n
489d5531 1315 - 2.*(pow(dReQ1n,2.)+pow(dImQ1n,2.))
1316 - (pow(dReQ2n,2.)+pow(dImQ2n,2.))+2.*dMult)
d66d46f7 1317 / (dMult*(dMult-1.)*(dMult-2.));
e999459d 1318
1319
1320 /*
1321 Double_t reQ2nQ1nstarQ1nstar = pow(dReQ1n,2.)*dReQ2n + 2.*dReQ1n*dImQ1n*dImQ2n - pow(dImQ1n,2.)*dReQ2n;
1322 Double_t three2n1n1nQC = (reQ2nQ1nstarQ1nstar-2.*(pow(dReQ1n,2.)+pow(dImQ1n,2.))
1323 - (pow(dReQ2n,2.)+pow(dImQ2n,2.))+2.*dMult)
1324 / (dMult*(dMult-1.)*(dMult-2.));
1325 */
1326
e558a733 1327 // Fill all-events profile:
e999459d 1328 f3pCorrelatorPro->Fill(0.5,three1n1n2n,dMult*(dMult-1.)*(dMult-2.));
5c078806 1329
e999459d 1330 // 3-particle azimuthal correlator <cos(n*(phi1+phi2-2phi3))> vs multiplicity:
5c078806 1331 if(dMult<fMinMultiplicity)
1332 {
e999459d 1333 f3pCorrelatorVsMPro->Fill(0.5,three1n1n2n,dMult*(dMult-1.)*(dMult-2.));
5c078806 1334 } else if(dMult>=fMinMultiplicity+fNoOfMultipicityBins*fMultipicityBinWidth)
1335 {
e999459d 1336 f3pCorrelatorVsMPro->Fill(0.5+fNoOfMultipicityBins+1,three1n1n2n,dMult*(dMult-1.)*(dMult-2.));
5c078806 1337 } else
1338 {
e999459d 1339 f3pCorrelatorVsMPro->Fill(1.5+(Int_t)((dMult-fMinMultiplicity)/fMultipicityBinWidth),three1n1n2n,dMult*(dMult-1.)*(dMult-2.));
5c078806 1340 }
1341
d66d46f7 1342 } // end of if(!(fUsePhiWeights || fUsePtWeights || fUseEtaWeights))
1343
e999459d 1344 // b) Calculate 3-p correlator with using particle weights:
d66d46f7 1345 if(fUsePhiWeights || fUsePtWeights || fUseEtaWeights)
1346 {
e999459d 1347 // ...
d66d46f7 1348 } // end of if(fUsePhiWeights || fUsePtWeights || fUseEtaWeights)
1349
1350} // end of void AliFlowAnalysisWithMixedHarmonics::Calculate3pCorrelator()
1351
1352//================================================================================================================
1353
1354void AliFlowAnalysisWithMixedHarmonics::CalculateNonIsotropicTerms()
1355{
e999459d 1356 // Calculate non-isotropic terms which appear in the decomposition of 3-p correlator <cos[n(phi1+phi2-2phi3)]>.
d66d46f7 1357
5c078806 1358 // For detector with uniform acceptance all these terms vanish. These non-isotropic terms are stored in fNonIsotropicTermsPro.
1359 // Binning of fNonIsotropicTermsPro is organized as follows:
489d5531 1360 // 1st bin: <<cos(n*phi1)>>
d66d46f7 1361 // 2nd bin: <<sin(n*phi1)>>
489d5531 1362 // 3rd bin: <<cos(2n*phi1)>>
d66d46f7 1363 // 4th bin: <<sin(2n*phi1)>>
489d5531 1364 // 5th bin: <<cos(n*(phi1+phi2)>>
d66d46f7 1365 // 6th bin: <<sin(n*(phi1+phi2)>>
489d5531 1366 // 7th bin: <<cos(n*(2phi1-phi2)>>
d66d46f7 1367 // 8th bin: <<sin(n*(2phi1-phi2)>>
e999459d 1368 // 9th bin: <<cos(n*(phi1-phi2-phi3)>> // not needed
1369 // 10th bin: <<sin(n*(phi1-phi2-phi3)>> // not needed
d66d46f7 1370
5c078806 1371 // a) Calculate using particle weights;
1372 // b) Calculate without using particle weights.
1373
d66d46f7 1374 // a) Calculate using particle weights:
1375 if(!(fUsePhiWeights || fUsePtWeights || fUseEtaWeights))
1376 {
489d5531 1377 // Multiplicity (number of RPs):
1378 Double_t dMult = (*fSpk)(0,0);
1379 // Real and imaginary parts of non-weighted Q-vectors (Q_{n,0}) evaluated in harmonics n and 2n:
1380 Double_t dReQ1n = (*fReQnk)(0,0);
1381 Double_t dReQ2n = (*fReQnk)(1,0);
1382 Double_t dImQ1n = (*fImQnk)(0,0);
d66d46f7 1383 Double_t dImQ2n = (*fImQnk)(1,0);
489d5531 1384 // 1-particle terms:
d66d46f7 1385 Double_t cosP1n = 0.; // <cos(n*(phi1))>
1386 Double_t sinP1n = 0.; // <sin(n*(phi1))>
1387 Double_t cosP2n = 0.; // <cos(2n*(phi1))>
489d5531 1388 Double_t sinP2n = 0.; // <sin(2n*(phi1))>
1389 if(dMult>0)
1390 {
d66d46f7 1391 cosP1n = dReQ1n/dMult;
1392 sinP1n = dImQ1n/dMult;
1393 cosP2n = dReQ2n/dMult;
1394 sinP2n = dImQ2n/dMult;
e999459d 1395 // All-event avarages:
d66d46f7 1396 fNonIsotropicTermsPro->Fill(0.5,cosP1n,dMult); // <<cos(n*(phi1))>>
1397 fNonIsotropicTermsPro->Fill(1.5,sinP1n,dMult); // <<sin(n*(phi1))>>
1398 fNonIsotropicTermsPro->Fill(2.5,cosP2n,dMult); // <<cos(2n*(phi1))>>
489d5531 1399 fNonIsotropicTermsPro->Fill(3.5,sinP2n,dMult); // <<sin(2n*(phi1))>>
e999459d 1400 // All-event avarages vs M:
e558a733 1401 if(dMult<fMinMultiplicity)
1402 {
1403 fNonIsotropicTermsVsMPro->Fill(0.5,0.5,cosP1n,dMult); // <<cos(n*(phi1))>>
1404 fNonIsotropicTermsVsMPro->Fill(1.5,0.5,sinP1n,dMult); // <<sin(n*(phi1))>>
1405 fNonIsotropicTermsVsMPro->Fill(2.5,0.5,cosP2n,dMult); // <<cos(2n*(phi1))>>
1406 fNonIsotropicTermsVsMPro->Fill(3.5,0.5,sinP2n,dMult); // <<sin(2n*(phi1))>>
1407 } else if(dMult>=fMinMultiplicity+fNoOfMultipicityBins*fMultipicityBinWidth)
1408 {
1409 fNonIsotropicTermsVsMPro->Fill(0.5,0.5+fNoOfMultipicityBins+1,cosP1n,dMult); // <<cos(n*(phi1))>>
1410 fNonIsotropicTermsVsMPro->Fill(1.5,0.5+fNoOfMultipicityBins+1,sinP1n,dMult); // <<sin(n*(phi1))>>
1411 fNonIsotropicTermsVsMPro->Fill(2.5,0.5+fNoOfMultipicityBins+1,cosP2n,dMult); // <<cos(2n*(phi1))>>
1412 fNonIsotropicTermsVsMPro->Fill(3.5,0.5+fNoOfMultipicityBins+1,sinP2n,dMult); // <<sin(2n*(phi1))>>
1413 } else
1414 {
1415 fNonIsotropicTermsVsMPro->Fill(0.5,1.5+(Int_t)((dMult-fMinMultiplicity)/fMultipicityBinWidth),cosP1n,dMult); // <<cos(n*(phi1))>>
1416 fNonIsotropicTermsVsMPro->Fill(1.5,1.5+(Int_t)((dMult-fMinMultiplicity)/fMultipicityBinWidth),sinP1n,dMult); // <<sin(n*(phi1))>>
1417 fNonIsotropicTermsVsMPro->Fill(2.5,1.5+(Int_t)((dMult-fMinMultiplicity)/fMultipicityBinWidth),cosP2n,dMult); // <<cos(2n*(phi1))>>
1418 fNonIsotropicTermsVsMPro->Fill(3.5,1.5+(Int_t)((dMult-fMinMultiplicity)/fMultipicityBinWidth),sinP2n,dMult); // <<sin(2n*(phi1))>>
1419 }
1420 } // end of if(dMult>0)
489d5531 1421 // 2-particle terms:
d66d46f7 1422 Double_t cosP1nP1n = 0.; // <cos(n*(phi1+phi2))>
1423 Double_t sinP1nP1n = 0.; // <sin(n*(phi1+phi2))>
1424 Double_t cosP2nM1n = 0.; // <cos(n*(2phi1-phi2))>
489d5531 1425 Double_t sinP2nM1n = 0.; // <sin(n*(2phi1-phi2))>
1426 if(dMult>1)
1427 {
d66d46f7 1428 cosP1nP1n = (pow(dReQ1n,2)-pow(dImQ1n,2)-dReQ2n)/(dMult*(dMult-1));
489d5531 1429 sinP1nP1n = (2.*dReQ1n*dImQ1n-dImQ2n)/(dMult*(dMult-1));
d66d46f7 1430 cosP2nM1n = (dReQ2n*dReQ1n+dImQ2n*dImQ1n-dReQ1n)/(dMult*(dMult-1));
489d5531 1431 sinP2nM1n = (dImQ2n*dReQ1n-dReQ2n*dImQ1n-dImQ1n)/(dMult*(dMult-1));
e999459d 1432 // All-event avarages:
d66d46f7 1433 fNonIsotropicTermsPro->Fill(4.5,cosP1nP1n,dMult*(dMult-1.)); // <<cos(n*(phi1+phi2))>>
489d5531 1434 fNonIsotropicTermsPro->Fill(5.5,sinP1nP1n,dMult*(dMult-1.)); // <<sin(n*(phi1+phi2))>>
d66d46f7 1435 fNonIsotropicTermsPro->Fill(6.5,cosP2nM1n,dMult*(dMult-1.)); // <<cos(n*(2phi1-phi2))>>
489d5531 1436 fNonIsotropicTermsPro->Fill(7.5,sinP2nM1n,dMult*(dMult-1.)); // <<sin(n*(2phi1-phi2))>>
e999459d 1437 // All-event avarages vs M:
e558a733 1438 if(dMult<fMinMultiplicity)
1439 {
1440 fNonIsotropicTermsVsMPro->Fill(4.5,0.5,cosP1nP1n,dMult*(dMult-1.)); // <<cos(n*(phi1+phi2))>>
1441 fNonIsotropicTermsVsMPro->Fill(5.5,0.5,sinP1nP1n,dMult*(dMult-1.)); // <<sin(n*(phi1+phi2))>>
1442 fNonIsotropicTermsVsMPro->Fill(6.5,0.5,cosP2nM1n,dMult*(dMult-1.)); // <<cos(n*(2phi1-phi2))>>
1443 fNonIsotropicTermsVsMPro->Fill(7.5,0.5,sinP2nM1n,dMult*(dMult-1.)); // <<sin(n*(2phi1-phi2))>>
1444 } else if(dMult>=fMinMultiplicity+fNoOfMultipicityBins*fMultipicityBinWidth)
1445 {
1446 fNonIsotropicTermsVsMPro->Fill(4.5,0.5+fNoOfMultipicityBins+1,cosP1nP1n,dMult*(dMult-1.)); // <<cos(n*(phi1+phi2))>>
1447 fNonIsotropicTermsVsMPro->Fill(5.5,0.5+fNoOfMultipicityBins+1,sinP1nP1n,dMult*(dMult-1.)); // <<sin(n*(phi1+phi2))>>
1448 fNonIsotropicTermsVsMPro->Fill(6.5,0.5+fNoOfMultipicityBins+1,cosP2nM1n,dMult*(dMult-1.)); // <<cos(n*(2phi1-phi2))>>
1449 fNonIsotropicTermsVsMPro->Fill(7.5,0.5+fNoOfMultipicityBins+1,sinP2nM1n,dMult*(dMult-1.)); // <<sin(n*(2phi1-phi2))>>
1450 } else
1451 {
1452 fNonIsotropicTermsVsMPro->Fill(4.5,1.5+(Int_t)((dMult-fMinMultiplicity)/fMultipicityBinWidth),cosP1nP1n,dMult*(dMult-1.)); // <<cos(n*(phi1+phi2))>>
1453 fNonIsotropicTermsVsMPro->Fill(5.5,1.5+(Int_t)((dMult-fMinMultiplicity)/fMultipicityBinWidth),sinP1nP1n,dMult*(dMult-1.)); // <<sin(n*(phi1+phi2))>>
1454 fNonIsotropicTermsVsMPro->Fill(6.5,1.5+(Int_t)((dMult-fMinMultiplicity)/fMultipicityBinWidth),cosP2nM1n,dMult*(dMult-1.)); // <<cos(n*(2phi1-phi2))>>
1455 fNonIsotropicTermsVsMPro->Fill(7.5,1.5+(Int_t)((dMult-fMinMultiplicity)/fMultipicityBinWidth),sinP2nM1n,dMult*(dMult-1.)); // <<sin(n*(2phi1-phi2))>>
1456 }
1457 } // end of if(dMult>1)
e999459d 1458 // 3-particle: correct and ready but not needed, hence commented out.
1459 /*
489d5531 1460 Double_t cosP1nM1nM1n = 0.; // <cos(n*(phi1-phi2-phi3))>
1461 Double_t sinP1nM1nM1n = 0.; // <sin(n*(phi1-phi2-phi3))>
1462 if(dMult>2)
1463 {
1464 cosP1nM1nM1n = (dReQ1n*(pow(dReQ1n,2)+pow(dImQ1n,2))-dReQ1n*dReQ2n-dImQ1n*dImQ2n-2.*(dMult-1)*dReQ1n)
d66d46f7 1465 / (dMult*(dMult-1)*(dMult-2));
489d5531 1466 sinP1nM1nM1n = (-dImQ1n*(pow(dReQ1n,2)+pow(dImQ1n,2))+dReQ1n*dImQ2n-dImQ1n*dReQ2n+2.*(dMult-1)*dImQ1n)
1467 / (dMult*(dMult-1)*(dMult-2));
489d5531 1468 // All-events avarages:
d66d46f7 1469 fNonIsotropicTermsPro->Fill(8.5,cosP1nM1nM1n,dMult*(dMult-1.)*(dMult-2.)); // <<cos(n*(phi1-phi2-phi3))>>
e558a733 1470 fNonIsotropicTermsPro->Fill(9.5,sinP1nM1nM1n,dMult*(dMult-1.)*(dMult-2.)); // <<sin(n*(phi1-phi2-phi3))>>
1471 // All-events avarages vs M:
1472 if(dMult<fMinMultiplicity)
1473 {
1474 fNonIsotropicTermsVsMPro->Fill(8.5,0.5,cosP1nM1nM1n,dMult*(dMult-1.)*(dMult-2.)); // <<cos(n*(phi1-phi2-phi3))>>
1475 fNonIsotropicTermsVsMPro->Fill(9.5,0.5,sinP1nM1nM1n,dMult*(dMult-1.)*(dMult-2.)); // <<sin(n*(phi1-phi2-phi3))>>
1476 } else if(dMult>=fMinMultiplicity+fNoOfMultipicityBins*fMultipicityBinWidth)
1477 {
1478 fNonIsotropicTermsVsMPro->Fill(8.5,0.5+fNoOfMultipicityBins+1,cosP1nM1nM1n,dMult*(dMult-1.)*(dMult-2.)); // <<cos(n*(phi1-phi2-phi3))>>
1479 fNonIsotropicTermsVsMPro->Fill(9.5,0.5+fNoOfMultipicityBins+1,sinP1nM1nM1n,dMult*(dMult-1.)*(dMult-2.)); // <<sin(n*(phi1-phi2-phi3))>>
1480 } else
1481 {
1482 fNonIsotropicTermsVsMPro->Fill(8.5,1.5+(Int_t)((dMult-fMinMultiplicity)/fMultipicityBinWidth),
1483 cosP1nM1nM1n,dMult*(dMult-1.)*(dMult-2.)); // <<cos(n*(phi1-phi2-phi3))>>
1484 fNonIsotropicTermsVsMPro->Fill(9.5,1.5+(Int_t)((dMult-fMinMultiplicity)/fMultipicityBinWidth),
1485 sinP1nM1nM1n,dMult*(dMult-1.)*(dMult-2.)); // <<sin(n*(phi1-phi2-phi3))>>
1486 }
1487 } // end of if(dMult>2)
e999459d 1488 */
d66d46f7 1489 } // end of if(!(fUsePhiWeights || fUsePtWeights || fUseEtaWeights))
1490
e999459d 1491 // b) Calculate with using particle weights:
d66d46f7 1492 if(fUsePhiWeights || fUsePtWeights || fUseEtaWeights)
1493 {
e999459d 1494 // ...
d66d46f7 1495 } // end of if(fUsePhiWeights || fUsePtWeights || fUseEtaWeights)
1496
1497} // end of void AliFlowAnalysisWithMixedHarmonics::CalculateNonIsotropicTerms()
1498
5c078806 1499//================================================================================================================
1500
1501void AliFlowAnalysisWithMixedHarmonics::CalculateDifferential3pCorrelator()
1502{
e999459d 1503 // Calculate differential 3-p azimuthal correlator cos[n(psi1+psi2-2phi3)] in terms of Q_{2n}, p_{n}, q1_{n} and q2_{n}.
5c078806 1504
1505 // a) Calculate differential 3-p correlator without using particle weights;
1506 // b) Calculate differential 3-p correlator with using particle weights.
1507
1508 // a) Calculate differential 3-p correlator without using particle weights:
1509 if(!(fUsePhiWeights || fUsePtWeights || fUseEtaWeights))
1510 {
1511 // Multiplicity (number of RPs):
1512 Double_t dMult = (*fSpk)(0,0);
1513 // Real and imaginary parts of non-weighted Q-vectors (Q_{n,0}) evaluated in harmonic 2n:
1514 Double_t dReQ2n = (*fReQnk)(1,0);
1515 Double_t dImQ2n = (*fImQnk)(1,0);
1516 for(Int_t sd=0;sd<2;sd++) // [(p1+p2)/2,|p1-p2|]
1517 {
1518 // looping over all bins and calculating reduced correlations:
1519 for(Int_t b=1;b<=fnBinsPt;b++)
1520 {
1521 // real and imaginary parts of p_{n}:
1522 Double_t p1nRe = fRePEBE[sd]->GetBinContent(b)*fRePEBE[sd]->GetBinEntries(b);
1523 Double_t p1nIm = fImPEBE[sd]->GetBinContent(b)*fImPEBE[sd]->GetBinEntries(b);
e999459d 1524 // overlap 1: to be improved (terminology)
1525 Double_t overlap1 = fOverlapEBE[0][sd]->GetBinContent(b)*fOverlapEBE[0][sd]->GetBinEntries(b);
1526 // overlap 2: to be improved (terminology)
1527 Double_t overlap2 = fOverlapEBE[1][sd]->GetBinContent(b)*fOverlapEBE[1][sd]->GetBinEntries(b);
5c078806 1528 // number of pairs of POIs in particular (p1+p2)/2 or |p1-p2| bin:
1529 Double_t mp = fRePEBE[sd]->GetBinEntries(b);
e999459d 1530 // number of pairs of POI1/RP and POI2 in particular (p1+p2)/2 or |p1-p2| bin:
1531 Double_t mOverlap1 = fOverlapEBE[0][sd]->GetBinEntries(b);
1532 // number of pairs of POI2/RP and POI1 in particular (p1+p2)/2 or |p1-p2| bin:
1533 Double_t mOverlap2 = fOverlapEBE[1][sd]->GetBinEntries(b);
1534 // e-b-e weight for cos[n(psi1+psi2-2phi3)]:
1535 Double_t weight = mp*dMult-mOverlap1-mOverlap2;
1536
1537 Double_t cosP2nphi1M1npsi2M1npsi2 = 0; // cos[n(psi1+psi2-2phi3)]
1538 if(weight>0.)
5c078806 1539 {
e999459d 1540 cosP2nphi1M1npsi2M1npsi2 = (p1nRe*dReQ2n+p1nIm*dImQ2n-overlap1-overlap2)/(weight);
5c078806 1541 }
e999459d 1542 f3pCorrelatorVsPtSumDiffPro[sd]->Fill(fPtMin+(b-1)*fPtBinWidth,cosP2nphi1M1npsi2M1npsi2,weight);
5c078806 1543 } // end of for(Int_t b=1;b<=fnBinsPt;b++)
1544 } // end of for(Int_t sd=0;sd<2;sd++)
1545 } // end of if(!(fUsePhiWeights || fUsePtWeights || fUseEtaWeights))
1546
1547 // b) Calculate differential 3-p correlator by using particle weights:
1548 if(fUsePhiWeights || fUsePtWeights || fUseEtaWeights)
1549 {
e999459d 1550 // ...
5c078806 1551 } // end of if(fUsePhiWeights || fUsePtWeights || fUseEtaWeights)
1552
1553} // end of void AliFlowAnalysisWithMixedHarmonics::CalculateDifferential3pCorrelator()
1554
1555//================================================================================================================
1556