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 fData = new Float_t[fNDim];
60 memcpy(fData, ref.fData, fNDim*sizeof(Float_t));
61 fVal[0] = ref.fVal[0];
62 fVal[1] = ref.fVal[1];
63 if(ref.fNpar&&ref.fPar){
65 fPar=new Double_t[fNpar];
66 memcpy(fPar, ref.fPar, fNpar*sizeof(Double_t));
68 if(ref.fNcov && ref.fCov){
70 fCov=new Double_t[fNcov];
71 memcpy(fCov, ref.fCov, fNcov*sizeof(Double_t));
76 //_________________________________________________________________
77 TKDNodeInfo::~TKDNodeInfo()
80 if(fData) delete [] fData;
81 if(fPar) delete [] fPar;
82 if(fCov) delete [] fCov;
85 //_________________________________________________________________
86 TKDNodeInfo& TKDNodeInfo::operator=(const TKDNodeInfo & ref)
88 // Info("operator==()", "...");
90 if(this == &ref) return *this;
92 if(fNDim != ref.fNDim){
96 memcpy(fData, ref.fData, fNDim*sizeof(Float_t));
97 fVal[0] = ref.fVal[0];
98 fVal[1] = ref.fVal[1];
99 if(ref.fNpar&&ref.fPar){
101 fPar=new Double_t[fNpar];
102 memcpy(fPar, ref.fPar, fNpar*sizeof(Double_t));
104 if(ref.fNcov && ref.fCov){
106 fCov=new Double_t[fNcov];
107 memcpy(fCov, ref.fCov, fNcov*sizeof(Double_t));
112 //_________________________________________________________________
113 void TKDNodeInfo::Build(Int_t dim)
115 // Allocate/Reallocate space for this node.
117 // Info("Build()", "...");
121 if(fData) delete [] fData;
122 fData = new Float_t[fNDim];
126 //_________________________________________________________________
127 void TKDNodeInfo::Bootstrap()
129 if(!fNpar || !fPar) return;
131 Int_t ndim = Int_t(.5*(TMath::Sqrt(1.+8.*fNpar)-1.))-1;
135 //_________________________________________________________________
136 void TKDNodeInfo::SetNode(Int_t ndim, Float_t *data, Float_t *pdf)
139 memcpy(fData, data, fNDim*sizeof(Float_t));
140 fVal[0]=pdf[0]; fVal[1]=pdf[1];
144 //_________________________________________________________________
145 void TKDNodeInfo::Print(const Option_t *opt) const
147 // Print the content of the node
148 Int_t dim = Int_t(fNDim/3.);
150 for(int idim=0; idim<dim; idim++) printf("%f ", fData?fData[idim]:0.);
151 printf("] f = [%f +- %f]\n", fVal[0], fVal[1]);
154 Float_t *bounds = &fData[dim];
156 for(int idim=0; idim<dim; idim++) printf("(%f %f) ", bounds[2*idim], bounds[2*idim+1]);
159 if(strcmp(opt, "a")!=0) return;
162 printf("Fit parameters : \n");
163 for(int ip=0; ip<fNpar; ip++) printf("p%d[%f] ", ip, fPar[ip]);
167 for(int ip(0), n(0); ip<fNpar; ip++){
168 for(int jp(ip); jp<fNpar; jp++) printf("c(%d %d)[%f] ", ip, jp, fCov[n++]);
173 //_________________________________________________________________
174 void TKDNodeInfo::Store(TVectorD const *par, TMatrixD const *cov)
176 // Store the parameters and the covariance in the node
178 if(!fPar){SetNpar(); fPar = new Double_t[fNpar];}
179 for(int ip=0; ip<fNpar; ip++) fPar[ip] = (*par)[ip];
182 if(!fCov){SetNcov(); fCov = new Double_t[fNcov];}
183 for(int ip(0), np(0); ip<fNpar; ip++)
184 for(int jp=ip; jp<fNpar; jp++) fCov[np++] = (*cov)(ip,jp);
187 //_________________________________________________________________
188 Bool_t TKDNodeInfo::CookPDF(const Double_t *point, Double_t &result, Double_t &error) const
190 // Recalculate the PDF for one node from the results of interpolation (parameters and covariance matrix)
192 Int_t ndim = Int_t(fNDim/3.);
193 if(ndim>10) return kFALSE; // support only up to 10 dimensions
194 //printf("ndim[%d] npar[%d] ncov[%d]\n", ndim, fNpar, fNcov);
196 Double_t fdfdp[66]; memset(fdfdp, 0, ndim*sizeof(Double_t));
199 for(int idim=0; idim<ndim; idim++){
200 fdfdp[ipar++] = point[idim];
201 for(int jdim=idim; jdim<ndim; jdim++) fdfdp[ipar++] = point[idim]*point[jdim];
204 // calculate estimation
205 result =0.; error = 0.;
206 for(int i=0; i<fNpar; i++) result += fdfdp[i]*fPar[i];
207 if(!fNcov) return kTRUE;
209 for(int i(0), n(0); i<fNpar; i++){
210 error += fdfdp[i]*fdfdp[i]*fCov[n++];
211 for(int j(i+1); j<fNpar; j++) error += 2.*fdfdp[i]*fdfdp[j]*fCov[n++];
213 error = TMath::Sqrt(error);
215 //printf("TKDNodeInfo::CookPDF() : %6.3f +- %6.3f\n", result, error);
222 //_________________________________________________________________
223 TKDNodeInfo::TKDNodeDraw::TKDNodeDraw()
229 SetFillColor(50+Int_t(gRandom->Uniform()*50.));
231 fCOG.SetMarkerStyle(3);
232 fCOG.SetMarkerSize(.7);
233 fCOG.SetMarkerColor(2);
237 //_________________________________________________________________
238 void TKDNodeInfo::TKDNodeDraw::Draw(Option_t* option)
244 //_________________________________________________________________
245 void TKDNodeInfo::TKDNodeDraw::SetNode(TKDNodeInfo *node, UChar_t size, UChar_t ax1, UChar_t ax2)
248 const Float_t kBorder = 0.;//1.E-4;
249 Float_t *bounds = &(node->Data()[size]);
250 fX1=bounds[2*ax1]+kBorder;
251 fX2=bounds[2*ax1+1]-kBorder;
252 fY1=bounds[2*ax2]+kBorder;
253 fY2=bounds[2*ax2+1]-kBorder;
255 Float_t x(node->Data()[ax1]), y(node->Data()[ax2]);
256 fCOG.SetX(x); fCOG.SetY(y);
260 //_________________________________________________________________
261 void TKDNodeInfo::TKDNodeDraw::Print(const Option_t* option) const
264 fNode->Print(option);