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.
23 // --- Added bits of code for checking V0s
24 // (from AliAnalysisTaskCheckStrange.cxx)
26 // --- Algorithm Description
27 // 1. Perform Physics Selection
28 // 2. Perform Primary Vertex |z|<10cm selection
29 // 3. Perform Primary Vertex NoTPCOnly vertexing selection (>0 contrib.)
30 // 4. Perform Pileup Rejection
32 // 5a. Fill TTree object with V0 information, candidates
34 // Please Report Any Bugs!
36 // --- David Dobrigkeit Chinellato
37 // (david.chinellato@gmail.com)
39 // +-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+
45 //class AliMCEventHandler;
54 #include <Riostream.h>
59 #include "THnSparse.h"
65 #include "AliCentrality.h"
66 #include "AliESDEvent.h"
67 #include "AliAODEvent.h"
68 #include "AliV0vertexer.h"
69 #include "AliCascadeVertexer.h"
70 #include "AliESDpid.h"
71 #include "AliESDtrack.h"
72 #include "AliESDtrackCuts.h"
73 #include "AliInputEventHandler.h"
74 #include "AliAnalysisManager.h"
75 #include "AliMCEventHandler.h"
77 #include "AliCFContainer.h"
78 #include "AliMultiplicity.h"
80 #include "AliESDcascade.h"
81 #include "AliAODcascade.h"
82 #include "AliESDUtils.h"
83 #include "AliESDHeader.h"
85 #include "AliAnalysisTaskExtractV0.h"
88 #include "TObjectTable.h"
90 ClassImp(AliAnalysisTaskExtractV0)
92 AliAnalysisTaskExtractV0::AliAnalysisTaskExtractV0()
93 : AliAnalysisTaskSE(), fListHistV0(0), fTree(0), fPIDResponse(0),fESDtrackCuts(0),
94 fkIsNuclear ( kFALSE ),
95 fkLowEnergyPP ( kFALSE ),
96 fkUseOnTheFly ( kFALSE ),
98 //------------------------------------------------
100 // --- Filled on an Event-by-event basis
101 //------------------------------------------------
102 fHistV0MultiplicityBeforeTrigSel(0),
103 fHistV0MultiplicityForTrigEvt(0),
104 fHistV0MultiplicityForSelEvt(0),
105 fHistV0MultiplicityForSelEvtNoTPCOnly(0),
106 fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup(0),
107 fHistMultiplicityBeforeTrigSel(0),
108 fHistMultiplicityForTrigEvt(0),
109 fHistMultiplicity(0),
110 fHistMultiplicityNoTPCOnly(0),
111 fHistMultiplicityNoTPCOnlyNoPileup(0),
113 //Raw Data for Vertex Z position estimator change
114 f2dHistMultiplicityVsVertexZBeforeTrigSel(0),
115 f2dHistMultiplicityVsVertexZForTrigEvt(0),
116 f2dHistMultiplicityVsVertexZ(0),
117 f2dHistMultiplicityVsVertexZNoTPCOnly(0),
118 f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup(0),
126 fHistSwappedV0Counter(0)
131 AliAnalysisTaskExtractV0::AliAnalysisTaskExtractV0(const char *name)
132 : AliAnalysisTaskSE(name), fListHistV0(0), fTree(0), fPIDResponse(0),fESDtrackCuts(0),
133 fkIsNuclear ( kFALSE ),
134 fkLowEnergyPP ( kFALSE ),
135 fkUseOnTheFly ( kFALSE ),
137 //------------------------------------------------
139 // --- Filled on an Event-by-event basis
140 //------------------------------------------------
141 fHistV0MultiplicityBeforeTrigSel(0),
142 fHistV0MultiplicityForTrigEvt(0),
143 fHistV0MultiplicityForSelEvt(0),
144 fHistV0MultiplicityForSelEvtNoTPCOnly(0),
145 fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup(0),
146 fHistMultiplicityBeforeTrigSel(0),
147 fHistMultiplicityForTrigEvt(0),
148 fHistMultiplicity(0),
149 fHistMultiplicityNoTPCOnly(0),
150 fHistMultiplicityNoTPCOnlyNoPileup(0),
152 //Raw Data for Vertex Z position estimator change
153 f2dHistMultiplicityVsVertexZBeforeTrigSel(0),
154 f2dHistMultiplicityVsVertexZForTrigEvt(0),
155 f2dHistMultiplicityVsVertexZ(0),
156 f2dHistMultiplicityVsVertexZNoTPCOnly(0),
157 f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup(0),
165 fHistSwappedV0Counter(0)
168 // Output slot #0 writes into a TList container (Lambda Histos and fTree)
169 DefineOutput(1, TList::Class());
170 DefineOutput(2, TTree::Class());
174 AliAnalysisTaskExtractV0::~AliAnalysisTaskExtractV0()
176 //------------------------------------------------
178 //------------------------------------------------
188 //cleanup esd track cuts object too...
190 delete fESDtrackCuts;
199 //________________________________________________________________________
200 void AliAnalysisTaskExtractV0::UserCreateOutputObjects()
203 //Create File-resident Tree, please.
206 fTree = new TTree("fTree","V0Candidates");
208 //------------------------------------------------
209 // fTree Branch definitions
210 //------------------------------------------------
212 //-----------BASIC-INFO---------------------------
213 /* 1*/ fTree->Branch("fTreeVariableChi2V0",&fTreeVariableChi2V0,"fTreeVariableChi2V0/F");
214 /* 2*/ fTree->Branch("fTreeVariableDcaV0Daughters",&fTreeVariableDcaV0Daughters,"fTreeVariableDcaV0Daughters/F");
215 /* 3*/ fTree->Branch("fTreeVariableDcaPosToPrimVertex",&fTreeVariableDcaPosToPrimVertex,"fTreeVariableDcaPosToPrimVertex/F");
216 /* 4*/ fTree->Branch("fTreeVariableDcaNegToPrimVertex",&fTreeVariableDcaNegToPrimVertex,"fTreeVariableDcaNegToPrimVertex/F");
217 /* 5*/ fTree->Branch("fTreeVariableV0Radius",&fTreeVariableV0Radius,"fTreeVariableV0Radius/F");
218 /* 6*/ fTree->Branch("fTreeVariablePt",&fTreeVariablePt,"fTreeVariablePt/F");
219 /* 7*/ fTree->Branch("fTreeVariableRapK0Short",&fTreeVariableRapK0Short,"fTreeVariableRapK0Short/F");
220 /* 8*/ fTree->Branch("fTreeVariableRapLambda",&fTreeVariableRapLambda,"fTreeVariableRapLambda/F");
221 /* 9*/ fTree->Branch("fTreeVariableInvMassK0s",&fTreeVariableInvMassK0s,"fTreeVariableInvMassK0s/F");
222 /*10*/ fTree->Branch("fTreeVariableInvMassLambda",&fTreeVariableInvMassLambda,"fTreeVariableInvMassLambda/F");
223 /*11*/ fTree->Branch("fTreeVariableInvMassAntiLambda",&fTreeVariableInvMassAntiLambda,"fTreeVariableInvMassAntiLambda/F");
224 /*12*/ fTree->Branch("fTreeVariableV0CosineOfPointingAngle",&fTreeVariableV0CosineOfPointingAngle,"fTreeVariableV0CosineOfPointingAngle/F");
225 /*13*/ fTree->Branch("fTreeVariableAlphaV0",&fTreeVariableAlphaV0,"fTreeVariableAlphaV0/F");
226 /*14*/ fTree->Branch("fTreeVariablePtArmV0",&fTreeVariablePtArmV0,"fTreeVariablePtArmV0/F");
227 /*15*/ fTree->Branch("fTreeVariableLeastNbrCrossedRows",&fTreeVariableLeastNbrCrossedRows,"fTreeVariableLeastNbrCrossedRows/I");
228 /*16*/ fTree->Branch("fTreeVariableLeastRatioCrossedRowsOverFindable",&fTreeVariableLeastRatioCrossedRowsOverFindable,"fTreeVariableLeastRatioCrossedRowsOverFindable/F");
229 //-----------MULTIPLICITY-INFO--------------------
230 /*17*/ fTree->Branch("fTreeVariableMultiplicity",&fTreeVariableMultiplicity,"fTreeVariableMultiplicity/I");
231 //------------------------------------------------
232 /*18*/ fTree->Branch("fTreeVariableDistOverTotMom",&fTreeVariableDistOverTotMom,"fTreeVariableDistOverTotMom/F");
233 /*19*/ fTree->Branch("fTreeVariableNSigmasPosProton",&fTreeVariableNSigmasPosProton,"fTreeVariableNSigmasPosProton/F");
234 /*20*/ fTree->Branch("fTreeVariableNSigmasPosPion",&fTreeVariableNSigmasPosPion,"fTreeVariableNSigmasPosPion/F");
235 /*21*/ fTree->Branch("fTreeVariableNSigmasNegProton",&fTreeVariableNSigmasNegProton,"fTreeVariableNSigmasNegProton/F");
236 /*22*/ fTree->Branch("fTreeVariableNSigmasNegPion",&fTreeVariableNSigmasNegPion,"fTreeVariableNSigmasNegPion/F");
237 /*23*/ fTree->Branch("fTreeVariableNegEta",&fTreeVariableNegEta,"fTreeVariableNegEta/F");
238 /*24*/ fTree->Branch("fTreeVariablePosEta",&fTreeVariablePosEta,"fTreeVariablePosEta/F");
239 /*25*/ fTree->Branch("fTreeVariableRunNumber",&fTreeVariableRunNumber,"fTreeVariableRunNumber/I");
240 /*26*/ fTree->Branch("fTreeVariableEventNumber",&fTreeVariableEventNumber,"fTreeVariableEventNumber/l");
241 //-----------FOR CTAU DEBUGGING: Full Phase Space + Decay Position Info
242 fTree->Branch("fTreeVariableV0x",&fTreeVariableV0x,"fTreeVariableV0x/F");
243 fTree->Branch("fTreeVariableV0y",&fTreeVariableV0y,"fTreeVariableV0y/F");
244 fTree->Branch("fTreeVariableV0z",&fTreeVariableV0z,"fTreeVariableV0z/F");
246 fTree->Branch("fTreeVariableV0Px",&fTreeVariableV0Px,"fTreeVariableV0Px/F");
247 fTree->Branch("fTreeVariableV0Py",&fTreeVariableV0Py,"fTreeVariableV0Py/F");
248 fTree->Branch("fTreeVariableV0Pz",&fTreeVariableV0Pz,"fTreeVariableV0Pz/F");
249 //------------------------------------------------
250 // Particle Identification Setup
251 //------------------------------------------------
253 AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
254 AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
255 fPIDResponse = inputHandler->GetPIDResponse();
259 if(! fESDtrackCuts ){
260 fESDtrackCuts = new AliESDtrackCuts();
264 //------------------------------------------------
265 // V0 Multiplicity Histograms
266 //------------------------------------------------
269 //Create File-resident Tree, please.
272 fListHistV0 = new TList();
273 fListHistV0->SetOwner(); // See http://root.cern.ch/root/html/TCollection.html#TCollection:SetOwner
275 if(! fHistV0MultiplicityBeforeTrigSel) {
276 fHistV0MultiplicityBeforeTrigSel = new TH1F("fHistV0MultiplicityBeforeTrigSel",
277 "V0s per event (before Trig. Sel.);Nbr of V0s/Evt;Events",
279 fListHistV0->Add(fHistV0MultiplicityBeforeTrigSel);
282 if(! fHistV0MultiplicityForTrigEvt) {
283 fHistV0MultiplicityForTrigEvt = new TH1F("fHistV0MultiplicityForTrigEvt",
284 "V0s per event (for triggered evt);Nbr of V0s/Evt;Events",
286 fListHistV0->Add(fHistV0MultiplicityForTrigEvt);
289 if(! fHistV0MultiplicityForSelEvt) {
290 fHistV0MultiplicityForSelEvt = new TH1F("fHistV0MultiplicityForSelEvt",
291 "V0s per event;Nbr of V0s/Evt;Events",
293 fListHistV0->Add(fHistV0MultiplicityForSelEvt);
296 if(! fHistV0MultiplicityForSelEvtNoTPCOnly) {
297 fHistV0MultiplicityForSelEvtNoTPCOnly = new TH1F("fHistV0MultiplicityForSelEvtNoTPCOnly",
298 "V0s per event;Nbr of V0s/Evt;Events",
300 fListHistV0->Add(fHistV0MultiplicityForSelEvtNoTPCOnly);
303 if(! fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup) {
304 fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup = new TH1F("fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup",
305 "V0s per event;Nbr of V0s/Evt;Events",
307 fListHistV0->Add(fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup);
310 //------------------------------------------------
311 // Track Multiplicity Histograms
312 //------------------------------------------------
314 if(! fHistMultiplicityBeforeTrigSel) {
315 fHistMultiplicityBeforeTrigSel = new TH1F("fHistMultiplicityBeforeTrigSel",
316 "Tracks per event;Nbr of Tracks;Events",
318 fListHistV0->Add(fHistMultiplicityBeforeTrigSel);
320 if(! fHistMultiplicityForTrigEvt) {
321 fHistMultiplicityForTrigEvt = new TH1F("fHistMultiplicityForTrigEvt",
322 "Tracks per event;Nbr of Tracks;Events",
324 fListHistV0->Add(fHistMultiplicityForTrigEvt);
326 if(! fHistMultiplicity) {
327 fHistMultiplicity = new TH1F("fHistMultiplicity",
328 "Tracks per event;Nbr of Tracks;Events",
330 fListHistV0->Add(fHistMultiplicity);
332 if(! fHistMultiplicityNoTPCOnly) {
333 fHistMultiplicityNoTPCOnly = new TH1F("fHistMultiplicityNoTPCOnly",
334 "Tracks per event;Nbr of Tracks;Events",
336 fListHistV0->Add(fHistMultiplicityNoTPCOnly);
338 if(! fHistMultiplicityNoTPCOnlyNoPileup) {
339 fHistMultiplicityNoTPCOnlyNoPileup = new TH1F("fHistMultiplicityNoTPCOnlyNoPileup",
340 "Tracks per event;Nbr of Tracks;Events",
342 fListHistV0->Add(fHistMultiplicityNoTPCOnlyNoPileup);
347 //Raw Data for Vertex Z position estimator change
348 //TH2F *f2dHistMultiplicityVsVertexZBeforeTrigSel; //! multiplicity distribution
349 //TH2F *f2dHistMultiplicityVsVertexZForTrigEvt; //! multiplicity distribution
350 //TH2F *f2dHistMultiplicityVsVertexZ; //! multiplicity distribution
351 //TH2F *f2dHistMultiplicityVsVertexZNoTPCOnly; //! multiplicity distribution
352 //TH2F *f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup; //! multiplicity distribution
354 if(! f2dHistMultiplicityVsVertexZBeforeTrigSel) {
355 f2dHistMultiplicityVsVertexZBeforeTrigSel = new TH2F("f2dHistMultiplicityVsVertexZBeforeTrigSel",
356 "Tracks per event", 200, 0, 200,400, -20, 20);
357 fListHistV0->Add(f2dHistMultiplicityVsVertexZBeforeTrigSel);
359 if(! f2dHistMultiplicityVsVertexZForTrigEvt) {
360 f2dHistMultiplicityVsVertexZForTrigEvt = new TH2F("f2dHistMultiplicityVsVertexZForTrigEvt",
361 "Tracks per event", 200, 0, 200, 400, -20, 20);
362 fListHistV0->Add(f2dHistMultiplicityVsVertexZForTrigEvt);
364 if(! f2dHistMultiplicityVsVertexZ) {
365 f2dHistMultiplicityVsVertexZ = new TH2F("f2dHistMultiplicityVsVertexZ",
366 "Tracks per event", 200, 0, 200, 400, -20, 20);
367 fListHistV0->Add(f2dHistMultiplicityVsVertexZ);
369 if(! f2dHistMultiplicityVsVertexZNoTPCOnly) {
370 f2dHistMultiplicityVsVertexZNoTPCOnly = new TH2F("f2dHistMultiplicityVsVertexZNoTPCOnly",
371 "Tracks per event", 200, 0, 200, 400, -20, 20);
372 fListHistV0->Add(f2dHistMultiplicityVsVertexZNoTPCOnly);
374 if(! f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup) {
375 f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup = new TH2F("f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup",
376 "Tracks per event", 200, 0, 200, 400, -20, 20);
377 fListHistV0->Add(f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup);
380 //-----------------------------------------------------
383 fHistPVx = new TH1F("fHistPVx",
384 "PV x position;Nbr of Evts;x",
386 fListHistV0->Add(fHistPVx);
389 fHistPVy = new TH1F("fHistPVy",
390 "PV y position;Nbr of Evts;y",
392 fListHistV0->Add(fHistPVy);
395 fHistPVz = new TH1F("fHistPVz",
396 "PV z position;Nbr of Evts;z",
398 fListHistV0->Add(fHistPVz);
400 if(! fHistPVxAnalysis) {
401 fHistPVxAnalysis = new TH1F("fHistPVxAnalysis",
402 "PV x position;Nbr of Evts;x",
404 fListHistV0->Add(fHistPVxAnalysis);
406 if(! fHistPVyAnalysis) {
407 fHistPVyAnalysis = new TH1F("fHistPVyAnalysis",
408 "PV y position;Nbr of Evts;y",
410 fListHistV0->Add(fHistPVyAnalysis);
412 if(! fHistPVzAnalysis) {
413 fHistPVzAnalysis = new TH1F("fHistPVzAnalysis",
414 "PV z position;Nbr of Evts;z",
416 fListHistV0->Add(fHistPVzAnalysis);
418 if(! fHistSwappedV0Counter) {
419 fHistSwappedV0Counter = new TH1F("fHistSwappedV0Counter",
420 "Swap or not histo;Swapped (1) or not (0); count",
422 fListHistV0->Add(fHistSwappedV0Counter);
424 //Regular output: Histograms
425 PostData(1, fListHistV0);
426 //TTree Object: Saved to base directory. Should cache to disk while saving.
427 //(Important to avoid excessive memory usage, particularly when merging)
430 }// end UserCreateOutputObjects
433 //________________________________________________________________________
434 void AliAnalysisTaskExtractV0::UserExec(Option_t *)
437 // Called for each event
438 //gObjectTable->Print();
439 AliESDEvent *lESDevent = 0x0;
441 //AliAODEvent *lAODevent = 0x0;
444 Double_t lTrkgPrimaryVtxPos[3] = {-100.0, -100.0, -100.0};
445 Double_t lBestPrimaryVtxPos[3] = {-100.0, -100.0, -100.0};
446 Double_t lMagneticField = -10.;
448 // Connect to the InputEvent
449 // After these lines, we should have an ESD/AOD event + the number of cascades in it.
451 // Appropriate for ESD analysis!
453 lESDevent = dynamic_cast<AliESDEvent*>( InputEvent() );
455 AliWarning("ERROR: lESDevent not available \n");
458 fTreeVariableRunNumber = lESDevent->GetRunNumber();
459 fTreeVariableEventNumber =
460 ( ( ((ULong64_t)lESDevent->GetPeriodNumber() ) << 36 ) |
461 ( ((ULong64_t)lESDevent->GetOrbitNumber () ) << 12 ) |
462 ((ULong64_t)lESDevent->GetBunchCrossNumber() ) );
464 //REVISED multiplicity estimator after 'multiplicity day' (2011)
465 Int_t lMultiplicity = -100;
467 if(fkIsNuclear == kFALSE) lMultiplicity = fESDtrackCuts->GetReferenceMultiplicity(lESDevent, AliESDtrackCuts::kTrackletsITSTPC,0.5);
469 //---> If this is a nuclear collision, then go nuclear on "multiplicity" variable...
470 //---> Warning: Experimental
471 if(fkIsNuclear == kTRUE){
472 AliCentrality* centrality;
473 centrality = lESDevent->GetCentrality();
474 lMultiplicity = ( ( Int_t ) ( centrality->GetCentralityPercentile( "V0M" ) ) );
475 if (centrality->GetQuality()>1) {
476 PostData(1, fListHistV0);
482 //Set variable for filling tree afterwards!
483 //---> pp case......: GetReferenceMultiplicity
484 //---> Pb-Pb case...: Centrality by V0M
485 fTreeVariableMultiplicity = lMultiplicity;
486 fHistMultiplicityBeforeTrigSel->Fill ( lMultiplicity );
487 fHistV0MultiplicityBeforeTrigSel->Fill ( lESDevent->GetNumberOfV0s() );
489 //------------------------------------------------
491 //------------------------------------------------
494 UInt_t maskIsSelected = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
495 Bool_t isSelected = 0;
496 isSelected = (maskIsSelected & AliVEvent::kMB) == AliVEvent::kMB;
498 //pp at 2.76TeV: special case, ignore FastOnly
499 if ( (fkLowEnergyPP == kTRUE) && (maskIsSelected& AliVEvent::kFastOnly) ){
500 PostData(1, fListHistV0);
504 //Standard Min-Bias Selection
505 if ( ! isSelected ) {
506 PostData(1, fListHistV0);
511 //------------------------------------------------
512 // After Trigger Selection
513 //------------------------------------------------
515 nV0s = lESDevent->GetNumberOfV0s();
517 fHistV0MultiplicityForTrigEvt ->Fill( nV0s );
518 fHistMultiplicityForTrigEvt ->Fill( lMultiplicity );
520 //------------------------------------------------
521 // Getting: Primary Vertex + MagField Info
522 //------------------------------------------------
524 const AliESDVertex *lPrimaryTrackingESDVtx = lESDevent->GetPrimaryVertexTracks();
525 // get the vtx stored in ESD found with tracks
526 lPrimaryTrackingESDVtx->GetXYZ( lTrkgPrimaryVtxPos );
528 const AliESDVertex *lPrimaryBestESDVtx = lESDevent->GetPrimaryVertex();
529 // get the best primary vertex available for the event
530 // As done in AliCascadeVertexer, we keep the one which is the best one available.
531 // between : Tracking vertex > SPD vertex > TPC vertex > default SPD vertex
532 // This one will be used for next calculations (DCA essentially)
533 lPrimaryBestESDVtx->GetXYZ( lBestPrimaryVtxPos );
535 Double_t tPrimaryVtxPosition[3];
536 const AliVVertex *primaryVtx = lESDevent->GetPrimaryVertex();
537 tPrimaryVtxPosition[0] = primaryVtx->GetX();
538 tPrimaryVtxPosition[1] = primaryVtx->GetY();
539 tPrimaryVtxPosition[2] = primaryVtx->GetZ();
540 fHistPVx->Fill( tPrimaryVtxPosition[0] );
541 fHistPVy->Fill( tPrimaryVtxPosition[1] );
542 fHistPVz->Fill( tPrimaryVtxPosition[2] );
544 f2dHistMultiplicityVsVertexZForTrigEvt->Fill( lMultiplicity, tPrimaryVtxPosition[2] );
546 //------------------------------------------------
547 // Primary Vertex Z position: SKIP
548 //------------------------------------------------
550 if(TMath::Abs(lBestPrimaryVtxPos[2]) > 10.0 ) {
551 AliWarning("Pb / | Z position of Best Prim Vtx | > 10.0 cm ... return !");
552 PostData(1, fListHistV0);
557 f2dHistMultiplicityVsVertexZ->Fill( lMultiplicity, tPrimaryVtxPosition[2] );
559 lMagneticField = lESDevent->GetMagneticField( );
560 fHistV0MultiplicityForSelEvt ->Fill( nV0s );
561 fHistMultiplicity->Fill(lMultiplicity);
563 //------------------------------------------------
564 // Only look at events with well-established PV
565 //------------------------------------------------
567 const AliESDVertex *lPrimaryTrackingESDVtxCheck = lESDevent->GetPrimaryVertexTracks();
568 const AliESDVertex *lPrimarySPDVtx = lESDevent->GetPrimaryVertexSPD();
569 if (!lPrimarySPDVtx->GetStatus() && !lPrimaryTrackingESDVtxCheck->GetStatus() ){
570 AliWarning("Pb / No SPD prim. vertex nor prim. Tracking vertex ... return !");
571 PostData(1, fListHistV0);
577 f2dHistMultiplicityVsVertexZNoTPCOnly->Fill( lMultiplicity, tPrimaryVtxPosition[2] );
578 fHistV0MultiplicityForSelEvtNoTPCOnly ->Fill( nV0s );
579 fHistMultiplicityNoTPCOnly->Fill(lMultiplicity);
581 //------------------------------------------------
583 //------------------------------------------------
585 // FIXME : quality selection regarding pile-up rejection
586 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
587 PostData(1, fListHistV0);
592 f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup->Fill( lMultiplicity, tPrimaryVtxPosition[2] );
593 fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup ->Fill( nV0s );
594 fHistMultiplicityNoTPCOnlyNoPileup->Fill(lMultiplicity);
596 //------------------------------------------------
597 // MAIN LAMBDA LOOP STARTS HERE
598 //------------------------------------------------
600 if( ! (lESDevent->GetPrimaryVertex()->IsFromVertexerZ() ) ){
601 fHistPVxAnalysis->Fill( tPrimaryVtxPosition[0] );
602 fHistPVyAnalysis->Fill( tPrimaryVtxPosition[1] );
603 fHistPVzAnalysis->Fill( tPrimaryVtxPosition[2] );
606 //Variable definition
607 Int_t lOnFlyStatus = 0;// nv0sOn = 0, nv0sOff = 0;
608 Double_t lChi2V0 = 0;
609 Double_t lDcaV0Daughters = 0, lDcaV0ToPrimVertex = 0;
610 Double_t lDcaPosToPrimVertex = 0, lDcaNegToPrimVertex = 0;
611 Double_t lV0CosineOfPointingAngle = 0;
612 Double_t lV0Radius = 0, lPt = 0;
613 Double_t lRapK0Short = 0, lRapLambda = 0;
614 Double_t lInvMassK0s = 0, lInvMassLambda = 0, lInvMassAntiLambda = 0;
615 Double_t lAlphaV0 = 0, lPtArmV0 = 0;
617 Double_t fMinV0Pt = 0;
618 Double_t fMaxV0Pt = 100;
621 nv0s = lESDevent->GetNumberOfV0s();
623 //for (Int_t iV0 = 0; iV0 < nv0s; iV0++)
624 for (Int_t iV0 = 0; iV0 < nv0s; iV0++) //extra-crazy test
625 {// This is the begining of the V0 loop
626 AliESDv0 *v0 = ((AliESDEvent*)lESDevent)->GetV0(iV0);
629 //---> Fix On-the-Fly candidates, count how many swapped
630 if( v0->GetParamN()->Charge() > 0 && v0->GetParamP()->Charge() < 0 ){
631 fHistSwappedV0Counter -> Fill( 1 );
633 fHistSwappedV0Counter -> Fill( 0 );
635 if ( fkUseOnTheFly ) CheckChargeV0(v0);
637 Double_t tDecayVertexV0[3]; v0->GetXYZ(tDecayVertexV0[0],tDecayVertexV0[1],tDecayVertexV0[2]);
640 v0->GetPxPyPz( tV0mom[0],tV0mom[1],tV0mom[2] );
641 Double_t lV0TotalMomentum = TMath::Sqrt(
642 tV0mom[0]*tV0mom[0]+tV0mom[1]*tV0mom[1]+tV0mom[2]*tV0mom[2] );
644 lV0Radius = TMath::Sqrt(tDecayVertexV0[0]*tDecayVertexV0[0]+tDecayVertexV0[1]*tDecayVertexV0[1]);
646 //Set Variables for later filling
647 fTreeVariableV0x = tDecayVertexV0[0];
648 fTreeVariableV0y = tDecayVertexV0[1];
649 fTreeVariableV0z = tDecayVertexV0[2];
651 //Set Variables for later filling
652 fTreeVariableV0Px = tV0mom[0];
653 fTreeVariableV0Py = tV0mom[1];
654 fTreeVariableV0Pz = tV0mom[2];
657 lRapK0Short = v0->RapK0Short();
658 lRapLambda = v0->RapLambda();
659 if ((lPt<fMinV0Pt)||(fMaxV0Pt<lPt)) continue;
661 UInt_t lKeyPos = (UInt_t)TMath::Abs(v0->GetPindex());
662 UInt_t lKeyNeg = (UInt_t)TMath::Abs(v0->GetNindex());
664 Double_t lMomPos[3]; v0->GetPPxPyPz(lMomPos[0],lMomPos[1],lMomPos[2]);
665 Double_t lMomNeg[3]; v0->GetNPxPyPz(lMomNeg[0],lMomNeg[1],lMomNeg[2]);
667 AliESDtrack *pTrack=((AliESDEvent*)lESDevent)->GetTrack(lKeyPos);
668 AliESDtrack *nTrack=((AliESDEvent*)lESDevent)->GetTrack(lKeyNeg);
669 if (!pTrack || !nTrack) {
670 Printf("ERROR: Could not retreive one of the daughter track");
674 //Daughter Eta for Eta selection, afterwards
675 fTreeVariableNegEta = nTrack->Eta();
676 fTreeVariablePosEta = pTrack->Eta();
678 // Filter like-sign V0 (next: add counter and distribution)
679 if ( pTrack->GetSign() == nTrack->GetSign()){
683 //________________________________________________________________________
684 // Track quality cuts
685 Float_t lPosTrackCrossedRows = pTrack->GetTPCClusterInfo(2,1);
686 Float_t lNegTrackCrossedRows = nTrack->GetTPCClusterInfo(2,1);
687 fTreeVariableLeastNbrCrossedRows = (Int_t) lPosTrackCrossedRows;
688 if( lNegTrackCrossedRows < fTreeVariableLeastNbrCrossedRows )
689 fTreeVariableLeastNbrCrossedRows = (Int_t) lNegTrackCrossedRows;
691 // TPC refit condition (done during reconstruction for Offline but not for On-the-fly)
692 if( !(pTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;
693 if( !(nTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;
695 if ( ( ( pTrack->GetTPCClusterInfo(2,1) ) < 70 ) || ( ( nTrack->GetTPCClusterInfo(2,1) ) < 70 ) ) continue;
697 //GetKinkIndex condition
698 if( pTrack->GetKinkIndex(0)>0 || nTrack->GetKinkIndex(0)>0 ) continue;
700 //Findable clusters > 0 condition
701 if( pTrack->GetTPCNclsF()<=0 || nTrack->GetTPCNclsF()<=0 ) continue;
703 //Compute ratio Crossed Rows / Findable clusters
704 //Note: above test avoids division by zero!
705 Float_t lPosTrackCrossedRowsOverFindable = lPosTrackCrossedRows / ((double)(pTrack->GetTPCNclsF()));
706 Float_t lNegTrackCrossedRowsOverFindable = lNegTrackCrossedRows / ((double)(nTrack->GetTPCNclsF()));
708 fTreeVariableLeastRatioCrossedRowsOverFindable = lPosTrackCrossedRowsOverFindable;
709 if( lNegTrackCrossedRowsOverFindable < fTreeVariableLeastRatioCrossedRowsOverFindable )
710 fTreeVariableLeastRatioCrossedRowsOverFindable = lNegTrackCrossedRowsOverFindable;
712 //Lowest Cut Level for Ratio Crossed Rows / Findable = 0.8, set here
713 if ( fTreeVariableLeastRatioCrossedRowsOverFindable < 0.8 ) continue;
715 //End track Quality Cuts
716 //________________________________________________________________________
718 lDcaPosToPrimVertex = TMath::Abs(pTrack->GetD(tPrimaryVtxPosition[0],
719 tPrimaryVtxPosition[1],
722 lDcaNegToPrimVertex = TMath::Abs(nTrack->GetD(tPrimaryVtxPosition[0],
723 tPrimaryVtxPosition[1],
726 lOnFlyStatus = v0->GetOnFlyStatus();
727 lChi2V0 = v0->GetChi2V0();
728 lDcaV0Daughters = v0->GetDcaV0Daughters();
729 lDcaV0ToPrimVertex = v0->GetD(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],tPrimaryVtxPosition[2]);
730 lV0CosineOfPointingAngle = v0->GetV0CosineOfPointingAngle(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],tPrimaryVtxPosition[2]);
731 fTreeVariableV0CosineOfPointingAngle=lV0CosineOfPointingAngle;
733 // Getting invariant mass infos directly from ESD
734 v0->ChangeMassHypothesis(310);
735 lInvMassK0s = v0->GetEffMass();
736 v0->ChangeMassHypothesis(3122);
737 lInvMassLambda = v0->GetEffMass();
738 v0->ChangeMassHypothesis(-3122);
739 lInvMassAntiLambda = v0->GetEffMass();
740 lAlphaV0 = v0->AlphaV0();
741 lPtArmV0 = v0->PtArmV0();
743 fTreeVariablePt = v0->Pt();
744 fTreeVariableChi2V0 = lChi2V0;
745 fTreeVariableDcaV0ToPrimVertex = lDcaV0ToPrimVertex;
746 fTreeVariableDcaV0Daughters = lDcaV0Daughters;
747 fTreeVariableV0CosineOfPointingAngle = lV0CosineOfPointingAngle;
748 fTreeVariableV0Radius = lV0Radius;
749 fTreeVariableDcaPosToPrimVertex = lDcaPosToPrimVertex;
750 fTreeVariableDcaNegToPrimVertex = lDcaNegToPrimVertex;
751 fTreeVariableInvMassK0s = lInvMassK0s;
752 fTreeVariableInvMassLambda = lInvMassLambda;
753 fTreeVariableInvMassAntiLambda = lInvMassAntiLambda;
754 fTreeVariableRapK0Short = lRapK0Short;
755 fTreeVariableRapLambda = lRapLambda;
756 fTreeVariableAlphaV0 = lAlphaV0;
757 fTreeVariablePtArmV0 = lPtArmV0;
759 //Official means of acquiring N-sigmas
760 fTreeVariableNSigmasPosProton = fPIDResponse->NumberOfSigmasTPC( pTrack, AliPID::kProton );
761 fTreeVariableNSigmasPosPion = fPIDResponse->NumberOfSigmasTPC( pTrack, AliPID::kPion );
762 fTreeVariableNSigmasNegProton = fPIDResponse->NumberOfSigmasTPC( nTrack, AliPID::kProton );
763 fTreeVariableNSigmasNegPion = fPIDResponse->NumberOfSigmasTPC( nTrack, AliPID::kPion );
765 //This requires an Invariant Mass Hypothesis afterwards
766 fTreeVariableDistOverTotMom = TMath::Sqrt(
767 TMath::Power( tDecayVertexV0[0] - lBestPrimaryVtxPos[0] , 2) +
768 TMath::Power( tDecayVertexV0[1] - lBestPrimaryVtxPos[1] , 2) +
769 TMath::Power( tDecayVertexV0[2] - lBestPrimaryVtxPos[2] , 2)
771 fTreeVariableDistOverTotMom /= (lV0TotalMomentum+1e-10); //avoid division by zero, to be sure
773 //------------------------------------------------
775 //------------------------------------------------
777 // The conditionals are meant to decrease excessive
780 //First Selection: Reject OnFly
781 if( (lOnFlyStatus == 0 && fkUseOnTheFly == kFALSE) || (lOnFlyStatus != 0 && fkUseOnTheFly == kTRUE ) ){
782 //Second Selection: rough 20-sigma band, parametric.
783 //K0Short: Enough to parametrize peak broadening with linear function.
784 Double_t lUpperLimitK0Short = (5.63707e-01) + (1.14979e-02)*fTreeVariablePt;
785 Double_t lLowerLimitK0Short = (4.30006e-01) - (1.10029e-02)*fTreeVariablePt;
786 //Lambda: Linear (for higher pt) plus exponential (for low-pt broadening)
787 //[0]+[1]*x+[2]*TMath::Exp(-[3]*x)
788 Double_t lUpperLimitLambda = (1.13688e+00) + (5.27838e-03)*fTreeVariablePt + (8.42220e-02)*TMath::Exp(-(3.80595e+00)*fTreeVariablePt);
789 Double_t lLowerLimitLambda = (1.09501e+00) - (5.23272e-03)*fTreeVariablePt - (7.52690e-02)*TMath::Exp(-(3.46339e+00)*fTreeVariablePt);
791 if( (fTreeVariableInvMassLambda < lUpperLimitLambda && fTreeVariableInvMassLambda > lLowerLimitLambda ) ||
792 (fTreeVariableInvMassAntiLambda < lUpperLimitLambda && fTreeVariableInvMassAntiLambda > lLowerLimitLambda ) ||
793 (fTreeVariableInvMassK0s < lUpperLimitK0Short && fTreeVariableInvMassK0s > lLowerLimitK0Short ) ){
794 //Pre-selection in case this is AA...
795 if( fkIsNuclear == kFALSE ) fTree->Fill();
796 if( fkIsNuclear == kTRUE){
797 //If this is a nuclear collision___________________
798 // ... pre-filter with TPC, daughter eta selection
799 if( (fTreeVariableInvMassLambda < lUpperLimitLambda && fTreeVariableInvMassLambda > lLowerLimitLambda
800 && TMath::Abs(fTreeVariableNSigmasPosProton) < 6.0 && TMath::Abs(fTreeVariableNSigmasNegPion) < 6.0 ) ||
801 (fTreeVariableInvMassAntiLambda < lUpperLimitLambda && fTreeVariableInvMassAntiLambda > lLowerLimitLambda
802 && TMath::Abs(fTreeVariableNSigmasNegProton) < 6.0 && TMath::Abs(fTreeVariableNSigmasPosPion) < 6.0 ) ||
803 (fTreeVariableInvMassK0s < lUpperLimitK0Short && fTreeVariableInvMassK0s > lLowerLimitK0Short
804 && TMath::Abs(fTreeVariableNSigmasNegPion) < 6.0 && TMath::Abs(fTreeVariableNSigmasPosPion) < 6.0 ) ){
806 if ( TMath::Abs(fTreeVariableNegEta)<0.8 && TMath::Abs(fTreeVariablePosEta)<0.8 ) fTree->Fill();
808 }//end nuclear_____________________________________
812 //------------------------------------------------
814 //------------------------------------------------
816 }// This is the end of the V0 loop
819 PostData(1, fListHistV0);
825 //________________________________________________________________________
826 void AliAnalysisTaskExtractV0::Terminate(Option_t *)
828 // Draw result to the screen
829 // Called once at the end of the query
830 // This will draw the V0 candidate multiplicity, whose
831 // number of entries corresponds to the number of triggered events.
832 TList *cRetrievedList = 0x0;
833 cRetrievedList = (TList*)GetOutputData(1);
835 Printf("ERROR - AliAnalysisTaskExtractV0 : ouput data container list not available\n");
838 fHistV0MultiplicityForTrigEvt = dynamic_cast<TH1F*> ( cRetrievedList->FindObject("fHistV0MultiplicityForTrigEvt") );
839 if (!fHistV0MultiplicityForTrigEvt) {
840 Printf("ERROR - AliAnalysisTaskExtractV0 : fHistV0MultiplicityForTrigEvt not available");
843 TCanvas *canCheck = new TCanvas("AliAnalysisTaskExtractV0","V0 Multiplicity",10,10,510,510);
844 canCheck->cd(1)->SetLogy();
845 fHistV0MultiplicityForTrigEvt->SetMarkerStyle(22);
846 fHistV0MultiplicityForTrigEvt->DrawCopy("E");
849 //________________________________________________________________________
850 void AliAnalysisTaskExtractV0::CheckChargeV0(AliESDv0 *v0)
852 // This function checks charge of negative and positive daughter tracks.
853 // If incorrectly defined (onfly vertexer), swaps out.
854 if( v0->GetParamN()->Charge() > 0 && v0->GetParamP()->Charge() < 0 ){
855 //V0 daughter track swapping is required! Note: everything is swapped here... P->N, N->P
856 Long_t lCorrectNidx = v0->GetPindex();
857 Long_t lCorrectPidx = v0->GetNindex();
858 Double32_t lCorrectNmom[3];
859 Double32_t lCorrectPmom[3];
860 v0->GetPPxPyPz( lCorrectNmom[0], lCorrectNmom[1], lCorrectNmom[2] );
861 v0->GetNPxPyPz( lCorrectPmom[0], lCorrectPmom[1], lCorrectPmom[2] );
863 AliExternalTrackParam lCorrectParamN(
864 v0->GetParamP()->GetX() ,
865 v0->GetParamP()->GetAlpha() ,
866 v0->GetParamP()->GetParameter() ,
867 v0->GetParamP()->GetCovariance()
869 AliExternalTrackParam lCorrectParamP(
870 v0->GetParamN()->GetX() ,
871 v0->GetParamN()->GetAlpha() ,
872 v0->GetParamN()->GetParameter() ,
873 v0->GetParamN()->GetCovariance()
875 lCorrectParamN.SetMostProbablePt( v0->GetParamP()->GetMostProbablePt() );
876 lCorrectParamP.SetMostProbablePt( v0->GetParamN()->GetMostProbablePt() );
878 //Get Variables___________________________________________________
879 Double_t lDcaV0Daughters = v0 -> GetDcaV0Daughters();
880 Double_t lCosPALocal = v0 -> GetV0CosineOfPointingAngle();
881 Bool_t lOnFlyStatusLocal = v0 -> GetOnFlyStatus();
883 //Create Replacement Object_______________________________________
884 AliESDv0 *v0correct = new AliESDv0(lCorrectParamN,lCorrectNidx,lCorrectParamP,lCorrectPidx);
885 v0correct->SetDcaV0Daughters ( lDcaV0Daughters );
886 v0correct->SetV0CosineOfPointingAngle ( lCosPALocal );
887 v0correct->ChangeMassHypothesis ( kK0Short );
888 v0correct->SetOnFlyStatus ( lOnFlyStatusLocal );
890 //Reverse Cluster info..._________________________________________
891 v0correct->SetClusters( v0->GetClusters( 1 ), v0->GetClusters ( 0 ) );
894 //Proper cleanup..._______________________________________________
898 //Just another cross-check and output_____________________________
899 if( v0->GetParamN()->Charge() > 0 && v0->GetParamP()->Charge() < 0 ) {
900 AliWarning("Found Swapped Charges, tried to correct but something FAILED!");
902 //AliWarning("Found Swapped Charges and fixed.");
904 //________________________________________________________________
906 //Don't touch it! ---
907 //Printf("Ah, nice. Charges are already ordered...");