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