]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibAlign.cxx
25936ec5d54c20690cee86b11917eac2f096b342
[u/mrichter/AliRoot.git] / TPC / AliTPCcalibAlign.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16
17 ///////////////////////////////////////////////////////////////////////////////
18 //                                                                           //
19 //     Class to make a internal alignemnt of TPC chambers                    //
20 //
21 //     Different linear tranformation investigated
22
23 //     12 parameters - arbitrary linear transformation 
24 //                     a00  a01 a02  a03     p[0]   p[1]  p[2]  p[9]
25 //                     a10  a11 a12  a13 ==> p[3]   p[4]  p[5]  p[10]
26 //                     a20  a21 a22  a23     p[6]   p[7]  p[8]  p[11] 
27 //
28 //      9 parameters - scaling fixed to 1
29 //                     a00  a01  a02 a03     1      p[0]  p[1]   p[6]
30 //                     a10  a11  a12 a13 ==> p[2]   1     p[3]   p[7]
31 //                     a20  a21  a22 a23     p[4]   p[5]  1      p[8] 
32 //
33 //      6 parameters - x-y rotation x-z, y-z tiliting
34 //                     a00  a01  a02 a03     1     -p[0]  0     p[3]
35 //                     a10  a11  a12 a13 ==> p[0]   1     0     p[4]
36 //                     a20  a21  a22 a23     p[1]   p[2]  1     p[5] 
37 //
38 ////
39 //// 
40
41 #include "TLinearFitter.h"
42 #include "AliTPCcalibAlign.h"
43 #include "AliExternalTrackParam.h"
44 #include "AliTPCTracklet.h"
45 #include "TH1D.h"
46 #include "TVectorD.h"
47 #include "TTreeStream.h"
48 #include "TFile.h"
49 #include "TF1.h"
50
51 #include <iostream>
52 #include <sstream>
53 using namespace std;
54
55 ClassImp(AliTPCcalibAlign)
56
57 AliTPCcalibAlign::AliTPCcalibAlign()
58   :  fDphiHistArray(72*72),
59      fDthetaHistArray(72*72),
60      fDyHistArray(72*72),
61      fDzHistArray(72*72),
62      fFitterArray12(72*72),
63      fFitterArray9(72*72),
64      fFitterArray6(72*72)
65 {
66   //
67   // Constructor
68   //
69   for (Int_t i=0;i<72*72;++i) {
70     fPoints[i]=0;
71   }
72 }
73
74 AliTPCcalibAlign::AliTPCcalibAlign(const Text_t *name, const Text_t *title)
75   :AliTPCcalibBase(),  
76    fDphiHistArray(72*72),
77    fDthetaHistArray(72*72),
78    fDyHistArray(72*72),
79    fDzHistArray(72*72),
80    fFitterArray12(72*72),
81    fFitterArray9(72*72),
82    fFitterArray6(72*72)
83 {
84   //
85   // Constructor
86   //
87   SetName(name);
88   SetTitle(title);
89   for (Int_t i=0;i<72*72;++i) {
90     fPoints[i]=0;
91   }
92 }
93
94 AliTPCcalibAlign::~AliTPCcalibAlign() {
95   //
96   // destructor
97   //
98 }
99
100 void AliTPCcalibAlign::Process(AliTPCseed *seed) {
101   TObjArray tracklets=
102     AliTPCTracklet::CreateTracklets(seed,AliTPCTracklet::kKalman,
103                                     kFALSE,20,2);
104  //  TObjArray trackletsL=
105 //     AliTPCTracklet::CreateTracklets(seed,AliTPCTracklet::kLinear,
106 //                                  kFALSE,20,2);
107 //   TObjArray trackletsQ=
108 //     AliTPCTracklet::CreateTracklets(seed,AliTPCTracklet::kQuadratic,
109 //                                  kFALSE,20,2);
110 //   TObjArray trackletsR=
111 //     AliTPCTracklet::CreateTracklets(seed,AliTPCTracklet::kRiemann,
112 //                                  kFALSE,20,2);
113   tracklets.SetOwner();
114  //  trackletsL.SetOwner();
115 //   trackletsQ.SetOwner();
116 //   trackletsR.SetOwner();
117   if (tracklets.GetEntries()==2) {
118     AliTPCTracklet *t1=static_cast<AliTPCTracklet*>(tracklets[0]);
119     AliTPCTracklet *t2=static_cast<AliTPCTracklet*>(tracklets[1]);
120     if (t1->GetSector()>t2->GetSector()) {
121       AliTPCTracklet* tmp=t1;
122       t1=t2;
123       t2=tmp;
124     }
125     AliExternalTrackParam *common1=0,*common2=0;
126     if (AliTPCTracklet::PropagateToMeanX(*t1,*t2,common1,common2))
127       ProcessTracklets(*common1,*common2,t1->GetSector(),t2->GetSector());
128     delete common1;
129     delete common2;
130   }
131   
132 }
133
134 void AliTPCcalibAlign::Analyze(){
135   //
136   // Analyze function 
137   //
138   EvalFitters();
139 }
140
141
142 void AliTPCcalibAlign::Terminate(){
143   //
144   // Terminate function
145   // call base terminate + Eval of fitters
146   //
147   if (GetDebugLevel()>0) Info("AliTPCcalibAlign","Teminate");
148   EvalFitters();
149   AliTPCcalibBase::Terminate();
150 }
151
152
153
154
155 void AliTPCcalibAlign::ProcessTracklets(const AliExternalTrackParam &tp1,
156                                         const AliExternalTrackParam &tp2,
157                                         Int_t s1,Int_t s2) {
158
159   //
160   //
161   //
162   FillHisto(tp1,tp2,s1,s2);
163   //
164   // Process function to fill fitters
165   //
166   Double_t t1[5],t2[5];
167   Double_t &x1=t1[0], &y1=t1[1], &z1=t1[2], &dydx1=t1[3], &dzdx1=t1[4];
168   Double_t &x2=t2[0], &y2=t2[1], &z2=t2[2], &dydx2=t2[3], &dzdx2=t2[4];
169   x1   =tp1.GetX();
170   y1   =tp1.GetY();
171   z1   =tp1.GetZ();
172   Double_t snp1=tp1.GetSnp();
173   dydx1=snp1/TMath::Sqrt(1.-snp1*snp1);
174   Double_t tgl1=tp1.GetTgl();
175   // dz/dx = 1/(cos(theta)*cos(phi))
176   dzdx1=tgl1/TMath::Sqrt(1.-snp1*snp1);
177   x2   =tp2.GetX();
178   y2   =tp2.GetY();
179   z2   =tp2.GetZ();
180   Double_t snp2=tp2.GetSnp();
181   dydx2=snp2/TMath::Sqrt(1.-snp2*snp2);
182   Double_t tgl2=tp2.GetTgl();
183   dzdx2=tgl2/TMath::Sqrt(1.-snp2*snp2);
184   //
185   //
186   //
187   if (fStreamLevel>1){
188     TTreeSRedirector *cstream = GetDebugStreamer();
189     if (cstream){
190       static TVectorD vec1(5);
191       static TVectorD vec2(5);
192       vec1.SetElements(t1);
193       vec2.SetElements(t2);
194       AliExternalTrackParam *p1 = &((AliExternalTrackParam&)tp1);
195       AliExternalTrackParam *p2 = &((AliExternalTrackParam&)tp2);
196       (*cstream)<<"Tracklet"<<
197         "tp1.="<<p1<<
198         "tp2.="<<p2<<
199         "v1.="<<&vec1<<
200         "v2.="<<&vec2<<
201         "s1="<<s1<<
202         "s2="<<s2<<
203         "\n";
204     }
205   }
206   Process12(t1,t2,GetOrMakeFitter12(s1,s2));
207   Process9(t1,t2,GetOrMakeFitter9(s1,s2));
208   Process6(t1,t2,GetOrMakeFitter6(s1,s2));
209   ++fPoints[GetIndex(s1,s2)];
210 }
211
212 void AliTPCcalibAlign::Process12(const Double_t *t1,
213                                  const Double_t *t2,
214                                  TLinearFitter *fitter) {
215   // x2    =  a00*x1 + a01*y1 + a02*z1 + a03
216   // y2    =  a10*x1 + a11*y1 + a12*z1 + a13
217   // z2    =  a20*x1 + a21*y1 + a22*z1 + a23
218   // dydx2 = (a10    + a11*dydx1 + a12*dzdx1)/(a00    + a01*dydx1 + a02*dzdx1)
219   // dzdx2 = (a20    + a21*dydx1 + a22*dzdx1)/(a00    + a01*dydx1 + a02*dzdx1)
220   //
221   //                     a00  a01 a02  a03     p[0]   p[1]  p[2]  p[9]
222   //                     a10  a11 a12  a13 ==> p[3]   p[4]  p[5]  p[10]
223   //                     a20  a21 a22  a23     p[6]   p[7]  p[8]  p[11] 
224
225
226
227   const Double_t &x1=t1[0], &y1=t1[1], &z1=t1[2], &dydx1=t1[3], &dzdx1=t1[4];
228   const Double_t /*&x2=t2[0],*/ &y2=t2[1], &z2=t2[2], &dydx2=t2[3], &dzdx2=t2[4];
229
230   // TODO:
231   Double_t sy    = 0.1;
232   Double_t sz    = 0.1;
233   Double_t sdydx = 0.001;
234   Double_t sdzdx = 0.001;
235
236   Double_t p[12];
237   Double_t value;
238
239   // x2  =  a00*x1 + a01*y1 + a02*z1 + a03
240   // y2  =  a10*x1 + a11*y1 + a12*z1 + a13
241   // y2' =  a10*x1 + a11*y1 + a12*z1 + a13 + (a01*y1 + a02*z1 + a03)*dydx2
242   for (Int_t i=0; i<12;i++) p[i]=0.;
243   p[3+0] = x1;          // a10
244   p[3+1] = y1;          // a11
245   p[3+2] = z1;          // a12
246   p[9+1] = 1.;          // a13
247   p[0+1] = y1*dydx2;    // a01
248   p[0+2] = z1*dydx2;    // a02
249   p[9+0] = dydx2;       // a03
250   value  = y2;
251   fitter->AddPoint(p,value,sy);
252
253   // x2  =  a00*x1 + a01*y1 + a02*z1 + a03
254   // z2  =  a20*x1 + a21*y1 + a22*z1 + a23
255   // z2' =  a20*x1 + a21*y1 + a22*z1 + a23 + (a01*y1 + a02*z1 + a03)*dzdx2;
256   for (Int_t i=0; i<12;i++) p[i]=0.;
257   p[6+0] = x1;           // a20 
258   p[6+1] = y1;           // a21
259   p[6+2] = z1;           // a22
260   p[9+2] = 1.;           // a23
261   p[0+1] = y1*dzdx2;     // a01
262   p[0+2] = z1*dzdx2;     // a02
263   p[9+0] = dzdx2;        // a03
264   value  = z2;
265   fitter->AddPoint(p,value,sz);
266
267   // dydx2 = (a10 + a11*dydx1 + a12*dzdx1)/( a00 + a01*dydx1 + a02*dzdx1)
268   // (a10 + a11*dydx1 + a12*dzdx1) - (a00 + a01*dydx1 + a02*dzdx1)*dydx2 = 0
269   for (Int_t i=0; i<12;i++) p[i]=0.;
270   p[3+0] = 1.;           // a10
271   p[3+1] = dydx1;        // a11
272   p[3+2] = dzdx1;        // a12
273   p[0+0] = -dydx2;       // a00
274   p[0+1] = -dydx1*dydx2; // a01
275   p[0+2] = -dzdx1*dydx2; // a02
276   value  = 0.;
277   fitter->AddPoint(p,value,sdydx);
278
279   // dzdx2 = (a20 + a21*dydx1 + a22*dzdx1)/( a00 + a01*dydx1 + a02*dzdx1)
280   // (a20 + a21*dydx1 + a22*dzdx1) - (a00 + a01*dydx1 + a02*dzdx1)*dzdx2 = 0
281   for (Int_t i=0; i<12;i++) p[i]=0.;
282   p[6+0] = 1;            // a20
283   p[6+1] = dydx1;        // a21
284   p[6+2] = dzdx1;        // a22
285   p[0+0] = -dzdx2;       // a00
286   p[0+1] = -dydx1*dzdx2; // a01
287   p[0+2] = -dzdx1*dzdx2; // a02
288   value  = 0.;
289   fitter->AddPoint(p,value,sdzdx);
290 }
291
292 void AliTPCcalibAlign::Process9(Double_t *t1,
293                                 Double_t *t2,
294                                 TLinearFitter *fitter) {
295   // x2    =  a00*x1 + a01*y1 + a02*z1 + a03
296   // y2    =  a10*x1 + a11*y1 + a12*z1 + a13
297   // z2    =  a20*x1 + a21*y1 + a22*z1 + a23
298   // dydx2 = (a10    + a11*dydx1 + a12*dzdx1)/(a00    + a01*dydx1 + a02*dzdx1)
299   // dzdx2 = (a20    + a21*dydx1 + a22*dzdx1)/(a00    + a01*dydx1 + a02*dzdx1)
300   //
301   //                     a00  a01  a02 a03     1      p[0]  p[1]   p[6]
302   //                     a10  a11  a12 a13 ==> p[2]   1     p[3]   p[7]
303   //                     a20  a21  a21 a23     p[4]   p[5]  1      p[8] 
304
305
306   Double_t &x1=t1[0], &y1=t1[1], &z1=t1[2], &dydx1=t1[3], &dzdx1=t1[4];
307   Double_t /*&x2=t2[0],*/ &y2=t2[1], &z2=t2[2], &dydx2=t2[3], &dzdx2=t2[4];
308
309   // TODO:
310   Double_t sy    = 0.1;
311   Double_t sz    = 0.1;
312   Double_t sdydx = 0.001;
313   Double_t sdzdx = 0.001;
314   //
315   Double_t p[12];
316   Double_t value;
317
318   // x2  =  a00*x1 + a01*y1 + a02*z1 + a03
319   // y2  =  a10*x1 + a11*y1 + a12*z1 + a13
320   // y2' =  a10*x1 + a11*y1 + a12*z1 + a13 + (a01*y1 + a02*z1 + a03)*dydx2
321   for (Int_t i=0; i<12;i++) p[i]=0.;
322   p[2]   += x1;           // a10
323   //p[]  +=1;             // a11
324   p[3]   += z1;           // a12    
325   p[7]   += 1;            // a13
326   p[0]   += y1*dydx2;     // a01
327   p[1]   += z1*dydx2;     // a02
328   p[6]   += dydx2;        // a03
329   value   = y2-y1;        //-a11
330   fitter->AddPoint(p,value,sy);
331   //
332   // x2  =  a00*x1 + a01*y1 + a02*z1 + a03
333   // z2  =  a20*x1 + a21*y1 + a22*z1 + a23
334   // z2' =  a20*x1 + a21*y1 + a22*z1 + a23 + (a01*y1 + a02*z1 + a03)*dzdx2;
335   for (Int_t i=0; i<12;i++) p[i]=0.;
336   p[4]   += x1;           // a20 
337   p[5]   += y1;           // a21
338   //p[]  += z1;           // a22
339   p[8]   += 1.;           // a23
340   p[0]   += y1*dzdx2;     // a01
341   p[1]   += z1*dzdx2;     // a02
342   p[6]   += dzdx2;        // a03
343   value  = z2-z1;         //-a22
344   fitter->AddPoint(p,value,sz);
345
346   // dydx2 = (a10 + a11*dydx1 + a12*dzdx1)/( a00 + a01*dydx1 + a02*dzdx1)
347   // (a10 + a11*dydx1 + a12*dzdx1) - (a00 + a01*dydx1 + a02*dzdx1)*dydx2 = 0
348   for (Int_t i=0; i<12;i++) p[i]=0.;
349   p[2]   += 1.;           // a10
350   //p[]  += dydx1;      // a11
351   p[3]   += dzdx1;        // a12
352   //p[]  += -dydx2;       // a00
353   p[0]   += -dydx1*dydx2; // a01
354   p[1]   += -dzdx1*dydx2; // a02
355   value  = -dydx1+dydx2;  // -a11 + a00
356   fitter->AddPoint(p,value,sdydx);
357
358   // dzdx2 = (a20 + a21*dydx1 + a22*dzdx1)/( a00 + a01*dydx1 + a02*dzdx1)
359   // (a20 + a21*dydx1 + a22*dzdx1) - (a00 + a01*dydx1 + a02*dzdx1)*dzdx2 = 0
360   for (Int_t i=0; i<12;i++) p[i]=0.;
361   p[4]   += 1;            // a20
362   p[5]   += dydx1;        // a21
363   //p[]  += dzdx1;        // a22
364   //p[]  += -dzdx2;       // a00
365   p[0]   += -dydx1*dzdx2; // a01
366   p[1]   += -dzdx1*dzdx2; // a02
367   value  = -dzdx1+dzdx2;  // -a22 + a00
368   fitter->AddPoint(p,value,sdzdx);
369 }
370
371 void AliTPCcalibAlign::Process6(Double_t *t1,
372                                 Double_t *t2,
373                                 TLinearFitter *fitter) {
374   // x2    =  1  *x1 +-a01*y1 + 0      +a03
375   // y2    =  a01*x1 + 1  *y1 + 0      +a13
376   // z2    =  a20*x1 + a21*y1 + 1  *z1 +a23
377   // dydx2 = (a10    + a11*dydx1 + a12*dzdx1)/(a00    + a01*dydx1 + a02*dzdx1)
378   // dzdx2 = (a20    + a21*dydx1 + a22*dzdx1)/(a00    + a01*dydx1 + a02*dzdx1)
379   //
380   //                     a00  a01  a02 a03     1     -p[0]  0     p[3]
381   //                     a10  a11  a12 a13 ==> p[0]   1     0     p[4]
382   //                     a20  a21  a21 a23     p[1]   p[2]  1     p[5] 
383
384   Double_t &x1=t1[0], &y1=t1[1], &z1=t1[2], &dydx1=t1[3], &dzdx1=t1[4];
385   Double_t /*&x2=t2[0],*/ &y2=t2[1], &z2=t2[2], &dydx2=t2[3], &dzdx2=t2[4];
386
387   // TODO:
388   Double_t sy    = 0.1;
389   Double_t sz    = 0.1;
390   Double_t sdydx = 0.001;
391   Double_t sdzdx = 0.001;
392
393   Double_t p[12];
394   Double_t value;
395   // x2  =  a00*x1 + a01*y1 + a02*z1 + a03
396   // y2  =  a10*x1 + a11*y1 + a12*z1 + a13
397   // y2' =  a10*x1 + a11*y1 + a12*z1 + a13 + (a01*y1 + a02*z1 + a03)*dydx2
398   for (Int_t i=0; i<12;i++) p[i]=0.;
399   p[0]   += x1;           // a10
400   //p[]  +=1;             // a11
401   //p[]  += z1;           // a12    
402   p[4]   += 1;            // a13
403   p[0]   += -y1*dydx2;    // a01
404   //p[]  += z1*dydx2;     // a02
405   p[3]   += dydx2;        // a03
406   value   = y2-y1;        //-a11
407   fitter->AddPoint(p,value,sy);
408   //
409   // x2  =  a00*x1 + a01*y1 + a02*z1 + a03
410   // z2  =  a20*x1 + a21*y1 + a22*z1 + a23
411   // z2' =  a20*x1 + a21*y1 + a22*z1 + a23 + (a01*y1 + a02*z1 + a03)*dzdx2;
412   for (Int_t i=0; i<12;i++) p[i]=0.;
413   p[1]   += x1;           // a20 
414   p[2]   += y1;           // a21
415   //p[]  += z1;           // a22
416   p[5]   += 1.;           // a23
417   p[0]   += -y1*dzdx2;    // a01
418   //p[]   += z1*dzdx2;     // a02
419   p[3]   += dzdx2;        // a03
420   value  = z2-z1;         //-a22
421   fitter->AddPoint(p,value,sz);
422
423   // dydx2 = (a10 + a11*dydx1 + a12*dzdx1)/( a00 + a01*dydx1 + a02*dzdx1)
424   // (a10 + a11*dydx1 + a12*dzdx1) - (a00 + a01*dydx1 + a02*dzdx1)*dydx2 = 0
425   for (Int_t i=0; i<12;i++) p[i]=0.;
426   p[0]   += 1.;           // a10
427   //p[]  += dydx1;      // a11
428   //p[]   += dzdx1;        // a12
429   //p[]  += -dydx2;       // a00
430   p[0]   +=  dydx1*dydx2; // a01
431   //p[]   += -dzdx1*dydx2; // a02
432   value  = -dydx1+dydx2;  // -a11 + a00
433   fitter->AddPoint(p,value,sdydx);
434
435   // dzdx2 = (a20 + a21*dydx1 + a22*dzdx1)/( a00 + a01*dydx1 + a02*dzdx1)
436   // (a20 + a21*dydx1 + a22*dzdx1) - (a00 + a01*dydx1 + a02*dzdx1)*dzdx2 = 0
437   for (Int_t i=0; i<12;i++) p[i]=0.;
438   p[1]   += 1;            // a20
439   p[2]   += dydx1;        // a21
440   //p[]  += dzdx1;        // a22
441   //p[]  += -dzdx2;       // a00
442   p[0]   +=  dydx1*dzdx2; // a01
443   //p[]  += -dzdx1*dzdx2; // a02
444   value  = -dzdx1+dzdx2;  // -a22 + a00
445   fitter->AddPoint(p,value,sdzdx);
446 }
447
448
449
450
451 void AliTPCcalibAlign::EvalFitters() {
452   //
453   // Analyze function 
454   // 
455   // Perform the fitting using linear fitters
456   //
457   Int_t kMinPoints =50;
458   TLinearFitter *f;
459   TFile fff("alignDebug.root","recreate");
460   for (Int_t s1=0;s1<72;++s1)
461     for (Int_t s2=0;s2<72;++s2){
462       if ((f=GetFitter12(s1,s2))&&fPoints[GetIndex(s1,s2)]>kMinPoints) {
463         //      cerr<<s1<<","<<s2<<": "<<fPoints[GetIndex(s1,s2)]<<endl;
464         if (f->Eval()!=0) {
465           cerr<<"Evaluation failed for "<<s1<<","<<s2<<endl;
466           f->Write(Form("f12_%d_%d",s1,s2));
467         }else{
468           f->Write(Form("f12_%d_%d",s1,s2));
469         }
470       }
471       if ((f=GetFitter9(s1,s2))&&fPoints[GetIndex(s1,s2)]>kMinPoints) {
472         //      cerr<<s1<<","<<s2<<": "<<fPoints[GetIndex(s1,s2)]<<endl;
473         if (f->Eval()!=0) {
474           cerr<<"Evaluation failed for "<<s1<<","<<s2<<endl;
475         }else{
476           f->Write(Form("f9_%d_%d",s1,s2));
477         }
478       }
479       if ((f=GetFitter6(s1,s2))&&fPoints[GetIndex(s1,s2)]>kMinPoints) {
480         //      cerr<<s1<<","<<s2<<": "<<fPoints[GetIndex(s1,s2)]<<endl;
481         if (f->Eval()!=0) {
482           cerr<<"Evaluation failed for "<<s1<<","<<s2<<endl;
483         }else{
484           f->Write(Form("f6_%d_%d",s1,s2));
485         }
486       }
487     }
488   this->Write("align");
489   /*
490                     
491   fitter->Eval();
492   fitter->Eval();
493   chi212 = align->GetChisquare()/(4.*entries);
494
495   TMatrixD mat(13,13);
496   TVectorD par(13);
497   align->GetParameters(par);
498   align->GetCovarianceMatrix(mat);
499
500   //
501   //
502   for (Int_t i=0; i<12;i++){
503     palign12(i)= par(i+1);
504     for (Int_t j=0; j<12;j++){
505       pcovar12(i,j)   = mat(i+1,j+1);
506       pcovar12(i,j) *= chi212;
507     }
508   }
509   //
510   for (Int_t i=0; i<12;i++){
511     psigma12(i)  = TMath::Sqrt(pcovar12(i,i));
512     palignR12(i) = palign12(i)/TMath::Sqrt(pcovar12(i,i));
513     for (Int_t j=0; j<12;j++){
514       pcovarN12(i,j) = pcovar12(i,j)/TMath::Sqrt(pcovar12(i,i)*pcovar12(j,j));
515     }
516   }
517   */
518 }
519
520 TLinearFitter* AliTPCcalibAlign::GetOrMakeFitter12(Int_t s1,Int_t s2) {
521   //
522   // get or make fitter - general linear transformation
523   //
524   static Int_t counter12=0;
525   static TF1 f12("f12","x[0]++x[1]++x[2]++x[3]++x[4]++x[5]++x[6]++x[7]++x[8]++x[9]++x[10]++x[11]");
526   TLinearFitter * fitter = GetFitter12(s1,s2);
527   if (fitter) return fitter;
528   //  fitter =new TLinearFitter(12,"x[0]++x[1]++x[2]++x[3]++x[4]++x[5]++x[6]++x[7]++x[8]++x[9]++x[10]++x[11]");
529   fitter =new TLinearFitter(&f12,"");
530   fitter->StoreData(kFALSE);
531   fFitterArray12.AddAt(fitter,GetIndex(s1,s2)); 
532   counter12++;
533   if (GetDebugLevel()>0) cerr<<"Creating fitter12 "<<s1<<","<<s2<<"  :  "<<counter12<<endl;
534   return fitter;
535 }
536
537 TLinearFitter* AliTPCcalibAlign::GetOrMakeFitter9(Int_t s1,Int_t s2) {
538   //
539   //get or make fitter - general linear transformation - no scaling
540   // 
541   static Int_t counter9=0;
542   static TF1 f9("f9","x[0]++x[1]++x[2]++x[3]++x[4]++x[5]++x[6]++x[7]++x[8]");
543   TLinearFitter * fitter = GetFitter9(s1,s2);
544   if (fitter) return fitter;
545   //  fitter =new TLinearFitter(9,"x[0]++x[1]++x[2]++x[3]++x[4]++x[5]++x[6]++x[7]++x[8]");
546   fitter =new TLinearFitter(&f9,"");
547   fitter->StoreData(kFALSE);
548   fFitterArray9.AddAt(fitter,GetIndex(s1,s2));
549   counter9++;
550   if (GetDebugLevel()>0) cerr<<"Creating fitter12 "<<s1<<","<<s2<<"  :  "<<counter9<<endl;
551   return fitter;
552 }
553
554 TLinearFitter* AliTPCcalibAlign::GetOrMakeFitter6(Int_t s1,Int_t s2) {
555   //
556   // get or make fitter  - 6 paramater linear tranformation
557   //                     - no scaling
558   //                     - rotation x-y
559   //                     - tilting x-z, y-z
560   static Int_t counter6=0;
561   static TF1 f6("f6","x[0]++x[1]++x[2]++x[3]++x[4]++x[5]");
562   TLinearFitter * fitter = GetFitter6(s1,s2);
563   if (fitter) return fitter;
564   //  fitter=new TLinearFitter(6,"x[0]++x[1]++x[2]++x[3]++x[4]++x[5]");
565   fitter=new TLinearFitter(&f6,"");
566   fitter->StoreData(kFALSE);
567   fFitterArray6.AddAt(fitter,GetIndex(s1,s2));
568   counter6++;
569   if (GetDebugLevel()>0) cerr<<"Creating fitter6 "<<s1<<","<<s2<<"  :  "<<counter6<<endl;
570   return fitter;
571 }
572
573
574
575
576
577 Bool_t AliTPCcalibAlign::GetTransformation12(Int_t s1,Int_t s2,TMatrixD &a) {
578   //
579   // GetTransformation matrix - 12 paramaters - generael linear transformation
580   //
581   if (!GetFitter12(s1,s2))
582     return false;
583   else {
584     TVectorD p(12);
585     GetFitter12(s1,s2)->GetParameters(p);
586     a.ResizeTo(4,4);
587     a[0][0]=p[0]; a[0][1]=p[1]; a[0][2]=p[2]; a[0][3]=p[9];
588     a[1][0]=p[3]; a[1][1]=p[4]; a[1][2]=p[5]; a[1][3]=p[10];
589     a[2][0]=p[6]; a[2][1]=p[7]; a[2][2]=p[8]; a[2][3]=p[11];
590     a[3][0]=0.;   a[3][1]=0.;   a[3][2]=0.;   a[3][3]=1.;
591     return true;
592   } 
593 }
594
595 Bool_t AliTPCcalibAlign::GetTransformation9(Int_t s1,Int_t s2,TMatrixD &a) {
596   //
597   // GetTransformation matrix - 9 paramaters - general linear transformation
598   //                            No scaling
599   //
600   if (!GetFitter9(s1,s2))
601     return false;
602   else {
603     TVectorD p(9);
604     GetFitter9(s1,s2)->GetParameters(p);
605     a.ResizeTo(4,4);
606     a[0][0]=1;    a[0][1]=p[0]; a[0][2]=p[1]; a[0][3]=p[6];
607     a[1][0]=p[2]; a[1][1]=1;    a[1][2]=p[3]; a[1][3]=p[7];
608     a[2][0]=p[4]; a[2][1]=p[5]; a[2][2]=1;    a[2][3]=p[8];
609     a[3][0]=0.;   a[3][1]=0.;   a[3][2]=0.;   a[3][3]=1.;
610     return true;
611   } 
612 }
613
614 Bool_t AliTPCcalibAlign::GetTransformation6(Int_t s1,Int_t s2,TMatrixD &a) {
615   //
616   // GetTransformation matrix - 6  paramaters
617   //                            3  translation
618   //                            1  rotation -x-y  
619   //                            2  tilting x-z y-z
620   if (!GetFitter6(s1,s2))
621     return false;
622   else {
623     TVectorD p(6);
624     GetFitter6(s1,s2)->GetParameters(p);
625     a.ResizeTo(4,4);
626     a[0][0]=1;       a[0][1]=-p[0];a[0][2]=0;    a[0][3]=p[3];
627     a[1][0]=p[0];    a[1][1]=1;    a[1][2]=0;    a[1][3]=p[4];
628     a[2][0]=p[1];    a[2][1]=p[2]; a[2][2]=1;    a[2][3]=p[5];
629     a[3][0]=0.;      a[3][1]=0.;   a[3][2]=0.;   a[3][3]=1.;
630     return true;
631   } 
632 }
633
634 void AliTPCcalibAlign::FillHisto(const AliExternalTrackParam &tp1,
635                                         const AliExternalTrackParam &tp2,
636                                         Int_t s1,Int_t s2) {
637   //
638   // Fill residual histograms
639   //
640   if (s2-s1==36) {//only inner-outer
641     GetHisto(kPhi,s1,s2,kTRUE)->Fill(TMath::ASin(tp1.GetSnp())-TMath::ASin(tp2.GetSnp()));    
642     GetHisto(kTheta,s1,s2,kTRUE)->Fill(TMath::ATan(tp1.GetTgl())-TMath::ATan(tp2.GetTgl()));
643     GetHisto(kY,s1,s2,kTRUE)->Fill(tp1.GetY()-tp2.GetY());
644     GetHisto(kZ,s1,s2,kTRUE)->Fill(tp1.GetZ()-tp2.GetZ());
645   }  
646 }
647
648
649
650 TH1 * AliTPCcalibAlign::GetHisto(HistoType type, Int_t s1, Int_t s2, Bool_t force)
651 {
652   //
653   // return specified residual histogram - it is only QA
654   // if force specified the histogram and given histogram is not existing 
655   //  new histogram is created
656   //
657   if (GetIndex(s1,s2)>=72*72) return 0;
658   TObjArray *histoArray=0;
659   switch (type) {
660   case kY:
661     histoArray = &fDyHistArray; break;
662   case kZ:
663     histoArray = &fDzHistArray; break;
664   case kPhi:
665     histoArray = &fDphiHistArray; break;
666   case kTheta:
667     histoArray = &fDthetaHistArray; break;
668   }
669   TH1 * histo= (TH1*)histoArray->At(GetIndex(s1,s2));
670   if (histo) return histo;
671   if (force==kFALSE) return 0; 
672   //
673   stringstream name;
674   stringstream title;
675   switch (type) {    
676   case kY:
677     name<<"hist_y_"<<s1<<"_"<<s2;
678     title<<"Y Missalignment for sectors "<<s1<<" and "<<s2;
679     histo =new TH1D(name.str().c_str(),title.str().c_str(),512,-0.3,0.3); // +/- 3 mm
680     break;
681   case kZ:
682     name<<"hist_z_"<<s1<<"_"<<s2;
683     title<<"Z Missalignment for sectors "<<s1<<" and "<<s2;
684     histo = new TH1D(name.str().c_str(),title.str().c_str(),512,-0.3,0.3); // +/- 3 mm
685     break;
686   case kPhi:
687     name<<"hist_phi_"<<s1<<"_"<<s2;
688     title<<"Phi Missalignment for sectors "<<s1<<" and "<<s2;
689     histo =new TH1D(name.str().c_str(),title.str().c_str(),512,-0.01,0.01); // +/- 10 mrad
690     break;
691   case kTheta:
692     name<<"hist_theta_"<<s1<<"_"<<s2;
693     title<<"Theta Missalignment for sectors "<<s1<<" and "<<s2;
694     histo =new TH1D(name.str().c_str(),title.str().c_str(),512,-0.01,0.01); // +/- 10 mrad
695     break;
696   }
697   histo->SetDirectory(0);
698   histoArray->AddAt(histo,GetIndex(s1,s2));
699   return histo;
700 }