]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibAlign.h
AliTriggerAnalisys:
[u/mrichter/AliRoot.git] / TPC / AliTPCcalibAlign.h
1 #ifndef ALITPCCALIBALIGN_H
2 #define ALITPCCALIBALIGN_H
3
4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5  * See cxx source for full Copyright notice                               */
6
7 ////
8 ////
9 ////
10
11 class TFile;
12 class TGraphErrors;
13 class TH1;
14 class THnSparse;
15 #include <TLinearFitter.h>
16 #include <TMatrixDfwd.h>
17 class TObjArray;
18 class TTree;
19
20 #include "AliTPCcalibBase.h"
21 class AliExternalComparison;
22 class AliExternalTrackParam;
23 class AliTPCPointCorrection;
24 class AliTPCseed;
25
26 class AliTPCcalibAlign:public AliTPCcalibBase {
27 public:
28   enum HistoType {kY=0, kZ =1, kPhi=2, kTheta=3, 
29                   kYPhi=4, kZTheta=5, 
30                   kYz=6,kZz=7,kPhiZ=8,kThetaZ=9};
31   enum FitType{ k6=0, k9=1, k12=2};
32   AliTPCcalibAlign();
33   AliTPCcalibAlign(const Text_t *name, const Text_t *title);
34   AliTPCcalibAlign(const AliTPCcalibAlign &align);
35   //
36   virtual ~AliTPCcalibAlign();
37   void     Process(AliESDEvent *event);
38   virtual void ProcessSeed(AliTPCseed *track);
39   virtual void Process(AliTPCseed */*track*/){ return ;}
40   virtual void Analyze();
41   virtual void Terminate();  
42   virtual Long64_t Merge(TCollection* const list);
43   void ExportTrackPoints(AliESDEvent *event);
44   //
45   //
46   void MakeReportDy(TFile *output); 
47   void MakeReportDyPhi(TFile *output);
48   //
49   void UpdatePointCorrection(AliTPCPointCorrection * correction);
50   //
51   virtual void EvalFitters(Int_t minPoints=20);
52   TH1 * GetHisto(HistoType type, Int_t s1, Int_t s2, Bool_t force=kFALSE);
53   void  MakeTree(const char *fname="alignTree.root", Int_t minPoints=20);
54   TGraphErrors * MakeGraph(Int_t sec0, Int_t sec1, Int_t dsec, 
55                            Int_t i0, Int_t i1, FitType type); 
56   Int_t  RefitLinear(const AliTPCseed * seed, Int_t isec, Double_t *fitParam, Int_t refSector, TMatrixD &param, TMatrixD&covar, Double_t xRef, Bool_t both=kFALSE);
57   
58   void ProcessTracklets(const AliExternalTrackParam &t1,
59                         const AliExternalTrackParam &t2,
60                         const AliTPCseed * seed,
61                         Int_t s1,Int_t s2);
62   
63   void UpdateClusterDeltaField(const AliTPCseed * seed);
64   void UpdateAlignSector(const AliTPCseed * seed,Int_t isec); 
65   Int_t GetIndex(Int_t s1,Int_t s2) const {return 72*s1+s2;}
66   //
67   inline const TMatrixD     * GetTransformation(Int_t s1,Int_t s2, Int_t fitType);
68   //
69   inline TLinearFitter* GetFitter12(Int_t s1,Int_t s2);
70   inline TLinearFitter* GetFitter9(Int_t s1,Int_t s2);
71   inline TLinearFitter* GetFitter6(Int_t s1,Int_t s2);
72   //
73   Bool_t GetTransformation12(Int_t s1,Int_t s2,TMatrixD &a);
74   Bool_t GetTransformation9(Int_t s1,Int_t s2,TMatrixD &a);
75   Bool_t GetTransformation6(Int_t s1,Int_t s2,TMatrixD &a);
76   Int_t  AcceptTracklet(const AliExternalTrackParam &tp1,
77                         const AliExternalTrackParam &tp2) const;
78   Int_t  AcceptTracklet(const Double_t *t1,
79                         const Double_t *t2) const;
80
81   void ProcessDiff(const AliExternalTrackParam &t1,
82                    const AliExternalTrackParam &t2,
83                    const AliTPCseed *seed,
84                    Int_t s1,Int_t s2);
85   void ProcessAlign(Double_t * t1, Double_t * t2, Int_t s1,Int_t s2);
86
87 //   Bool_t GetTransformationCovar12(Int_t s1,Int_t s2,TMatrixD &a, Bool_t norm=kFALSE);
88 //   Bool_t GetTransformationCovar9(Int_t s1,Int_t s2,TMatrixD &a, Bool_t norm=kFALSE);
89 //   Bool_t GetTransformationCovar6(Int_t s1,Int_t s2,TMatrixD &a, Bool_t norm=kFALSE);
90   void Add(AliTPCcalibAlign * align);
91   const Int_t *GetPoints() const {return fPoints;}
92   void     Process(AliESDtrack *const track, Int_t runNo=-1){AliTPCcalibBase::Process(track,runNo);};
93   TLinearFitter* GetOrMakeFitter12(Int_t s1,Int_t s2);
94   TLinearFitter* GetOrMakeFitter9(Int_t s1,Int_t s2);
95   TLinearFitter* GetOrMakeFitter6(Int_t s1,Int_t s2);
96   void Process12(const Double_t *t1, const Double_t *t2,
97                  TLinearFitter *fitter) const;
98   void Process9(const Double_t *const t1, const Double_t *const t2, TLinearFitter *fitter) const;
99   void Process6(const Double_t *const t1, const Double_t *const t2, TLinearFitter *fitter) const;
100   void ProcessTree(TTree * tree, AliExternalComparison *comp=0);
101   void GlobalAlign6(Int_t minPoints, Float_t sysError, Int_t niter);
102   //
103   // Cluster comparison Part
104   //
105   //
106   // For visualization and test purposes
107   //
108   Double_t Correct(Int_t type, Int_t value, Int_t s1, Int_t s2, Double_t x, Double_t y, Double_t z, Double_t phi,Double_t theta); 
109   static Double_t SCorrect(Int_t type, Int_t value, Int_t s1, Int_t s2, Double_t x, Double_t y, Double_t z, Double_t phi,Double_t theta){return Instance()->Correct(type,value,s1,s2,x,y,z,phi,theta);}
110   static AliTPCcalibAlign* Instance();
111   void SetInstance(AliTPCcalibAlign* const param){fgInstance = param;}
112   static void Constrain1Pt(AliExternalTrackParam &t1, const AliExternalTrackParam &t2, Bool_t noField);
113   void SetNoField(Bool_t noField){ fNoField=noField;}
114
115   //
116   // Kalman fileter for sectors
117   //
118   void MakeSectorKalman();
119   void UpdateSectorKalman(Int_t sector, Int_t quadrant0, Int_t quadrant1,  TMatrixD *const p0, TMatrixD *const c0, TMatrixD *const p1, TMatrixD *const c1);
120   void UpdateSectorKalman(TMatrixD &par0, TMatrixD &cov0, TMatrixD &para1, TMatrixD &cov1);
121   Double_t GetCorrectionSector(Int_t coord, Int_t sector, Double_t lx, Double_t ly, Double_t lz); 
122   static Double_t SGetCorrectionSector(Int_t coord, Int_t sector, Double_t lx, Double_t ly, Double_t lz); 
123
124   //
125   // Kalman filter for full TPC
126   //
127   void MakeKalman();
128   void UpdateKalman(Int_t sector0, Int_t sector1,  TMatrixD &p0, TMatrixD &c0, TMatrixD &p1, TMatrixD &c1);
129   void UpdateKalman(TMatrixD &par0, TMatrixD &cov0, TMatrixD &para1, TMatrixD &cov1);
130   //
131   //private:
132   static Int_t CheckCovariance(TMatrixD &covar);
133   //
134   //
135   void MakeResidualHistos();
136   void MakeResidualHistosTracklet();
137   THnSparse * GetClusterDelta(Int_t index) const  { return fClusterDelta[index];}
138   THnSparse * GetTrackletDelta(Int_t index) const  { return fTrackletDelta[index];}
139 public:
140   
141   void FillHisto(const Double_t *t1,
142                  const Double_t *t2,
143                  Int_t s1,Int_t s2);
144   void FillHisto(AliExternalTrackParam *tp1,
145                  AliExternalTrackParam *tp2,
146                  Int_t s1,Int_t s2);
147
148 protected:
149   THnSparse *fClusterDelta[2];  //clusters residuals
150   THnSparse *fTrackletDelta[4]; //track residuals
151
152   TObjArray fDphiHistArray;    // array of residual histograms  phi      -kPhi
153   TObjArray fDthetaHistArray;  // array of residual histograms  theta    -kTheta
154   TObjArray fDyHistArray;      // array of residual histograms  y        -kY
155   TObjArray fDzHistArray;      // array of residual histograms  z        -kZ
156   //
157   TObjArray fDyPhiHistArray;      // array of residual histograms  y     -kYPhi
158   TObjArray fDzThetaHistArray;    // array of residual histograms  z-z   -kZTheta
159   //
160   TObjArray fDphiZHistArray;      // array of residual histograms  phi   -kPhiz
161   TObjArray fDthetaZHistArray;    // array of residual histograms  theta -kThetaz
162   TObjArray fDyZHistArray;        // array of residual histograms  y     -kYz
163   TObjArray fDzZHistArray;        // array of residual histograms  z     -kZz
164   //
165   //
166   TObjArray fFitterArray12;    // array of fitters
167   TObjArray fFitterArray9;     // array of fitters
168   TObjArray fFitterArray6;     // array of fitters
169   //
170   TObjArray fMatrixArray12;    // array of transnformtation matrix
171   TObjArray fMatrixArray9;     // array of transnformtation matrix
172   TObjArray fMatrixArray6;     // array of transnformtation matrix 
173   //
174   //
175   //
176   //
177   TObjArray fCombinedMatrixArray6;      // array  combeined transformation matrix
178   //
179   //
180   Int_t fPoints[72*72];        // number of points in the fitter 
181   Bool_t fNoField;             // flag - no field data
182   // refernce x
183   Double_t fXIO;               // OROC-IROC boundary
184   Double_t fXmiddle;           // center of the TPC sector local X
185   Double_t fXquadrant;         // x quadrant
186   //
187   // Kalman filter for sector internal  alignemnt
188   //
189   TObjArray fArraySectorIntParam; // array of sector alignment parameters
190   TObjArray fArraySectorIntCovar; // array of sector alignment covariances 
191   //
192   // Kalman filter for global alignment
193   //
194   TMatrixD  *fSectorParamA;     // Kalman parameter   for A side
195   TMatrixD  *fSectorCovarA;     // Kalman covariance  for A side 
196   TMatrixD  *fSectorParamC;     // Kalman parameter   for A side
197   TMatrixD  *fSectorCovarC;     // Kalman covariance  for A side 
198   //
199   //
200   //
201   Bool_t    fUseInnerOuter;         // flag- use Inner Outer sector for left righ alignment
202   
203   static AliTPCcalibAlign*   fgInstance; //! Instance of this class (singleton implementation)
204 private:
205   AliTPCcalibAlign&  operator=(const AliTPCcalibAlign&);// not implemented
206
207   ClassDef(AliTPCcalibAlign,6)
208 };
209
210
211 TLinearFitter* AliTPCcalibAlign::GetFitter12(Int_t s1,Int_t s2) {
212   return static_cast<TLinearFitter*>(fFitterArray12[GetIndex(s1,s2)]);
213 }
214 TLinearFitter* AliTPCcalibAlign::GetFitter9(Int_t s1,Int_t s2) {
215   return static_cast<TLinearFitter*>(fFitterArray9[GetIndex(s1,s2)]);
216 }
217 TLinearFitter* AliTPCcalibAlign::GetFitter6(Int_t s1,Int_t s2) {
218   return static_cast<TLinearFitter*>(fFitterArray6[GetIndex(s1,s2)]);
219 }
220
221 const TMatrixD * AliTPCcalibAlign::GetTransformation(Int_t s1,Int_t s2, Int_t fitType){
222   if (fitType==0) return static_cast<TMatrixD*>(fMatrixArray6[GetIndex(s1,s2)]);
223   if (fitType==1) return static_cast<TMatrixD*>(fMatrixArray9[GetIndex(s1,s2)]);
224   if (fitType==2) return static_cast<TMatrixD*>(fMatrixArray12[GetIndex(s1,s2)]);
225   return 0;
226 }
227
228
229
230 #endif