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