]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
o use AliPIDResponse (Markus Fasel)
authorwiechula <wiechula@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 27 Sep 2012 10:27:12 +0000 (10:27 +0000)
committerwiechula <wiechula@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 27 Sep 2012 10:27:12 +0000 (10:27 +0000)
ANALYSIS/AliESDpidCuts.cxx
ANALYSIS/AliESDpidCuts.h

index 2970b5a815b278305ed857de02c2b4be3554007e..c46dc83a4f2605a3a5320151b5ee6f3087ceb20a 100644 (file)
 #include <TString.h>\r
 #include <TList.h>\r
 \r
+#include "AliAnalysisManager.h"\r
+#include "AliInputEventHandler.h"\r
 #include "AliESDtrack.h"\r
 #include "AliESDEvent.h"\r
 #include "AliLog.h"\r
-#include "AliESDpid.h"\r
+#include "AliPIDResponse.h"\r
 \r
 #include "AliESDpidCuts.h"\r
 \r
@@ -38,7 +40,7 @@ const Int_t AliESDpidCuts::kNcuts = 3;
 //_____________________________________________________________________\r
 AliESDpidCuts::AliESDpidCuts(const Char_t *name, const Char_t *title):\r
     AliAnalysisCuts(name, title)\r
-  , fESDpid(NULL)\r
+  , fPIDresponse(NULL)\r
   , fTPCsigmaCutRequired(0)\r
   , fTOFsigmaCutRequired(0)\r
   , fCutTPCclusterRatio(0.)\r
@@ -50,7 +52,6 @@ AliESDpidCuts::AliESDpidCuts(const Char_t *name, const Char_t *title):
   // Default constructor\r
   //\r
   \r
-  fESDpid = new AliESDpid;\r
   memset(fCutTPCnSigma, 0, sizeof(Float_t)* AliPID::kSPECIES * 2);\r
   memset(fCutTOFnSigma, 0, sizeof(Float_t)* AliPID::kSPECIES * 2);\r
 \r
@@ -62,7 +63,7 @@ AliESDpidCuts::AliESDpidCuts(const Char_t *name, const Char_t *title):
 //_____________________________________________________________________\r
 AliESDpidCuts::AliESDpidCuts(const AliESDpidCuts &ref):\r
     AliAnalysisCuts(ref)\r
-  , fESDpid(NULL)\r
+  , fPIDresponse(NULL)\r
   , fTPCsigmaCutRequired(ref.fTPCsigmaCutRequired)\r
   , fTOFsigmaCutRequired(ref.fTOFsigmaCutRequired)\r
   , fCutTPCclusterRatio(ref.fCutTPCclusterRatio)\r
@@ -73,7 +74,7 @@ AliESDpidCuts::AliESDpidCuts(const AliESDpidCuts &ref):
   //\r
   // Copy constructor\r
   //\r
-  fESDpid = new AliESDpid(*ref.fESDpid);\r
+  fPIDresponse = ref.fPIDresponse;\r
   memcpy(fCutTPCnSigma, ref.fCutTPCnSigma, sizeof(Float_t) * AliPID::kSPECIES * 2);\r
   memcpy(fCutTOFnSigma, ref.fCutTOFnSigma, sizeof(Float_t) * AliPID::kSPECIES * 2);\r
   \r
@@ -103,9 +104,7 @@ AliESDpidCuts::~AliESDpidCuts(){
   //\r
   // Destructor\r
   //\r
-  delete fESDpid;\r
 \r
-  delete fHcutStatistics;\r
   delete fHcutCorrelation;\r
   for(Int_t imode = 0; imode < 2; imode++){\r
     delete fHclusterRatio[imode];\r
@@ -116,6 +115,20 @@ AliESDpidCuts::~AliESDpidCuts(){
   }\r
 }\r
 \r
+//_____________________________________________________________________\r
+void AliESDpidCuts::Init(){\r
+  //\r
+  // Init function, get PID response from the Analysis Manager\r
+  //\r
+  AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();\r
+  if(mgr){\r
+    AliInputEventHandler *handler = dynamic_cast<AliInputEventHandler*>(mgr->GetInputEventHandler());\r
+    if(handler){\r
+      fPIDresponse = handler->GetPIDResponse();\r
+    }\r
+  }\r
+}\r
+\r
 //_____________________________________________________________________\r
 Bool_t AliESDpidCuts::IsSelected(TObject *obj){\r
   //\r
@@ -141,7 +154,7 @@ void AliESDpidCuts::Copy(TObject &c) const {
   //\r
   AliESDpidCuts &target = dynamic_cast<AliESDpidCuts &>(c);\r
 \r
-  target.fESDpid = new AliESDpid(*fESDpid);\r
+  target.fPIDresponse = fPIDresponse;\r
 \r
   target.fCutTPCclusterRatio = fCutTPCclusterRatio;\r
   target.fMinMomentumTOF = fMinMomentumTOF;\r
@@ -233,6 +246,10 @@ Bool_t AliESDpidCuts::AcceptTrack(const AliESDtrack *track, const AliESDEvent *e
   //\r
   // Check whether the tracks survived the cuts\r
   //\r
+  if(!fPIDresponse){\r
+    AliError("PID Response not available");\r
+    return 0;\r
+  }\r
   enum{\r
     kCutClusterRatioTPC,\r
     kCutNsigmaTPC,\r
@@ -252,7 +269,7 @@ Bool_t AliESDpidCuts::AcceptTrack(const AliESDtrack *track, const AliESDEvent *e
   // check TPC nSigma cut\r
   Float_t nsigmaTPC[AliPID::kSPECIES], nsigma;   // need all sigmas for QA plotting\r
   for(Int_t ispec = 0; ispec < AliPID::kSPECIES; ispec++){\r
-    nsigmaTPC[ispec] = nsigma = fESDpid->NumberOfSigmasTPC(track,(AliPID::EParticleType)ispec);\r
+    nsigmaTPC[ispec] = nsigma = fPIDresponse->NumberOfSigmasTPC(track,(AliPID::EParticleType)ispec);\r
     if(!(fTPCsigmaCutRequired & 1 << ispec)) continue;\r
     SETBIT(cutRequired, kCutNsigmaTPC); // We found at least one species where the n-Sigma Cut is required\r
     if(nsigma >= fCutTPCnSigma[2*ispec] && nsigma <= fCutTPCnSigma[2*ispec+1]) SETBIT(cutFullfiled, kCutNsigmaTPC);    // Fullfiled for at least one species\r
@@ -264,7 +281,7 @@ Bool_t AliESDpidCuts::AcceptTrack(const AliESDtrack *track, const AliESDEvent *e
   track->GetIntegratedTimes(times);\r
   for(Int_t ispec = 0; ispec < AliPID::kSPECIES; ispec++){\r
     \r
-    if(hasTOFpid && event) nsigmaTOF[ispec] = nsigma = fESDpid->NumberOfSigmasTOF(track,(AliPID::EParticleType)ispec, event->GetT0());\r
+    if(hasTOFpid && event) nsigmaTOF[ispec] = nsigma = fPIDresponse->NumberOfSigmasTOF(track,(AliPID::EParticleType)ispec);\r
     if(!(fTOFsigmaCutRequired && 1 << ispec)) continue;\r
     SETBIT(cutRequired, kCutNsigmaTOF);\r
     if(track->GetOuterParam()->P() >= fMinMomentumTOF){\r
index d7c9e52d512af90e77028322258779b2b290557d..667ca1992bf975a7c25cc45cb56a80529bf03beb 100644 (file)
@@ -35,6 +35,7 @@ class AliESDpidCuts : public AliAnalysisCuts{
     AliESDpidCuts &operator=(const AliESDpidCuts &ref);\r
     virtual ~AliESDpidCuts();\r
 \r
+    virtual void Init();\r
     virtual void Copy(TObject &c) const;\r
     virtual Long64_t Merge(TCollection *coll);\r
 \r
@@ -46,7 +47,8 @@ class AliESDpidCuts : public AliAnalysisCuts{
     virtual Bool_t IsSelected(TList * /*lst*/) {return kTRUE; }\r
     virtual Bool_t AcceptTrack(const AliESDtrack *track, const AliESDEvent *event);\r
 \r
-    AliESDpid *GetESDpid() { return fESDpid; };\r
+    void SetPIDResponse(AliPIDResponse * pidresponse) { fPIDresponse = pidresponse; }\r
+    AliPIDResponse *GetPIDresponse() { return fPIDresponse; };\r
     \r
     void SetTPCclusterRatioCut(Float_t clr) { fCutTPCclusterRatio = clr; }\r
     inline void SetTPCnSigmaCut(AliPID::EParticleType itype, Float_t nSigma);\r
@@ -57,7 +59,7 @@ class AliESDpidCuts : public AliAnalysisCuts{
   \r
   protected:\r
     static const Int_t kNcuts;                      // Number of Cuts\r
-    AliESDpid *fESDpid;                             //! PID helper (n-sigma-cut)\r
+    AliPIDResponse *fPIDresponse;                   //! PID helper (n-sigma-cut)\r
     Char_t  fTPCsigmaCutRequired;                   // Sigma cut Requirement for TPC and Particle Species\r
     Char_t  fTOFsigmaCutRequired;                   // Sigma cut Requirement for TOF and Particle Species\r
     Float_t fCutTPCnSigma[AliPID::kSPECIES * 2];    // Species dependent cut on the distance to the TPC dE/dx line\r