]> git.uio.no Git - u/mrichter/AliRoot.git/blame - ITS/AliTPCComparison4ITS.C
The access to several data members was changed from public to protected. The digitisa...
[u/mrichter/AliRoot.git] / ITS / AliTPCComparison4ITS.C
CommitLineData
c6f60ee1 1#ifndef __CINT__
2 #include "alles.h"
3 #include "AliTPCtracker.h"
4#endif
5static const kNGoodTracks = 120000;
6
7struct GoodTrack {
8 Int_t fEventN; //event number
9 Int_t lab;
10 Int_t code;
11 Float_t px,py,pz;
12 Float_t x,y,z;
13 Float_t pxg,pyg,pzg,ptg;
14 Bool_t flag;
15};
16Int_t good_tracks(GoodTrack *gt, Int_t max, Int_t eventn=1);
17
18Int_t AliTPCComparison4ITS(Int_t eventn=1) {
19 Int_t i;
20
21 gBenchmark->Start("AliTPCComparison");
22
23 TFile *cf=TFile::Open("AliTPCclusters.root");
24 if (!cf->IsOpen()) {cerr<<"Can't open AliTPCclusters.root !\n"; return 1;}
25 AliTPCParam *digp= (AliTPCParam*)cf->Get("75x40_100x60");
26 if (!digp) { cerr<<"TPC parameters have not been found !\n"; return 2; }
27
28 ///////////
29 AliTPCtracker *tracker =0;
30 TObjArray tarray(5000);
31 AliTPCtrack *iotrack=0;
32 Int_t nentr= 0;
33 Int_t eventptr[1000];
34 TFile *tf=TFile::Open("AliTPCtracksSorted.root");
35 TTree *tracktree=0;
36 // Load clusters
37 eventptr[0]=0;
38 eventptr[1]=0;
39 for (Int_t event=0;event<eventn; event++){
40 cf->cd();
41 tracker = new AliTPCtracker(digp,event);
42 tracker->LoadInnerSectors();
43 tracker->LoadOuterSectors();
44
45
46 char tname[100];
47 cout<<"eventn = "<<eventn<<endl;
48 if (eventn==-1) {
49 sprintf(tname,"TreeT_TPC");
50 }
51 else {
52 sprintf(tname,"TreeT_TPC_%d",event);
53 }
54
55 // Load tracks
56 if (!tf->IsOpen()) {cerr<<"Can't open AliTPCtracks.root !\n"; return 3;}
57 tracktree=(TTree*)tf->Get(tname);
58 if (!tracktree) {cerr<<"Can't get a tree with TPC tracks !\n"; return 4;}
59
60 TBranch *tbranch=tracktree->GetBranch("tracks");
61 Int_t nentr0=(Int_t)tracktree->GetEntries();
62 nentr+=nentr0;
63 cout<<"Got "<<nentr0<<" entries\n";
64 for (i=0; i<nentr0; i++) {
65 iotrack=new AliTPCtrack;
66 tbranch->SetAddress(&iotrack);
67 tracktree->GetEvent(i);
68 //iotrack->Dump();
69 tracker->CookLabel(iotrack,0.1);
70 tarray.AddLast(iotrack);
71 }
72 eventptr[event+1] = nentr; //store end of the event
73 delete tracker;
74 }
75 tf->Close();
76 cf->Close();
77
78 //MI change for a moment
79
80
81/////////////////////////////////////////////////////////////////////////
82 // GoodTrack gt[15000];
83 GoodTrack gt[kNGoodTracks];
84
85 Int_t ngood=0;
86 ifstream in("good_tracks_tpc");
87 /*********************************************************/
88 if (in) //Case good tracks already exists
89 {
90 cerr<<"Reading good tracks...\n";
91 while (in>>gt[ngood].fEventN>>gt[ngood].lab>>gt[ngood].code>>
92 gt[ngood].px>>gt[ngood].py>>gt[ngood].pz>>
93 gt[ngood].x >>gt[ngood].y >>gt[ngood].z>>
94 gt[ngood].pxg>>gt[ngood].pyg>>gt[ngood].pzg>>
95 gt[ngood].ptg>>gt[ngood].flag) {
96 ngood++;
97 cerr<<ngood<<'\r';
98 if (ngood==kNGoodTracks) {
99 cerr<<"Too many good tracks !\n";
100 break;
101 }
102 }
103 if (!in.eof()) cerr<<"Read error (good_tracks_tpc) !\n";
104 }
105 /*********************************************************/
106 else //Case good tracks do not yet exist
107 {
108 cerr<<"Marking good tracks (this will take a while)...\n";
109 ngood=good_tracks(gt,kNGoodTracks,eventn); //mi change
110 ofstream out("good_tracks_tpc");
111 if (out) {
112 for (Int_t ngd=0; ngd<ngood; ngd++)
113 out<<gt[ngd].fEventN<<' '<<gt[ngd].lab<<' '<<gt[ngd].code<<' '<<
114 gt[ngd].px<<' '<<gt[ngd].py<<' '<<gt[ngd].pz<<' '<<
115 gt[ngd].x <<' '<<gt[ngd].y <<' '<<gt[ngd].z<<' '<<
116 gt[ngd].pxg<<' '<<gt[ngd].pyg<<' '<<gt[ngd].pzg<<' '<<' '
117 <<gt[ngd].ptg<<gt[ngd].flag<<endl;
118 } else cerr<<"Can not open file (good_tracks_tpc) !\n";
119 out.close();
120
121 /*
122 cerr<<"Preparing tracks for matching with the ITS...\n";
123 tarray.Sort();
124 tf=TFile::Open("AliTPCtracks.root","recreate");
125 tracktree=new TTree("TPCf","Tree with TPC tracks");
126 tracktree->Branch("tracks","AliTPCtrack",&iotrack,32000,0);
127 for (i=0; i<nentr; i++) {
128 iotrack=(AliTPCtrack*)tarray.UncheckedAt(i);
129 tracktree->Fill();
130 }
131 tracktree->Write();
132 tf->Close();
133 */
134 }
135 /*********************************************************/
136 cerr<<"Number of good tracks : "<<ngood<<endl;
137
138 TH1F *hp=new TH1F("hp","PHI resolution",50,-20.,20.); hp->SetFillColor(4);
139 TH1F *hl=new TH1F("hl","LAMBDA resolution",50,-20,20);hl->SetFillColor(4);
140 TH1F *hpt=new TH1F("hpt","Relative Pt resolution",30,-10.,10.);
141 hpt->SetFillColor(2);
142 TH1F *hmpt=new TH1F("hmpt","Relative Pt resolution (pt>4GeV/c)",30,-60,60);
143 hmpt->SetFillColor(6);
144
145 AliTPCtrack *trk=(AliTPCtrack*)tarray.UncheckedAt(0);
146 Double_t pmin=0.1*(100/0.299792458/0.2/trk->GetConvConst());
147 Double_t pmax=6.0+pmin;
148
149 TH1F *hgood=new TH1F("hgood","Good tracks",30,pmin,pmax);
150 TH1F *hfound=new TH1F("hfound","Found tracks",30,pmin,pmax);
151 TH1F *hfake=new TH1F("hfake","Fake tracks",30,pmin,pmax);
152 TH1F *hg=new TH1F("hg","",30,pmin,pmax); //efficiency for good tracks
153 hg->SetLineColor(4); hg->SetLineWidth(2);
154 TH1F *hf=new TH1F("hf","Efficiency for fake tracks",30,pmin,pmax);
155 hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
156
157 TH1F *he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,100.);
158 TH2F *hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,250.);
159 hep->SetMarkerStyle(8);
160 hep->SetMarkerSize(0.4);
161
162 Int_t mingood = ngood; //MI change
163 Int_t * track_notfound = new Int_t[ngood];
164 Int_t itrack_notfound =0;
165 Int_t * track_fake = new Int_t[ngood];
166 Int_t itrack_fake = 0;
167 Int_t * track_multifound = new Int_t[ngood];
168 Int_t * track_multifound_n = new Int_t[ngood];
169 Int_t itrack_multifound =0;
170
171 while (ngood--) {
172 Int_t lab=gt[ngood].lab,tlab=-1;
173 Float_t ptg=
174 TMath::Sqrt(gt[ngood].px*gt[ngood].px + gt[ngood].py*gt[ngood].py);
175
176 if (ptg<pmin) continue;
177
178 hgood->Fill(ptg);
179 Int_t ievent = gt[ngood].fEventN;
180
181 AliTPCtrack *track=0;
182 // for (i=0; i<nentr; i++) {
183 for (i=eventptr[ievent]; i<eventptr[ievent+1]; i++) {
184
185 track=(AliTPCtrack*)tarray.UncheckedAt(i);
186 tlab=track->GetLabel();
187 if (lab==TMath::Abs(tlab)) break;
188 }
189 // if (i==nentr) {
190 if (i==eventptr[ievent+1]) {
191
192 cerr<<"Track "<<lab<<" was not found !\n";
193 track_notfound[itrack_notfound++]=lab;
194 continue;
195 }
196
197 //MI change - addition
198 Int_t micount=0;
199 Int_t mi;
200 AliTPCtrack * mitrack;
201 // for (mi=0; mi<nentr; mi++) {
202 for (mi=eventptr[ievent]; mi<eventptr[ievent+1]; mi++) {
203
204 mitrack=(AliTPCtrack*)tarray.UncheckedAt(mi);
205 if (lab==TMath::Abs(mitrack->GetLabel())) micount++;
206 }
207 if (micount>1) {
208 //cout<<"Track no. "<<lab<<" found "<<micount<<" times\n";
209 track_multifound[itrack_multifound]=lab;
210 track_multifound_n[itrack_multifound]=micount;
211 itrack_multifound++;
212 }
213
214
215 //
216 Double_t xk=gt[ngood].x;//digp->GetPadRowRadii(0,0);
217 track->PropagateTo(xk);
218
219 if (lab==tlab) hfound->Fill(ptg);
220 else {
221 track_fake[itrack_fake++]=lab;
222 hfake->Fill(ptg);
223 //cerr<<lab<<" fake\n";
224 }
225
226 Double_t par[5]; track->GetExternalParameters(xk,par);
227 Float_t phi=TMath::ASin(par[2]) + track->GetAlpha();
228 if (phi<-TMath::Pi()) phi+=2*TMath::Pi();
229 if (phi>=TMath::Pi()) phi-=2*TMath::Pi();
230 Float_t lam=TMath::ATan(par[3]);
231 Float_t pt_1=TMath::Abs(par[4]);
232
233 if (TMath::Abs(gt[ngood].code)==11 && ptg>4.) {
234 hmpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
235 } else {
236 Float_t phig=TMath::ATan2(gt[ngood].py,gt[ngood].px);
237 hp->Fill((phi - phig)*1000.);
238
239 Float_t lamg=TMath::ATan2(gt[ngood].pz,ptg);
240 hl->Fill((lam - lamg)*1000.);
241
242 hpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
243 }
244
245 Float_t mom=1./(pt_1*TMath::Cos(lam));
246 Float_t dedx=track->GetdEdx();
247 hep->Fill(mom,dedx,1.);
248 if (TMath::Abs(gt[ngood].code)==211)
249 if (mom>0.4 && mom<0.5) {
250 he->Fill(dedx,1.);
251 }
252
253 }
254
255
256 Stat_t ng=hgood->GetEntries(); cerr<<"Good tracks "<<ng<<endl;
257 Stat_t nf=hfound->GetEntries();
258 if (ng!=0)
259 cerr<<"\n\n\nIntegral efficiency is about "<<nf/ng*100.<<" %\n";
260
261 //MI change - addition
262 cout<<"Total number of found tracks ="<<nentr<<"\n";
263 cout<<"Total number of \"good\" tracks ="<<mingood<<"\n";
264
265
266 cout<<"\nList of Not found tracks :\n";
267 // Int_t i;
268 for ( i = 0; i< itrack_notfound; i++){
269 cout<<track_notfound[i]<<"\t";
270 if ((i%5)==4) cout<<"\n";
271 }
272 cout<<"\nList of fake tracks :\n";
273 for ( i = 0; i< itrack_fake; i++){
274 cout<<track_fake[i]<<"\t";
275 if ((i%5)==4) cout<<"\n";
276 }
277 cout<<"\nList of multiple found tracks :\n";
278 for ( i=0; i<itrack_multifound; i++) {
279 cout<<"id. "<<track_multifound[i]<<" found - "<<track_multifound_n[i]<<"times\n";
280 }
281 cout<<"\n\n\n";
282
283
284 gStyle->SetOptStat(111110);
285 gStyle->SetOptFit(1);
286
287 TCanvas *c1=new TCanvas("c1","",0,0,700,850);
288
289 TPad *p1=new TPad("p1","",0,0.3,.5,.6); p1->Draw();
290 p1->cd(); p1->SetFillColor(42); p1->SetFrameFillColor(10);
291 hp->SetFillColor(4); hp->SetXTitle("(mrad)"); hp->Fit("gaus"); c1->cd();
292
293 TPad *p2=new TPad("p2","",0.5,.3,1,.6); p2->Draw();
294 p2->cd(); p2->SetFillColor(42); p2->SetFrameFillColor(10);
295 hl->SetXTitle("(mrad)"); hl->Fit("gaus"); c1->cd();
296
297 TPad *p3=new TPad("p3","",0,0,0.5,0.3); p3->Draw();
298 p3->cd(); p3->SetFillColor(42); p3->SetFrameFillColor(10);
299 hpt->SetXTitle("(%)"); hpt->Fit("gaus"); c1->cd();
300
301 TPad *p4=new TPad("p4","",0.5,0,1,0.3); p4->Draw();
302 p4->cd(); p4->SetFillColor(42); p4->SetFrameFillColor(10);
303 hmpt->SetXTitle("(%)"); hmpt->Fit("gaus"); c1->cd();
304
305 TPad *p5=new TPad("p5","",0,0.6,1,1); p5->Draw(); p5->cd();
306 p5->SetFillColor(41); p5->SetFrameFillColor(10);
307 hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
308 hg->Divide(hfound,hgood,1,1.,"b");
309 hf->Divide(hfake,hgood,1,1.,"b");
310 hg->SetMaximum(1.4);
311 hg->SetYTitle("Tracking efficiency");
312 hg->SetXTitle("Pt (GeV/c)");
313 hg->Draw();
314
315 TLine *line1 = new TLine(pmin,1.0,pmax,1.0); line1->SetLineStyle(4);
316 line1->Draw("same");
317 TLine *line2 = new TLine(pmin,0.9,pmax,0.9); line2->SetLineStyle(4);
318 line2->Draw("same");
319
320 hf->SetFillColor(1);
321 hf->SetFillStyle(3013);
322 hf->SetLineColor(2);
323 hf->SetLineWidth(2);
324 hf->Draw("histsame");
325 TText *text = new TText(0.461176,0.248448,"Fake tracks");
326 text->SetTextSize(0.05);
327 text->Draw();
328 text = new TText(0.453919,1.11408,"Good tracks");
329 text->SetTextSize(0.05);
330 text->Draw();
331
332
333
334 TCanvas *c2=new TCanvas("c2","",320,32,530,590);
335
336 TPad *p6=new TPad("p6","",0.,0.,1.,.5); p6->Draw();
337 p6->cd(); p6->SetFillColor(42); p6->SetFrameFillColor(10);
338 he->SetFillColor(2); he->SetFillStyle(3005);
339 he->SetXTitle("Arbitrary Units");
340 he->Fit("gaus"); c2->cd();
341
342 TPad *p7=new TPad("p7","",0.,0.5,1.,1.); p7->Draw();
343 p7->cd(); p7->SetFillColor(42); p7->SetFrameFillColor(10);
344 hep->SetXTitle("p (Gev/c)"); hep->SetYTitle("dE/dX (Arb. Units)");
345 hep->Draw(); c1->cd();
346
347 gBenchmark->Stop("AliTPCComparison");
348 gBenchmark->Show("AliTPCComparison");
349
350
351 return 0;
352}
353
354
355/*********************************************************/
356/*********************************************************/
357/*********************************************************/
358/*********************************************************/
359
360Int_t good_tracks(GoodTrack *gt, Int_t max, Int_t eventn) {
361 //eventn - number of events in file
362
363 TFile *file=TFile::Open("rfio:galice.root");
364 if (!file->IsOpen()) {cerr<<"Can't open galice.root !\n"; exit(4);}
365
366 if (!(gAlice=(AliRun*)file->Get("gAlice"))) {
367 cerr<<"gAlice have not been found on galice.root !\n";
368 exit(5);
369 }
370
371 // Int_t np=gAlice->GetEvent(0); //MI change
372
373 AliTPC *TPC=(AliTPC*)gAlice->GetDetector("TPC");
374 Int_t ver = TPC->IsVersion();
375 cerr<<"TPC version "<<ver<<" has been found !\n";
376
377 AliTPCParam *digp=(AliTPCParam*)file->Get("75x40_100x60");
378 if (!digp) { cerr<<"TPC parameters have not been found !\n"; exit(6); }
379 TPC->SetParam(digp);
380
381 Int_t nrow_up=digp->GetNRowUp();
382 Int_t nrows=digp->GetNRowLow()+nrow_up;
383 Int_t zero=digp->GetZeroSup();
384 Int_t gap=Int_t(0.125*nrows);
385 Int_t good_number=Int_t(0.4*nrows);
386
387
388 //MI change
389 Int_t nt=0; //reset counter
390 char treeName[100]; //declare event identifier
391
392 for (Int_t event=0;event<eventn;event++){
393
394 Int_t np=gAlice->GetEvent(event);
395 Int_t *good=new Int_t[np];
396 for (Int_t ii=0; ii<np; ii++) good[ii]=0;
397
398
399 //MI change to be possible compile macro
400 //definition out of the swith statemnet
401 Int_t sectors_by_rows=0;
402 TTree *TD=0;
403 AliSimDigits da, *digits=&da;
404 Int_t *count=0;
405 switch (ver) {
406 case 1:
407 {
408 TFile *cf=TFile::Open("AliTPCclusters.root");
409 if (!cf->IsOpen()){cerr<<"Can't open AliTPCclusters.root !\n";exit(5);}
410 AliTPCClustersArray *ca=new AliTPCClustersArray;
411 ca->Setup(digp);
412 ca->SetClusterType("AliTPCcluster");
413 ca->ConnectTree("Segment Tree");
414 Int_t nrows=Int_t(ca->GetTree()->GetEntries());
415 for (Int_t n=0; n<nrows; n++) {
416 AliSegmentID *s=ca->LoadEntry(n);
417 Int_t sec,row;
418 digp->AdjustSectorRow(s->GetID(),sec,row);
419 AliTPCClustersRow &clrow = *ca->GetRow(sec,row);
420 Int_t ncl=clrow.GetArray()->GetEntriesFast();
421 while (ncl--) {
422 AliTPCcluster *c=(AliTPCcluster*)clrow[ncl];
423 Int_t lab=c->GetLabel(0);
424 if (lab<0) continue; //noise cluster
425 lab=TMath::Abs(lab);
426 if (sec>=digp->GetNInnerSector())
427 if (row==nrow_up-1 ) good[lab]|=0x1000;
428 if (sec>=digp->GetNInnerSector())
429 if (row==nrow_up-1-gap) good[lab]|=0x800;
430 good[lab]++;
431 }
432 ca->ClearRow(sec,row);
433 }
434 cf->Close();
435 }
436 break;
437 case 2:
438
439 sprintf(treeName,"TreeD_75x40_100x60_%d",event);
440 TD=(TTree*)gDirectory->Get(treeName);
441 TD->GetBranch("Segment")->SetAddress(&digits);
442 count = new Int_t[np];
443 Int_t i;
444 for (i=0; i<np; i++) count[i]=0;
445 sectors_by_rows=(Int_t)TD->GetEntries();
446 for (i=0; i<sectors_by_rows; i++) {
447 if (!TD->GetEvent(i)) continue;
448 Int_t sec,row;
449 digp->AdjustSectorRow(digits->GetID(),sec,row);
450 cerr<<sec<<' '<<row<<" \r";
451 digits->First();
452 do { //Many thanks to J.Chudoba who noticed this
453 Int_t it=digits->CurrentRow(), ip=digits->CurrentColumn();
454 Short_t dig = digits->GetDigit(it,ip);
455 Int_t idx0=digits->GetTrackID(it,ip,0);
456 Int_t idx1=digits->GetTrackID(it,ip,1);
457 Int_t idx2=digits->GetTrackID(it,ip,2);
458 if (idx0>=0 && dig>=zero) count[idx0]+=1;
459 if (idx1>=0 && dig>=zero) count[idx1]+=1;
460 if (idx2>=0 && dig>=zero) count[idx2]+=1;
461 } while (digits->Next());
462 for (Int_t j=0; j<np; j++) {
463 if (count[j]>1) {
464 if (sec>=digp->GetNInnerSector())
465 if (row==nrow_up-1 ) good[j]|=0x1000;
466 if (sec>=digp->GetNInnerSector())
467 if (row==nrow_up-1-gap) good[j]|=0x800;
468 good[j]++;
469 }
470 count[j]=0;
471 }
472 }
473 delete[] count;
474 break;
475 default:
476 cerr<<"Invalid TPC version !\n";
477 file->Close();
478 exit(7);
479 }
480
481 TTree *TH=gAlice->TreeH();
482 Int_t npart=(Int_t)TH->GetEntries();
483
484 while (npart--) {
485 AliTPChit *hit0=0;
486
487 TPC->ResetHits();
488 TH->GetEvent(npart);
489 AliTPChit * hit = (AliTPChit*) TPC->FirstHit(-1);
490 while (hit){
491 if (hit->fQ==0.) break;
492 hit = (AliTPChit*) TPC->NextHit();
493 }
494 if (hit) {
495 hit0 = new AliTPChit(*hit); //Make copy of hit
496 hit = hit0;
497 }
498 else continue;
499 AliTPChit *hit1=(AliTPChit*)TPC->NextHit();
500 if (hit1==0) continue;
501 if (hit1->fQ != 0.) continue;
502 Int_t i=hit->Track();
503 TParticle *p = (TParticle*)gAlice->Particle(i);
504
505 if (p->GetFirstMother()>=0) continue; //secondary particle
506 if (good[i] < 0x1000+0x800+2+good_number) continue;
507 if (p->Pt()<0.100) continue;
508 if (TMath::Abs(p->Pz()/p->Pt())>0.999) continue;
509
510 gt[nt].lab=i;
511 gt[nt].code=p->GetPdgCode();
512 //**** px py pz - in global coordinate system, x y z - in local !
513 gt[nt].px=hit->X(); gt[nt].py=hit->Y(); gt[nt].pz=hit->Z();
514 Float_t cs,sn; digp->AdjustCosSin(hit1->fSector,cs,sn);
515 gt[nt].fEventN=event; //MI change
516 gt[nt].x = hit1->X()*cs + hit1->Y()*sn;
517 gt[nt].y =-hit1->X()*sn + hit1->Y()*cs;
518 gt[nt].z = hit1->Z();
519 gt[nt].pxg = p->Px();
520 gt[nt].pyg = p->Py();
521 gt[nt].pzg = p->Pz();
522 gt[nt].ptg = p->Pt();
523 gt[nt].flag = 0;
524 nt++;
525 if (hit0) delete hit0;
526 cerr<<i<<" \r";
527 if (nt==max) {cerr<<"Too many good tracks !\n"; break;}
528 }
529 delete[] good;
530 }
531 //delete gAlice; gAlice=0;
532
533 file->Close();
534 gBenchmark->Stop("AliTPCComparison");
535 gBenchmark->Show("AliTPCComparison");
536 return nt;
537}
538
539
540