1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
16 // +-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+
18 // Modified version of AliAnalysisTaskCheckCascade.cxx.
19 // This is a 'hybrid' output version, in that it uses a classic TTree
20 // ROOT object to store the candidates, plus a couple of histograms filled on
21 // a per-event basis for storing variables too numerous to put in a tree.
24 // --- Algorithm Description
25 // 1. Loop over primaries in stack to acquire generated charged Xi
26 // 2. Loop over stack to find Cascades, fill TH3Fs "PrimRawPt"s for Efficiency
27 // 3. Perform Physics Selection
28 // 4. Perform Primary Vertex |z|<10cm selection
29 // 5. Perform Primary Vertex NoTPCOnly vertexing selection
30 // 6. Perform Pileup Rejection
32 // 7a. Fill TH3Fs "PrimAnalysisPt" for control purposes only
34 // Please Report Any Bugs!
36 // --- David Dobrigkeit Chinellato
37 // (david.chinellato@gmail.com)
39 // +-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+
45 //class AliMCEventHandler;
54 #include <Riostream.h>
60 #include "THnSparse.h"
67 #include "AliESDEvent.h"
68 #include "AliAODEvent.h"
69 #include "AliV0vertexer.h"
70 #include "AliCascadeVertexer.h"
71 #include "AliESDpid.h"
72 #include "AliESDtrack.h"
73 #include "AliESDtrackCuts.h"
74 #include "AliInputEventHandler.h"
75 #include "AliAnalysisManager.h"
76 #include "AliMCEventHandler.h"
77 #include "AliMCEvent.h"
80 #include "AliV0vertexer.h"
81 #include "AliCascadeVertexer.h"
83 #include "AliCFContainer.h"
84 #include "AliMultiplicity.h"
85 #include "AliAODMCParticle.h"
86 #include "AliESDcascade.h"
87 #include "AliAODcascade.h"
88 #include "AliESDUtils.h"
89 #include "AliGenEventHeader.h"
90 #include "AliAnalysisUtils.h"
92 #include "AliAnalysisTaskExtractCascade.h"
97 ClassImp(AliAnalysisTaskExtractCascade)
99 AliAnalysisTaskExtractCascade::AliAnalysisTaskExtractCascade()
100 : AliAnalysisTaskSE(), fListHist(0), fTreeCascade(0), fPIDResponse(0), fESDtrackCuts(0), fUtils(0),
101 fkIsNuclear ( kFALSE ),
102 fkSwitchINT7 ( kFALSE ),
103 fCentralityEstimator("V0M"),
104 fkpAVertexSelection( kFALSE ),
106 fkRunVertexers ( kFALSE ),
107 //------------------------------------------------
109 //------------------------------------------------
111 fTreeCascVarCharge(0),
112 fTreeCascVarMassAsXi(0),
113 fTreeCascVarMassAsOmega(0),
116 fTreeCascVarRapMC(0),
117 fTreeCascVarRapXi(0),
118 fTreeCascVarRapOmega(0),
119 fTreeCascVarNegEta(0),
120 fTreeCascVarPosEta(0),
121 fTreeCascVarBachEta(0),
122 fTreeCascVarDCACascDaughters(0),
123 fTreeCascVarDCABachToPrimVtx(0),
124 fTreeCascVarDCAV0Daughters(0),
125 fTreeCascVarDCAV0ToPrimVtx(0),
126 fTreeCascVarDCAPosToPrimVtx(0),
127 fTreeCascVarDCANegToPrimVtx(0),
128 fTreeCascVarCascCosPointingAngle(0),
129 fTreeCascVarCascRadius(0),
130 fTreeCascVarV0Mass(0),
131 fTreeCascVarV0CosPointingAngle(0),
132 fTreeCascVarV0CosPointingAngleSpecial(0),
133 fTreeCascVarV0Radius(0),
134 fTreeCascVarLeastNbrClusters(0),
135 fTreeCascVarMultiplicity(0),
136 fTreeCascVarMultiplicityV0A(0),
137 fTreeCascVarMultiplicityZNA(0),
138 fTreeCascVarMultiplicityTRK(0),
139 fTreeCascVarMultiplicitySPD(0),
140 fTreeCascVarDistOverTotMom(0),
142 fTreeCascVarPIDBachelor(0),
143 fTreeCascVarPIDNegative(0),
144 fTreeCascVarPIDPositive(0),
145 fTreeCascVarBachTransMom(0),
146 fTreeCascVarPosTransMom(0),
147 fTreeCascVarNegTransMom(0),
148 fTreeCascVarPosTransMomMC(0),
149 fTreeCascVarNegTransMomMC(0),
150 fTreeCascVarNegNSigmaPion(0),
151 fTreeCascVarNegNSigmaProton(0),
152 fTreeCascVarPosNSigmaPion(0),
153 fTreeCascVarPosNSigmaProton(0),
154 fTreeCascVarBachNSigmaPion(0),
155 fTreeCascVarBachNSigmaKaon(0),
157 fTreeCascVarkITSRefitBachelor(0),
158 fTreeCascVarkITSRefitNegative(0),
159 fTreeCascVarkITSRefitPositive(0),
161 //------------------------------------------------
163 // --- Filled on an Event-by-event basis
164 //------------------------------------------------
165 fHistV0MultiplicityBeforeTrigSel(0),
166 fHistV0MultiplicityForTrigEvt(0),
167 fHistV0MultiplicityForSelEvt(0),
168 fHistV0MultiplicityForSelEvtNoTPCOnly(0),
169 fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup(0),
170 fHistMultiplicityBeforeTrigSel(0),
171 fHistMultiplicityForTrigEvt(0),
172 fHistMultiplicity(0),
173 fHistMultiplicityNoTPCOnly(0),
174 fHistMultiplicityNoTPCOnlyNoPileup(0),
177 fHistMultiplicityV0ABeforeTrigSel(0),
178 fHistMultiplicityV0AForTrigEvt(0),
179 fHistMultiplicityV0A(0),
180 fHistMultiplicityV0ANoTPCOnly(0),
181 fHistMultiplicityV0ANoTPCOnlyNoPileup(0),
184 fHistMultiplicityZNABeforeTrigSel(0),
185 fHistMultiplicityZNAForTrigEvt(0),
186 fHistMultiplicityZNA(0),
187 fHistMultiplicityZNANoTPCOnly(0),
188 fHistMultiplicityZNANoTPCOnlyNoPileup(0),
191 fHistMultiplicityTRKBeforeTrigSel(0),
192 fHistMultiplicityTRKForTrigEvt(0),
193 fHistMultiplicityTRK(0),
194 fHistMultiplicityTRKNoTPCOnly(0),
195 fHistMultiplicityTRKNoTPCOnlyNoPileup(0),
198 fHistMultiplicitySPDBeforeTrigSel(0),
199 fHistMultiplicitySPDForTrigEvt(0),
200 fHistMultiplicitySPD(0),
201 fHistMultiplicitySPDNoTPCOnly(0),
202 fHistMultiplicitySPDNoTPCOnlyNoPileup(0),
214 AliAnalysisTaskExtractCascade::AliAnalysisTaskExtractCascade(const char *name)
215 : AliAnalysisTaskSE(name), fListHist(0), fTreeCascade(0), fPIDResponse(0), fESDtrackCuts(0), fUtils(0),
216 fkIsNuclear ( kFALSE ),
217 fkSwitchINT7 ( kFALSE ),
218 fCentralityEstimator("V0M"),
219 fkpAVertexSelection( kFALSE ),
221 fkRunVertexers ( kFALSE ),
222 //------------------------------------------------
224 //------------------------------------------------
226 fTreeCascVarCharge(0),
227 fTreeCascVarMassAsXi(0),
228 fTreeCascVarMassAsOmega(0),
231 fTreeCascVarRapMC(0),
232 fTreeCascVarRapXi(0),
233 fTreeCascVarRapOmega(0),
234 fTreeCascVarNegEta(0),
235 fTreeCascVarPosEta(0),
236 fTreeCascVarBachEta(0),
237 fTreeCascVarDCACascDaughters(0),
238 fTreeCascVarDCABachToPrimVtx(0),
239 fTreeCascVarDCAV0Daughters(0),
240 fTreeCascVarDCAV0ToPrimVtx(0),
241 fTreeCascVarDCAPosToPrimVtx(0),
242 fTreeCascVarDCANegToPrimVtx(0),
243 fTreeCascVarCascCosPointingAngle(0),
244 fTreeCascVarCascRadius(0),
245 fTreeCascVarV0Mass(0),
246 fTreeCascVarV0CosPointingAngle(0),
247 fTreeCascVarV0CosPointingAngleSpecial(0),
248 fTreeCascVarV0Radius(0),
249 fTreeCascVarLeastNbrClusters(0),
250 fTreeCascVarMultiplicity(0),
251 fTreeCascVarMultiplicityV0A(0),
252 fTreeCascVarMultiplicityZNA(0),
253 fTreeCascVarMultiplicityTRK(0),
254 fTreeCascVarMultiplicitySPD(0),
255 fTreeCascVarDistOverTotMom(0),
257 fTreeCascVarPIDBachelor(0),
258 fTreeCascVarPIDNegative(0),
259 fTreeCascVarPIDPositive(0),
260 fTreeCascVarBachTransMom(0),
261 fTreeCascVarPosTransMom(0),
262 fTreeCascVarNegTransMom(0),
263 fTreeCascVarPosTransMomMC(0),
264 fTreeCascVarNegTransMomMC(0),
265 fTreeCascVarNegNSigmaPion(0),
266 fTreeCascVarNegNSigmaProton(0),
267 fTreeCascVarPosNSigmaPion(0),
268 fTreeCascVarPosNSigmaProton(0),
269 fTreeCascVarBachNSigmaPion(0),
270 fTreeCascVarBachNSigmaKaon(0),
272 fTreeCascVarkITSRefitBachelor(0),
273 fTreeCascVarkITSRefitNegative(0),
274 fTreeCascVarkITSRefitPositive(0),
276 //------------------------------------------------
278 // --- Filled on an Event-by-event basis
279 //------------------------------------------------
280 fHistV0MultiplicityBeforeTrigSel(0),
281 fHistV0MultiplicityForTrigEvt(0),
282 fHistV0MultiplicityForSelEvt(0),
283 fHistV0MultiplicityForSelEvtNoTPCOnly(0),
284 fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup(0),
285 fHistMultiplicityBeforeTrigSel(0),
286 fHistMultiplicityForTrigEvt(0),
287 fHistMultiplicity(0),
288 fHistMultiplicityNoTPCOnly(0),
289 fHistMultiplicityNoTPCOnlyNoPileup(0),
292 fHistMultiplicityV0ABeforeTrigSel(0),
293 fHistMultiplicityV0AForTrigEvt(0),
294 fHistMultiplicityV0A(0),
295 fHistMultiplicityV0ANoTPCOnly(0),
296 fHistMultiplicityV0ANoTPCOnlyNoPileup(0),
299 fHistMultiplicityZNABeforeTrigSel(0),
300 fHistMultiplicityZNAForTrigEvt(0),
301 fHistMultiplicityZNA(0),
302 fHistMultiplicityZNANoTPCOnly(0),
303 fHistMultiplicityZNANoTPCOnlyNoPileup(0),
306 fHistMultiplicityTRKBeforeTrigSel(0),
307 fHistMultiplicityTRKForTrigEvt(0),
308 fHistMultiplicityTRK(0),
309 fHistMultiplicityTRKNoTPCOnly(0),
310 fHistMultiplicityTRKNoTPCOnlyNoPileup(0),
313 fHistMultiplicitySPDBeforeTrigSel(0),
314 fHistMultiplicitySPDForTrigEvt(0),
315 fHistMultiplicitySPD(0),
316 fHistMultiplicitySPDNoTPCOnly(0),
317 fHistMultiplicitySPDNoTPCOnlyNoPileup(0),
328 //Set Variables for re-running the cascade vertexers (as done for MS paper)
330 // New Loose : 1st step for the 7 TeV pp analysis
332 fV0VertexerSels[0] = 33. ; // max allowed chi2
333 fV0VertexerSels[1] = 0.02; // min allowed impact parameter for the 1st daughter (LHC09a4 : 0.05)
334 fV0VertexerSels[2] = 0.02; // min allowed impact parameter for the 2nd daughter (LHC09a4 : 0.05)
335 fV0VertexerSels[3] = 2.0 ; // max allowed DCA between the daughter tracks (LHC09a4 : 0.5)
336 fV0VertexerSels[4] = 0.95; // min allowed cosine of V0's pointing angle (LHC09a4 : 0.99)
337 fV0VertexerSels[5] = 1.0 ; // min radius of the fiducial volume (LHC09a4 : 0.2)
338 fV0VertexerSels[6] = 200. ; // max radius of the fiducial volume (LHC09a4 : 100.0)
340 fCascadeVertexerSels[0] = 33. ; // max allowed chi2 (same as PDC07)
341 fCascadeVertexerSels[1] = 0.05 ; // min allowed V0 impact parameter (PDC07 : 0.05 / LHC09a4 : 0.025 )
342 fCascadeVertexerSels[2] = 0.010; // "window" around the Lambda mass (PDC07 : 0.008 / LHC09a4 : 0.010 )
343 fCascadeVertexerSels[3] = 0.03 ; // min allowed bachelor's impact parameter (PDC07 : 0.035 / LHC09a4 : 0.025 )
344 fCascadeVertexerSels[4] = 2.0 ; // max allowed DCA between the V0 and the bachelor (PDC07 : 0.1 / LHC09a4 : 0.2 )
345 fCascadeVertexerSels[5] = 0.95 ; // min allowed cosine of the cascade pointing angle (PDC07 : 0.9985 / LHC09a4 : 0.998 )
346 fCascadeVertexerSels[6] = 0.4 ; // min radius of the fiducial volume (PDC07 : 0.9 / LHC09a4 : 0.2 )
347 fCascadeVertexerSels[7] = 100. ; // max radius of the fiducial volume (PDC07 : 100 / LHC09a4 : 100 )
349 // Output slot #0 writes into a TList container (Cascade)
350 DefineOutput(1, TList::Class());
351 DefineOutput(2, TTree::Class());
355 AliAnalysisTaskExtractCascade::~AliAnalysisTaskExtractCascade()
357 //------------------------------------------------
359 //------------------------------------------------
369 //cleanup esd track cuts object too...
371 delete fESDtrackCuts;
381 //________________________________________________________________________
382 void AliAnalysisTaskExtractCascade::UserCreateOutputObjects()
387 //------------------------------------------------
389 fTreeCascade = new TTree("fTreeCascade","CascadeCandidates");
391 //------------------------------------------------
392 // fTreeCascade Branch definitions - Cascade Tree
393 //------------------------------------------------
395 //------------------------------------------------
396 // fTreeCascade Branch definitions
397 //------------------------------------------------
399 //-----------BASIC-INFO---------------------------
400 /* 1*/ fTreeCascade->Branch("fTreeCascVarCharge",&fTreeCascVarCharge,"fTreeCascVarCharge/I");
401 /* 2*/ fTreeCascade->Branch("fTreeCascVarMassAsXi",&fTreeCascVarMassAsXi,"fTreeCascVarMassAsXi/F");
402 /* 3*/ fTreeCascade->Branch("fTreeCascVarMassAsOmega",&fTreeCascVarMassAsOmega,"fTreeCascVarMassAsOmega/F");
403 /* 4*/ fTreeCascade->Branch("fTreeCascVarPt",&fTreeCascVarPt,"fTreeCascVarPt/F");
404 /* 5*/ fTreeCascade->Branch("fTreeCascVarRapXi",&fTreeCascVarRapXi,"fTreeCascVarRapXi/F");
405 /* 6*/ fTreeCascade->Branch("fTreeCascVarRapOmega",&fTreeCascVarRapOmega,"fTreeCascVarRapOmega/F");
406 /* 7*/ fTreeCascade->Branch("fTreeCascVarNegEta",&fTreeCascVarNegEta,"fTreeCascVarNegEta/F");
407 /* 8*/ fTreeCascade->Branch("fTreeCascVarPosEta",&fTreeCascVarPosEta,"fTreeCascVarPosEta/F");
408 /* 9*/ fTreeCascade->Branch("fTreeCascVarBachEta",&fTreeCascVarBachEta,"fTreeCascVarBachEta/F");
409 //-----------INFO-FOR-CUTS------------------------
410 /*10*/ fTreeCascade->Branch("fTreeCascVarDCACascDaughters",&fTreeCascVarDCACascDaughters,"fTreeCascVarDCACascDaughters/F");
411 /*11*/ fTreeCascade->Branch("fTreeCascVarDCABachToPrimVtx",&fTreeCascVarDCABachToPrimVtx,"fTreeCascVarDCABachToPrimVtx/F");
412 /*12*/ fTreeCascade->Branch("fTreeCascVarDCAV0Daughters",&fTreeCascVarDCAV0Daughters,"fTreeCascVarDCAV0Daughters/F");
413 /*13*/ fTreeCascade->Branch("fTreeCascVarDCAV0ToPrimVtx",&fTreeCascVarDCAV0ToPrimVtx,"fTreeCascVarDCAV0ToPrimVtx/F");
414 /*14*/ fTreeCascade->Branch("fTreeCascVarDCAPosToPrimVtx",&fTreeCascVarDCAPosToPrimVtx,"fTreeCascVarDCAPosToPrimVtx/F");
415 /*15*/ fTreeCascade->Branch("fTreeCascVarDCANegToPrimVtx",&fTreeCascVarDCANegToPrimVtx,"fTreeCascVarDCANegToPrimVtx/F");
416 /*16*/ fTreeCascade->Branch("fTreeCascVarCascCosPointingAngle",&fTreeCascVarCascCosPointingAngle,"fTreeCascVarCascCosPointingAngle/F");
417 /*17*/ fTreeCascade->Branch("fTreeCascVarCascRadius",&fTreeCascVarCascRadius,"fTreeCascVarCascRadius/F");
418 /*18*/ fTreeCascade->Branch("fTreeCascVarV0Mass",&fTreeCascVarV0Mass,"fTreeCascVarV0Mass/F");
419 /*19*/ fTreeCascade->Branch("fTreeCascVarV0CosPointingAngle",&fTreeCascVarV0CosPointingAngle,"fTreeCascVarV0CosPointingAngle/F");
420 /*19*/ fTreeCascade->Branch("fTreeCascVarV0CosPointingAngleSpecial",&fTreeCascVarV0CosPointingAngleSpecial,"fTreeCascVarV0CosPointingAngleSpecial/F");
421 /*20*/ fTreeCascade->Branch("fTreeCascVarV0Radius",&fTreeCascVarV0Radius,"fTreeCascVarV0Radius/F");
422 /*21*/ fTreeCascade->Branch("fTreeCascVarLeastNbrClusters",&fTreeCascVarLeastNbrClusters,"fTreeCascVarLeastNbrClusters/I");
423 //-----------MULTIPLICITY-INFO--------------------
424 /*22*/ fTreeCascade->Branch("fTreeCascVarMultiplicity",&fTreeCascVarMultiplicity,"fTreeCascVarMultiplicity/I");
425 /*22*/ fTreeCascade->Branch("fTreeCascVarMultiplicityV0A",&fTreeCascVarMultiplicityV0A,"fTreeCascVarMultiplicityV0A/I");
426 /*22*/ fTreeCascade->Branch("fTreeCascVarMultiplicityZNA",&fTreeCascVarMultiplicityZNA,"fTreeCascVarMultiplicityZNA/I");
427 /*22*/ fTreeCascade->Branch("fTreeCascVarMultiplicityTRK",&fTreeCascVarMultiplicityTRK,"fTreeCascVarMultiplicityTRK/I");
428 /*22*/ fTreeCascade->Branch("fTreeCascVarMultiplicitySPD",&fTreeCascVarMultiplicitySPD,"fTreeCascVarMultiplicitySPD/I");
429 //-----------DECAY-LENGTH-INFO--------------------
430 /*23*/ fTreeCascade->Branch("fTreeCascVarDistOverTotMom",&fTreeCascVarDistOverTotMom,"fTreeCascVarDistOverTotMom/F");
431 //------------------------------------------------
432 /*24*/ fTreeCascade->Branch("fTreeCascVarNegNSigmaPion",&fTreeCascVarNegNSigmaPion,"fTreeCascVarNegNSigmaPion/F");
433 /*25*/ fTreeCascade->Branch("fTreeCascVarNegNSigmaProton",&fTreeCascVarNegNSigmaProton,"fTreeCascVarNegNSigmaProton/F");
434 /*26*/ fTreeCascade->Branch("fTreeCascVarPosNSigmaPion",&fTreeCascVarPosNSigmaPion,"fTreeCascVarPosNSigmaPion/F");
435 /*27*/ fTreeCascade->Branch("fTreeCascVarPosNSigmaProton",&fTreeCascVarPosNSigmaProton,"fTreeCascVarPosNSigmaProton/F");
436 /*28*/ fTreeCascade->Branch("fTreeCascVarBachNSigmaPion",&fTreeCascVarBachNSigmaPion,"fTreeCascVarBachNSigmaPion/F");
437 /*29*/ fTreeCascade->Branch("fTreeCascVarBachNSigmaKaon",&fTreeCascVarBachNSigmaKaon,"fTreeCascVarBachNSigmaKaon/F");
439 /*30*/ fTreeCascade->Branch("fTreeCascVarBachTransMom",&fTreeCascVarBachTransMom,"fTreeCascVarBachTransMom/F");
440 /*30*/ fTreeCascade->Branch("fTreeCascVarPosTransMom",&fTreeCascVarPosTransMom,"fTreeCascVarPosTransMom/F");
441 /*31*/ fTreeCascade->Branch("fTreeCascVarNegTransMom",&fTreeCascVarNegTransMom,"fTreeCascVarNegTransMom/F");
443 /*29*/ fTreeCascade->Branch("fTreeCascVarkITSRefitBachelor",&fTreeCascVarkITSRefitBachelor,"fTreeCascVarkITSRefitBachelor/O");
444 /*29*/ fTreeCascade->Branch("fTreeCascVarkITSRefitNegative",&fTreeCascVarkITSRefitNegative,"fTreeCascVarkITSRefitNegative/O");
445 /*29*/ fTreeCascade->Branch("fTreeCascVarkITSRefitPositive",&fTreeCascVarkITSRefitPositive,"fTreeCascVarkITSRefitPositive/O");
447 //------------------------------------------------
448 // Particle Identification Setup
449 //------------------------------------------------
451 AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
452 AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
453 fPIDResponse = inputHandler->GetPIDResponse();
457 if(! fESDtrackCuts ){
458 fESDtrackCuts = new AliESDtrackCuts();
461 fUtils = new AliAnalysisUtils();
464 //------------------------------------------------
465 // V0 Multiplicity Histograms
466 //------------------------------------------------
470 fListHist = new TList();
471 fListHist->SetOwner(); // See http://root.cern.ch/root/html/TCollection.html#TCollection:SetOwner
474 if(! fHistV0MultiplicityBeforeTrigSel) {
475 fHistV0MultiplicityBeforeTrigSel = new TH1F("fHistV0MultiplicityBeforeTrigSel",
476 "V0s per event (before Trig. Sel.);Nbr of V0s/Evt;Events",
478 fListHist->Add(fHistV0MultiplicityBeforeTrigSel);
481 if(! fHistV0MultiplicityForTrigEvt) {
482 fHistV0MultiplicityForTrigEvt = new TH1F("fHistV0MultiplicityForTrigEvt",
483 "V0s per event (for triggered evt);Nbr of V0s/Evt;Events",
485 fListHist->Add(fHistV0MultiplicityForTrigEvt);
488 if(! fHistV0MultiplicityForSelEvt) {
489 fHistV0MultiplicityForSelEvt = new TH1F("fHistV0MultiplicityForSelEvt",
490 "V0s per event;Nbr of V0s/Evt;Events",
492 fListHist->Add(fHistV0MultiplicityForSelEvt);
495 if(! fHistV0MultiplicityForSelEvtNoTPCOnly) {
496 fHistV0MultiplicityForSelEvtNoTPCOnly = new TH1F("fHistV0MultiplicityForSelEvtNoTPCOnly",
497 "V0s per event;Nbr of V0s/Evt;Events",
499 fListHist->Add(fHistV0MultiplicityForSelEvtNoTPCOnly);
501 if(! fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup) {
502 fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup = new TH1F("fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup",
503 "V0s per event;Nbr of V0s/Evt;Events",
505 fListHist->Add(fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup);
508 //------------------------------------------------
509 // Track Multiplicity Histograms
510 //------------------------------------------------
512 if(! fHistMultiplicityBeforeTrigSel) {
513 fHistMultiplicityBeforeTrigSel = new TH1F("fHistMultiplicityBeforeTrigSel",
514 "Tracks per event;Nbr of Tracks;Events",
516 fListHist->Add(fHistMultiplicityBeforeTrigSel);
518 if(! fHistMultiplicityForTrigEvt) {
519 fHistMultiplicityForTrigEvt = new TH1F("fHistMultiplicityForTrigEvt",
520 "Tracks per event;Nbr of Tracks;Events",
522 fListHist->Add(fHistMultiplicityForTrigEvt);
524 if(! fHistMultiplicity) {
525 fHistMultiplicity = new TH1F("fHistMultiplicity",
526 "Tracks per event;Nbr of Tracks;Events",
528 fListHist->Add(fHistMultiplicity);
530 if(! fHistMultiplicityNoTPCOnly) {
531 fHistMultiplicityNoTPCOnly = new TH1F("fHistMultiplicityNoTPCOnly",
532 "Tracks per event;Nbr of Tracks;Events",
534 fListHist->Add(fHistMultiplicityNoTPCOnly);
536 if(! fHistMultiplicityNoTPCOnlyNoPileup) {
537 fHistMultiplicityNoTPCOnlyNoPileup = new TH1F("fHistMultiplicityNoTPCOnlyNoPileup",
538 "Tracks per event;Nbr of Tracks;Events",
540 fListHist->Add(fHistMultiplicityNoTPCOnlyNoPileup);
544 //V0A Centrality (if PbPb / pPb)
545 if(! fHistMultiplicityV0ABeforeTrigSel) {
546 fHistMultiplicityV0ABeforeTrigSel = new TH1F("fHistMultiplicityV0ABeforeTrigSel",
547 "Centrality Distribution: V0A;V0A Centrality;Events",
549 fListHist->Add(fHistMultiplicityV0ABeforeTrigSel);
551 if(! fHistMultiplicityV0AForTrigEvt) {
552 fHistMultiplicityV0AForTrigEvt = new TH1F("fHistMultiplicityV0AForTrigEvt",
553 "Centrality Distribution: V0A;V0A Centrality;Events",
555 fListHist->Add(fHistMultiplicityV0AForTrigEvt);
557 if(! fHistMultiplicityV0A) {
558 fHistMultiplicityV0A = new TH1F("fHistMultiplicityV0A",
559 "Centrality Distribution: V0A;V0A Centrality;Events",
561 fListHist->Add(fHistMultiplicityV0A);
563 if(! fHistMultiplicityV0ANoTPCOnly) {
564 fHistMultiplicityV0ANoTPCOnly = new TH1F("fHistMultiplicityV0ANoTPCOnly",
565 "Centrality Distribution: V0A;V0A Centrality;Events",
567 fListHist->Add(fHistMultiplicityV0ANoTPCOnly);
569 if(! fHistMultiplicityV0ANoTPCOnlyNoPileup) {
570 fHistMultiplicityV0ANoTPCOnlyNoPileup = new TH1F("fHistMultiplicityV0ANoTPCOnlyNoPileup",
571 "Centrality Distribution: V0A;V0A Centrality;Events",
573 fListHist->Add(fHistMultiplicityV0ANoTPCOnlyNoPileup);
576 //ZNA Centrality (if PbPb / pPb)
577 if(! fHistMultiplicityZNABeforeTrigSel) {
578 fHistMultiplicityZNABeforeTrigSel = new TH1F("fHistMultiplicityZNABeforeTrigSel",
579 "Centrality Distribution: ZNA;ZNA Centrality;Events",
581 fListHist->Add(fHistMultiplicityZNABeforeTrigSel);
583 if(! fHistMultiplicityZNAForTrigEvt) {
584 fHistMultiplicityZNAForTrigEvt = new TH1F("fHistMultiplicityZNAForTrigEvt",
585 "Centrality Distribution: ZNA;ZNA Centrality;Events",
587 fListHist->Add(fHistMultiplicityZNAForTrigEvt);
589 if(! fHistMultiplicityZNA) {
590 fHistMultiplicityZNA = new TH1F("fHistMultiplicityZNA",
591 "Centrality Distribution: ZNA;ZNA Centrality;Events",
593 fListHist->Add(fHistMultiplicityZNA);
595 if(! fHistMultiplicityZNANoTPCOnly) {
596 fHistMultiplicityZNANoTPCOnly = new TH1F("fHistMultiplicityZNANoTPCOnly",
597 "Centrality Distribution: ZNA;ZNA Centrality;Events",
599 fListHist->Add(fHistMultiplicityZNANoTPCOnly);
601 if(! fHistMultiplicityZNANoTPCOnlyNoPileup) {
602 fHistMultiplicityZNANoTPCOnlyNoPileup = new TH1F("fHistMultiplicityZNANoTPCOnlyNoPileup",
603 "Centrality Distribution: ZNA;ZNA Centrality;Events",
605 fListHist->Add(fHistMultiplicityZNANoTPCOnlyNoPileup);
608 //TRK Centrality (if PbPb / pPb)
609 if(! fHistMultiplicityTRKBeforeTrigSel) {
610 fHistMultiplicityTRKBeforeTrigSel = new TH1F("fHistMultiplicityTRKBeforeTrigSel",
611 "Centrality Distribution: TRK;TRK Centrality;Events",
613 fListHist->Add(fHistMultiplicityTRKBeforeTrigSel);
615 if(! fHistMultiplicityTRKForTrigEvt) {
616 fHistMultiplicityTRKForTrigEvt = new TH1F("fHistMultiplicityTRKForTrigEvt",
617 "Centrality Distribution: TRK;TRK Centrality;Events",
619 fListHist->Add(fHistMultiplicityTRKForTrigEvt);
621 if(! fHistMultiplicityTRK) {
622 fHistMultiplicityTRK = new TH1F("fHistMultiplicityTRK",
623 "Centrality Distribution: TRK;TRK Centrality;Events",
625 fListHist->Add(fHistMultiplicityTRK);
627 if(! fHistMultiplicityTRKNoTPCOnly) {
628 fHistMultiplicityTRKNoTPCOnly = new TH1F("fHistMultiplicityTRKNoTPCOnly",
629 "Centrality Distribution: TRK;TRK Centrality;Events",
631 fListHist->Add(fHistMultiplicityTRKNoTPCOnly);
633 if(! fHistMultiplicityTRKNoTPCOnlyNoPileup) {
634 fHistMultiplicityTRKNoTPCOnlyNoPileup = new TH1F("fHistMultiplicityTRKNoTPCOnlyNoPileup",
635 "Centrality Distribution: TRK;TRK Centrality;Events",
637 fListHist->Add(fHistMultiplicityTRKNoTPCOnlyNoPileup);
640 //SPD Centrality (if PbPb / pPb)
641 if(! fHistMultiplicitySPDBeforeTrigSel) {
642 fHistMultiplicitySPDBeforeTrigSel = new TH1F("fHistMultiplicitySPDBeforeTrigSel",
643 "Centrality Distribution: SPD;SPD Centrality;Events",
645 fListHist->Add(fHistMultiplicitySPDBeforeTrigSel);
647 if(! fHistMultiplicitySPDForTrigEvt) {
648 fHistMultiplicitySPDForTrigEvt = new TH1F("fHistMultiplicitySPDForTrigEvt",
649 "Centrality Distribution: SPD;SPD Centrality;Events",
651 fListHist->Add(fHistMultiplicitySPDForTrigEvt);
653 if(! fHistMultiplicitySPD) {
654 fHistMultiplicitySPD = new TH1F("fHistMultiplicitySPD",
655 "Centrality Distribution: SPD;SPD Centrality;Events",
657 fListHist->Add(fHistMultiplicitySPD);
659 if(! fHistMultiplicitySPDNoTPCOnly) {
660 fHistMultiplicitySPDNoTPCOnly = new TH1F("fHistMultiplicitySPDNoTPCOnly",
661 "Centrality Distribution: SPD;SPD Centrality;Events",
663 fListHist->Add(fHistMultiplicitySPDNoTPCOnly);
665 if(! fHistMultiplicitySPDNoTPCOnlyNoPileup) {
666 fHistMultiplicitySPDNoTPCOnlyNoPileup = new TH1F("fHistMultiplicitySPDNoTPCOnlyNoPileup",
667 "Centrality Distribution: SPD;SPD Centrality;Events",
669 fListHist->Add(fHistMultiplicitySPDNoTPCOnlyNoPileup);
673 //----------------------------------
674 // Primary Vertex Position Histos
675 //----------------------------------
678 fHistPVx = new TH1F("fHistPVx",
679 "PV x position;Nbr of Evts;x",
681 fListHist->Add(fHistPVx);
684 fHistPVy = new TH1F("fHistPVy",
685 "PV y position;Nbr of Evts;y",
687 fListHist->Add(fHistPVy);
690 fHistPVz = new TH1F("fHistPVz",
691 "PV z position;Nbr of Evts;z",
693 fListHist->Add(fHistPVz);
696 if(! fHistPVxAnalysis) {
697 fHistPVxAnalysis = new TH1F("fHistPVxAnalysis",
698 "PV x position;Nbr of Evts;x",
700 fListHist->Add(fHistPVxAnalysis);
702 if(! fHistPVyAnalysis) {
703 fHistPVyAnalysis = new TH1F("fHistPVyAnalysis",
704 "PV y position;Nbr of Evts;y",
706 fListHist->Add(fHistPVyAnalysis);
708 if(! fHistPVzAnalysis) {
709 fHistPVzAnalysis = new TH1F("fHistPVzAnalysis",
710 "PV z position;Nbr of Evts;z",
712 fListHist->Add(fHistPVzAnalysis);
715 //List of Histograms: Normal
716 PostData(1, fListHist);
718 //TTree Object: Saved to base directory. Should cache to disk while saving.
719 //(Important to avoid excessive memory usage, particularly when merging)
720 PostData(2, fTreeCascade);
722 }// end UserCreateOutputObjects
725 //________________________________________________________________________
726 void AliAnalysisTaskExtractCascade::UserExec(Option_t *)
729 // Called for each event
731 AliESDEvent *lESDevent = 0x0;
733 Int_t lNumberOfV0s = -1;
734 Double_t lTrkgPrimaryVtxPos[3] = {-100.0, -100.0, -100.0};
735 Double_t lBestPrimaryVtxPos[3] = {-100.0, -100.0, -100.0};
736 Double_t lMagneticField = -10.;
738 // Connect to the InputEvent
739 // After these lines, we should have an ESD/AOD event + the number of V0s in it.
741 // Appropriate for ESD analysis!
743 lESDevent = dynamic_cast<AliESDEvent*>( InputEvent() );
745 AliWarning("ERROR: lESDevent not available \n");
749 /* --- Acquisition of exact event ID
750 fTreeVariableRunNumber = lESDevent->GetRunNumber();
751 fTreeVariableEventNumber =
752 ( ( ((ULong64_t)lESDevent->GetPeriodNumber() ) << 36 ) |
753 ( ((ULong64_t)lESDevent->GetOrbitNumber () ) << 12 ) |
754 ((ULong64_t)lESDevent->GetBunchCrossNumber() ) );
757 //------------------------------------------------
758 // Multiplicity Information Acquistion
759 //------------------------------------------------
761 //REVISED multiplicity estimator after 'multiplicity day' (2011)
762 Int_t lMultiplicity = -100;
763 Int_t lMultiplicityV0A = -100;
764 Int_t lMultiplicityZNA = -100;
765 Int_t lMultiplicityTRK = -100;
766 Int_t lMultiplicitySPD = -100;
769 if(fkIsNuclear == kFALSE) lMultiplicity = fESDtrackCuts->GetReferenceMultiplicity(lESDevent, AliESDtrackCuts::kTrackletsITSTPC, fEtaRefMult );
771 //---> If this is a nuclear collision, then go nuclear on "multiplicity" variable...
772 //---> Warning: Experimental
773 if(fkIsNuclear == kTRUE){
774 AliCentrality* centrality;
775 centrality = lESDevent->GetCentrality();
776 lMultiplicity = ( ( Int_t ) ( centrality->GetCentralityPercentile( fCentralityEstimator.Data() ) ) );
777 lMultiplicityV0A = ( ( Int_t ) ( centrality->GetCentralityPercentile( "V0A" ) ) );
778 lMultiplicityZNA = ( ( Int_t ) ( centrality->GetCentralityPercentile( "ZNA" ) ) );
779 lMultiplicityTRK = ( ( Int_t ) ( centrality->GetCentralityPercentile( "TRK" ) ) );
780 lMultiplicitySPD = ( ( Int_t ) ( centrality->GetCentralityPercentile( "CL1" ) ) );
781 if (centrality->GetQuality()>1) {
782 PostData(1, fListHist);
783 PostData(2, fTreeCascade);
788 //Set variable for filling tree afterwards!
789 //---> pp case......: GetReferenceMultiplicity
790 //---> Pb-Pb case...: Centrality by V0M
792 fTreeCascVarMultiplicity = lMultiplicity;
793 fTreeCascVarMultiplicityV0A = lMultiplicityV0A;
794 fTreeCascVarMultiplicityZNA = lMultiplicityZNA;
795 fTreeCascVarMultiplicityTRK = lMultiplicityTRK;
796 fTreeCascVarMultiplicitySPD = lMultiplicitySPD;
798 fHistV0MultiplicityBeforeTrigSel->Fill ( lESDevent->GetNumberOfV0s() );
799 fHistMultiplicityBeforeTrigSel->Fill ( lMultiplicity );
800 fHistMultiplicityV0ABeforeTrigSel->Fill ( lMultiplicityV0A );
801 fHistMultiplicityZNABeforeTrigSel->Fill ( lMultiplicityZNA );
802 fHistMultiplicityTRKBeforeTrigSel->Fill ( lMultiplicityTRK );
803 fHistMultiplicitySPDBeforeTrigSel->Fill ( lMultiplicitySPD );
805 //------------------------------------------------
807 //------------------------------------------------
809 UInt_t maskIsSelected = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
810 Bool_t isSelected = 0;
811 isSelected = (maskIsSelected & AliVEvent::kMB) == AliVEvent::kMB;
813 //pA triggering: CINT7
814 if( fkSwitchINT7 ) isSelected = (maskIsSelected & AliVEvent::kINT7) == AliVEvent::kINT7;
816 //Standard Min-Bias Selection
817 if ( ! isSelected ) {
818 PostData(1, fListHist);
819 PostData(2, fTreeCascade);
823 //------------------------------------------------
824 // Rerun cascade vertexer!
825 //------------------------------------------------
827 if( fkRunVertexers ){
828 lESDevent->ResetCascades();
829 lESDevent->ResetV0s();
831 AliV0vertexer lV0vtxer;
832 AliCascadeVertexer lCascVtxer;
834 lV0vtxer.SetDefaultCuts(fV0VertexerSels);
835 lCascVtxer.SetDefaultCuts(fCascadeVertexerSels);
837 lV0vtxer.Tracks2V0vertices(lESDevent);
838 lCascVtxer.V0sTracks2CascadeVertices(lESDevent);
840 //------------------------------------------------
841 // After Trigger Selection
842 //------------------------------------------------
844 lNumberOfV0s = lESDevent->GetNumberOfV0s();
846 //Set variable for filling tree afterwards!
847 fHistV0MultiplicityForTrigEvt->Fill(lNumberOfV0s);
848 fHistMultiplicityForTrigEvt->Fill ( lMultiplicity );
849 fHistMultiplicityV0AForTrigEvt ->Fill( lMultiplicityV0A );
850 fHistMultiplicityZNAForTrigEvt ->Fill( lMultiplicityZNA );
851 fHistMultiplicityTRKForTrigEvt ->Fill( lMultiplicityTRK );
852 fHistMultiplicitySPDForTrigEvt ->Fill( lMultiplicitySPD );
854 //------------------------------------------------
855 // Getting: Primary Vertex + MagField Info
856 //------------------------------------------------
858 const AliESDVertex *lPrimaryTrackingESDVtx = lESDevent->GetPrimaryVertexTracks();
859 // get the vtx stored in ESD found with tracks
860 lPrimaryTrackingESDVtx->GetXYZ( lTrkgPrimaryVtxPos );
862 const AliESDVertex *lPrimaryBestESDVtx = lESDevent->GetPrimaryVertex();
863 // get the best primary vertex available for the event
864 // As done in AliCascadeVertexer, we keep the one which is the best one available.
865 // between : Tracking vertex > SPD vertex > TPC vertex > default SPD vertex
866 // This one will be used for next calculations (DCA essentially)
867 lPrimaryBestESDVtx->GetXYZ( lBestPrimaryVtxPos );
869 Double_t lPrimaryVtxPosition[3];
870 const AliVVertex *primaryVtx = lESDevent->GetPrimaryVertex();
871 lPrimaryVtxPosition[0] = primaryVtx->GetX();
872 lPrimaryVtxPosition[1] = primaryVtx->GetY();
873 lPrimaryVtxPosition[2] = primaryVtx->GetZ();
874 fHistPVx->Fill( lPrimaryVtxPosition[0] );
875 fHistPVy->Fill( lPrimaryVtxPosition[1] );
876 fHistPVz->Fill( lPrimaryVtxPosition[2] );
878 //------------------------------------------------
879 // Primary Vertex Requirements Section:
880 // ---> pp and PbPb: Only requires |z|<10cm
881 // ---> pPb: all requirements checked at this stage
882 //------------------------------------------------
884 //Roberto's PV selection criteria, implemented 17th April 2013
886 /* vertex selection */
887 Bool_t fHasVertex = kFALSE;
889 const AliESDVertex *vertex = lESDevent->GetPrimaryVertexTracks();
890 if (vertex->GetNContributors() < 1) {
891 vertex = lESDevent->GetPrimaryVertexSPD();
892 if (vertex->GetNContributors() < 1) fHasVertex = kFALSE;
893 else fHasVertex = kTRUE;
894 TString vtxTyp = vertex->GetTitle();
896 vertex->GetCovarianceMatrix(cov);
897 Double_t zRes = TMath::Sqrt(cov[5]);
898 if (vtxTyp.Contains("vertexer:Z") && (zRes>0.25)) fHasVertex = kFALSE;
900 else fHasVertex = kTRUE;
902 //Is First event in chunk rejection: Still present!
903 if(fkpAVertexSelection==kTRUE && fHasVertex == kFALSE) {
904 AliWarning("Pb / | PV does not satisfy selection criteria!");
905 PostData(1, fListHist);
906 PostData(2, fTreeCascade);
910 //Is First event in chunk rejection: Still present!
911 if(fkpAVertexSelection==kTRUE && fUtils->IsFirstEventInChunk(lESDevent)) {
912 AliWarning("Pb / | This is the first event in the chunk!");
913 PostData(1, fListHist);
914 PostData(2, fTreeCascade);
918 //17 April Fix: Always do primary vertex Z selection, after pA vertex selection from Roberto
919 if(TMath::Abs(lBestPrimaryVtxPos[2]) > 10.0) {
920 AliWarning("Pb / | Z position of Best Prim Vtx | > 10.0 cm ... return !");
921 PostData(1, fListHist);
922 PostData(2, fTreeCascade);
927 lMagneticField = lESDevent->GetMagneticField( );
928 fHistV0MultiplicityForSelEvt ->Fill( lNumberOfV0s );
929 fHistMultiplicity->Fill(lMultiplicity);
930 fHistMultiplicityV0A->Fill(lMultiplicityV0A);
931 fHistMultiplicityZNA->Fill(lMultiplicityZNA);
932 fHistMultiplicityTRK->Fill(lMultiplicityTRK);
933 fHistMultiplicitySPD->Fill(lMultiplicitySPD);
935 //------------------------------------------------
936 // SKIP: Events with well-established PVtx
937 //------------------------------------------------
939 const AliESDVertex *lPrimaryTrackingESDVtxCheck = lESDevent->GetPrimaryVertexTracks();
940 const AliESDVertex *lPrimarySPDVtx = lESDevent->GetPrimaryVertexSPD();
941 if (!lPrimarySPDVtx->GetStatus() && !lPrimaryTrackingESDVtxCheck->GetStatus() && fkpAVertexSelection==kFALSE ){
942 AliWarning("Pb / No SPD prim. vertex nor prim. Tracking vertex ... return !");
943 PostData(1, fListHist);
944 PostData(2, fTreeCascade);
947 fHistV0MultiplicityForSelEvtNoTPCOnly ->Fill( lNumberOfV0s );
948 fHistMultiplicityNoTPCOnly->Fill(lMultiplicity);
949 fHistMultiplicityV0ANoTPCOnly->Fill(lMultiplicityV0A);
950 fHistMultiplicityZNANoTPCOnly->Fill(lMultiplicityZNA);
951 fHistMultiplicityTRKNoTPCOnly->Fill(lMultiplicityTRK);
952 fHistMultiplicitySPDNoTPCOnly->Fill(lMultiplicitySPD);
954 //------------------------------------------------
955 // Pileup Rejection Studies
956 //------------------------------------------------
958 // FIXME : quality selection regarding pile-up rejection
959 if(lESDevent->IsPileupFromSPD() && !fkIsNuclear ){// minContributors=3, minZdist=0.8, nSigmaZdist=3., nSigmaDiamXY=2., nSigmaDiamZ=5. -> see http://alisoft.cern.ch/viewvc/trunk/STEER/AliESDEvent.h?root=AliRoot&r1=41914&r2=42199&pathrev=42199
960 AliWarning("Pb / This is tagged as Pileup from SPD... return !");
961 PostData(1, fListHist);
962 PostData(2, fTreeCascade);
965 fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup ->Fill( lNumberOfV0s );
966 fHistMultiplicityNoTPCOnlyNoPileup->Fill(lMultiplicity);
967 fHistMultiplicityV0ANoTPCOnlyNoPileup->Fill(lMultiplicityV0A);
968 fHistMultiplicityZNANoTPCOnlyNoPileup->Fill(lMultiplicityZNA);
969 fHistMultiplicityTRKNoTPCOnlyNoPileup->Fill(lMultiplicityTRK);
970 fHistMultiplicitySPDNoTPCOnlyNoPileup->Fill(lMultiplicitySPD);
972 //Do control histograms without the IsFromVertexerZ events, but consider them in analysis...
973 if( ! (lESDevent->GetPrimaryVertex()->IsFromVertexerZ() ) ){
974 fHistPVxAnalysis->Fill( lPrimaryVtxPosition[0] );
975 fHistPVyAnalysis->Fill( lPrimaryVtxPosition[1] );
976 fHistPVzAnalysis->Fill( lPrimaryVtxPosition[2] );
979 //------------------------------------------------
980 // MAIN CASCADE LOOP STARTS HERE
981 //------------------------------------------------
982 // Code Credit: Antonin Maire (thanks^100)
983 // ---> This is an adaptation
985 Long_t ncascades = 0;
986 ncascades = lESDevent->GetNumberOfCascades();
988 for (Int_t iXi = 0; iXi < ncascades; iXi++){
989 //------------------------------------------------
991 //------------------------------------------------
992 //Double_t lTrkgPrimaryVtxRadius3D = -500.0;
993 //Double_t lBestPrimaryVtxRadius3D = -500.0;
995 // - 1st part of initialisation : variables needed to store AliESDCascade data members
996 Double_t lEffMassXi = 0. ;
997 //Double_t lChi2Xi = -1. ;
998 Double_t lDcaXiDaughters = -1. ;
999 Double_t lXiCosineOfPointingAngle = -1. ;
1000 Double_t lPosXi[3] = { -1000.0, -1000.0, -1000.0 };
1001 Double_t lXiRadius = -1000. ;
1003 // - 2nd part of initialisation : Nbr of clusters within TPC for the 3 daughter cascade tracks
1004 Int_t lPosTPCClusters = -1; // For ESD only ...//FIXME : wait for availability in AOD
1005 Int_t lNegTPCClusters = -1; // For ESD only ...
1006 Int_t lBachTPCClusters = -1; // For ESD only ...
1008 // - 3rd part of initialisation : about V0 part in cascades
1009 Double_t lInvMassLambdaAsCascDghter = 0.;
1010 //Double_t lV0Chi2Xi = -1. ;
1011 Double_t lDcaV0DaughtersXi = -1.;
1013 Double_t lDcaBachToPrimVertexXi = -1., lDcaV0ToPrimVertexXi = -1.;
1014 Double_t lDcaPosToPrimVertexXi = -1.;
1015 Double_t lDcaNegToPrimVertexXi = -1.;
1016 Double_t lV0CosineOfPointingAngleXi = -1. ;
1017 Double_t lV0CosineOfPointingAngleXiSpecial = -1. ;
1018 Double_t lPosV0Xi[3] = { -1000. , -1000., -1000. }; // Position of VO coming from cascade
1019 Double_t lV0RadiusXi = -1000.0;
1020 Double_t lV0quality = 0.;
1022 // - 4th part of initialisation : Effective masses
1023 Double_t lInvMassXiMinus = 0.;
1024 Double_t lInvMassXiPlus = 0.;
1025 Double_t lInvMassOmegaMinus = 0.;
1026 Double_t lInvMassOmegaPlus = 0.;
1028 // - 6th part of initialisation : extra info for QA
1029 Double_t lXiMomX = 0. , lXiMomY = 0., lXiMomZ = 0.;
1030 Double_t lXiTransvMom = 0. ;
1031 Double_t lXiTransvMomMC= 0. ;
1032 Double_t lXiTotMom = 0. ;
1034 Double_t lBachMomX = 0., lBachMomY = 0., lBachMomZ = 0.;
1035 //Double_t lBachTransvMom = 0.;
1036 //Double_t lBachTotMom = 0.;
1038 fTreeCascVarNegNSigmaPion = -100;
1039 fTreeCascVarNegNSigmaProton = -100;
1040 fTreeCascVarPosNSigmaPion = -100;
1041 fTreeCascVarPosNSigmaProton = -100;
1042 fTreeCascVarBachNSigmaPion = -100;
1043 fTreeCascVarBachNSigmaKaon = -100;
1045 Short_t lChargeXi = -2;
1046 //Double_t lV0toXiCosineOfPointingAngle = 0. ;
1048 Double_t lRapXi = -20.0, lRapOmega = -20.0, lRapMC = -20.0; // lEta = -20.0, lTheta = 360., lPhi = 720. ;
1049 //Double_t lAlphaXi = -200., lPtArmXi = -200.0;
1051 // -------------------------------------
1052 // II.ESD - Calculation Part dedicated to Xi vertices (ESD)
1054 AliESDcascade *xi = lESDevent->GetCascade(iXi);
1058 // - II.Step 1 : around primary vertex
1060 //lTrkgPrimaryVtxRadius3D = TMath::Sqrt( lTrkgPrimaryVtxPos[0] * lTrkgPrimaryVtxPos[0] +
1061 // lTrkgPrimaryVtxPos[1] * lTrkgPrimaryVtxPos[1] +
1062 // lTrkgPrimaryVtxPos[2] * lTrkgPrimaryVtxPos[2] );
1064 //lBestPrimaryVtxRadius3D = TMath::Sqrt( lBestPrimaryVtxPos[0] * lBestPrimaryVtxPos[0] +
1065 // lBestPrimaryVtxPos[1] * lBestPrimaryVtxPos[1] +
1066 // lBestPrimaryVtxPos[2] * lBestPrimaryVtxPos[2] );
1068 // - II.Step 2 : Assigning the necessary variables for specific AliESDcascade data members (ESD)
1071 xi->ChangeMassHypothesis(lV0quality , 3312); // default working hypothesis : cascade = Xi- decay
1073 lEffMassXi = xi->GetEffMassXi();
1074 //lChi2Xi = xi->GetChi2Xi();
1075 lDcaXiDaughters = xi->GetDcaXiDaughters();
1076 lXiCosineOfPointingAngle = xi->GetCascadeCosineOfPointingAngle( lBestPrimaryVtxPos[0],
1077 lBestPrimaryVtxPos[1],
1078 lBestPrimaryVtxPos[2] );
1079 // Take care : the best available vertex should be used (like in AliCascadeVertexer)
1081 xi->GetXYZcascade( lPosXi[0], lPosXi[1], lPosXi[2] );
1082 lXiRadius = TMath::Sqrt( lPosXi[0]*lPosXi[0] + lPosXi[1]*lPosXi[1] );
1084 // - II.Step 3 : around the tracks : Bach + V0 (ESD)
1085 // ~ Necessary variables for ESDcascade data members coming from the ESDv0 part (inheritance)
1088 UInt_t lIdxPosXi = (UInt_t) TMath::Abs( xi->GetPindex() );
1089 UInt_t lIdxNegXi = (UInt_t) TMath::Abs( xi->GetNindex() );
1090 UInt_t lBachIdx = (UInt_t) TMath::Abs( xi->GetBindex() );
1091 // Care track label can be negative in MC production (linked with the track quality)
1092 // However = normally, not the case for track index ...
1094 // FIXME : rejection of a double use of a daughter track (nothing but just a crosscheck of what is done in the cascade vertexer)
1095 if(lBachIdx == lIdxNegXi) {
1096 AliWarning("Pb / Idx(Bach. track) = Idx(Neg. track) ... continue!"); continue;
1098 if(lBachIdx == lIdxPosXi) {
1099 AliWarning("Pb / Idx(Bach. track) = Idx(Pos. track) ... continue!"); continue;
1102 AliESDtrack *pTrackXi = lESDevent->GetTrack( lIdxPosXi );
1103 AliESDtrack *nTrackXi = lESDevent->GetTrack( lIdxNegXi );
1104 AliESDtrack *bachTrackXi = lESDevent->GetTrack( lBachIdx );
1106 if (!pTrackXi || !nTrackXi || !bachTrackXi ) {
1107 AliWarning("ERROR: Could not retrieve one of the 3 ESD daughter tracks of the cascade ...");
1111 fTreeCascVarPosEta = pTrackXi->Eta();
1112 fTreeCascVarNegEta = nTrackXi->Eta();
1113 fTreeCascVarBachEta = bachTrackXi->Eta();
1115 Double_t lBMom[3], lNMom[3], lPMom[3];
1116 xi->GetBPxPyPz( lBMom[0], lBMom[1], lBMom[2] );
1117 xi->GetPPxPyPz( lPMom[0], lPMom[1], lPMom[2] );
1118 xi->GetNPxPyPz( lNMom[0], lNMom[1], lNMom[2] );
1120 fTreeCascVarBachTransMom = TMath::Sqrt( lBMom[0]*lBMom[0] + lBMom[1]*lBMom[1] );
1121 fTreeCascVarPosTransMom = TMath::Sqrt( lPMom[0]*lPMom[0] + lPMom[1]*lPMom[1] );
1122 fTreeCascVarNegTransMom = TMath::Sqrt( lNMom[0]*lNMom[0] + lNMom[1]*lNMom[1] );
1124 //------------------------------------------------
1125 // TPC dEdx information
1126 //------------------------------------------------
1127 fTreeCascVarNegNSigmaPion = fPIDResponse->NumberOfSigmasTPC( nTrackXi, AliPID::kPion );
1128 fTreeCascVarNegNSigmaProton = fPIDResponse->NumberOfSigmasTPC( nTrackXi, AliPID::kProton );
1129 fTreeCascVarPosNSigmaPion = fPIDResponse->NumberOfSigmasTPC( pTrackXi, AliPID::kPion );
1130 fTreeCascVarPosNSigmaProton = fPIDResponse->NumberOfSigmasTPC( pTrackXi, AliPID::kProton );
1131 fTreeCascVarBachNSigmaPion = fPIDResponse->NumberOfSigmasTPC( bachTrackXi, AliPID::kPion );
1132 fTreeCascVarBachNSigmaKaon = fPIDResponse->NumberOfSigmasTPC( bachTrackXi, AliPID::kKaon );
1134 //------------------------------------------------
1135 // TPC Number of clusters info
1136 // --- modified to save the smallest number
1137 // --- of TPC clusters for the 3 tracks
1138 //------------------------------------------------
1140 lPosTPCClusters = pTrackXi->GetTPCNcls();
1141 lNegTPCClusters = nTrackXi->GetTPCNcls();
1142 lBachTPCClusters = bachTrackXi->GetTPCNcls();
1144 // 1 - Poor quality related to TPCrefit
1145 ULong_t pStatus = pTrackXi->GetStatus();
1146 ULong_t nStatus = nTrackXi->GetStatus();
1147 ULong_t bachStatus = bachTrackXi->GetStatus();
1149 fTreeCascVarkITSRefitBachelor = kTRUE;
1150 fTreeCascVarkITSRefitNegative = kTRUE;
1151 fTreeCascVarkITSRefitPositive = kTRUE;
1153 if ((pStatus&AliESDtrack::kTPCrefit) == 0) { AliWarning("Pb / V0 Pos. track has no TPCrefit ... continue!"); continue; }
1154 if ((nStatus&AliESDtrack::kTPCrefit) == 0) { AliWarning("Pb / V0 Neg. track has no TPCrefit ... continue!"); continue; }
1155 if ((bachStatus&AliESDtrack::kTPCrefit) == 0) { AliWarning("Pb / Bach. track has no TPCrefit ... continue!"); continue; }
1157 //Extra Debug Information: booleans for ITS refit
1158 if ((pStatus&AliESDtrack::kITSrefit) == 0) { fTreeCascVarkITSRefitPositive = kFALSE; }
1159 if ((nStatus&AliESDtrack::kITSrefit) == 0) { fTreeCascVarkITSRefitNegative = kFALSE; }
1160 if ((bachStatus&AliESDtrack::kITSrefit) == 0) { fTreeCascVarkITSRefitBachelor = kFALSE; }
1162 // 2 - Poor quality related to TPC clusters: lowest cut of 70 clusters
1163 if(lPosTPCClusters < 70) { AliWarning("Pb / V0 Pos. track has less than 70 TPC clusters ... continue!"); continue; }
1164 if(lNegTPCClusters < 70) { AliWarning("Pb / V0 Neg. track has less than 70 TPC clusters ... continue!"); continue; }
1165 if(lBachTPCClusters < 70) { AliWarning("Pb / Bach. track has less than 70 TPC clusters ... continue!"); continue; }
1166 Int_t leastnumberofclusters = 1000;
1167 if( lPosTPCClusters < leastnumberofclusters ) leastnumberofclusters = lPosTPCClusters;
1168 if( lNegTPCClusters < leastnumberofclusters ) leastnumberofclusters = lNegTPCClusters;
1169 if( lBachTPCClusters < leastnumberofclusters ) leastnumberofclusters = lBachTPCClusters;
1171 lInvMassLambdaAsCascDghter = xi->GetEffMass();
1172 // This value shouldn't change, whatever the working hyp. is : Xi-, Xi+, Omega-, Omega+
1173 lDcaV0DaughtersXi = xi->GetDcaV0Daughters();
1174 //lV0Chi2Xi = xi->GetChi2V0();
1176 lV0CosineOfPointingAngleXi = xi->GetV0CosineOfPointingAngle( lBestPrimaryVtxPos[0],
1177 lBestPrimaryVtxPos[1],
1178 lBestPrimaryVtxPos[2] );
1179 //Modification: V0 CosPA wrt to Cascade decay vertex
1180 lV0CosineOfPointingAngleXiSpecial = xi->GetV0CosineOfPointingAngle( lPosXi[0],
1184 lDcaV0ToPrimVertexXi = xi->GetD( lBestPrimaryVtxPos[0],
1185 lBestPrimaryVtxPos[1],
1186 lBestPrimaryVtxPos[2] );
1188 lDcaBachToPrimVertexXi = TMath::Abs( bachTrackXi->GetD( lBestPrimaryVtxPos[0],
1189 lBestPrimaryVtxPos[1],
1191 // Note : AliExternalTrackParam::GetD returns an algebraic value ...
1193 xi->GetXYZ( lPosV0Xi[0], lPosV0Xi[1], lPosV0Xi[2] );
1194 lV0RadiusXi = TMath::Sqrt( lPosV0Xi[0]*lPosV0Xi[0] + lPosV0Xi[1]*lPosV0Xi[1] );
1196 lDcaPosToPrimVertexXi = TMath::Abs( pTrackXi ->GetD( lBestPrimaryVtxPos[0],
1197 lBestPrimaryVtxPos[1],
1200 lDcaNegToPrimVertexXi = TMath::Abs( nTrackXi ->GetD( lBestPrimaryVtxPos[0],
1201 lBestPrimaryVtxPos[1],
1204 // - II.Step 4 : around effective masses (ESD)
1205 // ~ change mass hypotheses to cover all the possibilities : Xi-/+, Omega -/+
1207 if( bachTrackXi->Charge() < 0 ) {
1209 xi->ChangeMassHypothesis(lV0quality , 3312);
1210 // Calculate the effective mass of the Xi- candidate.
1211 // pdg code 3312 = Xi-
1212 lInvMassXiMinus = xi->GetEffMassXi();
1215 xi->ChangeMassHypothesis(lV0quality , 3334);
1216 // Calculate the effective mass of the Xi- candidate.
1217 // pdg code 3334 = Omega-
1218 lInvMassOmegaMinus = xi->GetEffMassXi();
1221 xi->ChangeMassHypothesis(lV0quality , 3312); // Back to default hyp.
1222 }// end if negative bachelor
1225 if( bachTrackXi->Charge() > 0 ){
1227 xi->ChangeMassHypothesis(lV0quality , -3312);
1228 // Calculate the effective mass of the Xi+ candidate.
1229 // pdg code -3312 = Xi+
1230 lInvMassXiPlus = xi->GetEffMassXi();
1233 xi->ChangeMassHypothesis(lV0quality , -3334);
1234 // Calculate the effective mass of the Xi+ candidate.
1235 // pdg code -3334 = Omega+
1236 lInvMassOmegaPlus = xi->GetEffMassXi();
1239 xi->ChangeMassHypothesis(lV0quality , -3312); // Back to "default" hyp.
1240 }// end if positive bachelor
1241 // - II.Step 6 : extra info for QA (ESD)
1242 // miscellaneous pieces of info that may help regarding data quality assessment.
1245 xi->GetPxPyPz( lXiMomX, lXiMomY, lXiMomZ );
1246 lXiTransvMom = TMath::Sqrt( lXiMomX*lXiMomX + lXiMomY*lXiMomY );
1247 lXiTotMom = TMath::Sqrt( lXiMomX*lXiMomX + lXiMomY*lXiMomY + lXiMomZ*lXiMomZ );
1249 xi->GetBPxPyPz( lBachMomX, lBachMomY, lBachMomZ );
1250 //lBachTransvMom = TMath::Sqrt( lBachMomX*lBachMomX + lBachMomY*lBachMomY );
1251 //lBachTotMom = TMath::Sqrt( lBachMomX*lBachMomX + lBachMomY*lBachMomY + lBachMomZ*lBachMomZ );
1253 lChargeXi = xi->Charge();
1255 //lV0toXiCosineOfPointingAngle = xi->GetV0CosineOfPointingAngle( lPosXi[0], lPosXi[1], lPosXi[2] );
1257 lRapXi = xi->RapXi();
1258 lRapOmega = xi->RapOmega();
1260 //lTheta = xi->Theta() *180.0/TMath::Pi();
1261 //lPhi = xi->Phi() *180.0/TMath::Pi();
1262 //lAlphaXi = xi->AlphaXi();
1263 //lPtArmXi = xi->PtArmXi();
1265 //------------------------------------------------
1266 // Set Variables for adding to tree
1267 //------------------------------------------------
1269 /* 1*/ fTreeCascVarCharge = lChargeXi;
1270 /* 2*/ if(lInvMassXiMinus!=0) fTreeCascVarMassAsXi = lInvMassXiMinus;
1271 /* 2*/ if(lInvMassXiPlus!=0) fTreeCascVarMassAsXi = lInvMassXiPlus;
1272 /* 3*/ if(lInvMassOmegaMinus!=0) fTreeCascVarMassAsOmega = lInvMassOmegaMinus;
1273 /* 3*/ if(lInvMassOmegaPlus!=0) fTreeCascVarMassAsOmega = lInvMassOmegaPlus;
1274 /* 4*/ fTreeCascVarPt = lXiTransvMom;
1275 /* 4*/ fTreeCascVarPtMC = lXiTransvMomMC;
1276 /* 5*/ fTreeCascVarRapXi = lRapXi ;
1277 /* 5*/ fTreeCascVarRapMC = lRapMC ;
1278 /* 6*/ fTreeCascVarRapOmega = lRapOmega ;
1279 /* 7*/ fTreeCascVarDCACascDaughters = lDcaXiDaughters;
1280 /* 8*/ fTreeCascVarDCABachToPrimVtx = lDcaBachToPrimVertexXi;
1281 /* 9*/ fTreeCascVarDCAV0Daughters = lDcaV0DaughtersXi;
1282 /*10*/ fTreeCascVarDCAV0ToPrimVtx = lDcaV0ToPrimVertexXi;
1283 /*11*/ fTreeCascVarDCAPosToPrimVtx = lDcaPosToPrimVertexXi;
1284 /*12*/ fTreeCascVarDCANegToPrimVtx = lDcaNegToPrimVertexXi;
1285 /*13*/ fTreeCascVarCascCosPointingAngle = lXiCosineOfPointingAngle;
1286 /*14*/ fTreeCascVarCascRadius = lXiRadius;
1287 /*15*/ fTreeCascVarV0Mass = lInvMassLambdaAsCascDghter;
1288 /*16*/ fTreeCascVarV0CosPointingAngle = lV0CosineOfPointingAngleXi;
1289 /*16*/ fTreeCascVarV0CosPointingAngleSpecial = lV0CosineOfPointingAngleXiSpecial;
1290 /*17*/ fTreeCascVarV0Radius = lV0RadiusXi;
1291 /*20*/ fTreeCascVarLeastNbrClusters = leastnumberofclusters;
1292 /*21*/ fTreeCascVarMultiplicity = lMultiplicity; //multiplicity, whatever that may be
1294 /*23*/ fTreeCascVarDistOverTotMom = TMath::Sqrt(
1295 TMath::Power( lPosXi[0] - lBestPrimaryVtxPos[0] , 2) +
1296 TMath::Power( lPosXi[1] - lBestPrimaryVtxPos[1] , 2) +
1297 TMath::Power( lPosXi[2] - lBestPrimaryVtxPos[2] , 2)
1299 /*23*/ fTreeCascVarDistOverTotMom /= (lXiTotMom+1e-13);
1301 //All vars not specified here: specified elsewhere!
1303 //------------------------------------------------
1305 //------------------------------------------------
1307 // The conditional is meant to decrease excessive
1308 // memory usage! Be careful when loosening the
1311 //Xi Mass window: 150MeV wide
1312 //Omega mass window: 150MeV wide
1314 if( (fTreeCascVarMassAsXi<1.32+0.075&&fTreeCascVarMassAsXi>1.32-0.075) ||
1315 (fTreeCascVarMassAsOmega<1.68+0.075&&fTreeCascVarMassAsOmega>1.68-0.075) ){
1316 fTreeCascade->Fill();
1319 //------------------------------------------------
1321 //------------------------------------------------
1323 }// end of the Cascade loop (ESD or AOD)
1325 // Post output data.
1326 PostData(1, fListHist);
1327 PostData(2, fTreeCascade);
1330 //________________________________________________________________________
1331 void AliAnalysisTaskExtractCascade::Terminate(Option_t *)
1333 // Draw result to the screen
1334 // Called once at the end of the query
1336 TList *cRetrievedList = 0x0;
1337 cRetrievedList = (TList*)GetOutputData(1);
1338 if(!cRetrievedList){
1339 Printf("ERROR - AliAnalysisTaskExtractCascade : ouput data container list not available\n");
1343 fHistV0MultiplicityForTrigEvt = dynamic_cast<TH1F*> ( cRetrievedList->FindObject("fHistV0MultiplicityForTrigEvt") );
1344 if (!fHistV0MultiplicityForTrigEvt) {
1345 Printf("ERROR - AliAnalysisTaskExtractCascade : fHistV0MultiplicityForTrigEvt not available");
1349 TCanvas *canCheck = new TCanvas("AliAnalysisTaskExtractCascade","V0 Multiplicity",10,10,510,510);
1350 canCheck->cd(1)->SetLogy();
1352 fHistV0MultiplicityForTrigEvt->SetMarkerStyle(22);
1353 fHistV0MultiplicityForTrigEvt->DrawCopy("E");
1356 //----------------------------------------------------------------------------
1358 Double_t AliAnalysisTaskExtractCascade::MyRapidity(Double_t rE, Double_t rPz) const
1360 // Local calculation for rapidity
1361 Double_t ReturnValue = -100;
1362 if( (rE-rPz+1.e-13) != 0 && (rE+rPz) != 0 ){
1363 ReturnValue = 0.5*TMath::Log((rE+rPz)/(rE-rPz+1.e-13));