]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCtrackerParam.h
3e527251c1f0d7976a497cf20739ff528a0af6cf
[u/mrichter/AliRoot.git] / TPC / AliTPCtrackerParam.h
1 #ifndef ALITPCTRACKERPARAM_H
2 #define ALITPCTRACKERPARAM_H
3 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. */
4 /* See cxx source for full Copyright notice                               */
5
6 /* $Id$ */
7
8
9 /////////////////////////////////////////////////////////////////////////
10 //                                                                        
11 // This class builds AliTPCtrack objects from generated tracks to feed    
12 // ITS tracking (V2). The AliTPCtrack is built from its first hit in      
13 // the TPC. The track is assigned a Kalman-like covariance matrix         
14 // depending on its pT and pseudorapidity and track parameters are        
15 // smeared according to this covariance matrix.                           
16 // Output file contains sorted tracks, ready for matching with ITS.        
17 //                                                                        
18 // See implementation file for more details.  
19 //                                  
20 //                                                                        
21 //  Origin: Andrea Dainese, Padova - e-mail: andrea.dainese@pd.infn.it     
22 //                                                                        
23 /////////////////////////////////////////////////////////////////////////
24
25
26 //----- Root headers ---------
27 class TTree;
28 #include <TMatrixD.h>
29 //---- AliRoot headers -------
30 #include "AliConfig.h"
31 #include "AliTPCkineGrid.h"
32 #include "AliTPCtrack.h"
33 //----------------------------
34
35 class AliTPC;
36
37 class AliTPCtrackerParam:
38   public TObject
39 {
40    
41  public:
42   AliTPCtrackerParam(Int_t coll=0, Double_t Bz=0.4, Int_t n=1,
43                      const char* evfoldname = AliConfig::GetDefaultEventFolderName());
44   virtual ~AliTPCtrackerParam();
45
46   // this function performs the parameterized tracking
47   //
48   AliTPCtrackerParam(const AliTPCtrackerParam& p);
49   //
50   Int_t BuildTPCtracks(const TFile *inp, TFile *out);
51   // these functions are used to create a DB of cov. matrices,
52   // including regularization, efficiencies and dE/dx
53   void  AllGeantTracks() { fSelAndSmear=kFALSE; return; }
54   void  AnalyzedEdx(const Char_t *outName,Int_t pdg);
55   void  AnalyzePulls(const Char_t *outName);
56   void  AnalyzeResolutions(Int_t pdg);
57   void  CompareTPCtracks(const Char_t *galiceName="galice.root",
58                          const Char_t *trkGeaName="AliTPCtracksGeant.root",
59                          const Char_t *trkKalName="AliTPCtracksSorted.root",
60                          const Char_t *covmatName="CovMatrix.root",
61                          const Char_t *tpceffasciiName="TPCeff.dat",
62                          const Char_t *tpceffrootName="TPCeff.root");
63   void  DrawEffs(const Char_t *inName,Int_t pdg=211);
64   void  DrawPulls(const Char_t *inName,Int_t pdg=211,Int_t par=0);
65   void  MakeDataBase();
66   void  MergeEvents(Int_t evFirst=1,Int_t evLast=1);
67   void  RegularizeCovMatrix(const Char_t *outName,Int_t pdg);
68
69
70   //********* Internal class definition *******
71   class AliTPCtrackParam : public AliTPCtrack {
72   public:
73     AliTPCtrackParam():AliTPCtrack(){}
74     AliTPCtrackParam(const AliTPCtrack &t):AliTPCtrack(t){}
75
76     void AssignMass(Double_t mass) {SetMass(mass); return;}
77     
78   private:
79   };
80   //********* end of internal class ***********
81
82   //********* Internal class definition *******
83   class AliTPCseedGeant : public TObject {
84   public:
85     AliTPCseedGeant(Double_t x=0.,Double_t y=0.,Double_t z=0.,
86                     Double_t px=0.,Double_t py=0.,Double_t pz=0.,
87                     Int_t lab=0);
88     Int_t    GetLabel() const { return fLabel; }
89     Double_t GetAlpha() const { return fAlpha; }      
90     Double_t GetXL() const 
91       { return fXg*TMath::Cos(fAlpha)+fYg*TMath::Sin(fAlpha); }
92     Double_t GetYL() const 
93       { return -fXg*TMath::Sin(fAlpha)+fYg*TMath::Cos(fAlpha); }
94     Double_t GetZL() const { return fZg; }
95     Double_t GetPx() const { return fPx; } 
96     Double_t GetPy() const { return fPy; } 
97     Double_t GetPz() const { return fPz; } 
98     Double_t GetPt() const { return TMath::Sqrt(fPx*fPx+fPy*fPy); }
99     Double_t GetEta() const 
100       { return -TMath::Log(TMath::Tan(0.25*TMath::Pi()-0.5*TMath::ATan(fPz/GetPt()))); }
101     void     SetLabel(Int_t lab) { fLabel=lab; return; }
102     Bool_t   InTPCAcceptance() const {
103       if(TMath::Abs(GetZL()+(244.-GetXL())*fPz/GetPt())>252.) return kFALSE;
104       return kTRUE;
105     }
106
107   private:
108     Double_t fXg;     // global x of seed 
109     Double_t fYg;     // global y of seed
110     Double_t fZg;     // global z of seed
111     Double_t fPx;     // global px
112     Double_t fPy;     // global py
113     Double_t fPz;     // global pz
114     Double_t fAlpha;  // alpha angle
115     Int_t    fLabel;  // track label
116     Int_t    fSector; // TPC sector
117   };
118   //******* end of internal class ****************
119   
120  private:
121   AliTPCtrackerParam & operator=(const AliTPCtrackerParam & );
122
123   TString fEvFolderName;//! name of data folder
124
125   Int_t           fNevents;     // number of events in the file to be processed
126   Double_t        fBz;          // value of the z component of L3 field (Tesla)
127   Int_t           fColl;        // collision code (0: PbPb6000; 1: pp)
128   Bool_t          fSelAndSmear; // if kFALSE returns GEANT tracks 
129                                 // at TPC first hit 
130   TString         fDBfileName;  // DataBase file name
131   
132   AliTPCtrack     fTrack;    // current track
133
134   TTree          *fCovTree;  // tree with regularized cov matrices 
135                              // for the current track
136
137   AliTPCkineGrid *fDBgrid;   // grid for the cov matrix look-up table  
138   AliTPCkineGrid  fDBgridPi; //               "                  for pions  
139   AliTPCkineGrid  fDBgridKa; //               "                  for kaons
140   AliTPCkineGrid  fDBgridPr; //               "                  for protons
141   AliTPCkineGrid  fDBgridEl; //               "                  for electrons
142   AliTPCkineGrid  fDBgridMu; //               "                  for muons
143
144   AliTPCkineGrid *fEff;   // TPC efficiencies for the current track
145   AliTPCkineGrid  fEffPi; //           "        pions 
146   AliTPCkineGrid  fEffKa; //           "        kaons 
147   AliTPCkineGrid  fEffPr; //           "        protons 
148   AliTPCkineGrid  fEffEl; //           "        electrons 
149   AliTPCkineGrid  fEffMu; //           "        muons 
150
151   AliTPCkineGrid *fPulls;      // pulls for the current track
152   AliTPCkineGrid  fPullsPi[5]; //        "      pions
153   AliTPCkineGrid  fPullsKa[5]; //        "      muons
154   AliTPCkineGrid  fPullsPr[5]; //        "      protons
155   AliTPCkineGrid  fPullsEl[5]; //        "      electrons
156   AliTPCkineGrid  fPullsMu[5]; //        "      muons
157
158   TMatrixD       *fRegPar;     // regularization parameters for the curr. track
159   TMatrixD        fRegParPi;   //                  "        for pions         
160   TMatrixD        fRegParKa;   //                  "        for kaons
161   TMatrixD        fRegParPr;   //                  "        for protons
162   TMatrixD        fRegParEl;   //                  "        for electrons
163   TMatrixD        fRegParMu;   //                  "        for muons
164
165   AliTPCkineGrid *fdEdxMean;   // dEdx mean for the current track
166   AliTPCkineGrid  fdEdxMeanPi; //                "     pions
167   AliTPCkineGrid  fdEdxMeanKa; //                "     kaons
168   AliTPCkineGrid  fdEdxMeanPr; //                "     protons    
169   AliTPCkineGrid  fdEdxMeanEl; //                "     electrons
170   AliTPCkineGrid  fdEdxMeanMu; //                "     muons
171
172   AliTPCkineGrid *fdEdxRMS;    // dEdx RMS for the current track
173   AliTPCkineGrid  fdEdxRMSPi;  //                "     pions
174   AliTPCkineGrid  fdEdxRMSKa;  //                "     kaons
175   AliTPCkineGrid  fdEdxRMSPr;  //                "     protons    
176   AliTPCkineGrid  fdEdxRMSEl;  //                "     electrons
177   AliTPCkineGrid  fdEdxRMSMu;  //                "     muons
178
179
180   void     BuildTrack(AliTPCseedGeant *s,Int_t ch);
181   Int_t    CheckLabel(AliTPCseedGeant *s,Int_t nPart,
182                       Double_t *ptkine,Double_t *pzkine) const;
183   void     CookdEdx(Double_t pt,Double_t eta);
184   void     CookTrack(Double_t pt,Double_t eta);
185   TMatrixD GetSmearingMatrix(Double_t *cc, Double_t pt,Double_t eta) const;
186   void     InitializeKineGrid(Option_t *which);    
187   void     MakeSeedsFromHits(AliTPC *tpc,TTree *th,TObjArray &seedArray) const;
188   void     MakeSeedsFromRefs(TTree *ttr,
189                              TObjArray &seedArray) const;
190   Int_t    ReadAllData(const Char_t *inName);
191   Int_t    ReadDBgrid(const Char_t *inName);
192   Int_t    ReaddEdx(const Char_t *inName,Int_t pdg);
193   Int_t    ReadEffs(const Char_t *inName);
194   Int_t    ReadPulls(const Char_t *inName);
195   Int_t    ReadRegParams(const Char_t *inName,Int_t pdg);
196   Bool_t   SelectedTrack(Double_t pt, Double_t eta) const;
197   void     SetParticle(Int_t pdg);  
198   void     SmearTrack(Double_t *xx,Double_t *xxsm,TMatrixD cov) const;
199   Int_t    WritedEdx(const Char_t *outName,Int_t pdg);
200   Int_t    WriteEffs(const Char_t *outName);
201   Int_t    WritePulls(const Char_t *outName);
202   Int_t    WriteRegParams(const Char_t *outName,Int_t pdg);
203   
204   
205   ClassDef(AliTPCtrackerParam,1) // TPC tracking parameterization class
206 };
207
208 #endif
209
210
211
212
213
214
215
216
217
218
219
220
221