]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCCalibGlobalMisalignment.cxx
add filter bit to output name if it is not the default one
[u/mrichter/AliRoot.git] / TPC / AliTPCCalibGlobalMisalignment.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 // AliTPCCalibGlobalMisalignment class                                        //
20 // The class calculates the space point distortions due to simple         // 
21 // misalignments like shifts in caresian coordinates or a rotation        //
22 // of the TPC read out planes (A and C side)                              //
23 // Optionaly possible to use it for visualization of the alignemnt form the Alignment OCDB //
24 // fUseGeoManager has to be set to kTRUE to enable this option
25 //                                                                        //
26 // date: 06/05/2010                                                       //
27 // Authors: Stefan Rossegger, Jim Thomas, Magnus Mager                    //
28 ////////////////////////////////////////////////////////////////////////////
29
30 #include "AliTPCCalibGlobalMisalignment.h"
31 #include "TMath.h"
32 #include "TGeoMatrix.h"
33 #include "AliTPCROC.h"
34 #include "AliTPCcalibDB.h"
35 #include "AliTPCParam.h"
36 #include <TGeoPhysicalNode.h>
37
38
39 AliTPCCalibGlobalMisalignment::AliTPCCalibGlobalMisalignment()
40   : AliTPCCorrection("mialign","Misalignment"),
41     fXShift(0.),fYShift(0.),fZShift(0.),
42     fRotPhiA(0.),fRotPhiC(0.),
43     fdRPhiOffsetA(0.), fdRPhiOffsetC(0.), 
44     fQuadrantDQ1(0), fQuadrantDQ2(0), fQuadrantQ2(0), fMatrixGlobal(0),
45     fMatrixASide(0), fMatrixCSide(0),
46     fUseGeomanager(kFALSE)
47 {
48   //
49   // default constructor
50   //
51 }
52
53 AliTPCCalibGlobalMisalignment::~AliTPCCalibGlobalMisalignment() {
54   //
55   // default destructor
56   //
57   delete fQuadrantDQ1;   //OROC medium pads delta ly+ - ly-
58   delete fQuadrantDQ2;   //OROC long   pads delta ly+ - ly-
59   delete fQuadrantQ2;    //OROC long   pads - OROC medium pads
60
61 }
62  
63 void AliTPCCalibGlobalMisalignment::SetQuadranAlign(const TVectorD *dq1, const TVectorD *dq2, const TVectorD *q2){
64   //
65   // Set quadrant alignment
66   // 3 vectors for 36 (super) sectors
67   //
68   fQuadrantDQ1 = new TVectorD(*dq1);    //OROC medium pads delta ly+ - ly-
69   fQuadrantDQ2 = new TVectorD(*dq2);;   //OROC long   pads delta ly+ - ly-
70   fQuadrantQ2  = new TVectorD(*q2);;    //OROC long   pads - OROC medium pads
71 }
72
73 void AliTPCCalibGlobalMisalignment::SetGlobalAlign(const TGeoMatrix * matrixGlobal, const TGeoMatrix *matrixA, const TGeoMatrix *matrixC ){
74   //
75   // Set global misalignment as TGeoMatrix
76   // 
77   if (matrixGlobal) fMatrixGlobal = new TGeoHMatrix(*matrixGlobal);
78   if (matrixA) fMatrixASide = new TGeoHMatrix(*matrixA);
79   if (matrixC) fMatrixCSide = new TGeoHMatrix(*matrixC);
80 }
81
82
83 //void AliTPCCalibGlobalMisalignment::Init() {
84 //  //
85 // // Initialization funtion
86 //  //
87
88 //  // nothing to be initialized, results of this calibration class will go to the global aligment structure
89
90 //}
91
92 //void AliTPCCalibGlobalMisalignment::Update(const TTimeStamp &/*timeStamp*/) {
93 //  //
94 //  // Update function 
95 //  //
96 //
97 //  // nothing to be updated, results of this calibration class will go to the global aligment structure
98 //
99 //}
100
101
102
103 void AliTPCCalibGlobalMisalignment::GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]) {
104   //
105   // Calculates the simple correction due to a shift (in x,y,z) or an rotation of the TPC (around z)
106   //  
107   static AliTPCROC *tpcRoc =AliTPCROC::Instance();  
108   Double_t r=0, phi=0;
109   r   = TMath::Sqrt( x[0]*x[0] + x[1]*x[1] );
110   phi = TMath::ATan2(x[1],x[0]);
111   // Getsector number
112   Double_t sec=TMath::Nint(-0.5+(phi*9./TMath::Pi()));
113   if (sec<0) sec+=18;
114   Int_t isec = TMath::Nint(sec);    
115   if  (roc%36>=18) isec+=18;
116   //
117   // Get the point on the local coordiante frame
118   //
119   Double_t alpha=(sec+0.5)*TMath::Pi()/9;
120   Double_t pos[3]={0,0,x[2]};
121   pos[0]=  TMath::Cos(alpha)*x[0]+TMath::Sin(alpha)*x[1];
122   pos[1]= -TMath::Sin(alpha)*x[0]+TMath::Cos(alpha)*x[1];
123   if (pos[0]>tpcRoc->GetPadRowRadiiUp(0))  isec+=36;
124
125   //
126   // apply quadrant alignment if available - in local coordinate frame
127   //
128   Double_t posQG[3]={x[0],x[1],x[2]};
129   if (fQuadrantDQ1){
130     Double_t dly=0;
131     Bool_t isQ1 = TMath::Abs(pos[0]-161)<28;
132     Bool_t isQ2 = (pos[0]>189);
133     if (isQ1){
134       if (pos[1]>0.) dly+=(*fQuadrantDQ1)[isec];
135       if (pos[1]<0.) dly-=(*fQuadrantDQ1)[isec];      
136     }
137     if (isQ2){
138       dly+=(*fQuadrantQ2)[isec];
139       if (pos[1]>0.) dly+=(*fQuadrantDQ2)[isec];
140       if (pos[1]<0.) dly-=(*fQuadrantDQ2)[isec];
141     }
142     // Tranform the corrected point to the global frame
143     posQG[0]=  TMath::Cos(alpha)*pos[0]-TMath::Sin(alpha)*(pos[1]+dly);
144     posQG[1]=  TMath::Sin(alpha)*pos[0]+TMath::Cos(alpha)*(pos[1]+dly);
145   }
146   //
147   // rotation of the read-out planes
148   if  (roc%36<18) // A side
149     phi += fRotPhiA;
150   else         // C side
151     phi += fRotPhiC;
152   
153   // Simply adding a constant dRPHi residual. PURELY FOR CALIBRATION PURPOSES
154   if  (roc%36<18) // A side
155     phi += fdRPhiOffsetA/r;
156   else         // C side
157     phi += fdRPhiOffsetC/r;
158   
159   dx[0] = r * TMath::Cos(phi) - x[0];
160   dx[1] = r * TMath::Sin(phi) - x[1]; 
161   dx[2] = 0.; 
162
163   // Simple shifts
164   dx[0] -= fXShift;
165   dx[1] -= fYShift;
166   dx[2] -= fZShift;
167   // quadrant shifts
168   dx[0] += (posQG[0]-x[0]);
169   dx[1] += (posQG[1]-x[1]);
170   //
171   // alignment matrix in local frame
172   //
173   if (fUseGeomanager){ //loading from the OCDB
174     Double_t posC[3] ={pos[0],pos[1],pos[2]};
175     //
176     //2. correct the point in the local frame
177     AliTPCParam *param = AliTPCcalibDB::Instance()->GetParameters();
178     if (!param){
179       //AliFatal("OCDB not initialized");
180     }
181     TGeoHMatrix  *mat = param->GetClusterMatrix(isec);
182     //
183     if (mat) mat->LocalToMaster(pos,posC);
184     Double_t posCG[3]={posC[0],posC[1],posC[2]};
185     //3. tranform the corrected point to the global frame
186     posCG[0]=  TMath::Cos(alpha)*posC[0]-TMath::Sin(alpha)*posC[1];
187     posCG[1]=  TMath::Sin(alpha)*posC[0]+TMath::Cos(alpha)*posC[1];
188     posCG[2]=  posC[2];
189     //4. Add delta
190     dx[0]+=posCG[0]-x[0];
191     dx[1]+=posCG[1]-x[1];
192     dx[2]+=posCG[2]-x[2];
193   }
194   if (fMatrixGlobal){
195     // apply global alignment matrix
196     Double_t ppos[3]={x[0],x[1],x[2]};
197     Double_t pposC[3]={x[0],x[1],x[2]};
198     fMatrixGlobal->LocalToMaster(ppos,pposC);
199     dx[0]+=pposC[0]-ppos[0];
200     dx[1]+=pposC[1]-ppos[1];
201     dx[2]+=pposC[2]-ppos[2];
202   }
203
204   if (fMatrixASide && roc%36<18){
205     // apply global alignment matrix
206     Double_t ppos[3]={x[0],x[1],x[2]};
207     Double_t pposC[3]={x[0],x[1],x[2]};
208     fMatrixASide->LocalToMaster(ppos,pposC);
209     dx[0]+=pposC[0]-ppos[0];
210     dx[1]+=pposC[1]-ppos[1];
211     dx[2]+=pposC[2]-ppos[2];
212   }
213   if (fMatrixCSide && roc%36>=18){
214     // apply global alignment matrix
215     Double_t ppos[3]={x[0],x[1],x[2]};
216     Double_t pposC[3]={x[0],x[1],x[2]};
217     fMatrixCSide->LocalToMaster(ppos,pposC);
218     dx[0]+=pposC[0]-ppos[0];
219     dx[1]+=pposC[1]-ppos[1];
220     dx[2]+=pposC[2]-ppos[2];
221   }
222 }
223
224 void AliTPCCalibGlobalMisalignment::Print(Option_t* /*option*/ ) const {
225   //
226   // Print function to check the settings 
227   //
228   printf("%s",GetTitle());  
229   printf(" - Trivial Misalignments for calibration purposes: \n");
230   printf(" - X-Shift: %1.3f cm, Y-Shift: %1.3f cm, Z-Shift: %1.3f cm \n",fXShift,fYShift,fZShift);
231   printf(" - Phi-Rotations: A side: %1.5f rad, C side: %1.5f rad\n",fRotPhiA,fRotPhiC);
232   printf(" - dRPhi offsets: A side: %1.5f cm, C side: %1.5f cm\n",fdRPhiOffsetA,fdRPhiOffsetC);
233  
234  
235 }