]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCTransformation.cxx
20d3b2eb308ae5555c1166b5ee0d393beaf76e16
[u/mrichter/AliRoot.git] / TPC / AliTPCTransformation.cxx
1 /*
2   
3  Class AliTPCtransformation:
4  Should represent general non linear transformation. Currently tune for TPConly.
5  To be used:
6  1. Simulation-Digitization
7  2. Reconstruction - AliTPCTransform
8  3. Calibration/Alignment (KalmanFilter, Milipedde)
9  4. Set of transformation to be stored/retrieved as OCDB entry 
10
11  Base functionality:
12  
13  1. Double_t GetDeltaXYZ(Int_t coord, Int_t volID, Double_t param, Double_t x, Double_t y, Double_t z)
14  Get correction - return the delta of coordinate coord dx or dy or dz for given volID for point at point (x,y,z)
15  All coordinates are global
16
17  2. The transformation should work only for given volIDs and detector IDs
18     Currently Bitmask is used for filtering 
19
20
21  Transformation - naming convention:
22  //
23  XXX(local)YYYZZZ
24  TPClocaldLxdGX
25  XXX   - detector if detector specific
26  local - if local transforamtion
27  YYY   - type of transformation
28  ZZZ   - return type of transformation
29
30 */
31
32 #include <string.h>
33 #include "TRandom.h"
34 #include "TMath.h"
35 #include "TBits.h"
36 #include "TFormula.h"
37 #include "TF1.h"
38 #include "TLinearFitter.h"
39 #include "TFile.h"
40 #include "TObjString.h"
41
42 #include "TTreeStream.h"
43 #include "AliTrackPointArray.h"
44 #include "AliLog.h"
45 #include "AliTPCTransformation.h"
46
47 ClassImp(AliTPCTransformation)
48
49
50 AliTPCTransformation::GenFuncG    AliTPCTransformation::fgFormulas[10000];
51 TObjArray*  AliTPCTransformation::fgFormulasName = new TObjArray(10000);
52
53
54 void AliTPCTransformation::RegisterFormula(const char * name, GenFuncG formula){
55   //
56   // Add Formula to the list of formulas
57   //
58   Int_t last= fgFormulasName->GetEntries();
59   fgFormulasName->AddLast(new TObjString(name));
60   fgFormulas[last]=formula;
61 }
62
63 Int_t  AliTPCTransformation::BuildBasicFormulas(){
64   //
65
66   //
67   //build list of basic TPC formulas - corrections
68   //
69   RegisterFormula("TPCscalingRPol",(GenFuncG)(AliTPCTransformation::TPCscalingRPol));
70   RegisterFormula("TPCscalingRIFC",(GenFuncG)(AliTPCTransformation::TPCscalingRIFC));
71   RegisterFormula("TPCscalingROFC",(GenFuncG)(AliTPCTransformation::TPCscalingROFC));
72   //
73   RegisterFormula("TPCscalingZDr",(GenFuncG)(AliTPCTransformation::TPCscalingZDr));
74   RegisterFormula("TPCscalingPhiLocal",(GenFuncG)(AliTPCTransformation::TPCscalingPhiLocal));
75   //
76   // TPC Local X and Y misalignment + rotation 
77   //
78   RegisterFormula("TPClocaldLxdGX",(GenFuncG)(AliTPCTransformation::TPClocaldLxdGX));
79   RegisterFormula("TPClocaldLxdGY",(GenFuncG)(AliTPCTransformation::TPClocaldLxdGY));
80   RegisterFormula("TPClocaldLydGX",(GenFuncG)(AliTPCTransformation::TPClocaldLydGX));
81   RegisterFormula("TPClocaldLydGY",(GenFuncG)(AliTPCTransformation::TPClocaldLydGY));
82   RegisterFormula("TPClocaldRzdGX",(GenFuncG)(AliTPCTransformation::TPClocaldRzdGX));
83   RegisterFormula("TPClocaldRzdGY",(GenFuncG)(AliTPCTransformation::TPClocaldRzdGY));
84
85   //
86   // Z offset
87   //
88   RegisterFormula("TPCDeltaZ",(GenFuncG)(AliTPCTransformation::TPCDeltaZ));
89   RegisterFormula("TPCDeltaZMediumLong",(GenFuncG)(AliTPCTransformation::TPCDeltaZMediumLong));
90   RegisterFormula("TPCTiltingZ",(GenFuncG)(AliTPCTransformation::TPCTiltingZ));
91   return 0;
92 }
93
94 AliTPCTransformation::GenFuncG  AliTPCTransformation::FindFormula(const char * name){
95   //
96   // find formula - if registered
97   //
98   if (fgFormulasName->FindObject(name)==0) return 0;
99   Int_t entries = fgFormulasName->GetEntries();
100   for (Int_t i=0;i<entries;i++){
101     if (strcmp(fgFormulasName->At(i)->GetName(), name)==0){
102       return fgFormulas[i];
103     }
104   }
105   return 0;
106 }
107
108 Double_t AliTPCTransformation::Eval(const char * name, const Double_t*x,const Double_t*par){
109   //
110   // Only for test purposes - very slow
111   //
112   GenFuncG fun = FindFormula(name);
113   if (!fun) return 0;
114   return fun(x,par);
115 }
116
117
118 AliTPCTransformation::AliTPCTransformation():
119   TNamed(),
120   fNameX(0),          // x formula name
121   fNameY(0),          // y formula name
122   fNameZ(0),          // z formula name 
123   //
124   fBitMask(0),        // bitmaps - transformation only for specified volID
125   fCoordSystem(0),    // coord system of  output deltas 
126   fParam(0),          // free parameter of transformation 
127   fSigma(0),          // uncertainty of the parameter
128   fSigma2Time(0),     // change of the error in time - (For kalman filter) 
129   fFixedParam(0),     // fixed parameters of tranformation  
130   fIsActive(kTRUE),   // swith - On/Off
131   //
132   fInit(kFALSE),      // initialization flag - set to kTRUE if corresponding formulas found
133   fFormulaX(0),       // x formula - pointer to the function
134   fFormulaY(0),       // y formula - pointer to the function
135   fFormulaZ(0)       // z formula - pointer to the function
136   //
137 {
138   //
139   // default constructor
140   //
141 }
142
143
144
145 AliTPCTransformation::AliTPCTransformation(const char *name, TBits *mask, const char *fx, const char *fy, const char *fz, Int_t coordSystem):
146   TNamed(name,name),
147   fNameX(0),       // x formula name
148   fNameY(0),       // y formula name
149   fNameZ(0),       // z formula name
150   fBitMask(mask),   // bitmaps - transformation only for specified volID
151   fCoordSystem(coordSystem), // coordinate system of output deltas
152   fParam(0),          // free parameter of transformation 
153   fSigma(0),
154   fSigma2Time(0),     // change of sigma in time
155   fFixedParam(0),     // fixed parameters of tranformation  
156   fIsActive(kTRUE),   // swith - On/Off
157   //
158   fInit(kFALSE),      // initialization flag - set to kTRUE if corresponding formulas found
159   fFormulaX(0),       // x formula - pointer to the function
160   fFormulaY(0),       // y formula - pointer to the function
161   fFormulaZ(0)       // z formula - pointer to the function
162 {
163   //
164   // non default constructor
165   //
166   if (fx) fNameX= new TString(fx);
167   if (fy) fNameY= new TString(fy);
168   if (fz) fNameZ= new TString(fz);
169   Init();
170 }
171
172 AliTPCTransformation::AliTPCTransformation(const AliTPCTransformation&trafo):
173   TNamed(trafo),
174   fNameX(0),          // x formula name
175   fNameY(0),          // y formula name
176   fNameZ(0),          // z formula name 
177                      //
178   fBitMask(0),        // bitmaps - transformation only for specified volID
179   fCoordSystem(0),    // coord system of  output deltas 
180   fParam(trafo.fParam),          // free parameter of transformation 
181   fSigma(trafo.fSigma),          // uncertainty of the parameter
182   fSigma2Time(trafo.fSigma2Time),     // change of the error in time - (For kalman filter) 
183   fFixedParam(0),     // fixed parameters of tranformation
184   fIsActive(trafo.fIsActive),   // swith - On/Off
185   //
186   fInit(kFALSE),      // initialization flag - set to kTRUE if corresponding formulas found
187   fFormulaX(0),       // x formula - pointer to the function
188   fFormulaY(0),       // y formula - pointer to the function
189   fFormulaZ(0)       // z formula - pointer to the function
190 {
191   if (trafo.fNameX) fNameX = new TString(*(trafo.fNameX)); 
192   if (trafo.fNameY) fNameY = new TString(*(trafo.fNameY)); 
193   if (trafo.fNameZ) fNameZ = new TString(*(trafo.fNameZ)); 
194   if (trafo.fBitMask)  fBitMask = new TBits(*(trafo.fBitMask)); 
195 }
196
197 AliTPCTransformation::~AliTPCTransformation(){
198   //
199   // destructor
200   //
201   delete fNameX;
202   delete fNameY;
203   delete fNameZ;
204   delete fBitMask;
205   delete fFixedParam;
206 }
207 void AliTPCTransformation::SetParams(Double_t param, Double_t sigma, Double_t sigma2Time, TVectorD* fixedParams){
208   //
209   // Set parameters of transformation
210   //
211   fParam = param;
212   fSigma = sigma;
213   fSigma2Time = sigma2Time;
214   if (fFixedParam) delete fFixedParam;
215   fFixedParam = new TVectorD(*fixedParams);
216   Init();
217 }
218
219
220 Bool_t AliTPCTransformation::Init(){
221   //
222   // associate formulas with pointer to the function
223   //
224   Bool_t isOK=kTRUE;
225   if (fNameX) {
226     fFormulaX=FindFormula(fNameX->Data());
227     if (fFormulaX==0) isOK=kFALSE;
228   }
229   if (fNameY) {
230     fFormulaY=FindFormula(fNameY->Data());
231     if (fFormulaY==0) isOK=kFALSE;
232   }
233   if (fNameZ) {
234     fFormulaZ=FindFormula(fNameZ->Data());
235     if (!fFormulaZ) isOK=kFALSE;
236   }
237   return isOK;
238 }
239
240
241
242 TBits * AliTPCTransformation::BitsSide(Bool_t aside){
243   //
244   TBits * bits = new TBits(72);
245   for (Int_t i=0; i<72;i++){
246     if (i%36<18 && aside) (*bits)[i]=kTRUE;
247     if (i%36<18 && (!aside)) (*bits)[i]=kFALSE;
248     if (i%36>=18 && aside) (*bits)[i]=kFALSE;
249     if (i%36>=18 && (!aside)) (*bits)[i]=kTRUE;
250   }
251   return bits;
252 }
253
254 TBits * AliTPCTransformation::BitsAll(){
255   //
256   //
257   //
258   TBits * bits = new TBits(72);
259   for (Int_t i=0; i<72;i++){
260     (*bits)[i]=kTRUE;
261   }
262   return bits;
263 }
264
265 Double_t AliTPCTransformation::GetDeltaXYZ(Int_t coord, Int_t volID, Double_t param, Double_t x, Double_t y, Double_t z){
266   //
267   //
268   // coord - type of coordinate
269   //       - 0 -X
270   //         1 -Y
271   //         2 -Z
272   //         3 -R
273   //         4 -RPhi
274   if (!fIsActive) return 0;
275   if (fBitMask && (!(*fBitMask)[volID])) return 0;
276   Double_t xyz[5]={x,y,z, param,volID};
277   if (fCoordSystem==0){
278     // cartezian system
279     if (coord==0 && fFormulaX) return fFormulaX(xyz,fFixedParam->GetMatrixArray()); 
280     if (coord==1 && fFormulaY) return fFormulaY(xyz,fFixedParam->GetMatrixArray()); 
281     if (coord==2 && fFormulaZ) return fFormulaZ(xyz,fFixedParam->GetMatrixArray());
282   }
283   if (fCoordSystem==1){  
284     // cylindrical system
285     if (coord==2) {
286       if (fFormulaZ==0) return 0;
287       return fFormulaZ(xyz,fFixedParam->GetMatrixArray());
288     }
289     Double_t rrphiz[3]={0,0,0};
290     if (fFormulaX) rrphiz[0] = fFormulaX(xyz,fFixedParam->GetMatrixArray());
291     if (fFormulaY) rrphiz[1] = fFormulaY(xyz,fFixedParam->GetMatrixArray());
292     Double_t alpha = TMath::ATan2(y,x);
293     Double_t ca    = TMath::Cos(alpha);
294     Double_t sa    = TMath::Sin(alpha);
295     if (coord==0) return ca*rrphiz[0]-sa*rrphiz[1];
296     if (coord==1) return sa*rrphiz[0]+ca*rrphiz[1];
297   }
298   return 0;
299 }
300
301 Double_t  AliTPCTransformation::TPCscalingRPol(Double_t *xyz, Double_t * param){
302   //
303   // Scaling and shift of TPC radius
304   // xyz[0..2] - global xyz of point 
305   // xyz[3]    - scale parameter
306   // param[0]  - radial scaling power
307   // param[1]  - drift  scaling power
308   // radius  from -1(at rInner)   to 1 (rOuter)
309   // driftM  from -1(at 0 drift)  to 1 (250 cm drift)
310
311   Double_t rInner=78.8;
312   Double_t rOuter=258.0; 
313   Double_t deltaR  = rOuter-rInner;  
314   Double_t radius  = (TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1])-rInner)*2./deltaR; 
315   Double_t driftM  = (0.5 - TMath::Abs(xyz[2]/250.))*2.0;
316   Double_t delta   = TMath::Power(radius,param[0])*TMath::Power(driftM,param[1]);
317   return delta*xyz[3];
318 }
319
320
321 Double_t  AliTPCTransformation::TPCscalingZDr(Double_t *xyz, Double_t * param){
322   //
323   //
324   // Scaling and shift of TPC radius
325   // xyz[0..2] - global xyz of point 
326   // xyz[3]    - scale parameter
327   Double_t driftP  = TMath::Power(1. - TMath::Abs(xyz[2]/250.), param[0]);
328   Double_t deltaZ  = (xyz[2]>0) ? -driftP : driftP;
329   return deltaZ*xyz[3];
330 }
331
332 Double_t  AliTPCTransformation::TPCscalingPhiLocal(Double_t *xyz, Double_t * param){
333   //
334   //
335   // Scaling if the local y -phi
336   // xyz[0..2] - global xyz of point 
337   // xyz[3]    - scale parameter
338   // value = 1 for ful drift length and parameter 1
339   Double_t alpha       = TMath::ATan2(xyz[1],xyz[0]);
340   Double_t sector      = TMath::Nint(9*alpha/TMath::Pi()-0.5);
341   Double_t localAlpha  = (alpha-(sector+0.5)*TMath::Pi()/9.);
342   Double_t radius      = TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1])/250.;
343   //
344   Double_t deltaAlpha  = radius*TMath::Power(2.*9.*localAlpha/TMath::Pi(),param[0]);
345   return deltaAlpha*xyz[3];
346 }
347
348
349 Double_t       AliTPCTransformation::TPCscalingRIFC(Double_t *xyz, Double_t * param){
350   //
351   // inner field cage r distorion - proportinal to 1 over distance to the IFC
352   // param[0] - drift polynom order
353   // distortion at first pad row - is normalized to 
354   Double_t rInner=78.8;
355   Double_t rFirst=85.2; 
356   Double_t deltaR  = rFirst-rInner;
357   Double_t ndistR  = (TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1])-rInner)/deltaR;
358   Double_t driftM  = (0.5 - TMath::Abs(xyz[2]/250.))*2.;
359   Double_t value   = TMath::Power(driftM,param[0])/ndistR;
360   return xyz[3]*value;
361 }
362
363 Double_t       AliTPCTransformation::TPCscalingROFC(Double_t *xyz, Double_t * param){
364   //
365   // outer field cage r distorion - proportinal to 1 over distance to the OFC
366   // param[0] - drift polynom order
367   // driftM   - from -1 to 1 
368   //
369   Double_t rLast=245.8;
370   Double_t rOuter=258.0;  
371   Double_t deltaR  = rOuter-rLast;
372   Double_t ndistR  = (rOuter-TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]))/deltaR;
373   Double_t driftM  = (0.5 - TMath::Abs(xyz[2]/250.))*2.;
374   Double_t value   = TMath::Power(driftM,param[0])/ndistR;
375   return xyz[3]*value;
376 }
377
378
379 //
380 // TPC sector local misalignment 
381 //
382 //
383 //
384 Double_t AliTPCTransformation:: TPClocaldLxdGX(Double_t *xyz, Double_t * param){
385   //
386   // xyz - [0..2] - position 
387   //       [3]    - scale parameter
388   //       [4]    - volID
389   // param[0]= n  - cos(n *alpha)
390   // param[1]= n  - sin(n *alpha)
391   // param[2]     - indication - 0 - the same for IROC OROC 1 - opposite
392   // return delta in global coordiante system
393   //
394   Int_t    sector = TMath::Nint(xyz[4]);
395   Double_t alpha  = TMath::Pi()*(sector+0.5)/9;
396   Double_t ca     = TMath::Cos(alpha);
397   //  Double_t sa     = TMath::Sin(alpha); 
398   const Double_t xIROCOROC = 133.4;
399   Double_t factor = xyz[3];
400   if (param[0]>0)  factor*=TMath::Cos(alpha*param[0]);
401   if (param[1]>0)  factor*=TMath::Sin(alpha*param[1]);
402   if (param[2]>0.5 && TMath::Sqrt(xyz[1]*xyz[1]+xyz[0]*xyz[0])>xIROCOROC) factor*=-1;
403   return ca*factor;    
404 }
405
406 Double_t AliTPCTransformation::TPClocaldLxdGY(Double_t *xyz, Double_t * param){
407   //
408   // xyz - [0..2] - position 
409   //       [3]    - scale parameter
410   //       [4]    - volID
411   // param[0]= n  - cos(n *alpha)
412   // param[1]= n  - sin(n *alpha)
413   // param[2]     - indication - 0 - the same for IROC OROC 1 - opposite
414   // return delta in global coordiante system
415   //
416   Int_t    sector = TMath::Nint(xyz[4]);
417   Double_t alpha  = TMath::Pi()*(sector+0.5)/9;
418   //Double_t ca     = TMath::Cos(alpha);
419   Double_t sa     = TMath::Sin(alpha);
420   const Double_t xIROCOROC = 133.4;
421   Double_t factor = xyz[3];
422   if (param[0]>0)  factor*=TMath::Cos(alpha*param[0]);
423   if (param[1]>0)  factor*=TMath::Sin(alpha*param[1]);
424   if (param[2]>0.5 && TMath::Sqrt(xyz[1]*xyz[1]+xyz[0]*xyz[0])>xIROCOROC) factor*=-1;  
425   return   sa*factor;  
426 }
427
428 Double_t AliTPCTransformation:: TPClocaldLydGX(Double_t *xyz, Double_t * param){
429   //
430   // xyz - [0..2] - position 
431   //       [3]    - scale parameter
432   //       [4]    - volID
433   // param[0]= n  - cos(n *alpha)
434   // param[1]= n  - sin(n *alpha)
435   // param[2]     - indication - 0 - the same for IROC OROC 1 - opposite
436   // return delta in global coordiante system
437   //
438   Int_t    sector = TMath::Nint(xyz[4]);
439   Double_t alpha  = TMath::Pi()*(sector+0.5)/9;
440   //Double_t ca     = TMath::Cos(alpha);
441   Double_t sa     = TMath::Sin(alpha);
442   const Double_t xIROCOROC = 133.4;
443   Double_t factor = xyz[3];
444   if (param[0]>0)  factor*=TMath::Cos(alpha*param[0]);
445   if (param[1]>0)  factor*=TMath::Sin(alpha*param[1]);
446   if (param[2]>0.5 && TMath::Sqrt(xyz[1]*xyz[1]+xyz[0]*xyz[0])>xIROCOROC) factor*=-1;  
447   return            -sa*factor;  
448 }
449
450 Double_t AliTPCTransformation::TPClocaldLydGY(Double_t *xyz, Double_t * param){
451   //
452   // xyz - [0..2] - position 
453   //       [3]    - scale parameter
454   //       [4]    - volID
455   // param[0]= n  - cos(n *alpha)
456   // param[1]= n  - sin(n *alpha)
457   // param[2]     - indication - 0 - the same for IROC OROC 1 - opposite
458   // return delta in global coordiante system
459   //
460   Int_t    sector = TMath::Nint(xyz[4]);
461   Double_t alpha  = TMath::Pi()*(sector+0.5)/9;
462   Double_t ca     = TMath::Cos(alpha);
463   //Double_t sa     = TMath::Sin(alpha);
464   const Double_t xIROCOROC = 133.4;
465   Double_t factor = xyz[3];
466   if (param[0]>0)  factor*=TMath::Cos(alpha*param[0]);
467   if (param[1]>0)  factor*=TMath::Sin(alpha*param[1]);
468   if (param[2]>0.5 && TMath::Sqrt(xyz[1]*xyz[1]+xyz[0]*xyz[0])>xIROCOROC) factor*=-1;  
469   return   ca*factor;  
470 }
471
472
473 Double_t AliTPCTransformation::TPClocaldRzdGX(Double_t *xyz, Double_t * param){
474   //
475   // xyz - [0..2] - position 
476   //       [3]    - scale parameter - rotation angle in mrad
477   //       [4]    - volID
478   // param[0]= n  - cos(n *alpha)
479   // param[1]= n  - sin(n *alpha)
480   // param[2]     - indication - 0 - the same for IROC OROC 1 - opposite  
481   // return delta in global coordiante system
482   //
483   Int_t    sector = TMath::Nint(xyz[4]);
484   Double_t alpha  = TMath::Pi()*(sector+0.5)/9;
485   Double_t ca     = TMath::Cos(alpha);  
486   Double_t sa     = TMath::Sin(alpha);
487   Double_t lx     =  xyz[0]*ca + xyz[1]*sa;
488   Double_t ly     = -xyz[0]*sa + xyz[1]*ca;
489   //
490   const Double_t xIROCOROC = 133.4;  
491   lx-=xIROCOROC;
492   Double_t rot      =  xyz[3]*0.001;    // rotation in mrad
493   if (param[0]>0)  rot*=TMath::Cos(alpha*param[0]);
494   if (param[1]>0)  rot*=TMath::Sin(alpha*param[1]);
495   if (param[2]>0.5 && lx>0) rot*=-1;  
496   //
497   Double_t dlxR     =  - ly*rot; 
498   Double_t dlyR     =    lx*rot;
499   Double_t dgxR     =  dlxR*ca - dlyR*sa;
500   //Double_t dgyR     =  dlxR*sa + dlyR*ca;
501   return  dgxR;            
502 }
503
504 Double_t AliTPCTransformation::TPClocaldRzdGY(Double_t *xyz, Double_t * param){
505   //
506   // xyz - [0..2] - position 
507   //       [3]    - scale parameter - rotation angle in mrad
508   //       [4]    - volID
509   // param[0]= n  - cos(n *alpha)
510   // param[1]= n  - sin(n *alpha)
511   // param[2]     - indication - 0 - the same for IROC OROC 1 - opposite
512   // return delta in global coordiante system
513   //
514   Int_t    sector = TMath::Nint(xyz[4]);
515   Double_t alpha  = TMath::Pi()*(sector+0.5)/9;
516   Double_t ca     = TMath::Cos(alpha);  
517   Double_t sa     = TMath::Sin(alpha);
518   Double_t lx     =  xyz[0]*ca + xyz[1]*sa;
519   Double_t ly     = -xyz[0]*sa + xyz[1]*ca;
520   //
521   const Double_t xIROCOROC = 133.4;  
522   lx-=xIROCOROC;
523   Double_t rot      =  xyz[3]*0.001;    // rotation in mrad
524   if (param[0]>0)  rot*=TMath::Cos(alpha*param[0]);
525   if (param[1]>0)  rot*=TMath::Sin(alpha*param[1]);
526   if (param[2]>0.5 && lx>0) rot*=-1;  
527   Double_t dlxR     =  - ly*rot; 
528   Double_t dlyR     =    lx*rot;
529   //Double_t dgxR     =  dlxR*ca - dlyR*sa;
530   Double_t dgyR     =  dlxR*sa + dlyR*ca;
531   return  dgyR;            
532 }
533
534
535 Double_t        AliTPCTransformation::TPCDeltaZ(Double_t *xyz, Double_t * param){
536   //
537   // xyz - [0..2] - position 
538   //        [3]    - scale parameter
539   //        [4]    - volID
540   // return delta in global coordiante system
541   //
542   Int_t    sector = TMath::Nint(xyz[4]);
543   Double_t signZ  = (sector%36<18) ? 1: -1;  // drift direction
544   return signZ*xyz[3];     // IROC shift
545 }
546
547 Double_t        AliTPCTransformation::TPCDeltaZMediumLong(Double_t *xyz, Double_t * /*param*/){
548   //
549   // xyz - [0..2] - position 
550   //        [3]    - scale parameter
551   //        [4]    - volID
552   // return delta in global coordinate system 
553   //
554   Int_t    sector = TMath::Nint(xyz[4]);
555   Double_t signZ  = (sector%36<18) ? 1: -1;  // drift direction
556   if    (sector<36) return 0;     
557   //
558   Double_t radius  = (TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1])); 
559   const Double_t radiusLong = 198.1;
560   Double_t sign   = (radius<radiusLong) ? 1:-1;
561   return xyz[3]*sign*signZ;
562 }
563
564 Double_t       AliTPCTransformation::TPCTiltingZ(Double_t *xyz, Double_t * param){
565   // xyz - [0..2] - position 
566   //        [3]    - scale parameter
567   //        [4]    - volID
568   // param[0]      - n for cos
569   // param[1]      - n for sin
570   // return delta in global coordinate system 
571   Double_t radius = TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]);
572   const Double_t rFirst=85.2; 
573   const Double_t rLast =245.8;
574   Double_t radiusC  = (rFirst+rLast)*0.5;
575   Double_t deltaR = 2.0*(radius-radiusC)/(rLast-rFirst);
576   Double_t alpha       = TMath::ATan2(xyz[1],xyz[0]);
577   
578   if (param[0]>0) deltaR *= TMath::Cos(param[0]*alpha);
579   if (param[1]>0) deltaR *= TMath::Sin(param[1]*alpha);
580   return deltaR*xyz[3];
581 }
582