]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliV0Comparison.C
new AliSplineFit class. Performs recursive fits to data point arrays
[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           if (vertex->GetOnFlyStatus()) continue;
198
199           Int_t nidx=TMath::Abs(vertex->GetNindex());
200           Int_t pidx=TMath::Abs(vertex->GetPindex());
201
202           AliESDtrack *ntrack=event->GetTrack(nidx);
203           AliESDtrack *ptrack=event->GetTrack(pidx);
204
205           nlab=TMath::Abs(ntrack->GetLabel());
206           plab=TMath::Abs(ptrack->GetLabel());
207
208           /** Kinematical cuts **/
209           Double_t pxn,pyn,pzn; vertex->GetNPxPyPz(pxn,pyn,pzn); 
210           Double_t ptn=TMath::Sqrt(pxn*pxn + pyn*pyn);
211           if (ptn < ptncut) continue;
212           Double_t pxp,pyp,pzp; vertex->GetPPxPyPz(pxp,pyp,pzp); 
213           Double_t ptp=TMath::Sqrt(pxp*pxp + pyp*pyp);
214           if (ptp < ptpcut) continue;
215           Double_t kine=vertex->ChangeMassHypothesis(code);
216           //if (TMath::Abs(kine)>kinecut) continue;
217
218           Double_t mass=vertex->GetEffMass();
219           v0s->Fill(mass);
220           v0sf->Fill(mass);
221
222           AliTrackReference *nref=0, *pref=0;
223           Int_t j;
224           for (j=0; j<ngood; j++) {
225               nref=(AliTrackReference*)nrefs->UncheckedAt(j);
226               pref=(AliTrackReference*)prefs->UncheckedAt(j);
227               if (nref->Label() == nlab)
228               if (pref->Label() == plab) break;
229           }
230
231           if (TMath::Abs(mass-V0mass)>V0width) continue;
232
233           Double_t px,py,pz; vertex->GetPxPyPz(px,py,pz);
234           Double_t pt=TMath::Sqrt(px*px+py*py);
235
236           if (j==ngood) {
237              hfake->Fill(pt);
238              cout<<"Fake vertex: ("<<nlab<<","<<plab<<")\n";
239              continue;
240           }
241           v0sf->Fill(mass,-1);
242
243           pxg=nref->Px()+pref->Px(); pyg=nref->Py()+pref->Py();
244           pzg=nref->Pz()+pref->Pz(); 
245           ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
246           Double_t phig=TMath::ATan2(pyg,pxg), phi=TMath::ATan2(py,px);
247           Double_t lamg=TMath::ATan2(pzg,ptg), lam=TMath::ATan2(pz,pt);
248           hp->Fill((phi - phig)*1000.);
249           hl->Fill((lam - lamg)*1000.);
250           hpt->Fill((1/pt - 1/ptg)/(1/ptg)*100.);
251
252           Double_t x,y,z; vertex->GetXYZ(x,y,z);
253           hx->Fill((x - nref->X())*10);
254           hy->Fill((y - nref->Y())*10);
255           hz->Fill((z - nref->Z())*10);
256
257           hfound->Fill(ptg);
258           nref->SetLabel(-1);
259
260       }
261       for (i=0; i<ngood; i++) {
262          AliTrackReference *nref=(AliTrackReference*)nrefs->UncheckedAt(i);
263          AliTrackReference *pref=(AliTrackReference*)prefs->UncheckedAt(i);
264          Int_t pdg=(Int_t)nref->GetLength();//this is the mother's PDG !
265          if (code!=pdg) continue;
266          ng++;
267          pxg=nref->Px()+pref->Px(); pyg=nref->Py()+pref->Py(); 
268          ptg=TMath::Sqrt(pxg*pxg+pyg*pyg);
269          hgood->Fill(ptg);
270          nlab=nref->Label(); plab=pref->Label();
271          if (nlab < 0) continue;
272          cout<<"Vertex ("<<nlab<<','<<plab<<") has not been found !\n";
273       }
274       allgood+=ng;
275
276       cout<<"Number of found V0s : "<<nentr<<endl;
277       cout<<"Number of \"good\" V0s : "<<ng<<endl;
278
279       prefs->Clear();
280       nrefs->Clear();
281
282    } //**** End of the loop over events
283
284    delete event;
285    ef->Close();
286    
287    delete v0Tree;
288    refFile->Close();
289
290    Stat_t ng=hgood->GetEntries(), nf=hfound->GetEntries();
291    if (ng!=0) cout<<"Integral efficiency is about "<<nf/ng*100.<<" %\n";
292    cout<<"Total number of found V0s: "<<allfound<<" ("<<nf<<" in the peak)\n";
293    cout<<"Total number of \"good\" V0s: "<<allgood<<endl;
294
295    gStyle->SetOptStat(111110);
296    gStyle->SetOptFit(1);
297
298    TCanvas *c1=(TCanvas*)gROOT->FindObject("c1");
299    if (!c1) {
300       c1=new TCanvas("c1","",0,0,580,610);
301       c1->Divide(2,2);
302    }
303
304    Int_t minc=33; 
305
306    c1->cd(1);
307    gPad->SetFillColor(42); gPad->SetFrameFillColor(10); 
308    if (hp->GetEntries()<minc) hp->Draw(); else hp->Fit("gaus");
309    hl->Draw("same"); c1->cd();
310
311    c1->cd(2);
312    gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
313    if (hpt->GetEntries()<minc) hpt->Draw(); else hpt->Fit("gaus");
314    c1->cd();
315
316    c1->cd(3);
317    gPad->SetFillColor(42); gPad->SetFrameFillColor(10); 
318    if (hx->GetEntries()<minc) hx->Draw(); else hx->Fit("gaus"); 
319    hy->Draw("same"); c1->cd();
320
321    c1->cd(4);
322    gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
323    if (hz->GetEntries()<minc) hz->Draw(); else hz->Fit("gaus");
324
325    c1->Update();   
326
327    TCanvas *c2=(TCanvas*)gROOT->FindObject("c2");
328    if (!c2) {
329       c2=new TCanvas("c2","",600,0,580,610);
330       c2->Divide(1,2);
331    }
332
333    c2->cd(1);
334    gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
335    hfound->Sumw2(); hgood->Sumw2(); hfake->Sumw2();
336    hg->Divide(hfound,hgood,1,1.,"b");
337    hf->Divide(hfake,hgood,1,1.,"b");
338    hg->SetMaximum(1.4);
339    hg->SetYTitle("V0 reconstruction efficiency");
340    hg->SetXTitle("Pt (GeV/c)");
341    hg->Draw();
342
343    TLine *line1 = new TLine(pmin,1.0,pmax,1.0); line1->SetLineStyle(4);
344    line1->Draw("same");
345    TLine *line2 = new TLine(pmin,0.9,pmax,0.9); line2->SetLineStyle(4);
346    line2->Draw("same");
347
348    hf->SetFillColor(1);
349    hf->SetFillStyle(3013);
350    hf->SetLineColor(2);
351    hf->SetLineWidth(2);
352    hf->Draw("histsame");
353    TText *text = new TText(0.461176,0.248448,"Fake vertices");
354    text->SetTextSize(0.05);
355    text->Draw();
356    text = new TText(0.453919,1.11408,"Good vertices");
357    text->SetTextSize(0.05);
358    text->Draw();
359
360
361    c2->cd(2);
362    gPad->SetFillColor(42); gPad->SetFrameFillColor(10);
363    if (v0s->GetEntries()<minc) v0s->Draw();
364    else v0s->Fit("gaus","","",V0mass-V0width,V0mass+V0width);
365    v0sf->Draw("same");
366    Double_t max=v0s->GetMaximum();
367    TLine *line3 = new TLine(V0mass-V0width,0.,V0mass-V0width,max);
368    line3->Draw("same");
369    TLine *line4 = new TLine(V0mass+V0width,0.,V0mass+V0width,max);
370    line4->Draw("same");
371
372    c2->Update();
373
374    TFile fc("AliV0Comparison.root","RECREATE");
375    c1->Write();
376    c2->Write();
377    fc.Close();
378
379    gBenchmark->Stop("AliV0Comparison");
380    gBenchmark->Show("AliV0Comparison");
381
382    return 0;
383 }
384
385
386 Int_t GoodV0s(const Char_t *dir) {
387    if (gAlice) { 
388        delete gAlice->GetRunLoader();
389        delete gAlice;//if everything was OK here it is already NULL
390        gAlice = 0x0;
391    }
392
393    Char_t fname[100];
394    sprintf(fname,"%s/galice.root",dir);
395
396    AliRunLoader *rl = AliRunLoader::Open(fname,"COMPARISON");
397    if (!rl) {
398       ::Error("GoodTracksITS","Can't start session !");
399       return 1;
400    }
401
402    rl->LoadgAlice();
403    rl->LoadHeader();
404    rl->LoadKinematics();
405
406
407    Int_t nev=rl->GetNumberOfEvents();
408    ::Info("GoodV0s","Number of events : %d\n",nev);  
409
410    sprintf(fname,"%s/GoodTracksITS.root",dir);
411    TFile *itsFile=TFile::Open(fname);
412    if ((!itsFile)||(!itsFile->IsOpen())) {
413        ::Error("GoodV0s","Can't open the GoodTracksITS.root !");
414        delete rl;
415        return 5; 
416    }
417    TClonesArray dum("AliTrackReference",1000), *itsRefs=&dum;
418    TTree *itsTree=(TTree*)itsFile->Get("itsTree");
419    if (!itsTree) {
420        ::Error("GoodV0s","Can't get the ITS reference tree !");
421        delete rl;
422        return 6;
423    }
424    TBranch *itsBranch=itsTree->GetBranch("ITS");
425    if (!itsBranch) {
426       ::Error("GoodV0s","Can't get the ITS reference branch !");
427       delete rl;
428       return 7;
429    }
430    itsBranch->SetAddress(&itsRefs);
431
432
433    sprintf(fname,"%s/GoodV0s.root",dir);
434    TFile *v0File=TFile::Open(fname,"recreate");
435    TClonesArray dummy("AliTrackReference",1000), *nrefs=&dummy;
436    TClonesArray dumm("AliTrackReference",1000), *prefs=&dumm;
437    TTree v0Tree("v0Tree","Tree with info about the reconstructable V0s");
438    v0Tree.Branch("negative",&nrefs);
439    v0Tree.Branch("positive",&prefs);
440
441
442    /*** Some information about the cuts ***/
443    Double_t r2min=0.9*0.9;
444    Double_t r2max=2.9*2.9;
445
446
447
448    //********  Loop over generated events
449    for (Int_t e=0; e<nev; e++) {
450       rl->GetEvent(e);  v0File->cd();
451
452       Int_t np = rl->GetHeader()->GetNtrack();
453       cout<<"Event "<<e<<" Number of particles: "<<np<<endl;
454
455       itsTree->GetEvent(e);
456       Int_t nk=itsRefs->GetEntriesFast();
457
458       AliStack *stack=rl->Stack();
459
460       AliTrackReference *nref=0, *pref=0;
461
462       Int_t nv=0;
463       while (np--) {
464         //cerr<<np<<'\r';
465          TParticle *p0=stack->Particle(np);
466
467          /*** only these V0s are "good" ***/
468          Int_t code=p0->GetPdgCode();
469          if (code!=kK0Short)
470          if (code!=kLambda0)
471          if (code!=kLambda0Bar) continue; 
472
473          /*** daughter tracks should be "good" ***/
474          Int_t plab=p0->GetFirstDaughter(), nlab=p0->GetLastDaughter();
475          if (nlab==plab) continue;
476          if (nlab<0) continue;
477          if (plab<0) continue;
478          Int_t i;
479          TParticle * part = stack->Particle(plab);
480          if (part) {
481            TParticlePDG * partPDG = part->GetPDG();
482              if (partPDG && partPDG->Charge() < 0.) {
483                i=plab; plab=nlab; nlab=i;
484              }
485          }
486
487          for (i=0; i<nk; i++) {
488              nref=(AliTrackReference*)itsRefs->UncheckedAt(i);
489              if (nref->Label()==nlab) break;
490          }
491          if (i==nk) continue;
492          for (i=0; i<nk; i++) {
493              pref=(AliTrackReference*)itsRefs->UncheckedAt(i);
494              if (pref->Label()==plab) break;
495          }
496          if (i==nk) continue;
497
498          /*** fiducial volume ***/
499          TParticle *p=stack->Particle(nlab);
500          Double_t x=p->Vx(), y=p->Vy(), z=p->Vz(), r2=x*x+y*y;
501          if (r2<r2min) continue;
502          if (r2>r2max) continue;
503        
504          Int_t pdg=p0->GetPdgCode();
505          nref->SetLength(pdg);  //This will the V0's PDG !
506
507          new((*nrefs)[nv]) AliTrackReference(*nref);
508          new((*prefs)[nv]) AliTrackReference(*pref);
509
510          nv++;
511       }
512       itsRefs->Clear();
513
514       v0Tree.Fill();
515       nrefs->Clear(); prefs->Clear();
516
517    }  //**** end of the loop over generated events
518
519    v0Tree.Write();
520    v0File->Close();
521
522    delete itsTree;
523    itsFile->Close();
524
525    delete rl;
526    return 0;
527 }