1 #include "AliHBTReaderPPprod.h"
9 #include <AliTPCtrack.h>
10 #include <AliTPCParam.h>
11 #include <AliTPCtracker.h>
14 #include "AliHBTRun.h"
15 #include "AliHBTEvent.h"
16 #include "AliHBTParticle.h"
17 #include "AliHBTParticleCut.h"
20 ClassImp(AliHBTReaderPPprod)
23 AliHBTReaderPPprod(const Char_t* trackfilename,const Char_t* clusterfilename,
24 const Char_t* goodtracksfilename,const Char_t* galicefilename):
25 fTrackFileName(trackfilename),fClusterFileName(clusterfilename),
26 fGAliceFileName(galicefilename),
27 fGoodTPCTracksFileName(goodtracksfilename)
29 //constructor, only file names are set
31 // trackfilename = "AliTPCtracks.root"
32 // clusterfilename = "AliTPCclusters.root"
33 // goodtracksfilename = "good_tracks_tpc"
34 // galicefilename = "" - this means: Do not open gAlice file -
35 // just leave the global pointer untached
37 fParticles = new AliHBTRun();
38 fTracks = new AliHBTRun();
40 fTracksFile = 0x0; //files are opened during reading only
45 /********************************************************************/
47 AliHBTReaderPPprod::~AliHBTReaderPPprod()
52 /********************************************************************/
54 AliHBTEvent* AliHBTReaderPPprod::GetParticleEvent(Int_t n)
56 //returns Nth event with simulated particles
57 if (!fIsRead) Read(fParticles,fTracks);
58 return fParticles->GetEvent(n);
60 /********************************************************************/
61 AliHBTEvent* AliHBTReaderPPprod::GetTrackEvent(Int_t n)
63 //returns Nth event with reconstructed tracks
64 if (!fIsRead) Read(fParticles,fTracks);
65 return fTracks->GetEvent(n);
67 /********************************************************************/
69 Int_t AliHBTReaderPPprod::GetNumberOfPartEvents()
71 //returns number of events of particles
72 if (!fIsRead) Read(fParticles,fTracks);
73 return fParticles->GetNumberOfEvents();
76 /********************************************************************/
77 Int_t AliHBTReaderPPprod::GetNumberOfTrackEvents()
79 //returns number of events of tracks
80 if (!fIsRead) Read(fParticles,fTracks);
81 return fTracks->GetNumberOfEvents();
83 /********************************************************************/
86 Int_t AliHBTReaderPPprod::Read(AliHBTRun* particles, AliHBTRun *tracks)
88 //reads data and puts put to the particles and tracks objects
89 //reurns 0 if everything is OK
91 Int_t i; //iterator and some temprary values
93 if (!particles) //check if an object is instatiated
95 Error("AliHBTReaderPPprod::Read"," particles object must instatiated before passing it to the reader");
97 if (!tracks) //check if an object is instatiated
99 Error("AliHBTReaderPPprod::Read"," tracks object must instatiated before passing it to the reader");
101 particles->Reset();//clear runs == delete all old events
104 if( (i=OpenFiles()) )
106 Error("AliHBTReaderPPprod::Read","Exiting due to problems with opening files. Errorcode %d",i);
110 AliGoodTracksPP *goodTPCTracks = new AliGoodTracksPP(fGoodTPCTracksFileName);
113 Error("AliHBTReaderPPprod::Read","Exiting due to problems with opening files. Errorcode %d",i);
120 AliTPCParam *TPCParam= (AliTPCParam*)fClustersFile->Get("75x40_100x60");
123 Error("AliHBTReaderPPprod::Read","TPC parameters have not been found !\n");
127 TObjArray *tarray = new TObjArray(5000);
128 tarray->SetOwner();// this causes memory leak, but in some cases deleting is infinite loop
131 for(Int_t currentEvent =0; currentEvent<Nevents;currentEvent++)
133 cout<<"Reading Event "<<currentEvent<<endl;
134 /**************************************/
135 /**************************************/
136 /**************************************/
138 Char_t treename[100];
139 sprintf(treename,"TreeT_TPC_%d",currentEvent);
143 tracktree=(TTree*)fTracksFile->Get(treename);
146 Error("AliHBTReaderPPprod::Read","Can't get a tree with TPC tracks !\n");
150 TBranch *trackbranch=tracktree->GetBranch("tracks");
153 Error("AliHBTReaderPPprod::Read","Can't get a branch with TPC tracks !\n");
156 Int_t NTPCtracks=(Int_t)tracktree->GetEntries();
157 cout<<"Found "<<NTPCtracks<<" TPC tracks.\n";
158 //Copy tracks to array
160 AliTPCtrack *iotrack=0;
163 //AliTPCtracker *tracker=new AliTPCtracker(TPCParam,currentEvent);//I.B.
164 AliTPCtracker *tracker=new AliTPCtracker(TPCParam); //I.B.
165 tracker->SetEventNumber(currentEvent); //I.B.
168 Error("AliHBTReaderPPprod::Read","Can't get a tracker !\n");
171 //tracker->LoadInnerSectors(); //I.B.
172 //tracker->LoadOuterSectors(); //I.B.
173 tracker->LoadClusters(); //I.B.
175 for (i=0; i<NTPCtracks; i++)
177 iotrack=new AliTPCtrack;
178 trackbranch->SetAddress(&iotrack);
179 tracktree->GetEvent(i);
180 tracker->CookLabel(iotrack,0.1);
181 tarray->AddLast(iotrack);
185 fTracksFile->Delete(treename);//delete tree from memmory (and leave untached on disk)- we do not need it any more
186 fTracksFile->Delete("tracks");
194 Int_t & ngood = goodTPCTracks->fGoodInEvent[currentEvent];
198 Float_t phi, lam, pt;
202 for (i=0; i<ngood; i++)
204 const struct GoodTrack & gt = goodTPCTracks->GetTrack(currentEvent,i);
206 if(Pass(gt.code)) continue;
209 found = kFALSE; //guard in case we don't find track with such a label
210 for (Int_t j=0;j<NTPCtracks;j++)
212 iotrack = (AliTPCtrack*)tarray->At(j);
213 if (iotrack->GetLabel() == label)
222 "Sth is going wrong with tracks - there is no TPC track corresponding to goodtrack.\nGood tack label %d",label);
223 continue; //put comunicate on the screen and continue loop
226 Double_t mass = TDatabasePDG::Instance()->GetParticle(gt.code)->Mass();
227 Double_t pEtot = TMath::Sqrt(gt.px*gt.px + gt.py*gt.py + gt.pz*gt.pz + mass*mass);
229 AliHBTParticle* part = new AliHBTParticle(gt.code, gt.px, gt.py, gt.pz, pEtot, gt.x, gt.y, gt.z, 0.0);
230 if(Pass(part)) continue;
233 iotrack->PropagateTo(gt.x);
234 iotrack->GetExternalParameters(xk,par);
235 phi=TMath::ASin(par[2]) + iotrack->GetAlpha();
236 if (phi<-TMath::Pi()) phi+=2*TMath::Pi();
237 if (phi>=TMath::Pi()) phi-=2*TMath::Pi();
239 pt=1.0/TMath::Abs(par[4]);
241 Double_t tpx = pt * TMath::Cos(phi); //track x coordinate of momentum
242 Double_t tpy = pt * TMath::Sin(phi); //track x coordinate of momentum
243 Double_t tpz = pt * lam;
245 Double_t tEtot = TMath::Sqrt( tpx*tpx + tpy*tpy + tpz*tpz + mass*mass);
247 AliHBTParticle* track = new AliHBTParticle(gt.code, tpx, tpy , tpz, tEtot, 0., 0., 0., 0.);
248 if(Pass(track)) continue;
250 particles->AddParticle(currentEvent,part);
251 tracks->AddParticle(currentEvent,track);
255 /**************************************/
256 /**************************************/
257 /**************************************/
262 delete goodTPCTracks;
267 /********************************************************************/
268 Int_t AliHBTReaderPPprod::OpenFiles()
272 fTracksFile=TFile::Open("AliTPCtracks.root");
273 if (!fTracksFile->IsOpen())
275 Error("AliHBTReaderPPprod::OpenFiles","Can't open AliTPCtracks.root");
281 fClustersFile=TFile::Open("AliTPCclusters.root");
282 if (!fClustersFile->IsOpen())
284 Error("AliHBTReaderPPprod::OpenFiles","Can't open AliTPCclusters.root");
295 /********************************************************************/
297 void AliHBTReaderPPprod::CloseFiles()
299 fTracksFile->Close();
300 fClustersFile->Close();
303 /********************************************************************/
305 /********************************************************************/
306 /********************************************************************/
307 /********************************************************************/
310 AliGoodTracksPP::~AliGoodTracksPP()
312 delete [] fGoodInEvent;
313 for (Int_t i = 0;i<fNevents;i++)
317 /********************************************************************/
318 AliGoodTracksPP::AliGoodTracksPP(const TString& infilename)
322 ifstream in(infilename.Data());
326 cerr<<"Can not open file with Good TPC Tracks named:"<<infilename.Data()<<endl;
332 fGoodInEvent = new Int_t[fNevents];
333 fData = new struct GoodTrack* [fNevents];
336 for( i = 0;i<fNevents;i++)
339 fData[i] = new struct GoodTrack[500];
345 if(fGoodInEvent[evno]>=500)
347 cerr<<"AliGoodTracksPP::AliGoodTracksPP() : Not enough place in the array\n";
350 in>>fData[evno][fGoodInEvent[evno]].lab;
351 in>>fData[evno][fGoodInEvent[evno]].code;
352 in>>fData[evno][fGoodInEvent[evno]].px;
353 in>>fData[evno][fGoodInEvent[evno]].py;
354 in>>fData[evno][fGoodInEvent[evno]].pz;
355 in>>fData[evno][fGoodInEvent[evno]].x;
356 in>>fData[evno][fGoodInEvent[evno]].y;
357 in>>fData[evno][fGoodInEvent[evno]].z;
364 cout<<fData[evno][fGoodInEvent[evno]].lab;
365 cout<<" ";cout<<fData[evno][fGoodInEvent[evno]].code;
366 cout<<" ";cout<<fData[evno][fGoodInEvent[evno]].px;
367 cout<<" ";cout<<fData[evno][fGoodInEvent[evno]].py;
368 cout<<" ";cout<<fData[evno][fGoodInEvent[evno]].pz;
369 cout<<" ";cout<<fData[evno][fGoodInEvent[evno]].x;
370 cout<<" ";cout<<fData[evno][fGoodInEvent[evno]].y;
371 cout<<" ";cout<<fData[evno][fGoodInEvent[evno]].z;
374 fGoodInEvent[evno]++;
377 cout<<"AliGoodTracksPP::AliGoodTracksPP() .... Done\n";
382 const GoodTrack& AliGoodTracksPP::GetTrack(Int_t event, Int_t n) const
385 if( (event>fNevents) || (event<0))
387 gAlice->Fatal("AliGoodTracksPP::GetTrack","No such Event %d",event);
389 if( (n>fGoodInEvent[event]) || (n<0))
391 gAlice->Fatal("AliGoodTracksPP::GetTrack","No such Good TPC Track %d",n);
393 return fData[event][n];