]> git.uio.no Git - u/mrichter/AliRoot.git/commitdiff
Class for HMPID track inherited from AliKalmanTrack
authordibari <dibari@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 1 Aug 2008 13:31:30 +0000 (13:31 +0000)
committerdibari <dibari@f7af4fe6-9843-0410-8265-dc069ae4e863>
Fri, 1 Aug 2008 13:31:30 +0000 (13:31 +0000)
HMPID/AliHMPIDtrack.cxx [new file with mode: 0644]
HMPID/AliHMPIDtrack.h [new file with mode: 0644]

diff --git a/HMPID/AliHMPIDtrack.cxx b/HMPID/AliHMPIDtrack.cxx
new file mode 100644 (file)
index 0000000..22ccd55
--- /dev/null
@@ -0,0 +1,326 @@
+/**************************************************************************
+ * 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.                  *
+ **************************************************************************/
+
+#include "AliESDtrack.h" 
+#include "AliTracker.h" 
+#include "AliHMPIDtrack.h" 
+
+ClassImp(AliHMPIDtrack)
+
+//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
+AliHMPIDtrack::AliHMPIDtrack():AliKalmanTrack()
+{
+  //
+  // def. ctor
+  //
+}                                
+//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
+AliHMPIDtrack::AliHMPIDtrack(const AliHMPIDtrack& t):AliKalmanTrack(t)
+{
+  //
+  // cctor.
+  //
+ }                                
+//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
+AliHMPIDtrack::AliHMPIDtrack(const AliESDtrack& t):AliKalmanTrack()
+{
+  //
+  // Constructor from AliESDtrack
+  //
+  SetLabel(t.GetLabel());
+  SetChi2(0.);
+  SetMass(t.GetMass());
+
+  Set(t.GetX(),t.GetAlpha(),t.GetParameter(),t.GetCovariance());
+
+  if ((t.GetStatus()&AliESDtrack::kTIME) == 0) return;
+  StartTimeIntegral();
+  Double_t times[10]; t.GetIntegratedTimes(times); SetIntegratedTimes(times);
+  SetIntegratedLength(t.GetIntegratedLength());
+}              
+//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
+AliHMPIDtrack& AliHMPIDtrack::operator=(const AliHMPIDtrack &/*source*/)
+{
+  // ass. op.
+  
+  return *this;
+}
+//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
+Bool_t AliHMPIDtrack::PropagateTo(Double_t xk, Double_t xx0, Double_t xrho)
+{
+  //
+  // Propagates this track to a reference plane defined by "xk" [cm] 
+  // correcting for the mean crossed material.
+  // Arguments:
+  // "xx0"  - thickness/rad.length [units of the radiation length] 
+  // "xrho" - thickness*density    [g/cm^2] 
+  //  Returns: kTRUE if the track propagates to plane, else kFALSE
+  if (xk == GetX()) {
+    return kTRUE;
+  }
+  Double_t bz   = GetBz();
+  if (!AliExternalTrackParam::PropagateTo(xk,bz)) {
+    return kFALSE;
+  }
+  if (!AliExternalTrackParam::CorrectForMeanMaterial(xx0,xrho,GetMass())) { 
+    return kFALSE;
+  }
+  return kTRUE;            
+}//PropagateTo()  
+//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
+Bool_t AliHMPIDtrack::Rotate(Double_t alpha, Bool_t absolute)
+{
+  //
+  // Rotates track parameters in R*phi plane
+  // if absolute rotation alpha is in global system
+  // otherwise alpha rotation is relative to the current rotation angle
+  //  
+
+  if (absolute) {
+    alpha -= GetAlpha();
+  }
+
+  return AliExternalTrackParam::Rotate(GetAlpha()+alpha);
+
+}   
+//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++  
+Int_t AliHMPIDtrack::GetProlongation(Double_t xk, Double_t &y, Double_t &z)
+{
+  //
+  // Find a prolongation at given x
+  // Return 0 if it does not exist
+  //  
+
+  Double_t bz = GetBz();
+
+  if (!AliExternalTrackParam::GetYAt(xk,bz,y)) {
+    return 0;
+  }
+  if (!AliExternalTrackParam::GetZAt(xk,bz,z)) {
+    return 0;
+  }
+
+  return 1;  
+
+}
+//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
+Int_t   AliHMPIDtrack::PropagateToR(Double_t r,Double_t step)
+{
+  //
+  // Propagate track to the radial position
+  // Rotation always connected to the last track position
+  //
+
+  Double_t xyz0[3];
+  Double_t xyz1[3];
+  Double_t y;
+  Double_t z; 
+
+  Double_t radius = TMath::Sqrt(GetX()*GetX() + GetY()*GetY());
+  // Direction +-
+  Double_t dir    = (radius > r) ? -1.0 : 1.0;   
+
+  for (Double_t x = radius+dir*step; dir*x < dir*r; x += dir*step) {
+
+    GetXYZ(xyz0);      
+    Double_t alpha = TMath::ATan2(xyz0[1],xyz0[0]);
+    Rotate(alpha,kTRUE);
+    GetXYZ(xyz0);      
+    GetProlongation(x,y,z);
+    xyz1[0] = x * TMath::Cos(alpha) + y * TMath::Sin(alpha); 
+    xyz1[1] = x * TMath::Sin(alpha) - y * TMath::Cos(alpha);
+    xyz1[2] = z;
+    Double_t param[7];
+    AliTracker::MeanMaterialBudget(xyz0,xyz1,param);
+    if (param[1] <= 0) {
+      param[1] = 100000000;
+    }
+    PropagateTo(x,param[1],param[0]*param[4]);
+
+  } 
+
+  GetXYZ(xyz0);        
+  Double_t alpha = TMath::ATan2(xyz0[1],xyz0[0]);
+  Rotate(alpha,kTRUE);
+  GetXYZ(xyz0);        
+  GetProlongation(r,y,z);
+  xyz1[0] = r * TMath::Cos(alpha) + y * TMath::Sin(alpha); 
+  xyz1[1] = r * TMath::Sin(alpha) - y * TMath::Cos(alpha);
+  xyz1[2] = z;
+  Double_t param[7];
+  AliTracker::MeanMaterialBudget(xyz0,xyz1,param);
+
+  if (param[1] <= 0) {
+    param[1] = 100000000;
+  }
+  PropagateTo(r,param[1],param[0]*param[4]);
+
+  return 0;
+
+}
+//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
+Double_t AliHMPIDtrack::GetPredictedChi2(const AliCluster3D *c) const {
+  //
+  // Arguments: AliCluster3D
+  // Returns:   Chi2 of track for the cluster
+  Double_t      p[3]={c->GetX(),       c->GetY(),       c->GetZ()};
+  Double_t  covyz[3]={c->GetSigmaY2(), c->GetSigmaYZ(), c->GetSigmaZ2()};
+  Double_t covxyz[3]={c->GetSigmaX2(), c->GetSigmaXY(), c->GetSigmaXZ()};
+  return AliExternalTrackParam::GetPredictedChi2(p, covyz, covxyz);
+}//GetPredictedChi2()
+//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
+Bool_t AliHMPIDtrack::PropagateTo(const AliCluster3D *c) {
+  //
+  // Arguments: AliCluster3D
+  // Returns: kTRUE if the track propagates to the plane of the cluster
+  Double_t      oldX=GetX(),   oldY=GetY(), oldZ=GetZ();
+  Double_t      p[3]={c->GetX(), c->GetY(), c->GetZ()};
+  Double_t  covyz[3]={c->GetSigmaY2(), c->GetSigmaYZ(), c->GetSigmaZ2()};
+  Double_t covxyz[3]={c->GetSigmaX2(), c->GetSigmaXY(), c->GetSigmaXZ()};
+  Double_t bz=GetBz();
+  
+  if(!AliExternalTrackParam::PropagateTo(p, covyz, covxyz, bz)) return kFALSE;
+  if(IsStartedTimeIntegral()) 
+    {
+      Double_t d = TMath::Sqrt((GetX()-oldX)*(GetX()-oldX) + (GetY()-oldY)*(GetY()-oldY) + (GetZ()-oldZ)*(GetZ()-oldZ));
+      if (GetX()<oldX) d=-d;
+      AddTimeStep(d);
+  }
+  return kTRUE;
+}//PropagateTo()
+//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
+Double_t AliHMPIDtrack::GetBz() const {
+  // 
+  // Arguments: none
+  // Returns: Bz component of the magnetic field (kG)
+  //
+  if (AliTracker::UniformField()) return AliTracker::GetBz();
+  Double_t r[3]; GetXYZ(r);
+  return AliTracker::GetBz(r);
+}//GetBz()
+//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
+Bool_t AliHMPIDtrack::Intersect(Double_t pnt[3], Double_t norm[3], Double_t bz) const {
+  //+++++++++++++++++++++++++++++++++++++++++    
+  // Origin: K. Shileev (Kirill.Shileev@cern.ch)
+  // Finds point of intersection (if exists) of the helix with the plane. 
+  // Stores result in fX and fP.   
+  // Arguments: planePoint,planeNorm - the plane defined by any plane's point 
+  // and vector, normal to the plane
+  // Returns: kTrue if helix intersects the plane, kFALSE otherwise.
+  //+++++++++++++++++++++++++++++++++++++++++    
+  Double_t x0[3]; GetXYZ(x0); //get track position in MARS
+  //Printf("xxx Intersect OLD: bz: %lf xxx",bz);
+  
+  //estimates initial helix length up to plane
+  Double_t s=(pnt[0]-x0[0])*norm[0] + (pnt[1]-x0[1])*norm[1] + (pnt[2]-x0[2])*norm[2];
+  Double_t dist=99999,distPrev=dist;
+  Double_t x[3],p[3]; 
+  while(TMath::Abs(dist)>0.00001){
+    //calculates helix at the distance s from x0 ALONG the helix
+    Propagate(s,x,p,bz);
+    //distance between current helix position and plane
+    dist=(x[0]-pnt[0])*norm[0]+(x[1]-pnt[1])*norm[1]+(x[2]-pnt[2])*norm[2];  
+    if(TMath::Abs(dist) >= TMath::Abs(distPrev)) {return kFALSE;}
+    distPrev=dist;
+    s-=dist;
+  }
+  //on exit pnt is intersection point,norm is track vector at that point, 
+  //all in MARS
+  for (Int_t i=0; i<3; i++) {pnt[i]=x[i]; norm[i]=p[i];}
+  return kTRUE;
+}//Intersect()
+//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
+Bool_t AliHMPIDtrack::Intersect(AliHMPIDtrack *pTrk,Double_t pnt[3], Double_t norm[3]) {
+  // Finds point of intersection (if exists) of the helix with the plane. 
+  // Stores result in fX and fP.   
+  // Arguments: planePoint,planeNorm - the plane defined by any plane's point 
+  // and vector, normal to the plane
+  // Returns: kTrue if helix intersects the plane, kFALSE otherwise.
+  
+  //Printf(":::::: 111 :::::: pnt0: %lf pnt1: %lf pnt2: %lf norm0: %lf norm1: %lf norm2: %lf",pnt[0],pnt[1],pnt[2],norm[0],norm[1],norm[2]);
+  
+  Double_t bz=-1.0*GetBz();  
+  Double_t x0[3]; pTrk->GetXYZ(x0); //get track position in MARS
+  //Printf(":::::: 222 :::::: x0: %lf x1: %lf x2: %lf",x0[0],x0[1],x0[2]);
+
+  Double_t rad;
+  //estimates initial helix length up to plane
+  Double_t s=(pnt[0]-x0[0])*norm[0] + (pnt[1]-x0[1])*norm[1] + (pnt[2]-x0[2])*norm[2];
+  Double_t dist=99999,distPrev=dist;
+  Double_t x[3],p[3]; 
+  while(TMath::Abs(dist)>0.00001){
+    //calculates helix at the distance s from x0 ALONG the helix
+    Propagate(s,x,p,bz);
+    
+    //distance between current helix position and plane
+    dist=(x[0]-pnt[0])*norm[0]+(x[1]-pnt[1])*norm[1]+(x[2]-pnt[2])*norm[2];
+    if(TMath::Abs(dist) >= TMath::Abs(distPrev)) {return kFALSE;}
+    distPrev=dist;
+    s-=dist;
+  }
+   // Printf(":::::: 333 :::::: pnt0: %lf pnt1: %lf pnt2: %lf norm0: %lf norm1: %lf norm2: %lf",pnt[0],pnt[1],pnt[2],norm[0],norm[1],norm[2]);
+  //on exit pnt is intersection point,norm is track vector at that point, 
+  //all in MARS
+  rad=TMath::Sqrt(x[0]*x[0]+x[1]*x[1]); 
+  PropagateTo(rad);       //propagate with the average dens. rad. lenght
+  GetXYZAt(rad,bz,x); 
+  GetPxPyPz(p);                                                       //propagate to the radial distance of the PC or the Radiator    
+  for (Int_t i=0; i<3; i++) {pnt[i]=x[i]; norm[i]=p[i];}
+ // Printf(":::::: 444 :::::: pnt0: %lf pnt1: %lf pnt2: %lf norm0: %lf norm1: %lf norm2: %lf",pnt[0],pnt[1],pnt[2],norm[0],norm[1],norm[2]);
+ /*
+   pEsdTrk->SetHMPIDpx(p[0]);
+  pEsdTrk->SetHMPIDpy(p[1]);
+  pEsdTrk->SetHMPIDpz(p[2]);
+  */
+  
+  return kTRUE;
+}//Intersect()
+//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
+void AliHMPIDtrack::Propagate(Double_t len, Double_t x[3],Double_t p[3], Double_t bz) const {
+  //+++++++++++++++++++++++++++++++++++++++++    
+  // Origin: K. Shileev (Kirill.Shileev@cern.ch)
+  // Extrapolate track along simple helix in magnetic field
+  // Arguments: len -distance alogn helix, [cm]
+  //            bz  - mag field, [kGaus]   
+  // Returns: x and p contain extrapolated positon and momentum  
+  // The momentum returned for straight-line tracks is meaningless !
+  //+++++++++++++++++++++++++++++++++++++++++    
+  GetXYZ(x);    
+  if (OneOverPt() < kAlmost0 || TMath::Abs(bz) < kAlmost0Field ){ //straight-line tracks
+     Double_t unit[3]; GetDirection(unit);
+     x[0]+=unit[0]*len;   
+     x[1]+=unit[1]*len;   
+     x[2]+=unit[2]*len;
+
+     p[0]=unit[0]/kAlmost0;   
+     p[1]=unit[1]/kAlmost0;   
+     p[2]=unit[2]/kAlmost0;   
+  } else {
+     GetPxPyPz(p);
+     Double_t pp=GetP();
+     Double_t a = -kB2C*bz*GetSign(); ////////// what is kB2C
+     Double_t rho = a/pp;
+     x[0] += p[0]*TMath::Sin(rho*len)/a - p[1]*(1-TMath::Cos(rho*len))/a;
+     x[1] += p[1]*TMath::Sin(rho*len)/a + p[0]*(1-TMath::Cos(rho*len))/a;
+     x[2] += p[2]*len/pp;
+     Double_t p0=p[0];
+     p[0] = p0  *TMath::Cos(rho*len) - p[1]*TMath::Sin(rho*len);
+     p[1] = p[1]*TMath::Cos(rho*len) + p0  *TMath::Sin(rho*len);
+  }
+}//Propagate()
+//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
diff --git a/HMPID/AliHMPIDtrack.h b/HMPID/AliHMPIDtrack.h
new file mode 100644 (file)
index 0000000..e689fa8
--- /dev/null
@@ -0,0 +1,50 @@
+#ifndef ALIHMPIDTRACK_H
+#define ALIHMPIDTRACK_H
+
+/* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
+ * See cxx source for full Copyright notice                               */
+
+//+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
+// Class HMPID track matching in the common tracking framework
+//+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
+
+#include "TMath.h"
+#include "TVector2.h"
+
+#include "AliKalmanTrack.h"
+
+#include "AliCluster3D.h"
+#include "AliHMPIDCluster.h"
+
+class TObject;
+class AliESDtrack;
+
+class AliHMPIDtrack : public AliKalmanTrack {
+
+public:
+   AliHMPIDtrack();
+   AliHMPIDtrack(const AliHMPIDtrack& t);
+   AliHMPIDtrack(const AliESDtrack& t);
+   AliHMPIDtrack& operator=(const AliHMPIDtrack &/*source*/); // ass. op.
+   Double_t GetPredictedChi2(const AliCluster3D *c) const;                                
+   Bool_t   PropagateTo(const AliCluster3D *c);
+   Bool_t   PropagateTo(Double_t xr, Double_t x0 = 8.72, Double_t rho = 5.86e-3);                 //Use material definition as for TOF???
+   void     Propagate(Double_t len,Double_t x[3],Double_t p[3],Double_t bz) const;                //HMPID method moved from AliExternalTrackParam
+   Int_t    PropagateToR(Double_t r,Double_t step);
+   Bool_t   Rotate(Double_t alpha, Bool_t absolute);
+   Int_t    GetProlongation(Double_t xk, Double_t &y, Double_t &z);
+   Bool_t   Intersect(Double_t pnt[3], Double_t norm[3], Double_t bz) const;                      //HMPID method moved from AliExternalTrackParam
+   Bool_t   Intersect(AliHMPIDtrack *pTrk,Double_t pnt[3], Double_t norm[3]) ;                      //just for test 
+   Double_t GetBz() const;     
+
+              
+protected:
+   Bool_t   Update(const AliCluster */*c*/, Double_t /*chi2*/, Int_t /*idx*/) {return 0;}
+   Double_t GetPredictedChi2(const AliCluster */*c*/) const {return 0.;}
+   
+ private:
+   ClassDef(AliHMPIDtrack,0) // HMPID reconstructed tracks
+
+};                     
+
+#endif