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> AliTPCROCVoltError3D class </h2>
20 // The class calculates the space point distortions due to z offsets of the
21 // ROCs via the residual voltage technique (Poisson relaxation) in 3D.
22 // Since the GG (part of the ROCs) represents the closure of the FC in z direction,
23 // every misalignment in z produces not only dz distortions but also electrical
24 // field inhomogeneities throughout the volume, which produces additional dr and rd$\phi$ distortions.
26 // Each ROC can be misaligned (in z direction) in three ways. A general z0 offset,
27 // an inclination along the x and an inclination along the y axis. The z-misalignment's
28 // can be set via the function SetROCData(TMatrixD *mat) for each single chamber (ROC).
29 // The array size has to be (72,3) representing the 72 chambers according to the
30 // offline numbering scheme (IROC: roc$<$36; OROC: roc$\geq$36) and the three misalignment's
31 // (see the source code for further details).
33 // Internally, these z offsets (unit is cm) are recalculated into residual voltage
34 // equivalents in order to make use of the relaxation technique.
36 // One has two possibilities when calculating the $dz$ distortions. The resulting
37 // distortions are purely due to the change of the drift velocity (along with the
38 // change of the drift field) when the SetROCDisplacement is FALSE. <br>
39 // For this class, this is a rather unphysical setting and should be avoided. <br>
40 // When the flag is set to TRUE, the corresponding offset in z is added to the dz
41 // calculation of the outer ROCs. <br>
42 // For the Alice TPC gas, both effects are of similar magnitude. This means, if the
43 // drift length is sufficiently large, a z-offset of a chamber appears to have (approx.)
44 // twice the magnitude when one looks at the actual dz distortions.
46 // In addition, this class allows a correction regarding the different arrival times
47 // of the electrons due to the geometrical difference of the inner and outer chambers.
48 // The impact was simulated via Garfield. If the flag is set via the
49 // function SetElectronArrivalCorrection, the electron-arrival correction is added to the dz calculation.
52 // Begin_Macro(source)
54 // gROOT->SetStyle("Plain"); gStyle->SetPalette(1);
55 // TCanvas *c2 = new TCanvas("cAliTPCROCVoltError3D","cAliTPCROCVoltError3D",500,400);
56 // AliTPCROCVoltError3D roc;
57 // roc.SetROCDataFileName("$ALICE_ROOT/TPC/Calib/maps/TPCROCdzSurvey.root");
58 // roc.SetElectronArrivalCorrection(kFALSE); // Correction for electron arrival offset, IROC vs OROC
59 // roc.SetROCDisplacement(kTRUE); // include the chamber offset in z when calculating the dz
60 // roc.SetOmegaTauT1T2(0,1,1); // B=0
61 // roc.CreateHistoDZinXY(1.,300,300)->Draw("colz");
68 // Date: 08/08/2010 <br>
69 // Authors: Jim Thomas, Stefan Rossegger
71 // _________________________________________________________________
75 #include "TGeoGlobalMagField.h"
76 #include "AliTPCcalibDB.h"
77 #include "AliTPCParam.h"
83 #include "AliTPCROC.h"
84 #include "AliTPCROCVoltError3D.h"
86 ClassImp(AliTPCROCVoltError3D)
88 AliTPCROCVoltError3D::AliTPCROCVoltError3D()
89 : AliTPCCorrection("ROCVoltErrors","ROC z alignment Errors"),
91 fROCdisplacement(kTRUE),
92 fElectronArrivalCorrection(kTRUE),
98 // default constructor
101 // Array which will contain the solution according to the setted boundary conditions
102 // main input: z alignment of the Read Out chambers
103 // see InitROCVoltError3D() function
104 for ( Int_t k = 0 ; k < kNPhi ; k++ ) {
105 fLookUpErOverEz[k] = new TMatrixF(kNR,kNZ);
106 fLookUpEphiOverEz[k] = new TMatrixF(kNR,kNZ);
107 fLookUpDeltaEz[k] = new TMatrixF(kNR,kNZ);
109 // fROCDataFileName="$ALICE_ROOT/TPC/Calib/maps/TPCROCdzSurvey.root";
110 // SetROCDataFileName(fROCDataFileName.Data()); // initialization of fdzDataLinFit is included
114 AliTPCROCVoltError3D::~AliTPCROCVoltError3D() {
119 for ( Int_t k = 0 ; k < kNPhi ; k++ ) {
120 delete fLookUpErOverEz[k];
121 delete fLookUpEphiOverEz[k];
122 delete fLookUpDeltaEz[k];
125 delete fdzDataLinFit;
130 Bool_t AliTPCROCVoltError3D::AddCorrectionCompact(AliTPCCorrection* corr, Double_t weight){
132 // Add correction and make them compact
134 // - origin of distortion/correction are additive
135 // - only correction ot the same type supported ()
137 AliError("Zerro pointer - correction");
140 AliTPCROCVoltError3D * corrC = dynamic_cast<AliTPCROCVoltError3D *>(corr);
142 AliError(TString::Format("Inconsistent class types: %s\%s",IsA()->GetName(),corr->IsA()->GetName()).Data());
146 TMatrixD matrixDz=*(corrC->fdzDataLinFit);
148 if (fdzDataLinFit==NULL) {
149 fdzDataLinFit=new TMatrixD(matrixDz);
150 fROCdisplacement=corrC->fROCdisplacement;
151 fElectronArrivalCorrection=corrC->fElectronArrivalCorrection;
154 (*fdzDataLinFit)+=matrixDz;
163 void AliTPCROCVoltError3D::SetROCData(TMatrixD * matrix){
165 // Set a z alignment map of the chambers not via a file, but
166 // directly via a TMatrix(72,3), where dz = p0 + p1*(lx-133.4) + p2*ly (all in cm)
168 if (!fdzDataLinFit) fdzDataLinFit=new TMatrixD(*matrix);
169 else *fdzDataLinFit = *matrix;
173 void AliTPCROCVoltError3D::Init() {
175 // Initialization funtion
178 AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
179 if (!magF) AliError("Magneticd field - not initialized");
180 Double_t bzField = magF->SolenoidField()/10.; //field in T
181 AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
182 if (!param) AliError("Parameters - not initialized");
183 Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us] // From dataBase: to be updated: per second (ideally)
184 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
185 Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ;
186 // Correction Terms for effective omegaTau; obtained by a laser calibration run
187 SetOmegaTauT1T2(wt,fT1,fT2);
189 if (!fInitLookUp) InitROCVoltError3D();
192 void AliTPCROCVoltError3D::Update(const TTimeStamp &/*timeStamp*/) {
196 AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
197 if (!magF) AliError("Magneticd field - not initialized");
198 Double_t bzField = magF->SolenoidField()/10.; //field in T
199 AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
200 if (!param) AliError("Parameters - not initialized");
201 Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us] // From dataBase: to be updated: per second (ideally)
202 Double_t ezField = 400; // [V/cm] // to be updated: never (hopefully)
203 Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ;
204 // Correction Terms for effective omegaTau; obtained by a laser calibration run
205 SetOmegaTauT1T2(wt,fT1,fT2);
209 void AliTPCROCVoltError3D::SetROCDataFileName(const char * fname) {
211 // Set / load the ROC data (linear fit of ROC misalignments)
214 fROCDataFileName = fname;
216 TFile f(fROCDataFileName.Data(),"READ");
217 TMatrixD *m = (TMatrixD*) f.Get("dzSurveyLinFitData");
220 // prepare some space
222 if (fdzDataLinFit) delete fdzDataLinFit;
223 fdzDataLinFit = new TMatrixD(72,3);
224 TMatrixD &dataIntern = *fdzDataLinFit;
226 for (Int_t iroc=0;iroc<72;iroc++) {
227 dataIntern(iroc,0) = mf(iroc,0); // z0 offset
228 dataIntern(iroc,1) = mf(iroc,1); // slope in x
229 dataIntern(iroc,2) = mf(iroc,2); // slope in y
234 fInitLookUp = kFALSE;
238 void AliTPCROCVoltError3D::GetCorrection(const Float_t x[],const Short_t roc,Float_t dx[]) {
240 // Calculates the correction due e.g. residual voltage errors on the TPC boundaries
242 const Double_t kEpsilon=Double_t(FLT_MIN);
244 AliInfo("Lookup table was not initialized! Perform the inizialisation now ...");
245 InitROCVoltError3D();
247 static Bool_t forceInit=kTRUE; //temporary needed for back compatibility with old OCDB entries
248 if (forceInit&&fLookUpErOverEz[0]){
249 if (TMath::Abs(fLookUpErOverEz[0]->Sum())<kEpsilon){//temporary needed for back compatibility with old OCDB entries
250 ForceInitROCVoltError3D();
256 Int_t order = 1 ; // FIXME: hardcoded? Linear interpolation = 1, Quadratic = 2
258 Float_t intEr, intEphi, intDeltaEz;
262 r = TMath::Sqrt( x[0]*x[0] + x[1]*x[1] ) ;
263 phi = TMath::ATan2(x[1],x[0]) ;
264 if ( phi < 0 ) phi += TMath::TwoPi() ; // Table uses phi from 0 to 2*Pi
265 z = x[2] ; // Create temporary copy of x[2]
267 if ( (roc%36) < 18 ) {
268 sign = 1; // (TPC A side)
270 sign = -1; // (TPC C side)
273 if ( sign==1 && z < fgkZOffSet ) z = fgkZOffSet; // Protect against discontinuity at CE
274 if ( sign==-1 && z > -fgkZOffSet ) z = -fgkZOffSet; // Protect against discontinuity at CE
277 if ( (sign==1 && z<0) || (sign==-1 && z>0) ) // just a consistency check
278 AliError("ROC number does not correspond to z coordinate! Calculation of distortions is most likely wrong!");
280 // Get the Er and Ephi field integrals plus the integral over DeltaEz
281 intEr = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi,
282 fgkRList, fgkZList, fgkPhiList, fLookUpErOverEz );
283 intEphi = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi,
284 fgkRList, fgkZList, fgkPhiList, fLookUpEphiOverEz);
285 intDeltaEz = Interpolate3DTable(order, r, z, phi, kNR, kNZ, kNPhi,
286 fgkRList, fgkZList, fgkPhiList, fLookUpDeltaEz );
288 // printf("%lf %lf %lf\n",intEr,intEphi,intDeltaEz);
290 // Calculate distorted position
292 phi = phi + ( fC0*intEphi - fC1*intEr ) / r;
293 r = r + ( fC0*intEr + fC1*intEphi );
296 // Calculate correction in cartesian coordinates
297 dx[0] = r * TMath::Cos(phi) - x[0];
298 dx[1] = r * TMath::Sin(phi) - x[1];
299 dx[2] = intDeltaEz; // z distortion - (internally scaled with driftvelocity dependency
300 // on the Ez field plus the actual ROC misalignment (if set TRUE)
303 if (fElectronArrivalCorrection) {
305 // correction for the OROC (in average, a 0.014usec longer drift time
306 // due to different position of the anode wires) -> vd*dt -> 2.64*0.014 = 0.0369 cm
307 // FIXME: correction are token from Magboltz simulations
308 // should be loaded from a database
310 AliTPCROC * rocInfo = AliTPCROC::Instance();
311 Double_t rCrossingROC = (rocInfo->GetPadRowRadii(0,62)+rocInfo->GetPadRowRadii(36,0))/2;
313 if (r>rCrossingROC) {
315 dx[2] += 0.0369; // A side - negative correction
317 dx[2] -= 0.0369; // C side - positive correction
324 void AliTPCROCVoltError3D::InitROCVoltError3D() {
326 // Initialization of the Lookup table which contains the solutions of the
327 // Dirichlet boundary problem
328 // Calculation of the single 3D-Poisson solver is done just if needed
329 // (see basic lookup tables in header file)
332 const Int_t order = 1 ; // Linear interpolation = 1, Quadratic = 2
333 const Float_t gridSizeR = (fgkOFCRadius-fgkIFCRadius) / (kRows-1) ;
334 const Float_t gridSizeZ = fgkTPCZ0 / (kColumns-1) ;
335 const Float_t gridSizePhi = TMath::TwoPi() / ( 18.0 * kPhiSlicesPerSector);
337 // temporary arrays to create the boundary conditions
338 TMatrixD *arrayofArrayV[kPhiSlices], *arrayofCharge[kPhiSlices] ;
339 TMatrixD *arrayofEroverEz[kPhiSlices], *arrayofEphioverEz[kPhiSlices], *arrayofDeltaEz[kPhiSlices] ;
341 for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
342 arrayofArrayV[k] = new TMatrixD(kRows,kColumns) ;
343 arrayofCharge[k] = new TMatrixD(kRows,kColumns) ;
344 arrayofEroverEz[k] = new TMatrixD(kRows,kColumns) ;
345 arrayofEphioverEz[k] = new TMatrixD(kRows,kColumns) ;
346 arrayofDeltaEz[k] = new TMatrixD(kRows,kColumns) ;
349 // list of point as used in the poisson relation and the interpolation (during sum up)
350 Double_t rlist[kRows], zedlist[kColumns] , philist[kPhiSlices];
351 for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
352 philist[k] = gridSizePhi * k;
353 for ( Int_t i = 0 ; i < kRows ; i++ ) {
354 rlist[i] = fgkIFCRadius + i*gridSizeR ;
355 for ( Int_t j = 0 ; j < kColumns ; j++ ) { // Fill Vmatrix with Boundary Conditions
356 zedlist[j] = j * gridSizeZ ;
361 // ==========================================================================
362 // Solve Poisson's equation in 3D cylindrical coordinates by relaxation technique
363 // Allow for different size grid spacing in R and Z directions
365 const Int_t symmetry = 0;
367 // Set bondaries and solve Poisson's equation --------------------------
369 if ( !fInitLookUp ) {
371 AliInfo(Form("Solving the poisson equation (~ %d sec)",2*10*(int)(kPhiSlices/10)));
373 for ( Int_t side = 0 ; side < 2 ; side++ ) { // Solve Poisson3D twice; once for +Z and once for -Z
375 for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
376 TMatrixD &arrayV = *arrayofArrayV[k] ;
377 TMatrixD &charge = *arrayofCharge[k] ;
379 //Fill arrays with initial conditions. V on the boundary and Charge in the volume.
380 for ( Int_t i = 0 ; i < kRows ; i++ ) {
381 for ( Int_t j = 0 ; j < kColumns ; j++ ) { // Fill Vmatrix with Boundary Conditions
385 Float_t radius0 = rlist[i] ;
386 Float_t phi0 = gridSizePhi * k ;
388 // To avoid problems at sector boundaries, use an average of +- 1 degree from actual phi location
389 if ( j == (kColumns-1) ) {
390 arrayV(i,j) = 0.5* ( GetROCVoltOffset( side, radius0, phi0+0.02 ) + GetROCVoltOffset( side, radius0, phi0-0.02 ) ) ;
392 if (side==1) // C side
393 arrayV(i,j) = -arrayV(i,j); // minus sign on the C side to allow a consistent usage of global z when setting the boundaries
398 for ( Int_t i = 1 ; i < kRows-1 ; i++ ) {
399 for ( Int_t j = 1 ; j < kColumns-1 ; j++ ) {
405 // Solve Poisson's equation in 3D cylindrical coordinates by relaxation technique
406 // Allow for different size grid spacing in R and Z directions
408 PoissonRelaxation3D( arrayofArrayV, arrayofCharge,
409 arrayofEroverEz, arrayofEphioverEz, arrayofDeltaEz,
410 kRows, kColumns, kPhiSlices, gridSizePhi, kIterations,
411 symmetry, fROCdisplacement) ;
414 //Interpolate results onto a custom grid which is used just for these calculations.
416 for ( Int_t k = 0 ; k < kNPhi ; k++ ) {
417 phi = fgkPhiList[k] ;
419 TMatrixF &erOverEz = *fLookUpErOverEz[k] ;
420 TMatrixF &ephiOverEz = *fLookUpEphiOverEz[k];
421 TMatrixF &deltaEz = *fLookUpDeltaEz[k] ;
423 for ( Int_t j = 0 ; j < kNZ ; j++ ) {
425 z = TMath::Abs(fgkZList[j]) ; // Symmetric solution in Z that depends only on ABS(Z)
427 if ( side == 0 && fgkZList[j] < 0 ) continue; // Skip rest of this loop if on the wrong side
428 if ( side == 1 && fgkZList[j] > 0 ) continue; // Skip rest of this loop if on the wrong side
430 for ( Int_t i = 0 ; i < kNR ; i++ ) {
433 // Interpolate basicLookup tables; once for each rod, then sum the results
434 erOverEz(i,j) = Interpolate3DTable(order, r, z, phi, kRows, kColumns, kPhiSlices,
435 rlist, zedlist, philist, arrayofEroverEz );
436 ephiOverEz(i,j) = Interpolate3DTable(order, r, z, phi, kRows, kColumns, kPhiSlices,
437 rlist, zedlist, philist, arrayofEphioverEz);
438 deltaEz(i,j) = Interpolate3DTable(order, r, z, phi, kRows, kColumns, kPhiSlices,
439 rlist, zedlist, philist, arrayofDeltaEz );
441 if (side == 1) deltaEz(i,j) = - deltaEz(i,j); // negative coordinate system on C side
447 if ( side == 0 ) AliInfo(" A side done");
448 if ( side == 1 ) AliInfo(" C side done");
452 // clear the temporary arrays lists
453 for ( Int_t k = 0 ; k < kPhiSlices ; k++ ) {
454 delete arrayofArrayV[k];
455 delete arrayofCharge[k];
456 delete arrayofEroverEz[k];
457 delete arrayofEphioverEz[k];
458 delete arrayofDeltaEz[k];
467 Float_t AliTPCROCVoltError3D::GetROCVoltOffset(Int_t side, Float_t r0, Float_t phi0) {
469 // Returns the dz alignment data (in voltage equivalents) at
470 // the given position
473 Float_t xp = r0*TMath::Cos(phi0);
474 Float_t yp = r0*TMath::Sin(phi0);
476 // phi0 should be between 0 and 2pi
477 if (phi0<0) phi0+=TMath::TwoPi();
478 Int_t roc = (Int_t)TMath::Floor((TMath::RadToDeg()*phi0)/20);
479 if (side==1) roc+=18; // C side
480 if (r0>132) roc+=36; // OROC
482 // linear-plane data: z = z0 + kx*lx + ky*ly (rotation in local coordinates)
483 TMatrixD &fitData = *fdzDataLinFit;
486 Double_t secAlpha = TMath::DegToRad()*(10.+20.*(((Int_t)roc)%18));
487 Float_t lx = xp*TMath::Cos(secAlpha)+yp*TMath::Sin(secAlpha);
488 Float_t ly = yp*TMath::Cos(secAlpha)-xp*TMath::Sin(secAlpha);
490 // reference of rotation in lx is at the intersection between OROC and IROC
491 // necessary, since the Fitprozedure is otherwise useless
493 AliTPCROC * rocInfo = AliTPCROC::Instance();
494 Double_t lxRef = (rocInfo->GetPadRowRadii(0,62)+rocInfo->GetPadRowRadii(36,0))/2;
496 Float_t dz = fitData(roc,0)+fitData(roc,1)*(lx-lxRef) + fitData(roc,2)*ly; // value in cm
498 // aproximated Voltage-offset-aquivalent to the z misalignment
499 // (linearly scaled with the z position)
500 Double_t ezField = (fgkCathodeV-fgkGG)/fgkTPCZ0; // = ALICE Electric Field (V/cm) Magnitude ~ -400 V/cm;
501 Float_t voltOff = dz*ezField; // values in "Volt equivalents"
506 TH2F * AliTPCROCVoltError3D::CreateHistoOfZAlignment(Int_t side, Int_t nx, Int_t ny) {
508 // return a simple histogramm containing the input to the poisson solver
509 // (z positions of the Read-out chambers, linearly interpolated)
512 if (side==0) snprintf(hname,100,"survey_dz_Aside");
513 if (side==1) snprintf(hname,100,"survey_dz_Cside");
515 TH2F *h = new TH2F(hname,hname,nx,-250.,250.,ny,-250.,250.);
517 for (Int_t iy=1;iy<=ny;++iy) {
518 Double_t yp = h->GetYaxis()->GetBinCenter(iy);
519 for (Int_t ix=1;ix<=nx;++ix) {
520 Double_t xp = h->GetXaxis()->GetBinCenter(ix);
522 Float_t r0 = TMath::Sqrt(xp*xp+yp*yp);
523 Float_t phi0 = TMath::ATan2(yp,xp);
525 Float_t dz = GetROCVoltOffset(side,r0,phi0); // in [volt]
527 Double_t ezField = (fgkCathodeV-fgkGG)/fgkTPCZ0; // = ALICE Electric Field (V/cm) Magnitude ~ -400 V/cm;
528 dz = dz/ezField; // in [cm]
530 if (85.<=r0 && r0<=245.) {
531 h->SetBinContent(ix,iy,dz);
533 h->SetBinContent(ix,iy,0.);
538 h->GetXaxis()->SetTitle("x [cm]");
539 h->GetYaxis()->SetTitle("y [cm]");
540 h->GetZaxis()->SetTitle("dz [cm]");
542 // h->DrawCopy("colz");
547 void AliTPCROCVoltError3D::Print(const Option_t* option) const {
549 // Print function to check the settings of the Rod shifts and the rotated clips
550 // option=="a" prints the C0 and C1 coefficents for calibration purposes
553 TString opt = option; opt.ToLower();
554 printf("%s\n",GetTitle());
555 printf(" - z aligmnet of the TPC Read-Out chambers \n");
556 printf(" (linear interpolation within the chamber: dz = z0 + kx*(lx-133) + ky*ly [cm] ) \n");
557 printf(" Info: Check the following data-file for more details: %s \n",fROCDataFileName.Data());
558 printf(" fROCdisplacement\n",fROCdisplacement);
559 printf(" fElectronArrivalCorrection\n",fElectronArrivalCorrection);
561 if (opt.Contains("a")) { // Print all details
562 TMatrixD &fitData = *fdzDataLinFit;
563 printf(" A side: IROC ROCX=(z0,kx,ky): \n");
564 for (Int_t roc = 0; roc<18; roc++)
565 printf("ROC%d:(%.2e,%.2e,%.2e) ",roc,fitData(roc,0),fitData(roc,1),fitData(roc,2));
566 printf("\n A side: OROC ROCX=(z0,kx,ky): \n");
567 for (Int_t roc = 36; roc<54; roc++)
568 printf("ROC%d:(%.2e,%.2e,%.2e) ",roc,fitData(roc,0),fitData(roc,1),fitData(roc,2));
569 printf("\n C side: IROC ROCX=(z0,kx,ky): \n");
570 for (Int_t roc = 18; roc<36; roc++)
571 printf("ROC%d:(%.2e,%.2e,%.2e) ",roc,fitData(roc,0),fitData(roc,1),fitData(roc,2));
572 printf("\n C side: OROC ROCX=(z0,kx,ky): \n");
573 for (Int_t roc = 54; roc<72; roc++)
574 printf("ROC%d:(%.2e,%.2e,%.2e) ",roc,fitData(roc,0),fitData(roc,1),fitData(roc,2));
576 printf(" - T1: %1.4f, T2: %1.4f \n",fT1,fT2);
577 printf(" - C1: %1.4f, C0: %1.4f \n",fC1,fC0);
580 if (!fInitLookUp) AliError("Lookup table was not initialized! You should do InitROCVoltError3D() ...");