]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Added histos to count D+ and D- separately (Francesco)
authordainese <dainese@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 30 Sep 2010 17:05:00 +0000 (17:05 +0000)
committerdainese <dainese@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 30 Sep 2010 17:05:00 +0000 (17:05 +0000)
PWG3/vertexingHF/AliAnalysisTaskSEDplus.cxx
PWG3/vertexingHF/AliAnalysisTaskSEDplus.h

index eda6dea7efcd71dad54ef6e86c810b81b2f2f308..588b3eed75d96f83fa6f4346027d0416b5d8655d 100644 (file)
@@ -137,6 +137,8 @@ AliAnalysisTaskSEDplus::~AliAnalysisTaskSEDplus()
     if(fPtMaxHist[i]){ delete fPtMaxHist[i]; fPtMaxHist[i]=0;}
     if(fDCAHist[i]){ delete fDCAHist[i]; fDCAHist[i]=0;}
     if(fMassHistTC[i]){ delete fMassHistTC[i]; fMassHistTC[i]=0;}
+    if(fMassHistTCPlus[i]){ delete fMassHistTCPlus[i]; fMassHistTCPlus[i]=0;}
+    if(fMassHistTCMinus[i]){ delete fMassHistTCMinus[i]; fMassHistTCMinus[i]=0;}
 
     if(fMassHistLS[i]){ delete fMassHistLS[i]; fMassHistLS[i]=0;}
     if(fCosPHistLS[i]){ delete fCosPHistLS[i]; fCosPHistLS[i]=0;}
@@ -469,7 +471,12 @@ void AliAnalysisTaskSEDplus::UserCreateOutputObjects()
     hisname.Form("hMassPt%dTC",i);
     fMassHistTC[index]=new TH1F(hisname.Data(),hisname.Data(),nbins,fLowmasslimit,fUpmasslimit);
     fMassHistTC[index]->Sumw2();
-
+    hisname.Form("hMassPt%dTCPlus",i);
+    fMassHistTCPlus[index]=new TH1F(hisname.Data(),hisname.Data(),nbins,fLowmasslimit,fUpmasslimit);
+    fMassHistTCPlus[index]->Sumw2();
+    hisname.Form("hMassPt%dTCMinus",i);
+    fMassHistTCMinus[index]=new TH1F(hisname.Data(),hisname.Data(),nbins,fLowmasslimit,fUpmasslimit);
+    fMassHistTCMinus[index]->Sumw2();
 
 
 
@@ -533,6 +540,12 @@ void AliAnalysisTaskSEDplus::UserCreateOutputObjects()
     hisname.Form("hSigPt%dTC",i);
     fMassHistTC[index]=new TH1F(hisname.Data(),hisname.Data(),nbins,fLowmasslimit,fUpmasslimit);
     fMassHistTC[index]->Sumw2();
+    hisname.Form("hSigPt%dTCPlus",i);
+    fMassHistTCPlus[index]=new TH1F(hisname.Data(),hisname.Data(),nbins,fLowmasslimit,fUpmasslimit);
+    fMassHistTCPlus[index]->Sumw2();
+    hisname.Form("hSigPt%dTCMinus",i);
+    fMassHistTCMinus[index]=new TH1F(hisname.Data(),hisname.Data(),nbins,fLowmasslimit,fUpmasslimit);
+    fMassHistTCMinus[index]->Sumw2();
 
     hisname.Form("hLSPt%dLCnw",i);
     fMassHistLS[indexLS]=new TH1F(hisname.Data(),hisname.Data(),nbins,fLowmasslimit,fUpmasslimit);
@@ -594,6 +607,12 @@ void AliAnalysisTaskSEDplus::UserCreateOutputObjects()
     hisname.Form("hBkgPt%dTC",i);
     fMassHistTC[index]=new TH1F(hisname.Data(),hisname.Data(),nbins,fLowmasslimit,fUpmasslimit);
     fMassHistTC[index]->Sumw2();
+    hisname.Form("hBkgPt%dTCPlus",i);
+    fMassHistTCPlus[index]=new TH1F(hisname.Data(),hisname.Data(),nbins,fLowmasslimit,fUpmasslimit);
+    fMassHistTCPlus[index]->Sumw2();
+    hisname.Form("hBkgPt%dTCMinus",i);
+    fMassHistTCMinus[index]=new TH1F(hisname.Data(),hisname.Data(),nbins,fLowmasslimit,fUpmasslimit);
+    fMassHistTCMinus[index]->Sumw2();
 
     hisname.Form("hLSPt%dLCntrip",i);
     fMassHistLS[indexLS]=new TH1F(hisname.Data(),hisname.Data(),nbins,fLowmasslimit,fUpmasslimit);
@@ -649,6 +668,8 @@ void AliAnalysisTaskSEDplus::UserCreateOutputObjects()
     fOutput->Add(fPtMaxHist[i]);
     fOutput->Add(fDCAHist[i]);
     fOutput->Add(fMassHistTC[i]);
+    fOutput->Add(fMassHistTCPlus[i]);
+    fOutput->Add(fMassHistTCMinus[i]);
   }
   for(Int_t i=0; i<3*fNPtBins&&fDoLS; i++){
     fOutput->Add(fCosPHistLS[i]);
@@ -880,6 +901,8 @@ void AliAnalysisTaskSEDplus::UserExec(Option_t */*option*/)
          
          if(passTightCuts){
            fMassHistTC[index]->Fill(invMass);
+           if(d->GetCharge()>0) fMassHistTCPlus[index]->Fill(invMass);
+           else if(d->GetCharge()<0) fMassHistTCMinus[index]->Fill(invMass);
          }
        }
        
@@ -896,6 +919,8 @@ void AliAnalysisTaskSEDplus::UserExec(Option_t */*option*/)
              fDCAHist[index]->Fill(dca);
              if(passTightCuts){
                fMassHistTC[index]->Fill(invMass);            
+               if(d->GetCharge()>0) fMassHistTCPlus[index]->Fill(invMass);
+               else if(d->GetCharge()<0) fMassHistTCMinus[index]->Fill(invMass);
              }
            }       
            fYVsPtSig->Fill(ptCand,rapid);
@@ -912,6 +937,8 @@ void AliAnalysisTaskSEDplus::UserExec(Option_t */*option*/)
              fDCAHist[index]->Fill(dca);
              if(passTightCuts){
                fMassHistTC[index]->Fill(invMass);
+               if(d->GetCharge()>0) fMassHistTCPlus[index]->Fill(invMass);
+               else if(d->GetCharge()<0) fMassHistTCMinus[index]->Fill(invMass);
              }
            }   
          }
@@ -974,6 +1001,10 @@ void AliAnalysisTaskSEDplus::Terminate(Option_t */*option*/)
     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()));
@@ -1014,6 +1045,10 @@ void AliAnalysisTaskSEDplus::Terminate(Option_t */*option*/)
     
     hisname.Form("hSigPt%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()));
@@ -1052,6 +1087,10 @@ void AliAnalysisTaskSEDplus::Terminate(Option_t */*option*/)
     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()));
index 4c0f3614fb405d65b0f8e3bac615b7587da269b2..ddddc55aec78cf236687169a59910a2e0c252fb7 100644 (file)
@@ -77,6 +77,8 @@ class AliAnalysisTaskSEDplus : public AliAnalysisTaskSE
   TH1F*   fPtMaxHist[3*kMaxPtBins]; //!hist. for Pt Max (LC)
   TH1F*   fDCAHist[3*kMaxPtBins]; //!hist. for DCA (LC)
   TH1F *fMassHistTC[3*kMaxPtBins]; //!hist. for inv mass (TC)
+  TH1F *fMassHistTCPlus[3*kMaxPtBins]; //!hist. for D+ inv mass (TC)
+  TH1F *fMassHistTCMinus[3*kMaxPtBins]; //!hist. for D- inv mass (TC)
   TH1F *fMassHistLS[5*kMaxPtBins];//!hist. for LS inv mass (LC)
   TH1F *fCosPHistLS[3*kMaxPtBins];//!hist. for LS cuts variable 1 (LC)
   TH1F *fDLenHistLS[3*kMaxPtBins];//!hist. for LS cuts variable 2 (LC)
@@ -104,7 +106,7 @@ class AliAnalysisTaskSEDplus : public AliAnalysisTaskSE
   Bool_t fReadMC;    //flag for access to MC
   Bool_t fDoLS;      //flag to do LS analysis
   
-  ClassDef(AliAnalysisTaskSEDplus,7); // AliAnalysisTaskSE for the MC association of heavy-flavour decay candidates
+  ClassDef(AliAnalysisTaskSEDplus,8); // AliAnalysisTaskSE for the MC association of heavy-flavour decay candidates
 };
 
 #endif