1 /* Leading Charged Track+V0 Correlation.(Works for Real,Monte Carlo Data)
3 * University of Houston
4 * sandun.pahula.hewage@cern.ch
5 *****************************************************************************************/
16 #include <THnSparse.h>
18 #include <TDatabasePDG.h>
21 #include "AliAnalysisManager.h"
22 #include "AliAODTrack.h"
23 #include "AliAODEvent.h"
25 #include "AliAODVertex.h"
26 #include "AliAODPid.h"
27 #include "AliPIDResponse.h"
28 #include "AliEventPoolManager.h"
29 #include "AliCentrality.h"
30 #include "AliAODHandler.h"
31 #include "AliAODInputHandler.h"
32 #include "AliAODMCParticle.h"
33 #include "AliInputEventHandler.h"
34 #include "AliVParticle.h"
35 #include "AliMultiplicity.h"
36 #include "AliAODMCHeader.h"
38 #include "AliExternalTrackParam.h"
39 #include "AliAnalyseLeadingTrackUE.h"
41 #include "AliLeadingV0Correlation.h"
47 Double_t PI =TMath::Pi();
49 ClassImp(AliLeadingV0Correlation)
50 ClassImp(V0Correlationparticle)
52 //---------------------------------------------------------------------------------------
53 AliLeadingV0Correlation::AliLeadingV0Correlation()
54 : AliAnalysisTaskSE(),
89 fUseChargeHadrons (kTRUE),
92 fHist_Mult_B4_Trg_Sel (0),
93 fHist_Mult_Af_Trg_Sel (0),
94 fHist_Mult_PVz_Cut (0),
95 fHist_Mult_SPD_PVz (0),
96 fHist_Mult_SPD_PVz_Pileup (0),
100 fHistPVxAnalysis (0),
101 fHistPVyAnalysis (0),
102 fHistPVzAnalysis (0),
103 fHistEventViceGen (0),
104 fHistEventViceReconst (0),
108 fHistMCGenLAMXIPLS (0),
117 fHistMCAssoALAXiPlus (0),
120 fHistReconstSibGEN (0),
121 fHistReconstMixGEN (0),
122 fHistReconstSibASO (0),
123 fHistReconstMixASO (0),
124 fHistReconstSibFEED (0),
125 fHistReconstMixFEED (0),
128 fHistTriggerSibGEN (0),
129 fHistTriggerMixGEN (0),
130 fHistTriggerSibASO (0),
131 fHistTriggerMixASO (0)
134 for(Int_t iBin = 0; iBin < 100; iBin++){
135 fZvtxBins[iBin] = 0.;
136 fCentBins[iBin] = 0.;
139 //---------------------------------------------------------------------------------------
140 AliLeadingV0Correlation::AliLeadingV0Correlation(const char *name)
141 : AliAnalysisTaskSE(name),
176 fUseChargeHadrons (kTRUE),
179 fHist_Mult_B4_Trg_Sel (0),
180 fHist_Mult_Af_Trg_Sel (0),
181 fHist_Mult_PVz_Cut (0),
182 fHist_Mult_SPD_PVz (0),
183 fHist_Mult_SPD_PVz_Pileup (0),
187 fHistPVxAnalysis (0),
188 fHistPVyAnalysis (0),
189 fHistPVzAnalysis (0),
190 fHistEventViceGen (0),
191 fHistEventViceReconst (0),
195 fHistMCGenLAMXIPLS (0),
204 fHistMCAssoALAXiPlus (0),
207 fHistReconstSibGEN (0),
208 fHistReconstMixGEN (0),
209 fHistReconstSibASO (0),
210 fHistReconstMixASO (0),
211 fHistReconstSibFEED (0),
212 fHistReconstMixFEED (0),
215 fHistTriggerSibGEN (0),
216 fHistTriggerMixGEN (0),
217 fHistTriggerSibASO (0),
218 fHistTriggerMixASO (0)
222 for(Int_t iBin = 0; iBin < 100; iBin++){
223 fZvtxBins[iBin] = 0.;
224 fCentBins[iBin] = 0.;
226 DefineOutput(1, TList::Class());
229 //---------------------------------------------------------------------------------------
230 AliLeadingV0Correlation::~AliLeadingV0Correlation()
232 if (fOutputList && !AliAnalysisManager::GetAnalysisManager()->IsProofMode()) {
236 //---------------------------------------------------------------------------------------
237 void AliLeadingV0Correlation::UserCreateOutputObjects()
239 fAnalyseUE =new AliAnalyseLeadingTrackUE();
242 fAnalyseUE->SetParticleSelectionCriteria(fFilterBit,fUseChargeHadrons,fTrackEtaCut,fPtMin);
243 fAnalyseUE->DefineESDCuts(fFilterBit);
246 fOutputList = new TList();
247 fOutputList->SetOwner();
249 fHist_Mult_B4_Trg_Sel = new TH2F("fHist_Mult_B4_Trg_Sel","Tracks per event;Nbr of Tracks;Events", 1000, 0, 10000, 1000, 0, 10000);
250 fOutputList->Add(fHist_Mult_B4_Trg_Sel);
252 fHist_Mult_Af_Trg_Sel = new TH2F("fHist_Mult_Af_Trg_Sel","Tracks per event;Nbr of Tracks;Events",1000, 0, 10000, 1000, 0, 10000);
253 fOutputList->Add(fHist_Mult_Af_Trg_Sel);
255 fHist_Mult_PVz_Cut = new TH2F("fHist_Mult_PVz_Cut","Tracks per event;Nbr of Tracks;Events",1000, 0, 10000, 1000, 0, 10000);
256 fOutputList->Add(fHist_Mult_PVz_Cut);
258 fHist_Mult_SPD_PVz = new TH2F("fHist_Mult_SPD_PVz","Tracks per event;Nbr of Tracks;Events",1000, 0, 10000, 1000, 0, 10000);
259 fOutputList->Add(fHist_Mult_SPD_PVz);
261 fHist_Mult_SPD_PVz_Pileup = new TH2F("fHist_Mult_SPD_PVz_Pileup","Tracks per event;Nbr of Tracks;Events",1000, 0, 10000, 1000, 0, 10000);
262 fOutputList->Add(fHist_Mult_SPD_PVz_Pileup);
264 fHistPVx = new TH1F("fHistPVx","PV x position;Nbr of Evts;x", 200, -0.5, 0.5);
265 fOutputList->Add(fHistPVx);
267 fHistPVy = new TH1F("fHistPVy","PV y position;Nbr of Evts;y",200, -0.5, 0.5);
268 fOutputList->Add(fHistPVy);
270 fHistPVz = new TH1F("fHistPVz","PV z position;Nbr of Evts;z",400, -20, 20);
271 fOutputList->Add(fHistPVz);
273 fHistPVxAnalysis = new TH1F("fHistPVxAnalysis","PV x position;Nbr of Evts;x", 200, -0.5, 0.5);
274 fOutputList->Add(fHistPVxAnalysis);
276 fHistPVyAnalysis = new TH1F("fHistPVyAnalysis","PV y position;Nbr of Evts;y",200, -0.5, 0.5);
277 fOutputList->Add(fHistPVyAnalysis);
279 fHistPVzAnalysis = new TH1F("fHistPVzAnalysis","PV z position;Nbr of Evts;z",400, -20, 20);
280 fOutputList->Add(fHistPVzAnalysis);
282 //---------------------------------------------- Events histograms -----------------------------------------------------//
284 fHistEventViceGen= new TH2F("fHistEventViceGen", "fHistEventViceGen", 200, -20, 20, 10,0,1000);
285 fOutputList->Add(fHistEventViceGen);
287 fHistEventViceReconst= new TH2F("fHistEventViceReconst", "fHistEventViceReconst", 200, -20, 20, 10,0,1000);
288 fOutputList->Add(fHistEventViceReconst);
290 fHistMCGenLAM = new TH2F("fHistMCGenLAM" , "fHistMCGenLAM" ,140,1.06,1.2, 120, 0, fTriglow);
291 fOutputList->Add(fHistMCGenLAM);
293 fHistMCGenALAM = new TH2F("fHistMCGenALAM", "fHistMCGenALAM",140,1.06,1.2, 120, 0, fTriglow);
294 fOutputList->Add(fHistMCGenALAM);
296 fHistMCGenK0 = new TH2F("fHistMCGenK0" , "fHistMCGenK0" ,200,0.4,0.6, 120, 0, fTriglow);
297 fOutputList->Add(fHistMCGenK0);
299 fHistMCGenLAMXIPLS = new TH2F("fHistMCGenLAMXIPLS", "fHistMCGenLAMXIPLS",140,1.06,1.2, 120, 0, fTriglow);
300 fOutputList->Add(fHistMCGenLAMXIPLS);
302 fHistMCGenLAMXI = new TH2F("fHistMCGenLAMXI" , "fHistMCGenLAMXI" ,140,1.06,1.2, 120, 0, fTriglow);
303 fOutputList->Add(fHistMCGenLAMXI);
305 //New dimension for feed down corection
307 const Int_t ndimsK0 = 4;
308 Int_t binsK0[ndimsK0] = {200, 120,500,1000};
309 Double_t xminK0[ndimsK0] = {0.4, 0, 0,0.99};
310 Double_t xmaxK0[ndimsK0] = {0.6, fTriglow, 10, 1};
312 const Int_t ndimsLA = 4;
313 Int_t binsLA[ndimsLA] = { 140, 120,500,1000};
314 Double_t xminLA[ndimsLA] = {1.06, 0, 0,0.99};
315 Double_t xmaxLA[ndimsLA] = { 1.2, fTriglow, 10, 1};
317 fHistReconstK0= new THnSparseD("fHistReconstK0" , "fHistReconstK0",ndimsK0,binsK0,xminK0,xmaxK0);
318 fHistReconstK0->Sumw2();
319 fOutputList->Add(fHistReconstK0);
321 fHistReconstLA= new THnSparseD("fHistReconstLA" , "fHistReconstLA",ndimsLA,binsLA,xminLA,xmaxLA);
322 fHistReconstLA->Sumw2();
323 fOutputList->Add(fHistReconstLA);
325 fHistReconstALA= new THnSparseD("fHistReconstALA", "fHistReconstALA",ndimsLA,binsLA,xminLA,xmaxLA);
326 fHistReconstALA->Sumw2();
327 fOutputList->Add(fHistReconstALA);
329 fHistMCAssoK0= new THnSparseD("fHistMCAssoK0" , "fHistMCAssoK0" ,ndimsK0,binsK0,xminK0,xmaxK0);
330 fHistMCAssoK0->Sumw2();
331 fOutputList->Add(fHistMCAssoK0);
333 fHistMCAssoLA= new THnSparseD("fHistMCAssoLA" , "fHistMCAssoLA" ,ndimsLA,binsLA,xminLA,xmaxLA);
334 fHistMCAssoLA->Sumw2();
335 fOutputList->Add(fHistMCAssoLA);
337 fHistMCAssoALA= new THnSparseD("fHistMCAssoALA" , "fHistMCAssoALA" , ndimsLA,binsLA,xminLA,xmaxLA);
338 fHistMCAssoALA->Sumw2();
339 fOutputList->Add(fHistMCAssoALA);
341 fHistMCAssoLAXI= new THnSparseD("fHistMCAssoLAXI" , "fHistMCAssoLAXI" , ndimsLA,binsLA,xminLA,xmaxLA);
342 fHistMCAssoLAXI->Sumw2();
343 fOutputList->Add(fHistMCAssoLAXI);
345 fHistMCAssoALAXiPlus= new THnSparseD("fHistMCAssoALAXiPlus" , "fHistMCAssoALAXiPlus" , ndimsLA,binsLA,xminLA,xmaxLA);
346 fHistMCAssoALAXiPlus->Sumw2();
347 fOutputList->Add(fHistMCAssoALAXiPlus);
349 //--------------------------------------------Correlation Histos -----------------------------------------------------//
351 //0-pTK0,1-PhiK0,2-EtaK0,3-DPhiK0,4-DEtaK0,5-TYPE,6-CutSet
352 const Int_t ndimsv0CORR = 8;
353 Int_t binsv0CORR[ndimsv0CORR] = {120, 200, 200,CorrBinsX, CorrBinsY,4,500,1000};
355 Double_t xminv0CORR[ndimsv0CORR] = { 0, 0,-fTrackEtaCut, -PI/2,-2*fTrackEtaCut,0, 0,0.99};
357 Double_t xmaxv0CORR[ndimsv0CORR] = { fTriglow,2*PI, fTrackEtaCut, 3*PI/2, 2*fTrackEtaCut,4, 10, 1};
359 fHistReconstSib= new THnSparseD("fHistReconstSib", "fHistReconstSib", ndimsv0CORR, binsv0CORR, xminv0CORR, xmaxv0CORR);
360 fHistReconstSib->Sumw2();
361 fOutputList->Add(fHistReconstSib);
363 fHistReconstMix= new THnSparseD("fHistReconstMix", "fHistReconstMix", ndimsv0CORR, binsv0CORR, xminv0CORR, xmaxv0CORR);
364 fHistReconstMix->Sumw2();
365 fOutputList->Add(fHistReconstMix);
367 fHistReconstSibGEN= new THnSparseD("fHistReconstSibGEN", "fHistReconstSibGEN", ndimsv0CORR, binsv0CORR, xminv0CORR, xmaxv0CORR);
368 fHistReconstSibGEN->Sumw2();
369 fOutputList->Add(fHistReconstSibGEN);
371 fHistReconstMixGEN= new THnSparseD("fHistReconstMixGEN", "fHistReconstMixGEN", ndimsv0CORR, binsv0CORR, xminv0CORR, xmaxv0CORR);
372 fHistReconstMixGEN->Sumw2();
373 fOutputList->Add(fHistReconstMixGEN);
375 fHistReconstSibASO= new THnSparseD("fHistReconstSibASO", "fHistReconstSibASO", ndimsv0CORR, binsv0CORR, xminv0CORR, xmaxv0CORR);
376 fHistReconstSibASO->Sumw2();
377 fOutputList->Add(fHistReconstSibASO);
379 fHistReconstMixASO= new THnSparseD("fHistReconstMixASO", "fHistReconstMixASO", ndimsv0CORR, binsv0CORR, xminv0CORR, xmaxv0CORR);
380 fHistReconstMixASO->Sumw2();
381 fOutputList->Add(fHistReconstMixASO);
383 fHistReconstSibFEED= new THnSparseD("fHistReconstSibFEED", "fHistReconstSibFEED", ndimsv0CORR, binsv0CORR, xminv0CORR, xmaxv0CORR);
384 fHistReconstSibFEED->Sumw2();
385 fOutputList->Add(fHistReconstSibFEED);
387 fHistReconstMixFEED= new THnSparseD("fHistReconstMixFEED", "fHistReconstMixFEED", ndimsv0CORR, binsv0CORR, xminv0CORR, xmaxv0CORR);
388 fHistReconstMixFEED->Sumw2();
389 fOutputList->Add(fHistReconstMixFEED);
392 fHistTriggerSib= new TH3F("fHistTriggerSib", "fHistTriggerSib", 100, fTriglow, fTrighigh,200,0,2*PI,200,-fTrackEtaCut,fTrackEtaCut);
393 fHistTriggerSib->Sumw2();
394 fOutputList->Add(fHistTriggerSib);
396 fHistTriggerMix= new TH1F("fHistTriggerMix", "fHistTriggerMix", 100, fTriglow, fTrighigh);
397 fHistTriggerMix->Sumw2();
398 fOutputList->Add(fHistTriggerMix);
400 fHistTriggerSibGEN= new TH3F("fHistTriggerSibGEN", "fHistTriggerSibGEN", 100, fTriglow, fTrighigh,200,0,2*PI,200,-fTrackEtaCut,fTrackEtaCut);
401 fHistTriggerSibGEN->Sumw2();
402 fOutputList->Add(fHistTriggerSibGEN);
404 fHistTriggerMixGEN= new TH1F("fHistTriggerMixGEN", "fHistTriggerMixGEN", 100, fTriglow, fTrighigh);
405 fHistTriggerMixGEN->Sumw2();
406 fOutputList->Add(fHistTriggerMixGEN);
408 fHistTriggerSibASO= new TH3F("fHistTriggerSibASO", "fHistTriggerSibASO", 100, fTriglow, fTrighigh,200,0,2*PI,200,-fTrackEtaCut,fTrackEtaCut);
409 fHistTriggerSibASO->Sumw2();
410 fOutputList->Add(fHistTriggerSibASO);
412 fHistTriggerMixASO= new TH1F("fHistTriggerMixASO", "fHistTriggerMixASO", 100, fTriglow, fTrighigh);
413 fHistTriggerMixASO->Sumw2();
414 fOutputList->Add(fHistTriggerMixASO);
416 //----------------------------------------------Event Pool-----------------------------------------------------//
417 fPoolMgr = new AliEventPoolManager(fPoolMaxNEvents, fPoolMinNTracks, fNCentBins, fCentBins, fNzVtxBins, fZvtxBins);
418 if(!fPoolMgr) return;
420 PostData(1, fOutputList);
422 //---------------------------------------------------------------------------------------
423 void AliLeadingV0Correlation::UserExec(Option_t *)
426 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
427 AliInputEventHandler *inEvMain = (AliInputEventHandler*)(mgr->GetInputEventHandler());
428 if (!inEvMain) return;
430 // Pointers to PID Response objects.
431 fPIDResponse = inEvMain->GetPIDResponse();
432 if(!fPIDResponse) return;
434 fAODEvent = dynamic_cast<AliAODEvent*>(inEvMain->GetEvent());
435 if(!fAODEvent) return;
437 Int_t ltrackMultiplicity = 0;
438 Int_t lrefMultiplicity = 0;
440 //------------------------------------------------
441 // Before Physics Selection
442 //------------------------------------------------
443 ltrackMultiplicity = (InputEvent())->GetNumberOfTracks();
444 AliAODHeader * header = dynamic_cast<AliAODHeader*>(fAODEvent->GetHeader());
445 if(!header) AliFatal("Not a standard AOD");
446 lrefMultiplicity = header->GetRefMultiplicity();
448 fHist_Mult_B4_Trg_Sel->Fill(ltrackMultiplicity,lrefMultiplicity);
450 Double_t * CentBins = fCentBins;
451 Double_t poolmin = CentBins[0];
452 Double_t poolmax = CentBins[fNCentBins];
454 //----------------------------------------------------------
455 // Efficency denomenator comes before the physics selection
456 //----------------------------------------------------------
458 Double_t dimEventviceMC[2];
459 if(fAnalysisMC) //Efficency denomenator comes before the physics selection
461 AliAODMCHeader *aodMCheader = (AliAODMCHeader*)fAODEvent->FindListObject(AliAODMCHeader::StdBranchName());
462 if(!aodMCheader) return;
463 Float_t mcZv = aodMCheader->GetVtxZ();
465 if (TMath::Abs(mcZv) >= fpvzcut) return;
467 dimEventviceMC[0]=aodMCheader->GetVtxZ();
469 TClonesArray *mcArray = (TClonesArray*)fAODEvent->FindListObject(AliAODMCParticle::StdBranchName());
472 Int_t nMCTracks = mcArray->GetEntriesFast();
474 dimEventviceMC[1]=nMCTracks;
475 fHistEventViceGen->Fill(dimEventviceMC[0],dimEventviceMC[1]);
477 TObjArray *selectedTracksLeadingMC=fAnalyseUE->FindLeadingObjects(mcArray);
478 if(!selectedTracksLeadingMC) return;
479 selectedTracksLeadingMC->SetOwner(kTRUE);
481 TObjArray * selectedV0sMC =new TObjArray;
482 selectedV0sMC->SetOwner(kTRUE);
484 TObjArray * selectedV0sMCXI =new TObjArray;
485 selectedV0sMCXI->SetOwner(kTRUE);
487 for (Int_t iMC = 0; iMC<nMCTracks; iMC++)
489 AliAODMCParticle *mcTrack = (AliAODMCParticle*)mcArray->At(iMC);
490 if (!mcTrack) continue;
491 // Charged track Generated level
492 Double_t mcTrackPt = mcTrack->Pt();
493 if ((mcTrackPt<fPtMin)||(mcTrackPt>fTriglow)) continue;
495 Double_t mcTrackEta = mcTrack->Eta();
496 Double_t mcTrackPhi = mcTrack->Phi();
497 Bool_t TrIsPrime = mcTrack->IsPhysicalPrimary();
498 Bool_t TrPtMin = mcTrackPt>fPtMin;
499 Bool_t TrCharge = (mcTrack->Charge())!=0;
501 if (!TrIsPrime && !TrPtMin && !TrCharge) continue; //Check Point 1
503 // V0 Generated level
504 Int_t mcPartPdg = mcTrack->GetPdgCode();
506 Double_t mcRapidity = mcTrack->Y();
507 Bool_t V0RapMax = TMath::Abs(mcRapidity)<fRapidityCut;
508 Bool_t V0EtaMax = TMath::Abs(mcTrackEta)<fTrackEtaCut;
509 Double_t mcMass = mcTrack->M();
511 Double_t mcK0[3] = {mcMass,mcTrackPt,static_cast<Double_t>(nMCTracks)};
512 Double_t mcLa[3] = {mcMass,mcTrackPt,static_cast<Double_t>(nMCTracks)};
513 Double_t mcAl[3] = {mcMass,mcTrackPt,static_cast<Double_t>(nMCTracks)};
515 Int_t myTrackMotherLabel = mcTrack->GetMother();
517 AliAODMCParticle *mcMother = (AliAODMCParticle*)mcArray->At(myTrackMotherLabel);
518 if (!mcMother) continue;
519 Int_t MotherPdg = mcMother->GetPdgCode();
520 Bool_t IsK0 = mcPartPdg==310;
521 Bool_t IsLambda = mcPartPdg==3122;
522 Bool_t IsAntiLambda = mcPartPdg==-3122;
524 Bool_t IsXImin =MotherPdg== 3312;
525 Bool_t IsXIPlus =MotherPdg==-3312;
526 Bool_t IsXizero =MotherPdg== 3322;
527 Bool_t IsOmega =MotherPdg== 3334;
534 fHistMCGenK0->Fill(mcK0[0],mcK0[1]);
535 selectedV0sMC->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,1,0,0));
540 fHistMCGenLAM->Fill(mcLa[0],mcLa[1]);
541 selectedV0sMC->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,2,0,0));
546 fHistMCGenALAM->Fill(mcAl[0],mcAl[1]);
547 selectedV0sMC->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,3,0,0));
550 if (IsLambda && (IsXizero || IsXImin))
552 selectedV0sMCXI->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,1,0,0));
553 fHistMCGenLAMXI->Fill(mcLa[0],mcLa[1]);
556 if (IsLambda && IsOmega)
558 selectedV0sMCXI->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,2,0,0));
561 if (IsAntiLambda && IsXIPlus)
563 selectedV0sMCXI->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,3,0,0));
564 fHistMCGenLAMXIPLS->Fill(mcAl[0],mcAl[1]);
571 if (IsK0 && V0RapMax && TrIsPrime)
573 fHistMCGenK0->Fill(mcK0[0],mcK0[1]);
574 selectedV0sMC->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,1,0,0));
577 if (IsLambda && V0RapMax && TrIsPrime)
579 fHistMCGenLAM->Fill(mcLa[0],mcLa[1]);
580 selectedV0sMC->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,2,0,0));
583 if (IsAntiLambda && V0RapMax && TrIsPrime)
585 fHistMCGenALAM->Fill(mcAl[0],mcAl[1]);
586 selectedV0sMC->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,3,0,0));
589 if (IsLambda && V0RapMax && (IsXizero || IsXImin))
591 selectedV0sMCXI->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,1,0,0));
592 fHistMCGenLAMXI->Fill(mcLa[0],mcLa[1]);
595 if (IsLambda && V0RapMax && IsOmega)
597 selectedV0sMCXI->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,2,0,0));
600 if (IsAntiLambda && V0RapMax && IsXIPlus)
602 selectedV0sMCXI->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,3,0,0));
603 fHistMCGenLAMXIPLS->Fill(mcAl[0],mcAl[1]);
610 if (IsK0 && V0EtaMax && TrIsPrime)
612 fHistMCGenK0->Fill(mcK0[0],mcK0[1]);
613 selectedV0sMC->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,1,0,0));
616 if (IsLambda && V0EtaMax && TrIsPrime)
618 fHistMCGenLAM->Fill(mcLa[0],mcLa[1]);
619 selectedV0sMC->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,2,0,0));
622 if (IsAntiLambda && V0EtaMax && TrIsPrime)
624 fHistMCGenALAM->Fill(mcAl[0],mcAl[1]);
625 selectedV0sMC->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,3,0,0));
628 if (IsLambda && V0EtaMax && (IsXizero || IsXImin))
630 selectedV0sMCXI->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,1,0,0));
631 fHistMCGenLAMXI->Fill(mcLa[0],mcLa[1]);
634 if (IsLambda && V0EtaMax && IsOmega)
636 selectedV0sMCXI->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,2,0,0));
639 if (IsAntiLambda && V0EtaMax && IsXIPlus)
641 selectedV0sMCXI->Add(new V0Correlationparticle(mcTrackEta,mcTrackPhi,mcTrackPt,3,0,0));
642 fHistMCGenLAMXIPLS->Fill(mcAl[0],mcAl[1]);
647 AliInfo(Form("No case selected"));
652 FillCorrelationSibling(nMCTracks,selectedTracksLeadingMC,selectedV0sMC,fHistTriggerSibGEN,fHistReconstSibGEN);
653 FillCorrelationMixing(nMCTracks,mcZv,poolmax,poolmin,selectedTracksLeadingMC,selectedV0sMC,fHistTriggerMixGEN,fHistReconstMixGEN);
655 FillCorrelationSibling(nMCTracks,selectedTracksLeadingMC,selectedV0sMCXI,0,fHistReconstSibFEED);
656 FillCorrelationMixing(nMCTracks,mcZv,poolmax,poolmin,selectedTracksLeadingMC,selectedV0sMCXI,0,fHistReconstMixFEED);
659 // End Loop over MC condition
661 //------------------------------------------------
663 //------------------------------------------------
664 UInt_t maskIsSelected = inEvMain->IsEventSelected();
665 Bool_t isSelected = ((maskIsSelected & AliVEvent::kMB)== AliVEvent::kMB);
666 if (!isSelected) return;
668 //------------------------------------------------
669 // After Trigger Selection
670 //------------------------------------------------
672 fHist_Mult_Af_Trg_Sel->Fill(ltrackMultiplicity,lrefMultiplicity);
674 //------------------------------------------------
675 // Getting: Primary Vertex + MagField Info
676 //------------------------------------------------
677 Double_t dimEventviceReal[3];
678 Double_t lBestPrimaryVtxPos[3];
679 Double_t tPrimaryVtxPosition[3];
680 Double_t lV0Position[3];
682 AliAODVertex *lPrimaryBestAODVtx = fAODEvent->GetPrimaryVertex();
683 if (!lPrimaryBestAODVtx) return;
684 // get the best primary vertex available for the event
685 // As done in AliCascadeVertexer, we keep the one which is the best one available.
686 // between : Tracking vertex > SPD vertex > TPC vertex > default SPD vertex
687 // This one will be used for next calculations (DCA essentially)
688 lPrimaryBestAODVtx->GetXYZ(lBestPrimaryVtxPos);
690 const AliVVertex *primaryVtx = fAODEvent->GetPrimaryVertex();
691 if(!primaryVtx)return;
692 tPrimaryVtxPosition[0] = primaryVtx->GetX();
693 tPrimaryVtxPosition[1] = primaryVtx->GetY();
694 tPrimaryVtxPosition[2] = primaryVtx->GetZ();
695 fHistPVx->Fill( tPrimaryVtxPosition[0] );
696 fHistPVy->Fill( tPrimaryVtxPosition[1] );
697 fHistPVz->Fill( tPrimaryVtxPosition[2] );
699 //------------------------------------------------
700 // Primary Vertex Z position: SKIP
701 //------------------------------------------------
703 Double_t lPVx = lBestPrimaryVtxPos[0];
704 Double_t lPVy = lBestPrimaryVtxPos[1];
705 Double_t lPVz = lBestPrimaryVtxPos[2];
707 if ((TMath::Abs(lPVz)) >= fpvzcut) return ;
708 if (TMath::Abs(lPVx)<10e-5 && TMath::Abs(lPVy)<10e-5 && TMath::Abs(lPVz)<10e-5) return;
709 fHist_Mult_PVz_Cut->Fill(ltrackMultiplicity,lrefMultiplicity);
711 //------------------------------------------------
712 // Only look at events with well-established PV
713 //------------------------------------------------
715 const AliAODVertex *lPrimaryTrackingAODVtxCheck = fAODEvent->GetPrimaryVertex();
716 const AliAODVertex *lPrimarySPDVtx = fAODEvent->GetPrimaryVertexSPD();
717 if (!lPrimarySPDVtx && !lPrimaryTrackingAODVtxCheck )return;
719 fHist_Mult_SPD_PVz->Fill(ltrackMultiplicity,lrefMultiplicity);
720 //------------------------------------------------
722 //------------------------------------------------
724 // FIXME : quality selection regarding pile-up rejection
725 if(fAODEvent->IsPileupFromSPD()) return;
726 fHist_Mult_SPD_PVz_Pileup->Fill(ltrackMultiplicity,lrefMultiplicity);
728 fHistPVxAnalysis->Fill(tPrimaryVtxPosition[0]);
729 fHistPVyAnalysis->Fill(tPrimaryVtxPosition[1]);
730 fHistPVzAnalysis->Fill(tPrimaryVtxPosition[2]);
732 dimEventviceReal[0]=tPrimaryVtxPosition[2];
733 dimEventviceReal[1]=ltrackMultiplicity;
735 fHistEventViceReconst->Fill(dimEventviceReal[0],dimEventviceReal[1]);
737 //---------------------------------------------------------------------------------------------
739 Double_t lDcaPosToPrimVertex = 0;Double_t lDcaNegToPrimVertex = 0;Double_t lDcaV0Daughters = 0;
740 Double_t lV0cosPointAngle = 0;Double_t lV0DecayLength = 0;Double_t lV0Radius = 0;
741 Double_t lcTauLambda = 0;Double_t lcTauAntiLambda = 0;
742 Double_t lcTauK0s = 0;
743 Double_t lDCAV0PVz = 0;
745 Double_t lInvMassK0 = 0, lInvMassLambda = 0, lInvMassAntiLambda = 0;
746 Double_t lPtV0s = 0; Double_t lPhiV0s = 0; Double_t lEtaV0s = 0;
747 Double_t lRapK0 = 0, lRapLambda = 0, lRapAntiLambda = 0;
751 TObjArray *selectedTracksLeading=0;
752 selectedTracksLeading=fAnalyseUE->FindLeadingObjects(fAODEvent);
753 if(!selectedTracksLeading) return;
754 selectedTracksLeading->SetOwner(kTRUE);
756 TObjArray * selectedV0s = new TObjArray;
757 selectedV0s->SetOwner(kTRUE);
759 TObjArray * selectedV0sAssoc = new TObjArray;
760 selectedV0sAssoc->SetOwner(kTRUE);
762 Int_t nV0s = fAODEvent->GetNumberOfV0s();
764 for (Int_t i = 0; i < nV0s; i++)
765 { // start of V0 slection loop
766 AliAODv0* aodV0 = dynamic_cast<AliAODv0 *>(fAODEvent->GetV0(i));
767 if (!aodV0) continue;
769 if (((aodV0->Pt())<fPtMin)||((aodV0->Pt())>fTriglow)) continue;
772 AliAODTrack *myTrackPos=(AliAODTrack *)(aodV0->GetDaughter(0));
773 AliAODTrack *myTrackNeg=(AliAODTrack *)(aodV0->GetDaughter(1));
775 if (!myTrackPos || !myTrackNeg) continue;
777 if (!IsAcseptedV0(aodV0,myTrackPos,myTrackNeg)) continue;
779 // VO's main characteristics to check the reconstruction cuts
780 lDcaV0Daughters = aodV0->DcaV0Daughters();
781 lV0cosPointAngle = aodV0->CosPointingAngle(lBestPrimaryVtxPos);
783 aodV0->GetXYZ(lV0Position);
785 lV0Radius = TMath::Sqrt(lV0Position[0]*lV0Position[0]+lV0Position[1]*lV0Position[1]);
786 lV0DecayLength = TMath::Sqrt(TMath::Power(lV0Position[0] - tPrimaryVtxPosition[0],2) +
787 TMath::Power(lV0Position[1] - tPrimaryVtxPosition[1],2) +
788 TMath::Power(lV0Position[2] - tPrimaryVtxPosition[2],2));
790 // DCA between daughter and Primary Vertex:
791 if (myTrackPos) lDcaPosToPrimVertex = aodV0->DcaPosToPrimVertex();
792 if (myTrackNeg) lDcaNegToPrimVertex = aodV0->DcaNegToPrimVertex();
793 lDCAV0PVz = aodV0->DcaV0ToPrimVertex();
795 // Quality tracks cuts:
796 if ( !(IsAcseptedDaughterTrack(myTrackPos)) || !(IsAcseptedDaughterTrack(myTrackNeg)) ) { continue;}
799 lInvMassK0 = aodV0->MassK0Short();
800 lInvMassLambda = aodV0->MassLambda();
801 lInvMassAntiLambda = aodV0->MassAntiLambda();
803 lPtV0s = aodV0->Pt();
804 lPhiV0s= aodV0->Phi();
805 lEtaV0s= aodV0->Eta();
806 lPzV0s = aodV0->Pz();
809 lRapK0 = aodV0->RapK0Short();
810 lRapLambda = aodV0->RapLambda();
811 lRapAntiLambda = aodV0->Y(-3122);
813 if (lPtV0s==0) {continue;}
815 Float_t nSigmaPosPion = 0.;
816 Float_t nSigmaNegPion = 0.;
817 Float_t nSigmaPosProton = 0.;
818 Float_t nSigmaNegProton = 0.;
820 const AliAODPid *pPid = myTrackPos->GetDetPid();
821 const AliAODPid *nPid = myTrackNeg->GetDetPid();
825 Double_t pdMom = pPid->GetTPCmomentum();
826 if (pdMom<1.0 && (fcollidingSys=="PbPb"))
828 nSigmaPosPion = fPIDResponse->NumberOfSigmasTPC(myTrackPos, AliPID::kPion);
829 nSigmaPosProton = fPIDResponse->NumberOfSigmasTPC(myTrackPos, AliPID::kProton);
832 if (fcollidingSys=="PP")
834 nSigmaPosPion = fPIDResponse->NumberOfSigmasTPC(myTrackPos, AliPID::kPion);
835 nSigmaPosProton = fPIDResponse->NumberOfSigmasTPC(myTrackPos, AliPID::kProton);
841 Double_t ndMom = nPid->GetTPCmomentum();
842 if (ndMom<1.0 && (fcollidingSys=="PbPb"))
844 nSigmaNegPion = fPIDResponse->NumberOfSigmasTPC(myTrackNeg, AliPID::kPion);
845 nSigmaNegProton = fPIDResponse->NumberOfSigmasTPC(myTrackNeg, AliPID::kProton);
848 if (fcollidingSys=="PP")
850 nSigmaNegPion = fPIDResponse->NumberOfSigmasTPC(myTrackNeg, AliPID::kPion);
851 nSigmaNegProton = fPIDResponse->NumberOfSigmasTPC(myTrackNeg, AliPID::kProton);
854 Bool_t bpPion = TMath::Abs(nSigmaPosPion) <= fSigmaPID;
855 Bool_t bpProton = TMath::Abs(nSigmaPosProton) <= fSigmaPID;
856 Bool_t bnPion = TMath::Abs(nSigmaNegPion) <= fSigmaPID;
857 Bool_t bnProton = TMath::Abs(nSigmaNegProton) <= fSigmaPID;
859 Bool_t cutK0Pid = (bpPion && bnPion) ;
860 Bool_t cutLambdaPid = (bpProton && bnPion) ;
861 Bool_t cutAntiLambdaPid = (bpPion && bnProton);
862 //--------------------------------------------------
864 lPV0s = TMath::Sqrt(lPzV0s*lPzV0s + lPtV0s*lPtV0s);
866 if(lPV0s > 0) lcTauLambda = (lV0DecayLength*lInvMassLambda)/lPV0s;
867 if(lPV0s > 0) lcTauAntiLambda = (lV0DecayLength*lInvMassAntiLambda)/lPV0s;
868 if(lPV0s > 0) lcTauK0s = (lV0DecayLength*lInvMassK0)/lPV0s;
870 Bool_t k0ctcut = (lcTauK0s < fCutCTK0);
871 Bool_t lactcut = (lcTauLambda < fCutCTLa);
872 Bool_t alactcut= (lcTauAntiLambda < fCutCTLa);
874 Bool_t k0Rapcut = (TMath::Abs(lRapK0) < fRapidityCut);
875 Bool_t laRapcut = (TMath::Abs(lRapLambda) < fRapidityCut);
876 Bool_t alaRapcut= (TMath::Abs(lRapAntiLambda) < fRapidityCut);
878 Bool_t V0EtaMax= (TMath::Abs(lEtaV0s) < fTrackEtaCut);
880 Bool_t k0cutset = IsAcseptedK0(lV0Radius,lDcaPosToPrimVertex,lDcaNegToPrimVertex,lDcaV0Daughters,lV0cosPointAngle,lInvMassLambda,lInvMassAntiLambda);
881 Bool_t lacutset = IsAcseptedLA(lV0Radius,lDcaPosToPrimVertex,lDcaNegToPrimVertex,lDcaV0Daughters,lV0cosPointAngle,lInvMassK0);
882 Bool_t alacutset= IsAcseptedLA(lV0Radius,lDcaNegToPrimVertex,lDcaPosToPrimVertex,lDcaV0Daughters,lV0cosPointAngle,lInvMassK0);
884 Double_t spK0[4] = {lInvMassK0,lPtV0s,lDCAV0PVz,lV0cosPointAngle};
885 Double_t spLa[4] = {lInvMassLambda,lPtV0s,lDCAV0PVz,lV0cosPointAngle};
886 Double_t spAl[4] = {lInvMassAntiLambda,lPtV0s,lDCAV0PVz,lV0cosPointAngle};
890 fHistReconstK0->Fill(spK0);
891 if(IsK0InvMass(lInvMassK0))selectedV0s->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,1,lDCAV0PVz,lV0cosPointAngle));
893 fHistReconstLA->Fill(spLa);
894 if(IsLambdaInvMass(lInvMassLambda))selectedV0s->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,2,lDCAV0PVz,lV0cosPointAngle));
896 fHistReconstALA->Fill(spAl);
897 if(IsLambdaInvMass(lInvMassAntiLambda))selectedV0s->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,3,lDCAV0PVz,lV0cosPointAngle));
902 if(k0ctcut && k0Rapcut && k0cutset && cutK0Pid)
904 fHistReconstK0->Fill(spK0);
905 if(IsK0InvMass(lInvMassK0))selectedV0s->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,1,lDCAV0PVz,lV0cosPointAngle));
908 if (lactcut && laRapcut && lacutset && cutLambdaPid)
910 fHistReconstLA->Fill(spLa);
911 if(IsLambdaInvMass(lInvMassLambda))selectedV0s->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,2,lDCAV0PVz,lV0cosPointAngle));
914 if (alactcut && alaRapcut && alacutset && cutAntiLambdaPid)
916 fHistReconstALA->Fill(spAl);
917 if(IsLambdaInvMass(lInvMassAntiLambda))selectedV0s->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,3,lDCAV0PVz,lV0cosPointAngle));
923 if(k0ctcut && V0EtaMax && k0cutset && cutK0Pid)
925 fHistReconstK0->Fill(spK0);
926 if(IsK0InvMass(lInvMassK0))selectedV0s->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,1,lDCAV0PVz,lV0cosPointAngle));
929 if (lactcut && V0EtaMax && lacutset && cutLambdaPid)
931 fHistReconstLA->Fill(spLa);
932 if(IsLambdaInvMass(lInvMassLambda))selectedV0s->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,2,lDCAV0PVz,lV0cosPointAngle));
935 if (alactcut && V0EtaMax && alacutset && cutAntiLambdaPid)
937 fHistReconstALA->Fill(spAl);
938 if(IsLambdaInvMass(lInvMassAntiLambda))selectedV0s->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,3,lDCAV0PVz,lV0cosPointAngle));
943 AliInfo(Form("No case selected"));
949 TClonesArray *mcArray = (TClonesArray*)fAODEvent->FindListObject(AliAODMCParticle::StdBranchName());
952 Int_t myTrackPosLabel = TMath::Abs(myTrackPos->GetLabel());
953 Int_t myTrackNegLabel = TMath::Abs(myTrackNeg->GetLabel());
955 AliAODMCParticle *mcPosTrack = (AliAODMCParticle*)mcArray->At(myTrackPosLabel);
956 if(!mcPosTrack)continue;
957 AliAODMCParticle *mcNegTrack = (AliAODMCParticle*)mcArray->At(myTrackNegLabel);
958 if(!mcNegTrack)continue;
960 Int_t PosDaughterPdg = mcPosTrack->GetPdgCode();
961 Int_t NegDaughterPdg = mcNegTrack->GetPdgCode();
963 Int_t myTrackPosMotherLabel = mcPosTrack->GetMother();
964 Int_t myTrackNegMotherLabel = mcNegTrack->GetMother();
966 if ((myTrackPosMotherLabel==-1)||(myTrackNegMotherLabel==-1)) continue;
967 if (myTrackPosMotherLabel!=myTrackNegMotherLabel) continue;
969 AliAODMCParticle *mcPosMother = (AliAODMCParticle*)mcArray->At(myTrackPosMotherLabel);
970 if(!mcPosMother)continue;
971 Int_t MotherPdg = mcPosMother->GetPdgCode();
972 Bool_t IsPrime = mcPosMother->IsPhysicalPrimary();
974 Int_t myGrandMotherLabel = mcPosMother->GetMother();
975 AliAODMCParticle *mcGrandMother = (AliAODMCParticle*)mcArray->At(myGrandMotherLabel);
976 Int_t GrandMotherPdg = mcGrandMother->GetPdgCode();
978 Double_t rcK0[4] = {lInvMassK0,lPtV0s,lDCAV0PVz,lV0cosPointAngle};
979 Double_t rcLa[4] = {lInvMassLambda,lPtV0s,lDCAV0PVz,lV0cosPointAngle};
980 Double_t rcAl[4] = {lInvMassAntiLambda,lPtV0s,lDCAV0PVz,lV0cosPointAngle};
984 fHistMCAssoK0->Fill(rcK0);
985 if(IsK0InvMass(lInvMassK0))selectedV0sAssoc->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,1,lDCAV0PVz,lV0cosPointAngle));
987 fHistMCAssoLA->Fill(rcLa);
988 if(IsLambdaInvMass(lInvMassLambda))selectedV0sAssoc->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,2,lDCAV0PVz,lV0cosPointAngle));
990 fHistMCAssoALA->Fill(rcAl);
991 if(IsLambdaInvMass(lInvMassAntiLambda))selectedV0sAssoc->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,3,lDCAV0PVz,lV0cosPointAngle));
996 if ((k0ctcut && k0Rapcut && k0cutset)&&(MotherPdg == 310 &&
997 PosDaughterPdg== 211 &&
998 NegDaughterPdg== -211 &&
1001 fHistMCAssoK0->Fill(rcK0);
1002 if(IsK0InvMass(lInvMassK0))selectedV0sAssoc->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,1,lDCAV0PVz,lV0cosPointAngle));
1005 if ((lactcut && laRapcut && lacutset)&&(MotherPdg == 3122 &&
1006 PosDaughterPdg== 2212 &&
1007 NegDaughterPdg== -211 &&
1010 fHistMCAssoLA->Fill(rcLa);
1011 if(IsLambdaInvMass(lInvMassLambda))selectedV0sAssoc->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,2,lDCAV0PVz,lV0cosPointAngle));
1014 if ((alactcut && alaRapcut && alacutset)&&(MotherPdg == -3122 &&
1015 PosDaughterPdg== 211 &&
1016 NegDaughterPdg== -2212 &&
1019 fHistMCAssoALA->Fill(rcAl);
1020 if(IsLambdaInvMass(lInvMassAntiLambda))selectedV0sAssoc->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,3,lDCAV0PVz,lV0cosPointAngle));
1023 if ((lactcut && laRapcut && lacutset)&&(MotherPdg == 3122 &&
1024 PosDaughterPdg== 2212 &&
1025 NegDaughterPdg== -211 &&
1026 (GrandMotherPdg==3322 ||GrandMotherPdg==3312)))
1028 fHistMCAssoLAXI->Fill(rcLa);
1031 if ((alactcut && alaRapcut && alacutset)&&(MotherPdg == -3122 &&
1032 PosDaughterPdg== 211 &&
1033 NegDaughterPdg== -2212 &&
1034 GrandMotherPdg== -3312))
1036 fHistMCAssoALAXiPlus->Fill(rcAl);
1042 if ((k0ctcut && V0EtaMax && k0cutset)&&(MotherPdg == 310 &&
1043 PosDaughterPdg== 211 &&
1044 NegDaughterPdg== -211 &&
1047 fHistMCAssoK0->Fill(rcK0);
1048 if(IsK0InvMass(lInvMassK0))selectedV0sAssoc->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,1,lDCAV0PVz,lV0cosPointAngle));
1051 if ((lactcut && V0EtaMax && lacutset)&&(MotherPdg == 3122 &&
1052 PosDaughterPdg== 2212 &&
1053 NegDaughterPdg== -211 &&
1056 fHistMCAssoLA->Fill(rcLa);
1057 if(IsLambdaInvMass(lInvMassLambda))selectedV0sAssoc->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,2,lDCAV0PVz,lV0cosPointAngle));
1060 if ((alactcut && V0EtaMax && alacutset)&&(MotherPdg == -3122 &&
1061 PosDaughterPdg== 211 &&
1062 NegDaughterPdg== -2212 &&
1065 fHistMCAssoALA->Fill(rcAl);
1066 if(IsLambdaInvMass(lInvMassAntiLambda))selectedV0sAssoc->Add(new V0Correlationparticle(lEtaV0s,lPhiV0s,lPtV0s,3,lDCAV0PVz,lV0cosPointAngle));
1069 if ((lactcut && V0EtaMax && lacutset)&&(MotherPdg == 3122 &&
1070 PosDaughterPdg== 2212 &&
1071 NegDaughterPdg== -211 &&
1072 (GrandMotherPdg==3322 ||GrandMotherPdg==3312)))
1074 fHistMCAssoLAXI->Fill(rcLa);
1077 if ((alactcut && V0EtaMax && alacutset)&&(MotherPdg == -3122 &&
1078 PosDaughterPdg== 211 &&
1079 NegDaughterPdg== -2212 &&
1080 GrandMotherPdg== -3312))
1082 fHistMCAssoALAXiPlus->Fill(rcAl);
1087 AliInfo(Form("No case selected"));
1093 FillCorrelationSibling(ltrackMultiplicity,selectedTracksLeading,selectedV0s,fHistTriggerSib,fHistReconstSib);
1094 FillCorrelationMixing(ltrackMultiplicity,tPrimaryVtxPosition[2],poolmax,poolmin,selectedTracksLeading,selectedV0s,fHistTriggerMix,fHistReconstMix);
1096 FillCorrelationSibling(ltrackMultiplicity,selectedTracksLeading,selectedV0sAssoc,fHistTriggerSibASO,fHistReconstSibASO);
1097 FillCorrelationMixing(ltrackMultiplicity,lPVz,poolmax,poolmin,selectedTracksLeading,selectedV0sAssoc,fHistTriggerMixASO,fHistReconstMixASO);
1099 PostData(1,fOutputList);
1101 //---------------------------------------------------------------------------------------
1102 void AliLeadingV0Correlation::Terminate(Option_t *)
1104 //No need in the grid
1106 //---------------------------------------------------------------------------------------
1107 Bool_t AliLeadingV0Correlation::IsAcseptedDaughterTrack(const AliAODTrack *itrack)
1109 if(fCase==1 || fCase==2)
1110 if(TMath::Abs(itrack->Eta())>fTrackEtaCut)return kFALSE;
1112 if (!itrack->IsOn(AliAODTrack::kTPCrefit)) return kFALSE;
1114 Float_t nCrossedRowsTPC = itrack->GetTPCClusterInfo(2,1);
1115 if (nCrossedRowsTPC < fTPCClusters) return kFALSE;
1117 Int_t findable=itrack->GetTPCNclsF();
1118 if (findable <= 0) return kFALSE;
1120 if (nCrossedRowsTPC/findable < fTPCfindratio) return kFALSE;
1123 //---------------------------------------------------------------------------------------
1124 Bool_t AliLeadingV0Correlation::IsAcseptedV0(const AliAODv0* aodV0, const AliAODTrack* myTrackPos, const AliAODTrack* myTrackNeg)
1126 if (!aodV0) return kFALSE;
1128 // Offline reconstructed V0 only
1129 if (aodV0->GetOnFlyStatus()) return kFALSE;
1131 // Get daughters and check them
1132 myTrackPos=(AliAODTrack *)(aodV0->GetDaughter(0));
1133 myTrackNeg=(AliAODTrack *)(aodV0->GetDaughter(1));
1135 if (!myTrackPos||!myTrackNeg) return kFALSE;
1136 // Unlike signs of daughters
1137 if (myTrackPos->Charge() == myTrackNeg->Charge()) return kFALSE;
1139 // Track cuts for daughers
1140 if ( !(IsAcseptedDaughterTrack(myTrackPos)) || !(IsAcseptedDaughterTrack(myTrackNeg)) ) return kFALSE;
1142 // Minimum pt of daughters
1143 Double_t lPtPos = myTrackPos->Pt();
1144 Double_t lPtNeg = myTrackNeg->Pt();
1146 if (lPtPos<fPtMin || lPtNeg<fPtMin) return kFALSE;
1150 //---------------------------------------------------------------------------------------
1151 Bool_t AliLeadingV0Correlation::IsAcseptedK0(Double_t v0rad,
1159 if(v0rad >=fV0radius &&
1160 dcaptp >=fV0PostoPVz &&
1161 dcantp >=fV0NegtoPVz &&
1162 dcav0d <=fDCAV0Daughters &&
1164 TMath::Abs(massLa - 1.115683) > fRejectLamK0 &&
1165 TMath::Abs(massALa - 1.115683) > fRejectLamK0 )return kTRUE;
1168 //---------------------------------------------------------------------------------------
1169 Bool_t AliLeadingV0Correlation::IsAcseptedLA(Double_t v0rad,
1176 if(v0rad >=fV0radius &&
1177 dcaptp >=fV0PostoPVz &&
1178 dcantp >=fV0NegtoPVz &&
1179 dcav0d <=fDCAV0Daughters &&
1181 TMath::Abs(massK0 - 0.4976) > fRejectK0Lam &&
1182 TMath::Abs(massK0 - 0.4976) > fRejectK0Lam )return kTRUE;
1185 //---------------------------------------------------------------------------------------
1186 Bool_t AliLeadingV0Correlation::IsK0InvMass(const Double_t mass) const
1188 const Float_t massK0 = 0.497;
1190 return ((massK0-fMassCutK0)<=mass && mass<=(massK0 + fMassCutK0))?1:0;
1192 //---------------------------------------------------------------------------------------
1193 Bool_t AliLeadingV0Correlation::IsLambdaInvMass(const Double_t mass) const
1195 const Float_t massLambda = 1.116;
1197 return ((massLambda-fMassCutLa)<=mass && mass<=(massLambda + fMassCutLa))?1:0;
1199 //---------------------------------------------------------------------------------------
1200 Double_t AliLeadingV0Correlation::RangePhi(Double_t DPhi)
1202 if (DPhi < -TMath::Pi()/2) DPhi += 2*TMath::Pi();
1203 if (DPhi > 3*TMath::Pi()/2) DPhi -= 2*TMath::Pi();
1206 //---------------------------------------------------------------------------------------
1207 Bool_t AliLeadingV0Correlation::IsTrackFromV0(AliAODTrack* track)
1209 Int_t atrID = track->GetID();
1211 for(int i=0; i<fAODEvent->GetNumberOfV0s(); i++){ // loop over V0s
1212 AliAODv0* aodV0 = fAODEvent->GetV0(i);
1214 AliAODTrack *trackPos=(AliAODTrack *)(aodV0->GetDaughter(0));
1215 AliAODTrack *trackNeg=(AliAODTrack *)(aodV0->GetDaughter(1));
1217 if ( !(IsAcseptedDaughterTrack(trackPos)) || !(IsAcseptedDaughterTrack(trackNeg)) ) continue;
1218 //----------------------------------
1219 Int_t negID = trackNeg->GetID();
1220 Int_t posID = trackPos->GetID();
1222 if ((TMath::Abs(negID)+1)==(TMath::Abs(atrID))){ return kTRUE;}
1223 if ((TMath::Abs(posID)+1)==(TMath::Abs(atrID))){ return kTRUE;}
1224 //----------------------------------
1228 //---------------------------------------------------------------------------------------
1229 void AliLeadingV0Correlation::FillCorrelationSibling(Double_t MultipOrCent,
1230 TObjArray*triggerArray,
1231 TObjArray*selectedV0Array,
1233 THnSparse*associateHist)
1235 Double_t binsv0CORR[8];
1236 Double_t binsTrigSib[2];
1237 Int_t counterSibMCA=0;
1239 for(Int_t i=0;i<triggerArray->GetEntriesFast();i++)
1241 AliAODTrack* trigger = (AliAODTrack*)triggerArray->At(0);
1242 if(!trigger)continue;
1245 if(IsTrackFromV0(trigger))continue;
1247 Double_t triggerPt = trigger->Pt();
1248 Double_t triggerPhi = trigger->Phi();
1249 Double_t triggerEta = trigger->Eta();
1251 if(triggerPt<fTriglow||triggerPt>fTrighigh)continue;
1254 if(counterSibMCA==triggerArray->GetEntriesFast()){
1256 binsTrigSib[0]=triggerPt;
1257 binsTrigSib[1]=MultipOrCent;
1259 if(triggerHist)triggerHist->Fill(binsTrigSib[0],triggerPhi,triggerEta);
1261 for (Int_t j=0; j<selectedV0Array->GetEntriesFast(); j++){
1263 V0Correlationparticle* associate = (V0Correlationparticle*) selectedV0Array->At(j);
1264 if(!associate)continue;
1266 binsv0CORR[0]= associate->Pt();
1267 binsv0CORR[1]= associate->Phi();
1268 binsv0CORR[2]= associate->Eta();
1270 if(binsv0CORR[0]>triggerPt) continue;
1272 binsv0CORR[3]=RangePhi(triggerPhi-binsv0CORR[1]);
1273 binsv0CORR[4]=triggerEta-binsv0CORR[2];
1274 binsv0CORR[5]= associate->WhichCandidate();
1275 binsv0CORR[6]= associate->DCAPostoP();
1276 binsv0CORR[7]= associate->DCANegtoP();
1278 associateHist->Fill(binsv0CORR);
1283 //---------------------------------------------------------------------------------------
1284 void AliLeadingV0Correlation::FillCorrelationMixing(Double_t MultipOrCentMix,
1288 TObjArray*triggerArray,
1289 TObjArray*selectedV0Array,
1291 THnSparse*associateHist)
1293 if(TMath::Abs(pvxMix)>=fpvzcut || MultipOrCentMix>poolmax || MultipOrCentMix < poolmin)
1295 if(fcollidingSys=="PP")AliInfo(Form("pp Event with Zvertex = %.2f cm and multiplicity = %.0f out of pool bounds, SKIPPING",pvxMix,MultipOrCentMix));
1299 Double_t binsv0CORRMix[8];
1300 Double_t binsTrigMix[2];
1301 Double_t counterMix=0;
1303 AliEventPool* pool = fPoolMgr->GetEventPool(MultipOrCentMix, pvxMix);
1304 if (!pool) AliFatal(Form("No pool found for centrality = %f, zVtx = %f", MultipOrCentMix, pvxMix));
1306 if (pool->IsReady() || pool->NTracksInPool() > fPoolMinNTracks || pool->GetCurrentNEvents() > fMinEventsToMix)
1308 Int_t nMix = pool->GetCurrentNEvents();
1309 for (Int_t jMix=0; jMix<nMix; jMix++){
1311 TObjArray* mixEvents = pool->GetEvent(jMix);
1312 for (Int_t i=0; i<triggerArray->GetEntriesFast(); i++){
1314 AliAODTrack* trig = (AliAODTrack*)triggerArray->At(0);
1318 if(IsTrackFromV0(trig))continue;
1320 Double_t trigPhi = trig->Phi();
1321 Double_t trigEta = trig->Eta();
1322 Double_t trigPt = trig->Pt();
1324 if(trigPt<fTriglow||trigPt>fTrighigh)continue;
1327 if(counterMix==triggerArray->GetEntriesFast()){
1329 binsTrigMix[0]=trigPt;
1330 binsTrigMix[1]=MultipOrCentMix;
1332 if(triggerHist)triggerHist->Fill(binsTrigMix[0]);
1334 for (Int_t j=0; j<mixEvents->GetEntriesFast(); j++){
1336 V0Correlationparticle* associate = (V0Correlationparticle*) mixEvents->At(j);
1337 if(!associate)continue;
1339 binsv0CORRMix[0]= associate->Pt();
1340 binsv0CORRMix[1]= associate->Phi();
1341 binsv0CORRMix[2]= associate->Eta();
1343 if(binsv0CORRMix[0]>trigPt) continue;
1345 binsv0CORRMix[3]=RangePhi(trigPhi-binsv0CORRMix[1]);
1346 binsv0CORRMix[4]=trigEta-binsv0CORRMix[2];
1347 binsv0CORRMix[5]=associate->WhichCandidate();
1348 binsv0CORRMix[6]=associate->DCAPostoP();
1349 binsv0CORRMix[7]=associate->DCANegtoP();
1351 associateHist->Fill(binsv0CORRMix);
1358 TObjArray* tracksClone = new TObjArray;
1359 tracksClone->SetOwner(kTRUE);
1361 for (Int_t i=0; i<selectedV0Array->GetEntriesFast(); i++)
1363 V0Correlationparticle* particle = (V0Correlationparticle*) selectedV0Array->At(i);
1364 tracksClone->Add(new V0Correlationparticle(particle->Eta(),
1367 particle->WhichCandidate(),
1368 particle->DCAPostoP(),
1369 particle->DCANegtoP()));
1371 pool->UpdatePool(tracksClone);
1373 //---------------------------------------------------------------------------------------