]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCTransformation.cxx
Adding switch - usage of maximal or total charge for dEdx
[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::TPCscalingZDrift));
74   RegisterFormula("TPCscalingZDrGy",(GenFuncG)(AliTPCTransformation::TPCscalingZDriftGy));
75   RegisterFormula("TPCscalingZDriftT0",(GenFuncG)(AliTPCTransformation::TPCscalingZDriftT0));
76   RegisterFormula("TPCscalingPhiLocal",(GenFuncG)(AliTPCTransformation::TPCscalingPhiLocal));
77   RegisterFormula("TPClocalRPhiEdge",(GenFuncG)(AliTPCTransformation::TPClocalRPhiEdge));
78   //
79   // TPC Local X and Y misalignment + rotation 
80   //
81   RegisterFormula("TPClocaldLxdGX",(GenFuncG)(AliTPCTransformation::TPClocaldLxdGX));
82   RegisterFormula("TPClocaldLxdGY",(GenFuncG)(AliTPCTransformation::TPClocaldLxdGY));
83   RegisterFormula("TPClocaldLydGX",(GenFuncG)(AliTPCTransformation::TPClocaldLydGX));
84   RegisterFormula("TPClocaldLydGY",(GenFuncG)(AliTPCTransformation::TPClocaldLydGY));
85   RegisterFormula("TPClocaldRzdGX",(GenFuncG)(AliTPCTransformation::TPClocaldRzdGX));
86   RegisterFormula("TPClocaldRzdGY",(GenFuncG)(AliTPCTransformation::TPClocaldRzdGY));
87
88   //
89   // Z offset
90   //
91   RegisterFormula("TPCDeltaZ",(GenFuncG)(AliTPCTransformation::TPCDeltaZ));
92   RegisterFormula("TPCDeltaZMediumLong",(GenFuncG)(AliTPCTransformation::TPCDeltaZMediumLong));
93   RegisterFormula("TPCTiltingZ",(GenFuncG)(AliTPCTransformation::TPCTiltingZ));
94   return 0;
95 }
96
97 AliTPCTransformation::GenFuncG  AliTPCTransformation::FindFormula(const char * name){
98   //
99   // find formula - if registered
100   //
101   if (fgFormulasName->FindObject(name)==0) return 0;
102   Int_t entries = fgFormulasName->GetEntries();
103   for (Int_t i=0;i<entries;i++){
104     if (strcmp(fgFormulasName->At(i)->GetName(), name)==0){
105       return fgFormulas[i];
106     }
107   }
108   return 0;
109 }
110
111 Double_t AliTPCTransformation::Eval(const char * name, const Double_t*x,const Double_t*par){
112   //
113   // Only for test purposes - very slow
114   //
115   GenFuncG fun = FindFormula(name);
116   if (!fun) return 0;
117   return fun(x,par);
118 }
119
120
121 AliTPCTransformation::AliTPCTransformation():
122   TNamed(),
123   fNameX(0),          // x formula name
124   fNameY(0),          // y formula name
125   fNameZ(0),          // z formula name 
126   //
127   fBitMask(0),        // bitmaps - transformation only for specified volID
128   fCoordSystem(0),    // coord system of  output deltas 
129   fParam(0),          // free parameter of transformation 
130   fSigma(0),          // uncertainty of the parameter
131   fSigmaMax(0),       //maximal sigma (Not allowed to increase in propagate time by bigger factor)
132   fSigma2Time(0),     // change of the error in time - (For kalman filter) 
133   fFixedParam(0),     // fixed parameters of tranformation  
134   fIsActive(kTRUE),   // swith - On/Off
135   //
136   fInit(kFALSE),      // initialization flag - set to kTRUE if corresponding formulas found
137   fFormulaX(0),       // x formula - pointer to the function
138   fFormulaY(0),       // y formula - pointer to the function
139   fFormulaZ(0)       // z formula - pointer to the function
140   //
141 {
142   //
143   // default constructor
144   //
145 }
146
147
148
149 AliTPCTransformation::AliTPCTransformation(const char *name, TBits *mask, const char *fx, const char *fy, const char *fz, Int_t coordSystem):
150   TNamed(name,name),
151   fNameX(0),       // x formula name
152   fNameY(0),       // y formula name
153   fNameZ(0),       // z formula name
154   fBitMask(mask),   // bitmaps - transformation only for specified volID
155   fCoordSystem(coordSystem), // coordinate system of output deltas
156   fParam(0),          // free parameter of transformation 
157   fSigma(0),
158   fSigmaMax(0),       //maximal sigma (Not allowed to increase in propagate time by bigger factor)
159   fSigma2Time(0),     // change of sigma in time
160   fFixedParam(0),     // fixed parameters of tranformation  
161   fIsActive(kTRUE),   // swith - On/Off
162   //
163   fInit(kFALSE),      // initialization flag - set to kTRUE if corresponding formulas found
164   fFormulaX(0),       // x formula - pointer to the function
165   fFormulaY(0),       // y formula - pointer to the function
166   fFormulaZ(0)       // z formula - pointer to the function
167 {
168   //
169   // non default constructor
170   //
171   if (fx) fNameX= new TString(fx);
172   if (fy) fNameY= new TString(fy);
173   if (fz) fNameZ= new TString(fz);
174   Init();
175 }
176
177 AliTPCTransformation::AliTPCTransformation(const AliTPCTransformation&trafo):
178   TNamed(trafo),
179   fNameX(0),          // x formula name
180   fNameY(0),          // y formula name
181   fNameZ(0),          // z formula name 
182                      //
183   fBitMask(0),        // bitmaps - transformation only for specified volID
184   fCoordSystem(0),    // coord system of  output deltas 
185   fParam(trafo.fParam),          // free parameter of transformation 
186   fSigma(trafo.fSigma),          // uncertainty of the parameter
187   fSigmaMax(trafo.fSigma),       //maximal sigma (Not allowed to increase in propagate time by bigger factor)
188   fSigma2Time(trafo.fSigma2Time),     // change of the error in time - (For kalman filter) 
189   fFixedParam(0),     // fixed parameters of tranformation
190   fIsActive(trafo.fIsActive),   // swith - On/Off
191   //
192   fInit(kFALSE),      // initialization flag - set to kTRUE if corresponding formulas found
193   fFormulaX(0),       // x formula - pointer to the function
194   fFormulaY(0),       // y formula - pointer to the function
195   fFormulaZ(0)       // z formula - pointer to the function
196 {
197   if (trafo.fNameX) fNameX = new TString(*(trafo.fNameX)); 
198   if (trafo.fNameY) fNameY = new TString(*(trafo.fNameY)); 
199   if (trafo.fNameZ) fNameZ = new TString(*(trafo.fNameZ)); 
200   if (trafo.fBitMask)  fBitMask = new TBits(*(trafo.fBitMask)); 
201 }
202
203 AliTPCTransformation::~AliTPCTransformation(){
204   //
205   // destructor
206   //
207   delete fNameX;
208   delete fNameY;
209   delete fNameZ;
210   delete fBitMask;
211   delete fFixedParam;
212 }
213 void AliTPCTransformation::SetParams(Double_t param, Double_t sigma, Double_t sigma2Time, TVectorD* fixedParams){
214   //
215   // Set parameters of transformation
216   //
217   fParam = param;
218   fSigma = sigma;
219   fSigmaMax = sigma;
220   fSigma2Time = sigma2Time;
221   if (fFixedParam) delete fFixedParam;
222   fFixedParam = new TVectorD(*fixedParams);
223   Init();
224 }
225
226
227 Bool_t AliTPCTransformation::Init(){
228   //
229   // associate formulas with pointer to the function
230   //
231   Bool_t isOK=kTRUE;
232   if (fNameX) {
233     fFormulaX=FindFormula(fNameX->Data());
234     if (fFormulaX==0) isOK=kFALSE;
235   }
236   if (fNameY) {
237     fFormulaY=FindFormula(fNameY->Data());
238     if (fFormulaY==0) isOK=kFALSE;
239   }
240   if (fNameZ) {
241     fFormulaZ=FindFormula(fNameZ->Data());
242     if (!fFormulaZ) isOK=kFALSE;
243   }
244   return isOK;
245 }
246
247
248
249 TBits * AliTPCTransformation::BitsSide(Bool_t aside){
250   //
251   TBits * bits = new TBits(72);
252   for (Int_t i=0; i<72;i++){
253     if (i%36<18 && aside) (*bits)[i]=kTRUE;
254     if (i%36<18 && (!aside)) (*bits)[i]=kFALSE;
255     if (i%36>=18 && aside) (*bits)[i]=kFALSE;
256     if (i%36>=18 && (!aside)) (*bits)[i]=kTRUE;
257   }
258   return bits;
259 }
260
261 TBits * AliTPCTransformation::BitsAll(){
262   //
263   //
264   //
265   TBits * bits = new TBits(72);
266   for (Int_t i=0; i<72;i++){
267     (*bits)[i]=kTRUE;
268   }
269   return bits;
270 }
271
272 // Double_t AliTPCTransformation::GetDeltaXYZ(Int_t coord, Int_t volID, Double_t param, Double_t x, Double_t y, Double_t z){
273 //   //
274 //   //
275 //   // coord - type of coordinate
276 //   //       - 0 -X
277 //   //         1 -Y
278 //   //         2 -Z
279 //   //         3 -R
280 //   //         4 -RPhi
281 //   if (!fIsActive) return 0;
282 //   if (fBitMask && (!(*fBitMask)[volID])) return 0;
283 //   Double_t xyz[5]={x,y,z, param,volID};
284 //   if (fCoordSystem==0){
285 //     // cartezian system
286 //     if (coord==0 && fFormulaX) return fFormulaX(xyz,fFixedParam->GetMatrixArray()); 
287 //     if (coord==1 && fFormulaY) return fFormulaY(xyz,fFixedParam->GetMatrixArray()); 
288 //     if (coord==2 && fFormulaZ) return fFormulaZ(xyz,fFixedParam->GetMatrixArray());
289 //   }
290 //   if (fCoordSystem==1){  
291 //     // cylindrical system
292 //     if (coord==2) {
293 //       if (fFormulaZ==0) return 0;
294 //       return fFormulaZ(xyz,fFixedParam->GetMatrixArray());
295 //     }
296 //     Double_t rrphiz[3]={0,0,0};
297 //     if (fFormulaX) rrphiz[0] = fFormulaX(xyz,fFixedParam->GetMatrixArray());
298 //     if (fFormulaY) rrphiz[1] = fFormulaY(xyz,fFixedParam->GetMatrixArray());
299 //     Double_t alpha = TMath::ATan2(y,x);
300 //     Double_t ca    = TMath::Cos(alpha);
301 //     Double_t sa    = TMath::Sin(alpha);
302 //     if (coord==0) return ca*rrphiz[0]-sa*rrphiz[1];
303 //     if (coord==1) return sa*rrphiz[0]+ca*rrphiz[1];
304 //   }
305 //   return 0;
306 // }
307
308
309 Double_t AliTPCTransformation::GetDeltaXYZ(Int_t coord, Int_t volID, Double_t param, Double_t x, Double_t y, Double_t z){
310   //
311   //
312   // coord - type of coordinate
313   //       - 0 -X
314   //         1 -Y
315   //         2 -Z
316   //         3 -R
317   //         4 -RPhi
318   //         5 -Z
319   if (!fIsActive) return 0;
320   if (fBitMask && (!(*fBitMask)[volID])) return 0;
321   Double_t xyz[5]={x,y,z, param,volID};
322   Double_t alpha = TMath::ATan2(y,x);
323   Double_t ca    = TMath::Cos(alpha);
324   Double_t sa    = TMath::Sin(alpha);
325
326   if (fCoordSystem==0){
327     // cartezian system
328     Double_t dxdydz[3]={0,0,0};
329     if(fFormulaX) dxdydz[0]=fFormulaX(xyz,fFixedParam->GetMatrixArray());
330     if(fFormulaY) dxdydz[1]=fFormulaY(xyz,fFixedParam->GetMatrixArray());
331     if(fFormulaZ) dxdydz[2]=fFormulaZ(xyz,fFixedParam->GetMatrixArray());
332
333     if (coord==0) return dxdydz[0];
334     if (coord==1) return dxdydz[1];
335     if (coord==2) return dxdydz[2];
336     if (coord==3) return dxdydz[0]*ca+dxdydz[1]*sa;
337     if (coord==4) return -dxdydz[0]*sa+dxdydz[1]*ca;
338     if (coord==5) return dxdydz[2];
339   }
340   if (fCoordSystem==1){  
341     // cylindrical system
342     if (coord==2||coord==5) {
343       if (fFormulaZ==0) return 0;
344       return fFormulaZ(xyz,fFixedParam->GetMatrixArray());
345     }
346     Double_t rrphiz[3]={0,0,0};
347     if (fFormulaX) rrphiz[0] = fFormulaX(xyz,fFixedParam->GetMatrixArray());
348     if (fFormulaY) rrphiz[1] = fFormulaY(xyz,fFixedParam->GetMatrixArray());
349     Double_t alpha = TMath::ATan2(y,x);
350     Double_t ca    = TMath::Cos(alpha);
351     Double_t sa    = TMath::Sin(alpha);
352     if (coord==0) return ca*rrphiz[0]-sa*rrphiz[1];
353     if (coord==1) return sa*rrphiz[0]+ca*rrphiz[1];
354     if (coord==3) return rrphiz[0];
355     if (coord==4) return rrphiz[1];
356   }
357   return 0;
358 }
359
360
361
362 Double_t  AliTPCTransformation::TPCscalingRPol(Double_t *xyz, Double_t * param){
363   //
364   // Scaling and shift of TPC radius
365   // xyz[0..2] - global xyz of point 
366   // xyz[3]    - scale parameter
367   // param[0]  - radial scaling power
368   // param[1]  - drift  scaling power
369   // radius  from -1(at rInner)   to 1 (rOuter)
370   // driftM  from -1(at 0 drift)  to 1 (250 cm drift)
371
372   Double_t rInner=78.8;
373   Double_t rOuter=258.0; 
374   Double_t deltaR  = rOuter-rInner;  
375   Double_t radius  = (TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1])-rInner)*2./deltaR; 
376   Double_t driftM  = (0.5 - TMath::Abs(xyz[2]/250.))*2.0;
377   Double_t delta   = TMath::Power(radius,param[0])*TMath::Power(driftM,param[1]);
378   return delta*xyz[3];
379 }
380
381
382 Double_t  AliTPCTransformation::TPCscalingZDrift(Double_t *xyz, Double_t * param){
383   //
384   //
385   // Scaling and shift of TPC radius
386   // xyz[0..2] - global xyz of point 
387   // xyz[3]    - scale parameter
388   Double_t driftP  = TMath::Power(1. - TMath::Abs(xyz[2]/250.), param[0]);
389   Int_t    sector = TMath::Nint(xyz[4]);
390   Double_t deltaZ  = (sector%36<18) ? -driftP : driftP;
391   return deltaZ*xyz[3];
392 }
393
394 Double_t  AliTPCTransformation::TPCscalingZDriftT0(Double_t *xyz, Double_t * /*param*/){
395   //
396   //
397   // Z shift because time 0 offset
398   // opposite on A and C side
399   //
400   // xyz[0..2] - global xyz of point 
401   // xyz[3]    - scale parameter
402   Int_t    sector = TMath::Nint(xyz[4]);
403   Double_t sign  = (sector%36<18) ? -1 : 1;
404   return sign*xyz[3];
405 }
406
407
408 Double_t  AliTPCTransformation::TPCscalingZDriftGy(Double_t *xyz, Double_t * param){
409   //
410   //
411   // Scaling and shift of TPC radius
412   // xyz[0..2] - global xyz of point 
413   // xyz[3]    - scale parameter
414   Double_t driftP  = TMath::Power(1. - TMath::Abs(xyz[2]/250.), param[0]);
415   Double_t gy      = xyz[1]/250.;
416   Int_t    sector = TMath::Nint(xyz[4]);
417   Double_t deltaZ  = (sector%36<18) ? -driftP : driftP;
418   return deltaZ*xyz[3]*gy;
419 }
420
421
422
423 Double_t  AliTPCTransformation::TPCscalingPhiLocal(Double_t *xyz, Double_t * param){
424   //
425   //
426   // Scaling if the local y -phi
427   // xyz[0..2] - global xyz of point 
428   // xyz[3]    - scale parameter
429   // value = 1 for ful drift length and parameter 1
430   Double_t alpha       = TMath::ATan2(xyz[1],xyz[0]);
431   Double_t sector      = TMath::Nint(9*alpha/TMath::Pi()-0.5);
432   Double_t localAlpha  = (alpha-(sector+0.5)*TMath::Pi()/9.);
433   Double_t radius      = TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1])/250.;
434   //
435   Double_t deltaAlpha  = radius*TMath::Power(2.*9.*localAlpha/TMath::Pi(),param[0]);
436   return deltaAlpha*xyz[3];
437 }
438
439 Double_t  AliTPCTransformation::TPClocalRPhiEdge(Double_t *xyz, Double_t * param){
440   //
441   //
442   // Scaling if the local y -phi
443   // xyz[0..2] - global xyz of point 
444   // xyz[3]    - scale parameter
445   // param[0]  - dedge offset - should be around gap size/2.
446   // param[1]  - dedge factor - should be around gap size/2.
447   Double_t alpha       = TMath::ATan2(xyz[1],xyz[0]);
448   Double_t sector      = TMath::Nint(9*alpha/TMath::Pi()-0.5);
449   Double_t localAlpha  = (alpha-(sector+0.5)*TMath::Pi()/9.);
450   Double_t radius      = TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]);
451   Double_t deltaAlpha  = TMath::Pi()/18.-TMath::Abs(localAlpha);
452   Double_t distEdge    = (deltaAlpha*radius);
453   Double_t factor      = 1./(1.+(distEdge-param[0])/param[1]);
454   return factor*xyz[3]*((localAlpha>0)? -1.:1.);
455 }
456
457
458 Double_t       AliTPCTransformation::TPCscalingRIFC(Double_t *xyz, Double_t * param){
459   //
460   // inner field cage r distorion - proportinal to 1 over distance to the IFC
461   // param[0] - drift polynom order
462   // distortion at first pad row - is normalized to 
463   Double_t rInner=78.8;
464   Double_t rFirst=85.2; 
465   Double_t deltaR  = rFirst-rInner;
466   Double_t ndistR  = (TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1])-rInner)/deltaR;
467   Double_t driftM  = (0.5 - TMath::Abs(xyz[2]/250.))*2.;
468   Double_t value   = TMath::Power(driftM,param[0])/ndistR;
469   return xyz[3]*value;
470 }
471
472 Double_t       AliTPCTransformation::TPCscalingROFC(Double_t *xyz, Double_t * param){
473   //
474   // outer field cage r distorion - proportinal to 1 over distance to the OFC
475   // param[0] - drift polynom order
476   // driftM   - from -1 to 1 
477   //
478   Double_t rLast=245.8;
479   Double_t rOuter=258.0;  
480   Double_t deltaR  = rOuter-rLast;
481   Double_t ndistR  = (rOuter-TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]))/deltaR;
482   Double_t driftM  = (0.5 - TMath::Abs(xyz[2]/250.))*2.;
483   Double_t value   = TMath::Power(driftM,param[0])/ndistR;
484   return xyz[3]*value;
485 }
486
487
488 //
489 // TPC sector local misalignment 
490 //
491 //
492 //
493 Double_t AliTPCTransformation:: TPClocaldLxdGX(Double_t *xyz, Double_t * param){
494   //
495   // xyz - [0..2] - position 
496   //       [3]    - scale parameter
497   //       [4]    - volID
498   // param[0]= n  - cos(n *alpha)
499   // param[1]= n  - sin(n *alpha)
500   // param[2]     - indication - 0 - the same for IROC OROC 1 - opposite
501   // return delta in global coordiante system
502   //
503   Int_t    sector = TMath::Nint(xyz[4]);
504   Double_t alpha  = TMath::Pi()*(sector+0.5)/9;
505   Double_t ca     = TMath::Cos(alpha);
506   //  Double_t sa     = TMath::Sin(alpha); 
507   const Double_t xIROCOROC = 133.4;
508   Double_t factor = xyz[3];
509   if (param[0]>0)  factor*=TMath::Cos(alpha*param[0]);
510   if (param[1]>0)  factor*=TMath::Sin(alpha*param[1]);
511   if (param[2]>0.5 && TMath::Sqrt(xyz[1]*xyz[1]+xyz[0]*xyz[0])>xIROCOROC) factor*=-1;
512   return ca*factor;    
513 }
514
515 Double_t AliTPCTransformation::TPClocaldLxdGY(Double_t *xyz, Double_t * param){
516   //
517   // xyz - [0..2] - position 
518   //       [3]    - scale parameter
519   //       [4]    - volID
520   // param[0]= n  - cos(n *alpha)
521   // param[1]= n  - sin(n *alpha)
522   // param[2]     - indication - 0 - the same for IROC OROC 1 - opposite
523   // return delta in global coordiante system
524   //
525   Int_t    sector = TMath::Nint(xyz[4]);
526   Double_t alpha  = TMath::Pi()*(sector+0.5)/9;
527   //Double_t ca     = TMath::Cos(alpha);
528   Double_t sa     = TMath::Sin(alpha);
529   const Double_t xIROCOROC = 133.4;
530   Double_t factor = xyz[3];
531   if (param[0]>0)  factor*=TMath::Cos(alpha*param[0]);
532   if (param[1]>0)  factor*=TMath::Sin(alpha*param[1]);
533   if (param[2]>0.5 && TMath::Sqrt(xyz[1]*xyz[1]+xyz[0]*xyz[0])>xIROCOROC) factor*=-1;  
534   return   sa*factor;  
535 }
536
537 Double_t AliTPCTransformation:: TPClocaldLydGX(Double_t *xyz, Double_t * param){
538   //
539   // xyz - [0..2] - position 
540   //       [3]    - scale parameter
541   //       [4]    - volID
542   // param[0]= n  - cos(n *alpha)
543   // param[1]= n  - sin(n *alpha)
544   // param[2]     - indication - 0 - the same for IROC OROC 1 - opposite
545   // return delta in global coordiante system
546   //
547   Int_t    sector = TMath::Nint(xyz[4]);
548   Double_t alpha  = TMath::Pi()*(sector+0.5)/9;
549   //Double_t ca     = TMath::Cos(alpha);
550   Double_t sa     = TMath::Sin(alpha);
551   const Double_t xIROCOROC = 133.4;
552   Double_t factor = xyz[3];
553   if (param[0]>0)  factor*=TMath::Cos(alpha*param[0]);
554   if (param[1]>0)  factor*=TMath::Sin(alpha*param[1]);
555   if (param[2]>0.5 && TMath::Sqrt(xyz[1]*xyz[1]+xyz[0]*xyz[0])>xIROCOROC) factor*=-1;  
556   return            -sa*factor;  
557 }
558
559 Double_t AliTPCTransformation::TPClocaldLydGY(Double_t *xyz, Double_t * param){
560   //
561   // xyz - [0..2] - position 
562   //       [3]    - scale parameter
563   //       [4]    - volID
564   // param[0]= n  - cos(n *alpha)
565   // param[1]= n  - sin(n *alpha)
566   // param[2]     - indication - 0 - the same for IROC OROC 1 - opposite
567   // return delta in global coordiante system
568   //
569   Int_t    sector = TMath::Nint(xyz[4]);
570   Double_t alpha  = TMath::Pi()*(sector+0.5)/9;
571   Double_t ca     = TMath::Cos(alpha);
572   //Double_t sa     = TMath::Sin(alpha);
573   const Double_t xIROCOROC = 133.4;
574   Double_t factor = xyz[3];
575   if (param[0]>0)  factor*=TMath::Cos(alpha*param[0]);
576   if (param[1]>0)  factor*=TMath::Sin(alpha*param[1]);
577   if (param[2]>0.5 && TMath::Sqrt(xyz[1]*xyz[1]+xyz[0]*xyz[0])>xIROCOROC) factor*=-1;  
578   return   ca*factor;  
579 }
580
581
582 Double_t AliTPCTransformation::TPClocaldRzdGX(Double_t *xyz, Double_t * param){
583   //
584   // xyz - [0..2] - position 
585   //       [3]    - scale parameter - rotation angle in mrad
586   //       [4]    - volID
587   // param[0]= n  - cos(n *alpha)
588   // param[1]= n  - sin(n *alpha)
589   // param[2]     - indication - 0 - the same for IROC OROC 1 - opposite  
590   // return delta in global coordiante system
591   //
592   Int_t    sector = TMath::Nint(xyz[4]);
593   Double_t alpha  = TMath::Pi()*(sector+0.5)/9;
594   Double_t ca     = TMath::Cos(alpha);  
595   Double_t sa     = TMath::Sin(alpha);
596   Double_t lx     =  xyz[0]*ca + xyz[1]*sa;
597   Double_t ly     = -xyz[0]*sa + xyz[1]*ca;
598   //
599   const Double_t xIROCOROC = 133.4;  
600   lx-=xIROCOROC;
601   Double_t rot      =  xyz[3]*0.001;    // rotation in mrad
602   if (param[0]>0)  rot*=TMath::Cos(alpha*param[0]);
603   if (param[1]>0)  rot*=TMath::Sin(alpha*param[1]);
604   if (param[2]>0.5 && lx>0) rot*=-1;  
605   //
606   Double_t dlxR     =  - ly*rot; 
607   Double_t dlyR     =    lx*rot;
608   Double_t dgxR     =  dlxR*ca - dlyR*sa;
609   //Double_t dgyR     =  dlxR*sa + dlyR*ca;
610   return  dgxR;            
611 }
612
613 Double_t AliTPCTransformation::TPClocaldRzdGY(Double_t *xyz, Double_t * param){
614   //
615   // xyz - [0..2] - position 
616   //       [3]    - scale parameter - rotation angle in mrad
617   //       [4]    - volID
618   // param[0]= n  - cos(n *alpha)
619   // param[1]= n  - sin(n *alpha)
620   // param[2]     - indication - 0 - the same for IROC OROC 1 - opposite
621   // return delta in global coordiante system
622   //
623   Int_t    sector = TMath::Nint(xyz[4]);
624   Double_t alpha  = TMath::Pi()*(sector+0.5)/9;
625   Double_t ca     = TMath::Cos(alpha);  
626   Double_t sa     = TMath::Sin(alpha);
627   Double_t lx     =  xyz[0]*ca + xyz[1]*sa;
628   Double_t ly     = -xyz[0]*sa + xyz[1]*ca;
629   //
630   const Double_t xIROCOROC = 133.4;  
631   lx-=xIROCOROC;
632   Double_t rot      =  xyz[3]*0.001;    // rotation in mrad
633   if (param[0]>0)  rot*=TMath::Cos(alpha*param[0]);
634   if (param[1]>0)  rot*=TMath::Sin(alpha*param[1]);
635   if (param[2]>0.5 && lx>0) rot*=-1;  
636   Double_t dlxR     =  - ly*rot; 
637   Double_t dlyR     =    lx*rot;
638   //Double_t dgxR     =  dlxR*ca - dlyR*sa;
639   Double_t dgyR     =  dlxR*sa + dlyR*ca;
640   return  dgyR;            
641 }
642
643
644
645 Double_t        AliTPCTransformation::TPCDeltaZMediumLong(Double_t *xyz, Double_t * /*param*/){
646   //
647   // xyz - [0..2] - position 
648   //        [3]    - scale parameter
649   //        [4]    - volID
650   // return delta in global coordinate system 
651   //
652   Int_t    sector = TMath::Nint(xyz[4]);
653   Double_t signZ  = (sector%36<18) ? 1: -1;  // drift direction
654   if    (sector<36) return 0;     
655   //
656   const Double_t radiusLong = 198.1;
657   //
658   Double_t alpha  = TMath::Pi()*(sector+0.5)/9;
659   Double_t ca     = TMath::Cos(alpha);  
660   Double_t sa     = TMath::Sin(alpha);
661   Double_t lx     =  xyz[0]*ca + xyz[1]*sa;
662   Double_t sign   = (lx<radiusLong) ? 1:-1;
663   return xyz[3]*sign*signZ;
664 }
665
666 Double_t        AliTPCTransformation::TPCDeltaZ(Double_t *xyz, Double_t * param){
667   //
668   // xyz - [0..2] - position 
669   //        [3]    - scale parameter
670   //        [4]    - volID
671   // return delta in global coordiante system
672   //
673   Int_t    sector = TMath::Nint(xyz[4]);
674   Double_t delta  = (sector%36<18) ? 1: -1;  // drift direction
675   Double_t alpha  = TMath::Pi()*(sector+0.5)/9;
676   Double_t ca     = TMath::Cos(alpha);  
677   Double_t sa     = TMath::Sin(alpha);
678   Double_t lx     =  xyz[0]*ca + xyz[1]*sa;
679   //
680   const Double_t xIROCOROC = 133.4;  
681   if (param[0]>0) delta     *= TMath::Cos(param[0]*alpha);
682   if (param[1]>0) delta     *= TMath::Sin(param[1]*alpha);
683   if (param[2]>0.5 && lx >xIROCOROC) delta *=-1;
684   return delta*xyz[3];     // IROC shift
685 }
686
687
688 Double_t       AliTPCTransformation::TPCTiltingZ(Double_t *xyz, Double_t * param){
689   // xyz - [0..2] - position 
690   //        [3]    - scale parameter
691   //        [4]    - volID
692   // param[0]      - n for cos
693   // param[1]      - n for sin
694   // param[2]      - IROC-ORC relative (if >0.5 )
695   // return delta in global coordinate system 
696   const Double_t rFirst=85.2; 
697   const Double_t rLast =245.8;
698   const Double_t xIROCOROC = 133.4;  
699   //
700   Int_t    sector = TMath::Nint(xyz[4]);
701   Double_t alpha  = TMath::Pi()*(sector+0.5)/9;
702   Double_t ca     = TMath::Cos(alpha);  
703   Double_t sa     = TMath::Sin(alpha);
704   Double_t lx     =  xyz[0]*ca + xyz[1]*sa;
705   Double_t deltaR = 2.0*(lx-xIROCOROC)/(rLast-rFirst);  
706   if (param[0]>0) deltaR *= TMath::Cos(param[0]*alpha);
707   if (param[1]>0) deltaR *= TMath::Sin(param[1]*alpha);
708   if (param[2]>0.5 && lx >xIROCOROC) deltaR *=-1;
709   return deltaR*xyz[3];
710 }
711