]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/AliCascadeComparison.C
Fix for transient fSDigits, AliITSRawStream classes adapted to changed AliRawReader...
[u/mrichter/AliRoot.git] / ITS / AliCascadeComparison.C
1 /****************************************************************************
2  *           Very important, delicate and rather obscure macro.             *
3  *                                                                          *
4  *               Creates list of "findable" cascades,                       *
5  *             calculates efficiency, resolutions etc.                      *
6  *                                                                          *
7  *  Origin: Christian Kuhn, IReS, Strasbourg, christian.kuhn@ires.in2p3.fr  *
8  ****************************************************************************/
9
10 #ifndef __CINT__
11   #include <Riostream.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   #include "TPDGCode.h"
25
26   #include "AliRun.h"
27   #include "AliPDG.h"
28   #include "AliCascadeVertex.h"
29 #endif
30
31 struct GoodCascade {
32   Int_t nlab,plab;   // V0's daughter labels
33   Int_t blab;        // Bachelor label
34   Int_t code;
35   Float_t px,py,pz;
36   Float_t x,y,z;
37 };
38 Int_t good_cascades(GoodCascade *gt, Int_t max);
39
40 Int_t AliCascadeComparison(Int_t code=3312) {
41   //code= 3312; //kXiMinus 
42   //code=-3312; //kXiPlusBar
43   //code= 3334; //kOmegaMinus
44   //code=-3334; //kOmegaPlusBar
45
46    cerr<<"Doing comparison...\n";
47
48    const Double_t cascadeWindow=0.05, cascadeWidth=0.015; 
49    Double_t ptncut=0.12, ptpcut=0.33, kine0cut=0.003;
50    Double_t ptbcut=0.11, kinecut=0.002;
51    Double_t cascadeMass=1.32131;
52    switch (code) {
53    case kXiMinus:
54         break;
55    case kXiPlusBar:
56         ptncut=0.33; ptpcut=0.12;
57         break;
58    case kOmegaMinus: 
59         cascadeMass=1.67245;
60         kine0cut=0.001;
61         ptbcut=0.22; kinecut=0.006;
62         break; 
63    case kOmegaPlusBar:
64         cascadeMass=1.67245; 
65         kine0cut=0.001;
66         ptncut=0.33; ptpcut=0.12;
67         ptbcut=0.22; kinecut=0.006;
68         break;
69    default: cerr<<"Invalid PDG code !\n"; return 1;
70    }
71
72    TStopwatch timer;
73
74    /*** Load reconstructed cascades ***/
75    TFile *cf=TFile::Open("AliCascadeVertices.root");
76    if (!cf->IsOpen()){cerr<<"Can't open AliCascadeVertices.root !\n";return 2;}
77    TObjArray carray(1000);
78    TTree *cTree=(TTree*)cf->Get("TreeCasc");
79    TBranch *branch=cTree->GetBranch("cascades");
80    Int_t nentr=(Int_t)cTree->GetEntries();
81    for (Int_t i=0; i<nentr; i++) {
82        AliCascadeVertex *iovertex=new AliCascadeVertex; 
83        branch->SetAddress(&iovertex);
84        cTree->GetEvent(i);
85        carray.AddLast(iovertex);
86    }
87    delete cTree;
88    cf->Close();
89
90    /*** Check if the file with the "good" cascades exists ***/
91    GoodCascade gc[100];
92    Int_t ngood=0;
93    ifstream in("good_cascades");
94    if (in) {
95       cerr<<"Reading good cascades...\n";
96       while (in>>gc[ngood].nlab>>gc[ngood].plab>>
97                  gc[ngood].blab>>gc[ngood].code>>
98                  gc[ngood].px>>gc[ngood].py>>gc[ngood].pz>>
99                  gc[ngood].x >>gc[ngood].y >>gc[ngood].z) {
100          ngood++;
101          cerr<<ngood<<'\r';
102          if (ngood==100) {
103             cerr<<"Too many good cascades !\n";
104             break;
105          }
106       }
107       if (!in.eof()) cerr<<"Read error (good_cascades) !\n";
108    } else {
109      /*** generate a file with the "good" cascades ***/
110       cerr<<"Marking good cascades (this will take a while)...\n";
111       ngood=good_cascades(gc,100);
112       ofstream out("good_cascades");
113       if (out) {
114          for (Int_t ngd=0; ngd<ngood; ngd++)            
115             out<<gc[ngd].nlab<<' '<<gc[ngd].plab<<' '<<
116                  gc[ngd].blab<<' '<<gc[ngd].code<<' '<<
117                  gc[ngd].px<<' '<<gc[ngd].py<<' '<<gc[ngd].pz<<' '<<
118                  gc[ngd].x <<' '<<gc[ngd].y <<' '<<gc[ngd].z <<endl;
119       } else cerr<<"Can not open file (good_cascades) !\n";
120       out.close();
121    }
122
123    /*** create some histograms ***/
124    TH1F *hp=new TH1F("hp","Angular Resolution",30,-30.,30.); //phi resolution 
125    hp->SetXTitle("(mrad)"); hp->SetFillColor(2);
126    TH1F *hl=new TH1F("hl","Lambda Resolution",30,-30,30);
127    hl->SetXTitle("(mrad)"); hl->SetFillColor(1); hl->SetFillStyle(3013); 
128    TH1F *hpt=new TH1F("hpt","Relative Pt Resolution",30,-10.,10.); 
129    hpt->SetXTitle("(%)"); hpt->SetFillColor(2); 
130
131    TH1F *hx=new TH1F("hx","Position Resolution (X,Y)",30,-3.,3.); //x res. 
132    hx->SetXTitle("(mm)"); hx->SetFillColor(6);
133    TH1F *hy=new TH1F("hy","Position Resolution (Y)",30,-3.,3.);   //y res
134    hy->SetXTitle("(mm)"); hy->SetFillColor(1); hy->SetFillStyle(3013);
135    TH1F *hz=new TH1F("hz","Position Resolution (Z)",30,-3.,3.);   //z res. 
136    hz->SetXTitle("(mm)"); hz->SetFillColor(6);
137
138    Double_t pmin=0.2, pmax=4.2; Int_t nchan=20;
139    TH1F *hgood=new TH1F("hgood","Good Cascades",nchan,pmin,pmax);    
140    TH1F *hfound=new TH1F("hfound","Found Cascades",nchan,pmin,pmax);
141    TH1F *hfake=new TH1F("hfake","Fake Cascades",nchan,pmin,pmax);
142    TH1F *hg=new TH1F("hg","Efficiency for Good Cascades",nchan,pmin,pmax);
143    hg->SetLineColor(4); hg->SetLineWidth(2);
144    TH1F *hf=new TH1F("hf","Probability of Fake Cascades",nchan,pmin,pmax);
145    hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
146
147    Double_t mmin=cascadeMass-cascadeWindow, mmax=cascadeMass+cascadeWindow;
148    TH1F *cs =new TH1F("cs","Cascade Effective Mass",40, mmin, mmax);
149    cs->SetXTitle("(GeV)");
150    cs->SetLineColor(4); cs->SetLineWidth(4);
151    TH1F *csf =new TH1F("csf","Fake Cascade Effective Mass",40, mmin, mmax);
152    csf->SetXTitle("(GeV)"); csf->SetFillColor(6);
153
154    Double_t pxg=0.,pyg=0.,ptg=0.;
155    Int_t nlab=-1, plab=-1, blab=-1;
156    Int_t i;
157    for (i=0; i<nentr; i++) {
158        AliCascadeVertex *cascade=(AliCascadeVertex*)carray.UncheckedAt(i);
159        nlab=TMath::Abs(cascade->GetNlabel()); 
160        plab=TMath::Abs(cascade->GetPlabel());
161        blab=TMath::Abs(cascade->GetBlabel());
162
163        /** Kinematical cuts **/
164        Double_t pxn,pyn,pzn; cascade->GetNPxPyPz(pxn,pyn,pzn); 
165        Double_t ptn=TMath::Sqrt(pxn*pxn + pyn*pyn);
166        if (ptn < ptncut) continue;
167        Double_t pxp,pyp,pzp; cascade->GetPPxPyPz(pxp,pyp,pzp); 
168        Double_t ptp=TMath::Sqrt(pxp*pxp + pyp*pyp);
169        if (ptp < ptpcut) continue;
170        Double_t pxb,pyb,pzb; cascade->GetBPxPyPz(pxb,pyb,pzb); 
171        Double_t ptb=TMath::Sqrt(pxb*pxb + pyb*pyb);
172        if (ptb < ptbcut) continue;
173        Double_t kine0;
174        Double_t kine=cascade->ChangeMassHypothesis(kine0,code);
175        if (TMath::Abs(kine0)>kine0cut) continue;
176        //if (TMath::Abs(kine)>kinecut) continue;
177
178        Double_t mass=cascade->GetEffMass();
179        cs->Fill(mass);
180        csf->Fill(mass);
181
182        if (TMath::Abs(mass-cascadeMass)>cascadeWidth) continue;
183
184        Int_t j;
185        for (j=0; j<ngood; j++) {
186           if (gc[j].code != cascade->GetPdgCode()) continue;
187           if (gc[j].nlab == nlab)
188           if (gc[j].plab == plab)
189           if (gc[j].blab == blab) break;
190        }
191
192        Double_t px,py,pz; cascade->GetPxPyPz(px,py,pz);
193        Double_t pt=TMath::Sqrt(px*px+py*py);
194
195        if (j==ngood) {
196           hfake->Fill(pt);
197           cerr<<"Fake cascade: ("<<nlab<<","<<plab<<","<<blab<<")\n";
198           continue;
199        }
200        csf->Fill(mass,-1);
201
202        pxg=gc[j].px; pyg=gc[j].py; ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
203        Double_t phig=TMath::ATan2(pyg,pxg), phi=TMath::ATan2(py,px);
204        Double_t lamg=TMath::ATan2(gc[j].pz,ptg), lam=TMath::ATan2(pz,pt);
205        hp->Fill((phi - phig)*1000.);
206        hl->Fill((lam - lamg)*1000.);
207        hpt->Fill((1/pt - 1/ptg)/(1/ptg)*100.);
208
209        Double_t x,y,z; cascade->GetXYZ(x,y,z);
210        hx->Fill((x-gc[j].x)*10);
211        hy->Fill((y-gc[j].y)*10);
212        hz->Fill((z-gc[j].z)*10);
213
214        hfound->Fill(ptg);
215        gc[j].nlab=-1;
216
217    }
218    for (i=0; i<ngood; i++) {
219       if (gc[i].code != code) continue;
220       pxg=gc[i].px; pyg=gc[i].py; ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
221       hgood->Fill(ptg);
222       nlab=gc[i].nlab; plab=gc[i].plab; blab=gc[i].blab;
223       if (nlab < 0) continue;
224      cerr<<"Cascade ("<<nlab<<','<<plab<<","<<blab<<") has not been found !\n";
225    }
226
227    carray.Delete();
228
229    Stat_t ng=hgood->GetEntries();
230    Stat_t nf=hfound->GetEntries();
231
232    cerr<<"Number of found cascades: "<<nentr<<" ("<<nf<<" in the peak)\n";
233    cerr<<"Number of good cascades: "<<ng<<endl;
234
235    if (ng!=0) 
236       cerr<<"Integral efficiency is about "<<nf/ng*100.<<" %\n";
237
238    gStyle->SetOptStat(111110);
239    gStyle->SetOptFit(1);
240
241    TCanvas *c1=new TCanvas("c1","",0,0,580,610);
242    c1->Divide(2,2);
243
244    c1->cd(1);
245    gPad->SetFillColor(42); gPad->SetFrameFillColor(10); 
246    //hp->Fit("gaus");
247    hp->Draw();
248    hl->Draw("same"); c1->cd();
249
250    c1->cd(2);
251    gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
252    //hpt->Fit("gaus"); c1->cd();
253    hpt->Draw(); c1->cd();
254
255    c1->cd(3);
256    gPad->SetFillColor(42); gPad->SetFrameFillColor(10); 
257    //hx->Fit("gaus"); 
258    hx->Draw(); 
259    hy->Draw("same"); c1->cd();
260
261    c1->cd(4);
262    gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
263    //hz->Fit("gaus");
264    hz->Draw();
265
266
267    TCanvas *c2=new TCanvas("c2","",600,0,580,610);
268    c2->Divide(1,2);
269
270    c2->cd(1);
271    gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
272    hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
273    hg->Divide(hfound,hgood,1,1.,"b");
274    hf->Divide(hfake,hgood,1,1.,"b");
275    hg->SetMaximum(1.4);
276    hg->SetYTitle("Cascade reconstruction efficiency");
277    hg->SetXTitle("Pt (GeV/c)");
278    hg->Draw();
279
280    TLine *line1 = new TLine(pmin,1.0,pmax,1.0); line1->SetLineStyle(4);
281    line1->Draw("same");
282    TLine *line2 = new TLine(pmin,0.9,pmax,0.9); line2->SetLineStyle(4);
283    line2->Draw("same");
284
285    hf->SetFillColor(1);
286    hf->SetFillStyle(3013);
287    hf->SetLineColor(2);
288    hf->SetLineWidth(2);
289    hf->Draw("histsame");
290    TText *text = new TText(0.461176,0.248448,"Fake cascades");
291    text->SetTextSize(0.05);
292    text->Draw();
293    text = new TText(0.453919,1.11408,"Good cascades");
294    text->SetTextSize(0.05);
295    text->Draw();
296
297
298    c2->cd(2);
299    gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
300    //cs->Fit("gaus","","",cascadeMass-cascadeWidth,cascadeMass+cascadeWidth);
301    cs->Draw();
302    csf->Draw("same");
303    Double_t max=cs->GetMaximum();
304    TLine *line3 = 
305       new TLine(cascadeMass-cascadeWidth,0.,cascadeMass-cascadeWidth,max);
306    line3->Draw("same");
307    TLine *line4 = 
308       new TLine(cascadeMass+cascadeWidth,0.,cascadeMass+cascadeWidth,max);
309    line4->Draw("same");
310
311    timer.Stop(); timer.Print();
312
313    return 0;
314 }
315
316
317
318 Int_t good_cascades(GoodCascade *gc, Int_t max) {
319    Int_t nc=0;
320    /*** Get information about the cuts ***/
321    Double_t r2min=0.9*0.9;
322    Double_t r2max=2.9*2.9;
323
324    /*** Get labels of the "good" tracks ***/
325    Double_t dd; Int_t id, label[15000], ngt=0;
326    ifstream in("good_tracks_its");
327    if (!in) {
328      cerr<<"Can't open the file good_tracks_its \n";
329      return nc;
330    }
331    while (in>>label[ngt]>>id>>dd>>dd>>dd>>dd>>dd>>dd) {
332      ngt++;
333      if (ngt>=15000) {
334        cerr<<"Too many good ITS tracks !\n";
335        return nc;
336      }
337    }   
338    if (!in.eof()) {
339       cerr<<"Read error (good_tracks_its) !\n";
340       return nc;
341    }
342
343    /*** Get an access to the kinematics ***/
344    if (gAlice) {delete gAlice; gAlice=0;}
345
346    TFile *file=TFile::Open("galice.root");
347    if (!file->IsOpen()) {cerr<<"Can't open galice.root !\n"; exit(4);}
348    if (!(gAlice=(AliRun*)file->Get("gAlice"))) {
349      cerr<<"gAlice has not been found on galice.root !\n";
350      exit(5);
351    }
352
353    Int_t np=gAlice->GetEvent(0);
354    while (np--) {
355       cerr<<np<<'\r';
356       TParticle *cp=gAlice->Particle(np);
357
358       /*** only these cascades are "good" ***/
359       Int_t code=cp->GetPdgCode();
360       if (code!=kXiMinus)    if (code!=kXiPlusBar)
361       if (code!=kOmegaMinus) if (code!=kOmegaPlusBar) continue; 
362
363       /*** daughter tracks must be "good" ***/
364       Int_t v0lab=cp->GetFirstDaughter(), blab=cp->GetLastDaughter();
365       if (v0lab==blab) continue;
366       if (v0lab<0) continue;
367       if (blab<0) continue;
368
369       TParticle *p0=gAlice->Particle(v0lab);
370       TParticle *bp=gAlice->Particle(blab);
371       if ((p0->GetPdgCode()!=kLambda0) && (p0->GetPdgCode()!=kLambda0Bar)) {
372          TParticle *p=p0; p0=bp; bp=p;
373          Int_t i=v0lab; v0lab=blab; blab=i;         
374          if ((p0->GetPdgCode()!=kLambda0) && (p0->GetPdgCode()!=kLambda0Bar))
375                                                                    continue;
376       }
377
378       /** is the bachelor "good" ? **/
379       Int_t i;
380       for (i=0; i<ngt; i++) if (label[i]==blab) break;
381       if (i==ngt) continue; 
382
383       /** is the V0 "good" ? **/
384       Int_t plab=p0->GetFirstDaughter(), nlab=p0->GetLastDaughter();
385       if (nlab==plab) continue;
386       if (nlab<0) continue;
387       if (plab<0) continue;
388
389       for (i=0; i<ngt; i++) if (label[i]==nlab) break;
390       if (i==ngt) continue;
391       for (i=0; i<ngt; i++) if (label[i]==plab) break;
392       if (i==ngt) continue;
393
394       /*** fiducial volume ***/
395       Double_t x=bp->Vx(), y=bp->Vy(), r2=x*x+y*y; //bachelor
396       if (r2<r2min) continue;
397       if (r2>r2max) continue;
398       TParticle *pp=gAlice->Particle(plab);
399       x=pp->Vx(); y=pp->Vy(); r2=x*x+y*y;                      //V0
400       if (r2<r2min) continue;
401       if (r2>r2max) continue;
402
403       if (gAlice->Particle(plab)->GetPDG()->Charge() < 0.) {
404          i=plab; plab=nlab; nlab=i;
405       }
406
407       gc[nc].code=code;
408       gc[nc].plab=plab;   gc[nc].nlab=nlab; gc[nc].blab=blab;
409       gc[nc].px=cp->Px(); gc[nc].py=cp->Py(); gc[nc].pz=cp->Pz();
410       gc[nc].x=bp->Vx(); gc[nc].y=bp->Vy(); gc[nc].z=bp->Vz();
411       nc++;
412
413    }
414
415
416    return nc;
417 }