]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STAT/TKDNodeInfo.cxx
Adding missing object
[u/mrichter/AliRoot.git] / STAT / TKDNodeInfo.cxx
1 ////////////////////////////////////////////////////////
2 //
3 // Bucket representation for TKDInterpolator(Base)
4 //
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)
9 //   - boundaries
10 // and interpolated data like
11 //   - parameters of the local parabolic fit
12 //   - their covariance matrix
13 //  
14 // For drawing 2D projections the helper class TKDNodeInfo::TKDNodeDraw is used.
15 //
16 // Author Alexandru Bercuci <A.Bercuci@gsi.de>
17 //
18 ////////////////////////////////////////////////////////
19
20 #include "TKDNodeInfo.h"
21
22 #include "TVectorD.h"
23 #include "TMatrixD.h"
24 #include "TRandom.h"
25 #include "TMath.h"
26
27 ClassImp(TKDNodeInfo)
28 ClassImp(TKDNodeInfo::TKDNodeDraw)
29
30
31 //_________________________________________________________________
32 TKDNodeInfo::TKDNodeInfo(Int_t dim):
33   TObject()
34   ,fNDim(3*dim)
35   ,fData(NULL)
36   ,fNpar(0)
37   ,fNcov(0)
38   ,fPar(NULL)
39   ,fCov(NULL)
40 {
41   // Default constructor
42   fVal[0] = 0.; fVal[1] = 0.;
43   Build(dim);
44 }
45
46 //_________________________________________________________________
47 TKDNodeInfo::TKDNodeInfo(const TKDNodeInfo &ref):
48   TObject((TObject&) ref)
49   ,fNDim(fNDim)
50   ,fData(NULL)
51   ,fNpar(0)
52   ,fNcov(0)
53   ,fPar(NULL)
54   ,fCov(NULL)
55 {
56   // Copy constructor
57   Build(fNDim/3);
58
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){ 
63     fNpar = ref.fNpar;
64     fPar=new Double_t[fNpar];
65     memcpy(fPar, ref.fPar, fNpar*sizeof(Double_t));
66   }
67   if(ref.fNcov && ref.fCov){ 
68     fNcov = ref.fNcov;
69     fCov=new Double_t[fNcov];
70     memcpy(fCov, ref.fCov, fNcov*sizeof(Double_t));
71   }
72 }
73
74
75 //_________________________________________________________________
76 TKDNodeInfo::~TKDNodeInfo()
77 {
78   // Destructor
79   if(fData) delete [] fData;
80   if(fPar) delete [] fPar;
81   if(fCov) delete [] fCov;
82 }
83
84 //_________________________________________________________________
85 TKDNodeInfo& TKDNodeInfo::operator=(const TKDNodeInfo & ref)
86 {
87 //      Info("operator==()", "...");
88   
89   Int_t ndim = fNDim/3;
90   if(fNDim != ref.fNDim){
91     fNDim = ref.fNDim;
92     Build(ndim);
93   }
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){ 
98     fNpar = ref.fNpar;
99     fPar=new Double_t[fNpar];
100     memcpy(fPar, ref.fPar, fNpar*sizeof(Double_t));
101   }
102   if(ref.fNcov && ref.fCov){ 
103     fNcov = ref.fNcov;
104     fCov=new Double_t[fNcov];
105     memcpy(fCov, ref.fCov, fNcov*sizeof(Double_t));
106   }
107   return *this;
108 }
109
110 //_________________________________________________________________
111 void TKDNodeInfo::Build(Int_t dim)
112 {
113 // Allocate/Reallocate space for this node.
114
115 //      Info("Build()", "...");
116
117   if(!dim) return;
118   fNDim = 3*dim;
119   if(fData) delete [] fData;
120   fData = new Float_t[fNDim];
121   return;
122 }
123
124 //_________________________________________________________________
125 void TKDNodeInfo::Bootstrap()
126 {
127   if(!fNpar || !fPar) return;
128
129   Int_t ndim = Int_t(.5*(TMath::Sqrt(1.+8.*fNpar)-1.))-1;
130   fNDim = 3*ndim;
131 }
132
133 //_________________________________________________________________
134 void TKDNodeInfo::SetNode(Int_t ndim, Float_t *data, Float_t *pdf)
135 {
136   Build(ndim);
137   memcpy(fData, data, fNDim*sizeof(Float_t));
138   fVal[0]=pdf[0]; fVal[1]=pdf[1];
139 }
140
141
142 //_________________________________________________________________
143 void TKDNodeInfo::Print(const Option_t *opt) const
144 {
145   // Print the content of the node
146   Int_t dim = Int_t(fNDim/3.);
147   printf("x[");
148   for(int idim=0; idim<dim; idim++) printf("%f ", fData?fData[idim]:0.);
149   printf("] f = [%f +- %f]\n", fVal[0], fVal[1]);
150
151   if(fData){
152     Float_t *bounds = &fData[dim];
153     printf("range[");
154     for(int idim=0; idim<dim; idim++) printf("(%f %f) ", bounds[2*idim], bounds[2*idim+1]);
155     printf("]\n");
156   }
157   if(strcmp(opt, "a")!=0) return;
158
159   if(fNpar){ 
160     printf("Fit parameters : \n");
161     for(int ip=0; ip<fNpar; ip++) printf("p%d[%f] ", ip, fPar[ip]);
162     printf("\n");
163   }
164   if(!fNcov) return;
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++]);
167     printf("\n");
168   }
169 }
170
171 //_________________________________________________________________
172 void TKDNodeInfo::Store(TVectorD const *par, TMatrixD const *cov)
173 {
174 // Store the parameters and the covariance in the node
175
176   if(!fPar){SetNpar(); fPar = new Double_t[fNpar];}
177   for(int ip=0; ip<fNpar; ip++) fPar[ip] = (*par)[ip];
178
179   if(!cov) return;
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);
183 }
184
185 //_________________________________________________________________
186 Bool_t TKDNodeInfo::CookPDF(const Double_t *point, Double_t &result, Double_t &error) const
187 {
188 // Recalculate the PDF for one node from the results of interpolation (parameters and covariance matrix)
189
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);
193
194   Double_t fdfdp[66]; memset(fdfdp, 0, ndim*sizeof(Double_t));
195   Int_t ipar = 0;
196   fdfdp[ipar++] = 1.;
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];
200   }
201
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;
206
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++];
210   }     
211   error = TMath::Sqrt(error);
212   
213   //printf("TKDNodeInfo::CookPDF() : %6.3f +- %6.3f\n", result, error);
214
215   return kTRUE;
216 }
217
218
219
220 //_________________________________________________________________
221 TKDNodeInfo::TKDNodeDraw::TKDNodeDraw() 
222   :TBox()
223   ,fCOG()
224   ,fNode(NULL)
225 {
226   SetFillStyle(3002);
227   SetFillColor(50+Int_t(gRandom->Uniform()*50.));
228
229   fCOG.SetMarkerStyle(3);
230   fCOG.SetMarkerSize(.7);
231   fCOG.SetMarkerColor(2);
232 }
233
234
235 //_________________________________________________________________
236 void TKDNodeInfo::TKDNodeDraw::Draw(Option_t* option)
237 {
238   TBox::Draw(option);
239   fCOG.Draw("p");
240 }
241
242 //_________________________________________________________________
243 void TKDNodeInfo::TKDNodeDraw::SetNode(TKDNodeInfo *node, UChar_t size, UChar_t ax1, UChar_t ax2)
244 {
245   fNode=node;
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;
252   
253   Float_t x(node->Data()[ax1]), y(node->Data()[ax2]);
254   fCOG.SetX(x); fCOG.SetY(y);
255 }
256
257
258 //_________________________________________________________________
259 void TKDNodeInfo::TKDNodeDraw::Print(const Option_t* option) const
260 {
261   if(!fNode) return;
262   fNode->Print(option);
263 }