]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG4/PartCorrDep/AliAnaChargedParticles.cxx
Update and move to new raw reader (Sylwester)
[u/mrichter/AliRoot.git] / PWG4 / PartCorrDep / AliAnaChargedParticles.cxx
CommitLineData
477d6cee 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//_________________________________________________________________________
18//
19// Class for track selection and identification (not done now)
20// Tracks from the CTS are kept in the AOD.
21// Few histograms produced.
22//
23//-- Author: Gustavo Conesa (INFN-LNF)
24//_________________________________________________________________________
25
26
27// --- ROOT system ---
477d6cee 28#include "TParticle.h"
29#include "TH2F.h"
30
31//---- AliRoot system ----
32#include "AliAnaChargedParticles.h"
33#include "AliCaloTrackReader.h"
34#include "AliAODPWG4Particle.h"
35#include "AliStack.h"
477d6cee 36#include "AliFidutialCut.h"
37#include "AliAODTrack.h"
591cc579 38#include "AliAODMCParticle.h"
477d6cee 39
40ClassImp(AliAnaChargedParticles)
41
42//____________________________________________________________________________
43 AliAnaChargedParticles::AliAnaChargedParticles() :
44 AliAnaPartCorrBaseClass(),fPdg(0), fhPt(0),fhPhi(0),fhEta(0),
45 fhPtPion(0),fhPhiPion(0),fhEtaPion(0),
46 fhPtProton(0),fhPhiProton(0),fhEtaProton(0),
47 fhPtElectron(0),fhPhiElectron(0),fhEtaElectron(0),
48 fhPtKaon(0),fhPhiKaon(0),fhEtaKaon(0),
49 fhPtUnknown(0),fhPhiUnknown(0),fhEtaUnknown(0)
50
51{
52 //Default Ctor
53
54 //Initialize parameters
55 InitParameters();
56}
57
58//____________________________________________________________________________
59AliAnaChargedParticles::AliAnaChargedParticles(const AliAnaChargedParticles & ch) :
60 AliAnaPartCorrBaseClass(ch), fPdg(ch.fPdg), fhPt(ch.fhPt), fhPhi(ch.fhPhi),fhEta(ch.fhEta),
61 fhPtPion(ch. fhPtPion),fhPhiPion(ch.fhPhiPion),fhEtaPion(ch.fhEtaPion),
62 fhPtProton(ch.fhPtProton),fhPhiProton(ch.fhPhiProton),fhEtaProton(ch.fhEtaProton),
63 fhPtElectron(ch. fhPtElectron),fhPhiElectron(ch.fhPhiElectron),fhEtaElectron(ch.fhEtaElectron),
64 fhPtKaon(ch. fhPtKaon),fhPhiKaon(ch.fhPhiKaon),fhEtaKaon(ch.fhEtaKaon),
65 fhPtUnknown(ch.fhPtUnknown),fhPhiUnknown(ch.fhPhiUnknown),fhEtaUnknown(ch.fhEtaUnknown)
66
67{
68 // cpy ctor
69
70}
71
72//_________________________________________________________________________
73AliAnaChargedParticles & AliAnaChargedParticles::operator = (const AliAnaChargedParticles & ch)
74{
75 // assignment operator
76
77 if(this == &ch)return *this;
78 ((AliAnaPartCorrBaseClass *)this)->operator=(ch);
79
80 fPdg = ch.fPdg;
81 fhPt = ch.fhPt;
82 fhPhi = ch.fhPhi;
83 fhEta = ch.fhEta;
84
85 fhPtPion = ch. fhPtPion; fhPhiPion = ch.fhPhiPion; fhEtaPion = ch.fhEtaPion;
86 fhPtProton = ch.fhPtProton; fhPhiProton = ch.fhPhiProton; fhEtaProton = ch.fhEtaProton;
87 fhPtElectron = ch. fhPtElectron; fhPhiElectron = ch.fhPhiElectron; fhEtaElectron = ch.fhEtaElectron;
88 fhPtKaon = ch. fhPtKaon; fhPhiKaon = ch.fhPhiKaon; fhEtaKaon = ch.fhEtaKaon;
89 fhPtUnknown = ch.fhPtUnknown; fhPhiUnknown = ch.fhPhiUnknown; fhEtaUnknown = ch.fhEtaUnknown;
90
91 return *this;
92
93}
94
95//________________________________________________________________________
96TList * AliAnaChargedParticles::GetCreateOutputObjects()
97{
98 // Create histograms to be saved in output file and
99 // store them in fOutputContainer
100
101
102 TList * outputContainer = new TList() ;
103 outputContainer->SetName("ExampleHistos") ;
104
105 Int_t nptbins = GetHistoNPtBins();
106 Int_t nphibins = GetHistoNPhiBins();
107 Int_t netabins = GetHistoNEtaBins();
108 Float_t ptmax = GetHistoPtMax();
109 Float_t phimax = GetHistoPhiMax();
110 Float_t etamax = GetHistoEtaMax();
111 Float_t ptmin = GetHistoPtMin();
112 Float_t phimin = GetHistoPhiMin();
113 Float_t etamin = GetHistoEtaMin();
114
115 fhPt = new TH1F ("hPtCharged","p_T distribution", nptbins,ptmin,ptmax);
116 fhPt->SetXTitle("p_{T} (GeV/c)");
117 outputContainer->Add(fhPt);
118
119 fhPhi = new TH2F ("hPhiCharged","#phi distribution",nptbins,ptmin,ptmax, nphibins,phimin,phimax);
120 fhPhi->SetXTitle("#phi (rad)");
121 outputContainer->Add(fhPhi);
122
123 fhEta = new TH2F ("hEtaCharged","#eta distribution",nptbins,ptmin,ptmax, netabins,etamin,etamax);
124 fhEta->SetXTitle("#eta ");
125 outputContainer->Add(fhEta);
126
127
128 if(IsDataMC()){
129
a3aebfff 130 fhPtPion = new TH1F ("hPtMCPion","p_T distribution from #pi", nptbins,ptmin,ptmax);
477d6cee 131 fhPtPion->SetXTitle("p_{T} (GeV/c)");
132 outputContainer->Add(fhPtPion);
133
a3aebfff 134 fhPhiPion = new TH2F ("hPhiMCPion","#phi distribution from #pi",nptbins,ptmin,ptmax, nphibins,phimin,phimax);
477d6cee 135 fhPhiPion->SetXTitle("#phi (rad)");
136 outputContainer->Add(fhPhiPion);
137
a3aebfff 138 fhEtaPion = new TH2F ("hEtaMCPion","#eta distribution from #pi",nptbins,ptmin,ptmax, netabins,etamin,etamax);
477d6cee 139 fhEtaPion->SetXTitle("#eta ");
140 outputContainer->Add(fhEtaPion);
141
a3aebfff 142 fhPtProton = new TH1F ("hPtMCProton","p_T distribution from proton", nptbins,ptmin,ptmax);
477d6cee 143 fhPtProton->SetXTitle("p_{T} (GeV/c)");
144 outputContainer->Add(fhPtProton);
145
a3aebfff 146 fhPhiProton = new TH2F ("hPhiMCProton","#phi distribution from proton",nptbins,ptmin,ptmax, nphibins,phimin,phimax);
477d6cee 147 fhPhiProton->SetXTitle("#phi (rad)");
148 outputContainer->Add(fhPhiProton);
149
a3aebfff 150 fhEtaProton = new TH2F ("hEtaMCProton","#eta distribution from proton",nptbins,ptmin,ptmax, netabins,etamin,etamax);
477d6cee 151 fhEtaProton->SetXTitle("#eta ");
152 outputContainer->Add(fhEtaProton);
153
a3aebfff 154 fhPtKaon = new TH1F ("hPtMCKaon","p_T distribution from kaon", nptbins,ptmin,ptmax);
477d6cee 155 fhPtKaon->SetXTitle("p_{T} (GeV/c)");
156 outputContainer->Add(fhPtKaon);
157
a3aebfff 158 fhPhiKaon = new TH2F ("hPhiMCKaon","#phi distribution from kaon",nptbins,ptmin,ptmax, nphibins,phimin,phimax);
477d6cee 159 fhPhiKaon->SetXTitle("#phi (rad)");
160 outputContainer->Add(fhPhiKaon);
161
a3aebfff 162 fhEtaKaon = new TH2F ("hEtaMCKaon","#eta distribution from kaon",nptbins,ptmin,ptmax, netabins,etamin,etamax);
477d6cee 163 fhEtaKaon->SetXTitle("#eta ");
164 outputContainer->Add(fhEtaKaon);
165
a3aebfff 166 fhPtElectron = new TH1F ("hPtMCElectron","p_T distribution from electron", nptbins,ptmin,ptmax);
477d6cee 167 fhPtElectron->SetXTitle("p_{T} (GeV/c)");
168 outputContainer->Add(fhPtElectron);
169
a3aebfff 170 fhPhiElectron = new TH2F ("hPhiMCElectron","#phi distribution from electron",nptbins,ptmin,ptmax, nphibins,phimin,phimax);
477d6cee 171 fhPhiElectron->SetXTitle("#phi (rad)");
172 outputContainer->Add(fhPhiElectron);
173
a3aebfff 174 fhEtaElectron = new TH2F ("hEtaMCElectron","#eta distribution from electron",nptbins,ptmin,ptmax, netabins,etamin,etamax);
477d6cee 175 fhEtaElectron->SetXTitle("#eta ");
176 outputContainer->Add(fhEtaElectron);
177
a3aebfff 178 fhPtUnknown = new TH1F ("hPtMCUnknown","p_T distribution from unknown", nptbins,ptmin,ptmax);
477d6cee 179 fhPtUnknown->SetXTitle("p_{T} (GeV/c)");
180 outputContainer->Add(fhPtUnknown);
181
a3aebfff 182 fhPhiUnknown = new TH2F ("hPhiMCUnknown","#phi distribution from unknown",nptbins,ptmin,ptmax, nphibins,phimin,phimax);
477d6cee 183 fhPhiUnknown->SetXTitle("#phi (rad)");
184 outputContainer->Add(fhPhiUnknown);
185
a3aebfff 186 fhEtaUnknown = new TH2F ("hEtaMCUnknown","#eta distribution from unknown",nptbins,ptmin,ptmax, netabins,etamin,etamax);
477d6cee 187 fhEtaUnknown->SetXTitle("#eta ");
188 outputContainer->Add(fhEtaUnknown);
189
190 }
a3aebfff 191
477d6cee 192 return outputContainer;
a3aebfff 193
477d6cee 194}
195
196//__________________________________________________
197void AliAnaChargedParticles::InitParameters()
198{
199 //Initialize the parameters of the analysis.
200 SetOutputAODClassName("AliAODPWG4Particle");
a3aebfff 201 SetOutputAODName("PWG4Particle");
202
203 AddToHistogramsName("AnaCharged_");
204
477d6cee 205 fPdg = -1; //Select all tracks
206
207}
208
209//__________________________________________________________________
210void AliAnaChargedParticles::Print(const Option_t * opt) const
211{
212 //Print some relevant parameters set for the analysis
213 if(! opt)
214 return;
215
a3aebfff 216 printf("**** Print %s %s ****\n", GetName(), GetTitle() ) ;
217 AliAnaPartCorrBaseClass::Print(" ");
218
477d6cee 219 printf("Min Pt = %3.2f\n", GetMinPt());
220 printf("Max Pt = %3.2f\n", GetMaxPt());
221 printf("Select clusters with pdg %d \n",fPdg);
222
223}
224
225//____________________________________________________________________________
226void AliAnaChargedParticles::Init()
227{
228 //Init
229 //Do some checks
230 if(!GetReader()->IsCTSSwitchedOn()){
7cd4e982 231 printf("AliAnaChargedParticles::Init() - STOP!: You want to use CTS tracks in analysis but not read!! \n!!Check the configuration file!!\n");
477d6cee 232 abort();
233 }
234
235}
236
237//__________________________________________________________________
238void AliAnaChargedParticles::MakeAnalysisFillAOD()
239{
240 //Do analysis and fill aods
241 if(!GetAODCTS() || GetAODCTS()->GetEntriesFast() == 0) return ;
242 Int_t ntracks = GetAODCTS()->GetEntriesFast();
243
244 //Some prints
245 if(GetDebug() > 0)
a3aebfff 246 printf("AliAnaChargedParticles::MakeAnalysisFillAOD() - In CTS aod entries %d\n", ntracks);
477d6cee 247
248 //Fill AODParticle with CTS aods
249 TVector3 p3;
250 for(Int_t i = 0; i < ntracks; i++){
251
252 AliAODTrack * track = (AliAODTrack*) (GetAODCTS()->At(i));
253
254 //Fill AODParticle after some selection
255 Double_t mom[3] = {track->Px(),track->Py(),track->Pz()};
256 p3.SetXYZ(mom[0],mom[1],mom[2]);
257
258 //Acceptance selection
259 Bool_t in = GetFidutialCut()->IsInFidutialCut(mom,"CTS") ;
a3aebfff 260 if(GetDebug() > 1) printf("AliAnaChargedParticles::MakeAnalysisFillAOD() - Track pt %2.2f, phi %2.2f, in fidutial cut %d\n",p3.Pt(), p3.Phi(), in);
477d6cee 261 if(p3.Pt() > GetMinPt() && in) {
262 //Keep only particles identified with fPdg
263 //Selection not done for the moment
264 //Should be done here.
265
266 AliAODPWG4Particle tr = AliAODPWG4Particle(mom[0],mom[1],mom[2],0);
267 tr.SetDetector("CTS");
268 tr.SetLabel(track->GetLabel());
269 tr.SetTrackLabel(track->GetID(),-1);
7cd4e982 270 //Input from second AOD?
271 if(GetReader()->GetAODCTSNormalInputEntries() <= i) tr.SetInputFileIndex(1);
272
477d6cee 273 AddAODParticle(tr);
274 }//selection
275 }//loop
276
277 if(GetDebug() > 0)
a3aebfff 278 printf("AliAnaChargedParticles::MakeAnalysisFillAOD() - Final aod branch entries %d\n", GetOutputAODBranch()->GetEntriesFast());
477d6cee 279}
280
281//__________________________________________________________________
282void AliAnaChargedParticles::MakeAnalysisFillHistograms()
283{
284 //Do analysis and fill histograms
285
286 //Loop on stored AODParticles
287 Int_t naod = GetOutputAODBranch()->GetEntriesFast();
a3aebfff 288 if(GetDebug() > 0) printf("AliAnaChargedParticles::MakeAnalysisFillHistograms() - aod branch entries %d\n", naod);
477d6cee 289 for(Int_t iaod = 0; iaod < naod ; iaod++){
290 AliAODPWG4Particle* tr = (AliAODPWG4Particle*) (GetOutputAODBranch()->At(iaod));
291
292 fhPt->Fill(tr->Pt());
293 fhPhi->Fill(tr->Pt(), tr->Phi());
294 fhEta->Fill(tr->Pt(), tr->Eta());
295
296 if(IsDataMC()){
297 //Play with the MC stack if available
591cc579 298 Int_t mompdg = -1;
7cd4e982 299 Int_t label = tr->GetLabel();
591cc579 300 if(GetReader()->ReadStack()){
301 TParticle * mom = GetMCStack()->Particle(label);
302 mompdg =TMath::Abs(mom->GetPdgCode());
303 }
7cd4e982 304 else if(GetReader()->ReadAODMCParticles()){
591cc579 305 AliAODMCParticle * aodmom = 0;
306 //Get the list of MC particles
7cd4e982 307 aodmom = (AliAODMCParticle*) (GetReader()->GetAODMCParticles(tr->GetInputFileIndex()))->At(label);
308 mompdg =TMath::Abs(aodmom->GetPdgCode());
309 }
591cc579 310
7cd4e982 311 if(mompdg==211){
477d6cee 312 fhPtPion->Fill(tr->Pt());
313 fhPhiPion->Fill(tr->Pt(), tr->Phi());
314 fhEtaPion->Fill(tr->Pt(), tr->Eta());
315 }
316 else if(mompdg==2212){
317 fhPtProton->Fill(tr->Pt());
318 fhPhiProton->Fill(tr->Pt(), tr->Phi());
319 fhEtaProton->Fill(tr->Pt(), tr->Eta());
320 }
321 else if(mompdg==321){
322 fhPtKaon->Fill(tr->Pt());
323 fhPhiKaon->Fill(tr->Pt(), tr->Phi());
324 fhEtaKaon->Fill(tr->Pt(), tr->Eta());
325 }
326 else if(mompdg==11){
327 fhPtElectron->Fill(tr->Pt());
328 fhPhiElectron->Fill(tr->Pt(), tr->Phi());
329 fhEtaElectron->Fill(tr->Pt(), tr->Eta());
330 }
331 else {
332 //printf("unknown pdg %d\n",mompdg);
333 fhPtUnknown->Fill(tr->Pt());
334 fhPhiUnknown->Fill(tr->Pt(), tr->Phi());
335 fhEtaUnknown->Fill(tr->Pt(), tr->Eta());
336 }
337 }//Work with stack also
338 }// aod branch loop
339
340}