]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HBTAN/AliHBTReaderTPC.cxx
cuts on Q out, side, long added
[u/mrichter/AliRoot.git] / HBTAN / AliHBTReaderTPC.cxx
1 #include "AliHBTReaderTPC.h"
2
3 #include <iostream.h>
4 //#include <fstream.h>
5 #include <TTree.h>
6 #include <TFile.h>
7 #include <TParticle.h>
8
9 #include <AliRun.h>
10 #include <AliMagF.h>
11 #include <AliTPCtrack.h>
12 #include <AliTPCParam.h>
13 #include <AliTPCtracker.h>
14
15 #include "AliHBTRun.h"
16 #include "AliHBTEvent.h"
17 #include "AliHBTParticle.h"
18 #include "AliHBTParticleCut.h"
19
20
21 ClassImp(AliHBTReaderTPC)
22 //reader for TPC tracking
23 //needs galice.root, AliTPCtracks.root, AliTPCclusters.root, good_tracks_tpc 
24 //
25 //more info: http://alisoft.cern.ch/people/skowron/analyzer/index.html
26 //Piotr.Skowronski@cern.ch
27
28 AliHBTReaderTPC::
29  AliHBTReaderTPC(const Char_t* trackfilename,const Char_t* clusterfilename,
30                  const Char_t* galicefilename):
31                  fTrackFileName(trackfilename),fClusterFileName(clusterfilename),
32                  fGAliceFileName(galicefilename)
33 {
34   //constructor, 
35   //Defaults:
36   //  trackfilename = "AliTPCtracks.root"
37   //  clusterfilename = "AliTPCclusters.root"
38   //  galicefilename = ""  - this means: Do not open gAlice file - 
39   //                         just leave the global pointer untached
40   
41   fParticles = new AliHBTRun();
42   fTracks    = new AliHBTRun();
43   fIsRead = kFALSE;
44 }
45 /********************************************************************/
46 AliHBTReaderTPC::
47 AliHBTReaderTPC(TObjArray* dirs,
48                   const Char_t* trackfilename, const Char_t* clusterfilename,
49                   const Char_t* galicefilename):
50                   AliHBTReader(dirs), fTrackFileName(trackfilename),
51                   fClusterFileName(clusterfilename),fGAliceFileName(galicefilename)
52
53 {
54   //constructor, 
55   //Defaults:
56   //  trackfilename = "AliTPCtracks.root"
57   //  clusterfilename = "AliTPCclusters.root"
58   //  galicefilename = ""  - this means: Do not open gAlice file - 
59   //                         just leave the global pointer untached
60   
61   fParticles = new AliHBTRun();
62   fTracks    = new AliHBTRun();
63   fIsRead = kFALSE;
64   
65 }
66 /********************************************************************/
67
68 AliHBTReaderTPC::~AliHBTReaderTPC()
69  {
70  //desctructor
71    delete fParticles;
72    delete fTracks;
73  }
74 /********************************************************************/
75
76 AliHBTEvent* AliHBTReaderTPC::GetParticleEvent(Int_t n)
77  {
78  //returns Nth event with simulated particles
79    if (!fIsRead) 
80     if(Read(fParticles,fTracks))
81      {
82        Error("GetParticleEvent","Error in reading");
83        return 0x0;
84      }
85    return fParticles->GetEvent(n);
86  }
87 /********************************************************************/
88 AliHBTEvent* AliHBTReaderTPC::GetTrackEvent(Int_t n)
89  {
90  //returns Nth event with reconstructed tracks
91    if (!fIsRead) 
92     if(Read(fParticles,fTracks))
93      {
94        Error("GetTrackEvent","Error in reading");
95        return 0x0;
96      }
97    return fTracks->GetEvent(n);
98  }
99 /********************************************************************/
100
101 Int_t AliHBTReaderTPC::GetNumberOfPartEvents()
102  {
103  //returns number of events of particles
104    if (!fIsRead) 
105     if ( Read(fParticles,fTracks))
106      {
107        Error("GetNumberOfPartEvents","Error in reading");
108        return 0;
109      }
110    return fParticles->GetNumberOfEvents();
111  }
112
113 /********************************************************************/
114 Int_t AliHBTReaderTPC::GetNumberOfTrackEvents()
115  {
116  //returns number of events of tracks
117   if (!fIsRead)
118     if(Read(fParticles,fTracks))
119      {
120        Error("GetNumberOfTrackEvents","Error in reading");
121        return 0;
122      }
123   return fTracks->GetNumberOfEvents();
124  }
125 /********************************************************************/
126
127
128 Int_t AliHBTReaderTPC::Read(AliHBTRun* particles, AliHBTRun *tracks)
129  {
130  //reads data and puts put to the particles and tracks objects
131  //reurns 0 if everything is OK
132  //
133   cout<<"AliHBTReaderTPC::Read()"<<endl;
134   Int_t i; //iterator and some temprary values
135   Int_t Nevents = 0;
136   Int_t totalNevents = 0;
137   TFile *aTracksFile;//file with tracks
138   TFile *aClustersFile;//file with clusters
139   TFile *aGAliceFile;//!ile name with galice
140
141   if (!particles) //check if an object is instatiated
142    {
143      Error("Read"," particles object must instatiated before passing it to the reader");
144    }
145   if (!tracks)  //check if an object is instatiated
146    {
147      Error("Read"," tracks object must instatiated before passing it to the reader");
148    }
149   particles->Reset();//clear runs == delete all old events
150   tracks->Reset();
151  
152   TObjArray *tarray = new TObjArray(5000); //cotainer for tpc tracks
153   tarray->SetOwner(); //set the ownership of the objects it contains
154                       //when array is is deleted or cleared all objects 
155                       //that it contains are deleted
156   Int_t currentdir = 0;
157
158   Int_t Ndirs;
159   if (fDirs) //if array with directories is supplied by user
160    {
161      Ndirs = fDirs->GetEntries(); //get the number if directories
162    }
163   else
164    {
165      Ndirs = 0; //if the array is not supplied read only from current directory
166    }
167
168   do  //do{}while; is OK even if 0 dirs specified. In that case we try to read from "./"
169    {
170     
171     if( (i=OpenFiles(aTracksFile,aClustersFile,aGAliceFile,currentdir)) )
172      {
173        Error("Read","Exiting due to problems with opening files. Errorcode %d",i);
174        return i;
175      }
176   
177     
178     if (gAlice->TreeE())//check if tree E exists
179      {
180       Nevents = (Int_t)gAlice->TreeE()->GetEntries();//if yes get number of events in gAlice
181       cout<<"________________________________________________________\n";
182       cout<<"Found "<<Nevents<<" event(s) in directory "<<GetDirName(currentdir)<<endl;
183       cout<<"Setting Magnetic Field. Factor is "<<gAlice->Field()->Factor()<<endl;
184       AliKalmanTrack::SetConvConst(100/0.299792458/0.2/gAlice->Field()->Factor());  
185      }
186     else
187      {//if not return an error
188        Error("Read","Can not find Header tree (TreeE) in gAlice");
189        return 1;
190      }
191   
192     aClustersFile->cd();//set cluster file active 
193     AliTPCParam *TPCParam= (AliTPCParam*)aClustersFile->Get("75x40_100x60");
194     if (!TPCParam) 
195       { 
196        Error("Read","TPC parameters have not been found !\n");
197        return 1;
198       }
199
200   
201     for(Int_t currentEvent =0; currentEvent<Nevents;currentEvent++)//loop over all events
202      {
203        cout<<"Reading Event "<<currentEvent<<endl;
204        /**************************************/
205         /**************************************/
206          /**************************************/ 
207          
208          aTracksFile->cd();//set track file active
209           
210          Char_t  treename[100];
211          sprintf(treename,"TreeT_TPC_%d",currentEvent);//prepare name of the tree
212    
213          TTree *tracktree=0;
214          
215          tracktree=(TTree*)aTracksFile->Get(treename);//get the tree 
216          if (!tracktree) //check if we got the tree
217           {//if not return with error
218             Error("Read","Can't get a tree with TPC tracks !\n"); 
219             
220             return 1;
221           }
222    
223          TBranch *trackbranch=tracktree->GetBranch("tracks");//get the branch with tracks
224          if (!trackbranch) ////check if we got the branch
225           {//if not return with error
226             Error("Read","Can't get a branch with TPC tracks !\n"); 
227             return 2;
228           }
229          Int_t NTPCtracks=(Int_t)tracktree->GetEntries();//get number of TPC tracks 
230          cout<<"Found "<<NTPCtracks<<" TPC tracks.\n";
231          //Copy tracks to array
232          
233          AliTPCtrack *iotrack=0;
234          
235          aClustersFile->cd();//set cluster file active 
236          AliTPCtracker *tracker = new AliTPCtracker(TPCParam,currentEvent);//create the tacker for this event
237          if (!tracker) //check if it has created succeffuly
238           {//if not return with error
239             Error("Read","Can't get a tracker !\n"); 
240             return 3;
241           }
242          tracker->LoadInnerSectors();
243          tracker->LoadOuterSectors();
244    
245          for (i=0; i<NTPCtracks; i++) //loop over all tpc tracks
246           {
247             iotrack=new AliTPCtrack;   //create new tracks
248             trackbranch->SetAddress(&iotrack); //tell the branch ehere to put track data from tree(file)
249             tracktree->GetEvent(i); //stream track i to the iotrack
250             tracker->CookLabel(iotrack,0.1); //calculate (cook) the label of the tpc track
251                                              //which is the label of corresponding simulated particle 
252             tarray->AddLast(iotrack); //put the track in the array
253           }
254          
255          aTracksFile->Delete(treename);//delete tree from memmory (and leave untached on disk)- we do not need it any more
256          aTracksFile->Delete("tracks");//delete branch from memmory
257          delete tracker; //delete tracker
258          
259          tracker = 0x0;
260          trackbranch = 0x0;
261          tracktree = 0x0;
262    
263          Double_t xk;
264          Double_t par[5];
265          Float_t phi, lam, pt;//angles and transverse momentum
266          Int_t label; //label of the current track
267          
268          aGAliceFile->cd();
269          gAlice->GetEvent(currentEvent); 
270
271          gAlice->Particles();
272          
273          for (i=0; i<NTPCtracks; i++) //loop over all good tracks
274           { 
275             iotrack = (AliTPCtrack*)tarray->At(i);
276             label = iotrack->GetLabel();
277
278             if (label < 0) continue;
279             
280             TParticle *p = (TParticle*)gAlice->Particle(label);
281             
282             if(Pass(p->GetPdgCode())) continue; //check if we are intersted with particles of this type 
283                                         //if not take next partilce
284             
285             AliHBTParticle* part = new AliHBTParticle(*p);
286             if(Pass(part)) { delete part; continue;}//check if meets all criteria of any of our cuts
287                                                     //if it does not delete it and take next good track
288          
289          
290             iotrack->PropagateToVertex();
291
292             iotrack->GetExternalParameters(xk,par);     //get properties of the track
293             phi=TMath::ASin(par[2]) + iotrack->GetAlpha(); 
294             if (phi<-TMath::Pi()) phi+=2*TMath::Pi();
295             if (phi>=TMath::Pi()) phi-=2*TMath::Pi();
296             lam=par[3]; 
297             pt=1.0/TMath::Abs(par[4]);
298             
299             Double_t tpx = pt * TMath::Cos(phi); //track x coordinate of momentum
300             Double_t tpy = pt * TMath::Sin(phi); //track y coordinate of momentum
301             Double_t tpz = pt * lam; //track z coordinate of momentum
302             
303             Double_t mass = p->GetMass();
304             Double_t tEtot = TMath::Sqrt( tpx*tpx + tpy*tpy + tpz*tpz + mass*mass);//total energy of the track
305             
306             AliHBTParticle* track = new AliHBTParticle(p->GetPdgCode(), tpx, tpy , tpz, tEtot, 0., 0., 0., 0.);
307            if(Pass(track))//check if meets all criteria of any of our cuts
308                          //if it does not delete it and take next good track
309              { 
310                delete track;
311                delete part;
312                continue;
313              }
314             particles->AddParticle(totalNevents,part);//put track and particle on the run
315             tracks->AddParticle(totalNevents,track);
316
317           }
318          tarray->Clear(); //clear the array
319          
320         /**************************************/
321        /**************************************/
322       /**************************************/  
323      totalNevents++;
324     }
325   
326     //save environment (resouces) --
327     //clean your place after the work
328     CloseFiles(aTracksFile,aClustersFile,aGAliceFile); 
329     currentdir++;
330    }while(currentdir < Ndirs);
331
332
333   delete tarray;
334   fIsRead = kTRUE;
335   return 0;
336  }
337
338 /********************************************************************/
339 Int_t AliHBTReaderTPC::OpenFiles
340 (TFile*& aTracksFile, TFile*& aClustersFile, TFile*& agAliceFile,Int_t event)
341 {
342  //opens all the files
343    
344    
345    const TString& dirname = GetDirName(event); 
346    if (dirname == "")
347     {
348       Error("OpenFiles","Can not get directory name");
349       return 4;
350     }
351    
352    TString filename = dirname +"/"+ fTrackFileName;
353    aTracksFile = TFile::Open(filename.Data());
354    if ( aTracksFile  == 0x0 ) 
355      {
356        Error("OpenFiles","Can't open file with tacks named %s",filename.Data());
357        return 1;
358      }
359    if (!aTracksFile->IsOpen())
360      {
361        Error("OpenFiles","Can't open file with tacks named %s",filename.Data());
362        return 1;
363      }
364   
365    filename = dirname +"/"+ fClusterFileName;
366    aClustersFile = TFile::Open(filename.Data());
367    if ( aClustersFile == 0x0 )
368     {
369       Error("OpenFiles","Can't open file with TPC clusters named %s",filename.Data());
370       return 2;
371     }
372    if (!aClustersFile->IsOpen())
373     {
374       Error("OpenFiles","Can't open file with TPC clusters named %s",filename.Data());
375       return 2;
376     }
377
378    filename = dirname +"/"+ fGAliceFileName;
379    agAliceFile = TFile::Open(filename.Data());
380    if ( agAliceFile== 0x0)
381     {
382       Error("OpenFiles","Can't open file with TPC clusters named %s",filename.Data());
383       return 3;
384     }
385    if (!agAliceFile->IsOpen())
386     {
387       Error("OpenFiles","Can't open file with TPC clusters named %s",filename.Data());
388       return 3;
389     } 
390    
391    if (!(gAlice=(AliRun*)agAliceFile->Get("gAlice"))) 
392     {
393       Error("OpenFiles","gAlice have not been found on %s !\n",filename.Data());
394       return 5;
395     }
396
397    return 0; 
398 }
399 /********************************************************************/
400
401 /********************************************************************/
402   
403 void AliHBTReaderTPC::CloseFiles(TFile*& tracksFile, TFile*& clustersFile, TFile*& gAliceFile)
404 {
405   //closes the files
406   tracksFile->Close();
407   delete tracksFile;
408   tracksFile = 0x0;
409   clustersFile->Close();
410   delete clustersFile;
411   clustersFile = 0x0;
412
413   delete gAlice;
414   gAlice = 0;
415
416   if (gAliceFile) 
417    {
418      gAliceFile->Close();
419      delete gAliceFile;
420      gAliceFile = 0x0;
421    }
422 }
423
424 /********************************************************************/
425
426 /********************************************************************/
427 /********************************************************************/
428 /********************************************************************/
429