1 /**************************************************************************
\r
2 * Author : Benjamin Dönigus (benjamin.doenigus@cern.ch) *
\r
4 * Contributors are mentioned in the code where appropriate. *
\r
6 * Permission to use, copy, modify and distribute this software and its *
\r
7 * documentation strictly for non-commercial purposes is hereby granted *
\r
8 * without fee, provided that the above copyright notice appears in all *
\r
9 * copies and that both the copyright notice and this permission notice *
\r
10 * appear in the supporting documentation. The authors make no claims *
\r
11 * about the suitability of this software for any purpose. It is *
\r
12 * provided "as is" without express or implied warranty. *
\r
13 **************************************************************************/
\r
15 //-----------------------------------------------------------------
\r
16 // AliAnalysisTaskHdibaryonLPpi class
\r
17 // used to search for the H-Dibaryon in weak
\r
18 // (Lambda Proton Pion) and strong (Lambda Lambda) decays
\r
19 //-----------------------------------------------------------------
\r
21 #include "Riostream.h"
\r
25 #include "TSystem.h"
\r
31 #include "TCanvas.h"
\r
32 #include "TParticle.h"
\r
33 #include "TNtuple.h"
\r
34 #include "TObjString.h"
\r
35 #include "TLorentzVector.h"
\r
36 #include "TDatabasePDG.h"
\r
38 #include "AliAnalysisTask.h"
\r
39 #include "AliAnalysisManager.h"
\r
40 #include "AliMCEventHandler.h"
\r
41 #include "AliMCEvent.h"
\r
44 #include "AliESDpid.h"
\r
45 #include "AliESDEvent.h"
\r
46 #include "AliGenEventHeader.h"
\r
47 #include "AliESDInputHandler.h"
\r
49 #include "AliStack.h"
\r
50 #include "AliHeader.h"
\r
52 #include "AliKFParticle.h"
\r
53 #include "AliKFVertex.h"
\r
54 #include "AliESDtrackCuts.h"
\r
55 #include "AliESDv0Cuts.h"
\r
56 #include "AliAnalysisTaskHdibaryonLPpi.h"
\r
61 #include "TCanvas.h"
\r
63 #include "TParallelCoord.h"
\r
65 #include "AliMCParticle.h"
\r
66 #include "AliGenPythiaEventHeader.h"
\r
69 #include "AliESDtrack.h"
\r
70 #include "AliCentrality.h"
\r
72 #include "THnSparse.h"
\r
74 #include "AliVertexerTracks.h"
\r
76 using namespace std;
\r
77 ClassImp(AliAnalysisTaskHdibaryonLPpi)
\r
78 //________________________________________________________________________
\r
79 AliAnalysisTaskHdibaryonLPpi::AliAnalysisTaskHdibaryonLPpi() : AliAnalysisTaskSE()/*AliAnalysisTask(name, ""), fMCEvent(0)*/, fESD(0), fESDtrackCutsV0(0),
\r
87 fHistMassLambdaPi(0),
\r
89 fHistMassLambdaPPi(0),
\r
90 fHistMassLambdaP(0),
\r
91 fHistMassLambdaK(0),
\r
95 fHistMassK0offlC(0),
\r
99 fHistArmenterosPodolanski(0),
\r
100 fHistArmenterosPodolanskiCut(0),
\r
101 fHistHDibaryonInvaMassGen(0),
\r
102 fHistHDibaryonInvaMassGenRes(0),
\r
103 fHistAntiHDibaryonInvaMassGen(0),
\r
104 fHistHDibaryonInvaMassAso(0),
\r
105 fHistHDibaryonInvaMassAsoReso(0),
\r
106 fHistAntiHDibaryonInvaMassAso(0),
\r
108 fHistHPointingAngle(0),
\r
110 fHistMassLambdaFromH(0),
\r
111 fHistMassLambdaFromHtLorentz(0),
\r
113 fHistMassPpiReso(0),
\r
116 fHistProtonPIDBb(0),
\r
118 fHistProtonPIDLambda(0),
\r
119 fHistPionPIDLambda(0),
\r
120 fHistMCdcaPvtxDvtx(0),
\r
121 fHistMCdcaPvtxLvtx(0),
\r
122 fHistMCdcaDvtxLvtx(0),
\r
124 fHistMCdecayAngle(0),
\r
125 fHistMCpointingAngle(0),
\r
127 fHistMCdcaPvtxDvtxReso(0),
\r
128 fHistMCdcaPvtxLvtxReso(0),
\r
129 fHistMCdcaDvtxLvtxReso(0),
\r
130 fHistMCangleLHReso(0),
\r
131 fHistMCdecayAngleReso(0),
\r
132 fHistMCpointingAngleReso(0),
\r
134 fHistCentrality(0),
\r
135 fHistCentralityAC(0),
\r
136 fHistMultiplicity(0),
\r
143 fHistPtvsEtaGen(0),
\r
144 fHistPtvsEtaAso(0),
\r
150 fHistTriggerStat(0),
\r
151 fHistTriggerStatAfterEventSelection(0),
\r
152 fHistMassHcentMult(0),
\r
155 // DefaultConstructor
\r
159 //________________________________________________________________________
\r
160 AliAnalysisTaskHdibaryonLPpi::AliAnalysisTaskHdibaryonLPpi(const char *name) : AliAnalysisTaskSE(name)/*AliAnalysisTask(name, ""), fMCEvent(0)*/, fESD(0), fESDtrackCutsV0(0),
\r
168 fHistMassLambdaPi(0),
\r
169 fHistMassLambda(0),
\r
170 fHistMassLambdaPPi(0),
\r
171 fHistMassLambdaP(0),
\r
172 fHistMassLambdaK(0),
\r
174 fHistMassK0offl(0),
\r
175 fHistMassK0onlC(0),
\r
176 fHistMassK0offlC(0),
\r
177 fHistMassPQonl(0),
\r
178 fHistMassPQoffl(0),
\r
180 fHistArmenterosPodolanski(0),
\r
181 fHistArmenterosPodolanskiCut(0),
\r
182 fHistHDibaryonInvaMassGen(0),
\r
183 fHistHDibaryonInvaMassGenRes(0),
\r
184 fHistAntiHDibaryonInvaMassGen(0),
\r
185 fHistHDibaryonInvaMassAso(0),
\r
186 fHistHDibaryonInvaMassAsoReso(0),
\r
187 fHistAntiHDibaryonInvaMassAso(0),
\r
189 fHistHPointingAngle(0),
\r
191 fHistMassLambdaFromH(0),
\r
192 fHistMassLambdaFromHtLorentz(0),
\r
194 fHistMassPpiReso(0),
\r
197 fHistProtonPIDBb(0),
\r
199 fHistProtonPIDLambda(0),
\r
200 fHistPionPIDLambda(0),
\r
201 fHistMCdcaPvtxDvtx(0),
\r
202 fHistMCdcaPvtxLvtx(0),
\r
203 fHistMCdcaDvtxLvtx(0),
\r
205 fHistMCdecayAngle(0),
\r
206 fHistMCpointingAngle(0),
\r
208 fHistMCdcaPvtxDvtxReso(0),
\r
209 fHistMCdcaPvtxLvtxReso(0),
\r
210 fHistMCdcaDvtxLvtxReso(0),
\r
211 fHistMCangleLHReso(0),
\r
212 fHistMCdecayAngleReso(0),
\r
213 fHistMCpointingAngleReso(0),
\r
215 fHistCentrality(0),
\r
216 fHistCentralityAC(0),
\r
217 fHistMultiplicity(0),
\r
224 fHistPtvsEtaGen(0),
\r
225 fHistPtvsEtaAso(0),
\r
231 fHistTriggerStat(0),
\r
232 fHistTriggerStatAfterEventSelection(0),
\r
233 fHistMassHcentMult(0),
\r
239 // Define input and output slots here
\r
240 // Input from a TChain
\r
241 DefineInput(0, TChain::Class());
\r
242 // Output to TList container
\r
243 DefineOutput(1, TList::Class()); //full
\r
246 if (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())SetHasMC();
\r
250 fESDtrackCutsV0 = new AliESDtrackCuts("AliESDtrackCutsV0","AliESDtrackCutsV0");
\r
251 fESDtrackCutsV0->SetAcceptKinkDaughters(kFALSE);
\r
252 fESDtrackCutsV0->SetMinNClustersTPC(80);
\r
253 fESDtrackCutsV0->SetMaxChi2PerClusterTPC(5);
\r
254 fESDtrackCutsV0->SetRequireTPCRefit(kTRUE);
\r
255 fESDtrackCutsV0->SetEtaRange(-0.9,0.9);
\r
256 fESDtrackCutsV0->SetPtRange(0.2,1.5);
\r
257 fESDtrackCutsV0->SetMinDCAToVertexXY(3);
\r
258 fESDtrackCutsV0->SetMinDCAToVertexZ(3);
\r
260 fESDCutsV0 = new AliESDv0Cuts("AliESDCutsV0","AliESDCutsV0");
\r
261 fESDCutsV0->SetMaxDcaV0Daughters(1.0);
\r
262 fESDCutsV0->SetMinDcaNegToVertex(2.0);
\r
263 fESDCutsV0->SetMinDcaPosToVertex(2.0);
\r
266 fEsdTrackCuts = new AliESDtrackCuts("AliESDtrackCuts");
\r
267 fEsdTrackCuts->SetMinNClustersTPC(80);
\r
268 fEsdTrackCuts->SetAcceptKinkDaughters(kFALSE);
\r
269 fEsdTrackCuts->SetMaxChi2PerClusterTPC(5);
\r
270 fEsdTrackCuts->SetRequireTPCRefit(kTRUE);
\r
271 fEsdTrackCuts->SetEtaRange(-0.9,0.9);
\r
274 //____________________________________________________________
\r
275 AliAnalysisTaskHdibaryonLPpi::~AliAnalysisTaskHdibaryonLPpi(){
\r
280 fHistList->Clear();
\r
283 if(fEsdTrackCuts) delete fEsdTrackCuts;
\r
284 if(fESDtrackCutsV0) delete fESDtrackCutsV0;
\r
285 if(fESDCutsV0) delete fESDCutsV0;
\r
289 //________________________________________________________________________
\r
290 void AliAnalysisTaskHdibaryonLPpi::UserCreateOutputObjects()
\r
292 // Create histograms
\r
295 fHistList = new TList();
\r
297 fHistMassDPi = new TH1F("fHistMassDPi", "Invariant mass distribution p+#pi^{-} ", 500, 1.0, 1.25);
\r
298 fHistMassDPi->GetXaxis()->SetTitle("Invariant mass p+#pi^{-} (GeV/c^{2})");
\r
299 fHistMassDPi->GetYaxis()->SetTitle("Entries");
\r
300 fHistMassDPi->SetMarkerStyle(kFullCircle);
\r
302 fHistMassLPi = new TH1F("fHistMassLPi", "Offline Invariant mass distribution p+#pi^{-} ", 500, 1.0, 1.25);
\r
303 fHistMassLPi->GetXaxis()->SetTitle("Invariant mass p+#pi^{-} (GeV/c^{2})");
\r
304 fHistMassLPi->GetYaxis()->SetTitle("Entries");
\r
305 fHistMassLPi->SetMarkerStyle(kFullCircle);
\r
307 fHistMassLambdaPi = new TH1F("fHistMassLambdaPi", "Invariant mass distribution #Lambda+#pi^{-} ", 500, 1.2, 1.5);
\r
308 fHistMassLambdaPi->GetXaxis()->SetTitle("Invariant mass #Lambda+#pi^{-} (GeV/c^{2})");
\r
309 fHistMassLambdaPi->GetYaxis()->SetTitle("Entries");
\r
310 fHistMassLambdaPi->SetMarkerStyle(kFullCircle);
\r
312 fHistMassLambda = new TH1F("fHistMassLambda", "Invariant mass distribution of #Lambda for further analyis", 500, 1.0, 1.2);
\r
313 fHistMassLambda->GetXaxis()->SetTitle("Invariant mass p+#pi^{+} (GeV/c^{2})");
\r
314 fHistMassLambda->GetYaxis()->SetTitle("Entries");
\r
315 fHistMassLambda->SetMarkerStyle(kFullCircle);
\r
317 fHistMassLambdaPPi = new TH1F("fHistMassLambdaPPi", "Invariant mass distribution #Lambdap#pi^{-} ", 300, 2.2, 2.5);
\r
318 fHistMassLambdaPPi->GetXaxis()->SetTitle("Invariant mass #Lambdap#pi^{-} (GeV/c^{2})");
\r
319 fHistMassLambdaPPi->GetYaxis()->SetTitle("Entries");
\r
320 fHistMassLambdaPPi->SetMarkerStyle(kFullCircle);
\r
322 fHistMassLambdaP = new TH1F("fHistMassLambdaP", "Invariant mass distribution #Lambdap ", 300, 2.2, 2.5);
\r
323 fHistMassLambdaP->GetXaxis()->SetTitle("Invariant mass #Lambdap (GeV/c^{2})");
\r
324 fHistMassLambdaP->GetYaxis()->SetTitle("Entries");
\r
325 fHistMassLambdaP->SetMarkerStyle(kFullCircle);
\r
327 fHistMassLambdaK = new TH1F("fHistMassLambdaK", "Invariant mass distribution #LambdaK^{-} ", 300, 1.6, 1.9);
\r
328 fHistMassLambdaK->GetXaxis()->SetTitle("Invariant mass #LambdaK^{-} (GeV/c^{2})");
\r
329 fHistMassLambdaK->GetYaxis()->SetTitle("Entries");
\r
330 fHistMassLambdaK->SetMarkerStyle(kFullCircle);
\r
332 fHistMassK0onl = new TH1F("fHistMassK0onl", "Invariant mass distribution K_{s}^{0} online V0 finder ", 400, 0.2, 1.0);
\r
333 fHistMassK0onl->GetXaxis()->SetTitle("Invariant mass #pi^{+}#pi^{-} (GeV/c^{2})");
\r
334 fHistMassK0onl->GetYaxis()->SetTitle("Entries");
\r
335 fHistMassK0onl->SetMarkerStyle(kFullCircle);
\r
337 fHistMassK0offl = new TH1F("fHistMassK0offl", "Invariant mass distribution K_{s}^{0} offline V0 finder ", 400, 0.2, 1.0);
\r
338 fHistMassK0offl->GetXaxis()->SetTitle("Invariant mass #pi^{+}#pi^{-} (GeV/c^{2})");
\r
339 fHistMassK0offl->GetYaxis()->SetTitle("Entries");
\r
340 fHistMassK0offl->SetMarkerStyle(kFullCircle);
\r
342 fHistMassK0onlC = new TH1F("fHistMassK0onlC", "Invariant mass distribution K_{s}^{0} online V0 finder ", 400, 0.2, 1.0);
\r
343 fHistMassK0onlC->GetXaxis()->SetTitle("Invariant mass #pi^{+}#pi^{-} (GeV/c^{2})");
\r
344 fHistMassK0onlC->GetYaxis()->SetTitle("Entries");
\r
345 fHistMassK0onlC->SetMarkerStyle(kFullCircle);
\r
347 fHistMassK0offlC = new TH1F("fHistMassK0offlC", "Invariant mass distribution K_{s}^{0} offline V0 finder ", 400, 0.2, 1.0);
\r
348 fHistMassK0offlC->GetXaxis()->SetTitle("Invariant mass #pi^{+}#pi^{-} (GeV/c^{2})");
\r
349 fHistMassK0offlC->GetYaxis()->SetTitle("Entries");
\r
350 fHistMassK0offlC->SetMarkerStyle(kFullCircle);
\r
352 fHistMassPQonl = new TH1F("fHistMassPQonl", "Invariant mass distribution K_{s}^{0}p using online V0 finder ", 500, 1.3, 2.3);
\r
353 fHistMassPQonl->GetXaxis()->SetTitle("Invariant mass K_{s}^{0}p (GeV/c^{2})");
\r
354 fHistMassPQonl->GetYaxis()->SetTitle("Entries");
\r
355 fHistMassPQonl->SetMarkerStyle(kFullCircle);
\r
357 fHistMassPQoffl = new TH1F("fHistMassPQoffl", "Invariant mass distribution K_{s}^{0}p using offline V0 finder ", 500, 1.3, 2.3);
\r
358 fHistMassPQoffl->GetXaxis()->SetTitle("Invariant mass K_{s}^{0}p (GeV/c^{2})");
\r
359 fHistMassPQoffl->GetYaxis()->SetTitle("Entries");
\r
360 fHistMassPQoffl->SetMarkerStyle(kFullCircle);
\r
362 fHistDC = new TH1F("fHistDC", "Proper decay length", 500, 0.0, 25);
\r
363 fHistDC->GetXaxis()->SetTitle("c#tau (cm)");
\r
364 fHistDC->GetYaxis()->SetTitle("Entries");
\r
365 fHistDC->SetMarkerStyle(kFullCircle);
\r
367 fHistArmenterosPodolanski = new TH2F("fHistArmenterosPodolanski", "Armenteros-Podolanski plot", 200,-1.0,1.0, 500,0,1);
\r
368 fHistArmenterosPodolanski->GetXaxis()->SetTitle("#alpha");
\r
369 fHistArmenterosPodolanski->GetYaxis()->SetTitle("q_{t}");
\r
370 fHistArmenterosPodolanski->SetMarkerStyle(kFullCircle);
\r
372 fHistArmenterosPodolanskiCut = new TH2F("fHistArmenterosPodolanskiCut", "Armenteros-Podolanski plot after cut", 200,-1.0,1.0, 500,0,1);
\r
373 fHistArmenterosPodolanskiCut->GetXaxis()->SetTitle("#alpha");
\r
374 fHistArmenterosPodolanskiCut->GetYaxis()->SetTitle("q_{t}");
\r
375 fHistArmenterosPodolanskiCut->SetMarkerStyle(kFullCircle);
\r
377 fHistHDibaryonInvaMassGen = new TH1F("fHistHDibaryonInvaMassGen", "Generated #Lambda p #pi^{-}", 200, 2.1, 2.3);
\r
378 fHistHDibaryonInvaMassGen->GetYaxis()->SetTitle("Counts");
\r
379 fHistHDibaryonInvaMassGen->GetXaxis()->SetTitle("Invariant mass (GeV/c^{2})");
\r
381 fHistHDibaryonInvaMassGenRes = new TH1F("fHistHDibaryonInvaMassGenRes", "Generated #Lambda p #pi^{-} with particles in rapidity!!!", 200, 2.1, 2.3);
\r
382 fHistHDibaryonInvaMassGenRes->GetYaxis()->SetTitle("Counts");
\r
383 fHistHDibaryonInvaMassGenRes->GetXaxis()->SetTitle("Invariant mass (GeV/c^{2})");
\r
385 fHistAntiHDibaryonInvaMassGen = new TH1F("fHistAntiHDibaryonInvaMassGen", "Generated #bar{#Lambda} #bar{p} #pi^{+}", 200, 2.1, 2.3);
\r
386 fHistAntiHDibaryonInvaMassGen->GetYaxis()->SetTitle("Counts");
\r
387 fHistAntiHDibaryonInvaMassGen->GetXaxis()->SetTitle("Invariant mass (GeV/c^{2})");
\r
389 fHistHDibaryonInvaMassAso = new TH1F("fHistHDibaryonInvaMassAso", "Associated #Lambda p #pi^{-}", 200, 2.1, 2.3);
\r
390 fHistHDibaryonInvaMassAso->GetYaxis()->SetTitle("Counts");
\r
391 fHistHDibaryonInvaMassAso->GetXaxis()->SetTitle("Invariant mass (GeV/c^{2})");
\r
393 fHistHDibaryonInvaMassAsoReso = new TH1F("fHistHDibaryonInvaMassAsoReso", "Associated #Lambda p #pi^{-}", 200, 2.1, 2.3);
\r
394 fHistHDibaryonInvaMassAsoReso->GetYaxis()->SetTitle("Counts");
\r
395 fHistHDibaryonInvaMassAsoReso->GetXaxis()->SetTitle("Invariant mass (GeV/c^{2})");
\r
397 fHistAntiHDibaryonInvaMassAso = new TH1F("fHistAntiHDibaryonInvaMassAso", "Associated #bar{#Lambda} #bar{p} #pi^{+}", 200, 2.1, 2.3);
\r
398 fHistAntiHDibaryonInvaMassAso->GetYaxis()->SetTitle("Counts");
\r
399 fHistAntiHDibaryonInvaMassAso->GetXaxis()->SetTitle("Invariant mass (GeV/c^{2})");
\r
401 fHistCheck = new TH2F("fHistCheck", "Check online/offline", 200, -0.5, 1.5, 200, -0.5, 1.5);
\r
402 fHistCheck->GetXaxis()->SetTitle("offline");
\r
403 fHistCheck->GetYaxis()->SetTitle("online");
\r
404 fHistCheck->SetMarkerStyle(kFullCircle);
\r
406 fHistHPointingAngle= new TH1F("fHistHPointingAngle", "Pointing angle distribution for #Lambdap#pi^{-}", 200, 0., 2*TMath::Pi());
\r
407 fHistHPointingAngle->GetXaxis()->SetTitle("Pointing angle distribution for #Lambdap#pi^{-}");
\r
408 fHistHPointingAngle->GetYaxis()->SetTitle("Entries");
\r
409 fHistHPointingAngle->SetMarkerStyle(kFullCircle);
\r
411 fHistMassH= new TH1F("fHistMassH", "Invariant mass distribution #Lambdap#pi^{-} (GeV/c^{2}) after pointing angle cut", 3000, 2.2, 2.5);
\r
412 fHistMassH->GetXaxis()->SetTitle("Invariant mass #Lambdap#pi^{-} (GeV/c^{2})");
\r
413 fHistMassH->GetYaxis()->SetTitle("Entries");
\r
414 fHistMassH->SetMarkerStyle(kFullCircle);
\r
416 fHistMassLambdaFromH= new TH1F("fHistMassLambdaFromH", "Invariant mass distribution #Lambda (GeV/c^{2}) asking for Mother to be H", 300, 1.0, 1.3);
\r
417 fHistMassLambdaFromH->GetXaxis()->SetTitle("Invariant mass #Lambda (GeV/c^{2})");
\r
418 fHistMassLambdaFromH->GetYaxis()->SetTitle("Entries");
\r
419 fHistMassLambdaFromH->SetMarkerStyle(kFullCircle);
\r
421 fHistMassLambdaFromHtLorentz= new TH1F("fHistMassLambdaFromHtLorentz", "Invariant mass distribution #Lambda (GeV/c^{2}) asking for Mother to be H", 300, 1.0, 1.3);
\r
422 fHistMassLambdaFromHtLorentz->GetXaxis()->SetTitle("Invariant mass #Lambda (GeV/c^{2})");
\r
423 fHistMassLambdaFromHtLorentz->GetYaxis()->SetTitle("Entries");
\r
424 fHistMassLambdaFromHtLorentz->SetMarkerStyle(kFullCircle);
\r
426 fHistMassPpi= new TH1F("fHistMassPpi", "Invariant mass distribution of the p#pi^{-} used for combing with Lambda", 300, 1.0, 1.6);
\r
427 fHistMassPpi->GetXaxis()->SetTitle("Invariant mass p#pi^{-} (GeV/c^{2})");
\r
428 fHistMassPpi->GetYaxis()->SetTitle("Entries");
\r
429 fHistMassPpi->SetMarkerStyle(kFullCircle);
\r
431 fHistMassPpiReso= new TH1F("fHistMassPpiReso", "Invariant mass distribution of the p#pi^{-} used for combing with Lambda", 300, 1.0, 1.6);
\r
432 fHistMassPpiReso->GetXaxis()->SetTitle("Invariant mass p#pi^{-} (GeV/c^{2})");
\r
433 fHistMassPpiReso->GetYaxis()->SetTitle("Entries");
\r
434 fHistMassPpiReso->SetMarkerStyle(kFullCircle);
\r
436 fHistMassLpi= new TH1F("fHistMassLpi", "Invariant mass distribution of the #Lambda#pi^{-} used for combing with p", 300, 1.1, 1.7);
\r
437 fHistMassLpi->GetXaxis()->SetTitle("Invariant mass #Lambda#pi^{-} (GeV/c^{2})");
\r
438 fHistMassLpi->GetYaxis()->SetTitle("Entries");
\r
439 fHistMassLpi->SetMarkerStyle(kFullCircle);
\r
441 fHistMassLP= new TH1F("fHistMassLP", "Invariant mass distribution of the #Lambda p used for combing with #pi^{-}", 300, 2.0, 2.3);
\r
442 fHistMassLP->GetXaxis()->SetTitle("Invariant mass #Lambda p (GeV/c^{2})");
\r
443 fHistMassLP->GetYaxis()->SetTitle("Entries");
\r
444 fHistMassLP->SetMarkerStyle(kFullCircle);
\r
446 fHistProtonPIDBb = new TH2F("fHistProtonPIDBb", "dE/dx after p PID", 100, 0., 10, 100, 0, 100);
\r
447 fHistProtonPIDBb->GetYaxis()->SetTitle("TPC Signal");
\r
448 fHistProtonPIDBb->GetXaxis()->SetTitle("P (GeV/c)");
\r
449 fHistProtonPIDBb->SetOption("scat");
\r
450 fHistProtonPIDBb->SetMarkerStyle(kFullCircle);
\r
452 fHistPionPIDBb = new TH2F("fHistPionPIDBb", "dE/dx after K PID", 100, 0., 10, 100, 0, 100);
\r
453 fHistPionPIDBb->GetYaxis()->SetTitle("TPC Signal");
\r
454 fHistPionPIDBb->GetXaxis()->SetTitle("P (GeV/c)");
\r
455 fHistPionPIDBb->SetOption("scat");
\r
456 fHistPionPIDBb->SetMarkerStyle(kFullCircle);
\r
458 fHistProtonPIDLambda = new TH2F("fHistProtonPIDLambda", "dE/dx after p PID from V0", 100, 0., 10, 100, 0, 100);
\r
459 fHistProtonPIDLambda->GetYaxis()->SetTitle("TPC Signal");
\r
460 fHistProtonPIDLambda->GetXaxis()->SetTitle("P (GeV/c)");
\r
461 fHistProtonPIDLambda->SetOption("scat");
\r
462 fHistProtonPIDLambda->SetMarkerStyle(kFullCircle);
\r
464 fHistPionPIDLambda = new TH2F("fHistPionPIDLambda", "dE/dx after #pi PID from V0", 100, 0, 10, 100, 0, 100);
\r
465 fHistPionPIDLambda->GetYaxis()->SetTitle("TPC Signal");
\r
466 fHistPionPIDLambda->GetXaxis()->SetTitle("P (GeV/c)");
\r
467 fHistPionPIDLambda->SetOption("scat");
\r
468 fHistPionPIDLambda->SetMarkerStyle(kFullCircle);
\r
470 fHistMCdcaPvtxDvtx= new TH1F("fHistMCdcaPvtxDvtx", "MC True DCA Primary Vertex - H Decay Vertex", 300, -0.1, 11.9);
\r
471 fHistMCdcaPvtxDvtx->GetXaxis()->SetTitle("dca prim. vtx- decay vtx (cm)");
\r
472 fHistMCdcaPvtxDvtx->GetYaxis()->SetTitle("Entries");
\r
473 fHistMCdcaPvtxDvtx->SetMarkerStyle(kFullCircle);
\r
475 fHistMCdcaPvtxLvtx= new TH1F("fHistMCdcaPvtxLvtx", "MC True DCA Primary Vertex - Lambda Decay Vertex", 300, -0.1, 11.9);
\r
476 fHistMCdcaPvtxLvtx->GetXaxis()->SetTitle("dca prim. vtx-#Lambda decay vtx (cm)");
\r
477 fHistMCdcaPvtxLvtx->GetYaxis()->SetTitle("Entries");
\r
478 fHistMCdcaPvtxLvtx->SetMarkerStyle(kFullCircle);
\r
480 fHistMCdcaDvtxLvtx= new TH1F("fHistMCdcaDvtxLvtx", "MC True DCA H Decay Vertex - Lambda Decay Vertex", 300, -0.1, 11.9);
\r
481 fHistMCdcaDvtxLvtx->GetXaxis()->SetTitle("dca H deacy vtx-#Lambda decay vtx (cm)");
\r
482 fHistMCdcaDvtxLvtx->GetYaxis()->SetTitle("Entries");
\r
483 fHistMCdcaDvtxLvtx->SetMarkerStyle(kFullCircle);
\r
485 fHistMCangleLH= new TH1F("fHistMCangleLH", "MC True Angle between #Lambda and H", 300, 0., 2*TMath::Pi());
\r
486 fHistMCangleLH->GetXaxis()->SetTitle("Angle (#Lambda-H)");
\r
487 fHistMCangleLH->GetYaxis()->SetTitle("Entries");
\r
488 fHistMCangleLH->SetMarkerStyle(kFullCircle);
\r
490 fHistMCdecayAngle= new TH1F("fHistMCdecayAngle", "MC True Angle between decay products", 300, 0., 2*TMath::Pi());
\r
491 fHistMCdecayAngle->GetXaxis()->SetTitle("Angle (#Lambda-p#pi)");
\r
492 fHistMCdecayAngle->GetYaxis()->SetTitle("Entries");
\r
493 fHistMCdecayAngle->SetMarkerStyle(kFullCircle);
\r
495 fHistMCpointingAngle= new TH1F("fHistMCpointingAngle", "MC True Pointing Angle", 3000, 0., 2*TMath::Pi());
\r
496 fHistMCpointingAngle->GetXaxis()->SetTitle("Pointing Angle");
\r
497 fHistMCpointingAngle->GetYaxis()->SetTitle("Entries");
\r
498 fHistMCpointingAngle->SetMarkerStyle(kFullCircle);
\r
500 fHistMCap = new TH2F("fHistMCap", "True MC Armenteros-Podolanski", 200,-1.0,1.0, 500,0,1);
\r
501 fHistMCap->GetYaxis()->SetTitle("#alpha");
\r
502 fHistMCap->GetXaxis()->SetTitle("q_{t}");
\r
503 fHistMCap->SetOption("scat");
\r
504 fHistMCap->SetMarkerStyle(kFullCircle);
\r
506 fHistMCdcaPvtxDvtxReso= new TH1F("fHistMCdcaPvtxDvtxReso", "MC True DCA Primary Vertex - H Decay Vertex", 300, -0.1, 11.9);
\r
507 fHistMCdcaPvtxDvtxReso->GetXaxis()->SetTitle("dca prim. vtx- decay vtx (cm)");
\r
508 fHistMCdcaPvtxDvtxReso->GetYaxis()->SetTitle("Entries");
\r
509 fHistMCdcaPvtxDvtxReso->SetMarkerStyle(kFullCircle);
\r
511 fHistMCdcaPvtxLvtxReso= new TH1F("fHistMCdcaPvtxLvtxReso", "MC True DCA Primary Vertex - Lambda Decay Vertex", 300, -0.1, 11.9);
\r
512 fHistMCdcaPvtxLvtxReso->GetXaxis()->SetTitle("dca prim. vtx-#Lambda decay vtx (cm)");
\r
513 fHistMCdcaPvtxLvtxReso->GetYaxis()->SetTitle("Entries");
\r
514 fHistMCdcaPvtxLvtxReso->SetMarkerStyle(kFullCircle);
\r
516 fHistMCdcaDvtxLvtxReso= new TH1F("fHistMCdcaDvtxLvtxReso", "MC True DCA H Decay Vertex - Lambda Decay Vertex", 300, -0.1, 11.9);
\r
517 fHistMCdcaDvtxLvtxReso->GetXaxis()->SetTitle("dca H deacy vtx-#Lambda decay vtx (cm)");
\r
518 fHistMCdcaDvtxLvtxReso->GetYaxis()->SetTitle("Entries");
\r
519 fHistMCdcaDvtxLvtxReso->SetMarkerStyle(kFullCircle);
\r
521 fHistMCangleLHReso= new TH1F("fHistMCangleLHReso", "MC True Angle between #Lambda and H", 300, 0., 2*TMath::Pi());
\r
522 fHistMCangleLHReso->GetXaxis()->SetTitle("Angle (#Lambda-H)");
\r
523 fHistMCangleLHReso->GetYaxis()->SetTitle("Entries");
\r
524 fHistMCangleLHReso->SetMarkerStyle(kFullCircle);
\r
526 fHistMCdecayAngleReso= new TH1F("fHistMCdecayAngleReso", "MC True Angle between decay products", 300, 0., 2*TMath::Pi());
\r
527 fHistMCdecayAngleReso->GetXaxis()->SetTitle("Angle (#Lambda-p#pi)");
\r
528 fHistMCdecayAngleReso->GetYaxis()->SetTitle("Entries");
\r
529 fHistMCdecayAngleReso->SetMarkerStyle(kFullCircle);
\r
531 fHistMCpointingAngleReso= new TH1F("fHistMCpointingAngleReso", "MC True Pointing Angle", 300, 0., 2*TMath::Pi());
\r
532 fHistMCpointingAngleReso->GetXaxis()->SetTitle("Pointing Angle");
\r
533 fHistMCpointingAngleReso->GetYaxis()->SetTitle("Entries");
\r
534 fHistMCpointingAngleReso->SetMarkerStyle(kFullCircle);
\r
536 fHistMCapReso = new TH2F("fHistMCapReso", "True MC Armenteros-Podolanski", 200,-1.0,1.0, 500,0,1);
\r
537 fHistMCapReso->GetYaxis()->SetTitle("#alpha");
\r
538 fHistMCapReso->GetXaxis()->SetTitle("q_{t}");
\r
539 fHistMCapReso->SetOption("scat");
\r
540 fHistMCapReso->SetMarkerStyle(kFullCircle);
\r
542 fHistCentrality = new TH1F("Centrality ", "Centrality", 11, -0.5, 10.5);
\r
543 fHistCentrality ->GetXaxis()->SetTitle("Centrality");
\r
544 fHistCentrality ->GetYaxis()->SetTitle("Entries");
\r
546 fHistCentralityAC = new TH1F("CentralityAC ", "CentralityAC", 11, -0.5, 10.5);
\r
547 fHistCentralityAC ->GetXaxis()->SetTitle("Centrality");
\r
548 fHistCentralityAC ->GetYaxis()->SetTitle("Entries");
\r
550 fHistMultiplicity = new TH1F("Multiplicity ", "Multiplicity", 100, 0, 20000);
\r
551 fHistMultiplicity ->GetXaxis()->SetTitle("Centrality");
\r
552 fHistMultiplicity ->GetYaxis()->SetTitle("Entries");
\r
554 fHistList->Add(fHistMassDPi);
\r
555 fHistList->Add(fHistMassLPi);
\r
556 fHistList->Add(fHistMassLambdaPi);
\r
557 fHistList->Add(fHistMassLambda);
\r
558 fHistList->Add(fHistMassLambdaPPi);
\r
559 fHistList->Add(fHistMassLambdaP);
\r
560 fHistList->Add(fHistMassLambdaK);
\r
561 fHistList->Add(fHistMassK0onl);
\r
562 fHistList->Add(fHistMassK0offl);
\r
563 fHistList->Add(fHistMassK0onlC);
\r
564 fHistList->Add(fHistMassK0offlC);
\r
565 fHistList->Add(fHistMassPQonl);
\r
566 fHistList->Add(fHistMassPQoffl);
\r
567 fHistList->Add(fHistDC);
\r
568 fHistList->Add(fHistArmenterosPodolanski);
\r
569 fHistList->Add(fHistArmenterosPodolanskiCut);
\r
570 fHistList->Add(fHistHDibaryonInvaMassGen);
\r
571 fHistList->Add(fHistHDibaryonInvaMassGenRes);
\r
572 fHistList->Add(fHistAntiHDibaryonInvaMassGen);
\r
573 fHistList->Add(fHistHDibaryonInvaMassAso);
\r
574 fHistList->Add(fHistHDibaryonInvaMassAsoReso);
\r
575 fHistList->Add(fHistAntiHDibaryonInvaMassAso);
\r
576 fHistList->Add(fHistCheck);
\r
577 fHistList->Add(fHistHPointingAngle);
\r
578 fHistList->Add(fHistMassH);
\r
579 fHistList->Add(fHistMassPpi);
\r
580 fHistList->Add(fHistMassPpiReso);
\r
581 fHistList->Add(fHistMassLpi);
\r
582 fHistList->Add(fHistMassLP);
\r
583 fHistList->Add(fHistMassLambdaFromH);
\r
584 fHistList->Add(fHistMassLambdaFromHtLorentz);
\r
585 fHistList->Add(fHistProtonPIDBb);
\r
586 fHistList->Add(fHistPionPIDBb);
\r
587 fHistList->Add(fHistProtonPIDLambda);
\r
588 fHistList->Add(fHistPionPIDLambda);
\r
589 fHistList->Add(fHistMCdcaPvtxDvtx);
\r
590 fHistList->Add(fHistMCdcaPvtxLvtx);
\r
591 fHistList->Add(fHistMCdcaDvtxLvtx);
\r
592 fHistList->Add(fHistMCangleLH);
\r
593 fHistList->Add(fHistMCdecayAngle);
\r
594 fHistList->Add(fHistMCpointingAngle);
\r
595 fHistList->Add(fHistMCap);
\r
596 fHistList->Add(fHistMCdcaPvtxDvtxReso);
\r
597 fHistList->Add(fHistMCdcaPvtxLvtxReso);
\r
598 fHistList->Add(fHistMCdcaDvtxLvtxReso);
\r
599 fHistList->Add(fHistMCangleLHReso);
\r
600 fHistList->Add(fHistMCdecayAngleReso);
\r
601 fHistList->Add(fHistMCpointingAngleReso);
\r
602 fHistList->Add(fHistMCapReso);
\r
603 fHistList->Add(fHistCentrality);
\r
604 fHistList->Add(fHistCentralityAC);
\r
605 fHistList->Add(fHistMultiplicity);
\r
607 fHistHilf1= new TH1F("fHistHilf1", "HD", 300, 0., 10);
\r
608 fHistHilf1->GetXaxis()->SetTitle("");
\r
609 fHistHilf1->GetYaxis()->SetTitle("Entries");
\r
611 fHistHilf2= new TH1F("fHistHilf2", "HD", 300, 0., 10);
\r
612 fHistHilf2->GetXaxis()->SetTitle("");
\r
613 fHistHilf2->GetYaxis()->SetTitle("Entries");
\r
615 fHistHilf3= new TH1F("fHistHilf3", "HD", 300, 0., 10);
\r
616 fHistHilf3->GetXaxis()->SetTitle("");
\r
617 fHistHilf3->GetYaxis()->SetTitle("Entries");
\r
619 fHistHilf4= new TH1F("fHistHilf4", "HD", 300, 0., 10);
\r
620 fHistHilf4->GetXaxis()->SetTitle("");
\r
621 fHistHilf4->GetYaxis()->SetTitle("Entries");
\r
623 fHistHilf5= new TH1F("fHistHilf5", "HD", 300, 0., 10);
\r
624 fHistHilf5->GetXaxis()->SetTitle("");
\r
625 fHistHilf5->GetYaxis()->SetTitle("Entries");
\r
627 fHistHilf6= new TH1F("fHistHilf6", "HD", 300, 0., 10);
\r
628 fHistHilf6->GetXaxis()->SetTitle("");
\r
629 fHistHilf6->GetYaxis()->SetTitle("Entries");
\r
631 fHistPtvsEtaGen = new TH2F("fHistPtvsEtaGen", "p_{t} vs #eta from generated H", 200,0.0,10.0, 200,-1,1);
\r
632 fHistPtvsEtaGen->GetXaxis()->SetTitle("p_{t} (GeV/c)");
\r
633 fHistPtvsEtaGen->GetYaxis()->SetTitle("#eta");
\r
634 fHistPtvsEtaGen->SetOption("scat");
\r
635 fHistPtvsEtaGen->SetMarkerStyle(kFullCircle);
\r
637 fHistPtvsEtaAso = new TH2F("fHistPtvsEtaAso", "p_{t} vs #eta from associated H", 200,0.0,10.0, 200,-1,1);
\r
638 fHistPtvsEtaAso->GetYaxis()->SetTitle("p_{t} (GeV/c)");
\r
639 fHistPtvsEtaAso->GetXaxis()->SetTitle("#eta");
\r
640 fHistPtvsEtaAso->SetOption("scat");
\r
641 fHistPtvsEtaAso->SetMarkerStyle(kFullCircle);
\r
643 fHistPtvsYGen = new TH2F("fHistPtvsYGen", "p_{t} vs rapidity from generated H", 200,0.0,10.0, 200,-1,1);
\r
644 fHistPtvsYGen->GetXaxis()->SetTitle("p_{t} (GeV/c)");
\r
645 fHistPtvsYGen->GetYaxis()->SetTitle("y");
\r
646 fHistPtvsYGen->SetOption("scat");
\r
647 fHistPtvsYGen->SetMarkerStyle(kFullCircle);
\r
649 fHistPtvsYAso = new TH2F("fHistPtvsYAso", "p_{t} vs rapidity from associated H", 200,0.0,10.0, 200,-1,1);
\r
650 fHistPtvsYAso->GetXaxis()->SetTitle("p_{t} (GeV/c)");
\r
651 fHistPtvsYAso->GetYaxis()->SetTitle("y");
\r
652 fHistPtvsYAso->SetOption("scat");
\r
653 fHistPtvsYAso->SetMarkerStyle(kFullCircle);
\r
655 fHistRap= new TH1F("fHistRap", "Rapidity", 400, -2., 2);
\r
656 fHistRap->GetXaxis()->SetTitle("Y");
\r
657 fHistRap->GetYaxis()->SetTitle("Entries");
\r
658 fHistPtvsEtaAso->SetMarkerStyle(kFullCircle);
\r
660 fHistList->Add(fHistHilf1);
\r
661 fHistList->Add(fHistHilf2);
\r
662 fHistList->Add(fHistHilf3);
\r
663 fHistList->Add(fHistHilf4);
\r
664 fHistList->Add(fHistHilf5);
\r
665 fHistList->Add(fHistHilf6);
\r
666 fHistList->Add(fHistPtvsEtaGen);
\r
667 fHistList->Add(fHistPtvsEtaAso);
\r
668 fHistList->Add(fHistPtvsYGen);
\r
669 fHistList->Add(fHistPtvsYAso);
\r
670 fHistList->Add(fHistRap);
\r
672 fHistCount = new TH1F("fHistCount","test",17,0,17);
\r
673 fHistCount->GetXaxis()->SetBinLabel(1,"Events");
\r
674 fHistCount->GetXaxis()->SetBinLabel(2,"MC All");
\r
675 fHistCount->GetXaxis()->SetBinLabel(3,"MC from Primary Vtx");
\r
676 fHistCount->GetXaxis()->SetBinLabel(4,"Horst");
\r
677 fHistCount->GetXaxis()->SetBinLabel(5,"Lambda Candidates");
\r
678 fHistCount->GetXaxis()->SetBinLabel(6,"Sigma Candidates");
\r
679 fHistCount->GetXaxis()->SetBinLabel(7,"Horst");
\r
680 fHistCount->GetXaxis()->SetBinLabel(8,"Horst");
\r
681 fHistCount->GetXaxis()->SetBinLabel(9,"Horst");
\r
682 fHistCount->GetXaxis()->SetBinLabel(10,"MC All #bar{Lambda}(1520)s");
\r
683 fHistCount->GetXaxis()->SetBinLabel(11,"");
\r
684 fHistCount->GetXaxis()->SetBinLabel(12,"H-Dibaryon");
\r
685 fHistCount->GetXaxis()->SetBinLabel(13,"Hypertriton 2-Body");
\r
686 fHistCount->GetXaxis()->SetBinLabel(14,"Hypertriton 3-Body");
\r
687 fHistCount->GetXaxis()->SetBinLabel(15,"");
\r
688 fHistCount->GetXaxis()->SetBinLabel(16,"");
\r
689 fHistCount->GetXaxis()->SetBinLabel(17,"Lambdas!!!");
\r
690 fHistCount->SetStats(0);
\r
691 fHistCount->SetFillColor(38);
\r
692 fHistList->Add(fHistCount);
\r
694 //trigger statistics histogram
\r
695 fHistTriggerStat = new TH1F("fHistTriggerStat","Trigger statistics", 4,-0.5, 3.5);
\r
696 const Char_t* aTriggerNames[] = { "kMB", "kCentral", "kSemiCentral" };
\r
697 for ( Int_t ii=0; ii < 3; ii++ )
\r
698 fHistTriggerStat->GetXaxis()->SetBinLabel(ii+1, aTriggerNames[ii]);
\r
700 fHistTriggerStatAfterEventSelection = new TH1F("fHistTriggerStatAfterEventSelection","Trigger statistics after event selection", 4,-0.5, 3.5);
\r
701 for ( Int_t ii=0; ii < 3; ii++ )
\r
702 fHistTriggerStatAfterEventSelection->GetXaxis()->SetBinLabel(ii+1, aTriggerNames[ii]);
\r
703 fHistList->Add(fHistTriggerStat);
\r
704 fHistList->Add(fHistTriggerStatAfterEventSelection);
\r
706 fHistMassHcentMult = new TH3F("fHistMassHcentMult", "Inv. Mass vs. centrality vs. multiplicity", 100, 2.2, 2.3, 5, 0, 4, 300, 0, 6000);
\r
707 fHistMassHcentMult->GetXaxis()->SetTitle("Invariant mass #Lambdap#pi^{-} (GeV/c^{2})"); // inv. mass
\r
708 fHistMassHcentMult->GetYaxis()->SetTitle("Centrality"); // triggertype
\r
709 fHistMassHcentMult->GetZaxis()->SetTitle("Multiplicity"); // refTPC
\r
710 fHistList->Add(fHistMassHcentMult);
\r
713 const Double_t kz = 2*TMath::Pi();
\r
714 Int_t binsD01[16]={ 300, 200, 100, 100, 100, 100, 100, 100, 200, 200, 200, 200, 400, 200, 200, 3};
\r
715 Double_t xminD01[16]={2.0, 1.0, 0., -1, 0., 0., 0., 0., 0., 0., 0., 0., -1, 0., 0., 0};
\r
716 Double_t xmaxD01[16]={2.3, 1.2, kz, 1, 1, 10, 10, 5, 5, 5, 5, 100, 1, 100, 4000, 1};
\r
719 fHistNdim = new THnSparseF("fHistNdim","THnS;InvMass, InvMassLambda, pointingAngle, armPoAlpha, armPoQt, pTL, pTH, d0p, d0n, dcaHd, dca, decayL, cosPA, centr, multi, mcinf;InvMassH", 16,binsD01,xminD01,xmaxD01);
\r
720 fHistList->Add(fHistNdim);
\r
723 //________________________________________________________________________
\r
724 void AliAnalysisTaskHdibaryonLPpi::UserExec(Option_t *)
\r
727 // Called for each event
\r
729 //define improtant masses
\r
730 Double_t pionMass = 0.13957;
\r
731 Double_t protonMass = 0.93827;
\r
734 Long_t pdgPionPlus = 211;
\r
735 Long_t pdgPionMinus = -211;
\r
736 Long_t pdgProton = 2212;
\r
737 Long_t pdgAntiProton = -2212;
\r
738 Long_t pdgLambda = 3122;
\r
739 Long_t pdgAntiLambda = -3122;
\r
740 Long_t pdgHDibaryon = 1020000020;
\r
741 Long_t pdgAntiHDibaryon = -1020000020;
\r
743 AliStack* stack(NULL);
\r
746 if(!fMCEvent)return;
\r
748 AliHeader *head = fMCEvent->Header();
\r
750 AliGenPythiaEventHeader *pyheader = (AliGenPythiaEventHeader*)head->GenEventHeader();
\r
751 if(!pyheader)return;
\r
753 if(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler()){
\r
754 if(static_cast<AliMCEventHandler*>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())->MCEvent()){
\r
755 if(!(static_cast<AliMCEventHandler*>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())->TreeK()))return;
\r
756 if(!(static_cast<AliMCEventHandler*>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())->TreeTR()))return;
\r
757 if(!(static_cast<AliMCEventHandler*>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())->Init("local")))return;
\r
758 stack = static_cast<AliMCEventHandler*>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())->MCEvent()->Stack();
\r
759 if(!(static_cast<AliMCEventHandler*>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())->MCEvent()->Stack()->TreeK()))return;
\r
766 // -------------------------------------------------------
\r
767 // Loop for Inv. Mass via ESD tracks
\r
768 // -------------------------------------------------------
\r
770 fHistCount->Fill(0);
\r
772 fESD = dynamic_cast<AliESDEvent *>(fInputEvent);
\r
774 const AliESDVertex *vertex = fESD->GetPrimaryVertexTracks();
\r
775 if (vertex->GetNContributors()<1)
\r
778 vertex = fESD->GetPrimaryVertexSPD();
\r
779 if(vertex->GetNContributors()<1) {
\r
783 if (TMath::Abs(vertex->GetZv()) > 10) return;
\r
785 Int_t centrality = -5;
\r
786 Double_t centrPerc = -5;
\r
788 if (fESD->GetEventSpecie() == 4)
\r
790 AliCentrality *esdCentrality = fESD->GetCentrality();
\r
791 centrality = esdCentrality->GetCentralityClass10("V0M"); // centrality percentile determined with V0
\r
792 centrPerc = esdCentrality->GetCentralityPercentile("V0M");
\r
793 if (centrality < 0. || centrality > 8.) return; //0 bis 80 %
\r
794 // cout<<"Centrality: "<< centrality << endl;
\r
797 fHistCentrality->Fill(centrality);
\r
799 //*****************//
\r
801 //*****************//
\r
803 // Float_t percentile=centrality->GetCentralityPercentile("V0M");
\r
805 Bool_t isSelectedCentral = (((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected() & AliVEvent::kCentral);
\r
806 Bool_t isSelectedSemiCentral = (((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected() & AliVEvent::kSemiCentral);
\r
807 Bool_t isSelectedMB = (((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected() & AliVEvent::kMB);
\r
809 Int_t triggertype = 17;
\r
811 if(isSelectedCentral){
\r
812 fHistTriggerStat->Fill(1);
\r
816 if(isSelectedSemiCentral){
\r
817 fHistTriggerStat->Fill(2);
\r
822 fHistTriggerStat->Fill(0);
\r
826 // if(isSelectedCentral || isSelectedSemiCentral || isSelectedMB){
\r
828 //*******************************
\r
830 Int_t runNumber = 0;
\r
832 runNumber = fESD->GetRunNumber();
\r
834 if (!fPIDtpcESD) fPIDtpcESD = ((AliESDInputHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->GetESDpid();
\r
836 fPIDtpcESD = new AliESDpid(); // HACK FOR MC PBPB --> PLEASE REMOVE AS SOON AS POSSIBLE
\r
837 fPIDtpcESD->GetTPCResponse().SetBetheBlochParameters(1.28778e+00/50., 3.13539e+01, TMath::Exp(-3.16327e+01), 1.87901e+00, 6.41583e+00);
\r
843 TObjArray* listCrossV0 = fESDCutsV0->GetAcceptedV0s(fESD);
\r
844 Int_t nGoodV0s = listCrossV0->GetEntries();
\r
846 const AliESDVertex *esdVer = fESD->GetPrimaryVertex();
\r
847 AliESDVertex *esdVer1 = new AliESDVertex(*esdVer);
\r
849 AliVertexerTracks *vertexer = new AliVertexerTracks(fESD->GetMagneticField());
\r
850 TObjArray *trkArray = new TObjArray(2);
\r
851 AliVertexerTracks *vertexer1 = new AliVertexerTracks(fESD->GetMagneticField());
\r
852 TObjArray *trkArray1 = new TObjArray(2);
\r
854 AliKFParticle::SetField(fESD->GetMagneticField());
\r
856 AliKFVertex primVtx(*(fESD->GetPrimaryVertex()));
\r
858 Int_t refMultTpc = AliESDtrackCuts::GetReferenceMultiplicity(fESD, kTRUE);
\r
859 //cout<<"Multiplicity: "<< refMultTpc << endl;
\r
860 fHistMultiplicity->Fill(refMultTpc);
\r
862 Double_t mn[3] = {0,0,0};
\r
863 Double_t mp[3] = {0,0,0};
\r
864 Double_t mm[3] = {0,0,0};
\r
865 Double_t dd[3] = {0,0,0};
\r
866 Double_t dd1[3] = {0,0,0};
\r
867 const Double_t cProtonMass=TDatabasePDG::Instance()->GetParticle(2212)->Mass();
\r
868 const Double_t cPionMass=TDatabasePDG::Instance()->GetParticle(211)->Mass();
\r
869 const Double_t cElectronMass=TDatabasePDG::Instance()->GetParticle(11)->Mass();
\r
870 const Double_t cLambdaMass=TDatabasePDG::Instance()->GetParticle(3122)->Mass();
\r
871 Double_t decayLength=0;
\r
872 Double_t decayLengthH=0;
\r
874 //V0 Loop for Lambda and Anti-Lambda
\r
875 for(Int_t iV0MI = 0; iV0MI < nGoodV0s ; iV0MI++) {
\r
876 AliESDv0 * fV0MIs = fESD->GetV0(iV0MI);
\r
877 Int_t lOnFlyStatus = 0;
\r
879 lOnFlyStatus = fV0MIs->GetOnFlyStatus();
\r
880 Double_t lInvMassLambda=0;
\r
881 Double_t lInvMassLambdaPi=0;
\r
882 Double_t lPtLambda=0;
\r
883 Double_t lPzLambda=0;
\r
884 Double_t lPLambda=0;
\r
888 TLorentzVector posE;
\r
889 TLorentzVector negE;
\r
890 TLorentzVector photon;
\r
896 if (!lOnFlyStatus){
\r
901 // fHistMultiplicity->Fill(refMultTpc);
\r
902 fHistCentralityAC->Fill(centrality);
\r
903 fHistCheck->Fill(offl,onl);
\r
905 AliESDtrack* trackPosTest = fESD->GetTrack(fV0MIs->GetPindex());
\r
906 AliESDtrack* trackNegTest = fESD->GetTrack(fV0MIs->GetNindex());
\r
909 if (!fEsdTrackCuts->AcceptTrack(trackPosTest)) continue;
\r
910 if (!fESDtrackCutsV0->AcceptTrack(trackNegTest)) continue;
\r
913 //PID via specific energy loss in the TPC
\r
914 //define the arrays for the Bethe-Bloch-Parameters
\r
915 Double_t parProton[5] = {0,0,0,0,0};
\r
917 if(runNumber < 166500) //LHC10h
\r
919 parProton[0] = 1.45802; // ALEPH parameters for protons (pass2)
\r
920 parProton[1] = 27.4992;
\r
921 parProton[2] = 4.00313e-15;
\r
922 parProton[3] = 2.48485;
\r
923 parProton[4] = 8.31768;
\r
926 if(runNumber > 166500) //LHC11h
\r
928 parProton[0] = 1.11243; // ALEPH parameters for protons (pass2)
\r
929 parProton[1] = 26.1144;
\r
930 parProton[2] = 4.00313e-15;
\r
931 parProton[3] = 2.72969;
\r
932 parProton[4] = 9.15038;
\r
935 //Get the total momentum for each track at the inner readout of the TPC
\r
936 Double_t ptotN = trackNegTest->GetInnerParam()->GetP();
\r
937 Double_t ptotP = trackPosTest->GetInnerParam()->GetP();
\r
939 //define expected signals for the various species
\r
940 Double_t expSignalPionP = 0;
\r
941 Double_t expSignalPionN = 0;
\r
942 Double_t expSignalProtonN = 0;
\r
943 Double_t expSignalProtonP = 0;
\r
948 expSignalProtonN = AliExternalTrackParam::BetheBlochAleph(ptotN/(protonMass),parProton[0],parProton[1],parProton[2],parProton[3],parProton[4]);
\r
949 expSignalProtonP = AliExternalTrackParam::BetheBlochAleph(ptotP/(protonMass),parProton[0],parProton[1],parProton[2],parProton[3],parProton[4]);
\r
955 expSignalPionP = 0.7*AliExternalTrackParam::BetheBlochAleph(ptotP/(pionMass),parProton[0],parProton[1],parProton[2],parProton[3],parProton[4]);
\r
956 expSignalPionN = 0.7*AliExternalTrackParam::BetheBlochAleph(ptotN/(pionMass),parProton[0],parProton[1],parProton[2],parProton[3],parProton[4]);
\r
958 expSignalProtonN = 0.65*AliExternalTrackParam::BetheBlochAleph(ptotN/(protonMass),parProton[0],parProton[1],parProton[2],parProton[3],parProton[4]);
\r
959 expSignalProtonP = 0.65*AliExternalTrackParam::BetheBlochAleph(ptotP/(protonMass),parProton[0],parProton[1],parProton[2],parProton[3],parProton[4]);
\r
962 // PID cut on the nuclei (proton, deuteron, triton, helium3)
\r
963 Bool_t corrParticle = kFALSE;
\r
965 Bool_t posProton = kFALSE;
\r
969 if (TMath::Abs(fPIDtpcESD->NumberOfSigmasTPC(trackPosTest, AliPID::kProton)) < 3)
\r
972 corrParticle = kTRUE;
\r
978 if(//trackPosTest->GetTPCsignal() < 1200 &&
\r
979 TMath::Abs(trackPosTest->GetTPCsignal() - expSignalProtonP)/expSignalProtonP < 0.4)
\r
982 corrParticle = kTRUE;
\r
986 Bool_t negProton = kFALSE;
\r
990 if (TMath::Abs(fPIDtpcESD->NumberOfSigmasTPC(trackNegTest, AliPID::kProton)) < 3)
\r
993 corrParticle = kTRUE;
\r
999 if(//trackNegTest->GetTPCsignal() < 1200 &&
\r
1000 TMath::Abs(trackNegTest->GetTPCsignal() - expSignalProtonN)/expSignalProtonN < 0.4)
\r
1002 negProton = kTRUE;
\r
1003 corrParticle = kTRUE;
\r
1007 //PID cut for pions
\r
1008 //data: 3sigma cut on the pions
\r
1010 Bool_t negPion = kFALSE;
\r
1011 Bool_t posPion = kFALSE;
\r
1015 if (TMath::Abs(fPIDtpcESD->NumberOfSigmasTPC(trackPosTest, AliPID::kPion)) < 4) posPion=kTRUE; //pos daughter has to be a pion
\r
1016 if (TMath::Abs(fPIDtpcESD->NumberOfSigmasTPC(trackNegTest, AliPID::kPion)) < 4) negPion=kTRUE; // negative daughter has to be a pion
\r
1019 //MC: like the nuclei via the specific energyloss in the TPC
\r
1022 if (TMath::Abs(trackPosTest->GetTPCsignal() - expSignalPionP)/expSignalPionP < 0.4) posPion=kTRUE;
\r
1023 if (TMath::Abs(trackNegTest->GetTPCsignal() - expSignalPionN)/expSignalPionN < 0.4) negPion=kTRUE;
\r
1026 if (!(posProton==kTRUE)) continue;
\r
1027 if (!(negPion==kTRUE)) continue;
\r
1031 if( !(trackPosTest->GetStatus() & AliESDtrack::kTPCrefit)){
\r
1035 if( !(trackNegTest->GetStatus() & AliESDtrack::kTPCrefit)){
\r
1039 if( trackPosTest->GetSign() >0 && trackNegTest->GetSign() <0){
\r
1040 fV0MIs->GetNPxPyPz(mn[0],mn[1],mn[2]); //reconstructed cartesian momentum components of negative daughter
\r
1041 fV0MIs->GetPPxPyPz(mp[0],mp[1],mp[2]); //reconstructed cartesian momentum components of positive daughter
\r
1044 if( trackPosTest->GetSign() <0 && trackNegTest->GetSign() >0){
\r
1045 fV0MIs->GetPPxPyPz(mn[0],mn[1],mn[2]); //reconstructed cartesian momentum components of negative daughter
\r
1046 fV0MIs->GetNPxPyPz(mp[0],mp[1],mp[2]); //reconstructed cartesian momentum components of positive daughter
\r
1049 fV0MIs->GetPxPyPz(mm[0],mm[1],mm[2]); //reconstructed cartesian momentum components of mother
\r
1051 TVector3 vecN(mn[0],mn[1],mn[2]);
\r
1052 TVector3 vecP(mp[0],mp[1],mp[2]);
\r
1053 TVector3 vecM(mm[0],mm[1],mm[2]);
\r
1055 Double_t thetaP = acos((vecP * vecM)/(vecP.Mag() * vecM.Mag()));
\r
1056 Double_t thetaN = acos((vecN * vecM)/(vecN.Mag() * vecM.Mag()));
\r
1058 Double_t alfa = ((vecP.Mag())*cos(thetaP)-(vecN.Mag())*cos(thetaN))/
\r
1059 ((vecP.Mag())*cos(thetaP)+(vecN.Mag())*cos(thetaN)) ;
\r
1060 Double_t qt = vecP.Mag()*sin(thetaP);
\r
1062 fHistArmenterosPodolanski->Fill(alfa,qt); //Armenteros-Podolanski calculation
\r
1064 TLorentzVector k0;
\r
1065 TLorentzVector k0daugh1;
\r
1066 TLorentzVector k0daugh2;
\r
1067 TLorentzVector proton;
\r
1068 TLorentzVector pq;
\r
1070 k0daugh1.SetXYZM(mn[0],mn[1],mn[2],cPionMass);
\r
1071 k0daugh2.SetXYZM(mp[0],mp[1],mp[2],cPionMass);
\r
1072 k0=k0daugh1+k0daugh2;
\r
1074 fV0MIs->ChangeMassHypothesis(3122);
\r
1075 lInvMassLambda = fV0MIs->GetEffMass();
\r
1076 lPtLambda = fV0MIs->Pt();
\r
1077 lPzLambda = fV0MIs->Pz();
\r
1078 lPLambda = fV0MIs->P();
\r
1080 trkArray->AddAt(trackPosTest,0);
\r
1081 trkArray->AddAt(trackNegTest,1);
\r
1083 vertexer->SetVtxStart(esdVer1);
\r
1084 AliESDVertex *decayVertex = (AliESDVertex*)vertexer->VertexForSelectedESDTracks(trkArray);
\r
1086 dd[0]=fESD->GetPrimaryVertexSPD()->GetX()-decayVertex->GetX();
\r
1087 dd[1]=fESD->GetPrimaryVertexSPD()->GetY()-decayVertex->GetY();
\r
1088 dd[2]=fESD->GetPrimaryVertexSPD()->GetZ()-decayVertex->GetZ();
\r
1090 decayLength=sqrt(dd[0]*dd[0]+dd[1]*dd[1]+dd[2]*dd[2]);
\r
1092 if (decayVertex) delete decayVertex;
\r
1094 TLorentzVector negPio1;
\r
1095 TLorentzVector posProt1;
\r
1096 TLorentzVector posP;
\r
1097 TLorentzVector posProt;
\r
1098 TLorentzVector negK;
\r
1099 TLorentzVector negPio;
\r
1100 TLorentzVector negPi;
\r
1101 TLorentzVector omega;
\r
1102 TLorentzVector threeSum;
\r
1103 TLorentzVector fourSum;
\r
1104 TLorentzVector ppK;
\r
1105 TLorentzVector posPiK;
\r
1106 TLorentzVector negPiK;
\r
1107 TLorentzVector kaon;
\r
1108 TLorentzVector lambda;
\r
1109 TLorentzVector lambdaH;
\r
1110 TLorentzVector hDibaryon;
\r
1116 h.SetXYZ(-dd[0],-dd[1],-dd[2]);
\r
1118 if (onl==1)fHistMassDPi->Fill(lInvMassLambda);
\r
1120 if (offl==1||onl==1){
\r
1121 fHistMassLPi->Fill(lInvMassLambda);
\r
1123 negE.SetXYZM(mn[0],mn[1],mn[2],cElectronMass);
\r
1124 posE.SetXYZM(mp[0],mp[1],mp[2],cElectronMass);
\r
1127 negPiK.SetXYZM(mn[0],mn[1],mn[2],cPionMass);
\r
1128 posPiK.SetXYZM(mp[0],mp[1],mp[2],cPionMass);
\r
1129 kaon=posPiK+negPiK;
\r
1131 negPi.SetXYZM(mn[0],mn[1],mn[2],cPionMass);
\r
1132 posP.SetXYZM(mp[0],mp[1],mp[2],cProtonMass);
\r
1133 lambda=negPi+posP;
\r
1134 lambdaH.SetXYZM(mm[0],mm[1],mm[2],cLambdaMass);
\r
1136 if (lInvMassLambda>1.1113&&lInvMassLambda<1.1202){
\r
1139 if (qt<-2.21*alfa*alfa+2.945*alfa-0.887) continue;
\r
1140 if (qt>-2.21*alfa*alfa+2.945*alfa-0.873) continue;
\r
1141 if (photon.M()<0.005) continue;
\r
1142 if (kaon.M()>0.495 && kaon.M()<0.500 ) continue;
\r
1146 fHistMassLambda->Fill(lInvMassLambda);
\r
1148 Bool_t isCorrectlyAssociatedLambda = kFALSE;
\r
1149 Bool_t isPartialCorrectlyAssociatedLambda = kFALSE;
\r
1150 Int_t labelAssociatedH=-1;
\r
1151 Int_t labelLambda=-1;
\r
1154 Int_t labelPosTest = trackPosTest->GetLabel();
\r
1155 TParticle *tparticleDaughter = stack->Particle(TMath::Abs(labelPosTest));
\r
1156 Int_t labelMother = tparticleDaughter->GetFirstMother();
\r
1157 TParticle *tparticleMother = stack->Particle(TMath::Abs(labelMother));
\r
1159 Int_t labelOma = tparticleMother->GetFirstMother();
\r
1160 TParticle *tparticleOma = stack->Particle(TMath::Abs(labelOma));
\r
1162 if ((tparticleOma->GetPdgCode() < 0) && TMath::Abs(tparticleMother->GetPdgCode())==pdgLambda){// check mother to be Lambda
\r
1163 Int_t labelFirstDaughter = tparticleMother->GetDaughter(1);// Proton
\r
1164 Int_t labelThirdDaughter = tparticleMother->GetDaughter(0);// Pion
\r
1166 TParticle *tparticleFirstDaughter = stack->Particle(TMath::Abs(labelFirstDaughter));
\r
1167 TParticle *tparticleThirdDaughter = stack->Particle(TMath::Abs(labelThirdDaughter));
\r
1169 if((tparticleFirstDaughter->GetPdgCode() == pdgProton && tparticleThirdDaughter->GetPdgCode()== pdgPionMinus) ||
\r
1170 (tparticleFirstDaughter->GetPdgCode() == pdgPionMinus && tparticleThirdDaughter->GetPdgCode()== pdgProton)){ //daughter PDGs
\r
1171 isPartialCorrectlyAssociatedLambda = kTRUE;
\r
1172 labelLambda=labelMother;
\r
1177 if(tparticleOma->GetPdgCode() == pdgHDibaryon){ //check grandmother to be H PDG
\r
1178 if (TMath::Abs(tparticleMother->GetPdgCode())==pdgLambda){// check mother to be Lambda
\r
1179 Int_t labelFirstDaughter = tparticleMother->GetDaughter(1);// Proton
\r
1180 Int_t labelThirdDaughter = tparticleMother->GetDaughter(0);// Pion
\r
1182 TParticle *tparticleFirstDaughter = stack->Particle(TMath::Abs(labelFirstDaughter));
\r
1183 TParticle *tparticleThirdDaughter = stack->Particle(TMath::Abs(labelThirdDaughter));
\r
1185 if((tparticleFirstDaughter->GetPdgCode() == pdgProton && tparticleThirdDaughter->GetPdgCode()== pdgPionMinus) ||
\r
1186 (tparticleFirstDaughter->GetPdgCode() == pdgPionMinus && tparticleThirdDaughter->GetPdgCode()== pdgProton)){ //daughter PDGs
\r
1187 isCorrectlyAssociatedLambda = kTRUE;
\r
1188 labelAssociatedH=labelOma;
\r
1189 fHistMassLambdaFromH->Fill(lInvMassLambda);
\r
1190 fHistMassLambdaFromHtLorentz->Fill(lambda.M());
\r
1196 fHistProtonPIDLambda->Fill(trackPosTest->GetInnerParam()->GetP(), trackPosTest->GetTPCsignal());
\r
1197 fHistPionPIDLambda->Fill(trackNegTest->GetInnerParam()->GetP(), trackNegTest->GetTPCsignal());
\r
1199 //---------------------------------------------------------
\r
1200 // Proton track loop
\r
1201 //---------------------------------------------------------
\r
1202 fHistArmenterosPodolanskiCut->Fill(alfa,qt);
\r
1204 for (Int_t iTracksP = 0; iTracksP < fESD->GetNumberOfTracks(); iTracksP++) {
\r
1205 AliESDtrack* trackP = dynamic_cast<AliESDtrack*> (fESD->GetTrack(iTracksP));
\r
1206 if (trackP->GetSign()<0) continue;
\r
1208 if (iTracksP==fV0MIs->GetPindex())continue;
\r
1209 if (iTracksP==fV0MIs->GetNindex())continue;
\r
1211 if (!fEsdTrackCuts->AcceptTrack(trackP)) continue;
\r
1213 AliKFParticle protonKF( *(trackP), 2212);
\r
1215 if (!trackP->GetInnerParam()) continue;
\r
1222 if (TMath::Abs(fPIDtpcESD->NumberOfSigmasTPC(trackP, AliPID::kProton)) > 3) continue;
\r
1225 fHistProtonPIDBb->Fill(trackP->GetInnerParam()->GetP(), trackP->GetTPCsignal());
\r
1227 posProt.SetXYZM(trackP->Px(),trackP->Py(),trackP->Pz(),cProtonMass);
\r
1229 //!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
\r
1230 //Pion Track loop!!!!!!!!!!!!!!!!!!!!!
\r
1231 //!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
\r
1232 for (Int_t iTracksN = iTracksP+1; iTracksN < fESD->GetNumberOfTracks(); iTracksN++) {
\r
1233 AliESDtrack* trackN = dynamic_cast<AliESDtrack*> (fESD->GetTrack(iTracksN));
\r
1235 if (iTracksN==fV0MIs->GetPindex())continue;
\r
1236 if (iTracksN==fV0MIs->GetNindex())continue;
\r
1237 if (trackN->GetSign()>0) continue;
\r
1239 if (!fEsdTrackCuts->AcceptTrack(trackN)) continue;
\r
1240 if (!fESDtrackCutsV0->AcceptTrack(trackN)) continue;
\r
1242 negPi.SetXYZM(mn[0],mn[1],mn[2],cPionMass);
\r
1243 posP.SetXYZM(mp[0],mp[1],mp[2],cProtonMass);
\r
1244 negPio.SetXYZM(trackN->Px(),trackN->Py(),trackN->Pz(),cPionMass);
\r
1246 threeSum=negPi+posP+negPio;
\r
1247 lInvMassLambdaPi=threeSum.M();
\r
1249 fHistDC->Fill(decayLength*lInvMassLambda/lPLambda);
\r
1251 AliKFParticle posPionKF( *(trackN) ,-211);
\r
1253 if (!trackN->GetInnerParam()) continue;
\r
1259 if (TMath::Abs(fPIDtpcESD->NumberOfSigmasTPC(trackN, AliPID::kPion)) > 3) continue;
\r
1261 fHistPionPIDBb->Fill(trackN->GetInnerParam()->GetP(), trackN->GetTPCsignal());
\r
1263 trkArray1->AddAt(trackP,0);
\r
1264 trkArray1->AddAt(trackN,1);
\r
1266 vertexer1->SetVtxStart(esdVer1);
\r
1267 AliESDVertex *decayVertex1 = (AliESDVertex*)vertexer1->VertexForSelectedESDTracks(trkArray1);
\r
1269 dd1[0]=fESD->GetPrimaryVertexSPD()->GetX()-decayVertex1->GetX();
\r
1270 dd1[1]=fESD->GetPrimaryVertexSPD()->GetY()-decayVertex1->GetY();
\r
1271 dd1[2]=fESD->GetPrimaryVertexSPD()->GetZ()-decayVertex1->GetZ();
\r
1273 decayLengthH=sqrt(dd1[0]*dd1[0]+dd1[1]*dd1[1]+dd1[2]*dd1[2]);
\r
1275 Double_t bz = fESD->GetMagneticField();
\r
1277 trackP->PropagateToDCA(decayVertex1, bz, 10);
\r
1278 trackN->PropagateToDCA(decayVertex1, bz, 10);
\r
1280 Double_t xthiss(0.0);
\r
1281 Double_t xpp(0.0);
\r
1282 Double_t dca = trackN->GetDCA(trackP,bz,xthiss,xpp);
\r
1284 if (decayVertex1) delete decayVertex1;
\r
1285 h1.SetXYZ(-dd1[0],-dd1[1],-dd1[2]);
\r
1287 // if (dca>1) continue;
\r
1288 if (dca>0.1) continue;
\r
1290 fourSum=threeSum+posProt;
\r
1292 posProt1.SetXYZM(trackP->Px(),trackP->Py(),trackP->Pz(),cProtonMass);
\r
1293 negPio1.SetXYZM(trackN->Px(),trackN->Py(),trackN->Pz(),cPionMass);
\r
1294 hDibaryon=lambdaH+posProt1+negPio1;
\r
1295 Double_t hPointingAngle = hDibaryon.Angle(h);
\r
1296 Double_t pointingAngleH = hDibaryon.Angle(h1);
\r
1297 Double_t decayAngleH = h.Angle(h1);
\r
1298 TVector3 vecDist(dd[0]-dd1[0],dd[1]-dd1[1],dd[2]-dd1[2]);
\r
1299 fHistMassLambdaPPi->Fill(hDibaryon.M());
\r
1300 fHistHPointingAngle->Fill(pointingAngleH);
\r
1302 fHistMassHcentMult->Fill(hDibaryon.M(),triggertype,refMultTpc);
\r
1304 Double_t rapidity = hDibaryon.Rapidity();
\r
1305 if(rapidity > 1.0 || rapidity < -1.0) continue;
\r
1307 fHistRap->Fill(rapidity);
\r
1308 //if (pointingAngleH > 0.1) continue;
\r
1309 if (pointingAngleH > 0.05) continue;
\r
1311 ///////////////////////////
\r
1312 //MC part for Associated H
\r
1313 ///////////////////////////
\r
1315 if (HasMC() && isCorrectlyAssociatedLambda) {
\r
1316 Int_t labelP = trackP->GetLabel();
\r
1317 TParticle *tparticleDaughter = stack->Particle(TMath::Abs(labelP));
\r
1318 Int_t labelMother = tparticleDaughter->GetFirstMother();
\r
1319 TParticle *tparticleMother = stack->Particle(TMath::Abs(labelMother));
\r
1321 Int_t labelProton = trackP->GetLabel();
\r
1322 Int_t labelPion = trackN->GetLabel();
\r
1325 if(tparticleMother->GetPdgCode() == pdgHDibaryon && labelAssociatedH==labelMother){ //check mother PDG
\r
1326 Int_t labelFirstDaughter = tparticleMother->GetDaughter(0);
\r
1327 Int_t labelThirdDaughter = tparticleMother->GetDaughter(1);
\r
1328 Int_t labelSecondDaughter = labelFirstDaughter +1;
\r
1330 TParticle *tparticleSecondDaughter = stack->Particle(TMath::Abs(labelSecondDaughter));
\r
1331 TParticle *tparticleThirdDaughter = stack->Particle(TMath::Abs(labelThirdDaughter));
\r
1333 TLorentzVector ppi;
\r
1334 TLorentzVector lpi;
\r
1335 TLorentzVector lP;
\r
1337 ppi=posProt+negPio;
\r
1338 lpi=lambdaH+negPio;
\r
1339 lP=lambdaH+posProt;
\r
1341 if((tparticleThirdDaughter->GetPdgCode() == pdgPionMinus && labelPion==labelThirdDaughter)&&(tparticleSecondDaughter->GetPdgCode() == pdgProton||tparticleSecondDaughter->GetPdgCode() == pdgPionMinus) && labelProton==labelSecondDaughter) fHistMassPpi->Fill(ppi.M());
\r
1343 if(tparticleThirdDaughter->GetPdgCode() == pdgPionMinus && labelPion==labelThirdDaughter) fHistMassLpi->Fill(lpi.M());
\r
1345 if(tparticleSecondDaughter->GetPdgCode() == pdgProton && labelProton==labelSecondDaughter) fHistMassLP->Fill(lP.M());
\r
1347 if(tparticleSecondDaughter->GetPdgCode() == pdgProton && labelProton==labelSecondDaughter){//check second daughter PDG
\r
1348 if(tparticleThirdDaughter->GetPdgCode() == pdgPionMinus && labelPion==labelThirdDaughter){//check second daughter PDG
\r
1350 fHistHDibaryonInvaMassAso->Fill(hDibaryon.M());
\r
1352 Double_t distance01=vecDist.Mag();
\r
1353 fHistMCdcaPvtxDvtx->Fill(decayLengthH);
\r
1354 fHistMCdcaPvtxLvtx->Fill(decayLength);
\r
1355 fHistMCdcaDvtxLvtx->Fill(distance01);
\r
1356 fHistMCangleLH->Fill(hPointingAngle);
\r
1357 fHistMCdecayAngle->Fill(decayAngleH);
\r
1358 fHistMCpointingAngle->Fill(pointingAngleH);
\r
1359 fHistMCap->Fill(alfa,qt);
\r
1361 fHistHilf1->Fill(posPionKF.GetDistanceFromVertex(primVtx));
\r
1362 fHistHilf2->Fill(protonKF.GetDistanceFromVertex(primVtx));
\r
1363 fHistHilf3->Fill(protonKF.GetDistanceFromVertex(posPionKF));
\r
1364 fHistHilf6->Fill(dca);
\r
1365 fHistPtvsYAso->Fill(hDibaryon.Pt(),hDibaryon.Rapidity());
\r
1366 fHistPtvsEtaAso->Fill(hDibaryon.Pt(),hDibaryon.Eta());
\r
1368 }//end check for third daughter PDG
\r
1369 }//end check second daughter PDG
\r
1373 // cout<<"Trigger: "<<triggertype<<endl;
\r
1374 fHistMassH->Fill(hDibaryon.M());
\r
1375 fHistMassHcentMult->Fill(hDibaryon.M(),triggertype,refMultTpc);
\r
1376 ppK=lambdaH+posProt;
\r
1377 fHistMassLambdaP->Fill(ppK.M());
\r
1379 // fHistNdim = new THnSparseF("fHistNdim","THnS;InvMass, InvMassLambda, pointingAngle, armPoAlpha, armPoQt, pTL, pTH, d0p, d0n, dcaHd, dca, decayL, cosPA, centr, multi, mcinf;InvMassH", 16,binsD01,xminD01,xmaxD01);
\r
1381 Double_t vec[16]={hDibaryon.M(), lInvMassLambda, pointingAngleH, alfa, qt, lPtLambda, hDibaryon.Pt(), posPionKF.GetDistanceFromVertex(primVtx), protonKF.GetDistanceFromVertex(primVtx), dca, protonKF.GetDistanceFromVertex(posPionKF), TMath::Cos(pointingAngleH), centrPerc, refMultTpc, mcStatus};
\r
1382 fHistNdim->Fill(vec);
\r
1390 //!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
\r
1391 //Pure MC Part!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
\r
1392 //!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
\r
1394 // Monte Carlo for genenerated particles
\r
1395 if (HasMC()) //MC loop
\r
1400 Double_t momentumPionGen[3]={0,0,0};
\r
1401 Double_t momentumNucleonGen[3]={0,0,0};
\r
1402 Double_t momentumLambdaGen[3]={0,0,0};
\r
1404 Double_t energyPionGen = 0;
\r
1405 Double_t energyNucleonGen = 0;
\r
1406 Double_t energyLambdaGen = 0;
\r
1408 Double_t transversMomentumMotherGen = 0;
\r
1409 Double_t longitudinalMomentumMotherGen = 0;
\r
1410 Double_t totalEnergyMotherGen = 0;
\r
1412 Double_t rapidityGen = 2;
\r
1414 for(stackN = 0; stackN < stack->GetNtrack(); stackN++) //loop over stack
\r
1417 TParticle *tparticleMother = stack->Particle(stackN);
\r
1419 if(tparticleMother->GetPdgCode() == pdgLambda) fHistCount->Fill(16);
\r
1422 if(tparticleMother->GetPdgCode() == pdgHDibaryon) //check mother PDG
\r
1424 Int_t labelFirstDaughter = tparticleMother->GetDaughter(0);
\r
1425 Int_t labelThirdDaughter = tparticleMother->GetDaughter(1);
\r
1426 Int_t labelSecondDaughter = labelFirstDaughter +1;
\r
1428 TParticle *tparticleFirstDaughter = stack->Particle(TMath::Abs(labelFirstDaughter));
\r
1429 TParticle *tparticleSecondDaughter = stack->Particle(TMath::Abs(labelSecondDaughter));
\r
1430 TParticle *tparticleThirdDaughter = stack->Particle(TMath::Abs(labelThirdDaughter));
\r
1432 if(tparticleFirstDaughter->GetPdgCode() == pdgLambda) //check first daughter PDG
\r
1434 if(tparticleSecondDaughter->GetPdgCode() == pdgProton)//check second daughter PDG
\r
1436 if(tparticleThirdDaughter->GetPdgCode() == pdgPionMinus)//check second daughter PDG
\r
1438 momentumLambdaGen[0] = tparticleFirstDaughter->Px();
\r
1439 momentumLambdaGen[1] = tparticleFirstDaughter->Py();
\r
1440 momentumLambdaGen[2] = tparticleFirstDaughter->Pz();
\r
1442 momentumNucleonGen[0] = tparticleSecondDaughter->Px();
\r
1443 momentumNucleonGen[1] = tparticleSecondDaughter->Py();
\r
1444 momentumNucleonGen[2] = tparticleSecondDaughter->Pz();
\r
1446 momentumPionGen[0] = tparticleThirdDaughter->Px();
\r
1447 momentumPionGen[1] = tparticleThirdDaughter->Py();
\r
1448 momentumPionGen[2] = tparticleThirdDaughter->Pz();
\r
1450 TLorentzVector lorentzVectorLambda;
\r
1451 TLorentzVector lorentzVectorProton;
\r
1452 TLorentzVector lorentzVectorPion;
\r
1453 TLorentzVector lorentzVectorHDibaryon;
\r
1455 lorentzVectorLambda.SetXYZM(momentumLambdaGen[0],momentumLambdaGen[1],momentumLambdaGen[2],1.115);
\r
1456 lorentzVectorProton.SetXYZM(momentumNucleonGen[0],momentumNucleonGen[1],momentumNucleonGen[2],protonMass);
\r
1457 lorentzVectorPion.SetXYZM(momentumPionGen[0],momentumPionGen[1],momentumPionGen[2],pionMass);
\r
1459 lorentzVectorHDibaryon = lorentzVectorLambda + lorentzVectorProton + lorentzVectorPion;
\r
1460 rapidityGen=lorentzVectorHDibaryon.Rapidity();
\r
1461 transversMomentumMotherGen = lorentzVectorHDibaryon.Pt();
\r
1462 longitudinalMomentumMotherGen = lorentzVectorHDibaryon.Pz();
\r
1463 totalEnergyMotherGen = lorentzVectorHDibaryon.Energy();
\r
1465 if(rapidityGen > 1.0 || rapidityGen < -1 ) continue;
\r
1466 //lorentzVectorLambda
\r
1467 fHistHDibaryonInvaMassGen->Fill(lorentzVectorHDibaryon.M());
\r
1468 if (lorentzVectorLambda.Rapidity() > 1.0 || lorentzVectorLambda.Rapidity() < -1) continue;
\r
1469 if (lorentzVectorProton.Rapidity() > 1.0 || lorentzVectorProton.Rapidity() < -1) continue;
\r
1471 if (lorentzVectorPion.Rapidity() > 1.0 || lorentzVectorPion.Rapidity() < -1) continue;
\r
1472 fHistHDibaryonInvaMassGenRes->Fill(lorentzVectorHDibaryon.M());
\r
1473 fHistPtvsEtaGen->Fill(lorentzVectorHDibaryon.Pt(),lorentzVectorHDibaryon.Eta());
\r
1474 fHistPtvsYGen->Fill(lorentzVectorHDibaryon.Pt(),lorentzVectorHDibaryon.Rapidity());
\r
1475 fHistPtvsEtaGen->Fill(lorentzVectorHDibaryon.Pt(),lorentzVectorHDibaryon.Eta());
\r
1476 fHistCount->Fill(11);
\r
1477 }//end of check third daughter PDG
\r
1478 }//end of check second daughter PDG
\r
1479 }//end of check first daughter PDG
\r
1480 }//end of H-Dibaryon
\r
1483 if(tparticleMother->GetPdgCode() == pdgAntiHDibaryon) //check mother PDG
\r
1485 Int_t labelFirstDaughter = tparticleMother->GetDaughter(0);
\r
1486 Int_t labelThirdDaughter = tparticleMother->GetDaughter(1);
\r
1487 Int_t labelSecondDaughter = labelFirstDaughter +1;
\r
1489 TParticle *tparticleFirstDaughter = stack->Particle(TMath::Abs(labelFirstDaughter));
\r
1490 TParticle *tparticleSecondDaughter = stack->Particle(TMath::Abs(labelSecondDaughter));
\r
1491 TParticle *tparticleThirdDaughter = stack->Particle(TMath::Abs(labelThirdDaughter));
\r
1493 if(tparticleFirstDaughter->GetPdgCode() == pdgAntiLambda) //check first daughter PDG
\r
1495 if(tparticleSecondDaughter->GetPdgCode() == pdgAntiProton)//check second daughter PDG
\r
1497 if(tparticleThirdDaughter->GetPdgCode() == pdgPionPlus)//check second daughter PDG
\r
1499 momentumLambdaGen[0] = tparticleFirstDaughter->Px();
\r
1500 momentumLambdaGen[1] = tparticleFirstDaughter->Py();
\r
1501 momentumLambdaGen[2] = tparticleFirstDaughter->Pz();
\r
1503 momentumNucleonGen[0] = tparticleSecondDaughter->Px();
\r
1504 momentumNucleonGen[1] = tparticleSecondDaughter->Py();
\r
1505 momentumNucleonGen[2] = tparticleSecondDaughter->Pz();
\r
1507 momentumPionGen[0] = tparticleThirdDaughter->Px();
\r
1508 momentumPionGen[1] = tparticleThirdDaughter->Py();
\r
1509 momentumPionGen[2] = tparticleThirdDaughter->Pz();
\r
1511 energyLambdaGen = tparticleFirstDaughter->Energy();
\r
1512 energyNucleonGen = tparticleSecondDaughter->Energy();
\r
1513 energyPionGen = tparticleThirdDaughter->Energy();
\r
1515 TLorentzVector lorentzVectorLambda;
\r
1516 TLorentzVector lorentzVectorProton;
\r
1517 TLorentzVector lorentzVectorPion;
\r
1518 TLorentzVector lorentzVectorHDibaryon;
\r
1520 lorentzVectorLambda.SetXYZM(momentumLambdaGen[0],momentumLambdaGen[1],momentumLambdaGen[2],1.115);
\r
1521 lorentzVectorProton.SetXYZM(momentumNucleonGen[0],momentumNucleonGen[1],momentumNucleonGen[2],protonMass);
\r
1522 lorentzVectorPion.SetXYZM(momentumPionGen[0],momentumPionGen[1],momentumPionGen[2],pionMass);
\r
1524 lorentzVectorHDibaryon = lorentzVectorLambda + lorentzVectorProton + lorentzVectorPion;
\r
1526 rapidityGen=lorentzVectorHDibaryon.Rapidity();
\r
1527 if(rapidityGen > 1.0 || rapidityGen < -1 ) continue;
\r
1528 fHistAntiHDibaryonInvaMassGen->Fill(lorentzVectorHDibaryon.M());
\r
1529 }//end of check third daughter PDG
\r
1530 }//end of check second daughter PDG
\r
1531 }//end of check first daughter PDG
\r
1532 }//end of Anti-H-Dibaryon
\r
1536 // Post output data.
\r
1537 PostData(1,fHistList);
\r
1538 //PostData(0,fHistList);
\r
1540 if (listCrossV0) delete listCrossV0;
\r
1541 if (esdVer1) delete esdVer1;
\r
1542 if (vertexer) delete vertexer;
\r
1543 if (vertexer1) delete vertexer1;
\r
1544 if (trkArray) delete trkArray;
\r
1545 if (trkArray1) delete trkArray1;
\r
1548 //________________________________________________________________________
\r
1549 void AliAnalysisTaskHdibaryonLPpi::Terminate(Option_t *)
\r
1551 // Draw result to the screen
\r
1552 // Called once at the end of the query
\r