AliPhysics  2b88e80 (2b88e80)
AliFlowEventSimpleMaker.cxx
Go to the documentation of this file.
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  * *
4  * Author: The ALICE Off-line Project. *
5  * Contributors are mentioned in the code where appropriate. *
6  * *
7  * Permission to use, copy, modify and distribute this software and its *
8  * documentation strictly for non-commercial purposes is hereby granted *
9  * without fee, provided that the above copyright notice appears in all *
10  * copies and that both the copyright notice and this permission notice *
11  * appear in the supporting documentation. The authors make no claims *
12  * about the suitability of this software for any purpose. It is *
13  * provided "as is" without express or implied warranty. *
14  **************************************************************************/
16 // AliFlowEventSimpleMaker:
17 //
18 // Class to fill the AliFlowEventSimple
19 // with AliFlowTrackSimple objects
20 // Has fill methods for TTree, AliMCEvent, AliESDEvent and AliAODEvent
21 // author: N. van der Kolk (kolk@nikhef.nl)
23 
24 #include "Riostream.h"
25 #include "TTree.h"
26 #include "TParticle.h"
28 #include "AliFlowEventSimple.h"
29 #include "AliFlowTrackSimple.h"
30 #include "AliMCEvent.h"
31 #include "AliMCParticle.h"
32 #include "AliESDEvent.h"
33 #include "AliESDtrack.h"
34 #include "AliAODEvent.h"
35 #include "AliAODTrack.h"
36 #include "AliCFManager.h"
37 #include "AliFlowTrackSimpleCuts.h"
38 #include "assert.h"
39 
40 using std::endl;
41 using std::cout;
43 //-----------------------------------------------------------------------
45  fMCReactionPlaneAngle(0.),
46  fCount(0),
47  fNoOfLoops(1),
48  fEllipticFlowValue(0.),
49  fMultiplicityOfEvent(1000000000),
50  fMinMult(0),
51  fMaxMult(1000000000),
52  fEtaMinA(-1.0),
53  fEtaMaxA(-0.01),
54  fEtaMinB(0.01),
55  fEtaMaxB(1.0)
56 {
57  //constructor
58 }
59 
60 //-----------------------------------------------------------------------
62 {
63  //destructor
64 }
65 
66 //-----------------------------------------------------------------------
67 AliFlowEventSimple* AliFlowEventSimpleMaker::FillTracks( AliMCEvent* anInput, const AliCFManager* intCFManager, const AliCFManager* diffCFManager)
68 {
69  //Fills the event from the MC kinematic information
70 
71  Int_t iNumberOfInputTracks = anInput->GetNumberOfTracks() ;
72 
73  if (iNumberOfInputTracks==-1) {
74  cout<<"Skipping Event -- No MC information available for this event"<<endl;
75  return 0;
76  }
77 
78  Int_t iN = iNumberOfInputTracks; //maximum number of tracks in AliFlowEventSimple
79  Int_t iGoodTracks = 0; //number of good tracks
80  Int_t itrkN = 0; //track counter
81  Int_t iSelParticlesPOI = 0; //number of tracks selected for Diff
82  Int_t iSelParticlesRP = 0; //number of tracks selected for Int
83 
84  // cut on the multiplicity
85  if (intCFManager->CheckEventCuts(AliCFManager::kEvtGenCuts,anInput)) {
86  // cout<<"iNumberOfInputTracks = "<<iNumberOfInputTracks<<endl;
87  // create an AliFlowEventSimple
88  AliFlowEventSimple* pEvent = new AliFlowEventSimple(10);
89 
90  //loop over tracks
91  while (iGoodTracks < iN && itrkN < iNumberOfInputTracks) {
92  //get input particle
93  AliMCParticle* pParticle = (AliMCParticle*) anInput->GetTrack(itrkN);
94  //make new AliFlowTrackSimple
95  AliFlowTrackSimple* pTrack = new AliFlowTrackSimple();
96  pTrack->SetPt(pParticle->Pt() );
97  pTrack->SetEta(pParticle->Eta() );
98  pTrack->SetPhi(pParticle->Phi() );
99 
100  //check if pParticle passes the cuts
101  if (intCFManager->CheckParticleCuts(AliCFManager::kPartGenCuts,pParticle)) {
102  pTrack->SetForRPSelection(kTRUE);
103  //cout<<"integrated selection. PID = "<<pParticle->Particle()->GetPdgCode()<<endl;
104  }
105  if (diffCFManager->CheckParticleCuts(AliCFManager::kPartGenCuts,pParticle)) {
106  pTrack->SetForPOISelection(kTRUE);
107  //cout<<"differential selection. PID = "<<pParticle->Particle()->GetPdgCode()<<endl;
108  }
109 
110  //check if any bits are set
111  const TBits* bFlowBits = pTrack->GetFlowBits();
112  if (bFlowBits->CountBits() ==0) {
113  delete pTrack; } //track will not be used anymore
114  else {
115  pEvent->AddTrack(pTrack) ;
116  iGoodTracks++;
117 
118  if (pTrack->InRPSelection())
119  { iSelParticlesRP++; }
120  if (pTrack->InPOISelection())
121  { iSelParticlesPOI++; }
122  }
123 
124  itrkN++;
125  }
126 
127  pEvent-> SetEventNSelTracksRP(iSelParticlesRP);
129 
130  if (iSelParticlesRP >= fMinMult && iSelParticlesRP <= fMaxMult) {
131  if ( (++fCount % 100) == 0) {
132  cout<<" MC Reaction Plane Angle = "<< fMCReactionPlaneAngle << endl;
133  //
134  cout<<" iGoodTracks = "<<iGoodTracks<<endl;
135  cout<<" # of RP selected tracks = "<<iSelParticlesRP<<endl;
136  cout<<" # of POI selected tracks = "<<iSelParticlesPOI<<endl;
137  cout << "# " << fCount << " events processed" << endl;
138  }
139  return pEvent;
140  }
141  else {
142  cout<<"Not enough tracks in the FlowEventSimple"<<endl;
143  return 0;
144  }
145  }
146  else {
147  cout<<"Event does not pass multiplicity cuts"<<endl;
148  return 0;
149  }
150 
151 }
152 
153 //-----------------------------------------------------------------------
154 
155 AliFlowEventSimple* AliFlowEventSimpleMaker::FillTracks(AliESDEvent* anInput, const AliCFManager* intCFManager, const AliCFManager* diffCFManager)
156 {
157  //Fills the event from the ESD
158 
159  //flags for particles passing int. and diff. flow cuts
160  Bool_t bPassedRPFlowCuts = kFALSE;
161  Bool_t bPassedPOIFlowCuts = kFALSE;
162 
163  Int_t iNumberOfInputTracks = anInput->GetNumberOfTracks() ;
164 
165  Int_t iGoodTracks = 0; //number of good tracks
166  Int_t itrkN = 0; //track counter
167  Int_t iSelParticlesRP = 0; //number of tracks selected for Int
168  Int_t iSelParticlesPOI = 0; //number of tracks selected for Diff
169 
170  // cut on the multiplicity
171  if (intCFManager->CheckEventCuts(AliCFManager::kEvtRecCuts,anInput)) {
172  // cout<<"iNumberOfInputTracks = "<<iNumberOfInputTracks<<endl;
173  // create an AliFlowEventSimple
174  AliFlowEventSimple* pEvent = new AliFlowEventSimple(10);
175 
176  //loop over tracks
177  while (itrkN < iNumberOfInputTracks) {
178  AliESDtrack* pParticle = anInput->GetTrack(itrkN); //get input particle
179 
180  //check if pParticle passes the cuts
181  if (intCFManager->CheckParticleCuts(AliCFManager::kPartRecCuts,pParticle) &&
182  intCFManager->CheckParticleCuts(AliCFManager::kPartSelCuts,pParticle)) {
183  bPassedRPFlowCuts = kTRUE;
184  }
185  if (diffCFManager->CheckParticleCuts(AliCFManager::kPartRecCuts,pParticle) &&
186  diffCFManager->CheckParticleCuts(AliCFManager::kPartSelCuts,pParticle)) {
187  bPassedPOIFlowCuts = kTRUE;
188  }
189 
190  if (bPassedRPFlowCuts || bPassedPOIFlowCuts) {
191  //make new AliFLowTrackSimple
192  AliFlowTrackSimple* pTrack = new AliFlowTrackSimple();
193  pTrack->SetPt(pParticle->Pt() );
194  pTrack->SetEta(pParticle->Eta() );
195  if (fEllipticFlowValue>0.)
196  { pTrack->SetPhi(pParticle->Phi()-fEllipticFlowValue*TMath::Sin(2*(pParticle->Phi()-fMCReactionPlaneAngle))); cout<<"Added flow to particle"<<endl; }
197  else { pTrack->SetPhi(pParticle->Phi() ); }
198 
199  //marking the particles used for int. flow:
200  if(bPassedRPFlowCuts) {
201  pTrack->SetForRPSelection(kTRUE);
202  iSelParticlesRP++;
203  // assign particles to subevents
204  if (pTrack->Eta()>=fEtaMinA && pTrack->Eta()<=fEtaMaxA) {
205  pTrack->SetForSubevent(0);
206  }
207  if (pTrack->Eta()>=fEtaMinB && pTrack->Eta()<=fEtaMaxB) {
208  pTrack->SetForSubevent(1);
209  }
210 
211  }
212  //marking the particles used for diff. flow:
213  if(bPassedPOIFlowCuts) {
214  pTrack->SetForPOISelection(kTRUE);
215  iSelParticlesPOI++;
216  }
217  //adding particles which were used either for int. or diff. flow to the list
218  pEvent->AddTrack(pTrack);
219  iGoodTracks++;
220  }//end of if(bPassedIntFlowCuts || bPassedDiffFlowCuts)
221  itrkN++;
222  bPassedRPFlowCuts = kFALSE;
223  bPassedPOIFlowCuts = kFALSE;
224  }//end of while (itrkN < iNumberOfInputTracks)
225 
226  pEvent->SetEventNSelTracksRP(iSelParticlesRP);
228 
229 
230  if (iSelParticlesRP >= fMinMult && iSelParticlesRP <= fMaxMult) {
231  if ( (++fCount % 100) == 0) {
232  if (fMCReactionPlaneAngle != 0) cout<<" MC Reaction Plane Angle = "<< fMCReactionPlaneAngle << endl;
233  else cout<<" MC Reaction Plane Angle = unknown "<< endl;
234  cout<<" iGoodTracks = "<<iGoodTracks<<endl;
235  cout<<" # of RP selected tracks = "<<iSelParticlesRP<<endl;
236  cout<<" # of POI selected tracks = "<<iSelParticlesPOI<<endl;
237  cout << "# " << fCount << " events processed" << endl;
238  }
239  return pEvent;
240  }
241  else {
242  cout<<"Not enough tracks in the FlowEventSimple"<<endl;
243  return 0;
244  }
245  }
246  else {
247  cout<<"Event does not pass multiplicity cuts"<<endl;
248  return 0;
249  }
250 
251 }
252 
253 //-----------------------------------------------------------------------
254 AliFlowEventSimple* AliFlowEventSimpleMaker::FillTracks(AliAODEvent* anInput, const AliCFManager* intCFManager, const AliCFManager* diffCFManager)
255 {
256  //Fills the event from the AOD
257 
258  Int_t iNumberOfInputTracks = anInput->GetNumberOfTracks() ;
259 
260  Int_t iN = iNumberOfInputTracks; //maximum number of tracks in AliFlowEventSimple
261  Int_t iGoodTracks = 0; //number of good tracks
262  Int_t itrkN = 0; //track counter
263  Int_t iSelParticlesPOI = 0; //number of tracks selected for Diff
264  Int_t iSelParticlesRP = 0; //number of tracks selected for Int
265 
266  // cut on the multiplicity
267  if (intCFManager->CheckEventCuts(AliCFManager::kEvtRecCuts,anInput)) {
268  // cout<<"iNumberOfInputTracks = "<<iNumberOfInputTracks<<endl;
269  // create an AliFlowEventSimple
270  AliFlowEventSimple* pEvent = new AliFlowEventSimple(10);
271 
272  //loop over tracks
273  while (iGoodTracks < iN && itrkN < iNumberOfInputTracks) {
274  AliAODTrack* pParticle = dynamic_cast<AliAODTrack*>(anInput->GetTrack(itrkN));
275  assert((pParticle)&&"Not a standard AOD"); //get input particle
276  //make new AliFlowTrackSimple
277  AliFlowTrackSimple* pTrack = new AliFlowTrackSimple();
278  pTrack->SetPt(pParticle->Pt() );
279  pTrack->SetEta(pParticle->Eta() );
280  pTrack->SetPhi(pParticle->Phi() );
281 
282  //check if pParticle passes the cuts
283  if (intCFManager->CheckParticleCuts(AliCFManager::kPartRecCuts,pParticle) &&
284  intCFManager->CheckParticleCuts(AliCFManager::kPartSelCuts,pParticle)) {
285  pTrack->SetForRPSelection(kTRUE); }
286  if (diffCFManager->CheckParticleCuts(AliCFManager::kPartRecCuts,pParticle) &&
287  diffCFManager->CheckParticleCuts(AliCFManager::kPartSelCuts,pParticle)) {
288  pTrack->SetForPOISelection(kTRUE);}
289 
290 
291  //check if any bits are set
292  const TBits* bFlowBits = pTrack->GetFlowBits();
293  if (bFlowBits->CountBits() ==0) {
294  delete pTrack; } //track will not be used anymore
295  else {
296  pEvent->AddTrack(pTrack) ;
297  iGoodTracks++;
298 
299  if (pTrack->InRPSelection())
300  { iSelParticlesRP++; }
301  if (pTrack->InPOISelection())
302  { iSelParticlesPOI++; }
303 
304  }
305 
306  itrkN++;
307  }
308 
309  pEvent-> SetEventNSelTracksRP(iSelParticlesRP);
311 
312  if (iSelParticlesRP >= fMinMult && iSelParticlesRP <= fMaxMult) {
313  if ( (++fCount % 100) == 0) {
314  if (fMCReactionPlaneAngle != 0) cout<<" MC Reaction Plane Angle = "<< fMCReactionPlaneAngle << endl;
315  else cout<<" MC Reaction Plane Angle = unknown "<< endl;
316  cout<<" iGoodTracks = "<<iGoodTracks<<endl;
317  cout<<" # of RP selected tracks = "<<iSelParticlesRP<<endl;
318  cout<<" # of POI selected tracks = "<<iSelParticlesPOI<<endl;
319  cout << "# " << fCount << " events processed" << endl;
320  }
321  return pEvent;
322  }
323  else {
324  cout<<"Not enough tracks in the FlowEventSimple"<<endl;
325  return 0;
326  }
327  }
328  else {
329  cout<<"Event does not pass multiplicity cuts"<<endl;
330  return 0;
331  }
332 
333 }
334 
335 //-----------------------------------------------------------------------
336 AliFlowEventSimple* AliFlowEventSimpleMaker::FillTracks(AliESDEvent* anInput, const AliMCEvent* anInputMc, const AliCFManager* intCFManager, const AliCFManager* diffCFManager, Int_t anOption)
337 {
338  //fills the event with tracks from the ESD and kinematics from the MC info via the track label
339 
340 
341  if (!(anOption ==0 || anOption ==1)) {
342  cout<<"WRONG OPTION IN AliFlowEventSimpleMaker::FillTracks(AliESDEvent* anInput, AliMCEvent* anInputMc, Int_t anOption)"<<endl;
343  exit(1);
344  }
345 
346  Int_t iNumberOfInputTracks = anInput->GetNumberOfTracks() ;
347 
348  Int_t iNumberOfInputTracksMC = anInputMc->GetNumberOfTracks() ;
349  if (iNumberOfInputTracksMC==-1) {
350  cout<<"Skipping Event -- No MC information available for this event"<<endl;
351  return 0;
352  }
353 
354  Int_t iN = iNumberOfInputTracks; //maximum number of tracks in AliFlowEventSimple
355  Int_t iGoodTracks = 0; //number of good tracks
356  Int_t itrkN = 0; //track counter
357  Int_t iSelParticlesPOI = 0; //number of tracks selected for Diff
358  Int_t iSelParticlesRP = 0; //number of tracks selected for Int
359 
360  // cut on the multiplicity
361  if (intCFManager->CheckEventCuts(AliCFManager::kEvtRecCuts,anInput)) {
362  // cout<<"iNumberOfInputTracks = "<<iNumberOfInputTracks<<endl;
363  // create an AliFlowEventSimple
364  AliFlowEventSimple* pEvent = new AliFlowEventSimple(10);
365 
366  //loop over ESD tracks
367  while (iGoodTracks < iN && itrkN < iNumberOfInputTracks) {
368  AliESDtrack* pParticle = anInput->GetTrack(itrkN); //get input particle
369  //get Label
370  Int_t iLabel = pParticle->GetLabel();
371  //match to mc particle
372  AliMCParticle* pMcParticle = (AliMCParticle*) anInputMc->GetTrack(TMath::Abs(iLabel));
373 
374  //check
375  if (TMath::Abs(pParticle->GetLabel())!=pMcParticle->Label()) cout<<"pParticle->GetLabel()!=pMcParticle->Label() "<<pParticle->GetLabel()<<" "<<pMcParticle->Label()<<endl;
376 
377  //make new AliFlowTrackSimple
378  AliFlowTrackSimple* pTrack = new AliFlowTrackSimple();
379  if(anOption == 0) { //take the PID from the MC & the kinematics from the ESD
380  pTrack->SetPt(pParticle->Pt() );
381  pTrack->SetEta(pParticle->Eta() );
382  pTrack->SetPhi(pParticle->Phi() );
383  }
384  else if (anOption == 1) { //take the PID and kinematics from the MC
385  pTrack->SetPt(pMcParticle->Pt() );
386  pTrack->SetEta(pMcParticle->Eta() );
387  pTrack->SetPhi(pMcParticle->Phi() );
388  }
389  else { cout<<"Not a valid option"<<endl; }
390 
391  //check if pParticle passes the cuts
392  if(anOption == 0) {
393  //cout<<"take the PID from the MC & the kinematics from the ESD"<<endl;
394  if (intCFManager->CheckParticleCuts(AliCFManager::kPartGenCuts,pMcParticle,"mcGenCuts1") &&
395  intCFManager->CheckParticleCuts(AliCFManager::kPartRecCuts,pParticle)) {
396  pTrack->SetForRPSelection(kTRUE); }
397  if (diffCFManager->CheckParticleCuts(AliCFManager::kPartGenCuts,pMcParticle,"mcGenCuts2") &&
398  diffCFManager->CheckParticleCuts(AliCFManager::kPartRecCuts,pParticle)) {
399  pTrack->SetForPOISelection(kTRUE);}
400  }
401  else if (anOption == 1) {
402  //cout<<"take the PID and kinematics from the MC"<<endl;
403  if (intCFManager->CheckParticleCuts(AliCFManager::kPartGenCuts,pMcParticle)) {
404  pTrack->SetForRPSelection(kTRUE); }
405  if (diffCFManager->CheckParticleCuts(AliCFManager::kPartGenCuts,pMcParticle)) {
406  pTrack->SetForPOISelection(kTRUE);}
407  }
408  else { cout<<"Not a valid option"<<endl; }
409 
410  //check if any bits are set
411  const TBits* bFlowBits = pTrack->GetFlowBits();
412  if (bFlowBits->CountBits() ==0) {
413  delete pTrack; } //track will not be used anymore
414  else {
415  pEvent->AddTrack(pTrack) ;
416  iGoodTracks++;
417 
418  if (pTrack->InRPSelection())
419  { iSelParticlesRP++; }
420  if (pTrack->InPOISelection())
421  { iSelParticlesPOI++; }
422 
423  }
424 
425  itrkN++;
426  }
427 
428  pEvent->SetEventNSelTracksRP(iSelParticlesRP);
430 
431  if (iSelParticlesRP >= fMinMult && iSelParticlesRP <= fMaxMult) {
432  if ( (++fCount % 100) == 0) {
433  if (fMCReactionPlaneAngle != 0) cout<<" MC Reaction Plane Angle = "<< fMCReactionPlaneAngle << endl;
434  else cout<<" MC Reaction Plane Angle = unknown "<< endl;
435  cout << " Number of MC input tracks = " << iNumberOfInputTracksMC << endl;
436  cout<<" iGoodTracks = "<<iGoodTracks<<endl;
437  cout<<" # of RP selected tracks = "<<iSelParticlesRP<<endl;
438  cout<<" # of POI selected tracks = "<<iSelParticlesPOI<<endl;
439  cout << "# " << fCount << " events processed" << endl;
440  }
441  return pEvent;
442  }
443  else {
444  cout<<"Not enough tracks in the FlowEventSimple"<<endl;
445  return 0;
446  }
447  }
448  else {
449  cout<<"Event does not pass multiplicity cuts"<<endl;
450  return 0;
451  }
452 
453 }
454 
455 //local methods
456 //-----------------------------------------------------------------------
458 {
459  //fills the event from a TTree of kinematic.root files
460 
461  // number of times to use the same particle (trick to introduce nonflow)
462 
463  //flags for particles passing int. and diff. flow cuts
464  Bool_t bPassedRPFlowCuts = kFALSE;
465  Bool_t bPassedPOIFlowCuts = kFALSE;
466 
467  //track cut values
468  Double_t dPtMaxRP = rpCuts->GetPtMax();
469  Double_t dPtMinRP = rpCuts->GetPtMin();
470  Double_t dEtaMaxRP = rpCuts->GetEtaMax();
471  Double_t dEtaMinRP = rpCuts->GetEtaMin();
472  Double_t dPhiMaxRP = rpCuts->GetPhiMax();
473  Double_t dPhiMinRP = rpCuts->GetPhiMin();
474  Int_t iPIDRP = rpCuts->GetPID();
475 
476  Double_t dPtMaxPOI = poiCuts->GetPtMax();
477  Double_t dPtMinPOI = poiCuts->GetPtMin();
478  Double_t dEtaMaxPOI = poiCuts->GetEtaMax();
479  Double_t dEtaMinPOI = poiCuts->GetEtaMin();
480  Double_t dPhiMaxPOI = poiCuts->GetPhiMax();
481  Double_t dPhiMinPOI = poiCuts->GetPhiMin();
482  Int_t iPIDPOI = poiCuts->GetPID();
483 
484  Int_t iNumberOfInputTracks = anInput->GetEntries() ;
485 
486  TParticle* pParticle = new TParticle();
487  anInput->SetBranchAddress("Particles",&pParticle);
488  // AliFlowEventSimple* pEvent = new AliFlowEventSimple(iNumberOfInputTracks);
489  AliFlowEventSimple* pEvent = new AliFlowEventSimple(10);
490 
491  // Int_t fMultiplicityOfEvent = 576; //multiplicity for chi=1.5
492  // Int_t fMultiplicityOfEvent = 256; //multiplicity for chi=1
493  // Int_t fMultiplicityOfEvent = 164; //multiplicity for chi=0.8
494 
495  Int_t iGoodTracks = 0;
496  Int_t itrkN = 0;
497  Int_t iSelParticlesRP = 0;
498  Int_t iSelParticlesPOI = 0;
499 
500  while (itrkN < iNumberOfInputTracks) {
501  anInput->GetEntry(itrkN); //get input particle
502  if (pParticle->IsPrimary()) {
503  //checking the cuts for int. and diff. flow
504  if (pParticle->Pt() > dPtMinRP && pParticle->Pt() < dPtMaxRP &&
505  pParticle->Eta() > dEtaMinRP && pParticle->Eta() < dEtaMaxRP &&
506  pParticle->Phi() > dPhiMinRP && pParticle->Phi() < dPhiMaxRP &&
507  TMath::Abs(pParticle->GetPdgCode()) == iPIDRP) {
508  bPassedRPFlowCuts = kTRUE;
509  }
510 
511  if (pParticle->Pt() > dPtMinPOI && pParticle->Pt() < dPtMaxPOI &&
512  pParticle->Eta() > dEtaMinPOI && pParticle->Eta() < dEtaMaxPOI &&
513  pParticle->Phi() > dPhiMinPOI && pParticle->Phi() < dPhiMaxPOI &&
514  TMath::Abs(pParticle->GetPdgCode()) == iPIDPOI){
515  bPassedPOIFlowCuts = kTRUE;
516  }
517  }
518  if (bPassedRPFlowCuts || bPassedPOIFlowCuts) {
519  for(Int_t d=0;d<fNoOfLoops;d++) {
520  AliFlowTrackSimple* pTrack = new AliFlowTrackSimple();
521  pTrack->SetPt(pParticle->Pt());
522  pTrack->SetEta(pParticle->Eta());
523  pTrack->SetPhi(pParticle->Phi()-fEllipticFlowValue*TMath::Sin(2*(pParticle->Phi()-fMCReactionPlaneAngle)));
524 
525  //marking the particles used for int. flow:
526  if(bPassedRPFlowCuts && iSelParticlesRP < fMultiplicityOfEvent) {
527  pTrack->SetForRPSelection(kTRUE);
528  iSelParticlesRP++;
529  }
530  //marking the particles used for diff. flow:
531  if(bPassedPOIFlowCuts && iGoodTracks%fNoOfLoops==0) {
532  pTrack->SetForPOISelection(kTRUE);
533  iSelParticlesPOI++;
534  }
535  //adding a particles which were used either for int. or diff. flow to the list
536  pEvent->AddTrack(pTrack);
537  iGoodTracks++;
538  }//end of for(Int_t d=0;d<iLoops;d++)
539  }//end of if(bPassedIntFlowCuts || bPassedDiffFlowCuts)
540  itrkN++;
541  bPassedRPFlowCuts = kFALSE;
542  bPassedPOIFlowCuts = kFALSE;
543  }//end of while (itrkN < iNumberOfInputTracks)
544 
545  pEvent->SetEventNSelTracksRP(iSelParticlesRP);
547 
548  if ( (++fCount % 100) == 0) {
549  if (fMCReactionPlaneAngle != 0) cout<<" MC Reaction Plane Angle = "<< fMCReactionPlaneAngle << endl;
550  else cout<<" MC Reaction Plane Angle = unknown "<< endl;
551  cout<<" iGoodTracks = "<< iGoodTracks << endl;
552  cout<<" # of RP selected tracks = "<<iSelParticlesRP<<endl;
553  cout<<" # of POI selected tracks = "<<iSelParticlesPOI<<endl;
554  cout << "# " << fCount << " events processed" << endl;
555  }
556 
557  delete pParticle;
558  return pEvent;
559 }
560 
561 //-----------------------------------------------------------------------
563 {
564  //Fills the event from the MC kinematic information
565 
566  Int_t iNumberOfInputTracks = anInput->GetNumberOfTracks() ;
567 
568  AliFlowEventSimple* pEvent = new AliFlowEventSimple(10);
569 
570  //Int_t iN = 256; //multiplicity for chi=1
571  Int_t iN = iNumberOfInputTracks;
572  Int_t iGoodTracks = 0;
573  Int_t itrkN = 0;
574  Int_t iSelParticlesPOI = 0;
575  Int_t iSelParticlesRP = 0;
576 
577  //normal loop
578  while (iGoodTracks < iN && itrkN < iNumberOfInputTracks) {
579  AliMCParticle* pParticle = (AliMCParticle*) anInput->GetTrack(itrkN); //get input particle
580  //cut on tracks
581  if (TMath::Abs(pParticle->Eta()) < 0.9)
582  {
583  if(
584  TMath::Abs(pParticle->Particle()->GetPdgCode()) == 211
585  // TMath::Abs(pParticle->Particle()->GetPdgCode()) == 211 ||
586  // TMath::Abs(pParticle->Particle()->GetPdgCode()) == 321 ||
587  // TMath::Abs(pParticle->Particle()->GetPdgCode()) == 2212
588  )
589  {
590  AliFlowTrackSimple* pTrack = new AliFlowTrackSimple();
591  pTrack->SetPt(pParticle->Pt() );
592  pTrack->SetEta(pParticle->Eta() );
593  pTrack->SetPhi(pParticle->Phi() );
594  pTrack->SetForRPSelection(kTRUE);
595  pTrack->SetForPOISelection(kTRUE);
596 
597  if (pTrack->InRPSelection())
598  { iSelParticlesRP++; }
599  if (pTrack->InPOISelection())
600  { iSelParticlesPOI++; }
601  iGoodTracks++;
602  pEvent->AddTrack(pTrack) ;
603  }
604  /* else if(
605  TMath::Abs(pParticle->Particle()->GetPdgCode()) == 211
606  )
607  {
608  AliFlowTrackSimple* pTrack = new AliFlowTrackSimple();
609  pTrack->SetPt(pParticle->Pt() );
610  pTrack->SetEta(pParticle->Eta() );
611  pTrack->SetPhi(pParticle->Phi() );
612  pTrack->SetForRPSelection(kFALSE);
613  pTrack->SetForPOISelection(kTRUE);
614 
615  if (pTrack->InRPSelection())
616  { iSelParticlesRP++; }
617  if (pTrack->InPOISelection())
618  { iSelParticlesPOI++; }
619  iGoodTracks++;
620  pEvent->AddTrack(pTrack);
621  }
622  */
623  }
624 
625  itrkN++;
626  }
627 
628  pEvent-> SetEventNSelTracksRP(iSelParticlesRP);
630 
631  if ( (++fCount % 100) == 0) {
632  if (fMCReactionPlaneAngle != 0) cout<<" MC Reaction Plane Angle = "<< fMCReactionPlaneAngle << endl;
633  else cout<<" MC Reaction Plane Angle = unknown "<< endl;
634  cout<<" iGoodTracks = "<<iGoodTracks<<endl;
635  cout<<" # of RP selected tracks = "<<iSelParticlesRP<<endl;
636  cout<<" # of POI selected tracks = "<<iSelParticlesPOI<<endl;
637  cout << "# " << fCount << " events processed" << endl;
638  }
639 
640  return pEvent;
641 
642 }
643 
644 //-----------------------------------------------------------------------
646 {
647  //Fills the event from the ESD
648 
649  Int_t iNumberOfInputTracks = anInput->GetNumberOfTracks() ;
650 
651  AliFlowEventSimple* pEvent = new AliFlowEventSimple(10);
652 
653  //Int_t iN = 256; //multiplicity for chi=1
654  Int_t iN = iNumberOfInputTracks;
655  Int_t iGoodTracks = 0;
656  Int_t itrkN = 0;
657  Int_t iSelParticlesPOI = 0;
658  Int_t iSelParticlesRP = 0;
659 
660 
661 
662  //normal loop
663  while (iGoodTracks < iN && itrkN < iNumberOfInputTracks) {
664  AliESDtrack* pParticle = anInput->GetTrack(itrkN); //get input particle
665  //cut on tracks
666  if (TMath::Abs(pParticle->Eta()) < 0.9)
667 
668 
669 
670  {
671  AliFlowTrackSimple* pTrack = new AliFlowTrackSimple();
672  pTrack->SetPt(pParticle->Pt() );
673  pTrack->SetEta(pParticle->Eta() );
674  pTrack->SetPhi(pParticle->Phi() );
675  pTrack->SetForRPSelection(kTRUE);
676  pTrack->SetForPOISelection(kTRUE);
677 
678  if (pTrack->InRPSelection())
679  { iSelParticlesRP++; }
680  if (pTrack->InPOISelection())
681  { iSelParticlesPOI++; }
682  iGoodTracks++;
683  pEvent->AddTrack(pTrack) ;
684  }
685 
686  itrkN++;
687  }
688 
689  pEvent-> SetEventNSelTracksRP(iSelParticlesRP);
691 
692  if ( (++fCount % 100) == 0) {
693  if (fMCReactionPlaneAngle != 0) cout<<" MC Reaction Plane Angle = "<< fMCReactionPlaneAngle << endl;
694  else cout<<" MC Reaction Plane Angle = unknown "<< endl;
695  cout<<" iGoodTracks = "<<iGoodTracks<<endl;
696  cout<<" # of RP selected tracks = "<<iSelParticlesRP<<endl;
697  cout<<" # of POI selected tracks = "<<iSelParticlesPOI<<endl;
698  cout << "# " << fCount << " events processed" << endl;
699  }
700 
701  return pEvent;
702 }
703 
704 //-----------------------------------------------------------------------
705 
707 {
708  //Fills the event from the AOD
709 
710  Int_t iNumberOfInputTracks = anInput->GetNumberOfTracks() ;
711 
712  AliFlowEventSimple* pEvent = new AliFlowEventSimple(10);
713 
714  //Int_t iN = 256; //multiplicity for chi=1
715  Int_t iN = iNumberOfInputTracks;
716  Int_t iGoodTracks = 0;
717  Int_t itrkN = 0;
718  Int_t iSelParticlesPOI = 0;
719  Int_t iSelParticlesRP = 0;
720 
721  //normal loop
722  while (iGoodTracks < iN && itrkN < iNumberOfInputTracks) {
723  AliAODTrack* pParticle = dynamic_cast<AliAODTrack*>(anInput->GetTrack(itrkN));
724  assert((pParticle)&&"Not a standard AOD"); //get input particle
725  //cut on tracks
726  if (TMath::Abs(pParticle->Eta()) < 0.9)
727  {
728  AliFlowTrackSimple* pTrack = new AliFlowTrackSimple();
729  pTrack->SetPt(pParticle->Pt() );
730  pTrack->SetEta(pParticle->Eta() );
731  pTrack->SetPhi(pParticle->Phi() );
732  pTrack->SetForRPSelection(kTRUE);
733  pTrack->SetForPOISelection(kTRUE);
734 
735  if (pTrack->InRPSelection())
736  { iSelParticlesRP++; }
737  if (pTrack->InPOISelection())
738  { iSelParticlesPOI++; }
739  iGoodTracks++;
740  pEvent->AddTrack(pTrack) ;
741  }
742 
743  itrkN++;
744  }
745 
746  pEvent-> SetEventNSelTracksRP(iSelParticlesRP);
748 
749  if ( (++fCount % 100) == 0) {
750  if (fMCReactionPlaneAngle != 0) cout<<" MC Reaction Plane Angle = "<< fMCReactionPlaneAngle << endl;
751  else cout<<" MC Reaction Plane Angle = unknown "<< endl;
752  cout<<" iGoodTracks = "<<iGoodTracks<<endl;
753  cout<<" # of RP selected tracks = "<<iSelParticlesRP<<endl;
754  cout<<" # of POI selected tracks = "<<iSelParticlesPOI<<endl;
755  cout << "# " << fCount << " events processed" << endl;
756  }
757 
758  return pEvent;
759 }
760 
761 //-----------------------------------------------------------------------
762 AliFlowEventSimple* AliFlowEventSimpleMaker::FillTracks(AliESDEvent* anInput, const AliMCEvent* anInputMc, Int_t anOption)
763 {
764  //fills the event with tracks from the ESD and kinematics from the MC info via the track label
765 
766  if (!(anOption ==0 || anOption ==1)) {
767  cout<<"WRONG OPTION IN AliFlowEventSimpleMaker::FillTracks(AliESDEvent* anInput, AliMCEvent* anInputMc, Int_t anOption)"<<endl;
768  exit(1);
769  }
770 
771  Int_t iNumberOfInputTracks = anInput->GetNumberOfTracks() ;
772 
773  AliFlowEventSimple* pEvent = new AliFlowEventSimple(10);
774 
775  //Int_t iN = 256; //multiplicity for chi=1
776  Int_t iN = iNumberOfInputTracks;
777  Int_t iGoodTracks = 0;
778  Int_t itrkN = 0;
779  Int_t iSelParticlesPOI = 0;
780  Int_t iSelParticlesRP = 0;
781 
782  //normal loop
783  while (iGoodTracks < iN && itrkN < iNumberOfInputTracks) {
784  AliESDtrack* pParticle = anInput->GetTrack(itrkN); //get input particle
785  //get Label
786  Int_t iLabel = pParticle->GetLabel();
787  //match to mc particle
788  AliMCParticle* pMcParticle = (AliMCParticle*) anInputMc->GetTrack(TMath::Abs(iLabel));
789 
790  //check
791  if (TMath::Abs(pParticle->GetLabel())!=pMcParticle->Label()) cout<<"pParticle->GetLabel()!=pMcParticle->Label() "<<pParticle->GetLabel()<<" "<<pMcParticle->Label()<<endl;
792 
793  //cut on tracks
794  if (TMath::Abs(pParticle->Eta()) < 0.2)
795  {
796  if(
797  TMath::Abs(pMcParticle->Particle()->GetPdgCode()) == 211 //pions
798  // TMath::Abs(pMcParticle->Particle()->GetPdgCode()) == 211 ||
799  // TMath::Abs(pMcParticle->Particle()->GetPdgCode()) == 321 ||
800  // TMath::Abs(pMcParticle->Particle()->GetPdgCode()) == 2212
801  )
802  {
803  AliFlowTrackSimple* pTrack = new AliFlowTrackSimple();
804  if(anOption == 0) { //take the PID from the MC & the kinematics from the ESD
805  pTrack->SetPt(pParticle->Pt() );
806  pTrack->SetEta(pParticle->Eta() );
807  pTrack->SetPhi(pParticle->Phi() );
808  pTrack->SetForRPSelection(kTRUE);
809  pTrack->SetForPOISelection(kTRUE);
810  }
811  else if (anOption == 1) { //take the PID and kinematics from the MC
812  pTrack->SetPt(pMcParticle->Pt() );
813  pTrack->SetEta(pMcParticle->Eta() );
814  pTrack->SetPhi(pMcParticle->Phi() );
815  pTrack->SetForRPSelection(kTRUE);
816  pTrack->SetForPOISelection(kTRUE);
817  }
818  else { cout<<"Not a valid option"<<endl; }
819  if (pTrack->InRPSelection())
820  { iSelParticlesRP++; }
821  if (pTrack->InPOISelection())
822  { iSelParticlesPOI++; }
823  iGoodTracks++;
824  pEvent->AddTrack(pTrack) ;
825  }
826  }
827  itrkN++;
828  }
829 
830  pEvent-> SetEventNSelTracksRP(iSelParticlesRP);
832 
833  if ( (++fCount % 100) == 0) {
834  if (fMCReactionPlaneAngle != 0) cout<<" MC Reaction Plane Angle = "<< fMCReactionPlaneAngle << endl;
835  else cout<<" MC Reaction Plane Angle = unknown "<< endl;
836  cout<<" iGoodTracks = "<<iGoodTracks<<endl;
837  cout<<" # of RP selected tracks = "<<iSelParticlesRP<<endl;
838  cout<<" # of POI selected tracks = "<<iSelParticlesPOI<<endl;
839  cout << "# " << fCount << " events processed" << endl;
840  }
841 
842  return pEvent;
843 }
void SetForSubevent(Int_t i)
const TBits * GetFlowBits() const
double Double_t
Definition: External.C:58
void SetEta(Double_t eta)
void AddTrack(AliFlowTrackSimple *track)
Bool_t InRPSelection() const
void SetForRPSelection(Bool_t b=kTRUE)
AliFlowEventSimple * FillTracks(TTree *anInput, const AliFlowTrackSimpleCuts *rpCuts, const AliFlowTrackSimpleCuts *poiCuts)
int Int_t
Definition: External.C:63
void SetPhi(Double_t phi)
void SetPt(Double_t pt)
void SetEventNSelTracksRP(Int_t nr)
void SetForPOISelection(Bool_t b=kTRUE)
void SetMCReactionPlaneAngle(Double_t fPhiRP)
Bool_t InPOISelection(Int_t poiType=1) const
Double_t Eta() const
bool Bool_t
Definition: External.C:53