TRD included in the ESD chain (Yu.Belikov)
authorhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 13 Nov 2003 18:41:49 +0000 (18:41 +0000)
committerhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Thu, 13 Nov 2003 18:41:49 +0000 (18:41 +0000)
STEER/AliESDanalysis.C
STEER/AliESDtest.C
STEER/AliESDtrack.cxx
STEER/AliESDtrack.h

index 7355d2a..a29cf0f 100644 (file)
@@ -26,31 +26,31 @@ extern AliRun *gAlice;
 
 Int_t AliESDanalysis(Int_t nev=1) { 
   TH2F *tpcHist=
-     new TH2F("tpcHist","TPC dE/dX vs momentum",100,0.,3.,100,0.,500.);
+     new TH2F("tpcHist","TPC dE/dX vs momentum",100,0.,4.,100,0.,500.);
    tpcHist->SetXTitle("p (GeV/c)"); tpcHist->SetYTitle("dE/dx (Arb. Units)");
      tpcHist->SetMarkerStyle(8); 
      tpcHist->SetMarkerSize(0.3);
  
-  TH1F *piG=new TH1F("piG","",20,0.,3.);
-  TH1F *piR=new TH1F("piR","",20,0.,3.);
-  TH1F *piF=new TH1F("piF","",20,0.,3.);
-  TH1F *piGood=new TH1F("piGood","Combined PID for pions",20,0.,3.); 
+  TH1F *piG=new TH1F("piG","",20,0.,4.);
+  TH1F *piR=new TH1F("piR","",20,0.,4.);
+  TH1F *piF=new TH1F("piF","",20,0.,4.);
+  TH1F *piGood=new TH1F("piGood","Combined PID for pions",20,0.,4.); 
   piGood->SetLineColor(4); piGood->SetXTitle("p (GeV/c)");
-  TH1F *piFake=new TH1F("piFake","",20,0.,3.); piFake->SetLineColor(2);
+  TH1F *piFake=new TH1F("piFake","",20,0.,4.); piFake->SetLineColor(2);
 
-  TH1F *kaG=new TH1F("kaG","",20,0.,3.);
-  TH1F *kaR=new TH1F("kaR","",20,0.,3.);
-  TH1F *kaF=new TH1F("kaF","",20,0.,3.);
-  TH1F *kaGood=new TH1F("kaGood","Combined PID for kaons",20,0.,3.); 
+  TH1F *kaG=new TH1F("kaG","",20,0.,4.);
+  TH1F *kaR=new TH1F("kaR","",20,0.,4.);
+  TH1F *kaF=new TH1F("kaF","",20,0.,4.);
+  TH1F *kaGood=new TH1F("kaGood","Combined PID for kaons",20,0.,4.); 
   kaGood->SetLineColor(4); kaGood->SetXTitle("p (GeV/c)");
-  TH1F *kaFake=new TH1F("kaFake","",20,0.,3.); kaFake->SetLineColor(2);
+  TH1F *kaFake=new TH1F("kaFake","",20,0.,4.); kaFake->SetLineColor(2);
 
-  TH1F *prG=new TH1F("prG","",20,0.,3.);
-  TH1F *prR=new TH1F("prR","",20,0.,3.);
-  TH1F *prF=new TH1F("prF","",20,0.,3.);
-  TH1F *prGood=new TH1F("prGood","Combined PID for protons",20,0.,3.); 
+  TH1F *prG=new TH1F("prG","",20,0.,4.);
+  TH1F *prR=new TH1F("prR","",20,0.,4.);
+  TH1F *prF=new TH1F("prF","",20,0.,4.);
+  TH1F *prGood=new TH1F("prGood","Combined PID for protons",20,0.,4.); 
   prGood->SetLineColor(4); prGood->SetXTitle("p (GeV/c)");
-  TH1F *prFake=new TH1F("prFake","",20,0.,3.); prFake->SetLineColor(2);
+  TH1F *prFake=new TH1F("prFake","",20,0.,4.); prFake->SetLineColor(2);
 
    if (gAlice) {
       delete gAlice->GetRunLoader();
index 3ddf801..8683b03 100644 (file)
@@ -91,7 +91,6 @@ Int_t AliESDtest(Int_t nev=1) {
    );
 
 
-
 /**** The ITS corner ********************/
 
    AliITSLoader* itsl = (AliITSLoader*)rl->GetLoader("ITSLoader");
@@ -220,6 +219,8 @@ Int_t AliESDtest(Int_t nev=1) {
 
      rl->GetEvent(i);
  
+
+//***** Initial path towards the primary vertex
      TTree *tpcTree=tpcl->TreeR();
      if (!tpcTree) {
         cerr<<"Can't get the TPC cluster tree !\n";
@@ -228,43 +229,35 @@ Int_t AliESDtest(Int_t nev=1) {
      tpcTracker.LoadClusters(tpcTree);
      rc+=tpcTracker.Clusters2Tracks(event);
 
-
      TTree *itsTree=itsl->TreeR();
      if (!itsTree) {
-        cerr<<"Can't get the TPC cluster tree !\n";
+        cerr<<"Can't get the ITS cluster tree !\n";
         return 4;
      }     
      itsTracker.LoadClusters(itsTree);
      rc+=itsTracker.Clusters2Tracks(event);
 
-     rc+=vtxer.Tracks2V0vertices(event);            // V0 finding
-     rc+=cvtxer.V0sTracks2CascadeVertices(event);   // cascade finding
 
+//***** Back propagation towards the outer barrel detectors
      rc+=itsTracker.PropagateBack(event); 
-     itsTracker.UnloadClusters();
      //itsPID.MakePID(event);
      
      rc+=tpcTracker.PropagateBack(event);
-     tpcTracker.UnloadClusters();
      tpcPID.MakePID(event);
 
-
      TTree *trdTree=trdl->TreeR();
      if (!trdTree) {
-        cerr<<"Can't get the TPC cluster tree !\n";
+        cerr<<"Can't get the TRD cluster tree !\n";
         return 4;
      } 
      trdTracker.LoadClusters(trdTree);
      rc+=trdTracker.PropagateBack(event);
-     trdTracker.UnloadClusters();
-
 /*
      for (Int_t iTrack = 0; iTrack < event->GetNumberOfTracks(); iTrack++) {
        AliESDtrack* track = event->GetTrack(iTrack);
        trdPID->MakePID(track);
      }
 */
-
      TTree *tofTree=tofl->TreeD();
      if (!tofTree) {
         cerr<<"Can't get the TOF cluster tree !\n";
@@ -275,9 +268,31 @@ Int_t AliESDtest(Int_t nev=1) {
      tofPID.UnloadClusters();
 
 
-    //Here is the combined PID
+
+//***** Here is the combined PID
      AliESDpid::MakePID(event);
 
+
+
+//***** Now the final refit at the primary vertex...
+     rc+=trdTracker.RefitInward(event);
+     trdTracker.UnloadClusters();
+
+     rc+=tpcTracker.RefitInward(event);
+     tpcTracker.UnloadClusters();
+
+     rc+=itsTracker.RefitInward(event); 
+     itsTracker.UnloadClusters();
+
+
+
+//***** Hyperon reconstruction 
+     rc+=vtxer.Tracks2V0vertices(event);            // V0 finding
+     rc+=cvtxer.V0sTracks2CascadeVertices(event);   // cascade finding
+
+
+
+//***** Some final manipulations with this event 
      if (rc==0) {
         Char_t ename[100]; 
         sprintf(ename,"%d",i);
index e632663..1de3543 100644 (file)
@@ -25,7 +25,6 @@
 
 #include "AliESDtrack.h"
 #include "AliKalmanTrack.h"
-// #include "../ITS/AliITStrackV2.h"
 
 ClassImp(AliESDtrack)
 
@@ -55,7 +54,7 @@ fTRDncls(0),
 fTRDsignal(0),
 fTOFchi2(0),
 fTOFindex(0),
-fTOFsignal(0)
+fTOFsignal(-1)
 {
   //
   // The default ESD constructor 
@@ -68,10 +67,11 @@ fTOFsignal(0)
     fTRDr[i]=0;
     fTOFr[i]=0;
   }
-  for (Int_t i=0; i<5; fRp[i++]);
-  for (Int_t i=0; i<15; fRc[i++]);
-  for (Int_t i=0; i<6; fITSindex[i++]);
-  for (Int_t i=0; i<180; fTPCindex[i++]);
+  Int_t i;
+  for (i=0; i<5; i++)   fRp[i]=0.;
+  for (i=0; i<15; i++)  fRc[i]=0.;
+  for (i=0; i<6; i++)   fITSindex[i]=0;
+  for (i=0; i<180; i++) fTPCindex[i]=0;
 }
 
 //_______________________________________________________________________
@@ -120,6 +120,7 @@ Bool_t AliESDtrack::UpdateTrackParams(AliKalmanTrack *t, ULong_t flags) {
   case kTRDin: case kTRDout: case kTRDrefit:
     fTRDncls=t->GetNumberOfClusters();
     fTRDchi2=t->GetChi2();
+    for (Int_t i=0;i<fTRDncls;i++) fTRDindex[i]=t->GetClusterIndex(i);
     fTRDsignal=t->GetPIDsignal();
     break;
   default: 
@@ -142,7 +143,6 @@ Bool_t AliESDtrack::UpdateTrackParams(AliKalmanTrack *t, ULong_t flags) {
   
   if (flags == kITSin)
    {
-     //     AliITStrackV2* itstrack = dynamic_cast<AliITStrackV2*>(t);
      AliKalmanTrack *itstrack = t;
      if (itstrack)
       {
@@ -172,24 +172,6 @@ Bool_t AliESDtrack::UpdateTrackParams(AliKalmanTrack *t, ULong_t flags) {
   return kTRUE;
 }
 
-//_______________________________________________________________________
-void AliESDtrack::GetExternalParametersAt(Double_t x, Double_t p[5]) const {
-  //---------------------------------------------------------------------
-  // This function returns external representation of the track parameters
-  // at the plane x
-  //---------------------------------------------------------------------
-  Double_t dx=x-fRx;
-  Double_t c=fRp[4]/AliKalmanTrack::GetConvConst();
-  Double_t f1=fRp[2], f2=f1 + c*dx;
-  Double_t r1=sqrt(1.- f1*f1), r2=sqrt(1.- f2*f2);    
-
-  p[0]=fRp[0]+dx*(f1+f2)/(r1+r2);
-  p[1]=fRp[1]+dx*(f1+f2)/(f1*r2 + f2*r1)*fRp[3];
-  p[2]=fRp[2]+dx*c;
-  p[3]=fRp[3];
-  p[4]=fRp[4];
-}
-
 //_______________________________________________________________________
 void AliESDtrack::GetExternalParameters(Double_t &x, Double_t p[5]) const {
   //---------------------------------------------------------------------
@@ -221,7 +203,7 @@ void AliESDtrack::GetXYZ(Double_t *xyz) const {
   //---------------------------------------------------------------------
   // This function returns the global track position
   //---------------------------------------------------------------------
-  Double_t phi=TMath::ASin(fRp[2]) + fRalpha;
+  Double_t phi=TMath::ATan2(fRp[0],fRx) + fRalpha;
   Double_t r=TMath::Sqrt(fRx*fRx + fRp[0]*fRp[0]);
   xyz[0]=r*TMath::Cos(phi); xyz[1]=r*TMath::Sin(phi); xyz[2]=fRp[1]; 
 }
@@ -290,6 +272,15 @@ void AliESDtrack::GetTPCpid(Double_t *p) const {
   for (Int_t i=0; i<kSPECIES; i++) p[i]=fTPCr[i];
 }
 
+//_______________________________________________________________________
+Int_t AliESDtrack::GetTRDclusters(UInt_t *idx) const {
+  //---------------------------------------------------------------------
+  // This function returns indices of the assgined TRD clusters 
+  //---------------------------------------------------------------------
+  for (Int_t i=0; i<90; i++) idx[i]=fTRDindex[i];  // MI I prefer some constant
+  return fTRDncls;
+}
+
 //_______________________________________________________________________
 void AliESDtrack::SetTRDpid(const Double_t *p) {  
   // Sets values for the probability of each particle type (in TRD)
index 44ed3c7..798e579 100644 (file)
@@ -28,7 +28,6 @@ public:
   ULong_t GetStatus() const {return fFlags;}
   Int_t GetLabel() const {return fLabel;}
   Double_t GetAlpha() const {return fRalpha;}
-  void GetExternalParametersAt(Double_t x, Double_t p[5]) const;
   void GetExternalParameters(Double_t &x, Double_t p[5]) const;
   void GetExternalCovariance(Double_t cov[15]) const;
   Double_t GetIntegratedLength() const {return fTrackLength;}
@@ -51,6 +50,7 @@ public:
   void SetTRDpid(const Double_t *p);
   void GetTRDpid(Double_t *p) const;
   Float_t GetTRDsignal() const {return fTRDsignal;}
+  Int_t GetTRDclusters(UInt_t *idx) const;
   void    SetTRDpid(Int_t iSpecies, Float_t p);
   Float_t GetTRDpid(Int_t iSpecies) const;
 
@@ -124,12 +124,13 @@ protected:
   // TRD related track information
   Float_t fTRDchi2;        // chi2 in the TRD
   Int_t   fTRDncls;        // number of clusters assigned in the TRD
+  UInt_t  fTRDindex[90];   //! indices of the assigned TRD clusters
   Float_t fTRDsignal;      // detector's PID signal
   Float_t fTRDr[kSPECIES]; //! "detector response probabilities" (for the PID)
 
   // TOF related track information
   Float_t fTOFchi2;        // chi2 in the TOF
-  UInt_t  fTOFindex;       //! index of the assigned TOF cluster
+  UInt_t  fTOFindex;       // index of the assigned TOF cluster
   Float_t fTOFsignal;      // detector's PID signal
   Float_t fTOFr[kSPECIES]; // "detector response probabilities" (for the PID)