]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCComparison.C
Added documentation (Marian)
[u/mrichter/AliRoot.git] / TPC / AliTPCComparison.C
CommitLineData
7f6ddf58 1/****************************************************************************
2 * Very important, delicate and rather obscure macro. *
3 * *
4 * Creates list of "trackable" tracks, *
7f6ddf58 5 * calculates efficiency, resolutions etc. *
bf018fe4 6 * There is a possibility to run this macro over several events. *
7f6ddf58 7 * *
024a7fe9 8 * Origin: I.Belikov, CERN, Jouri.Belikov@cern.ch *
7f6ddf58 9 * with several nice improvements by: M.Ivanov, GSI, m.ivanov@gsi.de *
10 ****************************************************************************/
11
c630aafd 12#if !defined(__CINT__) || defined(__MAKECINT__)
5102bab6 13 #include <TMath.h>
14 #include <TError.h>
15 #include <Riostream.h>
16 #include <TH1F.h>
17 #include <TH2F.h>
18 #include <TTree.h>
19 #include <TParticle.h>
5102bab6 20 #include <TCanvas.h>
21 #include <TLine.h>
22 #include <TText.h>
23 #include <TBenchmark.h>
24 #include <TStyle.h>
8b462fd8 25 #include <TFile.h>
bf018fe4 26 #include <TROOT.h>
5102bab6 27
28 #include "AliStack.h"
29 #include "AliHeader.h"
30 #include "AliTrackReference.h"
88cb7938 31 #include "AliRunLoader.h"
c630aafd 32 #include "AliRun.h"
4514157e 33 #include "AliESDEvent.h"
fc023202 34 #include "AliESDtrack.h"
5102bab6 35
36 #include "AliSimDigits.h"
37 #include "AliTPC.h"
38 #include "AliTPCParamSR.h"
39 #include "AliTPCClustersArray.h"
40 #include "AliTPCClustersRow.h"
bf018fe4 41 #include "AliTPCcluster.h"
5102bab6 42 #include "AliTPCLoader.h"
be9c5115 43#endif
44
4fa1e629 45Int_t GoodTracksTPC(const Char_t *dir=".");
88cb7938 46
c630aafd 47extern AliRun *gAlice;
5102bab6 48extern TBenchmark *gBenchmark;
bf018fe4 49extern TROOT *gROOT;
fbb58996 50
4fa1e629 51static Int_t allgood=0;
52static Int_t allselected=0;
53static Int_t allfound=0;
54
55Int_t AliTPCComparison
56(Float_t ptcutl=0.2, Float_t ptcuth=10., const Char_t *dir=".") {
5102bab6 57 gBenchmark->Start("AliTPCComparison");
58
59 ::Info("AliTPCComparison.C","Doing comparison...");
c630aafd 60
88cb7938 61
bf018fe4 62
63 TH1F *hp=(TH1F*)gROOT->FindObject("hp");
64 if (!hp) hp=new TH1F("hp","PHI resolution",50,-20.,20.);
65 hp->SetFillColor(4);
66
67 TH1F *hl=(TH1F*)gROOT->FindObject("hl");
68 if (!hl) hl=new TH1F("hl","LAMBDA resolution",50,-20,20);
69 hl->SetFillColor(4);
70
71 TH1F *hpt=(TH1F*)gROOT->FindObject("hpt");
72 if (!hpt) hpt=new TH1F("hpt","Relative Pt resolution",30,-10.,10.);
73 hpt->SetFillColor(2);
74
75 TH1F *hmpt=(TH1F*)gROOT->FindObject("hmpt");
76 if (!hmpt)
77 hmpt=new TH1F("hmpt","Relative Pt resolution (pt>4GeV/c)",30,-60,60);
cc80f89e 78 hmpt->SetFillColor(6);
8c555625 79
bf018fe4 80
81
82 TH1F *hgood=(TH1F*)gROOT->FindObject("hgood");
83 if (!hgood) hgood=new TH1F("hgood","Good tracks",30,0.2,6.1);
84
85 TH1F *hfound=(TH1F*)gROOT->FindObject("hfound");
86 if (!hfound) hfound=new TH1F("hfound","Found tracks",30,0.2,6.1);
87
88 TH1F *hfake=(TH1F*)gROOT->FindObject("hfake");
89 if (!hfake) hfake=new TH1F("hfake","Fake tracks",30,0.2,6.1);
90
91 TH1F *hg=(TH1F*)gROOT->FindObject("hg");
92 if (!hg) hg=new TH1F("hg","Efficiency for good tracks",30,0.2,6.1);
8c555625 93 hg->SetLineColor(4); hg->SetLineWidth(2);
bf018fe4 94
95 TH1F *hf=(TH1F*)gROOT->FindObject("hf");
96 if (!hf) hf=new TH1F("hf","Efficiency for fake tracks",30,0.2,6.1);
8c555625 97 hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
98
bf018fe4 99 TH1F *he=(TH1F*)gROOT->FindObject("he");
100 if (!he)
101 he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,100.);
102
103 TH2F *hep=(TH2F*)gROOT->FindObject("hep");
104 if (!hep) hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,400.);
be9c5115 105 hep->SetMarkerStyle(8);
106 hep->SetMarkerSize(0.4);
8c555625 107
4fa1e629 108
109 Char_t fname[100];
110 sprintf(fname,"%s/GoodTracksTPC.root",dir);
111
112 TFile *refFile=TFile::Open(fname,"old");
113 if (!refFile || !refFile->IsOpen()) {
114 ::Info("AliTPCComparison.C","Marking good tracks (will take a while)...");
115 if (GoodTracksTPC(dir)) {
116 ::Error("AliTPCComparison.C","Can't generate the reference file !");
117 return 1;
118 }
119 }
120 refFile=TFile::Open(fname,"old");
121 if (!refFile || !refFile->IsOpen()) {
122 ::Error("AliTPCComparison.C","Can't open the reference file !");
123 return 2;
124 }
125
126 TTree *tpcTree=(TTree*)refFile->Get("tpcTree");
127 if (!tpcTree) {
128 ::Error("AliTPCComparison.C","Can't get the reference tree !");
129 return 3;
130 }
131 TBranch *branch=tpcTree->GetBranch("TPC");
132 if (!branch) {
133 ::Error("AliTPCComparison.C","Can't get the TPC branch !");
134 return 4;
135 }
136 TClonesArray dummy("AliTrackReference",1000), *refs=&dummy;
137 branch->SetAddress(&refs);
138
139
140 sprintf(fname,"%s/AliESDs.root",dir);
141 TFile *ef=TFile::Open(fname);
142 if ((!ef)||(!ef->IsOpen())) {
143 sprintf(fname,"%s/AliESDtpc.root",dir);
144 ef=TFile::Open(fname);
145 if ((!ef)||(!ef->IsOpen())) {
146 ::Error("AliTPCComparison.C","Can't open AliESDtpc.root !");
147 return 5;
148 }
149 }
4514157e 150 AliESDEvent* event = new AliESDEvent();
8b462fd8 151 TTree* esdTree = (TTree*) ef->Get("esdTree");
152 if (!esdTree) {
153 ::Error("AliTPCComparison.C", "no ESD tree found");
154 return 6;
155 }
4514157e 156 event->ReadFromTree(esdTree);
4fa1e629 157
158
159 //******* Loop over events *********
160 Int_t e=0;
8b462fd8 161 while (esdTree->GetEvent(e)) {
4fa1e629 162 cout<<endl<<endl<<"********* Processing event number: "<<e<<"*******\n";
163
4fa1e629 164 Int_t nentr=event->GetNumberOfTracks();
165 allfound+=nentr;
166
167 if (tpcTree->GetEvent(e++)==0) {
168 cerr<<"No reconstructable tracks !\n";
169 continue;
170 }
171
172 Int_t ngood=refs->GetEntriesFast();
173 allgood+=ngood;
174
175 const Int_t MAX=15000;
176 //MI change
177 Int_t track_notfound[MAX], itrack_notfound=0;
178 Int_t track_fake[MAX], itrack_fake=0;
179 Int_t track_multifound[MAX],track_multifound_n[MAX],itrack_multifound=0;
180
181 Int_t i;
182 for (Int_t k=0; k<ngood; k++) {
183 AliTrackReference *ref=(AliTrackReference*)refs->UncheckedAt(k);
184 Int_t lab=ref->Label(), tlab=-1;
185 Float_t ptg=TMath::Sqrt(ref->Px()*ref->Px() + ref->Py()*ref->Py());
186
187 if (ptg<1e-33) continue; // for those not crossing 0 pad row
188
189 if (ptg<ptcutl) continue;
190 if (ptg>ptcuth) continue;
191
192 allselected++;
193
194 hgood->Fill(ptg);
195
196 AliESDtrack *track=0;
197 for (i=0; i<nentr; i++) {
bf018fe4 198 track=event->GetTrack(i);
162085b7 199 tlab=track->GetTPCLabel();
cc80f89e 200 if (lab==TMath::Abs(tlab)) break;
4fa1e629 201 }
202 if (i==nentr) {
203 track_notfound[itrack_notfound++]=lab;
204 continue;
205 }
206
207 //MI change - addition
208 Int_t micount=0;
209 Int_t mi;
210 AliESDtrack * mitrack;
211 for (mi=0; mi<nentr; mi++) {
212 mitrack=event->GetTrack(mi);
213 if (lab==TMath::Abs(mitrack->GetTPCLabel())) micount++;
214 }
215 if (micount>1) {
216 track_multifound[itrack_multifound]=lab;
217 track_multifound_n[itrack_multifound]=micount;
218 itrack_multifound++;
219 }
220 if ((track->GetStatus()&AliESDtrack::kTPCrefit)==0) continue;
221 if (lab==tlab) hfound->Fill(ptg);
222 else {
223 track_fake[itrack_fake++]=lab;
224 hfake->Fill(ptg);
225 }
cf98c13f 226
4fa1e629 227 Double_t pxpypz[3]; track->GetInnerPxPyPz(pxpypz);
228 Float_t phi=TMath::ATan2(pxpypz[1],pxpypz[0]);
229 if (phi<-TMath::Pi()) phi+=2*TMath::Pi();
230 if (phi>=TMath::Pi()) phi-=2*TMath::Pi();
231 Double_t pt=TMath::Sqrt(pxpypz[0]*pxpypz[0]+pxpypz[1]*pxpypz[1]);
232 Float_t lam=TMath::ATan2(pxpypz[2],pt);
233 Float_t pt_1=1/pt;
234
235 Int_t pdg=(Int_t)ref->GetLength(); //this is particle's PDG !
236
237 if (TMath::Abs(pdg)==11 && ptg>4.) {//high momentum electrons
238 hmpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
239 } else {
240 Float_t phig=TMath::ATan2(ref->Py(),ref->Px());
241 hp->Fill((phi - phig)*1000.);
242
243 Float_t lamg=TMath::ATan2(ref->Pz(),ptg);
244 hl->Fill((lam - lamg)*1000.);
245
246 hpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
247 }
248
249 Float_t mom=pt/TMath::Cos(lam);
250 Float_t dedx=track->GetTPCsignal();
251 hep->Fill(mom,dedx,1.);
252 if (TMath::Abs(pdg)==211) //pions
253 if (mom>0.4 && mom<0.5) {
254 he->Fill(dedx,1.);
255 }
cf98c13f 256 }
4fa1e629 257
258 cout<<"\nList of Not found tracks :\n";
259 for ( i = 0; i< itrack_notfound; i++){
260 cout<<track_notfound[i]<<"\t";
261 if ((i%5)==4) cout<<"\n";
cf98c13f 262 }
4fa1e629 263 cout<<"\nList of fake tracks :\n";
264 for ( i = 0; i< itrack_fake; i++){
265 cout<<track_fake[i]<<"\t";
266 if ((i%5)==4) cout<<"\n";
cf98c13f 267 }
4fa1e629 268 cout<<"\nList of multiple found tracks :\n";
269 for ( i=0; i<itrack_multifound; i++) {
270 cout<<"id. "<<track_multifound[i]
271 <<" found - "<<track_multifound_n[i]<<"times\n";
cc80f89e 272 }
73042f01 273
4fa1e629 274 cout<<"Number of found tracks ="<<nentr<<endl;
275 cout<<"Number of \"good\" tracks ="<<ngood<<endl;
73042f01 276
4fa1e629 277 refs->Clear();
4fa1e629 278 }// ***** End of the loop over events
cc80f89e 279
8b462fd8 280 delete event;
4514157e 281 delete esdTree;
4fa1e629 282 ef->Close();
8c555625 283
4fa1e629 284 delete tpcTree;
285 refFile->Close();
7f6ddf58 286
4fa1e629 287 Stat_t ng=hgood->GetEntries(), nf=hfound->GetEntries();
288 if (ng!=0) cout<<"\n\nIntegral efficiency is about "<<nf/ng*100.<<" %\n";
289 cout<<"Total number selected of \"good\" tracks ="<<allselected<<endl<<endl;
290 cout<<"Total number of found tracks ="<<allfound<<endl;
291 cout<<"Total number of \"good\" tracks ="<<allgood<<endl;
292 cout<<endl;
7f6ddf58 293
8c555625 294 gStyle->SetOptStat(111110);
295 gStyle->SetOptFit(1);
296
297 TCanvas *c1=new TCanvas("c1","",0,0,700,850);
298
2d343a94 299 Int_t minc=33;
300
8c555625 301 TPad *p1=new TPad("p1","",0,0.3,.5,.6); p1->Draw();
302 p1->cd(); p1->SetFillColor(42); p1->SetFrameFillColor(10);
2d343a94 303 hp->SetFillColor(4); hp->SetXTitle("(mrad)");
304 if (hp->GetEntries()<minc) hp->Draw(); else hp->Fit("gaus"); c1->cd();
8c555625 305
306 TPad *p2=new TPad("p2","",0.5,.3,1,.6); p2->Draw();
307 p2->cd(); p2->SetFillColor(42); p2->SetFrameFillColor(10);
2d343a94 308 hl->SetXTitle("(mrad)");
309 if (hl->GetEntries()<minc) hl->Draw(); else hl->Fit("gaus"); c1->cd();
8c555625 310
311 TPad *p3=new TPad("p3","",0,0,0.5,0.3); p3->Draw();
312 p3->cd(); p3->SetFillColor(42); p3->SetFrameFillColor(10);
2d343a94 313 hpt->SetXTitle("(%)");
314 if (hpt->GetEntries()<minc) hpt->Draw(); else hpt->Fit("gaus"); c1->cd();
8c555625 315
316 TPad *p4=new TPad("p4","",0.5,0,1,0.3); p4->Draw();
317 p4->cd(); p4->SetFillColor(42); p4->SetFrameFillColor(10);
2d343a94 318 hmpt->SetXTitle("(%)");
319 if (hmpt->GetEntries()<minc) hmpt->Draw(); else hmpt->Fit("gaus"); c1->cd();
8c555625 320
321 TPad *p5=new TPad("p5","",0,0.6,1,1); p5->Draw(); p5->cd();
3c0f9266 322 p5->SetFillColor(41); p5->SetFrameFillColor(10);
8c555625 323 hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
324 hg->Divide(hfound,hgood,1,1.,"b");
325 hf->Divide(hfake,hgood,1,1.,"b");
326 hg->SetMaximum(1.4);
327 hg->SetYTitle("Tracking efficiency");
328 hg->SetXTitle("Pt (GeV/c)");
329 hg->Draw();
330
7f6ddf58 331 TLine *line1 = new TLine(0.1,1.0,6.1,1.0); line1->SetLineStyle(4);
8c555625 332 line1->Draw("same");
7f6ddf58 333 TLine *line2 = new TLine(0.1,0.9,6.1,0.9); line2->SetLineStyle(4);
8c555625 334 line2->Draw("same");
335
336 hf->SetFillColor(1);
337 hf->SetFillStyle(3013);
338 hf->SetLineColor(2);
339 hf->SetLineWidth(2);
340 hf->Draw("histsame");
341 TText *text = new TText(0.461176,0.248448,"Fake tracks");
342 text->SetTextSize(0.05);
343 text->Draw();
344 text = new TText(0.453919,1.11408,"Good tracks");
345 text->SetTextSize(0.05);
346 text->Draw();
3c0f9266 347
348
349
350 TCanvas *c2=new TCanvas("c2","",320,32,530,590);
351
352 TPad *p6=new TPad("p6","",0.,0.,1.,.5); p6->Draw();
353 p6->cd(); p6->SetFillColor(42); p6->SetFrameFillColor(10);
354 he->SetFillColor(2); he->SetFillStyle(3005);
355 he->SetXTitle("Arbitrary Units");
2d343a94 356 if (he->GetEntries()<minc) he->Draw(); else he->Fit("gaus"); c2->cd();
3c0f9266 357
358 TPad *p7=new TPad("p7","",0.,0.5,1.,1.); p7->Draw();
359 p7->cd(); p7->SetFillColor(42); p7->SetFrameFillColor(10);
360 hep->SetXTitle("p (Gev/c)"); hep->SetYTitle("dE/dX (Arb. Units)");
361 hep->Draw(); c1->cd();
362
2d343a94 363 TFile fc("AliTPCComparison.root","RECREATE");
364 c1->Write();
365 c2->Write();
366 fc.Close();
367
4fa1e629 368 gBenchmark->Stop("AliTPCComparison");
369 gBenchmark->Show("AliTPCComparison");
370
73042f01 371 return 0;
cc80f89e 372}
3c0f9266 373
cc80f89e 374
4fa1e629 375Int_t GoodTracksTPC(const Char_t *dir) {
376 if (gAlice) {
377 delete gAlice->GetRunLoader();
378 delete gAlice;//if everything was OK here it is already NULL
379 gAlice = 0x0;
380 }
381
382 Char_t fname[100];
383 sprintf(fname,"%s/galice.root",dir);
cc80f89e 384
4fa1e629 385 AliRunLoader *rl = AliRunLoader::Open(fname,"COMPARISON");
386 if (!rl) {
387 ::Error("GoodTracksTPC","Can't start session !");
388 return 1;
c630aafd 389 }
4fa1e629 390
391 rl->LoadgAlice();
392 rl->LoadHeader();
393 rl->LoadKinematics();
394 rl->LoadTrackRefs();
395
c630aafd 396 AliTPCLoader *tpcl = (AliTPCLoader *)rl->GetLoader("TPCLoader");
397 if (tpcl == 0x0) {
4fa1e629 398 ::Error("GoodTracksTPC","Can not find TPCLoader !");
399 delete rl;
400 return 2;
c630aafd 401 }
4fa1e629 402
88cb7938 403 AliTPC *TPC=(AliTPC*)rl->GetAliRun()->GetDetector("TPC");
73042f01 404 Int_t ver = TPC->IsVersion();
4fa1e629 405 cout<<"TPC version "<<ver<<" has been found !\n";
406 if (ver==1) tpcl->LoadRecPoints();
407 else if (ver==2) tpcl->LoadDigits();
408 else {
409 ::Error("GoodTracksTPC","Wrong TPC version !");
410 delete rl;
411 return 3;
412 }
73042f01 413
88cb7938 414 rl->CdGAFile();
415 AliTPCParamSR *digp=(AliTPCParamSR*)gDirectory->Get("75x40_100x60_150x60");
5102bab6 416 if (!digp) {
4fa1e629 417 ::Error("AliTPCHits2Digits.C","TPC parameters have not been found !");
418 delete rl;
419 return 4;
5102bab6 420 }
73042f01 421 TPC->SetParam(digp);
422
423 Int_t nrow_up=digp->GetNRowUp();
424 Int_t nrows=digp->GetNRowLow()+nrow_up;
425 Int_t zero=digp->GetZeroSup();
7f6ddf58 426 Int_t gap=Int_t(0.125*nrows), shift=Int_t(0.5*gap);
cc80f89e 427 Int_t good_number=Int_t(0.4*nrows);
428
4fa1e629 429 Int_t nev=rl->GetNumberOfEvents();
430 ::Info("GoodTracksTPC","Number of events : %d\n",nev);
431
432
433 sprintf(fname,"%s/GoodTracksTPC.root",dir);
434 TFile *file=TFile::Open(fname,"recreate");
435
436 TClonesArray dummy("AliTrackReference",1000), *refs=&dummy;
437 TTree tpcTree("tpcTree","Tree with info about the reconstructable TPC tracks");
438 tpcTree.Branch("TPC",&refs);
439
440 //******** Loop over generated events
441 for (Int_t e=0; e<nev; e++) {
442 Int_t nt=0;
443 refs->Clear();
444
445 Int_t i;
446
447 rl->GetEvent(e); file->cd();
448
449 Int_t np = rl->GetHeader()->GetNtrack();
450 cout<<"Event "<<e<<" Number of particles: "<<np<<endl;
451
452 //******** Fill the "good" masks
453 Int_t *good=new Int_t[np]; for (i=0; i<np; i++) good[i]=0;
454 switch (ver) {
455 case 1:
456 {
457 AliTPCClustersArray *pca=new AliTPCClustersArray, &ca=*pca;
458 ca.Setup(digp);
459 ca.SetClusterType("AliTPCcluster");
460 ca.ConnectTree(tpcl->TreeR());
461 Int_t nrows=Int_t(ca.GetTree()->GetEntries());
462 for (Int_t n=0; n<nrows; n++) {
7f6ddf58 463 AliSegmentID *s=ca.LoadEntry(n);
464 Int_t sec,row;
465 digp->AdjustSectorRow(s->GetID(),sec,row);
466 AliTPCClustersRow &clrow = *ca.GetRow(sec,row);
467 Int_t ncl=clrow.GetArray()->GetEntriesFast();
468 while (ncl--) {
469 AliTPCcluster *c=(AliTPCcluster*)clrow[ncl];
470 Int_t lab=c->GetLabel(0);
471 if (lab<0) continue; //noise cluster
472 lab=TMath::Abs(lab);
473
474 if (sec>=digp->GetNInnerSector())
475 if (row==nrow_up-1) good[lab]|=0x4000;
476 if (sec>=digp->GetNInnerSector())
477 if (row==nrow_up-1-gap) good[lab]|=0x1000;
478
479 if (sec>=digp->GetNInnerSector())
480 if (row==nrow_up-1-shift) good[lab]|=0x2000;
481 if (sec>=digp->GetNInnerSector())
482 if (row==nrow_up-1-gap-shift) good[lab]|=0x800;
483
484 good[lab]++;
485 }
486 ca.ClearRow(sec,row);
4fa1e629 487 }
488 delete pca;
489 }
490 break;
491 case 2:
492 {
493 TTree *TD=tpcl->TreeD();
88cb7938 494
4fa1e629 495 AliSimDigits da, *digits=&da;
496 TD->GetBranch("Segment")->SetAddress(&digits);
7f6ddf58 497
4fa1e629 498 Int_t *count = new Int_t[np];
499 Int_t i;
500 for (i=0; i<np; i++) count[i]=0;
7f6ddf58 501
4fa1e629 502 Int_t sectors_by_rows=(Int_t)TD->GetEntries();
503 for (i=0; i<sectors_by_rows; i++) {
7f6ddf58 504 if (!TD->GetEvent(i)) continue;
505 Int_t sec,row;
506 digp->AdjustSectorRow(digits->GetID(),sec,row);
4fa1e629 507 //cerr<<sec<<' '<<row<<'\r';
7f6ddf58 508 digits->First();
509 do { //Many thanks to J.Chudoba who noticed this
510 Int_t it=digits->CurrentRow(), ip=digits->CurrentColumn();
511 Short_t dig = digits->GetDigit(it,ip);
512 Int_t idx0=digits->GetTrackID(it,ip,0);
513 Int_t idx1=digits->GetTrackID(it,ip,1);
514 Int_t idx2=digits->GetTrackID(it,ip,2);
8ded1b7a 515 if (idx0>=0 && dig>=zero && idx0<np) count[idx0]+=1;
516 if (idx1>=0 && dig>=zero && idx1<np) count[idx1]+=1;
517 if (idx2>=0 && dig>=zero && idx2<np) count[idx2]+=1;
afc42102 518 } while (digits->Next());
7f6ddf58 519 for (Int_t j=0; j<np; j++) {
520 if (count[j]>1) {
521 if (sec>=digp->GetNInnerSector())
522 if (row==nrow_up-1 ) good[j]|=0x4000;
523 if (sec>=digp->GetNInnerSector())
524 if (row==nrow_up-1-gap) good[j]|=0x1000;
525
526 if (sec>=digp->GetNInnerSector())
527 if (row==nrow_up-1-shift) good[j]|=0x2000;
528 if (sec>=digp->GetNInnerSector())
529 if (row==nrow_up-1-gap-shift) good[j]|=0x800;
530 good[j]++;
531 }
532 count[j]=0;
533 }
4fa1e629 534 }
535 delete[] count;
536 }
537 break;
afc42102 538 }
7f6ddf58 539
4fa1e629 540
541 //****** select tracks which are "good" enough
542 AliStack* stack = rl->Stack();
543 for (i=0; i<np; i++) {
544 if ((good[i]&0x5000) != 0x5000)
545 if ((good[i]&0x2800) != 0x2800) continue;
546 if ((good[i]&0x7FF ) < good_number) continue;
547
548 TParticle *p = (TParticle*)stack->Particle(i);
549 if (p == 0x0) {
88cb7938 550 cerr<<"Can not get particle "<<i<<endl;
551 continue;
024a7fe9 552 }
4fa1e629 553 if (p->Pt()<0.100) continue;
554 if (TMath::Abs(p->Pz()/p->Pt())>0.999) continue;
bf018fe4 555
4fa1e629 556 Double_t vx=p->Vx(),vy=p->Vy();
557 if (TMath::Sqrt(vx*vx+vy*vy)>3.5) continue;
7f6ddf58 558
4fa1e629 559 AliTrackReference *ref=new((*refs)[nt]) AliTrackReference();
560
561 ref->SetLabel(i);
562 Int_t pdg=p->GetPdgCode();
563 ref->SetLength(pdg); //This will the particle's PDG !
564 ref->SetMomentum(0.,0.,0.);
565 ref->SetPosition(0.,0.,0.);
566
567 nt++;
568 }
569
570 //**** check if there is also information at the entrance of the TPC
571 TTree *TR=rl->TreeTR();
fc023202 572 TBranch *branch=TR->GetBranch("TrackReferences");
4fa1e629 573 if (branch==0) {
574 ::Error("GoodTracksTPC","No track references !");
575 delete rl;
576 return 5;
577 }
578 TClonesArray tpcdummy("AliTrackReference",1000), *tpcRefs=&tpcdummy;
579 branch->SetAddress(&tpcRefs);
580 Int_t nr=(Int_t)TR->GetEntries();
581 for (Int_t r=0; r<nr; r++) {
582 //cerr<<r<<' '<<nr<<'\r';
583 TR->GetEvent(r);
fc023202 584
585 Int_t nref = tpcRefs->GetEntriesFast();
586 if (!nref) continue;
587 AliTrackReference *tpcRef= 0x0;
588 for (Int_t iref=0; iref<nref; ++iref) {
589 tpcRef = (AliTrackReference*)tpcRefs->UncheckedAt(iref);
590 if (tpcRef->DetectorId() == AliTrackReference::kTPC) break;
591 tpcRef = 0x0;
592 }
593
594 if (!tpcRef) continue;
595
4fa1e629 596 Int_t j;
597 AliTrackReference *ref=0;
598 for (j=0; j<nt; j++) {
599 ref=(AliTrackReference *)refs->UncheckedAt(j);
600 if (ref->Label()==tpcRef->Label()) break;
601 ref=0;
602 }
603 if (ref==0) continue;
604
605 ref->SetMomentum(tpcRef->Px(),tpcRef->Py(),tpcRef->Pz());
606 ref->SetPosition(tpcRef->LocalX(),tpcRef->LocalY(),tpcRef->Z());
607
608 tpcRefs->Clear();
609 }
610
611 tpcTree.Fill();
612
613 delete[] good;
614 } //****** end of the loop over generated events
7f6ddf58 615
88cb7938 616
4fa1e629 617 tpcTree.Write();
618 file->Close();
7f6ddf58 619
4fa1e629 620 delete rl;
621 return 0;
8c555625 622}