]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/AliReaderESDTree.cxx
Set Probabilities to zero if there is no signal in any plane
[u/mrichter/AliRoot.git] / ANALYSIS / AliReaderESDTree.cxx
1 #include "AliReaderESDTree.h"
2 //_______________________________________________________________________
3 /////////////////////////////////////////////////////////////////////////
4 //
5 // class AliReaderESDTree
6 //
7 // Reader for MUON ESD Tree (only for rec)
8 //
9 // finck@subatech.in2p3.fr
10 //
11 /////////////////////////////////////////////////////////////////////////
12
13 #include <TString.h>
14 #include <TTree.h>
15 #include <TFile.h>
16
17
18 #include <AliRun.h>
19 #include <AliRunLoader.h>
20
21 #include <AliESD.h>
22 #include "AliAOD.h"
23
24 ClassImp(AliReaderESDTree)
25
26 AliReaderESDTree::AliReaderESDTree(const Char_t* esdfilename, const Char_t* galfilename):
27   AliReaderESD(esdfilename,galfilename),
28   fTree(0x0)
29 {
30 //ctor
31 }
32
33 /********************************************************************/
34 AliReaderESDTree::~AliReaderESDTree()
35 {
36 //dtor 
37  delete fTree;
38 }
39
40 /**********************************************************/
41 Int_t AliReaderESDTree::ReadNext()
42 {
43 //reads next event from fFile
44 //fRunLoader is for reading Kine
45   
46   if (AliVAODParticle::GetDebug())
47     Info("ReadNext","Entered");
48     
49   if (fEventSim == 0x0)  fEventSim = new AliAOD();
50   if (fEventRec == 0x0)  fEventRec = new AliAOD();
51   
52   fEventSim->Reset();
53   fEventRec->Reset();
54         
55   do  //do{}while; is OK even if 0 dirs specified. In that case we try to read from "./"
56     {
57       if (fFile == 0x0)
58         {
59           fFile = OpenFile(fCurrentDir);//rl is opened here
60           if (fFile == 0x0)
61             {
62               Error("ReadNext","Cannot get fFile for dir no. %d",fCurrentDir);
63               fCurrentDir++;
64               continue;
65             }
66           fCurrentEvent = 0;
67         }
68
69       static AliESD* esd = 0x0;
70       fTree->SetBranchAddress("ESD", &esd);
71       Int_t status = fTree->GetEvent(fCurrentEvent);
72
73       if (!status)
74         {
75           if (AliVAODParticle::GetDebug() > 2 )
76             {
77               Info("ReadNext","Can not find event# %d in Tree", fCurrentEvent);
78             }
79           fCurrentDir++;
80           delete fTree;
81           fTree = 0x0;
82           delete fFile;//we have to assume there is no more ESD objects in the fFile
83           fFile = 0x0;
84           delete fRunLoader;
85           fRunLoader = 0x0;
86           continue;
87         }
88
89       ReadESD(esd);
90       
91       fCurrentEvent++;
92       fNEventsRead++;
93       return 0;//success -> read one event
94     }while(fCurrentDir < GetNumberOfDirs());//end of loop over directories specified in fDirs Obj Array  
95    
96   return 1; //no more directories to read
97 }
98
99 /**********************************************************/
100 TFile* AliReaderESDTree::OpenFile(Int_t n)
101 {
102 //opens fFile with kine tree
103
104  const TString& dirname = GetDirName(n);
105  if (dirname == "")
106   {
107    Error("OpenFiles","Can not get directory name");
108    return 0x0;
109   }
110  TString filename = dirname +"/"+ fESDFileName;
111  TFile *ret = TFile::Open(filename.Data()); 
112
113  if (ret == 0x0)
114   {
115     Error("OpenFiles","Can't open fFile %s",filename.Data());
116     return 0x0;
117   }
118  if (!ret->IsOpen())
119   {
120     Error("OpenFiles","Can't open fFile  %s",filename.Data());
121     return 0x0;
122   }
123  
124  TString esdname = "esdTree";
125  fTree = dynamic_cast<TTree*> (ret->Get(esdname));
126
127  if (!fTree)
128   {
129     Error("OpenFiles","Can't open ESD Tree %s",esdname.Data());
130     delete ret;
131     return 0x0;
132
133   }
134  
135  if (fReadSim )
136   {
137    fRunLoader = AliRunLoader::Open(dirname +"/"+ fGAlFileName);
138    if (fRunLoader == 0x0)
139     {
140       Error("OpenFiles","Can't get RunLoader for directory %s",dirname.Data());
141       delete fTree;
142       delete ret;
143       return 0x0;
144     }
145     
146    fRunLoader->LoadHeader();
147    if (fRunLoader->LoadKinematics())
148     {
149       Error("Next","Error occured while loading kinematics.");
150       delete fRunLoader;
151       delete fTree;
152       delete ret;
153       return 0x0;
154     }
155   }
156    
157  return ret;
158 }