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