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