Applying the |eta|<0.8 cut for the daughter tracks
[u/mrichter/AliRoot.git] / PWGLF / STRANGENESS / LambdaK0PbPb / AliAnalysisTaskCTauPbPb.cxx
1 #include <TCanvas.h>
2 #include <TTree.h>
3 #include <TFile.h>
4 #include <TH1F.h>
5 #include <TH2F.h>
6 #include <TH3F.h>
7 #include <TPDGCode.h>
8 #include <TDatabasePDG.h>
9 #include <TParticlePDG.h>
10 #include <TParticle.h>
11 #include <TROOT.h>
12
13 #include "AliESDEvent.h"
14 #include "AliESDv0.h"
15 #include "AliESDcascade.h"
16
17 #include "AliCentrality.h"
18
19 #include "AliMCEvent.h"
20 #include "AliStack.h"
21
22 #include "AliPID.h"
23 #include "AliPIDResponse.h"
24
25 #include "AliInputEventHandler.h"
26 #include "AliAnalysisManager.h"
27
28 #include "AliAnalysisTaskCTauPbPb.h"
29
30 extern TROOT *gROOT;
31
32 ClassImp(AliAnalysisTaskCTauPbPb)
33
34 static Int_t    nbins=102;  // number of bins
35 static Double_t lMin=0.0, lMax=100.;
36 static Double_t pMin=0.0, pMax=10.;
37 static Double_t yMax=0.5;
38
39
40 //
41 //  This is a little task for checking the c*tau of the strange particles 
42 //
43
44 AliAnalysisTaskCTauPbPb::AliAnalysisTaskCTauPbPb(const char *name) :
45 AliAnalysisTaskSE(name),
46 fIsMC(kFALSE),
47 fCMin(0.),
48 fCMax(90.),
49 fOutput(0),
50 fMult(0),
51 fdEdx(0),
52 fdEdxPid(0),
53
54 fK0sM(0),
55 fK0sSi(0),
56 fK0sMC(0),
57 fK0sAs(0),
58
59 fLambdaM(0),
60 fLambdaSi(0),
61 fLambdaMC(0),
62 fLambdaAs(0),
63
64 fCPA(0),
65 fDCA(0),
66
67 fLambdaEff(0),
68 fLambdaPt(0),
69
70 fLambdaFromXi(0),
71 fXiM(0),
72 fXiSiP(0)
73 {
74   // Constructor. Initialization of pointers
75   DefineOutput(1, TList::Class());
76 }
77
78 void AliAnalysisTaskCTauPbPb::UserCreateOutputObjects()
79 {
80   fOutput = new TList(); 
81   fOutput->SetOwner();
82
83
84   fMult=new TH1F("fMult","Multiplicity",1100,0.,3300);
85   fMult->GetXaxis()->SetTitle("N tracks"); 
86   fOutput->Add(fMult);
87
88   fdEdx=new TH2F("fdEdx","dE/dx",50,0.2,3,50,0.,6.);
89   fOutput->Add(fdEdx);
90
91   fdEdxPid=new TH2F("fdEdxPid","dE/dx with PID",50,0.2,3,50,0.,6.);
92   fOutput->Add(fdEdxPid);
93
94   fK0sM = 
95   new TH2F("fK0sM", "Mass for K^{0}_{s}", nbins/2, 0.448, 0.548, 10,pMin,pMax);
96   fK0sM->GetXaxis()->SetTitle("Mass [GeV/c]"); 
97   fOutput->Add(fK0sM);
98
99   fK0sSi = 
100   new TH2F("fK0sSi","L_{T} vs p_{T} for K^{0}_{s}, side-band subtracted",
101   nbins,pMin,pMax,nbins,lMin,lMax);
102   fK0sSi->GetXaxis()->SetTitle("p_{T} [GeV/c]"); 
103   fK0sSi->GetYaxis()->SetTitle("L_{T} [cm]"); 
104   fOutput->Add(fK0sSi);
105
106   fK0sMC = 
107   new TH2F("fK0sMC","L_{T} vs p_{T} for K^{0}_{s}, from MC stack", 
108   nbins,pMin,pMax,nbins,lMin,lMax);
109   fK0sMC->GetXaxis()->SetTitle("p_{T} [GeV/c]"); 
110   fK0sMC->GetYaxis()->SetTitle("L_{T} [cm]"); 
111   fOutput->Add(fK0sMC);
112
113   fK0sAs = 
114   new TH2F("fK0sAs", "L_{T} vs p_{T} for K^{0}_{s}, associated", 
115   nbins,pMin,pMax,nbins,lMin,lMax);
116   fK0sAs->GetXaxis()->SetTitle("p_{T} [GeV/c]"); 
117   fK0sAs->GetYaxis()->SetTitle("L_{T} [cm]"); 
118   fOutput->Add(fK0sAs);
119
120   //----------------------
121
122   fLambdaM = 
123   new TH2F("fLambdaM","Mass for \\Lambda", nbins, 1.065, 1.165,nbins,pMin,pMax);
124   //new TH2F("fLambdaM","Mass for \\Lambda", nbins, 1.065, 1.165,10,0.1,1.1);
125   fLambdaM->GetXaxis()->SetTitle("Mass [GeV/c]"); 
126   fOutput->Add(fLambdaM);
127
128   fLambdaSi = 
129   new TH2F("fLambdaSi","L_{T} vs p_{T} for \\Lambda, side-band subtructed",
130   nbins,pMin,pMax,nbins,lMin,lMax);
131   fLambdaSi->GetXaxis()->SetTitle("p_{T} [GeV/c]"); 
132   fLambdaSi->GetYaxis()->SetTitle("L_{T} [cm]"); 
133   fOutput->Add(fLambdaSi);
134
135   fLambdaMC = 
136   new TH2F("fLambdaMC","c\\tau for \\Lambda, from MC stack", 
137   nbins,pMin,pMax,nbins,lMin,lMax);
138   fLambdaMC->GetXaxis()->SetTitle("p_{T} [GeV/c]"); 
139   fLambdaMC->GetYaxis()->SetTitle("L_{T} [cm]"); 
140   fOutput->Add(fLambdaMC);
141
142   fLambdaAs = 
143   new TH2F("fLambdaAs","c\\tau for \\Lambda, associated",
144   nbins,pMin,pMax,nbins,lMin,lMax);
145   fLambdaAs->GetXaxis()->SetTitle("p_{T} [GeV/c]"); 
146   fLambdaAs->GetYaxis()->SetTitle("L_{T} [cm]"); 
147   fOutput->Add(fLambdaAs);
148
149   fCPA=new TH1F("fCPA","Cosine of the pointing angle",30,0.9978,1.);
150   fOutput->Add(fCPA);
151   fDCA=new TH1F("fDCA","DCA between the daughters",30,0.,1.1);
152   fOutput->Add(fDCA);
153
154   fLambdaEff=fLambdaAs->ProjectionX();
155   fLambdaEff->SetName("fLambdaEff");
156   fLambdaEff->SetTitle("Efficiency for #Lambda");
157   fOutput->Add(fLambdaEff);
158
159   fLambdaPt=fLambdaAs->ProjectionX();
160   fLambdaPt->SetName("fLambdaPt");
161   fLambdaPt->SetTitle("Raw #Lambda pT spectrum");
162   fOutput->Add(fLambdaPt);
163
164   //----------------------
165
166   fLambdaFromXi=new TH3F("fLambdaFromXi","L_{T} vs p_{T} vs p_{T} of \\Xi for \\Lambda from Xi",
167   nbins,pMin,pMax,nbins,lMin,lMax,33,pMin,pMax+2);
168   fOutput->Add(fLambdaFromXi);
169
170   fXiM  = 
171   new TH2F("fXiM", "\\Xi mass distribution", 50, 1.271, 1.371,12,pMin,pMax+2);
172   fOutput->Add(fXiM);
173
174   fXiSiP  = new TH1F("fXiSiP", "Pt for \\Xi, side-band subracted",
175   33,pMin,pMax+2);
176   fOutput->Add(fXiSiP);
177
178
179   PostData(1, fOutput);
180 }
181
182 static Bool_t AcceptTrack(const AliESDtrack *t) {
183   if (!t->IsOn(AliESDtrack::kTPCrefit)) return kFALSE;
184   if (t->GetKinkIndex(0)>0) return kFALSE;
185
186   Float_t nCrossedRowsTPC = t->GetTPCClusterInfo(2,1); 
187   if (nCrossedRowsTPC < 70) return kFALSE;
188   Int_t findable=t->GetTPCNclsF();
189   if (findable <= 0) return kFALSE;
190   if (nCrossedRowsTPC/findable < 0.8) return kFALSE;
191
192   if (TMath::Abs(t->Eta()) > 0.8) return kFALSE;
193
194   return kTRUE;   
195 }
196
197 static Bool_t AcceptV0(const AliESDv0 *v0, const AliESDEvent *esd) {
198
199   if (v0->GetOnFlyStatus()) return kFALSE;
200
201   if (v0->Pt() < pMin) return kFALSE;
202
203   Int_t nidx=TMath::Abs(v0->GetNindex());
204   AliESDtrack *ntrack=esd->GetTrack(nidx);
205   if (!AcceptTrack(ntrack)) return kFALSE;
206
207   Int_t pidx=TMath::Abs(v0->GetPindex());
208   AliESDtrack *ptrack=esd->GetTrack(pidx);
209   if (!AcceptTrack(ptrack)) return kFALSE;
210
211   Float_t xy,z0;
212   ntrack->GetImpactParameters(xy,z0);
213   if (TMath::Abs(xy)<0.1) return kFALSE;
214   ptrack->GetImpactParameters(xy,z0);
215   if (TMath::Abs(xy)<0.1) return kFALSE;
216
217   Double_t dca=v0->GetDcaV0Daughters();
218   if (dca>1.0) return kFALSE;
219   //if (dca>0.7) return kFALSE;
220   //if (dca>0.4) return kFALSE;
221
222   Double_t cpa=v0->GetV0CosineOfPointingAngle();
223   if (cpa<0.998) return kFALSE;
224   //if (cpa<0.99875) return kFALSE;
225   //if (cpa<0.9995) return kFALSE;
226
227   Double_t xx,yy,zz; v0->GetXYZ(xx,yy,zz);
228   Double_t r2=xx*xx + yy*yy;
229   if (r2<0.9*0.9) return kFALSE;
230   if (r2>100*100) return kFALSE;
231
232   return kTRUE;
233 }
234
235 static Bool_t AcceptCascade(const AliESDcascade *cs, const AliESDEvent *esd) {
236
237   if (cs->Pt() < pMin) return kFALSE;
238
239   Int_t bidx=TMath::Abs(cs->GetBindex());
240   AliESDtrack *btrack=esd->GetTrack(bidx);
241   if (!AcceptTrack(btrack)) return kFALSE;
242
243   Float_t xy,z0; 
244   btrack->GetImpactParameters(xy,z0);
245   if (TMath::Abs(xy)<0.03) return kFALSE;
246
247   const AliESDVertex *vtx=esd->GetPrimaryVertexSPD();
248   if (!vtx->GetStatus()) {
249      vtx=esd->GetPrimaryVertexTracks();
250      if (!vtx->GetStatus()) return kFALSE;
251   }
252   Double_t xv=vtx->GetXv(), yv=vtx->GetYv(), zv=vtx->GetZv();
253   if (cs->GetCascadeCosineOfPointingAngle(xv,yv,zv) < 0.999) return kFALSE;
254
255   if (cs->GetDcaXiDaughters() > 0.3) return kFALSE;
256
257   return kTRUE;
258 }
259
260 void AliAnalysisTaskCTauPbPb::UserExec(Option_t *)
261 {
262
263   AliESDEvent *esd=(AliESDEvent *)InputEvent();
264
265   if (!esd) {
266     Printf("ERROR: esd not available");
267     return;
268   }
269
270   fMult->Fill(-100); //event counter  
271
272   // Physics selection
273   AliAnalysisManager *mgr= AliAnalysisManager::GetAnalysisManager();
274   AliInputEventHandler *hdr=(AliInputEventHandler*)mgr->GetInputEventHandler();
275   UInt_t maskIsSelected = hdr->IsEventSelected();
276   Bool_t isSelected = (maskIsSelected & AliVEvent::kMB);
277   if (!isSelected) return;
278
279   // Centrality selection
280   AliCentrality *cent=esd->GetCentrality();
281   if (!cent->IsEventInCentralityClass(fCMin,fCMax,"V0M")) return;
282
283   const AliESDVertex *vtx=esd->GetPrimaryVertexSPD();
284   if (!vtx->GetStatus()) {
285      vtx=esd->GetPrimaryVertexTracks();
286      if (!vtx->GetStatus()) return;
287   }
288   Double_t xv=vtx->GetXv(), yv=vtx->GetYv(), zv=vtx->GetZv();
289
290   if (TMath::Abs(zv) > 10.) return ;   
291  
292   AliPIDResponse *pidResponse = hdr->GetPIDResponse(); 
293  
294   //fMult->Fill(-100); //event counter  
295
296   //+++++++ MC
297   AliStack *stack = 0x0;
298   Double_t mcXv=0., mcYv=0., mcZv=0.;
299
300   if (fIsMC) {
301      AliMCEvent *mcEvent = MCEvent();
302      stack = mcEvent->Stack();
303      if (!stack) {
304         Printf("ERROR: stack not available");
305         return;
306      }
307
308      const AliVVertex *mcVtx=mcEvent->GetPrimaryVertex();
309
310      mcXv=mcVtx->GetX(); mcYv=mcVtx->GetY(); mcZv=mcVtx->GetZ();
311
312      Int_t ntrk=stack->GetNtrack(), ntrk0=ntrk;
313      while (ntrk--) {
314        TParticle *p0=stack->Particle(ntrk);
315        Int_t code=p0->GetPdgCode();
316        if (code != kK0Short)
317          if (code != kLambda0) continue;
318
319        Int_t plab=p0->GetFirstDaughter(), nlab=p0->GetLastDaughter();
320        if (nlab==plab) continue;
321        if (nlab<0) continue;
322        if (plab<0) continue;
323        if (nlab>=ntrk0) continue;
324        if (plab>=ntrk0) continue;
325        TParticle *part = stack->Particle(plab);
326        if (!part) continue;
327        TParticlePDG *partPDG = part->GetPDG();
328        if (!partPDG) continue;
329        Double_t charge=partPDG->Charge();
330        if (charge == 0.) continue;
331   
332        Double_t pt=p0->Pt();
333        if (pt<pMin) continue;
334        if (TMath::Abs(p0->Y())>yMax) continue;
335     
336        Double_t x=p0->Vx(), y=p0->Vy(), z=p0->Vz();
337        Double_t dx=mcXv-x, dy=mcYv-y, dz=mcZv-z;
338        Double_t l=TMath::Sqrt(dx*dx + dy*dy + dz*dz);
339
340        if (l > 0.01) continue; // secondary V0
341
342        x=part->Vx(); y=part->Vy();
343        dx=mcXv-x; dy=mcYv-y;
344        Double_t lt=TMath::Sqrt(dx*dx + dy*dy);
345
346        if (code == kK0Short) {
347           fK0sMC->Fill(pt,lt);
348        }
349        if (code == kLambda0) {
350           fLambdaMC->Fill(pt,lt);
351        }
352      }
353   }
354
355
356   Int_t ntrk=esd->GetNumberOfTracks();
357   Int_t mult=0;
358   Double_t nsig;
359   for (Int_t i=0; i<ntrk; i++) {
360     AliESDtrack *t=esd->GetTrack(i);
361     if (!t->IsOn(AliESDtrack::kTPCrefit)) continue;
362     Float_t xy,z0;
363     t->GetImpactParameters(xy,z0);
364     if (TMath::Abs(xy)>3.) continue;
365     if (TMath::Abs(z0)>3.) continue;
366     Double_t pt=t->Pt(),pz=t->Pz();
367     if (TMath::Abs(pz/pt)>0.8) continue;
368     mult++;
369
370     Double_t p=t->GetInnerParam()->GetP();
371     Double_t dedx=t->GetTPCsignal()/47.;
372     fdEdx->Fill(p,dedx,1);
373
374     nsig=pidResponse->NumberOfSigmasTPC(t,AliPID::kProton);
375     if (TMath::Abs(nsig) < 3.) fdEdxPid->Fill(p,dedx,1);
376
377   }
378   fMult->Fill(mult);
379
380
381   Int_t nv0 = esd->GetNumberOfV0s();
382   while (nv0--) {
383       AliESDv0 *v0=esd->GetV0(nv0);
384
385       if (!AcceptV0(v0,esd)) continue;
386
387       Int_t nidx=TMath::Abs(v0->GetNindex());
388       AliESDtrack *ntrack=esd->GetTrack(nidx);
389       Int_t pidx=TMath::Abs(v0->GetPindex());
390       AliESDtrack *ptrack=esd->GetTrack(pidx);
391
392       Double_t x,y,z; v0->GetXYZ(x,y,z);
393       Double_t dx=x-xv, dy=y-yv;
394       Double_t lt=TMath::Sqrt(dx*dx + dy*dy);
395
396       Double_t pt=v0->Pt();
397
398       Bool_t ctK=kTRUE; if (0.4977*lt/pt > 3*2.68) ctK=kFALSE;
399       Bool_t ctL=kTRUE; if (1.1157*lt/pt > 3*7.89) ctL=kFALSE;
400
401       //+++++++ MC
402       if (stack) {
403          Int_t ntrk=stack->GetNtrack();
404
405          Int_t nlab=TMath::Abs(ntrack->GetLabel());
406          Int_t plab=TMath::Abs(ptrack->GetLabel());
407
408          if (nlab<0) goto noas;      
409          if (nlab>=ntrk) goto noas;      
410          if (plab<0) goto noas;      
411          if (plab>=ntrk) goto noas;      
412
413          TParticle *np=stack->Particle(nlab);
414          TParticle *pp=stack->Particle(plab);
415          Int_t i0=pp->GetFirstMother();
416          //if (np->GetFirstMother() != i0) goto noas;
417
418          Int_t in0=np->GetFirstMother();
419          if (in0<0) goto noas;
420          if (in0>=ntrk) goto noas;
421          if (in0 != i0) { // did the negative daughter decay ?
422             TParticle *nnp=stack->Particle(in0);
423             if (nnp->GetFirstMother() != i0) goto noas;
424          }
425
426          if (i0<0) goto noas;
427          if (i0>=ntrk) goto noas;
428          TParticle *p0=stack->Particle(i0);
429
430          Int_t code=p0->GetPdgCode();
431          if (code != kK0Short)
432             if (code != kLambda0) goto noas;
433
434          if (p0->Pt()<pMin) goto noas;
435          if (TMath::Abs(p0->Y())>yMax ) goto noas;
436
437
438          Double_t dz=mcZv - p0->Vz(), dy=mcYv - p0->Vy(), dx=mcXv - p0->Vx();
439          Double_t l = TMath::Sqrt(dx*dx + dy*dy + dz*dz);
440
441          dx = mcXv - pp->Vx(); dy = mcYv - pp->Vy();
442          Double_t ltAs=TMath::Sqrt(dx*dx + dy*dy);
443          Double_t ptAs=p0->Pt();
444
445          if (l > 0.01) { // Secondary V0
446            if (code != kLambda0) goto noas;
447            Int_t nx=p0->GetFirstMother();
448            if (nx<0) goto noas;
449            if (nx>=ntrk) goto noas;
450            TParticle *xi=stack->Particle(nx);
451            Int_t xcode=xi->GetPdgCode();
452            if ( xcode != kXiMinus )
453              if( xcode != 3322 ) goto noas; 
454            fLambdaFromXi->Fill(ptAs,ltAs,xi->Pt());
455          } else {
456            if (code == kLambda0) {
457               if (ctL) fLambdaAs->Fill(ptAs,ltAs);
458            } else {
459               if (ctK)  fK0sAs->Fill(ptAs,ltAs);
460            } 
461          }
462  
463       }
464       //++++++++
465
466   noas:
467
468       Double_t dca=v0->GetDcaV0Daughters();
469       Double_t cpa=v0->GetV0CosineOfPointingAngle();
470
471       Double_t mass=0., m=0., s=0.;
472       if (ctK)
473       if (TMath::Abs(v0->RapK0Short())<yMax) {
474          v0->ChangeMassHypothesis(kK0Short);
475
476          mass=v0->GetEffMass();
477          fK0sM->Fill(mass,pt);
478
479          m=TDatabasePDG::Instance()->GetParticle(kK0Short)->Mass();
480          s=0.0044 + (0.008-0.0044)/(10-1)*(pt - 1.);
481          if (TMath::Abs(m-mass) < 3*s) {
482             fK0sSi->Fill(pt,lt);
483          }
484          if (TMath::Abs(m-mass + 4.5*s) < 1.5*s) {
485             fK0sSi->Fill(pt,lt,-1);
486          }
487          if (TMath::Abs(m-mass - 4.5*s) < 1.5*s) {
488             fK0sSi->Fill(pt,lt,-1);
489          }
490       }
491       
492       if (ctL)
493       if (TMath::Abs(v0->RapLambda())<yMax) {
494          Double_t p=ptrack->GetInnerParam()->GetP();
495          if (p<1.) {
496             nsig=pidResponse->NumberOfSigmasTPC(ptrack,AliPID::kProton);
497             if (TMath::Abs(nsig) > 3.) continue;
498          }
499          v0->ChangeMassHypothesis(kLambda0);
500
501          mass=v0->GetEffMass();
502          fLambdaM->Fill(mass,pt);
503
504          m=TDatabasePDG::Instance()->GetParticle(kLambda0)->Mass();
505          //s=0.0027 + (0.004-0.0027)/(10-1)*(pt-1);
506          //s=0.0015 + (0.002-0.0015)/(2.6-1)*(pt-1);
507          s=0.0023 + (0.004-0.0023)/(6-1)*(pt-1);
508          if (TMath::Abs(m-mass) < 3*s) {
509             fLambdaSi->Fill(pt,lt);
510             fCPA->Fill(cpa,1);
511             fDCA->Fill(dca,1);
512          }
513          if (TMath::Abs(m-mass + 4.5*s) < 1.5*s) {
514             fLambdaSi->Fill(pt,lt,-1);
515             fCPA->Fill(cpa,-1);
516             fDCA->Fill(dca,-1);
517          }
518          if (TMath::Abs(m-mass - 4.5*s) < 1.5*s) {
519             fLambdaSi->Fill(pt,lt,-1);
520             fCPA->Fill(cpa,-1);
521             fDCA->Fill(dca,-1);
522          }
523       }
524   }
525
526   Double_t kine0;
527   Int_t ncs=esd->GetNumberOfCascades();
528   for (Int_t i=0; i<ncs; i++) {
529       AliESDcascade *cs=esd->GetCascade(i);
530
531       if (TMath::Abs(cs->RapXi()) > yMax) continue;
532       if (!AcceptCascade(cs,esd)) continue;
533
534       AliESDv0 *v0 = (AliESDv0*)cs;
535       if (TMath::Abs(v0->RapLambda()) > yMax) continue;
536       if (!AcceptV0(v0,esd)) continue;
537
538       Double_t pt=cs->Pt();
539
540       Int_t charge=cs->Charge();      
541       if (charge < 0) {         
542          Int_t pidx=TMath::Abs(v0->GetPindex());
543          AliESDtrack *ptrack=esd->GetTrack(pidx);
544          Double_t p=ptrack->GetInnerParam()->GetP();
545          if (p<1.) {
546             nsig=pidResponse->NumberOfSigmasTPC(ptrack,AliPID::kProton);
547             if (TMath::Abs(nsig) > 3.) continue;
548          }
549          cs->ChangeMassHypothesis(kine0,kXiMinus);
550          Double_t mass=cs->GetEffMassXi();
551          pt=cs->Pt();       
552          fXiM->Fill(mass,pt);
553          Double_t m=TDatabasePDG::Instance()->GetParticle(kXiMinus)->Mass();
554          //Double_t s=0.0037;
555          Double_t s=0.002 + (0.0032-0.002)/(6-1.5)*(pt-1.5);
556          if (TMath::Abs(m-mass) < 3*s) {
557             fXiSiP->Fill(pt);
558          }
559          if (TMath::Abs(m-mass + 4.5*s) < 1.5*s) {
560             fXiSiP->Fill(pt,-1);
561          }
562          if (TMath::Abs(m-mass - 4.5*s) < 1.5*s) {
563             fXiSiP->Fill(pt,-1);
564          }
565       }
566   }
567
568 }
569
570 void AliAnalysisTaskCTauPbPb::Terminate(Option_t *)
571 {
572    // The Terminate() function is the last function to be called during
573    // a query. It always runs on the client, it can be used to present
574    // the results graphically or save the results to file.
575   
576   fOutput=(TList*)GetOutputData(1);
577   if (!fOutput) {
578      Printf("ERROR: fOutput not available");
579      return;
580   }
581  
582   fMult = dynamic_cast<TH1F*>(fOutput->FindObject("fMult")) ; 
583   if (!fMult) {
584      Printf("ERROR: fMult not available");
585      return;
586   }
587
588   fdEdx = dynamic_cast<TH2F*>(fOutput->FindObject("fdEdx")) ; 
589   if (!fdEdx) {
590      Printf("ERROR: fdEdx not available");
591      return;
592   }
593
594   fdEdxPid = dynamic_cast<TH2F*>(fOutput->FindObject("fdEdxPid")) ; 
595   if (!fdEdxPid) {
596      Printf("ERROR: fdEdxPid not available");
597      return;
598   }
599
600
601   fK0sMC = dynamic_cast<TH2F*>(fOutput->FindObject("fK0sMC")) ; 
602   if (!fK0sMC) {
603      Printf("ERROR: fK0sMC not available");
604      return;
605   }
606   TH1D *k0sMcPx=fK0sMC->ProjectionX(); k0sMcPx->Sumw2();
607   fK0sAs = dynamic_cast<TH2F*>(fOutput->FindObject("fK0sAs")) ; 
608   if (!fK0sAs) {
609      Printf("ERROR: fK0sAs not available");
610      return;
611   }
612   TH1D *k0sAsPx=fK0sAs->ProjectionX(); 
613   k0sAsPx->Sumw2(); //k0sAsPx->Scale(0.69);
614
615
616
617   fLambdaFromXi = dynamic_cast<TH3F*>(fOutput->FindObject("fLambdaFromXi")) ; 
618   if (!fLambdaFromXi) {
619      Printf("ERROR: fLambdaFromXi not available");
620      return;
621   }
622   TH1D *lambdaFromXiPx=fLambdaFromXi->ProjectionX(); lambdaFromXiPx->Sumw2();
623
624
625   fLambdaMC = dynamic_cast<TH2F*>(fOutput->FindObject("fLambdaMC")) ; 
626   if (!fLambdaMC) {
627      Printf("ERROR: fLambdaMC not available");
628      return;
629   }
630   TH1D *lambdaMcPx=fLambdaMC->ProjectionX(); lambdaMcPx->Sumw2();
631
632   fLambdaAs = dynamic_cast<TH2F*>(fOutput->FindObject("fLambdaAs")) ; 
633   if (!fLambdaAs) {
634      Printf("ERROR: fLambdaAs not available");
635      return;
636   }
637   TH1D *lambdaAsPx=fLambdaAs->ProjectionX(); 
638   lambdaAsPx->Sumw2(); //lambdaAsPx->Scale(0.64);
639
640   fLambdaSi = dynamic_cast<TH2F*>(fOutput->FindObject("fLambdaSi")) ; 
641   if (!fLambdaSi) {
642      Printf("ERROR: fLambdaSi not available");
643      return;
644   }
645   TH1D *lambdaSiPx=fLambdaSi->ProjectionX(); 
646   lambdaSiPx->SetName("fLambdaPt");
647   lambdaSiPx->Sumw2();
648
649   fLambdaEff = dynamic_cast<TH1D*>(fOutput->FindObject("fLambdaEff")) ; 
650   if (!fLambdaEff) {
651      Printf("ERROR: fLambdaEff not available");
652      return;
653   }
654   fLambdaPt = dynamic_cast<TH1D*>(fOutput->FindObject("fLambdaPt")) ; 
655   if (!fLambdaPt) {
656      Printf("ERROR: fLambdaPt not available");
657      return;
658   }
659
660
661   if (!gROOT->IsBatch()) {
662
663     TCanvas *c1 = new TCanvas("c1","Mulitplicity");
664     c1->SetLogy();
665     fMult->DrawCopy() ;
666
667     new TCanvas("c2","dE/dx");
668     fdEdx->DrawCopy() ;
669
670     new TCanvas("c3","dE/dx with PID");
671     fdEdxPid->DrawCopy() ;
672
673     if (fIsMC) {
674        /*
675        TH1D effK(*k0sAsPx); effK.SetTitle("Efficiency for K0s");
676        effK.Divide(k0sAsPx,k0sMcPx,1,1,"b");
677        new TCanvas("c4","Efficiency for K0s");
678        effK.DrawCopy("E") ;
679        */
680
681        fLambdaEff->Divide(lambdaAsPx,lambdaMcPx,1,1,"b");
682        new TCanvas("c5","Efficiency for #Lambda");
683        fLambdaEff->DrawCopy("E") ;
684
685        lambdaSiPx->Add(lambdaFromXiPx,-1);
686        lambdaSiPx->Divide(fLambdaEff);
687
688        new TCanvas("c6","Corrected #Lambda pt");
689        lambdaSiPx->SetTitle("Corrected #Lambda pt");
690       *fLambdaPt = *lambdaSiPx; 
691        fLambdaPt->SetLineColor(2);
692        fLambdaPt->DrawCopy("E");
693     
694        lambdaMcPx->DrawCopy("same");
695  
696     } else {
697        new TCanvas("c6","Raw #Lambda pt");
698        lambdaSiPx->SetTitle("Raw #Lambda pt");
699       *fLambdaPt = *lambdaSiPx; 
700        fLambdaPt->SetLineColor(2);
701        fLambdaPt->DrawCopy("E");
702     }
703   }
704 }