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
30 #include "TCanvas.h"
\r
31 #include "TParticle.h"
\r
32 #include "TNtuple.h"
\r
33 #include "TObjString.h"
\r
34 #include "TLorentzVector.h"
\r
35 #include "TDatabasePDG.h"
\r
37 #include "AliAnalysisTask.h"
\r
38 #include "AliAnalysisManager.h"
\r
39 #include "AliMCEventHandler.h"
\r
40 #include "AliMCEvent.h"
\r
43 #include "AliESDpid.h"
\r
44 #include "AliESDEvent.h"
\r
45 #include "AliGenEventHeader.h"
\r
46 #include "AliESDInputHandler.h"
\r
48 #include "AliStack.h"
\r
49 #include "AliHeader.h"
\r
51 #include "AliKFParticle.h"
\r
52 #include "AliKFVertex.h"
\r
53 #include "AliESDtrackCuts.h"
\r
54 #include "AliESDv0Cuts.h"
\r
55 #include "AliAnalysisTaskHdibaryonLPpi.h"
\r
60 #include "TCanvas.h"
\r
62 #include "TParallelCoord.h"
\r
64 #include "AliMCParticle.h"
\r
65 #include "AliGenPythiaEventHeader.h"
\r
68 #include "AliESDtrack.h"
\r
69 #include "AliCentrality.h"
\r
71 #include "THnSparse.h"
\r
73 #include "AliVertexerTracks.h"
\r
75 using namespace std;
\r
76 ClassImp(AliAnalysisTaskHdibaryonLPpi)
\r
77 //________________________________________________________________________
\r
78 AliAnalysisTaskHdibaryonLPpi::AliAnalysisTaskHdibaryonLPpi() : AliAnalysisTaskSE()/*AliAnalysisTask(name, ""), fMCEvent(0)*/, fESD(0), fESDtrackCutsV0(0),
\r
86 fHistMassLambdaPi(0),
\r
88 fHistMassLambdaPPi(0),
\r
89 fHistMassLambdaP(0),
\r
90 fHistMassLambdaK(0),
\r
94 fHistMassK0offlC(0),
\r
98 fHistArmenterosPodolanski(0),
\r
99 fHistArmenterosPodolanskiCut(0),
\r
100 fHistHDibaryonInvaMassGen(0),
\r
101 fHistHDibaryonInvaMassGenRes(0),
\r
102 fHistAntiHDibaryonInvaMassGen(0),
\r
103 fHistHDibaryonInvaMassAso(0),
\r
104 fHistHDibaryonInvaMassAsoReso(0),
\r
105 fHistAntiHDibaryonInvaMassAso(0),
\r
107 fHistHPointingAngle(0),
\r
109 fHistMassLambdaFromH(0),
\r
110 fHistMassLambdaFromHtLorentz(0),
\r
112 fHistMassPpiReso(0),
\r
115 fHistProtonPIDBb(0),
\r
117 fHistProtonPIDLambda(0),
\r
118 fHistPionPIDLambda(0),
\r
119 fHistMCdcaPvtxDvtx(0),
\r
120 fHistMCdcaPvtxLvtx(0),
\r
121 fHistMCdcaDvtxLvtx(0),
\r
123 fHistMCdecayAngle(0),
\r
124 fHistMCpointingAngle(0),
\r
126 fHistMCdcaPvtxDvtxReso(0),
\r
127 fHistMCdcaPvtxLvtxReso(0),
\r
128 fHistMCdcaDvtxLvtxReso(0),
\r
129 fHistMCangleLHReso(0),
\r
130 fHistMCdecayAngleReso(0),
\r
131 fHistMCpointingAngleReso(0),
\r
140 fHistPtvsEtaGen(0),
\r
141 fHistPtvsEtaAso(0),
\r
151 // Define input and output slots here
\r
152 // Input from a TChain
\r
153 DefineInput(0, TChain::Class());
\r
154 // Output to TList container
\r
155 DefineOutput(1, TList::Class()); //full
\r
158 if (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())SetHasMC();
\r
162 fESDtrackCutsV0 = new AliESDtrackCuts("AliESDtrackCutsV0","AliESDtrackCutsV0");
\r
163 fESDtrackCutsV0->SetAcceptKinkDaughters(kFALSE);
\r
164 fESDtrackCutsV0->SetMinNClustersTPC(80);
\r
165 fESDtrackCutsV0->SetMaxChi2PerClusterTPC(5);
\r
166 fESDtrackCutsV0->SetRequireTPCRefit(kTRUE);
\r
167 fESDtrackCutsV0->SetEtaRange(-0.9,0.9);
\r
169 fESDCutsV0 = new AliESDv0Cuts("AliESDCutsV0","AliESDCutsV0");
\r
170 fESDCutsV0->SetMaxDcaV0Daughters(2.0);
\r
171 fESDCutsV0->SetMinDcaNegToVertex(1.0);
\r
172 fESDCutsV0->SetMinDcaPosToVertex(1.0);
\r
175 fEsdTrackCuts = new AliESDtrackCuts("AliESDtrackCuts");
\r
176 fEsdTrackCuts->SetMinNClustersTPC(80);
\r
177 fEsdTrackCuts->SetAcceptKinkDaughters(kFALSE);
\r
178 fEsdTrackCuts->SetMaxChi2PerClusterTPC(5);
\r
179 fEsdTrackCuts->SetRequireTPCRefit(kTRUE);
\r
180 fEsdTrackCuts->SetEtaRange(-0.9,0.9);
\r
183 //________________________________________________________________________
\r
184 AliAnalysisTaskHdibaryonLPpi::AliAnalysisTaskHdibaryonLPpi(const char *name) : AliAnalysisTaskSE(name)/*AliAnalysisTask(name, ""), fMCEvent(0)*/, fESD(0), fESDtrackCutsV0(0),
\r
192 fHistMassLambdaPi(0),
\r
193 fHistMassLambda(0),
\r
194 fHistMassLambdaPPi(0),
\r
195 fHistMassLambdaP(0),
\r
196 fHistMassLambdaK(0),
\r
198 fHistMassK0offl(0),
\r
199 fHistMassK0onlC(0),
\r
200 fHistMassK0offlC(0),
\r
201 fHistMassPQonl(0),
\r
202 fHistMassPQoffl(0),
\r
204 fHistArmenterosPodolanski(0),
\r
205 fHistArmenterosPodolanskiCut(0),
\r
206 fHistHDibaryonInvaMassGen(0),
\r
207 fHistHDibaryonInvaMassGenRes(0),
\r
208 fHistAntiHDibaryonInvaMassGen(0),
\r
209 fHistHDibaryonInvaMassAso(0),
\r
210 fHistHDibaryonInvaMassAsoReso(0),
\r
211 fHistAntiHDibaryonInvaMassAso(0),
\r
213 fHistHPointingAngle(0),
\r
215 fHistMassLambdaFromH(0),
\r
216 fHistMassLambdaFromHtLorentz(0),
\r
218 fHistMassPpiReso(0),
\r
221 fHistProtonPIDBb(0),
\r
223 fHistProtonPIDLambda(0),
\r
224 fHistPionPIDLambda(0),
\r
225 fHistMCdcaPvtxDvtx(0),
\r
226 fHistMCdcaPvtxLvtx(0),
\r
227 fHistMCdcaDvtxLvtx(0),
\r
229 fHistMCdecayAngle(0),
\r
230 fHistMCpointingAngle(0),
\r
232 fHistMCdcaPvtxDvtxReso(0),
\r
233 fHistMCdcaPvtxLvtxReso(0),
\r
234 fHistMCdcaDvtxLvtxReso(0),
\r
235 fHistMCangleLHReso(0),
\r
236 fHistMCdecayAngleReso(0),
\r
237 fHistMCpointingAngleReso(0),
\r
246 fHistPtvsEtaGen(0),
\r
247 fHistPtvsEtaAso(0),
\r
257 // Define input and output slots here
\r
258 // Input from a TChain
\r
259 DefineInput(0, TChain::Class());
\r
260 // Output to TList container
\r
261 DefineOutput(1, TList::Class()); //full
\r
264 if (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())SetHasMC();
\r
268 fESDtrackCutsV0 = new AliESDtrackCuts("AliESDtrackCutsV0","AliESDtrackCutsV0");
\r
269 fESDtrackCutsV0->SetAcceptKinkDaughters(kFALSE);
\r
270 fESDtrackCutsV0->SetMinNClustersTPC(80);
\r
271 fESDtrackCutsV0->SetMaxChi2PerClusterTPC(5);
\r
272 fESDtrackCutsV0->SetRequireTPCRefit(kTRUE);
\r
273 fESDtrackCutsV0->SetEtaRange(-0.9,0.9);
\r
275 fESDCutsV0 = new AliESDv0Cuts("AliESDCutsV0","AliESDCutsV0");
\r
276 fESDCutsV0->SetMaxDcaV0Daughters(2.0);
\r
277 fESDCutsV0->SetMinDcaNegToVertex(1.0);
\r
278 fESDCutsV0->SetMinDcaPosToVertex(1.0);
\r
281 fEsdTrackCuts = new AliESDtrackCuts("AliESDtrackCuts");
\r
282 fEsdTrackCuts->SetMinNClustersTPC(80);
\r
283 fEsdTrackCuts->SetAcceptKinkDaughters(kFALSE);
\r
284 fEsdTrackCuts->SetMaxChi2PerClusterTPC(5);
\r
285 fEsdTrackCuts->SetRequireTPCRefit(kTRUE);
\r
286 fEsdTrackCuts->SetEtaRange(-0.9,0.9);
\r
289 //____________________________________________________________
\r
290 AliAnalysisTaskHdibaryonLPpi::~AliAnalysisTaskHdibaryonLPpi(){
\r
295 fHistList->Clear();
\r
298 if(fEsdTrackCuts) delete fEsdTrackCuts;
\r
299 if(fESDtrackCutsV0) delete fESDtrackCutsV0;
\r
300 if(fESDCutsV0) delete fESDCutsV0;
\r
304 //________________________________________________________________________
\r
305 void AliAnalysisTaskHdibaryonLPpi::UserCreateOutputObjects()
\r
307 // Create histograms
\r
310 fHistList = new TList();
\r
312 fHistMassDPi = new TH1F("fHistMassDPi", "Invariant mass distribution p+#pi^{-} ", 500, 1.0, 1.25);
\r
313 fHistMassDPi->GetXaxis()->SetTitle("Invariant mass p+#pi^{-} (GeV/c^{2})");
\r
314 fHistMassDPi->GetYaxis()->SetTitle("Entries");
\r
315 fHistMassDPi->SetMarkerStyle(kFullCircle);
\r
317 fHistMassLPi = new TH1F("fHistMassLPi", "Offline Invariant mass distribution p+#pi^{-} ", 500, 1.0, 1.25);
\r
318 fHistMassLPi->GetXaxis()->SetTitle("Invariant mass p+#pi^{-} (GeV/c^{2})");
\r
319 fHistMassLPi->GetYaxis()->SetTitle("Entries");
\r
320 fHistMassLPi->SetMarkerStyle(kFullCircle);
\r
322 fHistMassLambdaPi = new TH1F("fHistMassLambdaPi", "Invariant mass distribution #Lambda+#pi^{-} ", 500, 1.2, 1.5);
\r
323 fHistMassLambdaPi->GetXaxis()->SetTitle("Invariant mass #Lambda+#pi^{-} (GeV/c^{2})");
\r
324 fHistMassLambdaPi->GetYaxis()->SetTitle("Entries");
\r
325 fHistMassLambdaPi->SetMarkerStyle(kFullCircle);
\r
327 fHistMassLambda = new TH1F("fHistMassLambda", "Invariant mass distribution of #Lambda for further analyis", 500, 1.0, 1.2);
\r
328 fHistMassLambda->GetXaxis()->SetTitle("Invariant mass p+#pi^{+} (GeV/c^{2})");
\r
329 fHistMassLambda->GetYaxis()->SetTitle("Entries");
\r
330 fHistMassLambda->SetMarkerStyle(kFullCircle);
\r
332 fHistMassLambdaPPi = new TH1F("fHistMassLambdaPPi", "Invariant mass distribution #Lambdap#pi^{-} ", 300, 2.2, 2.5);
\r
333 fHistMassLambdaPPi->GetXaxis()->SetTitle("Invariant mass #Lambdap#pi^{-} (GeV/c^{2})");
\r
334 fHistMassLambdaPPi->GetYaxis()->SetTitle("Entries");
\r
335 fHistMassLambdaPPi->SetMarkerStyle(kFullCircle);
\r
337 fHistMassLambdaP = new TH1F("fHistMassLambdaP", "Invariant mass distribution #Lambdap ", 300, 2.2, 2.5);
\r
338 fHistMassLambdaP->GetXaxis()->SetTitle("Invariant mass #Lambdap (GeV/c^{2})");
\r
339 fHistMassLambdaP->GetYaxis()->SetTitle("Entries");
\r
340 fHistMassLambdaP->SetMarkerStyle(kFullCircle);
\r
342 fHistMassLambdaK = new TH1F("fHistMassLambdaK", "Invariant mass distribution #LambdaK^{-} ", 300, 1.6, 1.9);
\r
343 fHistMassLambdaK->GetXaxis()->SetTitle("Invariant mass #LambdaK^{-} (GeV/c^{2})");
\r
344 fHistMassLambdaK->GetYaxis()->SetTitle("Entries");
\r
345 fHistMassLambdaK->SetMarkerStyle(kFullCircle);
\r
347 fHistMassK0onl = new TH1F("fHistMassK0onl", "Invariant mass distribution K_{s}^{0} online V0 finder ", 400, 0.2, 1.0);
\r
348 fHistMassK0onl->GetXaxis()->SetTitle("Invariant mass #pi^{+}#pi^{-} (GeV/c^{2})");
\r
349 fHistMassK0onl->GetYaxis()->SetTitle("Entries");
\r
350 fHistMassK0onl->SetMarkerStyle(kFullCircle);
\r
352 fHistMassK0offl = new TH1F("fHistMassK0offl", "Invariant mass distribution K_{s}^{0} offline V0 finder ", 400, 0.2, 1.0);
\r
353 fHistMassK0offl->GetXaxis()->SetTitle("Invariant mass #pi^{+}#pi^{-} (GeV/c^{2})");
\r
354 fHistMassK0offl->GetYaxis()->SetTitle("Entries");
\r
355 fHistMassK0offl->SetMarkerStyle(kFullCircle);
\r
357 fHistMassK0onlC = new TH1F("fHistMassK0onlC", "Invariant mass distribution K_{s}^{0} online V0 finder ", 400, 0.2, 1.0);
\r
358 fHistMassK0onlC->GetXaxis()->SetTitle("Invariant mass #pi^{+}#pi^{-} (GeV/c^{2})");
\r
359 fHistMassK0onlC->GetYaxis()->SetTitle("Entries");
\r
360 fHistMassK0onlC->SetMarkerStyle(kFullCircle);
\r
362 fHistMassK0offlC = new TH1F("fHistMassK0offlC", "Invariant mass distribution K_{s}^{0} offline V0 finder ", 400, 0.2, 1.0);
\r
363 fHistMassK0offlC->GetXaxis()->SetTitle("Invariant mass #pi^{+}#pi^{-} (GeV/c^{2})");
\r
364 fHistMassK0offlC->GetYaxis()->SetTitle("Entries");
\r
365 fHistMassK0offlC->SetMarkerStyle(kFullCircle);
\r
367 fHistMassPQonl = new TH1F("fHistMassPQonl", "Invariant mass distribution K_{s}^{0}p using online V0 finder ", 500, 1.3, 2.3);
\r
368 fHistMassPQonl->GetXaxis()->SetTitle("Invariant mass K_{s}^{0}p (GeV/c^{2})");
\r
369 fHistMassPQonl->GetYaxis()->SetTitle("Entries");
\r
370 fHistMassPQonl->SetMarkerStyle(kFullCircle);
\r
372 fHistMassPQoffl = new TH1F("fHistMassPQoffl", "Invariant mass distribution K_{s}^{0}p using offline V0 finder ", 500, 1.3, 2.3);
\r
373 fHistMassPQoffl->GetXaxis()->SetTitle("Invariant mass K_{s}^{0}p (GeV/c^{2})");
\r
374 fHistMassPQoffl->GetYaxis()->SetTitle("Entries");
\r
375 fHistMassPQoffl->SetMarkerStyle(kFullCircle);
\r
377 fHistDC = new TH1F("fHistDC", "Proper decay length", 500, 0.0, 25);
\r
378 fHistDC->GetXaxis()->SetTitle("c#tau (cm)");
\r
379 fHistDC->GetYaxis()->SetTitle("Entries");
\r
380 fHistDC->SetMarkerStyle(kFullCircle);
\r
382 fHistArmenterosPodolanski = new TH2F("fHistArmenterosPodolanski", "Armenteros-Podolanski plot", 200,-1.0,1.0, 500,0,1);
\r
383 fHistArmenterosPodolanski->GetXaxis()->SetTitle("#alpha");
\r
384 fHistArmenterosPodolanski->GetYaxis()->SetTitle("q_{t}");
\r
385 fHistArmenterosPodolanski->SetMarkerStyle(kFullCircle);
\r
387 fHistArmenterosPodolanskiCut = new TH2F("fHistArmenterosPodolanskiCut", "Armenteros-Podolanski plot after cut", 200,-1.0,1.0, 500,0,1);
\r
388 fHistArmenterosPodolanskiCut->GetXaxis()->SetTitle("#alpha");
\r
389 fHistArmenterosPodolanskiCut->GetYaxis()->SetTitle("q_{t}");
\r
390 fHistArmenterosPodolanskiCut->SetMarkerStyle(kFullCircle);
\r
392 fHistHDibaryonInvaMassGen = new TH1F("fHistHDibaryonInvaMassGen", "Generated #Lambda p #pi^{-}", 200, 2.1, 2.3);
\r
393 fHistHDibaryonInvaMassGen->GetYaxis()->SetTitle("Counts");
\r
394 fHistHDibaryonInvaMassGen->GetXaxis()->SetTitle("Invariant mass (GeV/c^{2})");
\r
396 fHistHDibaryonInvaMassGenRes = new TH1F("fHistHDibaryonInvaMassGenRes", "Generated #Lambda p #pi^{-} with particles in rapidity!!!", 200, 2.1, 2.3);
\r
397 fHistHDibaryonInvaMassGenRes->GetYaxis()->SetTitle("Counts");
\r
398 fHistHDibaryonInvaMassGenRes->GetXaxis()->SetTitle("Invariant mass (GeV/c^{2})");
\r
400 fHistAntiHDibaryonInvaMassGen = new TH1F("fHistAntiHDibaryonInvaMassGen", "Generated #bar{#Lambda} #bar{p} #pi^{+}", 200, 2.1, 2.3);
\r
401 fHistAntiHDibaryonInvaMassGen->GetYaxis()->SetTitle("Counts");
\r
402 fHistAntiHDibaryonInvaMassGen->GetXaxis()->SetTitle("Invariant mass (GeV/c^{2})");
\r
404 fHistHDibaryonInvaMassAso = new TH1F("fHistHDibaryonInvaMassAso", "Associated #Lambda p #pi^{-}", 200, 2.1, 2.3);
\r
405 fHistHDibaryonInvaMassAso->GetYaxis()->SetTitle("Counts");
\r
406 fHistHDibaryonInvaMassAso->GetXaxis()->SetTitle("Invariant mass (GeV/c^{2})");
\r
408 fHistHDibaryonInvaMassAsoReso = new TH1F("fHistHDibaryonInvaMassAsoReso", "Associated #Lambda p #pi^{-}", 200, 2.1, 2.3);
\r
409 fHistHDibaryonInvaMassAsoReso->GetYaxis()->SetTitle("Counts");
\r
410 fHistHDibaryonInvaMassAsoReso->GetXaxis()->SetTitle("Invariant mass (GeV/c^{2})");
\r
412 fHistAntiHDibaryonInvaMassAso = new TH1F("fHistAntiHDibaryonInvaMassAso", "Associated #bar{#Lambda} #bar{p} #pi^{+}", 200, 2.1, 2.3);
\r
413 fHistAntiHDibaryonInvaMassAso->GetYaxis()->SetTitle("Counts");
\r
414 fHistAntiHDibaryonInvaMassAso->GetXaxis()->SetTitle("Invariant mass (GeV/c^{2})");
\r
416 fHistCheck = new TH2F("fHistCheck", "Check online/offline", 200, -0.5, 1.5, 200, -0.5, 1.5);
\r
417 fHistCheck->GetXaxis()->SetTitle("offline");
\r
418 fHistCheck->GetYaxis()->SetTitle("online");
\r
419 fHistCheck->SetMarkerStyle(kFullCircle);
\r
421 fHistHPointingAngle= new TH1F("fHistHPointingAngle", "Pointing angle distribution for #Lambdap#pi^{-}", 200, 0., 2*TMath::Pi());
\r
422 fHistHPointingAngle->GetXaxis()->SetTitle("Pointing angle distribution for #Lambdap#pi^{-}");
\r
423 fHistHPointingAngle->GetYaxis()->SetTitle("Entries");
\r
424 fHistHPointingAngle->SetMarkerStyle(kFullCircle);
\r
426 fHistMassH= new TH1F("fHistMassH", "Invariant mass distribution #Lambdap#pi^{-} (GeV/c^{2}) after pointing angle cut", 3000, 2.2, 2.5);
\r
427 fHistMassH->GetXaxis()->SetTitle("Invariant mass #Lambdap#pi^{-} (GeV/c^{2})");
\r
428 fHistMassH->GetYaxis()->SetTitle("Entries");
\r
429 fHistMassH->SetMarkerStyle(kFullCircle);
\r
431 fHistMassLambdaFromH= new TH1F("fHistMassLambdaFromH", "Invariant mass distribution #Lambda (GeV/c^{2}) asking for Mother to be H", 300, 1.0, 1.3);
\r
432 fHistMassLambdaFromH->GetXaxis()->SetTitle("Invariant mass #Lambda (GeV/c^{2})");
\r
433 fHistMassLambdaFromH->GetYaxis()->SetTitle("Entries");
\r
434 fHistMassLambdaFromH->SetMarkerStyle(kFullCircle);
\r
436 fHistMassLambdaFromHtLorentz= new TH1F("fHistMassLambdaFromHtLorentz", "Invariant mass distribution #Lambda (GeV/c^{2}) asking for Mother to be H", 300, 1.0, 1.3);
\r
437 fHistMassLambdaFromHtLorentz->GetXaxis()->SetTitle("Invariant mass #Lambda (GeV/c^{2})");
\r
438 fHistMassLambdaFromHtLorentz->GetYaxis()->SetTitle("Entries");
\r
439 fHistMassLambdaFromHtLorentz->SetMarkerStyle(kFullCircle);
\r
441 fHistMassPpi= new TH1F("fHistMassPpi", "Invariant mass distribution of the p#pi^{-} used for combing with Lambda", 300, 1.0, 1.6);
\r
442 fHistMassPpi->GetXaxis()->SetTitle("Invariant mass p#pi^{-} (GeV/c^{2})");
\r
443 fHistMassPpi->GetYaxis()->SetTitle("Entries");
\r
444 fHistMassPpi->SetMarkerStyle(kFullCircle);
\r
446 fHistMassPpiReso= new TH1F("fHistMassPpiReso", "Invariant mass distribution of the p#pi^{-} used for combing with Lambda", 300, 1.0, 1.6);
\r
447 fHistMassPpiReso->GetXaxis()->SetTitle("Invariant mass p#pi^{-} (GeV/c^{2})");
\r
448 fHistMassPpiReso->GetYaxis()->SetTitle("Entries");
\r
449 fHistMassPpiReso->SetMarkerStyle(kFullCircle);
\r
451 fHistMassLpi= new TH1F("fHistMassLpi", "Invariant mass distribution of the #Lambda#pi^{-} used for combing with p", 300, 1.1, 1.7);
\r
452 fHistMassLpi->GetXaxis()->SetTitle("Invariant mass #Lambda#pi^{-} (GeV/c^{2})");
\r
453 fHistMassLpi->GetYaxis()->SetTitle("Entries");
\r
454 fHistMassLpi->SetMarkerStyle(kFullCircle);
\r
456 fHistMassLP= new TH1F("fHistMassLP", "Invariant mass distribution of the #Lambda p used for combing with #pi^{-}", 300, 2.0, 2.3);
\r
457 fHistMassLP->GetXaxis()->SetTitle("Invariant mass #Lambda p (GeV/c^{2})");
\r
458 fHistMassLP->GetYaxis()->SetTitle("Entries");
\r
459 fHistMassLP->SetMarkerStyle(kFullCircle);
\r
461 fHistProtonPIDBb = new TH2F("fHistProtonPIDBb", "dE/dx after p PID", 100, 0., 10, 100, 0, 100);
\r
462 fHistProtonPIDBb->GetYaxis()->SetTitle("TPC Signal");
\r
463 fHistProtonPIDBb->GetXaxis()->SetTitle("P (GeV/c)");
\r
464 fHistProtonPIDBb->SetOption("scat");
\r
465 fHistProtonPIDBb->SetMarkerStyle(kFullCircle);
\r
467 fHistPionPIDBb = new TH2F("fHistPionPIDBb", "dE/dx after K PID", 100, 0., 10, 100, 0, 100);
\r
468 fHistPionPIDBb->GetYaxis()->SetTitle("TPC Signal");
\r
469 fHistPionPIDBb->GetXaxis()->SetTitle("P (GeV/c)");
\r
470 fHistPionPIDBb->SetOption("scat");
\r
471 fHistPionPIDBb->SetMarkerStyle(kFullCircle);
\r
473 fHistProtonPIDLambda = new TH2F("fHistProtonPIDLambda", "dE/dx after p PID from V0", 100, 0., 10, 100, 0, 100);
\r
474 fHistProtonPIDLambda->GetYaxis()->SetTitle("TPC Signal");
\r
475 fHistProtonPIDLambda->GetXaxis()->SetTitle("P (GeV/c)");
\r
476 fHistProtonPIDLambda->SetOption("scat");
\r
477 fHistProtonPIDLambda->SetMarkerStyle(kFullCircle);
\r
479 fHistPionPIDLambda = new TH2F("fHistPionPIDLambda", "dE/dx after #pi PID from V0", 100, 0, 10, 100, 0, 100);
\r
480 fHistPionPIDLambda->GetYaxis()->SetTitle("TPC Signal");
\r
481 fHistPionPIDLambda->GetXaxis()->SetTitle("P (GeV/c)");
\r
482 fHistPionPIDLambda->SetOption("scat");
\r
483 fHistPionPIDLambda->SetMarkerStyle(kFullCircle);
\r
485 fHistMCdcaPvtxDvtx= new TH1F("fHistMCdcaPvtxDvtx", "MC True DCA Primary Vertex - H Decay Vertex", 300, -0.1, 11.9);
\r
486 fHistMCdcaPvtxDvtx->GetXaxis()->SetTitle("dca prim. vtx- decay vtx (cm)");
\r
487 fHistMCdcaPvtxDvtx->GetYaxis()->SetTitle("Entries");
\r
488 fHistMCdcaPvtxDvtx->SetMarkerStyle(kFullCircle);
\r
490 fHistMCdcaPvtxLvtx= new TH1F("fHistMCdcaPvtxLvtx", "MC True DCA Primary Vertex - Lambda Decay Vertex", 300, -0.1, 11.9);
\r
491 fHistMCdcaPvtxLvtx->GetXaxis()->SetTitle("dca prim. vtx-#Lambda decay vtx (cm)");
\r
492 fHistMCdcaPvtxLvtx->GetYaxis()->SetTitle("Entries");
\r
493 fHistMCdcaPvtxLvtx->SetMarkerStyle(kFullCircle);
\r
495 fHistMCdcaDvtxLvtx= new TH1F("fHistMCdcaDvtxLvtx", "MC True DCA H Decay Vertex - Lambda Decay Vertex", 300, -0.1, 11.9);
\r
496 fHistMCdcaDvtxLvtx->GetXaxis()->SetTitle("dca H deacy vtx-#Lambda decay vtx (cm)");
\r
497 fHistMCdcaDvtxLvtx->GetYaxis()->SetTitle("Entries");
\r
498 fHistMCdcaDvtxLvtx->SetMarkerStyle(kFullCircle);
\r
500 fHistMCangleLH= new TH1F("fHistMCangleLH", "MC True Angle between #Lambda and H", 300, 0., 2*TMath::Pi());
\r
501 fHistMCangleLH->GetXaxis()->SetTitle("Angle (#Lambda-H)");
\r
502 fHistMCangleLH->GetYaxis()->SetTitle("Entries");
\r
503 fHistMCangleLH->SetMarkerStyle(kFullCircle);
\r
505 fHistMCdecayAngle= new TH1F("fHistMCdecayAngle", "MC True Angle between decay products", 300, 0., 2*TMath::Pi());
\r
506 fHistMCdecayAngle->GetXaxis()->SetTitle("Angle (#Lambda-p#pi)");
\r
507 fHistMCdecayAngle->GetYaxis()->SetTitle("Entries");
\r
508 fHistMCdecayAngle->SetMarkerStyle(kFullCircle);
\r
510 fHistMCpointingAngle= new TH1F("fHistMCpointingAngle", "MC True Pointing Angle", 3000, 0., 2*TMath::Pi());
\r
511 fHistMCpointingAngle->GetXaxis()->SetTitle("Pointing Angle");
\r
512 fHistMCpointingAngle->GetYaxis()->SetTitle("Entries");
\r
513 fHistMCpointingAngle->SetMarkerStyle(kFullCircle);
\r
515 fHistMCap = new TH2F("fHistMCap", "True MC Armenteros-Podolanski", 200,-1.0,1.0, 500,0,1);
\r
516 fHistMCap->GetYaxis()->SetTitle("#alpha");
\r
517 fHistMCap->GetXaxis()->SetTitle("q_{t}");
\r
518 fHistMCap->SetOption("scat");
\r
519 fHistMCap->SetMarkerStyle(kFullCircle);
\r
521 fHistMCdcaPvtxDvtxReso= new TH1F("fHistMCdcaPvtxDvtxReso", "MC True DCA Primary Vertex - H Decay Vertex", 300, -0.1, 11.9);
\r
522 fHistMCdcaPvtxDvtxReso->GetXaxis()->SetTitle("dca prim. vtx- decay vtx (cm)");
\r
523 fHistMCdcaPvtxDvtxReso->GetYaxis()->SetTitle("Entries");
\r
524 fHistMCdcaPvtxDvtxReso->SetMarkerStyle(kFullCircle);
\r
526 fHistMCdcaPvtxLvtxReso= new TH1F("fHistMCdcaPvtxLvtxReso", "MC True DCA Primary Vertex - Lambda Decay Vertex", 300, -0.1, 11.9);
\r
527 fHistMCdcaPvtxLvtxReso->GetXaxis()->SetTitle("dca prim. vtx-#Lambda decay vtx (cm)");
\r
528 fHistMCdcaPvtxLvtxReso->GetYaxis()->SetTitle("Entries");
\r
529 fHistMCdcaPvtxLvtxReso->SetMarkerStyle(kFullCircle);
\r
531 fHistMCdcaDvtxLvtxReso= new TH1F("fHistMCdcaDvtxLvtxReso", "MC True DCA H Decay Vertex - Lambda Decay Vertex", 300, -0.1, 11.9);
\r
532 fHistMCdcaDvtxLvtxReso->GetXaxis()->SetTitle("dca H deacy vtx-#Lambda decay vtx (cm)");
\r
533 fHistMCdcaDvtxLvtxReso->GetYaxis()->SetTitle("Entries");
\r
534 fHistMCdcaDvtxLvtxReso->SetMarkerStyle(kFullCircle);
\r
536 fHistMCangleLHReso= new TH1F("fHistMCangleLHReso", "MC True Angle between #Lambda and H", 300, 0., 2*TMath::Pi());
\r
537 fHistMCangleLHReso->GetXaxis()->SetTitle("Angle (#Lambda-H)");
\r
538 fHistMCangleLHReso->GetYaxis()->SetTitle("Entries");
\r
539 fHistMCangleLHReso->SetMarkerStyle(kFullCircle);
\r
541 fHistMCdecayAngleReso= new TH1F("fHistMCdecayAngleReso", "MC True Angle between decay products", 300, 0., 2*TMath::Pi());
\r
542 fHistMCdecayAngleReso->GetXaxis()->SetTitle("Angle (#Lambda-p#pi)");
\r
543 fHistMCdecayAngleReso->GetYaxis()->SetTitle("Entries");
\r
544 fHistMCdecayAngleReso->SetMarkerStyle(kFullCircle);
\r
546 fHistMCpointingAngleReso= new TH1F("fHistMCpointingAngleReso", "MC True Pointing Angle", 300, 0., 2*TMath::Pi());
\r
547 fHistMCpointingAngleReso->GetXaxis()->SetTitle("Pointing Angle");
\r
548 fHistMCpointingAngleReso->GetYaxis()->SetTitle("Entries");
\r
549 fHistMCpointingAngleReso->SetMarkerStyle(kFullCircle);
\r
551 fHistMCapReso = new TH2F("fHistMCapReso", "True MC Armenteros-Podolanski", 200,-1.0,1.0, 500,0,1);
\r
552 fHistMCapReso->GetYaxis()->SetTitle("#alpha");
\r
553 fHistMCapReso->GetXaxis()->SetTitle("q_{t}");
\r
554 fHistMCapReso->SetOption("scat");
\r
555 fHistMCapReso->SetMarkerStyle(kFullCircle);
\r
557 fHistList->Add(fHistMassDPi);
\r
558 fHistList->Add(fHistMassLPi);
\r
559 fHistList->Add(fHistMassLambdaPi);
\r
560 fHistList->Add(fHistMassLambda);
\r
561 fHistList->Add(fHistMassLambdaPPi);
\r
562 fHistList->Add(fHistMassLambdaP);
\r
563 fHistList->Add(fHistMassLambdaK);
\r
564 fHistList->Add(fHistMassK0onl);
\r
565 fHistList->Add(fHistMassK0offl);
\r
566 fHistList->Add(fHistMassK0onlC);
\r
567 fHistList->Add(fHistMassK0offlC);
\r
568 fHistList->Add(fHistMassPQonl);
\r
569 fHistList->Add(fHistMassPQoffl);
\r
570 fHistList->Add(fHistDC);
\r
571 fHistList->Add(fHistArmenterosPodolanski);
\r
572 fHistList->Add(fHistArmenterosPodolanskiCut);
\r
573 fHistList->Add(fHistHDibaryonInvaMassGen);
\r
574 fHistList->Add(fHistHDibaryonInvaMassGenRes);
\r
575 fHistList->Add(fHistAntiHDibaryonInvaMassGen);
\r
576 fHistList->Add(fHistHDibaryonInvaMassAso);
\r
577 fHistList->Add(fHistHDibaryonInvaMassAsoReso);
\r
578 fHistList->Add(fHistAntiHDibaryonInvaMassAso);
\r
579 fHistList->Add(fHistCheck);
\r
580 fHistList->Add(fHistHPointingAngle);
\r
581 fHistList->Add(fHistMassH);
\r
582 fHistList->Add(fHistMassPpi);
\r
583 fHistList->Add(fHistMassPpiReso);
\r
584 fHistList->Add(fHistMassLpi);
\r
585 fHistList->Add(fHistMassLP);
\r
586 fHistList->Add(fHistMassLambdaFromH);
\r
587 fHistList->Add(fHistMassLambdaFromHtLorentz);
\r
588 fHistList->Add(fHistProtonPIDBb);
\r
589 fHistList->Add(fHistPionPIDBb);
\r
590 fHistList->Add(fHistProtonPIDLambda);
\r
591 fHistList->Add(fHistPionPIDLambda);
\r
592 fHistList->Add(fHistMCdcaPvtxDvtx);
\r
593 fHistList->Add(fHistMCdcaPvtxLvtx);
\r
594 fHistList->Add(fHistMCdcaDvtxLvtx);
\r
595 fHistList->Add(fHistMCangleLH);
\r
596 fHistList->Add(fHistMCdecayAngle);
\r
597 fHistList->Add(fHistMCpointingAngle);
\r
598 fHistList->Add(fHistMCap);
\r
599 fHistList->Add(fHistMCdcaPvtxDvtxReso);
\r
600 fHistList->Add(fHistMCdcaPvtxLvtxReso);
\r
601 fHistList->Add(fHistMCdcaDvtxLvtxReso);
\r
602 fHistList->Add(fHistMCangleLHReso);
\r
603 fHistList->Add(fHistMCdecayAngleReso);
\r
604 fHistList->Add(fHistMCpointingAngleReso);
\r
605 fHistList->Add(fHistMCapReso);
\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
671 fHistList->Add(fHistPtvsEtaAso);
\r
673 fHistCount = new TH1F("fHistCount","test",17,0,17);
\r
674 fHistCount->GetXaxis()->SetBinLabel(1,"Events");
\r
675 fHistCount->GetXaxis()->SetBinLabel(2,"MC All");
\r
676 fHistCount->GetXaxis()->SetBinLabel(3,"MC from Primary Vtx");
\r
677 fHistCount->GetXaxis()->SetBinLabel(4,"Horst");
\r
678 fHistCount->GetXaxis()->SetBinLabel(5,"Lambda Candidates");
\r
679 fHistCount->GetXaxis()->SetBinLabel(6,"Sigma Candidates");
\r
680 fHistCount->GetXaxis()->SetBinLabel(7,"Horst");
\r
681 fHistCount->GetXaxis()->SetBinLabel(8,"Horst");
\r
682 fHistCount->GetXaxis()->SetBinLabel(9,"Horst");
\r
683 fHistCount->GetXaxis()->SetBinLabel(10,"MC All #bar{Lambda}(1520)s");
\r
684 fHistCount->GetXaxis()->SetBinLabel(11,"");
\r
685 fHistCount->GetXaxis()->SetBinLabel(12,"H-Dibaryon");
\r
686 fHistCount->GetXaxis()->SetBinLabel(13,"Hypertriton 2-Body");
\r
687 fHistCount->GetXaxis()->SetBinLabel(14,"Hypertriton 3-Body");
\r
688 fHistCount->GetXaxis()->SetBinLabel(15,"");
\r
689 fHistCount->GetXaxis()->SetBinLabel(16,"");
\r
690 fHistCount->GetXaxis()->SetBinLabel(17,"Lambdas!!!");
\r
691 fHistCount->SetStats(0);
\r
692 fHistCount->SetFillColor(38);
\r
693 fHistList->Add(fHistCount);
\r
696 //________________________________________________________________________
\r
697 void AliAnalysisTaskHdibaryonLPpi::UserExec(Option_t *)
\r
700 // Called for each event
\r
702 //define improtant masses
\r
703 Double_t pionMass = 0.13957;
\r
704 Double_t protonMass = 0.93827;
\r
707 Long_t pdgPionPlus = 211;
\r
708 Long_t pdgPionMinus = -211;
\r
709 Long_t pdgProton = 2212;
\r
710 Long_t pdgAntiProton = -2212;
\r
711 Long_t pdgLambda = 3122;
\r
712 Long_t pdgAntiLambda = -3122;
\r
713 Long_t pdgHDibaryon = 1020000020;
\r
714 Long_t pdgAntiHDibaryon = -1020000020;
\r
716 AliStack* stack(NULL);
\r
719 if(!fMCEvent)return;
\r
721 AliHeader *head = fMCEvent->Header();
\r
723 AliGenPythiaEventHeader *pyheader = (AliGenPythiaEventHeader*)head->GenEventHeader();
\r
724 if(!pyheader)return;
\r
726 if(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler()){
\r
727 if(static_cast<AliMCEventHandler*>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())->MCEvent()){
\r
728 if(!(static_cast<AliMCEventHandler*>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())->TreeK()))return;
\r
729 if(!(static_cast<AliMCEventHandler*>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())->TreeTR()))return;
\r
730 if(!(static_cast<AliMCEventHandler*>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())->Init("local")))return;
\r
731 stack = static_cast<AliMCEventHandler*>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())->MCEvent()->Stack();
\r
732 if(!(static_cast<AliMCEventHandler*>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())->MCEvent()->Stack()->TreeK()))return;
\r
739 // -------------------------------------------------------
\r
740 // Loop for Inv. Mass via ESD tracks
\r
741 // -------------------------------------------------------
\r
743 fHistCount->Fill(0);
\r
745 fESD = dynamic_cast<AliESDEvent *>(fInputEvent);
\r
747 const AliESDVertex *vertex = fESD->GetPrimaryVertexTracks();
\r
748 if (vertex->GetNContributors()<1)
\r
751 vertex = fESD->GetPrimaryVertexSPD();
\r
752 if(vertex->GetNContributors()<1) {
\r
756 if (TMath::Abs(vertex->GetZv()) > 10) return;
\r
758 Int_t centrality = 0;
\r
760 if (fESD->GetEventSpecie() == 4)
\r
762 AliCentrality *esdCentrality = fESD->GetCentrality();
\r
763 centrality = esdCentrality->GetCentralityClass10("V0M"); // centrality percentile determined with V0
\r
764 if (centrality > 8) return; //0 bis 80 %
\r
767 Int_t runNumber = 0;
\r
769 runNumber = fESD->GetRunNumber();
\r
771 if (!fPIDtpcESD) fPIDtpcESD = ((AliESDInputHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->GetESDpid();
\r
773 fPIDtpcESD = new AliESDpid(); // HACK FOR MC PBPB --> PLEASE REMOVE AS SOON AS POSSIBLE
\r
774 fPIDtpcESD->GetTPCResponse().SetBetheBlochParameters(1.28778e+00/50., 3.13539e+01, TMath::Exp(-3.16327e+01), 1.87901e+00, 6.41583e+00);
\r
780 TObjArray* listCrossV0 = fESDCutsV0->GetAcceptedV0s(fESD);
\r
781 Int_t nGoodV0s = listCrossV0->GetEntries();
\r
783 const AliESDVertex *esdVer = fESD->GetPrimaryVertex();
\r
784 AliESDVertex *esdVer1 = new AliESDVertex(*esdVer);
\r
786 AliVertexerTracks *vertexer = new AliVertexerTracks(fESD->GetMagneticField());
\r
787 TObjArray *trkArray = new TObjArray(2);
\r
788 AliVertexerTracks *vertexer1 = new AliVertexerTracks(fESD->GetMagneticField());
\r
789 TObjArray *trkArray1 = new TObjArray(2);
\r
791 AliKFParticle::SetField(fESD->GetMagneticField());
\r
793 AliKFVertex primVtx(*(fESD->GetPrimaryVertex()));
\r
795 Double_t mn[3] = {0,0,0};
\r
796 Double_t mp[3] = {0,0,0};
\r
797 Double_t mm[3] = {0,0,0};
\r
798 Double_t dd[3] = {0,0,0};
\r
799 Double_t dd1[3] = {0,0,0};
\r
800 const Double_t cProtonMass=TDatabasePDG::Instance()->GetParticle(2212)->Mass();
\r
801 const Double_t cPionMass=TDatabasePDG::Instance()->GetParticle(211)->Mass();
\r
802 const Double_t cElectronMass=TDatabasePDG::Instance()->GetParticle(11)->Mass();
\r
803 const Double_t cLambdaMass=TDatabasePDG::Instance()->GetParticle(3122)->Mass();
\r
804 Double_t decayLength=0;
\r
805 Double_t decayLengthH=0;
\r
807 //V0 Loop for Lambda and Anti-Lambda
\r
808 for(Int_t iV0MI = 0; iV0MI < nGoodV0s ; iV0MI++) {
\r
809 AliESDv0 * fV0MIs = fESD->GetV0(iV0MI);
\r
810 Int_t lOnFlyStatus = 0;
\r
812 lOnFlyStatus = fV0MIs->GetOnFlyStatus();
\r
813 Double_t lInvMassLambda=0;
\r
814 Double_t lInvMassLambdaPi=0;
\r
815 Double_t lPtLambda=0;
\r
816 Double_t lPzLambda=0;
\r
817 Double_t lPLambda=0;
\r
821 TLorentzVector posE;
\r
822 TLorentzVector negE;
\r
823 TLorentzVector photon;
\r
828 if (!lOnFlyStatus){
\r
832 fHistCheck->Fill(offl,onl);
\r
834 AliESDtrack* trackPosTest = fESD->GetTrack(fV0MIs->GetPindex());
\r
835 AliESDtrack* trackNegTest = fESD->GetTrack(fV0MIs->GetNindex());
\r
837 //PID via specific energy loss in the TPC
\r
838 //define the arrays for the Bethe-Bloch-Parameters
\r
839 Double_t parProton[5] = {0,0,0,0,0};
\r
841 if(runNumber < 166500) //LHC10h
\r
843 parProton[0] = 1.45802; // ALEPH parameters for protons (pass2)
\r
844 parProton[1] = 27.4992;
\r
845 parProton[2] = 4.00313e-15;
\r
846 parProton[3] = 2.48485;
\r
847 parProton[4] = 8.31768;
\r
850 if(runNumber > 166500) //LHC11h
\r
852 parProton[0] = 1.11243; // ALEPH parameters for protons (pass2)
\r
853 parProton[1] = 26.1144;
\r
854 parProton[2] = 4.00313e-15;
\r
855 parProton[3] = 2.72969;
\r
856 parProton[4] = 9.15038;
\r
859 //Get the total momentum for each track at the inner readout of the TPC
\r
860 Double_t ptotN = trackNegTest->GetInnerParam()->GetP();
\r
861 Double_t ptotP = trackPosTest->GetInnerParam()->GetP();
\r
863 //define expected signals for the various species
\r
864 Double_t expSignalPionP = 0;
\r
865 Double_t expSignalPionN = 0;
\r
866 Double_t expSignalProtonN = 0;
\r
867 Double_t expSignalProtonP = 0;
\r
872 expSignalProtonN = AliExternalTrackParam::BetheBlochAleph(ptotN/(protonMass),parProton[0],parProton[1],parProton[2],parProton[3],parProton[4]);
\r
873 expSignalProtonP = AliExternalTrackParam::BetheBlochAleph(ptotP/(protonMass),parProton[0],parProton[1],parProton[2],parProton[3],parProton[4]);
\r
879 expSignalPionP = 0.7*AliExternalTrackParam::BetheBlochAleph(ptotP/(pionMass),parProton[0],parProton[1],parProton[2],parProton[3],parProton[4]);
\r
880 expSignalPionN = 0.7*AliExternalTrackParam::BetheBlochAleph(ptotN/(pionMass),parProton[0],parProton[1],parProton[2],parProton[3],parProton[4]);
\r
882 expSignalProtonN = 0.65*AliExternalTrackParam::BetheBlochAleph(ptotN/(protonMass),parProton[0],parProton[1],parProton[2],parProton[3],parProton[4]);
\r
883 expSignalProtonP = 0.65*AliExternalTrackParam::BetheBlochAleph(ptotP/(protonMass),parProton[0],parProton[1],parProton[2],parProton[3],parProton[4]);
\r
886 // PID cut on the nuclei (proton, deuteron, triton, helium3)
\r
887 Bool_t corrParticle = kFALSE;
\r
889 Bool_t posProton = kFALSE;
\r
893 if (TMath::Abs(fPIDtpcESD->NumberOfSigmasTPC(trackPosTest, AliPID::kProton)) < 3)
\r
896 corrParticle = kTRUE;
\r
902 if(//trackPosTest->GetTPCsignal() < 1200 &&
\r
903 TMath::Abs(trackPosTest->GetTPCsignal() - expSignalProtonP)/expSignalProtonP < 0.4)
\r
906 corrParticle = kTRUE;
\r
910 Bool_t negProton = kFALSE;
\r
914 if (TMath::Abs(fPIDtpcESD->NumberOfSigmasTPC(trackNegTest, AliPID::kProton)) < 3)
\r
917 corrParticle = kTRUE;
\r
923 if(//trackNegTest->GetTPCsignal() < 1200 &&
\r
924 TMath::Abs(trackNegTest->GetTPCsignal() - expSignalProtonN)/expSignalProtonN < 0.4)
\r
927 corrParticle = kTRUE;
\r
931 //PID cut for pions
\r
932 //data: 3sigma cut on the pions
\r
934 Bool_t negPion = kFALSE;
\r
935 Bool_t posPion = kFALSE;
\r
939 if (TMath::Abs(fPIDtpcESD->NumberOfSigmasTPC(trackPosTest, AliPID::kPion)) < 4) posPion=kTRUE; //pos daughter has to be a pion
\r
940 if (TMath::Abs(fPIDtpcESD->NumberOfSigmasTPC(trackNegTest, AliPID::kPion)) < 4) negPion=kTRUE; // negative daughter has to be a pion
\r
943 //MC: like the nuclei via the specific energyloss in the TPC
\r
946 if (TMath::Abs(trackPosTest->GetTPCsignal() - expSignalPionP)/expSignalPionP < 0.4) posPion=kTRUE;
\r
947 if (TMath::Abs(trackNegTest->GetTPCsignal() - expSignalPionN)/expSignalPionN < 0.4) negPion=kTRUE;
\r
950 if (!(posProton==kTRUE)) continue;
\r
951 if (!(negPion==kTRUE)) continue;
\r
955 if( !(trackPosTest->GetStatus() & AliESDtrack::kTPCrefit)){
\r
959 if( !(trackNegTest->GetStatus() & AliESDtrack::kTPCrefit)){
\r
963 if( trackPosTest->GetSign() >0 && trackNegTest->GetSign() <0){
\r
964 fV0MIs->GetNPxPyPz(mn[0],mn[1],mn[2]); //reconstructed cartesian momentum components of negative daughter
\r
965 fV0MIs->GetPPxPyPz(mp[0],mp[1],mp[2]); //reconstructed cartesian momentum components of positive daughter
\r
968 if( trackPosTest->GetSign() <0 && trackNegTest->GetSign() >0){
\r
969 fV0MIs->GetPPxPyPz(mn[0],mn[1],mn[2]); //reconstructed cartesian momentum components of negative daughter
\r
970 fV0MIs->GetNPxPyPz(mp[0],mp[1],mp[2]); //reconstructed cartesian momentum components of positive daughter
\r
973 fV0MIs->GetPxPyPz(mm[0],mm[1],mm[2]); //reconstructed cartesian momentum components of mother
\r
975 TVector3 vecN(mn[0],mn[1],mn[2]);
\r
976 TVector3 vecP(mp[0],mp[1],mp[2]);
\r
977 TVector3 vecM(mm[0],mm[1],mm[2]);
\r
979 Double_t thetaP = acos((vecP * vecM)/(vecP.Mag() * vecM.Mag()));
\r
980 Double_t thetaN = acos((vecN * vecM)/(vecN.Mag() * vecM.Mag()));
\r
982 Double_t alfa = ((vecP.Mag())*cos(thetaP)-(vecN.Mag())*cos(thetaN))/
\r
983 ((vecP.Mag())*cos(thetaP)+(vecN.Mag())*cos(thetaN)) ;
\r
984 Double_t qt = vecP.Mag()*sin(thetaP);
\r
986 fHistArmenterosPodolanski->Fill(alfa,qt); //Armenteros-Podolanski calculation
\r
989 TLorentzVector k0daugh1;
\r
990 TLorentzVector k0daugh2;
\r
991 TLorentzVector proton;
\r
994 k0daugh1.SetXYZM(mn[0],mn[1],mn[2],cPionMass);
\r
995 k0daugh2.SetXYZM(mp[0],mp[1],mp[2],cPionMass);
\r
996 k0=k0daugh1+k0daugh2;
\r
998 fV0MIs->ChangeMassHypothesis(3122);
\r
999 lInvMassLambda = fV0MIs->GetEffMass();
\r
1000 lPtLambda = fV0MIs->Pt();
\r
1001 lPzLambda = fV0MIs->Pz();
\r
1002 lPLambda = fV0MIs->P();
\r
1004 trkArray->AddAt(trackPosTest,0);
\r
1005 trkArray->AddAt(trackNegTest,1);
\r
1007 vertexer->SetVtxStart(esdVer1);
\r
1008 AliESDVertex *decayVertex = (AliESDVertex*)vertexer->VertexForSelectedESDTracks(trkArray);
\r
1010 dd[0]=fESD->GetPrimaryVertexSPD()->GetX()-decayVertex->GetX();
\r
1011 dd[1]=fESD->GetPrimaryVertexSPD()->GetY()-decayVertex->GetY();
\r
1012 dd[2]=fESD->GetPrimaryVertexSPD()->GetZ()-decayVertex->GetZ();
\r
1014 decayLength=sqrt(dd[0]*dd[0]+dd[1]*dd[1]+dd[2]*dd[2]);
\r
1016 if (decayVertex) delete decayVertex;
\r
1018 TLorentzVector negPio1;
\r
1019 TLorentzVector posProt1;
\r
1020 TLorentzVector posP;
\r
1021 TLorentzVector posProt;
\r
1022 TLorentzVector negK;
\r
1023 TLorentzVector negPio;
\r
1024 TLorentzVector negPi;
\r
1025 TLorentzVector omega;
\r
1026 TLorentzVector threeSum;
\r
1027 TLorentzVector fourSum;
\r
1028 TLorentzVector ppK;
\r
1029 TLorentzVector posPiK;
\r
1030 TLorentzVector negPiK;
\r
1031 TLorentzVector kaon;
\r
1032 TLorentzVector lambda;
\r
1033 TLorentzVector lambdaH;
\r
1034 TLorentzVector hDibaryon;
\r
1038 h.SetXYZ(-dd[0],-dd[1],-dd[2]);
\r
1040 if (onl==1)fHistMassDPi->Fill(lInvMassLambda);
\r
1043 fHistMassLPi->Fill(lInvMassLambda);
\r
1045 negE.SetXYZM(mn[0],mn[1],mn[2],cElectronMass);
\r
1046 posE.SetXYZM(mp[0],mp[1],mp[2],cElectronMass);
\r
1049 negPiK.SetXYZM(mn[0],mn[1],mn[2],cPionMass);
\r
1050 posPiK.SetXYZM(mp[0],mp[1],mp[2],cPionMass);
\r
1051 kaon=posPiK+negPiK;
\r
1053 negPi.SetXYZM(mn[0],mn[1],mn[2],cPionMass);
\r
1054 posP.SetXYZM(mp[0],mp[1],mp[2],cProtonMass);
\r
1055 lambda=negPi+posP;
\r
1056 lambdaH.SetXYZM(mm[0],mm[1],mm[2],cLambdaMass);
\r
1058 if (lInvMassLambda>1.1113&&lInvMassLambda<1.1202){
\r
1061 if (qt<-2.21*alfa*alfa+2.945*alfa-0.887) continue;
\r
1062 if (qt>-2.21*alfa*alfa+2.945*alfa-0.873) continue;
\r
1063 if (photon.M()<0.005) continue;
\r
1064 if (kaon.M()>0.495 && kaon.M()<0.500 ) continue;
\r
1068 fHistMassLambda->Fill(lInvMassLambda);
\r
1070 Bool_t isCorrectlyAssociatedLambda = kFALSE;
\r
1071 Bool_t isPartialCorrectlyAssociatedLambda = kFALSE;
\r
1072 Int_t labelAssociatedH=-1;
\r
1073 Int_t labelLambda=-1;
\r
1076 Int_t labelPosTest = trackPosTest->GetLabel();
\r
1077 TParticle *tparticleDaughter = stack->Particle(TMath::Abs(labelPosTest));
\r
1078 Int_t labelMother = tparticleDaughter->GetFirstMother();
\r
1079 TParticle *tparticleMother = stack->Particle(TMath::Abs(labelMother));
\r
1081 Int_t labelOma = tparticleMother->GetFirstMother();
\r
1082 TParticle *tparticleOma = stack->Particle(TMath::Abs(labelOma));
\r
1084 if ((tparticleOma->GetPdgCode() < 0) && TMath::Abs(tparticleMother->GetPdgCode())==pdgLambda){// check mother to be Lambda
\r
1085 Int_t labelFirstDaughter = tparticleMother->GetDaughter(1);// Proton
\r
1086 Int_t labelThirdDaughter = tparticleMother->GetDaughter(0);// Pion
\r
1088 TParticle *tparticleFirstDaughter = stack->Particle(TMath::Abs(labelFirstDaughter));
\r
1089 TParticle *tparticleThirdDaughter = stack->Particle(TMath::Abs(labelThirdDaughter));
\r
1091 if((tparticleFirstDaughter->GetPdgCode() == pdgProton && tparticleThirdDaughter->GetPdgCode()== pdgPionMinus) ||
\r
1092 (tparticleFirstDaughter->GetPdgCode() == pdgPionMinus && tparticleThirdDaughter->GetPdgCode()== pdgProton)){ //daughter PDGs
\r
1093 isPartialCorrectlyAssociatedLambda = kTRUE;
\r
1094 labelLambda=labelMother;
\r
1099 if(tparticleOma->GetPdgCode() == pdgHDibaryon){ //check grandmother to be H PDG
\r
1100 if (TMath::Abs(tparticleMother->GetPdgCode())==pdgLambda){// check mother to be Lambda
\r
1101 Int_t labelFirstDaughter = tparticleMother->GetDaughter(1);// Proton
\r
1102 Int_t labelThirdDaughter = tparticleMother->GetDaughter(0);// Pion
\r
1104 TParticle *tparticleFirstDaughter = stack->Particle(TMath::Abs(labelFirstDaughter));
\r
1105 TParticle *tparticleThirdDaughter = stack->Particle(TMath::Abs(labelThirdDaughter));
\r
1107 if((tparticleFirstDaughter->GetPdgCode() == pdgProton && tparticleThirdDaughter->GetPdgCode()== pdgPionMinus) ||
\r
1108 (tparticleFirstDaughter->GetPdgCode() == pdgPionMinus && tparticleThirdDaughter->GetPdgCode()== pdgProton)){ //daughter PDGs
\r
1109 isCorrectlyAssociatedLambda = kTRUE;
\r
1110 labelAssociatedH=labelOma;
\r
1111 fHistMassLambdaFromH->Fill(lInvMassLambda);
\r
1112 fHistMassLambdaFromHtLorentz->Fill(lambda.M());
\r
1118 fHistProtonPIDLambda->Fill(trackPosTest->GetInnerParam()->GetP(), trackPosTest->GetTPCsignal());
\r
1119 fHistPionPIDLambda->Fill(trackNegTest->GetInnerParam()->GetP(), trackNegTest->GetTPCsignal());
\r
1121 //---------------------------------------------------------
\r
1122 // Proton track loop
\r
1123 //---------------------------------------------------------
\r
1124 fHistArmenterosPodolanskiCut->Fill(alfa,qt);
\r
1126 for (Int_t iTracksP = 0; iTracksP < fESD->GetNumberOfTracks(); iTracksP++) {
\r
1127 AliESDtrack* trackP = dynamic_cast<AliESDtrack*> (fESD->GetTrack(iTracksP));
\r
1128 if (trackP->GetSign()<0) continue;
\r
1130 if (iTracksP==fV0MIs->GetPindex())continue;
\r
1131 if (iTracksP==fV0MIs->GetNindex())continue;
\r
1133 if (!fEsdTrackCuts->AcceptTrack(trackP)) continue;
\r
1135 AliKFParticle protonKF( *(trackP), 2212);
\r
1137 if (!trackP->GetInnerParam()) continue;
\r
1144 if (TMath::Abs(fPIDtpcESD->NumberOfSigmasTPC(trackP, AliPID::kProton)) > 3) continue;
\r
1147 fHistProtonPIDBb->Fill(trackP->GetInnerParam()->GetP(), trackP->GetTPCsignal());
\r
1149 posProt.SetXYZM(trackP->Px(),trackP->Py(),trackP->Pz(),cProtonMass);
\r
1151 //!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
\r
1152 //Pion Track loop!!!!!!!!!!!!!!!!!!!!!
\r
1153 //!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
\r
1154 for (Int_t iTracksN = iTracksP+1; iTracksN < fESD->GetNumberOfTracks(); iTracksN++) {
\r
1155 AliESDtrack* trackN = dynamic_cast<AliESDtrack*> (fESD->GetTrack(iTracksN));
\r
1157 if (iTracksN==fV0MIs->GetPindex())continue;
\r
1158 if (iTracksN==fV0MIs->GetNindex())continue;
\r
1159 if (trackN->GetSign()>0) continue;
\r
1161 if (!fEsdTrackCuts->AcceptTrack(trackN)) continue;
\r
1163 negPi.SetXYZM(mn[0],mn[1],mn[2],cPionMass);
\r
1164 posP.SetXYZM(mp[0],mp[1],mp[2],cProtonMass);
\r
1165 negPio.SetXYZM(trackN->Px(),trackN->Py(),trackN->Pz(),cPionMass);
\r
1167 threeSum=negPi+posP+negPio;
\r
1168 lInvMassLambdaPi=threeSum.M();
\r
1170 fHistDC->Fill(decayLength*lInvMassLambda/lPLambda);
\r
1172 AliKFParticle posPionKF( *(trackN) ,-211);
\r
1174 if (!trackN->GetInnerParam()) continue;
\r
1180 if (TMath::Abs(fPIDtpcESD->NumberOfSigmasTPC(trackN, AliPID::kPion)) > 3) continue;
\r
1182 fHistPionPIDBb->Fill(trackN->GetInnerParam()->GetP(), trackN->GetTPCsignal());
\r
1184 trkArray1->AddAt(trackP,0);
\r
1185 trkArray1->AddAt(trackN,1);
\r
1187 vertexer1->SetVtxStart(esdVer1);
\r
1188 AliESDVertex *decayVertex1 = (AliESDVertex*)vertexer1->VertexForSelectedESDTracks(trkArray1);
\r
1190 dd1[0]=fESD->GetPrimaryVertexSPD()->GetX()-decayVertex1->GetX();
\r
1191 dd1[1]=fESD->GetPrimaryVertexSPD()->GetY()-decayVertex1->GetY();
\r
1192 dd1[2]=fESD->GetPrimaryVertexSPD()->GetZ()-decayVertex1->GetZ();
\r
1194 decayLengthH=sqrt(dd1[0]*dd1[0]+dd1[1]*dd1[1]+dd1[2]*dd1[2]);
\r
1196 Double_t bz = fESD->GetMagneticField();
\r
1198 trackP->PropagateToDCA(decayVertex1, bz, 10);
\r
1199 trackN->PropagateToDCA(decayVertex1, bz, 10);
\r
1201 Double_t xthiss(0.0);
\r
1202 Double_t xpp(0.0);
\r
1203 Double_t dca = trackN->GetDCA(trackP,bz,xthiss,xpp);
\r
1205 if (decayVertex1) delete decayVertex1;
\r
1206 h1.SetXYZ(-dd1[0],-dd1[1],-dd1[2]);
\r
1208 if (dca>1) continue;
\r
1210 fourSum=threeSum+posProt;
\r
1212 posProt1.SetXYZM(trackP->Px(),trackP->Py(),trackP->Pz(),cProtonMass);
\r
1213 negPio1.SetXYZM(trackN->Px(),trackN->Py(),trackN->Pz(),cPionMass);
\r
1214 hDibaryon=lambdaH+posProt1+negPio1;
\r
1215 Double_t hPointingAngle = hDibaryon.Angle(h);
\r
1216 Double_t pointingAngleH = hDibaryon.Angle(h1);
\r
1217 Double_t decayAngleH = h.Angle(h1);
\r
1218 TVector3 vecDist(dd[0]-dd1[0],dd[1]-dd1[1],dd[2]-dd1[2]);
\r
1219 fHistMassLambdaPPi->Fill(hDibaryon.M());
\r
1220 fHistHPointingAngle->Fill(pointingAngleH);
\r
1222 Double_t rapidity = hDibaryon.Rapidity();
\r
1223 if(rapidity > 1.0 || rapidity < -1.0) continue;
\r
1225 fHistRap->Fill(rapidity);
\r
1226 if (pointingAngleH > 0.1) continue;
\r
1228 ///////////////////////////
\r
1229 //MC part for Associated H
\r
1230 ///////////////////////////
\r
1232 if (HasMC() && isCorrectlyAssociatedLambda) {
\r
1233 Int_t labelP = trackP->GetLabel();
\r
1234 TParticle *tparticleDaughter = stack->Particle(TMath::Abs(labelP));
\r
1235 Int_t labelMother = tparticleDaughter->GetFirstMother();
\r
1236 TParticle *tparticleMother = stack->Particle(TMath::Abs(labelMother));
\r
1238 Int_t labelProton = trackP->GetLabel();
\r
1239 Int_t labelPion = trackN->GetLabel();
\r
1242 if(tparticleMother->GetPdgCode() == pdgHDibaryon && labelAssociatedH==labelMother){ //check mother PDG
\r
1243 Int_t labelFirstDaughter = tparticleMother->GetDaughter(0);
\r
1244 Int_t labelThirdDaughter = tparticleMother->GetDaughter(1);
\r
1245 Int_t labelSecondDaughter = labelFirstDaughter +1;
\r
1247 TParticle *tparticleSecondDaughter = stack->Particle(TMath::Abs(labelSecondDaughter));
\r
1248 TParticle *tparticleThirdDaughter = stack->Particle(TMath::Abs(labelThirdDaughter));
\r
1250 TLorentzVector ppi;
\r
1251 TLorentzVector lpi;
\r
1252 TLorentzVector lP;
\r
1254 ppi=posProt+negPio;
\r
1255 lpi=lambdaH+negPio;
\r
1256 lP=lambdaH+posProt;
\r
1258 if((tparticleThirdDaughter->GetPdgCode() == pdgPionMinus && labelPion==labelThirdDaughter)&&(tparticleSecondDaughter->GetPdgCode() == pdgProton||tparticleSecondDaughter->GetPdgCode() == pdgPionMinus) && labelProton==labelSecondDaughter) fHistMassPpi->Fill(ppi.M());
\r
1260 if(tparticleThirdDaughter->GetPdgCode() == pdgPionMinus && labelPion==labelThirdDaughter) fHistMassLpi->Fill(lpi.M());
\r
1262 if(tparticleSecondDaughter->GetPdgCode() == pdgProton && labelProton==labelSecondDaughter) fHistMassLP->Fill(lP.M());
\r
1264 if(tparticleSecondDaughter->GetPdgCode() == pdgProton && labelProton==labelSecondDaughter){//check second daughter PDG
\r
1265 if(tparticleThirdDaughter->GetPdgCode() == pdgPionMinus && labelPion==labelThirdDaughter){//check second daughter PDG
\r
1267 fHistHDibaryonInvaMassAso->Fill(hDibaryon.M());
\r
1269 Double_t distance01=vecDist.Mag();
\r
1270 fHistMCdcaPvtxDvtx->Fill(decayLengthH);
\r
1271 fHistMCdcaPvtxLvtx->Fill(decayLength);
\r
1272 fHistMCdcaDvtxLvtx->Fill(distance01);
\r
1273 fHistMCangleLH->Fill(hPointingAngle);
\r
1274 fHistMCdecayAngle->Fill(decayAngleH);
\r
1275 fHistMCpointingAngle->Fill(pointingAngleH);
\r
1276 fHistMCap->Fill(alfa,qt);
\r
1278 fHistHilf1->Fill(posPionKF.GetDistanceFromVertex(primVtx));
\r
1279 fHistHilf2->Fill(protonKF.GetDistanceFromVertex(primVtx));
\r
1280 fHistHilf3->Fill(protonKF.GetDistanceFromVertex(posPionKF));
\r
1281 fHistHilf6->Fill(dca);
\r
1282 fHistPtvsYAso->Fill(hDibaryon.Pt(),hDibaryon.Rapidity());
\r
1283 fHistPtvsEtaAso->Fill(hDibaryon.Pt(),hDibaryon.Eta());
\r
1284 }//end check for third daughter PDG
\r
1285 }//end check second daughter PDG
\r
1289 fHistMassH->Fill(hDibaryon.M());
\r
1290 ppK=lambdaH+posProt;
\r
1291 fHistMassLambdaP->Fill(ppK.M());
\r
1298 //!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
\r
1299 //Pure MC Part!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
\r
1300 //!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
\r
1302 // Monte Carlo for genenerated particles
\r
1303 if (HasMC()) //MC loop
\r
1308 Double_t momentumPionGen[3]={0,0,0};
\r
1309 Double_t momentumNucleonGen[3]={0,0,0};
\r
1310 Double_t momentumLambdaGen[3]={0,0,0};
\r
1312 Double_t energyPionGen = 0;
\r
1313 Double_t energyNucleonGen = 0;
\r
1314 Double_t energyLambdaGen = 0;
\r
1316 Double_t transversMomentumMotherGen = 0;
\r
1317 Double_t longitudinalMomentumMotherGen = 0;
\r
1318 Double_t totalEnergyMotherGen = 0;
\r
1320 Double_t rapidityGen = 2;
\r
1322 for(stackN = 0; stackN < stack->GetNtrack(); stackN++) //loop over stack
\r
1325 TParticle *tparticleMother = stack->Particle(stackN);
\r
1327 if(tparticleMother->GetPdgCode() == pdgLambda) fHistCount->Fill(16);
\r
1330 if(tparticleMother->GetPdgCode() == pdgHDibaryon) //check mother PDG
\r
1332 Int_t labelFirstDaughter = tparticleMother->GetDaughter(0);
\r
1333 Int_t labelThirdDaughter = tparticleMother->GetDaughter(1);
\r
1334 Int_t labelSecondDaughter = labelFirstDaughter +1;
\r
1336 TParticle *tparticleFirstDaughter = stack->Particle(TMath::Abs(labelFirstDaughter));
\r
1337 TParticle *tparticleSecondDaughter = stack->Particle(TMath::Abs(labelSecondDaughter));
\r
1338 TParticle *tparticleThirdDaughter = stack->Particle(TMath::Abs(labelThirdDaughter));
\r
1340 if(tparticleFirstDaughter->GetPdgCode() == pdgLambda) //check first daughter PDG
\r
1342 if(tparticleSecondDaughter->GetPdgCode() == pdgProton)//check second daughter PDG
\r
1344 if(tparticleThirdDaughter->GetPdgCode() == pdgPionMinus)//check second daughter PDG
\r
1346 momentumLambdaGen[0] = tparticleFirstDaughter->Px();
\r
1347 momentumLambdaGen[1] = tparticleFirstDaughter->Py();
\r
1348 momentumLambdaGen[2] = tparticleFirstDaughter->Pz();
\r
1350 momentumNucleonGen[0] = tparticleSecondDaughter->Px();
\r
1351 momentumNucleonGen[1] = tparticleSecondDaughter->Py();
\r
1352 momentumNucleonGen[2] = tparticleSecondDaughter->Pz();
\r
1354 momentumPionGen[0] = tparticleThirdDaughter->Px();
\r
1355 momentumPionGen[1] = tparticleThirdDaughter->Py();
\r
1356 momentumPionGen[2] = tparticleThirdDaughter->Pz();
\r
1358 TLorentzVector lorentzVectorLambda;
\r
1359 TLorentzVector lorentzVectorProton;
\r
1360 TLorentzVector lorentzVectorPion;
\r
1361 TLorentzVector lorentzVectorHDibaryon;
\r
1363 lorentzVectorLambda.SetXYZM(momentumLambdaGen[0],momentumLambdaGen[1],momentumLambdaGen[2],1.115);
\r
1364 lorentzVectorProton.SetXYZM(momentumNucleonGen[0],momentumNucleonGen[1],momentumNucleonGen[2],protonMass);
\r
1365 lorentzVectorPion.SetXYZM(momentumPionGen[0],momentumPionGen[1],momentumPionGen[2],pionMass);
\r
1367 lorentzVectorHDibaryon = lorentzVectorLambda + lorentzVectorProton + lorentzVectorPion;
\r
1368 rapidityGen=lorentzVectorHDibaryon.Rapidity();
\r
1369 transversMomentumMotherGen = lorentzVectorHDibaryon.Pt();
\r
1370 longitudinalMomentumMotherGen = lorentzVectorHDibaryon.Pz();
\r
1371 totalEnergyMotherGen = lorentzVectorHDibaryon.Energy();
\r
1373 if(rapidityGen > 1.0 || rapidityGen < -1 ) continue;
\r
1374 //lorentzVectorLambda
\r
1375 fHistHDibaryonInvaMassGen->Fill(lorentzVectorHDibaryon.M());
\r
1376 if (lorentzVectorLambda.Rapidity() > 1.0 || lorentzVectorLambda.Rapidity() < -1) continue;
\r
1377 if (lorentzVectorProton.Rapidity() > 1.0 || lorentzVectorProton.Rapidity() < -1) continue;
\r
1379 if (lorentzVectorPion.Rapidity() > 1.0 || lorentzVectorPion.Rapidity() < -1) continue;
\r
1380 fHistHDibaryonInvaMassGenRes->Fill(lorentzVectorHDibaryon.M());
\r
1381 fHistPtvsEtaGen->Fill(lorentzVectorHDibaryon.Pt(),lorentzVectorHDibaryon.Eta());
\r
1382 fHistPtvsYGen->Fill(lorentzVectorHDibaryon.Pt(),lorentzVectorHDibaryon.Rapidity());
\r
1383 fHistPtvsEtaGen->Fill(lorentzVectorHDibaryon.Pt(),lorentzVectorHDibaryon.Eta());
\r
1384 fHistCount->Fill(11);
\r
1385 }//end of check third daughter PDG
\r
1386 }//end of check second daughter PDG
\r
1387 }//end of check first daughter PDG
\r
1388 }//end of H-Dibaryon
\r
1391 if(tparticleMother->GetPdgCode() == pdgAntiHDibaryon) //check mother PDG
\r
1393 Int_t labelFirstDaughter = tparticleMother->GetDaughter(0);
\r
1394 Int_t labelThirdDaughter = tparticleMother->GetDaughter(1);
\r
1395 Int_t labelSecondDaughter = labelFirstDaughter +1;
\r
1397 TParticle *tparticleFirstDaughter = stack->Particle(TMath::Abs(labelFirstDaughter));
\r
1398 TParticle *tparticleSecondDaughter = stack->Particle(TMath::Abs(labelSecondDaughter));
\r
1399 TParticle *tparticleThirdDaughter = stack->Particle(TMath::Abs(labelThirdDaughter));
\r
1401 if(tparticleFirstDaughter->GetPdgCode() == pdgAntiLambda) //check first daughter PDG
\r
1403 if(tparticleSecondDaughter->GetPdgCode() == pdgAntiProton)//check second daughter PDG
\r
1405 if(tparticleThirdDaughter->GetPdgCode() == pdgPionPlus)//check second daughter PDG
\r
1407 momentumLambdaGen[0] = tparticleFirstDaughter->Px();
\r
1408 momentumLambdaGen[1] = tparticleFirstDaughter->Py();
\r
1409 momentumLambdaGen[2] = tparticleFirstDaughter->Pz();
\r
1411 momentumNucleonGen[0] = tparticleSecondDaughter->Px();
\r
1412 momentumNucleonGen[1] = tparticleSecondDaughter->Py();
\r
1413 momentumNucleonGen[2] = tparticleSecondDaughter->Pz();
\r
1415 momentumPionGen[0] = tparticleThirdDaughter->Px();
\r
1416 momentumPionGen[1] = tparticleThirdDaughter->Py();
\r
1417 momentumPionGen[2] = tparticleThirdDaughter->Pz();
\r
1419 energyLambdaGen = tparticleFirstDaughter->Energy();
\r
1420 energyNucleonGen = tparticleSecondDaughter->Energy();
\r
1421 energyPionGen = tparticleThirdDaughter->Energy();
\r
1423 TLorentzVector lorentzVectorLambda;
\r
1424 TLorentzVector lorentzVectorProton;
\r
1425 TLorentzVector lorentzVectorPion;
\r
1426 TLorentzVector lorentzVectorHDibaryon;
\r
1428 lorentzVectorLambda.SetXYZM(momentumLambdaGen[0],momentumLambdaGen[1],momentumLambdaGen[2],1.115);
\r
1429 lorentzVectorProton.SetXYZM(momentumNucleonGen[0],momentumNucleonGen[1],momentumNucleonGen[2],protonMass);
\r
1430 lorentzVectorPion.SetXYZM(momentumPionGen[0],momentumPionGen[1],momentumPionGen[2],pionMass);
\r
1432 lorentzVectorHDibaryon = lorentzVectorLambda + lorentzVectorProton + lorentzVectorPion;
\r
1434 rapidityGen=lorentzVectorHDibaryon.Rapidity();
\r
1435 if(rapidityGen > 1.0 || rapidityGen < -1 ) continue;
\r
1436 fHistAntiHDibaryonInvaMassGen->Fill(lorentzVectorHDibaryon.M());
\r
1437 }//end of check third daughter PDG
\r
1438 }//end of check second daughter PDG
\r
1439 }//end of check first daughter PDG
\r
1440 }//end of Anti-H-Dibaryon
\r
1444 // Post output data.
\r
1445 PostData(1,fHistList);
\r
1446 //PostData(0,fHistList);
\r
1448 if (listCrossV0) delete listCrossV0;
\r
1449 if (esdVer1) delete esdVer1;
\r
1450 if (vertexer) delete vertexer;
\r
1451 if (vertexer1) delete vertexer1;
\r
1452 if (trkArray) delete trkArray;
\r
1453 if (trkArray1) delete trkArray1;
\r
1456 //________________________________________________________________________
\r
1457 void AliAnalysisTaskHdibaryonLPpi::Terminate(Option_t *)
\r
1459 // Draw result to the screen
\r
1460 // Called once at the end of the query
\r