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