]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCBoundaryVoltError.cxx
Coverity fix
[u/mrichter/AliRoot.git] / TPC / AliTPCBoundaryVoltError.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> AliTPCBoundaryVoltError class   </h2>       
20 //   This class calculates the space point distortions due to residual voltage errors on 
21 //   the main boundaries of the TPC. For example, the inner vessel of the TPC is shifted 
22 //   by a certain amount, whereas the ROCs on the A side and the C side follow this mechanical 
23 //   shift (at the inner vessel) in z direction. This example can be named "conical deformation"
24 //   of the TPC field cage (see example below).                    
25 //   <p>
26 //   The boundary conditions can be set via two arrays (A and C side) which contain the 
27 //   residual voltage setting modeling a possible shift or an inhomogeneity on the TPC field 
28 //   cage. In order to avoid that the user splits the Central Electrode (CE), the settings for 
29 //   the C side is taken from the array on the A side (points: A.6 and A.7). The region betweem
30 //    the points is interpolated linearly onto the boundaries.
31 //   <p>
32 //   The class uses the PoissonRelaxation2D (see AliTPCCorrection) to calculate the resulting 
33 //   electrical field inhomogeneities in the (r,z)-plane. Then, the Langevin-integral formalism 
34 //   is used to calculate the space point distortions. <br>
35 //   Note: This class assumes a homogeneous magnetic field. 
36 //   <p>
37 //   One has two possibilities when calculating the $dz$ distortions. The resulting distortions 
38 //   are purely due to the change of the drift velocity (along with the change of the drift field) 
39 //   when the SetROCDisplacement is FALSE. This emulates for example a Gating-Grid Voltage offset 
40 //   without moving the ROCs. When the flag is set to TRUE, the ROCs are assumed to be misaligned 
41 //   and the corresponding offset in z is added.
42 // End_Html
43 //
44 // Begin_Macro(source) 
45 //   {
46 //   gROOT->SetStyle("Plain"); gStyle->SetPalette(1);
47 //   TCanvas *c2 = new TCanvas("c2","c2",500,300); 
48 //   AliTPCBoundaryVoltError bve;     
49 //   Float_t val = 40;// [Volt]; 40V corresponds to 1mm
50 //   /* IFC shift, CE follows, ROC follows by factor half */
51 //   Float_t boundA[8] = { val, val, val,0,0,0,0,val}; // voltages A-side
52 //   Float_t boundC[6] = {-val,-val,-val,0,0,0};       // voltages C-side  
53 //   bve.SetBoundariesA(boundA);                                        
54 //   bve.SetBoundariesC(boundC);                                        
55 //   bve.SetOmegaTauT1T2(-0.32,1,1);
56 //   bve.SetROCDisplacement(kTRUE); // include the chamber offset in z when calculating the dz distortions
57 //   bve.CreateHistoDRinZR(0)->Draw("surf2"); 
58 //   return c2;
59 //   } 
60 // End_Macro
61 //
62 // Begin_Html
63 //   <p>
64 // Date: 01/06/2010 <br>
65 // Authors: Jim Thomas, Stefan Rossegger      
66 // End_Html 
67 // _________________________________________________________________
68
69
70 #include "AliMagF.h"
71 #include "TGeoGlobalMagField.h"
72 #include "AliTPCcalibDB.h"
73 #include "AliTPCParam.h"
74 #include "AliLog.h"
75 #include "TMatrixD.h"
76
77 #include "TMath.h"
78 #include "AliTPCROC.h"
79 #include "AliTPCBoundaryVoltError.h"
80
81 ClassImp(AliTPCBoundaryVoltError)
82
83 AliTPCBoundaryVoltError::AliTPCBoundaryVoltError()
84   : AliTPCCorrection("BoundaryVoltError","Boundary Voltage Error"),
85     fC0(0.),fC1(0.),
86     fROCdisplacement(kTRUE),
87     fInitLookUp(kFALSE)
88 {
89   //
90   // default constructor
91   //
92   for (Int_t i=0; i<8; i++){
93     fBoundariesA[i]= 0;  
94     if (i<6) fBoundariesC[i]= 0;
95   }
96 }
97
98 AliTPCBoundaryVoltError::~AliTPCBoundaryVoltError() {
99   //
100   // default destructor
101   //
102 }
103
104
105
106 void AliTPCBoundaryVoltError::Init() {
107   //
108   // Initialization funtion
109   //
110   
111   AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
112   if (!magF) AliError("Magneticd field - not initialized");
113   Double_t bzField = magF->SolenoidField()/10.; //field in T
114   AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
115   if (!param) AliError("Parameters - not initialized");
116   Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us]   // From dataBase: to be updated: per second (ideally)
117   Double_t ezField = 400; // [V/cm]   // to be updated: never (hopefully)
118   Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ; 
119   // Correction Terms for effective omegaTau; obtained by a laser calibration run
120   SetOmegaTauT1T2(wt,fT1,fT2);
121
122   InitBoundaryVoltErrorDistortion();
123 }
124
125 void AliTPCBoundaryVoltError::Update(const TTimeStamp &/*timeStamp*/) {
126   //
127   // Update function 
128   //
129   AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
130   if (!magF) AliError("Magneticd field - not initialized");
131   Double_t bzField = magF->SolenoidField()/10.; //field in T
132   AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
133   if (!param) AliError("Parameters - not initialized");
134   Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us]  // From dataBase: to be updated: per second (ideally)
135   Double_t ezField = 400; // [V/cm]   // to be updated: never (hopefully)
136   Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ; 
137   // Correction Terms for effective omegaTau; obtained by a laser calibration run
138   SetOmegaTauT1T2(wt,fT1,fT2);
139
140 }
141
142
143
144 void AliTPCBoundaryVoltError::GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]) {
145   //
146   // Calculates the correction due e.g. residual voltage errors on the TPC boundaries
147   //   
148
149   if (!fInitLookUp) {
150     AliInfo("Lookup table was not initialized!  Perform the inizialisation now ...");
151     InitBoundaryVoltErrorDistortion();
152   }
153
154   Int_t   order     = 1 ;               // FIXME: hardcoded? Linear interpolation = 1, Quadratic = 2         
155                                         // note that the poisson solution was linearly mirroed on this grid!
156   Double_t intEr, intEphi, intdEz ;
157   Double_t r, phi, z ;
158   Int_t    sign;
159
160   r      =  TMath::Sqrt( x[0]*x[0] + x[1]*x[1] ) ;
161   phi    =  TMath::ATan2(x[1],x[0]) ;
162   if ( phi < 0 ) phi += TMath::TwoPi() ;                   // Table uses phi from 0 to 2*Pi
163   z      =  x[2] ;                                         // Create temporary copy of x[2]
164
165   if ( (roc%36) < 18 ) {
166     sign =  1;       // (TPC A side)
167   } else {
168     sign = -1;       // (TPC C side)
169   }
170   
171   if ( sign==1  && z <  fgkZOffSet ) z =  fgkZOffSet;    // Protect against discontinuity at CE
172   if ( sign==-1 && z > -fgkZOffSet ) z = -fgkZOffSet;    // Protect against discontinuity at CE
173   
174
175   intEphi = 0.0;  // Efield is symmetric in phi - 2D calculation
176
177   if ( (sign==1 && z<0) || (sign==-1 && z>0) ) // just a consistency check
178     AliError("ROC number does not correspond to z coordinate! Calculation of distortions is most likely wrong!");
179
180   // Get the E field integral
181   Interpolate2DEdistortion( order, r, z, fLookUpErOverEz, intEr );
182   // Get DeltaEz field integral
183   Interpolate2DEdistortion( order, r, z, fLookUpDeltaEz, intdEz );
184   
185   // Calculate distorted position
186   if ( r > 0.0 ) {
187     phi =  phi + ( fC0*intEphi - fC1*intEr ) / r;      
188     r   =  r   + ( fC0*intEr   + fC1*intEphi );  
189   }
190   
191   // Calculate correction in cartesian coordinates
192   dx[0] = r * TMath::Cos(phi) - x[0];
193   dx[1] = r * TMath::Sin(phi) - x[1]; 
194   dx[2] = intdEz;  // z distortion - (internally scaled with driftvelocity dependency 
195                    // on the Ez field plus the actual ROC misalignment (if set TRUE)
196
197
198 }
199
200 void AliTPCBoundaryVoltError::InitBoundaryVoltErrorDistortion() {
201   //
202   // Initialization of the Lookup table which contains the solutions of the 
203   // Dirichlet boundary problem
204   //
205
206   const Float_t  gridSizeR   =  (fgkOFCRadius-fgkIFCRadius) / (kRows-1) ;
207   const Float_t  gridSizeZ   =  fgkTPCZ0 / (kColumns-1) ;
208
209   TMatrixD voltArrayA(kRows,kColumns), voltArrayC(kRows,kColumns); // boundary vectors
210   TMatrixD chargeDensity(kRows,kColumns);                              // dummy charge
211   TMatrixD arrayErOverEzA(kRows,kColumns), arrayErOverEzC(kRows,kColumns); // solution
212   TMatrixD arrayDeltaEzA(kRows,kColumns),  arrayDeltaEzC(kRows,kColumns); // solution
213
214   Double_t  rList[kRows], zedList[kColumns] ;
215   
216   // Fill arrays with initial conditions.  V on the boundary and ChargeDensity in the volume.      
217   for ( Int_t j = 0 ; j < kColumns ; j++ ) {
218     Double_t zed = j*gridSizeZ ;
219     zedList[j] = zed ;
220     for ( Int_t i = 0 ; i < kRows ; i++ )  {
221       Double_t radius = fgkIFCRadius + i*gridSizeR ;
222       rList[i]           = radius ;
223       voltArrayA(i,j)        = 0;  // Initialize voltArrayA to zero
224       voltArrayC(i,j)        = 0;  // Initialize voltArrayC to zero
225       chargeDensity(i,j)     = 0;  // Initialize ChargeDensity to zero - not used in this class
226     }
227   }      
228
229
230   // check if boundary values are the same for the C side (for later, saving some calculation time)
231
232   Int_t symmetry = -1; // assume that  A and C side are identical (but anti-symmetric!) // e.g conical deformation
233   Int_t sVec[8];
234
235   // check if boundaries are different (regardless the sign)
236   for (Int_t i=0; i<8; i++) { 
237     if (TMath::Abs(TMath::Abs(fBoundariesA[i]) - TMath::Abs(fBoundariesC[i])) > 1e-5) 
238       symmetry = 0;  
239     sVec[i] = (Int_t)( TMath::Sign((Float_t)1.,fBoundariesA[i]) * TMath::Sign((Float_t)1.,fBoundariesC[i])); 
240   }
241   if (symmetry==-1) { // still the same values?
242     // check the kind of symmetry , if even ...
243     if (sVec[0]==1 && sVec[1]==1 && sVec[2]==1 && sVec[3]==1 && sVec[4]==1 && sVec[5]==1 && sVec[6]==1 && sVec[7]==1 ) 
244       symmetry =  1;
245     else if (sVec[0]==-1 && sVec[1]==-1 && sVec[2]==-1 && sVec[3]==-1 && sVec[4]==-1 && sVec[5]==-1 && sVec[6]==-1 && sVec[7]==-1 ) 
246       symmetry = -1;
247     else
248       symmetry =  0; // some of the values differ in the sign -> neither symmetric nor antisymmetric
249   }
250
251
252
253   // Solve the electrosatic problem in 2D 
254
255   // Fill the complete Boundary vectors
256   // Start at IFC at CE and work anti-clockwise through IFC, ROC, OFC, and CE (clockwise for C side)
257   for ( Int_t j = 0 ; j < kColumns ; j++ ) {
258     Double_t zed = j*gridSizeZ ;
259     for ( Int_t i = 0 ; i < kRows ; i++ ) { 
260       Double_t radius = fgkIFCRadius + i*gridSizeR ;
261
262       // A side boundary vectors
263       if ( i == 0 ) voltArrayA(i,j) += zed   *((fBoundariesA[1]-fBoundariesA[0])/((kColumns-1)*gridSizeZ))
264         + fBoundariesA[0] ; // IFC
265       if ( j == kColumns-1 ) voltArrayA(i,j) += (radius-fgkIFCRadius)*((fBoundariesA[3]-fBoundariesA[2])/((kRows-1)*gridSizeR))
266         + fBoundariesA[2] ; // ROC
267       if ( i == kRows-1 ) voltArrayA(i,j) += zed   *((fBoundariesA[4]-fBoundariesA[5])/((kColumns-1)*gridSizeZ))
268         + fBoundariesA[5] ; // OFC
269       if ( j == 0 ) voltArrayA(i,j) += (radius-fgkIFCRadius)*((fBoundariesA[6]-fBoundariesA[7])/((kRows-1)*gridSizeR))
270         + fBoundariesA[7] ; // CE
271       
272       if (symmetry==0) {
273         // C side boundary vectors
274         if ( i == 0 ) voltArrayC(i,j) += zed   *((fBoundariesC[1]-fBoundariesC[0])/((kColumns-1)*gridSizeZ))
275           + fBoundariesC[0] ; // IFC
276         if ( j == kColumns-1 ) voltArrayC(i,j) += (radius-fgkIFCRadius)*((fBoundariesC[3]-fBoundariesC[2])/((kRows-1)*gridSizeR))
277           + fBoundariesC[2] ; // ROC
278         if ( i == kRows-1 ) voltArrayC(i,j) += zed   *((fBoundariesC[4]-fBoundariesC[5])/((kColumns-1)*gridSizeZ))
279           + fBoundariesC[5] ; // OFC
280         if ( j == 0 ) voltArrayC(i,j) += (radius-fgkIFCRadius)*((fBoundariesC[6]-fBoundariesC[7])/((kRows-1)*gridSizeR))
281           + fBoundariesC[7] ; // CE
282       }
283
284     }
285   }
286
287   voltArrayA(0,0)               *= 0.5 ; // Use average boundary condition at corner
288   voltArrayA(kRows-1,0)         *= 0.5 ; // Use average boundary condition at corner
289   voltArrayA(0,kColumns-1)      *= 0.5 ; // Use average boundary condition at corner
290   voltArrayA(kRows-1,kColumns-1)*= 0.5 ; // Use average boundary condition at corner
291
292   if (symmetry==0) {
293     voltArrayC(0,0)               *= 0.5 ; // Use average boundary condition at corner
294     voltArrayC(kRows-1,0)         *= 0.5 ; // Use average boundary condition at corner
295     voltArrayC(0,kColumns-1)      *= 0.5 ; // Use average boundary condition at corner
296     voltArrayC(kRows-1,kColumns-1)*= 0.5 ; // Use average boundary condition at corner
297   }
298
299
300   // always solve the problem on the A side
301   PoissonRelaxation2D( voltArrayA, chargeDensity, arrayErOverEzA, arrayDeltaEzA, 
302                        kRows, kColumns, kIterations, fROCdisplacement ) ;
303
304   if (symmetry!=0) { // A and C side are the same ("anti-symmetric" or "symmetric")
305     for ( Int_t j = 0 ; j < kColumns ; j++ ) {
306       for ( Int_t i = 0 ; i < kRows ; i++ ) { 
307         arrayErOverEzC(i,j) = symmetry*arrayErOverEzA(i,j);
308         arrayDeltaEzC(i,j) = -symmetry*arrayDeltaEzA(i,j);
309       }
310     }
311   } else if (symmetry==0) { // A and C side are different - Solve the problem on the C side too
312     PoissonRelaxation2D( voltArrayC, chargeDensity, arrayErOverEzC, arrayDeltaEzC,
313                          kRows, kColumns, kIterations, fROCdisplacement ) ;
314     for ( Int_t j = 0 ; j < kColumns ; j++ ) {
315       for ( Int_t i = 0 ; i < kRows ; i++ ) { 
316         arrayDeltaEzC(i,j) = -arrayDeltaEzC(i,j); // negative z coordinate!
317       }
318     }
319   }
320
321   // Interpolate results onto standard grid for Electric Fields
322   Int_t ilow=0, jlow=0 ;
323   Double_t z,r;
324   Float_t saveEr[2] ;         
325   Float_t saveEz[2] ;         
326   for ( Int_t i = 0 ; i < kNZ ; ++i )  {
327     z = TMath::Abs( fgkZList[i] ) ;
328     for ( Int_t j = 0 ; j < kNR ; ++j ) {
329       // Linear interpolation !!
330       r = fgkRList[j] ;
331       Search( kRows,   rList, r, ilow ) ;          // Note switch - R in rows and Z in columns
332       Search( kColumns, zedList, z, jlow ) ;
333       if ( ilow < 0 ) ilow = 0 ;                   // check if out of range
334       if ( jlow < 0 ) jlow = 0 ;   
335       if ( ilow + 1  >=  kRows - 1 ) ilow =  kRows - 2 ;              
336       if ( jlow + 1  >=  kColumns - 1 ) jlow =  kColumns - 2 ; 
337       
338       if (fgkZList[i]>0) {         // A side solution
339         saveEr[0] = arrayErOverEzA(ilow,jlow) + 
340           (arrayErOverEzA(ilow,jlow+1)-arrayErOverEzA(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
341         saveEr[1] = arrayErOverEzA(ilow+1,jlow) + 
342           (arrayErOverEzA(ilow+1,jlow+1)-arrayErOverEzA(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;
343         saveEz[0] = arrayDeltaEzA(ilow,jlow) + 
344           (arrayDeltaEzA(ilow,jlow+1)-arrayDeltaEzA(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
345         saveEz[1] = arrayDeltaEzA(ilow+1,jlow) + 
346           (arrayDeltaEzA(ilow+1,jlow+1)-arrayDeltaEzA(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;
347
348       } else if (fgkZList[i]<0) {  // C side solution
349         saveEr[0] = arrayErOverEzC(ilow,jlow) + 
350           (arrayErOverEzC(ilow,jlow+1)-arrayErOverEzC(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
351         saveEr[1] = arrayErOverEzC(ilow+1,jlow) + 
352           (arrayErOverEzC(ilow+1,jlow+1)-arrayErOverEzC(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;
353         saveEz[0] = arrayDeltaEzC(ilow,jlow) + 
354           (arrayDeltaEzC(ilow,jlow+1)-arrayDeltaEzC(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
355         saveEz[1] = arrayDeltaEzC(ilow+1,jlow) + 
356           (arrayDeltaEzC(ilow+1,jlow+1)-arrayDeltaEzC(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;
357
358       } else {
359         AliWarning("Field calculation at z=0 (CE) is not allowed!");
360         saveEr[0]=0; saveEr[1]=0;
361         saveEz[0]=0; saveEz[1]=0;
362       }
363       fLookUpErOverEz[i][j] = saveEr[0] + (saveEr[1]-saveEr[0])*(r-rList[ilow])/gridSizeR ;
364       fLookUpDeltaEz[i][j]  = saveEz[0] + (saveEz[1]-saveEz[0])*(r-rList[ilow])/gridSizeR ;
365     }
366   }
367   
368   voltArrayA.Clear();
369   voltArrayC.Clear();
370   chargeDensity.Clear();
371   arrayErOverEzA.Clear();
372   arrayErOverEzC.Clear();
373   arrayDeltaEzA.Clear();
374   arrayDeltaEzC.Clear();
375   
376   fInitLookUp = kTRUE;
377
378 }
379
380 void AliTPCBoundaryVoltError::Print(const Option_t* option) const {
381   //
382   // Print function to check the settings of the boundary vectors
383   // option=="a" prints the C0 and C1 coefficents for calibration purposes
384   //
385
386   TString opt = option; opt.ToLower();
387   printf("%s\n",GetTitle());
388   printf(" - Voltage settings (on the TPC boundaries) - linearly interpolated\n");
389   printf("  : A-side (anti-clockwise)\n"); 
390   printf("     (0,1):\t IFC (CE) : %3.1f V \t IFC (ROC): %3.1f V \n",fBoundariesA[0],fBoundariesA[1]);
391   printf("     (2,3):\t ROC (IFC): %3.1f V \t ROC (OFC): %3.1f V \n",fBoundariesA[2],fBoundariesA[3]);
392   printf("     (4,5):\t OFC (ROC): %3.1f V \t OFC (CE) : %3.1f V \n",fBoundariesA[4],fBoundariesA[5]);
393   printf("     (6,7):\t CE  (OFC): %3.1f V \t CE  (IFC): %3.1f V \n",fBoundariesA[6],fBoundariesA[7]);
394   printf("  : C-side (clockwise)\n"); 
395   printf("     (0,1):\t IFC (CE) : %3.1f V \t IFC (ROC): %3.1f V \n",fBoundariesC[0],fBoundariesC[1]);
396   printf("     (2,3):\t ROC (IFC): %3.1f V \t ROC (OFC): %3.1f V \n",fBoundariesC[2],fBoundariesC[3]);
397   printf("     (4,5):\t OFC (ROC): %3.1f V \t OFC (CE) : %3.1f V \n",fBoundariesC[4],fBoundariesC[5]);
398   printf("     (6,7):\t CE  (OFC): %3.1f V \t CE  (IFC): %3.1f V \n",fBoundariesC[6],fBoundariesC[7]);
399
400   // Check wether the settings of the Central Electrode agree (on the A and C side)
401   // Note: they have to be antisymmetric!
402   if (( TMath::Abs(fBoundariesA[6]+fBoundariesC[6])>1e-5) || ( TMath::Abs(fBoundariesA[7]+fBoundariesC[7])>1e-5 ) ){
403     AliWarning("Boundary parameters for the Central Electrode (CE) are not anti-symmetric! HOW DID YOU MANAGE THAT?");
404     AliWarning("Congratulations, you just splitted the Central Electrode of the TPC!");
405     AliWarning("Non-physical settings of the boundary parameter at the Central Electrode");
406   }
407
408   if (opt.Contains("a")) { // Print all details
409     printf(" - T1: %1.4f, T2: %1.4f \n",fT1,fT2);
410     printf(" - C1: %1.4f, C0: %1.4f \n",fC1,fC0);
411   } 
412    
413   if (!fInitLookUp) 
414     AliError("Lookup table was not initialized! You should do InitBoundaryVoltErrorDistortion() ...");
415
416 }
417
418
419 void AliTPCBoundaryVoltError::SetBoundariesA(Float_t boundariesA[8]){
420   //
421   // set voltage errors on the TPC boundaries - A side 
422   //
423   // Start at IFC at the Central electrode and work anti-clockwise (clockwise for C side) through 
424   // IFC, ROC, OFC, and CE. The boundary conditions are currently defined to be a linear 
425   // interpolation between pairs of numbers in the Boundary (e.g. fBoundariesA) vector.  
426   // The first pair of numbers represent the beginning and end of the Inner Field cage, etc.
427   // The unit of the error potential vector is [Volt], whereas 1mm shift of the IFC would 
428   // correspond to ~ 40 V
429   // 
430   // Note: The setting for the CE will be passed to the C side!
431   
432   for (Int_t i=0; i<8; i++) {
433     fBoundariesA[i]= boundariesA[i];  
434     if (i>5) fBoundariesC[i]= -boundariesA[i]; // setting for the CE is passed to C side
435   }
436   fInitLookUp=kFALSE;
437 }
438 void AliTPCBoundaryVoltError::SetBoundariesC(Float_t boundariesC[6]){
439   //
440   // set voltage errors on the TPC boundaries - C side 
441   //
442   // Start at IFC at the Central electrode and work clockwise (for C side) through 
443   // IFC, ROC and OFC. The boundary conditions are currently defined to be a linear 
444   // interpolation between pairs of numbers in the Boundary (e.g. fBoundariesC) vector.  
445   // The first pair of numbers represent the beginning and end of the Inner Field cage, etc.
446   // The unit of the error potential vector is [Volt], whereas 1mm shift of the IFC would 
447   // correspond to ~ 40 V
448   // 
449   // Note: The setting for the CE will be taken from the A side (pos 6 and 7)!
450
451   for (Int_t i=0; i<6; i++) {
452     fBoundariesC[i]= boundariesC[i];  
453   }
454   fInitLookUp=kFALSE;
455 }