]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - JETAN/AliJetESDReader.cxx
EffC++ warnings corrected.
[u/mrichter/AliRoot.git] / JETAN / AliJetESDReader.cxx
... / ...
CommitLineData
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
16//-------------------------------------------------------------------------
17// Jet ESD Reader
18// ESD reader for jet analysis
19// Author: Mercedes Lopez Noriega (mercedes.lopez.noriega@cern.ch)
20//-------------------------------------------------------------------------
21
22
23#include <Riostream.h>
24#include <TSystem.h>
25#include <TLorentzVector.h>
26#include <TVector3.h>
27#include "AliJetESDReader.h"
28#include "AliJetESDReaderHeader.h"
29#include "AliESD.h"
30#include "AliESDtrack.h"
31
32ClassImp(AliJetESDReader)
33
34 AliJetESDReader::AliJetESDReader()
35{
36 // Constructor
37 fReaderHeader = 0x0;
38 fMass = 0;
39 fSign = 0;
40}
41
42//____________________________________________________________________________
43
44AliJetESDReader::~AliJetESDReader()
45{
46 // Destructor
47}
48
49//____________________________________________________________________________
50
51void AliJetESDReader::OpenInputFiles()
52
53{
54 // chain for the ESDs
55 fChain = new TChain("esdTree");
56
57 // get directory and pattern name from the header
58 const char* dirName=fReaderHeader->GetDirectory();
59 const char* pattern=fReaderHeader->GetPattern();
60
61 // Add files matching patters to the chain
62 void *dir = gSystem->OpenDirectory(dirName);
63 const char *name = 0x0;
64 int nesd = fReaderHeader->GetNesd();
65 int a = 0;
66 while ((name = gSystem->GetDirEntry(dir))){
67 if (a>=nesd) continue;
68 if (strstr(name,pattern)){
69 printf("Adding %s\n",name);
70 char path[256];
71 sprintf(path,"%s/%s",dirName,name);
72 fChain->AddFile(path,-1);
73 a++;
74 }
75 }
76 printf("%d ESDs added\n",a);
77
78 gSystem->FreeDirectory(dir);
79 fChain->SetBranchAddress("ESD", &fESD);
80
81 int nMax = fChain->GetEntries();
82 printf("\nTotal number of events in chain= %d",nMax);
83
84 // set number of events in header
85 if (fReaderHeader->GetLastEvent() == -1)
86 fReaderHeader->SetLastEvent(nMax);
87 else {
88 Int_t nUsr = fReaderHeader->GetLastEvent();
89 fReaderHeader->SetLastEvent(TMath::Min(nMax,nUsr));
90 }
91}
92
93//____________________________________________________________________________
94
95void AliJetESDReader::FillMomentumArray(Int_t event)
96{
97 // Fill momentum array
98
99 Int_t goodTrack = 0;
100 Int_t nt = 0;
101 Float_t pt, eta;
102 TVector3 p3;
103
104 // clear momentum array
105 ClearArray();
106
107 // get event from chain
108 fChain->GetEntry(event);
109
110 // get number of tracks in event (for the loop)
111 nt = fESD->GetNumberOfTracks();
112
113 // temporary storage of signal and pt cut flag
114 Int_t* sflag = new Int_t[nt];
115 Int_t* cflag = new Int_t[nt];
116
117 // get cuts set by user
118 Float_t ptMin = fReaderHeader->GetPtCut();
119 Float_t etaMin = fReaderHeader->GetFiducialEtaMin();
120 Float_t etaMax = fReaderHeader->GetFiducialEtaMax();
121
122 //loop over tracks
123 for (Int_t it = 0; it < nt; it++) {
124 AliESDtrack *track = fESD->GetTrack(it);
125 UInt_t status = track->GetStatus();
126
127 Double_t mom[3];
128 track->GetPxPyPz(mom);
129 p3.SetXYZ(mom[0],mom[1],mom[2]);
130 pt = p3.Pt();
131 if (((status & AliESDtrack::kITSrefit) == 0) ||
132 ((status & AliESDtrack::kTPCrefit) == 0)) continue; // quality check
133 if (((AliJetESDReaderHeader*) fReaderHeader)->ReadSignalOnly()
134 && TMath::Abs(track->GetLabel()) > 10000) continue; // quality check
135 if (((AliJetESDReaderHeader*) fReaderHeader)->ReadBkgdOnly()
136 && TMath::Abs(track->GetLabel()) < 10000) continue; // quality check
137 eta = p3.Eta();
138 if ( (eta > etaMax) || (eta < etaMin)) continue; // checking eta cut
139
140 new ((*fMomentumArray)[goodTrack]) TLorentzVector(p3,p3.Mag());
141 sflag[goodTrack]=0;
142 if (TMath::Abs(track->GetLabel()) < 10000) sflag[goodTrack]=1;
143 cflag[goodTrack]=0;
144 if (pt > ptMin) cflag[goodTrack]=1; // pt cut
145 goodTrack++;
146 }
147 // set the signal flags
148 fSignalFlag.Set(goodTrack,sflag);
149 fCutFlag.Set(goodTrack,cflag);
150}
151
152