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