/************************************************************************** * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * * * * Author: Satyajit Jena. * * Contributors are mentioned in the code where appropriate. * * * * Permission to use, copy, modify and distribute this software and its * * documentation strictly for non-commercial purposes is hereby granted * * without fee, provided that the above copyright notice appears in all * * copies and that both the copyright notice and this permission notice * * appear in the supporting documentation. The authors make no claims * * about the suitability of this software for any purpose. It is * * provided "as is" without express or implied warranty. * **************************************************************************/ //=========================================================================// // // // Analysis Task for Net-Charge Higher Moment Analysis // // Author: Satyajit Jena || Nirbhay K. Behera // // sjena@cern.ch || nbehera@cern.ch // // Last Modified: 30/01/2014: only for net-charge part // //=========================================================================// #include "TChain.h" #include "TList.h" #include "TFile.h" #include "TTree.h" #include "TH1D.h" #include "TH2D.h" #include "TH3D.h" #include "THnSparse.h" #include "TCanvas.h" #include "AliAnalysisTaskSE.h" #include "AliAnalysisManager.h" #include "AliAODEvent.h" #include "AliVTrack.h" #include "AliAODTrack.h" #include "AliAODInputHandler.h" #include "AliAODEvent.h" #include "AliStack.h" #include "AliAODMCHeader.h" #include "AliAODMCParticle.h" #include "AliMCEventHandler.h" #include "AliMCEvent.h" #include "AliStack.h" #include "AliGenHijingEventHeader.h" #include "AliGenEventHeader.h" #include "AliPID.h" #include "AliAODPid.h" #include "AliPIDResponse.h" #include "AliAODpidUtil.h" #include "AliPIDCombined.h" #include "AliHelperPID.h" #include "AliEbyEHigherMomentsTask.h" using std::cout; using std::endl; ClassImp(AliEbyEHigherMomentsTask) //----------------------------------------------------------------------- AliEbyEHigherMomentsTask::AliEbyEHigherMomentsTask( const char *name ) : AliAnalysisTaskSE( name ), fListOfHistos(0), fArrayMC(0), fAnalysisType(0), fCentralityEstimator("V0M"), fCentrality(0), fVxMax(3.), fVyMax(3.), fVzMax(10.), fPtLowerLimit(0.3), fPtHigherLimit(1.5), fNptBins(0), fBin(0), fEtaLowerLimit(-0.8), fEtaHigherLimit(0.8), fAODtrackCutBit(128), fEventCounter(0), fHistVxVy(0), fTHnCentNplusNminusCh(0), fTHnCentNplusNminusChTruth(0), fPtBinNplusNminusCh(0), fPtBinNplusNminusChTruth(0) { for ( Int_t i = 0; i < 4; i++) { fHistQA[i] = NULL; } DefineOutput(1, TList::Class()); // Outpput.... //DefineOutput(2, TList::Class()); } AliEbyEHigherMomentsTask::~AliEbyEHigherMomentsTask() { if(fListOfHistos) delete fListOfHistos; } //--------------------------------------------------------------------------------- void AliEbyEHigherMomentsTask::UserCreateOutputObjects() { fListOfHistos = new TList(); fListOfHistos->SetOwner(kTRUE); fEventCounter = new TH1D("fEventCounter","EventCounter", 10, 0.5,10.5); fEventCounter->GetXaxis()->SetBinLabel(1,"Event Accesed"); fEventCounter->GetXaxis()->SetBinLabel(2,"Within 0-80% centrality"); fEventCounter->GetXaxis()->SetBinLabel(5,"Have a vertex"); fEventCounter->GetXaxis()->SetBinLabel(6,"After vertex Cut"); fEventCounter->GetXaxis()->SetBinLabel(7,"Event Analyzed"); fListOfHistos->Add(fEventCounter); //For QA-Histograms fHistQA[0] = new TH1D("fHistQAvz", "Histo Vz After Cut", 500, -25., 25.); fHistQA[1] = new TH1D("fHistQAPt","p_{T} distribution",600,0,6); fHistQA[2] = new TH1D("fHistQAEta","#eta distribution",240,-1.2,1.2); fHistQA[3] = new TH1D("fHistQAPhi","#phi distribution",340,0,6.8); for(Int_t i = 0; i < 4; i++) { fListOfHistos->Add(fHistQA[i]); } fHistVxVy = new TH2D("fHistVxVy","Vertex-x Vs Vertex-y", 200, -1., 1., 200, -1., 1.); fListOfHistos->Add(fHistVxVy); const Int_t nDim = 3; Int_t fBinsCh[nDim] = {100, 1500, 1500}; Double_t fMinCh[nDim] = { -0.5, -0.5, -0.5 }; Double_t fMaxCh[nDim] = { 99.5, 1499.5, 1499.5}; const Int_t dim = fNptBins*2 + 1; //const Int_t dim = ; Int_t bin[dim]; Double_t min[dim]; Double_t max[dim]; bin[0] = 100; min[0] = -0.5; max[0] = 99.5; for(Int_t i = 1; i < dim; i++){ bin[i] = 900; min[i] = -0.5; max[i] = 899.5; } fTHnCentNplusNminusCh = new THnSparseD("fTHnCentNplusNminusCh","Cent-NplusChrg-NminusChrg", nDim, fBinsCh, fMinCh, fMaxCh); fTHnCentNplusNminusCh->GetAxis(0)->SetTitle("Centrality"); fTHnCentNplusNminusCh->GetAxis(1)->SetTitle("Nplus"); fTHnCentNplusNminusCh->GetAxis(2)->SetTitle("Nminus"); fListOfHistos->Add(fTHnCentNplusNminusCh); fPtBinNplusNminusCh = new THnSparseI("fPtBinNplusNminusCh","cent-nplus-nminus", dim, bin, min, max); fListOfHistos->Add(fPtBinNplusNminusCh); if( fAnalysisType == "MCAOD" ){ fTHnCentNplusNminusChTruth = new THnSparseD("fTHnCentNplusNminusChTruth","Cent-NplusChrg-NminusChrg", nDim, fBinsCh, fMinCh, fMaxCh); fTHnCentNplusNminusChTruth->GetAxis(0)->SetTitle("Centrality"); fTHnCentNplusNminusChTruth->GetAxis(1)->SetTitle("Nplus"); fTHnCentNplusNminusChTruth->GetAxis(2)->SetTitle("Nminus"); fListOfHistos->Add(fTHnCentNplusNminusChTruth); fPtBinNplusNminusChTruth = new THnSparseI("fPtBinNplusNminusChTruth","cent-nplus-nminus", dim, bin, min, max); fListOfHistos->Add(fPtBinNplusNminusChTruth); }//MCAOD--- //PostData(1, fListOfHistosQA); PostData(1, fListOfHistos); } //---------------------------------------------------------------------------------- void AliEbyEHigherMomentsTask::UserExec( Option_t * ){ fEventCounter->Fill(1); if(fAnalysisType == "AOD") { doAODEvent(); }//AOD--analysis----- else if(fAnalysisType == "MCAOD") { doMCAODEvent(); } else return; } //-------------------------------------------------------------------------------------- void AliEbyEHigherMomentsTask::doAODEvent(){ Double_t positiveSum = 0.; Double_t negativeSum = 0.; const Int_t dim = fNptBins*2; Double_t PtCh[dim]; for(Int_t idx = 0; idx < dim; idx++){ PtCh[idx] = 0.; } AliAnalysisManager* manager = AliAnalysisManager::GetAnalysisManager(); if (!manager) { cout<<"ERROR: Analysis manager not found."< (manager->GetInputEventHandler()); if (!inputHandler) { cout<<"ERROR: Input handler not found."<(InputEvent()); if (!fAOD) { cout<< "ERROR: AOD not found " <GetHeader(); fCentrality = (Int_t)aodHeader->GetCentralityP()->GetCentralityPercentile(fCentralityEstimator.Data()); if(fCentrality < 0 || fCentrality >= 81) return; fEventCounter->Fill(2); Int_t nTracks = fAOD->GetNumberOfTracks(); for(Int_t i = 0; i < nTracks; i++) { AliAODTrack* aodTrack = dynamic_cast(fAOD->GetTrack(i)); if(!aodTrack) { AliError(Form("ERROR: Could not retrieve AODtrack %d",i)); continue; } if(!AcceptTrack(aodTrack) ) continue; fBin = GetPtBin(aodTrack->Pt()); Short_t gCharge = aodTrack->Charge(); if(gCharge > 0){ positiveSum += 1.; PtCh[fBin] += 1; } if(gCharge < 0){ negativeSum += 1.; PtCh[fNptBins+fBin] += 1; } }//--------- Track Loop to select with filterbit Double_t fContainerCh[3] = { static_cast(fCentrality), positiveSum, negativeSum}; fTHnCentNplusNminusCh->Fill(fContainerCh); //cout << fCentrality << " " << positiveSum <<" " << negativeSum << endl; Double_t ptContainer[dim+1]; ptContainer[0] = fCentrality; for(Int_t i = 1; i <= dim; i++){ ptContainer[i] = PtCh[i-1]; //cout << PtCh[i-1] <<" ,"; } //cout << endl; fPtBinNplusNminusCh->Fill(ptContainer); fEventCounter->Fill(7); PostData(1, fListOfHistos); return; } //-------------------------------------------------------------------------------------- //-------------------------------------------------------------------------------------- void AliEbyEHigherMomentsTask::doMCAODEvent(){ //--------- Double_t positiveSumMCRec = 0.; Double_t negativeSumMCRec = 0.; Double_t positiveSumMCTruth = 0.; Double_t negativeSumMCTruth = 0.; const Int_t dim = fNptBins*2; Double_t PtCh[dim]; Double_t ptChMC[dim]; for(Int_t idx = 0; idx < dim; idx++){ PtCh[idx] = 0; ptChMC[idx] = 0; } //Connect to Anlaysis manager------ AliAnalysisManager* manager = AliAnalysisManager::GetAnalysisManager(); if (!manager) { cout<<"ERROR: Analysis manager not found."< (manager->GetInputEventHandler()); if (!inputHandler) { cout<<"ERROR: Input handler not found."<(InputEvent()); if (!fAOD) { cout<< "ERROR: AOD not found " <(fAOD->FindListObject(AliAODMCParticle::StdBranchName())); if (!fArrayMC) { AliFatal("Error: MC particles branch not found!\n"); return; } AliAODMCHeader *mcHdr=NULL; mcHdr=(AliAODMCHeader*)fAOD->GetList()->FindObject(AliAODMCHeader::StdBranchName()); if(!mcHdr) { printf("MC header branch not found!\n"); return; } if(!ProperVertex(fAOD)) return; AliAODHeader *aodHeader = fAOD->GetHeader(); fCentrality = (Int_t)aodHeader->GetCentralityP()->GetCentralityPercentile(fCentralityEstimator.Data()); if( fCentrality < 0 || fCentrality >= 81) return; fEventCounter->Fill(2); Int_t nTracks = fAOD->GetNumberOfTracks(); for(Int_t i = 0; i < nTracks; i++) { AliAODTrack* aodTrack = dynamic_cast(fAOD->GetTrack(i)); if(!aodTrack) { AliError(Form("ERROR: Could not retrieve AODtrack %d",i)); continue; } if(!AcceptTrack(aodTrack) ) continue; fBin = GetPtBin(aodTrack->Pt()); //cout << "Pt Bin " << fBin << endl; Short_t gCharge = aodTrack->Charge(); if(gCharge > 0){ positiveSumMCRec += 1.; PtCh[fBin] += 1; } if(gCharge < 0){ negativeSumMCRec += 1.; PtCh[fNptBins+fBin] += 1; } }//--------- Track Loop to select with filterbit //=============================================================================== //---------------------MC Truth-------------------------------------------------- //=============================================================================== Int_t nMCTrack = fArrayMC->GetEntriesFast(); for (Int_t iMC = 0; iMC < nMCTrack; iMC++) { AliAODMCParticle *partMC = (AliAODMCParticle*) fArrayMC->At(iMC); if(!partMC){ AliError(Form("ERROR: Could not retrieve AODMCtrack %d",iMC)); continue; } if(partMC->Charge() == 0) continue; if(!partMC->IsPhysicalPrimary()) continue; if (partMC->Eta() < fEtaLowerLimit || partMC->Eta() > fEtaHigherLimit) continue; if (partMC->Pt() < fPtLowerLimit || partMC->Pt() > fPtHigherLimit) continue; Short_t gCharge = partMC->Charge(); Int_t mcbin = GetPtBin(partMC->Pt()); if(gCharge > 0){ positiveSumMCTruth += 1.; ptChMC[mcbin] += 1.; } if(gCharge < 0){ negativeSumMCTruth += 1.; ptChMC[fNptBins+mcbin] += 1.; } }//MC-Truth Track loop-- Double_t fContainerCh[3] = { static_cast(fCentrality), positiveSumMCRec, negativeSumMCRec};//Reco. values ch. hadrons Double_t fContainerChTruth[3] = { static_cast(fCentrality), positiveSumMCTruth, negativeSumMCTruth}; fTHnCentNplusNminusCh->Fill(fContainerCh);//Fill the rec. ch. particles--- fTHnCentNplusNminusChTruth->Fill(fContainerChTruth);//MC -Truth ch. particles //cout << fCentrality << " " << positiveSumMCRec << " " << negativeSumMCRec <