]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Fixing algorithm for efficiency calculation (Diego)
authormartinez <martinez@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 16 Nov 2009 08:08:08 +0000 (08:08 +0000)
committermartinez <martinez@f7af4fe6-9843-0410-8265-dc069ae4e863>
Mon, 16 Nov 2009 08:08:08 +0000 (08:08 +0000)
PWG3/muon/AliAnalysisTaskTrigChEff.cxx
PWG3/muon/AliAnalysisTaskTrigChEff.h

index d529493093a3a9fe3acf5a2d43a5687bfbe7e230..3453b05e2858547269f21cb8dc400dcee3e50303 100644 (file)
@@ -1,7 +1,36 @@
+/**************************************************************************
+ * Copyright(c) 1998-2007, ALICE Experiment at CERN, All rights reserved. *
+ *                                                                        *
+ * Author: The ALICE Off-line Project.                                    *
+ * 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.                  *
+ **************************************************************************/
+
+//-----------------------------------------------------------------------------
+/// \class AliAnalysisTaskSingleMu
+/// Analysis task for single muons in the spectrometer.
+/// The output is a list of histograms.
+/// The macro class can run on AOD or ESDs.
+/// If Monte Carlo information is present, some basics checks are performed.
+///
+/// \author Diego Stocco
+//-----------------------------------------------------------------------------
+
+//----------------------------------------------------------------------------
+//    Implementation of the class for trigger chamber efficiency determinaltion
+//----------------------------------------------------------------------------
+
+
 #define AliAnalysisTaskTrigChEff_cxx
 
 // ROOT includes
-#include "TChain.h"
 #include "TH1.h"
 #include "TCanvas.h"
 #include "TROOT.h"
 
 // STEER includes
 #include "AliLog.h"
-
 #include "AliESDEvent.h"
 #include "AliESDMuonTrack.h"
-#include "AliESDInputHandler.h"
 
 // ANALYSIS includes
-#include "AliAnalysisTask.h"
-#include "AliAnalysisDataSlot.h"
-#include "AliAnalysisManager.h"
+#include "AliAnalysisTaskSE.h"
+
 #include "AliAnalysisTaskTrigChEff.h"
 
 ClassImp(AliAnalysisTaskTrigChEff)
 
 //________________________________________________________________________
 AliAnalysisTaskTrigChEff::AliAnalysisTaskTrigChEff(const char *name) :
-  AliAnalysisTask(name,""), 
-  fESD(0),
+  AliAnalysisTaskSE(name), 
   fUseGhosts(kFALSE),
   fList(0)
 {
   //
   /// Constructor.
   //
-  // Input slot #0 works with an Ntuple
-  DefineInput(0, TChain::Class());
-  // Output slot #0 writes into a TObjArray container
-  DefineOutput(0,  TList::Class());
+  // Output slot #1 writes into a TObjArray container
+  DefineOutput(1,  TList::Class());
 }
 
 //___________________________________________________________________________
-void AliAnalysisTaskTrigChEff::ConnectInputData(Option_t *) {
-  //
-  /// Connect ESD here
-  /// Called once
-  //
-
-  TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
-  if (!tree) {
-    Printf("ERROR: Could not read chain from input slot 0");
-  } else {
-    // Disable all branches and enable only the needed ones
-    // The next two lines are different when data produced as AliESDEvent is read
-    tree->SetBranchStatus("*", kFALSE);
-    tree->SetBranchStatus("MuonTracks.*", kTRUE);
-
-    AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
-      
-    if (!esdH) {
-      Printf("ERROR: Could not get ESDInputHandler");
-    } else
-      fESD = esdH->GetEvent();
-  }
-}
-
-//___________________________________________________________________________
-void AliAnalysisTaskTrigChEff::CreateOutputObjects() {
+void AliAnalysisTaskTrigChEff::UserCreateOutputObjects() {
   //
   /// Create histograms
   /// Called once
   //
-  Printf("   CreateOutputObjects of task %s\n", GetName());
 
   TString cathCode[2] = {"bendPlane", "nonBendPlane"};
   TString countTypeName[2] = {"CountInCh", "NonCountInCh"};
@@ -92,8 +89,8 @@ void AliAnalysisTaskTrigChEff::CreateOutputObjects() {
   Float_t boardLow = 1-0.5, boardHigh = kNboards+1.-0.5;
   const Char_t* boardName = "board";
 
-  Int_t angleBins = 140;
-  Float_t angleLow = -35., angleHigh = 35.;
+  Int_t angleBins = 280;
+  Float_t angleLow = -70., angleHigh = 70.;
   const Char_t* angleNameX = "#theta_{x} (deg)";
   const Char_t* angleNameY = "#theta_{y} (deg)";
 
@@ -117,7 +114,7 @@ void AliAnalysisTaskTrigChEff::CreateOutputObjects() {
   fList->AddAt(histo, kHtracksInBoard);
 
   for(Int_t hType=0; hType<kNcounts; hType++){
-    Int_t hindex = (hType==0) ? kHchamberAllEff : kHchamberNonEff;
+    Int_t hindex = (hType==0) ? kHchamberEff : kHchamberNonEff;
     for(Int_t cath=0; cath<kNcathodes; cath++){
       histoName = Form("%sChamber%s", cathCode[cath].Data(), countTypeName[hType].Data());
       histo = new TH1F(histoName, histoName,
@@ -130,7 +127,7 @@ void AliAnalysisTaskTrigChEff::CreateOutputObjects() {
   } // loop on counts
 
   for(Int_t hType=0; hType<kNcounts; hType++){
-    Int_t hindex = (hType==0) ? kHslatAllEff : kHslatNonEff;
+    Int_t hindex = (hType==0) ? kHslatEff : kHslatNonEff;
     for(Int_t cath=0; cath<kNcathodes; cath++){
       for(Int_t ch=0; ch<kNchambers; ch++){
        Int_t chCath = GetPlane(cath, ch);
@@ -146,7 +143,7 @@ void AliAnalysisTaskTrigChEff::CreateOutputObjects() {
   } // loop on counts
 
   for(Int_t hType=0; hType<kNcounts; hType++){
-    Int_t hindex = (hType==0) ? kHboardAllEff : kHboardNonEff;
+    Int_t hindex = (hType==0) ? kHboardEff : kHboardNonEff;
     for(Int_t cath=0; cath<kNcathodes; cath++){
       for(Int_t ch=0; ch<kNchambers; ch++){
        Int_t chCath = GetPlane(cath, ch);
@@ -176,33 +173,35 @@ void AliAnalysisTaskTrigChEff::CreateOutputObjects() {
 }
 
 //________________________________________________________________________
-void AliAnalysisTaskTrigChEff::Exec(Option_t *) {
+void AliAnalysisTaskTrigChEff::UserExec(Option_t *) {
   //
   /// Main loop
   /// Called for each event
   //
-  Int_t nTracks = 0, board = 0;
+  AliESDEvent* esdEvent = dynamic_cast<AliESDEvent*> (InputEvent());
+
+  if (!esdEvent) {
+    Printf("ERROR: esdEvent not available\n");
+    return;
+  }
+
+  Int_t slat = 0, board = 0;
   UShort_t pattern = 0;
   AliESDMuonTrack *esdTrack = 0x0;
 
   const Float_t kRadToDeg = 180./TMath::Pi();
+  Int_t nTracks = esdEvent->GetNumberOfMuonTracks();
 
-  if (!fESD) {
-    Printf("ERROR: fESD not available");
-    return;
-  }
-  nTracks = fESD->GetNumberOfMuonTracks();
-
-  // Object declaration
   const Int_t kFirstTrigCh = 11; //AliMpConstants::NofTrackingChambers()+1;
 
+  TArrayI othersEfficient(kNchambers);
+
   for (Int_t itrack = 0; itrack < nTracks; itrack++) {
-    esdTrack = fESD->GetMuonTrack(itrack);
-    pattern =  esdTrack->GetHitsPatternInTrigCh();
-    board = esdTrack->LoCircuit();
+    esdTrack = esdEvent->GetMuonTrack(itrack);
 
     if ( ! esdTrack->ContainTrackerData() && ! fUseGhosts ) continue;
 
+    pattern =  esdTrack->GetHitsPatternInTrigCh();
     Int_t effFlag = AliESDMuonTrack::GetEffFlag(pattern);
 
     if(effFlag < AliESDMuonTrack::kChEff) continue; // Track not good for efficiency calculation
@@ -210,48 +209,63 @@ void AliAnalysisTaskTrigChEff::Exec(Option_t *) {
     ((TH1F*)fList->At(kHthetaX))->Fill(esdTrack->GetThetaX() * kRadToDeg);
     ((TH1F*)fList->At(kHthetaY))->Fill(esdTrack->GetThetaY() * kRadToDeg);
 
-    Int_t slat = AliESDMuonTrack::GetSlatOrInfo(pattern);
+    othersEfficient.Reset(1);
+    for(Int_t cath=0; cath<kNcathodes; cath++){
+      for(Int_t ich=0; ich<kNchambers; ich++){
+       if( ! AliESDMuonTrack::IsChamberHit(pattern, cath, ich)){
+         for(Int_t jch=0; jch<kNchambers; jch++){
+           if ( jch != ich) {
+             othersEfficient[jch] = 0;
+             //AliInfo(Form("%s ch %i by New", baseOutString.Data(), jch));
+           }
+         } // loop on other chambers
+         break;
+       } // if chamber not efficient
+      } // loop on chambers
+    } // loop on cathodes
+
+    Bool_t rejectTrack = kTRUE;
+    for (Int_t ich=0; ich<kNchambers; ich++){
+      if ( othersEfficient[ich] > 0 ){
+       rejectTrack = kFALSE;
+       break;
+      }
+    }
+
+    if ( rejectTrack ) continue;
+
+    slat = AliESDMuonTrack::GetSlatOrInfo(pattern);
+    board = esdTrack->LoCircuit();
 
     if(effFlag >= AliESDMuonTrack::kSlatEff) ((TH1F*)fList->At(kHtracksInSlat))->Fill(slat);
     if(effFlag >= AliESDMuonTrack::kBoardEff) ((TH1F*)fList->At(kHtracksInBoard))->Fill(board);
 
     for(Int_t cath=0; cath<kNcathodes; cath++){
-      Int_t ineffCh = -1;
-      for(Int_t ich=0; ich<kNchambers; ich++){
-       if(!AliESDMuonTrack::IsChamberHit(pattern, cath, ich)){
-         ineffCh = ich;
-         break;
-       }
-      }
+      for(Int_t ch=0; ch<kNchambers; ch++){
+       if ( ! othersEfficient[ch] )
+         continue; // Reject track if the info of the chamber under study 
+                   // is necessary to create the track itself
+
+       Int_t whichType = AliESDMuonTrack::IsChamberHit(pattern, cath, ch) ? kChHit : kChNonHit;
 
-      Int_t nChambers = kNchambers;
-      for(Int_t ch=0; ch<nChambers; ch++){
-       Int_t whichType = kAllChEff;
-       Int_t currCh = ch;
-       if(ineffCh>=0){
-         whichType = kChNonEff;
-         currCh = ineffCh;
-         nChambers = -1;
-       }
-
-       Int_t iChamber = kFirstTrigCh + currCh;
-       Int_t hindex = (whichType==kAllChEff) ? kHchamberAllEff : kHchamberNonEff;
+       Int_t iChamber = kFirstTrigCh + ch;
+       Int_t hindex = ( whichType == kChHit ) ? kHchamberEff : kHchamberNonEff;
        ((TH1F*)fList->At(hindex + cath))->Fill(iChamber);
 
        if(effFlag < AliESDMuonTrack::kSlatEff) continue; // Track crossed different slats
-       Int_t chCath = GetPlane(cath, currCh);
-       hindex = (whichType==kAllChEff) ? kHslatAllEff : kHslatNonEff;
+       Int_t chCath = GetPlane(cath, ch);
+       hindex = ( whichType == kChHit ) ? kHslatEff : kHslatNonEff;
        ((TH1F*)fList->At(hindex + chCath))->Fill(slat);
 
        if(effFlag < AliESDMuonTrack::kBoardEff) continue; // Track crossed different boards
-       hindex = (whichType==kAllChEff) ? kHboardAllEff : kHboardNonEff;
+       hindex = ( whichType == kChHit ) ? kHboardEff : kHboardNonEff;
        ((TH1F*)fList->At(hindex + chCath))->Fill(board);
       } // loop on chambers
     } // loop on cathodes
-  }
+  } // loop on tracks
 
   // Post final data. It will be written to a file with option "RECREATE"
-  PostData(0, fList);
+  PostData(1, fList);
 }
 
 //________________________________________________________________________
@@ -272,7 +286,7 @@ void AliAnalysisTaskTrigChEff::Terminate(Option_t *) {
       can[cath]->Divide(2,2);
       for(Int_t ch=0; ch<kNchambers; ch++){
        Int_t chCath = GetPlane(cath, ch);
-       num = (TH1F*)(fList->At(kHboardAllEff + chCath)->Clone());
+       num = (TH1F*)(fList->At(kHboardEff + chCath)->Clone());
        den = (TH1F*)(fList->At(kHboardNonEff + chCath)->Clone());
        den->Add(num);
        num->Divide(den);
index 88a94fdce4bfa40e35799dd66c938249328e6da3..b75426969211004c24859209b9ab40ce788ecee3 100644 (file)
@@ -1,19 +1,21 @@
-#include "TH1.h"
-#include "TList.h"
+/// \ingroup "PWG3muon"
+/// \class AliAnalysisTaskTrigChEff
+/// \brief Analysis task for trigger chamber efficiency determination
+///
+//  Author Diego Stocco
 
-#include "AliESDEvent.h"
-#include "AliAODEvent.h"
+class TList;
 
-class AliAnalysisTaskTrigChEff : public AliAnalysisTask {
+class AliAnalysisTaskTrigChEff : public AliAnalysisTaskSE {
  public:
   AliAnalysisTaskTrigChEff(const char *name = "AliAnalysisTaskTrigChEff");
   virtual ~AliAnalysisTaskTrigChEff() {}
-  
-  virtual void   ConnectInputData(Option_t *);
-  virtual void   CreateOutputObjects();
-  virtual void   Exec(Option_t *option);
+
+  virtual void   UserCreateOutputObjects();
+  virtual void   UserExec(Option_t *option);
   virtual void   Terminate(Option_t *);
 
+  /// Use ghost tracks in calculations
   void SetUseGhostTracks(Bool_t useGhosts = kTRUE) { fUseGhosts = useGhosts; }
 
 protected:
@@ -25,10 +27,9 @@ private:
   /// Not implemented
   AliAnalysisTaskTrigChEff& operator = (const AliAnalysisTaskTrigChEff& rhs);
     
-  AliESDEvent* fESD; //!< ESDevent object
-  Bool_t fUseGhosts; //!< Flag to use also the trigger tracks not matching the tracker in eff. calculation
+  Bool_t fUseGhosts; ///< Flag to use also the trigger tracks not matching the tracker in eff. calculation
 
-  TList*  fList; //TList output object
+  TList*  fList; ///<TList output object
 
   enum {
     kNcathodes = 2,  ///< Number of cathodes
@@ -37,16 +38,16 @@ private:
     kNslats    = 18 ///< Number of slats
   };
 
-  enum {kAllChEff, kChNonEff, kNcounts};
+  enum {kChHit, kChNonHit, kNcounts};
 
   enum {
     kHtracksInSlat  = 0,  ///< Tracks in slat histogram index
     kHtracksInBoard = 1,  ///< Tracks in board histogram index
-    kHchamberAllEff = 2,  ///< N44 per cathode histogram index
+    kHchamberEff    = 2,  ///< N44 per cathode histogram index
     kHchamberNonEff = 4,  ///< N33 per cathode histogram index
-    kHslatAllEff    = 6,  ///< N44 per slat histogram index
+    kHslatEff       = 6,  ///< N44 per slat histogram index
     kHslatNonEff    = 14, ///< N33 per slat histogram index
-    kHboardAllEff   = 22, ///< N44 per board histogram index
+    kHboardEff      = 22, ///< N44 per board histogram index
     kHboardNonEff   = 30, ///< N33 per board histogram index
     kHthetaX        = 38, ///< Angular distribution theta_x
     kHthetaY        = 39  ///< Angular distribution theta_y
@@ -55,6 +56,6 @@ private:
   /// Given cathode and chamber, return plane number
   Int_t GetPlane(Int_t cathode, Int_t chamber) { return kNchambers*cathode + chamber; }
 
-  ClassDef(AliAnalysisTaskTrigChEff, 0); // Single muon analysis
+  ClassDef(AliAnalysisTaskTrigChEff, 1); // Trigger chamber efficiency analysis
 };