]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/Base/AliTPCROCVoltError3D.cxx
ATO-98 AliTPCRocVoltError3D::AddCorrectionCompact implmented and test included
[u/mrichter/AliRoot.git] / TPC / Base / AliTPCROCVoltError3D.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 // Begin_Html
19 //   <h2> AliTPCROCVoltError3D class   </h2>    
20 //   The class calculates the space point distortions due to z offsets of the 
21 //   ROCs via the residual voltage technique (Poisson relaxation) in 3D. 
22 //   Since the GG (part of the ROCs) represents the closure of the FC in z direction,
23 //   every misalignment in z produces not only dz distortions but also electrical 
24 //   field inhomogeneities throughout the volume, which produces additional dr and rd$\phi$ distortions.
25 //   <p>
26 //   Each ROC can be misaligned (in z direction) in three ways. A general z0 offset, 
27 //   an inclination along the x and an inclination along the y axis. The z-misalignment's
28 //   can be set via the function SetROCData(TMatrixD *mat) for each single chamber (ROC). 
29 //   The array size has to be (72,3) representing the 72 chambers according to the 
30 //   offline numbering scheme (IROC: roc$<$36; OROC: roc$\geq$36) and the three misalignment's
31 //   (see the source code for further details).
32 //   <p>
33 //   Internally, these z offsets (unit is cm)  are recalculated into residual voltage 
34 //   equivalents in order to make use of the relaxation technique. 
35 //   <p>
36 //   One has two possibilities when calculating the $dz$ distortions. The resulting 
37 //   distortions are purely due to the change of the drift velocity (along with the 
38 //   change of the drift field) when the SetROCDisplacement is FALSE. <br>
39 //   For this class, this is a rather unphysical setting and should be avoided. <br>
40 //   When the flag is set to TRUE, the corresponding offset in z is added to the dz 
41 //   calculation of the outer ROCs. <br>
42 //   For the Alice TPC gas, both effects are of similar magnitude. This means, if the 
43 //   drift length is sufficiently large, a z-offset of a chamber appears to have (approx.) 
44 //   twice the magnitude when one looks at the actual dz distortions.
45 //   <p>
46 //   In addition, this class allows a correction regarding the different arrival times 
47 //   of the electrons due to the geometrical difference of the inner and outer chambers.
48 //   The impact was simulated via Garfield. If the flag is set via the 
49 //   function SetElectronArrivalCorrection, the electron-arrival correction is added to the dz calculation.
50 // End_Html
51 //
52 // Begin_Macro(source)
53 //   {
54 //   gROOT->SetStyle("Plain"); gStyle->SetPalette(1);
55 //   TCanvas *c2 = new TCanvas("cAliTPCROCVoltError3D","cAliTPCROCVoltError3D",500,400); 
56 //   AliTPCROCVoltError3D roc; 
57 //   roc.SetROCDataFileName("$ALICE_ROOT/TPC/Calib/maps/TPCROCdzSurvey.root");
58 //   roc.SetElectronArrivalCorrection(kFALSE);  // Correction for electron arrival offset, IROC vs OROC
59 //   roc.SetROCDisplacement(kTRUE);   // include the chamber offset in z when calculating the dz 
60 //   roc.SetOmegaTauT1T2(0,1,1); // B=0
61 //   roc.CreateHistoDZinXY(1.,300,300)->Draw("colz"); 
62 //   return c2;
63 //   } 
64 // End_Macro
65 //
66 // Begin_Html
67 //   <p>
68 //   Date: 08/08/2010    <br>                                                 
69 //   Authors: Jim Thomas, Stefan Rossegger                                
70 // End_Html 
71 // _________________________________________________________________
72
73
74 #include "AliMagF.h"
75 #include "TGeoGlobalMagField.h"
76 #include "AliTPCcalibDB.h"
77 #include "AliTPCParam.h"
78 #include "AliLog.h"
79 #include "TMatrixD.h"
80 #include "TFile.h"
81
82 #include "TMath.h"
83 #include "AliTPCROC.h"
84 #include "AliTPCROCVoltError3D.h"
85
86 ClassImp(AliTPCROCVoltError3D)
87
88 AliTPCROCVoltError3D::AliTPCROCVoltError3D()
89   : AliTPCCorrection("ROCVoltErrors","ROC z alignment Errors"),
90     fC0(0.),fC1(0.),
91     fROCdisplacement(kTRUE),
92     fElectronArrivalCorrection(kTRUE),
93     fInitLookUp(kFALSE),
94     fROCDataFileName(""),  
95     fdzDataLinFit(0)
96 {
97   //
98   // default constructor
99   //
100
101   // Array which will contain the solution according to the setted boundary conditions
102   // main input: z alignment of the Read Out chambers
103   // see InitROCVoltError3D() function
104   for ( Int_t k = 0 ; k < kNPhi ; k++ ) {
105     fLookUpErOverEz[k]   =  new TMatrixF(kNR,kNZ);  
106     fLookUpEphiOverEz[k] =  new TMatrixF(kNR,kNZ);
107     fLookUpDeltaEz[k]    =  new TMatrixF(kNR,kNZ);   
108   }
109   //  fROCDataFileName="$ALICE_ROOT/TPC/Calib/maps/TPCROCdzSurvey.root";
110   //  SetROCDataFileName(fROCDataFileName.Data()); // initialization of fdzDataLinFit is included
111
112 }
113
114 AliTPCROCVoltError3D::~AliTPCROCVoltError3D() {
115   //
116   // destructor
117   //
118   
119   for ( Int_t k = 0 ; k < kNPhi ; k++ ) {
120     delete fLookUpErOverEz[k];
121     delete fLookUpEphiOverEz[k];
122     delete fLookUpDeltaEz[k];
123   }
124
125   delete fdzDataLinFit;
126 }
127
128
129
130 Bool_t AliTPCROCVoltError3D::AddCorrectionCompact(AliTPCCorrection* corr, Double_t weight){
131   //
132   // Add correction  and make them compact
133   // Assumptions:
134   //  - origin of distortion/correction are additive
135   //  - only correction ot the same type supported ()
136   if (corr==NULL) {
137     AliError("Zerro pointer - correction");
138     return kFALSE;
139   }  
140   AliTPCROCVoltError3D * corrC = dynamic_cast<AliTPCROCVoltError3D *>(corr);
141   if (corrC == NULL) {
142     AliError(TString::Format("Inconsistent class types: %s\%s",IsA()->GetName(),corr->IsA()->GetName()).Data());
143     return kFALSE;
144   }
145   //
146   TMatrixD matrixDz=*(corrC->fdzDataLinFit);
147   matrixDz*=weight;
148   if (fdzDataLinFit==NULL) {
149     fdzDataLinFit=new TMatrixD(matrixDz);
150   }
151   else{
152     (*fdzDataLinFit)+=matrixDz;
153   }
154   return kTRUE;
155 }
156
157
158
159
160
161 void AliTPCROCVoltError3D::SetROCData(TMatrixD * matrix){
162   //
163   // Set a z alignment map of the chambers not via a file, but
164   // directly via a TMatrix(72,3), where dz = p0 + p1*(lx-133.4) + p2*ly (all in cm)
165   //
166   if (!fdzDataLinFit) fdzDataLinFit=new TMatrixD(*matrix);
167   else *fdzDataLinFit = *matrix;
168 }
169
170
171 void AliTPCROCVoltError3D::Init() {
172   //
173   // Initialization funtion
174   //
175   
176   AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
177   if (!magF) AliError("Magneticd field - not initialized");
178   Double_t bzField = magF->SolenoidField()/10.; //field in T
179   AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
180   if (!param) AliError("Parameters - not initialized");
181   Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us]   // From dataBase: to be updated: per second (ideally)
182   Double_t ezField = 400; // [V/cm]   // to be updated: never (hopefully)
183   Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ; 
184   // Correction Terms for effective omegaTau; obtained by a laser calibration run
185   SetOmegaTauT1T2(wt,fT1,fT2);
186
187   if (!fInitLookUp) InitROCVoltError3D();
188 }
189
190 void AliTPCROCVoltError3D::Update(const TTimeStamp &/*timeStamp*/) {
191   //
192   // Update function 
193   //
194   AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
195   if (!magF) AliError("Magneticd field - not initialized");
196   Double_t bzField = magF->SolenoidField()/10.; //field in T
197   AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
198   if (!param) AliError("Parameters - not initialized");
199   Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us]  // From dataBase: to be updated: per second (ideally)
200   Double_t ezField = 400; // [V/cm]   // to be updated: never (hopefully)
201   Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ; 
202   // Correction Terms for effective omegaTau; obtained by a laser calibration run
203   SetOmegaTauT1T2(wt,fT1,fT2);
204
205 }
206
207 void  AliTPCROCVoltError3D::SetROCDataFileName(const char * fname) {
208   //
209   // Set / load the ROC data (linear fit of ROC misalignments)
210   //
211
212   fROCDataFileName = fname;
213   
214   TFile f(fROCDataFileName.Data(),"READ");
215   TMatrixD *m = (TMatrixD*) f.Get("dzSurveyLinFitData");
216   TMatrixD &mf = *m;
217
218   // prepare some space
219
220   if (fdzDataLinFit) delete fdzDataLinFit;
221   fdzDataLinFit = new TMatrixD(72,3);
222   TMatrixD &dataIntern = *fdzDataLinFit;
223   
224   for (Int_t iroc=0;iroc<72;iroc++) {
225     dataIntern(iroc,0) = mf(iroc,0);  // z0 offset
226     dataIntern(iroc,1) = mf(iroc,1);  // slope in x
227     dataIntern(iroc,2) = mf(iroc,2);  // slope in y
228   }
229
230   f.Close();
231
232   fInitLookUp = kFALSE;
233
234 }
235
236 void AliTPCROCVoltError3D::GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]) {
237   //
238   // Calculates the correction due e.g. residual voltage errors on the TPC boundaries
239   //   
240   const Double_t kEpsilon=Double_t(FLT_MIN);
241   if (!fInitLookUp) {
242     AliInfo("Lookup table was not initialized! Perform the inizialisation now ...");
243     InitROCVoltError3D();
244   }
245   static Bool_t forceInit=kTRUE; //temporary needed for back compatibility with old OCDB entries
246   if (forceInit&&fLookUpErOverEz[0]){
247     if (TMath::Abs(fLookUpErOverEz[0]->Sum())<kEpsilon){//temporary needed for back compatibility with old OCDB entries
248       ForceInitROCVoltError3D();
249     }
250     forceInit=kFALSE;
251   }
252
253   
254   Int_t   order     = 1 ;    // FIXME: hardcoded? Linear interpolation = 1, Quadratic = 2         
255
256   Float_t intEr, intEphi, intDeltaEz;
257   Double_t r, phi, z ;
258   Int_t    sign;
259
260   r      =  TMath::Sqrt( x[0]*x[0] + x[1]*x[1] ) ;
261   phi    =  TMath::ATan2(x[1],x[0]) ;
262   if ( phi < 0 ) phi += TMath::TwoPi() ;                   // Table uses phi from 0 to 2*Pi
263   z      =  x[2] ;                                         // Create temporary copy of x[2]
264
265   if ( (roc%36) < 18 ) {
266     sign =  1;       // (TPC A side)
267   } else {
268     sign = -1;       // (TPC C side)
269   }
270   
271   if ( sign==1  && z <  fgkZOffSet ) z =  fgkZOffSet;    // Protect against discontinuity at CE
272   if ( sign==-1 && z > -fgkZOffSet ) z = -fgkZOffSet;    // Protect against discontinuity at CE
273   
274
275   if ( (sign==1 && z<0) || (sign==-1 && z>0) ) // just a consistency check
276     AliError("ROC number does not correspond to z coordinate! Calculation of distortions is most likely wrong!");
277
278   // Get the Er and Ephi field integrals plus the integral over DeltaEz 
279   intEr      = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi, 
280                                   fgkRList, fgkZList, fgkPhiList, fLookUpErOverEz  );
281   intEphi    = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi, 
282                                   fgkRList, fgkZList, fgkPhiList, fLookUpEphiOverEz);
283   intDeltaEz = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi, 
284                                   fgkRList, fgkZList, fgkPhiList, fLookUpDeltaEz   );
285
286   //  printf("%lf %lf %lf\n",intEr,intEphi,intDeltaEz);
287
288   // Calculate distorted position
289   if ( r > 0.0 ) {
290     phi =  phi + ( fC0*intEphi - fC1*intEr ) / r;      
291     r   =  r   + ( fC0*intEr   + fC1*intEphi );  
292   }
293   
294   // Calculate correction in cartesian coordinates
295   dx[0] = r * TMath::Cos(phi) - x[0];
296   dx[1] = r * TMath::Sin(phi) - x[1]; 
297   dx[2] = intDeltaEz;  // z distortion - (internally scaled with driftvelocity dependency 
298                        // on the Ez field plus the actual ROC misalignment (if set TRUE)
299
300
301   if (fElectronArrivalCorrection) {
302
303     // correction for the OROC (in average, a 0.014usec longer drift time
304     // due to different position of the anode wires) -> vd*dt -> 2.64*0.014 = 0.0369 cm
305     // FIXME: correction are token from Magboltz simulations
306     //        should be loaded from a database
307  
308     AliTPCROC * rocInfo = AliTPCROC::Instance();
309     Double_t rCrossingROC  =  (rocInfo->GetPadRowRadii(0,62)+rocInfo->GetPadRowRadii(36,0))/2;
310   
311     if (r>rCrossingROC) {
312       if (sign==1)
313         dx[2] += 0.0369; // A side - negative correction
314       else
315         dx[2] -= 0.0369; // C side - positive correction
316     }
317     
318   }
319
320 }
321
322 void AliTPCROCVoltError3D::InitROCVoltError3D() {
323   //
324   // Initialization of the Lookup table which contains the solutions of the 
325   // Dirichlet boundary problem
326   // Calculation of the single 3D-Poisson solver is done just if needed
327   // (see basic lookup tables in header file)
328   //
329
330   const Int_t   order       =    1  ;  // Linear interpolation = 1, Quadratic = 2  
331   const Float_t gridSizeR   =  (fgkOFCRadius-fgkIFCRadius) / (kRows-1) ;
332   const Float_t gridSizeZ   =  fgkTPCZ0 / (kColumns-1) ;
333   const Float_t gridSizePhi =  TMath::TwoPi() / ( 18.0 * kPhiSlicesPerSector);
334
335   // temporary arrays to create the boundary conditions
336   TMatrixD *arrayofArrayV[kPhiSlices], *arrayofCharge[kPhiSlices] ; 
337   TMatrixD *arrayofEroverEz[kPhiSlices], *arrayofEphioverEz[kPhiSlices], *arrayofDeltaEz[kPhiSlices] ; 
338
339   for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
340     arrayofArrayV[k]     =   new TMatrixD(kRows,kColumns) ;
341     arrayofCharge[k]     =   new TMatrixD(kRows,kColumns) ;
342     arrayofEroverEz[k]   =   new TMatrixD(kRows,kColumns) ;
343     arrayofEphioverEz[k] =   new TMatrixD(kRows,kColumns) ;
344     arrayofDeltaEz[k]    =   new TMatrixD(kRows,kColumns) ;
345   }
346   
347   // list of point as used in the poisson relation and the interpolation (during sum up)
348   Double_t  rlist[kRows], zedlist[kColumns] , philist[kPhiSlices];
349   for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
350     philist[k] =  gridSizePhi * k;
351     for ( Int_t i = 0 ; i < kRows ; i++ )    {
352       rlist[i] = fgkIFCRadius + i*gridSizeR ;
353       for ( Int_t j = 0 ; j < kColumns ; j++ ) { // Fill Vmatrix with Boundary Conditions
354         zedlist[j]  = j * gridSizeZ ;
355       }
356     }
357   }
358
359   // ==========================================================================
360   // Solve Poisson's equation in 3D cylindrical coordinates by relaxation technique
361   // Allow for different size grid spacing in R and Z directions
362   
363   const Int_t   symmetry = 0;
364  
365   // Set bondaries and solve Poisson's equation --------------------------
366   
367   if ( !fInitLookUp ) {
368     
369     AliInfo(Form("Solving the poisson equation (~ %d sec)",2*10*(int)(kPhiSlices/10)));
370     
371     for ( Int_t side = 0 ; side < 2 ; side++ ) {  // Solve Poisson3D twice; once for +Z and once for -Z
372       
373       for ( Int_t k = 0 ; k < kPhiSlices ; k++ )  {
374         TMatrixD &arrayV    =  *arrayofArrayV[k] ;
375         TMatrixD &charge    =  *arrayofCharge[k] ;
376         
377         //Fill arrays with initial conditions.  V on the boundary and Charge in the volume.
378         for ( Int_t i = 0 ; i < kRows ; i++ ) {
379           for ( Int_t j = 0 ; j < kColumns ; j++ ) {  // Fill Vmatrix with Boundary Conditions
380             arrayV(i,j) = 0.0 ; 
381             charge(i,j) = 0.0 ;
382
383             Float_t radius0 = rlist[i] ;
384             Float_t phi0    = gridSizePhi * k ;
385             
386             // To avoid problems at sector boundaries, use an average of +- 1 degree from actual phi location
387             if ( j == (kColumns-1) ) {
388               arrayV(i,j) = 0.5*  ( GetROCVoltOffset( side, radius0, phi0+0.02 ) + GetROCVoltOffset( side, radius0, phi0-0.02 ) ) ;
389
390               if (side==1) // C side
391                 arrayV(i,j) = -arrayV(i,j); // minus sign on the C side to allow a consistent usage of global z when setting the boundaries
392             }
393           }
394         }      
395         
396         for ( Int_t i = 1 ; i < kRows-1 ; i++ ) { 
397           for ( Int_t j = 1 ; j < kColumns-1 ; j++ ) {
398             charge(i,j)  =  0.0 ;
399           }
400         }
401       }      
402       
403       // Solve Poisson's equation in 3D cylindrical coordinates by relaxation technique
404       // Allow for different size grid spacing in R and Z directions
405       
406       PoissonRelaxation3D( arrayofArrayV, arrayofCharge, 
407                            arrayofEroverEz, arrayofEphioverEz, arrayofDeltaEz,
408                            kRows, kColumns, kPhiSlices, gridSizePhi, kIterations, 
409                            symmetry, fROCdisplacement) ;
410       
411       
412       //Interpolate results onto a custom grid which is used just for these calculations.
413       Double_t  r, phi, z ;
414       for ( Int_t k = 0 ; k < kNPhi ; k++ ) {
415         phi = fgkPhiList[k] ;
416         
417         TMatrixF &erOverEz   =  *fLookUpErOverEz[k]  ;
418         TMatrixF &ephiOverEz =  *fLookUpEphiOverEz[k];
419         TMatrixF &deltaEz    =  *fLookUpDeltaEz[k]   ;
420         
421         for ( Int_t j = 0 ; j < kNZ ; j++ ) {
422
423           z = TMath::Abs(fgkZList[j]) ;  // Symmetric solution in Z that depends only on ABS(Z)
424   
425           if ( side == 0 &&  fgkZList[j] < 0 ) continue; // Skip rest of this loop if on the wrong side
426           if ( side == 1 &&  fgkZList[j] > 0 ) continue; // Skip rest of this loop if on the wrong side
427           
428           for ( Int_t i = 0 ; i < kNR ; i++ ) { 
429             r = fgkRList[i] ;
430
431             // Interpolate basicLookup tables; once for each rod, then sum the results
432             erOverEz(i,j)   = Interpolate3DTable(order, r, z, phi, kRows, kColumns, kPhiSlices, 
433                                                  rlist, zedlist, philist, arrayofEroverEz  );
434             ephiOverEz(i,j) = Interpolate3DTable(order, r, z, phi, kRows, kColumns, kPhiSlices,
435                                                  rlist, zedlist, philist, arrayofEphioverEz);
436             deltaEz(i,j)    = Interpolate3DTable(order, r, z, phi, kRows, kColumns, kPhiSlices,
437                                                  rlist, zedlist, philist, arrayofDeltaEz  );
438
439             if (side == 1)  deltaEz(i,j) = -  deltaEz(i,j); // negative coordinate system on C side
440
441           } // end r loop
442         }// end z loop
443       }// end phi loop
444
445       if ( side == 0 ) AliInfo(" A side done");
446       if ( side == 1 ) AliInfo(" C side done");
447     } // end side loop
448   }
449   
450   // clear the temporary arrays lists
451   for ( Int_t k = 0 ; k < kPhiSlices ; k++ )  {
452     delete arrayofArrayV[k];
453     delete arrayofCharge[k];
454     delete arrayofEroverEz[k];  
455     delete arrayofEphioverEz[k];
456     delete arrayofDeltaEz[k];
457   }
458  
459
460   fInitLookUp = kTRUE;
461
462 }
463
464
465 Float_t AliTPCROCVoltError3D::GetROCVoltOffset(Int_t side, Float_t r0, Float_t phi0) {
466   // 
467   // Returns the dz alignment data (in voltage equivalents) at 
468   // the given position
469   //
470
471   Float_t xp = r0*TMath::Cos(phi0);
472   Float_t yp = r0*TMath::Sin(phi0);
473   
474   // phi0 should be between 0 and 2pi 
475   if (phi0<0) phi0+=TMath::TwoPi();
476   Int_t roc = (Int_t)TMath::Floor((TMath::RadToDeg()*phi0)/20);
477   if (side==1) roc+=18; // C side
478   if (r0>132) roc+=36;  // OROC 
479   
480   // linear-plane data:  z = z0 + kx*lx + ky*ly (rotation in local coordinates)
481   TMatrixD &fitData = *fdzDataLinFit;
482
483   // local coordinates                                                          
484   Double_t secAlpha = TMath::DegToRad()*(10.+20.*(((Int_t)roc)%18));
485   Float_t lx = xp*TMath::Cos(secAlpha)+yp*TMath::Sin(secAlpha);
486   Float_t ly = yp*TMath::Cos(secAlpha)-xp*TMath::Sin(secAlpha);
487
488   // reference of rotation in lx is at the intersection between OROC and IROC
489   // necessary, since the Fitprozedure is otherwise useless
490   
491   AliTPCROC * rocInfo = AliTPCROC::Instance();
492   Double_t lxRef  = (rocInfo->GetPadRowRadii(0,62)+rocInfo->GetPadRowRadii(36,0))/2;
493   
494   Float_t dz = fitData(roc,0)+fitData(roc,1)*(lx-lxRef) + fitData(roc,2)*ly; // value in cm
495
496   // aproximated Voltage-offset-aquivalent to the z misalignment
497   // (linearly scaled with the z position)
498   Double_t ezField = (fgkCathodeV-fgkGG)/fgkTPCZ0; // = ALICE Electric Field (V/cm) Magnitude ~ -400 V/cm; 
499   Float_t voltOff = dz*ezField;            // values in "Volt equivalents"
500
501   return voltOff;
502 }
503
504 TH2F * AliTPCROCVoltError3D::CreateHistoOfZAlignment(Int_t side, Int_t nx, Int_t ny) {
505   //
506   // return a simple histogramm containing the input to the poisson solver
507   // (z positions of the Read-out chambers, linearly interpolated)
508
509   char hname[100];
510   if (side==0) snprintf(hname,100,"survey_dz_Aside");
511   if (side==1) snprintf(hname,100,"survey_dz_Cside");
512
513   TH2F *h = new TH2F(hname,hname,nx,-250.,250.,ny,-250.,250.);
514
515   for (Int_t iy=1;iy<=ny;++iy) {
516     Double_t yp = h->GetYaxis()->GetBinCenter(iy);
517     for (Int_t ix=1;ix<=nx;++ix) {
518       Double_t xp = h->GetXaxis()->GetBinCenter(ix);
519     
520       Float_t r0 = TMath::Sqrt(xp*xp+yp*yp);
521       Float_t phi0 = TMath::ATan2(yp,xp); 
522    
523       Float_t dz = GetROCVoltOffset(side,r0,phi0); // in [volt]
524
525       Double_t ezField = (fgkCathodeV-fgkGG)/fgkTPCZ0; // = ALICE Electric Field (V/cm) Magnitude ~ -400 V/cm; 
526       dz = dz/ezField;    // in [cm]
527
528       if (85.<=r0 && r0<=245.) {
529         h->SetBinContent(ix,iy,dz); 
530       } else {
531         h->SetBinContent(ix,iy,0.);
532       }
533     }
534   }
535   
536   h->GetXaxis()->SetTitle("x [cm]");
537   h->GetYaxis()->SetTitle("y [cm]");
538   h->GetZaxis()->SetTitle("dz [cm]");
539   h->SetStats(0);
540   //  h->DrawCopy("colz");
541
542   return h;
543
544
545 void AliTPCROCVoltError3D::Print(const Option_t* option) const {
546   //
547   // Print function to check the settings of the Rod shifts and the rotated clips
548   // option=="a" prints the C0 and C1 coefficents for calibration purposes
549   //
550
551   TString opt = option; opt.ToLower();
552   printf("%s\n",GetTitle());
553   printf(" - z aligmnet of the TPC Read-Out chambers \n");
554   printf("   (linear interpolation within the chamber:  dz = z0 + kx*(lx-133) + ky*ly [cm] ) \n");
555   printf("   Info: Check the following data-file for more details: %s \n",fROCDataFileName.Data());
556
557   if (opt.Contains("a")) { // Print all details
558     TMatrixD &fitData = *fdzDataLinFit;
559     printf(" A side:  IROC   ROCX=(z0,kx,ky): \n");
560     for (Int_t roc = 0; roc<18; roc++) 
561       printf("ROC%d:(%.2e,%.2e,%.2e) ",roc,fitData(roc,0),fitData(roc,1),fitData(roc,2));
562     printf("\n A side:  OROC   ROCX=(z0,kx,ky): \n");
563     for (Int_t roc = 36; roc<54; roc++) 
564       printf("ROC%d:(%.2e,%.2e,%.2e) ",roc,fitData(roc,0),fitData(roc,1),fitData(roc,2));
565     printf("\n C side:  IROC   ROCX=(z0,kx,ky): \n");
566     for (Int_t roc = 18; roc<36; roc++) 
567       printf("ROC%d:(%.2e,%.2e,%.2e) ",roc,fitData(roc,0),fitData(roc,1),fitData(roc,2));
568     printf("\n C side:  OROC   ROCX=(z0,kx,ky): \n");
569     for (Int_t roc = 54; roc<72; roc++) 
570       printf("ROC%d:(%.2e,%.2e,%.2e) ",roc,fitData(roc,0),fitData(roc,1),fitData(roc,2));
571     printf("\n\n");
572     printf(" - T1: %1.4f, T2: %1.4f \n",fT1,fT2);
573     printf(" - C1: %1.4f, C0: %1.4f \n",fC1,fC0);
574   }
575
576   if (!fInitLookUp) AliError("Lookup table was not initialized! You should do InitROCVoltError3D() ...");
577
578 }