New design of tracking classes (Yu.Belikov)
authorhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 13 Mar 2001 14:26:39 +0000 (14:26 +0000)
committerhristov <hristov@f7af4fe6-9843-0410-8265-dc069ae4e863>
Tue, 13 Mar 2001 14:26:39 +0000 (14:26 +0000)
13 files changed:
STEER/AliCluster.h
STEER/AliKalmanTrack.cxx
STEER/AliKalmanTrack.h
STEER/AliTracker.cxx [new file with mode: 0644]
STEER/AliTracker.h [new file with mode: 0644]
STEER/Makefile
STEER/STEERLinkDef.h
TPC/AliTPCComparison.C
TPC/AliTPCtrack.cxx
TPC/AliTPCtrack.h
TPC/AliTPCtracker.cxx
TPC/AliTPCtracker.h
TPC/TPCLinkDef.h

index c4fc46e..dcb94a0 100644 (file)
@@ -15,6 +15,7 @@
 class AliCluster : public TObject {
 public:
   AliCluster();
+  virtual ~AliCluster() {;}
   AliCluster(Int_t *lab, Float_t *hit);
   void SetLabel(Int_t lab, Int_t i) {fTracks[i]=lab;}
   void SetY(Float_t y)              {fY=y;}
@@ -28,6 +29,8 @@ public:
   Float_t GetSigmaY2()    const {return fSigmaY2;}
   Float_t GetSigmaZ2()    const {return fSigmaZ2;}
 
+  virtual void Use() = 0;
+
 protected:
   Int_t     fTracks[3];//labels of overlapped tracks
   Float_t   fY ;       //Y of cluster
index f552202..e5c4379 100644 (file)
 //-------------------------------------------------------------------------
 
 #include "AliKalmanTrack.h"
-#include "AliCluster.h"
-#include <TMath.h>
-#include <iostream.h>
 
 ClassImp(AliKalmanTrack)
 
-//_____________________________________________________________________________
-AliKalmanTrack::AliKalmanTrack(const AliKalmanTrack& t) {
-  //-----------------------------------------------------------------
-  // This is a copy constructor.
-  //-----------------------------------------------------------------
-  fLab=t.fLab;
-
-  fP0=t.fP0; fP1=t.fP1; fP2=t.fP2; fP3=t.fP3; fP4=t.fP4;
-
-  fC00=t.fC00;
-  fC10=t.fC10;  fC11=t.fC11;
-  fC20=t.fC20;  fC21=t.fC21;  fC22=t.fC22;
-  fC30=t.fC30;  fC31=t.fC31;  fC32=t.fC32;  fC33=t.fC33;
-  fC40=t.fC40;  fC41=t.fC41;  fC42=t.fC42;  fC43=t.fC43;  fC44=t.fC44;
-
-  fChi2=t.fChi2;
-  fN=t.fN;
-}
-
-//_____________________________________________________________________________
-Int_t AliKalmanTrack::Compare(const TObject *o) const {
-  //-----------------------------------------------------------------
-  // This function compares tracks according to the their curvature
-  //-----------------------------------------------------------------
-  AliKalmanTrack *t=(AliKalmanTrack*)o;
-  Double_t co=TMath::Abs(t->GetPt());
-  Double_t c =TMath::Abs(GetPt());
-  if (c<co) return 1;
-  else if (c>co) return -1;
-  return 0;
-}
-
-//_____________________________________________________________________________
-Double_t AliKalmanTrack::GetPredictedChi2(const AliCluster *c) const 
-{
-  //-----------------------------------------------------------------
-  // This function calculates a predicted chi2 increment.
-  //-----------------------------------------------------------------
-  Double_t r00=c->GetSigmaY2(), r01=0., r11=c->GetSigmaZ2();
-  r00+=fC00; r01+=fC10; r11+=fC11;
-
-  Double_t det=r00*r11 - r01*r01;
-  if (TMath::Abs(det) < 1.e-10) {
-    if (fN>4) cerr<<fN<<" AliKalmanTrack warning: Singular matrix !\n";
-    return 1e10;
-  }
-  Double_t tmp=r00; r00=r11; r11=tmp; r01=-r01;
-  
-  Double_t dy=c->GetY() - fP0, dz=c->GetZ() - fP1;
-  
-  return (dy*r00*dy + 2*r01*dy*dz + dz*r11*dz)/det;
-}
-
-//_____________________________________________________________________________
-void AliKalmanTrack::GetCovariance(Double_t cc[15]) const {
-  // return covariance maxtrix
-  cc[0 ]=fC00;
-  cc[1 ]=fC10;  cc[2 ]=fC11;
-  cc[3 ]=fC20;  cc[4 ]=fC21;  cc[5 ]=fC22;
-  cc[6 ]=fC30;  cc[7 ]=fC31;  cc[8 ]=fC32;  cc[9 ]=fC33;
-  cc[10]=fC40;  cc[11]=fC41;  cc[12]=fC42;  cc[13]=fC43;  cc[14]=fC44;
-}
-
-
-
index 1c637d5..a4da404 100644 (file)
@@ -16,42 +16,38 @@ class AliCluster;
 
 class AliKalmanTrack : public TObject {
 public:
-  AliKalmanTrack() {fN=0; fChi2=0; fLab=-3141593;}
-  AliKalmanTrack(const AliKalmanTrack& t);
-  virtual ~AliKalmanTrack() {}
-  Int_t Compare(const TObject *o) const;
-  void SetLabel(Int_t lab) {fLab=lab;} 
-
-  Double_t GetPredictedChi2(const AliCluster *cluster) const;
-  Bool_t IsSortable() const {return kTRUE;}
-  Int_t GetLabel() const {return fLab;}
-  void GetCovariance(Double_t cov[15]) const;
-  Double_t GetChi2() const {return fChi2;}
-  Int_t GetNumberOfClusters() const {return fN;}
-
- virtual Double_t GetPt() const=0;
- virtual Double_t GetP()  const=0;
- virtual void GetPxPyPz(Double_t &px, Double_t &py, Double_t &pz) const=0;
- virtual Int_t PropagateTo(Double_t xr,Double_t x0,Double_t rho,Double_t pm)=0;
- virtual void Update(const AliCluster* c, Double_t chi2, UInt_t i)=0;
-
-protected: 
+  AliKalmanTrack() {fLab=-3141593; fChi2=0; fN=0;}
+  AliKalmanTrack(const AliKalmanTrack &t) {fLab=t.fLab;fChi2=t.fChi2;fN=t.fN;}
+  virtual ~AliKalmanTrack(){};
+  void SetLabel(Int_t lab) {fLab=lab;}
+
+  Bool_t   IsSortable() const {return kTRUE;}
+  Int_t    GetLabel()   const {return fLab;}
+  Double_t GetChi2()    const {return fChi2;}
+  Int_t    GetNumberOfClusters() const {return fN;}
+  virtual Int_t GetClusterIndex(Int_t i) const { //reserved for AliTracker
+    printf("AliKalmanTrack::GetClusterIndex(Int_t i) must be overloaded !\n");
+    return 0;
+  } 
+
+  virtual Int_t Compare(const TObject *o) const=0; 
+
+  virtual void GetExternalParameters(Double_t &xr, Double_t x[5]) const=0;
+  virtual void GetExternalCovariance(Double_t cov[15]) const=0;
+
+  virtual Double_t GetPredictedChi2(const AliCluster *cluster) const=0;
+  virtual 
+  Int_t PropagateTo(Double_t xr,Double_t x0,Double_t rho,Double_t pm)=0;
+  virtual Int_t Update(const AliCluster* c, Double_t chi2, UInt_t i)=0;
+
+protected:
+  void SetChi2(Double_t chi2) {fChi2=chi2;} 
+  void SetNumberOfClusters(Int_t n) {fN=n;} 
+
+private: 
   Int_t fLab;             // track label
-
-  Double_t fP0;           // track parameter
-  Double_t fP1;           // track parameter
-  Double_t fP2;           // track parameter
-  Double_t fP3;           // track parameter
-  Double_t fP4;           // track parameter
-
-  Double_t fC00;                         // covariance
-  Double_t fC10, fC11;                   // matrix
-  Double_t fC20, fC21, fC22;             // of the
-  Double_t fC30, fC31, fC32, fC33;       // track
-  Double_t fC40, fC41, fC42, fC43, fC44; // parameters
-
   Double_t fChi2;         // total chi2 value for this track
-  Short_t fN;             // number of clusters 
+  Int_t fN;               // number of associated clusters
 
   ClassDef(AliKalmanTrack,1)    // Reconstructed track
 };
diff --git a/STEER/AliTracker.cxx b/STEER/AliTracker.cxx
new file mode 100644 (file)
index 0000000..c2d7f01
--- /dev/null
@@ -0,0 +1,86 @@
+/**************************************************************************
+ * 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 AliTracker class
+//
+//          Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
+//-------------------------------------------------------------------------
+#include <TMath.h>
+
+#include "AliTracker.h"
+#include "AliCluster.h"
+#include "AliKalmanTrack.h"
+
+ClassImp(AliTracker)
+
+//__________________________________________________________________________
+void AliTracker::CookLabel(AliKalmanTrack *t, Float_t wrong) const {
+  //--------------------------------------------------------------------
+  //This function "cooks" a track label. If label<0, this track is fake.
+  //--------------------------------------------------------------------
+  Int_t noc=t->GetNumberOfClusters();
+  Int_t *lb=new Int_t[noc];
+  Int_t *mx=new Int_t[noc];
+  AliCluster **clusters=new AliCluster*[noc];
+
+  Int_t i;
+  for (i=0; i<noc; i++) {
+     lb[i]=mx[i]=0;
+     Int_t index=t->GetClusterIndex(i);
+     clusters[i]=GetCluster(index);
+  }
+
+  Int_t lab=123456789;
+  for (i=0; i<noc; i++) {
+    AliCluster *c=clusters[i];
+    lab=TMath::Abs(c->GetLabel(0));
+    Int_t j;
+    for (j=0; j<noc; j++) if (lb[j]==lab || mx[j]==0) break;
+    lb[j]=lab;
+    (mx[j])++;
+  }
+
+  Int_t max=0;
+  for (i=0; i<noc; i++) if (mx[i]>max) {max=mx[i]; lab=lb[i];}
+    
+  for (i=0; i<noc; i++) {
+    AliCluster *c=clusters[i];
+    if (TMath::Abs(c->GetLabel(1)) == lab ||
+        TMath::Abs(c->GetLabel(2)) == lab ) max++;
+  }
+
+  if ((1.- Float_t(max)/noc) > wrong) lab=-lab;
+  t->SetLabel(lab);
+
+  delete[] lb;
+  delete[] mx;
+  delete[] clusters;
+}
+
+//____________________________________________________________________________
+void AliTracker::UseClusters(const AliKalmanTrack *t, Int_t from) const {
+  //------------------------------------------------------------------
+  //This function marks clusters associated with the track.
+  //------------------------------------------------------------------
+  Int_t noc=t->GetNumberOfClusters();
+  for (Int_t i=from; i<noc; i++) {
+     Int_t index=t->GetClusterIndex(i);
+     AliCluster *c=GetCluster(index); 
+     c->Use();   
+  }
+}
+
+
diff --git a/STEER/AliTracker.h b/STEER/AliTracker.h
new file mode 100644 (file)
index 0000000..0f13835
--- /dev/null
@@ -0,0 +1,34 @@
+#ifndef ALITRACKER_H
+#define ALITRACKER_H
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+//-------------------------------------------------------------------------
+//                          class AliTracker
+//
+//       Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch 
+//-------------------------------------------------------------------------
+#include <Rtypes.h>
+
+class AliKalmanTrack;
+class AliCluster;
+class TFile;
+
+class AliTracker {
+public:
+  AliTracker(){}
+  virtual ~AliTracker(){}
+  virtual Int_t Clusters2Tracks(const TFile *in, TFile *out)=0;
+  virtual Int_t PropagateBack(const TFile *in, TFile *out)=0;
+
+//protected:
+  virtual AliCluster *GetCluster(Int_t index) const=0;
+  virtual void  UseClusters(const AliKalmanTrack *t, Int_t from=0) const;
+  virtual void  CookLabel(AliKalmanTrack *t,Float_t wrong) const; 
+
+  ClassDef(AliTracker,1) //abstract tracker
+};
+
+#endif
+
+
index ee025b2..28878bc 100644 (file)
@@ -18,7 +18,7 @@ SRCS          = AliDetector.cxx       AliHeader.cxx   AliMagF.cxx \
                AliMagFDM.cxx   AliLegoGenerator.cxx AliLegoGeneratorXYZ.cxx\
                 AliLegoGeneratorPhiZ.cxx AliLegoGeneratorEta.cxx \
                 AliRndm.cxx \
-               AliKalmanTrack.cxx AliCluster.cxx \
+               AliKalmanTrack.cxx AliCluster.cxx AliTracker.cxx\
                AliMCQA.cxx AliPDG.cxx
 
 # C++ Headers
index 6c95ec3..281a323 100644 (file)
@@ -40,6 +40,7 @@
 #pragma link C++ class  AliHitMap+;
 #pragma link C++ class  AliCluster+;
 #pragma link C++ class  AliKalmanTrack+;
+#pragma link C++ class  AliTracker+;
 #pragma link C++ class  AliRndm+;
 #pragma link C++ class  AliMCQA+;
 #pragma link C++ class  AliPDG+;
index 589d2da..c595a55 100644 (file)
@@ -1,5 +1,7 @@
-//#include "alles.h"
-//#include "TParticle.h"
+#ifndef __CINT__
+  #include "alles.h"
+#endif
+
 struct GoodTrack {
   Int_t lab;
   Int_t code;
@@ -11,7 +13,6 @@ Int_t good_tracks(GoodTrack *gt, Int_t max);
 Int_t AliTPCComparison() {
   Int_t i;
    gBenchmark->Start("AliTPCComparison");
-   cerr<<"Doing comparison...\n";
 
    TFile *cf=TFile::Open("AliTPCclusters.root");
    if (!cf->IsOpen()) {cerr<<"Can't open AliTPCclusters.root !\n"; return 1;}
@@ -24,10 +25,7 @@ Int_t AliTPCComparison() {
    ca->SetClusterType("AliTPCcluster");
    ca->ConnectTree("Segment Tree");
    Int_t nentr=Int_t(ca->GetTree()->GetEntries());
-   for (i=0; i<nentr; i++) {
-     //AliSegmentID *s=;
-       ca->LoadEntry(i);
-   }
+   for (i=0; i<nentr; i++) ca->LoadEntry(i);
 
 // Load tracks
    TFile *tf=TFile::Open("AliTPCtracks.root");
@@ -37,8 +35,9 @@ Int_t AliTPCComparison() {
    TBranch *tbranch=tracktree->GetBranch("tracks");
    nentr=(Int_t)tracktree->GetEntries();
    
+   AliTPCtrack *iotrack=0;
    for (i=0; i<nentr; i++) {
-       AliTPCtrack *iotrack=new AliTPCtrack;
+       iotrack=new AliTPCtrack;
        tbranch->SetAddress(&iotrack);
        tracktree->GetEvent(i);
        iotrack->CookLabel(ca);
@@ -52,10 +51,12 @@ Int_t AliTPCComparison() {
 /////////////////////////////////////////////////////////////////////////
    GoodTrack gt[15000];
    Int_t ngood=0;
-   ifstream in("good_tracks");
+   ifstream in("good_tracks_tpc");
    if (in) {
       cerr<<"Reading good tracks...\n";
-      while (in>>gt[ngood].lab>>gt[ngood].code>>gt[ngood].px>>gt[ngood].py>>gt[ngood].pz>>gt[ngood].x>>gt[ngood].y>>gt[ngood].z) {
+      while (in>>gt[ngood].lab>>gt[ngood].code>>
+                 gt[ngood].px>>gt[ngood].py>>gt[ngood].pz>>
+                 gt[ngood].x >>gt[ngood].y >>gt[ngood].z) {
          ngood++;
          cerr<<ngood<<'\r';
          if (ngood==15000) {
@@ -63,16 +64,31 @@ Int_t AliTPCComparison() {
             break;
          }
       }
-      if (!in.eof()) cerr<<"Read error (good_tracks) !\n";
+      if (!in.eof()) cerr<<"Read error (good_tracks_tpc) !\n";
    } else {
       cerr<<"Marking good tracks (this will take a while)...\n";
       ngood=good_tracks(gt,15000);
-      ofstream out("good_tracks");
+      ofstream out("good_tracks_tpc");
       if (out) {
          for (Int_t ngd=0; ngd<ngood; ngd++)            
-           out<<gt[ngd].lab<<' '<<gt[ngd].code<<' '<<gt[ngd].px <<' '<<gt[ngd].py<<' '<<gt[ngd].pz<<' '<<gt[ngd].x  <<' '<<gt[ngd].y <<' '<<gt[ngd].z <<endl;
-      } else cerr<<"Can not open file (good_tracks) !\n";
+           out<<gt[ngd].lab<<' '<<gt[ngd].code<<' '<<
+                 gt[ngd].px<<' '<<gt[ngd].py<<' '<<gt[ngd].pz<<' '<<
+                 gt[ngd].x <<' '<<gt[ngd].y <<' '<<gt[ngd].z <<endl;
+      } else cerr<<"Can not open file (good_tracks_tpc) !\n";
       out.close();
+
+      cerr<<"Preparing tracks for matching with the ITS...\n";
+      tarray.Sort();
+      tf=TFile::Open("AliTPCtracks.root","recreate");
+      tracktree=new TTree("TreeT","Tree with TPC tracks");
+      tracktree->Branch("tracks","AliTPCtrack",&iotrack,32000,0);
+      for (i=0; i<nentr; i++) {
+          iotrack=(AliTPCtrack*)tarray.UncheckedAt(i);
+          tracktree->Fill();
+      }
+      tracktree->Write();
+      tf->Close();
+
    }
    cerr<<"Number of good tracks : "<<ngood<<endl;
 
@@ -84,16 +100,18 @@ Int_t AliTPCComparison() {
    TH1F *hmpt=new TH1F("hmpt","Relative Pt resolution (pt>4GeV/c)",30,-60,60); 
    hmpt->SetFillColor(6);
 
-   TH1F *hgood=new TH1F("hgood","Good tracks",20,0,7);    
-   TH1F *hfound=new TH1F("hfound","Found tracks",20,0,7);
-   TH1F *hfake=new TH1F("hfake","Fake tracks",20,0,7);
-   TH1F *hg=new TH1F("hg","",20,0,7); //efficiency for good tracks
+   TH1F *hgood=new TH1F("hgood","Good tracks",30,0.1,6.1);    
+   TH1F *hfound=new TH1F("hfound","Found tracks",30,0.1,6.1);
+   TH1F *hfake=new TH1F("hfake","Fake tracks",30,0.1,6.1);
+   TH1F *hg=new TH1F("hg","",30,0.1,6.1); //efficiency for good tracks
    hg->SetLineColor(4); hg->SetLineWidth(2);
-   TH1F *hf=new TH1F("hf","Efficiency for fake tracks",20,0,7);
+   TH1F *hf=new TH1F("hf","Efficiency for fake tracks",30,0.1,6.1);
    hf->SetFillColor(1); hf->SetFillStyle(3013); hf->SetLineWidth(2);
 
    TH1F *he =new TH1F("he","dE/dX for pions with 0.4<p<0.5 GeV/c",50,0.,100.);
    TH2F *hep=new TH2F("hep","dE/dX vs momentum",50,0.,2.,50,0.,200.);
+   hep->SetMarkerStyle(8);
+   hep->SetMarkerSize(0.4);
 
    Int_t mingood = ngood;  //MI change
    Int_t * track_notfound = new Int_t[ngood];
@@ -118,7 +136,7 @@ Int_t AliTPCComparison() {
           if (lab==TMath::Abs(tlab)) break;
       }
       if (i==nentr) {
-       //      cerr<<"Track "<<lab<<" was not found !\n";
+       //cerr<<"Track "<<lab<<" was not found !\n";
        track_notfound[itrack_notfound++]=lab;
         continue;
       }
@@ -150,23 +168,26 @@ Int_t AliTPCComparison() {
        //cerr<<lab<<" fake\n";
       }
 
-      Double_t px,py,pz,pt=fabs(track->GetPt());track->GetPxPyPz(px,py,pz);
+      Double_t par[5]; track->GetExternalParameters(xk,par);
+      Float_t phi=TMath::ASin(par[2]) + track->GetAlpha();
+      if (phi<-TMath::Pi()) phi+=2*TMath::Pi();
+      if (phi>=TMath::Pi()) phi-=2*TMath::Pi();
+      Float_t lam=TMath::ATan(par[3]); 
+      Float_t pt_1=TMath::Abs(par[4]);
 
       if (TMath::Abs(gt[ngood].code)==11 && ptg>4.) {
-         hmpt->Fill((1/pt - 1/ptg)/(1/ptg)*100.);
+         hmpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
       } else {
          Float_t phig=TMath::ATan2(gt[ngood].py,gt[ngood].px);
-         Float_t phi =TMath::ATan2(py,px );
          hp->Fill((phi - phig)*1000.);
 
          Float_t lamg=TMath::ATan2(gt[ngood].pz,ptg);
-         Float_t lam =TMath::ATan2(pz ,pt );
          hl->Fill((lam - lamg)*1000.);
 
-         hpt->Fill((1/pt - 1/ptg)/(1/ptg)*100.);
+         hpt->Fill((pt_1 - 1/ptg)/(1/ptg)*100.);
       }
 
-      Float_t mom=track->GetP();
+      Float_t mom=1./(pt_1*TMath::Cos(lam));
       Float_t dedx=track->GetdEdx();
       hep->Fill(mom,dedx,1.);
       if (TMath::Abs(gt[ngood].code)==211)
@@ -279,7 +300,7 @@ Int_t AliTPCComparison() {
 Int_t good_tracks(GoodTrack *gt, Int_t max) {
    Int_t nt=0;
 
-   TFile *file=TFile::Open("galice.root");
+   TFile *file=TFile::Open("rfio:galice.root");
    if (!file->IsOpen()) {cerr<<"Can't open galice.root !\n"; exit(4);}
 
    if (!(gAlice=(AliRun*)file->Get("gAlice"))) {
@@ -287,7 +308,7 @@ Int_t good_tracks(GoodTrack *gt, Int_t max) {
      exit(5);
    }
 
-   gAlice->GetEvent(0);   
+   Int_t np=gAlice->GetEvent(0);   
 
    AliTPC *TPC=(AliTPC*)gAlice->GetDetector("TPC");
    Int_t ver = TPC->IsVersion(); 
@@ -303,11 +324,6 @@ Int_t good_tracks(GoodTrack *gt, Int_t max) {
    Int_t gap=Int_t(0.125*nrows);
    Int_t good_number=Int_t(0.4*nrows);
 
-   TObjArray *particles=gAlice->Particles(); //MI change 
-   
-   //
-   //   Int_t np=particles->GetEntriesFast();
-   Int_t np=gAlice->GetNtrack(); //FCA correction
    Int_t *good=new Int_t[np];
    for (Int_t ii=0; ii<np; ii++) good[ii]=0;
 
@@ -393,7 +409,6 @@ Int_t good_tracks(GoodTrack *gt, Int_t max) {
    }
 
    TTree *TH=gAlice->TreeH();
-   //   TClonesArray *hits=TPC->Hits();
    Int_t npart=(Int_t)TH->GetEntries();
 
    while (npart--) {
@@ -413,20 +428,8 @@ Int_t good_tracks(GoodTrack *gt, Int_t max) {
       else continue;
       AliTPChit *hit1=(AliTPChit*)TPC->NextHit();      
       if (hit1==0) continue;
-      /*      Int_t nhits=hits->GetEntriesFast();
-      if (nhits==0) continue;
-      AliTPChit *hit;
-      Int_t j;
-      for (j=0; j<nhits-1; j++) {
-         hit=(AliTPChit*)hits->UncheckedAt(j);
-         if (hit->fQ==0.) break;
-      }
-      if (j==nhits-1) continue;
-      AliTPChit *hit1=(AliTPChit*)hits->UncheckedAt(j+1);
-      */
       if (hit1->fQ != 0.) continue;
       Int_t i=hit->Track();
-      //      TParticle *p = (TParticle*)particles->UncheckedAt(i);
       TParticle *p = (TParticle*)gAlice->Particle(i);
 
       if (p->GetFirstMother()>=0) continue;  //secondary particle
@@ -449,7 +452,7 @@ Int_t good_tracks(GoodTrack *gt, Int_t max) {
    }
    delete[] good;
 
-   //delete gAlice; gAlice=0;
+   delete gAlice; gAlice=0;
 
    file->Close();
    gBenchmark->Stop("AliTPCComparison");
index 2c8a8aa..dedeba1 100644 (file)
@@ -34,10 +34,11 @@ const Double_t cc[15], Double_t xref, Double_t alpha) : AliKalmanTrack() {
   //-----------------------------------------------------------------
   // This is the main track constructor.
   //-----------------------------------------------------------------
-  fdEdx=0.;
-
-  fAlpha=alpha;
   fX=xref;
+  fAlpha=alpha;
+  if (fAlpha<-TMath::Pi()) fAlpha += 2*TMath::Pi();
+  if (fAlpha>=TMath::Pi()) fAlpha -= 2*TMath::Pi();
+  fdEdx=0.;
 
   fP0=xx[0]; fP1=xx[1]; fP2=xx[2]; fP3=xx[3]; fP4=xx[4];
 
@@ -47,7 +48,8 @@ const Double_t cc[15], Double_t xref, Double_t alpha) : AliKalmanTrack() {
   fC30=cc[6];  fC31=cc[7];  fC32=cc[8];  fC33=cc[9];
   fC40=cc[10]; fC41=cc[11]; fC42=cc[12]; fC43=cc[13]; fC44=cc[14];
 
-  fIndex[fN++]=index;
+  fIndex[0]=index;
+  SetNumberOfClusters(1);
 }
 
 //_____________________________________________________________________________
@@ -55,12 +57,75 @@ AliTPCtrack::AliTPCtrack(const AliTPCtrack& t) : AliKalmanTrack(t) {
   //-----------------------------------------------------------------
   // This is a track copy constructor.
   //-----------------------------------------------------------------
-  for (Int_t i=0; i<fN; i++) fIndex[i]=t.fIndex[i];
-
+  fX=t.fX;
+  fAlpha=t.fAlpha;
   fdEdx=t.fdEdx;
 
-  fAlpha=t.fAlpha;
-  fX=t.fX;
+  fP0=t.fP0; fP1=t.fP1; fP2=t.fP2; fP3=t.fP3; fP4=t.fP4;
+
+  fC00=t.fC00;
+  fC10=t.fC10;  fC11=t.fC11;
+  fC20=t.fC20;  fC21=t.fC21;  fC22=t.fC22;
+  fC30=t.fC30;  fC31=t.fC31;  fC32=t.fC32;  fC33=t.fC33;
+  fC40=t.fC40;  fC41=t.fC41;  fC42=t.fC42;  fC43=t.fC43;  fC44=t.fC44;
+
+  Int_t n=GetNumberOfClusters();
+  for (Int_t i=0; i<n; i++) fIndex[i]=t.fIndex[i];
+}
+
+//_____________________________________________________________________________
+Int_t AliTPCtrack::Compare(const TObject *o) const {
+  //-----------------------------------------------------------------
+  // This function compares tracks according to the their curvature
+  //-----------------------------------------------------------------
+  AliTPCtrack *t=(AliTPCtrack*)o;
+  Double_t co=TMath::Abs(t->Get1Pt());
+  Double_t c =TMath::Abs(Get1Pt());
+  if (c>co) return 1;
+  else if (c<co) return -1;
+  return 0;
+}
+
+//_____________________________________________________________________________
+void AliTPCtrack::GetExternalCovariance(Double_t cc[15]) const {
+  //-------------------------------------------------------------------------
+  // This function returns an external representation of the covriance matrix.
+  //   (See comments in AliTPCtrack.h about external track representation)
+  //-------------------------------------------------------------------------
+  Double_t a=kConversionConstant;
+
+  Double_t c22=fX*fX*fC33-2*fX*fC32+fC22;
+  Double_t c42=fX*fC43-fC42;
+  Double_t c20=fX*fC30-fC20, c21=fX*fC31-fC21, c32=fX*fC33-fC32;
+  
+  cc[0 ]=fC00;
+  cc[1 ]=fC10;   cc[2 ]=fC11;
+  cc[3 ]=c20;    cc[4 ]=c21;    cc[5 ]=c22;
+  cc[6 ]=fC40;   cc[7 ]=fC41;   cc[8 ]=c42;   cc[9 ]=fC44; 
+  cc[10]=fC30*a; cc[11]=fC31*a; cc[12]=c32*a; cc[13]=fC43*a; cc[14]=fC33*a*a;
+
+}
+
+//_____________________________________________________________________________
+Double_t AliTPCtrack::GetPredictedChi2(const AliCluster *c) const 
+{
+  //-----------------------------------------------------------------
+  // This function calculates a predicted chi2 increment.
+  //-----------------------------------------------------------------
+  Double_t r00=c->GetSigmaY2(), r01=0., r11=c->GetSigmaZ2();
+  r00+=fC00; r01+=fC10; r11+=fC11;
+
+  Double_t det=r00*r11 - r01*r01;
+  if (TMath::Abs(det) < 1.e-10) {
+    Int_t n=GetNumberOfClusters();
+    if (n>4) cerr<<n<<" AliKalmanTrack warning: Singular matrix !\n";
+    return 1e10;
+  }
+  Double_t tmp=r00; r00=r11; r11=tmp; r01=-r01;
+  
+  Double_t dy=c->GetY() - fP0, dz=c->GetZ() - fP1;
+  
+  return (dy*r00*dy + 2*r01*dy*dz + dz*r11*dz)/det;
 }
 
 //_____________________________________________________________________________
@@ -70,7 +135,8 @@ Int_t AliTPCtrack::PropagateTo(Double_t xk,Double_t x0,Double_t rho,Double_t pm)
   // This function propagates a track to a reference plane x=xk.
   //-----------------------------------------------------------------
   if (TMath::Abs(fP3*xk - fP2) >= 0.99999) {
-    if (fN>4) cerr<<fN<<" AliTPCtrack warning: Propagation failed !\n";
+    Int_t n=GetNumberOfClusters();
+    if (n>4) cerr<<n<<" AliTPCtrack warning: Propagation failed !\n";
     return 0;
   }
 
@@ -115,9 +181,9 @@ Int_t AliTPCtrack::PropagateTo(Double_t xk,Double_t x0,Double_t rho,Double_t pm)
 
   //Multiple scattering******************
   Double_t d=sqrt((x1-fX)*(x1-fX)+(y1-fP0)*(y1-fP0)+(z1-fP1)*(z1-fP1));
-  Double_t p2=GetP()*GetP();
+  Double_t p2=(1.+ GetTgl()*GetTgl())/(Get1Pt()*Get1Pt());
   Double_t beta2=p2/(p2 + pm*pm);
-  Double_t theta2=14.1*14.1/(beta2*p2*1e6)*d/x0*rho/2.;
+  Double_t theta2=14.1*14.1/(beta2*p2*1e6)*d/x0*rho;
   //Double_t theta2=1.0259e-6*10*10/20/(beta2*p2)*d*rho;
 
   Double_t ey=fP3*fX - fP2, ez=fP4;
@@ -154,8 +220,7 @@ Int_t AliTPCtrack::PropagateToVertex(Double_t x0,Double_t rho,Double_t pm)
 }
 
 //_____________________________________________________________________________
-void AliTPCtrack::Update(const AliCluster *c, Double_t chisq, UInt_t index)
-{
+Int_t AliTPCtrack::Update(const AliCluster *c, Double_t chisq, UInt_t index) {
   //-----------------------------------------------------------------
   // This function associates a cluster with this track.
   //-----------------------------------------------------------------
@@ -173,8 +238,9 @@ void AliTPCtrack::Update(const AliCluster *c, Double_t chisq, UInt_t index)
   Double_t dy=c->GetY() - fP0, dz=c->GetZ() - fP1;
   Double_t cur=fP3 + k30*dy + k31*dz, eta=fP2 + k20*dy + k21*dz;
   if (TMath::Abs(cur*fX-eta) >= 0.99999) {
-    if (fN>4) cerr<<fN<<" AliTPCtrack warning: Filtering failed !\n";
-    return;
+    Int_t n=GetNumberOfClusters();
+    if (n>4) cerr<<n<<" AliTPCtrack warning: Filtering failed !\n";
+    return 0;
   }
 
   fP0 += k00*dy + k01*dz;
@@ -202,8 +268,12 @@ void AliTPCtrack::Update(const AliCluster *c, Double_t chisq, UInt_t index)
 
   fC44-=k40*c04+k41*c14; 
 
-  fIndex[fN++]=index;
-  fChi2 += chisq;
+  Int_t n=GetNumberOfClusters();
+  fIndex[n]=index;
+  SetNumberOfClusters(n+1);
+  SetChi2(GetChi2()+chisq);
+
+  return 1;
 }
 
 //_____________________________________________________________________________
@@ -213,6 +283,8 @@ Int_t AliTPCtrack::Rotate(Double_t alpha)
   // This function rotates this track.
   //-----------------------------------------------------------------
   fAlpha += alpha;
+  if (fAlpha<-TMath::Pi()) fAlpha += 2*TMath::Pi();
+  if (fAlpha>=TMath::Pi()) fAlpha -= 2*TMath::Pi();
   
   Double_t x1=fX, y1=fP0;
   Double_t ca=cos(alpha), sa=sin(alpha);
@@ -224,13 +296,15 @@ Int_t AliTPCtrack::Rotate(Double_t alpha)
   
   Double_t r2=fP3*fX - fP2;
   if (TMath::Abs(r2) >= 0.99999) {
-    if (fN>4) cerr<<fN<<" AliTPCtrack warning: Rotation failed !\n";
+    Int_t n=GetNumberOfClusters();
+    if (n>4) cerr<<n<<" AliTPCtrack warning: Rotation failed !\n";
     return 0;
   }
   
   Double_t y0=fP0 + sqrt(1.- r2*r2)/fP3;
   if ((fP0-y0)*fP3 >= 0.) {
-    if (fN>4) cerr<<fN<<" AliTPCtrack warning: Rotation failed !!!\n";
+    Int_t n=GetNumberOfClusters();
+    if (n>4) cerr<<n<<" AliTPCtrack warning: Rotation failed !!!\n";
     return 0;
   }
 
@@ -268,34 +342,18 @@ Int_t AliTPCtrack::Rotate(Double_t alpha)
 }
 
 //_____________________________________________________________________________
-void AliTPCtrack::GetPxPyPz(Double_t& px, Double_t& py, Double_t& pz) const 
-{
-  //-----------------------------------------------------------------
-  // This function returns reconstructed track momentum in the global system.
-  //-----------------------------------------------------------------
-  Double_t pt=TMath::Abs(GetPt()); // GeV/c
-  Double_t r=fP3*fX-fP2;
-  Double_t y0=fP0 + sqrt(1.- r*r)/fP3;
-  px=-pt*(fP0-y0)*fP3;    //cos(phi);
-  py=-pt*(fP2-fX *fP3);   //sin(phi);
-  pz=pt*fP4;
-  Double_t tmp=px*TMath::Cos(fAlpha) - py*TMath::Sin(fAlpha);
-  py=px*TMath::Sin(fAlpha) + py*TMath::Cos(fAlpha);
-  px=tmp;  
-}
-
-//_____________________________________________________________________________
 void AliTPCtrack::CookLabel(AliTPCClustersArray *ca) {
   //-----------------------------------------------------------------
   // This function cooks the track label. If label<0, this track is fake.
   //-----------------------------------------------------------------
-  Int_t *lb=new Int_t[fN];
-  Int_t *mx=new Int_t[fN];
-  AliTPCcluster **clusters=new AliTPCcluster*[fN];
+  Int_t n=GetNumberOfClusters();
+  Int_t *lb=new Int_t[n];
+  Int_t *mx=new Int_t[n];
+  AliTPCcluster **clusters=new AliTPCcluster*[n];
 
   Int_t i;
   Int_t sec,row,ncl;
-  for (i=0; i<fN; i++) {
+  for (i=0; i<n; i++) {
      lb[i]=mx[i]=0;
      GetCluster(i,sec,row,ncl);
      AliTPCClustersRow *clrow=ca->GetRow(sec,row);
@@ -303,40 +361,46 @@ void AliTPCtrack::CookLabel(AliTPCClustersArray *ca) {
   }
   
   Int_t lab=123456789;
-  for (i=0; i<fN; i++) {
+  for (i=0; i<n; i++) {
     AliTPCcluster *c=clusters[i];
     lab=TMath::Abs(c->GetLabel(0));
     Int_t j;
-    for (j=0; j<fN; j++)
+    for (j=0; j<n; j++)
       if (lb[j]==lab || mx[j]==0) break;
     lb[j]=lab;
     (mx[j])++;
   }
   
   Int_t max=0;
-  for (i=0; i<fN; i++) 
+  for (i=0; i<n; i++) 
     if (mx[i]>max) {max=mx[i]; lab=lb[i];}
     
-  for (i=0; i<fN; i++) {
+  for (i=0; i<n; i++) {
     AliTPCcluster *c=clusters[i];
     if (TMath::Abs(c->GetLabel(1)) == lab ||
         TMath::Abs(c->GetLabel(2)) == lab ) max++;
   }
   
   SetLabel(-lab);
-  if (1.-Float_t(max)/fN <= 0.10) {
+  if (1.-Float_t(max)/n <= 0.10) {
     //Int_t tail=Int_t(0.08*fN);
      Int_t tail=14;
      max=0;
      for (i=1; i<=tail; i++) {
-       AliTPCcluster *c=clusters[fN-i];
+       AliTPCcluster *c=clusters[n-i];
        if (lab == TMath::Abs(c->GetLabel(0)) ||
            lab == TMath::Abs(c->GetLabel(1)) ||
            lab == TMath::Abs(c->GetLabel(2))) max++;
      }
      if (max >= Int_t(0.5*tail)) SetLabel(lab);
   }
-
+/*
+if (lab==111)
+   for (i=0; i<n; i++) {
+     AliTPCcluster *c=clusters[i];
+     cerr<<i<<' '<<c->GetLabel(0)<<endl;
+   }
+*/
   delete[] lb;
   delete[] mx;
   delete[] clusters;
index bf083ad..36a2ea6 100644 (file)
@@ -9,9 +9,25 @@
 //   Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch 
 //-------------------------------------------------------
 
+/*****************************************************************************
+ *                          December 18, 2000                                *
+ *  Internal view of the TPC track parametrisation as well as the order of   *
+ *           track parameters are subject for possible changes !             *
+ *  Use GetExternalParameters() and GetExternalCovariance() to access TPC    *
+ *      track information regardless of its internal representation.         *
+ * This formation is now fixed in the following way:                         *
+ *      external param0:   local Y-coordinate of a track (cm)                *
+ *      external param1:   local Z-coordinate of a track (cm)                *
+ *      external param2:   local sine of the track momentum azimuth angle    *
+ *      external param3:   tangent of the track momentum dip angle           *
+ *      external param4:   1/pt (1/(GeV/c))                                  *
+ *****************************************************************************/
+
 #include <AliKalmanTrack.h>
 #include <TMath.h>
 
+const Double_t kConversionConstant=100/0.299792458/0.2; 
+
 class AliTPCClustersArray;
 class AliTPCcluster;
 
@@ -28,34 +44,63 @@ public:
   void CookLabel(AliTPCClustersArray *carray);
   void SetdEdx(Float_t dedx) {fdEdx=dedx;}
 
-  Float_t GetdEdx() const {return fdEdx;}
-  Double_t GetX()   const {return fX;}
+  Double_t GetX()     const {return fX;}
+  Double_t GetAlpha() const {return fAlpha;}
+  Double_t GetdEdx()  const {return fdEdx;}
+
   Double_t GetY()   const {return fP0;}
   Double_t GetZ()   const {return fP1;}
-  Double_t GetEta() const {return fP2;}
-  Double_t GetC()   const {return fP3;}
+  Double_t GetSnp() const {return fX*fP3 - fP2;}             
+  Double_t Get1Pt() const {return fP3*kConversionConstant;}             
   Double_t GetTgl() const {return fP4;}
+
   Double_t GetSigmaY2() const {return fC00;}
   Double_t GetSigmaZ2() const {return fC11;}
-  Double_t GetSigmaC2() const {return fC33;}
-  Double_t GetAlpha()   const {return fAlpha;}
+
+  Int_t Compare(const TObject *o) const;
+
+  void GetExternalParameters(Double_t& xr, Double_t x[5]) const ;
+  void GetExternalCovariance(Double_t cov[15]) const ;
+
+//******** To be removed next release !!! **************
+  Double_t GetEta() const {return fP2;}
+  Double_t GetC()   const {return fP3;}
+  void GetCovariance(Double_t cc[15]) const {
+    cc[0 ]=fC00;
+    cc[1 ]=fC10;  cc[2 ]=fC11;
+    cc[3 ]=fC20;  cc[4 ]=fC21;  cc[5 ]=fC22;
+    cc[6 ]=fC30;  cc[7 ]=fC31;  cc[8 ]=fC32;  cc[9 ]=fC33;
+    cc[10]=fC40;  cc[11]=fC41;  cc[12]=fC42;  cc[13]=fC43;  cc[14]=fC44;
+  }  
+//****************************************************** 
+
   void GetCluster(Int_t i, Int_t &sec, Int_t &row, Int_t &ncl) const;
 
-  Double_t GetPt() const {return 0.299792458*0.2/GetC()/100;}
-  Double_t GetP()  const {return TMath::Abs(GetPt())*
-    TMath::Sqrt(1.+GetTgl()*GetTgl());}
-  void GetPxPyPz(Double_t &px, Double_t &py, Double_t &pz) const;
+  virtual Double_t GetPredictedChi2(const AliCluster *cluster) const;
   Int_t PropagateTo(Double_t xr,
                     Double_t x0=28.94,Double_t rho=0.9e-3,Double_t pm=0.139);
-  void Update(const AliCluster* c, Double_t chi2, UInt_t i);
+  Int_t Update(const AliCluster* c, Double_t chi2, UInt_t i);
 
-protected: 
-  UInt_t fIndex[200];       // indices of associated clusters 
+private: 
+  Double_t fX;              // X-coordinate of this track (reference plane)
+  Double_t fAlpha;          // Rotation angle the local (TPC sector)
+                            // coordinate system and the global ALICE one.
 
-  Float_t fdEdx;            // dE/dx 
+  Double_t fdEdx;           // dE/dx 
 
-  Double_t fAlpha;          // rotation angle
-  Double_t fX;              // X-coordinate of this track (reference plane)
+  Double_t fP0;             // Y-coordinate of a track
+  Double_t fP1;             // Z-coordinate of a track
+  Double_t fP2;             // C*x0
+  Double_t fP3;             // track curvature
+  Double_t fP4;             // tangent of the track momentum dip angle
+
+  Double_t fC00;                         // covariance
+  Double_t fC10, fC11;                   // matrix
+  Double_t fC20, fC21, fC22;             // of the
+  Double_t fC30, fC31, fC32, fC33;       // track
+  Double_t fC40, fC41, fC42, fC43, fC44; // parameters
+  UInt_t fIndex[200];       // indices of associated clusters 
 
   ClassDef(AliTPCtrack,1)   // Time Projection Chamber reconstructed tracks
 };
@@ -70,6 +115,15 @@ void AliTPCtrack::GetCluster(Int_t i,Int_t &sec,Int_t &row,Int_t &ncl) const {
   ncl=(index&0x0000ffff)>>00;
 }
 
+inline 
+void AliTPCtrack::GetExternalParameters(Double_t& xr, Double_t x[5]) const {
+  //---------------------------------------------------------------------
+  // This function return external TPC track representation
+  //---------------------------------------------------------------------
+     xr=fX;          
+     x[0]=GetY(); x[1]=GetZ(); x[2]=GetSnp(); x[3]=GetTgl(); x[4]=Get1Pt();
+}
+
 #endif
 
 
index 14c9949..5cc4601 100644 (file)
 
 /*
 $Log$
+Revision 1.5  2000/12/20 07:51:59  kowal2
+Changes suggested by Alessandra and Paolo to avoid overlapped
+data fields in encapsulated classes.
+
 Revision 1.4  2000/11/02 07:27:16  kowal2
 code corrections
 
@@ -143,8 +147,8 @@ Int_t s, Int_t rf)
   //-----------------------------------------------------------------
   // This function tries to find a track prolongation.
   //-----------------------------------------------------------------
-  const Int_t kSKIP=(t.GetNumberOfClusters()<10) ? 
-                        10 : Int_t(0.5*sec->GetNRows());
+  const Int_t kSKIP=(t.GetNumberOfClusters()<10) ? 10 : 
+                    Int_t(0.5*sec->GetNRows());
   Int_t tryAgain=kSKIP;
   Double_t alpha=sec->GetAlpha();
   Int_t ns=Int_t(2*TMath::Pi()/alpha+0.5);
@@ -157,7 +161,7 @@ Int_t s, Int_t rf)
     UInt_t index=0;
     Double_t maxchi2=12.;
     const AliTPCRow &krow=sec[s][nr];
-    Double_t sy2=SigmaY2(t.GetX(),t.GetTgl(),t.GetPt());
+    Double_t sy2=SigmaY2(t.GetX(),t.GetTgl(),1./t.Get1Pt());
     Double_t sz2=SigmaZ2(t.GetX(),t.GetTgl());
     Double_t road=4.*sqrt(t.GetSigmaY2() + sy2), y=t.GetY(), z=t.GetZ();
 
@@ -184,8 +188,9 @@ Int_t s, Int_t rf)
     if (cl) {
       Float_t l=sec->GetPadPitchWidth();
       t.SetSampledEdx(cl->GetQ()/l,t.GetNumberOfClusters());
-      t.Update(cl,maxchi2,index);
-      tryAgain=kSKIP;
+      if (!t.Update(cl,maxchi2,index)) {
+         if (!tryAgain--) return 0;
+      } else tryAgain=kSKIP;
     } else {
       if (tryAgain==0) break;
       if (y > ymax) {
@@ -298,7 +303,7 @@ Int_t i1, Int_t i2)
         track->SetSampledEdx(kr1[is]->GetQ()/l,0);
 
         Int_t rc=FindProlongation(*track,sec,ns,i2);
-        if (rc<0 || track->GetNumberOfClusters()<(i1-i2)/2) delete track;
+        if (rc==0 || track->GetNumberOfClusters()<(i1-i2)/2) delete track;
         else seeds.AddLast(track); 
       }
     }
@@ -489,11 +494,12 @@ void AliTPCtracker::AliTPCseed::CookdEdx(Double_t low, Double_t up) {
   // This funtion calculates dE/dX within the "low" and "up" cuts.
   //-----------------------------------------------------------------
   Int_t i;
+  Int_t nc=GetNumberOfClusters();
 
   Int_t swap;//stupid sorting
   do {
     swap=0;
-    for (i=0; i<fN-1; i++) {
+    for (i=0; i<nc-1; i++) {
       if (fdEdxSample[i]<=fdEdxSample[i+1]) continue;
       Float_t tmp=fdEdxSample[i];
       fdEdxSample[i]=fdEdxSample[i+1]; fdEdxSample[i+1]=tmp;
@@ -501,7 +507,7 @@ void AliTPCtracker::AliTPCseed::CookdEdx(Double_t low, Double_t up) {
     }
   } while (swap);
 
-  Int_t nl=Int_t(low*fN), nu=Int_t(up*fN);
+  Int_t nl=Int_t(low*nc), nu=Int_t(up*nc);
   Float_t dedx=0;
   for (i=nl; i<=nu; i++) dedx += fdEdxSample[i];
   dedx /= (nu-nl+1);
@@ -513,8 +519,10 @@ void AliTPCtracker::AliTPCseed::UseClusters(AliTPCClustersArray *ca, Int_t n) {
   //-----------------------------------------------------------------
   // This function marks clusters associated with this track.
   //-----------------------------------------------------------------
+  Int_t nc=GetNumberOfClusters();
   Int_t sec,row,ncl;
-  for (Int_t i=n; i<fN; i++) {
+
+  for (Int_t i=n; i<nc; i++) {
      GetCluster(i,sec,row,ncl);
      AliTPCClustersRow *clrow=ca->GetRow(sec,row);
      AliTPCcluster *c=(AliTPCcluster*)(*clrow)[ncl]; 
index a361733..98d77de 100644 (file)
@@ -107,8 +107,8 @@ public:
                 const Double_t cc[15], Double_t xr, Double_t alpha): 
                 AliTPCtrack(index, xx, cc, xr, alpha) {}
      void SetSampledEdx(Float_t q, Int_t i) {
-        Double_t c=GetC(), e=GetEta(), t=GetTgl(), x=GetX();
-        q *= TMath::Sqrt((1-(c*x-e)*(c*x-e))/(1+t*t));
+        Double_t s=GetSnp(), t=GetTgl();
+        q *= TMath::Sqrt((1-s*s)/(1+t*t));
         fdEdxSample[i]=q;
      }
      void UseClusters(AliTPCClustersArray *ca, Int_t n=0);
index d7e139f..d91c9ef 100644 (file)
@@ -18,6 +18,7 @@
 #pragma link C++ class  AliTPCdigit+;
 #pragma link C++ class  AliTPCcluster+;
 #pragma link C++ class  AliTPCtrack+;
+#pragma link C++ class  AliTPCtracker+;
 
 #pragma link C++ class  AliTPCParam+;
 #pragma link C++ class  AliTPCParamSR-;