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