]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/Correlations/DPhi/AliLeadingV0Correlation.cxx
Missing using std::... declarations (Savannah bug 102892)
[u/mrichter/AliRoot.git] / PWGCF / Correlations / DPhi / AliLeadingV0Correlation.cxx
1 /* Leading Charged Track+V0 Correlation.(Works for Real,Monte Carlo Data)
2  *                            Sandun Jayarathna
3  *                          University of Houston
4  *                      sandun.pahula.hewage@cern.ch
5  *****************************************************************************************/
6 #include <TROOT.h>
7 #include <TList.h>
8 #include <TChain.h>
9 #include <TFile.h>
10 #include <TMath.h>
11 #include <TTree.h>
12 #include <TRandom.h>
13 #include <THnSparse.h>
14 #include <TPDGCode.h>
15 #include <TDatabasePDG.h>
16
17 #include "AliLog.h"
18 #include "AliAnalysisManager.h"
19 #include "AliAODTrack.h"
20 #include "AliAODEvent.h"
21 #include "AliAODv0.h"
22 #include "AliAODVertex.h"
23 #include "AliAODPid.h"
24 #include "AliPIDResponse.h"
25 #include "AliEventPoolManager.h"
26 #include "AliCentrality.h"
27 #include "AliAODHandler.h"
28 #include "AliAODInputHandler.h"
29 #include "AliAODMCParticle.h"
30 #include "AliInputEventHandler.h"
31 #include "AliVParticle.h"
32 #include "AliMultiplicity.h"
33 #include "AliAODMCHeader.h"
34 #include "AliPID.h"
35 #include "AliExternalTrackParam.h"
36 #include "AliAnalyseLeadingTrackUE.h"
37
38 #include "AliLeadingV0Correlation.h"
39
40 #define CorrBinsX 240
41 #define CorrBinsY 260
42 using std::endl;
43 using std::cout;
44
45
46 Double_t PI =TMath::Pi();
47
48 ClassImp(AliLeadingV0Correlation)
49 ClassImp(V0Correlationparticle)
50
51 //---------------------------------------------------------------------------------------
52 AliLeadingV0Correlation::AliLeadingV0Correlation()
53    : AliAnalysisTaskSE(),
54         fAODEvent(0x0),
55         fPoolMgr(0x0),
56         fPIDResponse(0x0),
57         fAnalyseUE(0x0),
58         fPoolMaxNEvents(0), 
59         fPoolMinNTracks(0), 
60         fMinEventsToMix(0),
61         fNzVtxBins(0), 
62         fNCentBins(0),
63         fcollidingSys(""),
64         fpvzcut(0),
65         fTrackEtaCut(0),
66         fFilterBit(128),
67         fAnalysisMC(0),
68         fCase(0),
69         fRemovePileUP(0),
70     fRemoveAutoCorr(0),
71         fRapidityCut(0),
72         fV0radius(0),
73         fV0PostoPVz(0),
74         fV0NegtoPVz(0),
75         fDCAV0Daughters(0),
76         fCPAK0(0),
77         fCPALam(0),
78         fRejectLamK0(0),
79         fRejectK0Lam(0),
80     fSigmaPID(0),
81         fCutCTK0(0),
82         fCutCTLa(0),
83         fMassCutK0(0),
84         fMassCutLa(0),
85         fUseChargeHadrons(kTRUE), 
86         fPtMin(0.15),
87         fOutputList(0),
88         fHistEventViceGen(0),
89         fHistEventViceReconst(0),
90         fHistMCGenK0(0),
91         fHistMCGenLAM(0),
92         fHistMCGenALAM(0),
93         fHistReconstK0(0),
94         fHistReconstLA(0),
95         fHistReconstALA(0),
96         fHistMCAssoK0(0),
97         fHistMCAssoLA(0),
98         fHistMCAssoALA(0),
99         fHistReconstSib(0),
100         fHistReconstMix(0),
101         fHistReconstSibMC(0),
102         fHistReconstMixMC(0),
103         fHistReconstSibMCAssoc(0),
104         fHistReconstMixMCAssoc(0),
105         fHistTriggerSib(0),
106         fHistTriggerMix(0),
107         fHistTriggerSibMC(0),
108         fHistTriggerMixMC(0)
109 {
110
111   for(Int_t iBin = 0; iBin < 100; iBin++){
112     fZvtxBins[iBin] = 0.;
113     fCentBins[iBin] = 0.;
114   }
115 }
116 //---------------------------------------------------------------------------------------
117 AliLeadingV0Correlation::AliLeadingV0Correlation(const char *name)
118    : AliAnalysisTaskSE(name),
119         fAODEvent(0x0),
120         fPoolMgr(0x0),
121     fPIDResponse(0x0),
122         fAnalyseUE(0x0),
123         fPoolMaxNEvents(0), 
124         fPoolMinNTracks(0), 
125         fMinEventsToMix(0),
126         fNzVtxBins(0), 
127         fNCentBins(0),
128         fcollidingSys(""),
129         fpvzcut(0),
130     fTrackEtaCut(0),
131     fFilterBit(128),
132         fAnalysisMC(0),
133         fCase(0),
134         fRemovePileUP(0),
135         fRemoveAutoCorr(0),
136     fRapidityCut(0),
137         fV0radius(0),
138         fV0PostoPVz(0),
139         fV0NegtoPVz(0),
140         fDCAV0Daughters(0),
141         fCPAK0(0),
142         fCPALam(0),
143         fRejectLamK0(0),
144         fRejectK0Lam(0),
145     fSigmaPID(0),
146         fCutCTK0(0),
147         fCutCTLa(0),
148         fMassCutK0(0),
149         fMassCutLa(0),
150         fUseChargeHadrons(kTRUE), 
151         fPtMin(0.15),
152         fOutputList(0),
153         fHistEventViceGen(0),
154         fHistEventViceReconst(0),
155         fHistMCGenK0(0),
156         fHistMCGenLAM(0),
157         fHistMCGenALAM(0),
158         fHistReconstK0(0),
159         fHistReconstLA(0),
160         fHistReconstALA(0),
161         fHistMCAssoK0(0),
162         fHistMCAssoLA(0),
163         fHistMCAssoALA(0),
164         fHistReconstSib(0),
165         fHistReconstMix(0),
166         fHistReconstSibMC(0),
167         fHistReconstMixMC(0),
168         fHistReconstSibMCAssoc(0),
169         fHistReconstMixMCAssoc(0),
170         fHistTriggerSib(0),
171         fHistTriggerMix(0),
172         fHistTriggerSibMC(0),
173         fHistTriggerMixMC(0)
174
175 {       
176   for(Int_t iBin = 0; iBin < 100; iBin++){
177     fZvtxBins[iBin] = 0.;
178     fCentBins[iBin] = 0.;
179   }
180   DefineOutput(1, TList::Class());                                            
181 }
182
183 //---------------------------------------------------------------------------------------
184 AliLeadingV0Correlation::~AliLeadingV0Correlation()
185 {
186    if (fOutputList && !AliAnalysisManager::GetAnalysisManager()->IsProofMode()) {
187       delete fOutputList;
188    }
189 }
190 //---------------------------------------------------------------------------------------
191 void AliLeadingV0Correlation::UserCreateOutputObjects()
192 {       
193         fAnalyseUE =new AliAnalyseLeadingTrackUE();
194         fAnalyseUE->SetParticleSelectionCriteria(fFilterBit,fUseChargeHadrons,fTrackEtaCut,fPtMin);
195         fAnalyseUE->DefineESDCuts(fFilterBit);
196         
197         fOutputList = new TList();
198         fOutputList->SetOwner();
199         
200         //---------------------------------------------- Events histograms -----------------------------------------------------//
201         //0-PVx,1-PVy,2-PVz,3-MULT,4-CENT
202         const Int_t ndimsEV = 5;
203         Int_t    binsEV[ndimsEV] = { 100,  100, 200, 1000,100};
204         Double_t xminEV[ndimsEV] = {-0.5, -0.5,-20 ,    0,  0};
205         Double_t xmaxEV[ndimsEV] = { 0.5,  0.5, 20 ,20000,100};
206         
207         fHistEventViceGen= new THnSparseD("fHistEventViceGen", "fHistEventViceGen", ndimsEV, binsEV, xminEV, xmaxEV);
208         fOutputList->Add(fHistEventViceGen);
209         
210         fHistEventViceReconst= new THnSparseD("fHistEventViceReconst", "fHistEventViceReconst", ndimsEV, binsEV, xminEV, xmaxEV);
211         fOutputList->Add(fHistEventViceReconst);
212         
213         //0-YK0,1-Pt
214         const Int_t ndimsGenMC = 4;
215         Int_t    binsGenMCLA[ndimsGenMC] = {160,240, 140,100};
216         Double_t xminGenMCLA[ndimsGenMC] = { -4,  0,1.06,  0};
217         Double_t xmaxGenMCLA[ndimsGenMC] = {  4, 12, 1.2,100};  
218         
219         Int_t    binsGenMCK0[ndimsGenMC] = {160,240, 200,100};
220         Double_t xminGenMCK0[ndimsGenMC] = { -4,  0, 0.4,  0};
221         Double_t xmaxGenMCK0[ndimsGenMC] = {  4, 12, 0.6,100};
222         
223         fHistMCGenLAM  = new THnSparseD("fHistMCGenLAM" , "fHistMCGenLAM" , ndimsGenMC, binsGenMCLA, xminGenMCLA, xmaxGenMCLA);
224         fOutputList->Add(fHistMCGenLAM);
225         
226         fHistMCGenALAM = new THnSparseD("fHistMCGenALAM", "fHistMCGenALAM", ndimsGenMC, binsGenMCLA, xminGenMCLA, xmaxGenMCLA);
227         fOutputList->Add(fHistMCGenALAM);
228         
229         fHistMCGenK0   = new THnSparseD("fHistMCGenK0"  , "fHistMCGenK0"  , ndimsGenMC, binsGenMCK0, xminGenMCK0, xmaxGenMCK0);
230         fOutputList->Add(fHistMCGenK0);
231         
232         const Int_t ndims=10;    //MK0  mLA  MALA DCAP DCAN  RV0 DVD   CPA  PT   cent
233         Int_t    bins[ndims] = {  200,  140,  140, 500, 500, 110,110,   200,240  ,100};
234         Double_t xmin[ndims] = {  0.4, 1.06, 1.06,   0,   0,   0,  0, 0.997,  0  ,  0};
235         Double_t xmax[ndims] = {  0.6,  1.2,  1.2,  10,  10, 110,1.1, 1.007, 12  ,100}; 
236         
237         
238         fHistReconstK0= new THnSparseD("fHistReconstK0"  , "fHistReconstK0",  ndims, bins, xmin, xmax);
239         fHistReconstK0->Sumw2();
240         fOutputList->Add(fHistReconstK0);
241         
242         fHistReconstLA= new THnSparseD("fHistReconstLA"  , "fHistReconstLA",  ndims, bins, xmin, xmax);
243         fHistReconstLA->Sumw2();
244         fOutputList->Add(fHistReconstLA);
245         
246         fHistReconstALA= new THnSparseD("fHistReconstALA", "fHistReconstALA", ndims, bins, xmin, xmax);
247         fHistReconstALA->Sumw2();
248         fOutputList->Add(fHistReconstALA);
249         
250         fHistMCAssoK0= new THnSparseD("fHistMCAssoK0"   , "fHistMCAssoK0"   , ndims, bins, xmin, xmax);
251         fHistMCAssoK0->Sumw2();
252         fOutputList->Add(fHistMCAssoK0);
253         
254         fHistMCAssoLA= new THnSparseD("fHistMCAssoLA"   , "fHistMCAssoLA"   , ndims, bins, xmin, xmax);
255         fHistMCAssoLA->Sumw2();
256         fOutputList->Add(fHistMCAssoLA);
257         
258         fHistMCAssoALA= new THnSparseD("fHistMCAssoALA" , "fHistMCAssoALA" ,  ndims, bins, xmin, xmax);
259         fHistMCAssoALA->Sumw2();
260         fOutputList->Add(fHistMCAssoALA);
261         
262         //--------------------------------------------Correlation Histos -----------------------------------------------------//
263         
264         //0-pTK0,1-PhiK0,2-EtaK0,3-DPhiK0,4-DEtaK0,5-TYPE,6-CutSet
265         const Int_t ndimsv0CORR = 7;       
266         Int_t    binsv0CORR[ndimsv0CORR] = {240, 200,          200,CorrBinsX,      CorrBinsY,4,100};
267         
268         Double_t xminv0CORR[ndimsv0CORR] = {  0,   0,-fTrackEtaCut,    -PI/2,-2*fTrackEtaCut,0,0};
269         
270         Double_t xmaxv0CORR[ndimsv0CORR] = { 12,2*PI, fTrackEtaCut,   3*PI/2, 2*fTrackEtaCut,4,100};
271         
272         fHistReconstSib= new THnSparseD("fHistReconstSib", "fHistReconstSib", ndimsv0CORR, binsv0CORR, xminv0CORR, xmaxv0CORR);
273         fHistReconstSib->Sumw2();
274         fOutputList->Add(fHistReconstSib);
275         
276         fHistReconstMix= new THnSparseD("fHistReconstMix", "fHistReconstMix", ndimsv0CORR, binsv0CORR, xminv0CORR, xmaxv0CORR);
277         fHistReconstMix->Sumw2();
278         fOutputList->Add(fHistReconstMix);
279         
280         fHistReconstSibMC= new THnSparseD("fHistReconstSibMC", "fHistReconstSibMC", ndimsv0CORR, binsv0CORR, xminv0CORR, xmaxv0CORR);
281         fHistReconstSibMC->Sumw2();
282         fOutputList->Add(fHistReconstSibMC);
283         
284         fHistReconstMixMC= new THnSparseD("fHistReconstMixMC", "fHistReconstMixMC", ndimsv0CORR, binsv0CORR, xminv0CORR, xmaxv0CORR);
285         fHistReconstMixMC->Sumw2();
286         fOutputList->Add(fHistReconstMixMC);
287         
288         fHistReconstSibMCAssoc= new THnSparseD("fHistReconstSibMCAssoc", "fHistReconstSibMCAssoc", ndimsv0CORR, binsv0CORR, xminv0CORR, xmaxv0CORR);
289         fHistReconstSibMCAssoc->Sumw2();
290         fOutputList->Add(fHistReconstSibMCAssoc);
291         
292         fHistReconstMixMCAssoc= new THnSparseD("fHistReconstMixMCAssoc", "fHistReconstMixMCAssoc", ndimsv0CORR, binsv0CORR, xminv0CORR, xmaxv0CORR);
293         fHistReconstMixMCAssoc->Sumw2();
294         fOutputList->Add(fHistReconstMixMCAssoc);
295         
296         //0-pt,1-PHI,2-Eta
297         const Int_t triggerdims       =4;
298         Int_t binsTrig[triggerdims]   ={240, 200,          200,100};
299         Double_t xminTrig[triggerdims]={  6,   0,-fTrackEtaCut,  0};
300         Double_t xmaxTrig[triggerdims]={ 12,2*PI, fTrackEtaCut,100};
301         
302         fHistTriggerSib= new THnSparseD("fHistTriggerSib", "fHistTriggerSib", triggerdims, binsTrig, xminTrig, xmaxTrig);
303         fHistTriggerSib->Sumw2();
304         fOutputList->Add(fHistTriggerSib);
305         
306         fHistTriggerMix= new THnSparseD("fHistTriggerMix", "fHistTriggerMix", triggerdims, binsTrig, xminTrig, xmaxTrig);
307         fHistTriggerMix->Sumw2();
308         fOutputList->Add(fHistTriggerMix);
309         
310         fHistTriggerSibMC= new THnSparseD("fHistTriggerSibMC", "fHistTriggerSibMC", triggerdims, binsTrig, xminTrig, xmaxTrig);
311         fHistTriggerSibMC->Sumw2();
312         fOutputList->Add(fHistTriggerSibMC);
313         
314         fHistTriggerMixMC= new THnSparseD("fHistTriggerMixMC", "fHistTriggerMixMC", triggerdims, binsTrig, xminTrig, xmaxTrig);
315         fHistTriggerMixMC->Sumw2();
316         fOutputList->Add(fHistTriggerMixMC);
317         
318         //----------------------------------------------Event Pool-----------------------------------------------------//
319         fPoolMgr = new AliEventPoolManager(fPoolMaxNEvents, fPoolMinNTracks, fNCentBins, fCentBins, fNzVtxBins, fZvtxBins);
320         if(!fPoolMgr) return;
321         
322         PostData(1, fOutputList);
323 }
324 //---------------------------------------------------------------------------------------
325 void AliLeadingV0Correlation::UserExec(Option_t *)
326 {
327         
328     AliAnalysisManager   *mgr      = AliAnalysisManager::GetAnalysisManager();
329     AliInputEventHandler *inEvMain = (AliInputEventHandler*)(mgr->GetInputEventHandler());
330         if (!inEvMain) return;
331         
332         // Pointers to PID Response objects.    
333         fPIDResponse = inEvMain->GetPIDResponse();
334         if(!fPIDResponse) return;
335         
336     fAODEvent = dynamic_cast<AliAODEvent*>(inEvMain->GetEvent());
337         if(!fAODEvent) return;
338         
339         Int_t multiplicity    = -1;
340         Int_t multiplicityMC  = -1;
341         Double_t MultipOrCent = -1; 
342         Double_t CentPecentMC = -1;
343         Double_t CentPecentAfterPhySel   = -1;
344         
345         if(fcollidingSys=="PbPb2010" || fcollidingSys=="PbPb2011"){   
346         AliCentrality *centralityObjMC = fAODEvent->GetHeader()->GetCentralityP();
347                 CentPecentMC  = centralityObjMC->GetCentralityPercentileUnchecked("V0M");
348                 if ((CentPecentMC < 0.)||(CentPecentMC > 90)) return;
349         }
350         
351         Double_t * CentBins = fCentBins;
352         Double_t poolmin    = CentBins[0];
353         Double_t poolmax    = CentBins[fNCentBins];
354         
355         
356         Double_t  dimEventviceMC[5];
357         if(fAnalysisMC)    //Efficency denomenator comes before the physics selection
358         {
359                 AliAODMCHeader *aodMCheader = (AliAODMCHeader*)fAODEvent->FindListObject(AliAODMCHeader::StdBranchName());
360                 Float_t mcZv = aodMCheader->GetVtxZ();
361                 
362                 if (TMath::Abs(mcZv) >= fpvzcut) return;
363                 
364                 dimEventviceMC[0]=aodMCheader->GetVtxX();
365                 dimEventviceMC[1]=aodMCheader->GetVtxY();
366                 dimEventviceMC[2]=aodMCheader->GetVtxZ();
367                 
368                 TClonesArray *mcArray = (TClonesArray*)fAODEvent->FindListObject(AliAODMCParticle::StdBranchName());
369                 if(!mcArray)return;
370                 
371                 TObjArray *selectedTracksLeadingMC=fAnalyseUE->FindLeadingObjects(mcArray);
372                 if(!selectedTracksLeadingMC) return;
373                 selectedTracksLeadingMC->SetOwner(kTRUE);
374                 
375                 TObjArray * selectedV0sMC =new TObjArray;
376                 selectedV0sMC->SetOwner(kTRUE);
377                 
378                 Int_t nMCTracks = mcArray->GetEntriesFast();
379                 
380                 if(fcollidingSys=="PbPb2010" || fcollidingSys=="PbPb2011") multiplicityMC=CentPecentMC;
381                 if(fcollidingSys=="PP")   multiplicityMC=nMCTracks;
382                 
383                 dimEventviceMC[3]=nMCTracks;
384                 dimEventviceMC[4]=CentPecentMC;
385                 fHistEventViceGen->Fill(dimEventviceMC);
386                 
387                 for (Int_t iMC = 0; iMC<nMCTracks; iMC++)
388                 {
389                         AliAODMCParticle *mcTrack = (AliAODMCParticle*)mcArray->At(iMC);
390                         if (!mcTrack) continue;
391                         // Charged track Generated level
392                         Double_t mcTrackPt  = mcTrack->Pt();
393                         if ((mcTrackPt<fPtMin)||(mcTrackPt>6.0)) continue;
394                         
395                         Double_t mcTrackEta = mcTrack->Eta();
396                         Double_t mcTrackPhi = mcTrack->Phi();
397                         Bool_t TrIsPrime    = mcTrack->IsPhysicalPrimary();
398                         Bool_t TrPtMin      = mcTrackPt>fPtMin;
399                         Bool_t TrCharge     = (mcTrack->Charge())!=0;
400                         
401                         if (!TrIsPrime  && !TrPtMin && !TrCharge) continue;  //Check Point 1
402                         
403                         // V0 Generated level
404                         Int_t mcPartPdg           = mcTrack->GetPdgCode();
405                         
406                         Double_t mcRapidity   = mcTrack->Y();
407                         Bool_t V0RapMax       = TMath::Abs(mcRapidity)<fRapidityCut;
408                         Double_t mcMass       = mcTrack->M();
409                         
410                         Double_t mcK0[4] = {mcRapidity,mcTrackPt,mcMass,multiplicityMC};
411                         Double_t mcLa[4] = {mcRapidity,mcTrackPt,mcMass,multiplicityMC};
412                         Double_t mcAl[4] = {mcRapidity,mcTrackPt,mcMass,multiplicityMC};
413                         
414                         
415                         Bool_t IsK0 = mcPartPdg==310;
416                         if (IsK0 && V0RapMax && TrIsPrime) 
417                         {
418                                 fHistMCGenK0->Fill(mcK0);
419                                 selectedV0sMC->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,1));
420                         } 
421                         
422                         Bool_t IsLambda = mcPartPdg==3122;
423                         if (IsLambda && V0RapMax && TrIsPrime) 
424                         {
425                                 fHistMCGenLAM->Fill(mcLa);
426                                 selectedV0sMC->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,2));
427                         }
428                         
429                         Bool_t IsAntiLambda = mcPartPdg==-3122;
430                         if (IsAntiLambda && V0RapMax && TrIsPrime) 
431                         {       
432                                 fHistMCGenALAM->Fill(mcAl);
433                                 selectedV0sMC->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,3));
434                         }                       
435                 }
436                 FillCorrelationSibling(multiplicityMC,selectedTracksLeadingMC,selectedV0sMC,fHistTriggerSibMC,fHistReconstSibMC);
437                 FillCorrelationMixing(multiplicityMC,mcZv,poolmax,poolmin,selectedTracksLeadingMC,selectedV0sMC,fHistTriggerMixMC,fHistReconstMixMC);
438     }
439         
440         // Physics selection  
441         UInt_t maskIsSelected = inEvMain->IsEventSelected();
442         Bool_t isSelected = ((maskIsSelected & AliVEvent::kMB)== AliVEvent::kMB 
443                                           || (maskIsSelected & AliVEvent::kCentral)== AliVEvent::kCentral 
444                                           || (maskIsSelected & AliVEvent::kSemiCentral)== AliVEvent::kSemiCentral);
445     if (!isSelected) return;
446         
447         // Remove pile up
448         if(fRemovePileUP)
449         if(fAODEvent->IsPileupFromSPD()) return;
450         
451         Double_t  dimEventviceReal[5];
452         Double_t  lPrimaryVtxPosition[3];
453         Double_t  lV0Position[3];
454         
455         
456         if(fcollidingSys=="PbPb2010" || fcollidingSys=="PbPb2011"){  //
457         AliCentrality *centralityObj = fAODEvent->GetHeader()->GetCentralityP();
458                 CentPecentAfterPhySel  = centralityObj->GetCentralityPercentileUnchecked("V0M");
459                 if ((CentPecentAfterPhySel < 0.)||(CentPecentAfterPhySel > 90)) return;
460         } //
461         
462         AliAODVertex *myPrimVertex = fAODEvent->GetPrimaryVertex();
463         if (!myPrimVertex) return;
464         myPrimVertex->GetXYZ(lPrimaryVtxPosition);
465         
466         Double_t lPVx = lPrimaryVtxPosition[0];
467         Double_t lPVy = lPrimaryVtxPosition[1];
468         Double_t lPVz = lPrimaryVtxPosition[2];
469         
470         if ((TMath::Abs(lPVz)) >= fpvzcut) return ;
471         if (TMath::Abs(lPVx)<10e-5 && TMath::Abs(lPVy)<10e-5 && TMath::Abs(lPVz)<10e-5) return;
472         
473         dimEventviceReal[0]=lPVx;
474     dimEventviceReal[1]=lPVy;
475     dimEventviceReal[2]=lPVz;
476         
477         multiplicity = fAODEvent->GetNTracks();
478         
479         dimEventviceReal[3]=multiplicity;
480         dimEventviceReal[4]=CentPecentAfterPhySel;
481         
482         fHistEventViceReconst->Fill(dimEventviceReal);
483         
484         if(fcollidingSys=="PP")MultipOrCent=multiplicity;
485         if(fcollidingSys=="PbPb2010" || fcollidingSys=="PbPb2011")MultipOrCent=CentPecentAfterPhySel;
486
487         //---------------------------------------------------------------------------------------------
488         
489         Double_t lDcaPosToPrimVertex = 0;Double_t lDcaNegToPrimVertex = 0;Double_t lDcaV0Daughters     = 0;
490         Double_t lV0cosPointAngle    = 0;Double_t lV0DecayLength      = 0;Double_t lV0Radius           = 0;
491         Double_t lcTauLambda         = 0;Double_t lcTauAntiLambda     = 0;   
492         Double_t lcTauK0s            = 0;   
493         
494         Double_t lInvMassK0   = 0, lInvMassLambda    = 0, lInvMassAntiLambda = 0;
495         Double_t lPtV0s       = 0; Double_t lPhiV0s  = 0; Double_t lEtaV0s   = 0;
496         Double_t lRapK0s      = 0, lRapLambda        = 0, lRapAntiLambda     = 0;
497         Double_t lPzV0s       = 0; Double_t lAlphaV0 = 0, lPtArmV0           = 0;
498         Double_t lPV0s        = 0;
499         
500         TObjArray *selectedTracksLeading=fAnalyseUE->FindLeadingObjects(fAODEvent);
501         if(!selectedTracksLeading) return;
502         selectedTracksLeading->SetOwner(kTRUE);
503         
504         TObjArray * selectedV0s = new TObjArray;
505         selectedV0s->SetOwner(kTRUE);
506         
507         TObjArray * selectedV0sAssoc = new TObjArray;
508         selectedV0sAssoc->SetOwner(kTRUE);
509         
510         Int_t nV0s = fAODEvent->GetNumberOfV0s();
511         
512         for (Int_t i = 0; i < nV0s; i++) 
513         { // start of V0 slection loop
514                 AliAODv0* aodV0 = dynamic_cast<AliAODv0 *>(fAODEvent->GetV0(i));
515                 if (!aodV0) continue;
516                 
517                 if (((aodV0->Pt())<fPtMin)||((aodV0->Pt())>6.0)) continue;
518                 
519                 // get daughters
520             AliAODTrack *myTrackPos=(AliAODTrack *)(aodV0->GetDaughter(0));
521         AliAODTrack *myTrackNeg=(AliAODTrack *)(aodV0->GetDaughter(1));
522                 
523                 if (!myTrackPos || !myTrackNeg) {Printf("ERROR: Could not retreive one of the daughter track");continue;}
524                 
525         if (!IsAcseptedV0(aodV0,myTrackPos,myTrackNeg)) continue;
526                 
527                 // VO's main characteristics to check the reconstruction cuts
528                 lDcaV0Daughters    = aodV0->DcaV0Daughters();
529                 lV0cosPointAngle   = aodV0->CosPointingAngle(lPrimaryVtxPosition);
530                 
531                 aodV0->GetXYZ(lV0Position);
532                 
533                 lV0Radius      = TMath::Sqrt(lV0Position[0]*lV0Position[0]+lV0Position[1]*lV0Position[1]);
534                 lV0DecayLength = TMath::Sqrt(TMath::Power(lV0Position[0] - lPrimaryVtxPosition[0],2) +
535                                                                          TMath::Power(lV0Position[1] - lPrimaryVtxPosition[1],2) +
536                                                                          TMath::Power(lV0Position[2] - lPrimaryVtxPosition[2],2));
537                 
538                 // DCA between daughter and Primary Vertex:
539                 if (myTrackPos) lDcaPosToPrimVertex = aodV0->DcaPosToPrimVertex();
540                 if (myTrackNeg) lDcaNegToPrimVertex = aodV0->DcaNegToPrimVertex();      
541                 
542                 // Quality tracks cuts:
543                 if ( !(IsAcseptedDaughterTrack(myTrackPos)) || !(IsAcseptedDaughterTrack(myTrackNeg)) ) { continue;}
544                 
545                 // Armenteros variables:
546                 lAlphaV0      =  aodV0->AlphaV0();
547                 lPtArmV0      =  aodV0->PtArmV0();
548                 
549                 // Invariant mass
550                 lInvMassK0         = aodV0->MassK0Short();
551                 lInvMassLambda     = aodV0->MassLambda();
552                 lInvMassAntiLambda = aodV0->MassAntiLambda();
553                 
554                 lPtV0s = aodV0->Pt();
555                 lPhiV0s= aodV0->Phi();
556                 lEtaV0s= aodV0->Eta();
557                 lPzV0s = aodV0->Pz();
558                 
559                 // Rapidity:
560                 lRapK0s    = aodV0->RapK0Short();
561                 lRapLambda = aodV0->RapLambda();
562                 lRapAntiLambda = aodV0->Y(-3122);               
563                 
564                 if (lPtV0s==0) {continue;}
565                 
566         Float_t nSigmaPosPion   = 0.;
567         Float_t nSigmaNegPion   = 0.;
568         Float_t nSigmaPosProton = 0.;
569         Float_t nSigmaNegProton = 0.;
570                 
571         const AliAODPid *pPid = myTrackPos->GetDetPid();
572         const AliAODPid *nPid = myTrackNeg->GetDetPid();
573                 
574         if (pPid)
575         {
576             Double_t pdMom = pPid->GetTPCmomentum();
577             if (pdMom<1.0 && (fcollidingSys=="PbPb2010"||fcollidingSys=="PbPb2011"))
578             {
579                 nSigmaPosPion   = fPIDResponse->NumberOfSigmasTPC(myTrackPos, AliPID::kPion);
580                 nSigmaPosProton = fPIDResponse->NumberOfSigmasTPC(myTrackPos, AliPID::kProton);
581             }
582                         
583                         if (fcollidingSys=="PP")
584             {
585                 nSigmaPosPion   = fPIDResponse->NumberOfSigmasTPC(myTrackPos, AliPID::kPion);
586                 nSigmaPosProton = fPIDResponse->NumberOfSigmasTPC(myTrackPos, AliPID::kProton);
587             }
588         }
589                 
590         if (nPid)
591         {
592             Double_t ndMom = nPid->GetTPCmomentum();
593             if (ndMom<1.0 && (fcollidingSys=="PbPb2010"||fcollidingSys=="PbPb2011"))
594             {
595                 nSigmaNegPion   = fPIDResponse->NumberOfSigmasTPC(myTrackNeg, AliPID::kPion);
596                 nSigmaNegProton = fPIDResponse->NumberOfSigmasTPC(myTrackNeg, AliPID::kProton);
597             }
598                         
599                         if (fcollidingSys=="PP")
600             {
601                 nSigmaNegPion   = fPIDResponse->NumberOfSigmasTPC(myTrackNeg, AliPID::kPion);
602                 nSigmaNegProton = fPIDResponse->NumberOfSigmasTPC(myTrackNeg, AliPID::kProton);
603             }
604         }
605                 Bool_t bpPion   = TMath::Abs(nSigmaPosPion)   <= fSigmaPID;
606         Bool_t bpProton = TMath::Abs(nSigmaPosProton) <= fSigmaPID;
607         Bool_t bnPion   = TMath::Abs(nSigmaNegPion)   <= fSigmaPID;
608         Bool_t bnProton = TMath::Abs(nSigmaNegProton) <= fSigmaPID;
609                 
610         Bool_t cutK0Pid         = (bpPion   && bnPion)  ;
611         Bool_t cutLambdaPid     = (bpProton && bnPion)  ;
612         Bool_t cutAntiLambdaPid = (bpPion   && bnProton);
613         //--------------------------------------------------
614                 
615                 lPV0s = TMath::Sqrt(lPzV0s*lPzV0s + lPtV0s*lPtV0s);
616                 
617                 if(lPV0s > 0) lcTauLambda     = (lV0DecayLength*lInvMassLambda)/lPV0s;
618                 if(lPV0s > 0) lcTauAntiLambda = (lV0DecayLength*lInvMassAntiLambda)/lPV0s; 
619                 if(lPV0s > 0) lcTauK0s        = (lV0DecayLength*lInvMassK0)/lPV0s;      
620                 
621                 Bool_t k0ctcut = (lcTauK0s        < fCutCTK0);
622                 Bool_t lactcut = (lcTauLambda     < fCutCTLa);
623                 Bool_t alactcut= (lcTauAntiLambda < fCutCTLa);
624
625                 Bool_t k0APcut = (lPtArmV0>(TMath::Abs(0.2*lAlphaV0)));
626                 
627                 Bool_t k0Rapcut = (TMath::Abs(lRapK0s)        < fRapidityCut);
628                 Bool_t laRapcut = (TMath::Abs(lRapLambda)     < fRapidityCut);
629                 Bool_t alaRapcut= (TMath::Abs(lRapAntiLambda) < fRapidityCut);
630                 
631                 
632                 Bool_t k0cutset = IsAcseptedK0(lV0Radius,lDcaPosToPrimVertex,lDcaNegToPrimVertex,lDcaV0Daughters,lV0cosPointAngle,lInvMassLambda,lInvMassAntiLambda);
633                 Bool_t lacutset = IsAcseptedLA(lV0Radius,lDcaPosToPrimVertex,lDcaNegToPrimVertex,lDcaV0Daughters,lV0cosPointAngle,lInvMassK0);
634                 Bool_t alacutset= IsAcseptedLA(lV0Radius,lDcaNegToPrimVertex,lDcaPosToPrimVertex,lDcaV0Daughters,lV0cosPointAngle,lInvMassK0);
635                 
636                 Double_t spK0[10] = {lInvMassK0, lInvMassLambda,lInvMassAntiLambda,lDcaPosToPrimVertex,lDcaNegToPrimVertex,lV0Radius,lDcaV0Daughters,lV0cosPointAngle,lPtV0s,MultipOrCent};
637                 Double_t spLa[10] = {lInvMassK0, lInvMassLambda,lInvMassAntiLambda,lDcaPosToPrimVertex,lDcaNegToPrimVertex,lV0Radius,lDcaV0Daughters,lV0cosPointAngle,lPtV0s,MultipOrCent};
638                 Double_t spAl[10] = {lInvMassK0, lInvMassLambda,lInvMassAntiLambda,lDcaPosToPrimVertex,lDcaNegToPrimVertex,lV0Radius,lDcaV0Daughters,lV0cosPointAngle,lPtV0s,MultipOrCent};
639         
640                 switch (fCase) {
641                         case 1:
642                                 fHistReconstK0->Fill(spK0); 
643                                 if(IsK0InvMass(lInvMassK0))selectedV0s->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,1));
644                                 
645                                 fHistReconstLA->Fill(spLa); 
646                                 if(IsLambdaInvMass(lInvMassLambda))selectedV0s->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,2));
647                                 
648                                 fHistReconstALA->Fill(spAl);
649                                 if(IsLambdaInvMass(lInvMassAntiLambda))selectedV0s->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,3));
650                                 
651                                 break;
652                                 
653                         case 2:
654                                 if(k0ctcut && k0Rapcut && k0cutset && cutK0Pid)
655                                 {
656                                         fHistReconstK0->Fill(spK0); 
657                                         if(IsK0InvMass(lInvMassK0))selectedV0s->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,1));
658                                 }
659                                 
660                                 if (lactcut && laRapcut && lacutset && cutLambdaPid)
661                                 {
662                                         fHistReconstLA->Fill(spLa); 
663                                         if(IsLambdaInvMass(lInvMassLambda))selectedV0s->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,2));
664                                 }
665                                 
666                                 if (alactcut && alaRapcut && alacutset && cutAntiLambdaPid)
667                                 {
668                                         fHistReconstALA->Fill(spAl);
669                                         if(IsLambdaInvMass(lInvMassAntiLambda))selectedV0s->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,3));
670                                 }
671
672                                 break;
673                                 
674                         case 3:
675                                 if(k0ctcut && k0Rapcut && k0cutset && cutK0Pid && k0APcut)
676                                 {
677                                         fHistReconstK0->Fill(spK0); 
678                                         if(IsK0InvMass(lInvMassK0))selectedV0s->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,1));
679                                 }
680                                 
681                                 if (lactcut && laRapcut && lacutset && cutLambdaPid)
682                                 {
683                                         fHistReconstLA->Fill(spLa); 
684                                         if(IsLambdaInvMass(lInvMassLambda))selectedV0s->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,2));
685                                 }
686                                 
687                                 if (alactcut && alaRapcut && alacutset && cutAntiLambdaPid)
688                                 {
689                                         fHistReconstALA->Fill(spAl);
690                                         if(IsLambdaInvMass(lInvMassAntiLambda))selectedV0s->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,3));
691                                 }
692                                 break;
693                                 
694                         default:
695                                 AliInfo(Form("No case selected"));
696                                 break;
697                 }
698                 
699         if (fAnalysisMC)
700         {
701                 TClonesArray *mcArray = (TClonesArray*)fAODEvent->FindListObject(AliAODMCParticle::StdBranchName());
702                 if(!mcArray)return;
703                 
704                         Int_t myTrackPosLabel        = TMath::Abs(myTrackPos->GetLabel());
705                         Int_t myTrackNegLabel        = TMath::Abs(myTrackNeg->GetLabel());
706                         
707                         AliAODMCParticle *mcPosTrack = (AliAODMCParticle*)mcArray->At(myTrackPosLabel);
708                         AliAODMCParticle *mcNegTrack = (AliAODMCParticle*)mcArray->At(myTrackNegLabel);
709                         
710                         Int_t PosDaughterPdg = mcPosTrack->GetPdgCode();
711                         Int_t NegDaughterPdg = mcNegTrack->GetPdgCode();
712                         
713                         Int_t myTrackPosMotherLabel = mcPosTrack->GetMother();
714                         Int_t myTrackNegMotherLabel = mcNegTrack->GetMother();
715                         
716                         if ((myTrackPosMotherLabel==-1)||(myTrackNegMotherLabel==-1)) continue;
717                         if (myTrackPosMotherLabel!=myTrackNegMotherLabel) continue;
718                         
719                         AliAODMCParticle *mcPosMother = (AliAODMCParticle*)mcArray->At(myTrackPosMotherLabel);
720                         Int_t MotherPdg  = mcPosMother->GetPdgCode();
721                         Bool_t IsPrime   = mcPosMother->IsPhysicalPrimary();
722                         
723                         Double_t rcK0[10] = {lInvMassK0, lInvMassLambda,lInvMassAntiLambda,lDcaPosToPrimVertex,lDcaNegToPrimVertex,lV0Radius,lDcaV0Daughters,lV0cosPointAngle,lPtV0s,MultipOrCent};
724                         Double_t rcLa[10] = {lInvMassK0, lInvMassLambda,lInvMassAntiLambda,lDcaPosToPrimVertex,lDcaNegToPrimVertex,lV0Radius,lDcaV0Daughters,lV0cosPointAngle,lPtV0s,MultipOrCent};
725                         Double_t rcAl[10] = {lInvMassK0, lInvMassLambda,lInvMassAntiLambda,lDcaPosToPrimVertex,lDcaNegToPrimVertex,lV0Radius,lDcaV0Daughters,lV0cosPointAngle,lPtV0s,MultipOrCent};
726                         
727                         switch (fCase) {
728                                 case 1:
729                                         fHistMCAssoK0->Fill(rcK0); 
730                                         if(IsK0InvMass(lInvMassK0))selectedV0sAssoc->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,1));
731                                         
732                                         fHistMCAssoLA->Fill(rcLa);
733                                         if(IsLambdaInvMass(lInvMassLambda))selectedV0sAssoc->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,2));
734                                         
735                                         fHistMCAssoALA->Fill(rcAl);
736                                         if(IsLambdaInvMass(lInvMassAntiLambda))selectedV0sAssoc->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,3));
737                                         
738                                         break;
739                                         
740                                 case 2:
741                                         if ((k0ctcut && k0Rapcut && k0cutset)&&(MotherPdg     ==  310 && 
742                                                                                                                         PosDaughterPdg==  211 && 
743                                                                                                                         NegDaughterPdg== -211 &&
744                                                                                                                         IsPrime))
745                                         {
746                                                 fHistMCAssoK0->Fill(rcK0); 
747                                                 if(IsK0InvMass(lInvMassK0))selectedV0sAssoc->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,1));
748                                         }
749                                         
750                                         if ((lactcut && laRapcut && lacutset)&&(MotherPdg     == 3122 && 
751                                                                                                                         PosDaughterPdg== 2212 && 
752                                                                                                                         NegDaughterPdg== -211 &&
753                                                                                                                         IsPrime)) 
754                                         {
755                                                 fHistMCAssoLA->Fill(rcLa);
756                                                 if(IsLambdaInvMass(lInvMassLambda))selectedV0sAssoc->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,2));
757                                         }
758                                         
759                                         if ((alactcut && alaRapcut && alacutset)&&(MotherPdg     == -3122 && 
760                                                                                                                            PosDaughterPdg==   211 && 
761                                                                                                                            NegDaughterPdg== -2212 &&
762                                                                                                                            IsPrime))
763                                         {
764                                                 fHistMCAssoALA->Fill(rcAl);
765                                                 if(IsLambdaInvMass(lInvMassAntiLambda))selectedV0sAssoc->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,3));
766                                         }
767                                         
768                                         break;
769                                         
770                                 case 3:
771                                         if ((k0ctcut && k0Rapcut && k0cutset && k0APcut)&&(MotherPdg     ==  310 && 
772                                                                                                                                            PosDaughterPdg==  211 && 
773                                                                                                                                            NegDaughterPdg== -211 &&
774                                                                                                                                            IsPrime))
775                                         {
776                                                 fHistMCAssoK0->Fill(rcK0); 
777                                                 if(IsK0InvMass(lInvMassK0))selectedV0sAssoc->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,1));
778                                         }
779                                         
780                                         if ((lactcut && laRapcut && lacutset)&&(MotherPdg     == 3122 && 
781                                                                                                                         PosDaughterPdg== 2212 && 
782                                                                                                                         NegDaughterPdg== -211 &&
783                                                                                                                         IsPrime)) 
784                                         {
785                                                 fHistMCAssoLA->Fill(rcLa);
786                                                 if(IsLambdaInvMass(lInvMassLambda))selectedV0sAssoc->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,2));
787                                         }
788                                         
789                                         if ((alactcut && alaRapcut && alacutset)&&(MotherPdg     == -3122 && 
790                                                                                                                            PosDaughterPdg==   211 && 
791                                                                                                                            NegDaughterPdg== -2212 &&
792                                                                                                                            IsPrime))
793                                         {
794                                                 fHistMCAssoALA->Fill(rcAl);
795                                                 if(IsLambdaInvMass(lInvMassAntiLambda))selectedV0sAssoc->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,3));
796                                         }
797                                         break;
798                                         
799                                 default:
800                                         AliInfo(Form("No case selected"));
801                                         break;
802                         }       
803         }
804         }       
805         
806         FillCorrelationSibling(MultipOrCent,selectedTracksLeading,selectedV0s,fHistTriggerSib,fHistReconstSib);
807         FillCorrelationMixing(MultipOrCent,lPVz,poolmax,poolmin,selectedTracksLeading,selectedV0s,fHistTriggerMix,fHistReconstMix);
808         
809         FillCorrelationSibling(MultipOrCent,selectedTracksLeading,selectedV0sAssoc,0,fHistReconstSibMCAssoc);
810         FillCorrelationMixing(MultipOrCent,lPVz,poolmax,poolmin,selectedTracksLeading,selectedV0sAssoc,0,fHistReconstMixMCAssoc);
811         
812         PostData(1,fOutputList);
813 }       
814 //---------------------------------------------------------------------------------------
815 void AliLeadingV0Correlation::Terminate(Option_t *)
816 {
817         //No need in the grid
818 }
819 //---------------------------------------------------------------------------------------
820 Bool_t AliLeadingV0Correlation::IsAcseptedDaughterTrack(const AliAODTrack *itrack)
821 {
822         if(TMath::Abs(itrack->Eta())>fTrackEtaCut)return kFALSE;
823         
824         if (!itrack->IsOn(AliAODTrack::kTPCrefit)) return kFALSE;
825         
826         Float_t nCrossedRowsTPC = itrack->GetTPCClusterInfo(2,1);
827         if (nCrossedRowsTPC < 70) return kFALSE;
828         
829         Int_t findable=itrack->GetTPCNclsF();
830         if (findable <= 0) return kFALSE;
831         
832         if (nCrossedRowsTPC/findable < 0.8) return kFALSE;
833         return kTRUE;
834 }
835 //---------------------------------------------------------------------------------------
836 Bool_t AliLeadingV0Correlation::IsAcseptedV0(const AliAODv0* aodV0, const AliAODTrack* myTrackPos, const AliAODTrack* myTrackNeg)
837 {
838         if (!aodV0) return kFALSE;
839         
840         // Offline reconstructed V0 only
841     if (aodV0->GetOnFlyStatus()) return kFALSE;
842         
843     // Get daughters and check them
844         myTrackPos=(AliAODTrack *)(aodV0->GetDaughter(0));
845         myTrackNeg=(AliAODTrack *)(aodV0->GetDaughter(1));
846         
847         if (!myTrackPos||!myTrackNeg) return kFALSE;
848         // Unlike signs of daughters
849     if (myTrackPos->Charge() == myTrackNeg->Charge()) return kFALSE;
850         
851         // Track cuts for daughers
852     if ( !(IsAcseptedDaughterTrack(myTrackPos)) || !(IsAcseptedDaughterTrack(myTrackNeg)) ) return kFALSE;
853         
854         // Minimum pt of daughters
855     Double_t lPtPos = myTrackPos->Pt();
856     Double_t lPtNeg = myTrackNeg->Pt();
857         
858         if (lPtPos<fPtMin || lPtNeg<fPtMin) return kFALSE;
859         
860         return kTRUE;
861 }
862 //---------------------------------------------------------------------------------------
863 Bool_t AliLeadingV0Correlation::IsAcseptedK0(Double_t v0rad,
864                                                         Double_t dcaptp,
865                                                         Double_t dcantp,
866                                                         Double_t dcav0d,
867                                                         Double_t cpa,
868                                                         Double_t massLa,
869                                                         Double_t massALa)
870 {       
871                         if(v0rad  >=fV0radius           &&
872                            dcaptp >=fV0PostoPVz         &&
873                            dcantp >=fV0NegtoPVz         &&
874                            dcav0d <=fDCAV0Daughters     &&
875                            cpa    >=fCPAK0                      &&
876                            TMath::Abs(massLa  - 1.115683) > fRejectLamK0 &&
877                            TMath::Abs(massALa - 1.115683) > fRejectLamK0 )return kTRUE;
878         return kFALSE;
879 }
880 //---------------------------------------------------------------------------------------
881 Bool_t AliLeadingV0Correlation::IsAcseptedLA(Double_t v0rad,
882                                                         Double_t dcaptp,
883                                                         Double_t dcantp,
884                                                         Double_t dcav0d,
885                                                         Double_t cpa,
886                                                         Double_t massK0)
887 {
888         if(v0rad  >=fV0radius           &&
889            dcaptp >=fV0PostoPVz         &&
890            dcantp >=fV0NegtoPVz         &&
891            dcav0d <=fDCAV0Daughters     &&
892            cpa    >=fCPALam                     &&
893            TMath::Abs(massK0  - 0.4976) > fRejectK0Lam &&
894            TMath::Abs(massK0  - 0.4976) > fRejectK0Lam )return kTRUE;
895         return kFALSE;
896 }
897 //---------------------------------------------------------------------------------------
898 Bool_t AliLeadingV0Correlation::IsK0InvMass(const Double_t mass) const 
899 {
900         const Float_t massK0            = 0.497; 
901         
902         return ((massK0-fMassCutK0)<=mass && mass<=(massK0 + fMassCutK0))?1:0;
903 }
904 //---------------------------------------------------------------------------------------
905 Bool_t AliLeadingV0Correlation::IsLambdaInvMass(const Double_t mass) const 
906 {
907         const Float_t massLambda        = 1.116; 
908         
909         return ((massLambda-fMassCutLa)<=mass && mass<=(massLambda + fMassCutLa))?1:0;
910 }
911 //---------------------------------------------------------------------------------------
912 Double_t AliLeadingV0Correlation::RangePhi(Double_t DPhi)
913 {
914         if (DPhi < -TMath::Pi()/2)  DPhi += 2*TMath::Pi();
915         if (DPhi > 3*TMath::Pi()/2) DPhi -= 2*TMath::Pi();      
916         return DPhi;    
917 }
918 //---------------------------------------------------------------------------------------
919 Bool_t AliLeadingV0Correlation::IsTrackFromV0(AliAODTrack* track)
920 {
921         Int_t atrID = track->GetID();
922
923         for(int i=0; i<fAODEvent->GetNumberOfV0s(); i++){ // loop over V0s
924                 AliAODv0* aodV0 = fAODEvent->GetV0(i);
925                 
926                 AliAODTrack *trackPos=(AliAODTrack *)(aodV0->GetDaughter(0));
927         AliAODTrack *trackNeg=(AliAODTrack *)(aodV0->GetDaughter(1));
928                         
929                 if ( !(IsAcseptedDaughterTrack(trackPos)) || !(IsAcseptedDaughterTrack(trackNeg)) ) continue;
930                 //----------------------------------
931                 Int_t negID = trackNeg->GetID();
932                 Int_t posID = trackPos->GetID();
933                 
934                 if ((TMath::Abs(negID)+1)==(TMath::Abs(atrID))){ return kTRUE;}
935                 if ((TMath::Abs(posID)+1)==(TMath::Abs(atrID))){ return kTRUE;}
936                 //----------------------------------
937         }
938         return kFALSE;
939 }
940 //---------------------------------------------------------------------------------------
941 void AliLeadingV0Correlation::FillCorrelationSibling(Double_t MultipOrCent,
942                                                                           TObjArray*triggerArray,
943                                                                           TObjArray*selectedV0Array,
944                                                                           THnSparse*triggerHist,
945                                                                           THnSparse*associateHist)
946 {
947         Double_t  binsv0CORR[7];
948         Double_t  binsTrigSib[4];
949         Int_t counterSibMCA=0;
950         
951     for(Int_t i=0;i<triggerArray->GetEntriesFast();i++)
952         {
953                 AliAODTrack* trigger = (AliAODTrack*)triggerArray->At(0);
954                 if(!trigger)continue;
955                 
956                 if(fRemoveAutoCorr) 
957                 if(IsTrackFromV0(trigger))continue;
958                         
959                 Double_t triggerPt  = trigger->Pt();
960                 Double_t triggerPhi = trigger->Phi();
961                 Double_t triggerEta = trigger->Eta();
962                 
963                 if(triggerPt<6.0||triggerPt>12.0)continue;
964                 counterSibMCA++;
965                 if(counterSibMCA==triggerArray->GetEntriesFast()){
966                         
967                         binsTrigSib[0]=triggerPt;
968                         binsTrigSib[1]=triggerPhi;
969                         binsTrigSib[2]=triggerEta;
970                         binsTrigSib[3]=MultipOrCent;
971                         
972                         if(triggerHist)triggerHist->Fill(binsTrigSib);
973                         
974                         for (Int_t j=0; j<selectedV0Array->GetEntriesFast(); j++){
975                                 
976                                 V0Correlationparticle* associate = (V0Correlationparticle*) selectedV0Array->At(j);
977                                 if(!associate)continue;
978                                 
979                                 binsv0CORR[0]= associate->Pt();
980                                 binsv0CORR[1]= associate->Phi();
981                                 binsv0CORR[2]= associate->Eta();
982                                 
983                                 if(binsv0CORR[0]>triggerPt) continue;
984                                 
985                                 binsv0CORR[3]=RangePhi(triggerPhi-binsv0CORR[1]);
986                                 binsv0CORR[4]=triggerEta-binsv0CORR[2];
987                                 binsv0CORR[5]= associate->WhichCandidate();
988                                 binsv0CORR[6]= MultipOrCent;
989                                 
990                                 associateHist->Fill(binsv0CORR);
991                         }
992                 }
993         }
994 }
995 //---------------------------------------------------------------------------------------
996 void AliLeadingV0Correlation::FillCorrelationMixing(Double_t MultipOrCentMix,
997                                                                    Double_t pvxMix,
998                                                                    Double_t poolmax,
999                                                                    Double_t poolmin,
1000                                                                    TObjArray*triggerArray,
1001                                                                    TObjArray*selectedV0Array,
1002                                                                    THnSparse*triggerHist,
1003                                                                    THnSparse*associateHist)
1004 {
1005         if(TMath::Abs(pvxMix)>=fpvzcut || MultipOrCentMix>poolmax || MultipOrCentMix < poolmin)
1006         {
1007                 if(fcollidingSys=="PP")AliInfo(Form("pp Event with Zvertex = %.2f cm and multiplicity = %.0f out of pool bounds, SKIPPING",pvxMix,MultipOrCentMix));
1008                 if(fcollidingSys=="PbPb2010" || fcollidingSys=="PbPb2011") AliInfo(Form("PbPb Event with Zvertex = %.2f cm and centrality = %.1f  out of pool bounds, SKIPPING",pvxMix,MultipOrCentMix));
1009                 return;
1010         }
1011         
1012         Double_t  binsv0CORRMix[7];
1013         Double_t  binsTrigMix[4];
1014         Double_t  counterMix=0;
1015         
1016         AliEventPool* pool = fPoolMgr->GetEventPool(MultipOrCentMix, pvxMix);
1017         if (!pool) AliFatal(Form("No pool found for centrality = %f, zVtx = %f", MultipOrCentMix, pvxMix));
1018         
1019     if (pool->IsReady() || pool->NTracksInPool() > fPoolMinNTracks  || pool->GetCurrentNEvents() > fMinEventsToMix)
1020         {
1021                 Int_t nMix = pool->GetCurrentNEvents();
1022                 for (Int_t jMix=0; jMix<nMix; jMix++){
1023                         
1024                         TObjArray* mixEvents = pool->GetEvent(jMix);
1025                         for (Int_t i=0; i<triggerArray->GetEntriesFast(); i++){
1026                                 
1027                                 AliAODTrack* trig = (AliAODTrack*)triggerArray->At(0);
1028                                 if(!trig)continue;
1029                                 
1030                                 if(fRemoveAutoCorr) 
1031                                 if(IsTrackFromV0(trig))continue;
1032                                 
1033                                 Double_t trigPhi  = trig->Phi();
1034                                 Double_t trigEta  = trig->Eta();
1035                                 Double_t trigPt   = trig->Pt();
1036                                 
1037                                 if(trigPt<6.0||trigPt>12.0)continue;
1038                                 counterMix++;
1039                                 
1040                                 if(counterMix==triggerArray->GetEntriesFast()){
1041                                         
1042                                         binsTrigMix[0]=trigPt;
1043                                         binsTrigMix[1]=trigPhi;
1044                                         binsTrigMix[2]=trigEta;
1045                                         binsTrigMix[3]=MultipOrCentMix;
1046                                         
1047                                         if(triggerHist)triggerHist->Fill(binsTrigMix);
1048                                         
1049                                         for (Int_t j=0; j<mixEvents->GetEntriesFast(); j++){
1050                                                 
1051                                                 V0Correlationparticle* associate = (V0Correlationparticle*) mixEvents->At(j);
1052                                                 if(!associate)continue;
1053                                                 
1054                                                 binsv0CORRMix[0]= associate->Pt();
1055                                                 binsv0CORRMix[1]= associate->Phi();
1056                                                 binsv0CORRMix[2]= associate->Eta();
1057                                                 
1058                                                 if(binsv0CORRMix[0]>trigPt) continue;
1059                                                 
1060                                                 binsv0CORRMix[3]=RangePhi(trigPhi-binsv0CORRMix[1]);
1061                                                 binsv0CORRMix[4]=trigEta-binsv0CORRMix[2];
1062                                                 binsv0CORRMix[5]=associate->WhichCandidate();
1063                                                 binsv0CORRMix[6]=MultipOrCentMix;
1064                                                 
1065                                                 associateHist->Fill(binsv0CORRMix);
1066                                         }
1067                                 }
1068                         }
1069                 }
1070         }
1071         
1072         TObjArray* tracksClone = new TObjArray;
1073         tracksClone->SetOwner(kTRUE);
1074         
1075         for (Int_t i=0; i<selectedV0Array->GetEntriesFast(); i++)
1076         {
1077                 V0Correlationparticle* particle = (V0Correlationparticle*) selectedV0Array->At(i);
1078                 tracksClone->Add(new V0Correlationparticle(particle->Eta(), 
1079                                                                                           particle->Phi(), 
1080                                                                                           particle->Pt(),
1081                                                                                           particle->WhichCandidate()));
1082         };
1083         pool->UpdatePool(tracksClone);
1084 }
1085 //---------------------------------------------------------------------------------------                                                                                       
1086                                                                 
1087                                         
1088                                                                                                 
1089