]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
adding event-by-event eta-phi histogram for general use
authorhdalsgaa <hdalsgaa@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 9 Jun 2010 11:21:31 +0000 (11:21 +0000)
committerhdalsgaa <hdalsgaa@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 9 Jun 2010 11:21:31 +0000 (11:21 +0000)
PWG2/FORWARD/analysis/AliFMDAnaParameters.cxx
PWG2/FORWARD/analysis/AliFMDAnaParameters.h
PWG2/FORWARD/analysis/AliFMDAnalysisTaskBackgroundCorrection.cxx
PWG2/FORWARD/analysis/AliFMDAnalysisTaskBackgroundCorrection.h

index bebb755a540a1a315d85e61a10f3a517d7ea5795..71dcc64fe511d396204c827d767c2aab81da943a 100644 (file)
@@ -88,7 +88,8 @@ AliFMDAnaParameters::AliFMDAnaParameters() :
   fSPDlowLimit(0),
   fSPDhighLimit(999999999),
   fCentralSelection(kFALSE),
-  fSharingObjectPresent(kTRUE)
+  fSharingObjectPresent(kTRUE),
+  fNumberOfEtaBinsToCut(1)
 {
   // Default constructor 
   fPhysicsSelection = new AliPhysicsSelection;
@@ -160,9 +161,12 @@ void AliFMDAnaParameters::Init(Bool_t forceReInit, UInt_t what)
   if (what & kEventSelectionEfficiency)   InitEventSelectionEff();
   if (what & kSharingEfficiency)          InitSharingEff();
   
-
-  
   fIsInit = kTRUE;
+  
+  if(fBackground)
+    FindEtaLimits();
+  
+  
 }
 //____________________________________________________________________
 
@@ -176,6 +180,8 @@ void AliFMDAnaParameters::InitBackground() {
   fBackground = dynamic_cast<AliFMDAnaCalibBackgroundCorrection*>(fin->Get(fgkBackgroundID));
   if (!fBackground) AliFatal("Invalid background object from CDB");
   
+  
+
 }
 
 //____________________________________________________________________
@@ -224,6 +230,24 @@ void AliFMDAnaParameters::InitSharingEff() {
     return; 
   }
   
+}
+//____________________________________________________________________
+void AliFMDAnaParameters::FindEtaLimits() {
+
+  fEtaLowBinLimits.SetBins(3,0,3,2,0,2,GetNvtxBins(),0,GetNvtxBins());
+  fEtaHighBinLimits.SetBins(3,0,3,2,0,2,GetNvtxBins(),0,GetNvtxBins());
+  for(Int_t det=1; det<=3;det++) {
+    Int_t nRings = (det==1 ? 1 : 2);
+    for (UShort_t ir = 0; ir < nRings; ir++) {
+      Char_t ringChar = (ir == 0 ? 'I' : 'O');
+      for(Int_t v =0; v<GetNvtxBins(); v++) {
+       fEtaLowBinLimits.SetBinContent(det,ir,v,GetFirstEtaBinFromMap(v, det, ringChar));
+       fEtaHighBinLimits.SetBinContent(det,ir,v,GetLastEtaBinFromMap(v, det, ringChar));
+       
+      }
+    }
+  }
+
 }
 //____________________________________________________________________
 
@@ -805,6 +829,68 @@ AliFMDAnaParameters::GetBaseStripLength(Char_t ring, UShort_t strip)
   return basearclength;
 }
 //____________________________________________________________________
+Int_t    AliFMDAnaParameters::GetFirstEtaBinFromMap(Int_t vtxbin, Int_t det, Char_t ring)
+{
+  TH2F* hBg = GetBackgroundCorrection(det,ring,vtxbin);
+  Int_t firstbin = -1;
+  Int_t nNonZeroFirst = 0;
+  
+  for(Int_t i=1;i<=hBg->GetNbinsX();i++) {
+    if(nNonZeroFirst == fNumberOfEtaBinsToCut && firstbin==-1) firstbin = i;
+    
+    for(Int_t j=1;j<=hBg->GetNbinsY();j++) {
+      
+      Float_t value = hBg->GetBinContent(i,j);
+      
+      if(value > 0.001 && nNonZeroFirst<fNumberOfEtaBinsToCut)
+       {nNonZeroFirst++; break;}
+      
+      
+    }
+  }
+  
+  return firstbin;
+
+}
+//____________________________________________________________________
+Int_t    AliFMDAnaParameters::GetLastEtaBinFromMap(Int_t vtxbin, Int_t det, Char_t ring)
+{
+  TH2F* hBg = GetBackgroundCorrection(det,ring,vtxbin);
+  Int_t lastbin=-1;
+  Int_t nNonZeroLast = 0;
+  for(Int_t i=hBg->GetNbinsX();i>0;i--) {
+    if(nNonZeroLast == fNumberOfEtaBinsToCut && lastbin==-1) lastbin = i;
+    
+    for(Int_t j=1;j<=hBg->GetNbinsY();j++) {
+      
+      Float_t value = hBg->GetBinContent(i,j);
+      
+      if(value > 0.001 && nNonZeroLast<fNumberOfEtaBinsToCut)
+       {nNonZeroLast++; break; }
+      
+      
+    }
+  }
+  
+  return lastbin;
+}
+
+//____________________________________________________________________
+Int_t    AliFMDAnaParameters::GetFirstEtaBinToInclude(Int_t vtxbin, Int_t det, Char_t ring)
+{
+  Int_t ringNumber = (ring == 'I' ? 0 : 1);
+  return fEtaLowBinLimits.GetBinContent(det,ringNumber,vtxbin);
+
+}
+
+//____________________________________________________________________
+Int_t    AliFMDAnaParameters::GetLastEtaBinToInclude(Int_t vtxbin, Int_t det, Char_t ring)
+{
+  Int_t ringNumber = (ring == 'I' ? 0 : 1);
+  return fEtaHighBinLimits.GetBinContent(det,ringNumber,vtxbin);
+  
+}
+//____________________________________________________________________
 //
 // EOF
 //
index 7747018453609b25dd75d507014004755d5c07c2..cf20363b7e79273f4a452e9c5f1a7e8777df7551 100644 (file)
@@ -35,6 +35,7 @@
 #include "TH2F.h"
 #include "TAxis.h"
 #include "TH1F.h"
+#include "TH3F.h"
 #include "AliFMDAnaCalibBackgroundCorrection.h"
 #include "AliFMDAnaCalibEnergyDistribution.h"
 #include "AliFMDAnaCalibEventSelectionEfficiency.h"
@@ -141,6 +142,13 @@ public:
   void     SetHighSPDLimit(Float_t cut) {fSPDhighLimit = cut;}
   void     SetCentralTriggerSelection(Bool_t selection) {fCentralSelection = selection;}
   Bool_t   SharingEffPresent() {return fSharingObjectPresent;}
+  Int_t    GetFirstEtaBinToInclude(Int_t vtxbin, Int_t det, Char_t ring) ;
+  Int_t    GetLastEtaBinToInclude(Int_t vtxbin, Int_t det, Char_t ring) ;
+  
+
+  void     SetNumberOfEtaBinsToCut(Int_t nbins) {fNumberOfEtaBinsToCut = nbins;}
+  Int_t    GetNumberOfEtaBinsToCut() const {return fNumberOfEtaBinsToCut;}
+  
 protected:
   
   AliFMDAnaParameters(const AliFMDAnaParameters& o) 
@@ -167,7 +175,8 @@ protected:
       fSPDlowLimit(o.fSPDlowLimit),
       fSPDhighLimit(o.fSPDhighLimit),   
       fCentralSelection(o.fCentralSelection),
-      fSharingObjectPresent(o.fSharingObjectPresent)
+      fSharingObjectPresent(o.fSharingObjectPresent),
+      fNumberOfEtaBinsToCut(o.fNumberOfEtaBinsToCut)
   {}
   AliFMDAnaParameters& operator=(const AliFMDAnaParameters&) { return *this; }
   virtual ~AliFMDAnaParameters() {}
@@ -180,7 +189,11 @@ protected:
   void InitEventSelectionEff();
   void InitSharingEff();
   
+  void     FindEtaLimits();
+  Int_t    GetFirstEtaBinFromMap(Int_t vtxbin, Int_t det, Char_t ring) ;
+  Int_t    GetLastEtaBinFromMap(Int_t vtxbin, Int_t det, Char_t ring) ;
 
+  
   TObjArray* GetBackgroundArray();
   
   TAxis* GetRefAxis();
@@ -218,6 +231,9 @@ protected:
   Float_t  fSPDhighLimit ;             // high limit of SPD tracklets
   Bool_t   fCentralSelection;         //if event selection is done centrally
   Bool_t   fSharingObjectPresent ;    //Do we have a sharing object ? 
+  Int_t    fNumberOfEtaBinsToCut;     //Number of eta bins to remove from edge effects
+  TH3F     fEtaLowBinLimits;
+  TH3F     fEtaHighBinLimits;
   ClassDef(AliFMDAnaParameters,1) // Manager of parameters
 };
 
index 312a318a3be897e98c3c685a84049858d089fff7..68672f697955282d74857c7aa7c4c1aea9526a04 100644 (file)
@@ -21,6 +21,7 @@
 #include "AliFMDAnaParameters.h"
 #include "AliESDInputHandler.h"
 #include "AliMultiplicity.h"
+#include "TProfile2D.h"
 //#include "AliFMDGeometry.h"
 
 ClassImp(AliFMDAnalysisTaskBackgroundCorrection)
@@ -148,6 +149,11 @@ void AliFMDAnalysisTaskBackgroundCorrection::CreateOutputObjects()
     
   }
   
+  TH2F* dNdetadphiHistogram = new TH2F("dNdetadphiHistogramTrVtx","dNdetadphiHistogramTrVtx;#eta;#Phi",pars->GetNetaBins(),-6,6,20,0,2*TMath::Pi());
+  
+  //dNdetadphiHistogram->SetErrorOption("g");
+  
+  fHitList->Add(dNdetadphiHistogram);
   
   
 }
@@ -298,6 +304,9 @@ void AliFMDAnalysisTaskBackgroundCorrection::Exec(Option_t */*option*/)
     AliWarning("No SPD background map found");
   
   //std::cout<<spdmult->GetNumberOfTracklets()<<"  "<<spdmult->GetNumberOfITSClusters(0)<<"    "<< spdmult->GetNumberOfSingleClusters()<<std::endl;
+  
+  CreatePerEventHistogram(vtxbin);
+  
   if(fStandalone) {
     PostData(0, fOutputList); 
   }
@@ -328,6 +337,69 @@ void AliFMDAnalysisTaskBackgroundCorrection::Terminate(Option_t */*option*/) {
     }
   }
 }
+
+//_____________________________________________________________________
+void AliFMDAnalysisTaskBackgroundCorrection::CreatePerEventHistogram(Int_t vtxbin) {
+  
+  AliFMDAnaParameters* pars = AliFMDAnaParameters::Instance();
+  TH2F* dNdetadphiHistogram = (TH2F*)fHitList->FindObject("dNdetadphiHistogramTrVtx");
+  dNdetadphiHistogram->Reset();
+  
+  for(Int_t det = 1; det<=3; det++) {
+    Int_t maxRing = (det == 1 ? 0 : 1);
+    for(Int_t ring = 0;ring<=maxRing;ring++) {
+      
+      Char_t ringChar = (ring == 0 ? 'I' : 'O');
+      
+      TH2F* multhistoriginal = (TH2F*)fOutputList->FindObject(Form("multTrVtx_FMD%d%c_vtxbin%d",det,ringChar,vtxbin));
+      TH2F* multhist = (TH2F*)multhistoriginal->Clone("tmp");
+      if(ringChar == 'O')
+       multhist->RebinY(2);
+      
+      for(Int_t i=pars->GetFirstEtaBinToInclude(vtxbin,det,ringChar); i<=pars->GetLastEtaBinToInclude(vtxbin,det,ringChar); i++) {
+       for(Int_t j=1; j<=multhist->GetNbinsY(); j++) {
+         if(multhist->GetBinContent(i,j) < 0.0001) continue;
+         
+         Bool_t overlap = kFALSE;
+       
+         if(det == 1 && ringChar =='I')
+           if(i<=pars->GetLastEtaBinToInclude(vtxbin,2,'I'))
+             overlap = kTRUE;
+                 
+         if(det == 2 && ringChar =='O')
+           if(i>=pars->GetFirstEtaBinToInclude(vtxbin,2,'I'))
+             overlap = kTRUE;
+         
+         if(det == 2 && ringChar =='I')
+           if(i<=pars->GetLastEtaBinToInclude(vtxbin,2,'O') || i>=pars->GetFirstEtaBinToInclude(vtxbin,1,'I'))
+             overlap = kTRUE;
+                 
+         if(det == 3 && ringChar =='I')
+           if(i>=pars->GetFirstEtaBinToInclude(vtxbin,3,'O'))
+             overlap = kTRUE;
+         
+         if(det == 3 && ringChar =='O')
+           if(i<=pars->GetLastEtaBinToInclude(vtxbin,3,'I'))
+             overlap = kTRUE;
+         
+         
+         
+
+         
+         if(overlap) {
+           dNdetadphiHistogram->SetBinContent(i,j,dNdetadphiHistogram->GetBinContent(i,j)+0.5*multhist->GetBinContent(i,j));
+             dNdetadphiHistogram->SetBinError(i,j,0.5*TMath::Sqrt(TMath::Power(dNdetadphiHistogram->GetBinError(i,j),2)+TMath::Power(multhist->GetBinError(i,j),2)));
+         }
+           else {
+             dNdetadphiHistogram->SetBinContent(i,j,multhist->GetBinContent(i,j));
+             dNdetadphiHistogram->SetBinError(i,j,multhist->GetBinError(i,j));
+           }
+       }
+      }
+      delete multhist;
+    }
+  }
+}
 //_____________________________________________________________________
 //
 //
index 20a01d7d16e7a1f16c8fabe885db3af514daf4c7..e161dfa6846c9e3fe55489e95ddba86d62ca9e24 100644 (file)
@@ -41,11 +41,12 @@ class AliFMDAnalysisTaskBackgroundCorrection : public AliAnalysisTask
     virtual void Exec(Option_t *option);
     virtual void Terminate(Option_t *option);
     virtual void SetDebugLevel(Int_t level) {fDebug = level;}
-    void SetInputList(TList* inputList) {fInputList = inputList;}
-    void SetOutputVertex(TObjString* vtxString) {fOutputVertexString = vtxString;}
-    //void SetInputVtx(TObjString* vtxString) {fVertexString = vtxString;}
-    void SetOutputList(TList* outputList) {fOutputList = outputList;}
-    void SetHitList(TList* hitList) {fHitList = hitList;}
+  void SetInputList(TList* inputList) {fInputList = inputList;}
+  void SetOutputVertex(TObjString* vtxString) {fOutputVertexString = vtxString;}
+  //void SetInputVtx(TObjString* vtxString) {fVertexString = vtxString;}
+  void SetOutputList(TList* outputList) {fOutputList = outputList;}
+  void SetHitList(TList* hitList) {fHitList = hitList;}
+  void         CreatePerEventHistogram(Int_t vtxbin);
  private:
     Int_t         fDebug;        //  Debug flag
     TList*        fOutputList;