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():
55 fAnalysisSettings(NULL),
56 fOppositeChargesPOI(kFALSE),
57 fEvaluateDifferential3pCorrelator(kFALSE),
58 fPrintOnTheScreen(kTRUE),
73 fUsePhiWeights(kFALSE),
74 fUsePtWeights(kFALSE),
75 fUseEtaWeights(kFALSE),
76 fUseParticleWeights(NULL),
81 fEvaluateNestedLoopsForRAD(kTRUE),
82 fRelativeAngleDistribution(NULL),
85 fEvaluateNestedLoopsForQC(kFALSE),
87 fEvaluateNestedLoopsForMH(kFALSE),
88 f3pCorrelatorPro(NULL)
92 // Base list to hold all output objects:
93 fHistList = new TList();
94 fHistListName = new TString("cobjNL");
95 fHistList->SetName(fHistListName->Data());
96 fHistList->SetOwner(kTRUE);
98 // List to hold histograms with phi, pt and eta weights:
99 fWeightsList = new TList();
101 // List to hold objects relevant for relative angle distributions:
102 fListRAD = new TList();
104 // List holding objects relevant for debugging and cross-checking of Q-cumulants class:
105 fListQC = new TList();
107 // List holding objects relevant for debugging and cross-checking of MH class:
108 fListMH = new TList();
110 // Initialize all arrays:
111 this->InitializeArraysForMH();
113 } // AliFlowAnalysisWithNestedLoops::AliFlowAnalysisWithNestedLoops()
115 //================================================================================================================
117 AliFlowAnalysisWithNestedLoops::~AliFlowAnalysisWithNestedLoops()
123 } // end of AliFlowAnalysisWithNestedLoops::~AliFlowAnalysisWithNestedLoops()
125 //================================================================================================================
127 void AliFlowAnalysisWithNestedLoops::Init()
129 // Initialize and book all objects.
131 // a) Cross check if the user settings make sense before starting;
132 // b) Access all common constants;
133 // c) Book and nest all lists in the base list fHistList;
134 // d) Book profile holding seetings for analysis with nested loops;
135 // e) Book common control histograms;
136 // f) Book all objects relevant for distributions;
137 // g) Book and fill histograms to hold phi, pt and eta weights;
138 // h) Store harmonic n.
141 //save old value and prevent histograms from being added to directory
142 //to avoid name clashes in case multiple analaysis objects are used
144 Bool_t oldHistAddStatus = TH1::AddDirectoryStatus();
145 TH1::AddDirectory(kFALSE);
147 TH1::SetDefaultSumw2();
149 this->CrossCheckSettings();
150 this->AccessConstants();
151 this->BookAndNestAllLists();
152 this->BookAndFillProfileHoldingSettings();
153 this->BookCommonHistograms();
154 this->BookEverythingForRAD();
155 this->BookEverythingForMH();
156 this->BookAndFillWeightsHistograms();
157 this->StoreHarmonic();
160 TH1::AddDirectory(oldHistAddStatus);
161 } // end of void AliFlowAnalysisWithNestedLoops::Init()
163 //================================================================================================================
165 void AliFlowAnalysisWithNestedLoops::Make(AliFlowEventSimple* anEvent)
167 // Running over data only in this method.
169 // a) Check all pointers used in this method;
170 // b) Fill common control histograms;
171 // c) Evaluate nested loops for relative angle distribution;
172 // d) Evaluate nested loops for mixed harmonics.
174 this->CheckPointersUsedInMake();
175 fCommonHists->FillControlHistograms(anEvent);
176 if(fEvaluateNestedLoopsForRAD) this->EvaluateNestedLoopsForRAD(anEvent);
177 if(fEvaluateNestedLoopsForMH) this->EvaluateNestedLoopsForMH(anEvent);
179 } // end of AliFlowAnalysisWithNestedLoops::Make(AliFlowEventSimple* anEvent)
181 //================================================================================================================
183 void AliFlowAnalysisWithNestedLoops::Finish()
185 // Calculate the final results.
187 // a) Access settings for analysis with mixed harmonics;
188 // b) Check all pointers used in this method;
189 // c) Print on the screen.
191 this->AccessSettings();
192 this->CheckPointersUsedInFinish();
193 if(fPrintOnTheScreen) this->PrintOnTheScreen();
195 } // end of AliFlowAnalysisWithNestedLoops::Finish()
197 //================================================================================================================
199 void AliFlowAnalysisWithNestedLoops::GetOutputHistograms(TList *outputListHistos)
201 // Get pointers to all objects saved in the output file.
203 // a) Get pointers for common control histograms.
206 this->SetHistList(outputListHistos);
210 cout<<" WARNING (NL): fHistList is NULL in AFAWNL::GOH() !!!!"<<endl;
214 this->GetPointersForBaseHistograms();
215 this->GetPointersForCommonHistograms();
216 this->AccessSettings();
217 if(fEvaluateNestedLoopsForRAD) this->GetPointersForRAD();
218 if(fEvaluateNestedLoopsForMH) this->GetPointersForMH();
222 cout<<" WARNING (NL): outputListHistos is NULL in AFAWNL::GOH() !!!!"<<endl;
227 } // end of void AliFlowAnalysisWithNestedLoops::GetOutputHistograms(TList *outputListHistos)
229 //================================================================================================================
231 void AliFlowAnalysisWithNestedLoops::GetPointersForBaseHistograms()
233 // Get pointers to base histograms.
235 TString analysisSettingsName = "fAnalysisSettings";
236 TProfile *analysisSettings = dynamic_cast<TProfile*>(fHistList->FindObject(analysisSettingsName.Data()));
239 this->SetAnalysisSettings(analysisSettings);
243 cout<<" WARNING (NL): analysisSettings is NULL in AFAWNL::GPFBH() !!!!"<<endl;
248 } // end of void AliFlowAnalysisWithNestedLoops::GetPointersForBaseHistograms()
250 //================================================================================================================
252 void AliFlowAnalysisWithNestedLoops::GetPointersForCommonHistograms()
254 // Get pointers to common control histograms.
256 TString commonHistsName = "AliFlowCommonHistNL";
257 AliFlowCommonHist *commonHist = dynamic_cast<AliFlowCommonHist*>(fHistList->FindObject(commonHistsName.Data()));
260 this->SetCommonHists(commonHist);
264 cout<<" WARNING (NL): commonHist is NULL in AFAWNL::GPFCH() !!!!"<<endl;
269 } // end of void AliFlowAnalysisWithNestedLoops::GetPointersForCommonHistograms()
271 //================================================================================================================
273 void AliFlowAnalysisWithNestedLoops::GetPointersForRAD()
275 // Get pointers to objects relevant for relative angle distributions.
277 TList *listRAD = NULL;
278 listRAD = dynamic_cast<TList*>(fHistList->FindObject("Relative Angle Distribution"));
281 cout<<"WARNING: listRAD is NULL in AFAWNL::GPFRAD() !!!!"<<endl;
285 TString relativeAngleDistributionName = "fRelativeAngleDistribution";
286 TH1D *relativeAngleDistribution = dynamic_cast<TH1D*>(listRAD->FindObject(relativeAngleDistributionName.Data()));
287 if(relativeAngleDistribution)
289 this->SetRelativeAngleDistribution(relativeAngleDistribution);
292 TString chargeName = "fCharge";
293 TH1D *charge = dynamic_cast<TH1D*>(listRAD->FindObject(chargeName.Data()));
296 this->SetCharge(charge);
299 } // end of void AliFlowAnalysisWithNestedLoops::GetPointersForRAD()
301 //================================================================================================================
303 void AliFlowAnalysisWithNestedLoops::GetPointersForMH()
305 // Get pointers to objects evaluated with nested loops.
307 TList *listMH = NULL;
308 listMH = dynamic_cast<TList*>(fHistList->FindObject("Mixed Harmonics"));
311 cout<<"WARNING: listMH is NULL in AFAWNL::GPFMH() !!!!"<<endl;
315 TString s3pCorrelatorProName = "f3pCorrelatorPro";
316 TProfile *p3pCorrelatorPro = dynamic_cast<TProfile*>(listMH->FindObject(s3pCorrelatorProName.Data()));
319 this->Set3pCorrelatorPro(p3pCorrelatorPro);
321 if(!fEvaluateDifferential3pCorrelator){return;}
322 TString psdFlag[2] = {"PtSum","PtDiff"};
323 for(Int_t sd=0;sd<2;sd++)
325 TProfile *p3pCorrelatorVsPtSumDiffDirectPro = dynamic_cast<TProfile*>(listMH->FindObject(Form("f3pCorrelatorDirectVs%s",psdFlag[sd].Data())));
326 if(p3pCorrelatorVsPtSumDiffDirectPro)
328 this->Set3pCorrelatorVsPtSumDiffDirectPro(p3pCorrelatorVsPtSumDiffDirectPro,sd);
330 } // end of for(Int_t sd=0;sd<2;sd++)
332 } // end of void AliFlowAnalysisWithNestedLoops::GetPointersForMH()
334 //================================================================================================================
336 void AliFlowAnalysisWithNestedLoops::WriteHistograms(TString outputFileName)
338 // Store the final results in output .root file.
339 TFile *output = new TFile(outputFileName.Data(),"RECREATE");
340 fHistList->Write(fHistList->GetName(),TObject::kSingleKey);
344 //================================================================================================================
346 void AliFlowAnalysisWithNestedLoops::WriteHistograms(TDirectoryFile *outputFileName)
348 // Store the final results in output .root file.
349 fHistList->SetName("cobjNL");
350 fHistList->SetOwner(kTRUE);
351 outputFileName->Add(fHistList);
352 outputFileName->Write(outputFileName->GetName(),TObject::kSingleKey);
355 //================================================================================================================
357 void AliFlowAnalysisWithNestedLoops::StoreHarmonic()
359 // Store harmonic n used in cos[n*(phi1+phi2-2phi3)] and cos[n*(psi1+psi2-2phi3)].
361 (fCommonHists->GetHarmonic())->Fill(0.5,fHarmonic);
363 } // end of void AliFlowAnalysisWithNestedLoops::StoreHarmonic()
365 //================================================================================================================
367 void AliFlowAnalysisWithNestedLoops::BookAndNestAllLists()
369 // Book and nest all list in base list fHistList.
372 fWeightsList->SetName("Weights");
373 fWeightsList->SetOwner(kTRUE);
374 fHistList->Add(fWeightsList);
375 // List for Relative Angle Distribution:
376 fListRAD->SetName("Relative Angle Distribution");
377 fListRAD->SetOwner(kTRUE);
378 if(fEvaluateNestedLoopsForRAD) fHistList->Add(fListRAD);
379 // List for Q-cumulants:
380 fListQC->SetName("Q-cumulants");
381 fListQC->SetOwner(kTRUE);
382 if(fEvaluateNestedLoopsForQC) fHistList->Add(fListQC);
383 // List for Mixed Harmonics:
384 fListMH->SetName("Mixed Harmonics");
385 fListMH->SetOwner(kTRUE);
386 if(fEvaluateNestedLoopsForMH) fHistList->Add(fListMH);
388 } // end of void AliFlowAnalysisWithNestedLoops::BookAndNestAllLists()
390 //================================================================================================================
392 void AliFlowAnalysisWithNestedLoops::BookAndFillProfileHoldingSettings()
394 // Book profile to hold all analysis settings.
396 TString analysisSettingsName = "fAnalysisSettings";
397 fAnalysisSettings = new TProfile(analysisSettingsName.Data(),"Settings for analysis with nested loops",7,0,7);
398 fAnalysisSettings->GetXaxis()->SetLabelSize(0.035);
399 fAnalysisSettings->GetXaxis()->SetBinLabel(1,"Nested loops for RAD?");
400 fAnalysisSettings->Fill(0.5,(Int_t)fEvaluateNestedLoopsForRAD);
401 fAnalysisSettings->GetXaxis()->SetBinLabel(2,"Nested loops for QC?");
402 fAnalysisSettings->Fill(1.5,(Int_t)fEvaluateNestedLoopsForQC);
403 fAnalysisSettings->GetXaxis()->SetBinLabel(3,"Nested loops for MH?");
404 fAnalysisSettings->Fill(2.5,(Int_t)fEvaluateNestedLoopsForMH);
405 fAnalysisSettings->GetXaxis()->SetBinLabel(4,"fHarmonic");
406 fAnalysisSettings->Fill(3.5,(Int_t)fHarmonic);
407 fAnalysisSettings->GetXaxis()->SetBinLabel(5,"Print on the screen?");
408 fAnalysisSettings->Fill(4.5,(Int_t)fPrintOnTheScreen);
409 fAnalysisSettings->GetXaxis()->SetBinLabel(6,"fOppositeChargesPOI");
410 fAnalysisSettings->Fill(5.5,fOppositeChargesPOI);
411 fAnalysisSettings->GetXaxis()->SetBinLabel(7,"fEvaluateDifferential3pCorrelator");
412 fAnalysisSettings->Fill(6.5,fEvaluateDifferential3pCorrelator);
413 fHistList->Add(fAnalysisSettings);
415 } // end of void AliFlowAnalysisWithNestedLoops::BookAndFillProfileHoldingSettings()
417 //================================================================================================================
419 void AliFlowAnalysisWithNestedLoops::BookCommonHistograms()
421 // Book common control histograms and common histograms for final results.
423 TString commonHistsName = "AliFlowCommonHistNL";
424 fCommonHists = new AliFlowCommonHist(commonHistsName.Data());
425 fHistList->Add(fCommonHists);
427 } // end of void AliFlowAnalysisWithNestedLoops::BookCommonHistograms()
429 //================================================================================================================
431 void AliFlowAnalysisWithNestedLoops::BookEverythingForRAD()
433 // Book all objects relevant calculation of relative angle distribution.
435 TString relativeAngleDistributionName = "fRelativeAngleDistribution";
436 fRelativeAngleDistribution = new TH1D(relativeAngleDistributionName.Data(),"Relative angle distribution",720,-TMath::TwoPi(),TMath::TwoPi());
437 fRelativeAngleDistribution->GetYaxis()->SetTitle("#frac{dN}{#Delta #phi}");
438 fRelativeAngleDistribution->GetXaxis()->SetTitle("#Delta #phi");
439 fListRAD->Add(fRelativeAngleDistribution);
441 TString chargeName = "fCharge";
442 fCharge = new TH1D(chargeName.Data(),"Charges",3,0,3);
444 fCharge->GetXaxis()->SetBinLabel(1,"+");
445 fCharge->GetXaxis()->SetBinLabel(2,"-");
446 fCharge->GetXaxis()->SetBinLabel(3,"rest");
447 fListRAD->Add(fCharge);
449 } // end fo void AliFlowAnalysisWithNestedLoops::BookEverythingForRAD()
451 //================================================================================================================
453 void AliFlowAnalysisWithNestedLoops::AccessConstants()
455 // Access needed common constants from AliFlowCommonConstants.
457 fnBinsPhi = AliFlowCommonConstants::GetMaster()->GetNbinsPhi();
458 fPhiMin = AliFlowCommonConstants::GetMaster()->GetPhiMin();
459 fPhiMax = AliFlowCommonConstants::GetMaster()->GetPhiMax();
460 if(fnBinsPhi) fPhiBinWidth = (fPhiMax-fPhiMin)/fnBinsPhi;
461 fnBinsPt = AliFlowCommonConstants::GetMaster()->GetNbinsPt();
462 fPtMin = AliFlowCommonConstants::GetMaster()->GetPtMin();
463 fPtMax = AliFlowCommonConstants::GetMaster()->GetPtMax();
464 if(fnBinsPt) fPtBinWidth = (fPtMax-fPtMin)/fnBinsPt;
465 fnBinsEta = AliFlowCommonConstants::GetMaster()->GetNbinsEta();
466 fEtaMin = AliFlowCommonConstants::GetMaster()->GetEtaMin();
467 fEtaMax = AliFlowCommonConstants::GetMaster()->GetEtaMax();
468 if(fnBinsEta) fEtaBinWidth = (fEtaMax-fEtaMin)/fnBinsEta;
470 } // end of void AliFlowAnalysisWithNestedLoops::AccessConstants()
472 //================================================================================================================
474 void AliFlowAnalysisWithNestedLoops::CrossCheckSettings()
476 // Cross-check if the user settings make sense.
480 } // end of void AliFlowAnalysisWithNestedLoops::CrossCheckSettings()
482 //================================================================================================================
484 void AliFlowAnalysisWithNestedLoops::BookAndFillWeightsHistograms()
486 // Book and fill (by accessing file "weights.root") histograms which hold phi, pt and eta weights.
490 cout<<"WARNING: fWeightsList is NULL in AFAWNL::BAFWH() !!!!"<<endl;
493 // Profile to hold flags for weights:
494 TString fUseParticleWeightsName = "fUseParticleWeightsNL";
495 fUseParticleWeights = new TProfile(fUseParticleWeightsName.Data(),"0 = particle weight not used, 1 = particle weight used ",3,0,3);
496 fUseParticleWeights->SetLabelSize(0.06);
497 (fUseParticleWeights->GetXaxis())->SetBinLabel(1,"w_{#phi}");
498 (fUseParticleWeights->GetXaxis())->SetBinLabel(2,"w_{p_{T}}");
499 (fUseParticleWeights->GetXaxis())->SetBinLabel(3,"w_{#eta}");
500 fUseParticleWeights->Fill(0.5,(Int_t)fUsePhiWeights);
501 fUseParticleWeights->Fill(1.5,(Int_t)fUsePtWeights);
502 fUseParticleWeights->Fill(2.5,(Int_t)fUseEtaWeights);
503 fWeightsList->Add(fUseParticleWeights);
507 if(fWeightsList->FindObject("phi_weights"))
509 fPhiWeights = dynamic_cast<TH1F*>(fWeightsList->FindObject("phi_weights"));
510 if(TMath::Abs(fPhiWeights->GetBinWidth(1)-fPhiBinWidth)>pow(10.,-6.))
513 cout<<"WARNING (NL): Inconsistent binning in histograms for phi-weights throughout the code."<<endl;
519 cout<<"WARNING (NL): fWeightsList->FindObject(\"phi_weights\") is NULL in AFAWNL::BAFWH() !!!!"<<endl;
522 } // end of if(fUsePhiWeights)
526 if(fWeightsList->FindObject("pt_weights"))
528 fPtWeights = dynamic_cast<TH1D*>(fWeightsList->FindObject("pt_weights"));
529 if(TMath::Abs(fPtWeights->GetBinWidth(1)-fPtBinWidth)>pow(10.,-6.))
532 cout<<"WARNING (NL): Inconsistent binning in histograms for pt-weights throughout the code."<<endl;
538 cout<<"WARNING (NL): fWeightsList->FindObject(\"pt_weights\") is NULL in AFAWNL::BAFWH() !!!!"<<endl;
541 } // end of if(fUsePtWeights)
545 if(fWeightsList->FindObject("eta_weights"))
547 fEtaWeights = dynamic_cast<TH1D*>(fWeightsList->FindObject("eta_weights"));
548 if(TMath::Abs(fEtaWeights->GetBinWidth(1)-fEtaBinWidth)>pow(10.,-6.))
551 cout<<"WARNING (NL): Inconsistent binning in histograms for eta-weights throughout the code."<<endl;
557 cout<<"WARNING: fUseEtaWeights && fWeightsList->FindObject(\"eta_weights\") is NULL in AFAWNL::BAFWH() !!!!"<<endl;
560 } // end of if(fUseEtaWeights)
562 } // end of AliFlowAnalysisWithNestedLoops::BookAndFillWeightsHistograms()
564 //================================================================================================================
566 void AliFlowAnalysisWithNestedLoops::CheckPointersUsedInMake()
568 // Check pointers used in method Make().
570 if(fEvaluateNestedLoopsForRAD) CheckPointersForRAD("Make");
571 if(fEvaluateNestedLoopsForMH) CheckPointersForMH("Make");
573 } // end of AliFlowAnalysisWithNestedLoops::CheckPointersUsedInMake()
575 //================================================================================================================
577 void AliFlowAnalysisWithNestedLoops::CheckPointersUsedInFinish()
579 // Check pointers used in method Finish().
581 if(fEvaluateNestedLoopsForRAD) CheckPointersForRAD("Finish");
582 if(fEvaluateNestedLoopsForMH) CheckPointersForMH("Finish");
584 } // end of AliFlowAnalysisWithNestedLoops::CheckPointersUsedInFinish()
586 //================================================================================================================
588 void AliFlowAnalysisWithNestedLoops::CheckPointersForRAD(TString where)
590 // Check pointers relevant for calculation of relative angle distribution.
592 if(!fRelativeAngleDistribution)
595 cout<<" WARNING (NL): fRelativeAngleDistribution is NULL in "<<where.Data()<<"() !!!!"<<endl;
600 if(strcmp(where.Data(),"Make") == 0)
602 // Check pointers used only in method Make():
605 else if(strcmp(where.Data(),"Finish") == 0)
607 // Check pointers used only in method Finish():
611 } // end of void AliFlowAnalysisWithNestedLoops::CheckPointersForRAD(TString where)
613 //================================================================================================================
615 void AliFlowAnalysisWithNestedLoops::CheckPointersForMH(TString where)
617 // Check pointers relevant for calculation of mixed harmonics.
619 if(strcmp(where.Data(),"Make") == 0)
621 // Check pointers used only in method Make():
622 if(!f3pCorrelatorPro)
625 cout<<" WARNING (NL): f3pCorrelatorPro is NULL in Make() !!!!"<<endl;
629 if(!fEvaluateDifferential3pCorrelator){return;}
630 for(Int_t sd=0;sd<2;sd++)
632 if(!(f3pCorrelatorVsPtSumDiffDirectPro[sd]))
635 cout<<" WARNING (NL): "<<Form("f3pCorrelatorVsPtSumDiffDirectPro[%d]",sd)<<" is NULL in Make() !!!!"<<endl;
639 } // end of for(Int_t sd=0;sd<2;sd++)
640 } // if(strcmp(where.Data(),"Make") == 0)
641 else if(strcmp(where.Data(),"Finish") == 0)
643 // Check pointers used only in method Finish():
644 if(!f3pCorrelatorPro)
647 cout<<" WARNING (NL): f3pCorrelatorPro is NULL in Finish() !!!!"<<endl;
651 } // else if(strcmp(where.Data(),"Finish") == 0)
653 } // end of void AliFlowAnalysisWithNestedLoops::CheckPointersForMH(TString where)
655 //================================================================================================================
657 void AliFlowAnalysisWithNestedLoops::AccessSettings()
659 // Access the settings for analysis.
661 fEvaluateNestedLoopsForRAD = (Bool_t)fAnalysisSettings->GetBinContent(1);
662 fEvaluateNestedLoopsForQC = (Bool_t)fAnalysisSettings->GetBinContent(2);
663 fEvaluateNestedLoopsForMH = (Bool_t)fAnalysisSettings->GetBinContent(3);
664 fHarmonic = (Int_t)fAnalysisSettings->GetBinContent(4);
665 fPrintOnTheScreen = (Bool_t)fAnalysisSettings->GetBinContent(5);
666 fOppositeChargesPOI = (Bool_t)fAnalysisSettings->GetBinContent(6);
667 fEvaluateDifferential3pCorrelator = (Bool_t)fAnalysisSettings->GetBinContent(7);
669 } // end of AliFlowAnalysisWithNestedLoops::AccessSettings()
671 //================================================================================================================
673 void AliFlowAnalysisWithNestedLoops::InitializeArraysForMH()
675 // Initialize arrays mixed harmonics calculations.
677 for(Int_t sd=0;sd<2;sd++) // sum or difference
679 f3pCorrelatorVsPtSumDiffDirectPro[sd] = NULL;
682 } // end of AliFlowAnalysisWithNestedLoops::InitializeArraysForMH()
684 //================================================================================================================
686 void AliFlowAnalysisWithNestedLoops::BookEverythingForMH()
688 // Book all objects relevant for mixed harmonics.
690 if(fEvaluateNestedLoopsForMH)
692 TString s3pCorrelatorProName = "f3pCorrelatorPro";
693 f3pCorrelatorPro = new TProfile(s3pCorrelatorProName.Data(),"",1,0,1);
694 f3pCorrelatorPro->SetStats(kFALSE);
695 f3pCorrelatorPro->GetXaxis()->SetLabelOffset(0.01);
696 f3pCorrelatorPro->GetXaxis()->SetLabelSize(0.05);
699 f3pCorrelatorPro->GetXaxis()->SetBinLabel(1,"#LT#LTcos(#phi_{1}+#phi_{2}-2#phi_{3})#GT#GT");
702 f3pCorrelatorPro->GetXaxis()->SetBinLabel(1,Form("#LT#LTcos[%i(#phi_{1}+#phi_{2}-2#phi_{3})]#GT#GT",fHarmonic));
704 fListMH->Add(f3pCorrelatorPro);
706 if(fEvaluateDifferential3pCorrelator)
708 TString psdFlag[2] = {"PtSum","PtDiff"};
709 TString psdTitleFlag[2] = {"(p_{t,1}+p_{t,2})/2","#left|p_{t,1}-p_{t,2}#right|"};
710 //TString s3pCorrelatorVsPtSumDiffDirectProName = "f3pCorrelatorVsPtSumDiffDirectPro";
711 for(Int_t sd=0;sd<2;sd++)
713 // to be improved: hardwired ,fnBinsPt,0.,fPtMax):
714 f3pCorrelatorVsPtSumDiffDirectPro[sd] = new TProfile(Form("f3pCorrelatorDirectVs%s",psdFlag[sd].Data()),"",fnBinsPt,0.,fPtMax);
715 //f3pCorrelatorVsPtSumDiffDirectPro[sd]->SetLabelSize(0.05);
716 //f3pCorrelatorVsPtSumDiffDirectPro[sd]->SetMarkerStyle(25);
717 f3pCorrelatorVsPtSumDiffDirectPro[sd]->GetXaxis()->SetTitle(psdTitleFlag[sd].Data());
718 fListMH->Add(f3pCorrelatorVsPtSumDiffDirectPro[sd]);
720 } // end of if(fEvaluateDifferential3pCorrelator)
721 } // end of if(fEvaluateNestedLoopsForMH)
723 } // end of AliFlowAnalysisWithNestedLoops::BookEverythingForMH()
725 //================================================================================================================
727 void AliFlowAnalysisWithNestedLoops::EvaluateNestedLoopsForRAD(AliFlowEventSimple *anEvent)
729 // Evaluate nested loops needed for calculation of relative angle distribution.
731 Double_t dPhi1=0., dPhi2=0.; // azimuthal angles in the laboratory frame
732 AliFlowTrackSimple *aftsTrack = NULL; // simple track
734 // Loop over data and store for each distinct pair phi1-phi2 in fRelativeAngleDistribution:
735 Int_t nPrim = anEvent->NumberOfTracks(); // nPrim = total number of primary tracks, i.e. nPrim = nRP + nPOI + rest, where:
736 // nRP = # of particles used to determine the reaction plane ("Reference Particles");
737 // nPOI = # of particles of interest for a detailed flow analysis ("Particles of Interest");
738 // rest = # of particles which are not niether RPs nor POIs.
739 // Start nested loops over data:
740 for(Int_t i=0;i<nPrim;i++)
742 aftsTrack=anEvent->GetTrack(i);
745 if(!aftsTrack->InRPSelection()) continue; // consider only tracks which are RPs
746 dPhi1 = aftsTrack->Phi();
747 for(Int_t j=0;j<nPrim;j++)
749 if(j==i) continue; // eliminating trivial contribution from autocorrelation
750 aftsTrack=anEvent->GetTrack(j);
753 if(!aftsTrack->InRPSelection()) continue; // consider only tracks which are RPs
754 dPhi2 = aftsTrack->Phi();
755 // Fill the histogram:
756 fRelativeAngleDistribution->Fill(dPhi1-dPhi2);
758 } // end of for(Int_t j=0;j<nPrim;j++)
759 } else // to if(aftsTrack)
762 cout<<" WARNING (NL): No particle! (i.e. aftsTrack is a NULL pointer in AFAWNL::Make().)"<<endl;
765 } // end of for(Int_t i=0;i<nPrim;i++)
767 } // end of void AliFlowAnalysisWithNestedLoops::EvaluateNestedLoopsForRAD(AliFlowEventSimple *anEvent)
769 //================================================================================================================
771 void AliFlowAnalysisWithNestedLoops::EvaluateNestedLoopsForMH(AliFlowEventSimple *anEvent)
773 // Evaluate nested loops needed for mixed harmonics.
774 // Remark: phi labels the azimuthal angle of RP particle and psi labels azimuthal angle of POI particle.
776 Int_t nPrim = anEvent->NumberOfTracks();
777 Int_t nRP = anEvent->GetEventNSelTracksRP();
778 AliFlowTrackSimple *aftsTrack = NULL;
779 Double_t phi1=0.,phi2=0.,phi3=0.; // azimuthal angles of RPs
780 Double_t psi1=0.,psi2=0.; // azimuthal angles of POIs
781 Int_t charge1=0,charge2=0; // charge of POIs
782 Double_t pt1=0.,pt2=0.; // transverse momenta of POIs
785 // Evaluting correlator cos[n(phi1+phi2-2*phi3)] with three nested loops:
788 for(Int_t i1=0;i1<nPrim;i1++)
790 aftsTrack=anEvent->GetTrack(i1);
791 if(!(aftsTrack->InRPSelection())) continue;
792 phi1 = aftsTrack->Phi();
793 for(Int_t i2=0;i2<nPrim;i2++)
796 aftsTrack = anEvent->GetTrack(i2);
797 if(!(aftsTrack->InRPSelection())) continue;
798 phi2 = aftsTrack->Phi();
799 for(Int_t i3=0;i3<nPrim;i3++)
801 if(i3==i1||i3==i2) continue;
802 aftsTrack=anEvent->GetTrack(i3);
803 if(!(aftsTrack->InRPSelection())) continue;
804 phi3 = aftsTrack->Phi();
805 f3pCorrelatorPro->Fill(0.5,cos(n*(phi1+phi2-2.*phi3)),1.);
806 } // end of for(Int_t i3=0;i3<nPrim;i3++)
807 } // end of for(Int_t i2=0;i2<nPrim;i2++)
808 } // end of for(Int_t i1=0;i1<nPrim;i1++)
811 // Evaluting correlator cos[n(psi1+psi2-2*phi3)] with three nested loops:
812 if(!fEvaluateDifferential3pCorrelator){return;}
813 for(Int_t i1=0;i1<nPrim;i1++)
815 aftsTrack=anEvent->GetTrack(i1);
816 if(!(aftsTrack->InPOISelection())){continue;}
817 psi1 = aftsTrack->Phi();
818 pt1 = aftsTrack->Pt();
819 charge1 = aftsTrack->Charge();
820 for(Int_t i2=0;i2<nPrim;i2++)
822 if(i2==i1){continue;}
823 aftsTrack = anEvent->GetTrack(i2);
824 if(!(aftsTrack->InPOISelection())){continue;}
825 psi2 = aftsTrack->Phi();
826 pt2 = aftsTrack->Pt();
827 charge2 = aftsTrack->Charge();
828 if(fOppositeChargesPOI && charge1 == charge2){continue;}
829 for(Int_t i3=0;i3<nPrim;i3++)
831 if(i3==i1||i3==i2){continue;}
832 aftsTrack=anEvent->GetTrack(i3);
833 if(!(aftsTrack->InRPSelection())){continue;}
834 phi3 = aftsTrack->Phi();
835 // Evaluate and store differential correlator cos[n(psi1+psi2-2*phi3)]:
836 Double_t ptSum = (pt1+pt2)/2.;
837 Double_t ptDiff = TMath::Abs(pt1-pt2);
838 Double_t diff3pCorrelator = TMath::Cos(n*(psi1+psi2-2.*phi3));
839 f3pCorrelatorVsPtSumDiffDirectPro[0]->Fill(ptSum,diff3pCorrelator,1.);
840 f3pCorrelatorVsPtSumDiffDirectPro[1]->Fill(ptDiff,diff3pCorrelator,1.);
841 } // end of for(Int_t i3=0;i3<nPrim;i3++)
842 } // end of for(Int_t i2=0;i2<nPrim;i2++)
843 } // end of for(Int_t i1=0;i1<nPrim;i1++)
845 } // end of void AliFlowAnalysisWithNestedLoops::EvaluateNestedLoopsForMH(AliFlowEventSimple *anEvent)
847 //================================================================================================================
849 void AliFlowAnalysisWithNestedLoops::PrintOnTheScreen()
851 // Print on the screen.
854 cout<<"****************************************************"<<endl;
855 cout<<"****************************************************"<<endl;
856 cout<<" Nested Loops "<<endl;
859 if(fEvaluateNestedLoopsForRAD)
861 cout<<" Evaluated for relative angle distribution."<<endl;
864 if(fEvaluateNestedLoopsForMH)
866 cout<<" Evaluated for mixed harmonics:"<<endl;
869 cout<<" cos["<<fHarmonic<<"(phi1+phi2-2phi3)] = "<<f3pCorrelatorPro->GetBinContent(1)<<" +/- " <<f3pCorrelatorPro->GetBinError(1)<<endl;
872 cout<<" cos(phi1+phi2-2phi3) = "<<f3pCorrelatorPro->GetBinContent(1)<<" +/- " <<f3pCorrelatorPro->GetBinError(1)<<endl;
874 if(fEvaluateDifferential3pCorrelator)
876 cout<<" Evaluated also for differential 3p correlator "<<endl;
879 cout<<" cos["<<fHarmonic<<"(psi1+psi2-2phi3)]. "<<endl;
882 cout<<" cos(psi1+psi2-2phi3). "<<endl;
884 } // end of if(fEvaluateDifferential3pCorrelator)
885 } // end of if(fEvaluateNestedLoopsForMH)
887 if(!fEvaluateNestedLoopsForRAD && !fEvaluateNestedLoopsForMH)
889 cout<<" Not evaluated."<<endl;
892 cout<<"****************************************************"<<endl;
893 cout<<"****************************************************"<<endl;
896 } // end of void AliFlowAnalysisWithNestedLoops::PrintOnTheScreen()