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);
256 fESDtrackCutsV0->SetPtRange(0.2,1.5);
\r
257 fESDtrackCutsV0->SetMinDCAToVertexXY(2);
\r //war inzwischen 1 & 3
258 fESDtrackCutsV0->SetMinDCAToVertexZ(2);
\r //war inzwischen 1 & 3
260 fESDCutsV0 = new AliESDv0Cuts("AliESDCutsV0","AliESDCutsV0");
\r
261 fESDCutsV0->SetMaxDcaV0Daughters(1.0);
262 fESDCutsV0->SetMinDcaNegToVertex(2);
\r //1.5
263 fESDCutsV0->SetMinDcaPosToVertex(2); //1.5
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
296 fHistList->SetOwner();
\r
298 fHistMassDPi = new TH1F("fHistMassDPi", "Invariant mass distribution p+#pi^{-} ", 500, 1.0, 1.25);
\r
299 fHistMassDPi->GetXaxis()->SetTitle("Invariant mass p+#pi^{-} (GeV/c^{2})");
\r
300 fHistMassDPi->GetYaxis()->SetTitle("Entries");
\r
301 fHistMassDPi->SetMarkerStyle(kFullCircle);
\r
303 fHistMassLPi = new TH1F("fHistMassLPi", "Offline Invariant mass distribution p+#pi^{-} ", 500, 1.0, 1.25);
\r
304 fHistMassLPi->GetXaxis()->SetTitle("Invariant mass p+#pi^{-} (GeV/c^{2})");
\r
305 fHistMassLPi->GetYaxis()->SetTitle("Entries");
\r
306 fHistMassLPi->SetMarkerStyle(kFullCircle);
\r
308 fHistMassLambdaPi = new TH1F("fHistMassLambdaPi", "Invariant mass distribution #Lambda+#pi^{-} ", 500, 1.2, 1.5);
\r
309 fHistMassLambdaPi->GetXaxis()->SetTitle("Invariant mass #Lambda+#pi^{-} (GeV/c^{2})");
\r
310 fHistMassLambdaPi->GetYaxis()->SetTitle("Entries");
\r
311 fHistMassLambdaPi->SetMarkerStyle(kFullCircle);
\r
313 fHistMassLambda = new TH1F("fHistMassLambda", "Invariant mass distribution of #Lambda for further analyis", 500, 1.0, 1.2);
\r
314 fHistMassLambda->GetXaxis()->SetTitle("Invariant mass p+#pi^{+} (GeV/c^{2})");
\r
315 fHistMassLambda->GetYaxis()->SetTitle("Entries");
\r
316 fHistMassLambda->SetMarkerStyle(kFullCircle);
\r
318 fHistMassLambdaPPi = new TH1F("fHistMassLambdaPPi", "Invariant mass distribution #Lambdap#pi^{-} ", 300, 2.2, 2.5);
\r
319 fHistMassLambdaPPi->GetXaxis()->SetTitle("Invariant mass #Lambdap#pi^{-} (GeV/c^{2})");
\r
320 fHistMassLambdaPPi->GetYaxis()->SetTitle("Entries");
\r
321 fHistMassLambdaPPi->SetMarkerStyle(kFullCircle);
\r
323 fHistMassLambdaP = new TH1F("fHistMassLambdaP", "Invariant mass distribution #Lambdap ", 300, 2.2, 2.5);
\r
324 fHistMassLambdaP->GetXaxis()->SetTitle("Invariant mass #Lambdap (GeV/c^{2})");
\r
325 fHistMassLambdaP->GetYaxis()->SetTitle("Entries");
\r
326 fHistMassLambdaP->SetMarkerStyle(kFullCircle);
\r
328 fHistMassLambdaK = new TH1F("fHistMassLambdaK", "Invariant mass distribution #LambdaK^{-} ", 300, 1.6, 1.9);
\r
329 fHistMassLambdaK->GetXaxis()->SetTitle("Invariant mass #LambdaK^{-} (GeV/c^{2})");
\r
330 fHistMassLambdaK->GetYaxis()->SetTitle("Entries");
\r
331 fHistMassLambdaK->SetMarkerStyle(kFullCircle);
\r
333 fHistMassK0onl = new TH1F("fHistMassK0onl", "Invariant mass distribution K_{s}^{0} online V0 finder ", 400, 0.2, 1.0);
\r
334 fHistMassK0onl->GetXaxis()->SetTitle("Invariant mass #pi^{+}#pi^{-} (GeV/c^{2})");
\r
335 fHistMassK0onl->GetYaxis()->SetTitle("Entries");
\r
336 fHistMassK0onl->SetMarkerStyle(kFullCircle);
\r
338 fHistMassK0offl = new TH1F("fHistMassK0offl", "Invariant mass distribution K_{s}^{0} offline V0 finder ", 400, 0.2, 1.0);
\r
339 fHistMassK0offl->GetXaxis()->SetTitle("Invariant mass #pi^{+}#pi^{-} (GeV/c^{2})");
\r
340 fHistMassK0offl->GetYaxis()->SetTitle("Entries");
\r
341 fHistMassK0offl->SetMarkerStyle(kFullCircle);
\r
343 fHistMassK0onlC = new TH1F("fHistMassK0onlC", "Invariant mass distribution K_{s}^{0} online V0 finder ", 400, 0.2, 1.0);
\r
344 fHistMassK0onlC->GetXaxis()->SetTitle("Invariant mass #pi^{+}#pi^{-} (GeV/c^{2})");
\r
345 fHistMassK0onlC->GetYaxis()->SetTitle("Entries");
\r
346 fHistMassK0onlC->SetMarkerStyle(kFullCircle);
\r
348 fHistMassK0offlC = new TH1F("fHistMassK0offlC", "Invariant mass distribution K_{s}^{0} offline V0 finder ", 400, 0.2, 1.0);
\r
349 fHistMassK0offlC->GetXaxis()->SetTitle("Invariant mass #pi^{+}#pi^{-} (GeV/c^{2})");
\r
350 fHistMassK0offlC->GetYaxis()->SetTitle("Entries");
\r
351 fHistMassK0offlC->SetMarkerStyle(kFullCircle);
\r
353 fHistMassPQonl = new TH1F("fHistMassPQonl", "Invariant mass distribution K_{s}^{0}p using online V0 finder ", 500, 1.3, 2.3);
\r
354 fHistMassPQonl->GetXaxis()->SetTitle("Invariant mass K_{s}^{0}p (GeV/c^{2})");
\r
355 fHistMassPQonl->GetYaxis()->SetTitle("Entries");
\r
356 fHistMassPQonl->SetMarkerStyle(kFullCircle);
\r
358 fHistMassPQoffl = new TH1F("fHistMassPQoffl", "Invariant mass distribution K_{s}^{0}p using offline V0 finder ", 500, 1.3, 2.3);
\r
359 fHistMassPQoffl->GetXaxis()->SetTitle("Invariant mass K_{s}^{0}p (GeV/c^{2})");
\r
360 fHistMassPQoffl->GetYaxis()->SetTitle("Entries");
\r
361 fHistMassPQoffl->SetMarkerStyle(kFullCircle);
\r
363 fHistDC = new TH1F("fHistDC", "Proper decay length", 500, 0.0, 25);
\r
364 fHistDC->GetXaxis()->SetTitle("c#tau (cm)");
\r
365 fHistDC->GetYaxis()->SetTitle("Entries");
\r
366 fHistDC->SetMarkerStyle(kFullCircle);
\r
368 fHistArmenterosPodolanski = new TH2F("fHistArmenterosPodolanski", "Armenteros-Podolanski plot", 200,-1.0,1.0, 500,0,1);
\r
369 fHistArmenterosPodolanski->GetXaxis()->SetTitle("#alpha");
\r
370 fHistArmenterosPodolanski->GetYaxis()->SetTitle("q_{t}");
\r
371 fHistArmenterosPodolanski->SetMarkerStyle(kFullCircle);
\r
373 fHistArmenterosPodolanskiCut = new TH2F("fHistArmenterosPodolanskiCut", "Armenteros-Podolanski plot after cut", 200,-1.0,1.0, 500,0,1);
\r
374 fHistArmenterosPodolanskiCut->GetXaxis()->SetTitle("#alpha");
\r
375 fHistArmenterosPodolanskiCut->GetYaxis()->SetTitle("q_{t}");
\r
376 fHistArmenterosPodolanskiCut->SetMarkerStyle(kFullCircle);
\r
378 fHistHDibaryonInvaMassGen = new TH1F("fHistHDibaryonInvaMassGen", "Generated #Lambda p #pi^{-}", 200, 2.1, 2.3);
\r
379 fHistHDibaryonInvaMassGen->GetYaxis()->SetTitle("Counts");
\r
380 fHistHDibaryonInvaMassGen->GetXaxis()->SetTitle("Invariant mass (GeV/c^{2})");
\r
382 fHistHDibaryonInvaMassGenRes = new TH1F("fHistHDibaryonInvaMassGenRes", "Generated #Lambda p #pi^{-} with particles in rapidity!!!", 200, 2.1, 2.3);
\r
383 fHistHDibaryonInvaMassGenRes->GetYaxis()->SetTitle("Counts");
\r
384 fHistHDibaryonInvaMassGenRes->GetXaxis()->SetTitle("Invariant mass (GeV/c^{2})");
\r
386 fHistAntiHDibaryonInvaMassGen = new TH1F("fHistAntiHDibaryonInvaMassGen", "Generated #bar{#Lambda} #bar{p} #pi^{+}", 200, 2.1, 2.3);
\r
387 fHistAntiHDibaryonInvaMassGen->GetYaxis()->SetTitle("Counts");
\r
388 fHistAntiHDibaryonInvaMassGen->GetXaxis()->SetTitle("Invariant mass (GeV/c^{2})");
\r
390 fHistHDibaryonInvaMassAso = new TH1F("fHistHDibaryonInvaMassAso", "Associated #Lambda p #pi^{-}", 200, 2.1, 2.3);
\r
391 fHistHDibaryonInvaMassAso->GetYaxis()->SetTitle("Counts");
\r
392 fHistHDibaryonInvaMassAso->GetXaxis()->SetTitle("Invariant mass (GeV/c^{2})");
\r
394 fHistHDibaryonInvaMassAsoReso = new TH1F("fHistHDibaryonInvaMassAsoReso", "Associated #Lambda p #pi^{-}", 200, 2.1, 2.3);
\r
395 fHistHDibaryonInvaMassAsoReso->GetYaxis()->SetTitle("Counts");
\r
396 fHistHDibaryonInvaMassAsoReso->GetXaxis()->SetTitle("Invariant mass (GeV/c^{2})");
\r
398 fHistAntiHDibaryonInvaMassAso = new TH1F("fHistAntiHDibaryonInvaMassAso", "Associated #bar{#Lambda} #bar{p} #pi^{+}", 200, 2.1, 2.3);
\r
399 fHistAntiHDibaryonInvaMassAso->GetYaxis()->SetTitle("Counts");
\r
400 fHistAntiHDibaryonInvaMassAso->GetXaxis()->SetTitle("Invariant mass (GeV/c^{2})");
\r
402 fHistCheck = new TH2F("fHistCheck", "Check online/offline", 200, -0.5, 1.5, 200, -0.5, 1.5);
\r
403 fHistCheck->GetXaxis()->SetTitle("offline");
\r
404 fHistCheck->GetYaxis()->SetTitle("online");
\r
405 fHistCheck->SetMarkerStyle(kFullCircle);
\r
407 fHistHPointingAngle= new TH1F("fHistHPointingAngle", "Pointing angle distribution for #Lambdap#pi^{-}", 200, 0., 2*TMath::Pi());
\r
408 fHistHPointingAngle->GetXaxis()->SetTitle("Pointing angle distribution for #Lambdap#pi^{-}");
\r
409 fHistHPointingAngle->GetYaxis()->SetTitle("Entries");
\r
410 fHistHPointingAngle->SetMarkerStyle(kFullCircle);
\r
412 fHistMassH= new TH1F("fHistMassH", "Invariant mass distribution #Lambdap#pi^{-} (GeV/c^{2}) after pointing angle cut", 3000, 2.2, 2.5);
\r
413 fHistMassH->GetXaxis()->SetTitle("Invariant mass #Lambdap#pi^{-} (GeV/c^{2})");
\r
414 fHistMassH->GetYaxis()->SetTitle("Entries");
\r
415 fHistMassH->SetMarkerStyle(kFullCircle);
\r
417 fHistMassLambdaFromH= new TH1F("fHistMassLambdaFromH", "Invariant mass distribution #Lambda (GeV/c^{2}) asking for Mother to be H", 300, 1.0, 1.3);
\r
418 fHistMassLambdaFromH->GetXaxis()->SetTitle("Invariant mass #Lambda (GeV/c^{2})");
\r
419 fHistMassLambdaFromH->GetYaxis()->SetTitle("Entries");
\r
420 fHistMassLambdaFromH->SetMarkerStyle(kFullCircle);
\r
422 fHistMassLambdaFromHtLorentz= new TH1F("fHistMassLambdaFromHtLorentz", "Invariant mass distribution #Lambda (GeV/c^{2}) asking for Mother to be H", 300, 1.0, 1.3);
\r
423 fHistMassLambdaFromHtLorentz->GetXaxis()->SetTitle("Invariant mass #Lambda (GeV/c^{2})");
\r
424 fHistMassLambdaFromHtLorentz->GetYaxis()->SetTitle("Entries");
\r
425 fHistMassLambdaFromHtLorentz->SetMarkerStyle(kFullCircle);
\r
427 fHistMassPpi= new TH1F("fHistMassPpi", "Invariant mass distribution of the p#pi^{-} used for combing with Lambda", 300, 1.0, 1.6);
\r
428 fHistMassPpi->GetXaxis()->SetTitle("Invariant mass p#pi^{-} (GeV/c^{2})");
\r
429 fHistMassPpi->GetYaxis()->SetTitle("Entries");
\r
430 fHistMassPpi->SetMarkerStyle(kFullCircle);
\r
432 fHistMassPpiReso= new TH1F("fHistMassPpiReso", "Invariant mass distribution of the p#pi^{-} used for combing with Lambda", 300, 1.0, 1.6);
\r
433 fHistMassPpiReso->GetXaxis()->SetTitle("Invariant mass p#pi^{-} (GeV/c^{2})");
\r
434 fHistMassPpiReso->GetYaxis()->SetTitle("Entries");
\r
435 fHistMassPpiReso->SetMarkerStyle(kFullCircle);
\r
437 fHistMassLpi= new TH1F("fHistMassLpi", "Invariant mass distribution of the #Lambda#pi^{-} used for combing with p", 300, 1.1, 1.7);
\r
438 fHistMassLpi->GetXaxis()->SetTitle("Invariant mass #Lambda#pi^{-} (GeV/c^{2})");
\r
439 fHistMassLpi->GetYaxis()->SetTitle("Entries");
\r
440 fHistMassLpi->SetMarkerStyle(kFullCircle);
\r
442 fHistMassLP= new TH1F("fHistMassLP", "Invariant mass distribution of the #Lambda p used for combing with #pi^{-}", 300, 2.0, 2.3);
\r
443 fHistMassLP->GetXaxis()->SetTitle("Invariant mass #Lambda p (GeV/c^{2})");
\r
444 fHistMassLP->GetYaxis()->SetTitle("Entries");
\r
445 fHistMassLP->SetMarkerStyle(kFullCircle);
\r
447 fHistProtonPIDBb = new TH2F("fHistProtonPIDBb", "dE/dx after p PID", 100, 0., 10, 100, 0, 100);
\r
448 fHistProtonPIDBb->GetYaxis()->SetTitle("TPC Signal");
\r
449 fHistProtonPIDBb->GetXaxis()->SetTitle("P (GeV/c)");
\r
450 fHistProtonPIDBb->SetOption("scat");
\r
451 fHistProtonPIDBb->SetMarkerStyle(kFullCircle);
\r
453 fHistPionPIDBb = new TH2F("fHistPionPIDBb", "dE/dx after K PID", 100, 0., 10, 100, 0, 100);
\r
454 fHistPionPIDBb->GetYaxis()->SetTitle("TPC Signal");
\r
455 fHistPionPIDBb->GetXaxis()->SetTitle("P (GeV/c)");
\r
456 fHistPionPIDBb->SetOption("scat");
\r
457 fHistPionPIDBb->SetMarkerStyle(kFullCircle);
\r
459 fHistProtonPIDLambda = new TH2F("fHistProtonPIDLambda", "dE/dx after p PID from V0", 100, 0., 10, 100, 0, 100);
\r
460 fHistProtonPIDLambda->GetYaxis()->SetTitle("TPC Signal");
\r
461 fHistProtonPIDLambda->GetXaxis()->SetTitle("P (GeV/c)");
\r
462 fHistProtonPIDLambda->SetOption("scat");
\r
463 fHistProtonPIDLambda->SetMarkerStyle(kFullCircle);
\r
465 fHistPionPIDLambda = new TH2F("fHistPionPIDLambda", "dE/dx after #pi PID from V0", 100, 0, 10, 100, 0, 100);
\r
466 fHistPionPIDLambda->GetYaxis()->SetTitle("TPC Signal");
\r
467 fHistPionPIDLambda->GetXaxis()->SetTitle("P (GeV/c)");
\r
468 fHistPionPIDLambda->SetOption("scat");
\r
469 fHistPionPIDLambda->SetMarkerStyle(kFullCircle);
\r
471 fHistMCdcaPvtxDvtx= new TH1F("fHistMCdcaPvtxDvtx", "MC True DCA Primary Vertex - H Decay Vertex", 300, -0.1, 11.9);
\r
472 fHistMCdcaPvtxDvtx->GetXaxis()->SetTitle("dca prim. vtx- decay vtx (cm)");
\r
473 fHistMCdcaPvtxDvtx->GetYaxis()->SetTitle("Entries");
\r
474 fHistMCdcaPvtxDvtx->SetMarkerStyle(kFullCircle);
\r
476 fHistMCdcaPvtxLvtx= new TH1F("fHistMCdcaPvtxLvtx", "MC True DCA Primary Vertex - Lambda Decay Vertex", 300, -0.1, 11.9);
\r
477 fHistMCdcaPvtxLvtx->GetXaxis()->SetTitle("dca prim. vtx-#Lambda decay vtx (cm)");
\r
478 fHistMCdcaPvtxLvtx->GetYaxis()->SetTitle("Entries");
\r
479 fHistMCdcaPvtxLvtx->SetMarkerStyle(kFullCircle);
\r
481 fHistMCdcaDvtxLvtx= new TH1F("fHistMCdcaDvtxLvtx", "MC True DCA H Decay Vertex - Lambda Decay Vertex", 300, -0.1, 11.9);
\r
482 fHistMCdcaDvtxLvtx->GetXaxis()->SetTitle("dca H deacy vtx-#Lambda decay vtx (cm)");
\r
483 fHistMCdcaDvtxLvtx->GetYaxis()->SetTitle("Entries");
\r
484 fHistMCdcaDvtxLvtx->SetMarkerStyle(kFullCircle);
\r
486 fHistMCangleLH= new TH1F("fHistMCangleLH", "MC True Angle between #Lambda and H", 300, 0., 2*TMath::Pi());
\r
487 fHistMCangleLH->GetXaxis()->SetTitle("Angle (#Lambda-H)");
\r
488 fHistMCangleLH->GetYaxis()->SetTitle("Entries");
\r
489 fHistMCangleLH->SetMarkerStyle(kFullCircle);
\r
491 fHistMCdecayAngle= new TH1F("fHistMCdecayAngle", "MC True Angle between decay products", 300, 0., 2*TMath::Pi());
\r
492 fHistMCdecayAngle->GetXaxis()->SetTitle("Angle (#Lambda-p#pi)");
\r
493 fHistMCdecayAngle->GetYaxis()->SetTitle("Entries");
\r
494 fHistMCdecayAngle->SetMarkerStyle(kFullCircle);
\r
496 fHistMCpointingAngle= new TH1F("fHistMCpointingAngle", "MC True Pointing Angle", 3000, 0., 2*TMath::Pi());
\r
497 fHistMCpointingAngle->GetXaxis()->SetTitle("Pointing Angle");
\r
498 fHistMCpointingAngle->GetYaxis()->SetTitle("Entries");
\r
499 fHistMCpointingAngle->SetMarkerStyle(kFullCircle);
\r
501 fHistMCap = new TH2F("fHistMCap", "True MC Armenteros-Podolanski", 200,-1.0,1.0, 500,0,1);
\r
502 fHistMCap->GetYaxis()->SetTitle("#alpha");
\r
503 fHistMCap->GetXaxis()->SetTitle("q_{t}");
\r
504 fHistMCap->SetOption("scat");
\r
505 fHistMCap->SetMarkerStyle(kFullCircle);
\r
507 fHistMCdcaPvtxDvtxReso= new TH1F("fHistMCdcaPvtxDvtxReso", "MC True DCA Primary Vertex - H Decay Vertex", 300, -0.1, 11.9);
\r
508 fHistMCdcaPvtxDvtxReso->GetXaxis()->SetTitle("dca prim. vtx- decay vtx (cm)");
\r
509 fHistMCdcaPvtxDvtxReso->GetYaxis()->SetTitle("Entries");
\r
510 fHistMCdcaPvtxDvtxReso->SetMarkerStyle(kFullCircle);
\r
512 fHistMCdcaPvtxLvtxReso= new TH1F("fHistMCdcaPvtxLvtxReso", "MC True DCA Primary Vertex - Lambda Decay Vertex", 300, -0.1, 11.9);
\r
513 fHistMCdcaPvtxLvtxReso->GetXaxis()->SetTitle("dca prim. vtx-#Lambda decay vtx (cm)");
\r
514 fHistMCdcaPvtxLvtxReso->GetYaxis()->SetTitle("Entries");
\r
515 fHistMCdcaPvtxLvtxReso->SetMarkerStyle(kFullCircle);
\r
517 fHistMCdcaDvtxLvtxReso= new TH1F("fHistMCdcaDvtxLvtxReso", "MC True DCA H Decay Vertex - Lambda Decay Vertex", 300, -0.1, 11.9);
\r
518 fHistMCdcaDvtxLvtxReso->GetXaxis()->SetTitle("dca H deacy vtx-#Lambda decay vtx (cm)");
\r
519 fHistMCdcaDvtxLvtxReso->GetYaxis()->SetTitle("Entries");
\r
520 fHistMCdcaDvtxLvtxReso->SetMarkerStyle(kFullCircle);
\r
522 fHistMCangleLHReso= new TH1F("fHistMCangleLHReso", "MC True Angle between #Lambda and H", 300, 0., 2*TMath::Pi());
\r
523 fHistMCangleLHReso->GetXaxis()->SetTitle("Angle (#Lambda-H)");
\r
524 fHistMCangleLHReso->GetYaxis()->SetTitle("Entries");
\r
525 fHistMCangleLHReso->SetMarkerStyle(kFullCircle);
\r
527 fHistMCdecayAngleReso= new TH1F("fHistMCdecayAngleReso", "MC True Angle between decay products", 300, 0., 2*TMath::Pi());
\r
528 fHistMCdecayAngleReso->GetXaxis()->SetTitle("Angle (#Lambda-p#pi)");
\r
529 fHistMCdecayAngleReso->GetYaxis()->SetTitle("Entries");
\r
530 fHistMCdecayAngleReso->SetMarkerStyle(kFullCircle);
\r
532 fHistMCpointingAngleReso= new TH1F("fHistMCpointingAngleReso", "MC True Pointing Angle", 300, 0., 2*TMath::Pi());
\r
533 fHistMCpointingAngleReso->GetXaxis()->SetTitle("Pointing Angle");
\r
534 fHistMCpointingAngleReso->GetYaxis()->SetTitle("Entries");
\r
535 fHistMCpointingAngleReso->SetMarkerStyle(kFullCircle);
\r
537 fHistMCapReso = new TH2F("fHistMCapReso", "True MC Armenteros-Podolanski", 200,-1.0,1.0, 500,0,1);
\r
538 fHistMCapReso->GetYaxis()->SetTitle("#alpha");
\r
539 fHistMCapReso->GetXaxis()->SetTitle("q_{t}");
\r
540 fHistMCapReso->SetOption("scat");
\r
541 fHistMCapReso->SetMarkerStyle(kFullCircle);
\r
543 fHistCentrality = new TH1F("Centrality ", "Centrality", 11, -0.5, 10.5);
\r
544 fHistCentrality ->GetXaxis()->SetTitle("Centrality");
\r
545 fHistCentrality ->GetYaxis()->SetTitle("Entries");
\r
547 fHistCentralityAC = new TH1F("CentralityAC ", "CentralityAC", 11, -0.5, 10.5);
\r
548 fHistCentralityAC ->GetXaxis()->SetTitle("Centrality");
\r
549 fHistCentralityAC ->GetYaxis()->SetTitle("Entries");
\r
551 fHistMultiplicity = new TH1F("Multiplicity ", "Multiplicity", 100, 0, 20000);
\r
552 fHistMultiplicity ->GetXaxis()->SetTitle("Centrality");
\r
553 fHistMultiplicity ->GetYaxis()->SetTitle("Entries");
\r
555 fHistList->Add(fHistMassDPi);
\r
556 fHistList->Add(fHistMassLPi);
\r
557 fHistList->Add(fHistMassLambdaPi);
\r
558 fHistList->Add(fHistMassLambda);
\r
559 fHistList->Add(fHistMassLambdaPPi);
\r
560 fHistList->Add(fHistMassLambdaP);
\r
561 fHistList->Add(fHistMassLambdaK);
\r
562 fHistList->Add(fHistMassK0onl);
\r
563 fHistList->Add(fHistMassK0offl);
\r
564 fHistList->Add(fHistMassK0onlC);
\r
565 fHistList->Add(fHistMassK0offlC);
\r
566 fHistList->Add(fHistMassPQonl);
\r
567 fHistList->Add(fHistMassPQoffl);
\r
568 fHistList->Add(fHistDC);
\r
569 fHistList->Add(fHistArmenterosPodolanski);
\r
570 fHistList->Add(fHistArmenterosPodolanskiCut);
\r
571 fHistList->Add(fHistHDibaryonInvaMassGen);
\r
572 fHistList->Add(fHistHDibaryonInvaMassGenRes);
\r
573 fHistList->Add(fHistAntiHDibaryonInvaMassGen);
\r
574 fHistList->Add(fHistHDibaryonInvaMassAso);
\r
575 fHistList->Add(fHistHDibaryonInvaMassAsoReso);
\r
576 fHistList->Add(fHistAntiHDibaryonInvaMassAso);
\r
577 fHistList->Add(fHistCheck);
\r
578 fHistList->Add(fHistHPointingAngle);
\r
579 fHistList->Add(fHistMassH);
\r
580 fHistList->Add(fHistMassPpi);
\r
581 fHistList->Add(fHistMassPpiReso);
\r
582 fHistList->Add(fHistMassLpi);
\r
583 fHistList->Add(fHistMassLP);
\r
584 fHistList->Add(fHistMassLambdaFromH);
\r
585 fHistList->Add(fHistMassLambdaFromHtLorentz);
\r
586 fHistList->Add(fHistProtonPIDBb);
\r
587 fHistList->Add(fHistPionPIDBb);
\r
588 fHistList->Add(fHistProtonPIDLambda);
\r
589 fHistList->Add(fHistPionPIDLambda);
\r
590 fHistList->Add(fHistMCdcaPvtxDvtx);
\r
591 fHistList->Add(fHistMCdcaPvtxLvtx);
\r
592 fHistList->Add(fHistMCdcaDvtxLvtx);
\r
593 fHistList->Add(fHistMCangleLH);
\r
594 fHistList->Add(fHistMCdecayAngle);
\r
595 fHistList->Add(fHistMCpointingAngle);
\r
596 fHistList->Add(fHistMCap);
\r
597 fHistList->Add(fHistMCdcaPvtxDvtxReso);
\r
598 fHistList->Add(fHistMCdcaPvtxLvtxReso);
\r
599 fHistList->Add(fHistMCdcaDvtxLvtxReso);
\r
600 fHistList->Add(fHistMCangleLHReso);
\r
601 fHistList->Add(fHistMCdecayAngleReso);
\r
602 fHistList->Add(fHistMCpointingAngleReso);
\r
603 fHistList->Add(fHistMCapReso);
\r
604 fHistList->Add(fHistCentrality);
\r
605 fHistList->Add(fHistCentralityAC);
\r
606 fHistList->Add(fHistMultiplicity);
\r
608 fHistHilf1= new TH1F("fHistHilf1", "HD", 300, 0., 10);
\r
609 fHistHilf1->GetXaxis()->SetTitle("");
\r
610 fHistHilf1->GetYaxis()->SetTitle("Entries");
\r
612 fHistHilf2= new TH1F("fHistHilf2", "HD", 300, 0., 10);
\r
613 fHistHilf2->GetXaxis()->SetTitle("");
\r
614 fHistHilf2->GetYaxis()->SetTitle("Entries");
\r
616 fHistHilf3= new TH1F("fHistHilf3", "HD", 300, 0., 10);
\r
617 fHistHilf3->GetXaxis()->SetTitle("");
\r
618 fHistHilf3->GetYaxis()->SetTitle("Entries");
\r
620 fHistHilf4= new TH1F("fHistHilf4", "HD", 300, 0., 10);
\r
621 fHistHilf4->GetXaxis()->SetTitle("");
\r
622 fHistHilf4->GetYaxis()->SetTitle("Entries");
\r
624 fHistHilf5= new TH1F("fHistHilf5", "HD", 300, 0., 10);
\r
625 fHistHilf5->GetXaxis()->SetTitle("");
\r
626 fHistHilf5->GetYaxis()->SetTitle("Entries");
\r
628 fHistHilf6= new TH1F("fHistHilf6", "HD", 300, 0., 10);
\r
629 fHistHilf6->GetXaxis()->SetTitle("");
\r
630 fHistHilf6->GetYaxis()->SetTitle("Entries");
\r
632 fHistPtvsEtaGen = new TH2F("fHistPtvsEtaGen", "p_{t} vs #eta from generated H", 200,0.0,10.0, 200,-1,1);
\r
633 fHistPtvsEtaGen->GetXaxis()->SetTitle("p_{t} (GeV/c)");
\r
634 fHistPtvsEtaGen->GetYaxis()->SetTitle("#eta");
\r
635 fHistPtvsEtaGen->SetOption("scat");
\r
636 fHistPtvsEtaGen->SetMarkerStyle(kFullCircle);
\r
638 fHistPtvsEtaAso = new TH2F("fHistPtvsEtaAso", "p_{t} vs #eta from associated H", 200,0.0,10.0, 200,-1,1);
\r
639 fHistPtvsEtaAso->GetYaxis()->SetTitle("p_{t} (GeV/c)");
\r
640 fHistPtvsEtaAso->GetXaxis()->SetTitle("#eta");
\r
641 fHistPtvsEtaAso->SetOption("scat");
\r
642 fHistPtvsEtaAso->SetMarkerStyle(kFullCircle);
\r
644 fHistPtvsYGen = new TH2F("fHistPtvsYGen", "p_{t} vs rapidity from generated H", 200,0.0,10.0, 200,-1,1);
\r
645 fHistPtvsYGen->GetXaxis()->SetTitle("p_{t} (GeV/c)");
\r
646 fHistPtvsYGen->GetYaxis()->SetTitle("y");
\r
647 fHistPtvsYGen->SetOption("scat");
\r
648 fHistPtvsYGen->SetMarkerStyle(kFullCircle);
\r
650 fHistPtvsYAso = new TH2F("fHistPtvsYAso", "p_{t} vs rapidity from associated H", 200,0.0,10.0, 200,-1,1);
\r
651 fHistPtvsYAso->GetXaxis()->SetTitle("p_{t} (GeV/c)");
\r
652 fHistPtvsYAso->GetYaxis()->SetTitle("y");
\r
653 fHistPtvsYAso->SetOption("scat");
\r
654 fHistPtvsYAso->SetMarkerStyle(kFullCircle);
\r
656 fHistRap= new TH1F("fHistRap", "Rapidity", 400, -2., 2);
\r
657 fHistRap->GetXaxis()->SetTitle("Y");
\r
658 fHistRap->GetYaxis()->SetTitle("Entries");
\r
659 fHistPtvsEtaAso->SetMarkerStyle(kFullCircle);
\r
661 fHistList->Add(fHistHilf1);
\r
662 fHistList->Add(fHistHilf2);
\r
663 fHistList->Add(fHistHilf3);
\r
664 fHistList->Add(fHistHilf4);
\r
665 fHistList->Add(fHistHilf5);
\r
666 fHistList->Add(fHistHilf6);
\r
667 fHistList->Add(fHistPtvsEtaGen);
\r
668 fHistList->Add(fHistPtvsEtaAso);
\r
669 fHistList->Add(fHistPtvsYGen);
\r
670 fHistList->Add(fHistPtvsYAso);
\r
671 fHistList->Add(fHistRap);
\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
695 //trigger statistics histogram
\r
696 fHistTriggerStat = new TH1F("fHistTriggerStat","Trigger statistics", 4,-0.5, 3.5);
\r
697 const Char_t* aTriggerNames[] = { "kMB", "kCentral", "kSemiCentral" };
\r
698 for ( Int_t ii=0; ii < 3; ii++ )
\r
699 fHistTriggerStat->GetXaxis()->SetBinLabel(ii+1, aTriggerNames[ii]);
\r
701 fHistTriggerStatAfterEventSelection = new TH1F("fHistTriggerStatAfterEventSelection","Trigger statistics after event selection", 4,-0.5, 3.5);
\r
702 for ( Int_t ii=0; ii < 3; ii++ )
\r
703 fHistTriggerStatAfterEventSelection->GetXaxis()->SetBinLabel(ii+1, aTriggerNames[ii]);
\r
704 fHistList->Add(fHistTriggerStat);
\r
705 fHistList->Add(fHistTriggerStatAfterEventSelection);
\r
707 fHistMassHcentMult = new TH3F("fHistMassHcentMult", "Inv. Mass vs. centrality vs. multiplicity", 100, 2.2, 2.3, 5, 0, 4, 300, 0, 6000);
\r
708 fHistMassHcentMult->GetXaxis()->SetTitle("Invariant mass #Lambdap#pi^{-} (GeV/c^{2})"); // inv. mass
\r
709 fHistMassHcentMult->GetYaxis()->SetTitle("Centrality"); // triggertype
\r
710 fHistMassHcentMult->GetZaxis()->SetTitle("Multiplicity"); // refTPC
\r
711 fHistList->Add(fHistMassHcentMult);
\r
714 const Double_t kz = 2*TMath::Pi();
\r
715 Int_t binsD01[16]={ 300, 200, 100, 100, 100, 100, 100, 100, 200, 200, 200, 200, 400, 200, 200, 3};
\r
716 Double_t xminD01[16]={2.0, 1.0, 0., -1, 0., 0., 0., 0., 0., 0., 0., 0., -1, 0., 0., 0};
\r
717 Double_t xmaxD01[16]={2.3, 1.2, kz, 1, 1, 10, 10, 5, 5, 5, 5, 100, 1, 100, 4000, 1};
\r
720 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
721 fHistList->Add(fHistNdim);
\r
723 AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
\r
724 AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
\r
725 fPIDtpcESD = inputHandler->GetPIDResponse();
\r
727 // Post output data (if histograms are not used later, PostData is at least called here)
\r
728 PostData(1, fHistList);
\r
732 //________________________________________________________________________
\r
733 void AliAnalysisTaskHdibaryonLPpi::UserExec(Option_t *)
\r
736 // Called for each event
\r
738 //define improtant masses
\r
739 Double_t pionMass = 0.13957;
\r
740 Double_t protonMass = 0.93827;
\r
743 Long_t pdgPionPlus = 211;
\r
744 Long_t pdgPionMinus = -211;
\r
745 Long_t pdgProton = 2212;
\r
746 Long_t pdgAntiProton = -2212;
\r
747 Long_t pdgLambda = 3122;
\r
748 Long_t pdgAntiLambda = -3122;
\r
749 Long_t pdgHDibaryon = 1020000020;
\r
750 Long_t pdgAntiHDibaryon = -1020000020;
\r
752 AliStack* stack(NULL);
\r
755 if(!fMCEvent)return;
\r
757 AliHeader *head = fMCEvent->Header();
\r
759 AliGenPythiaEventHeader *pyheader = (AliGenPythiaEventHeader*)head->GenEventHeader();
\r
760 if(!pyheader)return;
\r
762 if(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler()){
\r
763 if(static_cast<AliMCEventHandler*>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())->MCEvent()){
\r
764 if(!(static_cast<AliMCEventHandler*>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())->TreeK()))return;
\r
765 if(!(static_cast<AliMCEventHandler*>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())->TreeTR()))return;
\r
766 if(!(static_cast<AliMCEventHandler*>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())->Init("local")))return;
\r
767 stack = static_cast<AliMCEventHandler*>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())->MCEvent()->Stack();
\r
768 if(!(static_cast<AliMCEventHandler*>(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler())->MCEvent()->Stack()->TreeK()))return;
\r
775 // -------------------------------------------------------
\r
776 // Loop for Inv. Mass via ESD tracks
\r
777 // -------------------------------------------------------
\r
779 fHistCount->Fill(0);
\r
781 fESD = dynamic_cast<AliESDEvent *>(fInputEvent);
\r
784 //Printf("ERROR: fESD not available");
789 const AliESDVertex *vertex = fESD->GetPrimaryVertexTracks();
\r
790 if (vertex->GetNContributors()<1)
\r
793 vertex = fESD->GetPrimaryVertexSPD();
\r
794 if(vertex->GetNContributors()<1) {
\r
798 if (TMath::Abs(vertex->GetZv()) > 10) return;
\r
800 Int_t centrality = -5;
\r
801 Double_t centrPerc = -5;
\r
803 if (fESD->GetEventSpecie() == 4)
\r
805 AliCentrality *esdCentrality = fESD->GetCentrality();
\r
806 centrality = esdCentrality->GetCentralityClass10("V0M"); // centrality percentile determined with V0
\r
807 centrPerc = esdCentrality->GetCentralityPercentile("V0M");
\r
808 // if (centrality < 0. || centrality > 8.) return; //0 bis 80 %
\r
809 if (centrality > 8) return; //0 bis 80 %
\r
810 // cout<<"Centrality: "<< centrality << endl;
\r
814 fHistCentrality->Fill(centrality);
\r
816 //*****************//
\r
818 //*****************//
\r
820 // Float_t percentile=centrality->GetCentralityPercentile("V0M");
\r
822 Bool_t isSelectedCentral = (((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected() & AliVEvent::kCentral);
\r
823 Bool_t isSelectedSemiCentral = (((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected() & AliVEvent::kSemiCentral);
\r
824 Bool_t isSelectedMB = (((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected() & AliVEvent::kMB);
\r
826 Int_t triggertype = 17;
\r
828 if(isSelectedCentral){
\r
829 fHistTriggerStat->Fill(1);
\r
833 if(isSelectedSemiCentral){
\r
834 fHistTriggerStat->Fill(2);
\r
839 fHistTriggerStat->Fill(0);
\r
843 // if(isSelectedCentral || isSelectedSemiCentral || isSelectedMB){
\r
845 //*******************************
\r
847 Int_t runNumber = 0;
\r
849 runNumber = fESD->GetRunNumber();
\r
851 if (!fPIDtpcESD) fPIDtpcESD = ((AliESDInputHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->GetESDpid();
\r
853 fPIDtpcESD = new AliESDpid(); // HACK FOR MC PBPB --> PLEASE REMOVE AS SOON AS POSSIBLE
\r
854 fPIDtpcESD->GetTPCResponse().SetBetheBlochParameters(1.28778e+00/50., 3.13539e+01, TMath::Exp(-3.16327e+01), 1.87901e+00, 6.41583e+00);
\r
860 TObjArray* listCrossV0 = fESDCutsV0->GetAcceptedV0s(fESD);
\r
861 Int_t nGoodV0s = listCrossV0->GetEntries();
\r
863 const AliESDVertex *esdVer = fESD->GetPrimaryVertex();
\r
864 AliESDVertex *esdVer1 = new AliESDVertex(*esdVer);
\r
866 AliVertexerTracks *vertexer = new AliVertexerTracks(fESD->GetMagneticField());
\r
867 TObjArray *trkArray = new TObjArray(2);
\r
868 AliVertexerTracks *vertexer1 = new AliVertexerTracks(fESD->GetMagneticField());
\r
869 TObjArray *trkArray1 = new TObjArray(2);
\r
871 AliKFParticle::SetField(fESD->GetMagneticField());
\r
873 AliKFVertex primVtx(*(fESD->GetPrimaryVertex()));
\r
875 Int_t refMultTpc = AliESDtrackCuts::GetReferenceMultiplicity(fESD, kTRUE);
\r
876 //cout<<"Multiplicity: "<< refMultTpc << endl;
\r
877 fHistMultiplicity->Fill(refMultTpc);
\r
879 Double_t mn[3] = {0,0,0};
\r
880 Double_t mp[3] = {0,0,0};
\r
881 Double_t mm[3] = {0,0,0};
\r
882 Double_t dd[3] = {0,0,0};
\r
883 Double_t dd1[3] = {0,0,0};
\r
884 const Double_t cProtonMass=TDatabasePDG::Instance()->GetParticle(2212)->Mass();
\r
885 const Double_t cPionMass=TDatabasePDG::Instance()->GetParticle(211)->Mass();
\r
886 const Double_t cElectronMass=TDatabasePDG::Instance()->GetParticle(11)->Mass();
\r
887 const Double_t cLambdaMass=TDatabasePDG::Instance()->GetParticle(3122)->Mass();
\r
888 Double_t decayLength=0;
\r
889 Double_t decayLengthH=0;
\r
891 //V0 Loop for Lambda and Anti-Lambda
\r
892 for(Int_t iV0MI = 0; iV0MI < nGoodV0s ; iV0MI++) {
\r
893 AliESDv0 * fV0MIs = fESD->GetV0(iV0MI);
\r
894 Int_t lOnFlyStatus = 0;
\r
896 lOnFlyStatus = fV0MIs->GetOnFlyStatus();
\r
897 Double_t lInvMassLambda=0;
\r
898 Double_t lInvMassLambdaPi=0;
\r
899 Double_t lPtLambda=0;
\r
900 Double_t lPzLambda=0;
\r
901 Double_t lPLambda=0;
\r
905 TLorentzVector posE;
\r
906 TLorentzVector negE;
\r
907 TLorentzVector photon;
\r
913 if (!lOnFlyStatus){
\r
918 // fHistMultiplicity->Fill(refMultTpc);
\r
919 fHistCentralityAC->Fill(centrality);
\r
920 fHistCheck->Fill(offl,onl);
\r
922 AliESDtrack* trackPosTest = fESD->GetTrack(fV0MIs->GetPindex());
\r
923 AliESDtrack* trackNegTest = fESD->GetTrack(fV0MIs->GetNindex());
\r
926 if (!fEsdTrackCuts->AcceptTrack(trackPosTest)) continue;
\r
927 if (!fESDtrackCutsV0->AcceptTrack(trackNegTest)) continue;
\r
930 //PID via specific energy loss in the TPC
\r
931 //define the arrays for the Bethe-Bloch-Parameters
\r
932 Double_t parProton[5] = {0,0,0,0,0};
\r
934 if(runNumber < 166500) //LHC10h
\r
936 parProton[0] = 1.45802; // ALEPH parameters for protons (pass2)
\r
937 parProton[1] = 27.4992;
\r
938 parProton[2] = 4.00313e-15;
\r
939 parProton[3] = 2.48485;
\r
940 parProton[4] = 8.31768;
\r
943 if(runNumber > 166500) //LHC11h
\r
945 parProton[0] = 1.11243; // ALEPH parameters for protons (pass2)
\r
946 parProton[1] = 26.1144;
\r
947 parProton[2] = 4.00313e-15;
\r
948 parProton[3] = 2.72969;
\r
949 parProton[4] = 9.15038;
\r
952 //Get the total momentum for each track at the inner readout of the TPC
\r
953 Double_t ptotN = trackNegTest->GetInnerParam()->GetP();
\r
954 Double_t ptotP = trackPosTest->GetInnerParam()->GetP();
\r
956 //define expected signals for the various species
\r
957 Double_t expSignalPionP = 0;
\r
958 Double_t expSignalPionN = 0;
\r
959 Double_t expSignalProtonN = 0;
\r
960 Double_t expSignalProtonP = 0;
\r
965 expSignalProtonN = AliExternalTrackParam::BetheBlochAleph(ptotN/(protonMass),parProton[0],parProton[1],parProton[2],parProton[3],parProton[4]);
\r
966 expSignalProtonP = AliExternalTrackParam::BetheBlochAleph(ptotP/(protonMass),parProton[0],parProton[1],parProton[2],parProton[3],parProton[4]);
\r
972 expSignalPionP = 0.7*AliExternalTrackParam::BetheBlochAleph(ptotP/(pionMass),parProton[0],parProton[1],parProton[2],parProton[3],parProton[4]);
\r
973 expSignalPionN = 0.7*AliExternalTrackParam::BetheBlochAleph(ptotN/(pionMass),parProton[0],parProton[1],parProton[2],parProton[3],parProton[4]);
\r
975 expSignalProtonN = 0.65*AliExternalTrackParam::BetheBlochAleph(ptotN/(protonMass),parProton[0],parProton[1],parProton[2],parProton[3],parProton[4]);
\r
976 expSignalProtonP = 0.65*AliExternalTrackParam::BetheBlochAleph(ptotP/(protonMass),parProton[0],parProton[1],parProton[2],parProton[3],parProton[4]);
\r
979 // PID cut on the nuclei (proton, deuteron, triton, helium3)
\r
980 Bool_t corrParticle = kFALSE;
\r
982 Bool_t posProton = kFALSE;
\r
986 if (TMath::Abs(fPIDtpcESD->NumberOfSigmasTPC(trackPosTest, AliPID::kProton)) < 3)
\r
989 corrParticle = kTRUE;
\r
995 if(//trackPosTest->GetTPCsignal() < 1200 &&
\r
996 TMath::Abs(trackPosTest->GetTPCsignal() - expSignalProtonP)/expSignalProtonP < 0.4)
\r
999 corrParticle = kTRUE;
\r
1003 Bool_t negProton = kFALSE;
\r
1007 if (TMath::Abs(fPIDtpcESD->NumberOfSigmasTPC(trackNegTest, AliPID::kProton)) < 3)
\r
1009 negProton = kTRUE;
\r
1010 corrParticle = kTRUE;
\r
1016 if(//trackNegTest->GetTPCsignal() < 1200 &&
\r
1017 TMath::Abs(trackNegTest->GetTPCsignal() - expSignalProtonN)/expSignalProtonN < 0.4)
\r
1019 negProton = kTRUE;
\r
1020 corrParticle = kTRUE;
\r
1024 //PID cut for pions
\r
1025 //data: 3sigma cut on the pions
\r
1027 Bool_t negPion = kFALSE;
\r
1028 Bool_t posPion = kFALSE;
\r
1032 if (TMath::Abs(fPIDtpcESD->NumberOfSigmasTPC(trackPosTest, AliPID::kPion)) < 4) posPion=kTRUE; //pos daughter has to be a pion
\r
1033 if (TMath::Abs(fPIDtpcESD->NumberOfSigmasTPC(trackNegTest, AliPID::kPion)) < 4) negPion=kTRUE; // negative daughter has to be a pion
\r
1036 //MC: like the nuclei via the specific energyloss in the TPC
\r
1039 if (TMath::Abs(trackPosTest->GetTPCsignal() - expSignalPionP)/expSignalPionP < 0.4) posPion=kTRUE;
\r
1040 if (TMath::Abs(trackNegTest->GetTPCsignal() - expSignalPionN)/expSignalPionN < 0.4) negPion=kTRUE;
\r
1043 if (!(posProton==kTRUE)) continue;
\r
1044 if (!(negPion==kTRUE)) continue;
\r
1048 if( !(trackPosTest->GetStatus() & AliESDtrack::kTPCrefit)){
\r
1052 if( !(trackNegTest->GetStatus() & AliESDtrack::kTPCrefit)){
\r
1056 if( trackPosTest->GetSign() >0 && trackNegTest->GetSign() <0){
\r
1057 fV0MIs->GetNPxPyPz(mn[0],mn[1],mn[2]); //reconstructed cartesian momentum components of negative daughter
\r
1058 fV0MIs->GetPPxPyPz(mp[0],mp[1],mp[2]); //reconstructed cartesian momentum components of positive daughter
\r
1061 if( trackPosTest->GetSign() <0 && trackNegTest->GetSign() >0){
\r
1062 fV0MIs->GetPPxPyPz(mn[0],mn[1],mn[2]); //reconstructed cartesian momentum components of negative daughter
\r
1063 fV0MIs->GetNPxPyPz(mp[0],mp[1],mp[2]); //reconstructed cartesian momentum components of positive daughter
\r
1066 fV0MIs->GetPxPyPz(mm[0],mm[1],mm[2]); //reconstructed cartesian momentum components of mother
\r
1068 TVector3 vecN(mn[0],mn[1],mn[2]);
\r
1069 TVector3 vecP(mp[0],mp[1],mp[2]);
\r
1070 TVector3 vecM(mm[0],mm[1],mm[2]);
\r
1072 Double_t thetaP = acos((vecP * vecM)/(vecP.Mag() * vecM.Mag()));
\r
1073 Double_t thetaN = acos((vecN * vecM)/(vecN.Mag() * vecM.Mag()));
\r
1075 Double_t alfa = ((vecP.Mag())*cos(thetaP)-(vecN.Mag())*cos(thetaN))/
\r
1076 ((vecP.Mag())*cos(thetaP)+(vecN.Mag())*cos(thetaN)) ;
\r
1077 Double_t qt = vecP.Mag()*sin(thetaP);
\r
1079 fHistArmenterosPodolanski->Fill(alfa,qt); //Armenteros-Podolanski calculation
\r
1081 TLorentzVector k0;
\r
1082 TLorentzVector k0daugh1;
\r
1083 TLorentzVector k0daugh2;
\r
1084 TLorentzVector proton;
\r
1085 TLorentzVector pq;
\r
1087 k0daugh1.SetXYZM(mn[0],mn[1],mn[2],cPionMass);
\r
1088 k0daugh2.SetXYZM(mp[0],mp[1],mp[2],cPionMass);
\r
1089 k0=k0daugh1+k0daugh2;
\r
1091 fV0MIs->ChangeMassHypothesis(3122);
\r
1092 lInvMassLambda = fV0MIs->GetEffMass();
\r
1093 lPtLambda = fV0MIs->Pt();
\r
1094 lPzLambda = fV0MIs->Pz();
\r
1095 lPLambda = fV0MIs->P();
\r
1097 trkArray->AddAt(trackPosTest,0);
\r
1098 trkArray->AddAt(trackNegTest,1);
\r
1100 vertexer->SetVtxStart(esdVer1);
\r
1101 AliESDVertex *decayVertex = (AliESDVertex*)vertexer->VertexForSelectedESDTracks(trkArray);
\r
1103 dd[0]=fESD->GetPrimaryVertexSPD()->GetX()-decayVertex->GetX();
\r
1104 dd[1]=fESD->GetPrimaryVertexSPD()->GetY()-decayVertex->GetY();
\r
1105 dd[2]=fESD->GetPrimaryVertexSPD()->GetZ()-decayVertex->GetZ();
\r
1107 decayLength=sqrt(dd[0]*dd[0]+dd[1]*dd[1]+dd[2]*dd[2]);
\r
1109 if (decayVertex == NULL) cout << "Lambda decay vtx pointer NULL" << endl;
1110 if (decayVertex) delete decayVertex;
\r
1112 TLorentzVector negPio1;
\r
1113 TLorentzVector posProt1;
\r
1114 TLorentzVector posP;
\r
1115 TLorentzVector posProt;
\r
1116 TLorentzVector negK;
\r
1117 TLorentzVector negPio;
\r
1118 TLorentzVector negPi;
\r
1119 TLorentzVector omega;
\r
1120 TLorentzVector threeSum;
\r
1121 TLorentzVector fourSum;
\r
1122 TLorentzVector ppK;
\r
1123 TLorentzVector posPiK;
\r
1124 TLorentzVector negPiK;
\r
1125 TLorentzVector kaon;
\r
1126 TLorentzVector lambda;
\r
1127 TLorentzVector lambdaH;
\r
1128 TLorentzVector hDibaryon;
\r
1134 h.SetXYZ(-dd[0],-dd[1],-dd[2]);
\r
1137 if (offl==1)fHistMassDPi->Fill(lInvMassLambda);
\r
1140 fHistMassLPi->Fill(lInvMassLambda);
\r
1142 negE.SetXYZM(mn[0],mn[1],mn[2],cElectronMass);
\r
1143 posE.SetXYZM(mp[0],mp[1],mp[2],cElectronMass);
\r
1146 negPiK.SetXYZM(mn[0],mn[1],mn[2],cPionMass);
\r
1147 posPiK.SetXYZM(mp[0],mp[1],mp[2],cPionMass);
\r
1148 kaon=posPiK+negPiK;
\r
1150 negPi.SetXYZM(mn[0],mn[1],mn[2],cPionMass);
\r
1151 posP.SetXYZM(mp[0],mp[1],mp[2],cProtonMass);
\r
1152 lambda=negPi+posP;
\r
1153 lambdaH.SetXYZM(mm[0],mm[1],mm[2],cLambdaMass);
\r
1155 if (lInvMassLambda>1.1113&&lInvMassLambda<1.1202){
\r
1158 if (qt<-2.21*alfa*alfa+2.945*alfa-0.887) continue;
\r
1159 if (qt>-2.21*alfa*alfa+2.945*alfa-0.873) continue;
\r
1160 if (photon.M()<0.005) continue;
\r
1161 if (kaon.M()>0.495 && kaon.M()<0.500 ) continue;
\r
1165 fHistMassLambda->Fill(lInvMassLambda);
\r
1167 Bool_t isCorrectlyAssociatedLambda = kFALSE;
\r
1168 Bool_t isPartialCorrectlyAssociatedLambda = kFALSE;
\r
1169 Int_t labelAssociatedH=-1;
\r
1170 Int_t labelLambda=-1;
\r
1173 Int_t labelPosTest = trackPosTest->GetLabel();
\r
1174 TParticle *tparticleDaughter = stack->Particle(TMath::Abs(labelPosTest));
\r
1175 Int_t labelMother = tparticleDaughter->GetFirstMother();
\r
1176 TParticle *tparticleMother = stack->Particle(TMath::Abs(labelMother));
\r
1178 Int_t labelOma = tparticleMother->GetFirstMother();
\r
1179 TParticle *tparticleOma = stack->Particle(TMath::Abs(labelOma));
\r
1181 if ((tparticleOma->GetPdgCode() < 0) && TMath::Abs(tparticleMother->GetPdgCode())==pdgLambda){// check mother to be Lambda
\r
1182 Int_t labelFirstDaughter = tparticleMother->GetDaughter(1);// Proton
\r
1183 Int_t labelThirdDaughter = tparticleMother->GetDaughter(0);// Pion
\r
1185 TParticle *tparticleFirstDaughter = stack->Particle(TMath::Abs(labelFirstDaughter));
\r
1186 TParticle *tparticleThirdDaughter = stack->Particle(TMath::Abs(labelThirdDaughter));
\r
1188 if((tparticleFirstDaughter->GetPdgCode() == pdgProton && tparticleThirdDaughter->GetPdgCode()== pdgPionMinus) ||
\r
1189 (tparticleFirstDaughter->GetPdgCode() == pdgPionMinus && tparticleThirdDaughter->GetPdgCode()== pdgProton)){ //daughter PDGs
\r
1190 isPartialCorrectlyAssociatedLambda = kTRUE;
\r
1191 labelLambda=labelMother;
\r
1196 if(tparticleOma->GetPdgCode() == pdgHDibaryon){ //check grandmother to be H PDG
\r
1197 if (TMath::Abs(tparticleMother->GetPdgCode())==pdgLambda){// check mother to be Lambda
\r
1198 Int_t labelFirstDaughter = tparticleMother->GetDaughter(1);// Proton
\r
1199 Int_t labelThirdDaughter = tparticleMother->GetDaughter(0);// Pion
\r
1201 TParticle *tparticleFirstDaughter = stack->Particle(TMath::Abs(labelFirstDaughter));
\r
1202 TParticle *tparticleThirdDaughter = stack->Particle(TMath::Abs(labelThirdDaughter));
\r
1204 if((tparticleFirstDaughter->GetPdgCode() == pdgProton && tparticleThirdDaughter->GetPdgCode()== pdgPionMinus) ||
\r
1205 (tparticleFirstDaughter->GetPdgCode() == pdgPionMinus && tparticleThirdDaughter->GetPdgCode()== pdgProton)){ //daughter PDGs
\r
1206 isCorrectlyAssociatedLambda = kTRUE;
\r
1207 labelAssociatedH=labelOma;
\r
1208 fHistMassLambdaFromH->Fill(lInvMassLambda);
\r
1209 fHistMassLambdaFromHtLorentz->Fill(lambda.M());
\r
1215 fHistProtonPIDLambda->Fill(trackPosTest->GetInnerParam()->GetP(), trackPosTest->GetTPCsignal());
\r
1216 fHistPionPIDLambda->Fill(trackNegTest->GetInnerParam()->GetP(), trackNegTest->GetTPCsignal());
\r
1218 //---------------------------------------------------------
\r
1219 // Proton track loop
\r
1220 //---------------------------------------------------------
\r
1221 fHistArmenterosPodolanskiCut->Fill(alfa,qt);
\r
1223 for (Int_t iTracksP = 0; iTracksP < fESD->GetNumberOfTracks(); iTracksP++) {
\r
1224 AliESDtrack* trackP = dynamic_cast<AliESDtrack*> (fESD->GetTrack(iTracksP));
\r
1225 if (trackP->GetSign()<0) continue;
\r
1227 if (iTracksP==fV0MIs->GetPindex())continue;
\r
1228 if (iTracksP==fV0MIs->GetNindex())continue;
\r
1230 if (!fEsdTrackCuts->AcceptTrack(trackP)) continue;
\r
1232 AliKFParticle protonKF( *(trackP), 2212);
\r
1234 if (!trackP->GetInnerParam()) continue;
\r
1241 if (TMath::Abs(fPIDtpcESD->NumberOfSigmasTPC(trackP, AliPID::kProton)) > 3) continue;
\r
1244 fHistProtonPIDBb->Fill(trackP->GetInnerParam()->GetP(), trackP->GetTPCsignal());
\r
1246 posProt.SetXYZM(trackP->Px(),trackP->Py(),trackP->Pz(),cProtonMass);
\r
1248 //!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
\r
1249 //Pion Track loop!!!!!!!!!!!!!!!!!!!!!
\r
1250 //!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
\r
1251 for (Int_t iTracksN = iTracksP+1; iTracksN < fESD->GetNumberOfTracks(); iTracksN++) {
\r
1252 AliESDtrack* trackN = dynamic_cast<AliESDtrack*> (fESD->GetTrack(iTracksN));
\r
1254 if (iTracksN==fV0MIs->GetPindex())continue;
\r
1255 if (iTracksN==fV0MIs->GetNindex())continue;
\r
1256 if (trackN->GetSign()>0) continue;
\r
1258 if (!fEsdTrackCuts->AcceptTrack(trackN)) continue;
\r
1259 if (!fESDtrackCutsV0->AcceptTrack(trackN)) continue;
\r
1261 Double_t bz = fESD->GetMagneticField();
1263 Double_t xthiss(0.0);
1265 Double_t dca = trackN->GetDCA(trackP,bz,xthiss,xpp);
1266 if (dca>0.5) continue;
1268 negPi.SetXYZM(mn[0],mn[1],mn[2],cPionMass);
\r
1269 posP.SetXYZM(mp[0],mp[1],mp[2],cProtonMass);
\r
1270 negPio.SetXYZM(trackN->Px(),trackN->Py(),trackN->Pz(),cPionMass);
\r
1272 threeSum=negPi+posP+negPio;
\r
1273 lInvMassLambdaPi=threeSum.M();
\r
1275 fHistDC->Fill(decayLength*lInvMassLambda/lPLambda);
\r
1277 AliKFParticle posPionKF( *(trackN) ,-211);
\r
1279 if (!trackN->GetInnerParam()) continue;
\r
1285 if (TMath::Abs(fPIDtpcESD->NumberOfSigmasTPC(trackN, AliPID::kPion)) > 3) continue;
\r
1287 fHistPionPIDBb->Fill(trackN->GetInnerParam()->GetP(), trackN->GetTPCsignal());
\r
1289 trkArray1->AddAt(trackP,0);
\r
1290 trkArray1->AddAt(trackN,1);
\r
1292 vertexer1->SetVtxStart(esdVer1);
\r
1293 AliESDVertex *decayVertex1 = (AliESDVertex*)vertexer1->VertexForSelectedESDTracks(trkArray1);
\r
1295 dd1[0]=fESD->GetPrimaryVertexSPD()->GetX()-decayVertex1->GetX();
\r
1296 dd1[1]=fESD->GetPrimaryVertexSPD()->GetY()-decayVertex1->GetY();
\r
1297 dd1[2]=fESD->GetPrimaryVertexSPD()->GetZ()-decayVertex1->GetZ();
\r
1299 decayLengthH=sqrt(dd1[0]*dd1[0]+dd1[1]*dd1[1]+dd1[2]*dd1[2]);
\r
1301 // Double_t bz = fESD->GetMagneticField();
\r
1303 trackP->PropagateToDCA(decayVertex1, bz, 10);
\r
1304 trackN->PropagateToDCA(decayVertex1, bz, 10);
\r
1306 // Double_t xthiss(0.0);
\r
1307 // Double_t xpp(0.0);
\r
1308 // Double_t dca = trackN->GetDCA(trackP,bz,xthiss,xpp);
\r
1310 if (decayVertex1 == NULL) cout << "Secondary decay vtx pointer NULL" << endl;
1311 if (decayVertex1) delete decayVertex1;
\r
1312 h1.SetXYZ(-dd1[0],-dd1[1],-dd1[2]);
\r
1314 // if (dca>1) continue;
\r
1315 // if (dca>0.1) continue;
\r
1317 fourSum=threeSum+posProt;
\r
1319 posProt1.SetXYZM(trackP->Px(),trackP->Py(),trackP->Pz(),cProtonMass);
\r
1320 negPio1.SetXYZM(trackN->Px(),trackN->Py(),trackN->Pz(),cPionMass);
\r
1321 hDibaryon=lambdaH+posProt1+negPio1;
\r
1322 Double_t hPointingAngle = hDibaryon.Angle(h);
\r
1323 Double_t pointingAngleH = hDibaryon.Angle(h1);
\r
1324 Double_t decayAngleH = h.Angle(h1);
\r
1325 TVector3 vecDist(dd[0]-dd1[0],dd[1]-dd1[1],dd[2]-dd1[2]);
\r
1326 fHistMassLambdaPPi->Fill(hDibaryon.M());
\r
1327 fHistHPointingAngle->Fill(pointingAngleH);
\r
1329 fHistMassHcentMult->Fill(hDibaryon.M(),triggertype,refMultTpc);
\r
1331 Double_t rapidity = hDibaryon.Rapidity();
\r
1332 if(rapidity > 1.0 || rapidity < -1.0) continue;
\r
1334 //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
1335 //fHistNdim->Fill(vec);
\r
1337 fHistRap->Fill(rapidity);
\r
1338 //if (pointingAngleH > 0.1) continue;
\r
1339 if (pointingAngleH > 0.05) continue;
\r
1341 ///////////////////////////
\r
1342 //MC part for Associated H
\r
1343 ///////////////////////////
\r
1345 if (HasMC() && isCorrectlyAssociatedLambda) {
\r
1346 Int_t labelP = trackP->GetLabel();
\r
1347 TParticle *tparticleDaughter = stack->Particle(TMath::Abs(labelP));
\r
1348 Int_t labelMother = tparticleDaughter->GetFirstMother();
\r
1349 TParticle *tparticleMother = stack->Particle(TMath::Abs(labelMother));
\r
1351 Int_t labelProton = trackP->GetLabel();
\r
1352 Int_t labelPion = trackN->GetLabel();
\r
1355 if(tparticleMother->GetPdgCode() == pdgHDibaryon && labelAssociatedH==labelMother){ //check mother PDG
\r
1356 Int_t labelFirstDaughter = tparticleMother->GetDaughter(0);
\r
1357 Int_t labelThirdDaughter = tparticleMother->GetDaughter(1);
\r
1358 Int_t labelSecondDaughter = labelFirstDaughter +1;
\r
1360 TParticle *tparticleSecondDaughter = stack->Particle(TMath::Abs(labelSecondDaughter));
\r
1361 TParticle *tparticleThirdDaughter = stack->Particle(TMath::Abs(labelThirdDaughter));
\r
1363 TLorentzVector ppi;
\r
1364 TLorentzVector lpi;
\r
1365 TLorentzVector lP;
\r
1367 ppi=posProt+negPio;
\r
1368 lpi=lambdaH+negPio;
\r
1369 lP=lambdaH+posProt;
\r
1371 if((tparticleThirdDaughter->GetPdgCode() == pdgPionMinus && labelPion==labelThirdDaughter)&&(tparticleSecondDaughter->GetPdgCode() == pdgProton||tparticleSecondDaughter->GetPdgCode() == pdgPionMinus) && labelProton==labelSecondDaughter) fHistMassPpi->Fill(ppi.M());
\r
1373 if(tparticleThirdDaughter->GetPdgCode() == pdgPionMinus && labelPion==labelThirdDaughter) fHistMassLpi->Fill(lpi.M());
\r
1375 if(tparticleSecondDaughter->GetPdgCode() == pdgProton && labelProton==labelSecondDaughter) fHistMassLP->Fill(lP.M());
\r
1377 if(tparticleSecondDaughter->GetPdgCode() == pdgProton && labelProton==labelSecondDaughter){//check second daughter PDG
\r
1378 if(tparticleThirdDaughter->GetPdgCode() == pdgPionMinus && labelPion==labelThirdDaughter){//check second daughter PDG
\r
1380 fHistHDibaryonInvaMassAso->Fill(hDibaryon.M());
\r
1382 Double_t distance01=vecDist.Mag();
\r
1383 fHistMCdcaPvtxDvtx->Fill(decayLengthH);
\r
1384 fHistMCdcaPvtxLvtx->Fill(decayLength);
\r
1385 fHistMCdcaDvtxLvtx->Fill(distance01);
\r
1386 fHistMCangleLH->Fill(hPointingAngle);
\r
1387 fHistMCdecayAngle->Fill(decayAngleH);
\r
1388 fHistMCpointingAngle->Fill(pointingAngleH);
\r
1389 fHistMCap->Fill(alfa,qt);
\r
1391 fHistHilf1->Fill(posPionKF.GetDistanceFromVertex(primVtx));
\r
1392 fHistHilf2->Fill(protonKF.GetDistanceFromVertex(primVtx));
\r
1393 fHistHilf3->Fill(protonKF.GetDistanceFromVertex(posPionKF));
\r
1394 fHistHilf6->Fill(dca);
\r
1395 fHistPtvsYAso->Fill(hDibaryon.Pt(),hDibaryon.Rapidity());
\r
1396 fHistPtvsEtaAso->Fill(hDibaryon.Pt(),hDibaryon.Eta());
\r
1398 }//end check for third daughter PDG
\r
1399 }//end check second daughter PDG
\r
1403 // cout<<"Trigger: "<<triggertype<<endl;
\r
1404 fHistMassH->Fill(hDibaryon.M());
\r
1405 fHistMassHcentMult->Fill(hDibaryon.M(),triggertype,refMultTpc);
\r
1406 ppK=lambdaH+posProt;
\r
1407 fHistMassLambdaP->Fill(ppK.M());
\r
1409 //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
1411 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, static_cast<Double_t>(refMultTpc), static_cast<Double_t>(mcStatus)};
\r
1412 fHistNdim->Fill(vec);
\r
1420 //!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
\r
1421 //Pure MC Part!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
\r
1422 //!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
\r
1424 // Monte Carlo for genenerated particles
\r
1425 if (HasMC()) //MC loop
\r
1430 Double_t momentumPionGen[3]={0,0,0};
\r
1431 Double_t momentumNucleonGen[3]={0,0,0};
\r
1432 Double_t momentumLambdaGen[3]={0,0,0};
\r
1434 Double_t energyPionGen = 0;
\r
1435 Double_t energyNucleonGen = 0;
\r
1436 Double_t energyLambdaGen = 0;
\r
1438 Double_t transversMomentumMotherGen = 0;
\r
1439 Double_t longitudinalMomentumMotherGen = 0;
\r
1440 Double_t totalEnergyMotherGen = 0;
\r
1442 Double_t rapidityGen = 2;
\r
1444 for(stackN = 0; stackN < stack->GetNtrack(); stackN++) //loop over stack
\r
1447 TParticle *tparticleMother = stack->Particle(stackN);
\r
1449 if(tparticleMother->GetPdgCode() == pdgLambda) fHistCount->Fill(16);
\r
1452 if(tparticleMother->GetPdgCode() == pdgHDibaryon) //check mother PDG
\r
1454 Int_t labelFirstDaughter = tparticleMother->GetDaughter(0);
\r
1455 Int_t labelThirdDaughter = tparticleMother->GetDaughter(1);
\r
1456 Int_t labelSecondDaughter = labelFirstDaughter +1;
\r
1458 TParticle *tparticleFirstDaughter = stack->Particle(TMath::Abs(labelFirstDaughter));
\r
1459 TParticle *tparticleSecondDaughter = stack->Particle(TMath::Abs(labelSecondDaughter));
\r
1460 TParticle *tparticleThirdDaughter = stack->Particle(TMath::Abs(labelThirdDaughter));
\r
1462 if(tparticleFirstDaughter->GetPdgCode() == pdgLambda) //check first daughter PDG
\r
1464 if(tparticleSecondDaughter->GetPdgCode() == pdgProton)//check second daughter PDG
\r
1466 if(tparticleThirdDaughter->GetPdgCode() == pdgPionMinus)//check second daughter PDG
\r
1468 momentumLambdaGen[0] = tparticleFirstDaughter->Px();
\r
1469 momentumLambdaGen[1] = tparticleFirstDaughter->Py();
\r
1470 momentumLambdaGen[2] = tparticleFirstDaughter->Pz();
\r
1472 momentumNucleonGen[0] = tparticleSecondDaughter->Px();
\r
1473 momentumNucleonGen[1] = tparticleSecondDaughter->Py();
\r
1474 momentumNucleonGen[2] = tparticleSecondDaughter->Pz();
\r
1476 momentumPionGen[0] = tparticleThirdDaughter->Px();
\r
1477 momentumPionGen[1] = tparticleThirdDaughter->Py();
\r
1478 momentumPionGen[2] = tparticleThirdDaughter->Pz();
\r
1480 TLorentzVector lorentzVectorLambda;
\r
1481 TLorentzVector lorentzVectorProton;
\r
1482 TLorentzVector lorentzVectorPion;
\r
1483 TLorentzVector lorentzVectorHDibaryon;
\r
1485 lorentzVectorLambda.SetXYZM(momentumLambdaGen[0],momentumLambdaGen[1],momentumLambdaGen[2],1.115);
\r
1486 lorentzVectorProton.SetXYZM(momentumNucleonGen[0],momentumNucleonGen[1],momentumNucleonGen[2],protonMass);
\r
1487 lorentzVectorPion.SetXYZM(momentumPionGen[0],momentumPionGen[1],momentumPionGen[2],pionMass);
\r
1489 lorentzVectorHDibaryon = lorentzVectorLambda + lorentzVectorProton + lorentzVectorPion;
\r
1490 rapidityGen=lorentzVectorHDibaryon.Rapidity();
\r
1491 transversMomentumMotherGen = lorentzVectorHDibaryon.Pt();
\r
1492 longitudinalMomentumMotherGen = lorentzVectorHDibaryon.Pz();
\r
1493 totalEnergyMotherGen = lorentzVectorHDibaryon.Energy();
\r
1495 if(rapidityGen > 1.0 || rapidityGen < -1 ) continue;
\r
1496 //lorentzVectorLambda
\r
1497 fHistHDibaryonInvaMassGen->Fill(lorentzVectorHDibaryon.M());
\r
1498 if (lorentzVectorLambda.Rapidity() > 1.0 || lorentzVectorLambda.Rapidity() < -1) continue;
\r
1499 if (lorentzVectorProton.Rapidity() > 1.0 || lorentzVectorProton.Rapidity() < -1) continue;
\r
1501 if (lorentzVectorPion.Rapidity() > 1.0 || lorentzVectorPion.Rapidity() < -1) continue;
\r
1502 fHistHDibaryonInvaMassGenRes->Fill(lorentzVectorHDibaryon.M());
\r
1503 fHistPtvsEtaGen->Fill(lorentzVectorHDibaryon.Pt(),lorentzVectorHDibaryon.Eta());
\r
1504 fHistPtvsYGen->Fill(lorentzVectorHDibaryon.Pt(),lorentzVectorHDibaryon.Rapidity());
\r
1505 fHistPtvsEtaGen->Fill(lorentzVectorHDibaryon.Pt(),lorentzVectorHDibaryon.Eta());
\r
1506 fHistCount->Fill(11);
\r
1507 }//end of check third daughter PDG
\r
1508 }//end of check second daughter PDG
\r
1509 }//end of check first daughter PDG
\r
1510 }//end of H-Dibaryon
\r
1513 if(tparticleMother->GetPdgCode() == pdgAntiHDibaryon) //check mother PDG
\r
1515 Int_t labelFirstDaughter = tparticleMother->GetDaughter(0);
\r
1516 Int_t labelThirdDaughter = tparticleMother->GetDaughter(1);
\r
1517 Int_t labelSecondDaughter = labelFirstDaughter +1;
\r
1519 TParticle *tparticleFirstDaughter = stack->Particle(TMath::Abs(labelFirstDaughter));
\r
1520 TParticle *tparticleSecondDaughter = stack->Particle(TMath::Abs(labelSecondDaughter));
\r
1521 TParticle *tparticleThirdDaughter = stack->Particle(TMath::Abs(labelThirdDaughter));
\r
1523 if(tparticleFirstDaughter->GetPdgCode() == pdgAntiLambda) //check first daughter PDG
\r
1525 if(tparticleSecondDaughter->GetPdgCode() == pdgAntiProton)//check second daughter PDG
\r
1527 if(tparticleThirdDaughter->GetPdgCode() == pdgPionPlus)//check second daughter PDG
\r
1529 momentumLambdaGen[0] = tparticleFirstDaughter->Px();
\r
1530 momentumLambdaGen[1] = tparticleFirstDaughter->Py();
\r
1531 momentumLambdaGen[2] = tparticleFirstDaughter->Pz();
\r
1533 momentumNucleonGen[0] = tparticleSecondDaughter->Px();
\r
1534 momentumNucleonGen[1] = tparticleSecondDaughter->Py();
\r
1535 momentumNucleonGen[2] = tparticleSecondDaughter->Pz();
\r
1537 momentumPionGen[0] = tparticleThirdDaughter->Px();
\r
1538 momentumPionGen[1] = tparticleThirdDaughter->Py();
\r
1539 momentumPionGen[2] = tparticleThirdDaughter->Pz();
\r
1541 energyLambdaGen = tparticleFirstDaughter->Energy();
\r
1542 energyNucleonGen = tparticleSecondDaughter->Energy();
\r
1543 energyPionGen = tparticleThirdDaughter->Energy();
\r
1545 TLorentzVector lorentzVectorLambda;
\r
1546 TLorentzVector lorentzVectorProton;
\r
1547 TLorentzVector lorentzVectorPion;
\r
1548 TLorentzVector lorentzVectorHDibaryon;
\r
1550 lorentzVectorLambda.SetXYZM(momentumLambdaGen[0],momentumLambdaGen[1],momentumLambdaGen[2],1.115);
\r
1551 lorentzVectorProton.SetXYZM(momentumNucleonGen[0],momentumNucleonGen[1],momentumNucleonGen[2],protonMass);
\r
1552 lorentzVectorPion.SetXYZM(momentumPionGen[0],momentumPionGen[1],momentumPionGen[2],pionMass);
\r
1554 lorentzVectorHDibaryon = lorentzVectorLambda + lorentzVectorProton + lorentzVectorPion;
\r
1556 rapidityGen=lorentzVectorHDibaryon.Rapidity();
\r
1557 if(rapidityGen > 1.0 || rapidityGen < -1 ) continue;
\r
1558 fHistAntiHDibaryonInvaMassGen->Fill(lorentzVectorHDibaryon.M());
\r
1559 }//end of check third daughter PDG
\r
1560 }//end of check second daughter PDG
\r
1561 }//end of check first daughter PDG
\r
1562 }//end of Anti-H-Dibaryon
\r
1566 // Post output data.
\r
1567 PostData(1,fHistList);
\r
1568 //PostData(0,fHistList);
\r
1570 if (listCrossV0 == NULL) return;
1572 if (listCrossV0) delete listCrossV0;
\r
1573 if (esdVer1) delete esdVer1;
\r
1574 if (vertexer) delete vertexer;
\r
1575 if (vertexer1) delete vertexer1;
\r
1576 if (trkArray) delete trkArray;
\r
1577 if (trkArray1) delete trkArray1;
\r
1580 //________________________________________________________________________
\r
1581 void AliAnalysisTaskHdibaryonLPpi::Terminate(Option_t *)
\r
1583 // Draw result to the screen
\r
1584 // Called once at the end of the query
\r