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