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 **************************************************************************/
18 /***************************************************************
19 * Only in this class nested loops are used for flow analysis. *
20 * Nested loops are used to evaluate: *
22 * a) Distribution of relative angle difference (phi1-phi2); *
23 * b) Cross-check the results for mixed harmonics. *
25 * Author: Ante Bilandzic (abilandzic@gmail.com) *
26 ***************************************************************/
28 #define AliFlowAnalysisWithNestedLoops_cxx
30 #include "Riostream.h"
31 #include "AliFlowCommonConstants.h"
32 #include "AliFlowCommonHist.h"
33 #include "AliFlowCommonHistResults.h"
40 #include "AliFlowEventSimple.h"
41 #include "AliFlowTrackSimple.h"
42 #include "AliFlowAnalysisWithNestedLoops.h"
46 ClassImp(AliFlowAnalysisWithNestedLoops)
48 //================================================================================================================
50 AliFlowAnalysisWithNestedLoops::AliFlowAnalysisWithNestedLoops():
54 fAnalysisSettings(NULL),
55 fPrintOnTheScreen(kTRUE),
70 fUsePhiWeights(kFALSE),
71 fUsePtWeights(kFALSE),
72 fUseEtaWeights(kFALSE),
73 fUseParticleWeights(NULL),
78 fEvaluateNestedLoopsForRAD(kTRUE),
79 fRelativeAngleDistribution(NULL),
81 fEvaluateNestedLoopsForMH(kFALSE),
82 fCorrelatorIntegerMH(1),
83 fCrossCheckInPtSumBinNo(4),
84 fCrossCheckInPtDiffBinNo(5)
88 // Base list to hold all output objects:
89 fHistList = new TList();
90 fHistListName = new TString("cobjNL");
91 fHistList->SetName(fHistListName->Data());
92 fHistList->SetOwner(kTRUE);
94 // List to hold histograms with phi, pt and eta weights:
95 fWeightsList = new TList();
97 // List to hold objects relevant for relative angle distributions:
98 fListRAD = new TList();
100 // List holding objects relevant for debugging and cross-checking of MH class:
101 fListMH = new TList();
103 // Initialize all arrays:
104 this->InitializeArraysForMH();
106 } // AliFlowAnalysisWithNestedLoops::AliFlowAnalysisWithNestedLoops()
108 //================================================================================================================
110 AliFlowAnalysisWithNestedLoops::~AliFlowAnalysisWithNestedLoops()
116 } // end of AliFlowAnalysisWithNestedLoops::~AliFlowAnalysisWithNestedLoops()
118 //================================================================================================================
120 void AliFlowAnalysisWithNestedLoops::Init()
122 // Initialize and book all objects.
124 // a) Cross check if the user settings make sense before starting;
125 // b) Access all common constants;
126 // c) Book and nest all lists in the base list fHistList;
127 // d) Book profile holding seetings for analysis with nested loops;
128 // e) Book common control histograms;
129 // f) Book all objects relevant for distributions;
130 // g) Book and fill histograms to hold phi, pt and eta weights.
132 //save old value and prevent histograms from being added to directory
133 //to avoid name clashes in case multiple analaysis objects are used
135 Bool_t oldHistAddStatus = TH1::AddDirectoryStatus();
136 TH1::AddDirectory(kFALSE);
138 TH1::SetDefaultSumw2();
140 this->CrossCheckSettings();
141 this->AccessConstants();
142 this->BookAndNestAllLists();
143 this->BookAndFillProfileHoldingSettings();
144 this->BookCommonHistograms();
145 this->BookEverythingForRAD();
146 this->BookEverythingForMH();
147 this->BookAndFillWeightsHistograms();
150 TH1::AddDirectory(oldHistAddStatus);
151 } // end of void AliFlowAnalysisWithNestedLoops::Init()
153 //================================================================================================================
155 void AliFlowAnalysisWithNestedLoops::Make(AliFlowEventSimple* anEvent)
157 // Running over data only in this method.
159 // a) Check all pointers used in this method;
160 // b) Fill common control histograms;
161 // c) Evaluate nested loops for relative angle distribution;
162 // d) Evaluate nested loops for mixed harmonics.
164 this->CheckPointersUsedInMake();
165 fCommonHists->FillControlHistograms(anEvent);
166 if(fEvaluateNestedLoopsForRAD) this->EvaluateNestedLoopsForRAD(anEvent);
167 if(fEvaluateNestedLoopsForMH) this->EvaluateNestedLoopsForMH(anEvent);
169 } // end of AliFlowAnalysisWithNestedLoops::Make(AliFlowEventSimple* anEvent)
171 //================================================================================================================
173 void AliFlowAnalysisWithNestedLoops::Finish()
175 // Calculate the final results.
177 // a) Check all pointers used in this method;
178 // b) Access settings for analysis with mixed harmonics;
179 // c) Print on the screen.
181 this->CheckPointersUsedInFinish();
182 this->AccessSettings();
183 if(fPrintOnTheScreen) this->PrintOnTheScreen();
185 } // end of AliFlowAnalysisWithNestedLoops::Finish()
187 //================================================================================================================
189 void AliFlowAnalysisWithNestedLoops::GetOutputHistograms(TList *outputListHistos)
191 // Get pointers to all objects saved in the output file.
193 // a) Get pointers for common control histograms.
196 this->SetHistList(outputListHistos);
200 cout<<" WARNING (NL): fHistList is NULL in AFAWNL::GOH() !!!!"<<endl;
204 this->GetPointersForBaseHistograms();
205 this->GetPointersForCommonHistograms();
206 Bool_t bEvaluateNestedLoopsForRAD = (Bool_t) fAnalysisSettings->GetBinContent(1); // to be improved (not needed here?)
207 Bool_t bEvaluateNestedLoopsForMH = (Bool_t) fAnalysisSettings->GetBinContent(2); // to be improved (not needed here?)
208 if(bEvaluateNestedLoopsForRAD) this->GetPointersForRAD();
209 if(bEvaluateNestedLoopsForMH) this->GetPointersForMH();
213 cout<<" WARNING (NL): outputListHistos is NULL in AFAWNL::GOH() !!!!"<<endl;
218 } // end of void AliFlowAnalysisWithNestedLoops::GetOutputHistograms(TList *outputListHistos)
220 //================================================================================================================
222 void AliFlowAnalysisWithNestedLoops::GetPointersForBaseHistograms()
224 // Get pointers to base histograms.
226 TString analysisSettingsName = "fAnalysisSettings";
227 TProfile *analysisSettings = dynamic_cast<TProfile*>(fHistList->FindObject(analysisSettingsName.Data()));
230 this->SetAnalysisSettings(analysisSettings);
234 cout<<" WARNING (NL): analysisSettings is NULL in AFAWNL::GPFBH() !!!!"<<endl;
239 } // end of void AliFlowAnalysisWithNestedLoops::GetPointersForBaseHistograms()
241 //================================================================================================================
243 void AliFlowAnalysisWithNestedLoops::GetPointersForCommonHistograms()
245 // Get pointers to common control histograms.
247 TString commonHistsName = "AliFlowCommonHistNL";
248 AliFlowCommonHist *commonHist = dynamic_cast<AliFlowCommonHist*>(fHistList->FindObject(commonHistsName.Data()));
251 this->SetCommonHists(commonHist);
255 cout<<" WARNING (NL): commonHist is NULL in AFAWNL::GPFCH() !!!!"<<endl;
260 } // end of void AliFlowAnalysisWithNestedLoops::GetPointersForCommonHistograms()
262 //================================================================================================================
264 void AliFlowAnalysisWithNestedLoops::GetPointersForRAD()
266 // Get pointers to objects relevant for relative angle distributions.
268 TList *listRAD = NULL;
269 listRAD = dynamic_cast<TList*>(fHistList->FindObject("Relative Angle Distribution"));
272 cout<<"WARNING: listRAD is NULL in AFAWNL::GPFRAD() !!!!"<<endl;
276 TString relativeAngleDistributionName = "fRelativeAngleDistribution";
277 TH1D *relativeAngleDistribution = dynamic_cast<TH1D*>(listRAD->FindObject(relativeAngleDistributionName.Data()));
278 if(relativeAngleDistribution)
280 this->SetRelativeAngleDistribution(relativeAngleDistribution);
283 } // end of void AliFlowAnalysisWithNestedLoops::GetPointersForRAD()
285 //================================================================================================================
287 void AliFlowAnalysisWithNestedLoops::GetPointersForMH()
289 // Get pointers to objects evaluated with nested loops.
291 TList *listMH = NULL;
292 listMH = dynamic_cast<TList*>(fHistList->FindObject("Mixed Harmonics"));
295 cout<<"WARNING: listMH is NULL in AFAWNL::GPFMH() !!!!"<<endl;
299 TString psdFlag[2] = {"PtSum","PtDiff"};
300 for(Int_t sd=0;sd<2;sd++)
302 TProfile *p3pCorrelatorVsPtSumDiffDirectPro = dynamic_cast<TProfile*>(listMH->FindObject(Form("f3pCorrelatorDirectVs%s",psdFlag[sd].Data())));
303 if(p3pCorrelatorVsPtSumDiffDirectPro)
305 this->Set3pCorrelatorVsPtSumDiffDirectPro(p3pCorrelatorVsPtSumDiffDirectPro,sd);
307 } // end of for(Int_t sd=0;sd<2;sd++)
309 } // end of void AliFlowAnalysisWithNestedLoops::GetPointersForMH()
311 //================================================================================================================
313 void AliFlowAnalysisWithNestedLoops::WriteHistograms(TString outputFileName)
315 // Store the final results in output .root file.
316 TFile *output = new TFile(outputFileName.Data(),"RECREATE");
317 fHistList->Write(fHistList->GetName(),TObject::kSingleKey);
321 //================================================================================================================
323 void AliFlowAnalysisWithNestedLoops::WriteHistograms(TDirectoryFile *outputFileName)
325 // Store the final results in output .root file.
326 fHistList->SetName("cobjNL");
327 fHistList->SetOwner(kTRUE);
328 outputFileName->Add(fHistList);
329 outputFileName->Write(outputFileName->GetName(),TObject::kSingleKey);
332 //================================================================================================================
334 void AliFlowAnalysisWithNestedLoops::BookAndNestAllLists()
336 // Book and nest all list in base list fHistList.
339 fWeightsList->SetName("Weights");
340 fWeightsList->SetOwner(kTRUE);
341 fHistList->Add(fWeightsList);
342 // List for Relative Angle Distribution:
343 fListRAD->SetName("Relative Angle Distribution");
344 fListRAD->SetOwner(kTRUE);
345 if(fEvaluateNestedLoopsForRAD) fHistList->Add(fListRAD);
346 // List for Mixed Harmonics:
347 fListMH->SetName("Mixed Harmonics");
348 fListMH->SetOwner(kTRUE);
349 if(fEvaluateNestedLoopsForMH) fHistList->Add(fListMH);
351 } // end of void AliFlowAnalysisWithNestedLoops::BookAndNestAllLists()
353 //================================================================================================================
355 void AliFlowAnalysisWithNestedLoops::BookAndFillProfileHoldingSettings()
357 // Book profile to hold all analysis settings.
359 TString analysisSettingsName = "fAnalysisSettings";
360 fAnalysisSettings = new TProfile(analysisSettingsName.Data(),"Settings for analysis with nested loops",6,0,6);
361 fAnalysisSettings->GetXaxis()->SetLabelSize(0.035);
362 fAnalysisSettings->GetXaxis()->SetBinLabel(1,"Nested loops for RAD?");
363 fAnalysisSettings->Fill(0.5,(Int_t)fEvaluateNestedLoopsForRAD);
364 fAnalysisSettings->GetXaxis()->SetBinLabel(2,"Nested loops for MH?");
365 fAnalysisSettings->Fill(1.5,(Int_t)fEvaluateNestedLoopsForMH);
366 fAnalysisSettings->GetXaxis()->SetBinLabel(3,"Integer n in cos(n(2#phi_{1}-#psi_{2}-#psi_{3}))");
367 fAnalysisSettings->Fill(2.5,(Int_t)fCorrelatorIntegerMH);
368 fAnalysisSettings->GetXaxis()->SetBinLabel(4,"Print on the screen?");
369 fAnalysisSettings->Fill(3.5,(Int_t)fPrintOnTheScreen);
370 fAnalysisSettings->GetXaxis()->SetBinLabel(5,"fCrossCheckInPtSumBinNo");
371 fAnalysisSettings->Fill(4.5,fCrossCheckInPtSumBinNo);
372 fAnalysisSettings->GetXaxis()->SetBinLabel(6,"fCrossCheckInPtDiffBinNo");
373 fAnalysisSettings->Fill(5.5,fCrossCheckInPtDiffBinNo);
374 fHistList->Add(fAnalysisSettings);
376 } // end of void AliFlowAnalysisWithNestedLoops::BookAndFillProfileHoldingSettings()
378 //================================================================================================================
380 void AliFlowAnalysisWithNestedLoops::BookCommonHistograms()
382 // Book common control histograms and common histograms for final results.
384 TString commonHistsName = "AliFlowCommonHistNL";
385 fCommonHists = new AliFlowCommonHist(commonHistsName.Data());
386 fHistList->Add(fCommonHists);
388 } // end of void AliFlowAnalysisWithNestedLoops::BookCommonHistograms()
390 //================================================================================================================
392 void AliFlowAnalysisWithNestedLoops::BookEverythingForRAD()
394 // Book all objects relevant calculation of relative angle distribution.
396 TString relativeAngleDistributionName = "fRelativeAngleDistribution";
397 fRelativeAngleDistribution = new TH1D(relativeAngleDistributionName.Data(),"Relative angle distribution",720,-TMath::TwoPi(),TMath::TwoPi());
398 fRelativeAngleDistribution->GetYaxis()->SetTitle("#frac{dN}{#Delta #phi}");
399 fRelativeAngleDistribution->GetXaxis()->SetTitle("#Delta #phi");
400 fListRAD->Add(fRelativeAngleDistribution);
402 } // end fo void AliFlowAnalysisWithNestedLoops::BookEverythingForRAD()
404 //================================================================================================================
406 void AliFlowAnalysisWithNestedLoops::AccessConstants()
408 // Access needed common constants from AliFlowCommonConstants.
410 fnBinsPhi = AliFlowCommonConstants::GetMaster()->GetNbinsPhi();
411 fPhiMin = AliFlowCommonConstants::GetMaster()->GetPhiMin();
412 fPhiMax = AliFlowCommonConstants::GetMaster()->GetPhiMax();
413 if(fnBinsPhi) fPhiBinWidth = (fPhiMax-fPhiMin)/fnBinsPhi;
414 fnBinsPt = AliFlowCommonConstants::GetMaster()->GetNbinsPt();
415 fPtMin = AliFlowCommonConstants::GetMaster()->GetPtMin();
416 fPtMax = AliFlowCommonConstants::GetMaster()->GetPtMax();
417 if(fnBinsPt) fPtBinWidth = (fPtMax-fPtMin)/fnBinsPt;
418 fnBinsEta = AliFlowCommonConstants::GetMaster()->GetNbinsEta();
419 fEtaMin = AliFlowCommonConstants::GetMaster()->GetEtaMin();
420 fEtaMax = AliFlowCommonConstants::GetMaster()->GetEtaMax();
421 if(fnBinsEta) fEtaBinWidth = (fEtaMax-fEtaMin)/fnBinsEta;
423 } // end of void AliFlowAnalysisWithNestedLoops::AccessConstants()
425 //================================================================================================================
427 void AliFlowAnalysisWithNestedLoops::CrossCheckSettings()
429 // Cross-check if the user settings make sense.
431 } // end of void AliFlowAnalysisWithNestedLoops::CrossCheckSettings()
433 //================================================================================================================
435 void AliFlowAnalysisWithNestedLoops::BookAndFillWeightsHistograms()
437 // Book and fill (by accessing file "weights.root") histograms which hold phi, pt and eta weights.
441 cout<<"WARNING: fWeightsList is NULL in AFAWNL::BAFWH() !!!!"<<endl;
444 // Profile to hold flags for weights:
445 TString fUseParticleWeightsName = "fUseParticleWeightsNL";
446 fUseParticleWeights = new TProfile(fUseParticleWeightsName.Data(),"0 = particle weight not used, 1 = particle weight used ",3,0,3);
447 fUseParticleWeights->SetLabelSize(0.06);
448 (fUseParticleWeights->GetXaxis())->SetBinLabel(1,"w_{#phi}");
449 (fUseParticleWeights->GetXaxis())->SetBinLabel(2,"w_{p_{T}}");
450 (fUseParticleWeights->GetXaxis())->SetBinLabel(3,"w_{#eta}");
451 fUseParticleWeights->Fill(0.5,(Int_t)fUsePhiWeights);
452 fUseParticleWeights->Fill(1.5,(Int_t)fUsePtWeights);
453 fUseParticleWeights->Fill(2.5,(Int_t)fUseEtaWeights);
454 fWeightsList->Add(fUseParticleWeights);
458 if(fWeightsList->FindObject("phi_weights"))
460 fPhiWeights = dynamic_cast<TH1F*>(fWeightsList->FindObject("phi_weights"));
461 if(TMath::Abs(fPhiWeights->GetBinWidth(1)-fPhiBinWidth)>pow(10.,-6.))
464 cout<<"WARNING (NL): Inconsistent binning in histograms for phi-weights throughout the code."<<endl;
470 cout<<"WARNING (NL): fWeightsList->FindObject(\"phi_weights\") is NULL in AFAWNL::BAFWH() !!!!"<<endl;
473 } // end of if(fUsePhiWeights)
477 if(fWeightsList->FindObject("pt_weights"))
479 fPtWeights = dynamic_cast<TH1D*>(fWeightsList->FindObject("pt_weights"));
480 if(TMath::Abs(fPtWeights->GetBinWidth(1)-fPtBinWidth)>pow(10.,-6.))
483 cout<<"WARNING (NL): Inconsistent binning in histograms for pt-weights throughout the code."<<endl;
489 cout<<"WARNING (NL): fWeightsList->FindObject(\"pt_weights\") is NULL in AFAWNL::BAFWH() !!!!"<<endl;
492 } // end of if(fUsePtWeights)
496 if(fWeightsList->FindObject("eta_weights"))
498 fEtaWeights = dynamic_cast<TH1D*>(fWeightsList->FindObject("eta_weights"));
499 if(TMath::Abs(fEtaWeights->GetBinWidth(1)-fEtaBinWidth)>pow(10.,-6.))
502 cout<<"WARNING (NL): Inconsistent binning in histograms for eta-weights throughout the code."<<endl;
508 cout<<"WARNING: fUseEtaWeights && fWeightsList->FindObject(\"eta_weights\") is NULL in AFAWNL::BAFWH() !!!!"<<endl;
511 } // end of if(fUseEtaWeights)
513 } // end of AliFlowAnalysisWithNestedLoops::BookAndFillWeightsHistograms()
515 //================================================================================================================
517 void AliFlowAnalysisWithNestedLoops::CheckPointersUsedInMake()
519 // Check pointers used in method Make().
521 if(fEvaluateNestedLoopsForRAD) CheckPointersForRAD("Make");
522 if(fEvaluateNestedLoopsForMH) CheckPointersForMH("Make");
524 } // end of AliFlowAnalysisWithNestedLoops::CheckPointersUsedInMake()
526 //================================================================================================================
528 void AliFlowAnalysisWithNestedLoops::CheckPointersUsedInFinish()
530 // Check pointers used in method Finish().
532 if(fEvaluateNestedLoopsForRAD) CheckPointersForRAD("Finish");
533 if(fEvaluateNestedLoopsForMH) CheckPointersForMH("Finish");
535 } // end of AliFlowAnalysisWithNestedLoops::CheckPointersUsedInFinish()
537 //================================================================================================================
539 void AliFlowAnalysisWithNestedLoops::CheckPointersForRAD(TString where)
541 // Check pointers relevant for calculation of relative angle distribution.
543 if(!fRelativeAngleDistribution)
546 cout<<" WARNING (NL): !fRelativeAngleDistribution is NULL in "<<where.Data()<<"() !!!!"<<endl;
551 if(strcmp(where.Data(),"Make") == 0)
553 // Check pointers used only in method Make():
556 else if(strcmp(where.Data(),"Finish") == 0)
558 // Check pointers used only in method Finish():
562 } // end of void AliFlowAnalysisWithNestedLoops::CheckPointersForRAD(TString where)
564 //================================================================================================================
566 void AliFlowAnalysisWithNestedLoops::CheckPointersForMH(TString where)
568 // Check pointers relevant for calculation of mixed harmonics.
570 for(Int_t sd=0;sd<2;sd++)
572 if(!(f3pCorrelatorVsPtSumDiffDirectPro[sd]))
575 cout<<" WARNING (NL): !"<<Form("f3pCorrelatorVsPtSumDiffDirectPro[%d]",sd)<<" is NULL in "<<where.Data()<<"() !!!!"<<endl;
581 if(strcmp(where.Data(),"Make") == 0)
583 // Check pointers used only in method Make():
586 else if(strcmp(where.Data(),"Finish") == 0)
588 // Check pointers used only in method Finish():
592 } // end of void AliFlowAnalysisWithNestedLoops::CheckPointersForMH(TString where)
594 //================================================================================================================
596 void AliFlowAnalysisWithNestedLoops::AccessSettings()
598 // Access the settings for analysis.
600 fEvaluateNestedLoopsForRAD = (Bool_t)fAnalysisSettings->GetBinContent(1);
601 fEvaluateNestedLoopsForMH = (Bool_t)fAnalysisSettings->GetBinContent(2);
602 fCorrelatorIntegerMH = (Int_t)fAnalysisSettings->GetBinContent(3);
603 fPrintOnTheScreen = (Bool_t)fAnalysisSettings->GetBinContent(4);
604 fCrossCheckInPtSumBinNo = (Int_t)fAnalysisSettings->GetBinContent(5);
605 fCrossCheckInPtDiffBinNo = (Int_t)fAnalysisSettings->GetBinContent(6);
607 } // end of AliFlowAnalysisWithNestedLoops::AccessSettings()
609 //================================================================================================================
611 void AliFlowAnalysisWithNestedLoops::InitializeArraysForMH()
613 // Initialize arrays mixed harmonics calculations.
615 for(Int_t sd=0;sd<2;sd++) // sum or difference
617 f3pCorrelatorVsPtSumDiffDirectPro[sd] = NULL;
620 } // end of AliFlowAnalysisWithNestedLoops::InitializeArraysForMH()
622 //================================================================================================================
624 void AliFlowAnalysisWithNestedLoops::BookEverythingForMH()
626 // Book all objects relevant for mixed harmonics.
628 if(fEvaluateNestedLoopsForMH)
630 TString psdFlag[2] = {"PtSum","PtDiff"};
631 TString psdTitleFlag[2] = {"(p_{t,1}+p_{t,2})/2","#left|p_{t,1}-p_{t,2}#right|"};
632 //TString s3pCorrelatorVsPtSumDiffDirectProName = "f3pCorrelatorVsPtSumDiffDirectPro";
633 for(Int_t sd=0;sd<2;sd++)
635 // to be improved: hardwired ,fnBinsPt,0.,fPtMax):
636 f3pCorrelatorVsPtSumDiffDirectPro[sd] = new TProfile(Form("f3pCorrelatorDirectVs%s",psdFlag[sd].Data()),"",fnBinsPt,0.,fPtMax);
637 //f3pCorrelatorVsPtSumDiffDirectPro[sd]->SetLabelSize(0.05);
638 //f3pCorrelatorVsPtSumDiffDirectPro[sd]->SetMarkerStyle(25);
639 f3pCorrelatorVsPtSumDiffDirectPro[sd]->GetXaxis()->SetTitle(psdTitleFlag[sd].Data());
640 fListMH->Add(f3pCorrelatorVsPtSumDiffDirectPro[sd]);
642 } // end of if(fEvaluateNestedLoopsForMH)
644 } // end of AliFlowAnalysisWithNestedLoops::BookEverythingForMH()
646 //================================================================================================================
648 void AliFlowAnalysisWithNestedLoops::EvaluateNestedLoopsForRAD(AliFlowEventSimple *anEvent)
650 // Evaluate nested loops needed for calculation of relative angle distribution.
652 Double_t dPhi1=0., dPhi2=0.; // azimuthal angles in the laboratory frame
653 AliFlowTrackSimple *aftsTrack = NULL; // simple track
655 // Loop over data and store for each distinct pair phi1-phi2 in fRelativeAngleDistribution:
656 Int_t nPrim = anEvent->NumberOfTracks(); // nPrim = total number of primary tracks, i.e. nPrim = nRP + nPOI + rest, where:
657 // nRP = # of particles used to determine the reaction plane ("Reference Particles");
658 // nPOI = # of particles of interest for a detailed flow analysis ("Particles of Interest");
659 // rest = # of particles which are not niether RPs nor POIs.
660 // Start nested loops over data:
661 for(Int_t i=0;i<nPrim;i++)
663 aftsTrack=anEvent->GetTrack(i);
666 if(!aftsTrack->InRPSelection()) continue; // consider only tracks which are RPs
667 dPhi1 = aftsTrack->Phi();
668 for(Int_t j=0;j<nPrim;j++)
670 if(j==i) continue; // eliminating trivial contribution from autocorrelation
671 aftsTrack=anEvent->GetTrack(j);
674 if(!aftsTrack->InRPSelection()) continue; // consider only tracks which are RPs
675 dPhi2 = aftsTrack->Phi();
676 // Fill the histogram:
677 fRelativeAngleDistribution->Fill(dPhi1-dPhi2);
679 } // end of for(Int_t j=0;j<nPrim;j++)
680 } else // to if(aftsTrack)
683 cout<<" WARNING (NL): No particle! (i.e. aftsTrack is a NULL pointer in AFAWNL::Make().)"<<endl;
686 } // end of for(Int_t i=0;i<nPrim;i++)
688 } // end of void AliFlowAnalysisWithNestedLoops::EvaluateNestedLoopsForRAD(AliFlowEventSimple *anEvent)
690 //================================================================================================================
692 void AliFlowAnalysisWithNestedLoops::EvaluateNestedLoopsForMH(AliFlowEventSimple *anEvent)
694 // Evaluate nested loops needed for calculation of mixed harmonics.
695 // Remark: phi label azimuthal angle of RP particle and psi label azimuthal angle of POI particle.
697 Int_t nPrim = anEvent->NumberOfTracks();
698 AliFlowTrackSimple *aftsTrack = NULL;
699 Double_t phi1=0.,psi2=0.,psi3=0.; // angles in the correlator cos[n(2phi1-psi2-psi3)]
700 Double_t pt2=0.,pt3=0.; // transverse momenta of psi2 and psi3
701 Int_t n = fCorrelatorIntegerMH;
702 // Evaluting differential correlator cos[n(2phi1-psi2-psi3)] with three nested loops:
703 for(Int_t i1=0;i1<nPrim;i1++)
705 aftsTrack=anEvent->GetTrack(i1);
706 // RP condition (first particle in the correlator must be RP):
707 if(!(aftsTrack->InRPSelection())) continue;
708 phi1 = aftsTrack->Phi();
709 for(Int_t i2=0;i2<nPrim;i2++)
712 aftsTrack = anEvent->GetTrack(i2);
713 // POI condition (second particle in the correlator must be POI):
714 if(!(aftsTrack->InPOISelection())) continue;
715 psi2 = aftsTrack->Phi();
716 pt2 = aftsTrack->Pt();
717 for(Int_t i3=0;i3<nPrim;i3++)
719 if(i3==i1||i3==i2) continue;
720 aftsTrack=anEvent->GetTrack(i3);
721 // POI condition (third particle in the correlator must be POI):
722 if(!(aftsTrack->InPOISelection())) continue;
723 psi3 = aftsTrack->Phi();
724 pt3 = aftsTrack->Pt();
725 // Evaluate and store differential correlator cos[n(2phi1-psi2-psi3)]:
726 Double_t ptSum = (pt2+pt3)/2.;
727 Double_t ptDiff = TMath::Abs(pt2-pt3);
728 Double_t diff3pCorrelator = TMath::Cos(n*(2.*phi1-psi2-psi3));
729 f3pCorrelatorVsPtSumDiffDirectPro[0]->Fill(ptSum,diff3pCorrelator,1.);
730 f3pCorrelatorVsPtSumDiffDirectPro[1]->Fill(ptDiff,diff3pCorrelator,1.);
731 } // end of for(Int_t i3=0;i3<nPrim;i3++)
732 } // end of for(Int_t i2=0;i2<nPrim;i2++)
733 } // end of for(Int_t i1=0;i1<nPrim;i1++)
735 } // end of void AliFlowAnalysisWithNestedLoops::EvaluateNestedLoopsForMH(AliFlowEventSimple *anEvent)
737 //================================================================================================================
739 void AliFlowAnalysisWithNestedLoops::PrintOnTheScreen()
741 // Print on the screen.
744 cout<<"****************************************************"<<endl;
745 cout<<"****************************************************"<<endl;
746 cout<<" Nested Loops "<<endl;
749 if(fEvaluateNestedLoopsForRAD)
751 cout<<" Evaluated for relative angle distribution."<<endl;
754 if(fEvaluateNestedLoopsForMH)
756 cout<<" Evaluated for mixed harmonics."<<endl;
757 if(fCorrelatorIntegerMH!=1)
759 cout<< " cos["<<fCorrelatorIntegerMH<<"(2phi1-psi2-psi3)] = "<<endl;
762 cout<< " cos(2phi1-psi2-psi3) = "<<endl;
764 cout<< " a) in pt sum bin "<<fCrossCheckInPtSumBinNo<<": "<<
765 f3pCorrelatorVsPtSumDiffDirectPro[0]->GetBinContent(fCrossCheckInPtSumBinNo)<<
766 " +/- "<<f3pCorrelatorVsPtSumDiffDirectPro[0]->GetBinError(fCrossCheckInPtSumBinNo)<<endl;
767 cout<< " b) in pt diff bin "<<fCrossCheckInPtDiffBinNo<<": "<<
768 f3pCorrelatorVsPtSumDiffDirectPro[1]->GetBinContent(fCrossCheckInPtDiffBinNo)<<
769 " +/- "<<f3pCorrelatorVsPtSumDiffDirectPro[1]->GetBinError(fCrossCheckInPtDiffBinNo)<<endl;
772 if(!fEvaluateNestedLoopsForRAD && !fEvaluateNestedLoopsForMH)
774 cout<<" Not evaluated."<<endl;
777 cout<<"****************************************************"<<endl;
778 cout<<"****************************************************"<<endl;
781 } // end of void AliFlowAnalysisWithNestedLoops::PrintOnTheScreen()