]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Updated histos with multiplicity and centrality
authorprino <prino@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 28 Jun 2012 13:48:42 +0000 (13:48 +0000)
committerprino <prino@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 28 Jun 2012 13:48:42 +0000 (13:48 +0000)
PWGHF/vertexingHF/AliAnalysisTaskSEDplus.cxx
PWGHF/vertexingHF/AliAnalysisTaskSEDplus.h

index 44500b5ee05f30026238b0973021e4b0caafb665..ae739b3b6270b12462d67387d0fae31d24629ac5 100644 (file)
@@ -422,9 +422,9 @@ void AliAnalysisTaskSEDplus::UserCreateOutputObjects()
   TString hisname;
   Int_t index=0;
   Int_t nbins=GetNBinsHistos();
-  fHistCentrality[0]=new TH1F("centrality","centrality",100,0.5,30000.5);
-  fHistCentrality[1]=new TH1F("centrality(selectedCent)","centrality(selectedCent)",100,0.5,30000.5);
-  fHistCentrality[2]=new TH1F("centrality(OutofCent)","centrality(OutofCent)",100,0.5,30000.5);
+  fHistCentrality[0]=new TH2F("hCentrMult","centrality",100,0.5,30000.5,40,0.,100.);
+  fHistCentrality[1]=new TH2F("hCentrMult(selectedCent)","centrality(selectedCent)",100,0.5,30000.5,40,0.,100.);
+  fHistCentrality[2]=new TH2F("hCentrMult(OutofCent)","centrality(OutofCent)",100,0.5,30000.5,40,0.,100.);
   for(Int_t i=0;i<3;i++){
     fHistCentrality[i]->Sumw2();
     fOutput->Add(fHistCentrality[i]);
@@ -741,11 +741,12 @@ void AliAnalysisTaskSEDplus::UserExec(Option_t */*option*/)
 
   //Event selection
   Bool_t isEvSel=fRDCutsAnalysis->IsEventSelected(aod);
-  Float_t ntracks=aod->GetNTracks();//fRDCutsAnalysis->GetCentrality(aod);
-  fHistCentrality[0]->Fill(ntracks);
+  Float_t ntracks=aod->GetNTracks();
+  Float_t evCentr=fRDCutsAnalysis->GetCentrality(aod);
+  fHistCentrality[0]->Fill(ntracks,evCentr);
   if(fRDCutsAnalysis->GetWhyRejection()==5) fHistNEvents->Fill(2);
   if(fRDCutsAnalysis->GetWhyRejection()==1) fHistNEvents->Fill(3); 
-  if(fRDCutsAnalysis->GetWhyRejection()==2){fHistNEvents->Fill(4);fHistCentrality[2]->Fill(ntracks);}
+  if(fRDCutsAnalysis->GetWhyRejection()==2){fHistNEvents->Fill(4);fHistCentrality[2]->Fill(ntracks,evCentr);}
   if(fRDCutsAnalysis->GetWhyRejection()==6)fHistNEvents->Fill(5);
   if(fRDCutsAnalysis->GetWhyRejection()==7)fHistNEvents->Fill(6);
 
@@ -758,7 +759,7 @@ void AliAnalysisTaskSEDplus::UserExec(Option_t */*option*/)
   // printf("ntracklet===%d\n",tracklets);
   fSPDMult->Fill(tracklets);
 
-  fHistCentrality[1]->Fill(ntracks);
+  fHistCentrality[1]->Fill(ntracks,evCentr);
   fHistNEvents->Fill(1);
 
   TClonesArray *arrayMC=0;
index f63d8ded7db8ffc6ad5dafee2361b5b390f9dadd..2467114c49241db8bddb3838f8209d292e1cf265 100644 (file)
@@ -116,7 +116,7 @@ class AliAnalysisTaskSEDplus : public AliAnalysisTaskSE
   TH1F *fDCAHistLS[3*kMaxPtBins];//!hist. for LS cuts variable 6 (LC)
   TH1F *fMassHistLSTC[5*kMaxPtBins];//!hist. for LS inv mass (TC)
   TH2F *fCorreld0Kd0pi[3]; //!hist. for d0k*d0pi vs. d0k*d0pi (LC)
-  TH1F *fHistCentrality[3];//!hist. for cent distr (all,sel ev, )
+  TH2F *fHistCentrality[3];//!hist. for cent distr (all,sel ev, )
   THnSparseF *fHistMassPtImpParTC[5];//! histograms for impact paramter studies
     TH2F *fPtVsMass;    //! hist. of pt vs. mass (prod. cuts)
   TH2F *fPtVsMassTC;  //! hist. of pt vs. mass (analysis cuts)
@@ -146,7 +146,7 @@ class AliAnalysisTaskSEDplus : public AliAnalysisTaskSE
   Int_t  fDoLS;        // flag to do LS analysis
   Int_t fSystem;   //0=pp,1=PbPb
   
-  ClassDef(AliAnalysisTaskSEDplus,18); // AliAnalysisTaskSE for the MC association of heavy-flavour decay candidates
+  ClassDef(AliAnalysisTaskSEDplus,19); // AliAnalysisTaskSE for the MC association of heavy-flavour decay candidates
 };
 
 #endif