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