]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Update to select in centrality (Giacomo)
authordainese <dainese@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 24 Jan 2011 22:58:59 +0000 (22:58 +0000)
committerdainese <dainese@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 24 Jan 2011 22:58:59 +0000 (22:58 +0000)
PWG3/vertexingHF/AliAnalysisTaskSEDplus.cxx
PWG3/vertexingHF/AliAnalysisTaskSEDplus.h

index 36dad4b714e0bffc9072881b2c39bfe04d8c2938..a9d0feae850d2a1b902e14570e16805375807066 100644 (file)
@@ -136,6 +136,9 @@ AliAnalysisTaskSEDplus::~AliAnalysisTaskSEDplus()
     delete fHistNEvents;
     fHistNEvents=0;
   }  
+  for(Int_t i=0;i<3;i++){
+    if(fHistCentrality[i]){delete fHistCentrality[i]; fHistCentrality[i]=0;}
+  }
   for(Int_t i=0;i<3*fNPtBins;i++){
     if(fMassHist[i]){ delete fMassHist[i]; fMassHist[i]=0;}
     if(fCosPHist[i]){ delete fCosPHist[i]; fCosPHist[i]=0;}
@@ -450,6 +453,13 @@ void AliAnalysisTaskSEDplus::UserCreateOutputObjects()
   Int_t index=0;
   Int_t indexLS=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);
+  for(Int_t i=0;i<3;i++){
+    fHistCentrality[i]->Sumw2();
+    fOutput->Add(fHistCentrality[i]);
+  }
   for(Int_t i=0;i<fNPtBins;i++){
 
     index=GetHistoIndex(i);
@@ -697,7 +707,7 @@ void AliAnalysisTaskSEDplus::UserCreateOutputObjects()
   }
   
   
-  fHistNEvents = new TH1F("fHistNEvents", "number of events ",7,-0.5,6.5);
+  fHistNEvents = new TH1F("fHistNEvents", "number of events ",8,-0.5,7.5);
   fHistNEvents->GetXaxis()->SetBinLabel(1,"nEventsAnal");
   fHistNEvents->GetXaxis()->SetBinLabel(2,"nEvents with good vertex");
   fHistNEvents->GetXaxis()->SetBinLabel(3,"nEvents with PbPb HM trigger");
@@ -705,7 +715,8 @@ void AliAnalysisTaskSEDplus::UserCreateOutputObjects()
   fHistNEvents->GetXaxis()->SetBinLabel(5,"no. of candidate");
   fHistNEvents->GetXaxis()->SetBinLabel(6,"no. of D+ after loose cuts");
   fHistNEvents->GetXaxis()->SetBinLabel(7,"no. of D+ after tight cuts");
+  fHistNEvents->GetXaxis()->SetBinLabel(8,"no. of out centrality events");
+
   fHistNEvents->GetXaxis()->SetNdivisions(1,kFALSE);
   
   fHistNEvents->Sumw2();
@@ -790,15 +801,23 @@ void AliAnalysisTaskSEDplus::UserExec(Option_t */*option*/)
   if(!aod->GetPrimaryVertex()||TMath::Abs(aod->GetMagneticField())<0.001) return;
   fCounter->StoreEvent(aod,fReadMC);
   fHistNEvents->Fill(0); // count event
-  // Post the data already here
-  if(fRDCutsAnalysis->IsEventSelected(aod))fHistNEvents->Fill(1);
+
+  Bool_t isEvSel=fRDCutsAnalysis->IsEventSelected(aod);
+  Float_t centrality=aod->GetNTracks();//fRDCutsAnalysis->GetCentrality(aod);
+  fHistCentrality[0]->Fill(centrality);
   // trigger class for PbPb C0SMH-B-NOPF-ALLNOTRD
   TString trigclass=aod->GetFiredTriggerClasses();
   if(trigclass.Contains("C0SMH-B-NOPF-ALLNOTRD")||trigclass.Contains("C0SMH-B-NOPF-ALL")) fHistNEvents->Fill(2);
   if(fRDCutsAnalysis->GetWhyRejection()==1)fHistNEvents->Fill(3); 
-  
+  if(fRDCutsAnalysis->GetWhyRejection()==2){fHistNEvents->Fill(7);fHistCentrality[2]->Fill(centrality);}
+
+  // Post the data already here  
   PostData(1,fOutput);
-  
+  if(!isEvSel)return;
+
+  fHistCentrality[1]->Fill(centrality);
+  fHistNEvents->Fill(1);
+
   TClonesArray *arrayMC=0;
   AliAODMCHeader *mcHeader=0;
 
@@ -1079,172 +1098,17 @@ void AliAnalysisTaskSEDplus::Terminate(Option_t */*option*/)
     return;
   }
   fHistNEvents = dynamic_cast<TH1F*>(fOutput->FindObject("fHistNEvents"));
-  fYVsPt = dynamic_cast<TH2F*>(fOutput->FindObject("hYVsPt"));
-  fYVsPtTC = dynamic_cast<TH2F*>(fOutput->FindObject("hYVsPtTC"));
-  fYVsPtSig = dynamic_cast<TH2F*>(fOutput->FindObject("hYVsPtSig"));
-  fYVsPtSigTC = dynamic_cast<TH2F*>(fOutput->FindObject("hYVsPtSigTC"));
-  fPtVsMass = dynamic_cast<TH2F*>(fOutput->FindObject("hPtVsMass"));
-  fPtVsMassTC = dynamic_cast<TH2F*>(fOutput->FindObject("hPtVsMassTC"));
 
   TString hisname;
   Int_t index=0;
 
-  Int_t indexLS=0;
   for(Int_t i=0;i<fNPtBins;i++){
     index=GetHistoIndex(i);
-    if(fDoLS)indexLS=GetLSHistoIndex(i);
-    hisname.Form("hMassPt%d",i);
-    fMassHist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hCosPAllPt%d",i);
-    fCosPHist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hDLenAllPt%d",i);
-    fDLenHist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hSumd02AllPt%d",i);
-    fSumd02Hist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hSigVertAllPt%d",i);
-    fSigVertHist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hPtMaxAllPt%d",i);
-    fPtMaxHist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hDCAAllPt%d",i);
-    fDCAHist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hMassPt%dTC",i);
-    fMassHistTC[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hMassPt%dTCPlus",i);
-    fMassHistTCPlus[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hMassPt%dTCMinus",i);
-    fMassHistTCMinus[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    if(fDoLS){
-      hisname.Form("hLSPt%dLC",i);
-      fMassHistLS[indexLS]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      hisname.Form("hCosPAllPt%dLS",i);
-      fCosPHistLS[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      hisname.Form("hDLenAllPt%dLS",i);
-      fDLenHistLS[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      hisname.Form("hSumd02AllPt%dLS",i);
-      fSumd02HistLS[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      hisname.Form("hSigVertAllPt%dLS",i);
-      fSigVertHistLS[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      hisname.Form("hPtMaxAllPt%dLS",i);
-      fPtMaxHistLS[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      hisname.Form("hDCAAllPt%dLS",i);
-      fDCAHistLS[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      
-      hisname.Form("hLSPt%dTC",i);
-      fMassHistLSTC[indexLS]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      
-    } 
-    
-    index=GetSignalHistoIndex(i);    
-    if(fDoLS)indexLS++;
-    hisname.Form("hSigPt%d",i);
-    fMassHist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hCosPSigPt%d",i);
-    fCosPHist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hDLenSigPt%d",i);
-    fDLenHist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hSumd02SigPt%d",i);
-    fSumd02Hist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hSigVertSigPt%d",i);
-    fSigVertHist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hPtMaxSigPt%d",i);
-    fPtMaxHist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hDCASigPt%d",i);
-    fDCAHist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
     
-    hisname.Form("hSigPt%dTC",i);
+    hisname.Form("hMassPt%dTC",i);
     fMassHistTC[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hSigPt%dTCPlus",i);
-    fMassHistTCPlus[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hSigPt%dTCMinus",i);
-    fMassHistTCMinus[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    if(fDoLS){
-      hisname.Form("hLSPt%dLCnw",i);
-      fMassHistLS[indexLS]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      hisname.Form("hCosPSigPt%dLS",i);
-      fCosPHistLS[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      hisname.Form("hDLenSigPt%dLS",i);
-      fDLenHistLS[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      hisname.Form("hSumd02SigPt%dLS",i);
-      fSumd02HistLS[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      hisname.Form("hSigVertSigPt%dLS",i);
-      fSigVertHistLS[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      hisname.Form("hPtMaxSigPt%dLS",i);
-      fPtMaxHistLS[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      hisname.Form("hDCASigPt%dLS",i);
-      fDCAHistLS[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-
-      hisname.Form("hLSPt%dTCnw",i);
-      fMassHistLSTC[indexLS]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    }
+  } 
     
-    index=GetBackgroundHistoIndex(i); 
-    if(fDoLS)indexLS++;
-    hisname.Form("hBkgPt%d",i);
-    fMassHist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hCosPBkgPt%d",i);
-    fCosPHist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hDLenBkgPt%d",i);
-    fDLenHist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hSumd02BkgPt%d",i);
-    fSumd02Hist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hSigVertBkgPt%d",i);
-    fSigVertHist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hPtMaxBkgPt%d",i);
-    fPtMaxHist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hDCABkgPt%d",i);
-    fDCAHist[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hBkgPt%dTC",i);
-    fMassHistTC[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hBkgPt%dTCPlus",i);
-    fMassHistTCPlus[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    hisname.Form("hBkgPt%dTCMinus",i);
-    fMassHistTCMinus[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    if(fDoLS){
-      hisname.Form("hLSPt%dLCntrip",i);
-      fMassHistLS[indexLS]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-
-      hisname.Form("hCosPBkgPt%dLS",i);
-      fCosPHistLS[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      hisname.Form("hDLenBkgPt%dLS",i);
-      fDLenHistLS[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      hisname.Form("hSumd02BkgPt%dLS",i);
-      fSumd02HistLS[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      hisname.Form("hSigVertBkgPt%dLS",i);
-      fSigVertHistLS[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      hisname.Form("hPtMaxBkgPt%dLS",i);
-      fPtMaxHistLS[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      hisname.Form("hDCABkgPt%dLS",i);
-      fDCAHistLS[index]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      
-      hisname.Form("hLSPt%dTCntrip",i);
-      fMassHistLSTC[indexLS]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      
-      
-      indexLS++;
-      hisname.Form("hLSPt%dLCntripsinglecut",i);
-      fMassHistLS[indexLS]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      
-      hisname.Form("hLSPt%dTCntripsinglecut",i);
-      fMassHistLSTC[indexLS]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      
-      
-      indexLS++;
-      hisname.Form("hLSPt%dLCspc",i);
-      fMassHistLS[indexLS]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-      
-      hisname.Form("hLSPt%dTCspc",i);
-      fMassHistLSTC[indexLS]=dynamic_cast<TH1F*>(fOutput->FindObject(hisname.Data()));
-    }
-  }
-  fCounter = dynamic_cast<AliNormalizationCounter*>(GetOutputData(3)); 
-
-  if(fFillNtuple){
-    fNtupleDplus = dynamic_cast<TNtuple*>(GetOutputData(4));
-  }
-
   TCanvas *c1=new TCanvas("c1","D+ invariant mass distribution",500,500);
   c1->cd();
   TH1F *hMassPt=(TH1F*)fOutput->FindObject("hMassPt3TC");
index 53ac3bea79a5e722b67381f8c7797c20abe0d1a4..8f8b0d054c8aac10d8c43bf105ba953d13cf5747 100644 (file)
@@ -89,6 +89,7 @@ class AliAnalysisTaskSEDplus : public AliAnalysisTaskSE
   TH1F *fPtMaxHistLS[3*kMaxPtBins];//!hist. for LS cuts variable 5 (LC)
   TH1F *fDCAHistLS[3*kMaxPtBins];//!hist. for LS cuts variable 6 (LC)
   TH1F *fMassHistLSTC[5*kMaxPtBins];//!hist. for LS inv mass (TC)
+  TH1F *fHistCentrality[3];//!hist. for cent distr (all,sel ev, )
   TH2F *fPtVsMass;    //! hist. of pt vs. mass (prod. cuts)
   TH2F *fPtVsMassTC;  //! hist. of pt vs. mass (analysis cuts)
   TH2F *fYVsPt;       //! hist. of Y vs. Pt (prod. cuts)
@@ -110,7 +111,7 @@ class AliAnalysisTaskSEDplus : public AliAnalysisTaskSE
   Bool_t fUseStrangeness;//flag to enhance strangeness in MC to fit to data
   Bool_t fDoLS;      //flag to do LS analysis
   
-  ClassDef(AliAnalysisTaskSEDplus,9); // AliAnalysisTaskSE for the MC association of heavy-flavour decay candidates
+  ClassDef(AliAnalysisTaskSEDplus,10); // AliAnalysisTaskSE for the MC association of heavy-flavour decay candidates
 };
 
 #endif