kNN algorithm improved. IO Defined
[u/mrichter/AliRoot.git] / STAT / TKDInterpolator.cxx
1 #include "TKDInterpolator.h"
2
3 #include "TLinearFitter.h"
4 #include "TVector.h"
5 #include "TTree.h"
6 #include "TH2.h"
7 #include "TObjArray.h"
8 #include "TObjString.h"
9 #include "TPad.h"
10 #include "TBox.h"
11 #include "TGraph.h"
12 #include "TMarker.h"
13 #include "TRandom.h"
14 #include "TROOT.h"
15
16 ClassImp(TKDInterpolator)
17 ClassImp(TKDInterpolator::TKDNodeInfo)
18
19 /////////////////////////////////////////////////////////////////////
20 // Memory setup of protected data memebers
21 // fRefPoints : evaluation point of PDF for each terminal node of underlying KD Tree.
22 // | 1st terminal node (fNDim point coordinates) | 2nd terminal node (fNDim point coordinates) | ...
23 //
24 // fRefValues : evaluation value/error of PDF for each terminal node of underlying KD Tree.
25 // | 1st terminal node (value) | 2nd terminal node (value) | ... | 1st terminal node (error) | 2nd terminal node (error) | ...
26 //
27 // status = |0|0|0|0|0|1(tri-cubic weights)|1(STORE)|1 INT(0 COG )|
28 /////////////////////////////////////////////////////////////////////
29
30 //_________________________________________________________________
31 TKDInterpolator::TKDInterpolator() : TKDTreeIF()
32         ,fNTNodes(0)
33         ,fRefPoints(0x0)
34         ,fRefValues(0x0)
35         ,fCov(0x0)
36         ,fPar(0x0)
37         ,fPDFstatus(0x0)
38         ,fStatus(4)
39         ,fLambda(0)
40         ,fDepth(-1)
41         ,fBuffer(0x0)
42         ,fKDhelper(0x0)
43         ,fFitter(0x0)
44 {
45 // Default constructor. To be used with care since in this case building
46 // of data structure is completly left to the user responsability.
47 }
48
49 //_________________________________________________________________
50 TKDInterpolator::TKDInterpolator(Int_t npoints, Int_t ndim, UInt_t bsize, Float_t **data) : TKDTreeIF(npoints, ndim, bsize, data)
51         ,fNTNodes(GetNTerminalNodes())
52         ,fRefPoints(0x0)
53         ,fRefValues(0x0)
54         ,fCov(0x0)
55         ,fPar(0x0)
56         ,fPDFstatus(0x0)
57         ,fStatus(4)
58         ,fLambda(0)
59         ,fDepth(-1)
60         ,fBuffer(0x0)
61         ,fKDhelper(0x0)
62         ,fFitter(0x0)
63 {
64 // Wrapper constructor for the similar TKDTree one.
65         
66         Build();
67 }
68
69
70 //_________________________________________________________________
71 TKDInterpolator::TKDInterpolator(TTree *t, const Char_t *var, const Char_t *cut, UInt_t bsize, Long64_t nentries, Long64_t firstentry) : TKDTreeIF()
72         ,fNTNodes(0)
73         ,fRefPoints(0x0)
74         ,fRefValues(0x0)
75         ,fCov(0x0)
76         ,fPar(0x0)
77         ,fPDFstatus(0x0)
78         ,fStatus(4)
79         ,fLambda(0)
80         ,fDepth(-1)
81         ,fBuffer(0x0)
82         ,fKDhelper(0x0)
83         ,fFitter(0x0)
84 {
85 // Alocate data from a tree. The variables which have to be analysed are
86 // defined in the "var" parameter as a colon separated list. The format should
87 // be identical to that used by TTree::Draw().
88 //
89 // 
90
91         TObjArray *vars = TString(var).Tokenize(":");
92         fNDim = vars->GetEntriesFast(); fNDimm = 2*fNDim;
93         if(fNDim > 6/*kDimMax*/) Warning("TKDInterpolator(TTree*, const Char_t, const Char_t, UInt_t)", Form("Variable number exceed maximum dimension %d. Results are unpredictable.", 6/*kDimMax*/));
94         fBucketSize = bsize;
95
96         Int_t np;
97         Double_t *v;
98         for(int idim=0; idim<fNDim; idim++){
99                 if(!(np = t->Draw(((TObjString*)(*vars)[idim])->GetName(), cut, "goff", nentries, firstentry))){
100                         Warning("TKDInterpolator(TTree*, const Char_t, const Char_t, UInt_t)", Form("Can not access data for keys %s. Key defined on tree :", ((TObjString*)(*vars)[idim])->GetName() ));
101                         TIterator *it = (t->GetListOfLeaves())->MakeIterator();
102                         TObject *o;
103                         while(o = (*it)()) printf("\t%s\n", o->GetName());
104                         continue;
105                 }
106                 if(!fNpoints){
107                         fNpoints = np;
108                         Info("TKDInterpolator(TTree*, const Char_t, const Char_t, UInt_t)", Form("Allocating %d data points in %d dimensions.", fNpoints, fNDim));
109                         fData = new Float_t*[fNDim];
110                         for(int idim=0; idim<fNDim; idim++) fData[idim] = new Float_t[fNpoints];
111                         kDataOwner = kTRUE;
112                 }
113                 v = t->GetV1();
114                 for(int ip=0; ip<fNpoints; ip++) fData[idim][ip] = (Float_t)v[ip];
115         }
116         TKDTreeIF::Build();
117         fNTNodes = GetNTerminalNodes();
118         Build();
119 }
120
121 //_________________________________________________________________
122 TKDInterpolator::~TKDInterpolator()
123 {
124         if(fCov){
125                 delete [] fPar;
126                 delete [] fCov;
127                 delete [] fPDFstatus;
128         }
129         
130         if(fFitter) delete fFitter;
131         if(fKDhelper) delete fKDhelper;
132         if(fBuffer) delete [] fBuffer;
133         
134         if(fRefPoints){
135                 for(int idim=0; idim<fNDim; idim++) delete [] fRefPoints[idim] ;
136                 delete [] fRefPoints;
137         }
138         if(fRefValues) delete [] fRefValues;
139 }
140
141 //_________________________________________________________________
142 void TKDInterpolator::Build()
143 {
144 // Fill interpolator's data array i.e.
145 //  - estimation points 
146 //  - corresponding PDF values
147
148         if(!fBoundaries) MakeBoundaries();
149         fLambda = 1 + fNDim + fNDim*(fNDim+1)/2;
150
151         // allocate memory for data
152         fRefValues = new Float_t[fNTNodes];
153         fRefPoints = new Float_t*[fNDim];
154         for(int id=0; id<fNDim; id++){
155                 fRefPoints[id] = new Float_t[fNTNodes];
156                 for(int in=0; in<fNTNodes; in++) fRefPoints[id][in] = 0.;
157         }
158
159         Float_t *bounds = 0x0;
160         Int_t *indexPoints;
161         for(int inode=0, tnode = fNnodes; inode<fNTNodes-1; inode++, tnode++){
162                 fRefValues[inode] =  Float_t(fBucketSize)/fNpoints;
163                 bounds = GetBoundary(tnode);
164                 for(int idim=0; idim<fNDim; idim++) fRefValues[inode] /= (bounds[2*idim+1] - bounds[2*idim]);
165
166                 indexPoints = GetPointsIndexes(tnode);
167                 // loop points in this terminal node
168                 for(int idim=0; idim<fNDim; idim++){
169                         for(int ip = 0; ip<fBucketSize; ip++) fRefPoints[idim][inode] += fData[idim][indexPoints[ip]];
170                         fRefPoints[idim][inode] /= fBucketSize;
171                 }
172         }
173
174         // analyze last (incomplete) terminal node
175         Int_t counts = fNpoints%fBucketSize;
176         counts = counts ? counts : fBucketSize;
177         Int_t inode = fNTNodes - 1, tnode = inode + fNnodes;
178         fRefValues[inode] =  Float_t(counts)/fNpoints;
179         bounds = GetBoundary(tnode);
180         for(int idim=0; idim<fNDim; idim++) fRefValues[inode] /= (bounds[2*idim+1] - bounds[2*idim]);
181
182         indexPoints = GetPointsIndexes(tnode);
183         // loop points in this terminal node
184         for(int idim=0; idim<fNDim; idim++){
185                 for(int ip = 0; ip<counts; ip++) fRefPoints[idim][inode] += fData[idim][indexPoints[ip]];
186                 fRefPoints[idim][inode] /= counts;
187         }
188 }
189
190 //__________________________________________________________________
191 void TKDInterpolator::GetStatus()
192 {
193         printf("Interpolator Status :\n");
194         printf("  Method : %s\n", fStatus&1 ? "INT" : "COG");
195         printf("  Store  : %s\n", fStatus&2 ? "YES" : "NO");
196         printf("  Weights: %s\n", fStatus&4 ? "YES" : "NO");
197
198         printf("nodes %d\n", fNTNodes);        //Number of evaluation data points
199         printf("refs 0x%x\n", fRefPoints);    //[fNDim][fNTNodes]
200         printf("vals 0x%x\n", fRefValues);     //[fNTNodes]
201         for(int i=0; i<fNTNodes; i++){
202                 for(int idim=0; idim<fNDim; idim++) printf("%f ", fRefPoints[idim][i]);
203                 printf("[%f]\n", fRefValues[i]);
204         }
205
206         printf("cov 0x%x\n", fCov);           //[fNTNodes] cov matrix array for nodes
207         printf("pars 0x%x\n", fPar);           //[fNTNodes] parameters array for nodes
208         for(int i=0; i<fNTNodes; i++){
209                 printf("%d ", i);
210                 for(int ip=0; ip<3; ip++) printf("p%d[%f] ", ip, fPar[i](ip));
211                 printf("\n");
212         }
213         printf("pdf 0x%x\n", fPDFstatus);     //[fNTNodes] status bit for node's PDF
214 }
215
216 //_________________________________________________________________
217 Double_t TKDInterpolator::Eval(const Double_t *point, Double_t &result, Double_t &error)
218 {
219 // Evaluate PDF for "point". The result is returned in "result" and error in "error". The function returns the chi2 of the fit.
220 //
221 // Observations:
222 //
223 // 1. The default method used for interpolation is kCOG.
224 // 2. The initial number of neighbors used for the estimation is set to Int(alpha*fLambda) (alpha = 1.5)
225                         
226         Float_t pointF[50]; // local Float_t conversion for "point"
227         for(int idim=0; idim<fNDim; idim++) pointF[idim] = (Float_t)point[idim];
228         Int_t node = FindNode(pointF) - fNnodes;
229         if(fPDFstatus && (fStatus&1) && fPDFstatus[node]) return CookPDF(point, node, result, error);
230
231         // Allocate memory
232         if(!fBuffer) fBuffer = new Double_t[2*fLambda];
233         if(!fKDhelper) fKDhelper = new TKDTreeIF(fNTNodes, fNDim, 30, fRefPoints);
234         if(!fFitter) SetIntInterpolation(kFALSE);
235         
236         // generate parabolic for nD
237         //Float_t alpha = Float_t(2*lambda + 1) / fNTNodes; // the bandwidth or smoothing parameter
238         //Int_t npoints = Int_t(alpha * fNTNodes);
239         //printf("Params : %d NPoints %d\n", lambda, npoints);
240         // prepare workers
241
242         Int_t *index,  // indexes of NN 
243               ipar,    // local looping variable
244                                 npoints = Int_t(1.5*fLambda); // number of data points used for interpolation
245         Float_t *dist, // distances of NN
246                                         d,     // NN normalized distance
247                                         w0,    // work
248                                         w;     // tri-cubic weight function
249         Double_t sig   // bucket error 
250                 = TMath::Sqrt(1./fBucketSize);
251         do{
252                 // find nearest neighbors
253                 for(int idim=0; idim<fNDim; idim++) pointF[idim] = (Float_t)point[idim];
254                 if(!fKDhelper->FindNearestNeighbors(pointF, npoints+1, index, dist)){
255                         Error("Eval()", Form("Failed retriving %d neighbours for point:", npoints));
256                         for(int idim=0; idim<fNDim; idim++) printf("%f ", point[idim]);
257                         printf("\n");
258                         return -1;
259                 }
260                 // add points to fitter
261                 fFitter->ClearPoints();
262                 for(int in=0; in<npoints; in++){
263                         if(fStatus&1){ // INT
264                                 for(int idim=0; idim<fNDim; idim++) pointF[idim] = fRefPoints[idim][index[in]];
265                                 Float_t *bounds = GetBoundary(FindNode(pointF));
266                                 
267                                 ipar = 0;
268                                 for(int idim=0; idim<fNDim; idim++){
269                                         fBuffer[ipar++] = .5*(bounds[2*idim] + bounds[2*idim+1]);
270                                         fBuffer[ipar++] = (bounds[2*idim]*bounds[2*idim] + bounds[2*idim] * bounds[2*idim+1] + bounds[2*idim+1] * bounds[2*idim+1])/3.;
271                                         for(int jdim=idim+1; jdim<fNDim; jdim++) fBuffer[ipar++] = (bounds[2*idim] + bounds[2*idim+1]) * (bounds[2*jdim] + bounds[2*jdim+1]) * .25;
272                                 }
273                         } else { // COG
274                                 for(int idim=0; idim<fNDim; idim++) fBuffer[idim] = fRefPoints[idim][index[in]];
275                         }
276
277                         // calculate tri-cubic weighting function
278                         if(fStatus&4){
279                                 d = dist[in]/ dist[npoints];
280                                 w0 = (1. - d*d*d); w = w0*w0*w0;
281                         } else w = 1.;
282                          
283                         //for(int idim=0; idim<fNDim; idim++) printf("%f ", fBuffer[idim]);
284                         //printf("\nd[%f] w[%f] sig[%f]\n", d, w, sig);
285                         fFitter->AddPoint(fBuffer, fRefValues[index[in]], fRefValues[index[in]]*sig/w);
286                 }
287                 npoints += 4;
288         } while(fFitter->Eval());
289
290         // retrive fitter results
291         TMatrixD cov(fLambda, fLambda);
292         TVectorD par(fLambda);
293         fFitter->GetCovarianceMatrix(cov);
294         fFitter->GetParameters(par);
295         Double_t chi2 = fFitter->GetChisquare()/(npoints - 4 - fLambda);
296
297         // store results
298         if(fStatus&2 && fStatus&1){
299                 fPar[node] = par;
300                 fCov[node] = cov;
301                 fPDFstatus[node] = 1;
302         }
303                 
304         // Build df/dpi|x values
305         Double_t *fdfdp = &fBuffer[fLambda];
306         ipar = 0;
307         fdfdp[ipar++] = 1.;
308         for(int idim=0; idim<fNDim; idim++){
309                 fdfdp[ipar++] = point[idim];
310                 for(int jdim=idim; jdim<fNDim; jdim++) fdfdp[ipar++] = point[idim]*point[jdim];
311         }
312
313         // calculate estimation
314         result =0.; error = 0.;
315         for(int i=0; i<fLambda; i++){
316                 result += fdfdp[i]*par(i);
317                 for(int j=0; j<fLambda; j++) error += fdfdp[i]*fdfdp[j]*cov(i,j);
318         }       
319         error = TMath::Sqrt(error);
320
321         return chi2;
322 }
323
324 // //_________________________________________________________________
325 // Double_t TKDInterpolator::Eval1(const Double_t *point, Int_t npoints, Double_t &result, Double_t &error)
326 // {
327 // // Evaluate PDF at k-dimensional position "point". The initial number of
328 // // neighbour estimation points is set to "npoints". The default method
329 // // used for interpolation is kCOG.
330 // 
331 //      // calculate number of parameters in the parabolic expresion
332 //      Int_t lambda = 1 + fNDim + fNDim*(fNDim+1)/2;
333 // 
334 //      if(!fBuffer) fBuffer = new Double_t[lambda-1];
335 //      if(!fKDhelper) fKDhelper = new TKDTreeIF(GetNTerminalNodes(), fNDim, npoints, fRefPoints);
336 // 
337 //      if(!fFitter) fFitter = new TLinearFitter(lambda, Form("hyp%d", fNDim+1));
338 //      else fFitter->SetFormula(Form("hyp%d", fNDim+1));
339 // 
340 // 
341 //      Float_t pointF[50];
342 //      for(int idim=0; idim<fNDim; idim++) pointF[idim] = point[idim];
343 //      Int_t istart = 0;
344 //      Int_t *index, ipar;
345 //      Float_t *bounds, *dist, *w = new Float_t[fNDim];
346 //      Double_t uncertainty = TMath::Sqrt(1./fBucketSize);
347 //      fFitter->ClearPoints();
348 //      do{
349 //              if(!fKDhelper->FindNearestNeighbors(pointF, npoints+1, index, dist)){
350 //                      Error("Eval()", Form("Failed retriving %d neighbours for point:", npoints));
351 //                      for(int idim=0; idim<fNDim; idim++) printf("%f ", point[idim]);
352 //                      printf("\n");
353 //                      return -1;
354 //              }
355 //              for(int in=istart; in<npoints; in++){
356 //                      for(int idim=0; idim<fNDim; idim++) w[idim] = fRefPoints[idim][index[in]];
357 //                      bounds = GetBoundary(FindNode(w));
358 // 
359 //                      ipar = 0;
360 //                      for(int idim=0; idim<fNDim; idim++){
361 //                              fBuffer[ipar++] = .5*(bounds[2*idim] + bounds[2*idim+1]);
362 //                              fBuffer[ipar++] = (bounds[2*idim]*bounds[2*idim] + bounds[2*idim] * bounds[2*idim+1] + bounds[2*idim+1] * bounds[2*idim+1])/3.;
363 //                              for(int jdim=idim+1; jdim<fNDim; jdim++) fBuffer[ipar++] = (bounds[2*idim] + bounds[2*idim+1]) * (bounds[2*jdim] + bounds[2*jdim+1]) * .25;
364 //                      }
365 // 
366 //                      fFitter->AddPoint(fBuffer, fRefValues[index[in]], fRefValues[index[in]]*uncertainty);
367 //              }
368 //              istart = npoints;
369 //              npoints += 4;
370 //      } while(fFitter->Eval());
371 //      delete [] w;
372 // 
373 //      // calculate evaluation
374 // //   fFitter->PrintResults(3);
375 //      TMatrixD cov(lambda, lambda);
376 //      TVectorD par(lambda);
377 //      fFitter->GetCovarianceMatrix(cov);
378 //      fFitter->GetParameters(par);
379 // 
380 //      // Build temporary array to keep values df/dpi|x
381 //      Double_t f[100];
382 //      ipar = 0;
383 //      f[ipar++] = 1.;
384 //      for(int idim=0; idim<fNDim; idim++){
385 //              f[ipar++] = point[idim];
386 //              for(int jdim=idim; jdim<fNDim; jdim++) f[ipar++] = point[idim]*point[jdim];
387 //      }
388 //      result =0.; error = 0.;
389 //      for(int i=0; i<lambda; i++){
390 //              result += f[i]*par[i];
391 //              for(int j=0; j<lambda; j++) error += f[i]*f[j]*cov(i,j);
392 //      }
393 //      error = TMath::Sqrt(error);
394 //      Double_t chi2 = fFitter->GetChisquare()/(npoints - 4 - lambda);
395 // 
396 //      for(int ipar=0; ipar<lambda; ipar++) printf("%d %8.6e %8.6e\n", ipar, par[ipar], TMath::Sqrt(cov(ipar, ipar)));
397 //      printf("result %6.3f +- %6.3f [%f]\n", result, error, chi2);
398 //      return chi2;
399 // }
400
401
402 //_________________________________________________________________
403 void TKDInterpolator::DrawNodes(UInt_t ax1, UInt_t ax2, Int_t depth)
404 {
405 // Draw nodes structure projected on plane "ax1:ax2". The parameter
406 // "depth" specifies the bucket size per node. If depth == -1 draw only
407 // terminal nodes and evaluation points (default -1 i.e. bucket size per node equal bucket size specified by the user)
408 //
409 // Observation:
410 // This function creates the nodes (TBox) array for the specified depth
411 // but don't delete it. Abusing this function may cause memory leaks !
412
413
414         if(!fBoundaries) MakeBoundaries();
415
416         // Count nodes in specific view
417         Int_t nnodes = 0;
418         for(int inode = 0; inode <= 2*fNnodes; inode++){
419                 if(depth == -1){
420                         if(!IsTerminal(inode)) continue;
421                 } else if((inode+1) >> depth != 1) continue;
422                 nnodes++;
423         }
424
425         //printf("depth %d nodes %d\n", depth, nnodes);
426         
427         TH2 *h2 = 0x0;
428         if(!(h2 = (TH2S*)gROOT->FindObject("hNodes"))) h2 = new TH2S("hNodes", "", 100, fRange[2*ax1], fRange[2*ax1+1], 100, fRange[2*ax2], fRange[2*ax2+1]);
429         h2->GetXaxis()->SetTitle(Form("x_{%d}", ax1));
430         h2->GetYaxis()->SetTitle(Form("x_{%d}", ax2));
431         h2->Draw();
432         
433         const Float_t border = 0.;//1.E-4;
434         TBox *node_array = new TBox[nnodes], *node;
435         Float_t *bounds = 0x0;
436         nnodes = 0;
437         for(int inode = 0; inode <= 2*fNnodes; inode++){
438                 if(depth == -1){
439                         if(!IsTerminal(inode)) continue;
440                 } else if((inode+1) >> depth != 1) continue;
441
442                 node = &node_array[nnodes++];
443                 //node = new TBox(bounds[2*ax1]+border, bounds[2*ax2]+border, bounds[2*ax1+1]-border, bounds[2*ax2+1]-border);
444                 node->SetFillStyle(3002);       
445                 node->SetFillColor(50+Int_t(gRandom->Uniform()*50.));
446                 bounds = GetBoundary(inode);
447                 node->DrawBox(bounds[2*ax1]+border, bounds[2*ax2]+border, bounds[2*ax1+1]-border, bounds[2*ax2+1]-border);
448         }
449         if(depth != -1) return;
450
451         // Draw reference points
452         TGraph *ref = new TGraph(GetNTerminalNodes());
453         ref->SetMarkerStyle(3);
454         ref->SetMarkerSize(.7);
455         ref->SetMarkerColor(2);
456         for(int inode = 0; inode < GetNTerminalNodes(); inode++) ref->SetPoint(inode, fRefPoints[ax1][inode], fRefPoints[ax2][inode]);
457         ref->Draw("p");
458         return;
459 }
460
461 //_________________________________________________________________
462 void TKDInterpolator::DrawNode(Int_t tnode, UInt_t ax1, UInt_t ax2)
463 {
464 // Draw node "node" and the data points within.
465 //
466 // Observation:
467 // This function creates some graphical objects
468 // but don't delete it. Abusing this function may cause memory leaks !
469
470         if(tnode < 0 || tnode >= GetNTerminalNodes()){
471                 Warning("DrawNode()", Form("Terminal node %d outside defined range.", tnode));
472                 return;
473         }
474
475         Int_t inode = tnode;
476         tnode += fNnodes;
477         // select zone of interest in the indexes array
478         Int_t *index = GetPointsIndexes(tnode);
479         Int_t nPoints = (tnode == 2*fNnodes) ? fNpoints%fBucketSize : fBucketSize;
480
481         // draw data points
482         TGraph *g = new TGraph(nPoints);
483         g->SetMarkerStyle(7);
484         for(int ip = 0; ip<nPoints; ip++) g->SetPoint(ip, fData[ax1][index[ip]], fData[ax2][index[ip]]);
485
486         // draw estimation point
487         TMarker *m=new TMarker(fRefPoints[ax1][inode], fRefPoints[ax2][inode], 20);
488         m->SetMarkerColor(2);
489         m->SetMarkerSize(1.7);
490         
491         // draw node contour
492         Float_t *bounds = GetBoundary(tnode);
493         TBox *n = new TBox(bounds[2*ax1], bounds[2*ax2], bounds[2*ax1+1], bounds[2*ax2+1]);
494         n->SetFillStyle(0);
495
496         if(gPad) gPad->Clear(); 
497         g->Draw("ap");
498         m->Draw();
499         n->Draw();
500         
501         return;
502 }
503
504
505 //__________________________________________________________________
506 void TKDInterpolator::SetIntInterpolation(const Bool_t on)
507 {
508 // Set interpolation bit to "on" and build/delete memory
509         
510         if(on) fStatus += fStatus&1 ? 0 : 1;
511         else fStatus += fStatus&1 ? -1 : 0;
512         TString formula;
513         if(on) formula = Form("hyp%d", fLambda-1);
514         else {
515                 formula = "1";
516                 for(int idim=0; idim<fNDim; idim++){
517                         formula += Form("++x[%d]", idim);
518                         for(int jdim=idim; jdim<fNDim; jdim++) formula += Form("++x[%d]*x[%d]", idim, jdim);
519                 }
520         }
521         if(!fFitter) fFitter = new TLinearFitter(fLambda, formula.Data());
522         else fFitter->SetFormula(formula.Data());
523 }
524
525
526 //_________________________________________________________________
527 void TKDInterpolator::SetSetStore(const Bool_t on)
528 {
529 // Set store bit to "on" and build/delete memory
530         
531         if(on){
532                 fStatus += fStatus&2 ? 0 : 2;
533                 if(!fCov){
534                         fPDFstatus = new Bool_t[fNTNodes];
535                         fCov = new TMatrixD[fNTNodes];
536                         fPar = new TVectorD[fNTNodes];
537                         for(int i=0; i<fNTNodes; i++){
538                                 fPDFstatus[i] = kFALSE;
539                                 fCov[i].ResizeTo(fLambda, fLambda);
540                                 fPar[i].ResizeTo(fLambda);
541                         }
542                 }
543         } else {
544                 fStatus += fStatus&2 ? -2 : 0;
545                 if(fCov){
546                         delete [] fPar;
547                         delete [] fCov;
548                         delete [] fPDFstatus;
549                 }
550         }
551 }
552
553 //_________________________________________________________________
554 void TKDInterpolator::SetUseWeights(const Bool_t on)
555 {
556         if(on) fStatus += fStatus&4 ? 0 : 4;
557         else fStatus += fStatus&4 ? -4 : 0;
558 }
559
560
561 //_________________________________________________________________
562 Double_t TKDInterpolator::CookPDF(const Double_t *point, const Int_t node, Double_t &result, Double_t &error)
563 {
564 // Recalculate the PDF for one node from the results of interpolation (parameters and covariance matrix)
565
566         Info("CookPDF()", Form("Called for node %d", node));
567
568         if(!fBuffer) fBuffer = new Double_t[2*fLambda];
569         Double_t *fdfdp = &fBuffer[fLambda];
570         Int_t ipar = 0;
571         fdfdp[ipar++] = 1.;
572         for(int idim=0; idim<fNDim; idim++){
573                 fdfdp[ipar++] = point[idim];
574                 for(int jdim=idim; jdim<fNDim; jdim++) fdfdp[ipar++] = point[idim]*point[jdim];
575         }
576
577         // calculate estimation
578         result =0.; error = 0.;
579         for(int i=0; i<fLambda; i++){
580                 result += fdfdp[i]*fPar[node](i);
581                 for(int j=0; j<fLambda; j++) error += fdfdp[i]*fdfdp[j]*fCov[node](i,j);
582         }       
583         error = TMath::Sqrt(error);
584         printf("result[CookPDF] %6.3f +- %6.3f\n", result, error);
585
586         return 0.;
587 }
588