1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 ///////////////////////////////////////////////////////////////////////////////
19 // Class to make a internal alignemnt of TPC chambers //
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
26 //// Only the 12 parameter version works so far...
31 #include "TLinearFitter.h"
32 #include "AliTPCcalibAlign.h"
33 #include "AliExternalTrackParam.h"
34 #include "AliTPCTracklet.h"
37 #include "TTreeStream.h"
44 ClassImp(AliTPCcalibAlign)
46 AliTPCcalibAlign::AliTPCcalibAlign()
47 : fDphiHistArray(72*72),
48 fDthetaHistArray(72*72),
51 fFitterArray12(72*72),
58 for (Int_t i=0;i<72*72;++i) {
63 AliTPCcalibAlign::AliTPCcalibAlign(const Text_t *name, const Text_t *title)
65 fDphiHistArray(72*72),
66 fDthetaHistArray(72*72),
69 fFitterArray12(72*72),
78 for (Int_t i=0;i<72*72;++i) {
83 AliTPCcalibAlign::~AliTPCcalibAlign() {
89 void AliTPCcalibAlign::Process(AliTPCseed *seed) {
91 AliTPCTracklet::CreateTracklets(seed,AliTPCTracklet::kKalman,
93 // TObjArray trackletsL=
94 // AliTPCTracklet::CreateTracklets(seed,AliTPCTracklet::kLinear,
96 // TObjArray trackletsQ=
97 // AliTPCTracklet::CreateTracklets(seed,AliTPCTracklet::kQuadratic,
99 // TObjArray trackletsR=
100 // AliTPCTracklet::CreateTracklets(seed,AliTPCTracklet::kRiemann,
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;
114 AliExternalTrackParam *common1=0,*common2=0;
115 if (AliTPCTracklet::PropagateToMeanX(*t1,*t2,common1,common2))
116 ProcessTracklets(*common1,*common2,t1->GetSector(),t2->GetSector());
123 void AliTPCcalibAlign::Analyze(){
131 void AliTPCcalibAlign::Terminate(){
133 // Terminate function
134 // call base terminate + Eval of fitters
136 if (GetDebugLevel()>0) Info("AliTPCcalibAlign","Teminate");
138 AliTPCcalibBase::Terminate();
144 void AliTPCcalibAlign::ProcessTracklets(const AliExternalTrackParam &tp1,
145 const AliExternalTrackParam &tp2,
151 FillHisto(tp1,tp2,s1,s2);
153 // Process function to fill fitters
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];
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);
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);
177 TTreeSRedirector *cstream = GetDebugStreamer();
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"<<
195 Process12(t1,t2,GetOrMakeFitter12(s1,s2));
196 Process9(t1,t2,GetOrMakeFitter9(s1,s2));
197 Process6(t1,t2,GetOrMakeFitter6(s1,s2));
201 void AliTPCcalibAlign::Process12(const Double_t *t1,
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)
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];
219 Double_t sdydx = 0.001;
220 Double_t sdzdx = 0.001;
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.;
233 p[0+1] = y1*dydx2; // a01
234 p[0+2] = z1*dydx2; // a02
235 p[9+0] = dydx2; // a03
237 fitter->AddPoint(p,value,sy);
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.;
247 p[0+1] = y1*dzdx2; // a01
248 p[0+2] = z1*dzdx2; // a02
249 p[9+0] = dzdx2; // a03
251 fitter->AddPoint(p,value,sz);
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.;
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
263 fitter->AddPoint(p,value,sdydx);
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.;
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
275 fitter->AddPoint(p,value,sdzdx);
278 void AliTPCcalibAlign::Process9(Double_t *t1,
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)
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];
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.;
314 fitter->AddPoint(p,value,sy);
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.;
328 fitter->AddPoint(p,value,sz);
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.;
337 p[0+1] = -dydx1*dydx2;
338 p[0+2] = -dzdx1*dydx2;
340 fitter->AddPoint(p,value,sdydx);
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.;
349 p[0+1] = -dydx1*dzdx2;
350 p[0+2] = -dzdx1*dzdx2;
352 fitter->AddPoint(p,value,sdzdx);
355 void AliTPCcalibAlign::Process6(Double_t *t1,
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)
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];
373 Double_t sdydx = 0.001;
374 Double_t sdzdx = 0.001;
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.;
391 fitter->AddPoint(p,value,sy);
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.;
405 fitter->AddPoint(p,value,sz);
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.;
414 p[0+1] = -dydx1*dydx2;
415 p[0+2] = -dzdx1*dydx2;
417 fitter->AddPoint(p,value,sdydx);
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.;
426 p[0+1] = -dydx1*dzdx2;
427 p[0+2] = -dzdx1*dzdx2;
429 fitter->AddPoint(p,value,sdzdx);
435 void AliTPCcalibAlign::EvalFitters() {
439 // Perform the fitting using linear fitters
441 Int_t kMinPoints =50;
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;
449 cerr<<"Evaluation failed for "<<s1<<","<<s2<<endl;
450 f->Write(Form("f12_%d_%d",s1,s2));
452 f->Write(Form("f12_%d_%d",s1,s2));
455 if ((f=GetFitter9(s1,s2))&&fPoints[72*s1+s2]>kMinPoints) {
456 // cerr<<s1<<","<<s2<<": "<<fPoints[72*s1+s2]<<endl;
458 cerr<<"Evaluation failed for "<<s1<<","<<s2<<endl;
460 f->Write(Form("f9_%d_%d",s1,s2));
463 if ((f=GetFitter6(s1,s2))&&fPoints[72*s1+s2]>kMinPoints) {
464 // cerr<<s1<<","<<s2<<": "<<fPoints[72*s1+s2]<<endl;
466 cerr<<"Evaluation failed for "<<s1<<","<<s2<<endl;
468 f->Write(Form("f6_%d_%d",s1,s2));
472 this->Write("align");
477 chi212 = align->GetChisquare()/(4.*entries);
481 align->GetParameters(par);
482 align->GetCovarianceMatrix(mat);
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;
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));
504 TLinearFitter* AliTPCcalibAlign::GetOrMakeFitter12(Int_t s1,Int_t s2) {
506 // get or make fitter - general linear transformation
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));
516 TLinearFitter* AliTPCcalibAlign::GetOrMakeFitter9(Int_t s1,Int_t s2) {
518 //get or make fitter - general linear transformation - no scaling
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));
528 TLinearFitter* AliTPCcalibAlign::GetOrMakeFitter6(Int_t s1,Int_t s2) {
530 // get or make fitter - 6 paramater linear tranformation
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));
546 Bool_t AliTPCcalibAlign::GetTransformation12(Int_t s1,Int_t s2,TMatrixD &a) {
548 // GetTransformation matrix - 12 paramaters - generael linear transformation
550 if (!GetFitter12(s1,s2))
554 GetFitter12(s1,s2)->GetParameters(p);
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.;
564 Bool_t AliTPCcalibAlign::GetTransformation9(Int_t s1,Int_t s2,TMatrixD &a) {
566 // GetTransformation matrix - 9 paramaters - general linear transformation
569 if (!GetFitter9(s1,s2))
573 GetFitter9(s1,s2)->GetParameters(p);
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.;
583 Bool_t AliTPCcalibAlign::GetTransformation6(Int_t s1,Int_t s2,TMatrixD &a) {
585 // GetTransformation matrix - 6 paramaters
589 if (!GetFitter6(s1,s2))
593 GetFitter6(s1,s2)->GetParameters(p);
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.;
603 void AliTPCcalibAlign::FillHisto(const AliExternalTrackParam &tp1,
604 const AliExternalTrackParam &tp2,
607 // Fill residual histograms
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());
619 TH1 * AliTPCcalibAlign::GetHisto(HistoType type, Int_t s1, Int_t s2, Bool_t force)
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
626 if (GetIndex(s1,s2)>=72*72) return 0;
627 TObjArray *histoArray=0;
630 histoArray = &fDyHistArray; break;
632 histoArray = &fDzHistArray; break;
634 histoArray = &fDphiHistArray; break;
636 histoArray = &fDthetaHistArray; break;
638 TH1 * histo= (TH1*)histoArray->At(GetIndex(s1,s2));
639 if (histo) return histo;
640 if (force==kFALSE) return 0;
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
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
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
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
666 histo->SetDirectory(0);
667 histoArray->AddAt(histo,GetIndex(s1,s2));