]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGLF/STRANGENESS/LambdaK0PbPb/AliAnalysisTaskCTauPbPb.cxx
A new version of UserExec()
[u/mrichter/AliRoot.git] / PWGLF / STRANGENESS / LambdaK0PbPb / AliAnalysisTaskCTauPbPb.cxx
CommitLineData
e26966dc 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
30extern TROOT *gROOT;
31
32ClassImp(AliAnalysisTaskCTauPbPb)
33
e26966dc 34
35//
36// This is a little task for checking the c*tau of the strange particles
37//
38
39AliAnalysisTaskCTauPbPb::AliAnalysisTaskCTauPbPb(const char *name) :
40AliAnalysisTaskSE(name),
41fIsMC(kFALSE),
42fCMin(0.),
43fCMax(90.),
44fOutput(0),
45fMult(0),
46fdEdx(0),
47fdEdxPid(0),
48
49fK0sM(0),
50fK0sSi(0),
51fK0sMC(0),
52fK0sAs(0),
53
54fLambdaM(0),
55fLambdaSi(0),
56fLambdaMC(0),
57fLambdaAs(0),
58
c13a47f5 59fLambdaEff(0),
60fLambdaPt(0),
61
0b384e23 62fLambdaBarM(0),
63fLambdaBarSi(0),
64fLambdaBarMC(0),
65fLambdaBarAs(0),
66
c13a47f5 67fLambdaBarEff(0),
68fLambdaBarPt(0),
e26966dc 69
70fLambdaFromXi(0),
71fXiM(0),
0b384e23 72fXiSiP(0),
73
74fLambdaBarFromXiBar(0),
75fXiBarM(0),
76fXiBarSiP(0)
e26966dc 77{
78 // Constructor. Initialization of pointers
79 DefineOutput(1, TList::Class());
80}
81
82void AliAnalysisTaskCTauPbPb::UserCreateOutputObjects()
83{
c13a47f5 84 Int_t nbins=100; // number of bins
85 Double_t ltMax=100.;
86 Double_t ptMax=10.;
87
e26966dc 88 fOutput = new TList();
89 fOutput->SetOwner();
90
91
92 fMult=new TH1F("fMult","Multiplicity",1100,0.,3300);
93 fMult->GetXaxis()->SetTitle("N tracks");
94 fOutput->Add(fMult);
95
96 fdEdx=new TH2F("fdEdx","dE/dx",50,0.2,3,50,0.,6.);
97 fOutput->Add(fdEdx);
98
99 fdEdxPid=new TH2F("fdEdxPid","dE/dx with PID",50,0.2,3,50,0.,6.);
100 fOutput->Add(fdEdxPid);
101
102 fK0sM =
c13a47f5 103 new TH2F("fK0sM", "Mass for K^{0}_{s}", nbins/2, 0.448, 0.548, nbins,0.,ptMax);
e26966dc 104 fK0sM->GetXaxis()->SetTitle("Mass [GeV/c]");
105 fOutput->Add(fK0sM);
106
107 fK0sSi =
108 new TH2F("fK0sSi","L_{T} vs p_{T} for K^{0}_{s}, side-band subtracted",
c13a47f5 109 nbins,0.,ptMax,nbins,0.,ltMax);
e26966dc 110 fK0sSi->GetXaxis()->SetTitle("p_{T} [GeV/c]");
111 fK0sSi->GetYaxis()->SetTitle("L_{T} [cm]");
112 fOutput->Add(fK0sSi);
113
114 fK0sMC =
115 new TH2F("fK0sMC","L_{T} vs p_{T} for K^{0}_{s}, from MC stack",
c13a47f5 116 nbins,0.,ptMax,nbins,0.,ltMax);
e26966dc 117 fK0sMC->GetXaxis()->SetTitle("p_{T} [GeV/c]");
118 fK0sMC->GetYaxis()->SetTitle("L_{T} [cm]");
119 fOutput->Add(fK0sMC);
120
121 fK0sAs =
122 new TH2F("fK0sAs", "L_{T} vs p_{T} for K^{0}_{s}, associated",
c13a47f5 123 nbins,0.,ptMax,nbins,0.,ltMax);
e26966dc 124 fK0sAs->GetXaxis()->SetTitle("p_{T} [GeV/c]");
125 fK0sAs->GetYaxis()->SetTitle("L_{T} [cm]");
126 fOutput->Add(fK0sAs);
127
128 //----------------------
129
130 fLambdaM =
c13a47f5 131 new TH2F("fLambdaM","Mass for \\Lambda", nbins, 1.065, 1.165,nbins,0.,ptMax);
e26966dc 132 fLambdaM->GetXaxis()->SetTitle("Mass [GeV/c]");
133 fOutput->Add(fLambdaM);
134
135 fLambdaSi =
136 new TH2F("fLambdaSi","L_{T} vs p_{T} for \\Lambda, side-band subtructed",
c13a47f5 137 nbins,0.,ptMax,nbins,0.,ltMax);
e26966dc 138 fLambdaSi->GetXaxis()->SetTitle("p_{T} [GeV/c]");
139 fLambdaSi->GetYaxis()->SetTitle("L_{T} [cm]");
140 fOutput->Add(fLambdaSi);
141
142 fLambdaMC =
143 new TH2F("fLambdaMC","c\\tau for \\Lambda, from MC stack",
c13a47f5 144 nbins,0.,ptMax,nbins,0.,ltMax);
e26966dc 145 fLambdaMC->GetXaxis()->SetTitle("p_{T} [GeV/c]");
146 fLambdaMC->GetYaxis()->SetTitle("L_{T} [cm]");
147 fOutput->Add(fLambdaMC);
148
149 fLambdaAs =
150 new TH2F("fLambdaAs","c\\tau for \\Lambda, associated",
c13a47f5 151 nbins,0.,ptMax,nbins,0.,ltMax);
e26966dc 152 fLambdaAs->GetXaxis()->SetTitle("p_{T} [GeV/c]");
153 fLambdaAs->GetYaxis()->SetTitle("L_{T} [cm]");
154 fOutput->Add(fLambdaAs);
155
c13a47f5 156 //----------------------
157
158 fLambdaEff=fLambdaAs->ProjectionX();
159 fLambdaEff->SetName("fLambdaEff");
160 fLambdaEff->SetTitle("Efficiency for #Lambda");
161 fOutput->Add(fLambdaEff);
162
163 fLambdaPt=fLambdaAs->ProjectionX();
164 fLambdaPt->SetName("fLambdaPt");
165 fLambdaPt->SetTitle("Raw #Lambda pT spectrum");
166 fOutput->Add(fLambdaPt);
167
168 //----------------------
0b384e23 169
170 fLambdaBarM =
c13a47f5 171 new TH2F("fLambdaBarM","Mass for anti-\\Lambda", nbins, 1.065, 1.165,nbins,0.,ptMax);
0b384e23 172 fLambdaBarM->GetXaxis()->SetTitle("Mass [GeV/c]");
173 fOutput->Add(fLambdaBarM);
174
175 fLambdaBarSi =
176 new TH2F("fLambdaBarSi","L_{T} vs p_{T} for anti-\\Lambda, side-band subtructed",
c13a47f5 177 nbins,0.,ptMax,nbins,0.,ltMax);
0b384e23 178 fLambdaBarSi->GetXaxis()->SetTitle("p_{T} [GeV/c]");
179 fLambdaBarSi->GetYaxis()->SetTitle("L_{T} [cm]");
180 fOutput->Add(fLambdaBarSi);
181
182 fLambdaBarMC =
183 new TH2F("fLambdaBarMC","c\\tau for anti-\\Lambda, from MC stack",
c13a47f5 184 nbins,0.,ptMax,nbins,0.,ltMax);
0b384e23 185 fLambdaBarMC->GetXaxis()->SetTitle("p_{T} [GeV/c]");
186 fLambdaBarMC->GetYaxis()->SetTitle("L_{T} [cm]");
187 fOutput->Add(fLambdaBarMC);
188
189 fLambdaBarAs =
190 new TH2F("fLambdaBarAs","c\\tau for anti-\\Lambda, associated",
c13a47f5 191 nbins,0.,ptMax,nbins,0.,ltMax);
0b384e23 192 fLambdaBarAs->GetXaxis()->SetTitle("p_{T} [GeV/c]");
193 fLambdaBarAs->GetYaxis()->SetTitle("L_{T} [cm]");
194 fOutput->Add(fLambdaBarAs);
195
196
c13a47f5 197 //----------------------
0b384e23 198
c13a47f5 199 fLambdaBarEff=fLambdaBarAs->ProjectionX();
200 fLambdaBarEff->SetName("fLambdaBarEff");
201 fLambdaBarEff->SetTitle("Efficiency for anti-#Lambda");
202 fOutput->Add(fLambdaBarEff);
e26966dc 203
c13a47f5 204 fLambdaBarPt=fLambdaBarAs->ProjectionX();
205 fLambdaBarPt->SetName("fLambdaBarPt");
206 fLambdaBarPt->SetTitle("Raw anti-#Lambda pT spectrum");
207 fOutput->Add(fLambdaBarPt);
e26966dc 208
209 //----------------------
210
0b384e23 211 fLambdaFromXi=new TH3F("fLambdaFromXi","L_{T} vs p_{T} vs p_{T} of \\Xi for \\Lambda from \\Xi",
c13a47f5 212 nbins,0.,ptMax,nbins,0.,ltMax,33,0.,ptMax+2);
e26966dc 213 fOutput->Add(fLambdaFromXi);
214
215 fXiM =
c13a47f5 216 new TH2F("fXiM", "\\Xi mass distribution", 50, 1.271, 1.371,12,0.,ptMax+2);
e26966dc 217 fOutput->Add(fXiM);
218
219 fXiSiP = new TH1F("fXiSiP", "Pt for \\Xi, side-band subracted",
c13a47f5 220 33,0.,ptMax+2);
e26966dc 221 fOutput->Add(fXiSiP);
222
223
0b384e23 224 fLambdaBarFromXiBar=new TH3F("fLambdaBarFromXiBar","L_{T} vs p_{T} vs p_{T} of anti-\\Xi for anti-\\Lambda from anti-\\Xi",
c13a47f5 225 nbins,0.,ptMax,nbins,0.,ltMax,33,0.,ptMax+2);
0b384e23 226 fOutput->Add(fLambdaBarFromXiBar);
227
228 fXiBarM =
c13a47f5 229 new TH2F("fXiBarM", "anti-\\Xi mass distribution", 50, 1.271, 1.371,12,0.,ptMax+2);
0b384e23 230 fOutput->Add(fXiBarM);
231
232 fXiBarSiP = new TH1F("fXiBarSiP", "Pt for anti-\\Xi, side-band subracted",
c13a47f5 233 33,0.,ptMax+2);
0b384e23 234 fOutput->Add(fXiBarSiP);
235
236
e26966dc 237 PostData(1, fOutput);
238}
239
240static Bool_t AcceptTrack(const AliESDtrack *t) {
241 if (!t->IsOn(AliESDtrack::kTPCrefit)) return kFALSE;
242 if (t->GetKinkIndex(0)>0) return kFALSE;
243
244 Float_t nCrossedRowsTPC = t->GetTPCClusterInfo(2,1);
245 if (nCrossedRowsTPC < 70) return kFALSE;
246 Int_t findable=t->GetTPCNclsF();
247 if (findable <= 0) return kFALSE;
248 if (nCrossedRowsTPC/findable < 0.8) return kFALSE;
249
ff909752 250 if (TMath::Abs(t->Eta()) > 0.8) return kFALSE;
251
e26966dc 252 return kTRUE;
253}
254
255static Bool_t AcceptV0(const AliESDv0 *v0, const AliESDEvent *esd) {
256
257 if (v0->GetOnFlyStatus()) return kFALSE;
258
e26966dc 259 Int_t nidx=TMath::Abs(v0->GetNindex());
260 AliESDtrack *ntrack=esd->GetTrack(nidx);
261 if (!AcceptTrack(ntrack)) return kFALSE;
262
263 Int_t pidx=TMath::Abs(v0->GetPindex());
264 AliESDtrack *ptrack=esd->GetTrack(pidx);
265 if (!AcceptTrack(ptrack)) return kFALSE;
266
267 Float_t xy,z0;
268 ntrack->GetImpactParameters(xy,z0);
269 if (TMath::Abs(xy)<0.1) return kFALSE;
270 ptrack->GetImpactParameters(xy,z0);
271 if (TMath::Abs(xy)<0.1) return kFALSE;
272
273 Double_t dca=v0->GetDcaV0Daughters();
274 if (dca>1.0) return kFALSE;
275 //if (dca>0.7) return kFALSE;
276 //if (dca>0.4) return kFALSE;
277
278 Double_t cpa=v0->GetV0CosineOfPointingAngle();
279 if (cpa<0.998) return kFALSE;
280 //if (cpa<0.99875) return kFALSE;
281 //if (cpa<0.9995) return kFALSE;
282
283 Double_t xx,yy,zz; v0->GetXYZ(xx,yy,zz);
284 Double_t r2=xx*xx + yy*yy;
285 if (r2<0.9*0.9) return kFALSE;
286 if (r2>100*100) return kFALSE;
287
288 return kTRUE;
289}
290
291static Bool_t AcceptCascade(const AliESDcascade *cs, const AliESDEvent *esd) {
292
e26966dc 293 Int_t bidx=TMath::Abs(cs->GetBindex());
294 AliESDtrack *btrack=esd->GetTrack(bidx);
295 if (!AcceptTrack(btrack)) return kFALSE;
296
297 Float_t xy,z0;
298 btrack->GetImpactParameters(xy,z0);
299 if (TMath::Abs(xy)<0.03) return kFALSE;
300
301 const AliESDVertex *vtx=esd->GetPrimaryVertexSPD();
302 if (!vtx->GetStatus()) {
303 vtx=esd->GetPrimaryVertexTracks();
304 if (!vtx->GetStatus()) return kFALSE;
305 }
306 Double_t xv=vtx->GetXv(), yv=vtx->GetYv(), zv=vtx->GetZv();
307 if (cs->GetCascadeCosineOfPointingAngle(xv,yv,zv) < 0.999) return kFALSE;
308
309 if (cs->GetDcaXiDaughters() > 0.3) return kFALSE;
310
311 return kTRUE;
312}
313
1c102c9a 314static Bool_t AcceptPID(const AliPIDResponse *pidResponse,
315 const AliESDtrack *ptrack, AliStack *stack) {
316
f6357875 317 const AliExternalTrackParam *par=ptrack->GetInnerParam();
318 if (!par) return kTRUE;
319 if (par->GetP() > 1.) return kTRUE;
1c102c9a 320
c13a47f5 321
1c102c9a 322 if (stack) {
323 // MC PID
1c102c9a 324 Int_t plab=TMath::Abs(ptrack->GetLabel());
c13a47f5 325 TParticle *pp=stack->Particle(plab);
326 if (!pp) return kTRUE;
327 if (pp->GetPDG()->Charge() > 0) {
328 if (pp->GetPdgCode() == kProton) return kTRUE;
329 } else {
330 if (pp->GetPdgCode() == kProtonBar) return kTRUE;
331 }
1c102c9a 332 } else {
333 // Real PID
f6357875 334 Double_t nsig=pidResponse->NumberOfSigmasTPC(ptrack,AliPID::kProton);
335 if (TMath::Abs(nsig) < 3.) return kTRUE;
1c102c9a 336 }
337
f6357875 338 return kFALSE;
1c102c9a 339}
340
c13a47f5 341static TParticle*
342Associate(const AliESDtrack *ptrack,const AliESDtrack *ntrack,AliStack *stack,
343TParticle *&mcp) {
344 //
345 // Try to associate the V0 with the daughters ptrack and ntrack
346 // with the Monte Carlo
347 //
348 if (!stack) return 0;
349
350 Int_t nlab=TMath::Abs(ntrack->GetLabel());
351 TParticle *n=stack->Particle(nlab);
352 if (!n) return 0;
353
354 Int_t plab=TMath::Abs(ptrack->GetLabel());
355 TParticle *p=stack->Particle(plab);
356 if (!p) return 0;
357
358 Int_t imp=p->GetFirstMother();
359 if (imp<0) return 0;
360 if (imp>=stack->GetNtrack()) return 0;
361 TParticle *p0=stack->Particle(imp); // V0 particle, mother of the pos. track
362 if (!p0) return 0;
363
364 Int_t imn=n->GetFirstMother();
365 if (imp != imn) { // Check decays of the daughters
366 return 0; // Fixme
367 }
368
369 Int_t code=p0->GetPdgCode();
370 if (code != kK0Short)
371 if (code != kLambda0)
372 if (code != kLambda0Bar) return 0;
373
374 mcp=p;
375
376 return p0;
377
378}
379
380
e26966dc 381void AliAnalysisTaskCTauPbPb::UserExec(Option_t *)
382{
c13a47f5 383 const Double_t yMax=0.5;
384 const Double_t pMin=0.0;
385 const Double_t lMax=0.001;
e26966dc 386
387 AliESDEvent *esd=(AliESDEvent *)InputEvent();
388
389 if (!esd) {
390 Printf("ERROR: esd not available");
391 return;
392 }
393
0b384e23 394 // Vertex selection
395 const AliESDVertex *vtx=esd->GetPrimaryVertexSPD();
396 if (!vtx->GetStatus()) {
397 vtx=esd->GetPrimaryVertexTracks();
398 if (!vtx->GetStatus()) return;
399 }
400 Double_t xv=vtx->GetXv(), yv=vtx->GetYv(), zv=vtx->GetZv();
401
402 if (TMath::Abs(zv) > 10.) return ;
403
e26966dc 404
405 // Physics selection
406 AliAnalysisManager *mgr= AliAnalysisManager::GetAnalysisManager();
407 AliInputEventHandler *hdr=(AliInputEventHandler*)mgr->GetInputEventHandler();
408 UInt_t maskIsSelected = hdr->IsEventSelected();
409 Bool_t isSelected = (maskIsSelected & AliVEvent::kMB);
410 if (!isSelected) return;
411
0b384e23 412
413 fMult->Fill(-100); //event counter
414
415
e26966dc 416 // Centrality selection
417 AliCentrality *cent=esd->GetCentrality();
418 if (!cent->IsEventInCentralityClass(fCMin,fCMax,"V0M")) return;
419
e26966dc 420 AliPIDResponse *pidResponse = hdr->GetPIDResponse();
421
e26966dc 422
423 //+++++++ MC
424 AliStack *stack = 0x0;
425 Double_t mcXv=0., mcYv=0., mcZv=0.;
426
427 if (fIsMC) {
428 AliMCEvent *mcEvent = MCEvent();
429 stack = mcEvent->Stack();
430 if (!stack) {
431 Printf("ERROR: stack not available");
432 return;
433 }
434
435 const AliVVertex *mcVtx=mcEvent->GetPrimaryVertex();
436
437 mcXv=mcVtx->GetX(); mcYv=mcVtx->GetY(); mcZv=mcVtx->GetZ();
438
c13a47f5 439 Int_t ntrk=stack->GetNtrack();
e26966dc 440 while (ntrk--) {
441 TParticle *p0=stack->Particle(ntrk);
442 Int_t code=p0->GetPdgCode();
443 if (code != kK0Short)
0b384e23 444 if (code != kLambda0)
445 if (code != kLambda0Bar) continue;
e26966dc 446
447 Int_t plab=p0->GetFirstDaughter(), nlab=p0->GetLastDaughter();
448 if (nlab==plab) continue;
e26966dc 449 TParticle *part = stack->Particle(plab);
450 if (!part) continue;
451 TParticlePDG *partPDG = part->GetPDG();
452 if (!partPDG) continue;
453 Double_t charge=partPDG->Charge();
454 if (charge == 0.) continue;
455
456 Double_t pt=p0->Pt();
457 if (pt<pMin) continue;
458 if (TMath::Abs(p0->Y())>yMax) continue;
459
460 Double_t x=p0->Vx(), y=p0->Vy(), z=p0->Vz();
461 Double_t dx=mcXv-x, dy=mcYv-y, dz=mcZv-z;
462 Double_t l=TMath::Sqrt(dx*dx + dy*dy + dz*dz);
463
c13a47f5 464 if (l > lMax) continue; // secondary V0
e26966dc 465
466 x=part->Vx(); y=part->Vy();
467 dx=mcXv-x; dy=mcYv-y;
468 Double_t lt=TMath::Sqrt(dx*dx + dy*dy);
469
c13a47f5 470 switch (code) {
471 case kK0Short:
e26966dc 472 fK0sMC->Fill(pt,lt);
c13a47f5 473 break;
474 case kLambda0:
e26966dc 475 fLambdaMC->Fill(pt,lt);
c13a47f5 476 break;
477 case kLambda0Bar:
0b384e23 478 fLambdaBarMC->Fill(pt,lt);
c13a47f5 479 break;
480 default: break;
0b384e23 481 }
e26966dc 482 }
483 }
0b384e23 484 //+++++++
e26966dc 485
486
076bd7b1 487 Int_t ntrk1=esd->GetNumberOfTracks();
e26966dc 488 Int_t mult=0;
076bd7b1 489 for (Int_t i=0; i<ntrk1; i++) {
e26966dc 490 AliESDtrack *t=esd->GetTrack(i);
491 if (!t->IsOn(AliESDtrack::kTPCrefit)) continue;
492 Float_t xy,z0;
493 t->GetImpactParameters(xy,z0);
494 if (TMath::Abs(xy)>3.) continue;
495 if (TMath::Abs(z0)>3.) continue;
496 Double_t pt=t->Pt(),pz=t->Pz();
497 if (TMath::Abs(pz/pt)>0.8) continue;
498 mult++;
499
500 Double_t p=t->GetInnerParam()->GetP();
501 Double_t dedx=t->GetTPCsignal()/47.;
502 fdEdx->Fill(p,dedx,1);
503
1c102c9a 504 Double_t nsig=pidResponse->NumberOfSigmasTPC(t,AliPID::kProton);
e26966dc 505 if (TMath::Abs(nsig) < 3.) fdEdxPid->Fill(p,dedx,1);
506
507 }
508 fMult->Fill(mult);
509
510
511 Int_t nv0 = esd->GetNumberOfV0s();
512 while (nv0--) {
c13a47f5 513
514 Bool_t isK0s=kTRUE;
515 Bool_t isLambda=kTRUE;
516 Bool_t isLambdaBar=kTRUE;
517
e26966dc 518 AliESDv0 *v0=esd->GetV0(nv0);
519
c13a47f5 520 Double_t pt=v0->Pt();
521 if (pt < pMin) continue;
522
e26966dc 523 if (!AcceptV0(v0,esd)) continue;
524
525 Int_t nidx=TMath::Abs(v0->GetNindex());
526 AliESDtrack *ntrack=esd->GetTrack(nidx);
527 Int_t pidx=TMath::Abs(v0->GetPindex());
528 AliESDtrack *ptrack=esd->GetTrack(pidx);
529
530 Double_t x,y,z; v0->GetXYZ(x,y,z);
076bd7b1 531 Double_t dx1=x-xv, dy1=y-yv;
532 Double_t lt=TMath::Sqrt(dx1*dx1 + dy1*dy1);
e26966dc 533
c13a47f5 534 if (lt/pt > 3*7.89/1.1157) continue;
e26966dc 535
c13a47f5 536 if (0.4977*lt/pt > 3*2.68) isK0s=kFALSE;
537 if (1.1157*lt/pt > 3*7.89) isLambdaBar=isLambda=kFALSE;
e26966dc 538
c13a47f5 539 isLambda = isLambda && AcceptPID(pidResponse, ptrack, stack);
540 isLambdaBar = isLambdaBar && AcceptPID(pidResponse, ntrack, stack);
1c102c9a 541
c13a47f5 542 Double_t yK0s=TMath::Abs(v0->RapK0Short());
543 Double_t yLam=TMath::Abs(v0->RapLambda());
544 isK0s = isK0s && (yK0s < yMax);
545 isLambda = isLambda && (yLam < yMax);
546 isLambdaBar = isLambdaBar && (yLam < yMax);
e26966dc 547
548 Double_t mass=0., m=0., s=0.;
c13a47f5 549 if (isK0s) {
e26966dc 550 v0->ChangeMassHypothesis(kK0Short);
551
552 mass=v0->GetEffMass();
553 fK0sM->Fill(mass,pt);
554
555 m=TDatabasePDG::Instance()->GetParticle(kK0Short)->Mass();
556 s=0.0044 + (0.008-0.0044)/(10-1)*(pt - 1.);
557 if (TMath::Abs(m-mass) < 3*s) {
558 fK0sSi->Fill(pt,lt);
c13a47f5 559 } else {
560 isK0s=kFALSE;
e26966dc 561 }
562 if (TMath::Abs(m-mass + 4.5*s) < 1.5*s) {
563 fK0sSi->Fill(pt,lt,-1);
564 }
565 if (TMath::Abs(m-mass - 4.5*s) < 1.5*s) {
566 fK0sSi->Fill(pt,lt,-1);
567 }
568 }
569
c13a47f5 570 if (isLambda) {
e26966dc 571 v0->ChangeMassHypothesis(kLambda0);
572
573 mass=v0->GetEffMass();
574 fLambdaM->Fill(mass,pt);
575
576 m=TDatabasePDG::Instance()->GetParticle(kLambda0)->Mass();
577 //s=0.0027 + (0.004-0.0027)/(10-1)*(pt-1);
578 //s=0.0015 + (0.002-0.0015)/(2.6-1)*(pt-1);
579 s=0.0023 + (0.004-0.0023)/(6-1)*(pt-1);
580 if (TMath::Abs(m-mass) < 3*s) {
581 fLambdaSi->Fill(pt,lt);
c13a47f5 582 } else {
583 isLambda=kFALSE;
584 }
e26966dc 585 if (TMath::Abs(m-mass + 4.5*s) < 1.5*s) {
586 fLambdaSi->Fill(pt,lt,-1);
e26966dc 587 }
588 if (TMath::Abs(m-mass - 4.5*s) < 1.5*s) {
589 fLambdaSi->Fill(pt,lt,-1);
e26966dc 590 }
591 }
0b384e23 592
c13a47f5 593 if (isLambdaBar) {
0b384e23 594 v0->ChangeMassHypothesis(kLambda0Bar);
595
596 mass=v0->GetEffMass();
597 fLambdaBarM->Fill(mass,pt);
598
599 m=TDatabasePDG::Instance()->GetParticle(kLambda0Bar)->Mass();
600 //s=0.0027 + (0.004-0.0027)/(10-1)*(pt-1);
601 //s=0.0015 + (0.002-0.0015)/(2.6-1)*(pt-1);
602 s=0.0023 + (0.004-0.0023)/(6-1)*(pt-1);
603 if (TMath::Abs(m-mass) < 3*s) {
604 fLambdaBarSi->Fill(pt,lt);
c13a47f5 605 } else {
606 isLambdaBar=kFALSE;
0b384e23 607 }
608 if (TMath::Abs(m-mass + 4.5*s) < 1.5*s) {
609 fLambdaBarSi->Fill(pt,lt,-1);
0b384e23 610 }
611 if (TMath::Abs(m-mass - 4.5*s) < 1.5*s) {
612 fLambdaBarSi->Fill(pt,lt,-1);
0b384e23 613 }
614 }
c13a47f5 615
616 if (!fIsMC) continue;
617
618 //++++++ MC
619 if (!isK0s)
620 if (!isLambda)
621 if (!isLambdaBar) continue;//check MC only for the accepted V0s
622
623 TParticle *mcp=0;
624 TParticle *mc0=Associate(ptrack,ntrack,stack,mcp);
625 if (!mc0) continue;
626
627 Double_t ptAs=mc0->Pt();
628 if (ptAs < pMin) continue;
629 Double_t yAs=mc0->Y();
630 if (TMath::Abs(yAs) > yMax) continue;
631
632 Int_t code=mc0->GetPdgCode();
633
634 Double_t dx = mcXv - mcp->Vx(), dy = mcYv - mcp->Vy();
635 Double_t ltAs=TMath::Sqrt(dx*dx + dy*dy);
636
637 Double_t dz=mcZv - mc0->Vz(); dy=mcYv - mc0->Vy(); dx=mcXv - mc0->Vx();
638 Double_t l = TMath::Sqrt(dx*dx + dy*dy + dz*dz);
639 if (l<lMax) { // Primary V0s
640 switch (code) {
641 case kK0Short:
642 if (isK0s) fK0sAs->Fill(ptAs,ltAs);
643 break;
644 case kLambda0:
645 if (isLambda) fLambdaAs->Fill(ptAs,ltAs);
646 break;
647 case kLambda0Bar:
648 if (isLambdaBar) fLambdaBarAs->Fill(ptAs,ltAs);
649 break;
650 default: break;
651 }
652 } else {
653 if (code==kK0Short) continue;
654
655 Int_t nx=mc0->GetFirstMother();
656 if (nx<0) continue;
657 if (nx>=stack->GetNtrack()) continue;
658 TParticle *xi=stack->Particle(nx);
659 if (!xi) continue;
660 Int_t xcode=xi->GetPdgCode();
661
662 switch (code) {
663 case kLambda0:
664 if (isLambda)
665 if ((xcode==kXiMinus) || (xcode==3322))
666 fLambdaFromXi->Fill(ptAs,ltAs,xi->Pt());
667 break;
668 case kLambda0Bar:
669 if (isLambdaBar)
670 if ((xcode==kXiPlusBar)||(xcode==-3322))
671 fLambdaBarFromXiBar->Fill(ptAs,ltAs,xi->Pt());
672 break;
673 default: break;
674 }
675 }
676 //++++++
677
e26966dc 678 }
679
680 Double_t kine0;
681 Int_t ncs=esd->GetNumberOfCascades();
682 for (Int_t i=0; i<ncs; i++) {
683 AliESDcascade *cs=esd->GetCascade(i);
684
c13a47f5 685 if (cs->Pt() < pMin) continue;
e26966dc 686 if (TMath::Abs(cs->RapXi()) > yMax) continue;
687 if (!AcceptCascade(cs,esd)) continue;
688
689 AliESDv0 *v0 = (AliESDv0*)cs;
c13a47f5 690 if (v0->Pt() < pMin) continue;
e26966dc 691 if (TMath::Abs(v0->RapLambda()) > yMax) continue;
692 if (!AcceptV0(v0,esd)) continue;
693
694 Double_t pt=cs->Pt();
695
1c102c9a 696 Int_t pidx=TMath::Abs(v0->GetPindex());
697 AliESDtrack *ptrack=esd->GetTrack(pidx);
0b384e23 698 Bool_t isProton =AcceptPID(pidResponse, ptrack, stack);
1c102c9a 699
0b384e23 700 Int_t nidx=TMath::Abs(v0->GetNindex());
701 AliESDtrack *ntrack=esd->GetTrack(nidx);
702 Bool_t isProtonBar=AcceptPID(pidResponse, ntrack, stack);
1c102c9a 703
e26966dc 704 Int_t charge=cs->Charge();
1c102c9a 705 if (isProton)
e26966dc 706 if (charge < 0) {
e26966dc 707 cs->ChangeMassHypothesis(kine0,kXiMinus);
708 Double_t mass=cs->GetEffMassXi();
709 pt=cs->Pt();
710 fXiM->Fill(mass,pt);
711 Double_t m=TDatabasePDG::Instance()->GetParticle(kXiMinus)->Mass();
712 //Double_t s=0.0037;
713 Double_t s=0.002 + (0.0032-0.002)/(6-1.5)*(pt-1.5);
714 if (TMath::Abs(m-mass) < 3*s) {
715 fXiSiP->Fill(pt);
716 }
717 if (TMath::Abs(m-mass + 4.5*s) < 1.5*s) {
718 fXiSiP->Fill(pt,-1);
719 }
720 if (TMath::Abs(m-mass - 4.5*s) < 1.5*s) {
721 fXiSiP->Fill(pt,-1);
722 }
723 }
0b384e23 724 if (isProtonBar)
725 if (charge > 0) {
726 cs->ChangeMassHypothesis(kine0,kXiPlusBar);
727 Double_t mass=cs->GetEffMassXi();
728 pt=cs->Pt();
729 fXiBarM->Fill(mass,pt);
730 Double_t m=TDatabasePDG::Instance()->GetParticle(kXiPlusBar)->Mass();
731 //Double_t s=0.0037;
732 Double_t s=0.002 + (0.0032-0.002)/(6-1.5)*(pt-1.5);
733 if (TMath::Abs(m-mass) < 3*s) {
734 fXiBarSiP->Fill(pt);
735 }
736 if (TMath::Abs(m-mass + 4.5*s) < 1.5*s) {
737 fXiBarSiP->Fill(pt,-1);
738 }
739 if (TMath::Abs(m-mass - 4.5*s) < 1.5*s) {
740 fXiBarSiP->Fill(pt,-1);
741 }
742 }
e26966dc 743 }
744
745}
746
747void AliAnalysisTaskCTauPbPb::Terminate(Option_t *)
748{
749 // The Terminate() function is the last function to be called during
750 // a query. It always runs on the client, it can be used to present
751 // the results graphically or save the results to file.
752
753 fOutput=(TList*)GetOutputData(1);
754 if (!fOutput) {
755 Printf("ERROR: fOutput not available");
756 return;
757 }
758
759 fMult = dynamic_cast<TH1F*>(fOutput->FindObject("fMult")) ;
760 if (!fMult) {
761 Printf("ERROR: fMult not available");
762 return;
763 }
764
765 fdEdx = dynamic_cast<TH2F*>(fOutput->FindObject("fdEdx")) ;
766 if (!fdEdx) {
767 Printf("ERROR: fdEdx not available");
768 return;
769 }
770
771 fdEdxPid = dynamic_cast<TH2F*>(fOutput->FindObject("fdEdxPid")) ;
772 if (!fdEdxPid) {
773 Printf("ERROR: fdEdxPid not available");
774 return;
775 }
776
777
778 fK0sMC = dynamic_cast<TH2F*>(fOutput->FindObject("fK0sMC")) ;
779 if (!fK0sMC) {
780 Printf("ERROR: fK0sMC not available");
781 return;
782 }
783 TH1D *k0sMcPx=fK0sMC->ProjectionX(); k0sMcPx->Sumw2();
784 fK0sAs = dynamic_cast<TH2F*>(fOutput->FindObject("fK0sAs")) ;
785 if (!fK0sAs) {
786 Printf("ERROR: fK0sAs not available");
787 return;
788 }
789 TH1D *k0sAsPx=fK0sAs->ProjectionX();
790 k0sAsPx->Sumw2(); //k0sAsPx->Scale(0.69);
791
792
793
c13a47f5 794 // Lambda histograms
e26966dc 795 fLambdaFromXi = dynamic_cast<TH3F*>(fOutput->FindObject("fLambdaFromXi")) ;
796 if (!fLambdaFromXi) {
797 Printf("ERROR: fLambdaFromXi not available");
798 return;
799 }
800 TH1D *lambdaFromXiPx=fLambdaFromXi->ProjectionX(); lambdaFromXiPx->Sumw2();
801
e26966dc 802 fLambdaMC = dynamic_cast<TH2F*>(fOutput->FindObject("fLambdaMC")) ;
803 if (!fLambdaMC) {
804 Printf("ERROR: fLambdaMC not available");
805 return;
806 }
807 TH1D *lambdaMcPx=fLambdaMC->ProjectionX(); lambdaMcPx->Sumw2();
808
809 fLambdaAs = dynamic_cast<TH2F*>(fOutput->FindObject("fLambdaAs")) ;
810 if (!fLambdaAs) {
811 Printf("ERROR: fLambdaAs not available");
812 return;
813 }
814 TH1D *lambdaAsPx=fLambdaAs->ProjectionX();
815 lambdaAsPx->Sumw2(); //lambdaAsPx->Scale(0.64);
816
817 fLambdaSi = dynamic_cast<TH2F*>(fOutput->FindObject("fLambdaSi")) ;
818 if (!fLambdaSi) {
819 Printf("ERROR: fLambdaSi not available");
820 return;
821 }
822 TH1D *lambdaSiPx=fLambdaSi->ProjectionX();
823 lambdaSiPx->SetName("fLambdaPt");
824 lambdaSiPx->Sumw2();
825
826 fLambdaEff = dynamic_cast<TH1D*>(fOutput->FindObject("fLambdaEff")) ;
827 if (!fLambdaEff) {
828 Printf("ERROR: fLambdaEff not available");
829 return;
830 }
831 fLambdaPt = dynamic_cast<TH1D*>(fOutput->FindObject("fLambdaPt")) ;
832 if (!fLambdaPt) {
833 Printf("ERROR: fLambdaPt not available");
834 return;
835 }
836
837
c13a47f5 838 // anti-Lambda histograms
839 fLambdaBarFromXiBar =
840 dynamic_cast<TH3F*>(fOutput->FindObject("fLambdaBarFromXiBar")) ;
841 if (!fLambdaBarFromXiBar) {
842 Printf("ERROR: fLambdaBarFromXiBar not available");
843 return;
844 }
845 TH1D *lambdaBarFromXiBarPx=
846 fLambdaBarFromXiBar->ProjectionX("Bar"); lambdaBarFromXiBarPx->Sumw2();
847
848 fLambdaBarMC = dynamic_cast<TH2F*>(fOutput->FindObject("fLambdaBarMC")) ;
849 if (!fLambdaBarMC) {
850 Printf("ERROR: fLambdaBarMC not available");
851 return;
852 }
853 TH1D *lambdaBarMcPx=fLambdaBarMC->ProjectionX(); lambdaBarMcPx->Sumw2();
854
855 fLambdaBarAs = dynamic_cast<TH2F*>(fOutput->FindObject("fLambdaBarAs")) ;
856 if (!fLambdaBarAs) {
857 Printf("ERROR: fLambdaBarAs not available");
858 return;
859 }
860 TH1D *lambdaBarAsPx=fLambdaBarAs->ProjectionX();
861 lambdaBarAsPx->Sumw2(); //lambdaBarAsPx->Scale(0.64);
862
863 fLambdaBarSi = dynamic_cast<TH2F*>(fOutput->FindObject("fLambdaBarSi")) ;
864 if (!fLambdaBarSi) {
865 Printf("ERROR: fLambdaBarSi not available");
866 return;
867 }
868 TH1D *lambdaBarSiPx=fLambdaBarSi->ProjectionX();
869 lambdaBarSiPx->SetName("fLambdaBarPt");
870 lambdaBarSiPx->Sumw2();
871
872 fLambdaBarEff = dynamic_cast<TH1D*>(fOutput->FindObject("fLambdaBarEff")) ;
873 if (!fLambdaBarEff) {
874 Printf("ERROR: fLambdaBarEff not available");
875 return;
876 }
877 fLambdaBarPt = dynamic_cast<TH1D*>(fOutput->FindObject("fLambdaBarPt")) ;
878 if (!fLambdaBarPt) {
879 Printf("ERROR: fLambdaBarPt not available");
880 return;
881 }
882
883
e26966dc 884 if (!gROOT->IsBatch()) {
885
886 TCanvas *c1 = new TCanvas("c1","Mulitplicity");
887 c1->SetLogy();
888 fMult->DrawCopy() ;
889
890 new TCanvas("c2","dE/dx");
891 fdEdx->DrawCopy() ;
892
893 new TCanvas("c3","dE/dx with PID");
894 fdEdxPid->DrawCopy() ;
895
896 if (fIsMC) {
897 /*
898 TH1D effK(*k0sAsPx); effK.SetTitle("Efficiency for K0s");
899 effK.Divide(k0sAsPx,k0sMcPx,1,1,"b");
900 new TCanvas("c4","Efficiency for K0s");
901 effK.DrawCopy("E") ;
902 */
903
c13a47f5 904 //+++ Lambdas
e26966dc 905 fLambdaEff->Divide(lambdaAsPx,lambdaMcPx,1,1,"b");
906 new TCanvas("c5","Efficiency for #Lambda");
907 fLambdaEff->DrawCopy("E") ;
908
909 lambdaSiPx->Add(lambdaFromXiPx,-1);
910 lambdaSiPx->Divide(fLambdaEff);
911
912 new TCanvas("c6","Corrected #Lambda pt");
913 lambdaSiPx->SetTitle("Corrected #Lambda pt");
914 *fLambdaPt = *lambdaSiPx;
915 fLambdaPt->SetLineColor(2);
916 fLambdaPt->DrawCopy("E");
917
918 lambdaMcPx->DrawCopy("same");
919
c13a47f5 920
921 //+++ anti-Lambdas
922 fLambdaBarEff->Divide(lambdaBarAsPx,lambdaBarMcPx,1,1,"b");
923 new TCanvas("c7","Efficiency for anti-#Lambda");
924 fLambdaBarEff->DrawCopy("E") ;
925
926 lambdaBarSiPx->Add(lambdaBarFromXiBarPx,-1);
927 lambdaBarSiPx->Divide(fLambdaBarEff);
928
929 new TCanvas("c8","Corrected anti-#Lambda pt");
930 lambdaBarSiPx->SetTitle("Corrected anti-#Lambda pt");
931 *fLambdaBarPt = *lambdaBarSiPx;
932 fLambdaBarPt->SetLineColor(2);
933 fLambdaBarPt->DrawCopy("E");
934
935 lambdaBarMcPx->DrawCopy("same");
e26966dc 936 } else {
937 new TCanvas("c6","Raw #Lambda pt");
938 lambdaSiPx->SetTitle("Raw #Lambda pt");
939 *fLambdaPt = *lambdaSiPx;
940 fLambdaPt->SetLineColor(2);
941 fLambdaPt->DrawCopy("E");
c13a47f5 942
943
944 new TCanvas("c7","Raw anti-#Lambda pt");
945 lambdaBarSiPx->SetTitle("Raw anti-#Lambda pt");
946 *fLambdaBarPt = *lambdaBarSiPx;
947 fLambdaBarPt->SetLineColor(2);
948 fLambdaBarPt->DrawCopy("E");
e26966dc 949 }
950 }
951}