]>
Commit | Line | Data |
---|---|---|
e23730c7 | 1 | //******************************************************************** |
2 | // Example (very naive for the moment) of the data analysis | |
3 | // using the ESD classes. | |
4 | // It demonstrates the idea of the "combined PID" | |
5 | // applied to the Lambda0 reconstruction. | |
6 | // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch | |
7 | //******************************************************************** | |
8 | ||
9 | #if !defined( __CINT__) || defined(__MAKECINT__) | |
10 | #include <Riostream.h> | |
8b462fd8 | 11 | #include <TTree.h> |
9d0d65be | 12 | #include <TFile.h> |
13 | #include <TH1F.h> | |
14 | #include <TCanvas.h> | |
e23730c7 | 15 | |
16 | #include "AliESD.h" | |
17 | ||
18 | #endif | |
19 | ||
9d0d65be | 20 | Int_t AliESDv0Analysis(const Char_t *dir=".") { |
21 | TH1F *hm=(TH1F*)gROOT->FindObject("hm"); | |
22 | if (!hm) { | |
23 | hm=new TH1F("hm","Lambda+LambdaBar Effective Mass",60,1.065,1.165); | |
24 | hm->SetXTitle("Mass (GeV/c**2)"); | |
25 | } | |
26 | Char_t fname[100]; | |
27 | sprintf(fname,"%s/AliESDs.root",dir); | |
28 | TFile *ef=TFile::Open(fname); | |
29 | if (!ef||!ef->IsOpen()) {cerr<<"Can't AliESDs.root !\n"; return 1;} | |
30 | cerr<<"\n****** "<<fname<<" ******\n"; | |
e23730c7 | 31 | |
8b462fd8 | 32 | AliESD* event = new AliESD; |
9d0d65be | 33 | |
34 | ||
8b462fd8 | 35 | TTree* tree = (TTree*) ef->Get("esdTree"); |
36 | if (!tree) {cerr<<"no ESD tree found\n"; return 1;}; | |
37 | tree->SetBranchAddress("ESD", &event); | |
e23730c7 | 38 | |
e23730c7 | 39 | Int_t rc=0,n=0; |
e23730c7 | 40 | |
41 | //****** Tentative particle type "concentrations" | |
9d0d65be | 42 | Double_t c[5]={0.0, 0.0, 1, 0, 1}; |
4ca848eb | 43 | AliPID::SetPriors(c); |
e23730c7 | 44 | |
45 | //******* The loop over events | |
4ca848eb | 46 | while (tree->GetEvent(n)) { |
e23730c7 | 47 | |
48 | cerr<<"Processing event number : "<<n++<<endl; | |
49 | ||
e23730c7 | 50 | Int_t nv0=event->GetNumberOfV0s(); |
51 | cerr<<"Number of ESD v0s : "<<nv0<<endl; | |
52 | ||
53 | while (nv0--) { | |
54 | AliESDv0 *v0=event->GetV0(nv0); | |
9d0d65be | 55 | |
56 | Int_t protonIdx=v0->GetPindex(); | |
57 | Int_t pionIdx =v0->GetNindex(); | |
58 | ||
59 | v0->ChangeMassHypothesis(3122); | |
60 | Double_t mass=v0->GetEffMass(); | |
4ca848eb | 61 | if (mass>1.17) { //check also the LambdaBar hypothesis |
9d0d65be | 62 | v0->ChangeMassHypothesis(-3122); |
63 | mass=v0->GetEffMass(); | |
64 | if (mass>1.17) continue; | |
65 | Int_t tmp=protonIdx; protonIdx=pionIdx; pionIdx=tmp; | |
66 | } | |
67 | ||
68 | AliESDtrack *protonTrk=event->GetTrack(protonIdx); | |
69 | AliESDtrack *pionTrk =event->GetTrack(pionIdx); | |
70 | ||
71 | if (protonTrk->GetP()<0.5) continue; | |
72 | ||
73 | // Check if the "proton track" is a proton | |
74 | if ((protonTrk->GetStatus()&AliESDtrack::kESDpid)!=0) { | |
75 | Double_t r[10]; protonTrk->GetESDpid(r); | |
4ca848eb | 76 | AliPID pid(r); |
77 | Double_t pp=pid.GetProbability(AliPID::kProton); | |
78 | if (pp < pid.GetProbability(AliPID::kElectron)) continue; | |
79 | if (pp < pid.GetProbability(AliPID::kMuon)) continue; | |
80 | if (pp < pid.GetProbability(AliPID::kPion)) continue; | |
81 | if (pp < pid.GetProbability(AliPID::kKaon)) continue; | |
e23730c7 | 82 | } |
9d0d65be | 83 | |
84 | //Check if the "pion track" is a pion | |
85 | if ((pionTrk->GetStatus()&AliESDtrack::kESDpid)!=0) { | |
86 | Double_t r[10]; pionTrk->GetESDpid(r); | |
4ca848eb | 87 | AliPID pid(r); |
88 | Double_t ppi=pid.GetProbability(AliPID::kPion); | |
89 | if (ppi < pid.GetProbability(AliPID::kElectron)) continue; | |
90 | if (ppi < pid.GetProbability(AliPID::kMuon)) continue; | |
91 | if (ppi < pid.GetProbability(AliPID::kKaon)) continue; | |
92 | if (ppi < pid.GetProbability(AliPID::kProton)) continue; | |
93 | } | |
9d0d65be | 94 | |
e23730c7 | 95 | hm->Fill(mass); |
96 | } | |
9d0d65be | 97 | |
e23730c7 | 98 | } |
99 | ||
8b462fd8 | 100 | delete event; |
101 | ef->Close(); | |
102 | ||
9d0d65be | 103 | TCanvas *c1=(TCanvas*)gROOT->FindObject("c1"); |
104 | if (!c1) { | |
105 | c1=new TCanvas(); | |
106 | } | |
e23730c7 | 107 | |
9d0d65be | 108 | if (hm->GetEntries()>100) hm->Fit("gaus","","",1.11,1.12); |
109 | else hm->Draw(); | |
e23730c7 | 110 | |
9d0d65be | 111 | c1->Update(); |
e23730c7 | 112 | |
113 | return rc; | |
114 | } |