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