]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCFCVoltError3D.cxx
- Inclusion of AliTPCSpaceCharge3D in TPCbase library
[u/mrichter/AliRoot.git] / TPC / AliTPCFCVoltError3D.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>AliTPCFCVoltError3D class   </h2>       
20 //   The class calculates the space point distortions due to residual voltage errors 
21 //   on the Field Cage (FC) boundaries (rods and strips) of the TPC in 3D. It uses 
22 //   the Poisson relaxation technique in three dimension as implemented in the parent class. 
23 //   <p>
24 //   Although the calculation is performed in 3D, the calculation time was significantly 
25 //   reduced by using certain symmetry conditions within the calculation.
26 //   <p>
27 //   The input parameters can be set via the functions (e.g.) SetRodVoltShift(rod,dV) where 
28 //   rod is the number of the rod on which the voltage offset dV is set. The corresponding 
29 //   shift in z direction would be $dz=dV/40$ with an opposite sign for the C side. The 
30 //   rods are numbered in anti-clock-wise direction starting at $\phi=0$. Rods in the IFC 
31 //   are from 0 to 17, rods on the OFC go from 18 to 35. <br>
32 //   This convention is following the offline numbering scheme of the ROCs.
33 //   <p>
34 //   Different misalignment scenarios can be modeled: 
35 //   <ul>
36 //   <li> A rotated clip scenario is only possible at two positions (Rod 11 at IFC, rod 3(+18) at OFC) 
37 //        and can be set via SetRotatedClipVolt. The key feature is that at the mentioned positions,
38 //        the strip-ends were combined. At these positions, a systematic offset of one strip-end in
39 //        respect to the other is possible. 
40 //   <li> A normal rod offset, where the strips follow the movement of the rod, can be set via 
41 //        SetRodVoltShift. It has a anti-mirrored signature in terms of distortions when compared 
42 //        to the rotated clip. This misalignment is possible at each single rod of the FC.
43 //   <li> A simple rod offset, where the strips do not follow the shift, results in an even more 
44 //        localized distortion close to the rod. The difference to a rod shift, where the strips follow,
45 //        is more dominant on the OFC. This effect can be set via the function SetCopperRodShift.
46 //   </ul>
47 // End_Html
48 //
49 // Begin_Macro(source) 
50 //   {
51 //   gROOT->SetStyle("Plain"); gStyle->SetPalette(1);
52 //   TCanvas *c2 = new TCanvas("c2","c2",500,450); 
53 //   AliTPCFCVoltError3D fc;
54 //   fc.SetOmegaTauT1T2(0,1,1); 
55 //   fc.SetRotatedClipVoltA(0,40);
56 //   fc.SetRodVoltShiftA(3,40); 
57 //   fc.SetCopperRodShiftA(7+18,40);
58 //   fc.SetRodVoltShiftA(15+18,40); 
59 //   fc.CreateHistoDRPhiinXY(10)->Draw("cont4z");
60 //   return c2;
61 //   } 
62 // End_Macro
63 //
64 // Begin_Html
65 //   <p>
66 //   Date: 08/08/2010  <br>
67 //   Authors: Jim Thomas, Stefan Rossegger  
68 // End_Html 
69 // _________________________________________________________________
70
71
72 #include "AliMagF.h"
73 #include "TGeoGlobalMagField.h"
74 #include "AliTPCcalibDB.h"
75 #include "AliTPCParam.h"
76 #include "AliLog.h"
77 #include "TMatrixD.h"
78
79 #include "TMath.h"
80 #include "AliTPCROC.h"
81 #include "AliTPCFCVoltError3D.h"
82
83 ClassImp(AliTPCFCVoltError3D)
84
85 AliTPCFCVoltError3D::AliTPCFCVoltError3D()
86   : AliTPCCorrection("FieldCageVoltErrors","FieldCage (Rods) Voltage Errors"),
87     fC0(0.),fC1(0.),
88     fInitLookUp(kFALSE)
89 {
90   //
91   // default constructor
92   //
93
94   // flags for filled 'basic lookup tables'
95   for (Int_t i=0; i<6; i++){
96     fInitLookUpBasic[i]= kFALSE;  
97   }
98
99   // Boundary settings 
100   for (Int_t i=0; i<36; i++){
101     fRodVoltShiftA[i] = 0;  
102     fRodVoltShiftC[i] = 0;  
103   }
104   for (Int_t i=0; i<2; i++){
105     fRotatedClipVoltA[i] = 0;  
106     fRotatedClipVoltC[i] = 0;  
107   }
108   // 
109   for (Int_t i=0; i<36; i++){
110     fCopperRodShiftA[i] = 0;  
111     fCopperRodShiftC[i] = 0;  
112   }
113
114   // Array which will contain the solution according to the setted boundary conditions
115   // it represents a sum up of the 4 basic look up tables (created individually)
116   // see InitFCVoltError3D() function
117   for ( Int_t k = 0 ; k < kNPhi ; k++ ) {
118     fLookUpErOverEz[k]   =  new TMatrixD(kNR,kNZ);  
119     fLookUpEphiOverEz[k] =  new TMatrixD(kNR,kNZ);
120     fLookUpDeltaEz[k]    =  new TMatrixD(kNR,kNZ);   
121   }
122   
123   for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
124     fLookUpBasic1ErOverEz[k]   = 0;
125     fLookUpBasic1EphiOverEz[k] = 0; 
126     fLookUpBasic1DeltaEz[k]    = 0;
127
128     fLookUpBasic2ErOverEz[k]   = 0;
129     fLookUpBasic2EphiOverEz[k] = 0; 
130     fLookUpBasic2DeltaEz[k]    = 0;
131
132     fLookUpBasic3ErOverEz[k]   = 0;
133     fLookUpBasic3EphiOverEz[k] = 0; 
134     fLookUpBasic3DeltaEz[k]    = 0;
135
136     fLookUpBasic4ErOverEz[k]   = 0;
137     fLookUpBasic4EphiOverEz[k] = 0; 
138     fLookUpBasic4DeltaEz[k]    = 0;
139     
140     fLookUpBasic5ErOverEz[k]   = 0;
141     fLookUpBasic5EphiOverEz[k] = 0; 
142     fLookUpBasic5DeltaEz[k]    = 0;
143
144     fLookUpBasic6ErOverEz[k]   = 0;
145     fLookUpBasic6EphiOverEz[k] = 0; 
146     fLookUpBasic6DeltaEz[k]    = 0;
147   }
148
149 }
150
151 AliTPCFCVoltError3D::~AliTPCFCVoltError3D() {
152   //
153   // destructor
154   //
155   
156   for ( Int_t k = 0 ; k < kNPhi ; k++ ) {
157     delete fLookUpErOverEz[k];
158     delete fLookUpEphiOverEz[k];
159     delete fLookUpDeltaEz[k];
160   }
161
162   for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
163     delete fLookUpBasic1ErOverEz[k];  // does nothing if pointer is zero!
164     delete fLookUpBasic1EphiOverEz[k]; 
165     delete fLookUpBasic1DeltaEz[k]; 
166
167     delete fLookUpBasic2ErOverEz[k];  // does nothing if pointer is zero!
168     delete fLookUpBasic2EphiOverEz[k]; 
169     delete fLookUpBasic2DeltaEz[k]; 
170     
171     delete fLookUpBasic3ErOverEz[k];  // does nothing if pointer is zero!
172     delete fLookUpBasic3EphiOverEz[k]; 
173     delete fLookUpBasic3DeltaEz[k]; 
174
175     delete fLookUpBasic4ErOverEz[k];  // does nothing if pointer is zero!
176     delete fLookUpBasic4EphiOverEz[k]; 
177     delete fLookUpBasic4DeltaEz[k]; 
178
179     delete fLookUpBasic5ErOverEz[k];  // does nothing if pointer is zero!
180     delete fLookUpBasic5EphiOverEz[k]; 
181     delete fLookUpBasic5DeltaEz[k]; 
182
183     delete fLookUpBasic6ErOverEz[k];  // does nothing if pointer is zero!
184     delete fLookUpBasic6EphiOverEz[k]; 
185     delete fLookUpBasic6DeltaEz[k]; 
186
187   }
188 }
189
190 void AliTPCFCVoltError3D::Init() {
191   //
192   // Initialization funtion
193   //
194   
195   AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
196   if (!magF) AliError("Magneticd field - not initialized");
197   Double_t bzField = magF->SolenoidField()/10.; //field in T
198   AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
199   if (!param) AliError("Parameters - not initialized");
200   Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us]   // From dataBase: to be updated: per second (ideally)
201   Double_t ezField = 400; // [V/cm]   // to be updated: never (hopefully)
202   Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ; 
203   // Correction Terms for effective omegaTau; obtained by a laser calibration run
204   SetOmegaTauT1T2(wt,fT1,fT2);
205
206   if (!fInitLookUp) InitFCVoltError3D();
207 }
208
209 void AliTPCFCVoltError3D::Update(const TTimeStamp &/*timeStamp*/) {
210   //
211   // Update function 
212   //
213   AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
214   if (!magF) AliError("Magneticd field - not initialized");
215   Double_t bzField = magF->SolenoidField()/10.; //field in T
216   AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
217   if (!param) AliError("Parameters - not initialized");
218   Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us]  // From dataBase: to be updated: per second (ideally)
219   Double_t ezField = 400; // [V/cm]   // to be updated: never (hopefully)
220   Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ; 
221   // Correction Terms for effective omegaTau; obtained by a laser calibration run
222   SetOmegaTauT1T2(wt,fT1,fT2);
223
224
225 }
226
227
228
229 void AliTPCFCVoltError3D::GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]) {
230   //
231   // Calculates the correction due e.g. residual voltage errors on the TPC boundaries
232   //   
233
234   if (!fInitLookUp) {
235     AliInfo("Lookup table was not initialized! Perform the inizialisation now ...");
236     InitFCVoltError3D();
237   }
238
239   Int_t   order     = 1 ;               // FIXME: hardcoded? Linear interpolation = 1, Quadratic = 2         
240                                         // note that the poisson solution was linearly mirroed on this grid!
241   Double_t intEr, intEphi, intDeltaEz;
242   Double_t r, phi, z ;
243   Int_t    sign;
244
245   r      =  TMath::Sqrt( x[0]*x[0] + x[1]*x[1] ) ;
246   phi    =  TMath::ATan2(x[1],x[0]) ;
247   if ( phi < 0 ) phi += TMath::TwoPi() ;                   // Table uses phi from 0 to 2*Pi
248   z      =  x[2] ;                                         // Create temporary copy of x[2]
249
250   if ( (roc%36) < 18 ) {
251     sign =  1;       // (TPC A side)
252   } else {
253     sign = -1;       // (TPC C side)
254   }
255   
256   if ( sign==1  && z <  fgkZOffSet ) z =  fgkZOffSet;    // Protect against discontinuity at CE
257   if ( sign==-1 && z > -fgkZOffSet ) z = -fgkZOffSet;    // Protect against discontinuity at CE
258   
259
260   if ( (sign==1 && z<0) || (sign==-1 && z>0) ) // just a consistency check
261     AliError("ROC number does not correspond to z coordinate! Calculation of distortions is most likely wrong!");
262
263   // Get the Er and Ephi field integrals plus the integral over DeltaEz 
264   intEr      = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi, 
265                                   fgkRList, fgkZList, fgkPhiList, fLookUpErOverEz  );
266   intEphi    = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi, 
267                                   fgkRList, fgkZList, fgkPhiList, fLookUpEphiOverEz  );
268   intDeltaEz = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi, 
269                                   fgkRList, fgkZList, fgkPhiList, fLookUpDeltaEz  );
270
271   //  printf("%lf %lf %lf\n",intEr,intEphi,intDeltaEz);
272
273   // Calculate distorted position
274   if ( r > 0.0 ) {
275     phi =  phi + ( fC0*intEphi - fC1*intEr ) / r;      
276     r   =  r   + ( fC0*intEr   + fC1*intEphi );  
277   }
278   
279   // Calculate correction in cartesian coordinates
280   dx[0] = r * TMath::Cos(phi) - x[0];
281   dx[1] = r * TMath::Sin(phi) - x[1]; 
282   dx[2] = intDeltaEz;  // z distortion - (internally scaled with driftvelocity dependency 
283                        // on the Ez field plus the actual ROC misalignment (if set TRUE)
284
285 }
286
287 void AliTPCFCVoltError3D::InitFCVoltError3D() {
288   //
289   // Initialization of the Lookup table which contains the solutions of the 
290   // Dirichlet boundary problem
291   // Calculation of the single 3D-Poisson solver is done just if needed
292   // (see basic lookup tables in header file)
293   //
294
295   const Int_t   order       =    1  ;  // Linear interpolation = 1, Quadratic = 2  
296   const Float_t gridSizeR   =  (fgkOFCRadius-fgkIFCRadius) / (kRows-1) ;
297   const Float_t gridSizeZ   =  fgkTPCZ0 / (kColumns-1) ;
298   const Float_t gridSizePhi =  TMath::TwoPi() / ( 18.0 * kPhiSlicesPerSector);
299
300   // temporary arrays to create the boundary conditions
301   TMatrixD *arrayofArrayV[kPhiSlices], *arrayofCharge[kPhiSlices] ; 
302   for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
303     arrayofArrayV[k]     =   new TMatrixD(kRows,kColumns) ;
304     arrayofCharge[k]     =   new TMatrixD(kRows,kColumns) ;
305   }
306   Double_t  innerList[kPhiSlices], outerList[kPhiSlices] ;
307   
308   // list of point as used in the poisson relation and the interpolation (during sum up)
309   Double_t  rlist[kRows], zedlist[kColumns] , philist[kPhiSlices];
310   for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
311     philist[k] =  gridSizePhi * k;
312     for ( Int_t i = 0 ; i < kRows ; i++ )    {
313       rlist[i] = fgkIFCRadius + i*gridSizeR ;
314       for ( Int_t j = 0 ; j < kColumns ; j++ ) { // Fill Vmatrix with Boundary Conditions
315         zedlist[j]  = j * gridSizeZ ;
316       }
317     }
318   }
319
320   // ==========================================================================
321   // Solve Poisson's equation in 3D cylindrical coordinates by relaxation technique
322   // Allow for different size grid spacing in R and Z directions
323   
324   const Int_t   symmetry[6] =    {1,1,-1,-1,1,1}; // shifted rod (2x), rotated clip (2x), only rod shift on OFC (1x)
325
326   // check which lookup tables are needed ---------------------------------
327
328   Bool_t needTable[6] ={kFALSE,kFALSE,kFALSE,kFALSE,kFALSE,kFALSE};
329
330   // Shifted rods & strips
331   for ( Int_t rod = 0 ; rod < 18 ; rod++ ) {
332     if (fRodVoltShiftA[rod]!=0) needTable[0]=kTRUE;
333     if (fRodVoltShiftC[rod]!=0) needTable[0]=kTRUE;
334   }
335   for ( Int_t rod = 18 ; rod < 36 ; rod++ ) {
336     if (fRodVoltShiftA[rod]!=0) needTable[1]=kTRUE;
337     if (fRodVoltShiftC[rod]!=0) needTable[1]=kTRUE;
338   }
339   // Rotated clips
340   if (fRotatedClipVoltA[0]!=0 || fRotatedClipVoltC[0]!=0) needTable[2]=kTRUE;
341   if (fRotatedClipVoltA[1]!=0 || fRotatedClipVoltC[1]!=0) needTable[3]=kTRUE;
342  
343   // shifted Copper rods 
344   for ( Int_t rod = 0 ; rod < 18 ; rod++ ) {
345     if (fCopperRodShiftA[rod]!=0) needTable[4]=kTRUE;
346     if (fCopperRodShiftC[rod]!=0) needTable[4]=kTRUE;
347   }
348   // shifted Copper rods 
349   for ( Int_t rod = 18; rod < 36 ; rod++ ) {
350     if (fCopperRodShiftA[rod]!=0) needTable[5]=kTRUE;
351     if (fCopperRodShiftC[rod]!=0) needTable[5]=kTRUE;
352   }
353
354
355   // reserve the arrays for the basic lookup tables ----------------------
356   if (needTable[0] && !fInitLookUpBasic[0]) {
357     for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {   // Possibly make an extra table to be used for 0 == 360
358       fLookUpBasic1ErOverEz[k]   =   new TMatrixD(kRows,kColumns);
359       fLookUpBasic1EphiOverEz[k] =   new TMatrixD(kRows,kColumns);
360       fLookUpBasic1DeltaEz[k]    =   new TMatrixD(kRows,kColumns);
361       // will be deleted in destructor
362     }
363   }
364   if (needTable[1] && !fInitLookUpBasic[1]) {
365     for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {   // Possibly make an extra table to be used for 0 == 360
366       fLookUpBasic2ErOverEz[k]   =   new TMatrixD(kRows,kColumns);
367       fLookUpBasic2EphiOverEz[k] =   new TMatrixD(kRows,kColumns);
368       fLookUpBasic2DeltaEz[k]    =   new TMatrixD(kRows,kColumns);
369       // will be deleted in destructor
370     }
371   }
372   if (needTable[2] && !fInitLookUpBasic[2]) {
373     for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {   // Possibly make an extra table to be used for 0 == 360
374       fLookUpBasic3ErOverEz[k]   =   new TMatrixD(kRows,kColumns);
375       fLookUpBasic3EphiOverEz[k] =   new TMatrixD(kRows,kColumns);
376       fLookUpBasic3DeltaEz[k]    =   new TMatrixD(kRows,kColumns);
377       // will be deleted in destructor
378     }
379   }
380   if (needTable[3] && !fInitLookUpBasic[3]) {
381     for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {   // Possibly make an extra table to be used for 0 == 360
382       fLookUpBasic4ErOverEz[k]   =   new TMatrixD(kRows,kColumns);
383       fLookUpBasic4EphiOverEz[k] =   new TMatrixD(kRows,kColumns);
384       fLookUpBasic4DeltaEz[k]    =   new TMatrixD(kRows,kColumns);
385       // will be deleted in destructor
386     }
387   }
388   if (needTable[4] && !fInitLookUpBasic[4]) {
389     for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {   // Possibly make an extra table to be used for 0 == 360
390       fLookUpBasic5ErOverEz[k]   =   new TMatrixD(kRows,kColumns);
391       fLookUpBasic5EphiOverEz[k] =   new TMatrixD(kRows,kColumns);
392       fLookUpBasic5DeltaEz[k]    =   new TMatrixD(kRows,kColumns);
393       // will be deleted in destructor
394     }
395   }
396   if (needTable[5] && !fInitLookUpBasic[5]) {
397     for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {   // Possibly make an extra table to be used for 0 == 360
398       fLookUpBasic6ErOverEz[k]   =   new TMatrixD(kRows,kColumns);
399       fLookUpBasic6EphiOverEz[k] =   new TMatrixD(kRows,kColumns);
400       fLookUpBasic6DeltaEz[k]    =   new TMatrixD(kRows,kColumns);
401       // will be deleted in destructor
402     }
403   }
404  
405   // Set bondaries and solve Poisson's equation --------------------------
406  
407   for (Int_t look=0; look<6; look++) {
408    
409     Float_t inner18[18] = { 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 } ;  
410     Float_t outer18[18] = { 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 } ; 
411   
412     if (needTable[look] && !fInitLookUpBasic[look]) {
413
414       // Specify which rod is the reference; other distortions calculated by summing over multiple rotations of refrence
415       // Note: the interpolation later on depends on it!! Do not change if not really needed!
416       if (look==0) {
417         AliInfo(Form("IFC - ROD&Strip shift : Filling table (~ %d sec)",3*(int)(kPhiSlices/5)));
418         inner18[0] = 1;  
419       } else if (look==1) {
420         AliInfo(Form("OFC - ROD&Strip shift : Filling table (~ %d sec)",3*(int)(kPhiSlices/5)));
421         outer18[0] = 1;  
422       } else if (look==2) {
423         AliInfo(Form("IFC - Clip rot. : Filling table (~ %d sec)",3*(int)(kPhiSlices/5)));
424         inner18[0] = 1; 
425       } else if (look==3) {
426         AliInfo(Form("OFC - Clip rot. : Filling table (~ %d sec)",3*(int)(kPhiSlices/5)));
427         outer18[0] = 1;  
428       } else if (look==4) {
429         AliInfo(Form("IFC - CopperRod shift : Filling table (~ %d sec)",3*(int)(kPhiSlices/5)));
430         inner18[0] = 1;  
431       } else if (look==5) {
432         AliInfo(Form("OFC - CopperRod shift : Filling table (~ %d sec)",3*(int)(kPhiSlices/5)));
433         outer18[0] = 1;  
434       }
435       // Build potentials on boundary stripes for a rotated clip (SYMMETRY==-1) or a shifted rod (SYMMETRY==1)
436       if (look<4) {
437         // in these cases, the strips follow the actual rod shift (linear interpolation between the rods)
438         for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
439           Int_t slices = kPhiSlicesPerSector ;
440           Int_t ipoint = k/slices  ;  
441           innerList[k] = (((ipoint+1)*slices-k)*inner18[ipoint]-(k-ipoint*slices)*inner18[ipoint+1])/slices ;
442           outerList[k] = (((ipoint+1)*slices-k)*outer18[ipoint]-(k-ipoint*slices)*outer18[ipoint+1])/slices ;
443           if ( (k%slices) == 0 && symmetry[look] == -1 ) { innerList[k] = 0.0 ; outerList[k] = 0.0 ; } 
444           // above, force Zero if Anti-Sym
445         } 
446       } else if (look==4){
447         // in this case, we assume that the strips stay at the correct position, but the rods move
448         // the distortion is then much more localized around the rod (indicated by real data)
449         for ( Int_t k = 0 ; k < kPhiSlices ; k++ )
450           innerList[k] = outerList[k] = 0;
451         innerList[0] = inner18[0];     // point at rod
452         innerList[0] = inner18[0]/4*3;  // point close to rod (educated guess)
453       } else if (look==5){
454         // in this case, we assume that the strips stay at the correct position, but the copper plated OFC-rods move
455         // the distortion is then much more localized around the rod (indicated by real data)
456         for ( Int_t k = 0 ; k < kPhiSlices ; k++ )
457           innerList[k] = outerList[k] = 0;
458         outerList[0] = outer18[0];   // point at rod
459         outerList[1] = outer18[0]/4; // point close to rod (educated-`guessed` scaling)
460       }
461
462       // Fill arrays with initial conditions.  V on the boundary and Charge in the volume.
463       for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
464         TMatrixD &arrayV    =  *arrayofArrayV[k] ;
465         TMatrixD &charge    =  *arrayofCharge[k] ;
466         for ( Int_t i = 0 ; i < kRows ; i++ )    {
467           for ( Int_t j = 0 ; j < kColumns ; j++ ) { // Fill Vmatrix with Boundary Conditions
468             arrayV(i,j) = 0.0 ; 
469             charge(i,j) = 0.0 ;
470             if ( i == 0 )       arrayV(i,j) = innerList[k] ; 
471             if ( i == kRows-1 ) arrayV(i,j) = outerList[k] ; 
472           }
473         }      
474         // no charge in the volume
475         for ( Int_t i = 1 ; i < kRows-1 ; i++ )  { 
476           for ( Int_t j = 1 ; j < kColumns-1 ; j++ ) {
477             charge(i,j)  =  0.0 ;
478           }
479         }
480       }      
481      
482       // Solve Poisson's equation in 3D cylindrical coordinates by relaxation technique
483       // Allow for different size grid spacing in R and Z directions
484       if (look==0) {
485         PoissonRelaxation3D( arrayofArrayV, arrayofCharge, 
486                              fLookUpBasic1ErOverEz, fLookUpBasic1EphiOverEz, fLookUpBasic1DeltaEz,
487                              kRows, kColumns, kPhiSlices, gridSizePhi, kIterations, symmetry[0]) ;
488         AliInfo("IFC - ROD&Strip shift : done ");
489       } else if (look==1) {
490         PoissonRelaxation3D( arrayofArrayV, arrayofCharge, 
491                              fLookUpBasic2ErOverEz, fLookUpBasic2EphiOverEz, fLookUpBasic2DeltaEz,
492                              kRows, kColumns, kPhiSlices, gridSizePhi, kIterations, symmetry[1]) ;
493         
494         AliInfo("OFC - ROD&Strip shift : done ");
495       } else if (look==2) {
496         PoissonRelaxation3D( arrayofArrayV, arrayofCharge, 
497                              fLookUpBasic3ErOverEz, fLookUpBasic3EphiOverEz, fLookUpBasic3DeltaEz,
498                              kRows, kColumns, kPhiSlices, gridSizePhi, kIterations, symmetry[2]) ;
499         AliInfo("IFC - Clip rot. : done ");
500       } else if (look==3) {
501         PoissonRelaxation3D( arrayofArrayV, arrayofCharge, 
502                              fLookUpBasic4ErOverEz, fLookUpBasic4EphiOverEz, fLookUpBasic4DeltaEz,
503                              kRows, kColumns, kPhiSlices, gridSizePhi, kIterations, symmetry[3]) ;
504         AliInfo("OFC - Clip rot. : done ");
505       } else if (look==4) {
506         PoissonRelaxation3D( arrayofArrayV, arrayofCharge, 
507                              fLookUpBasic5ErOverEz, fLookUpBasic5EphiOverEz, fLookUpBasic5DeltaEz,
508                              kRows, kColumns, kPhiSlices, gridSizePhi, kIterations, symmetry[4]) ;
509         AliInfo("IFC - CopperRod shift : done ");
510       } else if (look==5) {
511         PoissonRelaxation3D( arrayofArrayV, arrayofCharge, 
512                              fLookUpBasic6ErOverEz, fLookUpBasic6EphiOverEz, fLookUpBasic6DeltaEz,
513                              kRows, kColumns, kPhiSlices, gridSizePhi, kIterations, symmetry[5]) ;
514         AliInfo("OFC - CopperRod shift : done ");
515       }
516       
517       fInitLookUpBasic[look] = kTRUE;
518     }
519   }
520   
521
522   // =============================================================================
523   // Create the final lookup table with corresponding ROD shifts or clip rotations
524
525   Float_t erIntegralSum   = 0.0 ;
526   Float_t ephiIntegralSum = 0.0 ;
527   Float_t ezIntegralSum   = 0.0 ;
528
529   Double_t phiPrime     = 0. ;
530   Double_t erIntegral   = 0. ;
531   Double_t ephiIntegral = 0. ;
532   Double_t ezIntegral   = 0. ;
533
534
535   AliInfo("Calculation of combined Look-up Table");
536
537   // Interpolate and sum the Basic lookup tables onto the standard grid
538   Double_t  r, phi, z ;
539   for ( Int_t k = 0 ; k < kNPhi ; k++ ) {
540     phi = fgkPhiList[k] ;
541
542     TMatrixD &erOverEz   =  *fLookUpErOverEz[k]  ;
543     TMatrixD &ephiOverEz =  *fLookUpEphiOverEz[k];
544     TMatrixD &deltaEz    =  *fLookUpDeltaEz[k]   ;
545
546     for ( Int_t i = 0 ; i < kNZ ; i++ ) {
547       z = TMath::Abs(fgkZList[i]) ;  // Symmetric solution in Z that depends only on ABS(Z)
548       for ( Int_t j = 0 ; j < kNR ; j++ ) { 
549         r = fgkRList[j] ;
550         // Interpolate basicLookup tables; once for each rod, then sum the results
551         
552         erIntegralSum   = 0.0 ;
553         ephiIntegralSum = 0.0 ;
554         ezIntegralSum   = 0.0 ;
555  
556         // SHIFTED RODS =========================================================
557
558         // IFC ROD SHIFTS +++++++++++++++++++++++++++++
559         for ( Int_t rod = 0 ; rod < 18 ; rod++ ) {
560           
561           erIntegral = ephiIntegral = ezIntegral = 0.0 ;
562           
563           if ( fRodVoltShiftA[rod] == 0 && fgkZList[i] > 0) continue ;
564           if ( fRodVoltShiftC[rod] == 0 && fgkZList[i] < 0) continue ;
565
566           // Rotate to a position where we have maps and prepare to sum
567           phiPrime =  phi - rod*kPhiSlicesPerSector*gridSizePhi ;  
568
569           if ( phiPrime < 0 ) phiPrime += TMath::TwoPi() ;   // Stay in range from 0 to TwoPi    
570
571           if ( (phiPrime >= 0) && (phiPrime <= kPhiSlices*gridSizePhi) ) {
572             
573             erIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices, 
574                                               rlist, zedlist, philist, fLookUpBasic1ErOverEz  );
575             ephiIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
576                                               rlist, zedlist, philist, fLookUpBasic1EphiOverEz);
577             ezIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
578                                               rlist, zedlist, philist, fLookUpBasic1DeltaEz   );
579             
580           }  else if ( (phiPrime <= TMath::TwoPi()) && (phiPrime >= (TMath::TwoPi()-kPhiSlices*gridSizePhi)) ){
581             
582             phiPrime     = TMath::TwoPi() - phiPrime ;
583             
584             erIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices, 
585                                               rlist, zedlist, philist, fLookUpBasic1ErOverEz  );
586             ephiIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
587                                               rlist, zedlist, philist, fLookUpBasic1EphiOverEz);
588             ezIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
589                                               rlist, zedlist, philist, fLookUpBasic1DeltaEz   );
590             
591             // Flip symmetry of phi integral for symmetric boundary conditions
592             if ( symmetry[0] ==  1 ) ephiIntegral *= -1  ;    
593             // Flip symmetry of r integral if anti-symmetric boundary conditions 
594             if ( symmetry[0] == -1 ) erIntegral   *= -1  ;    
595
596           }
597
598           if ( fgkZList[i] > 0 ) {
599             erIntegralSum   += fRodVoltShiftA[rod]*erIntegral   ;
600             ephiIntegralSum += fRodVoltShiftA[rod]*ephiIntegral ;
601             ezIntegralSum   += fRodVoltShiftA[rod]*ezIntegral;
602           } else if ( fgkZList[i] < 0 ) {
603             erIntegralSum   += fRodVoltShiftC[rod]*erIntegral   ;
604             ephiIntegralSum += fRodVoltShiftC[rod]*ephiIntegral ;
605             ezIntegralSum   -= fRodVoltShiftC[rod]*ezIntegral;
606           }
607         }
608
609         // OFC ROD SHIFTS +++++++++++++++++++++++++++++
610         for ( Int_t rod = 18 ; rod < 36 ; rod++ ) {
611
612           erIntegral = ephiIntegral = ezIntegral = 0.0 ;
613           
614           if ( fRodVoltShiftA[rod] == 0 && fgkZList[i] > 0) continue ;
615           if ( fRodVoltShiftC[rod] == 0 && fgkZList[i] < 0) continue ;
616
617           // Rotate to a position where we have maps and prepare to sum
618           phiPrime =  phi - (rod-18)*kPhiSlicesPerSector*gridSizePhi ;  
619                   
620           if ( phiPrime < 0 ) phiPrime += TMath::TwoPi() ;   // Stay in range from 0 to TwoPi    
621
622           if ( (phiPrime >= 0) && (phiPrime <= kPhiSlices*gridSizePhi) ) {
623             
624             erIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices, 
625                                               rlist, zedlist, philist, fLookUpBasic2ErOverEz  );
626             ephiIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
627                                               rlist, zedlist, philist, fLookUpBasic2EphiOverEz);
628             ezIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
629                                               rlist, zedlist, philist, fLookUpBasic2DeltaEz   );
630             
631           }  else if ( (phiPrime <= TMath::TwoPi()) && (phiPrime >= (TMath::TwoPi()-kPhiSlices*gridSizePhi)) ){
632             
633             phiPrime     = TMath::TwoPi() - phiPrime ;
634             
635             erIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices, 
636                                               rlist, zedlist, philist, fLookUpBasic2ErOverEz  );
637             ephiIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
638                                               rlist, zedlist, philist, fLookUpBasic2EphiOverEz);
639             ezIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
640                                               rlist, zedlist, philist, fLookUpBasic2DeltaEz   );
641             
642             // Flip symmetry of phi integral for symmetric boundary conditions
643             if ( symmetry[1] ==  1 ) ephiIntegral *= -1  ;    
644             // Flip symmetry of r integral if anti-symmetric boundary conditions 
645             if ( symmetry[1] == -1 ) erIntegral   *= -1  ;    
646
647           }
648
649           if ( fgkZList[i] > 0 ) {
650             erIntegralSum   += fRodVoltShiftA[rod]*erIntegral   ;
651             ephiIntegralSum += fRodVoltShiftA[rod]*ephiIntegral ;
652             ezIntegralSum   += fRodVoltShiftA[rod]*ezIntegral;
653           } else if ( fgkZList[i] < 0 ) {
654             erIntegralSum   += fRodVoltShiftC[rod]*erIntegral   ;
655             ephiIntegralSum += fRodVoltShiftC[rod]*ephiIntegral ;
656             ezIntegralSum   -= fRodVoltShiftC[rod]*ezIntegral;
657           }
658
659         } // rod loop - shited ROD
660
661
662         // Rotated clips =========================================================
663
664         Int_t rodIFC = 11; // resistor rod positions, IFC
665         Int_t rodOFC =  3; // resistor rod positions, OFC
666         // just one rod on IFC and OFC
667
668         // IFC ROTATED CLIP +++++++++++++++++++++++++++++
669         for ( Int_t rod = rodIFC ; rod < rodIFC+1 ; rod++ ) { // loop over 1 to keep the "ignore"-functionality
670
671           erIntegral = ephiIntegral = ezIntegral = 0.0 ;
672           if ( fRotatedClipVoltA[0] == 0 && fgkZList[i] > 0) continue ;
673           if ( fRotatedClipVoltC[0] == 0 && fgkZList[i] < 0) continue ;
674
675           // Rotate to a position where we have maps and prepare to sum
676           phiPrime =  phi - rod*kPhiSlicesPerSector*gridSizePhi ;  
677           
678           if ( phiPrime < 0 ) phiPrime += TMath::TwoPi() ;   // Stay in range from 0 to TwoPi    
679           
680           if ( (phiPrime >= 0) && (phiPrime <= kPhiSlices*gridSizePhi) ) {
681             
682             erIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices, 
683                                               rlist, zedlist, philist, fLookUpBasic3ErOverEz  );
684             ephiIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
685                                               rlist, zedlist, philist, fLookUpBasic3EphiOverEz);
686             ezIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
687                                               rlist, zedlist, philist, fLookUpBasic3DeltaEz   );
688             
689           }  else if ( (phiPrime <= TMath::TwoPi()) && (phiPrime >= (TMath::TwoPi()-kPhiSlices*gridSizePhi)) ){
690             
691             phiPrime     = TMath::TwoPi() - phiPrime ;
692             
693             erIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices, 
694                                               rlist, zedlist, philist, fLookUpBasic3ErOverEz  );
695             ephiIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
696                                               rlist, zedlist, philist, fLookUpBasic3EphiOverEz);
697             ezIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
698                                               rlist, zedlist, philist, fLookUpBasic3DeltaEz   );
699             
700             // Flip symmetry of phi integral for symmetric boundary conditions
701             if ( symmetry[2] ==  1 ) ephiIntegral *= -1  ;    
702             // Flip symmetry of r integral if anti-symmetric boundary conditions 
703             if ( symmetry[2] == -1 ) erIntegral   *= -1  ;    
704             
705           }
706           
707           if ( fgkZList[i] > 0 ) {
708             erIntegralSum   += fRotatedClipVoltA[0]*erIntegral   ;
709             ephiIntegralSum += fRotatedClipVoltA[0]*ephiIntegral ;
710             ezIntegralSum   += fRotatedClipVoltA[0]*ezIntegral;
711           } else if ( fgkZList[i] < 0 ) {
712             erIntegralSum   += fRotatedClipVoltC[0]*erIntegral   ;
713             ephiIntegralSum += fRotatedClipVoltC[0]*ephiIntegral ;
714             ezIntegralSum   -= fRotatedClipVoltC[0]*ezIntegral;
715           }
716         }
717
718         // OFC: ROTATED CLIP  +++++++++++++++++++++++++++++
719         for ( Int_t rod = rodOFC ; rod < rodOFC+1 ; rod++ ) { // loop over 1 to keep the "ignore"-functionality
720           
721           erIntegral = ephiIntegral = ezIntegral = 0.0 ;
722           
723           if ( fRotatedClipVoltA[1] == 0 && fgkZList[i] > 0) continue ;
724           if ( fRotatedClipVoltC[1] == 0 && fgkZList[i] < 0) continue ;
725
726           // Rotate to a position where we have maps and prepare to sum
727           phiPrime =  phi - rod*kPhiSlicesPerSector*gridSizePhi ;  
728           
729           
730           if ( phiPrime < 0 ) phiPrime += TMath::TwoPi() ;   // Stay in range from 0 to TwoPi    
731           
732           if ( (phiPrime >= 0) && (phiPrime <= kPhiSlices*gridSizePhi) ) {
733             
734             erIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices, 
735                                               rlist, zedlist, philist, fLookUpBasic4ErOverEz  );
736             ephiIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
737                                               rlist, zedlist, philist, fLookUpBasic4EphiOverEz);
738             ezIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
739                                               rlist, zedlist, philist, fLookUpBasic4DeltaEz   );
740             
741           }  else if ( (phiPrime <= TMath::TwoPi()) && (phiPrime >= (TMath::TwoPi()-kPhiSlices*gridSizePhi)) ){
742             
743             phiPrime     = TMath::TwoPi() - phiPrime ;
744             
745             erIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices, 
746                                               rlist, zedlist, philist, fLookUpBasic4ErOverEz  );
747             ephiIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
748                                               rlist, zedlist, philist, fLookUpBasic4EphiOverEz);
749             ezIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
750                                               rlist, zedlist, philist, fLookUpBasic4DeltaEz   );
751             
752             // Flip symmetry of phi integral for symmetric boundary conditions
753             if ( symmetry[3] ==  1 ) ephiIntegral *= -1  ;    
754             // Flip symmetry of r integral if anti-symmetric boundary conditions 
755             if ( symmetry[3] == -1 ) erIntegral   *= -1  ;    
756             
757           }
758           
759           if ( fgkZList[i] > 0 ) {
760             erIntegralSum   += fRotatedClipVoltA[1]*erIntegral   ;
761             ephiIntegralSum += fRotatedClipVoltA[1]*ephiIntegral ;
762             ezIntegralSum   += fRotatedClipVoltA[1]*ezIntegral;
763           } else if ( fgkZList[i] < 0 ) {
764             erIntegralSum   += fRotatedClipVoltC[1]*erIntegral   ;
765             ephiIntegralSum += fRotatedClipVoltC[1]*ephiIntegral ;
766             ezIntegralSum   -= fRotatedClipVoltC[1]*ezIntegral;
767           }
768         }
769
770         // IFC Cooper rod shift  +++++++++++++++++++++++++++++
771         for ( Int_t rod = 0 ; rod < 18 ; rod++ ) {
772           
773           erIntegral = ephiIntegral = ezIntegral = 0.0 ;
774           
775           if ( fCopperRodShiftA[rod] == 0 && fgkZList[i] > 0) continue ;
776           if ( fCopperRodShiftC[rod] == 0 && fgkZList[i] < 0) continue ;
777
778           // Rotate to a position where we have maps and prepare to sum
779           phiPrime =  phi - rod*kPhiSlicesPerSector*gridSizePhi ;  
780
781           if ( phiPrime < 0 ) phiPrime += TMath::TwoPi() ;   // Stay in range from 0 to TwoPi    
782
783           if ( (phiPrime >= 0) && (phiPrime <= kPhiSlices*gridSizePhi) ) {
784             
785             erIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices, 
786                                               rlist, zedlist, philist, fLookUpBasic5ErOverEz  );
787             ephiIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
788                                               rlist, zedlist, philist, fLookUpBasic5EphiOverEz);
789             ezIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
790                                               rlist, zedlist, philist, fLookUpBasic5DeltaEz   );
791             
792           }  else if ( (phiPrime <= TMath::TwoPi()) && (phiPrime >= (TMath::TwoPi()-kPhiSlices*gridSizePhi)) ){
793             
794             phiPrime     = TMath::TwoPi() - phiPrime ;
795             
796             erIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices, 
797                                               rlist, zedlist, philist, fLookUpBasic5ErOverEz  );
798             ephiIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
799                                               rlist, zedlist, philist, fLookUpBasic5EphiOverEz);
800             ezIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
801                                               rlist, zedlist, philist, fLookUpBasic5DeltaEz   );
802             
803             // Flip symmetry of phi integral for symmetric boundary conditions
804             if ( symmetry[4] ==  1 ) ephiIntegral *= -1  ;    
805             // Flip symmetry of r integral if anti-symmetric boundary conditions 
806             if ( symmetry[4] == -1 ) erIntegral   *= -1  ;    
807
808           }
809
810           if ( fgkZList[i] > 0 ) {
811             erIntegralSum   += fCopperRodShiftA[rod]*erIntegral   ;
812             ephiIntegralSum += fCopperRodShiftA[rod]*ephiIntegral ;
813             ezIntegralSum   += fCopperRodShiftA[rod]*ezIntegral;
814           } else if ( fgkZList[i] < 0 ) {
815             erIntegralSum   += fCopperRodShiftC[rod]*erIntegral   ;
816             ephiIntegralSum += fCopperRodShiftC[rod]*ephiIntegral ;
817             ezIntegralSum   -= fCopperRodShiftC[rod]*ezIntegral;
818           }
819         }
820
821         // OFC Cooper rod shift  +++++++++++++++++++++++++++++
822         for ( Int_t rod = 18 ; rod < 36 ; rod++ ) {
823           
824           erIntegral = ephiIntegral = ezIntegral = 0.0 ;
825           
826           if ( fCopperRodShiftA[rod] == 0 && fgkZList[i] > 0) continue ;
827           if ( fCopperRodShiftC[rod] == 0 && fgkZList[i] < 0) continue ;
828
829           // Rotate to a position where we have maps and prepare to sum
830           phiPrime =  phi - (rod-18)*kPhiSlicesPerSector*gridSizePhi ;  
831
832           if ( phiPrime < 0 ) phiPrime += TMath::TwoPi() ;   // Stay in range from 0 to TwoPi    
833
834           if ( (phiPrime >= 0) && (phiPrime <= kPhiSlices*gridSizePhi) ) {
835             
836             erIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices, 
837                                               rlist, zedlist, philist, fLookUpBasic6ErOverEz  );
838             ephiIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
839                                               rlist, zedlist, philist, fLookUpBasic6EphiOverEz);
840             ezIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
841                                               rlist, zedlist, philist, fLookUpBasic6DeltaEz   );
842             
843           }  else if ( (phiPrime <= TMath::TwoPi()) && (phiPrime >= (TMath::TwoPi()-kPhiSlices*gridSizePhi)) ){
844             
845             phiPrime     = TMath::TwoPi() - phiPrime ;
846             
847             erIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices, 
848                                               rlist, zedlist, philist, fLookUpBasic6ErOverEz  );
849             ephiIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
850                                               rlist, zedlist, philist, fLookUpBasic6EphiOverEz);
851             ezIntegral   = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
852                                               rlist, zedlist, philist, fLookUpBasic6DeltaEz   );
853             
854             // Flip symmetry of phi integral for symmetric boundary conditions
855             if ( symmetry[5] ==  1 ) ephiIntegral *= -1  ;    
856             // Flip symmetry of r integral if anti-symmetric boundary conditions 
857             if ( symmetry[5] == -1 ) erIntegral   *= -1  ;    
858
859           }
860
861           if ( fgkZList[i] > 0 ) {
862             erIntegralSum   += fCopperRodShiftA[rod]*erIntegral   ;
863             ephiIntegralSum += fCopperRodShiftA[rod]*ephiIntegral ;
864             ezIntegralSum   += fCopperRodShiftA[rod]*ezIntegral;
865           } else if ( fgkZList[i] < 0 ) {
866             erIntegralSum   += fCopperRodShiftC[rod]*erIntegral   ;
867             ephiIntegralSum += fCopperRodShiftC[rod]*ephiIntegral ;
868             ezIntegralSum   -= fCopperRodShiftC[rod]*ezIntegral;
869           }
870         }
871
872         // put the sum into the final lookup table
873         erOverEz(j,i)   =  erIntegralSum;
874         ephiOverEz(j,i) =  ephiIntegralSum;
875         deltaEz(j,i)    =  ezIntegralSum;
876         
877         //      if (j==1) printf("%lf %lf %lf: %lf %lf %lf\n",r, phi, z, erIntegralSum,ephiIntegralSum,ezIntegralSum);
878  
879       } // endo r loop
880     } // end of z loop
881   } // end of phi loop
882
883
884   // clear the temporary arrays lists
885   for ( Int_t k = 0 ; k < kPhiSlices ; k++ )  {
886     delete arrayofArrayV[k];
887     delete arrayofCharge[k];
888   }
889  
890   AliInfo(" done");
891   fInitLookUp = kTRUE;
892
893 }
894
895 void AliTPCFCVoltError3D::Print(const Option_t* option) const {
896   //
897   // Print function to check the settings of the Rod shifts and the rotated clips
898   // option=="a" prints the C0 and C1 coefficents for calibration purposes
899   //
900
901   TString opt = option; opt.ToLower();
902   printf("%s\n",GetTitle());
903   printf(" - ROD shifts  (residual voltage settings): 40V correspond to 1mm of shift\n");
904   Int_t count = 0;
905   printf("  : A-side - (Rod & Strips) \n     "); 
906   for (Int_t i=0; i<36;i++) {
907     if (fRodVoltShiftA[i]!=0) {
908       printf("Rod%2d:%3.1fV ",i,fRodVoltShiftA[i]);
909       count++;
910     }
911     if (count==0&&i==35) 
912       printf("-> all at 0 V\n");
913     else if (i==35){
914       printf("\n");
915       count=0;
916     }
917   } 
918   printf("  : C-side - (Rod & Strips) \n     "); 
919   for (Int_t i=0; i<36;i++) {
920     if (fRodVoltShiftC[i]!=0) {
921       printf("Rod%2d:%3.1fV ",i,fRodVoltShiftC[i]);
922       count++;
923     }
924     if (count==0&&i==35) 
925       printf("-> all at 0 V\n");
926     else if (i==35){
927       printf("\n");
928       count=0;
929     }
930   } 
931  
932   printf(" - Rotated clips (residual voltage settings): 40V correspond to 1mm of 'offset'\n");
933   if (fRotatedClipVoltA[0]!=0) {   printf("     A side (IFC): HV Rod: %3.1f V \n",fRotatedClipVoltA[0]); count++; }
934   if (fRotatedClipVoltA[1]!=0) {   printf("     A side (OFC): HV Rod: %3.1f V \n",fRotatedClipVoltA[1]); count++; }
935   if (fRotatedClipVoltC[0]!=0) {   printf("     C side (IFC): HV Rod: %3.1f V \n",fRotatedClipVoltC[0]); count++; }
936   if (fRotatedClipVoltC[1]!=0) {   printf("     C side (OFC): HV Rod: %3.1f V \n",fRotatedClipVoltC[1]); count++; }
937   if (count==0) 
938     printf("     -> no rotated clips \n");
939
940   count=0;
941   printf(" - Copper ROD shifts (without strips):\n");
942   printf("  : A-side - (Copper Rod shift) \n     "); 
943   for (Int_t i=0; i<36;i++) {
944     if (fCopperRodShiftA[i]!=0) {
945       printf("Rod%2d:%3.1fV ",i,fCopperRodShiftA[i]);
946       count++;
947     }
948     if (count==0&&i==35) 
949       printf("-> all at 0 V\n");
950     else if (i==35){
951       printf("\n");
952       count=0;
953     }
954   } 
955   printf("  : C-side - (Copper Rod shift) \n     "); 
956   for (Int_t i=0; i<36;i++) {
957     if (fCopperRodShiftC[i]!=0) {
958       printf("Rod%2d:%3.1fV ",i,fCopperRodShiftC[i]);
959       count++;
960     }
961     if (count==0&&i==35) 
962       printf("-> all at 0 V\n");
963     else if (i==35){
964       printf("\n");
965       count=0;
966     }
967   } 
968
969   if (opt.Contains("a")) { // Print all details
970     printf(" - T1: %1.4f, T2: %1.4f \n",fT1,fT2);
971     printf(" - C1: %1.4f, C0: %1.4f \n",fC1,fC0);
972   }
973
974   if (!fInitLookUp) AliError("Lookup table was not initialized! You should do InitFCVoltError3D() ...");
975
976 }