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