1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: Satyajit Jena. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
17 //=========================================================================//
19 // Analysis Task for Net-Charge Higher Moment Analysis //
20 // Author: Satyajit Jena || Nirbhay K. Behera //
21 // sjena@cern.ch || nbehera@cern.ch //
23 //=========================================================================//
31 #include "THnSparse.h"
34 #include "AliAnalysisTaskSE.h"
35 #include "AliAnalysisManager.h"
36 #include "AliAODEvent.h"
37 #include "AliVTrack.h"
38 #include "AliAODTrack.h"
39 #include "AliAODInputHandler.h"
40 #include "AliAODEvent.h"
42 #include "AliAODMCHeader.h"
43 #include "AliAODMCParticle.h"
44 #include "AliMCEventHandler.h"
45 #include "AliMCEvent.h"
47 #include "AliGenHijingEventHeader.h"
48 #include "AliGenEventHeader.h"
50 #include "AliAODPid.h"
51 #include "AliPIDResponse.h"
52 #include "AliAODpidUtil.h"
53 #include "AliPIDCombined.h"
54 #include "AliHelperPID.h"
56 #include "AliEbyEHigherMomentsTaskPID.h"
61 ClassImp(AliEbyEHigherMomentsTaskPID)
64 //-----------------------------------------------------------------------
65 AliEbyEHigherMomentsTaskPID::AliEbyEHigherMomentsTaskPID( const char *name )
66 : AliAnalysisTaskSE( name ),
70 fParticleSpecies(AliPID::kProton),
72 fCentralityEstimator("V0M"),
89 fTHnCentNplusNminusPid(0),
90 fTHnCentNplusNminusPidTruth(0),
91 fPtBinNplusNminusPid(0),
92 fPtBinNplusNminusPidTruth(0)
95 for ( Int_t i = 0; i < 4; i++) {
99 DefineOutput(1, TList::Class()); // Outpput...
102 AliEbyEHigherMomentsTaskPID::~AliEbyEHigherMomentsTaskPID(){
104 if(fListOfHistos) delete fListOfHistos;
105 if(fHelperPID) delete fHelperPID;
108 //---------------------------------------------------------------------------------
109 void AliEbyEHigherMomentsTaskPID::UserCreateOutputObjects() {
111 fListOfHistos = new TList();
112 fListOfHistos->SetOwner(kTRUE);
114 fEventCounter = new TH1D("fEventCounter","EventCounter", 10, 0.5,10.5);
115 fEventCounter->GetXaxis()->SetBinLabel(1,"Event Accesed");
116 fEventCounter->GetXaxis()->SetBinLabel(2,"Within 0-80% centrality");
117 fEventCounter->GetXaxis()->SetBinLabel(5,"Have a vertex");
118 fEventCounter->GetXaxis()->SetBinLabel(6,"After vertex Cut");
119 fEventCounter->GetXaxis()->SetBinLabel(7,"Event Analyzed");
120 fListOfHistos->Add(fEventCounter);
123 fHistQA[0] = new TH1D("fHistQAvz", "Histo Vz After Cut", 500, -25., 25.);
124 fHistQA[1] = new TH1D("fHistQAPt","p_{T} distribution",600,0,6);
125 fHistQA[2] = new TH1D("fHistQAEta","#eta distribution",240,-1.2,1.2);
126 fHistQA[3] = new TH1D("fHistQAPhi","#phi distribution",340,0,6.8);
128 for(Int_t i = 0; i < 4; i++)
130 fListOfHistos->Add(fHistQA[i]);
133 fHistVxVy = new TH2D("fHistVxVy","Vertex-x Vs Vertex-y", 200, -1., 1., 200, -1., 1.);
134 fListOfHistos->Add(fHistVxVy);
137 const Int_t nDim = 3;
138 const Int_t nPid = 5;
139 TString Species[nPid] = {"Electron","Muon","Pion","Kaon","Proton"};
141 Int_t fBins[nPid][nDim];
142 Double_t fMin[nPid][nDim];
143 Double_t fMax[nPid][nDim];
145 for( Int_t i = 0; i < nPid; i++ ){
146 for( Int_t j = 0; j < nDim; j++ ){
154 const Int_t dim = fNptBins*2 + 1;
155 Int_t bin[nPid][dim];
156 Double_t min[nPid][dim];
157 Double_t max[nPid][dim];
159 bin[2][0] = 100; min[2][0] = -0.5; max[2][0] = 99.5;
160 bin[3][0] = 100; min[3][0] = -0.5; max[3][0] = 99.5;
161 bin[4][0] = 100; min[4][0] = -0.5; max[4][0] = 99.5;
163 for(Int_t i = 1; i < dim; i++){
184 TString hname1, hname11;
185 TString htitle1, axisTitle1, axisTitle2;
187 Int_t gPid = (Int_t)fParticleSpecies;
190 fBins[2][0] = 100; fBins[2][1] = 1000; fBins[2][2] = 1000;
191 fMin[2][0] = -0.5; fMin[2][1] = -0.5; fMin[2][2] = -0.5;
192 fMax[2][0] = 99.5; fMax[2][1] = 999.5; fMax[2][2] = 999.5;
194 fBins[3][0] = 100; fBins[3][1] = 700; fBins[3][2] = 700;
195 fMin[3][0] = -0.5; fMin[3][1] = -0.5; fMin[3][2] = -0.5;
196 fMax[3][0] = 99.5; fMax[3][1] = 699.5; fMax[3][2] = 699.5;
198 fBins[4][0] = 100; fBins[4][1] = 400; fBins[4][2] = 400;
199 fMin[4][0] = -0.5; fMin[4][1] = -0.5; fMin[4][2] = -0.5;
200 fMax[4][0] = 99.5; fMax[4][1] = 399.5; fMax[4][2] = 399.5;
202 hname1 = "fCentNplusNminusPid"; hname1 +=gPid;
203 htitle1 = Species[gPid]; htitle1 +=" And Neg-"; htitle1 +=Species[gPid];
204 axisTitle1 ="Pos-"; axisTitle1 += Species[gPid];
205 axisTitle2 ="Neg-"; axisTitle2 += Species[gPid];
207 fTHnCentNplusNminusPid = new THnSparseD(hname1.Data(),htitle1.Data(), nDim, fBins[gPid], fMin[gPid],fMax[gPid]);
209 fTHnCentNplusNminusPid->GetAxis(0)->SetTitle("Centrality");
210 fTHnCentNplusNminusPid->GetAxis(1)->SetTitle(axisTitle1.Data());
211 fTHnCentNplusNminusPid->GetAxis(2)->SetTitle(axisTitle2.Data());
212 fListOfHistos->Add(fTHnCentNplusNminusPid);
214 TString hname2 = "fPtBinNplusNminusPid";
215 TString htitle2 = "cent-nplus-nminus-ptbinwise";
218 fPtBinNplusNminusPid = new THnSparseI(hname2.Data(),htitle2.Data(), dim, bin[gPid], min[gPid], max[gPid]);
219 fListOfHistos->Add(fPtBinNplusNminusPid);
224 if( fAnalysisType == "MCAOD" ){
226 hname11 = "fCentNplusNminusPidTruth"; hname11 +=gPid;
227 fTHnCentNplusNminusPidTruth = new THnSparseD(hname11.Data(),htitle1.Data(), nDim, fBins[gPid], fMin[gPid],fMax[gPid]);
229 fTHnCentNplusNminusPidTruth->GetAxis(0)->SetTitle("Centrality");
230 fTHnCentNplusNminusPidTruth->GetAxis(1)->SetTitle(axisTitle1.Data());
231 fTHnCentNplusNminusPidTruth->GetAxis(2)->SetTitle(axisTitle2.Data());
232 fListOfHistos->Add(fTHnCentNplusNminusPidTruth);
234 TString hname22 = "fPtBinNplusNminusPidTruth";
235 TString htitle22 = "cent-nplus-nminus-ptbinwise";
239 fPtBinNplusNminusPidTruth = new THnSparseI(hname22.Data(),htitle22.Data(), dim, bin[gPid], min[gPid], max[gPid]);
240 fListOfHistos->Add(fPtBinNplusNminusPidTruth);
244 PostData(1, fListOfHistos);
250 //----------------------------------------------------------------------------------
251 void AliEbyEHigherMomentsTaskPID::UserExec( Option_t * ){
253 fEventCounter->Fill(1);
255 if(fAnalysisType == "AOD") {
259 }//AOD--analysis-----
261 else if(fAnalysisType == "MCAOD") {
273 //--------------------------------------------------------------------------------------
274 void AliEbyEHigherMomentsTaskPID::doAODEvent(){
276 Double_t posPidSum = 0.;
277 Double_t negPidSum = 0.;
279 const Int_t dim = fNptBins*2;
282 for(Int_t idx = 0; idx < dim; idx++){
287 Int_t gPid = (Int_t)fParticleSpecies;
289 AliAnalysisManager* manager = AliAnalysisManager::GetAnalysisManager();
291 cout<<"ERROR: Analysis manager not found."<<endl;
294 //coneect to the inputHandler------------
295 AliAODInputHandler* inputHandler = dynamic_cast<AliAODInputHandler*> (manager->GetInputEventHandler());
297 cout<<"ERROR: Input handler not found."<<endl;
301 AliAODEvent *fAOD = dynamic_cast<AliAODEvent*>(InputEvent());
305 cout<< "ERROR: AOD not found " <<endl;
309 fPIDResponse = inputHandler->GetPIDResponse();
313 AliFatal("This Task needs the PID response attached to the inputHandler");
317 if(!ProperVertex(fAOD)) return;
319 AliAODHeader *aodHeader = fAOD->GetHeader();
321 fCentrality = (Int_t)aodHeader->GetCentralityP()->GetCentralityPercentile(fCentralityEstimator.Data());
323 if(fCentrality < 0 || fCentrality >= 81) return;
325 fEventCounter->Fill(2);
327 Int_t nTracks = fAOD->GetNumberOfTracks();
329 for(Int_t i = 0; i < nTracks; i++) {
331 AliAODTrack* aodTrack = dynamic_cast<AliAODTrack *>(fAOD->GetTrack(i));
334 AliError(Form("ERROR: Could not retrieve AODtrack %d",i));
338 if(!AcceptTrack(aodTrack) ) continue;
341 fBin = GetPtBin(aodTrack->Pt());
343 Short_t gCharge = aodTrack->Charge();
346 Double_t rap = aodTrack->Y(AliPID::ParticleMass(fParticleSpecies));
347 if( fabs(rap ) > fRapidityCut ) continue;//Rapidity cut
349 Int_t PID = fHelperPID->GetParticleSpecies((AliVTrack*)aodTrack, kFALSE);
351 if( (PID+2) == gPid ){
359 PtCh[fNptBins+fBin] += 1;
365 }//--------- Track Loop to select with filterbit
368 Double_t fContainerPid[3] = { static_cast<Double_t>(fCentrality), posPidSum, negPidSum};
369 fTHnCentNplusNminusPid->Fill(fContainerPid);
371 //cout << fCentrality <<" "<< posPidSum <<" " << negPidSum << endl;
373 Double_t ptContainer[dim+1];
375 ptContainer[0] = fCentrality;
377 for(Int_t i = 1; i <= dim; i++){
378 ptContainer[i] = PtCh[i-1];
379 //cout << PtCh[i-1] <<" ,";
383 fPtBinNplusNminusPid->Fill(ptContainer);
386 fEventCounter->Fill(7);
387 PostData(1, fListOfHistos);
391 //--------------------------------------------------------------------------------------
392 //--------------------------------------------------------------------------------------
393 void AliEbyEHigherMomentsTaskPID::doMCAODEvent(){
397 Double_t posPidSumMCRec = 0.;
398 Double_t negPidSumMCRec = 0.;
400 Double_t posPidSumMCTruth = 0.;
401 Double_t negPidSumMCTruth = 0.;
403 const Int_t dim = fNptBins*2;
405 Double_t ptChMC[dim];
407 for(Int_t idx = 0; idx < dim; idx++){
413 Int_t gPid = (Int_t)fParticleSpecies;
415 Int_t gPdgCode = AliPID::ParticleCode(fParticleSpecies);
417 //Connect to Anlaysis manager------
418 AliAnalysisManager* manager = AliAnalysisManager::GetAnalysisManager();
420 cout<<"ERROR: Analysis manager not found."<<endl;
424 AliAODInputHandler* inputHandler = dynamic_cast<AliAODInputHandler*> (manager->GetInputEventHandler());
426 cout<<"ERROR: Input handler not found."<<endl;
431 AliAODEvent *fAOD = dynamic_cast<AliAODEvent*>(InputEvent());
435 cout<< "ERROR: AOD not found " <<endl;
439 fPIDResponse =inputHandler->GetPIDResponse();
443 AliFatal("This Task needs the PID response attached to the inputHandler");
448 // ------------------------------------------------------------------
450 fArrayMC = dynamic_cast<TClonesArray*>(fAOD->FindListObject(AliAODMCParticle::StdBranchName()));
452 AliFatal("Error: MC particles branch not found!\n");
456 AliAODMCHeader *mcHdr=NULL;
457 mcHdr=(AliAODMCHeader*)fAOD->GetList()->FindObject(AliAODMCHeader::StdBranchName());
459 printf("MC header branch not found!\n");
463 if(!ProperVertex(fAOD)) return;
465 AliAODHeader *aodHeader = fAOD->GetHeader();
467 fCentrality = (Int_t)aodHeader->GetCentralityP()->GetCentralityPercentile(fCentralityEstimator.Data());
469 if( fCentrality < 0 || fCentrality >= 81) return;
471 fEventCounter->Fill(2);
473 Int_t nTracks = fAOD->GetNumberOfTracks();
475 for(Int_t i = 0; i < nTracks; i++) {
477 AliAODTrack* aodTrack = dynamic_cast<AliAODTrack *>(fAOD->GetTrack(i));
480 AliError(Form("ERROR: Could not retrieve AODtrack %d",i));
484 if(!AcceptTrack(aodTrack) ) continue;
486 fBin = GetPtBin(aodTrack->Pt());
488 //cout << "Pt Bin " << fBin << endl;
490 Short_t gCharge = aodTrack->Charge();
492 Double_t rap = aodTrack->Y(AliPID::ParticleMass(fParticleSpecies));
493 if( fabs(rap ) > fRapidityCut ) continue;//Rapidity cut
495 Int_t PID = fHelperPID->GetParticleSpecies((AliVTrack*)aodTrack, kFALSE);
497 if( (PID+2) == gPid ){
503 negPidSumMCRec += 1.;
504 PtCh[fNptBins+fBin] += 1;
508 }//--------- Track Loop to select with filterbit
510 //===============================================================================
511 //---------------------MC Truth--------------------------------------------------
512 //===============================================================================
514 Int_t nMCTrack = fArrayMC->GetEntriesFast();
516 for (Int_t iMC = 0; iMC < nMCTrack; iMC++) {
518 AliAODMCParticle *partMC = (AliAODMCParticle*) fArrayMC->At(iMC);
521 AliError(Form("ERROR: Could not retrieve AODMCtrack %d",iMC));
525 if(partMC->Charge() == 0) continue;
526 if(!partMC->IsPhysicalPrimary()) continue;
528 if (partMC->Eta() < fEtaLowerLimit || partMC->Eta() > fEtaHigherLimit) continue;
529 if (partMC->Pt() < fPtLowerLimit || partMC->Pt() > fPtHigherLimit) continue;
531 Short_t gCharge = partMC->Charge();
532 Int_t mcbin = GetPtBin(partMC->Pt());
534 Double_t rap = partMC->Y();
535 if( fabs(rap ) > fRapidityCut ) continue;//Rapidity cut
536 if(TMath::Abs(partMC->GetPdgCode()) != gPdgCode) continue;
539 posPidSumMCTruth += 1.;
543 negPidSumMCTruth += 1.;
544 ptChMC[fNptBins+mcbin] += 1;
548 }//MC-Truth Track loop--
551 //cout << fCentrality << " " << positiveSumMCRec << " " << negativeSumMCRec <<endl;
552 //cout <<" " << positiveSumMCTruth << " " << negativeSumMCTruth << endl;
553 //cout <<fCentrality<<" " << posPidSumMCRec << " " << negPidSumMCRec << endl;
554 //cout <<fCentrality<<" " << posPidSumMCTruth << " " << negPidSumMCTruth << endl;
555 //cout <<"---------------------------------" << endl;
557 Double_t fContainerPid[3] = { static_cast<Double_t>(fCentrality), posPidSumMCRec, negPidSumMCRec};//Reco.
558 Double_t fContainerPidTruth[3] = { static_cast<Double_t>(fCentrality), posPidSumMCTruth, negPidSumMCTruth};
560 fTHnCentNplusNminusPid->Fill(fContainerPid);//Fill the rec. pid tracks
561 fTHnCentNplusNminusPidTruth->Fill(fContainerPidTruth);//MC-Truth pid
563 Double_t ptContainer[dim+1];
564 Double_t ptContainerMC[dim+1];
565 ptContainer[0] = fCentrality;
566 ptContainerMC[0] = fCentrality;
568 for(Int_t i = 1; i <= dim; i++){
569 ptContainer[i] = PtCh[i-1];
570 ptContainerMC[i] = ptChMC[i-1];
571 //cout <<" "<< PtCh[i-1]<<endl;
572 //cout<< " Rec=" << PtCh[i-1];
577 fPtBinNplusNminusPid->Fill(ptContainer);
578 fPtBinNplusNminusPidTruth->Fill(ptContainerMC);
580 fEventCounter->Fill(7);
581 PostData(1, fListOfHistos);
586 //---------------------------------------------------------------------------------------
587 Bool_t AliEbyEHigherMomentsTaskPID::ProperVertex(AliAODEvent *fAOD) const{
589 Bool_t IsVtx = kFALSE;
591 const AliAODVertex *vertex = fAOD->GetPrimaryVertex();
595 vertex->GetCovarianceMatrix(fCov);
596 if(vertex->GetNContributors() > 0) {
599 Double_t lvx = vertex->GetX();
600 Double_t lvy = vertex->GetY();
601 Double_t lvz = vertex->GetZ();
603 fEventCounter->Fill(5);
605 if(TMath::Abs(lvx) < fVxMax) {
606 if(TMath::Abs(lvy) < fVyMax) {
607 if(TMath::Abs(lvz) < fVzMax) {
609 fEventCounter->Fill(6);
610 fHistQA[0]->Fill(lvz);
611 fHistVxVy->Fill(lvx,lvy);
618 }//Contributors check---
619 }//If vertex-----------
621 AliCentrality *centrality = fAOD->GetCentrality();
622 if (centrality->GetQuality() != 0) IsVtx = kFALSE;
626 //---------------------------------------------------------------------------------------
628 //---------------------------------------------------------------------------------------
629 Bool_t AliEbyEHigherMomentsTaskPID::AcceptTrack(AliAODTrack* track) const{
631 if(!track) return kFALSE;
632 if( track->Charge() == 0 ) return kFALSE;
634 Double_t pt = track->Pt();
635 Double_t eta = track->Eta();
636 if(!track->TestFilterBit(fAODtrackCutBit)) return kFALSE;
637 if( pt < fPtLowerLimit || pt > fPtHigherLimit ) return kFALSE;
638 if( eta < fEtaLowerLimit || eta > fEtaHigherLimit ) return kFALSE;
640 fHistQA[1]->Fill(pt);
641 fHistQA[2]->Fill(eta);
642 fHistQA[3]->Fill(track->Phi());
646 //------------------------------------------------------------------------
648 //------------------------------------------------------------------------
649 Int_t AliEbyEHigherMomentsTaskPID::GetPtBin(Double_t pt){
653 Double_t BinSize = (fPtHigherLimit - fPtLowerLimit)/fNptBins;
655 for(Int_t iBin = 0; iBin < fNptBins; iBin++){
657 Double_t xLow = fPtLowerLimit + iBin*BinSize;
658 Double_t xHigh = fPtLowerLimit + (iBin+1)*BinSize;
660 if( pt >= xLow && pt < xHigh){
662 //cout << "Pt Bin #" << bin <<"pt=" << pt << endl;
671 //------------------------------------------------------------------------
672 //------------------------------------------------------------------------
673 void AliEbyEHigherMomentsTaskPID::Terminate( Option_t * ){
675 Info("AliEbyEHigherMomentTask"," Task Successfully finished");
679 //------------------------------------------------------------------------