#include "TKDInterpolatorBase.h"
#include "TKDNodeInfo.h"
#include "TKDTree.h"
#include "TROOT.h"
#include "TClonesArray.h"
#include "TLinearFitter.h"
#include "TTree.h"
#include "TH2.h"
#include "TObjArray.h"
#include "TObjString.h"
#include "TMath.h"
#include "TVectorD.h"
#include "TMatrixD.h"
ClassImp(TKDInterpolatorBase)
TKDInterpolatorBase::TKDInterpolatorBase(Int_t dim) :
fNSize(dim)
,fNodes(NULL)
,fNodesDraw(NULL)
,fStatus(0)
,fLambda(1 + dim + (dim*(dim+1)>>1))
,fDepth(-1)
,fAlpha(.5)
,fRefPoints(NULL)
,fBuffer(NULL)
,fKDhelper(NULL)
,fFitter(NULL)
{
UseWeights();
}
Bool_t TKDInterpolatorBase::Build(Int_t n)
{
if(Int_t((1.+fAlpha)*fLambda) > n){
Error("TKDInterpolatorBase::Build()", "Minimum number of points [%d] needed for interpolation exceeds number of evaluation points [%d]. Please increase granularity.", Int_t((1.+fAlpha)*fLambda), n);
return kFALSE;
}
if(fNodes){
Warning("TKDInterpolatorBase::Build()", "Data already allocated.");
fNodes->Delete();
} else {
fNodes = new TClonesArray("TKDNodeInfo", n);
fNodes->SetOwner();
}
for(int in=0; in<n; in++) new ((*fNodes)[in]) TKDNodeInfo(fNSize);
return kTRUE;
}
Bool_t TKDInterpolatorBase::Bootstrap()
{
if(!fNodes){
Error("TKDInterpolatorBase::Bootstrap()", "Nodes missing. Nothing to bootstrap from.");
return kFALSE;
}
Int_t in = GetNTNodes(); TKDNodeInfo *n(NULL);
while(in--){
if(!(n=(TKDNodeInfo*)(*fNodes)[in])){
Error("TKDInterpolatorBase::Bootstrap()", "Node @ %d missing.", in);
return kFALSE;
}
n->Bootstrap();
if(!fNSize) fNSize = n->GetDimension();
}
fLambda = n->GetNpar();
return kTRUE;
}
TKDInterpolatorBase::~TKDInterpolatorBase()
{
if(fFitter) delete fFitter;
if(fKDhelper) delete fKDhelper;
if(fBuffer) delete [] fBuffer;
if(fRefPoints){
for(int idim=0; idim<fNSize; idim++) delete [] fRefPoints[idim] ;
delete [] fRefPoints;
}
if(fNodes){
fNodes->Delete();
delete fNodes;
}
if(fNodesDraw) delete [] fNodesDraw;
TH2 *h2=NULL;
if((h2 = (TH2S*)gROOT->FindObject("hKDnodes"))) delete h2;
}
Bool_t TKDInterpolatorBase::GetCOGPoint(Int_t inode, Float_t *&coord, Float_t &val, Float_t &err) const
{
if(inode < 0 || inode > GetNTNodes()) return kFALSE;
TKDNodeInfo *node = (TKDNodeInfo*)(*fNodes)[inode];
coord = &(node->Data()[0]);
val = node->Val()[0];
err = node->Val()[1];
return kTRUE;
}
TKDNodeInfo* TKDInterpolatorBase::GetNodeInfo(Int_t inode) const
{
if(!fNodes || inode >= GetNTNodes()) return NULL;
return (TKDNodeInfo*)(*fNodes)[inode];
}
Int_t TKDInterpolatorBase::GetNTNodes() const
{
return fNodes?fNodes->GetEntriesFast():0;
}
Bool_t TKDInterpolatorBase::GetRange(Int_t ax, Float_t &min, Float_t &max) const
{
if(!fNodes) return kFALSE;
Int_t ndim = ((TKDNodeInfo*)(*fNodes)[0])->GetDimension();
if(ax<0 || ax>=ndim){
min=0.; max=0.;
return kFALSE;
}
min=1.e10; max=-1.e10;
Float_t axmin, axmax;
for(Int_t in=GetNTNodes(); in--; ){
TKDNodeInfo *node = (TKDNodeInfo*)((*fNodes)[in]);
node->GetBoundary(ax, axmin, axmax);
if(axmin<min) min = axmin;
if(axmax>max) max = axmax;
}
return kTRUE;
}
void TKDInterpolatorBase::GetStatus(Option_t *opt)
{
printf("Interpolator Status[%d] :\n", fStatus);
printf(" Dim : %d [%d]\n", fNSize, fLambda);
printf(" Method : %s\n", UseCOG() ? "COG" : "INT");
printf(" Store : %s\n", HasStore() ? "YES" : "NO");
printf(" Weights: %s\n", UseWeights() ? "YES" : "NO");
if(strcmp(opt, "all") != 0 ) return;
printf("GetNTNodes() %d\n", GetNTNodes());
for(int i=0; i<GetNTNodes(); i++){
TKDNodeInfo *node = (TKDNodeInfo*)(*fNodes)[i];
printf("%d ", i); node->Print();
}
}
Double_t TKDInterpolatorBase::Eval(const Double_t *point, Double_t &result, Double_t &error, Bool_t force)
{
Float_t pointF[50];
for(int idim=0; idim<fNSize; idim++) pointF[idim] = (Float_t)point[idim];
Int_t nodeIndex = GetNodeIndex(pointF);
if(nodeIndex<0){
Error("TKDInterpolatorBase::Eval()", "Can not retrive node for data point.");
result = 0.;
error = 1.E10;
return 0.;
}
TKDNodeInfo *node = (TKDNodeInfo*)(*fNodes)[nodeIndex];
if(node->Par() && !force){
return node->CookPDF(point, result, error);
}
if(!fBuffer) fBuffer = new Double_t[2*fLambda];
if(!fKDhelper){
fRefPoints = new Float_t*[fNSize];
for(int id=0; id<fNSize; id++){
fRefPoints[id] = new Float_t[GetNTNodes()];
for(int in=0; in<GetNTNodes(); in++) fRefPoints[id][in] = ((TKDNodeInfo*)(*fNodes)[in])->Data()[id];
}
Info("TKDInterpolatorBase::Eval()", "Build TKDTree(%d, %d, %d)", GetNTNodes(), fNSize, kNhelper);
fKDhelper = new TKDTreeIF(GetNTNodes(), fNSize, kNhelper, fRefPoints);
fKDhelper->Build();
fKDhelper->MakeBoundariesExact();
}
if(!fFitter) fFitter = new TLinearFitter(fLambda, Form("hyp%d", fLambda-1));
Int_t ipar,
npoints_new = Int_t((1.+fAlpha)*fLambda),
npoints(0);
Int_t *index = new Int_t[2*npoints_new];
Float_t *dist = new Float_t[2*npoints_new],
d,
w0,
w;
Bool_t kDOWN = kFALSE;
do{
if(npoints){
Info("TKDInterpolatorBase::Eval()", "Interpolation failed. Trying to increase the number of interpolation points from %d to %d.", npoints, npoints_new);
}
if(npoints == npoints_new){
Error("TKDInterpolatorBase::Eval()", "Interpolation failed and number of interpolation points (%d) Can not be increased further.", npoints);
result = 0.;
error = 1.E10;
delete [] dist; delete [] index;
return 0.;
} else npoints = npoints_new;
if(npoints > GetNTNodes()){
Warning("TKDInterpolatorBase::Eval()", "The number of interpolation points requested (%d) exceeds number of PDF values (%d). Downscale.", npoints, GetNTNodes());
npoints = GetNTNodes();
kDOWN = kTRUE;
}
for(int idim=0; idim<fNSize; idim++) pointF[idim] = (Float_t)point[idim];
fKDhelper->FindNearestNeighbors(pointF, npoints+1, index, dist);
fFitter->ClearPoints();
TKDNodeInfo *tnode = NULL;
for(int in=0; in<npoints; in++){
tnode = (TKDNodeInfo*)(*fNodes)[index[in]];
if(UseCOG()){
Float_t *p = &(tnode->Data()[0]);
ipar = 0;
for(int idim=0; idim<fNSize; idim++){
fBuffer[ipar++] = p[idim];
for(int jdim=idim; jdim<fNSize; jdim++) fBuffer[ipar++] = p[idim]*p[jdim];
}
} else {
Float_t *bounds = &(tnode->Data()[fNSize]);
ipar = 0;
for(int idim=0; idim<fNSize; idim++){
fBuffer[ipar++] = .5*(bounds[2*idim] + bounds[2*idim+1]);
fBuffer[ipar++] = (bounds[2*idim]*bounds[2*idim] + bounds[2*idim] * bounds[2*idim+1] + bounds[2*idim+1] * bounds[2*idim+1])/3.;
for(int jdim=idim+1; jdim<fNSize; jdim++) fBuffer[ipar++] = (bounds[2*idim] + bounds[2*idim+1]) * (bounds[2*jdim] + bounds[2*jdim+1]) * .25;
}
}
if(UseWeights()){
d = dist[in]/dist[npoints];
w0 = (1. - d*d*d); w = w0*w0*w0;
if(w<1.e-30) continue;
} else w = 1.;
fFitter->AddPoint(fBuffer, tnode->Val()[0], tnode->Val()[1]/w);
}
npoints_new = npoints+ (kDOWN ? 0 : kdN);
} while(fFitter->Eval());
delete [] index;
delete [] dist;
TMatrixD cov(fLambda, fLambda);
TVectorD par(fLambda);
fFitter->GetCovarianceMatrix(cov);
fFitter->GetParameters(par);
Double_t chi2 = fFitter->GetChisquare()/(npoints - 4 - fLambda);
node->Store(&par, HasStore()?&cov:NULL);
Double_t *fdfdp = &fBuffer[fLambda];
ipar = 0;
fdfdp[ipar++] = 1.;
for(int idim=0; idim<fNSize; idim++){
fdfdp[ipar++] = point[idim];
for(int jdim=idim; jdim<fNSize; jdim++) fdfdp[ipar++] = point[idim]*point[jdim];
}
result =0.; error = 0.;
for(int i=0; i<fLambda; i++){
result += fdfdp[i]*par(i);
for(int j=0; j<fLambda; j++) error += fdfdp[i]*fdfdp[j]*cov(i,j);
}
error = TMath::Sqrt(error);
return chi2;
}
void TKDInterpolatorBase::DrawProjection(UInt_t ax1, UInt_t ax2)
{
Float_t ax1min, ax1max, ax2min, ax2max;
GetRange(ax1, ax1min, ax1max);
GetRange(ax2, ax2min, ax2max);
TH2 *h2 = NULL;
if(!(h2 = (TH2S*)gROOT->FindObject("hKDnodes"))){
h2 = new TH2S("hKDnodes", "", 100, ax1min, ax1max, 100, ax2min, ax2max);
}
h2->GetXaxis()->SetRangeUser(ax1min, ax1max);
h2->GetXaxis()->SetTitle(Form("x_{%d}", ax1));
h2->GetYaxis()->SetRangeUser(ax2min, ax2max);
h2->GetYaxis()->SetTitle(Form("x_{%d}", ax2));
h2->Draw();
if(!fNodesDraw) fNodesDraw = new TKDNodeInfo::TKDNodeDraw[GetNTNodes()];
TKDNodeInfo::TKDNodeDraw *box = NULL;
for(Int_t in=GetNTNodes(); in--; ){
box = &(fNodesDraw[in]);
box->SetNode((TKDNodeInfo*)((*fNodes)[in]), fNSize, ax1, ax2);
box->Draw();
}
return;
}
void TKDInterpolatorBase::SetAlpha(Float_t a)
{
if(a<0.5){
Warning("TKDInterpolatorBase::SetAlpha()", "The scale parameter has to be larger than 0.5");
fAlpha = 0.5;
return;
}
if(Int_t((a+1.)*fLambda) > GetNTNodes()){
fAlpha = TMath::Max(0.5, Float_t(GetNTNodes())/fLambda-1.);
Warning("TKDInterpolatorBase::SetAlpha()", "Interpolation neighborhood exceeds number of evaluation points. Downscale alpha to %f", fAlpha);
return;
}
fAlpha = a;
return;
}
TKDInterpolatorBase.cxx:1 TKDInterpolatorBase.cxx:2 TKDInterpolatorBase.cxx:3 TKDInterpolatorBase.cxx:4 TKDInterpolatorBase.cxx:5 TKDInterpolatorBase.cxx:6 TKDInterpolatorBase.cxx:7 TKDInterpolatorBase.cxx:8 TKDInterpolatorBase.cxx:9 TKDInterpolatorBase.cxx:10 TKDInterpolatorBase.cxx:11 TKDInterpolatorBase.cxx:12 TKDInterpolatorBase.cxx:13 TKDInterpolatorBase.cxx:14 TKDInterpolatorBase.cxx:15 TKDInterpolatorBase.cxx:16 TKDInterpolatorBase.cxx:17 TKDInterpolatorBase.cxx:18 TKDInterpolatorBase.cxx:19 TKDInterpolatorBase.cxx:20 TKDInterpolatorBase.cxx:21 TKDInterpolatorBase.cxx:22 TKDInterpolatorBase.cxx:23 TKDInterpolatorBase.cxx:24 TKDInterpolatorBase.cxx:25 TKDInterpolatorBase.cxx:26 TKDInterpolatorBase.cxx:27 TKDInterpolatorBase.cxx:28 TKDInterpolatorBase.cxx:29 TKDInterpolatorBase.cxx:30 TKDInterpolatorBase.cxx:31 TKDInterpolatorBase.cxx:32 TKDInterpolatorBase.cxx:33 TKDInterpolatorBase.cxx:34 TKDInterpolatorBase.cxx:35 TKDInterpolatorBase.cxx:36 TKDInterpolatorBase.cxx:37 TKDInterpolatorBase.cxx:38 TKDInterpolatorBase.cxx:39 TKDInterpolatorBase.cxx:40 TKDInterpolatorBase.cxx:41 TKDInterpolatorBase.cxx:42 TKDInterpolatorBase.cxx:43 TKDInterpolatorBase.cxx:44 TKDInterpolatorBase.cxx:45 TKDInterpolatorBase.cxx:46 TKDInterpolatorBase.cxx:47 TKDInterpolatorBase.cxx:48 TKDInterpolatorBase.cxx:49 TKDInterpolatorBase.cxx:50 TKDInterpolatorBase.cxx:51 TKDInterpolatorBase.cxx:52 TKDInterpolatorBase.cxx:53 TKDInterpolatorBase.cxx:54 TKDInterpolatorBase.cxx:55 TKDInterpolatorBase.cxx:56 TKDInterpolatorBase.cxx:57 TKDInterpolatorBase.cxx:58 TKDInterpolatorBase.cxx:59 TKDInterpolatorBase.cxx:60 TKDInterpolatorBase.cxx:61 TKDInterpolatorBase.cxx:62 TKDInterpolatorBase.cxx:63 TKDInterpolatorBase.cxx:64 TKDInterpolatorBase.cxx:65 TKDInterpolatorBase.cxx:66 TKDInterpolatorBase.cxx:67 TKDInterpolatorBase.cxx:68 TKDInterpolatorBase.cxx:69 TKDInterpolatorBase.cxx:70 TKDInterpolatorBase.cxx:71 TKDInterpolatorBase.cxx:72 TKDInterpolatorBase.cxx:73 TKDInterpolatorBase.cxx:74 TKDInterpolatorBase.cxx:75 TKDInterpolatorBase.cxx:76 TKDInterpolatorBase.cxx:77 TKDInterpolatorBase.cxx:78 TKDInterpolatorBase.cxx:79 TKDInterpolatorBase.cxx:80 TKDInterpolatorBase.cxx:81 TKDInterpolatorBase.cxx:82 TKDInterpolatorBase.cxx:83 TKDInterpolatorBase.cxx:84 TKDInterpolatorBase.cxx:85 TKDInterpolatorBase.cxx:86 TKDInterpolatorBase.cxx:87 TKDInterpolatorBase.cxx:88 TKDInterpolatorBase.cxx:89 TKDInterpolatorBase.cxx:90 TKDInterpolatorBase.cxx:91 TKDInterpolatorBase.cxx:92 TKDInterpolatorBase.cxx:93 TKDInterpolatorBase.cxx:94 TKDInterpolatorBase.cxx:95 TKDInterpolatorBase.cxx:96 TKDInterpolatorBase.cxx:97 TKDInterpolatorBase.cxx:98 TKDInterpolatorBase.cxx:99 TKDInterpolatorBase.cxx:100 TKDInterpolatorBase.cxx:101 TKDInterpolatorBase.cxx:102 TKDInterpolatorBase.cxx:103 TKDInterpolatorBase.cxx:104 TKDInterpolatorBase.cxx:105 TKDInterpolatorBase.cxx:106 TKDInterpolatorBase.cxx:107 TKDInterpolatorBase.cxx:108 TKDInterpolatorBase.cxx:109 TKDInterpolatorBase.cxx:110 TKDInterpolatorBase.cxx:111 TKDInterpolatorBase.cxx:112 TKDInterpolatorBase.cxx:113 TKDInterpolatorBase.cxx:114 TKDInterpolatorBase.cxx:115 TKDInterpolatorBase.cxx:116 TKDInterpolatorBase.cxx:117 TKDInterpolatorBase.cxx:118 TKDInterpolatorBase.cxx:119 TKDInterpolatorBase.cxx:120 TKDInterpolatorBase.cxx:121 TKDInterpolatorBase.cxx:122 TKDInterpolatorBase.cxx:123 TKDInterpolatorBase.cxx:124 TKDInterpolatorBase.cxx:125 TKDInterpolatorBase.cxx:126 TKDInterpolatorBase.cxx:127 TKDInterpolatorBase.cxx:128 TKDInterpolatorBase.cxx:129 TKDInterpolatorBase.cxx:130 TKDInterpolatorBase.cxx:131 TKDInterpolatorBase.cxx:132 TKDInterpolatorBase.cxx:133 TKDInterpolatorBase.cxx:134 TKDInterpolatorBase.cxx:135 TKDInterpolatorBase.cxx:136 TKDInterpolatorBase.cxx:137 TKDInterpolatorBase.cxx:138 TKDInterpolatorBase.cxx:139 TKDInterpolatorBase.cxx:140 TKDInterpolatorBase.cxx:141 TKDInterpolatorBase.cxx:142 TKDInterpolatorBase.cxx:143 TKDInterpolatorBase.cxx:144 TKDInterpolatorBase.cxx:145 TKDInterpolatorBase.cxx:146 TKDInterpolatorBase.cxx:147 TKDInterpolatorBase.cxx:148 TKDInterpolatorBase.cxx:149 TKDInterpolatorBase.cxx:150 TKDInterpolatorBase.cxx:151 TKDInterpolatorBase.cxx:152 TKDInterpolatorBase.cxx:153 TKDInterpolatorBase.cxx:154 TKDInterpolatorBase.cxx:155 TKDInterpolatorBase.cxx:156 TKDInterpolatorBase.cxx:157 TKDInterpolatorBase.cxx:158 TKDInterpolatorBase.cxx:159 TKDInterpolatorBase.cxx:160 TKDInterpolatorBase.cxx:161 TKDInterpolatorBase.cxx:162 TKDInterpolatorBase.cxx:163 TKDInterpolatorBase.cxx:164 TKDInterpolatorBase.cxx:165 TKDInterpolatorBase.cxx:166 TKDInterpolatorBase.cxx:167 TKDInterpolatorBase.cxx:168 TKDInterpolatorBase.cxx:169 TKDInterpolatorBase.cxx:170 TKDInterpolatorBase.cxx:171 TKDInterpolatorBase.cxx:172 TKDInterpolatorBase.cxx:173 TKDInterpolatorBase.cxx:174 TKDInterpolatorBase.cxx:175 TKDInterpolatorBase.cxx:176 TKDInterpolatorBase.cxx:177 TKDInterpolatorBase.cxx:178 TKDInterpolatorBase.cxx:179 TKDInterpolatorBase.cxx:180 TKDInterpolatorBase.cxx:181 TKDInterpolatorBase.cxx:182 TKDInterpolatorBase.cxx:183 TKDInterpolatorBase.cxx:184 TKDInterpolatorBase.cxx:185 TKDInterpolatorBase.cxx:186 TKDInterpolatorBase.cxx:187 TKDInterpolatorBase.cxx:188 TKDInterpolatorBase.cxx:189 TKDInterpolatorBase.cxx:190 TKDInterpolatorBase.cxx:191 TKDInterpolatorBase.cxx:192 TKDInterpolatorBase.cxx:193 TKDInterpolatorBase.cxx:194 TKDInterpolatorBase.cxx:195 TKDInterpolatorBase.cxx:196 TKDInterpolatorBase.cxx:197 TKDInterpolatorBase.cxx:198 TKDInterpolatorBase.cxx:199 TKDInterpolatorBase.cxx:200 TKDInterpolatorBase.cxx:201 TKDInterpolatorBase.cxx:202 TKDInterpolatorBase.cxx:203 TKDInterpolatorBase.cxx:204 TKDInterpolatorBase.cxx:205 TKDInterpolatorBase.cxx:206 TKDInterpolatorBase.cxx:207 TKDInterpolatorBase.cxx:208 TKDInterpolatorBase.cxx:209 TKDInterpolatorBase.cxx:210 TKDInterpolatorBase.cxx:211 TKDInterpolatorBase.cxx:212 TKDInterpolatorBase.cxx:213 TKDInterpolatorBase.cxx:214 TKDInterpolatorBase.cxx:215 TKDInterpolatorBase.cxx:216 TKDInterpolatorBase.cxx:217 TKDInterpolatorBase.cxx:218 TKDInterpolatorBase.cxx:219 TKDInterpolatorBase.cxx:220 TKDInterpolatorBase.cxx:221 TKDInterpolatorBase.cxx:222 TKDInterpolatorBase.cxx:223 TKDInterpolatorBase.cxx:224 TKDInterpolatorBase.cxx:225 TKDInterpolatorBase.cxx:226 TKDInterpolatorBase.cxx:227 TKDInterpolatorBase.cxx:228 TKDInterpolatorBase.cxx:229 TKDInterpolatorBase.cxx:230 TKDInterpolatorBase.cxx:231 TKDInterpolatorBase.cxx:232 TKDInterpolatorBase.cxx:233 TKDInterpolatorBase.cxx:234 TKDInterpolatorBase.cxx:235 TKDInterpolatorBase.cxx:236 TKDInterpolatorBase.cxx:237 TKDInterpolatorBase.cxx:238 TKDInterpolatorBase.cxx:239 TKDInterpolatorBase.cxx:240 TKDInterpolatorBase.cxx:241 TKDInterpolatorBase.cxx:242 TKDInterpolatorBase.cxx:243 TKDInterpolatorBase.cxx:244 TKDInterpolatorBase.cxx:245 TKDInterpolatorBase.cxx:246 TKDInterpolatorBase.cxx:247 TKDInterpolatorBase.cxx:248 TKDInterpolatorBase.cxx:249 TKDInterpolatorBase.cxx:250 TKDInterpolatorBase.cxx:251 TKDInterpolatorBase.cxx:252 TKDInterpolatorBase.cxx:253 TKDInterpolatorBase.cxx:254 TKDInterpolatorBase.cxx:255 TKDInterpolatorBase.cxx:256 TKDInterpolatorBase.cxx:257 TKDInterpolatorBase.cxx:258 TKDInterpolatorBase.cxx:259 TKDInterpolatorBase.cxx:260 TKDInterpolatorBase.cxx:261 TKDInterpolatorBase.cxx:262 TKDInterpolatorBase.cxx:263 TKDInterpolatorBase.cxx:264 TKDInterpolatorBase.cxx:265 TKDInterpolatorBase.cxx:266 TKDInterpolatorBase.cxx:267 TKDInterpolatorBase.cxx:268 TKDInterpolatorBase.cxx:269 TKDInterpolatorBase.cxx:270 TKDInterpolatorBase.cxx:271 TKDInterpolatorBase.cxx:272 TKDInterpolatorBase.cxx:273 TKDInterpolatorBase.cxx:274 TKDInterpolatorBase.cxx:275 TKDInterpolatorBase.cxx:276 TKDInterpolatorBase.cxx:277 TKDInterpolatorBase.cxx:278 TKDInterpolatorBase.cxx:279 TKDInterpolatorBase.cxx:280 TKDInterpolatorBase.cxx:281 TKDInterpolatorBase.cxx:282 TKDInterpolatorBase.cxx:283 TKDInterpolatorBase.cxx:284 TKDInterpolatorBase.cxx:285 TKDInterpolatorBase.cxx:286 TKDInterpolatorBase.cxx:287 TKDInterpolatorBase.cxx:288 TKDInterpolatorBase.cxx:289 TKDInterpolatorBase.cxx:290 TKDInterpolatorBase.cxx:291 TKDInterpolatorBase.cxx:292 TKDInterpolatorBase.cxx:293 TKDInterpolatorBase.cxx:294 TKDInterpolatorBase.cxx:295 TKDInterpolatorBase.cxx:296 TKDInterpolatorBase.cxx:297 TKDInterpolatorBase.cxx:298 TKDInterpolatorBase.cxx:299 TKDInterpolatorBase.cxx:300 TKDInterpolatorBase.cxx:301 TKDInterpolatorBase.cxx:302 TKDInterpolatorBase.cxx:303 TKDInterpolatorBase.cxx:304 TKDInterpolatorBase.cxx:305 TKDInterpolatorBase.cxx:306 TKDInterpolatorBase.cxx:307 TKDInterpolatorBase.cxx:308 TKDInterpolatorBase.cxx:309 TKDInterpolatorBase.cxx:310 TKDInterpolatorBase.cxx:311 TKDInterpolatorBase.cxx:312 TKDInterpolatorBase.cxx:313 TKDInterpolatorBase.cxx:314 TKDInterpolatorBase.cxx:315 TKDInterpolatorBase.cxx:316 TKDInterpolatorBase.cxx:317 TKDInterpolatorBase.cxx:318 TKDInterpolatorBase.cxx:319 TKDInterpolatorBase.cxx:320 TKDInterpolatorBase.cxx:321 TKDInterpolatorBase.cxx:322 TKDInterpolatorBase.cxx:323 TKDInterpolatorBase.cxx:324 TKDInterpolatorBase.cxx:325 TKDInterpolatorBase.cxx:326 TKDInterpolatorBase.cxx:327 TKDInterpolatorBase.cxx:328 TKDInterpolatorBase.cxx:329 TKDInterpolatorBase.cxx:330 TKDInterpolatorBase.cxx:331 TKDInterpolatorBase.cxx:332 TKDInterpolatorBase.cxx:333 TKDInterpolatorBase.cxx:334 TKDInterpolatorBase.cxx:335 TKDInterpolatorBase.cxx:336 TKDInterpolatorBase.cxx:337 TKDInterpolatorBase.cxx:338 TKDInterpolatorBase.cxx:339 TKDInterpolatorBase.cxx:340 TKDInterpolatorBase.cxx:341 TKDInterpolatorBase.cxx:342 TKDInterpolatorBase.cxx:343 TKDInterpolatorBase.cxx:344 TKDInterpolatorBase.cxx:345 TKDInterpolatorBase.cxx:346 TKDInterpolatorBase.cxx:347 TKDInterpolatorBase.cxx:348 TKDInterpolatorBase.cxx:349 TKDInterpolatorBase.cxx:350 TKDInterpolatorBase.cxx:351 TKDInterpolatorBase.cxx:352 TKDInterpolatorBase.cxx:353 TKDInterpolatorBase.cxx:354 TKDInterpolatorBase.cxx:355 TKDInterpolatorBase.cxx:356 TKDInterpolatorBase.cxx:357 TKDInterpolatorBase.cxx:358 TKDInterpolatorBase.cxx:359 TKDInterpolatorBase.cxx:360 TKDInterpolatorBase.cxx:361 TKDInterpolatorBase.cxx:362 TKDInterpolatorBase.cxx:363 TKDInterpolatorBase.cxx:364 TKDInterpolatorBase.cxx:365 TKDInterpolatorBase.cxx:366 TKDInterpolatorBase.cxx:367 TKDInterpolatorBase.cxx:368 TKDInterpolatorBase.cxx:369 TKDInterpolatorBase.cxx:370 TKDInterpolatorBase.cxx:371 TKDInterpolatorBase.cxx:372 TKDInterpolatorBase.cxx:373 TKDInterpolatorBase.cxx:374 TKDInterpolatorBase.cxx:375 TKDInterpolatorBase.cxx:376 TKDInterpolatorBase.cxx:377