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