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