]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliV0Comparison.C
The default thickness of the chips is set to 150 mkm (D.Elia). Removing some obsolete...
[u/mrichter/AliRoot.git] / STEER / 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 #if !defined(__CINT__) || defined(__MAKECINT__)
11   #include <TMath.h>
12   #include <TError.h>
13   #include <Riostream.h>
14   #include <TH1F.h>
15   #include <TTree.h>
16   #include <TParticle.h>
17   #include <TCanvas.h>
18   #include <TLine.h>
19   #include <TText.h>
20   #include <TBenchmark.h>
21   #include <TStyle.h>
22   #include <TFile.h>
23   #include <TROOT.h>
24
25   #include "AliStack.h"
26   #include "AliHeader.h"
27   #include "AliTrackReference.h"
28   #include "AliRunLoader.h"
29   #include "AliRun.h"
30   #include "AliESD.h"
31 #endif
32
33 Int_t GoodV0s(const Char_t *dir=".");
34
35 extern AliRun *gAlice;
36 extern TBenchmark *gBenchmark;
37 extern TROOT *gROOT;
38
39 static Int_t allgood=0;
40 static Int_t allfound=0;
41
42 Int_t AliV0Comparison(Int_t code=310, const Char_t *dir=".") { 
43    //Lambda=3122, LambdaBar=-3122
44    gBenchmark->Start("AliV0Comparison");
45
46    ::Info("AliV0Comparison.C","Doing comparison...");
47
48
49    const Double_t V0window=0.05;
50    Double_t ptncut=0.13, ptpcut=0.13, kinecut=0.03; 
51    Double_t V0mass=0.497672, V0width=0.020;
52    switch (code) {
53    case kK0Short: 
54         break;
55    case kLambda0:    
56         V0mass=1.115683; V0width=0.015; ptpcut=0.50; kinecut=0.002; 
57         break;
58    case kLambda0Bar: 
59         V0mass=1.115683; V0width=0.015; ptncut=0.50; kinecut=0.002;
60         break;
61    default: cerr<<"Invalid PDG code !\n"; return 1;
62    }
63
64    TH1F *hp=(TH1F*)gROOT->FindObject("hp");
65    if (!hp) hp=new TH1F("hp","PHI Resolution",30,-30.,30.);
66    hp->SetXTitle("(mrad)"); hp->SetFillColor(2);
67
68    TH1F *hl=(TH1F*)gROOT->FindObject("hl");
69    if (!hl) hl=new TH1F("hl","LAMBDA Resolution",30,-30,30);
70    hl->SetXTitle("(mrad)"); hl->SetFillColor(1); hl->SetFillStyle(3013);
71  
72    TH1F *hpt=(TH1F*)gROOT->FindObject("hpt");
73    if (!hpt) hpt=new TH1F("hpt","Relative Pt Resolution",30,-10.,10.); 
74    hpt->SetXTitle("(%)"); hpt->SetFillColor(2); 
75
76    TH1F *hx=(TH1F*)gROOT->FindObject("hx");
77    if (!hx) hx=new TH1F("hx","Position Resolution (X)",30,-3.,3.);
78    hx->SetXTitle("(mm)"); hx->SetFillColor(6);
79
80    TH1F *hy=(TH1F*)gROOT->FindObject("hy");
81    if (!hy) hy=new TH1F("hy","Position Resolution (Y)",30,-3.,3.); 
82    hy->SetXTitle("(mm)"); hy->SetFillColor(1); hy->SetFillStyle(3013);
83
84    TH1F *hz=(TH1F*)gROOT->FindObject("hz");
85    if (!hz) hz=new TH1F("hz","Position Resolution (Z)",30,-3.,3.);
86    hz->SetXTitle("(mm)"); hz->SetFillColor(6);
87
88
89    Double_t pmin=0.2, pmax=4.2; Int_t nchan=20;
90    TH1F *hgood=(TH1F*)gROOT->FindObject("hgood");    
91    if (!hgood) hgood=new TH1F("hgood","Good Vertices",nchan,pmin,pmax);
92     
93    TH1F *hfound=(TH1F*)gROOT->FindObject("hfound");
94    if (!hfound) hfound=new TH1F("hfound","Found Vertices",nchan,pmin,pmax);
95
96    TH1F *hfake=(TH1F*)gROOT->FindObject("hfake");
97    if (!hfake) hfake=new TH1F("hfake","Fake Vertices",nchan,pmin,pmax);
98
99    TH1F *hg=(TH1F*)gROOT->FindObject("hg");
100    if (!hg) hg=new TH1F("hg","Efficiency for Good Vertices",nchan,pmin,pmax);
101    hg->SetLineColor(4); hg->SetLineWidth(2);
102
103    TH1F *hf=(TH1F*)gROOT->FindObject("hf");
104    if (!hf) hf=new TH1F("hf","Probability of Fake Vertices",nchan,pmin,pmax);
105    hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
106
107
108    Double_t mmin=V0mass-V0window, mmax=V0mass+V0window;
109    TH1F *v0s=(TH1F*)gROOT->FindObject("v0s");
110    if (!v0s) v0s=new TH1F("v0s","V0s Effective Mass",40, mmin, mmax);
111    v0s->SetXTitle("(GeV)"); v0s->SetLineColor(4); v0s->SetLineWidth(4);
112
113    TH1F *v0sf=(TH1F*)gROOT->FindObject("v0sf");
114    if (!v0sf) v0sf=new TH1F("v0sf","Fake V0s Effective Mass",40, mmin, mmax);
115    v0sf->SetXTitle("(GeV)"); v0sf->SetFillColor(6);
116
117
118    Char_t fname[100];
119    sprintf(fname,"%s/GoodV0s.root",dir);
120
121    TFile *refFile=TFile::Open(fname,"old");
122    if (!refFile || !refFile->IsOpen()) {
123    ::Info("AliV0Comparison.C","Marking good V0s (will take a while)...");
124      if (GoodV0s(dir)) {
125         ::Error("AliV0Comparison.C","Can't generate the reference file !");
126         return 1;
127      }
128    }
129    refFile=TFile::Open(fname,"old");
130    if (!refFile || !refFile->IsOpen()) {
131      ::Error("AliV0Comparison.C","Can't open the reference file !");
132      return 1;
133    }   
134   
135    TTree *v0Tree=(TTree*)refFile->Get("v0Tree");
136    if (!v0Tree) {
137      ::Error("AliV0Comparison.C","Can't get the reference tree !");
138      return 2;
139    }
140    TBranch *pbranch=v0Tree->GetBranch("positive");
141    if (!pbranch) {
142      ::Error("AliV0Comparison.C","Can't get the positive daughter branch !");
143      return 3;
144    }
145    TClonesArray dummy("AliTrackReference",1000), *prefs=&dummy;
146    pbranch->SetAddress(&prefs);
147
148    TBranch *nbranch=v0Tree->GetBranch("negative");
149    if (!nbranch) {
150      ::Error("AliV0Comparison.C","Can't get the negative daughter branch !");
151      return 4;
152    }
153    TClonesArray dumm("AliTrackReference",1000), *nrefs=&dumm;
154    nbranch->SetAddress(&nrefs);
155
156    
157    sprintf(fname,"%s/AliESDs.root",dir);
158    TFile *ef=TFile::Open(fname);
159    if ((!ef)||(!ef->IsOpen())) {
160       sprintf(fname,"%s/AliESDv0.root",dir);
161       ef=TFile::Open(fname);
162       if ((!ef)||(!ef->IsOpen())) {
163          ::Error("AliV0Comparison.C","Can't open AliESDv0.root !");
164          return 5;
165       }
166    }
167    AliESD* event = new AliESD;
168    TTree* esdTree = (TTree*) ef->Get("esdTree");
169    if (!esdTree) {
170       ::Error("AliV0Comparison.C", "no ESD tree found");
171       return 6;
172    }
173    esdTree->SetBranchAddress("ESD", &event);
174
175
176    //******* Loop over events *********
177    Int_t e=0;
178    while (esdTree->GetEvent(e)) {
179       cout<<endl<<endl<<"********* Processing event number: "<<e<<"*******\n";
180  
181       Int_t nentr=event->GetNumberOfV0s();
182       allfound+=nentr;
183
184       if (v0Tree->GetEvent(e++)==0) {
185          cerr<<"No reconstructable V0s !\n";
186          continue;
187       }
188
189       Int_t ngood=prefs->GetEntriesFast(),ng=0; 
190
191       Double_t pxg=0.,pyg=0.,pzg=0.,ptg=0.;
192       Int_t nlab=-1, plab=-1;
193       Int_t i;
194       for (i=0; i<nentr; i++) {
195           AliESDv0 *vertex=event->GetV0(i);
196
197           Int_t nidx=TMath::Abs(vertex->GetNindex());
198           Int_t pidx=TMath::Abs(vertex->GetPindex());
199
200           AliESDtrack *ntrack=event->GetTrack(nidx);
201           AliESDtrack *ptrack=event->GetTrack(pidx);
202
203           nlab=TMath::Abs(ntrack->GetLabel());
204           plab=TMath::Abs(ptrack->GetLabel());
205
206           /** Kinematical cuts **/
207           Double_t pxn,pyn,pzn; vertex->GetNPxPyPz(pxn,pyn,pzn); 
208           Double_t ptn=TMath::Sqrt(pxn*pxn + pyn*pyn);
209           if (ptn < ptncut) continue;
210           Double_t pxp,pyp,pzp; vertex->GetPPxPyPz(pxp,pyp,pzp); 
211           Double_t ptp=TMath::Sqrt(pxp*pxp + pyp*pyp);
212           if (ptp < ptpcut) continue;
213           Double_t kine=vertex->ChangeMassHypothesis(code);
214           //if (TMath::Abs(kine)>kinecut) continue;
215
216           Double_t mass=vertex->GetEffMass();
217           v0s->Fill(mass);
218           v0sf->Fill(mass);
219
220           AliTrackReference *nref=0, *pref=0;
221           Int_t j;
222           for (j=0; j<ngood; j++) {
223               nref=(AliTrackReference*)nrefs->UncheckedAt(j);
224               pref=(AliTrackReference*)prefs->UncheckedAt(j);
225               if (nref->Label() == nlab)
226               if (pref->Label() == plab) break;
227           }
228
229           if (TMath::Abs(mass-V0mass)>V0width) continue;
230
231           Double_t px,py,pz; vertex->GetPxPyPz(px,py,pz);
232           Double_t pt=TMath::Sqrt(px*px+py*py);
233
234           if (j==ngood) {
235              hfake->Fill(pt);
236              cout<<"Fake vertex: ("<<nlab<<","<<plab<<")\n";
237              continue;
238           }
239           v0sf->Fill(mass,-1);
240
241           pxg=nref->Px()+pref->Px(); pyg=nref->Py()+pref->Py();
242           pzg=nref->Pz()+pref->Pz(); 
243           ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
244           Double_t phig=TMath::ATan2(pyg,pxg), phi=TMath::ATan2(py,px);
245           Double_t lamg=TMath::ATan2(pzg,ptg), lam=TMath::ATan2(pz,pt);
246           hp->Fill((phi - phig)*1000.);
247           hl->Fill((lam - lamg)*1000.);
248           hpt->Fill((1/pt - 1/ptg)/(1/ptg)*100.);
249
250           Double_t x,y,z; vertex->GetXYZ(x,y,z);
251           hx->Fill((x - nref->X())*10);
252           hy->Fill((y - nref->Y())*10);
253           hz->Fill((z - nref->Z())*10);
254
255           hfound->Fill(ptg);
256           nref->SetLabel(-1);
257
258       }
259       for (i=0; i<ngood; i++) {
260          AliTrackReference *nref=(AliTrackReference*)nrefs->UncheckedAt(i);
261          AliTrackReference *pref=(AliTrackReference*)prefs->UncheckedAt(i);
262          Int_t pdg=(Int_t)nref->GetLength();//this is the mother's PDG !
263          if (code!=pdg) continue;
264          ng++;
265          pxg=nref->Px()+pref->Px(); pyg=nref->Py()+pref->Py(); 
266          ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
267          hgood->Fill(ptg);
268          nlab=nref->Label(); plab=pref->Label();
269          if (nlab < 0) continue;
270          cout<<"Vertex ("<<nlab<<','<<plab<<") has not been found !\n";
271       }
272       allgood+=ng;
273
274       cout<<"Number of found V0s : "<<nentr<<endl;
275       cout<<"Number of \"good\" V0s : "<<ng<<endl;
276
277       prefs->Clear();
278       nrefs->Clear();
279
280    } //**** End of the loop over events
281
282    delete event;
283    ef->Close();
284    
285    delete v0Tree;
286    refFile->Close();
287
288    Stat_t ng=hgood->GetEntries(), nf=hfound->GetEntries();
289    if (ng!=0) cout<<"Integral efficiency is about "<<nf/ng*100.<<" %\n";
290    cout<<"Total number of found V0s: "<<allfound<<" ("<<nf<<" in the peak)\n";
291    cout<<"Total number of \"good\" V0s: "<<allgood<<endl;
292
293    gStyle->SetOptStat(111110);
294    gStyle->SetOptFit(1);
295
296    TCanvas *c1=(TCanvas*)gROOT->FindObject("c1");
297    if (!c1) {
298       c1=new TCanvas("c1","",0,0,580,610);
299       c1->Divide(2,2);
300    }
301
302    Int_t minc=33; 
303
304    c1->cd(1);
305    gPad->SetFillColor(42); gPad->SetFrameFillColor(10); 
306    if (hp->GetEntries()<minc) hp->Draw(); else hp->Fit("gaus");
307    hl->Draw("same"); c1->cd();
308
309    c1->cd(2);
310    gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
311    if (hpt->GetEntries()<minc) hpt->Draw(); else hpt->Fit("gaus");
312    c1->cd();
313
314    c1->cd(3);
315    gPad->SetFillColor(42); gPad->SetFrameFillColor(10); 
316    if (hx->GetEntries()<minc) hx->Draw(); else hx->Fit("gaus"); 
317    hy->Draw("same"); c1->cd();
318
319    c1->cd(4);
320    gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
321    if (hz->GetEntries()<minc) hz->Draw(); else hz->Fit("gaus");
322
323    c1->Update();   
324
325    TCanvas *c2=(TCanvas*)gROOT->FindObject("c2");
326    if (!c2) {
327       c2=new TCanvas("c2","",600,0,580,610);
328       c2->Divide(1,2);
329    }
330
331    c2->cd(1);
332    gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
333    hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
334    hg->Divide(hfound,hgood,1,1.,"b");
335    hf->Divide(hfake,hgood,1,1.,"b");
336    hg->SetMaximum(1.4);
337    hg->SetYTitle("V0 reconstruction efficiency");
338    hg->SetXTitle("Pt (GeV/c)");
339    hg->Draw();
340
341    TLine *line1 = new TLine(pmin,1.0,pmax,1.0); line1->SetLineStyle(4);
342    line1->Draw("same");
343    TLine *line2 = new TLine(pmin,0.9,pmax,0.9); line2->SetLineStyle(4);
344    line2->Draw("same");
345
346    hf->SetFillColor(1);
347    hf->SetFillStyle(3013);
348    hf->SetLineColor(2);
349    hf->SetLineWidth(2);
350    hf->Draw("histsame");
351    TText *text = new TText(0.461176,0.248448,"Fake vertices");
352    text->SetTextSize(0.05);
353    text->Draw();
354    text = new TText(0.453919,1.11408,"Good vertices");
355    text->SetTextSize(0.05);
356    text->Draw();
357
358
359    c2->cd(2);
360    gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
361    if (v0s->GetEntries()<minc) v0s->Draw();
362    else v0s->Fit("gaus","","",V0mass-V0width,V0mass+V0width);
363    v0sf->Draw("same");
364    Double_t max=v0s->GetMaximum();
365    TLine *line3 = new TLine(V0mass-V0width,0.,V0mass-V0width,max);
366    line3->Draw("same");
367    TLine *line4 = new TLine(V0mass+V0width,0.,V0mass+V0width,max);
368    line4->Draw("same");
369
370    c2->Update();
371
372    TFile fc("AliV0Comparison.root","RECREATE");
373    c1->Write();
374    c2->Write();
375    fc.Close();
376
377    gBenchmark->Stop("AliV0Comparison");
378    gBenchmark->Show("AliV0Comparison");
379
380    return 0;
381 }
382
383
384 Int_t GoodV0s(const Char_t *dir) {
385    if (gAlice) { 
386        delete gAlice->GetRunLoader();
387        delete gAlice;//if everything was OK here it is already NULL
388        gAlice = 0x0;
389    }
390
391    Char_t fname[100];
392    sprintf(fname,"%s/galice.root",dir);
393
394    AliRunLoader *rl = AliRunLoader::Open(fname,"COMPARISON");
395    if (!rl) {
396       ::Error("GoodTracksITS","Can't start session !");
397       return 1;
398    }
399
400    rl->LoadgAlice();
401    rl->LoadHeader();
402    rl->LoadKinematics();
403
404
405    Int_t nev=rl->GetNumberOfEvents();
406    ::Info("GoodV0s","Number of events : %d\n",nev);  
407
408    sprintf(fname,"%s/GoodTracksITS.root",dir);
409    TFile *itsFile=TFile::Open(fname);
410    if ((!itsFile)||(!itsFile->IsOpen())) {
411        ::Error("GoodV0s","Can't open the GoodTracksITS.root !");
412        delete rl;
413        return 5; 
414    }
415    TClonesArray dum("AliTrackReference",1000), *itsRefs=&dum;
416    TTree *itsTree=(TTree*)itsFile->Get("itsTree");
417    if (!itsTree) {
418        ::Error("GoodV0s","Can't get the ITS reference tree !");
419        delete rl;
420        return 6;
421    }
422    TBranch *itsBranch=itsTree->GetBranch("ITS");
423    if (!itsBranch) {
424       ::Error("GoodV0s","Can't get the ITS reference branch !");
425       delete rl;
426       return 7;
427    }
428    itsBranch->SetAddress(&itsRefs);
429
430
431    sprintf(fname,"%s/GoodV0s.root",dir);
432    TFile *v0File=TFile::Open(fname,"recreate");
433    TClonesArray dummy("AliTrackReference",1000), *nrefs=&dummy;
434    TClonesArray dumm("AliTrackReference",1000), *prefs=&dumm;
435    TTree v0Tree("v0Tree","Tree with info about the reconstructable V0s");
436    v0Tree.Branch("negative",&nrefs);
437    v0Tree.Branch("positive",&prefs);
438
439
440    /*** Some information about the cuts ***/
441    Double_t r2min=0.9*0.9;
442    Double_t r2max=2.9*2.9;
443
444
445
446    //********  Loop over generated events
447    for (Int_t e=0; e<nev; e++) {
448       rl->GetEvent(e);  v0File->cd();
449
450       Int_t np = rl->GetHeader()->GetNtrack();
451       cout<<"Event "<<e<<" Number of particles: "<<np<<endl;
452
453       itsTree->GetEvent(e);
454       Int_t nk=itsRefs->GetEntriesFast();
455
456       AliStack *stack=rl->Stack();
457
458       AliTrackReference *nref=0, *pref=0;
459
460       Int_t nv=0;
461       while (np--) {
462         //cerr<<np<<'\r';
463          TParticle *p0=stack->Particle(np);
464
465          /*** only these V0s are "good" ***/
466          Int_t code=p0->GetPdgCode();
467          if (code!=kK0Short)
468          if (code!=kLambda0)
469          if (code!=kLambda0Bar) continue; 
470
471          /*** daughter tracks should be "good" ***/
472          Int_t plab=p0->GetFirstDaughter(), nlab=p0->GetLastDaughter();
473          if (nlab==plab) continue;
474          if (nlab<0) continue;
475          if (plab<0) continue;
476          Int_t i;
477          if (stack->Particle(plab)->GetPDG()->Charge() < 0.) {
478             i=plab; plab=nlab; nlab=i;
479          }
480
481          for (i=0; i<nk; i++) {
482              nref=(AliTrackReference*)itsRefs->UncheckedAt(i);
483              if (nref->Label()==nlab) break;
484          }
485          if (i==nk) continue;
486          for (i=0; i<nk; i++) {
487              pref=(AliTrackReference*)itsRefs->UncheckedAt(i);
488              if (pref->Label()==plab) break;
489          }
490          if (i==nk) continue;
491
492          /*** fiducial volume ***/
493          TParticle *p=stack->Particle(nlab);
494          Double_t x=p->Vx(), y=p->Vy(), z=p->Vz(), r2=x*x+y*y;
495          if (r2<r2min) continue;
496          if (r2>r2max) continue;
497        
498          Int_t pdg=p0->GetPdgCode();
499          nref->SetLength(pdg);  //This will the V0's PDG !
500
501          new((*nrefs)[nv]) AliTrackReference(*nref);
502          new((*prefs)[nv]) AliTrackReference(*pref);
503
504          nv++;
505       }
506       itsRefs->Clear();
507
508       v0Tree.Fill();
509       nrefs->Clear(); prefs->Clear();
510
511    }  //**** end of the loop over generated events
512
513    v0Tree.Write();
514    v0File->Close();
515
516    delete itsTree;
517    itsFile->Close();
518
519    delete rl;
520    return 0;
521 }