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 "AliAnalysisUtils.h"
86 #include "AliAnalysisTaskExtractV0pPb.h"
89 #include "TObjectTable.h"
91 ClassImp(AliAnalysisTaskExtractV0pPb)
93 AliAnalysisTaskExtractV0pPb::AliAnalysisTaskExtractV0pPb()
94 : AliAnalysisTaskSE(), fListHistV0(0), fTree(0), fTreeEvents(0), fPIDResponse(0),
95 fTPCdEdxSelection ( kTRUE ),
96 //------------------------------------------------
98 fTreeVariableChi2V0(0),
99 fTreeVariableDcaV0Daughters(0),
100 fTreeVariableDcaV0ToPrimVertex(0),
101 fTreeVariableDcaPosToPrimVertex(0),
102 fTreeVariableDcaNegToPrimVertex(0),
103 fTreeVariableV0CosineOfPointingAngle(0),
104 fTreeVariableV0Radius(0),
106 fTreeVariableRapK0Short(0),
107 fTreeVariableRapLambda(0),
108 fTreeVariableInvMassK0s(0),
109 fTreeVariableInvMassLambda(0),
110 fTreeVariableInvMassAntiLambda(0),
111 fTreeVariableAlphaV0(0),
112 fTreeVariablePtArmV0(0),
113 fTreeVariableNegTotMomentum(0),
114 fTreeVariablePosTotMomentum(0),
115 fTreeVariableNegdEdxSig(0),
116 fTreeVariablePosdEdxSig(0),
117 fTreeVariableNegEta(0),
118 fTreeVariablePosEta(0),
120 fTreeVariableNSigmasPosProton(0),
121 fTreeVariableNSigmasPosPion(0),
122 fTreeVariableNSigmasNegProton(0),
123 fTreeVariableNSigmasNegPion(0),
125 fTreeVariableDistOverTotMom(0),
126 fTreeVariableLeastNbrCrossedRows(0),
127 fTreeVariableLeastRatioCrossedRowsOverFindable(0),
128 fTreeVariableCentrality(0),
131 fTreeEventsCentrality(0),
133 //------------------------------------------------
135 // --- Filled on an Event-by-event basis
136 //------------------------------------------------
138 fHistCentralityProcessed(0),
139 fHistCentralityTrigEvt(0),
140 fHistCentralityHasVtx(0),
141 fHistCentralityVtxZ(0)
146 AliAnalysisTaskExtractV0pPb::AliAnalysisTaskExtractV0pPb(const char *name)
147 : AliAnalysisTaskSE(name), fListHistV0(0), fTree(0), fTreeEvents(0), fPIDResponse(0),
148 fTPCdEdxSelection ( kTRUE ),
149 //------------------------------------------------
151 fTreeVariableChi2V0(0),
152 fTreeVariableDcaV0Daughters(0),
153 fTreeVariableDcaV0ToPrimVertex(0),
154 fTreeVariableDcaPosToPrimVertex(0),
155 fTreeVariableDcaNegToPrimVertex(0),
156 fTreeVariableV0CosineOfPointingAngle(0),
157 fTreeVariableV0Radius(0),
159 fTreeVariableRapK0Short(0),
160 fTreeVariableRapLambda(0),
161 fTreeVariableInvMassK0s(0),
162 fTreeVariableInvMassLambda(0),
163 fTreeVariableInvMassAntiLambda(0),
164 fTreeVariableAlphaV0(0),
165 fTreeVariablePtArmV0(0),
166 fTreeVariableNegTotMomentum(0),
167 fTreeVariablePosTotMomentum(0),
168 fTreeVariableNegdEdxSig(0),
169 fTreeVariablePosdEdxSig(0),
170 fTreeVariableNegEta(0),
171 fTreeVariablePosEta(0),
173 fTreeVariableNSigmasPosProton(0),
174 fTreeVariableNSigmasPosPion(0),
175 fTreeVariableNSigmasNegProton(0),
176 fTreeVariableNSigmasNegPion(0),
178 fTreeVariableDistOverTotMom(0),
179 fTreeVariableLeastNbrCrossedRows(0),
180 fTreeVariableLeastRatioCrossedRowsOverFindable(0),
181 fTreeVariableCentrality(0),
184 fTreeEventsCentrality(0),
186 //------------------------------------------------
188 // --- Filled on an Event-by-event basis
189 //------------------------------------------------
191 fHistCentralityProcessed(0),
192 fHistCentralityTrigEvt(0),
193 fHistCentralityHasVtx(0),
194 fHistCentralityVtxZ(0)
196 // Output slot #0 writes into a TList container (Lambda Histos and fTree)
197 DefineOutput(1, TList::Class());
198 DefineOutput(2, TTree::Class());
199 DefineOutput(3, TTree::Class());
203 AliAnalysisTaskExtractV0pPb::~AliAnalysisTaskExtractV0pPb()
205 //------------------------------------------------
207 //------------------------------------------------
225 //________________________________________________________________________
226 void AliAnalysisTaskExtractV0pPb::UserCreateOutputObjects()
229 //Create File-resident Tree, please.
232 fTree = new TTree("fTree","V0Candidates");
234 //------------------------------------------------
235 // fTree Branch definitions
236 //------------------------------------------------
238 //-----------BASIC-INFO---------------------------
239 /* 1*/ fTree->Branch("fTreeVariableChi2V0",&fTreeVariableChi2V0,"fTreeVariableChi2V0/F");
240 /* 2*/ fTree->Branch("fTreeVariableDcaV0Daughters",&fTreeVariableDcaV0Daughters,"fTreeVariableDcaV0Daughters/F");
241 /* 3*/ fTree->Branch("fTreeVariableDcaPosToPrimVertex",&fTreeVariableDcaPosToPrimVertex,"fTreeVariableDcaPosToPrimVertex/F");
242 /* 4*/ fTree->Branch("fTreeVariableDcaNegToPrimVertex",&fTreeVariableDcaNegToPrimVertex,"fTreeVariableDcaNegToPrimVertex/F");
243 /* 5*/ fTree->Branch("fTreeVariableV0Radius",&fTreeVariableV0Radius,"fTreeVariableV0Radius/F");
244 /* 6*/ fTree->Branch("fTreeVariablePt",&fTreeVariablePt,"fTreeVariablePt/F");
245 /* 7*/ fTree->Branch("fTreeVariableRapK0Short",&fTreeVariableRapK0Short,"fTreeVariableRapK0Short/F");
246 /* 8*/ fTree->Branch("fTreeVariableRapLambda",&fTreeVariableRapLambda,"fTreeVariableRapLambda/F");
247 /* 9*/ fTree->Branch("fTreeVariableInvMassK0s",&fTreeVariableInvMassK0s,"fTreeVariableInvMassK0s/F");
248 /*10*/ fTree->Branch("fTreeVariableInvMassLambda",&fTreeVariableInvMassLambda,"fTreeVariableInvMassLambda/F");
249 /*11*/ fTree->Branch("fTreeVariableInvMassAntiLambda",&fTreeVariableInvMassAntiLambda,"fTreeVariableInvMassAntiLambda/F");
250 /*12*/ fTree->Branch("fTreeVariableV0CosineOfPointingAngle",&fTreeVariableV0CosineOfPointingAngle,"fTreeVariableV0CosineOfPointingAngle/F");
251 /*13*/ fTree->Branch("fTreeVariableAlphaV0",&fTreeVariableAlphaV0,"fTreeVariableAlphaV0/F");
252 /*14*/ fTree->Branch("fTreeVariablePtArmV0",&fTreeVariablePtArmV0,"fTreeVariablePtArmV0/F");
253 /*15*/ fTree->Branch("fTreeVariableLeastNbrCrossedRows",&fTreeVariableLeastNbrCrossedRows,"fTreeVariableLeastNbrCrossedRows/I");
254 /*16*/ fTree->Branch("fTreeVariableLeastRatioCrossedRowsOverFindable",&fTreeVariableLeastRatioCrossedRowsOverFindable,"fTreeVariableLeastRatioCrossedRowsOverFindable/F");
255 //-----------MULTIPLICITY-INFO--------------------
256 /*17*/ fTree->Branch("fTreeVariableCentrality",&fTreeVariableCentrality,"fTreeVariableCentrality/F");
257 //------------------------------------------------
258 /*18*/ fTree->Branch("fTreeVariableDistOverTotMom",&fTreeVariableDistOverTotMom,"fTreeVariableDistOverTotMom/F");
259 /*19*/ fTree->Branch("fTreeVariableNSigmasPosProton",&fTreeVariableNSigmasPosProton,"fTreeVariableNSigmasPosProton/F");
260 /*20*/ fTree->Branch("fTreeVariableNSigmasPosPion",&fTreeVariableNSigmasPosPion,"fTreeVariableNSigmasPosPion/F");
261 /*21*/ fTree->Branch("fTreeVariableNSigmasNegProton",&fTreeVariableNSigmasNegProton,"fTreeVariableNSigmasNegProton/F");
262 /*22*/ fTree->Branch("fTreeVariableNSigmasNegPion",&fTreeVariableNSigmasNegPion,"fTreeVariableNSigmasNegPion/F");
263 /*23*/ fTree->Branch("fTreeVariableNegEta",&fTreeVariableNegEta,"fTreeVariableNegEta/F");
264 /*24*/ fTree->Branch("fTreeVariablePosEta",&fTreeVariablePosEta,"fTreeVariablePosEta/F");
266 //------------------------------------------------
267 // fTreeEvents Branch definitions
268 //------------------------------------------------
270 fTreeEvents = new TTree("fTreeEvents","Events");
271 /*01*/ fTreeEvents->Branch("fTreeEventsCentrality",&fTreeEventsCentrality,"fTreeEventsCentrality/F");
273 //------------------------------------------------
274 // Particle Identification Setup
275 //------------------------------------------------
277 AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
278 AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
279 fPIDResponse = inputHandler->GetPIDResponse();
281 //------------------------------------------------
282 // V0 Multiplicity Histograms
283 //------------------------------------------------
286 //Create File-resident Tree, please.
289 fListHistV0 = new TList();
290 fListHistV0->SetOwner(); // See http://root.cern.ch/root/html/TCollection.html#TCollection:SetOwner
292 //------------------------------------------------
293 // V0A Centrality Histograms
294 //------------------------------------------------
296 //Default V0A Centrality (if PbPb)
297 if(! fHistCentralityProcessed) {
298 fHistCentralityProcessed = new TH1F("fHistCentralityProcessed",
299 "All processed Events;V0A Centrality;Events",
301 fListHistV0->Add(fHistCentralityProcessed);
303 if(! fHistCentralityTrigEvt) {
304 fHistCentralityTrigEvt = new TH1F("fHistCentralityTrigEvt",
305 "PS selected Events;V0A Centrality;Events",
307 fListHistV0->Add(fHistCentralityTrigEvt);
309 if(! fHistCentralityHasVtx) {
310 fHistCentralityHasVtx = new TH1F("fHistCentralityHasVtx",
311 "Events having Vertex;V0A Centrality;Events",
313 fListHistV0->Add(fHistCentralityHasVtx);
315 if(! fHistCentralityVtxZ) {
316 fHistCentralityVtxZ = new TH1F("fHistCentralityVtxZ",
317 "Vertex |z|<10cm;V0A Centrality;Events",
319 fListHistV0->Add(fHistCentralityVtxZ);
322 //Regular output: Histograms
323 PostData(1, fListHistV0);
324 //TTree Object: Saved to base directory. Should cache to disk while saving.
325 //(Important to avoid excessive memory usage, particularly when merging)
327 PostData(3, fTreeEvents);
329 }// end UserCreateOutputObjects
332 //________________________________________________________________________
333 void AliAnalysisTaskExtractV0pPb::UserExec(Option_t *)
336 // Called for each event
337 //gObjectTable->Print();
338 AliESDEvent *lESDevent = 0x0;
340 //AliAODEvent *lAODevent = 0x0;
342 Double_t lTrkgPrimaryVtxPos[3] = {-100.0, -100.0, -100.0};
343 Double_t lBestPrimaryVtxPos[3] = {-100.0, -100.0, -100.0};
344 Double_t lMagneticField = -10.;
346 // Connect to the InputEvent
347 // After these lines, we should have an ESD/AOD event + the number of cascades in it.
349 // Appropriate for ESD analysis!
351 lESDevent = dynamic_cast<AliESDEvent*>( InputEvent() );
353 AliWarning("ERROR: lESDevent not available \n");
357 //REVISED multiplicity estimator after 'multiplicity day' (2011)
358 Float_t lMultiplicity = -100;
359 fTreeVariableCentrality = -100;
361 //---> If this is a nuclear collision, then go nuclear on "multiplicity" variable...
362 //---> Warning: Experimental
363 AliCentrality* centrality;
364 centrality = lESDevent->GetCentrality();
365 if (centrality->GetQuality()>1) {
366 PostData(1, fListHistV0);
368 PostData(3, fTreeEvents);
371 lMultiplicity = centrality->GetCentralityPercentile( "V0A" );
374 //Set variable for filling tree afterwards!
376 fTreeVariableCentrality = lMultiplicity;
377 fTreeEventsCentrality = lMultiplicity;
379 fHistCentralityProcessed->Fill ( fTreeVariableCentrality );
382 //------------------------------------------------
384 //------------------------------------------------
387 UInt_t maskIsSelected = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
388 Bool_t isSelected = 0;
389 isSelected = (maskIsSelected & AliVEvent::kINT7) == AliVEvent::kINT7;
391 //Trigger Selection for kINT7
393 PostData(1, fListHistV0);
395 PostData(3, fTreeEvents);
399 //------------------------------------------------
400 // After Trigger Selection
401 //------------------------------------------------
403 fHistCentralityTrigEvt -> Fill( fTreeVariableCentrality );
405 //------------------------------------------------
406 // Getting: Primary Vertex + MagField Info
407 //------------------------------------------------
409 const AliESDVertex *lPrimaryTrackingESDVtx = lESDevent->GetPrimaryVertexTracks();
410 // get the vtx stored in ESD found with tracks
411 lPrimaryTrackingESDVtx->GetXYZ( lTrkgPrimaryVtxPos );
413 const AliESDVertex *lPrimaryBestESDVtx = lESDevent->GetPrimaryVertex();
414 // get the best primary vertex available for the event
415 // As done in AliCascadeVertexer, we keep the one which is the best one available.
416 // between : Tracking vertex > SPD vertex > TPC vertex > default SPD vertex
417 // This one will be used for next calculations (DCA essentially)
418 lPrimaryBestESDVtx->GetXYZ( lBestPrimaryVtxPos );
420 Double_t tPrimaryVtxPosition[3];
421 const AliVVertex *primaryVtx = lESDevent->GetPrimaryVertex();
422 tPrimaryVtxPosition[0] = primaryVtx->GetX();
423 tPrimaryVtxPosition[1] = primaryVtx->GetY();
424 tPrimaryVtxPosition[2] = primaryVtx->GetZ();
426 //------------------------------------------------
427 // Primary Vertex Requirements Section:
428 // ---> pp and PbPb: Only requires |z|<10cm
429 // ---> pPb: all requirements checked at this stage
430 //------------------------------------------------
432 //Roberto's PV selection criteria, implemented 17th April 2013
434 /* vertex selection */
435 Bool_t fHasVertex = kFALSE;
437 const AliESDVertex *vertex = lESDevent->GetPrimaryVertexTracks();
438 if (vertex->GetNContributors() < 1) {
439 vertex = lESDevent->GetPrimaryVertexSPD();
440 if (vertex->GetNContributors() < 1) fHasVertex = kFALSE;
441 else fHasVertex = kTRUE;
442 TString vtxTyp = vertex->GetTitle();
444 vertex->GetCovarianceMatrix(cov);
445 Double_t zRes = TMath::Sqrt(cov[5]);
446 if (vtxTyp.Contains("vertexer:Z") && (zRes>0.25)) fHasVertex = kFALSE;
448 else fHasVertex = kTRUE;
450 //Is First event in chunk rejection: Still present!
451 if(fHasVertex == kFALSE) {
452 AliWarning("Pb / | PV does not satisfy selection criteria!");
453 PostData(1, fListHistV0);
455 PostData(3, fTreeEvents);
459 fHistCentralityHasVtx -> Fill ( fTreeVariableCentrality );
461 //17 April Fix: Always do primary vertex Z selection, after pA vertex selection from Roberto
462 if(TMath::Abs(lBestPrimaryVtxPos[2]) > 10.0) {
463 AliWarning("Pb / | pPb case | Z position of Best Prim Vtx | > 10.0 cm ... return !");
464 PostData(1, fListHistV0);
466 PostData(3, fTreeEvents);
470 fHistCentralityVtxZ -> Fill ( fTreeVariableCentrality );
471 lMagneticField = lESDevent->GetMagneticField( );
473 //Fill Event Tree: Analysis Selection Level
476 //------------------------------------------------
477 // MAIN LAMBDA LOOP STARTS HERE
478 //------------------------------------------------
480 //Variable definition
481 Int_t lOnFlyStatus = 0;// nv0sOn = 0, nv0sOff = 0;
482 Double_t lChi2V0 = 0;
483 Double_t lDcaV0Daughters = 0, lDcaV0ToPrimVertex = 0;
484 Double_t lDcaPosToPrimVertex = 0, lDcaNegToPrimVertex = 0;
485 Double_t lV0CosineOfPointingAngle = 0;
486 Double_t lV0Radius = 0, lPt = 0;
487 Double_t lRapK0Short = 0, lRapLambda = 0;
488 Double_t lInvMassK0s = 0, lInvMassLambda = 0, lInvMassAntiLambda = 0;
489 Double_t lAlphaV0 = 0, lPtArmV0 = 0;
491 Double_t fMinV0Pt = 0;
492 Double_t fMaxV0Pt = 100;
495 nv0s = lESDevent->GetNumberOfV0s();
497 //for (Int_t iV0 = 0; iV0 < nv0s; iV0++)
498 for (Int_t iV0 = 0; iV0 < nv0s; iV0++) //extra-crazy test
499 {// This is the begining of the V0 loop
500 AliESDv0 *v0 = ((AliESDEvent*)lESDevent)->GetV0(iV0);
503 Double_t tDecayVertexV0[3]; v0->GetXYZ(tDecayVertexV0[0],tDecayVertexV0[1],tDecayVertexV0[2]);
506 v0->GetPxPyPz( tV0mom[0],tV0mom[1],tV0mom[2] );
507 Double_t lV0TotalMomentum = TMath::Sqrt(
508 tV0mom[0]*tV0mom[0]+tV0mom[1]*tV0mom[1]+tV0mom[2]*tV0mom[2] );
510 lV0Radius = TMath::Sqrt(tDecayVertexV0[0]*tDecayVertexV0[0]+tDecayVertexV0[1]*tDecayVertexV0[1]);
513 lRapK0Short = v0->RapK0Short();
514 lRapLambda = v0->RapLambda();
515 if ((lPt<fMinV0Pt)||(fMaxV0Pt<lPt)) continue;
517 UInt_t lKeyPos = (UInt_t)TMath::Abs(v0->GetPindex());
518 UInt_t lKeyNeg = (UInt_t)TMath::Abs(v0->GetNindex());
520 Double_t lMomPos[3]; v0->GetPPxPyPz(lMomPos[0],lMomPos[1],lMomPos[2]);
521 Double_t lMomNeg[3]; v0->GetNPxPyPz(lMomNeg[0],lMomNeg[1],lMomNeg[2]);
523 AliESDtrack *pTrack=((AliESDEvent*)lESDevent)->GetTrack(lKeyPos);
524 AliESDtrack *nTrack=((AliESDEvent*)lESDevent)->GetTrack(lKeyNeg);
525 if (!pTrack || !nTrack) {
526 Printf("ERROR: Could not retreive one of the daughter track");
530 //Daughter Eta for Eta selection, afterwards
531 fTreeVariableNegEta = nTrack->Eta();
532 fTreeVariablePosEta = pTrack->Eta();
534 // Filter like-sign V0 (next: add counter and distribution)
535 if ( pTrack->GetSign() == nTrack->GetSign()){
539 //________________________________________________________________________
540 // Track quality cuts
541 Float_t lPosTrackCrossedRows = pTrack->GetTPCClusterInfo(2,1);
542 Float_t lNegTrackCrossedRows = nTrack->GetTPCClusterInfo(2,1);
543 fTreeVariableLeastNbrCrossedRows = (Int_t) lPosTrackCrossedRows;
544 if( lNegTrackCrossedRows < fTreeVariableLeastNbrCrossedRows )
545 fTreeVariableLeastNbrCrossedRows = (Int_t) lNegTrackCrossedRows;
547 // TPC refit condition (done during reconstruction for Offline but not for On-the-fly)
548 if( !(pTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;
549 if( !(nTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;
552 //fTreeVariablePosTrackStatus = pTrack->GetStatus();
553 //fTreeVariableNegTrackStatus = nTrack->GetStatus();
555 if ( ( ( pTrack->GetTPCClusterInfo(2,1) ) < 70 ) || ( ( nTrack->GetTPCClusterInfo(2,1) ) < 70 ) ) continue;
557 //GetKinkIndex condition
558 if( pTrack->GetKinkIndex(0)>0 || nTrack->GetKinkIndex(0)>0 ) continue;
560 //Findable clusters > 0 condition
561 if( pTrack->GetTPCNclsF()<=0 || nTrack->GetTPCNclsF()<=0 ) continue;
563 //Compute ratio Crossed Rows / Findable clusters
564 //Note: above test avoids division by zero!
565 Float_t lPosTrackCrossedRowsOverFindable = lPosTrackCrossedRows / ((double)(pTrack->GetTPCNclsF()));
566 Float_t lNegTrackCrossedRowsOverFindable = lNegTrackCrossedRows / ((double)(nTrack->GetTPCNclsF()));
568 fTreeVariableLeastRatioCrossedRowsOverFindable = lPosTrackCrossedRowsOverFindable;
569 if( lNegTrackCrossedRowsOverFindable < fTreeVariableLeastRatioCrossedRowsOverFindable )
570 fTreeVariableLeastRatioCrossedRowsOverFindable = lNegTrackCrossedRowsOverFindable;
572 //Lowest Cut Level for Ratio Crossed Rows / Findable = 0.8, set here
573 //if ( fTreeVariableLeastRatioCrossedRowsOverFindable < 0.8 ) continue;
575 //End track Quality Cuts
576 //________________________________________________________________________
578 lDcaPosToPrimVertex = TMath::Abs(pTrack->GetD(tPrimaryVtxPosition[0],
579 tPrimaryVtxPosition[1],
582 lDcaNegToPrimVertex = TMath::Abs(nTrack->GetD(tPrimaryVtxPosition[0],
583 tPrimaryVtxPosition[1],
586 lOnFlyStatus = v0->GetOnFlyStatus();
587 lChi2V0 = v0->GetChi2V0();
588 lDcaV0Daughters = v0->GetDcaV0Daughters();
589 lDcaV0ToPrimVertex = v0->GetD(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],tPrimaryVtxPosition[2]);
590 lV0CosineOfPointingAngle = v0->GetV0CosineOfPointingAngle(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],tPrimaryVtxPosition[2]);
591 fTreeVariableV0CosineOfPointingAngle=lV0CosineOfPointingAngle;
593 // Getting invariant mass infos directly from ESD
594 v0->ChangeMassHypothesis(310);
595 lInvMassK0s = v0->GetEffMass();
596 v0->ChangeMassHypothesis(3122);
597 lInvMassLambda = v0->GetEffMass();
598 v0->ChangeMassHypothesis(-3122);
599 lInvMassAntiLambda = v0->GetEffMass();
600 lAlphaV0 = v0->AlphaV0();
601 lPtArmV0 = v0->PtArmV0();
603 fTreeVariablePt = v0->Pt();
604 fTreeVariableChi2V0 = lChi2V0;
605 fTreeVariableDcaV0ToPrimVertex = lDcaV0ToPrimVertex;
606 fTreeVariableDcaV0Daughters = lDcaV0Daughters;
607 fTreeVariableV0CosineOfPointingAngle = lV0CosineOfPointingAngle;
608 fTreeVariableV0Radius = lV0Radius;
609 fTreeVariableDcaPosToPrimVertex = lDcaPosToPrimVertex;
610 fTreeVariableDcaNegToPrimVertex = lDcaNegToPrimVertex;
611 fTreeVariableInvMassK0s = lInvMassK0s;
612 fTreeVariableInvMassLambda = lInvMassLambda;
613 fTreeVariableInvMassAntiLambda = lInvMassAntiLambda;
614 fTreeVariableRapK0Short = lRapK0Short;
615 fTreeVariableRapLambda = lRapLambda;
616 fTreeVariableAlphaV0 = lAlphaV0;
617 fTreeVariablePtArmV0 = lPtArmV0;
619 //Official means of acquiring N-sigmas
620 fTreeVariableNSigmasPosProton = fPIDResponse->NumberOfSigmasTPC( pTrack, AliPID::kProton );
621 fTreeVariableNSigmasPosPion = fPIDResponse->NumberOfSigmasTPC( pTrack, AliPID::kPion );
622 fTreeVariableNSigmasNegProton = fPIDResponse->NumberOfSigmasTPC( nTrack, AliPID::kProton );
623 fTreeVariableNSigmasNegPion = fPIDResponse->NumberOfSigmasTPC( nTrack, AliPID::kPion );
625 //This requires an Invariant Mass Hypothesis afterwards
626 fTreeVariableDistOverTotMom = TMath::Sqrt(
627 TMath::Power( tDecayVertexV0[0] - lBestPrimaryVtxPos[0] , 2) +
628 TMath::Power( tDecayVertexV0[1] - lBestPrimaryVtxPos[1] , 2) +
629 TMath::Power( tDecayVertexV0[2] - lBestPrimaryVtxPos[2] , 2)
631 fTreeVariableDistOverTotMom /= (lV0TotalMomentum+1e-10); //avoid division by zero, to be sure
633 //------------------------------------------------
635 //------------------------------------------------
637 // The conditionals are meant to decrease excessive
640 //First Selection: Reject OnFly
641 if( lOnFlyStatus == 0 ){
642 //Second Selection: rough 20-sigma band, parametric.
643 //K0Short: Enough to parametrize peak broadening with linear function.
644 Double_t lUpperLimitK0Short = (5.63707e-01) + (1.14979e-02)*fTreeVariablePt;
645 Double_t lLowerLimitK0Short = (4.30006e-01) - (1.10029e-02)*fTreeVariablePt;
646 //Lambda: Linear (for higher pt) plus exponential (for low-pt broadening)
647 //[0]+[1]*x+[2]*TMath::Exp(-[3]*x)
648 Double_t lUpperLimitLambda = (1.13688e+00) + (5.27838e-03)*fTreeVariablePt + (8.42220e-02)*TMath::Exp(-(3.80595e+00)*fTreeVariablePt);
649 Double_t lLowerLimitLambda = (1.09501e+00) - (5.23272e-03)*fTreeVariablePt - (7.52690e-02)*TMath::Exp(-(3.46339e+00)*fTreeVariablePt);
651 if( (fTreeVariableInvMassLambda < lUpperLimitLambda && fTreeVariableInvMassLambda > lLowerLimitLambda ) ||
652 (fTreeVariableInvMassAntiLambda < lUpperLimitLambda && fTreeVariableInvMassAntiLambda > lLowerLimitLambda ) ||
653 (fTreeVariableInvMassK0s < lUpperLimitK0Short && fTreeVariableInvMassK0s > lLowerLimitK0Short ) ){
654 //Pre-selection in case this is AA...
655 //If this is a nuclear collision___________________
656 // ... pre-filter with TPC, daughter eta selection
657 if( (fTreeVariableInvMassLambda < lUpperLimitLambda && fTreeVariableInvMassLambda > lLowerLimitLambda
658 && ( fTPCdEdxSelection==kFALSE ||( fTPCdEdxSelection==kTRUE && TMath::Abs(fTreeVariableNSigmasPosProton) < 6.0 && TMath::Abs(fTreeVariableNSigmasNegPion) < 6.0 ) ) ) ||
659 (fTreeVariableInvMassAntiLambda < lUpperLimitLambda && fTreeVariableInvMassAntiLambda > lLowerLimitLambda
660 && ( fTPCdEdxSelection==kFALSE ||( fTPCdEdxSelection==kTRUE && TMath::Abs(fTreeVariableNSigmasNegProton) < 6.0 && TMath::Abs(fTreeVariableNSigmasPosPion) < 6.0 ) ) ) ||
661 (fTreeVariableInvMassK0s < lUpperLimitK0Short && fTreeVariableInvMassK0s > lLowerLimitK0Short
662 && ( fTPCdEdxSelection==kFALSE ||( fTPCdEdxSelection==kTRUE && TMath::Abs(fTreeVariableNSigmasNegPion) < 6.0 && TMath::Abs(fTreeVariableNSigmasPosPion) < 6.0 ) ) ) ){
664 if ( TMath::Abs(fTreeVariableNegEta)<0.8 && TMath::Abs(fTreeVariablePosEta)<0.8) fTree->Fill();
666 }//end nuclear_____________________________________
669 //------------------------------------------------
671 //------------------------------------------------
673 }// This is the end of the V0 loop
676 PostData(1, fListHistV0);
678 PostData(3, fTreeEvents);
681 //________________________________________________________________________
682 void AliAnalysisTaskExtractV0pPb::Terminate(Option_t *)
684 // Draw result to the screen
685 // Called once at the end of the query
686 // This will draw the V0 candidate multiplicity, whose
687 // number of entries corresponds to the number of triggered events.
688 TList *cRetrievedList = 0x0;
689 cRetrievedList = (TList*)GetOutputData(1);
691 Printf("ERROR - AliAnalysisTaskExtractV0pPb : ouput data container list not available\n");
694 fHistCentralityProcessed = dynamic_cast<TH1F*> ( cRetrievedList->FindObject("fHistCentralityProcessed") );
695 if (!fHistCentralityProcessed) {
696 Printf("ERROR - AliAnalysisTaskExtractV0pPb : fHistV0MultiplicityForTrigEvt not available");
699 TCanvas *canCheck = new TCanvas("AliAnalysisTaskExtractV0pPb","V0 Multiplicity",10,10,510,510);
700 canCheck->cd(1)->SetLogy();
701 fHistCentralityProcessed->SetMarkerStyle(22);
702 fHistCentralityProcessed->DrawCopy("E");