1 ////////////////////////////////////////////////////////
3 // Bucket representation for TKDInterpolator(Base)
5 // The class store data and provides the interface to draw and print.
6 // The bucket - generalized histogram bin in N dimensions is represented by unprocessed data like
7 // - experimental PDF value and statistical error
8 // - COG position (n-tuplu)
10 // and interpolated data like
11 // - parameters of the local parabolic fit
12 // - their covariance matrix
14 // For drawing 2D projections the helper class TKDNodeInfo::TKDNodeDraw is used.
16 // Author Alexandru Bercuci <A.Bercuci@gsi.de>
18 ////////////////////////////////////////////////////////
20 #include "TKDNodeInfo.h"
28 ClassImp(TKDNodeInfo::TKDNodeDraw)
31 //_________________________________________________________________
32 TKDNodeInfo::TKDNodeInfo(Int_t dim):
41 // Default constructor
42 fVal[0] = 0.; fVal[1] = 0.;
46 //_________________________________________________________________
47 TKDNodeInfo::TKDNodeInfo(const TKDNodeInfo &ref):
48 TObject((TObject&) ref)
59 memcpy(fData, ref.fData, fNDim*sizeof(Float_t));
60 fVal[0] = ref.fVal[0];
61 fVal[1] = ref.fVal[1];
62 if(ref.fNpar&&ref.fPar){
64 fPar=new Double_t[fNpar];
65 memcpy(fPar, ref.fPar, fNpar*sizeof(Double_t));
67 if(ref.fNcov && ref.fCov){
69 fCov=new Double_t[fNcov];
70 memcpy(fCov, ref.fCov, fNcov*sizeof(Double_t));
75 //_________________________________________________________________
76 TKDNodeInfo::~TKDNodeInfo()
79 if(fData) delete [] fData;
80 if(fPar) delete [] fPar;
81 if(fCov) delete [] fCov;
84 //_________________________________________________________________
85 TKDNodeInfo& TKDNodeInfo::operator=(const TKDNodeInfo & ref)
87 // Info("operator==()", "...");
90 if(fNDim != ref.fNDim){
94 memcpy(fData, ref.fData, fNDim*sizeof(Float_t));
95 fVal[0] = ref.fVal[0];
96 fVal[1] = ref.fVal[1];
97 if(ref.fNpar&&ref.fPar){
99 fPar=new Double_t[fNpar];
100 memcpy(fPar, ref.fPar, fNpar*sizeof(Double_t));
102 if(ref.fNcov && ref.fCov){
104 fCov=new Double_t[fNcov];
105 memcpy(fCov, ref.fCov, fNcov*sizeof(Double_t));
110 //_________________________________________________________________
111 void TKDNodeInfo::Build(Int_t dim)
113 // Allocate/Reallocate space for this node.
115 // Info("Build()", "...");
119 if(fData) delete [] fData;
120 fData = new Float_t[fNDim];
124 //_________________________________________________________________
125 void TKDNodeInfo::Bootstrap()
127 if(!fNpar || !fPar) return;
129 Int_t ndim = Int_t(.5*(TMath::Sqrt(1.+8.*fNpar)-1.))-1;
133 //_________________________________________________________________
134 void TKDNodeInfo::SetNode(Int_t ndim, Float_t *data, Float_t *pdf)
137 memcpy(fData, data, fNDim*sizeof(Float_t));
138 fVal[0]=pdf[0]; fVal[1]=pdf[1];
142 //_________________________________________________________________
143 void TKDNodeInfo::Print(const Option_t *opt) const
145 // Print the content of the node
146 Int_t dim = Int_t(fNDim/3.);
148 for(int idim=0; idim<dim; idim++) printf("%f ", fData?fData[idim]:0.);
149 printf("] f = [%f +- %f]\n", fVal[0], fVal[1]);
152 Float_t *bounds = &fData[dim];
154 for(int idim=0; idim<dim; idim++) printf("(%f %f) ", bounds[2*idim], bounds[2*idim+1]);
157 if(strcmp(opt, "a")!=0) return;
160 printf("Fit parameters : \n");
161 for(int ip=0; ip<fNpar; ip++) printf("p%d[%f] ", ip, fPar[ip]);
165 for(int ip(0), n(0); ip<fNpar; ip++){
166 for(int jp(ip); jp<fNpar; jp++) printf("c(%d %d)[%f] ", ip, jp, fCov[n++]);
171 //_________________________________________________________________
172 void TKDNodeInfo::Store(TVectorD const *par, TMatrixD const *cov)
174 // Store the parameters and the covariance in the node
176 if(!fPar){SetNpar(); fPar = new Double_t[fNpar];}
177 for(int ip=0; ip<fNpar; ip++) fPar[ip] = (*par)[ip];
180 if(!fCov){SetNcov(); fCov = new Double_t[fNcov];}
181 for(int ip(0), np(0); ip<fNpar; ip++)
182 for(int jp=ip; jp<fNpar; jp++) fCov[np++] = (*cov)(ip,jp);
185 //_________________________________________________________________
186 Bool_t TKDNodeInfo::CookPDF(const Double_t *point, Double_t &result, Double_t &error) const
188 // Recalculate the PDF for one node from the results of interpolation (parameters and covariance matrix)
190 Int_t ndim = Int_t(fNDim/3.);
191 if(ndim>10) return kFALSE; // support only up to 10 dimensions
192 //printf("ndim[%d] npar[%d] ncov[%d]\n", ndim, fNpar, fNcov);
194 Double_t fdfdp[66]; memset(fdfdp, 0, ndim*sizeof(Double_t));
197 for(int idim=0; idim<ndim; idim++){
198 fdfdp[ipar++] = point[idim];
199 for(int jdim=idim; jdim<ndim; jdim++) fdfdp[ipar++] = point[idim]*point[jdim];
202 // calculate estimation
203 result =0.; error = 0.;
204 for(int i=0; i<fNpar; i++) result += fdfdp[i]*fPar[i];
205 if(!fNcov) return kTRUE;
207 for(int i(0), n(0); i<fNpar; i++){
208 error += fdfdp[i]*fdfdp[i]*fCov[n++];
209 for(int j(i+1); j<fNpar; j++) error += 2.*fdfdp[i]*fdfdp[j]*fCov[n++];
211 error = TMath::Sqrt(error);
213 //printf("TKDNodeInfo::CookPDF() : %6.3f +- %6.3f\n", result, error);
220 //_________________________________________________________________
221 TKDNodeInfo::TKDNodeDraw::TKDNodeDraw()
227 SetFillColor(50+Int_t(gRandom->Uniform()*50.));
229 fCOG.SetMarkerStyle(3);
230 fCOG.SetMarkerSize(.7);
231 fCOG.SetMarkerColor(2);
235 //_________________________________________________________________
236 void TKDNodeInfo::TKDNodeDraw::Draw(Option_t* option)
242 //_________________________________________________________________
243 void TKDNodeInfo::TKDNodeDraw::SetNode(TKDNodeInfo *node, UChar_t size, UChar_t ax1, UChar_t ax2)
246 const Float_t kBorder = 0.;//1.E-4;
247 Float_t *bounds = &(node->Data()[size]);
248 fX1=bounds[2*ax1]+kBorder;
249 fX2=bounds[2*ax1+1]-kBorder;
250 fY1=bounds[2*ax2]+kBorder;
251 fY2=bounds[2*ax2+1]-kBorder;
253 Float_t x(node->Data()[ax1]), y(node->Data()[ax2]);
254 fCOG.SetX(x); fCOG.SetY(y);
258 //_________________________________________________________________
259 void TKDNodeInfo::TKDNodeDraw::Print(const Option_t* option) const
262 fNode->Print(option);