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