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