1 #ifndef ALITPCCALIBALIGN_H
2 #define ALITPCCALIBALIGN_H
4 /* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5 * See cxx source for full Copyright notice */
12 #include "TObjArray.h"
13 #include "TLinearFitter.h"
14 #include "AliTPCcalibBase.h"
17 class AliExternalTrackParam;
18 class AliExternalComparison;
23 class AliTPCPointCorrection;
26 class AliTPCcalibAlign:public AliTPCcalibBase {
28 enum HistoType {kY=0, kZ =1, kPhi=2, kTheta=3,
30 kYz=6,kZz=7,kPhiZ=8,kThetaZ=9};
31 enum FitType{ k6=0, k9=1, k12=2};
33 AliTPCcalibAlign(const Text_t *name, const Text_t *title);
34 AliTPCcalibAlign(const AliTPCcalibAlign &align);
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 void ExportTrackPoints(AliESDEvent *event);
45 void MakeReportDy(TFile *output);
46 void MakeReportDyPhi(TFile *output);
48 void UpdatePointCorrection(AliTPCPointCorrection * correction);
50 virtual void EvalFitters(Int_t minPoints=20);
51 TH1 * GetHisto(HistoType type, Int_t s1, Int_t s2, Bool_t force=kFALSE);
52 void MakeTree(const char *fname="alignTree.root", Int_t minPoints=20);
53 TGraphErrors * MakeGraph(Int_t sec0, Int_t sec1, Int_t dsec,
54 Int_t i0, Int_t i1, FitType type);
55 Int_t RefitLinear(const AliTPCseed * seed, Int_t isec, Double_t *fitParam, Int_t refSector, TMatrixD ¶m, TMatrixD&covar, Double_t xRef, Bool_t both=kFALSE);
57 void ProcessTracklets(const AliExternalTrackParam &t1,
58 const AliExternalTrackParam &t2,
59 const AliTPCseed * seed,
62 void UpdateAlignSector(const AliTPCseed * seed,Int_t isec);
63 inline Int_t GetIndex(Int_t s1,Int_t s2){return 72*s1+s2;}
65 inline const TMatrixD * GetTransformation(Int_t s1,Int_t s2, Int_t fitType);
67 inline TLinearFitter* GetFitter12(Int_t s1,Int_t s2);
68 inline TLinearFitter* GetFitter9(Int_t s1,Int_t s2);
69 inline TLinearFitter* GetFitter6(Int_t s1,Int_t s2);
71 Bool_t GetTransformation12(Int_t s1,Int_t s2,TMatrixD &a);
72 Bool_t GetTransformation9(Int_t s1,Int_t s2,TMatrixD &a);
73 Bool_t GetTransformation6(Int_t s1,Int_t s2,TMatrixD &a);
74 Int_t AcceptTracklet(const AliExternalTrackParam &tp1,
75 const AliExternalTrackParam &tp2);
76 Int_t AcceptTracklet(const Double_t *t1,
79 void ProcessDiff(const AliExternalTrackParam &t1,
80 const AliExternalTrackParam &t2,
81 const AliTPCseed *seed,
83 void ProcessAlign(Double_t * t1, Double_t * t2, Int_t s1,Int_t s2);
85 // Bool_t GetTransformationCovar12(Int_t s1,Int_t s2,TMatrixD &a, Bool_t norm=kFALSE);
86 // Bool_t GetTransformationCovar9(Int_t s1,Int_t s2,TMatrixD &a, Bool_t norm=kFALSE);
87 // Bool_t GetTransformationCovar6(Int_t s1,Int_t s2,TMatrixD &a, Bool_t norm=kFALSE);
88 void Add(AliTPCcalibAlign * align);
89 Int_t *GetPoints() {return fPoints;}
90 void Process(AliESDtrack *track, Int_t runNo=-1){AliTPCcalibBase::Process(track,runNo);};
91 TLinearFitter* GetOrMakeFitter12(Int_t s1,Int_t s2);
92 TLinearFitter* GetOrMakeFitter9(Int_t s1,Int_t s2);
93 TLinearFitter* GetOrMakeFitter6(Int_t s1,Int_t s2);
94 void Process12(const Double_t *t1, const Double_t *t2,
95 TLinearFitter *fitter);
96 void Process9(Double_t *t1, Double_t *t2, TLinearFitter *fitter);
97 void Process6(Double_t *t1, Double_t *t2, TLinearFitter *fitter);
98 void ProcessTree(TTree * tree, AliExternalComparison *comp=0);
99 void GlobalAlign6(Int_t minPoints, Float_t sysError, Int_t niter);
101 // Cluster comparison Part
104 // For visualization and test purposes
106 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);
107 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);}
108 static AliTPCcalibAlign* Instance();
109 void SetInstance(AliTPCcalibAlign*param){fgInstance = param;}
110 static void Constrain1Pt(AliExternalTrackParam &t1, const AliExternalTrackParam &t2, Bool_t noField);
111 void SetNoField(Bool_t noField){ fNoField=noField;}
114 // Kalman fileter for sectors
116 void MakeSectorKalman();
117 void UpdateSectorKalman(Int_t sector, Int_t quadrant0, Int_t quadrant1, TMatrixD *p0, TMatrixD *c0, TMatrixD *p1, TMatrixD *c1);
118 void UpdateSectorKalman(TMatrixD &par0, TMatrixD &cov0, TMatrixD ¶1, TMatrixD &cov1);
119 Double_t GetCorrectionSector(Int_t coord, Int_t sector, Double_t lx, Double_t ly, Double_t lz);
120 static Double_t SGetCorrectionSector(Int_t coord, Int_t sector, Double_t lx, Double_t ly, Double_t lz);
123 // Kalman filter for full TPC
126 void UpdateKalman(Int_t sector0, Int_t sector1, TMatrixD &p0, TMatrixD &c0, TMatrixD &p1, TMatrixD &c1);
127 void UpdateKalman(TMatrixD &par0, TMatrixD &cov0, TMatrixD ¶1, TMatrixD &cov1);
130 static Int_t CheckCovariance(TMatrixD &covar);
133 void FillHisto(const Double_t *t1,
137 TObjArray fDphiHistArray; // array of residual histograms phi -kPhi
138 TObjArray fDthetaHistArray; // array of residual histograms theta -kTheta
139 TObjArray fDyHistArray; // array of residual histograms y -kY
140 TObjArray fDzHistArray; // array of residual histograms z -kZ
142 TObjArray fDyPhiHistArray; // array of residual histograms y -kYPhi
143 TObjArray fDzThetaHistArray; // array of residual histograms z-z -kZTheta
145 TObjArray fDphiZHistArray; // array of residual histograms phi -kPhiz
146 TObjArray fDthetaZHistArray; // array of residual histograms theta -kThetaz
147 TObjArray fDyZHistArray; // array of residual histograms y -kYz
148 TObjArray fDzZHistArray; // array of residual histograms z -kZz
151 TObjArray fFitterArray12; // array of fitters
152 TObjArray fFitterArray9; // array of fitters
153 TObjArray fFitterArray6; // array of fitters
155 TObjArray fMatrixArray12; // array of transnformtation matrix
156 TObjArray fMatrixArray9; // array of transnformtation matrix
157 TObjArray fMatrixArray6; // array of transnformtation matrix
160 TObjArray fCombinedMatrixArray6; // array combeined transformation matrix
162 AliExternalComparison *fCompTracklet; //tracklet comparison
164 Int_t fPoints[72*72]; // number of points in the fitter
165 Bool_t fNoField; // flag - no field data
167 Double_t fXIO; // OROC-IROC boundary
168 Double_t fXmiddle; // center of the TPC sector local X
169 Double_t fXquadrant; // x quadrant
171 // Kalman filter for sector internal alignemnt
173 TObjArray fArraySectorIntParam; // array of sector alignment parameters
174 TObjArray fArraySectorIntCovar; // array of sector alignment covariances
176 // Kalman filter for global alignment
178 TMatrixD *fSectorParamA; // Kalman parameter for A side
179 TMatrixD *fSectorCovarA; // Kalman covariance for A side
180 TMatrixD *fSectorParamC; // Kalman parameter for A side
181 TMatrixD *fSectorCovarC; // Kalman covariance for A side
185 Bool_t fUseInnerOuter; // flag- use Inner Outer sector for left righ alignment
187 static AliTPCcalibAlign* fgInstance; //! Instance of this class (singleton implementation)
189 AliTPCcalibAlign& operator=(const AliTPCcalibAlign&);// not implemented
191 ClassDef(AliTPCcalibAlign,3)
195 TLinearFitter* AliTPCcalibAlign::GetFitter12(Int_t s1,Int_t s2) {
196 return static_cast<TLinearFitter*>(fFitterArray12[GetIndex(s1,s2)]);
198 TLinearFitter* AliTPCcalibAlign::GetFitter9(Int_t s1,Int_t s2) {
199 return static_cast<TLinearFitter*>(fFitterArray9[GetIndex(s1,s2)]);
201 TLinearFitter* AliTPCcalibAlign::GetFitter6(Int_t s1,Int_t s2) {
202 return static_cast<TLinearFitter*>(fFitterArray6[GetIndex(s1,s2)]);
205 const TMatrixD * AliTPCcalibAlign::GetTransformation(Int_t s1,Int_t s2, Int_t fitType){
206 if (fitType==0) return static_cast<TMatrixD*>(fMatrixArray6[GetIndex(s1,s2)]);
207 if (fitType==1) return static_cast<TMatrixD*>(fMatrixArray9[GetIndex(s1,s2)]);
208 if (fitType==2) return static_cast<TMatrixD*>(fMatrixArray12[GetIndex(s1,s2)]);