]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STAT/TKDNodeInfo.cxx
1a1fa7494f0d5f523ed9744c8afa3fe9e3b7d8d2
[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(ref.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   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){ 
64     fNpar = ref.fNpar;
65     fPar=new Double_t[fNpar];
66     memcpy(fPar, ref.fPar, fNpar*sizeof(Double_t));
67   }
68   if(ref.fNcov && ref.fCov){ 
69     fNcov = ref.fNcov;
70     fCov=new Double_t[fNcov];
71     memcpy(fCov, ref.fCov, fNcov*sizeof(Double_t));
72   }
73 }
74
75
76 //_________________________________________________________________
77 TKDNodeInfo::~TKDNodeInfo()
78 {
79   // Destructor
80   if(fData) delete [] fData;
81   if(fPar) delete [] fPar;
82   if(fCov) delete [] fCov;
83 }
84
85 //_________________________________________________________________
86 TKDNodeInfo& TKDNodeInfo::operator=(const TKDNodeInfo & ref)
87 {
88 //      Info("operator==()", "...");
89   
90   Int_t ndim = fNDim/3;
91   if(fNDim != ref.fNDim){
92     fNDim = ref.fNDim;
93     Build(ndim);
94   }
95   memcpy(fData, ref.fData, fNDim*sizeof(Float_t));
96   fVal[0] = ref.fVal[0];
97   fVal[1] = ref.fVal[1];
98   if(ref.fNpar&&ref.fPar){ 
99     fNpar = ref.fNpar;
100     fPar=new Double_t[fNpar];
101     memcpy(fPar, ref.fPar, fNpar*sizeof(Double_t));
102   }
103   if(ref.fNcov && ref.fCov){ 
104     fNcov = ref.fNcov;
105     fCov=new Double_t[fNcov];
106     memcpy(fCov, ref.fCov, fNcov*sizeof(Double_t));
107   }
108   return *this;
109 }
110
111 //_________________________________________________________________
112 void TKDNodeInfo::Build(Int_t dim)
113 {
114 // Allocate/Reallocate space for this node.
115
116 //      Info("Build()", "...");
117
118   if(!dim) return;
119   fNDim = 3*dim;
120   if(fData) delete [] fData;
121   fData = new Float_t[fNDim];
122   return;
123 }
124
125 //_________________________________________________________________
126 void TKDNodeInfo::Bootstrap()
127 {
128   if(!fNpar || !fPar) return;
129
130   Int_t ndim = Int_t(.5*(TMath::Sqrt(1.+8.*fNpar)-1.))-1;
131   fNDim = 3*ndim;
132 }
133
134 //_________________________________________________________________
135 void TKDNodeInfo::SetNode(Int_t ndim, Float_t *data, Float_t *pdf)
136 {
137   Build(ndim);
138   memcpy(fData, data, fNDim*sizeof(Float_t));
139   fVal[0]=pdf[0]; fVal[1]=pdf[1];
140 }
141
142
143 //_________________________________________________________________
144 void TKDNodeInfo::Print(const Option_t *opt) const
145 {
146   // Print the content of the node
147   Int_t dim = Int_t(fNDim/3.);
148   printf("x[");
149   for(int idim=0; idim<dim; idim++) printf("%f ", fData?fData[idim]:0.);
150   printf("] f = [%f +- %f]\n", fVal[0], fVal[1]);
151
152   if(fData){
153     Float_t *bounds = &fData[dim];
154     printf("range[");
155     for(int idim=0; idim<dim; idim++) printf("(%f %f) ", bounds[2*idim], bounds[2*idim+1]);
156     printf("]\n");
157   }
158   if(strcmp(opt, "a")!=0) return;
159
160   if(fNpar){ 
161     printf("Fit parameters : \n");
162     for(int ip=0; ip<fNpar; ip++) printf("p%d[%f] ", ip, fPar[ip]);
163     printf("\n");
164   }
165   if(!fNcov) return;
166   for(int ip(0), n(0); ip<fNpar; ip++){
167     for(int jp(ip); jp<fNpar; jp++) printf("c(%d %d)[%f] ", ip, jp, fCov[n++]);
168     printf("\n");
169   }
170 }
171
172 //_________________________________________________________________
173 void TKDNodeInfo::Store(TVectorD const *par, TMatrixD const *cov)
174 {
175 // Store the parameters and the covariance in the node
176
177   if(!fPar){SetNpar(); fPar = new Double_t[fNpar];}
178   for(int ip=0; ip<fNpar; ip++) fPar[ip] = (*par)[ip];
179
180   if(!cov) return;
181   if(!fCov){SetNcov(); fCov = new Double_t[fNcov];}
182   for(int ip(0), np(0); ip<fNpar; ip++)
183     for(int jp=ip; jp<fNpar; jp++) fCov[np++] = (*cov)(ip,jp);
184 }
185
186 //_________________________________________________________________
187 Bool_t TKDNodeInfo::CookPDF(const Double_t *point, Double_t &result, Double_t &error) const
188 {
189 // Recalculate the PDF for one node from the results of interpolation (parameters and covariance matrix)
190
191   Int_t ndim = Int_t(fNDim/3.);
192   if(ndim>10) return kFALSE; // support only up to 10 dimensions
193   //printf("ndim[%d] npar[%d] ncov[%d]\n", ndim, fNpar, fNcov);
194
195   Double_t fdfdp[66]; memset(fdfdp, 0, ndim*sizeof(Double_t));
196   Int_t ipar = 0;
197   fdfdp[ipar++] = 1.;
198   for(int idim=0; idim<ndim; idim++){
199     fdfdp[ipar++] = point[idim];
200     for(int jdim=idim; jdim<ndim; jdim++) fdfdp[ipar++] = point[idim]*point[jdim];
201   }
202
203   // calculate estimation
204   result =0.; error = 0.;
205   for(int i=0; i<fNpar; i++) result += fdfdp[i]*fPar[i];
206   if(!fNcov) return kTRUE;
207
208   for(int i(0), n(0); i<fNpar; i++){
209     error += fdfdp[i]*fdfdp[i]*fCov[n++];
210     for(int j(i+1); j<fNpar; j++) error += 2.*fdfdp[i]*fdfdp[j]*fCov[n++];
211   }     
212   error = TMath::Sqrt(error);
213   
214   //printf("TKDNodeInfo::CookPDF() : %6.3f +- %6.3f\n", result, error);
215
216   return kTRUE;
217 }
218
219
220
221 //_________________________________________________________________
222 TKDNodeInfo::TKDNodeDraw::TKDNodeDraw() 
223   :TBox()
224   ,fCOG()
225   ,fNode(NULL)
226 {
227   SetFillStyle(3002);
228   SetFillColor(50+Int_t(gRandom->Uniform()*50.));
229
230   fCOG.SetMarkerStyle(3);
231   fCOG.SetMarkerSize(.7);
232   fCOG.SetMarkerColor(2);
233 }
234
235
236 //_________________________________________________________________
237 void TKDNodeInfo::TKDNodeDraw::Draw(Option_t* option)
238 {
239   TBox::Draw(option);
240   fCOG.Draw("p");
241 }
242
243 //_________________________________________________________________
244 void TKDNodeInfo::TKDNodeDraw::SetNode(TKDNodeInfo *node, UChar_t size, UChar_t ax1, UChar_t ax2)
245 {
246   fNode=node;
247   const Float_t kBorder = 0.;//1.E-4;
248   Float_t *bounds = &(node->Data()[size]);
249   fX1=bounds[2*ax1]+kBorder;
250   fX2=bounds[2*ax1+1]-kBorder;
251   fY1=bounds[2*ax2]+kBorder;
252   fY2=bounds[2*ax2+1]-kBorder;
253   
254   Float_t x(node->Data()[ax1]), y(node->Data()[ax2]);
255   fCOG.SetX(x); fCOG.SetY(y);
256 }
257
258
259 //_________________________________________________________________
260 void TKDNodeInfo::TKDNodeDraw::Print(const Option_t* option) const
261 {
262   if(!fNode) return;
263   fNode->Print(option);
264 }