]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCBoundaryVoltError.cxx
A bug leading to the break segmentation fixed by Jens Wiechula.
[u/mrichter/AliRoot.git] / TPC / AliTPCBoundaryVoltError.cxx
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 capable of calculating 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
58 ClassImp(AliTPCBoundaryVoltError)
59
60 AliTPCBoundaryVoltError::AliTPCBoundaryVoltError()
61   : AliTPCCorrection("BoundaryVoltError","Boundary Voltage Error"),
62     fC0(0.),fC1(0.),
63     fROCdisplacement(kTRUE),
64     fInitLookUp(kFALSE)
65 {
66   //
67   // default constructor
68   //
69   for (Int_t i=0; i<8; i++){
70     fBoundariesA[i]= 0;  
71     if (i<6) fBoundariesC[i]= 0;
72   }
73 }
74
75 AliTPCBoundaryVoltError::~AliTPCBoundaryVoltError() {
76   //
77   // default destructor
78   //
79 }
80
81
82
83 void AliTPCBoundaryVoltError::Init() {
84   //
85   // Initialization funtion
86   //
87   
88   AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
89   if (!magF) AliError("Magneticd field - not initialized");
90   Double_t bzField = magF->SolenoidField()/10.; //field in T
91   AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
92   if (!param) AliError("Parameters - not initialized");
93   Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us]   // From dataBase: to be updated: per second (ideally)
94   Double_t ezField = 400; // [V/cm]   // to be updated: never (hopefully)
95   Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ; 
96   // Correction Terms for effective omegaTau; obtained by a laser calibration run
97   SetOmegaTauT1T2(wt,fT1,fT2);
98
99   InitBoundaryVoltErrorDistortion();
100 }
101
102 void AliTPCBoundaryVoltError::Update(const TTimeStamp &/*timeStamp*/) {
103   //
104   // Update function 
105   //
106   AliMagF* magF= (AliMagF*)TGeoGlobalMagField::Instance()->GetField();
107   if (!magF) AliError("Magneticd field - not initialized");
108   Double_t bzField = magF->SolenoidField()/10.; //field in T
109   AliTPCParam *param= AliTPCcalibDB::Instance()->GetParameters();
110   if (!param) AliError("Parameters - not initialized");
111   Double_t vdrift = param->GetDriftV()/1000000.; // [cm/us]  // From dataBase: to be updated: per second (ideally)
112   Double_t ezField = 400; // [V/cm]   // to be updated: never (hopefully)
113   Double_t wt = -10.0 * (bzField*10) * vdrift / ezField ; 
114   // Correction Terms for effective omegaTau; obtained by a laser calibration run
115   SetOmegaTauT1T2(wt,fT1,fT2);
116
117 }
118
119
120
121 void 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) {
127     AliInfo("Lookup table was not initialized!  Perform the inizialisation now ...");
128     InitBoundaryVoltErrorDistortion();
129   }
130
131   Int_t   order     = 1 ;               // FIXME: hardcoded? Linear interpolation = 1, Quadratic = 2         
132                                         // note that the poisson solution was linearly mirroed on this grid!
133   Double_t intEr, intEphi, intdEz ;
134   Double_t r, phi, z ;
135   Int_t    sign;
136
137   r      =  TMath::Sqrt( x[0]*x[0] + x[1]*x[1] ) ;
138   phi    =  TMath::ATan2(x[1],x[0]) ;
139   if ( phi < 0 ) phi += TMath::TwoPi() ;                   // Table uses phi from 0 to 2*Pi
140   z      =  x[2] ;                                         // Create temporary copy of x[2]
141
142   if ( (roc%36) < 18 ) {
143     sign =  1;       // (TPC A side)
144   } else {
145     sign = -1;       // (TPC C side)
146   }
147   
148   if ( sign==1  && z <  fgkZOffSet ) z =  fgkZOffSet;    // Protect against discontinuity at CE
149   if ( sign==-1 && z > -fgkZOffSet ) z = -fgkZOffSet;    // Protect against discontinuity at CE
150   
151
152   intEphi = 0.0;  // Efield is symmetric in phi - 2D calculation
153
154   if ( (sign==1 && z<0) || (sign==-1 && z>0) ) // just a consistency check
155     AliError("ROC number does not correspond to z coordinate! Calculation of distortions is most likely wrong!");
156
157   // Get the E field integral
158   Interpolate2DEdistortion( order, r, z, fLookUpErOverEz, intEr );
159   // Get DeltaEz field integral
160   Interpolate2DEdistortion( order, r, z, fLookUpDeltaEz, intdEz );
161   
162   // Calculate distorted position
163   if ( r > 0.0 ) {
164     phi =  phi + ( fC0*intEphi - fC1*intEr ) / r;      
165     r   =  r   + ( fC0*intEr   + fC1*intEphi );  
166   }
167   
168   // Calculate correction in cartesian coordinates
169   dx[0] = r * TMath::Cos(phi) - x[0];
170   dx[1] = r * TMath::Sin(phi) - x[1]; 
171   dx[2] = intdEz;  // z distortion - (internally scaled with driftvelocity dependency 
172                    // on the Ez field plus the actual ROC misalignment (if set TRUE)
173
174
175 }
176
177 void AliTPCBoundaryVoltError::InitBoundaryVoltErrorDistortion() {
178   //
179   // Initialization of the Lookup table which contains the solutions of the 
180   // Dirichlet boundary problem
181   //
182
183   const Float_t  gridSizeR   =  (fgkOFCRadius-fgkIFCRadius) / (kRows-1) ;
184   const Float_t  gridSizeZ   =  fgkTPCZ0 / (kColumns-1) ;
185
186   TMatrixD voltArrayA(kRows,kColumns), voltArrayC(kRows,kColumns); // boundary vectors
187   TMatrixD chargeDensity(kRows,kColumns);                              // dummy charge
188   TMatrixD arrayErOverEzA(kRows,kColumns), arrayErOverEzC(kRows,kColumns); // solution
189   TMatrixD arrayDeltaEzA(kRows,kColumns),  arrayDeltaEzC(kRows,kColumns); // solution
190
191   Double_t  rList[kRows], zedList[kColumns] ;
192   
193   // Fill arrays with initial conditions.  V on the boundary and ChargeDensity in the volume.      
194   for ( Int_t j = 0 ; j < kColumns ; j++ ) {
195     Double_t zed = j*gridSizeZ ;
196     zedList[j] = zed ;
197     for ( Int_t i = 0 ; i < kRows ; i++ )  {
198       Double_t radius = fgkIFCRadius + i*gridSizeR ;
199       rList[i]           = radius ;
200       voltArrayA(i,j)        = 0;  // Initialize voltArrayA to zero
201       voltArrayC(i,j)        = 0;  // Initialize voltArrayC to zero
202       chargeDensity(i,j)     = 0;  // Initialize ChargeDensity to zero - not used in this class
203     }
204   }      
205
206
207   // check if boundary values are the same for the C side (for later, saving some calculation time)
208
209   Int_t symmetry = -1; // assume that  A and C side are identical (but anti-symmetric!) // e.g conical deformation
210   Int_t sVec[8];
211
212   // check if boundaries are different (regardless the sign)
213   for (Int_t i=0; i<8; i++) { 
214     if (TMath::Abs(TMath::Abs(fBoundariesA[i]) - TMath::Abs(fBoundariesC[i])) > 1e-5) 
215       symmetry = 0;  
216     sVec[i] = (Int_t)( TMath::Sign((Float_t)1.,fBoundariesA[i]) * TMath::Sign((Float_t)1.,fBoundariesC[i])); 
217   }
218   if (symmetry==-1) { // still the same values?
219     // check the kind of symmetry , if even ...
220     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 ) 
221       symmetry =  1;
222     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 ) 
223       symmetry = -1;
224     else
225       symmetry =  0; // some of the values differ in the sign -> neither symmetric nor antisymmetric
226   }
227
228
229
230   // Solve the electrosatic problem in 2D 
231
232   // Fill the complete Boundary vectors
233   // Start at IFC at CE and work anti-clockwise through IFC, ROC, OFC, and CE (clockwise for C side)
234   for ( Int_t j = 0 ; j < kColumns ; j++ ) {
235     Double_t zed = j*gridSizeZ ;
236     for ( Int_t i = 0 ; i < kRows ; i++ ) { 
237       Double_t radius = fgkIFCRadius + i*gridSizeR ;
238
239       // A side boundary vectors
240       if ( i == 0 ) voltArrayA(i,j) += zed   *((fBoundariesA[1]-fBoundariesA[0])/((kColumns-1)*gridSizeZ))
241         + fBoundariesA[0] ; // IFC
242       if ( j == kColumns-1 ) voltArrayA(i,j) += (radius-fgkIFCRadius)*((fBoundariesA[3]-fBoundariesA[2])/((kRows-1)*gridSizeR))
243         + fBoundariesA[2] ; // ROC
244       if ( i == kRows-1 ) voltArrayA(i,j) += zed   *((fBoundariesA[4]-fBoundariesA[5])/((kColumns-1)*gridSizeZ))
245         + fBoundariesA[5] ; // OFC
246       if ( j == 0 ) voltArrayA(i,j) += (radius-fgkIFCRadius)*((fBoundariesA[6]-fBoundariesA[7])/((kRows-1)*gridSizeR))
247         + fBoundariesA[7] ; // CE
248       
249       if (symmetry==0) {
250         // C side boundary vectors
251         if ( i == 0 ) voltArrayC(i,j) += zed   *((fBoundariesC[1]-fBoundariesC[0])/((kColumns-1)*gridSizeZ))
252           + fBoundariesC[0] ; // IFC
253         if ( j == kColumns-1 ) voltArrayC(i,j) += (radius-fgkIFCRadius)*((fBoundariesC[3]-fBoundariesC[2])/((kRows-1)*gridSizeR))
254           + fBoundariesC[2] ; // ROC
255         if ( i == kRows-1 ) voltArrayC(i,j) += zed   *((fBoundariesC[4]-fBoundariesC[5])/((kColumns-1)*gridSizeZ))
256           + fBoundariesC[5] ; // OFC
257         if ( j == 0 ) voltArrayC(i,j) += (radius-fgkIFCRadius)*((fBoundariesC[6]-fBoundariesC[7])/((kRows-1)*gridSizeR))
258           + fBoundariesC[7] ; // CE
259       }
260
261     }
262   }
263
264   voltArrayA(0,0)               *= 0.5 ; // Use average boundary condition at corner
265   voltArrayA(kRows-1,0)         *= 0.5 ; // Use average boundary condition at corner
266   voltArrayA(0,kColumns-1)      *= 0.5 ; // Use average boundary condition at corner
267   voltArrayA(kRows-1,kColumns-1)*= 0.5 ; // Use average boundary condition at corner
268
269   if (symmetry==0) {
270     voltArrayC(0,0)               *= 0.5 ; // Use average boundary condition at corner
271     voltArrayC(kRows-1,0)         *= 0.5 ; // Use average boundary condition at corner
272     voltArrayC(0,kColumns-1)      *= 0.5 ; // Use average boundary condition at corner
273     voltArrayC(kRows-1,kColumns-1)*= 0.5 ; // Use average boundary condition at corner
274   }
275
276
277   // always solve the problem on the A side
278   PoissonRelaxation2D( voltArrayA, chargeDensity, arrayErOverEzA, arrayDeltaEzA, 
279                        kRows, kColumns, kIterations, fROCdisplacement ) ;
280
281   if (symmetry!=0) { // A and C side are the same ("anti-symmetric" or "symmetric")
282     for ( Int_t j = 0 ; j < kColumns ; j++ ) {
283       for ( Int_t i = 0 ; i < kRows ; i++ ) { 
284         arrayErOverEzC(i,j) = symmetry*arrayErOverEzA(i,j);
285         arrayDeltaEzC(i,j) = -symmetry*arrayDeltaEzA(i,j);
286       }
287     }
288   } else if (symmetry==0) { // A and C side are different - Solve the problem on the C side too
289     PoissonRelaxation2D( voltArrayC, chargeDensity, arrayErOverEzC, arrayDeltaEzC,
290                          kRows, kColumns, kIterations, fROCdisplacement ) ;
291     for ( Int_t j = 0 ; j < kColumns ; j++ ) {
292       for ( Int_t i = 0 ; i < kRows ; i++ ) { 
293         arrayDeltaEzC(i,j) = -arrayDeltaEzC(i,j); // negative z coordinate!
294       }
295     }
296   }
297
298   // Interpolate results onto standard grid for Electric Fields
299   Int_t ilow=0, jlow=0 ;
300   Double_t z,r;
301   Float_t saveEr[2] ;         
302   Float_t saveEz[2] ;         
303   for ( Int_t i = 0 ; i < kNZ ; ++i )  {
304     z = TMath::Abs( fgkZList[i] ) ;
305     for ( Int_t j = 0 ; j < kNR ; ++j ) {
306       // Linear interpolation !!
307       r = fgkRList[j] ;
308       Search( kRows,   rList, r, ilow ) ;          // Note switch - R in rows and Z in columns
309       Search( kColumns, zedList, z, jlow ) ;
310       if ( ilow < 0 ) ilow = 0 ;                   // check if out of range
311       if ( jlow < 0 ) jlow = 0 ;   
312       if ( ilow + 1  >=  kRows - 1 ) ilow =  kRows - 2 ;              
313       if ( jlow + 1  >=  kColumns - 1 ) jlow =  kColumns - 2 ; 
314       
315       if (fgkZList[i]>0) {         // A side solution
316         saveEr[0] = arrayErOverEzA(ilow,jlow) + 
317           (arrayErOverEzA(ilow,jlow+1)-arrayErOverEzA(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
318         saveEr[1] = arrayErOverEzA(ilow+1,jlow) + 
319           (arrayErOverEzA(ilow+1,jlow+1)-arrayErOverEzA(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;
320         saveEz[0] = arrayDeltaEzA(ilow,jlow) + 
321           (arrayDeltaEzA(ilow,jlow+1)-arrayDeltaEzA(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
322         saveEz[1] = arrayDeltaEzA(ilow+1,jlow) + 
323           (arrayDeltaEzA(ilow+1,jlow+1)-arrayDeltaEzA(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;
324
325       } else if (fgkZList[i]<0) {  // C side solution
326         saveEr[0] = arrayErOverEzC(ilow,jlow) + 
327           (arrayErOverEzC(ilow,jlow+1)-arrayErOverEzC(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
328         saveEr[1] = arrayErOverEzC(ilow+1,jlow) + 
329           (arrayErOverEzC(ilow+1,jlow+1)-arrayErOverEzC(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;
330         saveEz[0] = arrayDeltaEzC(ilow,jlow) + 
331           (arrayDeltaEzC(ilow,jlow+1)-arrayDeltaEzC(ilow,jlow))*(z-zedList[jlow])/gridSizeZ ;
332         saveEz[1] = arrayDeltaEzC(ilow+1,jlow) + 
333           (arrayDeltaEzC(ilow+1,jlow+1)-arrayDeltaEzC(ilow+1,jlow))*(z-zedList[jlow])/gridSizeZ ;
334
335       } else {
336         AliWarning("Field calculation at z=0 (CE) is not allowed!");
337         saveEr[0]=0; saveEr[1]=0;
338         saveEz[0]=0; saveEz[1]=0;
339       }
340       fLookUpErOverEz[i][j] = saveEr[0] + (saveEr[1]-saveEr[0])*(r-rList[ilow])/gridSizeR ;
341       fLookUpDeltaEz[i][j]  = saveEz[0] + (saveEz[1]-saveEz[0])*(r-rList[ilow])/gridSizeR ;
342     }
343   }
344   
345   voltArrayA.Clear();
346   voltArrayC.Clear();
347   chargeDensity.Clear();
348   arrayErOverEzA.Clear();
349   arrayErOverEzC.Clear();
350   arrayDeltaEzA.Clear();
351   arrayDeltaEzC.Clear();
352   
353   fInitLookUp = kTRUE;
354
355 }
356
357 void AliTPCBoundaryVoltError::Print(const Option_t* option) const {
358   //
359   // Print function to check the settings of the boundary vectors
360   // option=="a" prints the C0 and C1 coefficents for calibration purposes
361   //
362
363   TString opt = option; opt.ToLower();
364   printf("%s\n",GetTitle());
365   printf(" - Voltage settings (on the TPC boundaries) - linearly interpolated\n");
366   printf("  : A-side (anti-clockwise)\n"); 
367   printf("     (0,1):\t IFC (CE) : %3.1f V \t IFC (ROC): %3.1f V \n",fBoundariesA[0],fBoundariesA[1]);
368   printf("     (2,3):\t ROC (IFC): %3.1f V \t ROC (OFC): %3.1f V \n",fBoundariesA[2],fBoundariesA[3]);
369   printf("     (4,5):\t OFC (ROC): %3.1f V \t OFC (CE) : %3.1f V \n",fBoundariesA[4],fBoundariesA[5]);
370   printf("     (6,7):\t CE  (OFC): %3.1f V \t CE  (IFC): %3.1f V \n",fBoundariesA[6],fBoundariesA[7]);
371   printf("  : C-side (clockwise)\n"); 
372   printf("     (0,1):\t IFC (CE) : %3.1f V \t IFC (ROC): %3.1f V \n",fBoundariesC[0],fBoundariesC[1]);
373   printf("     (2,3):\t ROC (IFC): %3.1f V \t ROC (OFC): %3.1f V \n",fBoundariesC[2],fBoundariesC[3]);
374   printf("     (4,5):\t OFC (ROC): %3.1f V \t OFC (CE) : %3.1f V \n",fBoundariesC[4],fBoundariesC[5]);
375   printf("     (6,7):\t CE  (OFC): %3.1f V \t CE  (IFC): %3.1f V \n",fBoundariesC[6],fBoundariesC[7]);
376
377   // Check wether the settings of the Central Electrode agree (on the A and C side)
378   // Note: they have to be antisymmetric!
379   if (( TMath::Abs(fBoundariesA[6]+fBoundariesC[6])>1e-5) || ( TMath::Abs(fBoundariesA[7]+fBoundariesC[7])>1e-5 ) ){
380     AliWarning("Boundary parameters for the Central Electrode (CE) are not anti-symmetric! HOW DID YOU MANAGE THAT?");
381     AliWarning("Congratulations, you just splitted the Central Electrode of the TPC!");
382     AliWarning("Non-physical settings of the boundary parameter at the Central Electrode");
383   }
384
385   if (opt.Contains("a")) { // Print all details
386     printf(" - T1: %1.4f, T2: %1.4f \n",fT1,fT2);
387     printf(" - C1: %1.4f, C0: %1.4f \n",fC1,fC0);
388   } 
389    
390   if (!fInitLookUp) 
391     AliError("Lookup table was not initialized! You should do InitBoundaryVoltErrorDistortion() ...");
392
393 }
394
395
396 void AliTPCBoundaryVoltError::SetBoundariesA(Float_t boundariesA[8]){
397   //
398   // set voltage errors on the TPC boundaries - A side 
399   //
400   // Start at IFC at the Central electrode and work anti-clockwise (clockwise for C side) through 
401   // IFC, ROC, OFC, and CE. The boundary conditions are currently defined to be a linear 
402   // interpolation between pairs of numbers in the Boundary (e.g. fBoundariesA) vector.  
403   // The first pair of numbers represent the beginning and end of the Inner Field cage, etc.
404   // The unit of the error potential vector is [Volt], whereas 1mm shift of the IFC would 
405   // correspond to ~ 40 V
406   // 
407   // Note: The setting for the CE will be passed to the C side!
408   
409   for (Int_t i=0; i<8; i++) {
410     fBoundariesA[i]= boundariesA[i];  
411     if (i>5) fBoundariesC[i]= -boundariesA[i]; // setting for the CE is passed to C side
412   }
413   fInitLookUp=kFALSE;
414 }
415 void AliTPCBoundaryVoltError::SetBoundariesC(Float_t boundariesC[6]){
416   //
417   // set voltage errors on the TPC boundaries - C side 
418   //
419   // Start at IFC at the Central electrode and work clockwise (for C side) through 
420   // IFC, ROC and OFC. The boundary conditions are currently defined to be a linear 
421   // interpolation between pairs of numbers in the Boundary (e.g. fBoundariesC) vector.  
422   // The first pair of numbers represent the beginning and end of the Inner Field cage, etc.
423   // The unit of the error potential vector is [Volt], whereas 1mm shift of the IFC would 
424   // correspond to ~ 40 V
425   // 
426   // Note: The setting for the CE will be taken from the A side (pos 6 and 7)!
427
428   for (Int_t i=0; i<6; i++) {
429     fBoundariesC[i]= boundariesC[i];  
430   }
431   fInitLookUp=kFALSE;
432 }