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