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