]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibAlign.cxx
Additional protection: no update of V0 momentum without mag. field
[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 //     12 parameters - arbitrary transformation 
23 //      9 parameters - scaling fixed to 1
24 //      6 parameters - 
25 ////
26 //// Only the 12 parameter version works so far...
27 ////
28 ////
29 //// 
30
31 #include "TLinearFitter.h"
32 #include "AliTPCcalibAlign.h"
33 #include "AliExternalTrackParam.h"
34 #include "AliTPCTracklet.h"
35 #include "TH1D.h"
36
37 #include <iostream>
38 #include <sstream>
39 using namespace std;
40
41 ClassImp(AliTPCcalibAlign)
42
43 AliTPCcalibAlign::AliTPCcalibAlign()
44   :  fDphiHistArray(72*72),
45      fDthetaHistArray(72*72),
46      fDyHistArray(72*72),
47      fDzHistArray(72*72),
48      fFitterArray12(72*72),
49      fFitterArray9(72*72),
50      fFitterArray6(72*72)
51 {
52   //
53   // Constructor
54   //
55   for (Int_t i=0;i<72*72;++i) {
56     fPoints[i]=0;
57   }
58 }
59
60 AliTPCcalibAlign::~AliTPCcalibAlign() {
61   //
62   // destructor
63   //
64 }
65
66 void AliTPCcalibAlign::Process(AliTPCseed *seed) {
67   TObjArray tracklets=
68     AliTPCTracklet::CreateTracklets(seed,AliTPCTracklet::kKalman,
69                                     kFALSE,20,2);
70  //  TObjArray trackletsL=
71 //     AliTPCTracklet::CreateTracklets(seed,AliTPCTracklet::kLinear,
72 //                                  kFALSE,20,2);
73 //   TObjArray trackletsQ=
74 //     AliTPCTracklet::CreateTracklets(seed,AliTPCTracklet::kQuadratic,
75 //                                  kFALSE,20,2);
76 //   TObjArray trackletsR=
77 //     AliTPCTracklet::CreateTracklets(seed,AliTPCTracklet::kRiemann,
78 //                                  kFALSE,20,2);
79   tracklets.SetOwner();
80  //  trackletsL.SetOwner();
81 //   trackletsQ.SetOwner();
82 //   trackletsR.SetOwner();
83   if (tracklets.GetEntries()==2) {
84     AliTPCTracklet *t1=static_cast<AliTPCTracklet*>(tracklets[0]);
85     AliTPCTracklet *t2=static_cast<AliTPCTracklet*>(tracklets[1]);
86     if (t1->GetSector()>t2->GetSector()) {
87       AliTPCTracklet* tmp=t1;
88       t1=t2;
89       t2=tmp;
90     }
91     AliExternalTrackParam *common1=0,*common2=0;
92     if (AliTPCTracklet::PropagateToMeanX(*t1,*t2,common1,common2))
93       ProcessTracklets(*common1,*common2,t1->GetSector(),t2->GetSector());
94     delete common1;
95     delete common2;
96   }
97   
98 }
99
100 void AliTPCcalibAlign::ProcessTracklets(const AliExternalTrackParam &tp1,
101                                         const AliExternalTrackParam &tp2,
102                                         Int_t s1,Int_t s2) {
103
104   if (s2-s1==36) {//only inner-outer
105     if (!fDphiHistArray[s1*72+s2]) {
106       stringstream name;
107       stringstream title;
108       name<<"hist_phi_"<<s1<<"_"<<s2;
109       title<<"Phi Missalignment for sectors "<<s1<<" and "<<s2;
110       fDphiHistArray[s1*72+s2]=new TH1D(name.str().c_str(),title.str().c_str(),1024,-0.01,0.01); // +/- 10 mrad
111       ((TH1D*)fDphiHistArray[s1*72+s2])->SetDirectory(0);
112     }
113     if (!fDthetaHistArray[s1*72+s2]) {
114       stringstream name;
115       stringstream title;
116       name<<"hist_theta_"<<s1<<"_"<<s2;
117       title<<"Theta Missalignment for sectors "<<s1<<" and "<<s2;
118       fDthetaHistArray[s1*72+s2]=new TH1D(name.str().c_str(),title.str().c_str(),1024,-0.01,0.01); // +/- 10 mrad
119       ((TH1D*)fDthetaHistArray[s1*72+s2])->SetDirectory(0);
120     }
121     if (!fDyHistArray[s1*72+s2]) {
122       stringstream name;
123       stringstream title;
124       name<<"hist_y_"<<s1<<"_"<<s2;
125       title<<"Y Missalignment for sectors "<<s1<<" and "<<s2;
126       fDyHistArray[s1*72+s2]=new TH1D(name.str().c_str(),title.str().c_str(),1024,-0.3,0.3); // +/- 3 mm
127       ((TH1D*)fDyHistArray[s1*72+s2])->SetDirectory(0);
128     }
129     if (!fDzHistArray[s1*72+s2]) {
130       stringstream name;
131       stringstream title;
132       name<<"hist_z_"<<s1<<"_"<<s2;
133       title<<"Z Missalignment for sectors "<<s1<<" and "<<s2;
134       fDzHistArray[s1*72+s2]=new TH1D(name.str().c_str(),title.str().c_str(),1024,-0.3,0.3); // +/- 3 mm
135       ((TH1D*)fDzHistArray[s1*72+s2])->SetDirectory(0);
136     }
137     static_cast<TH1D*>(fDphiHistArray[s1*72+s2])
138       ->Fill(TMath::ASin(tp1.GetSnp())
139              -TMath::ASin(tp2.GetSnp()));
140     static_cast<TH1D*>(fDthetaHistArray[s1*72+s2])
141       ->Fill(TMath::ATan(tp1.GetTgl())
142              -TMath::ATan(tp2.GetTgl()));
143     static_cast<TH1D*>(fDyHistArray[s1*72+s2])
144       ->Fill(tp1.GetY()
145              -tp2.GetY());
146     static_cast<TH1D*>(fDzHistArray[s1*72+s2])
147       ->Fill(tp1.GetZ()
148              -tp2.GetZ());
149   }
150   return;
151
152
153
154   //
155   // Process function to fill fitters
156   //
157   Double_t t1[5],t2[5];
158   Double_t &x1=t1[0], &y1=t1[1], &z1=t1[2], &dydx1=t1[3], &dzdx1=t1[4];
159   Double_t /*&x2=t2[0],*/ &y2=t2[1], &z2=t2[2], &dydx2=t2[3]/*, &dzdx2=t2[4]*/;
160   x1   =tp1.GetX();
161   y1   =tp1.GetY();
162   z1   =tp1.GetZ();
163   Double_t snp1=tp1.GetSnp();
164   dydx1=snp1/TMath::Sqrt(1.-snp1*snp1);
165   Double_t tgl1=tp1.GetTgl();
166   // dz/dx = 1/(cos(theta)*cos(phi))
167   dzdx1=1./TMath::Sqrt((1.+tgl1*tgl1)*(1.-snp1*snp1));
168   //  x2  =tp1->GetX();
169   y2   =tp2.GetY();
170   z2   =tp2.GetZ();
171   Double_t snp2=tp2.GetSnp();
172   dydx2=snp2/TMath::Sqrt(1.-snp2*snp2);
173   Double_t tgl2=tp2.GetTgl();
174   dzdx1=1./TMath::Sqrt((1.+tgl2*tgl2)*(1.-snp2*snp2));
175
176   Process12(t1,t2,GetOrMakeFitter12(s1,s2));
177   Process9(t1,t2,GetOrMakeFitter12(s1,s2));
178   Process6(t1,t2,GetOrMakeFitter12(s1,s2));
179   ++fPoints[72*s1+s2];
180 }
181
182 void AliTPCcalibAlign::Process12(Double_t *t1,
183                                  Double_t *t2,
184                                  TLinearFitter *fitter) {
185   // x2    =  a00*x1 + a01*y1 + a02*z1 + a03
186   // y2    =  a10*x1 + a11*y1 + a12*z1 + a13
187   // z2    =  a20*x1 + a21*y1 + a22*z1 + a23
188   // dydx2 = (a10    + a11*dydx1 + a12*dzdx1)/(a00    + a01*dydx1 + a02*dzdx1)
189   // dzdx2 = (a20    + a21*dydx1 + a22*dzdx1)/(a00    + a01*dydx1 + a02*dzdx1)
190   //
191   //       a00  a01 a02  a03     p[0]   p[1]  p[2]  p[9]
192   //       a10  a11 a12  a13 ==> p[3]   p[4]  p[5]  p[10]
193   //       a20  a21 a22  a23     p[6]   p[7]  p[8]  p[11] 
194   Double_t &x1=t1[0], &y1=t1[1], &z1=t1[2], &dydx1=t1[3], &dzdx1=t1[4];
195   Double_t /*&x2=t2[0],*/ &y2=t2[1], &z2=t2[2], &dydx2=t2[3], &dzdx2=t2[4];
196
197   // TODO:
198   Double_t sy    = 1.;
199   Double_t sz    = 1.;
200   Double_t sdydx = 1.;
201   Double_t sdzdx = 1.;
202
203   Double_t p[12];
204   Double_t value;
205
206   // x2  =  a00*x1 + a01*y1 + a02*z1 + a03
207   // y2  =  a10*x1 + a11*y1 + a12*z1 + a13
208   // y2' =  a10*x1 + a11*y1 + a12*z1 + a13 + (a01*y1 + a02*z1 + a03)*dydx2
209   for (Int_t i=0; i<12;i++) p[i]=0.;
210   p[3+0] = x1;          // a10
211   p[3+1] = y1;          // a11
212   p[3+2] = z1;          // a12
213   p[9+1] = 1.;          // a13
214   p[0+1] = y1*dydx2;    // a01
215   p[0+2] = z1*dydx2;    // a02
216   p[9+0] = dydx2;       // a03
217   value  = y2;
218   fitter->AddPoint(p,value,sy);
219
220   // x2  =  a00*x1 + a01*y1 + a02*z1 + a03
221   // z2  =  a20*x1 + a21*y1 + a22*z1 + a23
222   // z2' =  a20*x1 + a21*y1 + a22*z1 + a23 + (a01*y1 + a02*z1 + a03)*dzdx2;
223   for (Int_t i=0; i<12;i++) p[i]=0.;
224   p[6+0] = x1;           // a20 
225   p[6+1] = y1;           // a21
226   p[6+2] = z1;           // a22
227   p[9+2] = 1.;           // a23
228   p[0+1] = y1*dzdx2;     // a01
229   p[0+2] = z1*dzdx2;     // a02
230   p[9+0] = dzdx2;        // a03
231   value  = z2;
232   fitter->AddPoint(p,value,sz);
233
234   // dydx2 = (a10 + a11*dydx1 + a12*dzdx1)/( a00 + a01*dydx1 + a02*dzdx1)
235   // (a10 + a11*dydx1 + a12*dzdx1) - (a00 + a01*dydx1 + a02*dzdx1)*dydx2 = 0
236   for (Int_t i=0; i<12;i++) p[i]=0.;
237   p[3+0] = 1.;           // a10
238   p[3+1] = dydx1;        // a11
239   p[3+2] = dzdx1;        // a12
240   p[0+0] = -dydx2;       // a00
241   p[0+1] = -dydx1*dydx2; // a01
242   p[0+2] = -dzdx1*dydx2; // a02
243   value  = 0.;
244   fitter->AddPoint(p,value,sdydx);
245
246   // dzdx2 = (a20 + a21*dydx1 + a22*dzdx1)/( a00 + a01*dydx1 + a02*dzdx1)
247   // (a20 + a21*dydx1 + a22*dzdx1) - (a00 + a01*dydx1 + a02*dzdx1)*dzdx2 = 0
248   for (Int_t i=0; i<12;i++) p[i]=0.;
249   p[6+0] = 1;            // a20
250   p[6+1] = dydx1;        // a21
251   p[6+2] = dzdx1;        // a22
252   p[0+0] = -dzdx2;       // a00
253   p[0+1] = -dydx1*dzdx2; // a01
254   p[0+2] = -dzdx1*dzdx2; // a02
255   value  = 0.;
256   fitter->AddPoint(p,value,sdzdx);
257 }
258
259 void AliTPCcalibAlign::Process9(Double_t *t1,
260                                 Double_t *t2,
261                                 TLinearFitter *fitter) {
262   // x2    =  a00*x1 + a01*y1 + a02*z1 + a03
263   // y2    =  a10*x1 + a11*y1 + a12*z1 + a13
264   // z2    =  a20*x1 + a21*y1 + a22*z1 + a23
265   // dydx2 = (a10    + a11*dydx1 + a12*dzdx1)/(a00    + a01*dydx1 + a02*dzdx1)
266   // dzdx2 = (a20    + a21*dydx1 + a22*dzdx1)/(a00    + a01*dydx1 + a02*dzdx1)
267   //
268   //       a00  a01 a02  a03     p[0]   p[1]  p[2]  p[9]
269   //       a10  a11 a12  a13 ==> p[3]   p[4]  p[5]  p[10]
270   //       a20  a21 a22  a23     p[6]   p[7]  p[8]  p[11] 
271   Double_t &x1=t1[0], &y1=t1[1], &z1=t1[2], &dydx1=t1[3], &dzdx1=t1[4];
272   Double_t /*&x2=t2[0],*/ &y2=t2[1], &z2=t2[2], &dydx2=t2[3], &dzdx2=t2[4];
273
274   // TODO:
275   Double_t sy    = 1.;
276   Double_t sz    = 1.;
277   Double_t sdydx = 1.;
278   Double_t sdzdx = 1.;
279
280   Double_t p[12];
281   Double_t value;
282
283   // x2  =  a00*x1 + a01*y1 + a02*z1 + a03
284   // y2  =  a10*x1 + a11*y1 + a12*z1 + a13
285   // y2' =  a10*x1 + a11*y1 + a12*z1 + a13 + (a01*y1 + a03)*dydx2
286   for (Int_t i=0; i<12;i++) p[i]=0.;
287   p[3+0] = x1;
288   p[3+1] = y1;
289   p[3+2] = z1;
290   p[9+1] = 1.;
291   p[0+1] = y1*dydx2;
292   p[0+2] = z1*dydx2;
293   p[9+0] = dydx2;
294   value  = y2;
295   fitter->AddPoint(p,value,sy);
296
297   // x2  =  a00*x1 + a01*y1 + a02*z1 + a03
298   // z2  =  a20*x1 + a21*y1 + a22*z1 + a23
299   // z2' =  a20*x1 + a21*y1 + a22*z1 + a23 + (a01*y1 + a03)*dzdx2;
300   for (Int_t i=0; i<12;i++) p[i]=0.;
301   p[6+0] = x1;
302   p[6+1] = y1;
303   p[6+2] = z1;
304   p[9+2] = 1.;
305   p[0+1] = y1*dzdx2;
306   p[0+2] = z1*dzdx2;
307   p[9+0] = dzdx2;
308   value  = z2;
309   fitter->AddPoint(p,value,sz);
310
311   // dydx2 = (a10 + a11*dydx1 + a12*dzdx1)/( a00 + a01*dydx1 + a02*dzdx1)
312   // (a10 + a11*dydx1 + a12*dzdx1) - (a00 + a01*dydx1 + a02*dzdx1)*dydx2 = 0
313   for (Int_t i=0; i<12;i++) p[i]=0.;
314   p[3+0] = 1.;
315   p[3+1] = dydx1;
316   p[3+2] = dzdx1;
317   p[0+0] = -dydx2;
318   p[0+1] = -dydx1*dydx2;
319   p[0+2] = -dzdx1*dydx2;
320   value  = 0.;
321   fitter->AddPoint(p,value,sdydx);
322
323   // dzdx2 = (a20 + a21*dydx1 + a22*dzdx1)/( a00 + a01*dydx1 + a02*dzdx1)
324   // (a20 + a21*dydx1 + a22*dzdx1) - (a00 + a01*dydx1 + a02*dzdx1)*dzdx2 = 0
325   for (Int_t i=0; i<12;i++) p[i]=0.;
326   p[6+0] = 1;
327   p[6+1] = dydx1;
328   p[6+2] = dzdx1;
329   p[0+0] = -dzdx2;
330   p[0+1] = -dydx1*dzdx2;
331   p[0+2] = -dzdx1*dzdx2;
332   value  = 0.;
333   fitter->AddPoint(p,value,sdzdx);
334 }
335
336 void AliTPCcalibAlign::Process6(Double_t *t1,
337                                 Double_t *t2,
338                                 TLinearFitter *fitter) {
339   // x2    =  1  *x1 +-a01*y1 + 0      +a03
340   // y2    =  a01*x1 + 1  *y1 + 0      +a13
341   // z2    =  a20*x1 + a21*y1 + 1  *z1 +a23
342   // dydx2 = (a10    + a11*dydx1 + a12*dzdx1)/(a00    + a01*dydx1 + a02*dzdx1)
343   // dzdx2 = (a20    + a21*dydx1 + a22*dzdx1)/(a00    + a01*dydx1 + a02*dzdx1)
344   //
345   //       1   -a01 0    a03     x     -p[0]  x     p[3]
346   //       a10  1   0    a13 ==> p[0]   x     x     p[4]
347   //       a20  a21 1    a23     p[1]   p[2]  x     p[5] 
348   Double_t &x1=t1[0], &y1=t1[1], &z1=t1[2], &dydx1=t1[3], &dzdx1=t1[4];
349   Double_t /*&x2=t2[0],*/ &y2=t2[1], &z2=t2[2], &dydx2=t2[3], &dzdx2=t2[4];
350
351   // TODO:
352   Double_t sy    = 1.;
353   Double_t sz    = 1.;
354   Double_t sdydx = 1.;
355   Double_t sdzdx = 1.;
356
357   Double_t p[12];
358   Double_t value;
359
360   // x2  =  1  *x1 +-a01*y1 + 0      +a03
361   // y2  =  a01*x1 + 1  *y1 + 0      +a13
362   // y2' =  a10*x1 + a11*y1 + a12*z1 + a13 + (a01*y1 + a02*z1 + a03)*dydx2
363   for (Int_t i=0; i<12;i++) p[i]=0.;
364   p[3+0] = x1;
365   p[3+1] = y1;
366   p[3+2] = z1;
367   p[9+1] = 1.;
368   p[0+1] = y1*dydx2;
369   p[0+2] = z1*dydx2;
370   p[9+0] = dydx2;
371   value  = y2;
372   fitter->AddPoint(p,value,sy);
373
374   // x2  =  1  *x1 +-a01*y1 + 0      + a03
375   // z2  =  a20*x1 + a21*y1 + 1  *z1 + a23
376   // z2' =  a20*x1 + a21*y1 + a22*z1 + a23 + (a01*y1 + a02*z1 +a03)*dzdx2;
377   for (Int_t i=0; i<12;i++) p[i]=0.;
378   p[6+0] = x1;
379   p[6+1] = y1;
380   p[6+2] = z1;
381   p[9+2] = 1.;
382   p[0+1] = y1*dzdx2;
383   p[0+2] = z1*dzdx2;
384   p[9+0] = dzdx2;
385   value  = z2;
386   fitter->AddPoint(p,value,sz);
387
388   // dydx2 = (a10 + a11*dydx1 + a12*dzdx1)/( a00 + a01*dydx1 + a02*dzdx1)
389   // (a10 + a11*dydx1 + a12*dzdx1) - (a00 + a01*dydx1 + a02*dzdx1)*dydx2 = 0
390   for (Int_t i=0; i<12;i++) p[i]=0.;
391   p[3+0] = 1.;
392   p[3+1] = dydx1;
393   p[3+2] = dzdx1;
394   p[0+0] = -dydx2;
395   p[0+1] = -dydx1*dydx2;
396   p[0+2] = -dzdx1*dydx2;
397   value  = 0.;
398   fitter->AddPoint(p,value,sdydx);
399
400   // dzdx2 = (a20 + a21*dydx1 + a22*dzdx1)/( a00 + a01*dydx1 + a02*dzdx1)
401   // (a20 + a21*dydx1 + a22*dzdx1) - (a00 + a01*dydx1 + a02*dzdx1)*dzdx2 = 0
402   for (Int_t i=0; i<12;i++) p[i]=0.;
403   p[6+0] = 1;
404   p[6+1] = dydx1;
405   p[6+2] = dzdx1;
406   p[0+0] = -dzdx2;
407   p[0+1] = -dydx1*dzdx2;
408   p[0+2] = -dzdx1*dzdx2;
409   value  = 0.;
410   fitter->AddPoint(p,value,sdzdx);
411 }
412
413 void AliTPCcalibAlign::Eval() {
414   TLinearFitter *f;
415   for (Int_t s1=0;s1<72;++s1)
416     for (Int_t s2=0;s2<72;++s2)
417       if ((f=GetFitter12(s1,s2))&&fPoints[72*s1+s2]>12) {
418         //      cerr<<s1<<","<<s2<<": "<<fPoints[72*s1+s2]<<endl;
419         if (!f->Eval()) {
420           cerr<<"Evaluation failed for "<<s1<<","<<s2<<endl;
421         }
422       }
423   /*
424                     
425   fitter->Eval();
426   fitter->Eval();
427   chi212 = align->GetChisquare()/(4.*entries);
428
429   TMatrixD mat(13,13);
430   TVectorD par(13);
431   align->GetParameters(par);
432   align->GetCovarianceMatrix(mat);
433
434   //
435   //
436   for (Int_t i=0; i<12;i++){
437     palign12(i)= par(i+1);
438     for (Int_t j=0; j<12;j++){
439       pcovar12(i,j)   = mat(i+1,j+1);
440       pcovar12(i,j) *= chi212;
441     }
442   }
443   //
444   for (Int_t i=0; i<12;i++){
445     psigma12(i)  = TMath::Sqrt(pcovar12(i,i));
446     palignR12(i) = palign12(i)/TMath::Sqrt(pcovar12(i,i));
447     for (Int_t j=0; j<12;j++){
448       pcovarN12(i,j) = pcovar12(i,j)/TMath::Sqrt(pcovar12(i,i)*pcovar12(j,j));
449     }
450   }
451   */
452 }
453
454 Bool_t AliTPCcalibAlign::GetTransformation12(Int_t s1,Int_t s2,TMatrixD &a) {
455   if (!GetFitter12(s1,s2))
456     return false;
457   else {
458     TVectorD p(12);
459     cerr<<"foo"<<endl;
460     GetFitter12(s1,s2)->GetParameters(p);
461     cerr<<"bar"<<endl;
462     a.ResizeTo(4,4);
463     a[0][0]=p[0];
464     a[0][1]=p[1];
465     a[0][2]=p[2];
466     a[0][3]=p[9];
467     a[1][0]=p[3];
468     a[1][1]=p[4];
469     a[1][2]=p[5];
470     a[1][3]=p[10];
471     a[2][0]=p[6];
472     a[2][1]=p[7];
473     a[2][2]=p[8];
474     a[2][3]=p[11];
475     a[3][0]=0.;
476     a[3][1]=0.;
477     a[3][2]=0.;
478     a[3][3]=1.;
479     return true;
480   } 
481 }
482
483 Bool_t AliTPCcalibAlign::GetTransformation9(Int_t s1,Int_t s2,TMatrixD &a) {
484   if (!GetFitter9(s1,s2))
485     return false;
486   else {
487     TVectorD p(9);
488     GetFitter9(s1,s2)->GetParameters(p);
489     a.ResizeTo(4,4);
490     a[0][0]=p[0];
491     a[0][1]=p[1];
492     a[0][2]=p[2];
493     a[0][3]=p[9];
494     a[1][0]=p[3];
495     a[1][1]=p[4];
496     a[1][2]=p[5];
497     a[1][3]=p[10];
498     a[2][0]=p[6];
499     a[2][1]=p[7];
500     a[2][2]=p[8];
501     a[2][3]=p[11];
502     a[3][0]=0.;
503     a[3][1]=0.;
504     a[3][2]=0.;
505     a[3][3]=1.;
506     return true;
507   } 
508 }
509
510 Bool_t AliTPCcalibAlign::GetTransformation6(Int_t s1,Int_t s2,TMatrixD &a) {
511   if (!GetFitter6(s1,s2))
512     return false;
513   else {
514     TVectorD p(6);
515     cerr<<"foo"<<endl;
516     GetFitter6(s1,s2)->GetParameters(p);
517     cerr<<"bar"<<endl;
518     a.ResizeTo(4,4);
519     a[0][0]=p[0];
520     a[0][1]=p[1];
521     a[0][2]=p[2];
522     a[0][3]=p[9];
523     a[1][0]=p[3];
524     a[1][1]=p[4];
525     a[1][2]=p[5];
526     a[1][3]=p[10];
527     a[2][0]=p[6];
528     a[2][1]=p[7];
529     a[2][2]=p[8];
530     a[2][3]=p[11];
531     a[3][0]=0.;
532     a[3][1]=0.;
533     a[3][2]=0.;
534     a[3][3]=1.;
535     return true;
536   } 
537 }