]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliV0Comparison.C
Updated V2 stream of tracking (Yu.Belikov). The new long waited features are: 1)...
[u/mrichter/AliRoot.git] / ITS / AliV0Comparison.C
1 #ifndef __CINT__
2   #include <iostream.h>
3   #include <fstream.h>
4
5   #include "TH1.h"
6   #include "TFile.h"
7   #include "TTree.h"
8   #include "TObjArray.h"
9   #include "TStyle.h"
10   #include "TCanvas.h"
11   #include "TLine.h"
12   #include "TText.h"
13   #include "TParticle.h"
14   #include "TStopwatch.h"
15
16   #include "AliRun.h"
17   #include "AliPDG.h"
18   #include "AliV0vertex.h"
19 #endif
20
21 struct GoodVertex {
22   Int_t nlab,plab;
23   Int_t code;
24   Float_t px,py,pz;
25   Float_t x,y,z;
26 };
27 Int_t good_vertices(GoodVertex *gt, Int_t max);
28
29 Int_t AliV0Comparison(Int_t code=310) { //Lambda=3122, LambdaBar=-3122
30    cerr<<"Doing comparison...\n";
31
32    const Double_t V0window=0.05, V0width=0.015; 
33    Double_t V0mass=0.5;
34    switch (code) {
35    case kK0Short:    V0mass=0.497672; break;
36    case kLambda0:    V0mass=1.115683; break;
37    case kLambda0Bar: V0mass=1.115683; break;
38    default: cerr<<"Invalid PDG code !\n"; return 1;
39    }
40
41    TStopwatch timer;
42
43    /*** Load reconstructed vertices ***/
44    TFile *vf=TFile::Open("AliV0vertices.root");
45    if (!vf->IsOpen()) {cerr<<"Can't open AliV0vertices.root !\n"; return 2;}
46    TObjArray varray(1000);
47    TTree *vTree=(TTree*)vf->Get("TreeV");
48    TBranch *branch=vTree->GetBranch("vertices");
49    Int_t nentr=(Int_t)vTree->GetEntries();
50    for (Int_t i=0; i<nentr; i++) {
51        AliV0vertex *iovertex=new AliV0vertex; branch->SetAddress(&iovertex);
52        vTree->GetEvent(i);
53        varray.AddLast(iovertex);
54    }
55
56    /*** Check if the file with the "good" vertices exists ***/
57    GoodVertex gv[1000];
58    Int_t ngood=0;
59    ifstream in("good_vertices");
60    if (in) {
61       cerr<<"Reading good vertices...\n";
62       while (in>>gv[ngood].nlab>>gv[ngood].plab>>gv[ngood].code>>
63                  gv[ngood].px>>gv[ngood].py>>gv[ngood].pz>>
64                  gv[ngood].x >>gv[ngood].y >>gv[ngood].z) {
65          ngood++;
66          cerr<<ngood<<'\r';
67          if (ngood==1000) {
68             cerr<<"Too many good vertices !\n";
69             break;
70          }
71       }
72       if (!in.eof()) cerr<<"Read error (good_vertices) !\n";
73    } else {
74      /*** generate a file with the "good" vertices ***/
75       cerr<<"Marking good vertices (this will take a while)...\n";
76       ngood=good_vertices(gv,1000);
77       ofstream out("good_vertices");
78       if (out) {
79          for (Int_t ngd=0; ngd<ngood; ngd++)            
80             out<<gv[ngd].nlab<<' '<<gv[ngd].plab<<' '<<gv[ngd].code<<' '<<
81                  gv[ngd].px<<' '<<gv[ngd].py<<' '<<gv[ngd].pz<<' '<<
82                  gv[ngd].x <<' '<<gv[ngd].y <<' '<<gv[ngd].z <<endl;
83       } else cerr<<"Can not open file (good_vertices) !\n";
84       out.close();
85    }
86
87    vf->Close();
88
89
90    TH1F *hp=new TH1F("hp","Angular Resolution",30,-30.,30.); //phi resolution 
91    hp->SetXTitle("(mrad)"); hp->SetFillColor(2);
92    TH1F *hl=new TH1F("hl","Lambda Resolution",30,-30,30);
93    hl->SetXTitle("(mrad)"); hl->SetFillColor(1); hl->SetFillStyle(3013); 
94    TH1F *hpt=new TH1F("hpt","Relative Pt Resolution",30,-10.,10.); 
95    hpt->SetXTitle("(%)"); hpt->SetFillColor(2); 
96
97    TH1F *hx=new TH1F("hx","Position Resolution (X,Y)",30,-3.,3.); //x res. 
98    hx->SetXTitle("(mm)"); hx->SetFillColor(6);
99    TH1F *hy=new TH1F("hy","Position Resolution (Y)",30,-3.,3.);   //y res
100    hy->SetXTitle("(mm)"); hy->SetFillColor(1); hy->SetFillStyle(3013);
101    TH1F *hz=new TH1F("hz","Position Resolution (Z)",30,-3.,3.);   //z res. 
102    hz->SetXTitle("(mm)"); hz->SetFillColor(6);
103
104
105    Double_t pmin=0.2, pmax=4.2; Int_t nchan=20;
106    TH1F *hgood=new TH1F("hgood","Good Vertices",nchan,pmin,pmax);    
107    TH1F *hfound=new TH1F("hfound","Found Vertices",nchan,pmin,pmax);
108    TH1F *hfake=new TH1F("hfake","Fake Vertices",nchan,pmin,pmax);
109    TH1F *hg=new TH1F("hg","Efficiency for Good Vertices",nchan,pmin,pmax);
110    hg->SetLineColor(4); hg->SetLineWidth(2);
111    TH1F *hf=new TH1F("hf","Probability of Fake Vertices",nchan,pmin,pmax);
112    hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
113
114    Double_t mmin=V0mass-V0window, mmax=V0mass+V0window;
115    TH1F *v0s =new TH1F("v0s","V0s Effective Mass",40, mmin, mmax);
116    v0s->SetXTitle("(GeV/c)"); v0s->SetFillColor(6);
117
118
119    Double_t pxg=0.,pyg=0.,ptg=0.;
120    Int_t nlab=-1, plab=-1;
121    Int_t i;
122    for (i=0; i<nentr; i++) {
123        AliV0vertex *vertex=(AliV0vertex*)varray.UncheckedAt(i);
124        nlab=TMath::Abs(vertex->GetNlabel());
125        plab=TMath::Abs(vertex->GetPlabel());
126
127        vertex->ChangeMassHypothesis(code);
128
129        Double_t mass=vertex->GetEffMass();
130        v0s->Fill(mass);
131
132        Int_t j;
133        for (j=0; j<ngood; j++) {
134           if (gv[j].code != vertex->GetPdgCode()) continue;
135           if (gv[j].nlab == nlab)
136           if (gv[j].plab == plab) break;
137        }
138
139        if (TMath::Abs(mass-V0mass)>V0width) continue;
140
141        Double_t px,py,pz; vertex->GetPxPyPz(px,py,pz);
142        Double_t pt=TMath::Sqrt(px*px+py*py);
143
144        if (j==ngood) {
145           hfake->Fill(pt);
146           cerr<<"Fake vertex: ("<<nlab<<","<<plab<<")\n";
147           continue;
148        }
149
150        pxg=gv[j].px; pyg=gv[j].py; ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
151        Double_t phig=TMath::ATan2(pyg,pxg), phi=TMath::ATan2(py,px);
152        Double_t lamg=TMath::ATan2(gv[j].pz,ptg), lam=TMath::ATan2(pz,pt);
153        hp->Fill((phi - phig)*1000.);
154        hl->Fill((lam - lamg)*1000.);
155        hpt->Fill((1/pt - 1/ptg)/(1/ptg)*100.);
156
157        Double_t x,y,z; vertex->GetXYZ(x,y,z);
158        hx->Fill((x-gv[j].x)*10);
159        hy->Fill((y-gv[j].y)*10);
160        hz->Fill((z-gv[j].z)*10);
161
162        hfound->Fill(ptg);
163        gv[j].nlab=-1;
164
165    }
166    for (i=0; i<ngood; i++) {
167       if (gv[i].code != code) continue;
168       pxg=gv[i].px; pyg=gv[i].py; ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
169       hgood->Fill(ptg);
170       nlab=gv[i].nlab; plab=gv[i].plab;
171       if (nlab < 0) continue;
172       cerr<<"Vertex ("<<nlab<<','<<plab<<") has not been found !\n";
173    }
174
175    varray.Delete();
176
177    Stat_t ng=hgood->GetEntries();
178    Stat_t nf=hfound->GetEntries();
179
180    cerr<<"Number of found vertices: "<<nentr<<" ("<<nf<<" in the peak)\n";
181    cerr<<"Number of good vertices: "<<ng<<endl;
182
183    if (ng!=0) 
184       cerr<<"Integral efficiency is about "<<nf/ng*100.<<" %\n";
185
186    gStyle->SetOptStat(111110);
187    gStyle->SetOptFit(1);
188
189    TCanvas *c1=new TCanvas("c1","",0,0,580,610);
190    c1->Divide(2,2);
191
192    c1->cd(1);
193    gPad->SetFillColor(42); gPad->SetFrameFillColor(10); 
194    //hp->Fit("gaus");
195    hp->Draw();
196    hl->Draw("same"); c1->cd();
197
198    c1->cd(2);
199    gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
200    //hpt->Fit("gaus"); c1->cd();
201    hpt->Draw(); c1->cd();
202
203    c1->cd(3);
204    gPad->SetFillColor(42); gPad->SetFrameFillColor(10); 
205    //hx->Fit("gaus"); 
206    hx->Draw(); 
207    hy->Draw("same"); c1->cd();
208
209    c1->cd(4);
210    gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
211    //hz->Fit("gaus");
212    hz->Draw();
213
214
215    TCanvas *c2=new TCanvas("c2","",600,0,580,610);
216    c2->Divide(1,2);
217
218    c2->cd(1);
219    gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
220    hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
221    hg->Divide(hfound,hgood,1,1.,"b");
222    hf->Divide(hfake,hgood,1,1.,"b");
223    hg->SetMaximum(1.4);
224    hg->SetYTitle("V0 reconstruction efficiency");
225    hg->SetXTitle("Pt (GeV/c)");
226    hg->Draw();
227
228    TLine *line1 = new TLine(pmin,1.0,pmax,1.0); line1->SetLineStyle(4);
229    line1->Draw("same");
230    TLine *line2 = new TLine(pmin,0.9,pmax,0.9); line2->SetLineStyle(4);
231    line2->Draw("same");
232
233    hf->SetFillColor(1);
234    hf->SetFillStyle(3013);
235    hf->SetLineColor(2);
236    hf->SetLineWidth(2);
237    hf->Draw("histsame");
238    TText *text = new TText(0.461176,0.248448,"Fake vertices");
239    text->SetTextSize(0.05);
240    text->Draw();
241    text = new TText(0.453919,1.11408,"Good vertices");
242    text->SetTextSize(0.05);
243    text->Draw();
244
245
246    c2->cd(2);
247    gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
248    v0s->SetXTitle("(GeV/c)"); 
249    v0s->Fit("gaus","","",V0mass-V0width,V0mass+V0width);
250    Double_t max=v0s->GetMaximum();
251    TLine *line3 = new TLine(V0mass-V0width,0.,V0mass-V0width,max);
252    line3->Draw("same");
253    TLine *line4 = new TLine(V0mass+V0width,0.,V0mass+V0width,max);
254    line4->Draw("same");
255
256    timer.Stop(); timer.Print();
257
258    return 0;
259 }
260
261 Int_t good_vertices(GoodVertex *gv, Int_t max) {
262    Int_t nv=0;
263    /*** Get information about the cuts ***/
264    Double_t r2min=0.9*0.9;
265    Double_t r2max=2.9*2.9;
266
267    /*** Get labels of the "good" tracks ***/
268    Double_t dd; Int_t id, label[15000], ngt=0;
269    ifstream in("good_tracks_its");
270    if (!in) {
271      cerr<<"Can't open the file good_tracks_its \n";
272      return nv;
273    }
274    while (in>>label[ngt]>>id>>dd>>dd>>dd>>dd>>dd>>dd) {
275      ngt++;
276      if (ngt>=15000) {
277        cerr<<"Too many good ITS tracks !\n";
278        return nv;
279      }
280    }   
281    if (!in.eof()) {
282       cerr<<"Read error (good_tracks_its) !\n";
283       return nv;
284    }
285
286    /*** Get an access to the kinematics ***/
287    if (gAlice) {delete gAlice; gAlice=0;}
288
289    TFile *file=TFile::Open("galice.root");
290    if (!file->IsOpen()) {cerr<<"Can't open galice.root !\n"; exit(4);}
291    if (!(gAlice=(AliRun*)file->Get("gAlice"))) {
292      cerr<<"gAlice has not been found on galice.root !\n";
293      exit(5);
294    }
295
296    Int_t np=gAlice->GetEvent(0);
297    while (np--) {
298       cerr<<np<<'\r';
299       TParticle *p0=gAlice->Particle(np);
300
301       /*** only these V0s are "good" ***/
302       Int_t code=p0->GetPdgCode();
303       if (code!=kK0Short)
304       if (code!=kLambda0)
305       if (code!=kLambda0Bar) continue; 
306
307       /*** daughter tracks should be "good" ***/
308       Int_t plab=p0->GetFirstDaughter(), nlab=p0->GetLastDaughter();
309       if (nlab==plab) continue;
310       if (nlab<0) continue;
311       if (plab<0) continue;
312       Int_t i;
313       for (i=0; i<ngt; i++) if (label[i]==nlab) break;
314       if (i==ngt) continue;
315       for (i=0; i<ngt; i++) if (label[i]==plab) break;
316       if (i==ngt) continue;
317
318       /*** fiducial volume ***/
319       TParticle *p=gAlice->Particle(nlab);
320       Double_t x=p->Vx(), y=p->Vy(), z=p->Vz(), r2=x*x+y*y;
321       if (r2<r2min) continue;
322       if (r2>r2max) continue;
323  
324       gv[nv].code=code;
325       gv[nv].plab=plab; gv[nv].nlab=nlab;
326       gv[nv].px=p0->Px(); gv[nv].py=p0->Py(); gv[nv].pz=p0->Pz();
327       gv[nv].x=x; gv[nv].y=y; gv[nv].z=z;
328       nv++;
329    }
330  
331    delete gAlice; gAlice=0;
332
333    file->Close();  
334
335    return nv;
336 }