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