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 fkSwitchINT7 ( kFALSE ),
96 fkUseOnTheFly ( kFALSE ),
97 fkTakeAllTracks ( kFALSE ),
98 fCentralityEstimator("V0M"),
99 fkLightWeight ( kFALSE ),
100 //------------------------------------------------
102 fTreeVariableChi2V0(0),
103 fTreeVariableDcaV0Daughters(0),
104 fTreeVariableDcaV0ToPrimVertex(0),
105 fTreeVariableDcaPosToPrimVertex(0),
106 fTreeVariableDcaNegToPrimVertex(0),
107 fTreeVariableV0CosineOfPointingAngle(0),
108 fTreeVariableV0Radius(0),
110 fTreeVariableRapK0Short(0),
111 fTreeVariableRapLambda(0),
112 fTreeVariableInvMassK0s(0),
113 fTreeVariableInvMassLambda(0),
114 fTreeVariableInvMassAntiLambda(0),
115 fTreeVariableAlphaV0(0),
116 fTreeVariablePtArmV0(0),
117 fTreeVariableNegTotMomentum(0),
118 fTreeVariablePosTotMomentum(0),
119 fTreeVariableNegdEdxSig(0),
120 fTreeVariablePosdEdxSig(0),
121 fTreeVariableNegEta(0),
122 fTreeVariablePosEta(0),
124 fTreeVariableNSigmasPosProton(0),
125 fTreeVariableNSigmasPosPion(0),
126 fTreeVariableNSigmasNegProton(0),
127 fTreeVariableNSigmasNegPion(0),
129 fTreeVariableDistOverTotMom(0),
130 fTreeVariableLeastNbrCrossedRows(0),
131 fTreeVariableLeastRatioCrossedRowsOverFindable(0),
132 fTreeVariableMultiplicity(0),
134 fTreeVariableRunNumber(0),
135 fTreeVariableEventNumber(0),
141 fTreeVariableV0Px(0),
142 fTreeVariableV0Py(0),
143 fTreeVariableV0Pz(0),
149 fTreeVariableNegTrackStatus(0),
150 fTreeVariablePosTrackStatus(0),
152 //------------------------------------------------
154 // --- Filled on an Event-by-event basis
155 //------------------------------------------------
156 fHistV0MultiplicityBeforeTrigSel(0),
157 fHistV0MultiplicityForTrigEvt(0),
158 fHistV0MultiplicityForSelEvt(0),
159 fHistV0MultiplicityForSelEvtNoTPCOnly(0),
160 fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup(0),
161 fHistMultiplicityBeforeTrigSel(0),
162 fHistMultiplicityForTrigEvt(0),
163 fHistMultiplicity(0),
164 fHistMultiplicityNoTPCOnly(0),
165 fHistMultiplicityNoTPCOnlyNoPileup(0),
167 //Raw Data for Vertex Z position estimator change
168 f2dHistMultiplicityVsVertexZBeforeTrigSel(0),
169 f2dHistMultiplicityVsVertexZForTrigEvt(0),
170 f2dHistMultiplicityVsVertexZ(0),
171 f2dHistMultiplicityVsVertexZNoTPCOnly(0),
172 f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup(0),
180 fHistSwappedV0Counter(0)
185 AliAnalysisTaskExtractV0::AliAnalysisTaskExtractV0(const char *name)
186 : AliAnalysisTaskSE(name), fListHistV0(0), fTree(0), fPIDResponse(0),fESDtrackCuts(0),
187 fkIsNuclear ( kFALSE ),
188 fkSwitchINT7 ( kFALSE ),
189 fkUseOnTheFly ( kFALSE ),
190 fkTakeAllTracks ( kFALSE ),
191 fCentralityEstimator("V0M"),
192 fkLightWeight ( kFALSE ),
193 //------------------------------------------------
195 fTreeVariableChi2V0(0),
196 fTreeVariableDcaV0Daughters(0),
197 fTreeVariableDcaV0ToPrimVertex(0),
198 fTreeVariableDcaPosToPrimVertex(0),
199 fTreeVariableDcaNegToPrimVertex(0),
200 fTreeVariableV0CosineOfPointingAngle(0),
201 fTreeVariableV0Radius(0),
203 fTreeVariableRapK0Short(0),
204 fTreeVariableRapLambda(0),
205 fTreeVariableInvMassK0s(0),
206 fTreeVariableInvMassLambda(0),
207 fTreeVariableInvMassAntiLambda(0),
208 fTreeVariableAlphaV0(0),
209 fTreeVariablePtArmV0(0),
210 fTreeVariableNegTotMomentum(0),
211 fTreeVariablePosTotMomentum(0),
212 fTreeVariableNegdEdxSig(0),
213 fTreeVariablePosdEdxSig(0),
214 fTreeVariableNegEta(0),
215 fTreeVariablePosEta(0),
217 fTreeVariableNSigmasPosProton(0),
218 fTreeVariableNSigmasPosPion(0),
219 fTreeVariableNSigmasNegProton(0),
220 fTreeVariableNSigmasNegPion(0),
222 fTreeVariableDistOverTotMom(0),
223 fTreeVariableLeastNbrCrossedRows(0),
224 fTreeVariableLeastRatioCrossedRowsOverFindable(0),
225 fTreeVariableMultiplicity(0),
227 fTreeVariableRunNumber(0),
228 fTreeVariableEventNumber(0),
234 fTreeVariableV0Px(0),
235 fTreeVariableV0Py(0),
236 fTreeVariableV0Pz(0),
242 fTreeVariableNegTrackStatus(0),
243 fTreeVariablePosTrackStatus(0),
245 //------------------------------------------------
247 // --- Filled on an Event-by-event basis
248 //------------------------------------------------
249 fHistV0MultiplicityBeforeTrigSel(0),
250 fHistV0MultiplicityForTrigEvt(0),
251 fHistV0MultiplicityForSelEvt(0),
252 fHistV0MultiplicityForSelEvtNoTPCOnly(0),
253 fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup(0),
254 fHistMultiplicityBeforeTrigSel(0),
255 fHistMultiplicityForTrigEvt(0),
256 fHistMultiplicity(0),
257 fHistMultiplicityNoTPCOnly(0),
258 fHistMultiplicityNoTPCOnlyNoPileup(0),
260 //Raw Data for Vertex Z position estimator change
261 f2dHistMultiplicityVsVertexZBeforeTrigSel(0),
262 f2dHistMultiplicityVsVertexZForTrigEvt(0),
263 f2dHistMultiplicityVsVertexZ(0),
264 f2dHistMultiplicityVsVertexZNoTPCOnly(0),
265 f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup(0),
273 fHistSwappedV0Counter(0)
276 // Output slot #0 writes into a TList container (Lambda Histos and fTree)
277 DefineOutput(1, TList::Class());
278 DefineOutput(2, TTree::Class());
282 AliAnalysisTaskExtractV0::~AliAnalysisTaskExtractV0()
284 //------------------------------------------------
286 //------------------------------------------------
296 //cleanup esd track cuts object too...
298 delete fESDtrackCuts;
307 //________________________________________________________________________
308 void AliAnalysisTaskExtractV0::UserCreateOutputObjects()
311 //Create File-resident Tree, please.
314 fTree = new TTree("fTree","V0Candidates");
316 //------------------------------------------------
317 // fTree Branch definitions
318 //------------------------------------------------
320 //-----------BASIC-INFO---------------------------
321 /* 1*/ fTree->Branch("fTreeVariableChi2V0",&fTreeVariableChi2V0,"fTreeVariableChi2V0/F");
322 /* 2*/ fTree->Branch("fTreeVariableDcaV0Daughters",&fTreeVariableDcaV0Daughters,"fTreeVariableDcaV0Daughters/F");
323 /* 3*/ fTree->Branch("fTreeVariableDcaPosToPrimVertex",&fTreeVariableDcaPosToPrimVertex,"fTreeVariableDcaPosToPrimVertex/F");
324 /* 4*/ fTree->Branch("fTreeVariableDcaNegToPrimVertex",&fTreeVariableDcaNegToPrimVertex,"fTreeVariableDcaNegToPrimVertex/F");
325 /* 5*/ fTree->Branch("fTreeVariableV0Radius",&fTreeVariableV0Radius,"fTreeVariableV0Radius/F");
326 /* 6*/ fTree->Branch("fTreeVariablePt",&fTreeVariablePt,"fTreeVariablePt/F");
327 /* 7*/ fTree->Branch("fTreeVariableRapK0Short",&fTreeVariableRapK0Short,"fTreeVariableRapK0Short/F");
328 /* 8*/ fTree->Branch("fTreeVariableRapLambda",&fTreeVariableRapLambda,"fTreeVariableRapLambda/F");
329 /* 9*/ fTree->Branch("fTreeVariableInvMassK0s",&fTreeVariableInvMassK0s,"fTreeVariableInvMassK0s/F");
330 /*10*/ fTree->Branch("fTreeVariableInvMassLambda",&fTreeVariableInvMassLambda,"fTreeVariableInvMassLambda/F");
331 /*11*/ fTree->Branch("fTreeVariableInvMassAntiLambda",&fTreeVariableInvMassAntiLambda,"fTreeVariableInvMassAntiLambda/F");
332 /*12*/ fTree->Branch("fTreeVariableV0CosineOfPointingAngle",&fTreeVariableV0CosineOfPointingAngle,"fTreeVariableV0CosineOfPointingAngle/F");
333 /*13*/ fTree->Branch("fTreeVariableAlphaV0",&fTreeVariableAlphaV0,"fTreeVariableAlphaV0/F");
334 /*14*/ fTree->Branch("fTreeVariablePtArmV0",&fTreeVariablePtArmV0,"fTreeVariablePtArmV0/F");
335 /*15*/ fTree->Branch("fTreeVariableLeastNbrCrossedRows",&fTreeVariableLeastNbrCrossedRows,"fTreeVariableLeastNbrCrossedRows/I");
336 /*16*/ fTree->Branch("fTreeVariableLeastRatioCrossedRowsOverFindable",&fTreeVariableLeastRatioCrossedRowsOverFindable,"fTreeVariableLeastRatioCrossedRowsOverFindable/F");
337 //-----------MULTIPLICITY-INFO--------------------
338 /*17*/ fTree->Branch("fTreeVariableMultiplicity",&fTreeVariableMultiplicity,"fTreeVariableMultiplicity/I");
339 //------------------------------------------------
340 /*18*/ fTree->Branch("fTreeVariableDistOverTotMom",&fTreeVariableDistOverTotMom,"fTreeVariableDistOverTotMom/F");
341 /*19*/ fTree->Branch("fTreeVariableNSigmasPosProton",&fTreeVariableNSigmasPosProton,"fTreeVariableNSigmasPosProton/F");
342 /*20*/ fTree->Branch("fTreeVariableNSigmasPosPion",&fTreeVariableNSigmasPosPion,"fTreeVariableNSigmasPosPion/F");
343 /*21*/ fTree->Branch("fTreeVariableNSigmasNegProton",&fTreeVariableNSigmasNegProton,"fTreeVariableNSigmasNegProton/F");
344 /*22*/ fTree->Branch("fTreeVariableNSigmasNegPion",&fTreeVariableNSigmasNegPion,"fTreeVariableNSigmasNegPion/F");
345 /*23*/ fTree->Branch("fTreeVariableNegEta",&fTreeVariableNegEta,"fTreeVariableNegEta/F");
346 /*24*/ fTree->Branch("fTreeVariablePosEta",&fTreeVariablePosEta,"fTreeVariablePosEta/F");
347 /*25*/ fTree->Branch("fTreeVariableRunNumber",&fTreeVariableRunNumber,"fTreeVariableRunNumber/I");
348 /*26*/ fTree->Branch("fTreeVariableEventNumber",&fTreeVariableEventNumber,"fTreeVariableEventNumber/l");
350 if( fkLightWeight == kFALSE ){
351 //-----------FOR CTAU DEBUGGING: Full Phase Space + Decay Position Info
352 fTree->Branch("fTreeVariablePVx",&fTreeVariablePVx,"fTreeVariablePVx/F");
353 fTree->Branch("fTreeVariablePVy",&fTreeVariablePVy,"fTreeVariablePVy/F");
354 fTree->Branch("fTreeVariablePVz",&fTreeVariablePVz,"fTreeVariablePVz/F");
356 fTree->Branch("fTreeVariableV0x",&fTreeVariableV0x,"fTreeVariableV0x/F");
357 fTree->Branch("fTreeVariableV0y",&fTreeVariableV0y,"fTreeVariableV0y/F");
358 fTree->Branch("fTreeVariableV0z",&fTreeVariableV0z,"fTreeVariableV0z/F");
360 fTree->Branch("fTreeVariableV0Px",&fTreeVariableV0Px,"fTreeVariableV0Px/F");
361 fTree->Branch("fTreeVariableV0Py",&fTreeVariableV0Py,"fTreeVariableV0Py/F");
362 fTree->Branch("fTreeVariableV0Pz",&fTreeVariableV0Pz,"fTreeVariableV0Pz/F");
364 fTree->Branch("fTreeVariableNegTrackStatus",&fTreeVariableNegTrackStatus,"fTreeVariableNegTrackStatus/l");
365 fTree->Branch("fTreeVariablePosTrackStatus",&fTreeVariablePosTrackStatus,"fTreeVariablePosTrackStatus/l");
368 //------------------------------------------------
369 // Particle Identification Setup
370 //------------------------------------------------
372 AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager();
373 AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler());
374 fPIDResponse = inputHandler->GetPIDResponse();
378 if(! fESDtrackCuts ){
379 fESDtrackCuts = new AliESDtrackCuts();
383 //------------------------------------------------
384 // V0 Multiplicity Histograms
385 //------------------------------------------------
388 //Create File-resident Tree, please.
391 fListHistV0 = new TList();
392 fListHistV0->SetOwner(); // See http://root.cern.ch/root/html/TCollection.html#TCollection:SetOwner
394 if(! fHistV0MultiplicityBeforeTrigSel) {
395 fHistV0MultiplicityBeforeTrigSel = new TH1F("fHistV0MultiplicityBeforeTrigSel",
396 "V0s per event (before Trig. Sel.);Nbr of V0s/Evt;Events",
398 fListHistV0->Add(fHistV0MultiplicityBeforeTrigSel);
401 if(! fHistV0MultiplicityForTrigEvt) {
402 fHistV0MultiplicityForTrigEvt = new TH1F("fHistV0MultiplicityForTrigEvt",
403 "V0s per event (for triggered evt);Nbr of V0s/Evt;Events",
405 fListHistV0->Add(fHistV0MultiplicityForTrigEvt);
408 if(! fHistV0MultiplicityForSelEvt) {
409 fHistV0MultiplicityForSelEvt = new TH1F("fHistV0MultiplicityForSelEvt",
410 "V0s per event;Nbr of V0s/Evt;Events",
412 fListHistV0->Add(fHistV0MultiplicityForSelEvt);
415 if(! fHistV0MultiplicityForSelEvtNoTPCOnly) {
416 fHistV0MultiplicityForSelEvtNoTPCOnly = new TH1F("fHistV0MultiplicityForSelEvtNoTPCOnly",
417 "V0s per event;Nbr of V0s/Evt;Events",
419 fListHistV0->Add(fHistV0MultiplicityForSelEvtNoTPCOnly);
422 if(! fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup) {
423 fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup = new TH1F("fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup",
424 "V0s per event;Nbr of V0s/Evt;Events",
426 fListHistV0->Add(fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup);
429 //------------------------------------------------
430 // Track Multiplicity Histograms
431 //------------------------------------------------
433 if(! fHistMultiplicityBeforeTrigSel) {
434 fHistMultiplicityBeforeTrigSel = new TH1F("fHistMultiplicityBeforeTrigSel",
435 "Tracks per event;Nbr of Tracks;Events",
437 fListHistV0->Add(fHistMultiplicityBeforeTrigSel);
439 if(! fHistMultiplicityForTrigEvt) {
440 fHistMultiplicityForTrigEvt = new TH1F("fHistMultiplicityForTrigEvt",
441 "Tracks per event;Nbr of Tracks;Events",
443 fListHistV0->Add(fHistMultiplicityForTrigEvt);
445 if(! fHistMultiplicity) {
446 fHistMultiplicity = new TH1F("fHistMultiplicity",
447 "Tracks per event;Nbr of Tracks;Events",
449 fListHistV0->Add(fHistMultiplicity);
451 if(! fHistMultiplicityNoTPCOnly) {
452 fHistMultiplicityNoTPCOnly = new TH1F("fHistMultiplicityNoTPCOnly",
453 "Tracks per event;Nbr of Tracks;Events",
455 fListHistV0->Add(fHistMultiplicityNoTPCOnly);
457 if(! fHistMultiplicityNoTPCOnlyNoPileup) {
458 fHistMultiplicityNoTPCOnlyNoPileup = new TH1F("fHistMultiplicityNoTPCOnlyNoPileup",
459 "Tracks per event;Nbr of Tracks;Events",
461 fListHistV0->Add(fHistMultiplicityNoTPCOnlyNoPileup);
466 //Raw Data for Vertex Z position estimator change
467 //TH2F *f2dHistMultiplicityVsVertexZBeforeTrigSel; //! multiplicity distribution
468 //TH2F *f2dHistMultiplicityVsVertexZForTrigEvt; //! multiplicity distribution
469 //TH2F *f2dHistMultiplicityVsVertexZ; //! multiplicity distribution
470 //TH2F *f2dHistMultiplicityVsVertexZNoTPCOnly; //! multiplicity distribution
471 //TH2F *f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup; //! multiplicity distribution
473 if(! f2dHistMultiplicityVsVertexZBeforeTrigSel) {
474 f2dHistMultiplicityVsVertexZBeforeTrigSel = new TH2F("f2dHistMultiplicityVsVertexZBeforeTrigSel",
475 "Tracks per event", 200, 0, 200,400, -20, 20);
476 fListHistV0->Add(f2dHistMultiplicityVsVertexZBeforeTrigSel);
478 if(! f2dHistMultiplicityVsVertexZForTrigEvt) {
479 f2dHistMultiplicityVsVertexZForTrigEvt = new TH2F("f2dHistMultiplicityVsVertexZForTrigEvt",
480 "Tracks per event", 200, 0, 200, 400, -20, 20);
481 fListHistV0->Add(f2dHistMultiplicityVsVertexZForTrigEvt);
483 if(! f2dHistMultiplicityVsVertexZ) {
484 f2dHistMultiplicityVsVertexZ = new TH2F("f2dHistMultiplicityVsVertexZ",
485 "Tracks per event", 200, 0, 200, 400, -20, 20);
486 fListHistV0->Add(f2dHistMultiplicityVsVertexZ);
488 if(! f2dHistMultiplicityVsVertexZNoTPCOnly) {
489 f2dHistMultiplicityVsVertexZNoTPCOnly = new TH2F("f2dHistMultiplicityVsVertexZNoTPCOnly",
490 "Tracks per event", 200, 0, 200, 400, -20, 20);
491 fListHistV0->Add(f2dHistMultiplicityVsVertexZNoTPCOnly);
493 if(! f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup) {
494 f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup = new TH2F("f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup",
495 "Tracks per event", 200, 0, 200, 400, -20, 20);
496 fListHistV0->Add(f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup);
499 //-----------------------------------------------------
502 fHistPVx = new TH1F("fHistPVx",
503 "PV x position;Nbr of Evts;x",
505 fListHistV0->Add(fHistPVx);
508 fHistPVy = new TH1F("fHistPVy",
509 "PV y position;Nbr of Evts;y",
511 fListHistV0->Add(fHistPVy);
514 fHistPVz = new TH1F("fHistPVz",
515 "PV z position;Nbr of Evts;z",
517 fListHistV0->Add(fHistPVz);
519 if(! fHistPVxAnalysis) {
520 fHistPVxAnalysis = new TH1F("fHistPVxAnalysis",
521 "PV x position;Nbr of Evts;x",
523 fListHistV0->Add(fHistPVxAnalysis);
525 if(! fHistPVyAnalysis) {
526 fHistPVyAnalysis = new TH1F("fHistPVyAnalysis",
527 "PV y position;Nbr of Evts;y",
529 fListHistV0->Add(fHistPVyAnalysis);
531 if(! fHistPVzAnalysis) {
532 fHistPVzAnalysis = new TH1F("fHistPVzAnalysis",
533 "PV z position;Nbr of Evts;z",
535 fListHistV0->Add(fHistPVzAnalysis);
537 if(! fHistSwappedV0Counter) {
538 fHistSwappedV0Counter = new TH1F("fHistSwappedV0Counter",
539 "Swap or not histo;Swapped (1) or not (0); count",
541 fListHistV0->Add(fHistSwappedV0Counter);
543 //Regular output: Histograms
544 PostData(1, fListHistV0);
545 //TTree Object: Saved to base directory. Should cache to disk while saving.
546 //(Important to avoid excessive memory usage, particularly when merging)
549 }// end UserCreateOutputObjects
552 //________________________________________________________________________
553 void AliAnalysisTaskExtractV0::UserExec(Option_t *)
556 // Called for each event
557 //gObjectTable->Print();
558 AliESDEvent *lESDevent = 0x0;
560 //AliAODEvent *lAODevent = 0x0;
563 Double_t lTrkgPrimaryVtxPos[3] = {-100.0, -100.0, -100.0};
564 Double_t lBestPrimaryVtxPos[3] = {-100.0, -100.0, -100.0};
565 Double_t lMagneticField = -10.;
567 // Connect to the InputEvent
568 // After these lines, we should have an ESD/AOD event + the number of cascades in it.
570 // Appropriate for ESD analysis!
572 lESDevent = dynamic_cast<AliESDEvent*>( InputEvent() );
574 AliWarning("ERROR: lESDevent not available \n");
577 fTreeVariableRunNumber = lESDevent->GetRunNumber();
578 fTreeVariableEventNumber =
579 ( ( ((ULong64_t)lESDevent->GetPeriodNumber() ) << 36 ) |
580 ( ((ULong64_t)lESDevent->GetOrbitNumber () ) << 12 ) |
581 ((ULong64_t)lESDevent->GetBunchCrossNumber() ) );
583 //REVISED multiplicity estimator after 'multiplicity day' (2011)
584 Int_t lMultiplicity = -100;
586 if(fkIsNuclear == kFALSE) lMultiplicity = fESDtrackCuts->GetReferenceMultiplicity(lESDevent, AliESDtrackCuts::kTrackletsITSTPC,0.5);
588 //---> If this is a nuclear collision, then go nuclear on "multiplicity" variable...
589 //---> Warning: Experimental
590 if(fkIsNuclear == kTRUE){
591 AliCentrality* centrality;
592 centrality = lESDevent->GetCentrality();
593 lMultiplicity = ( ( Int_t ) ( centrality->GetCentralityPercentile( fCentralityEstimator.Data() ) ) );
594 if (centrality->GetQuality()>1) {
595 PostData(1, fListHistV0);
601 //Set variable for filling tree afterwards!
602 //---> pp case......: GetReferenceMultiplicity
603 //---> Pb-Pb case...: Centrality by V0M
604 fTreeVariableMultiplicity = lMultiplicity;
605 fHistMultiplicityBeforeTrigSel->Fill ( lMultiplicity );
606 fHistV0MultiplicityBeforeTrigSel->Fill ( lESDevent->GetNumberOfV0s() );
608 //------------------------------------------------
610 //------------------------------------------------
613 UInt_t maskIsSelected = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
614 Bool_t isSelected = 0;
615 isSelected = (maskIsSelected & AliVEvent::kMB) == AliVEvent::kMB;
617 //pA triggering: CINT7
618 if( fkSwitchINT7 ) isSelected = (maskIsSelected & AliVEvent::kINT7) == AliVEvent::kINT7;
620 //Standard Min-Bias Selection
621 if ( ! isSelected ) {
622 PostData(1, fListHistV0);
627 //------------------------------------------------
628 // After Trigger Selection
629 //------------------------------------------------
631 nV0s = lESDevent->GetNumberOfV0s();
633 fHistV0MultiplicityForTrigEvt ->Fill( nV0s );
634 fHistMultiplicityForTrigEvt ->Fill( lMultiplicity );
636 //------------------------------------------------
637 // Getting: Primary Vertex + MagField Info
638 //------------------------------------------------
640 const AliESDVertex *lPrimaryTrackingESDVtx = lESDevent->GetPrimaryVertexTracks();
641 // get the vtx stored in ESD found with tracks
642 lPrimaryTrackingESDVtx->GetXYZ( lTrkgPrimaryVtxPos );
644 const AliESDVertex *lPrimaryBestESDVtx = lESDevent->GetPrimaryVertex();
645 // get the best primary vertex available for the event
646 // As done in AliCascadeVertexer, we keep the one which is the best one available.
647 // between : Tracking vertex > SPD vertex > TPC vertex > default SPD vertex
648 // This one will be used for next calculations (DCA essentially)
649 lPrimaryBestESDVtx->GetXYZ( lBestPrimaryVtxPos );
651 Double_t tPrimaryVtxPosition[3];
652 const AliVVertex *primaryVtx = lESDevent->GetPrimaryVertex();
653 tPrimaryVtxPosition[0] = primaryVtx->GetX();
654 tPrimaryVtxPosition[1] = primaryVtx->GetY();
655 tPrimaryVtxPosition[2] = primaryVtx->GetZ();
656 fHistPVx->Fill( tPrimaryVtxPosition[0] );
657 fHistPVy->Fill( tPrimaryVtxPosition[1] );
658 fHistPVz->Fill( tPrimaryVtxPosition[2] );
660 f2dHistMultiplicityVsVertexZForTrigEvt->Fill( lMultiplicity, tPrimaryVtxPosition[2] );
662 //------------------------------------------------
663 // Primary Vertex Z position: SKIP
664 //------------------------------------------------
666 if(TMath::Abs(lBestPrimaryVtxPos[2]) > 10.0 ) {
667 AliWarning("Pb / | Z position of Best Prim Vtx | > 10.0 cm ... return !");
668 PostData(1, fListHistV0);
673 f2dHistMultiplicityVsVertexZ->Fill( lMultiplicity, tPrimaryVtxPosition[2] );
675 lMagneticField = lESDevent->GetMagneticField( );
676 fHistV0MultiplicityForSelEvt ->Fill( nV0s );
677 fHistMultiplicity->Fill(lMultiplicity);
679 //------------------------------------------------
680 // Only look at events with well-established PV
681 //------------------------------------------------
683 const AliESDVertex *lPrimaryTrackingESDVtxCheck = lESDevent->GetPrimaryVertexTracks();
684 const AliESDVertex *lPrimarySPDVtx = lESDevent->GetPrimaryVertexSPD();
685 if (!lPrimarySPDVtx->GetStatus() && !lPrimaryTrackingESDVtxCheck->GetStatus() ){
686 AliWarning("Pb / No SPD prim. vertex nor prim. Tracking vertex ... return !");
687 PostData(1, fListHistV0);
693 f2dHistMultiplicityVsVertexZNoTPCOnly->Fill( lMultiplicity, tPrimaryVtxPosition[2] );
694 fHistV0MultiplicityForSelEvtNoTPCOnly ->Fill( nV0s );
695 fHistMultiplicityNoTPCOnly->Fill(lMultiplicity);
697 //------------------------------------------------
699 //------------------------------------------------
701 // FIXME : quality selection regarding pile-up rejection
702 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
703 PostData(1, fListHistV0);
708 f2dHistMultiplicityVsVertexZNoTPCOnlyNoPileup->Fill( lMultiplicity, tPrimaryVtxPosition[2] );
709 fHistV0MultiplicityForSelEvtNoTPCOnlyNoPileup ->Fill( nV0s );
710 fHistMultiplicityNoTPCOnlyNoPileup->Fill(lMultiplicity);
712 //------------------------------------------------
713 // MAIN LAMBDA LOOP STARTS HERE
714 //------------------------------------------------
716 if( ! (lESDevent->GetPrimaryVertex()->IsFromVertexerZ() ) ){
717 fHistPVxAnalysis->Fill( tPrimaryVtxPosition[0] );
718 fHistPVyAnalysis->Fill( tPrimaryVtxPosition[1] );
719 fHistPVzAnalysis->Fill( tPrimaryVtxPosition[2] );
722 fTreeVariablePVx = tPrimaryVtxPosition[0];
723 fTreeVariablePVy = tPrimaryVtxPosition[1];
724 fTreeVariablePVz = tPrimaryVtxPosition[2];
726 //Variable definition
727 Int_t lOnFlyStatus = 0;// nv0sOn = 0, nv0sOff = 0;
728 Double_t lChi2V0 = 0;
729 Double_t lDcaV0Daughters = 0, lDcaV0ToPrimVertex = 0;
730 Double_t lDcaPosToPrimVertex = 0, lDcaNegToPrimVertex = 0;
731 Double_t lV0CosineOfPointingAngle = 0;
732 Double_t lV0Radius = 0, lPt = 0;
733 Double_t lRapK0Short = 0, lRapLambda = 0;
734 Double_t lInvMassK0s = 0, lInvMassLambda = 0, lInvMassAntiLambda = 0;
735 Double_t lAlphaV0 = 0, lPtArmV0 = 0;
737 Double_t fMinV0Pt = 0;
738 Double_t fMaxV0Pt = 100;
741 nv0s = lESDevent->GetNumberOfV0s();
743 //for (Int_t iV0 = 0; iV0 < nv0s; iV0++)
744 for (Int_t iV0 = 0; iV0 < nv0s; iV0++) //extra-crazy test
745 {// This is the begining of the V0 loop
746 AliESDv0 *v0 = ((AliESDEvent*)lESDevent)->GetV0(iV0);
749 //---> Fix On-the-Fly candidates, count how many swapped
750 if( v0->GetParamN()->Charge() > 0 && v0->GetParamP()->Charge() < 0 ){
751 fHistSwappedV0Counter -> Fill( 1 );
753 fHistSwappedV0Counter -> Fill( 0 );
755 if ( fkUseOnTheFly ) CheckChargeV0(v0);
757 Double_t tDecayVertexV0[3]; v0->GetXYZ(tDecayVertexV0[0],tDecayVertexV0[1],tDecayVertexV0[2]);
760 v0->GetPxPyPz( tV0mom[0],tV0mom[1],tV0mom[2] );
761 Double_t lV0TotalMomentum = TMath::Sqrt(
762 tV0mom[0]*tV0mom[0]+tV0mom[1]*tV0mom[1]+tV0mom[2]*tV0mom[2] );
764 lV0Radius = TMath::Sqrt(tDecayVertexV0[0]*tDecayVertexV0[0]+tDecayVertexV0[1]*tDecayVertexV0[1]);
766 //Set Variables for later filling
767 fTreeVariableV0x = tDecayVertexV0[0];
768 fTreeVariableV0y = tDecayVertexV0[1];
769 fTreeVariableV0z = tDecayVertexV0[2];
771 //Set Variables for later filling
772 fTreeVariableV0Px = tV0mom[0];
773 fTreeVariableV0Py = tV0mom[1];
774 fTreeVariableV0Pz = tV0mom[2];
777 lRapK0Short = v0->RapK0Short();
778 lRapLambda = v0->RapLambda();
779 if ((lPt<fMinV0Pt)||(fMaxV0Pt<lPt)) continue;
781 UInt_t lKeyPos = (UInt_t)TMath::Abs(v0->GetPindex());
782 UInt_t lKeyNeg = (UInt_t)TMath::Abs(v0->GetNindex());
784 Double_t lMomPos[3]; v0->GetPPxPyPz(lMomPos[0],lMomPos[1],lMomPos[2]);
785 Double_t lMomNeg[3]; v0->GetNPxPyPz(lMomNeg[0],lMomNeg[1],lMomNeg[2]);
787 AliESDtrack *pTrack=((AliESDEvent*)lESDevent)->GetTrack(lKeyPos);
788 AliESDtrack *nTrack=((AliESDEvent*)lESDevent)->GetTrack(lKeyNeg);
789 if (!pTrack || !nTrack) {
790 Printf("ERROR: Could not retreive one of the daughter track");
794 //Daughter Eta for Eta selection, afterwards
795 fTreeVariableNegEta = nTrack->Eta();
796 fTreeVariablePosEta = pTrack->Eta();
798 // Filter like-sign V0 (next: add counter and distribution)
799 if ( pTrack->GetSign() == nTrack->GetSign()){
803 //________________________________________________________________________
804 // Track quality cuts
805 Float_t lPosTrackCrossedRows = pTrack->GetTPCClusterInfo(2,1);
806 Float_t lNegTrackCrossedRows = nTrack->GetTPCClusterInfo(2,1);
807 fTreeVariableLeastNbrCrossedRows = (Int_t) lPosTrackCrossedRows;
808 if( lNegTrackCrossedRows < fTreeVariableLeastNbrCrossedRows )
809 fTreeVariableLeastNbrCrossedRows = (Int_t) lNegTrackCrossedRows;
811 // TPC refit condition (done during reconstruction for Offline but not for On-the-fly)
812 if( !(pTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;
813 if( !(nTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;
816 fTreeVariablePosTrackStatus = pTrack->GetStatus();
817 fTreeVariableNegTrackStatus = nTrack->GetStatus();
819 if ( ( ( ( pTrack->GetTPCClusterInfo(2,1) ) < 70 ) || ( ( nTrack->GetTPCClusterInfo(2,1) ) < 70 ) )&&(fkTakeAllTracks==kFALSE) ) continue;
821 //GetKinkIndex condition
822 if( pTrack->GetKinkIndex(0)>0 || nTrack->GetKinkIndex(0)>0 ) continue;
824 //Findable clusters > 0 condition
825 if( pTrack->GetTPCNclsF()<=0 || nTrack->GetTPCNclsF()<=0 ) continue;
827 //Compute ratio Crossed Rows / Findable clusters
828 //Note: above test avoids division by zero!
829 Float_t lPosTrackCrossedRowsOverFindable = lPosTrackCrossedRows / ((double)(pTrack->GetTPCNclsF()));
830 Float_t lNegTrackCrossedRowsOverFindable = lNegTrackCrossedRows / ((double)(nTrack->GetTPCNclsF()));
832 fTreeVariableLeastRatioCrossedRowsOverFindable = lPosTrackCrossedRowsOverFindable;
833 if( lNegTrackCrossedRowsOverFindable < fTreeVariableLeastRatioCrossedRowsOverFindable )
834 fTreeVariableLeastRatioCrossedRowsOverFindable = lNegTrackCrossedRowsOverFindable;
836 //Lowest Cut Level for Ratio Crossed Rows / Findable = 0.8, set here
837 if ( (fTreeVariableLeastRatioCrossedRowsOverFindable < 0.8)&&(fkTakeAllTracks==kFALSE) ) continue;
839 //End track Quality Cuts
840 //________________________________________________________________________
842 lDcaPosToPrimVertex = TMath::Abs(pTrack->GetD(tPrimaryVtxPosition[0],
843 tPrimaryVtxPosition[1],
846 lDcaNegToPrimVertex = TMath::Abs(nTrack->GetD(tPrimaryVtxPosition[0],
847 tPrimaryVtxPosition[1],
850 lOnFlyStatus = v0->GetOnFlyStatus();
851 lChi2V0 = v0->GetChi2V0();
852 lDcaV0Daughters = v0->GetDcaV0Daughters();
853 lDcaV0ToPrimVertex = v0->GetD(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],tPrimaryVtxPosition[2]);
854 lV0CosineOfPointingAngle = v0->GetV0CosineOfPointingAngle(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],tPrimaryVtxPosition[2]);
855 fTreeVariableV0CosineOfPointingAngle=lV0CosineOfPointingAngle;
857 // Getting invariant mass infos directly from ESD
858 v0->ChangeMassHypothesis(310);
859 lInvMassK0s = v0->GetEffMass();
860 v0->ChangeMassHypothesis(3122);
861 lInvMassLambda = v0->GetEffMass();
862 v0->ChangeMassHypothesis(-3122);
863 lInvMassAntiLambda = v0->GetEffMass();
864 lAlphaV0 = v0->AlphaV0();
865 lPtArmV0 = v0->PtArmV0();
867 fTreeVariablePt = v0->Pt();
868 fTreeVariableChi2V0 = lChi2V0;
869 fTreeVariableDcaV0ToPrimVertex = lDcaV0ToPrimVertex;
870 fTreeVariableDcaV0Daughters = lDcaV0Daughters;
871 fTreeVariableV0CosineOfPointingAngle = lV0CosineOfPointingAngle;
872 fTreeVariableV0Radius = lV0Radius;
873 fTreeVariableDcaPosToPrimVertex = lDcaPosToPrimVertex;
874 fTreeVariableDcaNegToPrimVertex = lDcaNegToPrimVertex;
875 fTreeVariableInvMassK0s = lInvMassK0s;
876 fTreeVariableInvMassLambda = lInvMassLambda;
877 fTreeVariableInvMassAntiLambda = lInvMassAntiLambda;
878 fTreeVariableRapK0Short = lRapK0Short;
879 fTreeVariableRapLambda = lRapLambda;
880 fTreeVariableAlphaV0 = lAlphaV0;
881 fTreeVariablePtArmV0 = lPtArmV0;
883 //Official means of acquiring N-sigmas
884 fTreeVariableNSigmasPosProton = fPIDResponse->NumberOfSigmasTPC( pTrack, AliPID::kProton );
885 fTreeVariableNSigmasPosPion = fPIDResponse->NumberOfSigmasTPC( pTrack, AliPID::kPion );
886 fTreeVariableNSigmasNegProton = fPIDResponse->NumberOfSigmasTPC( nTrack, AliPID::kProton );
887 fTreeVariableNSigmasNegPion = fPIDResponse->NumberOfSigmasTPC( nTrack, AliPID::kPion );
889 //This requires an Invariant Mass Hypothesis afterwards
890 fTreeVariableDistOverTotMom = TMath::Sqrt(
891 TMath::Power( tDecayVertexV0[0] - lBestPrimaryVtxPos[0] , 2) +
892 TMath::Power( tDecayVertexV0[1] - lBestPrimaryVtxPos[1] , 2) +
893 TMath::Power( tDecayVertexV0[2] - lBestPrimaryVtxPos[2] , 2)
895 fTreeVariableDistOverTotMom /= (lV0TotalMomentum+1e-10); //avoid division by zero, to be sure
897 //------------------------------------------------
899 //------------------------------------------------
901 // The conditionals are meant to decrease excessive
904 //First Selection: Reject OnFly
905 if( (lOnFlyStatus == 0 && fkUseOnTheFly == kFALSE) || (lOnFlyStatus != 0 && fkUseOnTheFly == kTRUE ) ){
906 //Second Selection: rough 20-sigma band, parametric.
907 //K0Short: Enough to parametrize peak broadening with linear function.
908 Double_t lUpperLimitK0Short = (5.63707e-01) + (1.14979e-02)*fTreeVariablePt;
909 Double_t lLowerLimitK0Short = (4.30006e-01) - (1.10029e-02)*fTreeVariablePt;
910 //Lambda: Linear (for higher pt) plus exponential (for low-pt broadening)
911 //[0]+[1]*x+[2]*TMath::Exp(-[3]*x)
912 Double_t lUpperLimitLambda = (1.13688e+00) + (5.27838e-03)*fTreeVariablePt + (8.42220e-02)*TMath::Exp(-(3.80595e+00)*fTreeVariablePt);
913 Double_t lLowerLimitLambda = (1.09501e+00) - (5.23272e-03)*fTreeVariablePt - (7.52690e-02)*TMath::Exp(-(3.46339e+00)*fTreeVariablePt);
915 if( (fTreeVariableInvMassLambda < lUpperLimitLambda && fTreeVariableInvMassLambda > lLowerLimitLambda ) ||
916 (fTreeVariableInvMassAntiLambda < lUpperLimitLambda && fTreeVariableInvMassAntiLambda > lLowerLimitLambda ) ||
917 (fTreeVariableInvMassK0s < lUpperLimitK0Short && fTreeVariableInvMassK0s > lLowerLimitK0Short ) ){
918 //Pre-selection in case this is AA...
919 if( fkIsNuclear == kFALSE ) fTree->Fill();
920 if( fkIsNuclear == kTRUE){
921 //If this is a nuclear collision___________________
922 // ... pre-filter with TPC, daughter eta selection
923 if( (fTreeVariableInvMassLambda < lUpperLimitLambda && fTreeVariableInvMassLambda > lLowerLimitLambda
924 && TMath::Abs(fTreeVariableNSigmasPosProton) < 6.0 && TMath::Abs(fTreeVariableNSigmasNegPion) < 6.0 ) ||
925 (fTreeVariableInvMassAntiLambda < lUpperLimitLambda && fTreeVariableInvMassAntiLambda > lLowerLimitLambda
926 && TMath::Abs(fTreeVariableNSigmasNegProton) < 6.0 && TMath::Abs(fTreeVariableNSigmasPosPion) < 6.0 ) ||
927 (fTreeVariableInvMassK0s < lUpperLimitK0Short && fTreeVariableInvMassK0s > lLowerLimitK0Short
928 && TMath::Abs(fTreeVariableNSigmasNegPion) < 6.0 && TMath::Abs(fTreeVariableNSigmasPosPion) < 6.0 ) ){
930 if ( TMath::Abs(fTreeVariableNegEta)<0.8 && TMath::Abs(fTreeVariablePosEta)<0.8 ) fTree->Fill();
932 }//end nuclear_____________________________________
936 //------------------------------------------------
938 //------------------------------------------------
940 }// This is the end of the V0 loop
943 PostData(1, fListHistV0);
949 //________________________________________________________________________
950 void AliAnalysisTaskExtractV0::Terminate(Option_t *)
952 // Draw result to the screen
953 // Called once at the end of the query
954 // This will draw the V0 candidate multiplicity, whose
955 // number of entries corresponds to the number of triggered events.
956 TList *cRetrievedList = 0x0;
957 cRetrievedList = (TList*)GetOutputData(1);
959 Printf("ERROR - AliAnalysisTaskExtractV0 : ouput data container list not available\n");
962 fHistV0MultiplicityForTrigEvt = dynamic_cast<TH1F*> ( cRetrievedList->FindObject("fHistV0MultiplicityForTrigEvt") );
963 if (!fHistV0MultiplicityForTrigEvt) {
964 Printf("ERROR - AliAnalysisTaskExtractV0 : fHistV0MultiplicityForTrigEvt not available");
967 TCanvas *canCheck = new TCanvas("AliAnalysisTaskExtractV0","V0 Multiplicity",10,10,510,510);
968 canCheck->cd(1)->SetLogy();
969 fHistV0MultiplicityForTrigEvt->SetMarkerStyle(22);
970 fHistV0MultiplicityForTrigEvt->DrawCopy("E");
973 //________________________________________________________________________
974 void AliAnalysisTaskExtractV0::CheckChargeV0(AliESDv0 *v0)
976 // This function checks charge of negative and positive daughter tracks.
977 // If incorrectly defined (onfly vertexer), swaps out.
978 if( v0->GetParamN()->Charge() > 0 && v0->GetParamP()->Charge() < 0 ){
979 //V0 daughter track swapping is required! Note: everything is swapped here... P->N, N->P
980 Long_t lCorrectNidx = v0->GetPindex();
981 Long_t lCorrectPidx = v0->GetNindex();
982 Double32_t lCorrectNmom[3];
983 Double32_t lCorrectPmom[3];
984 v0->GetPPxPyPz( lCorrectNmom[0], lCorrectNmom[1], lCorrectNmom[2] );
985 v0->GetNPxPyPz( lCorrectPmom[0], lCorrectPmom[1], lCorrectPmom[2] );
987 AliExternalTrackParam lCorrectParamN(
988 v0->GetParamP()->GetX() ,
989 v0->GetParamP()->GetAlpha() ,
990 v0->GetParamP()->GetParameter() ,
991 v0->GetParamP()->GetCovariance()
993 AliExternalTrackParam lCorrectParamP(
994 v0->GetParamN()->GetX() ,
995 v0->GetParamN()->GetAlpha() ,
996 v0->GetParamN()->GetParameter() ,
997 v0->GetParamN()->GetCovariance()
999 lCorrectParamN.SetMostProbablePt( v0->GetParamP()->GetMostProbablePt() );
1000 lCorrectParamP.SetMostProbablePt( v0->GetParamN()->GetMostProbablePt() );
1002 //Get Variables___________________________________________________
1003 Double_t lDcaV0Daughters = v0 -> GetDcaV0Daughters();
1004 Double_t lCosPALocal = v0 -> GetV0CosineOfPointingAngle();
1005 Bool_t lOnFlyStatusLocal = v0 -> GetOnFlyStatus();
1007 //Create Replacement Object_______________________________________
1008 AliESDv0 *v0correct = new AliESDv0(lCorrectParamN,lCorrectNidx,lCorrectParamP,lCorrectPidx);
1009 v0correct->SetDcaV0Daughters ( lDcaV0Daughters );
1010 v0correct->SetV0CosineOfPointingAngle ( lCosPALocal );
1011 v0correct->ChangeMassHypothesis ( kK0Short );
1012 v0correct->SetOnFlyStatus ( lOnFlyStatusLocal );
1014 //Reverse Cluster info..._________________________________________
1015 v0correct->SetClusters( v0->GetClusters( 1 ), v0->GetClusters ( 0 ) );
1018 //Proper cleanup..._______________________________________________
1019 v0correct->Delete();
1022 //Just another cross-check and output_____________________________
1023 if( v0->GetParamN()->Charge() > 0 && v0->GetParamP()->Charge() < 0 ) {
1024 AliWarning("Found Swapped Charges, tried to correct but something FAILED!");
1026 //AliWarning("Found Swapped Charges and fixed.");
1028 //________________________________________________________________
1030 //Don't touch it! ---
1031 //Printf("Ah, nice. Charges are already ordered...");