]> git.uio.no Git - u/mrichter/AliRoot.git/blame - JETAN/AliJetESDReader.cxx
Removing useless const to avoid warnings on alphacxx6
[u/mrichter/AliRoot.git] / JETAN / AliJetESDReader.cxx
CommitLineData
99e5fe42 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
20// mercedes.lopez.noriega@cern.ch
21//
22
23
24#include <Riostream.h>
25#include <TSystem.h>
26#include <TLorentzVector.h>
27#include "AliJetESDReader.h"
28#include "AliJetESDReaderHeader.h"
29#include "AliESD.h"
30#include "AliESDtrack.h"
31
32ClassImp(AliJetESDReader)
33
34AliJetESDReader::AliJetESDReader()
35{
36 // Constructor
37
38 fReaderHeader = 0x0;
39 fMass = 0;
40 fSign = 0;
41}
42
43//____________________________________________________________________________
44
45AliJetESDReader::~AliJetESDReader()
46{
47 // Destructor
48 delete fReaderHeader;
49}
50
51//____________________________________________________________________________
52
53void AliJetESDReader::OpenInputFiles()
54
55{
56 // chain for the ESDs
57 fChain = new TChain("esdTree");
58 fChainMC = new TChain("mcStackTree");
59
60 // get directory and pattern name from the header
61 const char* dirName=fReaderHeader->GetDirectory();
62 const char* pattern=fReaderHeader->GetPattern();
63
64 // Add files matching patters to the chain
65 void *dir = gSystem->OpenDirectory(dirName);
66 const char *name = 0x0;
67 while ((name = gSystem->GetDirEntry(dir))){
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,"esdTree");
73 fChainMC->AddFile(path,-1,"mcStackTree");
74 }
75 }
76
77 gSystem ->FreeDirectory(dir);
78 fChain ->SetBranchAddress("ESD", &fESD);
79 fChainMC->SetBranchAddress("Header", &fAliHeader);
80 fChainMC->SetBranchAddress("Stack", &fArrayMC);
81
82 int nMax = fChain->GetEntries();
83 printf("\nTotal number of events in chain= %d",nMax);
84
85 // set number of events in header
86 if (fReaderHeader->GetLastEvent() == -1)
87 fReaderHeader->SetLastEvent(nMax);
88 else {
89 Int_t nUsr = fReaderHeader->GetLastEvent();
90 fReaderHeader->SetLastEvent(TMath::Min(nMax,nUsr));
91 }
92}
93
94//____________________________________________________________________________
95
96void AliJetESDReader::FillMomentumArray(Int_t event)
97{
98// Fills the momentum array
99 Int_t goodTrack = 0;
100 Int_t nt = 0;
101 Float_t pt;
102 Double_t p[3]; // track 3 momentum vector
103
104 // clear momentum array
105 ClearArray();
106 // get event from chain
107 fChain->GetEntry(event);
108 fChainMC->GetEntry(event);
109 // get number of tracks in event (for the loop)
110 nt = fESD->GetNumberOfTracks();
111 // get cuts set by user
112 Float_t ptMin = ((AliJetESDReaderHeader*) fReaderHeader)->GetPtCut();
113
114 //loop over tracks
115 for (Int_t it = 0; it < nt; it++) {
116 AliESDtrack *track = fESD->GetTrack(it);
117 UInt_t status = track->GetStatus();
118 if ((status & AliESDtrack::kITSrefit) == 0) continue; // quality check
119
120 if (((AliJetESDReaderHeader*) fReaderHeader)->ReadSignalOnly() && TMath::Abs(track->GetLabel()) > 10000) continue;
121
122 track->GetPxPyPz(p);
123 pt = TMath::Sqrt(p[0] * p[0] + p[1] * p[1]); // pt of the track
124
125 if (pt < ptMin) continue; //check cuts
126
127 new ((*fMomentumArray)[goodTrack])
128 TLorentzVector(p[0], p[1], p[2],
129 TMath::Sqrt(pt * pt +p[2] * p[2]));
130 goodTrack++;
131 }
132
133 printf("\nIn event %d, number of good tracks %d \n", event, goodTrack);
134}
135
136