]>
Commit | Line | Data |
---|---|---|
1 | #ifndef __CINT__ | |
2 | #include <Riostream.h> | |
3 | #include <fstream.h> | |
4 | ||
5 | #include "AliRun.h" | |
6 | #include "AliITS.h" | |
7 | #include "AliITSgeom.h" | |
8 | #include "AliITStrackerV2.h" | |
9 | #include "AliITStrackV2.h" | |
10 | #include "AliITSclusterV2.h" | |
11 | ||
12 | #include "TFile.h" | |
13 | #include "TTree.h" | |
14 | #include "TH1.h" | |
15 | #include "TH2.h" | |
16 | #include "TObjArray.h" | |
17 | #include "TStyle.h" | |
18 | #include "TCanvas.h" | |
19 | #include "TLine.h" | |
20 | #include "TText.h" | |
21 | #include "TParticle.h" | |
22 | #endif | |
23 | ||
24 | struct GoodTrackITS { | |
25 | Int_t lab; | |
26 | Int_t code; | |
27 | Float_t px,py,pz; | |
28 | Float_t x,y,z; | |
29 | }; | |
30 | ||
31 | Int_t AliITSComparisonV2(Int_t event=0) { | |
32 | cerr<<"Doing comparison...\n"; | |
33 | ||
34 | const Int_t MAX=15000; | |
35 | Int_t good_tracks_its(GoodTrackITS *gt, const Int_t max, const Int_t event); | |
36 | ||
37 | Int_t nentr=0; TObjArray tarray(2000); | |
38 | {/* Load tracks */ | |
39 | TFile *tf=TFile::Open("AliITStracksV2.root"); | |
40 | if (!tf->IsOpen()) {cerr<<"Can't open AliITStracksV2.root !\n"; return 3;} | |
41 | char tname[100]; sprintf(tname,"TreeT_ITS_%d",event); | |
42 | TTree *tracktree=(TTree*)tf->Get(tname); | |
43 | if (!tracktree) {cerr<<"Can't get a tree with ITS tracks !\n"; return 4;} | |
44 | TBranch *tbranch=tracktree->GetBranch("tracks"); | |
45 | nentr=(Int_t)tracktree->GetEntries(); | |
46 | ||
47 | for (Int_t i=0; i<nentr; i++) { | |
48 | AliITStrackV2 *iotrack=new AliITStrackV2; | |
49 | tbranch->SetAddress(&iotrack); | |
50 | tracktree->GetEvent(i); | |
51 | tarray.AddLast(iotrack); | |
52 | /*if (itsLabel != 1234) continue; | |
53 | Int_t nc=iotrack->GetNumberOfClusters(); | |
54 | for (Int_t k=0; k<nc; k++) { | |
55 | Int_t index=iotrack->GetClusterIndex(k); | |
56 | AliITSclusterV2 *c=tracker.GetCluster(index); | |
57 | cout<<c->GetLabel(0)<<' '<<c->GetLabel(1)<<' '<<c->GetLabel(2)<<endl; | |
58 | }*/ | |
59 | } | |
60 | delete tracktree; //Thanks to Mariana Bondila | |
61 | tf->Close(); | |
62 | } | |
63 | ||
64 | /* Generate a list of "good" tracks */ | |
65 | GoodTrackITS gt[MAX]; | |
66 | Int_t ngood=0; | |
67 | ifstream in("good_tracks_its"); | |
68 | if (in) { | |
69 | cerr<<"Reading good tracks...\n"; | |
70 | while (in>>gt[ngood].lab>>gt[ngood].code>> | |
71 | gt[ngood].px>>gt[ngood].py>>gt[ngood].pz>> | |
72 | gt[ngood].x >>gt[ngood].y >>gt[ngood].z) { | |
73 | ngood++; | |
74 | cerr<<ngood<<'\r'; | |
75 | if (ngood==MAX) { | |
76 | cerr<<"Too many good tracks !\n"; | |
77 | break; | |
78 | } | |
79 | } | |
80 | if (!in.eof()) cerr<<"Read error (good_tracks_its) !\n"; | |
81 | } else { | |
82 | cerr<<"Marking good tracks (this will take a while)...\n"; | |
83 | ngood=good_tracks_its(gt,MAX,event); | |
84 | ofstream out("good_tracks_its"); | |
85 | if (out) { | |
86 | for (Int_t ngd=0; ngd<ngood; ngd++) | |
87 | out<<gt[ngd].lab<<' '<<gt[ngd].code<<' ' | |
88 | <<gt[ngd].px<<' '<<gt[ngd].py<<' '<<gt[ngd].pz<<' ' | |
89 | <<gt[ngd].x <<' '<<gt[ngd].y <<' '<<gt[ngd].z <<endl; | |
90 | } else cerr<<"Can not open file (good_tracks_its) !\n"; | |
91 | out.close(); | |
92 | } | |
93 | ||
94 | TH1F *hp=new TH1F("hp","PHI resolution",50,-20.,20.); hp->SetFillColor(4); | |
95 | TH1F *hl=new TH1F("hl","LAMBDA resolution",50,-20,20);hl->SetFillColor(4); | |
96 | TH1F *hpt=new TH1F("hpt","Relative Pt resolution",30,-10.,10.); | |
97 | hpt->SetFillColor(2); | |
98 | TH1F *hmpt=new TH1F("hmpt","Transverse impact parameter",30,-300,300); | |
99 | hmpt->SetFillColor(6); | |
100 | TH1F *hz=new TH1F("hz","Longitudinal impact parameter",30,-300,300); | |
101 | ||
102 | AliITStrackV2 *trk=(AliITStrackV2*)tarray.UncheckedAt(0); | |
103 | Double_t pmin=0.1*(100/0.299792458/0.2/trk->GetConvConst()); | |
104 | Double_t pmax=6.0+pmin; | |
105 | ||
106 | TH1F *hgood=new TH1F("hgood","Good tracks",30,pmin,pmax); | |
107 | TH1F *hfound=new TH1F("hfound","Found tracks",30,pmin,pmax); | |
108 | TH1F *hfake=new TH1F("hfake","Fake tracks",30,pmin,pmax); | |
109 | TH1F *hg=new TH1F("hg","",30,pmin,pmax); //efficiency for good tracks | |
110 | hg->SetLineColor(4); hg->SetLineWidth(2); | |
111 | TH1F *hf=new TH1F("hf","Efficiency for fake tracks",30,pmin,pmax); | |
112 | hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2); | |
113 | ||
114 | //TH1F *hptw=new TH1F("hptw","Weghted pt",30,pmax,pmin); | |
115 | ||
116 | TH1F *he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,100.); | |
117 | TH2F *hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,200.); | |
118 | hep->SetMarkerStyle(8); | |
119 | hep->SetMarkerSize(0.4); | |
120 | ||
121 | ||
122 | Int_t notf[MAX], nnotf=0; | |
123 | Int_t fake[MAX], nfake=0; | |
124 | Int_t mult[MAX], numb[MAX], nmult=0; | |
125 | ||
126 | Int_t ng; | |
127 | for (ng=0; ng<ngood; ng++) { | |
128 | Int_t lab=gt[ng].lab, tlab=-1; | |
129 | Double_t pxg=gt[ng].px, pyg=gt[ng].py, pzg=gt[ng].pz; | |
130 | Double_t ptg=TMath::Sqrt(pxg*pxg+pyg*pyg); | |
131 | ||
132 | ||
133 | if (ptg>pmin) hgood->Fill(ptg); | |
134 | ||
135 | AliITStrackV2 *track=0; | |
136 | Int_t cnt=0; | |
137 | Int_t j; | |
138 | for (j=0; j<nentr; j++) { | |
139 | AliITStrackV2 *trk=(AliITStrackV2*)tarray.UncheckedAt(j); | |
140 | Int_t lbl=trk->GetLabel(); | |
141 | if (lab==TMath::Abs(lbl)) { | |
142 | if (cnt==0) {track=trk; tlab=lbl;} | |
143 | cnt++; | |
144 | } | |
145 | } | |
146 | if (cnt==0) { | |
147 | notf[nnotf++]=lab; | |
148 | continue; | |
149 | } else if (cnt>1){ | |
150 | mult[nmult]=lab; | |
151 | numb[nmult]=cnt; nmult++; | |
152 | } | |
153 | ||
154 | if (ptg>pmin) { | |
155 | if (lab==tlab) hfound->Fill(ptg); | |
156 | else { | |
157 | fake[nfake++]=lab; | |
158 | hfake->Fill(ptg); | |
159 | } | |
160 | } | |
161 | ||
162 | track->PropagateTo(3.,0.0028,65.19); | |
163 | track->PropagateToVertex(); | |
164 | ||
165 | Double_t xv,par[5]; track->GetExternalParameters(xv,par); | |
166 | Float_t phi=TMath::ASin(par[2]) + track->GetAlpha(); | |
167 | if (phi<-TMath::Pi()) phi+=2*TMath::Pi(); | |
168 | if (phi>=TMath::Pi()) phi-=2*TMath::Pi(); | |
169 | Float_t lam=TMath::ATan(par[3]); | |
170 | Float_t pt_1=TMath::Abs(par[4]); | |
171 | ||
172 | Double_t phig=TMath::ATan2(pyg,pxg); | |
173 | hp->Fill((phi - phig)*1000.); | |
174 | ||
175 | Double_t lamg=TMath::ATan2(pzg,ptg); | |
176 | hl->Fill((lam - lamg)*1000.); | |
177 | ||
178 | Double_t d=10000*track->GetD(); | |
179 | hmpt->Fill(d); | |
180 | ||
181 | //hptw->Fill(ptg,TMath::Abs(d)); | |
182 | ||
183 | Double_t z=10000*track->GetZ(); | |
184 | hz->Fill(z); | |
185 | ||
186 | hpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.); | |
187 | ||
188 | Float_t mom=1./(pt_1*TMath::Cos(lam)); | |
189 | Float_t dedx=track->GetdEdx(); | |
190 | hep->Fill(mom,dedx,1.); | |
191 | if (TMath::Abs(gt[ng].code)==211) | |
192 | if (mom>0.4 && mom<0.5) he->Fill(dedx,1.); | |
193 | } | |
194 | ||
195 | ||
196 | cout<<"\nList of Not found tracks :\n"; | |
197 | for (ng=0; ng<nnotf; ng++){ | |
198 | cout<<notf[ng]<<"\t"; | |
199 | if ((ng%9)==8) cout<<"\n"; | |
200 | } | |
201 | cout<<"\n\nList of fake tracks :\n"; | |
202 | for (ng=0; ng<nfake; ng++){ | |
203 | cout<<fake[ng]<<"\t"; | |
204 | if ((ng%9)==8) cout<<"\n"; | |
205 | } | |
206 | cout<<"\n\nList of multiple found tracks :\n"; | |
207 | for (ng=0; ng<nmult; ng++) { | |
208 | cout<<"id. "<<mult[ng] | |
209 | <<" found - "<<numb[ng]<<"times\n"; | |
210 | } | |
211 | cout<<endl; | |
212 | ||
213 | cerr<<"Number of good tracks : "<<ngood<<endl; | |
214 | cerr<<"Number of found tracks : "<<nentr<<endl; | |
215 | ||
216 | ng=(Int_t)hgood->GetEntries(); //cerr<<"Good tracks "<<ng<<endl; | |
217 | if (ng!=0) | |
218 | cerr<<"Integral efficiency is about "<<hfound->GetEntries()/ng*100.<<" %\n"; | |
219 | cerr<<endl; | |
220 | ||
221 | gStyle->SetOptStat(111110); | |
222 | gStyle->SetOptFit(1); | |
223 | ||
224 | TCanvas *c1=new TCanvas("c1","",0,0,700,850); | |
225 | ||
226 | TPad *p1=new TPad("p1","",0,0.3,.5,.6); p1->Draw(); | |
227 | p1->cd(); p1->SetFillColor(42); p1->SetFrameFillColor(10); | |
228 | hp->SetFillColor(4); hp->SetXTitle("(mrad)"); hp->Fit("gaus"); c1->cd(); | |
229 | ||
230 | TPad *p2=new TPad("p2","",0.5,.3,1,.6); p2->Draw(); | |
231 | p2->cd(); p2->SetFillColor(42); p2->SetFrameFillColor(10); | |
232 | hl->SetXTitle("(mrad)"); hl->Fit("gaus"); c1->cd(); | |
233 | ||
234 | TPad *p3=new TPad("p3","",0,0,0.5,0.3); p3->Draw(); | |
235 | p3->cd(); p3->SetFillColor(42); p3->SetFrameFillColor(10); | |
236 | hpt->SetXTitle("(%)"); hpt->Fit("gaus"); c1->cd(); | |
237 | ||
238 | TPad *p4=new TPad("p4","",0.5,0,1,0.3); p4->Draw(); | |
239 | p4->cd(); p4->SetFillColor(42); p4->SetFrameFillColor(10); | |
240 | hmpt->SetXTitle("(micron)"); hmpt->Fit("gaus"); hz->Draw("same"); c1->cd(); | |
241 | //hfound->Sumw2(); | |
242 | //hptw->Sumw2(); | |
243 | //hg->SetMaximum(333); | |
244 | //hg->SetYTitle("Impact Parameter Resolution (micron)"); | |
245 | //hg->SetXTitle("Pt (GeV/c)"); | |
246 | //hg->GetXaxis()->SetRange(0,10); | |
247 | //hg->Divide(hptw,hfound,1,1.); | |
248 | //hg->DrawCopy(); c1->cd(); | |
249 | ||
250 | ||
251 | TPad *p5=new TPad("p5","",0,0.6,1,1); p5->Draw(); p5->cd(); | |
252 | p5->SetFillColor(41); p5->SetFrameFillColor(10); | |
253 | hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2(); | |
254 | hg->Divide(hfound,hgood,1,1.,"b"); | |
255 | hf->Divide(hfake,hgood,1,1.,"b"); | |
256 | hg->SetMaximum(1.4); | |
257 | hg->SetYTitle("Tracking efficiency"); | |
258 | hg->SetXTitle("Pt (GeV/c)"); | |
259 | hg->Draw(); | |
260 | ||
261 | TLine *line1 = new TLine(pmin,1.0,pmax,1.0); line1->SetLineStyle(4); | |
262 | line1->Draw("same"); | |
263 | TLine *line2 = new TLine(pmin,0.9,pmax,0.9); line2->SetLineStyle(4); | |
264 | line2->Draw("same"); | |
265 | ||
266 | hf->SetFillColor(1); | |
267 | hf->SetFillStyle(3013); | |
268 | hf->SetLineColor(2); | |
269 | hf->SetLineWidth(2); | |
270 | hf->Draw("histsame"); | |
271 | TText *text = new TText(0.461176,0.248448,"Fake tracks"); | |
272 | text->SetTextSize(0.05); | |
273 | text->Draw(); | |
274 | text = new TText(0.453919,1.11408,"Good tracks"); | |
275 | text->SetTextSize(0.05); | |
276 | text->Draw(); | |
277 | ||
278 | ||
279 | ||
280 | TCanvas *c2=new TCanvas("c2","",320,32,530,590); | |
281 | ||
282 | TPad *p6=new TPad("p6","",0.,0.,1.,.5); p6->Draw(); | |
283 | p6->cd(); p6->SetFillColor(42); p6->SetFrameFillColor(10); | |
284 | he->SetFillColor(2); he->SetFillStyle(3005); | |
285 | he->SetXTitle("Arbitrary Units"); | |
286 | he->Fit("gaus"); c2->cd(); | |
287 | ||
288 | TPad *p7=new TPad("p7","",0.,0.5,1.,1.); p7->Draw(); | |
289 | p7->cd(); p7->SetFillColor(42); p7->SetFrameFillColor(10); | |
290 | hep->SetXTitle("p (Gev/c)"); hep->SetYTitle("dE/dX (Arb. Units)"); | |
291 | hep->Draw(); c1->cd(); | |
292 | ||
293 | ||
294 | return 0; | |
295 | } | |
296 | ||
297 | Int_t good_tracks_its(GoodTrackITS *gt, const Int_t max, const Int_t event) { | |
298 | if (gAlice) {delete gAlice; gAlice=0;} | |
299 | ||
300 | TFile *file=TFile::Open("galice.root"); | |
301 | if (!file->IsOpen()) {cerr<<"Can't open galice.root !\n"; exit(4);} | |
302 | if (!(gAlice=(AliRun*)file->Get("gAlice"))) { | |
303 | cerr<<"gAlice have not been found on galice.root !\n"; | |
304 | exit(5); | |
305 | } | |
306 | ||
307 | Int_t np=gAlice->GetEvent(event); | |
308 | ||
309 | Int_t *good=new Int_t[np]; | |
310 | Int_t k; | |
311 | for (k=0; k<np; k++) good[k]=0; | |
312 | ||
313 | AliITS *ITS=(AliITS*)gAlice->GetDetector("ITS"); | |
314 | if (!ITS) { | |
315 | cerr<<"can't get ITS !\n"; exit(8); | |
316 | } | |
317 | AliITSgeom *geom=ITS->GetITSgeom(); | |
318 | if (!geom) { | |
319 | cerr<<"can't get ITS geometry !\n"; exit(9); | |
320 | } | |
321 | ||
322 | TFile *cf=TFile::Open("AliITSclustersV2.root"); | |
323 | if (!cf->IsOpen()){ | |
324 | cerr<<"Can't open AliITSclustersV2.root !\n"; exit(6); | |
325 | } | |
326 | char cname[100]; sprintf(cname,"TreeC_ITS_%d",event); | |
327 | TTree *cTree=(TTree*)cf->Get(cname); | |
328 | if (!cTree) { | |
329 | cerr<<"Can't get cTree !\n"; exit(7); | |
330 | } | |
331 | TBranch *branch=cTree->GetBranch("Clusters"); | |
332 | if (!branch) { | |
333 | cerr<<"Can't get clusters branch !\n"; exit(8); | |
334 | } | |
335 | TClonesArray *clusters=new TClonesArray("AliITSclusterV2",10000); | |
336 | branch->SetAddress(&clusters); | |
337 | ||
338 | Int_t entr=(Int_t)cTree->GetEntries(); | |
339 | for (k=0; k<entr; k++) { | |
340 | cTree->GetEvent(k); | |
341 | Int_t ncl=clusters->GetEntriesFast(); if (ncl==0) continue; | |
342 | Int_t lay,lad,det; geom->GetModuleId(k,lay,lad,det); | |
343 | if (lay<1 || lay>6) { | |
344 | cerr<<"wrong layer !\n"; exit(10); | |
345 | } | |
346 | while (ncl--) { | |
347 | AliITSclusterV2 *pnt=(AliITSclusterV2*)clusters->UncheckedAt(ncl); | |
348 | Int_t l0=pnt->GetLabel(0); | |
349 | Int_t l1=pnt->GetLabel(1); | |
350 | Int_t l2=pnt->GetLabel(2); | |
351 | Int_t mask=1<<(lay-1); | |
352 | if (l0>=0) good[l0]|=mask; | |
353 | if (l1>=0) good[l1]|=mask; | |
354 | if (l2>=0) good[l2]|=mask; | |
355 | } | |
356 | } | |
357 | clusters->Delete(); delete clusters; | |
358 | delete cTree; //Thanks to Mariana Bondila | |
359 | cf->Close(); | |
360 | ||
361 | ifstream in("good_tracks_tpc"); | |
362 | if (!in) { | |
363 | cerr<<"can't get good_tracks_tpc !\n"; exit(11); | |
364 | } | |
365 | Int_t nt=0; | |
366 | Double_t px,py,pz,x,y,z; | |
367 | Int_t code,lab; | |
368 | while (in>>lab>>code>>px>>py>>pz>>x>>y>>z) { | |
369 | if (good[lab] != 0x3F) continue; | |
370 | TParticle *p = (TParticle*)gAlice->Particle(lab); | |
371 | gt[nt].lab=lab; | |
372 | gt[nt].code=p->GetPdgCode(); | |
373 | //**** px py pz - in global coordinate system | |
374 | gt[nt].px=p->Px(); gt[nt].py=p->Py(); gt[nt].pz=p->Pz(); | |
375 | gt[nt].x=gt[nt].y=gt[nt].z=0.; | |
376 | nt++; | |
377 | if (nt==max) {cerr<<"Too many good tracks !\n"; break;} | |
378 | } | |
379 | ||
380 | delete[] good; | |
381 | ||
382 | delete gAlice; gAlice=0; | |
383 | file->Close(); | |
384 | ||
385 | return nt; | |
386 | } | |
387 | ||
388 |