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 /// \class AliTPCROCVoltError3D
17 /// \brief The class calculates the space point distortions due to z offsets of the
19 /// ROCs via the residual voltage technique (Poisson relaxation) in 3D.
20 /// Since the GG (part of the ROCs) represents the closure of the FC in z direction,
21 /// every misalignment in z produces not only dz distortions but also electrical
22 /// field inhomogeneities throughout the volume, which produces additional dr and rd$\phi$ distortions.
24 /// Each ROC can be misaligned (in z direction) in three ways. A general z0 offset,
25 /// an inclination along the x and an inclination along the y axis. The z-misalignment's
26 /// can be set via the function SetROCData(TMatrixD *mat) for each single chamber (ROC).
27 /// The array size has to be (72,3) representing the 72 chambers according to the
28 /// offline numbering scheme (IROC: roc$<$36; OROC: roc$\geq$36) and the three misalignment's
29 /// (see the source code for further details).
31 /// Internally, these z offsets (unit is cm) are recalculated into residual voltage
32 /// equivalents in order to make use of the relaxation technique.
34 /// One has two possibilities when calculating the $dz$ distortions. The resulting
35 /// distortions are purely due to the change of the drift velocity (along with the
36 /// change of the drift field) when the SetROCDisplacement is FALSE.
37 /// For this class, this is a rather unphysical setting and should be avoided.
38 /// When the flag is set to TRUE, the corresponding offset in z is added to the dz
39 /// calculation of the outer ROCs.
40 /// For the Alice TPC gas, both effects are of similar magnitude. This means, if the
41 /// drift length is sufficiently large, a z-offset of a chamber appears to have (approx.)
42 /// twice the magnitude when one looks at the actual dz distortions.
44 /// In addition, this class allows a correction regarding the different arrival times
45 /// of the electrons due to the geometrical difference of the inner and outer chambers.
46 /// The impact was simulated via Garfield. If the flag is set via the
47 /// function SetElectronArrivalCorrection, the electron-arrival correction is added to the dz calculation.
48 /// ![Picture from ROOT macro](AliTPCROCVoltError3D_cxx_99b6733.png)
50 /// \author Jim Thomas, Stefan Rossegger
55 #include "TGeoGlobalMagField.h"
56 #include "AliTPCcalibDB.h"
57 #include "AliTPCParam.h"
63 #include "AliTPCROC.h"
64 #include "AliTPCROCVoltError3D.h"
67 ClassImp(AliTPCROCVoltError3D)
70 AliTPCROCVoltError3D::AliTPCROCVoltError3D()
71 : AliTPCCorrection("ROCVoltErrors","ROC z alignment Errors"),
73 fROCdisplacement(kTRUE),
74 fElectronArrivalCorrection(kTRUE),
80 // default constructor
83 // Array which will contain the solution according to the setted boundary conditions
84 // main input: z alignment of the Read Out chambers
85 // see InitROCVoltError3D() function
86 for ( Int_t k = 0 ; k < kNPhi ; k++ ) {
87 fLookUpErOverEz[k] = new TMatrixF(kNR,kNZ);
88 fLookUpEphiOverEz[k] = new TMatrixF(kNR,kNZ);
89 fLookUpDeltaEz[k] = new TMatrixF(kNR,kNZ);
91 // fROCDataFileName="$ALICE_ROOT/TPC/Calib/maps/TPCROCdzSurvey.root";
92 // SetROCDataFileName(fROCDataFileName.Data()); // initialization of fdzDataLinFit is included
96 AliTPCROCVoltError3D::~AliTPCROCVoltError3D() {
99 for ( Int_t k = 0 ; k < kNPhi ; k++ ) {
100 delete fLookUpErOverEz[k];
101 delete fLookUpEphiOverEz[k];
102 delete fLookUpDeltaEz[k];
105 delete fdzDataLinFit;
110 Bool_t AliTPCROCVoltError3D::AddCorrectionCompact(AliTPCCorrection* corr, Double_t weight){
111 /// Add correction and make them compact
113 /// - origin of distortion/correction are additive
114 /// - only correction ot the same type supported ()
117 AliError("Zerro pointer - correction");
120 AliTPCROCVoltError3D * corrC = dynamic_cast<AliTPCROCVoltError3D *>(corr);
122 AliError(TString::Format("Inconsistent class types: %s\%s",IsA()->GetName(),corr->IsA()->GetName()).Data());
126 TMatrixD matrixDz=*(corrC->fdzDataLinFit);
128 if (fdzDataLinFit==NULL) {
129 fdzDataLinFit=new TMatrixD(matrixDz);
130 fROCdisplacement=corrC->fROCdisplacement;
131 fElectronArrivalCorrection=corrC->fElectronArrivalCorrection;
134 (*fdzDataLinFit)+=matrixDz;
143 void AliTPCROCVoltError3D::SetROCData(TMatrixD * matrix){
144 /// Set a z alignment map of the chambers not via a file, but
145 /// directly via a TMatrix(72,3), where dz = p0 + p1*(lx-133.4) + p2*ly (all in cm)
147 if (!fdzDataLinFit) fdzDataLinFit=new TMatrixD(*matrix);
148 else *fdzDataLinFit = *matrix;
152 void AliTPCROCVoltError3D::Init() {
153 /// Initialization funtion
155 AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
156 if (!magF) AliError("Magneticd field - not initialized");
157 Double_t bzField = magF->SolenoidField()/10.; //field in T
158 AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
159 if (!param) AliError("Parameters - not initialized");
160 Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us] // From dataBase: to be updated: per second (ideally)
161 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
162 Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ;
163 // Correction Terms for effective omegaTau; obtained by a laser calibration run
164 SetOmegaTauT1T2(wt,fT1,fT2);
166 if (!fInitLookUp) InitROCVoltError3D();
169 void AliTPCROCVoltError3D::Update(const TTimeStamp &/*timeStamp*/) {
172 AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
173 if (!magF) AliError("Magneticd field - not initialized");
174 Double_t bzField = magF->SolenoidField()/10.; //field in T
175 AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
176 if (!param) AliError("Parameters - not initialized");
177 Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us] // From dataBase: to be updated: per second (ideally)
178 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
179 Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ;
180 // Correction Terms for effective omegaTau; obtained by a laser calibration run
181 SetOmegaTauT1T2(wt,fT1,fT2);
185 void AliTPCROCVoltError3D::SetROCDataFileName(const char * fname) {
186 /// Set / load the ROC data (linear fit of ROC misalignments)
188 fROCDataFileName = fname;
190 TFile f(fROCDataFileName.Data(),"READ");
191 TMatrixD *m = (TMatrixD*) f.Get("dzSurveyLinFitData");
194 // prepare some space
196 if (fdzDataLinFit) delete fdzDataLinFit;
197 fdzDataLinFit = new TMatrixD(72,3);
198 TMatrixD &dataIntern = *fdzDataLinFit;
200 for (Int_t iroc=0;iroc<72;iroc++) {
201 dataIntern(iroc,0) = mf(iroc,0); // z0 offset
202 dataIntern(iroc,1) = mf(iroc,1); // slope in x
203 dataIntern(iroc,2) = mf(iroc,2); // slope in y
208 fInitLookUp = kFALSE;
212 void AliTPCROCVoltError3D::GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]) {
213 /// Calculates the correction due e.g. residual voltage errors on the TPC boundaries
215 const Double_t kEpsilon=Double_t(FLT_MIN);
217 AliInfo("Lookup table was not initialized! Perform the inizialisation now ...");
218 InitROCVoltError3D();
220 static Bool_t forceInit=kTRUE; //temporary needed for back compatibility with old OCDB entries
221 if (forceInit&&fLookUpErOverEz[0]){
222 if (TMath::Abs(fLookUpErOverEz[0]->Sum())<kEpsilon){//temporary needed for back compatibility with old OCDB entries
223 ForceInitROCVoltError3D();
229 Int_t order = 1 ; // FIXME: hardcoded? Linear interpolation = 1, Quadratic = 2
231 Float_t intEr, intEphi, intDeltaEz;
235 r = TMath::Sqrt( x[0]*x[0] + x[1]*x[1] ) ;
236 phi = TMath::ATan2(x[1],x[0]) ;
237 if ( phi < 0 ) phi += TMath::TwoPi() ; // Table uses phi from 0 to 2*Pi
238 z = x[2] ; // Create temporary copy of x[2]
240 if ( (roc%36) < 18 ) {
241 sign = 1; // (TPC A side)
243 sign = -1; // (TPC C side)
246 if ( sign==1 && z < fgkZOffSet ) z = fgkZOffSet; // Protect against discontinuity at CE
247 if ( sign==-1 && z > -fgkZOffSet ) z = -fgkZOffSet; // Protect against discontinuity at CE
250 if ( (sign==1 && z<0) || (sign==-1 && z>0) ) // just a consistency check
251 AliError("ROC number does not correspond to z coordinate! Calculation of distortions is most likely wrong!");
253 // Get the Er and Ephi field integrals plus the integral over DeltaEz
254 intEr = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi,
255 fgkRList, fgkZList, fgkPhiList, fLookUpErOverEz );
256 intEphi = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi,
257 fgkRList, fgkZList, fgkPhiList, fLookUpEphiOverEz);
258 intDeltaEz = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi,
259 fgkRList, fgkZList, fgkPhiList, fLookUpDeltaEz );
261 // printf("%lf %lf %lf\n",intEr,intEphi,intDeltaEz);
263 // Calculate distorted position
265 phi = phi + ( fC0*intEphi - fC1*intEr ) / r;
266 r = r + ( fC0*intEr + fC1*intEphi );
269 // Calculate correction in cartesian coordinates
270 dx[0] = r * TMath::Cos(phi) - x[0];
271 dx[1] = r * TMath::Sin(phi) - x[1];
272 dx[2] = intDeltaEz; // z distortion - (internally scaled with driftvelocity dependency
273 // on the Ez field plus the actual ROC misalignment (if set TRUE)
276 if (fElectronArrivalCorrection) {
278 // correction for the OROC (in average, a 0.014usec longer drift time
279 // due to different position of the anode wires) -> vd*dt -> 2.64*0.014 = 0.0369 cm
280 // FIXME: correction are token from Magboltz simulations
281 // should be loaded from a database
283 AliTPCROC * rocInfo = AliTPCROC::Instance();
284 Double_t rCrossingROC = (rocInfo->GetPadRowRadii(0,62)+rocInfo->GetPadRowRadii(36,0))/2;
286 if (r>rCrossingROC) {
288 dx[2] += 0.0369; // A side - negative correction
290 dx[2] -= 0.0369; // C side - positive correction
297 void AliTPCROCVoltError3D::InitROCVoltError3D() {
298 /// Initialization of the Lookup table which contains the solutions of the
299 /// Dirichlet boundary problem
300 /// Calculation of the single 3D-Poisson solver is done just if needed
301 /// (see basic lookup tables in header file)
303 const Int_t order = 1 ; // Linear interpolation = 1, Quadratic = 2
304 const Float_t gridSizeR = (fgkOFCRadius-fgkIFCRadius) / (kRows-1) ;
305 const Float_t gridSizeZ = fgkTPCZ0 / (kColumns-1) ;
306 const Float_t gridSizePhi = TMath::TwoPi() / ( 18.0 * kPhiSlicesPerSector);
308 // temporary arrays to create the boundary conditions
309 TMatrixD *arrayofArrayV[kPhiSlices], *arrayofCharge[kPhiSlices] ;
310 TMatrixD *arrayofEroverEz[kPhiSlices], *arrayofEphioverEz[kPhiSlices], *arrayofDeltaEz[kPhiSlices] ;
312 for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
313 arrayofArrayV[k] = new TMatrixD(kRows,kColumns) ;
314 arrayofCharge[k] = new TMatrixD(kRows,kColumns) ;
315 arrayofEroverEz[k] = new TMatrixD(kRows,kColumns) ;
316 arrayofEphioverEz[k] = new TMatrixD(kRows,kColumns) ;
317 arrayofDeltaEz[k] = new TMatrixD(kRows,kColumns) ;
320 // list of point as used in the poisson relation and the interpolation (during sum up)
321 Double_t rlist[kRows], zedlist[kColumns] , philist[kPhiSlices];
322 for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
323 philist[k] = gridSizePhi * k;
324 for ( Int_t i = 0 ; i < kRows ; i++ ) {
325 rlist[i] = fgkIFCRadius + i*gridSizeR ;
326 for ( Int_t j = 0 ; j < kColumns ; j++ ) { // Fill Vmatrix with Boundary Conditions
327 zedlist[j] = j * gridSizeZ ;
332 // ==========================================================================
333 // Solve Poisson's equation in 3D cylindrical coordinates by relaxation technique
334 // Allow for different size grid spacing in R and Z directions
336 const Int_t symmetry = 0;
338 // Set bondaries and solve Poisson's equation --------------------------
340 if ( !fInitLookUp ) {
342 AliInfo(Form("Solving the poisson equation (~ %d sec)",2*10*(int)(kPhiSlices/10)));
344 for ( Int_t side = 0 ; side < 2 ; side++ ) { // Solve Poisson3D twice; once for +Z and once for -Z
346 for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
347 TMatrixD &arrayV = *arrayofArrayV[k] ;
348 TMatrixD &charge = *arrayofCharge[k] ;
350 //Fill arrays with initial conditions. V on the boundary and Charge in the volume.
351 for ( Int_t i = 0 ; i < kRows ; i++ ) {
352 for ( Int_t j = 0 ; j < kColumns ; j++ ) { // Fill Vmatrix with Boundary Conditions
356 Float_t radius0 = rlist[i] ;
357 Float_t phi0 = gridSizePhi * k ;
359 // To avoid problems at sector boundaries, use an average of +- 1 degree from actual phi location
360 if ( j == (kColumns-1) ) {
361 arrayV(i,j) = 0.5* ( GetROCVoltOffset( side, radius0, phi0+0.02 ) + GetROCVoltOffset( side, radius0, phi0-0.02 ) ) ;
363 if (side==1) // C side
364 arrayV(i,j) = -arrayV(i,j); // minus sign on the C side to allow a consistent usage of global z when setting the boundaries
369 for ( Int_t i = 1 ; i < kRows-1 ; i++ ) {
370 for ( Int_t j = 1 ; j < kColumns-1 ; j++ ) {
376 // Solve Poisson's equation in 3D cylindrical coordinates by relaxation technique
377 // Allow for different size grid spacing in R and Z directions
379 PoissonRelaxation3D( arrayofArrayV, arrayofCharge,
380 arrayofEroverEz, arrayofEphioverEz, arrayofDeltaEz,
381 kRows, kColumns, kPhiSlices, gridSizePhi, kIterations,
382 symmetry, fROCdisplacement) ;
385 //Interpolate results onto a custom grid which is used just for these calculations.
387 for ( Int_t k = 0 ; k < kNPhi ; k++ ) {
388 phi = fgkPhiList[k] ;
390 TMatrixF &erOverEz = *fLookUpErOverEz[k] ;
391 TMatrixF &ephiOverEz = *fLookUpEphiOverEz[k];
392 TMatrixF &deltaEz = *fLookUpDeltaEz[k] ;
394 for ( Int_t j = 0 ; j < kNZ ; j++ ) {
396 z = TMath::Abs(fgkZList[j]) ; // Symmetric solution in Z that depends only on ABS(Z)
398 if ( side == 0 && fgkZList[j] < 0 ) continue; // Skip rest of this loop if on the wrong side
399 if ( side == 1 && fgkZList[j] > 0 ) continue; // Skip rest of this loop if on the wrong side
401 for ( Int_t i = 0 ; i < kNR ; i++ ) {
404 // Interpolate basicLookup tables; once for each rod, then sum the results
405 erOverEz(i,j) = Interpolate3DTable(order, r, z, phi, kRows, kColumns, kPhiSlices,
406 rlist, zedlist, philist, arrayofEroverEz );
407 ephiOverEz(i,j) = Interpolate3DTable(order, r, z, phi, kRows, kColumns, kPhiSlices,
408 rlist, zedlist, philist, arrayofEphioverEz);
409 deltaEz(i,j) = Interpolate3DTable(order, r, z, phi, kRows, kColumns, kPhiSlices,
410 rlist, zedlist, philist, arrayofDeltaEz );
412 if (side == 1) deltaEz(i,j) = - deltaEz(i,j); // negative coordinate system on C side
418 if ( side == 0 ) AliInfo(" A side done");
419 if ( side == 1 ) AliInfo(" C side done");
423 // clear the temporary arrays lists
424 for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
425 delete arrayofArrayV[k];
426 delete arrayofCharge[k];
427 delete arrayofEroverEz[k];
428 delete arrayofEphioverEz[k];
429 delete arrayofDeltaEz[k];
438 Float_t AliTPCROCVoltError3D::GetROCVoltOffset(Int_t side, Float_t r0, Float_t phi0) {
439 /// Returns the dz alignment data (in voltage equivalents) at
440 /// the given position
442 Float_t xp = r0*TMath::Cos(phi0);
443 Float_t yp = r0*TMath::Sin(phi0);
445 // phi0 should be between 0 and 2pi
446 if (phi0<0) phi0+=TMath::TwoPi();
447 Int_t roc = (Int_t)TMath::Floor((TMath::RadToDeg()*phi0)/20);
448 if (side==1) roc+=18; // C side
449 if (r0>132) roc+=36; // OROC
451 // linear-plane data: z = z0 + kx*lx + ky*ly (rotation in local coordinates)
452 TMatrixD &fitData = *fdzDataLinFit;
455 Double_t secAlpha = TMath::DegToRad()*(10.+20.*(((Int_t)roc)%18));
456 Float_t lx = xp*TMath::Cos(secAlpha)+yp*TMath::Sin(secAlpha);
457 Float_t ly = yp*TMath::Cos(secAlpha)-xp*TMath::Sin(secAlpha);
459 // reference of rotation in lx is at the intersection between OROC and IROC
460 // necessary, since the Fitprozedure is otherwise useless
462 AliTPCROC * rocInfo = AliTPCROC::Instance();
463 Double_t lxRef = (rocInfo->GetPadRowRadii(0,62)+rocInfo->GetPadRowRadii(36,0))/2;
465 Float_t dz = fitData(roc,0)+fitData(roc,1)*(lx-lxRef) + fitData(roc,2)*ly; // value in cm
467 // aproximated Voltage-offset-aquivalent to the z misalignment
468 // (linearly scaled with the z position)
469 Double_t ezField = (fgkCathodeV-fgkGG)/fgkTPCZ0; // = ALICE Electric Field (V/cm) Magnitude ~ -400 V/cm;
470 Float_t voltOff = dz*ezField; // values in "Volt equivalents"
475 TH2F * AliTPCROCVoltError3D::CreateHistoOfZAlignment(Int_t side, Int_t nx, Int_t ny) {
476 /// return a simple histogramm containing the input to the poisson solver
477 /// (z positions of the Read-out chambers, linearly interpolated)
480 if (side==0) snprintf(hname,100,"survey_dz_Aside");
481 if (side==1) snprintf(hname,100,"survey_dz_Cside");
483 TH2F *h = new TH2F(hname,hname,nx,-250.,250.,ny,-250.,250.);
485 for (Int_t iy=1;iy<=ny;++iy) {
486 Double_t yp = h->GetYaxis()->GetBinCenter(iy);
487 for (Int_t ix=1;ix<=nx;++ix) {
488 Double_t xp = h->GetXaxis()->GetBinCenter(ix);
490 Float_t r0 = TMath::Sqrt(xp*xp+yp*yp);
491 Float_t phi0 = TMath::ATan2(yp,xp);
493 Float_t dz = GetROCVoltOffset(side,r0,phi0); // in [volt]
495 Double_t ezField = (fgkCathodeV-fgkGG)/fgkTPCZ0; // = ALICE Electric Field (V/cm) Magnitude ~ -400 V/cm;
496 dz = dz/ezField; // in [cm]
498 if (85.<=r0 && r0<=245.) {
499 h->SetBinContent(ix,iy,dz);
501 h->SetBinContent(ix,iy,0.);
506 h->GetXaxis()->SetTitle("x [cm]");
507 h->GetYaxis()->SetTitle("y [cm]");
508 h->GetZaxis()->SetTitle("dz [cm]");
510 // h->DrawCopy("colz");
515 void AliTPCROCVoltError3D::Print(const Option_t* option) const {
516 /// Print function to check the settings of the Rod shifts and the rotated clips
517 /// option=="a" prints the C0 and C1 coefficents for calibration purposes
519 TString opt = option; opt.ToLower();
520 printf("%s\n",GetTitle());
521 printf(" - z aligmnet of the TPC Read-Out chambers \n");
522 printf(" (linear interpolation within the chamber: dz = z0 + kx*(lx-133) + ky*ly [cm] ) \n");
523 printf(" Info: Check the following data-file for more details: %s \n",fROCDataFileName.Data());
524 printf(" fROCdisplacement\n",fROCdisplacement);
525 printf(" fElectronArrivalCorrection\n",fElectronArrivalCorrection);
527 if (opt.Contains("a")) { // Print all details
528 TMatrixD &fitData = *fdzDataLinFit;
529 printf(" A side: IROC ROCX=(z0,kx,ky): \n");
530 for (Int_t roc = 0; roc<18; roc++)
531 printf("ROC%d:(%.2e,%.2e,%.2e) ",roc,fitData(roc,0),fitData(roc,1),fitData(roc,2));
532 printf("\n A side: OROC ROCX=(z0,kx,ky): \n");
533 for (Int_t roc = 36; roc<54; roc++)
534 printf("ROC%d:(%.2e,%.2e,%.2e) ",roc,fitData(roc,0),fitData(roc,1),fitData(roc,2));
535 printf("\n C side: IROC ROCX=(z0,kx,ky): \n");
536 for (Int_t roc = 18; roc<36; roc++)
537 printf("ROC%d:(%.2e,%.2e,%.2e) ",roc,fitData(roc,0),fitData(roc,1),fitData(roc,2));
538 printf("\n C side: OROC ROCX=(z0,kx,ky): \n");
539 for (Int_t roc = 54; roc<72; roc++)
540 printf("ROC%d:(%.2e,%.2e,%.2e) ",roc,fitData(roc,0),fitData(roc,1),fitData(roc,2));
542 printf(" - T1: %1.4f, T2: %1.4f \n",fT1,fT2);
543 printf(" - C1: %1.4f, C0: %1.4f \n",fC1,fC0);
546 if (!fInitLookUp) AliError("Lookup table was not initialized! You should do InitROCVoltError3D() ...");