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