]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliCheb3D.cxx
be8981b914317efd8f3dd0a4cdb00ec664fe5c11
[u/mrichter/AliRoot.git] / STEER / AliCheb3D.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 #include <TString.h>
17 #include <TSystem.h>
18 #include <TROOT.h>
19 #include <TRandom.h>
20 #include <stdio.h>
21 #include <TMethodCall.h>
22 #include <TMath.h>
23 #include <TH1.h>
24 #include "AliCheb3D.h"
25 #include "AliCheb3DCalc.h"
26
27 ClassImp(AliCheb3D)
28
29 //__________________________________________________________________________________________
30 AliCheb3D::AliCheb3D() : 
31   fDimOut(0), 
32   fPrec(0), 
33   fChebCalc(1), 
34   fMaxCoefs(0), 
35   fResTmp(0), 
36   fGrid(0), 
37   fUsrFunName(""), 
38   fUsrMacro(0) 
39 {
40   for (int i=3;i--;) {
41     fBMin[i] = fBMax[i] = fBScale[i] = fBOffset[i] = 0;
42     fNPoints[i] = 0;
43     fGridOffs[i] = 0.;
44   }
45 }
46
47 //__________________________________________________________________________________________
48 AliCheb3D::AliCheb3D(const AliCheb3D& src) : 
49   TNamed(src),
50   fDimOut(src.fDimOut), 
51   fPrec(src.fPrec), 
52   fChebCalc(1), 
53   fMaxCoefs(src.fMaxCoefs), 
54   fResTmp(0),
55   fGrid(0), 
56   fUsrFunName(src.fUsrFunName), 
57   fUsrMacro(0)
58 {
59   // read coefs from text file
60   for (int i=3;i--;) {
61     fBMin[i]    = src.fBMin[i];
62     fBMax[i]    = src.fBMax[i];
63     fBScale[i]  = src.fBScale[i];
64     fBOffset[i] = src.fBOffset[i];
65     fNPoints[i] = src.fNPoints[i];
66     fGridOffs[i] = src.fGridOffs[i];
67   }
68   for (int i=0;i<fDimOut;i++) {
69     AliCheb3DCalc* cbc = src.GetChebCalc(i);
70     if (cbc) fChebCalc.AddAtAndExpand(new AliCheb3DCalc(*cbc),i);
71   }
72 }
73
74 //__________________________________________________________________________________________
75 AliCheb3D::AliCheb3D(const char* inpFile) : 
76   fDimOut(0), 
77   fPrec(0),  
78   fChebCalc(1),
79   fMaxCoefs(0),  
80   fResTmp(0),
81   fGrid(0), 
82   fUsrFunName(""), 
83   fUsrMacro(0)
84 {
85   // read coefs from text file
86   for (int i=3;i--;) {
87     fBMin[i] = fBMax[i] = fBScale[i] = fBOffset[i] = 0;
88     fNPoints[i] = 0;
89     fGridOffs[i] = 0.;
90   }
91   LoadData(inpFile);
92 }
93
94 //__________________________________________________________________________________________
95 AliCheb3D::AliCheb3D(FILE* stream) : 
96   fDimOut(0), 
97   fPrec(0), 
98   fChebCalc(1), 
99   fMaxCoefs(0),
100   fResTmp(0),
101   fGrid(0),
102   fUsrFunName(""),
103   fUsrMacro(0)
104 {
105   // read coefs from stream
106   for (int i=3;i--;) {
107     fBMin[i] = fBMax[i] = fBScale[i] = fBOffset[i] = 0;
108     fNPoints[i] = 0;
109     fGridOffs[i] = 0.;
110   }
111   LoadData(stream);
112 }
113
114 //__________________________________________________________________________________________
115 #ifdef _INC_CREATION_ALICHEB3D_
116 AliCheb3D::AliCheb3D(const char* funName, int DimOut, Float_t  *bmin,Float_t  *bmax, Int_t *npoints, Float_t prec) : 
117   TNamed(funName,funName), 
118   fDimOut(0), 
119   fPrec(TMath::Max(1.E-12f,prec)), 
120   fChebCalc(1), 
121   fMaxCoefs(0), 
122   fResTmp(0), 
123   fGrid(0), 
124   fUsrFunName("") ,
125   fUsrMacro(0)
126 {
127   // Construct the parameterization for the function
128   // funName : name of the file containing the function: void funName(Float_t * inp,Float_t * out)
129   // DimOut  : dimension of the vector computed by the user function
130   // bmin    : array of 3 elements with the lower boundaries of the region where the function is defined
131   // bmax    : array of 3 elements with the upper boundaries of the region where the function is defined
132   // npoints : array of 3 elements with the number of points to compute in each of 3 dimension
133   // prec    : max allowed absolute difference between the user function and computed parameterization on the requested grid
134   //
135   if (DimOut<1) {Error("AliCheb3D","Requested output dimension is %d\nStop\n",fDimOut); exit(1);}
136   for (int i=3;i--;) {
137     fBMin[i] = fBMax[i] = fBScale[i] = fBOffset[i] = 0;
138     fNPoints[i] = 0;
139     fGridOffs[i] = 0.;
140   }
141   SetDimOut(DimOut);
142   PrepareBoundaries(bmin,bmax);
143   DefineGrid(npoints);
144   SetUsrFunction(funName);
145   ChebFit();
146   //
147 }
148 #endif
149
150 //__________________________________________________________________________________________
151 #ifdef _INC_CREATION_ALICHEB3D_
152 AliCheb3D::AliCheb3D(void (*ptr)(float*,float*), int DimOut, Float_t  *bmin,Float_t  *bmax, Int_t *npoints, Float_t prec) : 
153   fDimOut(0), 
154   fPrec(TMath::Max(1.E-12f,prec)), 
155   fChebCalc(1), 
156   fMaxCoefs(0), 
157   fResTmp(0), 
158   fGrid(0), 
159   fUsrFunName(""),
160   fUsrMacro(0)
161 {
162   // Construct the parameterization for the function
163   // ptr     : pointer on the function: void fun(Float_t * inp,Float_t * out)
164   // DimOut  : dimension of the vector computed by the user function
165   // bmin    : array of 3 elements with the lower boundaries of the region where the function is defined
166   // bmax    : array of 3 elements with the upper boundaries of the region where the function is defined
167   // npoints : array of 3 elements with the number of points to compute in each of 3 dimension
168   // prec    : max allowed absolute difference between the user function and computed parameterization on the requested grid
169   //
170   if (DimOut<1) {Error("AliCheb3D","Requested output dimension is %d\nStop\n",fDimOut); exit(1);}
171   if (DimOut<1) {Error("AliCheb3D","Requested output dimension is %d\nStop\n",fDimOut); exit(1);}
172   for (int i=3;i--;) {
173     fBMin[i] = fBMax[i] = fBScale[i] = fBOffset[i] = 0;
174     fNPoints[i] = 0;
175     fGridOffs[i] = 0.;
176   }
177   SetDimOut(DimOut);
178   PrepareBoundaries(bmin,bmax);
179   DefineGrid(npoints);
180   SetUsrFunction(ptr);
181   ChebFit();
182   //
183 }
184 #endif
185
186 //__________________________________________________________________________________________
187 #ifdef _INC_CREATION_ALICHEB3D_
188 AliCheb3D::AliCheb3D(void (*ptr)(float*,float*), int DimOut, Float_t  *bmin,Float_t  *bmax, Int_t *npX,Int_t *npY,Int_t *npZ, Float_t prec) : 
189   fDimOut(0), 
190   fPrec(TMath::Max(1.E-12f,prec)), 
191   fChebCalc(1), 
192   fMaxCoefs(0), 
193   fResTmp(0), 
194   fGrid(0), 
195   fUsrFunName(""),
196   fUsrMacro(0)
197 {
198   // Construct very economic  parameterization for the function
199   // ptr     : pointer on the function: void fun(Float_t * inp,Float_t * out)
200   // DimOut  : dimension of the vector computed by the user function
201   // bmin    : array of 3 elements with the lower boundaries of the region where the function is defined
202   // bmax    : array of 3 elements with the upper boundaries of the region where the function is defined
203   // npX     : array of 3 elements with the number of points to compute in each dimension for 1st component 
204   // npY     : array of 3 elements with the number of points to compute in each dimension for 2nd component 
205   // npZ     : array of 3 elements with the number of points to compute in each dimension for 3d  component 
206   // prec    : max allowed absolute difference between the user function and computed parameterization on the requested grid
207   //
208   if (DimOut<1) {Error("AliCheb3D","Requested output dimension is %d\nStop\n",fDimOut); exit(1);}
209   if (DimOut<1) {Error("AliCheb3D","Requested output dimension is %d\nStop\n",fDimOut); exit(1);}
210   for (int i=3;i--;) {
211     fBMin[i] = fBMax[i] = fBScale[i] = fBOffset[i] = 0;
212     fNPoints[i] = 0;
213     fGridOffs[i] = 0.;
214   }
215   SetDimOut(DimOut);
216   PrepareBoundaries(bmin,bmax);
217   SetUsrFunction(ptr);
218   //
219   DefineGrid(npX);
220   ChebFit(0);
221   DefineGrid(npY);
222   ChebFit(1);
223   DefineGrid(npZ);
224   ChebFit(2);
225   //
226 }
227 #endif
228
229
230 //__________________________________________________________________________________________
231 #ifdef _INC_CREATION_ALICHEB3D_
232 AliCheb3D::AliCheb3D(void (*ptr)(float*,float*), int DimOut, Float_t  *bmin,Float_t  *bmax, Float_t prec, Bool_t run) : 
233   fDimOut(0), 
234   fPrec(TMath::Max(1.E-12f,prec)), 
235   fChebCalc(1), 
236   fMaxCoefs(0), 
237   fResTmp(0), 
238   fGrid(0), 
239   fUsrFunName(""),
240   fUsrMacro(0)
241 {
242   // Construct very economic  parameterization for the function with automatic calculation of the root's grid
243   // ptr     : pointer on the function: void fun(Float_t * inp,Float_t * out)
244   // DimOut  : dimension of the vector computed by the user function
245   // bmin    : array of 3 elements with the lower boundaries of the region where the function is defined
246   // bmax    : array of 3 elements with the upper boundaries of the region where the function is defined
247   // prec    : max allowed absolute difference between the user function and computed parameterization on the requested grid
248   //
249   if (DimOut!=3) {Error("AliCheb3D","This constructor works only for 3D fits, %dD fit was requested\n",fDimOut); exit(1);}
250   if (DimOut<1) {Error("AliCheb3D","Requested output dimension is %d\nStop\n",fDimOut); exit(1);}
251   for (int i=3;i--;) {
252     fBMin[i] = fBMax[i] = fBScale[i] = fBOffset[i] = 0;
253     fNPoints[i] = 0;
254     fGridOffs[i] = 0.;
255   }
256   SetDimOut(DimOut);
257   PrepareBoundaries(bmin,bmax);
258   SetUsrFunction(ptr);
259   //
260   if (run) {
261     int gridNC[3][3];
262     EstimateNPoints(prec,gridNC);
263     DefineGrid(gridNC[0]);
264     ChebFit(0);
265     DefineGrid(gridNC[1]);
266     ChebFit(1);
267     DefineGrid(gridNC[2]);
268     ChebFit(2);
269   }
270   //
271 }
272 #endif
273
274
275 //__________________________________________________________________________________________
276 AliCheb3D& AliCheb3D::operator=(const AliCheb3D& rhs)
277 {
278   // assignment operator
279   //
280   if (this != &rhs) {
281     Clear();
282     fDimOut   = rhs.fDimOut;
283     fPrec     = rhs.fPrec;
284     fMaxCoefs = rhs.fMaxCoefs;
285     fUsrFunName = rhs.fUsrFunName;
286     fUsrMacro   = 0;
287     for (int i=3;i--;) {
288       fBMin[i]    = rhs.fBMin[i];
289       fBMax[i]    = rhs.fBMax[i];
290       fBScale[i]  = rhs.fBScale[i];
291       fBOffset[i] = rhs.fBOffset[i];
292       fNPoints[i] = rhs.fNPoints[i];
293     } 
294     for (int i=0;i<fDimOut;i++) {
295       AliCheb3DCalc* cbc = rhs.GetChebCalc(i);
296       if (cbc) fChebCalc.AddAtAndExpand(new AliCheb3DCalc(*cbc),i);
297     }    
298   }
299   return *this;
300   //
301 }
302
303 //__________________________________________________________________________________________
304 void AliCheb3D::Clear(const Option_t*)
305 {
306   // clear all dynamic structures
307   //
308   if (fResTmp)        { delete[] fResTmp; fResTmp = 0; }
309   if (fGrid)          { delete[] fGrid;   fGrid   = 0; }
310   if (fUsrMacro)      { delete fUsrMacro; fUsrMacro = 0;}
311   fChebCalc.SetOwner(kTRUE);
312   fChebCalc.Delete();
313   //
314 }
315
316 //__________________________________________________________________________________________
317 void AliCheb3D::Print(const Option_t* opt) const
318 {
319   // print info
320   //
321   printf("%s: Chebyshev parameterization for 3D->%dD function. Precision: %e\n",GetName(),fDimOut,fPrec);
322   printf("Region of validity: [%+.5e:%+.5e] [%+.5e:%+.5e] [%+.5e:%+.5e]\n",fBMin[0],fBMax[0],fBMin[1],fBMax[1],fBMin[2],fBMax[2]);
323   TString opts = opt; opts.ToLower();
324   if (opts.Contains("l")) for (int i=0;i<fDimOut;i++) {printf("Output dimension %d:\n",i+1); GetChebCalc(i)->Print();}
325   //
326 }
327
328 //__________________________________________________________________________________________
329 void AliCheb3D::PrepareBoundaries(const Float_t  *bmin, const Float_t  *bmax)
330 {
331   // Set and check boundaries defined by user, prepare coefficients for their conversion to [-1:1] interval
332   //
333   for (int i=3;i--;) {
334     fBMin[i]   = bmin[i];
335     fBMax[i]   = bmax[i];
336     fBScale[i] = bmax[i]-bmin[i];
337     if (fBScale[i]<=0) { 
338       Error("PrepareBoundaries","Boundaries for %d-th dimension are not increasing: %+.4e %+.4e\nStop\n",i,fBMin[i],fBMax[i]);
339       exit(1);
340     }
341     fBOffset[i] = bmin[i] + fBScale[i]/2.0;
342     fBScale[i] = 2./fBScale[i];
343   }
344   //
345 }
346
347
348 //__________________________________________________________________________________________
349 #ifdef _INC_CREATION_ALICHEB3D_
350
351 // Pointer on user function (faster altrnative to TMethodCall)
352 void (*gUsrFunAliCheb3D) (float* ,float* );
353
354 void AliCheb3D::EvalUsrFunction() 
355 {
356   // call user supplied function
357   if   (gUsrFunAliCheb3D) gUsrFunAliCheb3D(fArgsTmp,fResTmp);
358   else fUsrMacro->Execute(); 
359 }
360
361 void AliCheb3D::SetUsrFunction(const char* name)
362 {
363   // load user macro with function definition and compile it
364   gUsrFunAliCheb3D = 0; 
365   fUsrFunName = name;
366   gSystem->ExpandPathName(fUsrFunName);
367   if (fUsrMacro) delete fUsrMacro;
368   TString tmpst = fUsrFunName;
369   tmpst += "+"; // prepare filename to compile
370   if (gROOT->LoadMacro(tmpst.Data())) {Error("SetUsrFunction","Failed to load user function from %s\nStop\n",name); exit(1);}
371   fUsrMacro = new TMethodCall();        
372   tmpst = tmpst.Data() + tmpst.Last('/')+1; //Strip away any path preceding the macro file name
373   int dot = tmpst.Last('.');
374   if (dot>0) tmpst.Resize(dot);
375   fUsrMacro->InitWithPrototype(tmpst.Data(),"Float_t *,Float_t *");
376   long args[2];
377   args[0] = (long)fArgsTmp;
378   args[1] = (long)fResTmp;
379   fUsrMacro->SetParamPtrs(args); 
380   //
381 }
382 #endif
383
384 //__________________________________________________________________________________________
385 #ifdef _INC_CREATION_ALICHEB3D_
386 void AliCheb3D::SetUsrFunction(void (*ptr)(float*,float*))
387 {
388   // assign user training function
389   //
390   if (fUsrMacro) delete fUsrMacro;
391   fUsrMacro = 0;
392   fUsrFunName = "";
393   gUsrFunAliCheb3D = ptr;
394 }
395 #endif
396
397 //__________________________________________________________________________________________
398 #ifdef _INC_CREATION_ALICHEB3D_
399 void AliCheb3D::EvalUsrFunction(const Float_t  *x, Float_t  *res) 
400 {
401   // evaluate user function value
402   //
403   for (int i=3;i--;) fArgsTmp[i] = x[i];
404   if   (gUsrFunAliCheb3D) gUsrFunAliCheb3D(fArgsTmp,fResTmp);
405   else fUsrMacro->Execute(); 
406   for (int i=fDimOut;i--;) res[i] = fResTmp[i];
407 }
408 #endif
409
410 //__________________________________________________________________________________________
411 #ifdef _INC_CREATION_ALICHEB3D_
412 Int_t AliCheb3D::CalcChebCoefs(const Float_t  *funval,int np, Float_t  *outCoefs, Float_t  prec)
413 {
414   // Calculate Chebyshev coeffs using precomputed function values at np roots.
415   // If prec>0, estimate the highest coeff number providing the needed precision
416   //
417   double sm;                 // do summations in double to minimize the roundoff error
418   for (int ic=0;ic<np;ic++) { // compute coeffs
419     sm = 0;          
420     for (int ir=0;ir<np;ir++) {
421       float  rt = TMath::Cos( ic*(ir+0.5)*TMath::Pi()/np);
422       sm += funval[ir]*rt;
423     }
424     outCoefs[ic] = Float_t( sm * ((ic==0) ? 1./np : 2./np) );
425   }
426   //
427   if (prec<=0) return np;
428   //
429   sm = 0;
430   int cfMax = 0;
431   for (cfMax=np;cfMax--;) {
432     sm += TMath::Abs(outCoefs[cfMax]);
433     if (sm>=prec) break;
434   }
435   if (++cfMax==0) cfMax=1;
436   return cfMax;
437   //
438 }
439 #endif
440
441 //__________________________________________________________________________________________
442 #ifdef _INC_CREATION_ALICHEB3D_
443 void AliCheb3D::DefineGrid(Int_t* npoints)
444 {
445   // prepare the grid of Chebyshev roots in each dimension
446   const int kMinPoints = 1;
447   int ntot = 0;
448   fMaxCoefs = 1;
449   for (int id=3;id--;) { 
450     fNPoints[id] = npoints[id];
451     if (fNPoints[id]<kMinPoints) {
452       Error("DefineGrid","at %d-th dimension %d point is requested, at least %d is needed\nStop\n",id,fNPoints[id],kMinPoints);
453       exit(1);
454     }
455     ntot += fNPoints[id];
456     fMaxCoefs *= fNPoints[id];
457   }
458   printf("Computing Chebyshev nodes on [%2d/%2d/%2d] grid\n",npoints[0],npoints[1],npoints[2]);
459   if (fGrid) delete[] fGrid;
460   fGrid = new Float_t [ntot];
461   //
462   int curp = 0;
463   for (int id=3;id--;) { 
464     int np = fNPoints[id];
465     fGridOffs[id] = curp;
466     for (int ip=0;ip<np;ip++) {
467       Float_t x = TMath::Cos( TMath::Pi()*(ip+0.5)/np );
468       fGrid[curp++] = MapToExternal(x,id);
469     }
470   }
471   //
472 }
473 #endif
474
475 //__________________________________________________________________________________________
476 #ifdef _INC_CREATION_ALICHEB3D_
477 Int_t AliCheb3D::ChebFit()
478 {
479   // prepare parameterization for all output dimensions
480   int ir=0; 
481   for (int i=fDimOut;i--;) ir+=ChebFit(i); 
482   return ir;
483 }
484 #endif
485
486 //__________________________________________________________________________________________
487 #ifdef _INC_CREATION_ALICHEB3D_
488 Int_t AliCheb3D::ChebFit(int dmOut)
489 {
490   // prepare paramaterization of 3D function for dmOut-th dimension 
491   int maxDim = 0;
492   for (int i=0;i<3;i++) if (maxDim<fNPoints[i]) maxDim = fNPoints[i];
493   Float_t  *fvals      = new Float_t [ fNPoints[0] ];
494   Float_t  *tmpCoef3D  = new Float_t [ fNPoints[0]*fNPoints[1]*fNPoints[2] ]; 
495   Float_t  *tmpCoef2D  = new Float_t [ fNPoints[0]*fNPoints[1] ]; 
496   Float_t  *tmpCoef1D  = new Float_t [ maxDim ];
497   //
498   Float_t rTiny = 0.1*fPrec/Float_t(maxDim); // neglect coefficient below this threshold
499   //
500   // 1D Cheb.fit for 0-th dimension at current steps of remaining dimensions
501   int ncmax = 0;
502   //
503   printf("Dim%d : 00.00%% Done",dmOut);fflush(stdout);
504   AliCheb3DCalc* cheb =  GetChebCalc(dmOut);
505   //
506   float ncals2count = fNPoints[2]*fNPoints[1]*fNPoints[0];
507   float ncals = 0;
508   float frac = 0;
509   float fracStep = 0.001;
510   //
511   for (int id2=fNPoints[2];id2--;) {
512     fArgsTmp[2] = fGrid[ fGridOffs[2]+id2 ];
513     //
514     for (int id1=fNPoints[1];id1--;) {
515       fArgsTmp[1] = fGrid[ fGridOffs[1]+id1 ];
516       //
517       for (int id0=fNPoints[0];id0--;) {
518         fArgsTmp[0] = fGrid[ fGridOffs[0]+id0 ];
519         EvalUsrFunction();     // compute function values at Chebyshev roots of 0-th dimension
520         fvals[id0] =  fResTmp[dmOut];
521         float fr = (++ncals)/ncals2count;
522         if (fr-frac>=fracStep) {
523           frac = fr;
524           printf("\b\b\b\b\b\b\b\b\b\b\b");
525           printf("%05.2f%% Done",fr*100);
526           fflush(stdout);
527         }
528         //
529       }
530       int nc = CalcChebCoefs(fvals,fNPoints[0], tmpCoef1D, fPrec);
531       for (int id0=fNPoints[0];id0--;) tmpCoef2D[id1 + id0*fNPoints[1]] = tmpCoef1D[id0];
532       if (ncmax<nc) ncmax = nc;              // max coefs to be kept in dim0 to guarantee needed precision
533     }
534     //
535     // once each 1d slice of given 2d slice is parametrized, parametrize the Cheb.coeffs
536     for (int id0=fNPoints[0];id0--;) {
537       CalcChebCoefs( tmpCoef2D+id0*fNPoints[1], fNPoints[1], tmpCoef1D, -1);
538       for (int id1=fNPoints[1];id1--;) tmpCoef3D[id2 + fNPoints[2]*(id1+id0*fNPoints[1])] = tmpCoef1D[id1];
539     }
540   }
541   //
542   // now fit the last dimensions Cheb.coefs
543   for (int id0=fNPoints[0];id0--;) {
544     for (int id1=fNPoints[1];id1--;) {
545       CalcChebCoefs( tmpCoef3D+ fNPoints[2]*(id1+id0*fNPoints[1]), fNPoints[2], tmpCoef1D, -1);
546       for (int id2=fNPoints[2];id2--;) tmpCoef3D[id2+ fNPoints[2]*(id1+id0*fNPoints[1])] = tmpCoef1D[id2]; // store on place
547     }
548   }
549   //
550   // now find 2D surface which separates significant coefficients of 3D matrix from nonsignificant ones (up to fPrec)
551   UShort_t *tmpCoefSurf = new UShort_t[ fNPoints[0]*fNPoints[1] ];
552   for (int id0=fNPoints[0];id0--;) for (int id1=fNPoints[1];id1--;) tmpCoefSurf[id1+id0*fNPoints[1]]=0;  
553   Double_t resid = 0;
554   for (int id0=fNPoints[0];id0--;) {
555     for (int id1=fNPoints[1];id1--;) {
556       for (int id2=fNPoints[2];id2--;) {
557         int id = id2 + fNPoints[2]*(id1+id0*fNPoints[1]);
558         Float_t  cfa = TMath::Abs(tmpCoef3D[id]);
559         if (cfa < rTiny) {tmpCoef3D[id] = 0; continue;} // neglect coefs below the threshold
560         resid += cfa;
561         if (resid<fPrec) continue; // this coeff is negligible
562         // otherwise go back 1 step
563         resid -= cfa;
564         tmpCoefSurf[id1+id0*fNPoints[1]] = id2+1; // how many coefs to keep
565         break;
566       }
567     }
568   }
569   /*
570   printf("\n\nCoeffs\n");  
571   int cnt = 0;
572   for (int id0=0;id0<fNPoints[0];id0++) {
573     for (int id1=0;id1<fNPoints[1];id1++) {
574       for (int id2=0;id2<fNPoints[2];id2++) {
575         printf("%2d%2d%2d %+.4e |",id0,id1,id2,tmpCoef3D[cnt++]);
576       }
577       printf("\n");
578     }
579     printf("\n");
580   }
581   */
582   // see if there are rows to reject, find max.significant column at each row
583   int nRows = fNPoints[0];
584   UShort_t *tmpCols = new UShort_t[nRows]; 
585   for (int id0=fNPoints[0];id0--;) {
586     int id1 = fNPoints[1];
587     while (id1>0 && tmpCoefSurf[(id1-1)+id0*fNPoints[1]]==0) id1--;
588     tmpCols[id0] = id1;
589   }
590   // find max significant row
591   for (int id0=nRows;id0--;) {if (tmpCols[id0]>0) break; nRows--;}
592   // find max significant column and fill the permanent storage for the max sigificant column of each row
593   cheb->InitRows(nRows);                  // create needed arrays;
594   UShort_t *nColsAtRow = cheb->GetNColsAtRow();
595   UShort_t *colAtRowBg = cheb->GetColAtRowBg();
596   int nCols = 0;
597   int NElemBound2D = 0;
598   for (int id0=0;id0<nRows;id0++) {
599     nColsAtRow[id0] = tmpCols[id0];     // number of columns to store for this row
600     colAtRowBg[id0] = NElemBound2D;     // begining of this row in 2D boundary surface
601     NElemBound2D += tmpCols[id0];
602     if (nCols<nColsAtRow[id0]) nCols = nColsAtRow[id0];
603   }
604   cheb->InitCols(nCols);
605   delete[] tmpCols;
606   //  
607   // create the 2D matrix defining the boundary of significance for 3D coeffs.matrix 
608   // and count the number of siginifacnt coefficients
609   //
610   cheb->InitElemBound2D(NElemBound2D);
611   UShort_t *coefBound2D0 = cheb->GetCoefBound2D0();
612   UShort_t *coefBound2D1 = cheb->GetCoefBound2D1();
613   fMaxCoefs = 0; // redefine number of coeffs
614   for (int id0=0;id0<nRows;id0++) {
615     int nCLoc = nColsAtRow[id0];
616     int col0  = colAtRowBg[id0];
617     for (int id1=0;id1<nCLoc;id1++) {
618       coefBound2D0[col0 + id1] = tmpCoefSurf[id1+id0*fNPoints[1]];  // number of coefs to store for 3-d dimension
619       coefBound2D1[col0 + id1] = fMaxCoefs;
620       fMaxCoefs += coefBound2D0[col0 + id1];
621     }
622   }
623   //
624   // create final compressed 3D matrix for significant coeffs
625   cheb->InitCoefs(fMaxCoefs);
626   Float_t  *coefs = cheb->GetCoefs();
627   int count = 0;
628   for (int id0=0;id0<nRows;id0++) {
629     int ncLoc = nColsAtRow[id0];
630     int col0  = colAtRowBg[id0];
631     for (int id1=0;id1<ncLoc;id1++) {
632       int ncf2 = coefBound2D0[col0 + id1];
633       for (int id2=0;id2<ncf2;id2++) {
634         coefs[count++] = tmpCoef3D[id2 + fNPoints[2]*(id1+id0*fNPoints[1])];
635       }
636     }
637   }
638   /*
639   printf("\n\nNewSurf\n");
640   for (int id0=0;id0<fNPoints[0];id0++) {
641     for (int id1=0;id1<fNPoints[1];id1++) {
642       printf("(%2d %2d) %2d |",id0,id1,tmpCoefSurf[id1+id0*fNPoints[1]]);  
643     }
644     printf("\n");
645   }
646   */
647   //
648   delete[] tmpCoefSurf;
649   delete[] tmpCoef1D;
650   delete[] tmpCoef2D;
651   delete[] tmpCoef3D;
652   delete[] fvals;
653   //
654   printf("\b\b\b\b\b\b\b\b\b\b\b\b");
655   printf("100.00%% Done\n");
656   return 1;
657 }
658 #endif
659
660 //_______________________________________________
661 #ifdef _INC_CREATION_ALICHEB3D_
662 void AliCheb3D::SaveData(const char* outfile,Bool_t append) const
663 {
664   // writes coefficients data to output text file, optionallt appending on the end of existing file
665   TString strf = outfile;
666   gSystem->ExpandPathName(strf);
667   FILE* stream = fopen(strf,append ? "a":"w");
668   SaveData(stream);
669   fclose(stream);
670   //
671 }
672 #endif
673
674 //_______________________________________________
675 #ifdef _INC_CREATION_ALICHEB3D_
676 void AliCheb3D::SaveData(FILE* stream) const
677 {
678   // writes coefficients data to existing output stream
679   //
680   fprintf(stream,"\n# These are automatically generated data for the Chebyshev interpolation of 3D->%dD function\n",fDimOut); 
681   fprintf(stream,"#\nSTART %s\n",GetName());
682   fprintf(stream,"# Dimensionality of the output\n%d\n",fDimOut);
683   fprintf(stream,"# Interpolation abs. precision\n%+.8e\n",fPrec);
684   //
685   fprintf(stream,"# Lower boundaries of interpolation region\n");
686   for (int i=0;i<3;i++) fprintf(stream,"%+.8e\n",fBMin[i]);
687   fprintf(stream,"# Upper boundaries of interpolation region\n");
688   for (int i=0;i<3;i++) fprintf(stream,"%+.8e\n",fBMax[i]);
689   fprintf(stream,"# Parameterization for each output dimension follows:\n");
690   //
691   for (int i=0;i<fDimOut;i++) GetChebCalc(i)->SaveData(stream);
692   fprintf(stream,"#\nEND %s\n#\n",GetName());
693   //
694 }
695 #endif
696
697 //__________________________________________________________________________________________
698 #ifdef _INC_CREATION_ALICHEB3D_
699 void AliCheb3D::InvertSign()
700 {
701   // invert the sign of all parameterizations
702   for (int i=fDimOut;i--;) {
703     AliCheb3DCalc* par =  GetChebCalc(i);
704     int ncf = par->GetNCoefs();
705     float *coefs = par->GetCoefs();
706     for (int j=ncf;j--;) coefs[j] = -coefs[j];
707   }
708 }
709 #endif
710
711
712 //_______________________________________________
713 void AliCheb3D::LoadData(const char* inpFile)
714 {
715   // load coefficients data from txt file
716   //
717   TString strf = inpFile;
718   gSystem->ExpandPathName(strf);
719   FILE* stream = fopen(strf.Data(),"r");
720   LoadData(stream);
721   fclose(stream);
722   //
723 }
724
725 //_______________________________________________
726 void AliCheb3D::LoadData(FILE* stream)
727 {
728   // load coefficients data from stream
729   //
730   if (!stream) {Error("LoadData","No stream provided.\nStop"); exit(1);}
731   TString buffs;
732   Clear();
733   AliCheb3DCalc::ReadLine(buffs,stream);
734   if (!buffs.BeginsWith("START")) {Error("LoadData","Expected: \"START <fit_name>\", found \"%s\"\nStop\n",buffs.Data());exit(1);}
735   SetName(buffs.Data()+buffs.First(' ')+1);
736   //
737   AliCheb3DCalc::ReadLine(buffs,stream); // N output dimensions
738   fDimOut = buffs.Atoi(); 
739   if (fDimOut<1) {Error("LoadData","Expected: '<number_of_output_dimensions>', found \"%s\"\nStop\n",buffs.Data());exit(1);}
740   //
741   SetDimOut(fDimOut);
742   //
743   AliCheb3DCalc::ReadLine(buffs,stream); // Interpolation abs. precision
744   fPrec = buffs.Atof();
745   if (fPrec<=0) {Error("LoadData","Expected: '<abs.precision>', found \"%s\"\nStop\n",buffs.Data());exit(1);}
746   //
747   for (int i=0;i<3;i++) { // Lower boundaries of interpolation region
748     AliCheb3DCalc::ReadLine(buffs,stream);
749     fBMin[i] = buffs.Atof(); 
750   }
751   for (int i=0;i<3;i++) { // Upper boundaries of interpolation region
752     AliCheb3DCalc::ReadLine(buffs,stream);
753     fBMax[i] = buffs.Atof(); 
754   }
755   PrepareBoundaries(fBMin,fBMax);
756   //
757   // data for each output dimension
758   for (int i=0;i<fDimOut;i++) GetChebCalc(i)->LoadData(stream);
759   //
760   // check end_of_data record
761   AliCheb3DCalc::ReadLine(buffs,stream);
762   if (!buffs.BeginsWith("END") || !buffs.Contains(GetName())) {
763     Error("LoadData","Expected \"END %s\", found \"%s\".\nStop\n",GetName(),buffs.Data());
764     exit(1);
765   }
766   //
767 }
768
769 //_______________________________________________
770 void AliCheb3D::SetDimOut(const int d)
771 {
772   // init output dimensions
773   fDimOut = d;
774   if (fResTmp) delete fResTmp;
775   fResTmp = new Float_t[fDimOut];
776   fChebCalc.Delete();
777   for (int i=0;i<d;i++) fChebCalc.AddAtAndExpand(new AliCheb3DCalc(),i);
778 }
779
780 //_______________________________________________
781 void AliCheb3D::ShiftBound(int id,float dif)
782 {
783   // modify the bounds of the grid
784   //
785   if (id<0||id>2) {printf("Maximum 3 dimensions are supported\n"); return;}
786   fBMin[id] += dif;
787   fBMax[id] += dif;
788   fBOffset[id] += dif;
789 }
790
791 //_______________________________________________
792 #ifdef _INC_CREATION_ALICHEB3D_
793 TH1* AliCheb3D::TestRMS(int idim,int npoints,TH1* histo)
794 {
795   // fills the difference between the original function and parameterization (for idim-th component of the output)
796   // to supplied histogram. Calculations are done in npoints random points. 
797   // If the hostgram was not supplied, it will be created. It is up to the user to delete it! 
798   if (!fUsrMacro) {
799     printf("No user function is set\n");
800     return 0;
801   }
802   if (!histo) histo = new TH1D(GetName(),"Control: Function - Parametrization",100,-2*fPrec,2*fPrec);
803   for (int ip=npoints;ip--;) {
804     gRandom->RndmArray(3,(Float_t *)fArgsTmp);
805     for (int i=3;i--;) fArgsTmp[i] = fBMin[i] + fArgsTmp[i]*(fBMax[i]-fBMin[i]);
806     EvalUsrFunction();
807     Float_t valFun = fResTmp[idim];
808     Eval(fArgsTmp,fResTmp);
809     Float_t valPar = fResTmp[idim];
810     histo->Fill(valFun - valPar);
811   }
812   return histo;
813   //
814 }
815 #endif
816
817 //_______________________________________________
818 #ifdef _INC_CREATION_ALICHEB3D_
819
820 void AliCheb3D::EstimateNPoints(float Prec, int gridBC[3][3],Int_t npd1,Int_t npd2,Int_t npd3)
821 {
822   // Estimate number of points to generate a training data
823   //
824   const int kScp = 9;
825   const float kScl[9] = {0.1,0.2,0.3,0.4,0.5,0.6,0.7,0.8,0.9};
826   //
827   const float sclDim[2] = {0.001,0.999};
828   const int   compDim[3][2] = { {1,2}, {2,0}, {0,1} };
829   static float xyz[3];
830   Int_t npdTst[3] = {npd1,npd2,npd3};
831   //
832
833   for (int i=3;i--;)for (int j=3;j--;) gridBC[i][j] = -1;
834   //
835   for (int idim=0;idim<3;idim++) {
836     float dimMN = fBMin[idim] + sclDim[0]*(fBMax[idim]-fBMin[idim]);
837     float dimMX = fBMin[idim] + sclDim[1]*(fBMax[idim]-fBMin[idim]);
838     //
839     int id1 = compDim[idim][0]; // 1st fixed dim
840     int id2 = compDim[idim][1]; // 2nd fixed dim
841     for (int i1=0;i1<kScp;i1++) {
842       xyz[ id1 ] = fBMin[id1] + kScl[i1]*( fBMax[id1]-fBMin[id1] );
843       for (int i2=0;i2<kScp;i2++) {
844         xyz[ id2 ] = fBMin[id2] + kScl[i2]*( fBMax[id2]-fBMin[id2] );
845         int* npt = GetNCNeeded(xyz,idim, dimMN,dimMX, Prec, npdTst[idim]); // npoints for Bx,By,Bz
846         for (int ib=0;ib<3;ib++) if (npt[ib]>gridBC[ib][idim]) gridBC[ib][idim] = npt[ib];
847       }
848     }
849   }
850 }
851
852 /*
853 void AliCheb3D::EstimateNPoints(float Prec, int gridBC[3][3])
854 {
855   // Estimate number of points to generate a training data
856   //
857   const float sclA[9] = {0.1, 0.5, 0.9, 0.1, 0.5, 0.9, 0.1, 0.5, 0.9} ;
858   const float sclB[9] = {0.1, 0.1, 0.1, 0.5, 0.5, 0.5, 0.9, 0.9, 0.9} ;
859   const float sclDim[2] = {0.01,0.99};
860   const int   compDim[3][2] = { {1,2}, {2,0}, {0,1} };
861   static float xyz[3];
862   //
863   for (int i=3;i--;)for (int j=3;j--;) gridBC[i][j] = -1;
864   //
865   for (int idim=0;idim<3;idim++) {
866     float dimMN = fBMin[idim] + sclDim[0]*(fBMax[idim]-fBMin[idim]);
867     float dimMX = fBMin[idim] + sclDim[1]*(fBMax[idim]-fBMin[idim]);
868     //
869     for (int it=0;it<9;it++) { // test in 9 points
870       int id1 = compDim[idim][0]; // 1st fixed dim
871       int id2 = compDim[idim][1]; // 2nd fixed dim
872       xyz[ id1 ] = fBMin[id1] + sclA[it]*( fBMax[id1]-fBMin[id1] );
873       xyz[ id2 ] = fBMin[id2] + sclB[it]*( fBMax[id2]-fBMin[id2] );
874       //
875       int* npt = GetNCNeeded(xyz,idim, dimMN,dimMX, Prec); // npoints for Bx,By,Bz
876       for (int ib=0;ib<3;ib++) if (npt[ib]>gridBC[ib][idim]) gridBC[ib][idim] = npt[ib];//+2;
877       //
878     }
879   }
880 }
881
882
883 int* AliCheb3D::GetNCNeeded(float xyz[3],int DimVar, float mn,float mx, float prec)
884 {
885   // estimate needed number of chebyshev coefs for given function description in DimVar dimension
886   // The values for two other dimensions must be set beforehand
887   //
888   static int curNC[3];
889   static int retNC[3];
890   const int kMaxPoint = 400;
891   float* gridVal = new float[3*kMaxPoint];
892   float* coefs   = new float[3*kMaxPoint];
893   //
894   float scale = mx-mn;
895   float offs  = mn + scale/2.0;
896   scale = 2./scale;
897   // 
898   int curNP;
899   int maxNC=-1;
900   int maxNCPrev=-1;
901   for (int i=0;i<3;i++) retNC[i] = -1;
902   for (int i=0;i<3;i++) fArgsTmp[i] = xyz[i];
903   //
904   for (curNP=3; curNP<kMaxPoint; curNP+=3) { 
905     maxNCPrev = maxNC;
906     //
907     for (int i=0;i<curNP;i++) { // get function values on Cheb. nodes
908       float x = TMath::Cos( TMath::Pi()*(i+0.5)/curNP );
909       fArgsTmp[DimVar] =  x/scale+offs; // map to requested interval
910       EvalUsrFunction();
911       for (int ib=3;ib--;) gridVal[ib*kMaxPoint + i] = fResTmp[ib];
912     }
913     //
914     for (int ib=0;ib<3;ib++) {
915       curNC[ib] = AliCheb3D::CalcChebCoefs(&gridVal[ib*kMaxPoint], curNP, &coefs[ib*kMaxPoint],prec);
916       if (maxNC < curNC[ib]) maxNC = curNC[ib];
917       if (retNC[ib] < curNC[ib]) retNC[ib] = curNC[ib];
918     }
919     if ( (curNP-maxNC)>3 &&  (maxNC-maxNCPrev)<1 ) break;
920     maxNCPrev = maxNC;
921     //
922   }
923   delete[] gridVal;
924   delete[] coefs;
925   return retNC;
926   //
927 }
928 */
929
930
931 int* AliCheb3D::GetNCNeeded(float xyz[3],int DimVar, float mn,float mx, float prec, Int_t npCheck)
932 {
933   // estimate needed number of chebyshev coefs for given function description in DimVar dimension
934   // The values for two other dimensions must be set beforehand
935   //
936   static int retNC[3];
937   static int npChLast = 0;
938   static float *gridVal=0,*coefs=0;
939   if (npCheck<3) npCheck = 3;
940   if (npChLast<npCheck) {
941     if (gridVal) delete[] gridVal;
942     if (coefs)   delete[] coefs;
943     gridVal = new float[3*npCheck];
944     coefs   = new float[3*npCheck];
945     npChLast = npCheck;
946   }
947   //
948   float scale = mx-mn;
949   float offs  = mn + scale/2.0;
950   scale = 2./scale;
951   //
952   for (int i=0;i<3;i++) fArgsTmp[i] = xyz[i];
953   for (int i=0;i<npCheck;i++) {
954     fArgsTmp[DimVar] =  TMath::Cos( TMath::Pi()*(i+0.5)/npCheck)/scale+offs; // map to requested interval
955     EvalUsrFunction();
956     for (int ib=3;ib--;) gridVal[ib*npCheck + i] = fResTmp[ib];
957   } 
958   //
959   for (int ib=0;ib<3;ib++) retNC[ib] = AliCheb3D::CalcChebCoefs(&gridVal[ib*npCheck], npCheck, &coefs[ib*npCheck],prec);
960   return retNC;
961   //
962 }
963
964
965 #endif