]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCtrackerParam.h
delete the AliSurvey objs after use
[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 //  adapted to ESD output by: Marcello Lunardon, Padova
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 #include "AliESDEvent.h"
34 //----------------------------
35
36 class AliTPC;
37
38 class AliTPCtrackerParam:
39   public TObject
40 {   
41  public:
42   
43   AliTPCtrackerParam(Int_t coll=0, Double_t Bz=0.4,
44                      const char* evfoldname = AliConfig::GetDefaultEventFolderName());
45   virtual ~AliTPCtrackerParam();
46
47   // this function performs the parameterized tracking
48   //
49   AliTPCtrackerParam(const AliTPCtrackerParam& p);
50   //
51   Int_t Init();
52   Int_t BuildTPCtracks(AliESDEvent* event);
53   // these functions are used to create a DB of cov. matrices,
54   // including regularization, efficiencies and dE/dx
55   void  AllGeantTracks() { fSelAndSmear=kFALSE; return; }
56   void  AnalyzedEdx(const Char_t *outName,Int_t pdg);
57   void  AnalyzePulls(const Char_t *outName);
58   void  AnalyzeResolutions(Int_t pdg);
59   void  CompareTPCtracks(Int_t nEvents=1,
60                          const Char_t *galiceName="galice.root",
61                          const Char_t *trkGeaName="AliTPCtracksGeant.root",
62                          const Char_t *trkKalName="AliTPCtracksSorted.root",
63                          const Char_t *covmatName="CovMatrix.root",
64                          const Char_t *tpceffasciiName="TPCeff.dat",
65                          const Char_t *tpceffrootName="TPCeff.root");
66   void  DrawEffs(const Char_t *inName,Int_t pdg=211);
67   void  DrawPulls(const Char_t *inName,Int_t pdg=211,Int_t par=0);
68   void  MakeDataBase();
69   void  MergeEvents(Int_t evFirst=1,Int_t evLast=1);
70   void  RegularizeCovMatrix(const Char_t *outName,Int_t pdg);
71
72
73   //********* Internal class definition *******
74   class AliTPCtrackParam : public AliTPCtrack {
75   public:
76     AliTPCtrackParam():AliTPCtrack(){}
77     AliTPCtrackParam(const AliTPCtrack &t):AliTPCtrack(t){}
78
79     void AssignMass(Double_t mass) {SetMass(mass); return;}
80     
81   private:
82   };
83   //********* end of internal class ***********
84
85   //********* Internal class definition *******
86   class AliTPCseedGeant : public TObject {
87   public:
88     AliTPCseedGeant():TObject(),
89       fXg(0.),
90       fYg(0.),
91       fZg(0.),
92       fPx(0.),
93       fPy(0.),
94       fPz(0.),
95       fAlpha(0.),
96       fLabel(0),
97       fSector(0){}
98
99     AliTPCseedGeant(Double_t x=0.,Double_t y=0.,Double_t z=0.,
100                     Double_t px=0.,Double_t py=0.,Double_t pz=0.,
101                     Int_t lab=0);
102     Int_t    GetLabel() const { return fLabel; }
103     Double_t GetAlpha() const { return fAlpha; }      
104     Double_t GetXG() const { return fXg; }
105     Double_t GetYG() const { return fYg; }
106     Double_t GetZG() const { return fZg; }
107     Double_t GetXL() const 
108       { return fXg*TMath::Cos(fAlpha)+fYg*TMath::Sin(fAlpha); }
109     Double_t GetYL() const 
110       { return -fXg*TMath::Sin(fAlpha)+fYg*TMath::Cos(fAlpha); }
111     Double_t GetZL() const { return fZg; }
112     Double_t GetPx() const { return fPx; } 
113     Double_t GetPy() const { return fPy; } 
114     Double_t GetPz() const { return fPz; } 
115     Double_t GetPt() const { return TMath::Sqrt(fPx*fPx+fPy*fPy); }
116     Double_t GetEta() const 
117       { return -TMath::Log(TMath::Tan(0.25*TMath::Pi()-0.5*TMath::ATan(fPz/GetPt()))); }
118     void     SetLabel(Int_t lab) { fLabel=lab; return; }
119     Bool_t   InTPCAcceptance() const {
120       if(TMath::Abs(GetZL()+(244.-GetXL())*fPz/GetPt())>252.) return kFALSE;
121       return kTRUE;
122     }
123
124   private:
125     Double_t fXg;     // global x of seed 
126     Double_t fYg;     // global y of seed
127     Double_t fZg;     // global z of seed
128     Double_t fPx;     // global px
129     Double_t fPy;     // global py
130     Double_t fPz;     // global pz
131     Double_t fAlpha;  // alpha angle
132     Int_t    fLabel;  // track label
133     Int_t    fSector; // TPC sector
134   };
135   //******* end of internal class ****************
136   
137  private:
138   AliTPCtrackerParam & operator=(const AliTPCtrackerParam & );
139
140   TString fEvFolderName;//! name of data folder
141
142   Double_t        fBz;          // value of the z component of L3 field (Tesla)
143   Int_t           fColl;        // collision code (0: PbPb6000; 1: pp)
144   Bool_t          fSelAndSmear; // if kFALSE returns GEANT tracks 
145                                 // at TPC first hit 
146   TString         fDBfileName;  // DataBase file name
147   
148   AliTPCtrack     fTrack;    // current track
149
150   TTree          *fCovTree;  // tree with regularized cov matrices 
151                              // for the current track
152   TTree           *fCovTreePi[30];
153   TTree           *fCovTreeKa[30];
154   TTree           *fCovTreePr[30];
155   TTree           *fCovTreeEl[30];
156   TTree           *fCovTreeMu[30];
157
158   AliTPCkineGrid *fDBgrid;   // grid for the cov matrix look-up table  
159   AliTPCkineGrid  fDBgridPi; //               "                  for pions  
160   AliTPCkineGrid  fDBgridKa; //               "                  for kaons
161   AliTPCkineGrid  fDBgridPr; //               "                  for protons
162   AliTPCkineGrid  fDBgridEl; //               "                  for electrons
163   AliTPCkineGrid  fDBgridMu; //               "                  for muons
164
165   AliTPCkineGrid *fEff;   // TPC efficiencies for the current track
166   AliTPCkineGrid  fEffPi; //           "        pions 
167   AliTPCkineGrid  fEffKa; //           "        kaons 
168   AliTPCkineGrid  fEffPr; //           "        protons 
169   AliTPCkineGrid  fEffEl; //           "        electrons 
170   AliTPCkineGrid  fEffMu; //           "        muons 
171
172   AliTPCkineGrid *fPulls;      // pulls for the current track
173   AliTPCkineGrid  fPullsPi[5]; //        "      pions
174   AliTPCkineGrid  fPullsKa[5]; //        "      muons
175   AliTPCkineGrid  fPullsPr[5]; //        "      protons
176   AliTPCkineGrid  fPullsEl[5]; //        "      electrons
177   AliTPCkineGrid  fPullsMu[5]; //        "      muons
178
179   TMatrixD       *fRegPar;     // regularization parameters for the curr. track
180   TMatrixD        fRegParPi;   //                  "        for pions         
181   TMatrixD        fRegParKa;   //                  "        for kaons
182   TMatrixD        fRegParPr;   //                  "        for protons
183   TMatrixD        fRegParEl;   //                  "        for electrons
184   TMatrixD        fRegParMu;   //                  "        for muons
185
186   AliTPCkineGrid *fdEdxMean;   // dEdx mean for the current track
187   AliTPCkineGrid  fdEdxMeanPi; //                "     pions
188   AliTPCkineGrid  fdEdxMeanKa; //                "     kaons
189   AliTPCkineGrid  fdEdxMeanPr; //                "     protons    
190   AliTPCkineGrid  fdEdxMeanEl; //                "     electrons
191   AliTPCkineGrid  fdEdxMeanMu; //                "     muons
192
193   AliTPCkineGrid *fdEdxRMS;    // dEdx RMS for the current track
194   AliTPCkineGrid  fdEdxRMSPi;  //                "     pions
195   AliTPCkineGrid  fdEdxRMSKa;  //                "     kaons
196   AliTPCkineGrid  fdEdxRMSPr;  //                "     protons    
197   AliTPCkineGrid  fdEdxRMSEl;  //                "     electrons
198   AliTPCkineGrid  fdEdxRMSMu;  //                "     muons
199
200
201   void     BuildTrack(AliTPCseedGeant *s,Int_t ch);
202   void     CookdEdx(Double_t pt,Double_t eta);
203   void     CookTrack(Double_t pt,Double_t eta);
204   TMatrixD GetSmearingMatrix(Double_t *cc, Double_t pt,Double_t eta) const;
205   void     InitializeKineGrid(Option_t *which);    
206   void     MakeSeedsFromHits(AliTPC *tpc,TTree *th,TObjArray &seedArray) const;
207   void     MakeSeedsFromRefs(TTree *ttr,
208                              TObjArray &seedArray) const;
209   Int_t    ReadAllData(const Char_t *inName);
210   Int_t    ReadDBgrid(const Char_t *inName);
211   Int_t    ReaddEdx(const Char_t *inName,Int_t pdg);
212   Int_t    ReadEffs(const Char_t *inName);
213   Int_t    ReadPulls(const Char_t *inName);
214   Int_t    ReadRegParams(const Char_t *inName,Int_t pdg);
215   Int_t    ReadCovTrees(const Char_t* inName); 
216   Bool_t   SelectedTrack(Double_t pt, Double_t eta) const;
217   void     SetParticle(Int_t pdg);  
218   void     SmearTrack(Double_t *xx,Double_t *xxsm,TMatrixD cov) const;
219   Int_t    WritedEdx(const Char_t *outName,Int_t pdg);
220   Int_t    WriteEffs(const Char_t *outName);
221   Int_t    WritePulls(const Char_t *outName);
222   Int_t    WriteRegParams(const Char_t *outName,Int_t pdg);
223   
224   
225   ClassDef(AliTPCtrackerParam,1) // TPC tracking parameterization class
226 };
227
228 #endif
229
230
231
232
233
234
235
236
237
238
239
240
241