]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HBTAN/AliHBTReaderKineTree.cxx
Track Points Added. Class needed for Anti-Merging cut
[u/mrichter/AliRoot.git] / HBTAN / AliHBTReaderKineTree.cxx
1 #include "AliHBTReaderKineTree.h"
2
3 #include <TString.h>
4 #include <TObjString.h>
5 #include <TTree.h>
6 #include <TFile.h>
7 #include <TParticle.h>
8
9 #include <AliRun.h>
10 #include <AliRunLoader.h>
11 #include <AliStack.h>
12 #include <AliHeader.h>
13
14 #include "AliHBTRun.h"
15 #include "AliHBTEvent.h"
16 #include "AliHBTParticle.h"
17 #include "AliHBTParticleCut.h"
18
19 ClassImp(AliHBTReaderKineTree)
20 /**********************************************************/
21 const TString AliHBTReaderKineTree::fgkEventFolderName("HBTReaderKineTree");
22
23 AliHBTReaderKineTree::AliHBTReaderKineTree():
24  fFileName("galice.root"),
25  fRunLoader(0x0)
26 {
27   //ctor
28 }
29 /**********************************************************/
30
31 AliHBTReaderKineTree::AliHBTReaderKineTree(TString& fname):
32  fFileName(fname),
33  fRunLoader(0x0)
34 {
35   //ctor
36 }
37 /**********************************************************/
38
39 AliHBTReaderKineTree::AliHBTReaderKineTree(TObjArray* dirs,const Char_t *filename):
40  AliHBTReader(dirs),
41  fFileName(filename),
42  fRunLoader(0x0)
43 {
44   //ctor
45 }
46 /**********************************************************/
47
48 AliHBTReaderKineTree::~AliHBTReaderKineTree()
49 {
50   //dtor
51   delete fRunLoader;
52 }
53 /**********************************************************/
54
55 void AliHBTReaderKineTree::Rewind()
56 {
57   delete fRunLoader;
58   fRunLoader = 0x0;
59   fCurrentDir = 0;
60   fNEventsRead= 0;  
61 }
62 /**********************************************************/
63
64 Int_t AliHBTReaderKineTree::ReadNext()
65 {
66  //Reads Kinematics Tree
67   
68  Info("Read","");
69  if (fParticlesEvent == 0x0)  fParticlesEvent = new AliHBTEvent();
70  fParticlesEvent->Reset();
71
72  do  //do{}while; is OK even if 0 dirs specified. In that case we try to read from "./"
73   { 
74     if (fRunLoader == 0x0) 
75       if (OpenNextFile()) 
76         { 
77           fCurrentDir++;
78           continue;
79         }
80     
81     if (fCurrentEvent == fRunLoader->GetNumberOfEvents())
82      {
83        //read next directory
84        delete fRunLoader;//close current session
85        fRunLoader = 0x0;//assure pointer is null
86        fCurrentDir++;//go to next dir
87        continue;//directory counter is increased inside in case of error
88      }
89      
90     Info("ReadNext","Reading Event %d",fCurrentEvent);
91     
92     fRunLoader->GetEvent(fCurrentEvent);
93     
94     AliStack* stack = fRunLoader->Stack();
95     if (!stack)
96      {
97        Error("ReadNext","Can not get stack for event %d",fCurrentEvent);
98        continue;
99      }
100     Int_t npart = stack->GetNtrack();
101     for (Int_t i = 0;i<npart; i++)
102       {
103          TParticle * p = stack->Particle(i);
104 //         if (p->GetFirstMother() >= 0) continue; do not apply with pythia etc
105          
106          if(Pass(p->GetPdgCode())) continue; //check if we are intersted with particles of this type 
107                                              //if not take next partilce
108          
109          AliHBTParticle* part = new AliHBTParticle(*p,i);
110          if(Pass(part)) { delete part; continue;}//check if meets all criteria of any of our cuts
111                                                   //if it does not delete it and take next good track
112          fParticlesEvent->AddParticle(part);//put particle in event
113       }
114     Info("ReadNext","Read %d particles from event %d (event %d in dir %d).",
115                      fParticlesEvent->GetNumberOfParticles(),
116                      fNEventsRead,fCurrentEvent,fCurrentDir);
117       
118     fCurrentEvent++;
119     fNEventsRead++;
120     return 0;
121   }while(fCurrentDir < GetNumberOfDirs());//end of loop over directories specified in fDirs Obj Array
122   
123  return 1;
124 }
125 /**********************************************************/
126
127 Int_t AliHBTReaderKineTree::OpenNextFile()
128 {
129 //opens file with kine tree
130  Info("OpenNextFile","________________________________________________________");
131  
132  const TString& dirname = GetDirName(fCurrentEvent);
133  if (dirname == "")
134   {
135    Error("OpenNextFile","Can not get directory name");
136    return 1;
137   }
138  TString filename = dirname +"/"+ fFileName;
139
140  fRunLoader = AliRunLoader::Open(filename.Data(),fgkEventFolderName,"READ"); 
141
142  if ( fRunLoader == 0x0)
143   {
144     Error("OpenNextFile","Can't open session from file %s",filename.Data());
145     return 0x0;
146   }
147   
148  if (fRunLoader->GetNumberOfEvents() <= 0)
149   {
150     Error("OpenNextFile","There is no events in this directory.");
151     delete fRunLoader;
152     fRunLoader = 0x0;
153     return 1;
154   }
155   
156  if (fRunLoader->LoadKinematics())
157   {
158     Error("OpenNextFile","Error occured while loading kinematics.");
159     return 1;
160   }
161   
162  fCurrentEvent = 0;
163  return 0;
164 }