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