]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCROCVoltError3D.cxx
Additional protection
[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("c2","c2",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 TMatrixD(kNR,kNZ);  
105     fLookUpEphiOverEz[k] =  new TMatrixD(kNR,kNZ);
106     fLookUpDeltaEz[k]    =  new TMatrixD(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
207   if (!fInitLookUp) {
208     AliInfo("Lookup table was not initialized! Perform the inizialisation now ...");
209     InitROCVoltError3D();
210     return;
211   }
212
213   Int_t   order     = 1 ;    // FIXME: hardcoded? Linear interpolation = 1, Quadratic = 2         
214
215   Double_t intEr, intEphi, intDeltaEz;
216   Double_t r, phi, z ;
217   Int_t    sign;
218
219   r      =  TMath::Sqrt( x[0]*x[0] + x[1]*x[1] ) ;
220   phi    =  TMath::ATan2(x[1],x[0]) ;
221   if ( phi < 0 ) phi += TMath::TwoPi() ;                   // Table uses phi from 0 to 2*Pi
222   z      =  x[2] ;                                         // Create temporary copy of x[2]
223
224   if ( (roc%36) < 18 ) {
225     sign =  1;       // (TPC A side)
226   } else {
227     sign = -1;       // (TPC C side)
228   }
229   
230   if ( sign==1  && z <  fgkZOffSet ) z =  fgkZOffSet;    // Protect against discontinuity at CE
231   if ( sign==-1 && z > -fgkZOffSet ) z = -fgkZOffSet;    // Protect against discontinuity at CE
232   
233
234   if ( (sign==1 && z<0) || (sign==-1 && z>0) ) // just a consistency check
235     AliError("ROC number does not correspond to z coordinate! Calculation of distortions is most likely wrong!");
236
237   // Get the Er and Ephi field integrals plus the integral over DeltaEz 
238   intEr      = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi, 
239                                   fgkRList, fgkZList, fgkPhiList, fLookUpErOverEz  );
240   intEphi    = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi, 
241                                   fgkRList, fgkZList, fgkPhiList, fLookUpEphiOverEz);
242   intDeltaEz = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi, 
243                                   fgkRList, fgkZList, fgkPhiList, fLookUpDeltaEz   );
244
245   //  printf("%lf %lf %lf\n",intEr,intEphi,intDeltaEz);
246
247   // Calculate distorted position
248   if ( r > 0.0 ) {
249     phi =  phi + ( fC0*intEphi - fC1*intEr ) / r;      
250     r   =  r   + ( fC0*intEr   + fC1*intEphi );  
251   }
252   
253   // Calculate correction in cartesian coordinates
254   dx[0] = r * TMath::Cos(phi) - x[0];
255   dx[1] = r * TMath::Sin(phi) - x[1]; 
256   dx[2] = intDeltaEz;  // z distortion - (internally scaled with driftvelocity dependency 
257                        // on the Ez field plus the actual ROC misalignment (if set TRUE)
258
259
260   if (fElectronArrivalCorrection) {
261
262     // correction for the OROC (in average, a 0.014usec longer drift time
263     // due to different position of the anode wires) -> vd*dt -> 2.64*0.014 = 0.0369 cm
264     // FIXME: correction are token from Magboltz simulations
265     //        should be loaded from a database
266  
267     AliTPCROC * rocInfo = AliTPCROC::Instance();
268     Double_t rCrossingROC  =  (rocInfo->GetPadRowRadii(0,62)+rocInfo->GetPadRowRadii(36,0))/2;
269   
270     if (r>rCrossingROC) {
271       if (sign==1)
272         dx[2] += 0.0369; // A side - negative correction
273       else
274         dx[2] -= 0.0369; // C side - positive correction
275     }
276     
277   }
278
279 }
280
281 void AliTPCROCVoltError3D::InitROCVoltError3D() {
282   //
283   // Initialization of the Lookup table which contains the solutions of the 
284   // Dirichlet boundary problem
285   // Calculation of the single 3D-Poisson solver is done just if needed
286   // (see basic lookup tables in header file)
287   //
288
289   const Int_t   order       =    1  ;  // Linear interpolation = 1, Quadratic = 2  
290   const Float_t gridSizeR   =  (fgkOFCRadius-fgkIFCRadius) / (kRows-1) ;
291   const Float_t gridSizeZ   =  fgkTPCZ0 / (kColumns-1) ;
292   const Float_t gridSizePhi =  TMath::TwoPi() / ( 18.0 * kPhiSlicesPerSector);
293
294   // temporary arrays to create the boundary conditions
295   TMatrixD *arrayofArrayV[kPhiSlices], *arrayofCharge[kPhiSlices] ; 
296   TMatrixD *arrayofEroverEz[kPhiSlices], *arrayofEphioverEz[kPhiSlices], *arrayofDeltaEz[kPhiSlices] ; 
297
298   for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
299     arrayofArrayV[k]     =   new TMatrixD(kRows,kColumns) ;
300     arrayofCharge[k]     =   new TMatrixD(kRows,kColumns) ;
301     arrayofEroverEz[k]   =   new TMatrixD(kRows,kColumns) ;
302     arrayofEphioverEz[k] =   new TMatrixD(kRows,kColumns) ;
303     arrayofDeltaEz[k]    =   new TMatrixD(kRows,kColumns) ;
304   }
305   
306   // list of point as used in the poisson relation and the interpolation (during sum up)
307   Double_t  rlist[kRows], zedlist[kColumns] , philist[kPhiSlices];
308   for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
309     philist[k] =  gridSizePhi * k;
310     for ( Int_t i = 0 ; i < kRows ; i++ )    {
311       rlist[i] = fgkIFCRadius + i*gridSizeR ;
312       for ( Int_t j = 0 ; j < kColumns ; j++ ) { // Fill Vmatrix with Boundary Conditions
313         zedlist[j]  = j * gridSizeZ ;
314       }
315     }
316   }
317
318   // ==========================================================================
319   // Solve Poisson's equation in 3D cylindrical coordinates by relaxation technique
320   // Allow for different size grid spacing in R and Z directions
321   
322   const Int_t   symmetry = 0;
323  
324   // Set bondaries and solve Poisson's equation --------------------------
325   
326   if ( !fInitLookUp ) {
327     
328     AliInfo(Form("Solving the poisson equation (~ %d sec)",2*10*(int)(kPhiSlices/10)));
329     
330     for ( Int_t side = 0 ; side < 2 ; side++ ) {  // Solve Poisson3D twice; once for +Z and once for -Z
331       
332       for ( Int_t k = 0 ; k < kPhiSlices ; k++ )  {
333         TMatrixD &arrayV    =  *arrayofArrayV[k] ;
334         TMatrixD &charge    =  *arrayofCharge[k] ;
335         
336         //Fill arrays with initial conditions.  V on the boundary and Charge in the volume.
337         for ( Int_t i = 0 ; i < kRows ; i++ ) {
338           for ( Int_t j = 0 ; j < kColumns ; j++ ) {  // Fill Vmatrix with Boundary Conditions
339             arrayV(i,j) = 0.0 ; 
340             charge(i,j) = 0.0 ;
341
342             Float_t radius0 = rlist[i] ;
343             Float_t phi0    = gridSizePhi * k ;
344             
345             // To avoid problems at sector boundaries, use an average of +- 1 degree from actual phi location
346             if ( j == (kColumns-1) ) {
347               arrayV(i,j) = 0.5*  ( GetROCVoltOffset( side, radius0, phi0+0.02 ) + GetROCVoltOffset( side, radius0, phi0-0.02 ) ) ;
348
349               if (side==1) // C side
350                 arrayV(i,j) = -arrayV(i,j); // minus sign on the C side to allow a consistent usage of global z when setting the boundaries
351             }
352           }
353         }      
354         
355         for ( Int_t i = 1 ; i < kRows-1 ; i++ ) { 
356           for ( Int_t j = 1 ; j < kColumns-1 ; j++ ) {
357             charge(i,j)  =  0.0 ;
358           }
359         }
360       }      
361       
362       // Solve Poisson's equation in 3D cylindrical coordinates by relaxation technique
363       // Allow for different size grid spacing in R and Z directions
364       
365       PoissonRelaxation3D( arrayofArrayV, arrayofCharge, 
366                            arrayofEroverEz, arrayofEphioverEz, arrayofDeltaEz,
367                            kRows, kColumns, kPhiSlices, gridSizePhi, kIterations, 
368                            symmetry, fROCdisplacement) ;
369       
370       
371       //Interpolate results onto a custom grid which is used just for these calculations.
372       Double_t  r, phi, z ;
373       for ( Int_t k = 0 ; k < kNPhi ; k++ ) {
374         phi = fgkPhiList[k] ;
375         
376         TMatrixD &erOverEz   =  *fLookUpErOverEz[k]  ;
377         TMatrixD &ephiOverEz =  *fLookUpEphiOverEz[k];
378         TMatrixD &deltaEz    =  *fLookUpDeltaEz[k]   ;
379         
380         for ( Int_t j = 0 ; j < kNZ ; j++ ) {
381
382           z = TMath::Abs(fgkZList[j]) ;  // Symmetric solution in Z that depends only on ABS(Z)
383   
384           if ( side == 0 &&  fgkZList[j] < 0 ) continue; // Skip rest of this loop if on the wrong side
385           if ( side == 1 &&  fgkZList[j] > 0 ) continue; // Skip rest of this loop if on the wrong side
386           
387           for ( Int_t i = 0 ; i < kNR ; i++ ) { 
388             r = fgkRList[i] ;
389
390             // Interpolate basicLookup tables; once for each rod, then sum the results
391             erOverEz(i,j)   = Interpolate3DTable(order, r, z, phi, kRows, kColumns, kPhiSlices, 
392                                                  rlist, zedlist, philist, arrayofEroverEz  );
393             ephiOverEz(i,j) = Interpolate3DTable(order, r, z, phi, kRows, kColumns, kPhiSlices,
394                                                  rlist, zedlist, philist, arrayofEphioverEz);
395             deltaEz(i,j)    = Interpolate3DTable(order, r, z, phi, kRows, kColumns, kPhiSlices,
396                                                  rlist, zedlist, philist, arrayofDeltaEz  );
397
398             if (side == 1)  deltaEz(i,j) = -  deltaEz(i,j); // negative coordinate system on C side
399
400           } // end r loop
401         }// end z loop
402       }// end phi loop
403
404       if ( side == 0 ) AliInfo(" A side done");
405       if ( side == 1 ) AliInfo(" C side done");
406     } // end side loop
407   }
408   
409   // clear the temporary arrays lists
410   for ( Int_t k = 0 ; k < kPhiSlices ; k++ )  {
411     delete arrayofArrayV[k];
412     delete arrayofCharge[k];
413     delete arrayofEroverEz[k];  
414     delete arrayofEphioverEz[k];
415     delete arrayofDeltaEz[k];
416   }
417  
418
419   fInitLookUp = kTRUE;
420
421 }
422
423
424 Float_t AliTPCROCVoltError3D::GetROCVoltOffset(Int_t side, Float_t r0, Float_t phi0) {
425   // 
426   // Returns the dz alignment data (in voltage equivalents) at 
427   // the given position
428   //
429
430   Float_t xp = r0*TMath::Cos(phi0);
431   Float_t yp = r0*TMath::Sin(phi0);
432   
433   // phi0 should be between 0 and 2pi 
434   if (phi0<0) phi0+=TMath::TwoPi();
435   Int_t roc = (Int_t)TMath::Floor((TMath::RadToDeg()*phi0)/20);
436   if (side==1) roc+=18; // C side
437   if (r0>132) roc+=36;  // OROC 
438   
439   // linear-plane data:  z = z0 + kx*lx + ky*ly (rotation in local coordinates)
440   TMatrixD &fitData = *fdzDataLinFit;
441
442   // local coordinates                                                          
443   Double_t secAlpha = TMath::DegToRad()*(10.+20.*(((Int_t)roc)%18));
444   Float_t lx = xp*TMath::Cos(secAlpha)+yp*TMath::Sin(secAlpha);
445   Float_t ly = yp*TMath::Cos(secAlpha)-xp*TMath::Sin(secAlpha);
446
447   // reference of rotation in lx is at the intersection between OROC and IROC
448   // necessary, since the Fitprozedure is otherwise useless
449   
450   AliTPCROC * rocInfo = AliTPCROC::Instance();
451   Double_t lxRef  = (rocInfo->GetPadRowRadii(0,62)+rocInfo->GetPadRowRadii(36,0))/2;
452   
453   Float_t dz = fitData(roc,0)+fitData(roc,1)*(lx-lxRef) + fitData(roc,2)*ly; // value in cm
454
455   // aproximated Voltage-offset-aquivalent to the z misalignment
456   // (linearly scaled with the z position)
457   Double_t ezField = (fgkCathodeV-fgkGG)/fgkTPCZ0; // = ALICE Electric Field (V/cm) Magnitude ~ -400 V/cm; 
458   Float_t voltOff = dz*ezField;            // values in "Volt equivalents"
459
460   return voltOff;
461 }
462
463 TH2F * AliTPCROCVoltError3D::CreateHistoOfZAlignment(Int_t side, Int_t nx, Int_t ny) {
464   //
465   // return a simple histogramm containing the input to the poisson solver
466   // (z positions of the Read-out chambers, linearly interpolated)
467
468   char hname[100];
469   if (side==0) snprintf(hname,100,"survey_dz_Aside");
470   if (side==1) snprintf(hname,100,"survey_dz_Cside");
471
472   TH2F *h = new TH2F(hname,hname,nx,-250.,250.,ny,-250.,250.);
473
474   for (Int_t iy=1;iy<=ny;++iy) {
475     Double_t yp = h->GetYaxis()->GetBinCenter(iy);
476     for (Int_t ix=1;ix<=nx;++ix) {
477       Double_t xp = h->GetXaxis()->GetBinCenter(ix);
478     
479       Float_t r0 = TMath::Sqrt(xp*xp+yp*yp);
480       Float_t phi0 = TMath::ATan2(yp,xp); 
481    
482       Float_t dz = GetROCVoltOffset(side,r0,phi0); // in [volt]
483
484       Double_t ezField = (fgkCathodeV-fgkGG)/fgkTPCZ0; // = ALICE Electric Field (V/cm) Magnitude ~ -400 V/cm; 
485       dz = dz/ezField;    // in [cm]
486
487       if (85.<=r0 && r0<=245.) {
488         h->SetBinContent(ix,iy,dz); 
489       } else {
490         h->SetBinContent(ix,iy,0.);
491       }
492     }
493   }
494   
495   h->GetXaxis()->SetTitle("x [cm]");
496   h->GetYaxis()->SetTitle("y [cm]");
497   h->GetZaxis()->SetTitle("dz [cm]");
498   h->SetStats(0);
499   //  h->DrawCopy("colz");
500
501   return h;
502
503
504 void AliTPCROCVoltError3D::Print(const Option_t* option) const {
505   //
506   // Print function to check the settings of the Rod shifts and the rotated clips
507   // option=="a" prints the C0 and C1 coefficents for calibration purposes
508   //
509
510   TString opt = option; opt.ToLower();
511   printf("%s\n",GetTitle());
512   printf(" - z aligmnet of the TPC Read-Out chambers \n");
513   printf("   (linear interpolation within the chamber:  dz = z0 + kx*(lx-133) + ky*ly [cm] ) \n");
514   printf("   Info: Check the following data-file for more details: %s \n",fROCDataFileName.Data());
515
516   if (opt.Contains("a")) { // Print all details
517     TMatrixD &fitData = *fdzDataLinFit;
518     printf(" A side:  IROC   ROCX=(z0,kx,ky): \n");
519     for (Int_t roc = 0; roc<18; roc++) 
520       printf("ROC%d:(%.2e,%.2e,%.2e) ",roc,fitData(roc,0),fitData(roc,1),fitData(roc,2));
521     printf("\n A side:  OROC   ROCX=(z0,kx,ky): \n");
522     for (Int_t roc = 36; roc<54; roc++) 
523       printf("ROC%d:(%.2e,%.2e,%.2e) ",roc,fitData(roc,0),fitData(roc,1),fitData(roc,2));
524     printf("\n C side:  IROC   ROCX=(z0,kx,ky): \n");
525     for (Int_t roc = 18; roc<36; roc++) 
526       printf("ROC%d:(%.2e,%.2e,%.2e) ",roc,fitData(roc,0),fitData(roc,1),fitData(roc,2));
527     printf("\n C side:  OROC   ROCX=(z0,kx,ky): \n");
528     for (Int_t roc = 54; roc<72; roc++) 
529       printf("ROC%d:(%.2e,%.2e,%.2e) ",roc,fitData(roc,0),fitData(roc,1),fitData(roc,2));
530     printf("\n\n");
531     printf(" - T1: %1.4f, T2: %1.4f \n",fT1,fT2);
532     printf(" - C1: %1.4f, C0: %1.4f \n",fC1,fC0);
533   }
534
535   if (!fInitLookUp) AliError("Lookup table was not initialized! You should do InitROCVoltError3D() ...");
536
537 }