]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/FlavourJetTasks/AliAnalysisTaskSEPicoV0Maker.cxx
Coverity fixes
[u/mrichter/AliRoot.git] / PWGJE / FlavourJetTasks / AliAnalysisTaskSEPicoV0Maker.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2008, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /////////////////////////////////////////////////////////////
17 //
18 // AliAnalysisTaskSE for RecoDecay object (K0 short, Lambda,
19 // D mesons ...) filtering
20 //
21 // Author: X-M. Zhang, xmzhang@lbl.gov
22 /////////////////////////////////////////////////////////////
23
24 #include <iostream>
25
26 #include <TString.h>
27 #include <TH1D.h>
28 #include <TH2D.h>
29 #include <THnSparse.h>
30 #include <TMath.h>
31 #include <TClonesArray.h>
32 #include <TDatabasePDG.h>
33 #include <TParticle.h>
34
35 #include "AliAODEvent.h"
36 #include "AliESDEvent.h"
37 #include "AliMCEvent.h"
38 #include "AliStack.h"
39 #include "AliVVertex.h"
40 #include "AliAODv0.h"
41 #include "AliESDv0.h"
42 #include "AliAODTrack.h"
43 #include "AliESDtrack.h"
44 #include "AliStack.h"
45 #include "AliMCParticle.h"
46 #include "AliAODMCParticle.h"
47 #include "AliHeader.h"
48 #include "AliGenDPMjetEventHeader.h"
49 #include "AliV0vertexer.h"
50 #include "AliAnalysisUtils.h"
51 #include "AliInputEventHandler.h"
52 #include "AliCentrality.h"
53 #include "AliPIDResponse.h"
54
55 #include "AliPicoHeaderCJ.h"
56 #include "AliPicoV0RD.h"
57 #include "AliPicoV0MC.h"
58 #include "AliAnalysisTaskSEPicoV0Maker.h"
59
60 ClassImp(AliAnalysisTaskSEPicoV0Maker)
61 //=============================================================================
62
63 const Double_t AliAnalysisTaskSEPicoV0Maker::fgkMassPion   = 0.13957;
64 const Double_t AliAnalysisTaskSEPicoV0Maker::fgkMassKshort = 0.497614;
65 const Double_t AliAnalysisTaskSEPicoV0Maker::fgkMassProton = 0.938272;
66 const Double_t AliAnalysisTaskSEPicoV0Maker::fgkMassLambda = 1.11568;
67
68 //_____________________________________________________________________________
69 AliAnalysisTaskSEPicoV0Maker::AliAnalysisTaskSEPicoV0Maker() :
70 AliAnalysisTaskSE(),
71 fEventAOD(0),
72 fEventESD(0),
73 fCentInfo(0),
74 fRespoPID(0),
75 fAnaUtils(0),
76 fEventAcptMask(0),
77 fTriggerMask(0),
78 fCollisionType(0),
79 fCentEst(""),
80 fIsAnaInfoMC(kFALSE),
81 fIsRefitV0sESD(kFALSE),
82 fIsSkipFastOnly(kFALSE),
83 fIsDPMjetMC(kFALSE),
84 fRapidityShift(0.),
85 fCutMinEventVtxContr(0),
86 fCutMaxEventVzAbs(0.),
87 fCutMinV0Pt(0.),
88 fCutMaxV0Pt(0.),
89 fCutMinV0Rap(0.),
90 fCutMaxV0Rap(0.),
91 fCutMinDauPt(0.),
92 fCutMinDauEta(0.),
93 fCutMaxDauEta(0.),
94 fCutMaxV0Chi2(0.),
95 fCutMinV0Radius(0.),
96 fCutMaxV0Radius(0.),
97 fCutMaxDausDCA(0.),
98 fCutMinDauDCAtoPV(0.),
99 fCutMinDauXrowsTPC(0.),
100 fCutMinDauXrowsOverFindableClusTPC(0.),
101 fCutMaxKshortSigmaTPC(0.),
102 fCutMinKshortCosPA(0.),
103 fCutMaxKshortCtau(0.),
104 fCutMaxKshortArmFrac(0.),
105 fCutMinKshortDeltaM(0.),
106 fCutMaxLambdaSigmaTPC(0.),
107 fCutMinLambdaCosPA(0.),
108 fCutMaxLambdaCtau(0.),
109 fCutMaxLambdaArmFrac(0.),
110 fCutMinLambdaDeletaM(0.),
111 fPicoV0sClArr(0),
112 fOutputListEH(0),
113 fOutputListMC(0)
114 {
115 //
116 // Default constructor
117 //
118
119   for (Int_t i=3; i--;) fPrimaryVtx[i] = 0.;
120 }
121
122 //_____________________________________________________________________________
123 AliAnalysisTaskSEPicoV0Maker::AliAnalysisTaskSEPicoV0Maker(const char *name, Bool_t bIsMC) :
124 AliAnalysisTaskSE(name),
125 fEventAOD(0),
126 fEventESD(0),
127 fCentInfo(0),
128 fRespoPID(0),
129 fAnaUtils(0),
130 fEventAcptMask(0),
131 fTriggerMask(0),
132 fCollisionType(AliPicoHeaderCJ::kPP),
133 fCentEst("V0M"),
134 fIsAnaInfoMC(bIsMC),
135 fIsRefitV0sESD(kFALSE),
136 fIsSkipFastOnly(kFALSE),
137 fIsDPMjetMC(kFALSE),
138 fRapidityShift(0.),
139 fCutMinEventVtxContr(2),
140 fCutMaxEventVzAbs(10.),
141 fCutMinV0Pt(0.),
142 fCutMaxV0Pt(100.),
143 fCutMinV0Rap(-10.),
144 fCutMaxV0Rap(10.),
145 fCutMinDauPt(0.),
146 fCutMinDauEta(-10.),
147 fCutMaxDauEta(10.),
148 fCutMaxV0Chi2(33.),
149 fCutMinV0Radius(0.3),
150 fCutMaxV0Radius(200.),
151 fCutMaxDausDCA(1.5),
152 fCutMinDauDCAtoPV(0.05),
153 fCutMinDauXrowsTPC(70.),
154 fCutMinDauXrowsOverFindableClusTPC(0.8),
155 fCutMaxKshortSigmaTPC(-1.),
156 fCutMinKshortCosPA(0.95),
157 fCutMaxKshortCtau(30.),
158 fCutMaxKshortArmFrac(-1.),
159 fCutMinKshortDeltaM(0.003),
160 fCutMaxLambdaSigmaTPC(7.),
161 fCutMinLambdaCosPA(0.993),
162 fCutMaxLambdaCtau(40.),
163 fCutMaxLambdaArmFrac(-1.),
164 fCutMinLambdaDeletaM(-1.),
165 fPicoV0sClArr(0),
166 fOutputListEH(0),
167 fOutputListMC(0)
168 {
169 //
170 // Constructor
171 //
172
173   for (Int_t i=3; i--;) fPrimaryVtx[i] = 0.;
174
175   DefineOutput(1, TList::Class());
176   if (fIsAnaInfoMC) DefineOutput(2, TList::Class());
177 }
178
179 //_____________________________________________________________________________
180 AliAnalysisTaskSEPicoV0Maker::~AliAnalysisTaskSEPicoV0Maker()
181 {
182 //
183 // Default destructor
184 //
185
186
187   if (fEventAOD) { delete fEventAOD; fEventAOD = 0; }
188   if (fEventESD) { delete fEventESD; fEventESD = 0; }
189   if (fCentInfo) { delete fCentInfo; fCentInfo = 0; }
190   if (fRespoPID) { delete fRespoPID; fRespoPID = 0; }
191   if (fAnaUtils) { delete fAnaUtils; fAnaUtils = 0; }
192
193   if (fPicoV0sClArr) { delete fPicoV0sClArr; fPicoV0sClArr = 0; }
194   if (fOutputListEH) { delete fOutputListEH; fOutputListEH = 0; }
195   if (fOutputListMC) { delete fOutputListMC; fOutputListMC = 0; }
196 }
197
198 //_____________________________________________________________________________
199 void AliAnalysisTaskSEPicoV0Maker::Init()
200 {
201 //
202 //  AliAnalysisTaskSEPicoV0Maker::Init
203 //
204
205   return;
206 }
207
208 //_____________________________________________________________________________
209 void AliAnalysisTaskSEPicoV0Maker::UserCreateOutputObjects()
210 {
211 //
212 //  AliAnalysisTaskSEPicoV0Maker::UserCreateOutputObjects
213 //
214
215   InitAnalysis();
216 //=============================================================================
217
218   if (fIsAnaInfoMC) {
219     fPicoV0sClArr = new TClonesArray("AliPicoV0MC");
220     fPicoV0sClArr->SetName("PicoV0s");
221   } else {
222     fPicoV0sClArr = new TClonesArray("AliPicoV0RD");
223     fPicoV0sClArr->SetName("PicoV0s");
224   }
225 //=============================================================================
226
227   fOutputListEH = new TList();
228   fOutputListEH->SetOwner();
229   CreateHistogramsEH();
230   PostData(1, fOutputListEH);
231
232   if (fIsAnaInfoMC) {
233     fOutputListMC = new TList();
234     fOutputListMC->SetOwner();
235     CreateHistogramsMC();
236     PostData(2, fOutputListMC);
237   }
238
239   return;
240 }
241
242 //_____________________________________________________________________________
243 void AliAnalysisTaskSEPicoV0Maker::UserExec(Option_t */*opt*/)
244 {
245 //
246 //  AliAnalysisTaskSEPicoV0Maker::UserExec
247 //
248
249   fPicoV0sClArr->Delete();
250   if (IsEventNotAcpt()) return;
251   if (!(InputEvent()->FindListObject("PicoV0s"))) InputEvent()->AddObject(fPicoV0sClArr);
252 //=============================================================================
253
254   FillHistogramsEH();
255   if (IsEventNotINEL()) return;
256   if (fIsAnaInfoMC) FillHistogramsMC();
257 //=============================================================================
258
259   if (IsEventNotMBsa()) return;
260
261
262   FillPicoV0s();
263   return;
264 }
265
266 //_____________________________________________________________________________
267 void AliAnalysisTaskSEPicoV0Maker::Terminate(Option_t */*opt*/)
268 {
269 //
270 //  AliAnalysisTaskSEPicoV0Maker::Terminate
271 //
272
273   return;
274 }
275
276 //_____________________________________________________________________________
277 void AliAnalysisTaskSEPicoV0Maker::NotifyRun()
278 {
279 //
280 //  AliAnalysisTaskSEPicoV0Maker::NotifyRun
281 //
282
283   return;
284 }
285
286 //_____________________________________________________________________________
287 void AliAnalysisTaskSEPicoV0Maker::FillPicoV0s()
288 {
289 //
290 //  AliAnalysisTaskSEPicoV0Maker::FillPicoRecoV0s
291 //
292
293   Int_t nV0s = 0;
294   AliAODv0 *pV0AOD = 0;
295   AliESDv0 *pV0ESD = 0;
296   if (fEventAOD) nV0s = fEventAOD->GetNumberOfV0s();
297   if (fEventESD) nV0s = fEventESD->GetNumberOfV0s();
298   if (nV0s<=0) return;
299 //=============================================================================
300
301
302   TH2D *hKshortPtInvM = dynamic_cast<TH2D*>(fOutputListEH->FindObject("hKshortPtInvM"));
303   TH2D *hLambdaPtInvM = dynamic_cast<TH2D*>(fOutputListEH->FindObject("hLambdaPtInvM"));
304   TH2D *hAntiLaPtInvM = dynamic_cast<TH2D*>(fOutputListEH->FindObject("hAntiLaPtInvM"));
305 //=============================================================================
306
307   AliPicoV0RD *pV0RD = 0;
308   AliPicoV0MC *pV0MC = 0;
309   Int_t nAt = fPicoV0sClArr->GetEntriesFast();
310
311   for (Int_t iV0=0; iV0<nV0s; iV0++) {
312
313     if (fEventAOD) {
314       pV0AOD = fEventAOD->GetV0(iV0); if (!pV0AOD) continue;
315
316       if (fIsAnaInfoMC) {
317         pV0MC = SelectV0CandidateMC(pV0AOD); if (!pV0MC) { pV0AOD = 0; continue; }
318       } else {
319         pV0RD = SelectV0CandidateRD(pV0AOD); if (!pV0RD) { pV0AOD = 0; continue; }
320       }
321     }
322
323     if (fEventESD) {
324       pV0ESD = fEventESD->GetV0(iV0); if (!pV0ESD) continue;
325
326       if (fIsAnaInfoMC) {
327         pV0MC = SelectV0CandidateMC(pV0ESD); if (!pV0MC) { pV0ESD = 0; continue; }
328       } else {
329         pV0RD = SelectV0CandidateRD(pV0ESD); if (!pV0RD) { pV0ESD = 0; continue; }
330       }
331     }
332
333     if (pV0RD) {
334       pV0RD->FillKshortPtInvM(hKshortPtInvM);
335       pV0RD->FillLambdaPtInvM(hLambdaPtInvM);
336       pV0RD->FillAntiLaPtInvM(hAntiLaPtInvM);
337       new ((*fPicoV0sClArr)[nAt++]) AliPicoV0RD(*pV0RD);
338       delete pV0RD; pV0RD=0;
339     }
340
341     if (pV0MC) {
342       pV0MC->FillKshortPtInvM(hKshortPtInvM);
343       pV0MC->FillLambdaPtInvM(hLambdaPtInvM);
344       pV0MC->FillAntiLaPtInvM(hAntiLaPtInvM);
345       new ((*fPicoV0sClArr)[nAt++]) AliPicoV0MC(*pV0MC);
346       delete pV0MC; pV0MC=0;
347     }
348
349     if (pV0AOD) pV0AOD=0;
350     if (pV0ESD) pV0ESD=0;
351   }
352
353   return;
354 }
355
356 //_____________________________________________________________________________
357 AliPicoV0RD* AliAnalysisTaskSEPicoV0Maker::SelectV0CandidateRD(AliAODv0 const *pV0)
358 {
359 //
360 //  AliAnalysisTaskSEPicoV0Maker::SelectV0CandidateRD
361 //
362
363   if  (pV0->GetOnFlyStatus()) return 0x0;
364   if ((pV0->Chi2V0())>fCutMaxV0Chi2) return 0x0;
365 //=============================================================================
366
367   Double_t dV0Pt  = pV0->Pt(); if ((dV0Pt<fCutMinV0Pt) || (dV0Pt>fCutMaxV0Pt)) return 0x0;
368   Double_t dKaRap = pV0->RapK0Short(); if ((dKaRap<fCutMinV0Rap) || (dKaRap>fCutMaxV0Rap)) return 0x0;
369   Double_t dLaRap = pV0->RapLambda();  if ((dLaRap<fCutMinV0Rap) || (dLaRap>fCutMaxV0Rap)) return 0x0;
370 //=============================================================================
371
372   Double_t dV0Vtx[3]; pV0->GetXYZ(dV0Vtx);
373   Double_t dV0Radius = TMath::Sqrt(dV0Vtx[0]*dV0Vtx[0] + dV0Vtx[1]*dV0Vtx[1]);
374   if ((dV0Radius<fCutMinV0Radius) || (dV0Radius>fCutMaxV0Radius)) return 0x0;
375
376   Double_t dDausDCA = pV0->DcaV0Daughters(); if (dDausDCA>fCutMaxDausDCA) return 0x0;
377   Double_t dPosDCAtoPV = pV0->DcaPosToPrimVertex(); if (dPosDCAtoPV<fCutMinDauDCAtoPV) return 0x0;
378   Double_t dNegDCAtoPV = pV0->DcaNegToPrimVertex(); if (dNegDCAtoPV<fCutMinDauDCAtoPV) return 0x0;
379 //=============================================================================
380
381   AliAODTrack *pDauPos = (AliAODTrack*)pV0->GetDaughter(0); if (!pDauPos) return 0x0;
382   AliAODTrack *pDauNeg = (AliAODTrack*)pV0->GetDaughter(1); if (!pDauNeg) return 0x0;
383
384   if (!(pDauPos->GetStatus() & AliESDtrack::kTPCrefit)) return 0x0;
385   if (!(pDauNeg->GetStatus() & AliESDtrack::kTPCrefit)) return 0x0;
386
387   if ((pDauPos->GetProdVertex()->GetType())==((Char_t)AliAODVertex::kKink)) return 0x0;
388   if ((pDauNeg->GetProdVertex()->GetType())==((Char_t)AliAODVertex::kKink)) return 0x0;
389
390   Float_t dPosXrowsTPC = pDauPos->GetTPCClusterInfo(2,1);
391   Float_t dNegXrowsTPC = pDauNeg->GetTPCClusterInfo(2,1);
392   Float_t dDauXrowsTPC = dPosXrowsTPC; if (dDauXrowsTPC>dNegXrowsTPC) dDauXrowsTPC = dNegXrowsTPC;
393   if (dDauXrowsTPC<fCutMinDauXrowsTPC) return 0x0;
394
395   UShort_t wPosTPCNClsF = pDauPos->GetTPCNclsF(); if (wPosTPCNClsF<=0) return 0x0;
396   UShort_t wNegTPCNClsF = pDauNeg->GetTPCNclsF(); if (wNegTPCNClsF<=0) return 0x0;
397   Double_t dPosXrowsOverFindableClusTPC = ((Double_t)dPosXrowsTPC) / ((Double_t)wPosTPCNClsF);
398   Double_t dNegXrowsOverFindableClusTPC = ((Double_t)dNegXrowsTPC) / ((Double_t)wNegTPCNClsF);
399
400   Double_t dDauXrowsOverFindableClusTPC = dPosXrowsOverFindableClusTPC;
401   if (dDauXrowsOverFindableClusTPC>dNegXrowsOverFindableClusTPC) dDauXrowsOverFindableClusTPC = dNegXrowsOverFindableClusTPC;
402   if (dDauXrowsOverFindableClusTPC<fCutMinDauXrowsOverFindableClusTPC) return 0x0;
403 //=============================================================================
404
405   Short_t nPosCharge = pDauPos->Charge();
406   Short_t nNegCharge = pDauNeg->Charge();
407   if ((nPosCharge==0) || (nNegCharge==0) || (nPosCharge==nNegCharge)) return 0x0;
408
409   Double_t dPosPxPyPz[3] = { 0., 0., 0. };
410   Double_t dNegPxPyPz[3] = { 0., 0., 0. };
411   if ((nPosCharge<0) && (nNegCharge>0)) {
412     pDauPos = (AliAODTrack*)pV0->GetDaughter(1);
413     pDauNeg = (AliAODTrack*)pV0->GetDaughter(0);
414
415     dPosPxPyPz[0] = pV0->MomNegX(); dPosPxPyPz[1] = pV0->MomNegY(); dPosPxPyPz[2] = pV0->MomNegZ();
416     dNegPxPyPz[0] = pV0->MomPosX(); dNegPxPyPz[1] = pV0->MomPosY(); dNegPxPyPz[2] = pV0->MomPosZ();
417   } else {
418     dPosPxPyPz[0] = pV0->MomPosX(); dPosPxPyPz[1] = pV0->MomPosY(); dPosPxPyPz[2] = pV0->MomPosZ();
419     dNegPxPyPz[0] = pV0->MomNegX(); dNegPxPyPz[1] = pV0->MomNegY(); dNegPxPyPz[2] = pV0->MomNegZ();
420   }
421
422   TVector3 v3Pos(dPosPxPyPz);
423   TVector3 v3Neg(dNegPxPyPz);
424
425   if ((v3Pos.Pt()<fCutMinDauPt) || (v3Neg.Pt()<fCutMinDauPt)) return 0x0;
426   Double_t dPosEta = v3Pos.Eta(); if ((dPosEta<fCutMinDauEta) || (dPosEta>fCutMaxDauEta)) return 0x0;
427   Double_t dNegEta = v3Neg.Eta(); if ((dNegEta<fCutMinDauEta) || (dNegEta>fCutMaxDauEta)) return 0x0;
428 //=============================================================================
429
430   Bool_t bIsKshort = kFALSE;
431   Bool_t bIsLambda = kFALSE;
432   Bool_t bIsAntiLa = kFALSE;
433
434   Float_t dPosPionSigmaTPC = fRespoPID->NumberOfSigmasTPC(pDauPos,AliPID::kPion);
435   Float_t dNegPionSigmaTPC = fRespoPID->NumberOfSigmasTPC(pDauNeg,AliPID::kPion);
436
437   Float_t dPosProtonSigmaTPC = fRespoPID->NumberOfSigmasTPC(pDauPos,AliPID::kProton);
438   Float_t dNegProtonSigmaTPC = fRespoPID->NumberOfSigmasTPC(pDauNeg,AliPID::kProton);
439
440   if (fCutMaxKshortSigmaTPC>0.) {
441     bIsKshort = ((TMath::Abs(dPosPionSigmaTPC)<fCutMaxKshortSigmaTPC) &&
442                  (TMath::Abs(dNegPionSigmaTPC)<fCutMaxKshortSigmaTPC));
443   } else {
444     bIsKshort = kTRUE;
445   }
446
447   if (fCutMaxLambdaSigmaTPC>0.) {
448     bIsLambda = ((TMath::Abs(dPosProtonSigmaTPC)<fCutMaxLambdaSigmaTPC) &&
449                  (TMath::Abs(dNegPionSigmaTPC)  <fCutMaxLambdaSigmaTPC));
450
451     bIsAntiLa = ((TMath::Abs(dPosPionSigmaTPC)  <fCutMaxLambdaSigmaTPC) &&
452                  (TMath::Abs(dNegProtonSigmaTPC)<fCutMaxLambdaSigmaTPC));
453   } else {
454     bIsLambda = kTRUE;
455     bIsAntiLa = kTRUE;
456   }
457
458   if (!(bIsKshort || bIsLambda || bIsAntiLa)) return 0x0;
459 //=============================================================================
460
461   Double_t dV0CosPA = pV0->CosPointingAngle(fPrimaryVtx);
462
463   if (bIsKshort) if (dV0CosPA<fCutMinKshortCosPA) {
464     bIsKshort = kFALSE;
465   }
466
467   if (bIsLambda || bIsAntiLa) if (dV0CosPA<fCutMinLambdaCosPA) {
468     bIsLambda = kFALSE;
469     bIsAntiLa = kFALSE;
470   }
471
472   if (!(bIsKshort || bIsLambda || bIsAntiLa)) return 0x0;
473 //=============================================================================
474
475   Double_t dV0DistToPV = 0.;
476   for (Int_t i=0; i<3; i++) dV0DistToPV +=  ((dV0Vtx[i]-fPrimaryVtx[i]) * (dV0Vtx[i]-fPrimaryVtx[i]));
477   Double_t dV0DistToPVoverP = TMath::Sqrt(dV0DistToPV) / (pV0->P()+1e-10);
478
479   if (bIsKshort) if ((dV0DistToPVoverP*fgkMassKshort)>fCutMaxKshortCtau) {
480     bIsKshort = kFALSE;
481   }
482
483   if (bIsLambda || bIsAntiLa) if ((dV0DistToPVoverP*fgkMassLambda)>fCutMaxLambdaCtau) {
484     bIsLambda = kFALSE;
485     bIsAntiLa = kFALSE;
486   }
487
488   if (!(bIsKshort || bIsLambda || bIsAntiLa)) return 0x0;
489 //=============================================================================
490
491   Double_t dV0ArmFrac = pV0->PtArmV0() / (TMath::Abs(pV0->AlphaV0())+1e-12);
492
493   if (bIsKshort && (fCutMaxKshortArmFrac>0.)) if (dV0ArmFrac>fCutMaxKshortArmFrac) {
494     bIsKshort = kFALSE;
495   }
496
497   if ((bIsLambda || bIsAntiLa) && fCutMaxLambdaArmFrac>0.) if (dV0ArmFrac>fCutMaxLambdaArmFrac) {
498     bIsLambda = kFALSE;
499     bIsAntiLa = kFALSE;
500   }
501
502   if (!(bIsKshort || bIsLambda || bIsAntiLa)) return 0x0;
503 //=============================================================================
504
505   TLorentzVector vPosPion; vPosPion.SetVectM(v3Pos, fgkMassPion);
506   TLorentzVector vNegPion; vNegPion.SetVectM(v3Neg, fgkMassPion);
507
508   TLorentzVector vPosProton; vPosProton.SetVectM(v3Pos, fgkMassProton);
509   TLorentzVector vNegProton; vNegProton.SetVectM(v3Neg, fgkMassProton);
510
511   TLorentzVector vKshort = vPosPion   + vNegPion;
512   TLorentzVector vLamvda = vPosProton + vNegPion;
513   TLorentzVector vAntiLa = vNegProton + vPosPion;
514
515   Double_t dKshortInvM = vKshort.M();
516   Double_t dLambdaInvM = vLamvda.M();
517   Double_t dAntiLaInvM = vAntiLa.M();
518
519   if (bIsKshort) {
520     Double_t dLower = 0.430006 - 0.0110029*dV0Pt;
521     Double_t dUpper = 0.563707 + 0.0114979*dV0Pt;
522
523     if ((dKshortInvM<dLower) || (dKshortInvM>dUpper)) bIsKshort = kFALSE;
524   }
525
526   if (bIsLambda || bIsAntiLa) {
527     Double_t dLower = 1.09501 - 0.00523272*dV0Pt - 0.075269*TMath::Exp(-3.46339*dV0Pt);
528     Double_t dUpper = 1.13688 + 0.00527838*dV0Pt + 0.084222*TMath::Exp(-3.80595*dV0Pt);
529
530     if (bIsLambda) if ((dLambdaInvM<dLower) || (dLambdaInvM>dUpper)) bIsLambda = kFALSE;
531     if (bIsAntiLa) if ((dAntiLaInvM<dLower) || (dAntiLaInvM>dUpper)) bIsAntiLa = kFALSE;
532   }
533
534   if (!(bIsKshort || bIsLambda || bIsAntiLa)) return 0x0;
535 //=============================================================================
536
537   if (bIsKshort && (fCutMinKshortDeltaM>0.)) {
538     if ((TMath::Abs(dLambdaInvM-fgkMassLambda)<fCutMinKshortDeltaM) ||
539         (TMath::Abs(dAntiLaInvM-fgkMassLambda)<fCutMinKshortDeltaM)) bIsKshort = kFALSE;
540   }
541
542   if ((bIsLambda || bIsAntiLa) && (fCutMinLambdaDeletaM>0.)) {
543     if ((TMath::Abs(dKshortInvM-fgkMassKshort)<fCutMinLambdaDeletaM)) {
544       bIsLambda = kFALSE;
545       bIsAntiLa = kFALSE;
546     }
547   }
548
549   if (!(bIsKshort || bIsLambda || bIsAntiLa)) return 0x0;
550 //=============================================================================
551
552   UInt_t wMask = 0;
553   if (bIsKshort) wMask |= AliPicoHeaderCJ::kKshort;
554   if (bIsLambda) wMask |= AliPicoHeaderCJ::kLambda;
555   if (bIsAntiLa) wMask |= AliPicoHeaderCJ::kAntiLambda;
556
557   Bool_t bPosInJC = kFALSE;
558   Bool_t bNegInJC = kFALSE;
559
560 /*AliAODTrack *pTrkAOD = 0;
561   Int_t idPos = pDauPos->GetID();
562   Int_t idNeg = pDauNeg->GetID();
563   if (fJetContisClArr) for (Int_t i=0; i<fJetContisClArr->GetEntriesFast(); i++) {
564     pTrkAOD = static_cast<AliAODTrack*>(fJetContisClArr->At(i)); if (!pTrkAOD) continue;
565
566     Int_t id = pTrkAOD->GetID();
567     if (idPos==id) bPosInJC = kTRUE;
568     if (idNeg==id) bNegInJC = kTRUE;
569     if (bPosInJC && bNegInJC) { pTrkAOD = 0; break; }
570
571     pTrkAOD = 0;
572   }*/
573
574   AliPicoV0RD *pPicoV0 = new AliPicoV0RD(wMask,
575                                          dV0Radius,
576                                          dV0CosPA,
577                                          dV0DistToPVoverP,
578                                          dDausDCA,
579                                          dPosDCAtoPV,
580                                          dNegDCAtoPV,
581                                          dDauXrowsTPC,
582                                          dDauXrowsOverFindableClusTPC,
583                                          v3Pos.Px(), v3Pos.Py(), v3Pos.Pz(),
584                                          v3Neg.Px(), v3Neg.Py(), v3Neg.Pz(),
585                                          bPosInJC, bNegInJC,
586                                          dPosPionSigmaTPC, dPosProtonSigmaTPC,
587                                          dNegPionSigmaTPC, dNegProtonSigmaTPC);
588
589   return pPicoV0;
590 }
591
592 //_____________________________________________________________________________
593 AliPicoV0RD* AliAnalysisTaskSEPicoV0Maker::SelectV0CandidateRD(AliESDv0 const *pV0)
594 {
595 //
596 //  AliAnalysisTaskSEPicoV0Maker::SelectV0CandidateRD
597 //
598
599   if (pV0->GetOnFlyStatus()) return 0x0;
600   if (pV0->GetChi2V0()>fCutMaxV0Chi2) return 0x0;
601 //=============================================================================
602
603   Double_t dV0Pt = pV0->Pt(); if ((dV0Pt<fCutMinV0Pt) || (dV0Pt>fCutMaxV0Pt)) return 0x0;
604   Double_t dKaRap = pV0->RapK0Short(); if ((dKaRap<fCutMinV0Rap) || (dKaRap>fCutMaxV0Rap)) return 0x0;
605   Double_t dLaRap = pV0->RapLambda();  if ((dLaRap<fCutMinV0Rap) || (dLaRap>fCutMaxV0Rap)) return 0x0;
606 //=============================================================================
607
608   Double_t dV0Vtx[3];  pV0->GetXYZ(dV0Vtx[0], dV0Vtx[1], dV0Vtx[2]);
609   Double_t dV0Radius = TMath::Sqrt(dV0Vtx[0]*dV0Vtx[0] + dV0Vtx[1]*dV0Vtx[1]);
610   if ((dV0Radius<fCutMinV0Radius) || (dV0Radius>fCutMaxV0Radius)) return 0x0;
611
612   Double_t dDausDCA = pV0->GetDcaV0Daughters(); if (dDausDCA>fCutMaxDausDCA) return 0x0;
613 //=============================================================================
614
615   Int_t nPosIndex = pV0->GetPindex(); if (nPosIndex<0) return 0x0;
616   Int_t nNegIndex = pV0->GetNindex(); if (nNegIndex<0) return 0x0;
617
618   AliESDtrack *pDauPos = fEventESD->GetTrack(nPosIndex); if (!pDauPos) return 0x0;
619   AliESDtrack *pDauNeg = fEventESD->GetTrack(nNegIndex); if (!pDauNeg) return 0x0;
620
621   Double_t dMegField = fEventESD->GetMagneticField();
622   Double_t dPosDCAtoPV = TMath::Abs(pDauPos->GetD(fPrimaryVtx[0],fPrimaryVtx[1],dMegField)); if (dPosDCAtoPV<fCutMinDauDCAtoPV) return 0x0;
623   Double_t dNegDCAtoPV = TMath::Abs(pDauNeg->GetD(fPrimaryVtx[0],fPrimaryVtx[1],dMegField)); if (dNegDCAtoPV<fCutMinDauDCAtoPV) return 0x0;
624 //=============================================================================
625
626   if (!(pDauPos->GetStatus() & AliESDtrack::kTPCrefit)) return 0x0;
627   if (!(pDauNeg->GetStatus() & AliESDtrack::kTPCrefit)) return 0x0;
628   if ((pDauPos->GetKinkIndex(0)>0) || (pDauNeg->GetKinkIndex(0)>0)) return 0x0;
629
630   Float_t dPosXrowsTPC = pDauPos->GetTPCClusterInfo(2,1);
631   Float_t dNegXrowsTPC = pDauNeg->GetTPCClusterInfo(2,1);
632   Float_t dDauXrowsTPC = dPosXrowsTPC; if (dDauXrowsTPC>dNegXrowsTPC) dDauXrowsTPC = dNegXrowsTPC;
633   if (dDauXrowsTPC<fCutMinDauXrowsTPC) return 0x0;
634
635   UShort_t wPosTPCNClsF = pDauPos->GetTPCNclsF(); if (wPosTPCNClsF<=0) return 0x0;
636   UShort_t wNegTPCNClsF = pDauNeg->GetTPCNclsF(); if (wNegTPCNClsF<=0) return 0x0;
637   Double_t dPosXrowsOverFindableClusTPC = ((Double_t)dPosXrowsTPC) / ((Double_t)wPosTPCNClsF);
638   Double_t dNegXrowsOverFindableClusTPC = ((Double_t)dNegXrowsTPC) / ((Double_t)wNegTPCNClsF);
639
640   Double_t dDauXrowsOverFindableClusTPC = dPosXrowsOverFindableClusTPC;
641   if (dDauXrowsOverFindableClusTPC>dNegXrowsOverFindableClusTPC) dDauXrowsOverFindableClusTPC = dNegXrowsOverFindableClusTPC;
642   if (dDauXrowsOverFindableClusTPC<fCutMinDauXrowsOverFindableClusTPC) return 0x0;
643 //=============================================================================
644
645   Short_t nPosCharge = pDauPos->Charge();
646   Short_t nNegCharge = pDauNeg->Charge();
647   if ((nPosCharge==0) || (nNegCharge==0) || (nPosCharge==nNegCharge)) return 0x0;
648
649   Double_t dPosPxPyPz[3] = { 0., 0., 0. };
650   Double_t dNegPxPyPz[3] = { 0., 0., 0. };
651   if ((nPosCharge<0) && (nNegCharge>0)) {
652     pDauPos = fEventESD->GetTrack(nNegIndex);
653     pDauNeg = fEventESD->GetTrack(nPosIndex);
654
655     pV0->GetNPxPyPz(dPosPxPyPz[0], dPosPxPyPz[1], dPosPxPyPz[2]);
656     pV0->GetPPxPyPz(dNegPxPyPz[0], dNegPxPyPz[1], dNegPxPyPz[2]);
657   } else {
658     pV0->GetPPxPyPz(dPosPxPyPz[0], dPosPxPyPz[1], dPosPxPyPz[2]);
659     pV0->GetNPxPyPz(dNegPxPyPz[0], dNegPxPyPz[1], dNegPxPyPz[2]);
660   }
661
662   TVector3 v3Pos(dPosPxPyPz);
663   TVector3 v3Neg(dNegPxPyPz);
664
665   if ((v3Pos.Pt()<fCutMinDauPt) || (v3Neg.Pt()<fCutMinDauPt)) return 0x0;
666   Double_t dPosEta = v3Pos.Eta(); if ((dPosEta<fCutMinDauEta) || (dPosEta>fCutMaxDauEta)) return 0x0;
667   Double_t dNegEta = v3Neg.Eta(); if ((dNegEta<fCutMinDauEta) || (dNegEta>fCutMaxDauEta)) return 0x0;
668 //=============================================================================
669
670   Bool_t bIsKshort = kFALSE;
671   Bool_t bIsLambda = kFALSE;
672   Bool_t bIsAntiLa = kFALSE;
673
674   Float_t dPosPionSigmaTPC = fRespoPID->NumberOfSigmasTPC(pDauPos,AliPID::kPion);
675   Float_t dNegPionSigmaTPC = fRespoPID->NumberOfSigmasTPC(pDauNeg,AliPID::kPion);
676
677   Float_t dPosProtonSigmaTPC = fRespoPID->NumberOfSigmasTPC(pDauPos,AliPID::kProton);
678   Float_t dNegProtonSigmaTPC = fRespoPID->NumberOfSigmasTPC(pDauNeg,AliPID::kProton);
679
680   if (fCutMaxKshortSigmaTPC>0.) {
681     bIsKshort = ((TMath::Abs(dPosPionSigmaTPC)<fCutMaxKshortSigmaTPC) &&
682                  (TMath::Abs(dNegPionSigmaTPC)<fCutMaxKshortSigmaTPC));
683   } else {
684     bIsKshort = kTRUE;
685   }
686
687   if (fCutMaxLambdaSigmaTPC>0.) {
688     bIsLambda = ((TMath::Abs(dPosProtonSigmaTPC)<fCutMaxLambdaSigmaTPC) &&
689                  (TMath::Abs(dNegPionSigmaTPC)  <fCutMaxLambdaSigmaTPC));
690
691     bIsAntiLa = ((TMath::Abs(dPosPionSigmaTPC)  <fCutMaxLambdaSigmaTPC) &&
692                  (TMath::Abs(dNegProtonSigmaTPC)<fCutMaxLambdaSigmaTPC));
693   } else {
694     bIsLambda = kTRUE;
695     bIsAntiLa = kTRUE;
696   }
697
698   if (!(bIsKshort || bIsLambda || bIsAntiLa)) return 0x0;
699 //=============================================================================
700
701   Double_t dV0CosPA = pV0->GetV0CosineOfPointingAngle(fPrimaryVtx[0], fPrimaryVtx[1], fPrimaryVtx[2]);
702
703   if (bIsKshort) if (dV0CosPA<fCutMinKshortCosPA) {
704     bIsKshort = kFALSE;
705   }
706
707   if (bIsLambda || bIsAntiLa) if (dV0CosPA<fCutMinLambdaCosPA) {
708     bIsLambda = kFALSE;
709     bIsAntiLa = kFALSE;
710   }
711
712   if (!(bIsKshort || bIsLambda || bIsAntiLa)) return 0x0;
713 //=============================================================================
714
715   Double_t dV0DistToPV = 0.;
716   for (Int_t i=0; i<3; i++) dV0DistToPV +=  ((dV0Vtx[i]-fPrimaryVtx[i]) * (dV0Vtx[i]-fPrimaryVtx[i]));
717   Double_t dV0DistToPVoverP = TMath::Sqrt(dV0DistToPV) / (pV0->P()+1e-10);
718
719   if (bIsKshort) if ((dV0DistToPVoverP*fgkMassKshort)>fCutMaxKshortCtau) {
720     bIsKshort = kFALSE;
721   }
722
723   if (bIsLambda || bIsAntiLa) if ((dV0DistToPVoverP*fgkMassLambda)>fCutMaxLambdaCtau) {
724     bIsLambda = kFALSE;
725     bIsAntiLa = kFALSE;
726   }
727
728   if (!(bIsKshort || bIsLambda || bIsAntiLa)) return 0x0;
729 //=============================================================================
730
731   Double_t dV0ArmFrac = pV0->PtArmV0() / (TMath::Abs(pV0->AlphaV0())+1e-12);
732
733   if (bIsKshort && (fCutMaxKshortArmFrac>0.)) if (dV0ArmFrac>fCutMaxKshortArmFrac) {
734     bIsKshort = kFALSE;
735   }
736
737   if ((bIsLambda && bIsAntiLa) && (fCutMaxLambdaArmFrac>0.)) if (dV0ArmFrac>fCutMaxLambdaArmFrac) {
738     bIsLambda = kFALSE;
739     bIsAntiLa = kFALSE;
740   }
741
742   if (!(bIsKshort || bIsLambda || bIsAntiLa)) return 0x0;
743 //=============================================================================
744
745   TLorentzVector vPosPion; vPosPion.SetVectM(v3Pos, fgkMassPion);
746   TLorentzVector vNegPion; vNegPion.SetVectM(v3Neg, fgkMassPion);
747
748   TLorentzVector vPosProton; vPosProton.SetVectM(v3Pos, fgkMassProton);
749   TLorentzVector vNegProton; vNegProton.SetVectM(v3Neg, fgkMassProton);
750
751   TLorentzVector vKshort = vPosPion   + vNegPion;
752   TLorentzVector vLamvda = vPosProton + vNegPion;
753   TLorentzVector vAntiLa = vNegProton + vPosPion;
754
755   Double_t dKshortInvM = vKshort.M();
756   Double_t dLambdaInvM = vLamvda.M();
757   Double_t dAntiLaInvM = vAntiLa.M();
758
759   if (bIsKshort) {
760     Double_t dLower = 0.430006 - 0.0110029*dV0Pt;
761     Double_t dUpper = 0.563707 + 0.0114979*dV0Pt;
762
763     if ((dKshortInvM<dLower) || (dKshortInvM>dUpper)) bIsKshort = kFALSE;
764   }
765
766   if (bIsLambda || bIsAntiLa) {
767     Double_t dLower = 1.09501 - 0.00523272*dV0Pt - 0.075269*TMath::Exp(-3.46339*dV0Pt);
768     Double_t dUpper = 1.13688 + 0.00527838*dV0Pt + 0.084222*TMath::Exp(-3.80595*dV0Pt);
769
770     if (bIsLambda) if ((dLambdaInvM<dLower) || (dLambdaInvM>dUpper)) bIsLambda = kFALSE;
771     if (bIsAntiLa) if ((dAntiLaInvM<dLower) || (dAntiLaInvM>dUpper)) bIsAntiLa = kFALSE;
772   }
773
774   if (!(bIsKshort || bIsLambda || bIsAntiLa)) return 0x0;
775 //=============================================================================
776
777   if (bIsKshort && (fCutMinKshortDeltaM>0.)) {
778     if ((TMath::Abs(dLambdaInvM-fgkMassLambda)<fCutMinKshortDeltaM) ||
779         (TMath::Abs(dAntiLaInvM-fgkMassLambda)<fCutMinKshortDeltaM)) {
780       bIsKshort = kFALSE;
781     }
782   }
783
784   if ((bIsLambda || bIsAntiLa) && (fCutMinLambdaDeletaM>0.)) {
785     if ((TMath::Abs(dKshortInvM-fgkMassKshort)<fCutMinLambdaDeletaM)) {
786       bIsLambda = kFALSE;
787       bIsAntiLa = kFALSE;
788     }
789   }
790
791   if (!(bIsKshort || bIsLambda || bIsAntiLa)) return 0x0;
792 //=============================================================================
793
794   UInt_t wMask = 0;
795   if (bIsKshort) wMask |= AliPicoHeaderCJ::kKshort;
796   if (bIsLambda) wMask |= AliPicoHeaderCJ::kLambda;
797   if (bIsAntiLa) wMask |= AliPicoHeaderCJ::kAntiLambda;
798
799   Bool_t bPosInJC = kFALSE;
800   Bool_t bNegInJC = kFALSE;
801   AliPicoV0RD *pPicoV0 = new AliPicoV0RD(wMask,
802                                          dV0Radius,
803                                          dV0CosPA,
804                                          dV0DistToPVoverP,
805                                          dDausDCA,
806                                          dPosDCAtoPV,
807                                          dNegDCAtoPV,
808                                          dDauXrowsTPC,
809                                          dDauXrowsOverFindableClusTPC,
810                                          v3Pos.Px(), v3Pos.Py(), v3Pos.Pz(),
811                                          v3Neg.Px(), v3Neg.Py(), v3Neg.Pz(),
812                                          bPosInJC, bNegInJC,
813                                          dPosPionSigmaTPC, dPosProtonSigmaTPC,
814                                          dNegPionSigmaTPC, dNegProtonSigmaTPC);
815   
816   return pPicoV0;
817 }
818
819 //_____________________________________________________________________________
820 AliPicoV0MC* AliAnalysisTaskSEPicoV0Maker::SelectV0CandidateMC(AliAODv0 const *pV0RD)
821 {
822 //
823 //  AliAnalysisTaskSEPicoV0Maker::SelectV0CandidateMC
824 //
825
826   if (pV0RD->GetOnFlyStatus()) return 0x0;
827   if ((pV0RD->Chi2V0())>fCutMaxV0Chi2) return 0x0;
828
829   Double_t dV0Pt = pV0RD->Pt(); if ((dV0Pt<fCutMinV0Pt) || (dV0Pt>fCutMaxV0Pt)) return 0x0;
830 //=============================================================================
831
832   Double_t dV0Vtx[3]; pV0RD->GetXYZ(dV0Vtx);
833   Double_t dV0Radius = TMath::Sqrt(dV0Vtx[0]*dV0Vtx[0] + dV0Vtx[1]*dV0Vtx[1]);
834   if ((dV0Radius<fCutMinV0Radius) || (dV0Radius>fCutMaxV0Radius)) return 0x0;
835
836   Double_t dDausDCA = pV0RD->DcaV0Daughters(); if (dDausDCA>fCutMaxDausDCA) return 0x0;
837   Double_t dPosDCAtoPV = pV0RD->DcaPosToPrimVertex(); if (dPosDCAtoPV<fCutMinDauDCAtoPV) return 0x0;
838   Double_t dNegDCAtoPV = pV0RD->DcaNegToPrimVertex(); if (dNegDCAtoPV<fCutMinDauDCAtoPV) return 0x0;
839 //=============================================================================
840
841   AliAODTrack *pDauPosRD = (AliAODTrack*)pV0RD->GetDaughter(0); if (!pDauPosRD) return 0x0;
842   AliAODTrack *pDauNegRD = (AliAODTrack*)pV0RD->GetDaughter(1); if (!pDauNegRD) return 0x0;
843
844   if (!(pDauPosRD->GetStatus() & AliESDtrack::kTPCrefit)) return 0x0;
845   if (!(pDauNegRD->GetStatus() & AliESDtrack::kTPCrefit)) return 0x0;
846
847   if ((pDauPosRD->GetProdVertex()->GetType())==((Char_t)AliAODVertex::kKink)) return 0x0;
848   if ((pDauNegRD->GetProdVertex()->GetType())==((Char_t)AliAODVertex::kKink)) return 0x0;
849
850   Float_t dPosXrowsTPC = pDauPosRD->GetTPCClusterInfo(2,1);
851   Float_t dNegXrowsTPC = pDauNegRD->GetTPCClusterInfo(2,1);
852   Float_t dDauXrowsTPC = dPosXrowsTPC; if (dDauXrowsTPC>dNegXrowsTPC) dDauXrowsTPC = dNegXrowsTPC;
853   if (dDauXrowsTPC<fCutMinDauXrowsTPC) return 0x0;
854
855   UShort_t wPosTPCNClsF = pDauPosRD->GetTPCNclsF(); if (wPosTPCNClsF<=0) return 0x0;
856   UShort_t wNegTPCNClsF = pDauNegRD->GetTPCNclsF(); if (wNegTPCNClsF<=0) return 0x0;
857   Double_t dPosXrowsOverFindableClusTPC = ((Double_t)dPosXrowsTPC) / ((Double_t)wPosTPCNClsF);
858   Double_t dNegXrowsOverFindableClusTPC = ((Double_t)dNegXrowsTPC) / ((Double_t)wNegTPCNClsF);
859
860   Double_t dDauXrowsOverFindableClusTPC = dPosXrowsOverFindableClusTPC;
861   if (dDauXrowsOverFindableClusTPC>dNegXrowsOverFindableClusTPC) dDauXrowsOverFindableClusTPC = dNegXrowsOverFindableClusTPC;
862   if (dDauXrowsOverFindableClusTPC<fCutMinDauXrowsOverFindableClusTPC) return 0x0;
863 //=============================================================================
864
865   Short_t nPosCharge = pDauPosRD->Charge();
866   Short_t nNegCharge = pDauNegRD->Charge();
867   if ((nPosCharge==0) || (nNegCharge==0) || (nPosCharge==nNegCharge)) return 0x0;
868
869   Double_t dPosPxPyPz[3] = { 0., 0., 0. };
870   Double_t dNegPxPyPz[3] = { 0., 0., 0. };
871   if ((nPosCharge<0) && (nNegCharge>0)) {
872     pDauPosRD = (AliAODTrack*)pV0RD->GetDaughter(1);
873     pDauNegRD = (AliAODTrack*)pV0RD->GetDaughter(0);
874
875     dPosPxPyPz[0] = pV0RD->MomNegX(); dPosPxPyPz[1] = pV0RD->MomNegY(); dPosPxPyPz[2] = pV0RD->MomNegZ();
876     dNegPxPyPz[0] = pV0RD->MomPosX(); dNegPxPyPz[1] = pV0RD->MomPosY(); dNegPxPyPz[2] = pV0RD->MomPosZ();
877   } else {
878     dPosPxPyPz[0] = pV0RD->MomPosX(); dPosPxPyPz[1] = pV0RD->MomPosY(); dPosPxPyPz[2] = pV0RD->MomPosZ();
879     dNegPxPyPz[0] = pV0RD->MomNegX(); dNegPxPyPz[1] = pV0RD->MomNegY(); dNegPxPyPz[2] = pV0RD->MomNegZ();
880   }
881
882   TVector3 v3Pos(dPosPxPyPz);
883   TVector3 v3Neg(dNegPxPyPz);
884
885   if ((v3Pos.Pt()<fCutMinDauPt) || (v3Neg.Pt()<fCutMinDauPt)) return 0x0;
886   Double_t dPosEta = v3Pos.Eta(); if ((dPosEta<fCutMinDauEta) || (dPosEta>fCutMaxDauEta)) return 0x0;
887   Double_t dNegEta = v3Neg.Eta(); if ((dNegEta<fCutMinDauEta) || (dNegEta>fCutMaxDauEta)) return 0x0;
888 //=============================================================================
889
890   Int_t inp = TMath::Abs(pDauPosRD->GetLabel()); if (inp<0) return 0x0;
891   Int_t inn = TMath::Abs(pDauNegRD->GetLabel()); if (inn<0) return 0x0;
892   AliAODMCParticle *pDauPosMC = (AliAODMCParticle*)MCEvent()->GetTrack(inp); if (!pDauPosMC) return 0x0;
893   AliAODMCParticle *pDauNegMC = (AliAODMCParticle*)MCEvent()->GetTrack(inn); if (!pDauNegMC) return 0x0;
894
895   Int_t imp = pDauPosMC->GetMother(); if (imp<0) return 0x0;
896   Int_t imn = pDauNegMC->GetMother(); if (imn<0) return 0x0;
897   if (imp != imn) return 0x0;
898
899   AliAODMCParticle *pV0MC = (AliAODMCParticle*)MCEvent()->GetTrack(imp); if (!pV0MC) return 0x0;
900   if (((pV0MC->Y())<fCutMinV0Rap) || ((pV0MC->Y())>fCutMaxV0Rap)) return 0x0;
901
902   Int_t idvMC = pV0MC->GetPdgCode();
903   Int_t idp = pDauPosMC->GetPdgCode();
904   Int_t idn = pDauNegMC->GetPdgCode();
905   Bool_t bIsKshort = ((idp==211)  && (idn==-211)  && (idvMC== 310));
906   Bool_t bIsLambda = ((idp==2212) && (idn==-211)  && (idvMC== 3122));
907   Bool_t bIsAntiLa = ((idp==211)  && (idn==-2212) && (idvMC==-3122));
908   if (!(bIsKshort || bIsLambda || bIsAntiLa)) return 0x0;
909 //=============================================================================
910
911   UInt_t wsvMC = 0;
912   if (pV0MC->IsPrimary())                wsvMC |= AliPicoHeaderCJ::kPrimary;
913   if (pV0MC->IsPhysicalPrimary())        wsvMC |= AliPicoHeaderCJ::kPhysicalPrimary;
914   if (pV0MC->IsSecondaryFromWeakDecay()) wsvMC |= AliPicoHeaderCJ::kSecondaryFromWeakDecay;
915   if (pV0MC->IsSecondaryFromMaterial())  wsvMC |= AliPicoHeaderCJ::kSecondaryFromMaterial;
916
917   Int_t  idmMC = 0;
918   UInt_t wsmMC = 0;
919   Double_t dMotherPt  = 0.;
920   Double_t dMotherEta = 0.;
921   Double_t dMotherRap = 0.;
922   if (bIsLambda || bIsAntiLa) {
923     Int_t imv = pV0MC->GetMother(); if (imv>=0) {
924       AliAODMCParticle *pMother = (AliAODMCParticle*)MCEvent()->GetTrack(imv);
925
926       if (pMother) {
927         idmMC = pMother->GetPdgCode();
928         if ((bIsLambda && ((idmMC== 3312) || (idmMC== 3322))) ||
929             (bIsAntiLa && ((idmMC==-3312) || (idmMC==-3322)))) {
930           dMotherPt  = pMother->Pt();
931           dMotherEta = pMother->Eta();
932           dMotherRap = pMother->Y();
933
934           if (pMother->IsPrimary())                wsmMC |= AliPicoHeaderCJ::kPrimary;
935           if (pMother->IsPhysicalPrimary())        wsmMC |= AliPicoHeaderCJ::kPhysicalPrimary;
936           if (pMother->IsSecondaryFromWeakDecay()) wsmMC |= AliPicoHeaderCJ::kSecondaryFromWeakDecay;
937           if (pMother->IsSecondaryFromMaterial())  wsmMC |= AliPicoHeaderCJ::kSecondaryFromMaterial;
938         }
939       }
940     }
941   }
942 //=============================================================================
943
944   Double_t dV0CosPA = pV0RD->CosPointingAngle(fPrimaryVtx);
945
946   if (bIsKshort) if (dV0CosPA<fCutMinKshortCosPA) {
947     bIsKshort = kFALSE;
948   }
949
950   if (bIsLambda || bIsAntiLa) if (dV0CosPA<fCutMinLambdaCosPA) {
951     bIsLambda = kFALSE;
952     bIsAntiLa = kFALSE;
953   }
954
955   if (!(bIsKshort || bIsLambda || bIsAntiLa)) return 0x0;
956 //=============================================================================
957
958   Double_t dV0DistToPV = 0.;
959   for (Int_t i=0; i<3; i++) dV0DistToPV +=  ((dV0Vtx[i]-fPrimaryVtx[i]) * (dV0Vtx[i]-fPrimaryVtx[i]));
960   Double_t dV0DistToPVoverP = TMath::Sqrt(dV0DistToPV) / (pV0RD->P()+1e-10);
961
962   if (bIsKshort) if ((dV0DistToPVoverP*fgkMassKshort)>fCutMaxKshortCtau) {
963     bIsKshort = kFALSE;
964   }
965
966   if (bIsLambda || bIsAntiLa) if ((dV0DistToPVoverP*fgkMassLambda)>fCutMaxLambdaCtau) {
967     bIsLambda = kFALSE;
968     bIsAntiLa = kFALSE;
969   }
970
971   if (!(bIsKshort || bIsLambda || bIsAntiLa)) return 0x0;
972 //=============================================================================
973
974   Double_t dV0ArmFrac = pV0RD->PtArmV0() / (TMath::Abs(pV0RD->AlphaV0())+1e-12);
975
976   if (bIsKshort && (fCutMaxKshortArmFrac>0.)) if (dV0ArmFrac>fCutMaxKshortArmFrac) {
977     bIsKshort = kFALSE;
978   }
979
980   if ((bIsLambda || bIsAntiLa) && fCutMaxLambdaArmFrac>0.) if (dV0ArmFrac>fCutMaxLambdaArmFrac) {
981     bIsLambda = kFALSE;
982     bIsAntiLa = kFALSE;
983   }
984
985   if (!(bIsKshort || bIsLambda || bIsAntiLa)) return 0x0;
986 //=============================================================================
987
988   Int_t wMask = 0;
989   if (bIsKshort) {
990     TLorentzVector vPosPion; vPosPion.SetVectM(v3Pos, fgkMassPion);
991     TLorentzVector vNegPion; vNegPion.SetVectM(v3Neg, fgkMassPion);
992     TLorentzVector vKshort = vPosPion + vNegPion;
993
994     Double_t dKshortInvM = vKshort.M();
995     Double_t dLower = 0.430006 - 0.0110029*dV0Pt;
996     Double_t dUpper = 0.563707 + 0.0114979*dV0Pt;
997     if ((dKshortInvM<dLower) || (dKshortInvM>dUpper)) return 0x0;
998
999     if (fCutMinKshortDeltaM>0.) {
1000       TLorentzVector vPosProton; vPosProton.SetVectM(v3Pos, fgkMassProton);
1001       TLorentzVector vNegProton; vNegProton.SetVectM(v3Neg, fgkMassProton);
1002
1003       TLorentzVector vLamvda = vPosProton + vNegPion;
1004       TLorentzVector vAntiLa = vNegProton + vPosPion;
1005
1006       Double_t dLambdaInvM = vLamvda.M();
1007       Double_t dAntiLaInvM = vAntiLa.M();
1008       if ((TMath::Abs(dLambdaInvM-fgkMassLambda)<fCutMinKshortDeltaM) ||
1009           (TMath::Abs(dAntiLaInvM-fgkMassLambda)<fCutMinKshortDeltaM)) return 0x0;
1010     }
1011
1012     wMask = AliPicoHeaderCJ::kKshort;
1013   }
1014
1015   if (bIsLambda) {
1016     TLorentzVector vPosProton; vPosProton.SetVectM(v3Pos, fgkMassProton);
1017     TLorentzVector vNegPion;   vNegPion.SetVectM(v3Neg, fgkMassPion);
1018     TLorentzVector vLamvda = vPosProton + vNegPion;
1019
1020     Double_t dLambdaInvM = vLamvda.M();
1021     Double_t dLower = 1.09501 - 0.00523272*dV0Pt - 0.075269*TMath::Exp(-3.46339*dV0Pt);
1022     Double_t dUpper = 1.13688 + 0.00527838*dV0Pt + 0.084222*TMath::Exp(-3.80595*dV0Pt);
1023     if ((dLambdaInvM<dLower) || (dLambdaInvM>dUpper)) return 0x0;
1024
1025     if (fCutMinLambdaDeletaM>0.) {
1026       TLorentzVector vPosPion; vPosPion.SetVectM(v3Pos, fgkMassPion);
1027       TLorentzVector vKshort = vPosPion + vNegPion;
1028
1029       Double_t dKshortInvM = vKshort.M();
1030       if ((TMath::Abs(dKshortInvM-fgkMassKshort)<fCutMinLambdaDeletaM)) return 0x0;
1031     }
1032
1033     wMask = AliPicoHeaderCJ::kLambda;
1034   }
1035
1036   if (bIsAntiLa) {
1037     TLorentzVector vNegProton; vNegProton.SetVectM(v3Neg, fgkMassProton);
1038     TLorentzVector vPosPion;   vPosPion.SetVectM(v3Pos, fgkMassPion);
1039     TLorentzVector vAntiLa = vNegProton + vPosPion;
1040
1041     Double_t dAntiLaInvM = vAntiLa.M();
1042     Double_t dLower = 1.09501 - 0.00523272*dV0Pt - 0.075269*TMath::Exp(-3.46339*dV0Pt);
1043     Double_t dUpper = 1.13688 + 0.00527838*dV0Pt + 0.084222*TMath::Exp(-3.80595*dV0Pt);
1044     if ((dAntiLaInvM<dLower) || (dAntiLaInvM>dUpper)) return 0x0;
1045
1046     if (fCutMinLambdaDeletaM>0.) {
1047       TLorentzVector vNegPion; vNegPion.SetVectM(v3Neg, fgkMassPion);
1048       TLorentzVector vKshort = vPosPion + vNegPion;
1049
1050       Double_t dKshortInvM = vKshort.M();
1051       if ((TMath::Abs(dKshortInvM-fgkMassKshort)<fCutMinLambdaDeletaM)) return 0x0;
1052     }
1053
1054     wMask = AliPicoHeaderCJ::kAntiLambda;
1055   }
1056 //=============================================================================
1057
1058   Bool_t bPosInJC = kFALSE;
1059   Bool_t bNegInJC = kFALSE;
1060   AliPicoV0MC *pPicoV0 = new AliPicoV0MC(wMask,
1061                                          dV0Radius,
1062                                          dV0CosPA,
1063                                          dV0DistToPVoverP,
1064                                          dDausDCA,
1065                                          dPosDCAtoPV,
1066                                          dNegDCAtoPV,
1067                                          dDauXrowsTPC,
1068                                          dDauXrowsOverFindableClusTPC,
1069                                          v3Pos.Px(), v3Pos.Py(), v3Pos.Pz(),
1070                                          v3Neg.Px(), v3Neg.Py(), v3Neg.Pz(),
1071                                          bPosInJC, bNegInJC,
1072                                          idvMC, wsvMC, pV0MC->Px(), pV0MC->Py(), pV0MC->Pz(), pV0MC->E(),
1073                                          idmMC, wsmMC, dMotherPt, dMotherEta, dMotherRap);
1074
1075
1076   return pPicoV0;
1077 }
1078
1079 //_____________________________________________________________________________
1080 AliPicoV0MC* AliAnalysisTaskSEPicoV0Maker::SelectV0CandidateMC(AliESDv0 const *pV0RD)
1081 {
1082 //
1083 //  AliAnalysisTaskSEPicoV0Maker::SelectV0CandidateMC
1084 //
1085
1086   AliStack *pStack = MCEvent()->Stack(); if (!pStack) return 0x0;
1087   Int_t nPrimary = pStack->GetNprimary();
1088 //=============================================================================
1089
1090   if (pV0RD->GetOnFlyStatus()) return 0x0;
1091   if (pV0RD->GetChi2V0()>fCutMaxV0Chi2) return 0x0;
1092
1093   Double_t dV0Pt = pV0RD->Pt(); if ((dV0Pt<fCutMinV0Pt) || (dV0Pt>fCutMaxV0Pt)) return 0x0;
1094 //=============================================================================
1095   
1096   Double_t dV0Vtx[3];  pV0RD->GetXYZ(dV0Vtx[0], dV0Vtx[1], dV0Vtx[2]);
1097   Double_t dV0Radius = TMath::Sqrt(dV0Vtx[0]*dV0Vtx[0] + dV0Vtx[1]*dV0Vtx[1]);
1098   if ((dV0Radius<fCutMinV0Radius) || (dV0Radius>fCutMaxV0Radius)) return 0x0; 
1099
1100   Double_t dDausDCA = pV0RD->GetDcaV0Daughters(); if (dDausDCA>fCutMaxDausDCA) return 0x0;
1101 //=============================================================================
1102
1103   Int_t nPosIndex = TMath::Abs(pV0RD->GetPindex()); if (nPosIndex<0) return 0x0;
1104   Int_t nNegIndex = TMath::Abs(pV0RD->GetNindex()); if (nNegIndex<0) return 0x0;
1105
1106   AliESDtrack *pDauPosRD = fEventESD->GetTrack(nPosIndex); if (!pDauPosRD) return 0x0;
1107   AliESDtrack *pDauNegRD = fEventESD->GetTrack(nNegIndex); if (!pDauNegRD) return 0x0;
1108
1109   Double_t dMegField = fEventESD->GetMagneticField();
1110   Double_t dPosDCAtoPV = TMath::Abs(pDauPosRD->GetD(fPrimaryVtx[0],fPrimaryVtx[1],dMegField)); if (dPosDCAtoPV<fCutMinDauDCAtoPV) return 0x0;
1111   Double_t dNegDCAtoPV = TMath::Abs(pDauNegRD->GetD(fPrimaryVtx[0],fPrimaryVtx[1],dMegField)); if (dNegDCAtoPV<fCutMinDauDCAtoPV) return 0x0;
1112 //=============================================================================
1113
1114   if (!(pDauPosRD->GetStatus() & AliESDtrack::kTPCrefit)) return 0x0;
1115   if (!(pDauNegRD->GetStatus() & AliESDtrack::kTPCrefit)) return 0x0;
1116   if ((pDauPosRD->GetKinkIndex(0)>0) || (pDauNegRD->GetKinkIndex(0)>0)) return 0x0;
1117
1118   Float_t dPosXrowsTPC = pDauPosRD->GetTPCClusterInfo(2,1);
1119   Float_t dNegXrowsTPC = pDauNegRD->GetTPCClusterInfo(2,1);
1120   Float_t dDauXrowsTPC = dPosXrowsTPC; if (dDauXrowsTPC>dNegXrowsTPC) dDauXrowsTPC = dNegXrowsTPC;
1121   if (dDauXrowsTPC<fCutMinDauXrowsTPC) return 0x0;
1122
1123   UShort_t wPosTPCNClsF = pDauPosRD->GetTPCNclsF(); if (wPosTPCNClsF<=0) return 0x0;
1124   UShort_t wNegTPCNClsF = pDauNegRD->GetTPCNclsF(); if (wNegTPCNClsF<=0) return 0x0;
1125   Double_t dPosXrowsOverFindableClusTPC = ((Double_t)dPosXrowsTPC) / ((Double_t)wPosTPCNClsF);
1126   Double_t dNegXrowsOverFindableClusTPC = ((Double_t)dNegXrowsTPC) / ((Double_t)wNegTPCNClsF);
1127
1128   Double_t dDauXrowsOverFindableClusTPC = dPosXrowsOverFindableClusTPC;
1129   if (dDauXrowsOverFindableClusTPC>dNegXrowsOverFindableClusTPC) dDauXrowsOverFindableClusTPC = dNegXrowsOverFindableClusTPC;
1130   if (dDauXrowsOverFindableClusTPC<fCutMinDauXrowsOverFindableClusTPC) return 0x0;
1131 //=============================================================================
1132
1133   Short_t nPosCharge = pDauPosRD->Charge();
1134   Short_t nNegCharge = pDauNegRD->Charge();
1135   if ((nPosCharge==0) || (nNegCharge==0) || (nPosCharge==nNegCharge)) return 0x0;
1136
1137   Double_t dPosPxPyPz[3] = { 0., 0., 0. };
1138   Double_t dNegPxPyPz[3] = { 0., 0., 0. };
1139   if ((nPosCharge<0) && (nNegCharge>0)) {
1140     pDauPosRD = fEventESD->GetTrack(nNegIndex);
1141     pDauNegRD = fEventESD->GetTrack(nPosIndex);
1142
1143     pV0RD->GetNPxPyPz(dPosPxPyPz[0], dPosPxPyPz[1], dPosPxPyPz[2]);
1144     pV0RD->GetPPxPyPz(dNegPxPyPz[0], dNegPxPyPz[1], dNegPxPyPz[2]);
1145   } else {
1146     pV0RD->GetPPxPyPz(dPosPxPyPz[0], dPosPxPyPz[1], dPosPxPyPz[2]);
1147     pV0RD->GetNPxPyPz(dNegPxPyPz[0], dNegPxPyPz[1], dNegPxPyPz[2]);
1148   }
1149
1150   TVector3 v3Pos(dPosPxPyPz);
1151   TVector3 v3Neg(dNegPxPyPz);
1152
1153   if ((v3Pos.Pt()<fCutMinDauPt) || (v3Neg.Pt()<fCutMinDauPt)) return 0x0;
1154   Double_t dPosEta = v3Pos.Eta(); if ((dPosEta<fCutMinDauEta) || (dPosEta>fCutMaxDauEta)) return 0x0;
1155   Double_t dNegEta = v3Neg.Eta(); if ((dNegEta<fCutMinDauEta) || (dNegEta>fCutMaxDauEta)) return 0x0;
1156 //=============================================================================
1157
1158   Int_t inp = TMath::Abs(pDauPosRD->GetLabel()); if (inp<0) return 0x0;
1159   Int_t inn = TMath::Abs(pDauNegRD->GetLabel()); if (inn<0) return 0x0;
1160   TParticle *pDauPosMC = ((AliMCParticle*)MCEvent()->GetTrack(inp))->Particle(); if (!pDauPosMC) return 0x0;
1161   TParticle *pDauNegMC = ((AliMCParticle*)MCEvent()->GetTrack(inn))->Particle(); if (!pDauNegMC) return 0x0;
1162
1163   Int_t imp = pDauPosMC->GetFirstMother(); if (imp<0) return 0x0;
1164   Int_t imn = pDauNegMC->GetFirstMother(); if (imn<0) return 0x0;
1165   if (imp != imn) return 0x0;
1166
1167   TParticle *pV0MC = ((AliMCParticle*)MCEvent()->GetTrack(imp))->Particle(); if (!pV0MC) return 0x0;
1168   if (((pV0MC->Y())<fCutMinV0Rap) || ((pV0MC->Y())>fCutMaxV0Rap)) return 0x0;
1169
1170   Int_t idvMC = pV0MC->GetPdgCode();
1171   Int_t idp = pDauPosMC->GetPdgCode();
1172   Int_t idn = pDauNegMC->GetPdgCode();
1173   Bool_t bIsKshort = ((idp==211)  && (idn==-211)  && (idvMC== 310));
1174   Bool_t bIsLambda = ((idp==2212) && (idn==-211)  && (idvMC== 3122));
1175   Bool_t bIsAntiLa = ((idp==211)  && (idn==-2212) && (idvMC==-3122));
1176   if (!(bIsKshort || bIsLambda || bIsAntiLa)) return 0x0;
1177 //=============================================================================
1178
1179   UInt_t wsvMC = 0;
1180   if (imp<nPrimary)                          wsvMC |= AliPicoHeaderCJ::kPrimary;
1181   if (pStack->IsPhysicalPrimary(imp))        wsvMC |= AliPicoHeaderCJ::kPhysicalPrimary;
1182   if (pStack->IsSecondaryFromWeakDecay(imp)) wsvMC |= AliPicoHeaderCJ::kSecondaryFromWeakDecay;
1183   if (pStack->IsSecondaryFromMaterial(imp))  wsvMC |= AliPicoHeaderCJ::kSecondaryFromMaterial;
1184
1185   Int_t  idmMC = 0;
1186   UInt_t wsmMC = 0;
1187   Double_t dMotherPt  = 0.;
1188   Double_t dMotherEta = 0.;
1189   Double_t dMotherRap = 0.;
1190   if (bIsLambda || bIsAntiLa) {
1191     Int_t imv = pV0MC->GetFirstMother(); if (imv>=0) {
1192       TParticle *pMother = ((AliMCParticle*)MCEvent()->GetTrack(imv))->Particle();
1193
1194       if (pMother) {
1195         idmMC = pMother->GetPdgCode();
1196         if ((bIsLambda && ((idmMC== 3312) || (idmMC== 3322))) ||
1197             (bIsAntiLa && ((idmMC==-3312) || (idmMC==-3322)))) {
1198           dMotherPt  = pMother->Pt();
1199           dMotherEta = pMother->Eta();
1200           dMotherRap = pMother->Y();
1201
1202           if (imp<nPrimary)                          wsmMC |= AliPicoHeaderCJ::kPrimary;
1203           if (pStack->IsPhysicalPrimary(imv))        wsmMC |= AliPicoHeaderCJ::kPhysicalPrimary;
1204           if (pStack->IsSecondaryFromWeakDecay(imv)) wsmMC |= AliPicoHeaderCJ::kSecondaryFromWeakDecay;
1205           if (pStack->IsSecondaryFromMaterial(imv))  wsmMC |= AliPicoHeaderCJ::kSecondaryFromMaterial;
1206         }
1207       }
1208     }
1209   }
1210 //=============================================================================
1211
1212   Double_t dV0CosPA = pV0RD->GetV0CosineOfPointingAngle(fPrimaryVtx[0], fPrimaryVtx[1], fPrimaryVtx[2]);
1213
1214   if (bIsKshort) if (dV0CosPA<fCutMinKshortCosPA) {
1215     bIsKshort = kFALSE;
1216   }
1217
1218   if (bIsLambda || bIsAntiLa) if (dV0CosPA<fCutMinLambdaCosPA) {
1219     bIsLambda = kFALSE;
1220     bIsAntiLa = kFALSE;
1221   }
1222
1223   if (!(bIsKshort || bIsLambda || bIsAntiLa)) return 0x0;
1224 //=============================================================================
1225
1226   Double_t dV0DistToPV = 0.;
1227   for (Int_t i=0; i<3; i++) dV0DistToPV +=  ((dV0Vtx[i]-fPrimaryVtx[i]) * (dV0Vtx[i]-fPrimaryVtx[i]));
1228   Double_t dV0DistToPVoverP = TMath::Sqrt(dV0DistToPV) / (pV0RD->P()+1e-10);
1229
1230   if (bIsKshort) if ((dV0DistToPVoverP*fgkMassKshort)>fCutMaxKshortCtau) {
1231     bIsKshort = kFALSE;
1232   }
1233
1234   if (bIsLambda || bIsAntiLa) if ((dV0DistToPVoverP*fgkMassLambda)>fCutMaxLambdaCtau) {
1235     bIsLambda = kFALSE;
1236     bIsAntiLa = kFALSE;
1237   }
1238
1239   if (!(bIsKshort || bIsLambda || bIsAntiLa)) return 0x0;
1240 //=============================================================================
1241
1242   Double_t dV0ArmFrac = pV0RD->PtArmV0() / (TMath::Abs(pV0RD->AlphaV0())+1e-12);
1243
1244   if (bIsKshort && (fCutMaxKshortArmFrac>0.)) if (dV0ArmFrac>fCutMaxKshortArmFrac) {
1245     bIsKshort = kFALSE;
1246   }
1247
1248   if ((bIsLambda && bIsAntiLa) && (fCutMaxLambdaArmFrac>0.)) if (dV0ArmFrac>fCutMaxLambdaArmFrac) {
1249     bIsLambda = kFALSE;
1250     bIsAntiLa = kFALSE;
1251   }
1252
1253   if (!(bIsKshort || bIsLambda || bIsAntiLa)) return 0x0;
1254 //=============================================================================
1255
1256   Int_t wMask = 0;
1257
1258   if (bIsKshort) {
1259     TLorentzVector vPosPion; vPosPion.SetVectM(v3Pos, fgkMassPion);
1260     TLorentzVector vNegPion; vNegPion.SetVectM(v3Neg, fgkMassPion);
1261     TLorentzVector vKshort = vPosPion + vNegPion;
1262
1263     Double_t dKshortInvM = vKshort.M();
1264     Double_t dLower = 0.430006 - 0.0110029*dV0Pt;
1265     Double_t dUpper = 0.563707 + 0.0114979*dV0Pt;
1266     if ((dKshortInvM<dLower) || (dKshortInvM>dUpper)) return 0x0;
1267
1268     if (fCutMinKshortDeltaM>0.) {
1269       TLorentzVector vPosProton; vPosProton.SetVectM(v3Pos, fgkMassProton);
1270       TLorentzVector vNegProton; vNegProton.SetVectM(v3Neg, fgkMassProton);
1271
1272       TLorentzVector vLamvda = vPosProton + vNegPion;
1273       TLorentzVector vAntiLa = vNegProton + vPosPion;
1274
1275       Double_t dLambdaInvM = vLamvda.M();
1276       Double_t dAntiLaInvM = vAntiLa.M();
1277       if ((TMath::Abs(dLambdaInvM-fgkMassLambda)<fCutMinKshortDeltaM) ||
1278           (TMath::Abs(dAntiLaInvM-fgkMassLambda)<fCutMinKshortDeltaM)) return 0x0;
1279     }
1280
1281     wMask = AliPicoHeaderCJ::kKshort;
1282   }
1283
1284   if (bIsLambda) {
1285     TLorentzVector vPosProton; vPosProton.SetVectM(v3Pos, fgkMassProton);
1286     TLorentzVector vNegPion; vNegPion.SetVectM(v3Neg, fgkMassPion);
1287     TLorentzVector vLamvda = vPosProton + vNegPion;
1288
1289     Double_t dLambdaInvM = vLamvda.M();
1290     Double_t dLower = 1.09501 - 0.00523272*dV0Pt - 0.075269*TMath::Exp(-3.46339*dV0Pt);
1291     Double_t dUpper = 1.13688 + 0.00527838*dV0Pt + 0.084222*TMath::Exp(-3.80595*dV0Pt);
1292     if ((dLambdaInvM<dLower) || (dLambdaInvM>dUpper)) return 0x0;
1293
1294     if (fCutMinLambdaDeletaM>0.) {
1295       TLorentzVector vPosPion; vPosPion.SetVectM(v3Pos, fgkMassPion);
1296       TLorentzVector vKshort = vPosPion + vNegPion;
1297
1298       Double_t dKshortInvM = vKshort.M();
1299       if ((TMath::Abs(dKshortInvM-fgkMassKshort)<fCutMinLambdaDeletaM)) return 0x0;
1300     }
1301
1302     wMask = AliPicoHeaderCJ::kLambda;
1303   }
1304
1305   if (bIsAntiLa) {
1306     TLorentzVector vNegProton; vNegProton.SetVectM(v3Neg, fgkMassProton);
1307     TLorentzVector vPosPion; vPosPion.SetVectM(v3Pos, fgkMassPion);
1308     TLorentzVector vAntiLa = vNegProton + vPosPion;
1309   
1310     Double_t dAntiLaInvM = vAntiLa.M();
1311     Double_t dLower = 1.09501 - 0.00523272*dV0Pt - 0.075269*TMath::Exp(-3.46339*dV0Pt);
1312     Double_t dUpper = 1.13688 + 0.00527838*dV0Pt + 0.084222*TMath::Exp(-3.80595*dV0Pt);
1313     if ((dAntiLaInvM<dLower) || (dAntiLaInvM>dUpper)) return 0x0;
1314   
1315     if (fCutMinLambdaDeletaM>0.) {
1316       TLorentzVector vNegPion; vNegPion.SetVectM(v3Neg, fgkMassPion);
1317       TLorentzVector vKshort = vPosPion + vNegPion;
1318
1319       Double_t dKshortInvM = vKshort.M();
1320       if ((TMath::Abs(dKshortInvM-fgkMassKshort)<fCutMinLambdaDeletaM)) return 0x0;
1321     }
1322
1323     wMask = AliPicoHeaderCJ::kAntiLambda;
1324   }
1325 //=============================================================================
1326
1327   Bool_t bPosInJC = kFALSE;
1328   Bool_t bNegInJC = kFALSE;
1329   AliPicoV0MC *pPicoV0 = new AliPicoV0MC(wMask,
1330                                          dV0Radius,
1331                                          dV0CosPA,
1332                                          dV0DistToPVoverP,
1333                                          dDausDCA,
1334                                          dPosDCAtoPV,
1335                                          dNegDCAtoPV,
1336                                          dDauXrowsTPC,
1337                                          dDauXrowsOverFindableClusTPC,
1338                                          v3Pos.Px(), v3Pos.Py(), v3Pos.Pz(),
1339                                          v3Neg.Px(), v3Neg.Py(), v3Neg.Pz(),
1340                                          bPosInJC, bNegInJC,
1341                                          idvMC, wsvMC, pV0MC->Px(), pV0MC->Py(), pV0MC->Pz(), pV0MC->Energy(),
1342                                          idmMC, wsmMC, dMotherPt, dMotherEta, dMotherRap);
1343
1344
1345   return pPicoV0;
1346 }
1347
1348 //_____________________________________________________________________________
1349 Bool_t AliAnalysisTaskSEPicoV0Maker::IsEventNotAcpt()
1350 {
1351 //
1352 //  AliAnalysisTaskSEPicoV0Maker::IsEventNotAcpt
1353 //
1354
1355   fEventAcptMask = 0;
1356   if (!InputEvent())  return (fEventAcptMask==0);
1357   if (!fInputHandler) return (fEventAcptMask==0);
1358
1359   if (fCollisionType!=(AliPicoHeaderCJ::kPP)) {
1360     fCentInfo = InputEvent()->GetCentrality();
1361     if (!fCentInfo) return (fEventAcptMask==0);
1362   }
1363
1364   fEventAOD = dynamic_cast<AliAODEvent*>(InputEvent());
1365   fEventESD = dynamic_cast<AliESDEvent*>(InputEvent());
1366   if ((!fEventAOD) && (!fEventESD)) return (fEventAcptMask==0);
1367
1368   fRespoPID = fInputHandler->GetPIDResponse();
1369   if (!fRespoPID) return kTRUE;
1370
1371   if (fIsAnaInfoMC) {
1372      if (MCEvent()) {
1373       if (MCEvent()->GetNumberOfTracks()<=0) return (fEventAcptMask==0);
1374     } else return (fEventAcptMask==0);
1375
1376     AliHeader *pHeader = MCEvent()->Header(); if (!pHeader) return (fEventAcptMask==0);
1377
1378     if (fIsDPMjetMC) {
1379       AliGenDPMjetEventHeader *pDPMjetH = dynamic_cast<AliGenDPMjetEventHeader*>(pHeader->GenEventHeader());
1380
1381       if (pDPMjetH) {
1382         Int_t nd0=0, nd1=0, nd2=0; pDPMjetH->GetNDiffractive(nd1, nd2, nd0);
1383         if ((nd1+nd2) != (pDPMjetH->ProjectileParticipants() + pDPMjetH->TargetParticipants())) return (fEventAcptMask==0);
1384       }
1385     }
1386   }
1387
1388   fEventAcptMask |= AliPicoHeaderCJ::kEventAccCheck;
1389 //=============================================================================
1390
1391   if (fCollisionType==(AliPicoHeaderCJ::kPP)) {
1392     fEventAcptMask |= AliPicoHeaderCJ::kEventAccMult;
1393   } else {
1394     if (fCentInfo->GetQuality()==0)
1395       fEventAcptMask |= AliPicoHeaderCJ::kEventAccMult;
1396     else
1397       return (fEventAcptMask==0);
1398   }
1399 //=============================================================================
1400
1401   UInt_t wMask = fInputHandler->IsEventSelected();
1402   if ((wMask & fTriggerMask) != fTriggerMask) return (fEventAcptMask==0);
1403   if (fIsSkipFastOnly) if ((wMask & AliVEvent::kFastOnly) == AliVEvent::kFastOnly) return (fEventAcptMask==0);
1404
1405   fEventAcptMask |= AliPicoHeaderCJ::kEventAccTrigger;
1406 //=============================================================================
1407
1408   const AliVVertex *pVertex = InputEvent()->GetPrimaryVertex(); if (!pVertex) return (fEventAcptMask==0);
1409   pVertex->GetXYZ(fPrimaryVtx); if (TMath::Abs(fPrimaryVtx[2])>fCutMaxEventVzAbs) return (fEventAcptMask==0);
1410
1411   if ((fCollisionType==(AliPicoHeaderCJ::kPA)) || (fCollisionType==(AliPicoHeaderCJ::kAP))) {
1412     if ( fAnaUtils->IsFirstEventInChunk(InputEvent()))    return (fEventAcptMask==0);
1413     if (!fAnaUtils->IsVertexSelected2013pA(InputEvent())) return (fEventAcptMask==0);
1414
1415 /*  if (fEventAOD) {
1416       const AliAODVertex *pVtxSPD = fEventAOD->GetPrimaryVertexSPD();
1417       const AliAODVertex *pVtxTrk = fEventAOD->GetPrimaryVertex();
1418       if ((!pVtxSPD) && (!pVtxTrk)) return (fEventAcptMask==0);
1419     }*/
1420
1421 /*  if (fEventESD) {
1422       Bool_t fHasVertex = kFALSE;
1423       const AliESDVertex *pVtxESD = fEventESD->GetPrimaryVertexTracks();
1424       if (pVtxESD->GetNContributors()<1) {
1425         pVtxESD = fEventESD->GetPrimaryVertexSPD();
1426         if (pVtxESD->GetNContributors()<1) fHasVertex = kFALSE;
1427         else fHasVertex = kTRUE;
1428
1429         TString vtxTyp = pVtxESD->GetTitle();
1430         Double_t cov[6] = { 0., 0., 0., 0., 0., 0. };
1431         pVtxESD->GetCovarianceMatrix(cov);
1432         Double_t zRes = TMath::Sqrt(cov[5]);
1433         if (vtxTyp.Contains("vertexer:Z") && (zRes>0.25)) fHasVertex = kFALSE;
1434       } else fHasVertex = kTRUE;
1435
1436       if (!fHasVertex) return (fEventAcptMask==0);
1437     }*/
1438
1439   } else {
1440     if (fEventAOD) {
1441       const AliAODVertex *pVtxSPD = fEventAOD->GetPrimaryVertexSPD(); if (!pVtxSPD) return (fEventAcptMask==0);
1442       const AliAODVertex *pVtxTrk = fEventAOD->GetPrimaryVertex();    if (!pVtxTrk) return (fEventAcptMask==0);
1443     }
1444
1445     if (fEventESD) {
1446       const AliESDVertex *pVtxPri = fEventESD->GetPrimaryVertex();       if (!pVtxPri) return (fEventAcptMask==0);
1447       const AliESDVertex *pVtxSPD = fEventESD->GetPrimaryVertexSPD();    if (!pVtxSPD) return (fEventAcptMask==0);
1448       const AliESDVertex *pVtxTrk = fEventESD->GetPrimaryVertexTracks(); if (!pVtxTrk) return (fEventAcptMask==0);
1449       if ((!(pVtxPri->GetStatus())) && (!(pVtxSPD->GetStatus())) && (!(pVtxTrk->GetStatus()))) return (fEventAcptMask==0);
1450     }
1451   }
1452
1453   fEventAcptMask |= AliPicoHeaderCJ::kEventAccVertex;
1454 //=============================================================================
1455
1456   if ((fCollisionType==AliPicoHeaderCJ::kPP) ||
1457       (fCollisionType==AliPicoHeaderCJ::kPA) ||
1458       (fCollisionType==AliPicoHeaderCJ::kAP)) {
1459     if (fAnaUtils->IsPileUpEvent(InputEvent())) return (fEventAcptMask==0);
1460   }
1461
1462   fEventAcptMask |= AliPicoHeaderCJ::kEventAccPileup;
1463 //=============================================================================
1464
1465   if (fIsRefitV0sESD && fEventESD) {
1466     Double_t dCuts[7] = { fCutMaxV0Chi2,
1467                           fCutMinDauDCAtoPV,
1468                           fCutMinDauDCAtoPV,
1469                           fCutMaxDausDCA,
1470                           fCutMinKshortCosPA,
1471                           fCutMinV0Radius,
1472                           fCutMaxV0Radius };
1473
1474     fEventESD->ResetV0s();
1475     AliV0vertexer aV0vtxer;
1476     aV0vtxer.SetDefaultCuts(dCuts);
1477     aV0vtxer.Tracks2V0vertices(fEventESD);
1478   }
1479
1480   return (fEventAcptMask==0);
1481 }
1482
1483 //_____________________________________________________________________________
1484 Bool_t AliAnalysisTaskSEPicoV0Maker::IsEventNotINEL()
1485 {
1486 //
1487 //  AliAnalysisTaskSEPicoV0Maker::IsEventNotINEL
1488 //
1489
1490   if ((fEventAcptMask & AliPicoHeaderCJ::kEventAccCheck) != AliPicoHeaderCJ::kEventAccCheck) return kTRUE;
1491   if ((fEventAcptMask & AliPicoHeaderCJ::kEventAccMult)  != AliPicoHeaderCJ::kEventAccMult)  return kTRUE;
1492
1493   return kFALSE;
1494 }
1495
1496 //_____________________________________________________________________________
1497 Bool_t AliAnalysisTaskSEPicoV0Maker::IsEventNotMBsa()
1498 {
1499 //
1500 //  AliAnalysisTaskSEPicoV0Maker::IsEventNotMBsa
1501 //
1502
1503   if (IsEventNotINEL()) return kTRUE;
1504
1505   if ((fEventAcptMask & AliPicoHeaderCJ::kEventAccTrigger) != AliPicoHeaderCJ::kEventAccTrigger) return kTRUE;
1506   if ((fEventAcptMask & AliPicoHeaderCJ::kEventAccVertex)  != AliPicoHeaderCJ::kEventAccVertex)  return kTRUE;
1507   if ((fEventAcptMask & AliPicoHeaderCJ::kEventAccPileup)  != AliPicoHeaderCJ::kEventAccPileup)  return kTRUE;
1508
1509   return kFALSE;
1510 }
1511
1512 //_____________________________________________________________________________
1513 void AliAnalysisTaskSEPicoV0Maker::FillHistogramsEH()
1514 {
1515 //
1516 //  AliAnalysisTaskSEPicoV0Maker::FillHistogramsEH
1517 //
1518
1519   Float_t dV0M = fCentInfo->GetCentralityPercentile("V0M");
1520   Float_t dV0A = fCentInfo->GetCentralityPercentile("V0A");
1521   Float_t dCL1 = fCentInfo->GetCentralityPercentile("CL1");
1522   Float_t dZNA = fCentInfo->GetCentralityPercentile("ZNA");
1523
1524   if ((fEventAcptMask & AliPicoHeaderCJ::kEventAccCheck) == AliPicoHeaderCJ::kEventAccCheck) {
1525     ((TH1D*)fOutputListEH->FindObject("hEventAccCheck_V0M"))->Fill(dV0M);
1526     ((TH1D*)fOutputListEH->FindObject("hEventAccCheck_V0A"))->Fill(dV0A);
1527     ((TH1D*)fOutputListEH->FindObject("hEventAccCheck_CL1"))->Fill(dCL1);
1528     ((TH1D*)fOutputListEH->FindObject("hEventAccCheck_ZNA"))->Fill(dZNA);
1529   }
1530
1531   if ((fEventAcptMask & AliPicoHeaderCJ::kEventAccMult) == AliPicoHeaderCJ::kEventAccMult) {
1532     ((TH1D*)fOutputListEH->FindObject("hEventAccMult_V0M"))->Fill(dV0M);
1533     ((TH1D*)fOutputListEH->FindObject("hEventAccMult_V0A"))->Fill(dV0A);
1534     ((TH1D*)fOutputListEH->FindObject("hEventAccMult_CL1"))->Fill(dCL1);
1535     ((TH1D*)fOutputListEH->FindObject("hEventAccMult_ZNA"))->Fill(dZNA);
1536   }
1537
1538   if ((fEventAcptMask & AliPicoHeaderCJ::kEventAccTrigger) == AliPicoHeaderCJ::kEventAccTrigger) {
1539     ((TH1D*)fOutputListEH->FindObject("hEventAccTrigger_V0M"))->Fill(dV0M);
1540     ((TH1D*)fOutputListEH->FindObject("hEventAccTrigger_V0A"))->Fill(dV0A);
1541     ((TH1D*)fOutputListEH->FindObject("hEventAccTrigger_CL1"))->Fill(dCL1);
1542     ((TH1D*)fOutputListEH->FindObject("hEventAccTrigger_ZNA"))->Fill(dZNA);
1543   }
1544
1545   if ((fEventAcptMask & AliPicoHeaderCJ::kEventAccVertex) == AliPicoHeaderCJ::kEventAccVertex) {
1546     ((TH1D*)fOutputListEH->FindObject("hEventAccVertex_V0M"))->Fill(dV0M);
1547     ((TH1D*)fOutputListEH->FindObject("hEventAccVertex_V0A"))->Fill(dV0A);
1548     ((TH1D*)fOutputListEH->FindObject("hEventAccVertex_CL1"))->Fill(dCL1);
1549     ((TH1D*)fOutputListEH->FindObject("hEventAccVertex_ZNA"))->Fill(dZNA);
1550   }
1551
1552   if ((fEventAcptMask & AliPicoHeaderCJ::kEventAccPileup) == AliPicoHeaderCJ::kEventAccPileup) {
1553     ((TH1D*)fOutputListEH->FindObject("hEventAccPileup_V0M"))->Fill(dV0M);
1554     ((TH1D*)fOutputListEH->FindObject("hEventAccPileup_V0A"))->Fill(dV0A);
1555     ((TH1D*)fOutputListEH->FindObject("hEventAccPileup_CL1"))->Fill(dCL1);
1556     ((TH1D*)fOutputListEH->FindObject("hEventAccPileup_ZNA"))->Fill(dZNA);
1557   }
1558
1559   return;
1560 }
1561
1562 //_____________________________________________________________________________
1563 void AliAnalysisTaskSEPicoV0Maker::FillHistogramsMC()
1564 {
1565 //
1566 //  AliAnalysisTaskSEPicoV0Maker::FillHistogramsMC
1567 //
1568
1569   Int_t   nPrimary = 0;
1570   AliStack *pStack = 0;
1571
1572   if (fEventESD) {
1573     pStack   = MCEvent()->Stack(); if (!pStack) return;
1574     nPrimary = pStack->GetNprimary();
1575   }
1576 //=============================================================================
1577
1578   Double_t dEvType = -0.5;
1579   if (IsEventNotMBsa()) dEvType = 0.5;
1580
1581   Double_t dV0M = fCentInfo->GetCentralityPercentile("V0M");
1582   Double_t dV0A = fCentInfo->GetCentralityPercentile("V0A");
1583   Double_t dCL1 = fCentInfo->GetCentralityPercentile("CL1");
1584   Double_t dZNA = fCentInfo->GetCentralityPercentile("ZNA");
1585
1586   THnSparseD *hsV0 = dynamic_cast<THnSparseD*>(fOutputListMC->FindObject("hsV0"));
1587   THnSparseD *hsXi = dynamic_cast<THnSparseD*>(fOutputListMC->FindObject("hsXi"));
1588
1589   if (hsV0 == 0 || hsXi == 0) { // Keep Coverity happy
1590     AliFatal("Cannot find hsV0 or hsXi; should not happen");
1591     return;
1592   }  // Keep Coverity happy
1593 //=============================================================================
1594
1595   TParticle        *pESD = 0;
1596   AliAODMCParticle *pAOD = 0;
1597   for (Int_t i=0; i<MCEvent()->GetNumberOfTracks(); i++) {
1598     if (fEventAOD) { pAOD = (AliAODMCParticle*)MCEvent()->GetTrack(i);              if (!pAOD) continue; }
1599     if (fEventESD) { pESD =   ((AliMCParticle*)MCEvent()->GetTrack(i))->Particle(); if (!pESD) continue; }
1600
1601     Bool_t bPri = kFALSE;
1602     if (pAOD) bPri = pAOD->IsPrimary();
1603     if (pESD) bPri = (i<nPrimary);
1604
1605     Bool_t bPhy = kFALSE;
1606     if (pAOD) bPhy =   pAOD->IsPhysicalPrimary();
1607     if (pESD) bPhy = pStack->IsPhysicalPrimary(i);
1608     if ((!bPri) && (!bPhy)) { pAOD=0; pESD=0; continue; }
1609
1610     Int_t id = 0;
1611     if (pAOD) id = pAOD->GetPdgCode();
1612     if (pESD) id = pESD->GetPdgCode();
1613
1614     Bool_t bXi = (bPri && ((id==3312) || (id==-3312)));
1615     Bool_t bV0 = (bPhy && ((id==3122) || (id==-3122) || (id==310)));
1616     if (!(bXi || bV0)) { pAOD=0; pESD=0; continue; }
1617
1618     Double_t  dEta = 0.;
1619     if (pAOD) dEta = pAOD->Eta();
1620     if (pESD) dEta = pESD->Eta();
1621     if ((dEta<-5.) || (dEta>=5.)) { pAOD=0; pESD=0; continue; }
1622
1623     Double_t  dRapLab = 0.;
1624     if (pAOD) dRapLab = pAOD->Y();
1625     if (pESD) dRapLab = pESD->Y();
1626     if ((dRapLab<-5.) || (dRapLab>=5.)) { pAOD=0; pESD=0; continue; }
1627
1628     Double_t dRapCMS = dRapLab + fRapidityShift;
1629     if ((dRapCMS<-5.) || (dRapCMS>=5.)) { pAOD=0; pESD=0; continue; }
1630
1631     Double_t dVar[10];
1632     if (pAOD) dVar[9] = pAOD->Pt();
1633     if (pESD) dVar[9] = pESD->Pt();
1634
1635     dVar[1] = dEvType;
1636     dVar[2] = dV0M;
1637     dVar[3] = dV0A;
1638     dVar[4] = dCL1;
1639     dVar[5] = dZNA;
1640     dVar[6] = dEta;
1641     dVar[7] = dRapLab;
1642     dVar[8] = dRapCMS;
1643
1644     if (bXi) {
1645       if (id== 3312) dVar[0] = -0.5;
1646       if (id==-3312) dVar[0] =  0.5;
1647       hsXi->Fill(dVar);
1648     }
1649
1650     if (bV0) {
1651       if (id== 310 ) dVar[0] = 0.;
1652       if (id== 3122) dVar[0] = 1.;
1653       if (id==-3122) dVar[0] = 2.;
1654       hsV0->Fill(dVar);
1655     }
1656
1657     pAOD = 0;
1658     pESD = 0;
1659   }
1660
1661   return;
1662 }
1663
1664 //_____________________________________________________________________________
1665 void AliAnalysisTaskSEPicoV0Maker::CreateHistogramsEH()
1666 {
1667 //
1668 //  AliAnalysisTaskSEPicoV0Maker::CreateHistogramsEH
1669 //
1670
1671   Bool_t bStatusTmpH = TH1::AddDirectoryStatus();
1672   TH1::AddDirectory(kFALSE);
1673 //=============================================================================
1674
1675   TH1D *h1 = 0;
1676   h1 = new TH1D("hEventAccCheck_V0M", "", 210, -10., 200.); h1->Sumw2(); fOutputListEH->Add(h1);
1677   h1 = new TH1D("hEventAccCheck_V0A", "", 210, -10., 200.); h1->Sumw2(); fOutputListEH->Add(h1);
1678   h1 = new TH1D("hEventAccCheck_CL1", "", 210, -10., 200.); h1->Sumw2(); fOutputListEH->Add(h1);
1679   h1 = new TH1D("hEventAccCheck_ZNA", "", 210, -10., 200.); h1->Sumw2(); fOutputListEH->Add(h1);
1680
1681   h1 = new TH1D("hEventAccMult_V0M", "", 210, -10., 200.); h1->Sumw2(); fOutputListEH->Add(h1);
1682   h1 = new TH1D("hEventAccMult_V0A", "", 210, -10., 200.); h1->Sumw2(); fOutputListEH->Add(h1);
1683   h1 = new TH1D("hEventAccMult_CL1", "", 210, -10., 200.); h1->Sumw2(); fOutputListEH->Add(h1);
1684   h1 = new TH1D("hEventAccMult_ZNA", "", 210, -10., 200.); h1->Sumw2(); fOutputListEH->Add(h1);
1685
1686   h1 = new TH1D("hEventAccTrigger_V0M", "", 210, -10., 200.); h1->Sumw2(); fOutputListEH->Add(h1);
1687   h1 = new TH1D("hEventAccTrigger_V0A", "", 210, -10., 200.); h1->Sumw2(); fOutputListEH->Add(h1);
1688   h1 = new TH1D("hEventAccTrigger_CL1", "", 210, -10., 200.); h1->Sumw2(); fOutputListEH->Add(h1);
1689   h1 = new TH1D("hEventAccTrigger_ZNA", "", 210, -10., 200.); h1->Sumw2(); fOutputListEH->Add(h1);
1690
1691   h1 = new TH1D("hEventAccVertex_V0M", "", 210, -10., 200.); h1->Sumw2(); fOutputListEH->Add(h1);
1692   h1 = new TH1D("hEventAccVertex_V0A", "", 210, -10., 200.); h1->Sumw2(); fOutputListEH->Add(h1);
1693   h1 = new TH1D("hEventAccVertex_CL1", "", 210, -10., 200.); h1->Sumw2(); fOutputListEH->Add(h1);
1694   h1 = new TH1D("hEventAccVertex_ZNA", "", 210, -10., 200.); h1->Sumw2(); fOutputListEH->Add(h1);
1695
1696   h1 = new TH1D("hEventAccPileup_V0M", "", 210, -10., 200.); h1->Sumw2(); fOutputListEH->Add(h1);
1697   h1 = new TH1D("hEventAccPileup_V0A", "", 210, -10., 200.); h1->Sumw2(); fOutputListEH->Add(h1);
1698   h1 = new TH1D("hEventAccPileup_CL1", "", 210, -10., 200.); h1->Sumw2(); fOutputListEH->Add(h1);
1699   h1 = new TH1D("hEventAccPileup_ZNA", "", 210, -10., 200.); h1->Sumw2(); fOutputListEH->Add(h1);
1700
1701
1702   TH2D *h2 = 0;
1703   h2 = new TH2D("hKshortPtInvM", "", 1000, 0., 100., 300, fgkMassKshort-0.15, fgkMassKshort+0.15);
1704   h2->Sumw2(); fOutputListEH->Add(h2); h2=0;
1705
1706   h2 = new TH2D("hLambdaPtInvM", "", 1000, 0., 100., 200, fgkMassLambda-0.10, fgkMassLambda+0.10);
1707   h2->Sumw2(); fOutputListEH->Add(h2); h2=0;
1708
1709   h2 = new TH2D("hAntiLaPtInvM", "", 1000, 0., 100., 200, fgkMassLambda-0.10, fgkMassLambda+0.10);
1710   h2->Sumw2(); fOutputListEH->Add(h2); h2=0;
1711
1712   TH1::AddDirectory(bStatusTmpH);
1713   return;
1714 }
1715
1716 //_____________________________________________________________________________
1717 void AliAnalysisTaskSEPicoV0Maker::CreateHistogramsMC()
1718 {
1719 //
1720 //  AliAnalysisTaskSEPicoV0Maker::CreateHistogramsMC
1721 //
1722
1723 // TODO: add the multiplicity bins
1724
1725   Bool_t bStatusTmpH = TH1::AddDirectoryStatus();
1726   TH1::AddDirectory(kFALSE);
1727
1728   const Int_t nV0 = 10; // 0: particle type
1729                         //    ==0, Kshort
1730                         //    ==1, Lambda
1731                         //    ==2, AntiLa
1732                         // 1: Event type
1733                         //    ==-0.5, INEL
1734                         //    == 0.5, MB
1735                         // 2: V0M
1736                         // 3: V0A
1737                         // 4: CL1
1738                         // 5: ZNA
1739                         // 6: eta
1740                         // 7: rap in Lab
1741                         // 8: rap in CMS
1742                         // 9: Pt
1743   const Int_t    nV0Bin[nV0] = {  3,    2,  210,  210,  210,  210, 100, 100, 100, 1000  };
1744   const Double_t dV0Min[nV0] = { -0.5, -1., -10., -10., -10., -10., -5., -5., -5.,   0. };
1745   const Double_t dV0Max[nV0] = {  2.5,  1., 200., 200., 200., 200.,  5.,  5.,  5., 100. };
1746   THnSparseD *hsV0 = new THnSparseD("hsV0", "", nV0, nV0Bin, dV0Min, dV0Max); fOutputListMC->Add(hsV0);
1747 //=============================================================================
1748
1749   const Int_t nXi = 10; // 0: particle type
1750                         //    ==-0.5, XiNeg
1751                         //    == 0.5, XiPos
1752                         // 1: Event type
1753                         //    ==-0.5, INEL
1754                         //    == 0.5, MB
1755                         // 2: V0M
1756                         // 3: V0A
1757                         // 4: CL1
1758                         // 5: ZNA
1759                         // 6: eta
1760                         // 7: rap in Lab
1761                         // 8: rap in CMS
1762                         // 9: Pt
1763   const Int_t    nXiBin[nV0] = {  2,   2,  210,  210,  210,  210, 100, 100, 100, 1000  };
1764   const Double_t dXiMin[nV0] = { -1., -1., -10., -10., -10., -10., -5., -5., -5.,   0. };
1765   const Double_t dXiMax[nV0] = {  1.,  1., 200., 200., 200., 200.,  5.,  5.,  5., 100. };
1766   THnSparseD *hsXi = new THnSparseD("hsXi", "", nXi, nXiBin, dXiMin, dXiMax); fOutputListMC->Add(hsXi);
1767
1768   TH1::AddDirectory(bStatusTmpH);
1769   return;
1770 }
1771
1772 //_____________________________________________________________________________
1773 void AliAnalysisTaskSEPicoV0Maker::InitAnalysis()
1774 {
1775 //
1776 //  AliAnalysisTaskSEPicoV0Maker::InitAnalysis
1777 //
1778
1779   if (fCollisionType==(AliPicoHeaderCJ::kPP)) InitParamsPP();
1780   if (fCollisionType==(AliPicoHeaderCJ::kPA)) InitParamsPA();
1781   if (fCollisionType==(AliPicoHeaderCJ::kAP)) InitParamsAP();
1782   if (fCollisionType==(AliPicoHeaderCJ::kAA)) InitParamsAA();
1783
1784   fAnaUtils = new AliAnalysisUtils();
1785   fAnaUtils->SetMinVtxContr(fCutMinEventVtxContr);
1786   fAnaUtils->SetMaxVtxZ(fCutMaxEventVzAbs);
1787
1788   return;
1789 }
1790
1791 //_____________________________________________________________________________
1792 void AliAnalysisTaskSEPicoV0Maker::InitParamsPP()
1793 {
1794 //
1795 //  AliAnalysisTaskSEPicoV0Maker::InitParametersPP
1796 //
1797
1798   fRapidityShift = 0.;
1799
1800   fCutMaxEventVzAbs = 10.;
1801 //=============================================================================
1802
1803   fCutMaxV0Chi2   = 33.;
1804   fCutMinV0Radius = 0.3;   // default: >0.5; uncertainty: >0.3, 0.4, 0.6, 0.7;
1805   fCutMaxV0Radius = 200.;  // default: not applied in pp
1806 //=============================================================================
1807
1808   fCutMaxDausDCA                     = 1.5;  // default: <1;    uncertainty: <0.5, 0.75, 1.25, 1.5
1809   fCutMinDauDCAtoPV                  = 0.05; // default: >0.06; uncertainty: >0.05, 0.055, 0.07, 0.08
1810   fCutMinDauXrowsTPC                 = 70.;  // default: >70;   uncertainty: >75, 80
1811   fCutMinDauXrowsOverFindableClusTPC = 0.8;  // default: >0.8;  uncertainty: >0.95
1812 //=============================================================================
1813
1814   fCutMaxKshortSigmaTPC = -1.;   // default: <5;     uncertainty: w/o cut
1815   fCutMinKshortCosPA    = 0.95;  // default: >0.97;  uncertainty: >0.95, 0.96, 0.98, 0.99
1816   fCutMaxKshortCtau     = 30.;   // default: <20;    uncertainty: <12, 30
1817   fCutMaxKshortArmFrac  = -1.;   // default: not applied in pp
1818   fCutMinKshortDeltaM   = 0.003; // default: >0.005; uncertainty: >0.003, 0.006
1819 //=============================================================================
1820
1821   fCutMaxLambdaSigmaTPC  = 7.;    // default: <5;     uncertainty: 4, 6, 7
1822   fCutMinLambdaCosPA     = 0.993; // default: >0.995; uncertainty: >0.993, 0.994, 0.996, 0.997
1823   fCutMaxLambdaCtau      = 40.;   // default: <30;    uncertainty: <20, 40
1824   fCutMaxLambdaArmFrac   = -1.;   // default: not applied in pp
1825   fCutMinLambdaDeletaM   = -1.;   // default: >0.01;  uncertainty: w/o rejection
1826
1827   return;
1828 }
1829
1830 //_____________________________________________________________________________
1831 void AliAnalysisTaskSEPicoV0Maker::InitParamsPA()
1832 {
1833 //
1834 //  AliAnalysisTaskSEPicoV0Maker::InitParametersPA
1835 //
1836
1837   InitParamsPP();
1838
1839   fRapidityShift = 0.465;
1840
1841   fCutMaxLambdaSigmaTPC = 6;  // default: <5; uncertaity: <4, 6
1842
1843   return;
1844 }
1845
1846 //_____________________________________________________________________________
1847 void AliAnalysisTaskSEPicoV0Maker::InitParamsAP()
1848 {
1849 //
1850 //  AliAnalysisTaskSEPicoV0Maker::InitParametersAP
1851 //
1852
1853   InitParamsPP();
1854
1855   fRapidityShift = -0.465;
1856
1857   fCutMaxLambdaSigmaTPC = 6;  // default: <5; uncertaity: <4, 6
1858
1859   return;
1860 }
1861
1862 //_____________________________________________________________________________
1863 void AliAnalysisTaskSEPicoV0Maker::InitParamsAA()
1864 {
1865 //
1866 //  AliAnalysisTaskSEPicoV0Maker::InitParametersAA
1867 //
1868
1869   InitParamsPP();
1870
1871   fCutMinV0Radius = 0.9;  // default: 5; uncertainty: varying around
1872   fCutMaxV0Radius = 100.;
1873
1874 //fCutMinDauPt      = 0.16;  // default: 0.16
1875   fCutMinDauDCAtoPV = 0.08;  // default: 0.1; uncertainty: 0.08, 0.12
1876
1877   fCutMaxKshortSigmaTPC = 6;      // default: <5;     uncertainty: <4, 6;
1878   fCutMinKshortCosPA    = 0.997;  // default: >0.998; uncertainty: 0.997, 0.999
1879   fCutMaxKshortArmFrac  = 0.2;    // default: <0.2
1880
1881   fCutMaxLambdaSigmaTPC = 6;      // default: <5;     uncertaity: <4, 6
1882   fCutMinLambdaCosPA    = 0.997;  // default: >0.998; uncertainty: 0.997, 0.999
1883   fCutMinLambdaDeletaM  = 0.008;  // default: >0.01;  uncertainty: 0.008, 0.012
1884
1885   return;
1886 }