]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HBTAN/AliHBTReaderKineTree.cxx
24f844390a2cc0715324eaeedefc90faa180d9af
[u/mrichter/AliRoot.git] / HBTAN / AliHBTReaderKineTree.cxx
1 #include "AliHBTReaderKineTree.h"
2
3 #include <Riostream.h>
4 #include <TString.h>
5 #include <TObjString.h>
6 #include <TTree.h>
7 #include <TFile.h>
8 #include <TParticle.h>
9
10 #include <AliRun.h>
11 #include <AliRunLoader.h>
12 #include <AliStack.h>
13 #include <AliHeader.h>
14
15 #include "AliHBTRun.h"
16 #include "AliHBTEvent.h"
17 #include "AliHBTParticle.h"
18 #include "AliHBTParticleCut.h"
19
20 ClassImp(AliHBTReaderKineTree)
21 /**********************************************************/
22 const TString AliHBTReaderKineTree::fgkEventFolderName("HBTReaderKineTree");
23
24 AliHBTReaderKineTree::AliHBTReaderKineTree():fFileName("galice.root")
25 {
26   fParticles = new AliHBTRun();
27   fIsRead = kFALSE;
28 }
29
30 AliHBTReaderKineTree::
31 AliHBTReaderKineTree(TString& fname):fFileName(fname)
32 {
33   fParticles = new AliHBTRun();
34   fIsRead = kFALSE;
35 }
36
37
38 /**********************************************************/
39 AliHBTReaderKineTree::
40 AliHBTReaderKineTree(TObjArray* dirs,const Char_t *filename):AliHBTReader(dirs),fFileName(filename)
41 {
42   fParticles = new AliHBTRun();
43   fIsRead = kFALSE;
44 }
45 /**********************************************************/
46
47 AliHBTEvent* AliHBTReaderKineTree::GetParticleEvent(Int_t n)
48  {
49  //returns Nth event with simulated particles
50    if (!fIsRead) 
51     if(Read(fParticles,0x0))
52      {
53        Error("GetParticleEvent","Error in reading");
54        return 0x0;
55      }
56
57    return fParticles->GetEvent(n);
58  }
59
60 /********************************************************************/
61
62 Int_t AliHBTReaderKineTree::GetNumberOfPartEvents()
63  {
64  //returns number of events of particles
65    if (!fIsRead)
66     if(Read(fParticles,0x0))
67      {
68        Error("GetNumberOfPartEvents","Error in reading");
69        return 0;
70      }
71    return fParticles->GetNumberOfEvents();
72  }
73
74
75 /**********************************************************/
76 Int_t AliHBTReaderKineTree::
77 Read(AliHBTRun* particles, AliHBTRun *tracks)
78  {
79  cout<<"AliHBTReaderKineTree::Read()"<<endl;
80  if (!particles) //check if an object is instatiated
81    {
82      Error("Read"," particles object must instatiated before passing it to the reader");
83      return 1;
84    }  
85  particles->Reset();//clear runs == delete all old events
86  
87  Int_t Ndirs;//total number of directories specified in fDirs array
88  Int_t Nevents; //number of events read in current directory
89  Int_t totalNevents = 0; //total number of read events 
90  Int_t currentdir = 0; //number of current directory name is fDirs array
91
92  if (fDirs) //if array with directories is supplied by user
93   {
94     Ndirs = fDirs->GetEntries(); //get the number if directories
95   }
96  else
97   {
98     Ndirs = 0; //if the array is not supplied read only from current directory
99   }
100
101  do  //do{}while; is OK even if 0 dirs specified. In that case we try to read from "./"
102   { 
103     cout<<"________________________________________________________\n";
104     AliRunLoader * rl = OpenFile(currentdir);
105
106     if (rl == 0x0)
107      {
108       Error("Read","Cannot open session");
109       return 2;
110      }
111     rl->LoadHeader();
112     rl->LoadKinematics("READ");
113     Nevents = rl->GetNumberOfEvents();
114     
115     for(Int_t currentEvent =0; currentEvent<Nevents;currentEvent++)//loop over all events
116      {
117 //      cout<<"processing event "<<currentEvent<<" in current dir."<<endl;
118       rl->GetEvent(currentEvent);
119       AliStack* stack = rl->Stack();
120       if (!stack)
121        {
122          Error("Read","Can not get stack for event %d",currentEvent);
123          continue;
124        }
125       Int_t npart = stack->GetNtrack();
126       Int_t nnn = 0;
127       for (Int_t i = 0;i<npart; i++)
128        {
129          
130          TParticle * p = stack->Particle(i);
131 //         if (p->GetFirstMother() >= 0) continue; do not apply with pythia etc
132          
133          if(Pass(p->GetPdgCode())) continue; //check if we are intersted with particles of this type 
134                                              //if not take next partilce
135          
136          AliHBTParticle* part = new AliHBTParticle(*p,i);
137          if(Pass(part)) { delete part; continue;}//check if meets all criteria of any of our cuts
138                                                   //if it does not delete it and take next good track
139          particles->AddParticle(totalNevents,part);//put track and particle on the run
140
141          if ( (nnn%100) == 0) 
142           {
143             cout<<nnn<<"\r";
144             fflush(0);
145           }
146          nnn++;
147        }
148       cout<<"Total read "<<nnn<<endl;
149       totalNevents++;
150      }
151      delete rl;
152      currentdir++;
153   }while(currentdir < Ndirs);//end of loop over directories specified in fDirs Obj Array
154  fIsRead = kTRUE;
155  return 0;
156 }
157
158 AliRunLoader* AliHBTReaderKineTree::OpenFile(Int_t n)
159 {
160 //opens file with kine tree
161
162  const TString& dirname = GetDirName(n);
163  if (dirname == "")
164   {
165    Error("OpenFiles","Can not get directory name");
166    return 0x0;
167   }
168  TString filename = dirname +"/"+ fFileName;
169
170  AliRunLoader *ret = AliRunLoader::Open(filename.Data(),fgkEventFolderName,"READ"); 
171
172  if ( ret == 0x0)
173   {
174     Error("OpenFiles","Can't open session from file %s",filename.Data());
175     return 0x0;
176   }
177  
178  return ret;
179 }