]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCBoundaryVoltError.cxx
Adding the function:
[u/mrichter/AliRoot.git] / TPC / AliTPCBoundaryVoltError.cxx
CommitLineData
1b923461 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16//////////////////////////////////////////////////////////////////////////////
17// //
18// AliTPCBoundaryVoltError class //
19// The class calculates the space point distortions due to residual voltage //
20// errors on the main boundaries of the TPC. For example, the inner vessel //
21// of the TPC is shifted by a certain amount, whereas the ROCs on the A side//
22// and the ROCs on the C side follow this mechanical shift (at the inner //
23// vessel) in z direction (see example below). This example is commonly //
24// named "conical deformation" of the TPC field cage. //
25// //
26// The class allows "effective Omega Tau" corrections. //
27// //
28// NOTE: This class is not capable of calculation z distortions due to //
29// drift velocity change in dependence of the electric field!!! //
30// //
31// date: 01/06/2010 //
32// Authors: Jim Thomas, Stefan Rossegger //
33// //
34// Example usage (e.g +1mm shift of "conical deformation") //
35// AliTPCBoundaryVoltError bve; //
36// Float_t boundA[8] = {-40,-40,-40,0,0,0,0,-40}; // voltages A-side //
37// Float_t boundC[6] = { 40, 40, 40,0,0,0}; // voltages C-side //
38// bve.SetBoundariesA(boundA); //
39// bve.SetBoundariesC(boundC); //
40// bve.SetOmegaTauT1T2(0.32,1.,1.); // values ideally from OCDB //
41// // initialization of the look up //
42// bve.InitBoundaryVoltErrorDistortion(); //
43// // plot dRPhi distortions ... //
44// bve.CreateHistoDRPhiinZR(1.,100,100)->Draw("surf2"); //
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
58ClassImp(AliTPCBoundaryVoltError)
59
60AliTPCBoundaryVoltError::AliTPCBoundaryVoltError()
61 : AliTPCCorrection("BoundaryVoltError","Boundary Voltage Error"),
62 fC0(0.),fC1(0.),
63 fInitLookUp(kFALSE)
64{
65 //
66 // default constructor
67 //
68 for (Int_t i=0; i<8; i++){
69 fBoundariesA[i]= 0;
70 if (i<6) fBoundariesC[i]= 0;
71 }
72}
73
74AliTPCBoundaryVoltError::~AliTPCBoundaryVoltError() {
75 //
76 // default destructor
77 //
78}
79
80
81
82void AliTPCBoundaryVoltError::Init() {
83 //
84 // Initialization funtion
85 //
86
87 AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
88 if (!magF) AliError("Magneticd field - not initialized");
89 Double_t bzField = magF->SolenoidField()/10.; //field in T
90 AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
91 if (!param) AliError("Parameters - not initialized");
92 Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us] // From dataBase: to be updated: per second (ideally)
93 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
94 Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ;
95 // Correction Terms for effective omegaTau; obtained by a laser calibration run
96 SetOmegaTauT1T2(wt,fT1,fT2);
97
98 InitBoundaryVoltErrorDistortion();
99}
100
101void AliTPCBoundaryVoltError::Update(const TTimeStamp &/*timeStamp*/) {
102 //
103 // Update function
104 //
105 AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
106 if (!magF) AliError("Magneticd field - not initialized");
107 Double_t bzField = magF->SolenoidField()/10.; //field in T
108 AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
109 if (!param) AliError("Parameters - not initialized");
110 Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us] // From dataBase: to be updated: per second (ideally)
111 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
112 Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ;
113 // Correction Terms for effective omegaTau; obtained by a laser calibration run
114 SetOmegaTauT1T2(wt,fT1,fT2);
115
116
117}
118
119
120
121void AliTPCBoundaryVoltError::GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]) {
122 //
123 // Calculates the correction due e.g. residual voltage errors on the TPC boundaries
124 //
125
126 if (!fInitLookUp) AliError("Lookup table was not initialized! You should do InitBoundaryVoltErrorDistortion() ...");
127
128 Int_t order = 1 ; // FIXME: hardcoded? Linear interpolation = 1, Quadratic = 2
129 // note that the poisson solution was linearly mirroed on this grid!
130 Double_t intEr, intEphi ;
131 Double_t r, phi, z ;
132 Int_t sign;
133
134 r = TMath::Sqrt( x[0]*x[0] + x[1]*x[1] ) ;
135 phi = TMath::ATan2(x[1],x[0]) ;
136 if ( phi < 0 ) phi += TMath::TwoPi() ; // Table uses phi from 0 to 2*Pi
137 z = x[2] ; // Create temporary copy of x[2]
138
139 if ( (roc%36) < 18 ) {
140 sign = 1; // (TPC A side)
141 } else {
142 sign = -1; // (TPC C side)
143 }
144
145 if ( sign==1 && z < fgkZOffSet ) z = fgkZOffSet; // Protect against discontinuity at CE
146 if ( sign==-1 && z > -fgkZOffSet ) z = -fgkZOffSet; // Protect against discontinuity at CE
147
148
149 intEphi = 0.0; // Efield is symmetric in phi - 2D calculation
150
151 if ( (sign==1 && z<0) || (sign==-1 && z>0) ) // just a consistency check
152 AliError("ROC number does not correspond to z coordinate! Calculation of distortions is most likely wrong!");
153
154 // Get the E field integral
155 Interpolate2DEdistortion( order, r, z, fLookUpErOverEz, intEr );
156
157 // Calculate distorted position
158 if ( r > 0.0 ) {
159 phi = phi + ( fC0*intEphi - fC1*intEr ) / r;
160 r = r + ( fC0*intEr + fC1*intEphi );
161 }
162
163 // Calculate correction in cartesian coordinates
164 dx[0] = r * TMath::Cos(phi) - x[0];
165 dx[1] = r * TMath::Sin(phi) - x[1];
166 dx[2] = 0.; // z distortion not implemented (1st order distortions)
167
168}
169
170void AliTPCBoundaryVoltError::InitBoundaryVoltErrorDistortion() {
171 //
172 // Initialization of the Lookup table which contains the solutions of the
173 // Dirichlet boundary problem
174 //
175
176 const Float_t gridSizeR = (fgkOFCRadius-fgkIFCRadius) / (kRows-1) ;
177 const Float_t gridSizeZ = fgkTPCZ0 / (kColumns-1) ;
178
179 TMatrixD voltArrayA(kRows,kColumns), voltArrayC(kRows,kColumns); // boundary vectors
180 TMatrixD chargeDensity(kRows,kColumns); // dummy charge
181 TMatrixD arrayErOverEzA(kRows,kColumns), arrayErOverEzC(kRows,kColumns); // solution
182
183 Double_t rList[kRows], zedList[kColumns] ;
184
185 // Fill arrays with initial conditions. V on the boundary and ChargeDensity in the volume.
186 for ( Int_t j = 0 ; j < kColumns ; j++ ) {
187 Double_t zed = j*gridSizeZ ;
188 zedList[j] = zed ;
189 for ( Int_t i = 0 ; i < kRows ; i++ ) {
190 Double_t radius = fgkIFCRadius + i*gridSizeR ;
191 rList[i] = radius ;
192 voltArrayA(i,j) = 0; // Initialize voltArrayA to zero
193 voltArrayC(i,j) = 0; // Initialize voltArrayC to zero
194 chargeDensity(i,j) = 0; // Initialize ChargeDensity to zero - not used in this class
195 }
196 }
197
198
199 // check if boundary values are the same for the C side (for later, saving some calculation time)
200
201 Int_t symmetry = -1; // assume that A and C side are identical (but anti-symmetric!) // e.g conical deformation
202 Int_t sVec[8];
203
204 // check if boundaries are different (regardless the sign)
205 for (Int_t i=0; i<8; i++) {
206 if ((TMath::Abs(fBoundariesA[i]) - TMath::Abs(fBoundariesC[i])) > 1e-5) symmetry = 0;
207 sVec[i] = (Int_t) (TMath::Sign((Float_t)1.,fBoundariesA[i])*TMath::Sign((Float_t)1.,fBoundariesC[i])); // == -1 for anti-symmetry
208 }
209 if (symmetry==-1) { // still the same values?
210 // check the kind of symmetry , if even ...
211 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 )
212 symmetry = 1;
213 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 )
214 symmetry = -1;
215 else
216 symmetry = 0; // some of the values differ in the sign -> neither symmetric nor antisymmetric
217 }
218
219
220
221 // Solve the electrosatic problem in 2D
222
223 // Fill the complete Boundary vectors
224 // Start at IFC at CE and work anti-clockwise through IFC, ROC, OFC, and CE (clockwise for C side)
225 for ( Int_t j = 0 ; j < kColumns ; j++ ) {
226 Double_t zed = j*gridSizeZ ;
227 for ( Int_t i = 0 ; i < kRows ; i++ ) {
228 Double_t radius = fgkIFCRadius + i*gridSizeR ;
229
230 // A side boundary vectors
231 if ( i == 0 ) voltArrayA(i,j) += zed *((fBoundariesA[1]-fBoundariesA[0])/((kColumns-1)*gridSizeZ))
232 + fBoundariesA[0] ; // IFC
233 if ( j == kColumns-1 ) voltArrayA(i,j) += radius*((fBoundariesA[3]-fBoundariesA[2])/((kRows-1)*gridSizeR+fgkIFCRadius))
234 + fBoundariesA[2] ; // ROC
235 if ( i == kRows-1 ) voltArrayA(i,j) += zed *((fBoundariesA[4]-fBoundariesA[5])/((kColumns-1)*gridSizeZ))
236 + fBoundariesA[5] ; // OFC
237 if ( j == 0 ) voltArrayA(i,j) += radius*((fBoundariesA[6]-fBoundariesA[7])/((kRows-1)*gridSizeR+fgkIFCRadius))
238 + fBoundariesA[7] ; // CE
239
240 if (symmetry==0) {
241 // C side boundary vectors
242 if ( i == 0 ) voltArrayC(i,j) += zed *((fBoundariesC[1]-fBoundariesC[0])/((kColumns-1)*gridSizeZ))
243 + fBoundariesC[0] ; // IFC
244 if ( j == kColumns-1 ) voltArrayC(i,j) += radius*((fBoundariesC[3]-fBoundariesC[2])/((kRows-1)*gridSizeR+fgkIFCRadius))
245 + fBoundariesC[2] ; // ROC
246 if ( i == kRows-1 ) voltArrayC(i,j) += zed *((fBoundariesC[4]-fBoundariesC[5])/((kColumns-1)*gridSizeZ))
247 + fBoundariesC[5] ; // OFC
248 if ( j == 0 ) voltArrayC(i,j) += radius*((fBoundariesC[6]-fBoundariesC[7])/((kRows-1)*gridSizeR+fgkIFCRadius))
249 + fBoundariesC[7] ; // CE
250
251 }
252 }
253 }
254
255 voltArrayA(0,0) *= 0.5 ; // Use average boundary condition at corner
256 voltArrayA(kRows-1,0) *= 0.5 ; // Use average boundary condition at corner
257 voltArrayA(0,kColumns-1) *= 0.5 ; // Use average boundary condition at corner
258 voltArrayA(kRows-1,kColumns-1)*= 0.5 ; // Use average boundary condition at corner
259
260 if (symmetry==0) {
261 voltArrayC(0,0) *= 0.5 ; // Use average boundary condition at corner
262 voltArrayC(kRows-1,0) *= 0.5 ; // Use average boundary condition at corner
263 voltArrayC(0,kColumns-1) *= 0.5 ; // Use average boundary condition at corner
264 voltArrayC(kRows-1,kColumns-1)*= 0.5 ; // Use average boundary condition at corner
265 }
266
267
268 // always solve the problem on the A side
269 PoissonRelaxation2D( voltArrayA, chargeDensity, arrayErOverEzA, kRows, kColumns, kIterations ) ;
270
271 if (symmetry!=0) { // A and C side are the same ("anti-symmetric" or "symmetric")
272 for ( Int_t j = 0 ; j < kColumns ; j++ ) {
273 for ( Int_t i = 0 ; i < kRows ; i++ ) {
274 arrayErOverEzC(i,j) = symmetry*arrayErOverEzA(i,j);
275 }
276 }
277 } else if (symmetry==0) { // A and C side are different - Solve the problem on the C side too
278 PoissonRelaxation2D( voltArrayC, chargeDensity, arrayErOverEzC, kRows, kColumns, kIterations ) ;
279 }
280
281 //Interpolate results onto standard grid for Electric Fields
282 Int_t ilow=0, jlow=0 ;
283 Double_t z,r;
284 Float_t saveEr[2] ;
285 for ( Int_t i = 0 ; i < kNZ ; ++i ) {
286 z = TMath::Abs( fgkZList[i] ) ;
287 for ( Int_t j = 0 ; j < kNR ; ++j ) {
288 // Linear interpolation !!
289 r = fgkRList[j] ;
290 Search( kRows, rList, r, ilow ) ; // Note switch - R in rows and Z in columns
291 Search( kColumns, zedList, z, jlow ) ;
292 if ( ilow < 0 ) ilow = 0 ; // check if out of range
293 if ( jlow < 0 ) jlow = 0 ;
294 if ( ilow + 1 >= kRows - 1 ) ilow = kRows - 2 ;
295 if ( jlow + 1 >= kColumns - 1 ) jlow = kColumns - 2 ;
296
297 if (fgkZList[i]>0) { // A side solution
298 saveEr[0] = arrayErOverEzA(ilow,jlow) +
299 (arrayErOverEzA(ilow,jlow+1)-arrayErOverEzA(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
300 saveEr[1] = arrayErOverEzA(ilow+1,jlow) +
301 (arrayErOverEzA(ilow+1,jlow+1)-arrayErOverEzA(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;
302 } else if (fgkZList[i]<0) { // C side solution
303 saveEr[0] = arrayErOverEzC(ilow,jlow) +
304 (arrayErOverEzC(ilow,jlow+1)-arrayErOverEzC(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
305 saveEr[1] = arrayErOverEzC(ilow+1,jlow) +
306 (arrayErOverEzC(ilow+1,jlow+1)-arrayErOverEzC(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;
307 } else {
308 AliWarning("Field calculation at z=0 (CE) is not allowed!");
309 saveEr[0]=0; saveEr[1]=0;
310 }
311 fLookUpErOverEz[i][j] = saveEr[0] + (saveEr[1]-saveEr[0])*(r-rList[ilow])/gridSizeR ;
312 }
313 }
314
315 /* delete [] saveEr;
316 delete [] sVec;
317 delete [] rList;
318 delete [] zedList;
319 */
320
321 fInitLookUp = kTRUE;
322
323}
324
325void AliTPCBoundaryVoltError::Print(const Option_t* option) const {
326 //
327 // Print function to check the settings of the boundary vectors
328 // option=="a" prints the C0 and C1 coefficents for calibration purposes
329 //
330
331 TString opt = option; opt.ToLower();
332 printf("%s\n",GetTitle());
333 printf(" - Voltage settings (on the TPC boundaries) - linearly interpolated\n");
334 printf(" : A-side (anti-clockwise)\n");
335 printf(" (0,1):\t IFC (CE) : %3.1f V \t IFC (ROC): %3.1f V \n",fBoundariesA[0],fBoundariesA[1]);
336 printf(" (2,3):\t ROC (IFC): %3.1f V \t ROC (OFC): %3.1f V \n",fBoundariesA[2],fBoundariesA[3]);
337 printf(" (4,5):\t OFC (ROC): %3.1f V \t OFC (CE) : %3.1f V \n",fBoundariesA[4],fBoundariesA[5]);
338 printf(" (6,7):\t CE (OFC): %3.1f V \t CE (IFC): %3.1f V \n",fBoundariesA[6],fBoundariesA[7]);
339 printf(" : C-side (clockwise)\n");
340 printf(" (0,1):\t IFC (CE) : %3.1f V \t IFC (ROC): %3.1f V \n",fBoundariesC[0],fBoundariesC[1]);
341 printf(" (2,3):\t ROC (IFC): %3.1f V \t ROC (OFC): %3.1f V \n",fBoundariesC[2],fBoundariesC[3]);
342 printf(" (4,5):\t OFC (ROC): %3.1f V \t OFC (CE) : %3.1f V \n",fBoundariesC[4],fBoundariesC[5]);
343 printf(" (6,7):\t CE (OFC): %3.1f V \t CE (IFC): %3.1f V \n",fBoundariesC[6],fBoundariesC[7]);
344
345 // Check wether the settings of the Central Electrode agree (on the A and C side)
346 // Note: they have to be antisymmetric!
347 if (( TMath::Abs(fBoundariesA[6]+fBoundariesC[6])>1e-5) || ( TMath::Abs(fBoundariesA[7]+fBoundariesC[7])>1e-5 ) ){
348 AliWarning("Boundary parameters for the Central Electrode (CE) are not anti-symmetric! HOW DID YOU MANAGE THAT?");
349 AliWarning("Congratulations, you just splitted the Central Electrode of the TPC!");
350 AliWarning("Non-physical settings of the boundary parameter at the Central Electrode");
351 }
352
353 if (opt.Contains("a")) { // Print all details
354 printf(" - T1: %1.4f, T2: %1.4f \n",fT1,fT2);
355 printf(" - C1: %1.4f, C0: %1.4f \n",fC1,fC0);
356 }
357
358 if (!fInitLookUp) AliError("Lookup table was not initialized! You should do InitBoundaryVoltErrorDistortion() ...");
359
360}
361
362
363void AliTPCBoundaryVoltError::SetBoundariesA(Float_t boundariesA[8]){
364 //
365 // set voltage errors on the TPC boundaries - A side
366 //
367 // Start at IFC at the Central electrode and work anti-clockwise (clockwise for C side) through
368 // IFC, ROC, OFC, and CE. The boundary conditions are currently defined to be a linear
369 // interpolation between pairs of numbers in the Boundary (e.g. fBoundariesA) vector.
370 // The first pair of numbers represent the beginning and end of the Inner Field cage, etc.
371 // The unit of the error potential vector is [Volt], whereas 1mm shift of the IFC would
372 // correspond to ~ 40 V
373 //
374 // Note: The setting for the CE will be passed to the C side!
375
376 for (Int_t i=0; i<8; i++) {
377 fBoundariesA[i]= boundariesA[i];
378 if (i>5) fBoundariesC[i]= -boundariesA[i]; // setting for the CE is passed to C side
379 }
380
381}
382void AliTPCBoundaryVoltError::SetBoundariesC(Float_t boundariesC[6]){
383 //
384 // set voltage errors on the TPC boundaries - A side
385 //
386 // Start at IFC at the Central electrode and work clockwise (for C side) through
387 // IFC, ROC and OFC. The boundary conditions are currently defined to be a linear
388 // interpolation between pairs of numbers in the Boundary (e.g. fBoundariesC) vector.
389 // The first pair of numbers represent the beginning and end of the Inner Field cage, etc.
390 // The unit of the error potential vector is [Volt], whereas 1mm shift of the IFC would
391 // correspond to ~ 40 V
392 //
393 // Note: The setting for the CE will be taken from the A side (pos 6 and 7)!
394
395 for (Int_t i=0; i<6; i++) {
396 fBoundariesC[i]= boundariesC[i];
397 }
398
399}