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