]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibAlign.cxx
Bug fix - Teminate - as virtula function
[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 #include "TVectorD.h"
37 #include "TTreeStream.h"
38 #include "TFile.h"
39
40 #include <iostream>
41 #include <sstream>
42 using namespace std;
43
44 ClassImp(AliTPCcalibAlign)
45
46 AliTPCcalibAlign::AliTPCcalibAlign()
47   :  fDphiHistArray(72*72),
48      fDthetaHistArray(72*72),
49      fDyHistArray(72*72),
50      fDzHistArray(72*72),
51      fFitterArray12(72*72),
52      fFitterArray9(72*72),
53      fFitterArray6(72*72)
54 {
55   //
56   // Constructor
57   //
58   for (Int_t i=0;i<72*72;++i) {
59     fPoints[i]=0;
60   }
61 }
62
63 AliTPCcalibAlign::AliTPCcalibAlign(const Text_t *name, const Text_t *title)
64   :AliTPCcalibBase(),  
65    fDphiHistArray(72*72),
66    fDthetaHistArray(72*72),
67    fDyHistArray(72*72),
68    fDzHistArray(72*72),
69    fFitterArray12(72*72),
70    fFitterArray9(72*72),
71    fFitterArray6(72*72)
72 {
73   //
74   // Constructor
75   //
76   SetName(name);
77   SetTitle(title);
78   for (Int_t i=0;i<72*72;++i) {
79     fPoints[i]=0;
80   }
81 }
82
83 AliTPCcalibAlign::~AliTPCcalibAlign() {
84   //
85   // destructor
86   //
87 }
88
89 void AliTPCcalibAlign::Process(AliTPCseed *seed) {
90   TObjArray tracklets=
91     AliTPCTracklet::CreateTracklets(seed,AliTPCTracklet::kKalman,
92                                     kFALSE,20,2);
93  //  TObjArray trackletsL=
94 //     AliTPCTracklet::CreateTracklets(seed,AliTPCTracklet::kLinear,
95 //                                  kFALSE,20,2);
96 //   TObjArray trackletsQ=
97 //     AliTPCTracklet::CreateTracklets(seed,AliTPCTracklet::kQuadratic,
98 //                                  kFALSE,20,2);
99 //   TObjArray trackletsR=
100 //     AliTPCTracklet::CreateTracklets(seed,AliTPCTracklet::kRiemann,
101 //                                  kFALSE,20,2);
102   tracklets.SetOwner();
103  //  trackletsL.SetOwner();
104 //   trackletsQ.SetOwner();
105 //   trackletsR.SetOwner();
106   if (tracklets.GetEntries()==2) {
107     AliTPCTracklet *t1=static_cast<AliTPCTracklet*>(tracklets[0]);
108     AliTPCTracklet *t2=static_cast<AliTPCTracklet*>(tracklets[1]);
109     if (t1->GetSector()>t2->GetSector()) {
110       AliTPCTracklet* tmp=t1;
111       t1=t2;
112       t2=tmp;
113     }
114     AliExternalTrackParam *common1=0,*common2=0;
115     if (AliTPCTracklet::PropagateToMeanX(*t1,*t2,common1,common2))
116       ProcessTracklets(*common1,*common2,t1->GetSector(),t2->GetSector());
117     delete common1;
118     delete common2;
119   }
120   
121 }
122
123 void AliTPCcalibAlign::Analyze(){
124   //
125   // Analyze function 
126   //
127   EvalFitters();
128 }
129
130
131 void AliTPCcalibAlign::Terminate(){
132   //
133   // Terminate function
134   // call base terminate + Eval of fitters
135   //
136   if (GetDebugLevel()>0) Info("AliTPCcalibAlign","Teminate");
137   EvalFitters();
138   AliTPCcalibBase::Terminate();
139 }
140
141
142
143
144 void AliTPCcalibAlign::ProcessTracklets(const AliExternalTrackParam &tp1,
145                                         const AliExternalTrackParam &tp2,
146                                         Int_t s1,Int_t s2) {
147
148
149   if (s2-s1==36) {//only inner-outer
150     if (!fDphiHistArray[GetIndex(s1,s2)]) {
151       stringstream name;
152       stringstream title;
153       name<<"hist_phi_"<<s1<<"_"<<s2;
154       title<<"Phi Missalignment for sectors "<<s1<<" and "<<s2;
155       fDphiHistArray[GetIndex(s1,s2)]=new TH1D(name.str().c_str(),title.str().c_str(),512,-0.01,0.01); // +/- 10 mrad
156       ((TH1D*)fDphiHistArray[GetIndex(s1,s2)])->SetDirectory(0);
157     }
158     if (!fDthetaHistArray[GetIndex(s1,s2)]) {
159       stringstream name;
160       stringstream title;
161       name<<"hist_theta_"<<s1<<"_"<<s2;
162       title<<"Theta Missalignment for sectors "<<s1<<" and "<<s2;
163       fDthetaHistArray[GetIndex(s1,s2)]=new TH1D(name.str().c_str(),title.str().c_str(),512,-0.01,0.01); // +/- 10 mrad
164       ((TH1D*)fDthetaHistArray[GetIndex(s1,s2)])->SetDirectory(0);
165     }
166     if (!fDyHistArray[GetIndex(s1,s2)]) {
167       stringstream name;
168       stringstream title;
169       name<<"hist_y_"<<s1<<"_"<<s2;
170       title<<"Y Missalignment for sectors "<<s1<<" and "<<s2;
171       fDyHistArray[GetIndex(s1,s2)]=new TH1D(name.str().c_str(),title.str().c_str(),512,-0.3,0.3); // +/- 3 mm
172       ((TH1D*)fDyHistArray[GetIndex(s1,s2)])->SetDirectory(0);
173     }
174     if (!fDzHistArray[GetIndex(s1,s2)]) {
175       stringstream name;
176       stringstream title;
177       name<<"hist_z_"<<s1<<"_"<<s2;
178       title<<"Z Missalignment for sectors "<<s1<<" and "<<s2;
179       fDzHistArray[GetIndex(s1,s2)]=new TH1D(name.str().c_str(),title.str().c_str(),512,-0.3,0.3); // +/- 3 mm
180       ((TH1D*)fDzHistArray[GetIndex(s1,s2)])->SetDirectory(0);
181     }
182     static_cast<TH1D*>(fDphiHistArray[GetIndex(s1,s2)])
183       ->Fill(TMath::ASin(tp1.GetSnp())
184              -TMath::ASin(tp2.GetSnp()));
185     static_cast<TH1D*>(fDthetaHistArray[GetIndex(s1,s2)])
186       ->Fill(TMath::ATan(tp1.GetTgl())
187              -TMath::ATan(tp2.GetTgl()));
188     static_cast<TH1D*>(fDyHistArray[GetIndex(s1,s2)])
189       ->Fill(tp1.GetY()
190              -tp2.GetY());
191     static_cast<TH1D*>(fDzHistArray[GetIndex(s1,s2)])
192       ->Fill(tp1.GetZ()
193              -tp2.GetZ());
194   }
195   
196
197
198
199   //
200   // Process function to fill fitters
201   //
202   Double_t t1[5],t2[5];
203   Double_t &x1=t1[0], &y1=t1[1], &z1=t1[2], &dydx1=t1[3], &dzdx1=t1[4];
204   Double_t &x2=t2[0], &y2=t2[1], &z2=t2[2], &dydx2=t2[3], &dzdx2=t2[4];
205   x1   =tp1.GetX();
206   y1   =tp1.GetY();
207   z1   =tp1.GetZ();
208   Double_t snp1=tp1.GetSnp();
209   dydx1=snp1/TMath::Sqrt(1.-snp1*snp1);
210   Double_t tgl1=tp1.GetTgl();
211   // dz/dx = 1/(cos(theta)*cos(phi))
212   dzdx1=tgl1/TMath::Sqrt(1.-snp1*snp1);
213   x2   =tp2.GetX();
214   y2   =tp2.GetY();
215   z2   =tp2.GetZ();
216   Double_t snp2=tp2.GetSnp();
217   dydx2=snp2/TMath::Sqrt(1.-snp2*snp2);
218   Double_t tgl2=tp2.GetTgl();
219   dzdx2=tgl2/TMath::Sqrt(1.-snp2*snp2);
220
221   //
222   //
223   //
224   if (fStreamLevel>1){
225     TTreeSRedirector *cstream = GetDebugStreamer();
226     if (cstream){
227       static TVectorD vec1(5);
228       static TVectorD vec2(5);
229       vec1.SetElements(t1);
230       vec2.SetElements(t2);
231       AliExternalTrackParam *p1 = &((AliExternalTrackParam&)tp1);
232       AliExternalTrackParam *p2 = &((AliExternalTrackParam&)tp2);
233       (*cstream)<<"Tracklet"<<
234         "tp1.="<<p1<<
235         "tp2.="<<p2<<
236         "v1.="<<&vec1<<
237         "v2.="<<&vec2<<
238         "s1="<<s1<<
239         "s2="<<s2<<
240         "\n";
241     }
242   }
243   Process12(t1,t2,GetOrMakeFitter12(s1,s2));
244   Process9(t1,t2,GetOrMakeFitter9(s1,s2));
245   Process6(t1,t2,GetOrMakeFitter6(s1,s2));
246   ++fPoints[72*s1+s2];
247 }
248
249 void AliTPCcalibAlign::Process12(const Double_t *t1,
250                                  const Double_t *t2,
251                                  TLinearFitter *fitter) {
252   // x2    =  a00*x1 + a01*y1 + a02*z1 + a03
253   // y2    =  a10*x1 + a11*y1 + a12*z1 + a13
254   // z2    =  a20*x1 + a21*y1 + a22*z1 + a23
255   // dydx2 = (a10    + a11*dydx1 + a12*dzdx1)/(a00    + a01*dydx1 + a02*dzdx1)
256   // dzdx2 = (a20    + a21*dydx1 + a22*dzdx1)/(a00    + a01*dydx1 + a02*dzdx1)
257   //
258   //       a00  a01 a02  a03     p[0]   p[1]  p[2]  p[9]
259   //       a10  a11 a12  a13 ==> p[3]   p[4]  p[5]  p[10]
260   //       a20  a21 a22  a23     p[6]   p[7]  p[8]  p[11] 
261   const Double_t &x1=t1[0], &y1=t1[1], &z1=t1[2], &dydx1=t1[3], &dzdx1=t1[4];
262   const Double_t /*&x2=t2[0],*/ &y2=t2[1], &z2=t2[2], &dydx2=t2[3], &dzdx2=t2[4];
263
264   // TODO:
265   Double_t sy    = 0.1;
266   Double_t sz    = 0.1;
267   Double_t sdydx = 0.001;
268   Double_t sdzdx = 0.001;
269
270   Double_t p[12];
271   Double_t value;
272
273   // x2  =  a00*x1 + a01*y1 + a02*z1 + a03
274   // y2  =  a10*x1 + a11*y1 + a12*z1 + a13
275   // y2' =  a10*x1 + a11*y1 + a12*z1 + a13 + (a01*y1 + a02*z1 + a03)*dydx2
276   for (Int_t i=0; i<12;i++) p[i]=0.;
277   p[3+0] = x1;          // a10
278   p[3+1] = y1;          // a11
279   p[3+2] = z1;          // a12
280   p[9+1] = 1.;          // a13
281   p[0+1] = y1*dydx2;    // a01
282   p[0+2] = z1*dydx2;    // a02
283   p[9+0] = dydx2;       // a03
284   value  = y2;
285   fitter->AddPoint(p,value,sy);
286
287   // x2  =  a00*x1 + a01*y1 + a02*z1 + a03
288   // z2  =  a20*x1 + a21*y1 + a22*z1 + a23
289   // z2' =  a20*x1 + a21*y1 + a22*z1 + a23 + (a01*y1 + a02*z1 + a03)*dzdx2;
290   for (Int_t i=0; i<12;i++) p[i]=0.;
291   p[6+0] = x1;           // a20 
292   p[6+1] = y1;           // a21
293   p[6+2] = z1;           // a22
294   p[9+2] = 1.;           // a23
295   p[0+1] = y1*dzdx2;     // a01
296   p[0+2] = z1*dzdx2;     // a02
297   p[9+0] = dzdx2;        // a03
298   value  = z2;
299   fitter->AddPoint(p,value,sz);
300
301   // dydx2 = (a10 + a11*dydx1 + a12*dzdx1)/( a00 + a01*dydx1 + a02*dzdx1)
302   // (a10 + a11*dydx1 + a12*dzdx1) - (a00 + a01*dydx1 + a02*dzdx1)*dydx2 = 0
303   for (Int_t i=0; i<12;i++) p[i]=0.;
304   p[3+0] = 1.;           // a10
305   p[3+1] = dydx1;        // a11
306   p[3+2] = dzdx1;        // a12
307   p[0+0] = -dydx2;       // a00
308   p[0+1] = -dydx1*dydx2; // a01
309   p[0+2] = -dzdx1*dydx2; // a02
310   value  = 0.;
311   fitter->AddPoint(p,value,sdydx);
312
313   // dzdx2 = (a20 + a21*dydx1 + a22*dzdx1)/( a00 + a01*dydx1 + a02*dzdx1)
314   // (a20 + a21*dydx1 + a22*dzdx1) - (a00 + a01*dydx1 + a02*dzdx1)*dzdx2 = 0
315   for (Int_t i=0; i<12;i++) p[i]=0.;
316   p[6+0] = 1;            // a20
317   p[6+1] = dydx1;        // a21
318   p[6+2] = dzdx1;        // a22
319   p[0+0] = -dzdx2;       // a00
320   p[0+1] = -dydx1*dzdx2; // a01
321   p[0+2] = -dzdx1*dzdx2; // a02
322   value  = 0.;
323   fitter->AddPoint(p,value,sdzdx);
324 }
325
326 void AliTPCcalibAlign::Process9(Double_t *t1,
327                                 Double_t *t2,
328                                 TLinearFitter *fitter) {
329   // x2    =  a00*x1 + a01*y1 + a02*z1 + a03
330   // y2    =  a10*x1 + a11*y1 + a12*z1 + a13
331   // z2    =  a20*x1 + a21*y1 + a22*z1 + a23
332   // dydx2 = (a10    + a11*dydx1 + a12*dzdx1)/(a00    + a01*dydx1 + a02*dzdx1)
333   // dzdx2 = (a20    + a21*dydx1 + a22*dzdx1)/(a00    + a01*dydx1 + a02*dzdx1)
334   //
335   //       a00  a01 a02  a03     p[0]   p[1]  p[2]  p[9]
336   //       a10  a11 a12  a13 ==> p[3]   p[4]  p[5]  p[10]
337   //       a20  a21 a22  a23     p[6]   p[7]  p[8]  p[11] 
338   Double_t &x1=t1[0], &y1=t1[1], &z1=t1[2], &dydx1=t1[3], &dzdx1=t1[4];
339   Double_t /*&x2=t2[0],*/ &y2=t2[1], &z2=t2[2], &dydx2=t2[3], &dzdx2=t2[4];
340
341   // TODO:
342   Double_t sy    = 1.;
343   Double_t sz    = 1.;
344   Double_t sdydx = 1.;
345   Double_t sdzdx = 1.;
346
347   Double_t p[12];
348   Double_t value;
349
350   // x2  =  a00*x1 + a01*y1 + a02*z1 + a03
351   // y2  =  a10*x1 + a11*y1 + a12*z1 + a13
352   // y2' =  a10*x1 + a11*y1 + a12*z1 + a13 + (a01*y1 + a03)*dydx2
353   for (Int_t i=0; i<12;i++) p[i]=0.;
354   p[3+0] = x1;
355   p[3+1] = y1;
356   p[3+2] = z1;
357   p[9+1] = 1.;
358   p[0+1] = y1*dydx2;
359   p[0+2] = z1*dydx2;
360   p[9+0] = dydx2;
361   value  = y2;
362   fitter->AddPoint(p,value,sy);
363
364   // x2  =  a00*x1 + a01*y1 + a02*z1 + a03
365   // z2  =  a20*x1 + a21*y1 + a22*z1 + a23
366   // z2' =  a20*x1 + a21*y1 + a22*z1 + a23 + (a01*y1 + a03)*dzdx2;
367   for (Int_t i=0; i<12;i++) p[i]=0.;
368   p[6+0] = x1;
369   p[6+1] = y1;
370   p[6+2] = z1;
371   p[9+2] = 1.;
372   p[0+1] = y1*dzdx2;
373   p[0+2] = z1*dzdx2;
374   p[9+0] = dzdx2;
375   value  = z2;
376   fitter->AddPoint(p,value,sz);
377
378   // dydx2 = (a10 + a11*dydx1 + a12*dzdx1)/( a00 + a01*dydx1 + a02*dzdx1)
379   // (a10 + a11*dydx1 + a12*dzdx1) - (a00 + a01*dydx1 + a02*dzdx1)*dydx2 = 0
380   for (Int_t i=0; i<12;i++) p[i]=0.;
381   p[3+0] = 1.;
382   p[3+1] = dydx1;
383   p[3+2] = dzdx1;
384   p[0+0] = -dydx2;
385   p[0+1] = -dydx1*dydx2;
386   p[0+2] = -dzdx1*dydx2;
387   value  = 0.;
388   fitter->AddPoint(p,value,sdydx);
389
390   // dzdx2 = (a20 + a21*dydx1 + a22*dzdx1)/( a00 + a01*dydx1 + a02*dzdx1)
391   // (a20 + a21*dydx1 + a22*dzdx1) - (a00 + a01*dydx1 + a02*dzdx1)*dzdx2 = 0
392   for (Int_t i=0; i<12;i++) p[i]=0.;
393   p[6+0] = 1;
394   p[6+1] = dydx1;
395   p[6+2] = dzdx1;
396   p[0+0] = -dzdx2;
397   p[0+1] = -dydx1*dzdx2;
398   p[0+2] = -dzdx1*dzdx2;
399   value  = 0.;
400   fitter->AddPoint(p,value,sdzdx);
401 }
402
403 void AliTPCcalibAlign::Process6(Double_t *t1,
404                                 Double_t *t2,
405                                 TLinearFitter *fitter) {
406   // x2    =  1  *x1 +-a01*y1 + 0      +a03
407   // y2    =  a01*x1 + 1  *y1 + 0      +a13
408   // z2    =  a20*x1 + a21*y1 + 1  *z1 +a23
409   // dydx2 = (a10    + a11*dydx1 + a12*dzdx1)/(a00    + a01*dydx1 + a02*dzdx1)
410   // dzdx2 = (a20    + a21*dydx1 + a22*dzdx1)/(a00    + a01*dydx1 + a02*dzdx1)
411   //
412   //       1   -a01 0    a03     x     -p[0]  x     p[3]
413   //       a10  1   0    a13 ==> p[0]   x     x     p[4]
414   //       a20  a21 1    a23     p[1]   p[2]  x     p[5] 
415   Double_t &x1=t1[0], &y1=t1[1], &z1=t1[2], &dydx1=t1[3], &dzdx1=t1[4];
416   Double_t /*&x2=t2[0],*/ &y2=t2[1], &z2=t2[2], &dydx2=t2[3], &dzdx2=t2[4];
417
418   // TODO:
419   Double_t sy    = 0.1;
420   Double_t sz    = 0.1;
421   Double_t sdydx = 0.001;
422   Double_t sdzdx = 0.001;
423
424   Double_t p[12];
425   Double_t value;
426
427   // x2  =  1  *x1 +-a01*y1 + 0      +a03
428   // y2  =  a01*x1 + 1  *y1 + 0      +a13
429   // y2' =  a10*x1 + a11*y1 + a12*z1 + a13 + (a01*y1 + a02*z1 + a03)*dydx2
430   for (Int_t i=0; i<12;i++) p[i]=0.;
431   p[3+0] = x1;
432   p[3+1] = y1;
433   p[3+2] = z1;
434   p[9+1] = 1.;
435   p[0+1] = y1*dydx2;
436   p[0+2] = z1*dydx2;
437   p[9+0] = dydx2;
438   value  = y2;
439   fitter->AddPoint(p,value,sy);
440
441   // x2  =  1  *x1 +-a01*y1 + 0      + a03
442   // z2  =  a20*x1 + a21*y1 + 1  *z1 + a23
443   // z2' =  a20*x1 + a21*y1 + a22*z1 + a23 + (a01*y1 + a02*z1 +a03)*dzdx2;
444   for (Int_t i=0; i<12;i++) p[i]=0.;
445   p[6+0] = x1;
446   p[6+1] = y1;
447   p[6+2] = z1;
448   p[9+2] = 1.;
449   p[0+1] = y1*dzdx2;
450   p[0+2] = z1*dzdx2;
451   p[9+0] = dzdx2;
452   value  = z2;
453   fitter->AddPoint(p,value,sz);
454
455   // dydx2 = (a10 + a11*dydx1 + a12*dzdx1)/( a00 + a01*dydx1 + a02*dzdx1)
456   // (a10 + a11*dydx1 + a12*dzdx1) - (a00 + a01*dydx1 + a02*dzdx1)*dydx2 = 0
457   for (Int_t i=0; i<12;i++) p[i]=0.;
458   p[3+0] = 1.;
459   p[3+1] = dydx1;
460   p[3+2] = dzdx1;
461   p[0+0] = -dydx2;
462   p[0+1] = -dydx1*dydx2;
463   p[0+2] = -dzdx1*dydx2;
464   value  = 0.;
465   fitter->AddPoint(p,value,sdydx);
466
467   // dzdx2 = (a20 + a21*dydx1 + a22*dzdx1)/( a00 + a01*dydx1 + a02*dzdx1)
468   // (a20 + a21*dydx1 + a22*dzdx1) - (a00 + a01*dydx1 + a02*dzdx1)*dzdx2 = 0
469   for (Int_t i=0; i<12;i++) p[i]=0.;
470   p[6+0] = 1;
471   p[6+1] = dydx1;
472   p[6+2] = dzdx1;
473   p[0+0] = -dzdx2;
474   p[0+1] = -dydx1*dzdx2;
475   p[0+2] = -dzdx1*dzdx2;
476   value  = 0.;
477   fitter->AddPoint(p,value,sdzdx);
478 }
479
480
481
482
483 void AliTPCcalibAlign::EvalFitters() {
484   //
485   // Analyze function 
486   // 
487   // Perform the fitting using linear fitters
488   //
489   Int_t kMinPoints =50;
490   TLinearFitter *f;
491   TFile fff("alignDebug.root","recreate");
492   for (Int_t s1=0;s1<72;++s1)
493     for (Int_t s2=0;s2<72;++s2){
494       if ((f=GetFitter12(s1,s2))&&fPoints[72*s1+s2]>kMinPoints) {
495         //      cerr<<s1<<","<<s2<<": "<<fPoints[72*s1+s2]<<endl;
496         if (!f->Eval()) {
497           cerr<<"Evaluation failed for "<<s1<<","<<s2<<endl;
498           f->Write(Form("f12_%d_%d",s1,s2));
499         }else{
500           f->Write(Form("f12_%d_%d",s1,s2));
501         }
502       }
503       if ((f=GetFitter9(s1,s2))&&fPoints[72*s1+s2]>kMinPoints) {
504         //      cerr<<s1<<","<<s2<<": "<<fPoints[72*s1+s2]<<endl;
505         if (!f->Eval()) {
506           cerr<<"Evaluation failed for "<<s1<<","<<s2<<endl;
507         }else{
508           f->Write(Form("f9_%d_%d",s1,s2));
509         }
510       }
511       if ((f=GetFitter6(s1,s2))&&fPoints[72*s1+s2]>kMinPoints) {
512         //      cerr<<s1<<","<<s2<<": "<<fPoints[72*s1+s2]<<endl;
513         if (!f->Eval()) {
514           cerr<<"Evaluation failed for "<<s1<<","<<s2<<endl;
515         }else{
516           f->Write(Form("f6_%d_%d",s1,s2));
517         }
518       }
519     }
520   /*
521                     
522   fitter->Eval();
523   fitter->Eval();
524   chi212 = align->GetChisquare()/(4.*entries);
525
526   TMatrixD mat(13,13);
527   TVectorD par(13);
528   align->GetParameters(par);
529   align->GetCovarianceMatrix(mat);
530
531   //
532   //
533   for (Int_t i=0; i<12;i++){
534     palign12(i)= par(i+1);
535     for (Int_t j=0; j<12;j++){
536       pcovar12(i,j)   = mat(i+1,j+1);
537       pcovar12(i,j) *= chi212;
538     }
539   }
540   //
541   for (Int_t i=0; i<12;i++){
542     psigma12(i)  = TMath::Sqrt(pcovar12(i,i));
543     palignR12(i) = palign12(i)/TMath::Sqrt(pcovar12(i,i));
544     for (Int_t j=0; j<12;j++){
545       pcovarN12(i,j) = pcovar12(i,j)/TMath::Sqrt(pcovar12(i,i)*pcovar12(j,j));
546     }
547   }
548   */
549 }
550
551 Bool_t AliTPCcalibAlign::GetTransformation12(Int_t s1,Int_t s2,TMatrixD &a) {
552   if (!GetFitter12(s1,s2))
553     return false;
554   else {
555     TVectorD p(12);
556     cerr<<"foo"<<endl;
557     GetFitter12(s1,s2)->GetParameters(p);
558     cerr<<"bar"<<endl;
559     a.ResizeTo(4,4);
560     a[0][0]=p[0];
561     a[0][1]=p[1];
562     a[0][2]=p[2];
563     a[0][3]=p[9];
564     a[1][0]=p[3];
565     a[1][1]=p[4];
566     a[1][2]=p[5];
567     a[1][3]=p[10];
568     a[2][0]=p[6];
569     a[2][1]=p[7];
570     a[2][2]=p[8];
571     a[2][3]=p[11];
572     a[3][0]=0.;
573     a[3][1]=0.;
574     a[3][2]=0.;
575     a[3][3]=1.;
576     return true;
577   } 
578 }
579
580 Bool_t AliTPCcalibAlign::GetTransformation9(Int_t s1,Int_t s2,TMatrixD &a) {
581   if (!GetFitter9(s1,s2))
582     return false;
583   else {
584     TVectorD p(9);
585     GetFitter9(s1,s2)->GetParameters(p);
586     a.ResizeTo(4,4);
587     a[0][0]=p[0];
588     a[0][1]=p[1];
589     a[0][2]=p[2];
590     a[0][3]=p[9];
591     a[1][0]=p[3];
592     a[1][1]=p[4];
593     a[1][2]=p[5];
594     a[1][3]=p[10];
595     a[2][0]=p[6];
596     a[2][1]=p[7];
597     a[2][2]=p[8];
598     a[2][3]=p[11];
599     a[3][0]=0.;
600     a[3][1]=0.;
601     a[3][2]=0.;
602     a[3][3]=1.;
603     return true;
604   } 
605 }
606
607 Bool_t AliTPCcalibAlign::GetTransformation6(Int_t s1,Int_t s2,TMatrixD &a) {
608   if (!GetFitter6(s1,s2))
609     return false;
610   else {
611     TVectorD p(6);
612     cerr<<"foo"<<endl;
613     GetFitter6(s1,s2)->GetParameters(p);
614     cerr<<"bar"<<endl;
615     a.ResizeTo(4,4);
616     a[0][0]=p[0];
617     a[0][1]=p[1];
618     a[0][2]=p[2];
619     a[0][3]=p[9];
620     a[1][0]=p[3];
621     a[1][1]=p[4];
622     a[1][2]=p[5];
623     a[1][3]=p[10];
624     a[2][0]=p[6];
625     a[2][1]=p[7];
626     a[2][2]=p[8];
627     a[2][3]=p[11];
628     a[3][0]=0.;
629     a[3][1]=0.;
630     a[3][2]=0.;
631     a[3][3]=1.;
632     return true;
633   } 
634 }