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