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