AliRoot Core  edcc906 (edcc906)
AliTPCFCVoltError3D.cxx
Go to the documentation of this file.
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 
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 
65 ClassImp(AliTPCFCVoltError3D)
67 
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'
78  for (Int_t i=0; i<6; i++){
79  fInitLookUpBasic[i]= kFALSE;
80  }
81 
82  // Boundary settings
83  for (Int_t i=0; i<36; i++){
84  fRodVoltShiftA[i] = 0;
85  fRodVoltShiftC[i] = 0;
86  }
87  for (Int_t i=0; i<2; i++){
88  fRotatedClipVoltA[i] = 0;
89  fRotatedClipVoltC[i] = 0;
90  }
91  //
92  for (Int_t i=0; i<36; i++){
93  fCopperRodShiftA[i] = 0;
94  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  for ( Int_t k = 0 ; k < kNPhi ; k++ ) {
101  fLookUpErOverEz[k] = new TMatrixF(kNR,kNZ);
102  fLookUpEphiOverEz[k] = new TMatrixF(kNR,kNZ);
103  fLookUpDeltaEz[k] = new TMatrixF(kNR,kNZ);
104  }
105 
106  for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
107  fLookUpBasic1ErOverEz[k] = 0;
108  fLookUpBasic1EphiOverEz[k] = 0;
109  fLookUpBasic1DeltaEz[k] = 0;
110 
111  fLookUpBasic2ErOverEz[k] = 0;
112  fLookUpBasic2EphiOverEz[k] = 0;
113  fLookUpBasic2DeltaEz[k] = 0;
114 
115  fLookUpBasic3ErOverEz[k] = 0;
116  fLookUpBasic3EphiOverEz[k] = 0;
117  fLookUpBasic3DeltaEz[k] = 0;
118 
119  fLookUpBasic4ErOverEz[k] = 0;
120  fLookUpBasic4EphiOverEz[k] = 0;
121  fLookUpBasic4DeltaEz[k] = 0;
122 
123  fLookUpBasic5ErOverEz[k] = 0;
124  fLookUpBasic5EphiOverEz[k] = 0;
125  fLookUpBasic5DeltaEz[k] = 0;
126 
127  fLookUpBasic6ErOverEz[k] = 0;
128  fLookUpBasic6EphiOverEz[k] = 0;
129  fLookUpBasic6DeltaEz[k] = 0;
130  }
131 
132 }
133 
136 
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!
145  delete fLookUpBasic1EphiOverEz[k];
146  delete fLookUpBasic1DeltaEz[k];
147 
148  delete fLookUpBasic2ErOverEz[k]; // does nothing if pointer is zero!
149  delete fLookUpBasic2EphiOverEz[k];
150  delete fLookUpBasic2DeltaEz[k];
151 
152  delete fLookUpBasic3ErOverEz[k]; // does nothing if pointer is zero!
153  delete fLookUpBasic3EphiOverEz[k];
154  delete fLookUpBasic3DeltaEz[k];
155 
156  delete fLookUpBasic4ErOverEz[k]; // does nothing if pointer is zero!
157  delete fLookUpBasic4EphiOverEz[k];
158  delete fLookUpBasic4DeltaEz[k];
159 
160  delete fLookUpBasic5ErOverEz[k]; // does nothing if pointer is zero!
161  delete fLookUpBasic5EphiOverEz[k];
162  delete fLookUpBasic5DeltaEz[k];
163 
164  delete fLookUpBasic6ErOverEz[k]; // does nothing if pointer is zero!
165  delete fLookUpBasic6EphiOverEz[k];
166  delete fLookUpBasic6DeltaEz[k];
167 
168  }
169 }
170 
171 
177 
178  if (corr==NULL) {
179  AliError("Zerro pointer - correction");
180  return kFALSE;
181  }
182  AliTPCFCVoltError3D * corrC = dynamic_cast<AliTPCFCVoltError3D *>(corr);
183  if (corrC == NULL) {
184  AliError(TString::Format("Inconsistent class types: %s\%s",IsA()->GetName(),corr->IsA()->GetName()).Data());
185  return kFALSE;
186  }
187  //
188  for (Int_t isec=0; isec<36; isec++){
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
193  }
194  for (Int_t isec=0; isec<2; isec++){
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 
206 
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
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)
214  Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ;
215  // Correction Terms for effective omegaTau; obtained by a laser calibration run
216  SetOmegaTauT1T2(wt,fT1,fT2);
217 
219 }
220 
221 void AliTPCFCVoltError3D::Update(const TTimeStamp &/*timeStamp*/) {
223 
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
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)
231  Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ;
232  // Correction Terms for effective omegaTau; obtained by a laser calibration run
233  SetOmegaTauT1T2(wt,fT1,fT2);
234 
235 
236 }
237 
238 
239 
240 void AliTPCFCVoltError3D::GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]) {
242 
243  const Double_t kEpsilon=Double_t(FLT_MIN);
244 
245  if (!fInitLookUp) {
246  AliInfo("Lookup table was not initialized! Perform the inizialisation now ...");
248  }
249 
250  static Bool_t forceInit=kTRUE; //temporary needed for back compatibility with old OCDB entries
251  if (forceInit &&fLookUpErOverEz[0]){
252  if (TMath::Abs(fLookUpErOverEz[0]->Sum())<kEpsilon){//temporary needed for back compatibility with old OCDB entries
254  }
255  forceInit=kFALSE;
256  }
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  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  }
275 
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
278 
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 
283  // Get the Er and Ephi field integrals plus the integral over DeltaEz
284  intEr = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi,
286  intEphi = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi,
288  intDeltaEz = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi,
290 
291  // printf("%lf %lf %lf\n",intEr,intEphi,intDeltaEz);
292 
293  // Calculate distorted position
294  if ( r > 0.0 ) {
295  phi = phi + ( fC0*intEphi - fC1*intEr ) / r;
296  r = r + ( fC0*intEr + fC1*intEphi );
297  }
298 
299  // Calculate correction in cartesian coordinates
300  dx[0] = r * TMath::Cos(phi) - x[0];
301  dx[1] = r * TMath::Sin(phi) - x[1];
302  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 }
306 
312 
313  const Int_t order = 1 ; // Linear interpolation = 1, Quadratic = 2
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
319  TMatrixD *arrayofArrayV[kPhiSlices], *arrayofCharge[kPhiSlices] ;
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] ;
325 
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
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  Bool_t needTable[6] ={kFALSE,kFALSE,kFALSE,kFALSE,kFALSE,kFALSE};
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;
360 
361  // shifted Copper rods
362  for ( Int_t rod = 0 ; rod < 18 ; rod++ ) {
363  if (fCopperRodShiftA[rod]!=0) needTable[4]=kTRUE;
364  if (fCopperRodShiftC[rod]!=0) needTable[4]=kTRUE;
365  }
366  // shifted Copper rods
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;
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
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
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
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
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
411  // will be deleted in destructor
412  }
413  }
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
419  // will be deleted in destructor
420  }
421  }
422 
423  // Set bondaries and solve Poisson's equation --------------------------
424 
425  for (Int_t look=0; look<6; look++) {
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 
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)));
436  inner18[0] = 1;
437  } else if (look==1) {
438  AliInfo(Form("OFC - ROD&Strip shift : Filling table (~ %d sec)",3*(int)(kPhiSlices/5)));
439  outer18[0] = 1;
440  } else if (look==2) {
441  AliInfo(Form("IFC - Clip rot. : Filling table (~ %d sec)",3*(int)(kPhiSlices/5)));
442  inner18[0] = 1;
443  } else if (look==3) {
444  AliInfo(Form("OFC - Clip rot. : Filling table (~ %d sec)",3*(int)(kPhiSlices/5)));
445  outer18[0] = 1;
446  } else if (look==4) {
447  AliInfo(Form("IFC - CopperRod shift : Filling table (~ %d sec)",3*(int)(kPhiSlices/5)));
448  inner18[0] = 1;
449  } else if (look==5) {
450  AliInfo(Form("OFC - CopperRod shift : Filling table (~ %d sec)",3*(int)(kPhiSlices/5)));
451  outer18[0] = 1;
452  }
453  // Build potentials on boundary stripes for a rotated clip (SYMMETRY==-1) or a shifted rod (SYMMETRY==1)
454  if (look<4) {
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 ;
458  Int_t ipoint = k/slices ;
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 ;
461  if ( (k%slices) == 0 && symmetry[look] == -1 ) { innerList[k] = 0.0 ; outerList[k] = 0.0 ; }
462  // above, force Zero if Anti-Sym
463  }
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){
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
486  arrayV(i,j) = 0.0 ;
487  charge(i,j) = 0.0 ;
488  if ( i == 0 ) arrayV(i,j) = innerList[k] ;
489  if ( i == kRows-1 ) arrayV(i,j) = outerList[k] ;
490  }
491  }
492  // no charge in the volume
493  for ( Int_t i = 1 ; i < kRows-1 ; i++ ) {
494  for ( Int_t j = 1 ; j < kColumns-1 ; j++ ) {
495  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  if (look==0) {
503  PoissonRelaxation3D( arrayofArrayV, arrayofCharge,
505  kRows, kColumns, kPhiSlices, gridSizePhi, kIterations, symmetry[0]) ;
506  AliInfo("IFC - ROD&Strip shift : done ");
507  } else if (look==1) {
508  PoissonRelaxation3D( arrayofArrayV, arrayofCharge,
510  kRows, kColumns, kPhiSlices, gridSizePhi, kIterations, symmetry[1]) ;
511 
512  AliInfo("OFC - ROD&Strip shift : done ");
513  } else if (look==2) {
514  PoissonRelaxation3D( arrayofArrayV, arrayofCharge,
516  kRows, kColumns, kPhiSlices, gridSizePhi, kIterations, symmetry[2]) ;
517  AliInfo("IFC - Clip rot. : done ");
518  } else if (look==3) {
519  PoissonRelaxation3D( arrayofArrayV, arrayofCharge,
521  kRows, kColumns, kPhiSlices, gridSizePhi, kIterations, symmetry[3]) ;
522  AliInfo("OFC - Clip rot. : done ");
523  } else if (look==4) {
524  PoissonRelaxation3D( arrayofArrayV, arrayofCharge,
526  kRows, kColumns, kPhiSlices, gridSizePhi, kIterations, symmetry[4]) ;
527  AliInfo("IFC - CopperRod shift : done ");
528  } else if (look==5) {
529  PoissonRelaxation3D( arrayofArrayV, arrayofCharge,
531  kRows, kColumns, kPhiSlices, gridSizePhi, kIterations, symmetry[5]) ;
532  AliInfo("OFC - CopperRod shift : done ");
533  }
534 
535  fInitLookUpBasic[look] = kTRUE;
536  }
537  }
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  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 
560  TMatrixF &erOverEz = *fLookUpErOverEz[k] ;
561  TMatrixF &ephiOverEz = *fLookUpEphiOverEz[k];
562  TMatrixF &deltaEz = *fLookUpDeltaEz[k] ;
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)
566  for ( Int_t j = 0 ; j < kNR ; j++ ) {
567  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  for ( Int_t rod = 0 ; rod < 18 ; rod++ ) {
578 
579  erIntegral = ephiIntegral = ezIntegral = 0.0 ;
580 
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
585  phiPrime = phi - rod*kPhiSlicesPerSector*gridSizePhi ;
586 
587  if ( phiPrime < 0 ) phiPrime += TMath::TwoPi() ; // Stay in range from 0 to TwoPi
588 
589  if ( (phiPrime >= 0) && (phiPrime <= kPhiSlices*gridSizePhi) ) {
590 
591  erIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
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 );
597 
598  } else if ( (phiPrime <= TMath::TwoPi()) && (phiPrime >= (TMath::TwoPi()-kPhiSlices*gridSizePhi)) ){
599 
600  phiPrime = TMath::TwoPi() - phiPrime ;
601 
602  erIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
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 );
608 
609  // Flip symmetry of phi integral for symmetric boundary conditions
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 ;
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 ;
631 
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
636  phiPrime = phi - (rod-18)*kPhiSlicesPerSector*gridSizePhi ;
637 
638  if ( phiPrime < 0 ) phiPrime += TMath::TwoPi() ; // Stay in range from 0 to TwoPi
639 
640  if ( (phiPrime >= 0) && (phiPrime <= kPhiSlices*gridSizePhi) ) {
641 
642  erIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
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 );
648 
649  } else if ( (phiPrime <= TMath::TwoPi()) && (phiPrime >= (TMath::TwoPi()-kPhiSlices*gridSizePhi)) ){
650 
651  phiPrime = TMath::TwoPi() - phiPrime ;
652 
653  erIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
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 );
659 
660  // Flip symmetry of phi integral for symmetric boundary conditions
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 ;
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
694  phiPrime = phi - rod*kPhiSlicesPerSector*gridSizePhi ;
695 
696  if ( phiPrime < 0 ) phiPrime += TMath::TwoPi() ; // Stay in range from 0 to TwoPi
697 
698  if ( (phiPrime >= 0) && (phiPrime <= kPhiSlices*gridSizePhi) ) {
699 
700  erIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
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 );
706 
707  } else if ( (phiPrime <= TMath::TwoPi()) && (phiPrime >= (TMath::TwoPi()-kPhiSlices*gridSizePhi)) ){
708 
709  phiPrime = TMath::TwoPi() - phiPrime ;
710 
711  erIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
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 );
717 
718  // Flip symmetry of phi integral for symmetric boundary conditions
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 
723  }
724 
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
738 
739  erIntegral = ephiIntegral = ezIntegral = 0.0 ;
740 
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
745  phiPrime = phi - rod*kPhiSlicesPerSector*gridSizePhi ;
746 
747 
748  if ( phiPrime < 0 ) phiPrime += TMath::TwoPi() ; // Stay in range from 0 to TwoPi
749 
750  if ( (phiPrime >= 0) && (phiPrime <= kPhiSlices*gridSizePhi) ) {
751 
752  erIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
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 );
758 
759  } else if ( (phiPrime <= TMath::TwoPi()) && (phiPrime >= (TMath::TwoPi()-kPhiSlices*gridSizePhi)) ){
760 
761  phiPrime = TMath::TwoPi() - phiPrime ;
762 
763  erIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
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 );
769 
770  // Flip symmetry of phi integral for symmetric boundary conditions
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 
775  }
776 
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 
788  // IFC Cooper rod shift +++++++++++++++++++++++++++++
789  for ( Int_t rod = 0 ; rod < 18 ; rod++ ) {
790 
791  erIntegral = ephiIntegral = ezIntegral = 0.0 ;
792 
793  if ( fCopperRodShiftA[rod] == 0 && fgkZList[i] > 0) continue ;
794  if ( fCopperRodShiftC[rod] == 0 && fgkZList[i] < 0) continue ;
795 
796  // Rotate to a position where we have maps and prepare to sum
797  phiPrime = phi - rod*kPhiSlicesPerSector*gridSizePhi ;
798 
799  if ( phiPrime < 0 ) phiPrime += TMath::TwoPi() ; // Stay in range from 0 to TwoPi
800 
801  if ( (phiPrime >= 0) && (phiPrime <= kPhiSlices*gridSizePhi) ) {
802 
803  erIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
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 );
809 
810  } else if ( (phiPrime <= TMath::TwoPi()) && (phiPrime >= (TMath::TwoPi()-kPhiSlices*gridSizePhi)) ){
811 
812  phiPrime = TMath::TwoPi() - phiPrime ;
813 
814  erIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
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 );
820 
821  // Flip symmetry of phi integral for symmetric boundary conditions
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 ;
825 
826  }
827 
828  if ( fgkZList[i] > 0 ) {
829  erIntegralSum += fCopperRodShiftA[rod]*erIntegral ;
830  ephiIntegralSum += fCopperRodShiftA[rod]*ephiIntegral ;
831  ezIntegralSum += fCopperRodShiftA[rod]*ezIntegral;
832  } else if ( fgkZList[i] < 0 ) {
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++ ) {
841 
842  erIntegral = ephiIntegral = ezIntegral = 0.0 ;
843 
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
848  phiPrime = phi - (rod-18)*kPhiSlicesPerSector*gridSizePhi ;
849 
850  if ( phiPrime < 0 ) phiPrime += TMath::TwoPi() ; // Stay in range from 0 to TwoPi
851 
852  if ( (phiPrime >= 0) && (phiPrime <= kPhiSlices*gridSizePhi) ) {
853 
854  erIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
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 );
860 
861  } else if ( (phiPrime <= TMath::TwoPi()) && (phiPrime >= (TMath::TwoPi()-kPhiSlices*gridSizePhi)) ){
862 
863  phiPrime = TMath::TwoPi() - phiPrime ;
864 
865  erIntegral = Interpolate3DTable(order, r, z, phiPrime, kRows, kColumns, kPhiSlices,
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 );
871 
872  // Flip symmetry of phi integral for symmetric boundary conditions
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 ;
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;
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;
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  for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
904  delete arrayofArrayV[k];
905  delete arrayofCharge[k];
906  }
907 
908  AliInfo(" done");
909  fInitLookUp = kTRUE;
910 
911 }
912 
913 void AliTPCFCVoltError3D::Print(const Option_t* option) const {
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;
921  printf(" : A-side - (Rod & Strips) \n ");
922  for (Int_t i=0; i<36;i++) {
923  if (fRodVoltShiftA[i]!=0) {
924  printf("Rod%2d:%3.1fV ",i,fRodVoltShiftA[i]);
925  count++;
926  }
927  if (count==0&&i==35)
928  printf("-> all at 0 V\n");
929  else if (i==35){
930  printf("\n");
931  count=0;
932  }
933  }
934  printf(" : C-side - (Rod & Strips) \n ");
935  for (Int_t i=0; i<36;i++) {
936  if (fRodVoltShiftC[i]!=0) {
937  printf("Rod%2d:%3.1fV ",i,fRodVoltShiftC[i]);
938  count++;
939  }
940  if (count==0&&i==35)
941  printf("-> all at 0 V\n");
942  else if (i==35){
943  printf("\n");
944  count=0;
945  }
946  }
947 
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++; }
953  if (count==0)
954  printf(" -> no rotated clips \n");
955 
956  count=0;
957  printf(" - Copper ROD shifts (without strips):\n");
958  printf(" : A-side - (Copper Rod shift) \n ");
959  for (Int_t i=0; i<36;i++) {
960  if (fCopperRodShiftA[i]!=0) {
961  printf("Rod%2d:%3.1fV ",i,fCopperRodShiftA[i]);
962  count++;
963  }
964  if (count==0&&i==35)
965  printf("-> all at 0 V\n");
966  else if (i==35){
967  printf("\n");
968  count=0;
969  }
970  }
971  printf(" : C-side - (Copper Rod shift) \n ");
972  for (Int_t i=0; i<36;i++) {
973  if (fCopperRodShiftC[i]!=0) {
974  printf("Rod%2d:%3.1fV ",i,fCopperRodShiftC[i]);
975  count++;
976  }
977  if (count==0&&i==35)
978  printf("-> all at 0 V\n");
979  else if (i==35){
980  printf("\n");
981  count=0;
982  }
983  }
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 }
static AliTPCcalibDB * Instance()
printf("Chi2/npoints = %f\n", TMath::Sqrt(chi2/npoints))
const Double_t kEpsilon
TMatrixD * fLookUpBasic6EphiOverEz[kPhiSlices]
! Array to store electric field integral
Double_t fgkZList[kNZ]
points in the z direction (for the lookup table)
Float_t fRodVoltShiftC[36]
Rod (&strips) shift in Volt (40V~1mm)
AliTPCFCVoltError3D class.
TMatrixD * fLookUpBasic5DeltaEz[kPhiSlices]
! Array to store electric field integral
TMatrixD * fLookUpBasic3DeltaEz[kPhiSlices]
! Array to store electric field integral
virtual void SetOmegaTauT1T2(Float_t omegaTau, Float_t t1, Float_t t2)
Manager and of geomety classes for set: TPC.
Definition: AliTPCParam.h:18
Float_t GetDriftV() const
Definition: AliTPCParam.h:342
Double_t Interpolate3DTable(Int_t order, Double_t x, Double_t y, Double_t z, Int_t nx, Int_t ny, Int_t nz, const Double_t xv[], const Double_t yv[], const Double_t zv[], TMatrixD **arrayofArrays)
TMatrixD * fLookUpBasic2EphiOverEz[kPhiSlices]
! Array to store electric field integral
virtual void Print(const Option_t *option="") const
TMatrixD * fLookUpBasic6DeltaEz[kPhiSlices]
! Array to store electric field integral
Float_t fRotatedClipVoltA[2]
rotated clips at HV rod
static const Double_t fgkOFCRadius
Mean Radius of the Outer Field Cage (252.55 min, 256.45 max) (cm)
AliTPCParam * GetParameters() const
virtual Bool_t AddCorrectionCompact(AliTPCCorrection *corr, Double_t weight)
Double_t SolenoidField() const
Definition: AliMagF.h:67
Float_t fCopperRodShiftC[36]
only Rod shift
virtual void GetCorrection(const Float_t x[], const Short_t roc, Float_t dx[])
Double_t fT1
tensor term of wt - T1
static const Double_t fgkZOffSet
Offset from CE: calculate all distortions closer to CE as if at this point.
TMatrixD * fLookUpBasic3EphiOverEz[kPhiSlices]
! Array to store electric field integral
TMatrixD * fLookUpBasic3ErOverEz[kPhiSlices]
! Array to store electric field integral
TMatrixD * fLookUpBasic4EphiOverEz[kPhiSlices]
! Array to store electric field integral
AliTPCCorrection class.
virtual void Update(const TTimeStamp &timeStamp)
TMatrixD * fLookUpBasic6ErOverEz[kPhiSlices]
! Array to store electric field integral
static const Double_t fgkTPCZ0
nominal gating grid position
#define AliInfo(message)
Definition: AliLog.h:484
Float_t fRotatedClipVoltC[2]
rotated clips at HV rod
TMatrixD * fLookUpBasic4ErOverEz[kPhiSlices]
! Array to store electric field integral
Bool_t fInitLookUpBasic[6]
! flag if the basic lookup was created (shifted Rod (IFC,OFC) or rotated clip (IFC,OFC))
Bool_t fInitLookUp
flag to check if the Look Up table was created (SUM)
TMatrixF * fLookUpEphiOverEz[kNPhi]
Array to store electric field integral (int Er/Ez)
Float_t fCopperRodShiftA[36]
only Rod shift
TMatrixD * fLookUpBasic5ErOverEz[kPhiSlices]
! Array to store electric field integral
TMatrixD * fLookUpBasic2ErOverEz[kPhiSlices]
! Array to store electric field integral
TMatrixD * fLookUpBasic2DeltaEz[kPhiSlices]
! Array to store electric field integral
TMatrixD * fLookUpBasic1DeltaEz[kPhiSlices]
! Array to store electric field integral (int Ez)
TMatrixD * fLookUpBasic1ErOverEz[kPhiSlices]
! Array to store electric field integral (int Er/Ez)
Float_t fRodVoltShiftA[36]
Rod (&strips) shift in Volt (40V~1mm)
Float_t fC1
coefficient C1 (compare Jim Thomas&#39;s notes for definitions)
Double_t fgkPhiList[kNPhi]
points in the phi direction (for the lookup table)
Double_t fT2
tensor term of wt - T2
Double_t vdrift
Definition: DriftKalman.C:98
Double_t fgkRList[kNR]
points in the radial direction (for the lookup table)
TMatrixD * fLookUpBasic1EphiOverEz[kPhiSlices]
! Array to store electric field integral (int Ephi/Ez)
#define AliError(message)
Definition: AliLog.h:591
TMatrixF * fLookUpDeltaEz[kNPhi]
Array to store electric field integral (int Er/Ez)
TMatrixF * fLookUpErOverEz[kNPhi]
Array to store electric field integral (int Er/Ez)
TMatrixD * fLookUpBasic4DeltaEz[kPhiSlices]
! Array to store electric field integral
Float_t fC0
coefficient C0 (compare Jim Thomas&#39;s notes for definitions)
static const Double_t fgkIFCRadius
Mean Radius of the Inner Field Cage ( 82.43 min, 83.70 max) (cm)/hera/alice/wiechula/calib/guiTrees.
class TMatrixT< Double_t > TMatrixD
TMatrixD * fLookUpBasic5EphiOverEz[kPhiSlices]
! Array to store electric field integral
void PoissonRelaxation3D(TMatrixD **arrayofArrayV, TMatrixD **arrayofChargeDensities, TMatrixD **arrayofEroverEz, TMatrixD **arrayofEPhioverEz, TMatrixD **arrayofEz, Int_t rows, Int_t columns, Int_t phislices, Float_t deltaphi, Int_t iterations, Int_t summetry, Bool_t rocDisplacement=kTRUE, IntegrationType integrationType=kIntegral)