]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliITSComparisonV2.C
First commit of ITS tracking version 2 (Yu.Belikov)
[u/mrichter/AliRoot.git] / ITS / AliITSComparisonV2.C
1 #ifndef __CINT__
2   #include <iostream.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 "TObjArray.h"
16   #include "TStyle.h"
17   #include "TCanvas.h"
18   #include "TLine.h"
19   #include "TText.h"
20   #include "TParticle.h"
21 #endif
22
23 struct GoodTrack {
24   Int_t lab;
25   Int_t code;
26   Float_t px,py,pz;
27   Float_t x,y,z;
28 };
29 Int_t good_tracks(GoodTrack *gt, Int_t max);
30
31 Int_t AliITSComparisonV2() {
32    cerr<<"Doing comparison...\n";
33
34    TFile *cf=TFile::Open("AliITSclustersV2.root");
35    if (!cf->IsOpen()) {cerr<<"Can't open AliITSclustersV2.root !\n"; return 1;}
36    AliITSgeom *geom=(AliITSgeom*)cf->Get("AliITSgeom");
37    if (!geom) { cerr<<"Can't get the ITS geometry !\n"; return 2; }
38    AliITStrackerV2 tracker(geom);   
39
40 // Load tracks
41    TFile *tf=TFile::Open("AliITStracksV2.root");
42    if (!tf->IsOpen()) {cerr<<"Can't open AliITStracksV2.root !\n"; return 3;}
43    TObjArray tarray(2000);
44    TTree *tracktree=(TTree*)tf->Get("TreeT");
45    TBranch *tbranch=tracktree->GetBranch("tracks");
46    Int_t nentr=(Int_t)tracktree->GetEntries(),i;
47    for (i=0; i<nentr; i++) {
48        AliITStrackV2 *iotrack=new AliITStrackV2;
49        tbranch->SetAddress(&iotrack);
50        tracktree->GetEvent(i);
51
52        Int_t tpcLabel=iotrack->GetLabel();
53        tracker.CookLabel(iotrack,0.);
54        Int_t itsLabel=iotrack->GetLabel();
55        if (itsLabel != tpcLabel) iotrack->SetLabel(-TMath::Abs(itsLabel));
56        if (tpcLabel < 0)         iotrack->SetLabel(-TMath::Abs(itsLabel));
57        /*
58        if (itsLabel==1234) {
59          Int_t nc=iotrack->GetNumberOfClusters();
60          for (Int_t k=0; k<nc; k++) {
61            Int_t index=iotrack->GetClusterIndex(k);
62            AliITSclusterV2 *c=tracker.GetCluster(index);
63            cout<<c->GetLabel(0)<<' '<<c->GetLabel(1)<<' '<<c->GetLabel(2)<<endl;
64          }
65        }
66        */
67        tarray.AddLast(iotrack);
68    }   
69    tf->Close();
70    cf->Close();
71
72 /////////////////////////////////////////////////////////////////////////
73    GoodTrack gt[15000];
74    Int_t ngood=0;
75    ifstream in("good_tracks_its");
76    if (in) {
77       cerr<<"Reading good tracks...\n";
78       while (in>>gt[ngood].lab>>gt[ngood].code>>
79                  gt[ngood].px>>gt[ngood].py>>gt[ngood].pz>>
80                  gt[ngood].x >>gt[ngood].y >>gt[ngood].z) {
81          ngood++;
82          cerr<<ngood<<'\r';
83          if (ngood==15000) {
84             cerr<<"Too many good tracks !\n";
85             break;
86          }
87       }
88       if (!in.eof()) cerr<<"Read error (good_tracks_its) !\n";
89    } else {
90       cerr<<"Marking good tracks (this will take a while)...\n";
91       ngood=good_tracks(gt,15000);
92       ofstream out("good_tracks_its");
93       if (out) {
94          for (Int_t ngd=0; ngd<ngood; ngd++)            
95             out<<gt[ngd].lab<<' '<<gt[ngd].code<<' '<<
96                  gt[ngd].px<<' '<<gt[ngd].py<<' '<<gt[ngd].pz<<' '<<
97                  gt[ngd].x <<' '<<gt[ngd].y <<' '<<gt[ngd].z <<endl;
98       } else cerr<<"Can not open file (good_tracks_its) !\n";
99       out.close();
100    }
101    cerr<<"Number of good tracks : "<<ngood<<endl;
102
103    TH1F *hp=new TH1F("hp","PHI resolution",50,-20.,20.); hp->SetFillColor(4);
104    TH1F *hl=new TH1F("hl","LAMBDA resolution",50,-20,20);hl->SetFillColor(4);
105    TH1F *hpt=new TH1F("hpt","Relative Pt resolution",30,-10.,10.); 
106    hpt->SetFillColor(2); 
107    TH1F *hmpt=new TH1F("hmpt","Transverse impact parameter",30,-300,300); 
108    hmpt->SetFillColor(6);
109    TH1F *hz=new TH1F("hz","Longitudinal impact parameter",30,-300,300); 
110    //hmpt->SetFillColor(6);
111
112    TH1F *hgood=new TH1F("hgood","Good tracks",30,0.1,6.1);    
113    TH1F *hfound=new TH1F("hfound","Found tracks",30,0.1,6.1);
114    TH1F *hfake=new TH1F("hfake","Fake tracks",30,0.1,6.1);
115    TH1F *hg=new TH1F("hg","",30,0.1,6.1); //efficiency for good tracks
116    hg->SetLineColor(4); hg->SetLineWidth(2);
117    TH1F *hf=new TH1F("hf","Efficiency for fake tracks",30,0.1,6.1);
118    hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
119
120    TH1F *hptw=new TH1F("hptw","Weghted pt",30,0.1,6.1);
121
122    while (ngood--) {
123       Int_t lab=gt[ngood].lab, tlab=-1;
124       Double_t pxg=gt[ngood].px, pyg=gt[ngood].py, pzg=gt[ngood].pz;
125       Double_t ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
126
127       hgood->Fill(ptg);
128
129       AliITStrackV2 *track=0;
130       Int_t j;
131       for (j=0; j<nentr; j++) {
132           track=(AliITStrackV2*)tarray.UncheckedAt(j);
133           tlab=track->GetLabel();
134           if (lab==TMath::Abs(tlab)) break;
135       }
136       if (j==nentr) {
137         cerr<<"Track "<<lab<<" was not found !\n";
138         continue;
139       }
140       track->Propagate(track->GetAlpha(),3.,0.1/65.19*1.848,0.1*1.848);
141       track->PropagateToVertex();
142
143       if (lab==tlab) hfound->Fill(ptg);
144       else { hfake->Fill(ptg); cerr<<lab<<" fake\n";}
145
146       Double_t xv,par[5]; track->GetExternalParameters(xv,par);
147       Float_t phi=TMath::ASin(par[2]) + track->GetAlpha();
148       if (phi<-TMath::Pi()) phi+=2*TMath::Pi();
149       if (phi>=TMath::Pi()) phi-=2*TMath::Pi();
150       Float_t lam=TMath::ATan(par[3]); 
151       Float_t pt_1=TMath::Abs(par[4]);
152
153       Double_t phig=TMath::ATan2(pyg,pxg);
154       hp->Fill((phi - phig)*1000.);
155
156       Double_t lamg=TMath::ATan2(pzg,ptg);
157       hl->Fill((lam - lamg)*1000.);
158
159       Double_t d=10000*track->GetD();
160       hmpt->Fill(d);
161
162       hptw->Fill(ptg,TMath::Abs(d));
163
164       Double_t z=10000*track->GetZ();
165       hz->Fill(z);
166
167 //if (TMath::Abs(gt[ngood].code)==11 && ptg>4.)
168       hpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
169
170    }
171
172    Stat_t ng=hgood->GetEntries(); cerr<<"Good tracks "<<ng<<endl;
173    Stat_t nf=hfound->GetEntries();
174    if (ng!=0) 
175       cerr<<"Integral efficiency is about "<<nf/ng*100.<<" %\n";
176
177    gStyle->SetOptStat(111110);
178    gStyle->SetOptFit(1);
179
180    TCanvas *c1=new TCanvas("c1","",0,0,700,850);
181
182    TPad *p1=new TPad("p1","",0,0.3,.5,.6); p1->Draw();
183    p1->cd(); p1->SetFillColor(42); p1->SetFrameFillColor(10); 
184    hp->SetFillColor(4);  hp->SetXTitle("(mrad)"); hp->Fit("gaus"); c1->cd();
185
186    TPad *p2=new TPad("p2","",0.5,.3,1,.6); p2->Draw(); 
187    p2->cd(); p2->SetFillColor(42); p2->SetFrameFillColor(10);
188    hl->SetXTitle("(mrad)"); hl->Fit("gaus"); c1->cd();
189
190    TPad *p3=new TPad("p3","",0,0,0.5,0.3); p3->Draw();
191    p3->cd(); p3->SetFillColor(42); p3->SetFrameFillColor(10); 
192    hpt->SetXTitle("(%)"); hpt->Fit("gaus"); c1->cd();
193
194    TPad *p4=new TPad("p4","",0.5,0,1,0.3); p4->Draw();
195    p4->cd(); p4->SetFillColor(42); p4->SetFrameFillColor(10);
196    hmpt->SetXTitle("(micron)"); hmpt->Fit("gaus"); hz->Draw("same"); c1->cd();
197    //hfound->Sumw2();
198    //hptw->Sumw2(); 
199    //hg->SetMaximum(333);
200    //hg->SetYTitle("Impact Parameter Resolution (micron)");
201    //hg->SetXTitle("Pt (GeV/c)");
202    //hg->GetXaxis()->SetRange(0,10);
203    //hg->Divide(hptw,hfound,1,1.);
204    //hg->DrawCopy(); c1->cd();
205    
206
207    TPad *p5=new TPad("p5","",0,0.6,1,1); p5->Draw(); p5->cd(); 
208    p5->SetFillColor(41); p5->SetFrameFillColor(10);
209    hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
210    hg->Divide(hfound,hgood,1,1.,"b");
211    hf->Divide(hfake,hgood,1,1.,"b");
212    hg->SetMaximum(1.4);
213    hg->SetYTitle("Tracking efficiency");
214    hg->SetXTitle("Pt (GeV/c)");
215    hg->Draw();
216
217    TLine *line1 = new TLine(0,1.0,7,1.0); line1->SetLineStyle(4);
218    line1->Draw("same");
219    TLine *line2 = new TLine(0,0.9,7,0.9); line2->SetLineStyle(4);
220    line2->Draw("same");
221
222    hf->SetFillColor(1);
223    hf->SetFillStyle(3013);
224    hf->SetLineColor(2);
225    hf->SetLineWidth(2);
226    hf->Draw("histsame");
227    TText *text = new TText(0.461176,0.248448,"Fake tracks");
228    text->SetTextSize(0.05);
229    text->Draw();
230    text = new TText(0.453919,1.11408,"Good tracks");
231    text->SetTextSize(0.05);
232    text->Draw();
233
234    return 0;
235 }
236
237 Int_t good_tracks(GoodTrack *gt, Int_t max) {
238    if (gAlice) {delete gAlice; gAlice=0;}
239
240    TFile *file=TFile::Open("galice.root");
241    if (!file->IsOpen()) {cerr<<"Can't open galice.root !\n"; exit(4);}
242    if (!(gAlice=(AliRun*)file->Get("gAlice"))) {
243      cerr<<"gAlice have not been found on galice.root !\n";
244      exit(5);
245    }
246
247    Int_t np=gAlice->GetEvent(0);
248
249    Int_t *good=new Int_t[np];
250    Int_t k;
251    for (k=0; k<np; k++) good[k]=0;
252
253    AliITS *ITS=(AliITS*)gAlice->GetDetector("ITS");
254    if (!ITS) {
255       cerr<<"can't get ITS !\n"; exit(8);
256    }
257    AliITSgeom *geom=ITS->GetITSgeom();
258    if (!geom) {
259       cerr<<"cen't get ITS geometry !\n"; exit(9);
260    }
261
262    TFile *cf=TFile::Open("AliITSclustersV2.root");
263    if (!cf->IsOpen()){
264       cerr<<"Can't open AliITSclustersV2.root !\n"; exit(6);
265    }
266    TTree *cTree=(TTree*)cf->Get("cTree");
267    if (!cTree) {
268       cerr<<"Can't get cTree !\n"; exit(7);
269    }
270    TBranch *branch=cTree->GetBranch("Clusters");
271    if (!branch) {
272       cerr<<"Can't get clusters branch !\n"; exit(8);
273    }
274    TClonesArray *clusters=new TClonesArray("AliITSclusterV2",10000);
275    branch->SetAddress(&clusters);
276
277    Int_t entr=(Int_t)cTree->GetEntries();
278    for (k=0; k<entr; k++) {
279      if (!cTree->GetEvent(k)) continue;
280      Int_t lay,lad,det;  geom->GetModuleId(k-1,lay,lad,det);
281      if (lay<1 || lay>6) {
282         cerr<<"wrong layer !\n"; exit(10);
283      }
284      Int_t ncl=clusters->GetEntriesFast();
285      while (ncl--) {
286         AliITSclusterV2 *pnt=(AliITSclusterV2*)clusters->UncheckedAt(ncl);
287         Int_t l0=pnt->GetLabel(0);
288         Int_t l1=pnt->GetLabel(1);
289         Int_t l2=pnt->GetLabel(2);
290         Int_t mask=1<<(lay-1);
291         if (l0>=0) good[l0]|=mask; 
292         if (l1>=0) good[l1]|=mask; 
293         if (l2>=0) good[l2]|=mask;
294      }
295    }
296    clusters->Delete(); delete clusters;
297    cf->Close();
298
299    ifstream in("good_tracks_tpc");
300    if (!in) {
301      cerr<<"can't get good_tracks_tpc !\n"; exit(11);
302    }
303    Int_t nt=0;
304    Double_t px,py,pz,x,y,z;
305    Int_t code,lab;
306    while (in>>lab>>code>>px>>py>>pz>>x>>y>>z) {
307       if (good[lab] != 0x3F) continue;
308       TParticle *p = (TParticle*)gAlice->Particle(lab);
309       gt[nt].lab=lab;
310       gt[nt].code=p->GetPdgCode();
311 //**** px py pz - in global coordinate system
312       gt[nt].px=p->Px(); gt[nt].py=p->Py(); gt[nt].pz=p->Pz();
313       gt[nt].x=gt[nt].y=gt[nt].z=0.;
314       nt++;
315       if (nt==max) {cerr<<"Too many good tracks !\n"; break;}
316    }
317
318    delete[] good;
319
320    delete gAlice; gAlice=0;
321    file->Close();
322
323    return nt;
324 }
325
326