Comparison working with ESD (Yu.Belikov)
[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() {
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         if ((t->GetStatus()&AliESDtrack::kITSin)==0) continue;
124         AliITStrackV2 *iotrack=new AliITStrackV2(*t);
125         tarray.AddLast(iotrack);
126      }
127      delete event;
128    }
129    ef->Close();
130    }
131    Int_t nentr=tarray.GetEntriesFast();
132
133    TH1F *hp=new TH1F("hp","PHI resolution",50,-20.,20.); hp->SetFillColor(4);
134    TH1F *hl=new TH1F("hl","LAMBDA resolution",50,-20,20);hl->SetFillColor(4);
135    TH1F *hpt=new TH1F("hpt","Relative Pt resolution",30,-10.,10.); 
136    hpt->SetFillColor(2); 
137    TH1F *hmpt=new TH1F("hmpt","Transverse impact parameter",30,-300,300); 
138    hmpt->SetFillColor(6);
139    TH1F *hz=new TH1F("hz","Longitudinal impact parameter",30,-300,300); 
140
141    AliITStrackV2 *trk=(AliITStrackV2*)tarray.UncheckedAt(0);
142    Double_t pmin=0.1*(100/0.299792458/0.2/trk->GetConvConst());
143    Double_t pmax=6.0+pmin;
144
145    TH1F *hgood=new TH1F("hgood","Good tracks",30,pmin,pmax);    
146    TH1F *hfound=new TH1F("hfound","Found tracks",30,pmin,pmax);
147    TH1F *hfake=new TH1F("hfake","Fake tracks",30,pmin,pmax);
148    TH1F *hg=new TH1F("hg","",30,pmin,pmax); //efficiency for good tracks
149    hg->SetLineColor(4); hg->SetLineWidth(2);
150    TH1F *hf=new TH1F("hf","Efficiency for fake tracks",30,pmin,pmax);
151    hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
152
153    //TH1F *hptw=new TH1F("hptw","Weghted pt",30,pmax,pmin);
154
155    TH1F *he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,100.);
156    TH2F *hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,200.);
157    hep->SetMarkerStyle(8);
158    hep->SetMarkerSize(0.4);
159
160
161    Int_t notf[MAX], nnotf=0;
162    Int_t fake[MAX], nfake=0;
163    Int_t mult[MAX], numb[MAX], nmult=0;
164
165    Int_t ng;
166    for (ng=0; ng<ngood; ng++) {
167       Int_t lab=gt[ng].lab, tlab=-1;
168       Double_t pxg=gt[ng].px, pyg=gt[ng].py, pzg=gt[ng].pz;
169       Double_t ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
170
171
172       if (ptg>pmin) hgood->Fill(ptg);
173
174       AliITStrackV2 *track=0;
175       Int_t cnt=0;
176       Int_t j;
177       for (j=0; j<nentr; j++) {
178           AliITStrackV2 *trk=(AliITStrackV2*)tarray.UncheckedAt(j);
179           Int_t lbl=trk->GetLabel();
180           if (lab==TMath::Abs(lbl)) {
181             if (cnt==0) {track=trk; tlab=lbl;}
182             cnt++;
183           }
184       }
185       if (cnt==0) {
186         notf[nnotf++]=lab;
187         continue;
188       } else if (cnt>1){
189         mult[nmult]=lab;
190         numb[nmult]=cnt; nmult++;        
191       }
192
193       if (ptg>pmin) {
194         if (lab==tlab) hfound->Fill(ptg);
195         else {
196           fake[nfake++]=lab;
197           hfake->Fill(ptg); 
198         }
199       }
200
201       track->PropagateTo(3.,0.0028,65.19);
202       track->PropagateToVertex();
203
204       Double_t xv,par[5]; track->GetExternalParameters(xv,par);
205       Float_t phi=TMath::ASin(par[2]) + track->GetAlpha();
206       if (phi<-TMath::Pi()) phi+=2*TMath::Pi();
207       if (phi>=TMath::Pi()) phi-=2*TMath::Pi();
208       Float_t lam=TMath::ATan(par[3]); 
209       Float_t pt_1=TMath::Abs(par[4]);
210
211       Double_t phig=TMath::ATan2(pyg,pxg);
212       hp->Fill((phi - phig)*1000.);
213
214       Double_t lamg=TMath::ATan2(pzg,ptg);
215       hl->Fill((lam - lamg)*1000.);
216
217       Double_t d=10000*track->GetD();
218       hmpt->Fill(d);
219
220       //hptw->Fill(ptg,TMath::Abs(d));
221
222       Double_t z=10000*track->GetZ();
223       hz->Fill(z);
224
225       hpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
226
227       Float_t mom=1./(pt_1*TMath::Cos(lam));
228       Float_t dedx=track->GetdEdx();
229       hep->Fill(mom,dedx,1.);
230       if (TMath::Abs(gt[ng].code)==211)
231          if (mom>0.4 && mom<0.5) he->Fill(dedx,1.);
232    }
233
234
235    cout<<"\nList of Not found tracks :\n";
236    for (ng=0; ng<nnotf; ng++){
237      cout<<notf[ng]<<"\t";
238      if ((ng%9)==8) cout<<"\n";
239    }
240    cout<<"\n\nList of fake  tracks :\n";
241    for (ng=0; ng<nfake; ng++){
242      cout<<fake[ng]<<"\t";
243      if ((ng%9)==8) cout<<"\n";
244    }
245    cout<<"\n\nList of multiple found tracks :\n";
246    for (ng=0; ng<nmult; ng++) {
247        cout<<"id.   "<<mult[ng]
248            <<"     found - "<<numb[ng]<<"times\n";
249    }
250    cout<<endl;
251
252    cerr<<"Number of good tracks : "<<ngood<<endl;
253    cerr<<"Number of found tracks : "<<nentr<<endl;
254
255    ng=(Int_t)hgood->GetEntries(); //cerr<<"Good tracks "<<ng<<endl;
256    if (ng!=0)  
257    cerr<<"Integral efficiency is about "<<hfound->GetEntries()/ng*100.<<" %\n";
258    cerr<<endl;
259
260    gStyle->SetOptStat(111110);
261    gStyle->SetOptFit(1);
262
263    TCanvas *c1=new TCanvas("c1","",0,0,700,850);
264
265    TPad *p1=new TPad("p1","",0,0.3,.5,.6); p1->Draw();
266    p1->cd(); p1->SetFillColor(42); p1->SetFrameFillColor(10); 
267    hp->SetFillColor(4);  hp->SetXTitle("(mrad)"); hp->Fit("gaus"); c1->cd();
268
269    TPad *p2=new TPad("p2","",0.5,.3,1,.6); p2->Draw(); 
270    p2->cd(); p2->SetFillColor(42); p2->SetFrameFillColor(10);
271    hl->SetXTitle("(mrad)"); hl->Fit("gaus"); c1->cd();
272
273    TPad *p3=new TPad("p3","",0,0,0.5,0.3); p3->Draw();
274    p3->cd(); p3->SetFillColor(42); p3->SetFrameFillColor(10); 
275    hpt->SetXTitle("(%)"); hpt->Fit("gaus"); c1->cd();
276
277    TPad *p4=new TPad("p4","",0.5,0,1,0.3); p4->Draw();
278    p4->cd(); p4->SetFillColor(42); p4->SetFrameFillColor(10);
279    hmpt->SetXTitle("(micron)"); hmpt->Fit("gaus"); hz->Draw("same"); c1->cd();
280    //hfound->Sumw2();
281    //hptw->Sumw2(); 
282    //hg->SetMaximum(333);
283    //hg->SetYTitle("Impact Parameter Resolution (micron)");
284    //hg->SetXTitle("Pt (GeV/c)");
285    //hg->GetXaxis()->SetRange(0,10);
286    //hg->Divide(hptw,hfound,1,1.);
287    //hg->DrawCopy(); c1->cd();
288    
289
290    TPad *p5=new TPad("p5","",0,0.6,1,1); p5->Draw(); p5->cd(); 
291    p5->SetFillColor(41); p5->SetFrameFillColor(10);
292    hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
293    hg->Divide(hfound,hgood,1,1.,"b");
294    hf->Divide(hfake,hgood,1,1.,"b");
295    hg->SetMaximum(1.4);
296    hg->SetYTitle("Tracking efficiency");
297    hg->SetXTitle("Pt (GeV/c)");
298    hg->Draw();
299
300    TLine *line1 = new TLine(pmin,1.0,pmax,1.0); line1->SetLineStyle(4);
301    line1->Draw("same");
302    TLine *line2 = new TLine(pmin,0.9,pmax,0.9); line2->SetLineStyle(4);
303    line2->Draw("same");
304
305    hf->SetFillColor(1);
306    hf->SetFillStyle(3013);
307    hf->SetLineColor(2);
308    hf->SetLineWidth(2);
309    hf->Draw("histsame");
310    TText *text = new TText(0.461176,0.248448,"Fake tracks");
311    text->SetTextSize(0.05);
312    text->Draw();
313    text = new TText(0.453919,1.11408,"Good tracks");
314    text->SetTextSize(0.05);
315    text->Draw();
316
317    TCanvas *c2=new TCanvas("c2","",320,32,530,590);
318    TPad *p6=new TPad("p6","",0.,0.,1.,.5); p6->Draw();
319    p6->cd(); p6->SetFillColor(42); p6->SetFrameFillColor(10); 
320    he->SetFillColor(2); he->SetFillStyle(3005);  
321    he->SetXTitle("Arbitrary Units"); 
322    he->Fit("gaus"); c2->cd();
323
324    TPad *p7=new TPad("p7","",0.,0.5,1.,1.); p7->Draw(); 
325    p7->cd(); p7->SetFillColor(42); p7->SetFrameFillColor(10);
326    hep->SetXTitle("p (Gev/c)"); hep->SetYTitle("dE/dX (Arb. Units)"); 
327    hep->Draw(); c1->cd();
328    
329    return 0;
330 }
331
332 Int_t good_tracks_its(GoodTrackITS *gt, const Int_t max, const char* evfoldname) {
333    AliRunLoader* rl = AliRunLoader::GetRunLoader(evfoldname);
334    if (rl == 0x0) {
335       ::Fatal("AliITSComparisonV2.C::good_tracks_its",
336               "Can not find Run Loader in Folder Named %s",
337               evfoldname);
338    }
339
340    AliITSLoader* itsl = (AliITSLoader*)rl->GetLoader("ITSLoader");
341    if (itsl == 0x0) {
342        cerr<<"AliITSComparisonV2.C : Can not find ITSLoader\n";
343        delete rl;
344        return 3;
345    }
346    
347    rl->LoadgAlice();
348    rl->LoadHeader();
349    Int_t np = rl->GetHeader()->GetNtrack();
350
351    Int_t *good=new Int_t[np];
352    Int_t k;
353    for (k=0; k<np; k++) good[k]=0;
354
355    AliITS *ITS=(AliITS*)rl->GetAliRun()->GetDetector("ITS");
356    if (!ITS) {
357       cerr<<"can't get ITS !\n"; exit(8);
358    }
359    AliITSgeom *geom=ITS->GetITSgeom();
360    if (!geom) {
361       cerr<<"can't get ITS geometry !\n"; exit(9);
362    }
363
364    itsl->LoadRecPoints();
365    TTree *cTree=itsl->TreeR();
366    if (!cTree) {
367       cerr<<"Can't get cTree !\n"; exit(7);
368    }
369    TBranch *branch=cTree->GetBranch("Clusters");
370    if (!branch) {
371       cerr<<"Can't get clusters branch !\n"; exit(8);
372    }
373    TClonesArray *clusters=new TClonesArray("AliITSclusterV2",10000);
374    branch->SetAddress(&clusters);
375
376    Int_t entr=(Int_t)cTree->GetEntries();
377    for (k=0; k<entr; k++) {
378      cTree->GetEvent(k);
379      Int_t ncl=clusters->GetEntriesFast(); if (ncl==0) continue;
380      Int_t lay,lad,det;  geom->GetModuleId(k,lay,lad,det);
381      if (lay<1 || lay>6) {
382         cerr<<"wrong layer !\n"; exit(10);
383      }
384      while (ncl--) {
385         AliITSclusterV2 *pnt=(AliITSclusterV2*)clusters->UncheckedAt(ncl);
386         Int_t l0=pnt->GetLabel(0);
387           if (l0>=np) {cerr<<"Wrong label: "<<l0<<endl; continue;}
388         Int_t l1=pnt->GetLabel(1);
389           if (l1>=np) {cerr<<"Wrong label: "<<l1<<endl; continue;}
390         Int_t l2=pnt->GetLabel(2);
391           if (l2>=np) {cerr<<"Wrong label: "<<l2<<endl; continue;}
392         Int_t mask=1<<(lay-1);
393         if (l0>=0) good[l0]|=mask; 
394         if (l1>=0) good[l1]|=mask; 
395         if (l2>=0) good[l2]|=mask;
396      }
397    }
398    clusters->Delete(); delete clusters;
399    itsl->UnloadRawClusters();
400    
401    ifstream in("good_tracks_tpc");
402    if (!in) {
403      cerr<<"can't get good_tracks_tpc !\n"; exit(11);
404    }
405    
406    rl->LoadKinematics();
407    AliStack* stack = rl->Stack();
408    Int_t nt=0;
409    Double_t px,py,pz,x,y,z;
410    Int_t code,lab;
411    while (in>>lab>>code>>px>>py>>pz>>x>>y>>z) {
412       if (good[lab] != 0x3F) continue;
413       TParticle *p = (TParticle*)stack->Particle(lab);
414       if (p == 0x0) {
415          cerr<<"Can not get particle "<<lab<<endl;
416          nt++;
417          if (nt==max) {cerr<<"Too many good tracks !\n"; break;}
418          continue;
419       }
420       gt[nt].lab=lab;
421       gt[nt].code=p->GetPdgCode();
422 //**** px py pz - in global coordinate system
423       gt[nt].px=p->Px(); gt[nt].py=p->Py(); gt[nt].pz=p->Pz();
424       gt[nt].x=gt[nt].y=gt[nt].z=0.;
425       nt++;
426       if (nt==max) {cerr<<"Too many good tracks !\n"; break;}
427    }
428
429    delete[] good;
430
431    rl->UnloadKinematics();
432    rl->UnloadHeader();
433    rl->UnloadgAlice();
434
435    return nt;
436 }
437
438