]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Adding TRD in the combined PID (Yu.Belikov)
authorhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 8 Sep 2004 05:55:08 +0000 (05:55 +0000)
committerhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Wed, 8 Sep 2004 05:55:08 +0000 (05:55 +0000)
STEER/AliESDComparison.C
TRD/AliTRDReconstructor.cxx
TRD/AliTRDpidESD.cxx [new file with mode: 0644]
TRD/AliTRDpidESD.h [new file with mode: 0644]
TRD/TRDrecLinkDef.h
TRD/libTRDrec.pkg

index a5202d2debb335bbbb779aeff5959cdbaa24c5af..8fd54e9df97a6e5e466a2ea9a3860b3ff61d2af2 100644 (file)
@@ -151,6 +151,7 @@ Int_t AliESDComparison(const Char_t *dir=".") {
          UInt_t status=AliESDtrack::kESDpid;
          status|=AliESDtrack::kITSpid; 
          status|=AliESDtrack::kTPCpid; 
+         status|=AliESDtrack::kTRDpid; 
          status|=AliESDtrack::kTOFpid; 
 
         if ((t->GetStatus()&status) == status) {
@@ -165,6 +166,7 @@ Int_t AliESDComparison(const Char_t *dir=".") {
            Int_t code=part->GetPdgCode();
 
            Double_t r[10]; t->GetESDpid(r);
+           //t->GetTRDpid(r);
 
            Double_t rcc=0.;
            Int_t i;
index 5c3bc6e1bf6c10591dc98f81aa0c1df6fb265384..d48f229f9b1e7f20db09fb73eebcf2ff59e224bf 100644 (file)
@@ -27,6 +27,7 @@
 #include "AliTRDparameter.h"
 #include "AliTRDclusterizerV1.h"
 #include "AliTRDtracker.h"
+#include "AliTRDpidESD.h"
 #include <TFile.h>
 
 
@@ -73,8 +74,17 @@ AliTracker* AliTRDReconstructor::CreateTracker(AliRunLoader* runLoader) const
 
 //_____________________________________________________________________________
 void AliTRDReconstructor::FillESD(AliRunLoader* /*runLoader*/, 
-                                 AliESD* /*esd*/) const
+                                 AliESD* esd) const
 {
+// make PID
+
+  Double_t parTRD[] = {
+    187., // Min. Ionizing Particle signal.  Check it !!!
+    0.23, // relative resolution             Check it !!!
+    10.   // PID range (in sigmas)
+  };
+  AliTRDpidESD trdPID(parTRD);
+  trdPID.MakePID(esd);
 }
 
 
diff --git a/TRD/AliTRDpidESD.cxx b/TRD/AliTRDpidESD.cxx
new file mode 100644 (file)
index 0000000..77c6b8b
--- /dev/null
@@ -0,0 +1,84 @@
+/**************************************************************************
+ * Copyright(c) 1998-1999, 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.                  *
+ **************************************************************************/
+
+//-----------------------------------------------------------------
+//           Implementation of the TRD PID class
+// Very naive one... And the implementation is even poorer... 
+// Should be made better by the detector experts...
+//-----------------------------------------------------------------
+
+#include "AliTRDpidESD.h"
+#include "AliESD.h"
+#include "AliESDtrack.h"
+
+ClassImp(AliTRDpidESD)
+
+//_________________________________________________________________________
+AliTRDpidESD::AliTRDpidESD(Double_t *param)
+{
+  //
+  //  The main constructor
+  //
+  fMIP=param[0];   // MIP signal
+  fRes=param[1];   // relative resolution
+  fRange=param[2]; // PID "range" (in sigmas)
+}
+
+Double_t AliTRDpidESD::Bethe(Double_t bg) {
+  //
+  // This is the Bethe-Bloch function normalised to 1 at the minimum
+  //
+  Double_t bg2=bg*bg;
+  Double_t bethe;
+  if (bg<3.5e1) 
+      bethe=(1.+ bg2)/bg2*(log(5940*bg2) - bg2/(1.+ bg2));
+  else // Density effect ( approximately :) 
+      bethe=(1.+ bg2)/bg2*(log(3.5*5940*bg) - bg2/(1.+ bg2));
+  return bethe/11.091;
+}
+
+//_________________________________________________________________________
+Int_t AliTRDpidESD::MakePID(AliESD *event)
+{
+  //
+  //  This function calculates the "detector response" PID probabilities 
+  //
+  static const Double_t masses[]={
+    0.000511, 0.105658, 0.139570, 0.493677, 0.938272, 1.875613
+  };
+  Int_t ntrk=event->GetNumberOfTracks();
+  for (Int_t i=0; i<ntrk; i++) {
+    AliESDtrack *t=event->GetTrack(i);
+    if ((t->GetStatus()&AliESDtrack::kTRDin)==0)
+       if ((t->GetStatus()&AliESDtrack::kTRDout)==0)
+          if ((t->GetStatus()&AliESDtrack::kTRDrefit)==0) continue;
+    Int_t ns=AliESDtrack::kSPECIES;
+    Double_t p[10];
+    for (Int_t j=0; j<ns; j++) {
+      Double_t mass=masses[j];
+      Double_t mom=t->GetP();
+      Double_t dedx=t->GetTRDsignal()/fMIP;
+      Double_t bethe=Bethe(mom/mass); 
+      Double_t sigma=fRes*bethe;
+      if (TMath::Abs(dedx-bethe) > fRange*sigma) {
+       p[j]=TMath::Exp(-0.5*fRange*fRange)/sigma;
+        continue;
+      }
+      p[j]=TMath::Exp(-0.5*(dedx-bethe)*(dedx-bethe)/(sigma*sigma))/sigma;
+    }
+    t->SetTRDpid(p);
+  }
+  return 0;
+}
diff --git a/TRD/AliTRDpidESD.h b/TRD/AliTRDpidESD.h
new file mode 100644 (file)
index 0000000..2792045
--- /dev/null
@@ -0,0 +1,29 @@
+#ifndef ALITRDPIDESD_H
+#define ALITRDPIDESD_H
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+//-------------------------------------------------------
+//                    TRD PID class
+// A very naive design... And the implementation is even poorer... 
+// Should be made better by the detector experts...
+//-------------------------------------------------------
+#include <Rtypes.h>
+
+class AliESD;
+
+class AliTRDpidESD {
+public:
+  AliTRDpidESD(Double_t *param);
+  Int_t MakePID(AliESD *event);
+  static Double_t Bethe(Double_t bg);
+private:
+  Double_t fMIP;          // dEdx for MIP
+  Double_t fRes;          // relative dEdx resolution
+  Double_t fRange;        // one particle type PID range (in sigmas)
+  ClassDef(AliTRDpidESD,1)   // TRD PID class
+};
+
+#endif
+
+
index 913e21c80eb8abced2eb0efafe44b819fb754cdb..03f8f6f75873f8052c6b6e8cf3d550fbd0042f71 100644 (file)
@@ -27,6 +27,7 @@
 #pragma link C++ class  AliTRDpid+;
 #pragma link C++ class  AliTRDpidLQ+;
 #pragma link C++ class  AliTRDPartID+;
+#pragma link C++ class  AliTRDpidESD+;
 
 #pragma link C++ class  AliTRDReconstructor+;
 
index bbbeaba22f35a53a31606497ec614a85bf41cc48..80f12c91a00d773ef97a95fdeca7501aeb26157f 100644 (file)
@@ -12,6 +12,7 @@ SRCS= AliTRDpixel.cxx \
       AliTRDpid.cxx \
       AliTRDpidLQ.cxx \
       AliTRDPartID.cxx \
+      AliTRDpidESD.cxx \
       AliTRDReconstructor.cxx
 
 HDRS= $(SRCS:.cxx=.h)