]> git.uio.no Git - u/mrichter/AliRoot.git/blobdiff - TPC/AliTPCTransformation.cxx
Clean-up of include files
[u/mrichter/AliRoot.git] / TPC / AliTPCTransformation.cxx
index 9536c3a4a261547221f94e511e76bea9572feeb2..ba192319c96bc471c6b4e485b2b80d2297bc8ad4 100644 (file)
@@ -70,9 +70,14 @@ Int_t  AliTPCTransformation::BuildBasicFormulas(){
   RegisterFormula("TPCscalingRIFC",(GenFuncG)(AliTPCTransformation::TPCscalingRIFC));
   RegisterFormula("TPCscalingROFC",(GenFuncG)(AliTPCTransformation::TPCscalingROFC));
   //
-  RegisterFormula("TPCscalingZDr",(GenFuncG)(AliTPCTransformation::TPCscalingZDr));
-  RegisterFormula("TPCscalingZDrGy",(GenFuncG)(AliTPCTransformation::TPCscalingZDrGy));
+  RegisterFormula("TPCdeltaFCROC",(GenFuncG)(AliTPCTransformation::TPCdeltaFCROC));
+  RegisterFormula("TPCdeltaFCCE",(GenFuncG)(AliTPCTransformation::TPCdeltaFCCE));
+  //
+  RegisterFormula("TPCscalingZDr",(GenFuncG)(AliTPCTransformation::TPCscalingZDrift));
+  RegisterFormula("TPCscalingZDrGy",(GenFuncG)(AliTPCTransformation::TPCscalingZDriftGy));
+  RegisterFormula("TPCscalingZDriftT0",(GenFuncG)(AliTPCTransformation::TPCscalingZDriftT0));
   RegisterFormula("TPCscalingPhiLocal",(GenFuncG)(AliTPCTransformation::TPCscalingPhiLocal));
+  RegisterFormula("TPClocalRPhiEdge",(GenFuncG)(AliTPCTransformation::TPClocalRPhiEdge));
   //
   // TPC Local X and Y misalignment + rotation 
   //
@@ -126,6 +131,7 @@ AliTPCTransformation::AliTPCTransformation():
   fCoordSystem(0),    // coord system of  output deltas 
   fParam(0),          // free parameter of transformation 
   fSigma(0),          // uncertainty of the parameter
+  fSigmaMax(0),       //maximal sigma (Not allowed to increase in propagate time by bigger factor)
   fSigma2Time(0),     // change of the error in time - (For kalman filter) 
   fFixedParam(0),     // fixed parameters of tranformation  
   fIsActive(kTRUE),   // swith - On/Off
@@ -152,6 +158,7 @@ AliTPCTransformation::AliTPCTransformation(const char *name, TBits *mask, const
   fCoordSystem(coordSystem), // coordinate system of output deltas
   fParam(0),          // free parameter of transformation 
   fSigma(0),
+  fSigmaMax(0),       //maximal sigma (Not allowed to increase in propagate time by bigger factor)
   fSigma2Time(0),     // change of sigma in time
   fFixedParam(0),     // fixed parameters of tranformation  
   fIsActive(kTRUE),   // swith - On/Off
@@ -180,6 +187,7 @@ AliTPCTransformation::AliTPCTransformation(const AliTPCTransformation&trafo):
   fCoordSystem(0),    // coord system of  output deltas 
   fParam(trafo.fParam),          // free parameter of transformation 
   fSigma(trafo.fSigma),          // uncertainty of the parameter
+  fSigmaMax(trafo.fSigma),       //maximal sigma (Not allowed to increase in propagate time by bigger factor)
   fSigma2Time(trafo.fSigma2Time),     // change of the error in time - (For kalman filter) 
   fFixedParam(0),     // fixed parameters of tranformation
   fIsActive(trafo.fIsActive),   // swith - On/Off
@@ -189,6 +197,9 @@ AliTPCTransformation::AliTPCTransformation(const AliTPCTransformation&trafo):
   fFormulaY(0),       // y formula - pointer to the function
   fFormulaZ(0)       // z formula - pointer to the function
 {
+  //
+  // comment are above
+  //
   if (trafo.fNameX) fNameX = new TString(*(trafo.fNameX)); 
   if (trafo.fNameY) fNameY = new TString(*(trafo.fNameY)); 
   if (trafo.fNameZ) fNameZ = new TString(*(trafo.fNameZ)); 
@@ -205,12 +216,13 @@ AliTPCTransformation::~AliTPCTransformation(){
   delete fBitMask;
   delete fFixedParam;
 }
-void AliTPCTransformation::SetParams(Double_t param, Double_t sigma, Double_t sigma2Time, TVectorD* fixedParams){
+void AliTPCTransformation::SetParams(Double_t param, Double_t sigma, Double_t sigma2Time, const TVectorD *const fixedParams){
   //
   // Set parameters of transformation
   //
   fParam = param;
   fSigma = sigma;
+  fSigmaMax = sigma;
   fSigma2Time = sigma2Time;
   if (fFixedParam) delete fFixedParam;
   fFixedParam = new TVectorD(*fixedParams);
@@ -241,6 +253,8 @@ Bool_t AliTPCTransformation::Init(){
 
 
 TBits * AliTPCTransformation::BitsSide(Bool_t aside){
+  //
+  // Set bits for given side
   //
   TBits * bits = new TBits(72);
   for (Int_t i=0; i<72;i++){
@@ -254,7 +268,7 @@ TBits * AliTPCTransformation::BitsSide(Bool_t aside){
 
 TBits * AliTPCTransformation::BitsAll(){
   //
-  //
+  // Set all bits to kTRUE
   //
   TBits * bits = new TBits(72);
   for (Int_t i=0; i<72;i++){
@@ -263,6 +277,43 @@ TBits * AliTPCTransformation::BitsAll(){
   return bits;
 }
 
+// Double_t AliTPCTransformation::GetDeltaXYZ(Int_t coord, Int_t volID, Double_t param, Double_t x, Double_t y, Double_t z){
+//   //
+//   //
+//   // coord - type of coordinate
+//   //       - 0 -X
+//   //         1 -Y
+//   //         2 -Z
+//   //         3 -R
+//   //         4 -RPhi
+//   if (!fIsActive) return 0;
+//   if (fBitMask && (!(*fBitMask)[volID])) return 0;
+//   Double_t xyz[5]={x,y,z, param,volID};
+//   if (fCoordSystem==0){
+//     // cartezian system
+//     if (coord==0 && fFormulaX) return fFormulaX(xyz,fFixedParam->GetMatrixArray()); 
+//     if (coord==1 && fFormulaY) return fFormulaY(xyz,fFixedParam->GetMatrixArray()); 
+//     if (coord==2 && fFormulaZ) return fFormulaZ(xyz,fFixedParam->GetMatrixArray());
+//   }
+//   if (fCoordSystem==1){  
+//     // cylindrical system
+//     if (coord==2) {
+//       if (fFormulaZ==0) return 0;
+//       return fFormulaZ(xyz,fFixedParam->GetMatrixArray());
+//     }
+//     Double_t rrphiz[3]={0,0,0};
+//     if (fFormulaX) rrphiz[0] = fFormulaX(xyz,fFixedParam->GetMatrixArray());
+//     if (fFormulaY) rrphiz[1] = fFormulaY(xyz,fFixedParam->GetMatrixArray());
+//     Double_t alpha = TMath::ATan2(y,x);
+//     Double_t ca    = TMath::Cos(alpha);
+//     Double_t sa    = TMath::Sin(alpha);
+//     if (coord==0) return ca*rrphiz[0]-sa*rrphiz[1];
+//     if (coord==1) return sa*rrphiz[0]+ca*rrphiz[1];
+//   }
+//   return 0;
+// }
+
+
 Double_t AliTPCTransformation::GetDeltaXYZ(Int_t coord, Int_t volID, Double_t param, Double_t x, Double_t y, Double_t z){
   //
   //
@@ -272,34 +323,51 @@ Double_t AliTPCTransformation::GetDeltaXYZ(Int_t coord, Int_t volID, Double_t pa
   //         2 -Z
   //         3 -R
   //         4 -RPhi
+  //         5 -Z
   if (!fIsActive) return 0;
   if (fBitMask && (!(*fBitMask)[volID])) return 0;
   Double_t xyz[5]={x,y,z, param,volID};
+  Double_t alpha = TMath::ATan2(y,x);
+  Double_t ca    = TMath::Cos(alpha);
+  Double_t sa    = TMath::Sin(alpha);
+
   if (fCoordSystem==0){
     // cartezian system
-    if (coord==0 && fFormulaX) return fFormulaX(xyz,fFixedParam->GetMatrixArray()); 
-    if (coord==1 && fFormulaY) return fFormulaY(xyz,fFixedParam->GetMatrixArray()); 
-    if (coord==2 && fFormulaZ) return fFormulaZ(xyz,fFixedParam->GetMatrixArray());
+    Double_t dxdydz[3]={0,0,0};
+    if(fFormulaX) dxdydz[0]=fFormulaX(xyz,fFixedParam->GetMatrixArray());
+    if(fFormulaY) dxdydz[1]=fFormulaY(xyz,fFixedParam->GetMatrixArray());
+    if(fFormulaZ) dxdydz[2]=fFormulaZ(xyz,fFixedParam->GetMatrixArray());
+
+    if (coord==0) return dxdydz[0];
+    if (coord==1) return dxdydz[1];
+    if (coord==2) return dxdydz[2];
+    if (coord==3) return dxdydz[0]*ca+dxdydz[1]*sa;
+    if (coord==4) return -dxdydz[0]*sa+dxdydz[1]*ca;
+    if (coord==5) return dxdydz[2];
   }
   if (fCoordSystem==1){  
     // cylindrical system
-    if (coord==2) {
+    if (coord==2||coord==5) {
       if (fFormulaZ==0) return 0;
       return fFormulaZ(xyz,fFixedParam->GetMatrixArray());
     }
     Double_t rrphiz[3]={0,0,0};
     if (fFormulaX) rrphiz[0] = fFormulaX(xyz,fFixedParam->GetMatrixArray());
     if (fFormulaY) rrphiz[1] = fFormulaY(xyz,fFixedParam->GetMatrixArray());
-    Double_t alpha = TMath::ATan2(y,x);
-    Double_t ca    = TMath::Cos(alpha);
-    Double_t sa    = TMath::Sin(alpha);
+    alpha = TMath::ATan2(y,x);
+    ca    = TMath::Cos(alpha);
+    sa    = TMath::Sin(alpha);
     if (coord==0) return ca*rrphiz[0]-sa*rrphiz[1];
     if (coord==1) return sa*rrphiz[0]+ca*rrphiz[1];
+    if (coord==3) return rrphiz[0];
+    if (coord==4) return rrphiz[1];
   }
   return 0;
 }
 
-Double_t  AliTPCTransformation::TPCscalingRPol(Double_t *xyz, Double_t * param){
+
+
+Double_t  AliTPCTransformation::TPCscalingRPol(Double_t *xyz, const Double_t * const param){
   //
   // Scaling and shift of TPC radius
   // xyz[0..2] - global xyz of point 
@@ -319,19 +387,33 @@ Double_t  AliTPCTransformation::TPCscalingRPol(Double_t *xyz, Double_t * param){
 }
 
 
-Double_t  AliTPCTransformation::TPCscalingZDr(Double_t *xyz, Double_t * param){
+Double_t  AliTPCTransformation::TPCscalingZDrift(Double_t *xyz, const Double_t * const param){
   //
   //
   // Scaling and shift of TPC radius
   // xyz[0..2] - global xyz of point 
   // xyz[3]    - scale parameter
   Double_t driftP  = TMath::Power(1. - TMath::Abs(xyz[2]/250.), param[0]);
-  Double_t deltaZ  = (xyz[2]>0) ? -driftP : driftP;
+  Int_t    sector = TMath::Nint(xyz[4]);
+  Double_t deltaZ  = (sector%36<18) ? -driftP : driftP;
   return deltaZ*xyz[3];
 }
 
+Double_t  AliTPCTransformation::TPCscalingZDriftT0(Double_t *xyz, const Double_t * const /*param*/){
+  //
+  //
+  // Z shift because time 0 offset
+  // opposite on A and C side
+  //
+  // xyz[0..2] - global xyz of point 
+  // xyz[3]    - scale parameter
+  Int_t    sector = TMath::Nint(xyz[4]);
+  Double_t sign  = (sector%36<18) ? -1 : 1;
+  return sign*xyz[3];
+}
+
 
-Double_t  AliTPCTransformation::TPCscalingZDrGy(Double_t *xyz, Double_t * param){
+Double_t  AliTPCTransformation::TPCscalingZDriftGy(Double_t *xyz, const Double_t * const param){
   //
   //
   // Scaling and shift of TPC radius
@@ -339,13 +421,14 @@ Double_t  AliTPCTransformation::TPCscalingZDrGy(Double_t *xyz, Double_t * param)
   // xyz[3]    - scale parameter
   Double_t driftP  = TMath::Power(1. - TMath::Abs(xyz[2]/250.), param[0]);
   Double_t gy      = xyz[1]/250.;
-  Double_t deltaZ  = (xyz[2]>0) ? -driftP : driftP;
+  Int_t    sector = TMath::Nint(xyz[4]);
+  Double_t deltaZ  = (sector%36<18) ? -driftP : driftP;
   return deltaZ*xyz[3]*gy;
 }
 
 
 
-Double_t  AliTPCTransformation::TPCscalingPhiLocal(Double_t *xyz, Double_t * param){
+Double_t  AliTPCTransformation::TPCscalingPhiLocal(Double_t *xyz, const Double_t * const param){
   //
   //
   // Scaling if the local y -phi
@@ -361,8 +444,26 @@ Double_t  AliTPCTransformation::TPCscalingPhiLocal(Double_t *xyz, Double_t * par
   return deltaAlpha*xyz[3];
 }
 
+Double_t  AliTPCTransformation::TPClocalRPhiEdge(Double_t *xyz, const Double_t *const param){
+  //
+  //
+  // Scaling if the local y -phi
+  // xyz[0..2] - global xyz of point 
+  // xyz[3]    - scale parameter
+  // param[0]  - dedge offset - should be around gap size/2.
+  // param[1]  - dedge factor - should be around gap size/2.
+  Double_t alpha       = TMath::ATan2(xyz[1],xyz[0]);
+  Double_t sector      = TMath::Nint(9*alpha/TMath::Pi()-0.5);
+  Double_t localAlpha  = (alpha-(sector+0.5)*TMath::Pi()/9.);
+  Double_t radius      = TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]);
+  Double_t deltaAlpha  = TMath::Pi()/18.-TMath::Abs(localAlpha);
+  Double_t distEdge    = (deltaAlpha*radius);
+  Double_t factor      = 1./(1.+(distEdge-param[0])/param[1]);
+  return factor*xyz[3]*((localAlpha>0)? -1.:1.);
+}
+
 
-Double_t       AliTPCTransformation::TPCscalingRIFC(Double_t *xyz, Double_t * param){
+Double_t       AliTPCTransformation::TPCscalingRIFC(Double_t *xyz, const Double_t * const param){
   //
   // inner field cage r distorion - proportinal to 1 over distance to the IFC
   // param[0] - drift polynom order
@@ -376,7 +477,7 @@ Double_t       AliTPCTransformation::TPCscalingRIFC(Double_t *xyz, Double_t * pa
   return xyz[3]*value;
 }
 
-Double_t       AliTPCTransformation::TPCscalingROFC(Double_t *xyz, Double_t * param){
+Double_t       AliTPCTransformation::TPCscalingROFC(Double_t *xyz, const Double_t * const param){
   //
   // outer field cage r distorion - proportinal to 1 over distance to the OFC
   // param[0] - drift polynom order
@@ -392,12 +493,70 @@ Double_t       AliTPCTransformation::TPCscalingROFC(Double_t *xyz, Double_t * pa
 }
 
 
+Double_t       AliTPCTransformation::TPCdeltaFCROC(Double_t *xyz, const Double_t *const param){
+  // 
+  // delta R(Z) ROC induced
+  // param[0] - switch  0 - use distance to IFC - 1 - distance to IFC
+  // param[1] - kFC scaling factor  (multiplication factor  of (OFC-IFC))
+  // param[2] - kROC scaling factor 
+  // parameters [1] and [2] should be obtained from the electric field
+  //            simulation
+  //
+  Double_t rInner=78.8;
+  Double_t rFirst=85.2; 
+  Double_t rLast=245.8;
+  Double_t rOuter=258.0;  
+
+  Double_t radius  = TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]);
+  //calculate distance to the FC - inner or outer 
+  Double_t deltaFC = (param[0]<0.5)? TMath::Abs(radius-rFirst) : TMath::Abs(radius-rLast);
+  deltaFC/=(rOuter-rInner);
+  Double_t scalingFC = 1./(1.+deltaFC/(param[1]));
+  //
+  Double_t drift = 1.-TMath::Abs(xyz[2]/250.);  // normalized drift length
+  Double_t scalingROC = (1.-1./(1.+drift/param[2]));
+  //
+  return xyz[3]*scalingFC*scalingROC;
+}
+
+
+Double_t       AliTPCTransformation::TPCdeltaFCCE(Double_t *xyz, const Double_t *const param){
+  // 
+  // delta R(Z) CE (central electrode) induced
+  // param[0] - switch  0 - use distance to IFC - 1 - distance to IFC
+  // param[1] - kFC scaling factor  (multiplication factor  of (OFC-IFC))
+  // param[2] - kCE scaling factor 
+  // parameters [1] and [2] should be obtained from the electric field
+  //            simulation
+  Double_t rInner=78.8;
+  Double_t rFirst=85.2; 
+  Double_t rLast =245.8;
+  Double_t rOuter=258.0;  
+
+  Double_t radius  = TMath::Sqrt(xyz[0]*xyz[0]+xyz[1]*xyz[1]);
+  //calculate distance to the FC - inner or outer 
+  Double_t deltaFC = (param[0]<0.5)? TMath::Abs(radius-rFirst) : TMath::Abs(radius-rLast);
+  deltaFC/=(rOuter-rInner);
+  Double_t scalingFC = 1./(1.+deltaFC/(param[1]));
+  //
+  Double_t drift     = 1.-TMath::Abs(xyz[2]/250.);  // normalized drift length
+  Double_t scalingCE = 1/(1.+(1.-drift)/param[2]);  // 
+  //
+  return xyz[3]*scalingFC*scalingCE;
+}
+
+
+
+
+
+
+
 //
 // TPC sector local misalignment 
 //
 //
 //
-Double_t AliTPCTransformation:: TPClocaldLxdGX(Double_t *xyz, Double_t * param){
+Double_t AliTPCTransformation:: TPClocaldLxdGX(Double_t *xyz, const Double_t *const param){
   //
   // xyz - [0..2] - position 
   //       [3]    - scale parameter
@@ -419,7 +578,7 @@ Double_t AliTPCTransformation:: TPClocaldLxdGX(Double_t *xyz, Double_t * param){
   return ca*factor;    
 }
 
-Double_t AliTPCTransformation::TPClocaldLxdGY(Double_t *xyz, Double_t * param){
+Double_t AliTPCTransformation::TPClocaldLxdGY(Double_t *xyz, const Double_t *const param){
   //
   // xyz - [0..2] - position 
   //       [3]    - scale parameter
@@ -441,7 +600,7 @@ Double_t AliTPCTransformation::TPClocaldLxdGY(Double_t *xyz, Double_t * param){
   return   sa*factor;  
 }
 
-Double_t AliTPCTransformation:: TPClocaldLydGX(Double_t *xyz, Double_t * param){
+Double_t AliTPCTransformation:: TPClocaldLydGX(Double_t *xyz, const Double_t *const param){
   //
   // xyz - [0..2] - position 
   //       [3]    - scale parameter
@@ -463,7 +622,7 @@ Double_t AliTPCTransformation:: TPClocaldLydGX(Double_t *xyz, Double_t * param){
   return            -sa*factor;  
 }
 
-Double_t AliTPCTransformation::TPClocaldLydGY(Double_t *xyz, Double_t * param){
+Double_t AliTPCTransformation::TPClocaldLydGY(Double_t *xyz, const Double_t *const param){
   //
   // xyz - [0..2] - position 
   //       [3]    - scale parameter
@@ -486,7 +645,7 @@ Double_t AliTPCTransformation::TPClocaldLydGY(Double_t *xyz, Double_t * param){
 }
 
 
-Double_t AliTPCTransformation::TPClocaldRzdGX(Double_t *xyz, Double_t * param){
+Double_t AliTPCTransformation::TPClocaldRzdGX(Double_t *xyz, const Double_t *const param){
   //
   // xyz - [0..2] - position 
   //       [3]    - scale parameter - rotation angle in mrad
@@ -517,7 +676,7 @@ Double_t AliTPCTransformation::TPClocaldRzdGX(Double_t *xyz, Double_t * param){
   return  dgxR;            
 }
 
-Double_t AliTPCTransformation::TPClocaldRzdGY(Double_t *xyz, Double_t * param){
+Double_t AliTPCTransformation::TPClocaldRzdGY(Double_t *xyz, const Double_t *const param){
   //
   // xyz - [0..2] - position 
   //       [3]    - scale parameter - rotation angle in mrad
@@ -570,7 +729,7 @@ Double_t        AliTPCTransformation::TPCDeltaZMediumLong(Double_t *xyz, Double_
   return xyz[3]*sign*signZ;
 }
 
-Double_t        AliTPCTransformation::TPCDeltaZ(Double_t *xyz, Double_t * param){
+Double_t        AliTPCTransformation::TPCDeltaZ(Double_t *xyz, const Double_t *const param){
   //
   // xyz - [0..2] - position 
   //        [3]    - scale parameter
@@ -592,7 +751,7 @@ Double_t        AliTPCTransformation::TPCDeltaZ(Double_t *xyz, Double_t * param)
 }
 
 
-Double_t       AliTPCTransformation::TPCTiltingZ(Double_t *xyz, Double_t * param){
+Double_t       AliTPCTransformation::TPCTiltingZ(Double_t *xyz, const Double_t *const param){
   // xyz - [0..2] - position 
   //        [3]    - scale parameter
   //        [4]    - volID