]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FLOW/FlowEventMakers/FlowEventSimpleMaker.cxx
name change Int/Diff RP/POI
[u/mrichter/AliRoot.git] / PWG2 / FLOW / FlowEventMakers / FlowEventSimpleMaker.cxx
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  **************************************************************************/
15 /* $Id */
16
17 #include "Riostream.h"
18 #include "FlowEventSimpleMaker.h"
19 #include "AliFlowCommon/AliFlowEventSimple.h"
20 #include "AliFlowCommon/AliFlowTrackSimple.h"
21 #include "TTree.h"
22 #include "TParticle.h"
23 #include "AliFlowCommon/AliFlowTrackSimpleCuts.h"
24
25
26 // FlowEventSimpleMaker:
27 // Class to fill the AliFlowEventSimple
28 // with AliFlowTrackSimple objects
29 // ouside the AliRoot Framework
30 // Has fill methods for TTree, 
31
32 ClassImp(FlowEventSimpleMaker)
33 //----------------------------------------------------------------------- 
34 FlowEventSimpleMaker::FlowEventSimpleMaker()
35 {
36   //constructor
37 }
38
39 //-----------------------------------------------------------------------   
40 FlowEventSimpleMaker::~FlowEventSimpleMaker()
41 {
42   //destructor
43 }
44
45 //-----------------------------------------------------------------------   
46 AliFlowEventSimple* FlowEventSimpleMaker::FillTracks(TTree* anInput, AliFlowTrackSimpleCuts* RPCuts, AliFlowTrackSimpleCuts* POICuts)
47 {
48   //fills the event from a TTree of kinematic.root files
49   
50   // number of times to use the same particle (trick to introduce nonflow)
51   Int_t iLoops = 1;
52   
53   //flags for particles passing RP and POI cuts
54   Bool_t bPassedRPCuts  = kFALSE;
55   Bool_t bPassedPOICuts = kFALSE;
56   
57   //track cut values
58   Double_t dPtMaxRP  = RPCuts->GetPtMax();
59   Double_t dPtMinRP  = RPCuts->GetPtMin();
60   Double_t dEtaMaxRP = RPCuts->GetEtaMax();
61   Double_t dEtaMinRP = RPCuts->GetEtaMin();
62   Double_t dPhiMaxRP = RPCuts->GetPhiMax();
63   Double_t dPhiMinRP = RPCuts->GetPhiMin();
64   Int_t iPIDRP       = RPCuts->GetPID();
65   
66   Double_t dPtMaxPOI  = POICuts->GetPtMax();
67   Double_t dPtMinPOI  = POICuts->GetPtMin();
68   Double_t dEtaMaxPOI = POICuts->GetEtaMax();
69   Double_t dEtaMinPOI = POICuts->GetEtaMin();
70   Double_t dPhiMaxPOI = POICuts->GetPhiMax();
71   Double_t dPhiMinPOI = POICuts->GetPhiMin();
72   Int_t iPIDPOI       = POICuts->GetPID();
73   
74   Int_t iNumberOfInputTracks = anInput->GetEntries() ;
75   //cerr<<"iNumberOfInputTracks = "<<iNumberOfInputTracks<<endl;
76   TParticle* pParticle = new TParticle();
77   anInput->SetBranchAddress("Particles",&pParticle);  
78   //  AliFlowEventSimple* pEvent = new AliFlowEventSimple(iNumberOfInputTracks);
79   AliFlowEventSimple* pEvent = new AliFlowEventSimple(10);
80   //cerr<<pEvent<<" pEvent "<<endl;
81   
82   Int_t iN = iNumberOfInputTracks; // additional variable to artificially fix the number of tracks
83   //  Int_t iN = 576; //multiplicity for chi=1.5
84   //  Int_t iN = 256; //multiplicity for chi=1
85   //  Int_t iN = 164; //multiplicity for chi=0.8
86   
87   Int_t iGoodTracks = 0;
88   Int_t itrkN = 0;
89   Int_t iSelParticlesPOI = 0;
90   Int_t iSelParticlesRP  = 0;
91   
92   while (itrkN < iNumberOfInputTracks) {
93     anInput->GetEntry(itrkN);   //get input particle
94     //checking the cuts for int. and diff. flow
95     if (pParticle->Pt() > dPtMinRP && pParticle->Pt() < dPtMaxRP &&
96         pParticle->Eta() > dEtaMinRP && pParticle->Eta() < dEtaMaxRP &&
97         pParticle->Phi() > dPhiMinRP && pParticle->Phi() < dPhiMaxRP &&
98         TMath::Abs(pParticle->GetPdgCode()) == iPIDRP) { 
99       bPassedRPCuts = kTRUE; 
100     } 
101     
102     if (pParticle->Pt() > dPtMinPOI && pParticle->Pt() < dPtMaxPOI &&
103         pParticle->Eta() > dEtaMinPOI && pParticle->Eta() < dEtaMaxPOI &&
104         pParticle->Phi() > dPhiMinPOI && pParticle->Phi() < dPhiMaxPOI &&
105         TMath::Abs(pParticle->GetPdgCode()) == iPIDPOI){ 
106       bPassedPOICuts = kTRUE; 
107     }
108     
109     if (bPassedRPCuts || bPassedPOICuts) {
110       for(Int_t d=0;d<iLoops;d++) {
111         AliFlowTrackSimple* pTrack = new AliFlowTrackSimple();
112         pTrack->SetPt(pParticle->Pt());
113         pTrack->SetEta(pParticle->Eta());
114         pTrack->SetPhi(pParticle->Phi());
115         
116         //marking the particles used for int. flow:
117         if(bPassedRPCuts && iSelParticlesRP < iN*iLoops) {  
118           //pTrack->SetForIntegratedFlow(kTRUE);
119           pTrack->SetForPRSelection(kTRUE);
120           iSelParticlesRP++;
121         }
122         //marking the particles used for diff. flow:
123         if(bPassedPOICuts) {
124           pTrack->SetForPOISelection(kTRUE);
125           iSelParticlesPOI++;
126         }
127         //adding a particles which were used either for int. or diff. flow to the list
128         pEvent->TrackCollection()->Add(pTrack);
129         iGoodTracks++;
130       }//end of for(Int_t d=0;d<iLoops;d++)
131     }//end of if(bPassedIntFlowCuts || bPassedDiffFlowCuts) 
132     itrkN++;  
133     bPassedRPCuts  = kFALSE;
134     bPassedPOICuts = kFALSE;
135   }//end of while (itrkN < iNumberOfInputTracks)
136   
137   pEvent->SetEventNSelTracksRP(iSelParticlesRP);  
138   pEvent->SetNumberOfTracks(iGoodTracks);//tracks used either for int. or for diff. flow
139
140   cout<<" iGoodTracks = "<<iGoodTracks<<endl;
141   cout<<" # of selected tracks for RP  = "<<iSelParticlesRP<<endl;
142   cout<<" # of selected tracks for POI = "<<iSelParticlesPOI<<endl;  
143
144   delete pParticle;
145   return pEvent;
146 }
147
148
149
150
151