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 AliTPCBoundaryVoltError
17 : ///
18 : /// <h2> AliTPCBoundaryVoltError class </h2>
19 : /// This class calculates the space point distortions due to residual voltage errors on
20 : /// the main boundaries of the TPC. For example, the inner vessel of the TPC is shifted
21 : /// by a certain amount, whereas the ROCs on the A side and the C side follow this mechanical
22 : /// shift (at the inner vessel) in z direction. This example can be named "conical deformation"
23 : /// of the TPC field cage (see example below).
24 : ///
25 : /// The boundary conditions can be set via two arrays (A and C side) which contain the
26 : /// residual voltage setting modeling a possible shift or an inhomogeneity on the TPC field
27 : /// cage. In order to avoid that the user splits the Central Electrode (CE), the settings for
28 : /// the C side is taken from the array on the A side (points: A.6 and A.7). The region betweem
29 : /// the points is interpolated linearly onto the boundaries.
30 : ///
31 : /// The class uses the PoissonRelaxation2D (see AliTPCCorrection) to calculate the resulting
32 : /// electrical field inhomogeneities in the (r,z)-plane. Then, the Langevin-integral formalism
33 : /// is used to calculate the space point distortions.
34 : /// Note: This class assumes a homogeneous magnetic field.
35 : ///
36 : /// One has two possibilities when calculating the $dz$ distortions. The resulting distortions
37 : /// are purely due to the change of the drift velocity (along with the change of the drift field)
38 : /// when the SetROCDisplacement is FALSE. This emulates for example a Gating-Grid Voltage offset
39 : /// without moving the ROCs. When the flag is set to TRUE, the ROCs are assumed to be misaligned
40 : /// and the corresponding offset in z is added.
41 : /// 
42 : ///
43 : /// \author Jim Thomas, Stefan Rossegger
44 : /// \date 01/06/2010
45 :
46 :
47 : #include "AliMagF.h"
48 : #include "TGeoGlobalMagField.h"
49 : #include "AliTPCcalibDB.h"
50 : #include "AliTPCParam.h"
51 : #include "AliLog.h"
52 : #include "TMatrixD.h"
53 :
54 : #include "TMath.h"
55 : #include "AliTPCROC.h"
56 : #include "AliTPCBoundaryVoltError.h"
57 :
58 : /// \cond CLASSIMP
59 24 : ClassImp(AliTPCBoundaryVoltError)
60 : /// \endcond
61 :
62 : AliTPCBoundaryVoltError::AliTPCBoundaryVoltError()
63 0 : : AliTPCCorrection("BoundaryVoltError","Boundary Voltage Error"),
64 0 : fC0(0.),fC1(0.),
65 0 : fROCdisplacement(kTRUE),
66 0 : fInitLookUp(kFALSE)
67 0 : {
68 : //
69 : // default constructor
70 : //
71 0 : for (Int_t i=0; i<8; i++){
72 0 : fBoundariesA[i]= 0;
73 0 : if (i<6) fBoundariesC[i]= 0;
74 : }
75 0 : }
76 :
77 0 : AliTPCBoundaryVoltError::~AliTPCBoundaryVoltError() {
78 : /// default destructor
79 :
80 0 : }
81 :
82 :
83 :
84 :
85 : Bool_t AliTPCBoundaryVoltError::AddCorrectionCompact(AliTPCCorrection* corr, Double_t weight){
86 : /// Add correction and make them compact
87 : /// Assumptions:
88 : /// - origin of distortion/correction are additive
89 : /// - only correction ot the same type supported ()
90 :
91 0 : if (corr==NULL) {
92 0 : AliError("Zerro pointer - correction");
93 0 : return kFALSE;
94 : }
95 0 : AliTPCBoundaryVoltError* corrC = dynamic_cast<AliTPCBoundaryVoltError *>(corr);
96 0 : if (corrC == NULL) {
97 0 : AliError(TString::Format("Inconsistent class types: %s\%s",IsA()->GetName(),corr->IsA()->GetName()).Data());
98 0 : return kFALSE;
99 : }
100 0 : if (fROCdisplacement!=corrC->fROCdisplacement){
101 0 : AliError(TString::Format("Inconsistent fROCdisplacement : %s\%s",IsA()->GetName(),corr->IsA()->GetName()).Data());
102 0 : return kFALSE;
103 : }
104 0 : for (Int_t i=0;i <8; i++){
105 0 : fBoundariesA[i]+= corrC->fBoundariesA[i]*weight;
106 0 : fBoundariesC[i]+= corrC->fBoundariesC[i]*weight;
107 : }
108 : //
109 0 : return kTRUE;
110 0 : }
111 :
112 :
113 :
114 :
115 : void AliTPCBoundaryVoltError::Init() {
116 : /// Initialization funtion
117 :
118 0 : AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
119 0 : if (!magF) AliError("Magneticd field - not initialized");
120 0 : Double_t bzField = magF->SolenoidField()/10.; //field in T
121 0 : AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
122 0 : if (!param) AliError("Parameters - not initialized");
123 0 : Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us] // From dataBase: to be updated: per second (ideally)
124 : Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
125 0 : Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ;
126 : // Correction Terms for effective omegaTau; obtained by a laser calibration run
127 0 : SetOmegaTauT1T2(wt,fT1,fT2);
128 :
129 0 : InitBoundaryVoltErrorDistortion();
130 0 : }
131 :
132 : void AliTPCBoundaryVoltError::Update(const TTimeStamp &/*timeStamp*/) {
133 : /// Update function
134 :
135 0 : AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
136 0 : if (!magF) AliError("Magneticd field - not initialized");
137 0 : Double_t bzField = magF->SolenoidField()/10.; //field in T
138 0 : AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
139 0 : if (!param) AliError("Parameters - not initialized");
140 0 : Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us] // From dataBase: to be updated: per second (ideally)
141 : Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
142 0 : Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ;
143 : // Correction Terms for effective omegaTau; obtained by a laser calibration run
144 0 : SetOmegaTauT1T2(wt,fT1,fT2);
145 :
146 0 : }
147 :
148 :
149 :
150 : void AliTPCBoundaryVoltError::GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]) {
151 : /// Calculates the correction due e.g. residual voltage errors on the TPC boundaries
152 :
153 0 : if (!fInitLookUp) {
154 0 : AliInfo("Lookup table was not initialized! Perform the inizialisation now ...");
155 0 : InitBoundaryVoltErrorDistortion();
156 0 : }
157 :
158 : Int_t order = 1 ; // FIXME: hardcoded? Linear interpolation = 1, Quadratic = 2
159 : // note that the poisson solution was linearly mirroed on this grid!
160 0 : Double_t intEr, intEphi, intdEz ;
161 : Double_t r, phi, z ;
162 : Int_t sign;
163 :
164 0 : r = TMath::Sqrt( x[0]*x[0] + x[1]*x[1] ) ;
165 0 : phi = TMath::ATan2(x[1],x[0]) ;
166 0 : if ( phi < 0 ) phi += TMath::TwoPi() ; // Table uses phi from 0 to 2*Pi
167 0 : z = x[2] ; // Create temporary copy of x[2]
168 :
169 0 : if ( (roc%36) < 18 ) {
170 : sign = 1; // (TPC A side)
171 0 : } else {
172 : sign = -1; // (TPC C side)
173 : }
174 :
175 0 : if ( sign==1 && z < fgkZOffSet ) z = fgkZOffSet; // Protect against discontinuity at CE
176 0 : if ( sign==-1 && z > -fgkZOffSet ) z = -fgkZOffSet; // Protect against discontinuity at CE
177 :
178 :
179 : intEphi = 0.0; // Efield is symmetric in phi - 2D calculation
180 :
181 0 : if ( (sign==1 && z<0) || (sign==-1 && z>0) ) // just a consistency check
182 0 : AliError("ROC number does not correspond to z coordinate! Calculation of distortions is most likely wrong!");
183 :
184 : // Get the E field integral
185 0 : Interpolate2DEdistortion( order, r, z, fLookUpErOverEz, intEr );
186 : // Get DeltaEz field integral
187 0 : Interpolate2DEdistortion( order, r, z, fLookUpDeltaEz, intdEz );
188 :
189 : // Calculate distorted position
190 0 : if ( r > 0.0 ) {
191 0 : phi = phi + ( fC0*intEphi - fC1*intEr ) / r;
192 0 : r = r + ( fC0*intEr + fC1*intEphi );
193 0 : }
194 :
195 : // Calculate correction in cartesian coordinates
196 0 : dx[0] = r * TMath::Cos(phi) - x[0];
197 0 : dx[1] = r * TMath::Sin(phi) - x[1];
198 0 : dx[2] = intdEz; // z distortion - (internally scaled with driftvelocity dependency
199 : // on the Ez field plus the actual ROC misalignment (if set TRUE)
200 :
201 :
202 0 : }
203 :
204 : void AliTPCBoundaryVoltError::InitBoundaryVoltErrorDistortion() {
205 : /// Initialization of the Lookup table which contains the solutions of the
206 : /// Dirichlet boundary problem
207 :
208 0 : const Float_t gridSizeR = (fgkOFCRadius-fgkIFCRadius) / (kRows-1) ;
209 0 : const Float_t gridSizeZ = fgkTPCZ0 / (kColumns-1) ;
210 :
211 0 : TMatrixD voltArrayA(kRows,kColumns), voltArrayC(kRows,kColumns); // boundary vectors
212 0 : TMatrixD chargeDensity(kRows,kColumns); // dummy charge
213 0 : TMatrixD arrayErOverEzA(kRows,kColumns), arrayErOverEzC(kRows,kColumns); // solution
214 0 : TMatrixD arrayDeltaEzA(kRows,kColumns), arrayDeltaEzC(kRows,kColumns); // solution
215 :
216 0 : Double_t rList[kRows], zedList[kColumns] ;
217 :
218 : // Fill arrays with initial conditions. V on the boundary and ChargeDensity in the volume.
219 0 : for ( Int_t j = 0 ; j < kColumns ; j++ ) {
220 0 : Double_t zed = j*gridSizeZ ;
221 0 : zedList[j] = zed ;
222 0 : for ( Int_t i = 0 ; i < kRows ; i++ ) {
223 0 : Double_t radius = fgkIFCRadius + i*gridSizeR ;
224 0 : rList[i] = radius ;
225 0 : voltArrayA(i,j) = 0; // Initialize voltArrayA to zero
226 0 : voltArrayC(i,j) = 0; // Initialize voltArrayC to zero
227 0 : chargeDensity(i,j) = 0; // Initialize ChargeDensity to zero - not used in this class
228 : }
229 : }
230 :
231 :
232 : // check if boundary values are the same for the C side (for later, saving some calculation time)
233 :
234 : Int_t symmetry = -1; // assume that A and C side are identical (but anti-symmetric!) // e.g conical deformation
235 0 : Int_t sVec[8];
236 :
237 : // check if boundaries are different (regardless the sign)
238 0 : for (Int_t i=0; i<8; i++) {
239 0 : if (TMath::Abs(TMath::Abs(fBoundariesA[i]) - TMath::Abs(fBoundariesC[i])) > 1e-5)
240 0 : symmetry = 0;
241 0 : sVec[i] = (Int_t)( TMath::Sign((Float_t)1.,fBoundariesA[i]) * TMath::Sign((Float_t)1.,fBoundariesC[i]));
242 : }
243 0 : if (symmetry==-1) { // still the same values?
244 : // check the kind of symmetry , if even ...
245 0 : if (sVec[0]==1 && sVec[1]==1 && sVec[2]==1 && sVec[3]==1 && sVec[4]==1 && sVec[5]==1 && sVec[6]==1 && sVec[7]==1 )
246 0 : symmetry = 1;
247 0 : else if (sVec[0]==-1 && sVec[1]==-1 && sVec[2]==-1 && sVec[3]==-1 && sVec[4]==-1 && sVec[5]==-1 && sVec[6]==-1 && sVec[7]==-1 )
248 0 : symmetry = -1;
249 : else
250 : symmetry = 0; // some of the values differ in the sign -> neither symmetric nor antisymmetric
251 : }
252 :
253 :
254 :
255 : // Solve the electrosatic problem in 2D
256 :
257 : // Fill the complete Boundary vectors
258 : // Start at IFC at CE and work anti-clockwise through IFC, ROC, OFC, and CE (clockwise for C side)
259 0 : for ( Int_t j = 0 ; j < kColumns ; j++ ) {
260 0 : Double_t zed = j*gridSizeZ ;
261 0 : for ( Int_t i = 0 ; i < kRows ; i++ ) {
262 0 : Double_t radius = fgkIFCRadius + i*gridSizeR ;
263 :
264 : // A side boundary vectors
265 0 : if ( i == 0 ) voltArrayA(i,j) += zed *((fBoundariesA[1]-fBoundariesA[0])/((kColumns-1)*gridSizeZ))
266 0 : + fBoundariesA[0] ; // IFC
267 0 : if ( j == kColumns-1 ) voltArrayA(i,j) += (radius-fgkIFCRadius)*((fBoundariesA[3]-fBoundariesA[2])/((kRows-1)*gridSizeR))
268 0 : + fBoundariesA[2] ; // ROC
269 0 : if ( i == kRows-1 ) voltArrayA(i,j) += zed *((fBoundariesA[4]-fBoundariesA[5])/((kColumns-1)*gridSizeZ))
270 0 : + fBoundariesA[5] ; // OFC
271 0 : if ( j == 0 ) voltArrayA(i,j) += (radius-fgkIFCRadius)*((fBoundariesA[6]-fBoundariesA[7])/((kRows-1)*gridSizeR))
272 0 : + fBoundariesA[7] ; // CE
273 :
274 0 : if (symmetry==0) {
275 : // C side boundary vectors
276 0 : if ( i == 0 ) voltArrayC(i,j) += zed *((fBoundariesC[1]-fBoundariesC[0])/((kColumns-1)*gridSizeZ))
277 0 : + fBoundariesC[0] ; // IFC
278 0 : if ( j == kColumns-1 ) voltArrayC(i,j) += (radius-fgkIFCRadius)*((fBoundariesC[3]-fBoundariesC[2])/((kRows-1)*gridSizeR))
279 0 : + fBoundariesC[2] ; // ROC
280 0 : if ( i == kRows-1 ) voltArrayC(i,j) += zed *((fBoundariesC[4]-fBoundariesC[5])/((kColumns-1)*gridSizeZ))
281 0 : + fBoundariesC[5] ; // OFC
282 0 : if ( j == 0 ) voltArrayC(i,j) += (radius-fgkIFCRadius)*((fBoundariesC[6]-fBoundariesC[7])/((kRows-1)*gridSizeR))
283 0 : + fBoundariesC[7] ; // CE
284 : }
285 :
286 : }
287 : }
288 :
289 0 : voltArrayA(0,0) *= 0.5 ; // Use average boundary condition at corner
290 0 : voltArrayA(kRows-1,0) *= 0.5 ; // Use average boundary condition at corner
291 0 : voltArrayA(0,kColumns-1) *= 0.5 ; // Use average boundary condition at corner
292 0 : voltArrayA(kRows-1,kColumns-1)*= 0.5 ; // Use average boundary condition at corner
293 :
294 0 : if (symmetry==0) {
295 0 : voltArrayC(0,0) *= 0.5 ; // Use average boundary condition at corner
296 0 : voltArrayC(kRows-1,0) *= 0.5 ; // Use average boundary condition at corner
297 0 : voltArrayC(0,kColumns-1) *= 0.5 ; // Use average boundary condition at corner
298 0 : voltArrayC(kRows-1,kColumns-1)*= 0.5 ; // Use average boundary condition at corner
299 0 : }
300 :
301 :
302 : // always solve the problem on the A side
303 0 : PoissonRelaxation2D( voltArrayA, chargeDensity, arrayErOverEzA, arrayDeltaEzA,
304 0 : kRows, kColumns, kIterations, fROCdisplacement ) ;
305 :
306 0 : if (symmetry!=0) { // A and C side are the same ("anti-symmetric" or "symmetric")
307 0 : for ( Int_t j = 0 ; j < kColumns ; j++ ) {
308 0 : for ( Int_t i = 0 ; i < kRows ; i++ ) {
309 0 : arrayErOverEzC(i,j) = symmetry*arrayErOverEzA(i,j);
310 0 : arrayDeltaEzC(i,j) = -symmetry*arrayDeltaEzA(i,j);
311 : }
312 : }
313 0 : } else if (symmetry==0) { // A and C side are different - Solve the problem on the C side too
314 0 : PoissonRelaxation2D( voltArrayC, chargeDensity, arrayErOverEzC, arrayDeltaEzC,
315 0 : kRows, kColumns, kIterations, fROCdisplacement ) ;
316 0 : for ( Int_t j = 0 ; j < kColumns ; j++ ) {
317 0 : for ( Int_t i = 0 ; i < kRows ; i++ ) {
318 0 : arrayDeltaEzC(i,j) = -arrayDeltaEzC(i,j); // negative z coordinate!
319 : }
320 : }
321 0 : }
322 :
323 : // Interpolate results onto standard grid for Electric Fields
324 0 : Int_t ilow=0, jlow=0 ;
325 : Double_t z,r;
326 : Float_t saveEr[2] ;
327 : Float_t saveEz[2] ;
328 0 : for ( Int_t i = 0 ; i < kNZ ; ++i ) {
329 0 : z = TMath::Abs( fgkZList[i] ) ;
330 0 : for ( Int_t j = 0 ; j < kNR ; ++j ) {
331 : // Linear interpolation !!
332 0 : r = fgkRList[j] ;
333 0 : Search( kRows, rList, r, ilow ) ; // Note switch - R in rows and Z in columns
334 0 : Search( kColumns, zedList, z, jlow ) ;
335 0 : if ( ilow < 0 ) ilow = 0 ; // check if out of range
336 0 : if ( jlow < 0 ) jlow = 0 ;
337 0 : if ( ilow + 1 >= kRows - 1 ) ilow = kRows - 2 ;
338 0 : if ( jlow + 1 >= kColumns - 1 ) jlow = kColumns - 2 ;
339 :
340 0 : if (fgkZList[i]>0) { // A side solution
341 0 : saveEr[0] = arrayErOverEzA(ilow,jlow) +
342 0 : (arrayErOverEzA(ilow,jlow+1)-arrayErOverEzA(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
343 0 : saveEr[1] = arrayErOverEzA(ilow+1,jlow) +
344 0 : (arrayErOverEzA(ilow+1,jlow+1)-arrayErOverEzA(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;
345 0 : saveEz[0] = arrayDeltaEzA(ilow,jlow) +
346 0 : (arrayDeltaEzA(ilow,jlow+1)-arrayDeltaEzA(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
347 0 : saveEz[1] = arrayDeltaEzA(ilow+1,jlow) +
348 0 : (arrayDeltaEzA(ilow+1,jlow+1)-arrayDeltaEzA(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;
349 :
350 0 : } else if (fgkZList[i]<0) { // C side solution
351 0 : saveEr[0] = arrayErOverEzC(ilow,jlow) +
352 0 : (arrayErOverEzC(ilow,jlow+1)-arrayErOverEzC(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
353 0 : saveEr[1] = arrayErOverEzC(ilow+1,jlow) +
354 0 : (arrayErOverEzC(ilow+1,jlow+1)-arrayErOverEzC(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;
355 0 : saveEz[0] = arrayDeltaEzC(ilow,jlow) +
356 0 : (arrayDeltaEzC(ilow,jlow+1)-arrayDeltaEzC(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
357 0 : saveEz[1] = arrayDeltaEzC(ilow+1,jlow) +
358 0 : (arrayDeltaEzC(ilow+1,jlow+1)-arrayDeltaEzC(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;
359 :
360 0 : } else {
361 0 : AliWarning("Field calculation at z=0 (CE) is not allowed!");
362 : saveEr[0]=0; saveEr[1]=0;
363 : saveEz[0]=0; saveEz[1]=0;
364 : }
365 0 : fLookUpErOverEz[i][j] = saveEr[0] + (saveEr[1]-saveEr[0])*(r-rList[ilow])/gridSizeR ;
366 0 : fLookUpDeltaEz[i][j] = saveEz[0] + (saveEz[1]-saveEz[0])*(r-rList[ilow])/gridSizeR ;
367 : }
368 : }
369 :
370 0 : voltArrayA.Clear();
371 0 : voltArrayC.Clear();
372 0 : chargeDensity.Clear();
373 0 : arrayErOverEzA.Clear();
374 0 : arrayErOverEzC.Clear();
375 0 : arrayDeltaEzA.Clear();
376 0 : arrayDeltaEzC.Clear();
377 :
378 0 : fInitLookUp = kTRUE;
379 :
380 0 : }
381 :
382 : void AliTPCBoundaryVoltError::Print(const Option_t* option) const {
383 : /// Print function to check the settings of the boundary vectors
384 : /// option=="a" prints the C0 and C1 coefficents for calibration purposes
385 :
386 0 : TString opt = option; opt.ToLower();
387 0 : printf("%s\n",GetTitle());
388 0 : printf(" - Voltage settings (on the TPC boundaries) - linearly interpolated\n");
389 0 : printf(" : A-side (anti-clockwise)\n");
390 0 : printf(" (0,1):\t IFC (CE) : %3.1f V \t IFC (ROC): %3.1f V \n",fBoundariesA[0],fBoundariesA[1]);
391 0 : printf(" (2,3):\t ROC (IFC): %3.1f V \t ROC (OFC): %3.1f V \n",fBoundariesA[2],fBoundariesA[3]);
392 0 : printf(" (4,5):\t OFC (ROC): %3.1f V \t OFC (CE) : %3.1f V \n",fBoundariesA[4],fBoundariesA[5]);
393 0 : printf(" (6,7):\t CE (OFC): %3.1f V \t CE (IFC): %3.1f V \n",fBoundariesA[6],fBoundariesA[7]);
394 0 : printf(" : C-side (clockwise)\n");
395 0 : printf(" (0,1):\t IFC (CE) : %3.1f V \t IFC (ROC): %3.1f V \n",fBoundariesC[0],fBoundariesC[1]);
396 0 : printf(" (2,3):\t ROC (IFC): %3.1f V \t ROC (OFC): %3.1f V \n",fBoundariesC[2],fBoundariesC[3]);
397 0 : printf(" (4,5):\t OFC (ROC): %3.1f V \t OFC (CE) : %3.1f V \n",fBoundariesC[4],fBoundariesC[5]);
398 0 : printf(" (6,7):\t CE (OFC): %3.1f V \t CE (IFC): %3.1f V \n",fBoundariesC[6],fBoundariesC[7]);
399 :
400 : // Check wether the settings of the Central Electrode agree (on the A and C side)
401 : // Note: they have to be antisymmetric!
402 0 : if (( TMath::Abs(fBoundariesA[6]+fBoundariesC[6])>1e-5) || ( TMath::Abs(fBoundariesA[7]+fBoundariesC[7])>1e-5 ) ){
403 0 : AliWarning("Boundary parameters for the Central Electrode (CE) are not anti-symmetric! HOW DID YOU MANAGE THAT?");
404 0 : AliWarning("Congratulations, you just splitted the Central Electrode of the TPC!");
405 0 : AliWarning("Non-physical settings of the boundary parameter at the Central Electrode");
406 : }
407 :
408 0 : if (opt.Contains("a")) { // Print all details
409 0 : printf(" - T1: %1.4f, T2: %1.4f \n",fT1,fT2);
410 0 : printf(" - C1: %1.4f, C0: %1.4f \n",fC1,fC0);
411 : }
412 :
413 0 : if (!fInitLookUp)
414 0 : AliError("Lookup table was not initialized! You should do InitBoundaryVoltErrorDistortion() ...");
415 :
416 0 : }
417 :
418 :
419 : void AliTPCBoundaryVoltError::SetBoundariesA(Float_t boundariesA[8]){
420 : /// set voltage errors on the TPC boundaries - A side
421 : ///
422 : /// Start at IFC at the Central electrode and work anti-clockwise (clockwise for C side) through
423 : /// IFC, ROC, OFC, and CE. The boundary conditions are currently defined to be a linear
424 : /// interpolation between pairs of numbers in the Boundary (e.g. fBoundariesA) vector.
425 : /// The first pair of numbers represent the beginning and end of the Inner Field cage, etc.
426 : /// The unit of the error potential vector is [Volt], whereas 1mm shift of the IFC would
427 : /// correspond to ~ 40 V
428 : ///
429 : /// Note: The setting for the CE will be passed to the C side!
430 :
431 0 : for (Int_t i=0; i<8; i++) {
432 0 : fBoundariesA[i]= boundariesA[i];
433 0 : if (i>5) fBoundariesC[i]= -boundariesA[i]; // setting for the CE is passed to C side
434 : }
435 0 : fInitLookUp=kFALSE;
436 0 : }
437 : void AliTPCBoundaryVoltError::SetBoundariesC(Float_t boundariesC[6]){
438 : /// set voltage errors on the TPC boundaries - C side
439 : ///
440 : /// Start at IFC at the Central electrode and work clockwise (for C side) through
441 : /// IFC, ROC and OFC. The boundary conditions are currently defined to be a linear
442 : /// interpolation between pairs of numbers in the Boundary (e.g. fBoundariesC) vector.
443 : /// The first pair of numbers represent the beginning and end of the Inner Field cage, etc.
444 : /// The unit of the error potential vector is [Volt], whereas 1mm shift of the IFC would
445 : /// correspond to ~ 40 V
446 : ///
447 : /// Note: The setting for the CE will be taken from the A side (pos 6 and 7)!
448 :
449 0 : for (Int_t i=0; i<6; i++) {
450 0 : fBoundariesC[i]= boundariesC[i];
451 : }
452 0 : fInitLookUp=kFALSE;
453 0 : }
|