1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 // _________________________________________________________________
19 // <h2> AliTPCBoundaryVoltError class </h2>
20 // This class calculates the space point distortions due to residual voltage errors on
21 // the main boundaries of the TPC. For example, the inner vessel of the TPC is shifted
22 // by a certain amount, whereas the ROCs on the A side and the C side follow this mechanical
23 // shift (at the inner vessel) in z direction. This example can be named "conical deformation"
24 // of the TPC field cage (see example below).
26 // The boundary conditions can be set via two arrays (A and C side) which contain the
27 // residual voltage setting modeling a possible shift or an inhomogeneity on the TPC field
28 // cage. In order to avoid that the user splits the Central Electrode (CE), the settings for
29 // the C side is taken from the array on the A side (points: A.6 and A.7). The region betweem
30 // the points is interpolated linearly onto the boundaries.
32 // The class uses the PoissonRelaxation2D (see AliTPCCorrection) to calculate the resulting
33 // electrical field inhomogeneities in the (r,z)-plane. Then, the Langevin-integral formalism
34 // is used to calculate the space point distortions. <br>
35 // Note: This class assumes a homogeneous magnetic field.
37 // One has two possibilities when calculating the $dz$ distortions. The resulting distortions
38 // are purely due to the change of the drift velocity (along with the change of the drift field)
39 // when the SetROCDisplacement is FALSE. This emulates for example a Gating-Grid Voltage offset
40 // without moving the ROCs. When the flag is set to TRUE, the ROCs are assumed to be misaligned
41 // and the corresponding offset in z is added.
44 // Begin_Macro(source)
46 // gROOT->SetStyle("Plain"); gStyle->SetPalette(1);
47 // TCanvas *c2 = new TCanvas("cAliTPCBoundaryVoltError","cAliTPCBoundaryVoltError",500,300);
48 // AliTPCBoundaryVoltError bve;
49 // Float_t val = 40;// [Volt]; 40V corresponds to 1mm
50 // /* IFC shift, CE follows, ROC follows by factor half */
51 // Float_t boundA[8] = { val, val, val,0,0,0,0,val}; // voltages A-side
52 // Float_t boundC[6] = {-val,-val,-val,0,0,0}; // voltages C-side
53 // bve.SetBoundariesA(boundA);
54 // bve.SetBoundariesC(boundC);
55 // bve.SetOmegaTauT1T2(-0.32,1,1);
56 // bve.SetROCDisplacement(kTRUE); // include the chamber offset in z when calculating the dz distortions
57 // bve.CreateHistoDRinZR(0)->Draw("surf2");
64 // Date: 01/06/2010 <br>
65 // Authors: Jim Thomas, Stefan Rossegger
67 // _________________________________________________________________
71 #include "TGeoGlobalMagField.h"
72 #include "AliTPCcalibDB.h"
73 #include "AliTPCParam.h"
78 #include "AliTPCROC.h"
79 #include "AliTPCBoundaryVoltError.h"
81 ClassImp(AliTPCBoundaryVoltError)
83 AliTPCBoundaryVoltError::AliTPCBoundaryVoltError()
84 : AliTPCCorrection("BoundaryVoltError","Boundary Voltage Error"),
86 fROCdisplacement(kTRUE),
90 // default constructor
92 for (Int_t i=0; i<8; i++){
94 if (i<6) fBoundariesC[i]= 0;
98 AliTPCBoundaryVoltError::~AliTPCBoundaryVoltError() {
100 // default destructor
107 Bool_t AliTPCBoundaryVoltError::AddCorrectionCompact(AliTPCCorrection* corr, Double_t weight){
109 // Add correction and make them compact
111 // - origin of distortion/correction are additive
112 // - only correction ot the same type supported ()
114 AliError("Zerro pointer - correction");
117 AliTPCBoundaryVoltError* corrC = dynamic_cast<AliTPCBoundaryVoltError *>(corr);
119 AliError(TString::Format("Inconsistent class types: %s\%s",IsA()->GetName(),corr->IsA()->GetName()).Data());
122 if (fROCdisplacement!=corrC->fROCdisplacement){
123 AliError(TString::Format("Inconsistent fROCdisplacement : %s\%s",IsA()->GetName(),corr->IsA()->GetName()).Data());
126 for (Int_t i=0;i <8; i++){
127 fBoundariesA[i]+= corrC->fBoundariesA[i]*weight;
128 fBoundariesC[i]+= corrC->fBoundariesC[i]*weight;
137 void AliTPCBoundaryVoltError::Init() {
139 // Initialization funtion
142 AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
143 if (!magF) AliError("Magneticd field - not initialized");
144 Double_t bzField = magF->SolenoidField()/10.; //field in T
145 AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
146 if (!param) AliError("Parameters - not initialized");
147 Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us] // From dataBase: to be updated: per second (ideally)
148 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
149 Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ;
150 // Correction Terms for effective omegaTau; obtained by a laser calibration run
151 SetOmegaTauT1T2(wt,fT1,fT2);
153 InitBoundaryVoltErrorDistortion();
156 void AliTPCBoundaryVoltError::Update(const TTimeStamp &/*timeStamp*/) {
160 AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
161 if (!magF) AliError("Magneticd field - not initialized");
162 Double_t bzField = magF->SolenoidField()/10.; //field in T
163 AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
164 if (!param) AliError("Parameters - not initialized");
165 Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us] // From dataBase: to be updated: per second (ideally)
166 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
167 Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ;
168 // Correction Terms for effective omegaTau; obtained by a laser calibration run
169 SetOmegaTauT1T2(wt,fT1,fT2);
175 void AliTPCBoundaryVoltError::GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]) {
177 // Calculates the correction due e.g. residual voltage errors on the TPC boundaries
181 AliInfo("Lookup table was not initialized! Perform the inizialisation now ...");
182 InitBoundaryVoltErrorDistortion();
185 Int_t order = 1 ; // FIXME: hardcoded? Linear interpolation = 1, Quadratic = 2
186 // note that the poisson solution was linearly mirroed on this grid!
187 Double_t intEr, intEphi, intdEz ;
191 r = TMath::Sqrt( x[0]*x[0] + x[1]*x[1] ) ;
192 phi = TMath::ATan2(x[1],x[0]) ;
193 if ( phi < 0 ) phi += TMath::TwoPi() ; // Table uses phi from 0 to 2*Pi
194 z = x[2] ; // Create temporary copy of x[2]
196 if ( (roc%36) < 18 ) {
197 sign = 1; // (TPC A side)
199 sign = -1; // (TPC C side)
202 if ( sign==1 && z < fgkZOffSet ) z = fgkZOffSet; // Protect against discontinuity at CE
203 if ( sign==-1 && z > -fgkZOffSet ) z = -fgkZOffSet; // Protect against discontinuity at CE
206 intEphi = 0.0; // Efield is symmetric in phi - 2D calculation
208 if ( (sign==1 && z<0) || (sign==-1 && z>0) ) // just a consistency check
209 AliError("ROC number does not correspond to z coordinate! Calculation of distortions is most likely wrong!");
211 // Get the E field integral
212 Interpolate2DEdistortion( order, r, z, fLookUpErOverEz, intEr );
213 // Get DeltaEz field integral
214 Interpolate2DEdistortion( order, r, z, fLookUpDeltaEz, intdEz );
216 // Calculate distorted position
218 phi = phi + ( fC0*intEphi - fC1*intEr ) / r;
219 r = r + ( fC0*intEr + fC1*intEphi );
222 // Calculate correction in cartesian coordinates
223 dx[0] = r * TMath::Cos(phi) - x[0];
224 dx[1] = r * TMath::Sin(phi) - x[1];
225 dx[2] = intdEz; // z distortion - (internally scaled with driftvelocity dependency
226 // on the Ez field plus the actual ROC misalignment (if set TRUE)
231 void AliTPCBoundaryVoltError::InitBoundaryVoltErrorDistortion() {
233 // Initialization of the Lookup table which contains the solutions of the
234 // Dirichlet boundary problem
237 const Float_t gridSizeR = (fgkOFCRadius-fgkIFCRadius) / (kRows-1) ;
238 const Float_t gridSizeZ = fgkTPCZ0 / (kColumns-1) ;
240 TMatrixD voltArrayA(kRows,kColumns), voltArrayC(kRows,kColumns); // boundary vectors
241 TMatrixD chargeDensity(kRows,kColumns); // dummy charge
242 TMatrixD arrayErOverEzA(kRows,kColumns), arrayErOverEzC(kRows,kColumns); // solution
243 TMatrixD arrayDeltaEzA(kRows,kColumns), arrayDeltaEzC(kRows,kColumns); // solution
245 Double_t rList[kRows], zedList[kColumns] ;
247 // Fill arrays with initial conditions. V on the boundary and ChargeDensity in the volume.
248 for ( Int_t j = 0 ; j < kColumns ; j++ ) {
249 Double_t zed = j*gridSizeZ ;
251 for ( Int_t i = 0 ; i < kRows ; i++ ) {
252 Double_t radius = fgkIFCRadius + i*gridSizeR ;
254 voltArrayA(i,j) = 0; // Initialize voltArrayA to zero
255 voltArrayC(i,j) = 0; // Initialize voltArrayC to zero
256 chargeDensity(i,j) = 0; // Initialize ChargeDensity to zero - not used in this class
261 // check if boundary values are the same for the C side (for later, saving some calculation time)
263 Int_t symmetry = -1; // assume that A and C side are identical (but anti-symmetric!) // e.g conical deformation
266 // check if boundaries are different (regardless the sign)
267 for (Int_t i=0; i<8; i++) {
268 if (TMath::Abs(TMath::Abs(fBoundariesA[i]) - TMath::Abs(fBoundariesC[i])) > 1e-5)
270 sVec[i] = (Int_t)( TMath::Sign((Float_t)1.,fBoundariesA[i]) * TMath::Sign((Float_t)1.,fBoundariesC[i]));
272 if (symmetry==-1) { // still the same values?
273 // check the kind of symmetry , if even ...
274 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 )
276 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 )
279 symmetry = 0; // some of the values differ in the sign -> neither symmetric nor antisymmetric
284 // Solve the electrosatic problem in 2D
286 // Fill the complete Boundary vectors
287 // Start at IFC at CE and work anti-clockwise through IFC, ROC, OFC, and CE (clockwise for C side)
288 for ( Int_t j = 0 ; j < kColumns ; j++ ) {
289 Double_t zed = j*gridSizeZ ;
290 for ( Int_t i = 0 ; i < kRows ; i++ ) {
291 Double_t radius = fgkIFCRadius + i*gridSizeR ;
293 // A side boundary vectors
294 if ( i == 0 ) voltArrayA(i,j) += zed *((fBoundariesA[1]-fBoundariesA[0])/((kColumns-1)*gridSizeZ))
295 + fBoundariesA[0] ; // IFC
296 if ( j == kColumns-1 ) voltArrayA(i,j) += (radius-fgkIFCRadius)*((fBoundariesA[3]-fBoundariesA[2])/((kRows-1)*gridSizeR))
297 + fBoundariesA[2] ; // ROC
298 if ( i == kRows-1 ) voltArrayA(i,j) += zed *((fBoundariesA[4]-fBoundariesA[5])/((kColumns-1)*gridSizeZ))
299 + fBoundariesA[5] ; // OFC
300 if ( j == 0 ) voltArrayA(i,j) += (radius-fgkIFCRadius)*((fBoundariesA[6]-fBoundariesA[7])/((kRows-1)*gridSizeR))
301 + fBoundariesA[7] ; // CE
304 // C side boundary vectors
305 if ( i == 0 ) voltArrayC(i,j) += zed *((fBoundariesC[1]-fBoundariesC[0])/((kColumns-1)*gridSizeZ))
306 + fBoundariesC[0] ; // IFC
307 if ( j == kColumns-1 ) voltArrayC(i,j) += (radius-fgkIFCRadius)*((fBoundariesC[3]-fBoundariesC[2])/((kRows-1)*gridSizeR))
308 + fBoundariesC[2] ; // ROC
309 if ( i == kRows-1 ) voltArrayC(i,j) += zed *((fBoundariesC[4]-fBoundariesC[5])/((kColumns-1)*gridSizeZ))
310 + fBoundariesC[5] ; // OFC
311 if ( j == 0 ) voltArrayC(i,j) += (radius-fgkIFCRadius)*((fBoundariesC[6]-fBoundariesC[7])/((kRows-1)*gridSizeR))
312 + fBoundariesC[7] ; // CE
318 voltArrayA(0,0) *= 0.5 ; // Use average boundary condition at corner
319 voltArrayA(kRows-1,0) *= 0.5 ; // Use average boundary condition at corner
320 voltArrayA(0,kColumns-1) *= 0.5 ; // Use average boundary condition at corner
321 voltArrayA(kRows-1,kColumns-1)*= 0.5 ; // Use average boundary condition at corner
324 voltArrayC(0,0) *= 0.5 ; // Use average boundary condition at corner
325 voltArrayC(kRows-1,0) *= 0.5 ; // Use average boundary condition at corner
326 voltArrayC(0,kColumns-1) *= 0.5 ; // Use average boundary condition at corner
327 voltArrayC(kRows-1,kColumns-1)*= 0.5 ; // Use average boundary condition at corner
331 // always solve the problem on the A side
332 PoissonRelaxation2D( voltArrayA, chargeDensity, arrayErOverEzA, arrayDeltaEzA,
333 kRows, kColumns, kIterations, fROCdisplacement ) ;
335 if (symmetry!=0) { // A and C side are the same ("anti-symmetric" or "symmetric")
336 for ( Int_t j = 0 ; j < kColumns ; j++ ) {
337 for ( Int_t i = 0 ; i < kRows ; i++ ) {
338 arrayErOverEzC(i,j) = symmetry*arrayErOverEzA(i,j);
339 arrayDeltaEzC(i,j) = -symmetry*arrayDeltaEzA(i,j);
342 } else if (symmetry==0) { // A and C side are different - Solve the problem on the C side too
343 PoissonRelaxation2D( voltArrayC, chargeDensity, arrayErOverEzC, arrayDeltaEzC,
344 kRows, kColumns, kIterations, fROCdisplacement ) ;
345 for ( Int_t j = 0 ; j < kColumns ; j++ ) {
346 for ( Int_t i = 0 ; i < kRows ; i++ ) {
347 arrayDeltaEzC(i,j) = -arrayDeltaEzC(i,j); // negative z coordinate!
352 // Interpolate results onto standard grid for Electric Fields
353 Int_t ilow=0, jlow=0 ;
357 for ( Int_t i = 0 ; i < kNZ ; ++i ) {
358 z = TMath::Abs( fgkZList[i] ) ;
359 for ( Int_t j = 0 ; j < kNR ; ++j ) {
360 // Linear interpolation !!
362 Search( kRows, rList, r, ilow ) ; // Note switch - R in rows and Z in columns
363 Search( kColumns, zedList, z, jlow ) ;
364 if ( ilow < 0 ) ilow = 0 ; // check if out of range
365 if ( jlow < 0 ) jlow = 0 ;
366 if ( ilow + 1 >= kRows - 1 ) ilow = kRows - 2 ;
367 if ( jlow + 1 >= kColumns - 1 ) jlow = kColumns - 2 ;
369 if (fgkZList[i]>0) { // A side solution
370 saveEr[0] = arrayErOverEzA(ilow,jlow) +
371 (arrayErOverEzA(ilow,jlow+1)-arrayErOverEzA(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
372 saveEr[1] = arrayErOverEzA(ilow+1,jlow) +
373 (arrayErOverEzA(ilow+1,jlow+1)-arrayErOverEzA(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;
374 saveEz[0] = arrayDeltaEzA(ilow,jlow) +
375 (arrayDeltaEzA(ilow,jlow+1)-arrayDeltaEzA(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
376 saveEz[1] = arrayDeltaEzA(ilow+1,jlow) +
377 (arrayDeltaEzA(ilow+1,jlow+1)-arrayDeltaEzA(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;
379 } else if (fgkZList[i]<0) { // C side solution
380 saveEr[0] = arrayErOverEzC(ilow,jlow) +
381 (arrayErOverEzC(ilow,jlow+1)-arrayErOverEzC(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
382 saveEr[1] = arrayErOverEzC(ilow+1,jlow) +
383 (arrayErOverEzC(ilow+1,jlow+1)-arrayErOverEzC(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;
384 saveEz[0] = arrayDeltaEzC(ilow,jlow) +
385 (arrayDeltaEzC(ilow,jlow+1)-arrayDeltaEzC(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
386 saveEz[1] = arrayDeltaEzC(ilow+1,jlow) +
387 (arrayDeltaEzC(ilow+1,jlow+1)-arrayDeltaEzC(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;
390 AliWarning("Field calculation at z=0 (CE) is not allowed!");
391 saveEr[0]=0; saveEr[1]=0;
392 saveEz[0]=0; saveEz[1]=0;
394 fLookUpErOverEz[i][j] = saveEr[0] + (saveEr[1]-saveEr[0])*(r-rList[ilow])/gridSizeR ;
395 fLookUpDeltaEz[i][j] = saveEz[0] + (saveEz[1]-saveEz[0])*(r-rList[ilow])/gridSizeR ;
401 chargeDensity.Clear();
402 arrayErOverEzA.Clear();
403 arrayErOverEzC.Clear();
404 arrayDeltaEzA.Clear();
405 arrayDeltaEzC.Clear();
411 void AliTPCBoundaryVoltError::Print(const Option_t* option) const {
413 // Print function to check the settings of the boundary vectors
414 // option=="a" prints the C0 and C1 coefficents for calibration purposes
417 TString opt = option; opt.ToLower();
418 printf("%s\n",GetTitle());
419 printf(" - Voltage settings (on the TPC boundaries) - linearly interpolated\n");
420 printf(" : A-side (anti-clockwise)\n");
421 printf(" (0,1):\t IFC (CE) : %3.1f V \t IFC (ROC): %3.1f V \n",fBoundariesA[0],fBoundariesA[1]);
422 printf(" (2,3):\t ROC (IFC): %3.1f V \t ROC (OFC): %3.1f V \n",fBoundariesA[2],fBoundariesA[3]);
423 printf(" (4,5):\t OFC (ROC): %3.1f V \t OFC (CE) : %3.1f V \n",fBoundariesA[4],fBoundariesA[5]);
424 printf(" (6,7):\t CE (OFC): %3.1f V \t CE (IFC): %3.1f V \n",fBoundariesA[6],fBoundariesA[7]);
425 printf(" : C-side (clockwise)\n");
426 printf(" (0,1):\t IFC (CE) : %3.1f V \t IFC (ROC): %3.1f V \n",fBoundariesC[0],fBoundariesC[1]);
427 printf(" (2,3):\t ROC (IFC): %3.1f V \t ROC (OFC): %3.1f V \n",fBoundariesC[2],fBoundariesC[3]);
428 printf(" (4,5):\t OFC (ROC): %3.1f V \t OFC (CE) : %3.1f V \n",fBoundariesC[4],fBoundariesC[5]);
429 printf(" (6,7):\t CE (OFC): %3.1f V \t CE (IFC): %3.1f V \n",fBoundariesC[6],fBoundariesC[7]);
431 // Check wether the settings of the Central Electrode agree (on the A and C side)
432 // Note: they have to be antisymmetric!
433 if (( TMath::Abs(fBoundariesA[6]+fBoundariesC[6])>1e-5) || ( TMath::Abs(fBoundariesA[7]+fBoundariesC[7])>1e-5 ) ){
434 AliWarning("Boundary parameters for the Central Electrode (CE) are not anti-symmetric! HOW DID YOU MANAGE THAT?");
435 AliWarning("Congratulations, you just splitted the Central Electrode of the TPC!");
436 AliWarning("Non-physical settings of the boundary parameter at the Central Electrode");
439 if (opt.Contains("a")) { // Print all details
440 printf(" - T1: %1.4f, T2: %1.4f \n",fT1,fT2);
441 printf(" - C1: %1.4f, C0: %1.4f \n",fC1,fC0);
445 AliError("Lookup table was not initialized! You should do InitBoundaryVoltErrorDistortion() ...");
450 void AliTPCBoundaryVoltError::SetBoundariesA(Float_t boundariesA[8]){
452 // set voltage errors on the TPC boundaries - A side
454 // Start at IFC at the Central electrode and work anti-clockwise (clockwise for C side) through
455 // IFC, ROC, OFC, and CE. The boundary conditions are currently defined to be a linear
456 // interpolation between pairs of numbers in the Boundary (e.g. fBoundariesA) vector.
457 // The first pair of numbers represent the beginning and end of the Inner Field cage, etc.
458 // The unit of the error potential vector is [Volt], whereas 1mm shift of the IFC would
459 // correspond to ~ 40 V
461 // Note: The setting for the CE will be passed to the C side!
463 for (Int_t i=0; i<8; i++) {
464 fBoundariesA[i]= boundariesA[i];
465 if (i>5) fBoundariesC[i]= -boundariesA[i]; // setting for the CE is passed to C side
469 void AliTPCBoundaryVoltError::SetBoundariesC(Float_t boundariesC[6]){
471 // set voltage errors on the TPC boundaries - C side
473 // Start at IFC at the Central electrode and work clockwise (for C side) through
474 // IFC, ROC and OFC. The boundary conditions are currently defined to be a linear
475 // interpolation between pairs of numbers in the Boundary (e.g. fBoundariesC) vector.
476 // The first pair of numbers represent the beginning and end of the Inner Field cage, etc.
477 // The unit of the error potential vector is [Volt], whereas 1mm shift of the IFC would
478 // correspond to ~ 40 V
480 // Note: The setting for the CE will be taken from the A side (pos 6 and 7)!
482 for (Int_t i=0; i<6; i++) {
483 fBoundariesC[i]= boundariesC[i];