]>
Commit | Line | Data |
---|---|---|
21f3a443 | 1 | /************************************************************************** |
2 | * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * | |
3 | * * | |
4 | * Author: The ALICE Off-line Project. * | |
5 | * Contributors are mentioned in the code where appropriate. * | |
6 | * * | |
7 | * Permission to use, copy, modify and distribute this software and its * | |
8 | * documentation strictly for non-commercial purposes is hereby granted * | |
9 | * without fee, provided that the above copyright notice appears in all * | |
10 | * copies and that both the copyright notice and this permission notice * | |
11 | * appear in the supporting documentation. The authors make no claims * | |
12 | * about the suitability of this software for any purpose. It is * | |
13 | * provided "as is" without express or implied warranty. * | |
14 | **************************************************************************/ | |
15 | ||
16 | ||
17 | /////////////////////////////////////////////////////////////////////////// | |
18 | // Class TStatToolkit | |
19 | // | |
20 | // Subset of matheamtical functions not included in the TMath | |
21 | // | |
22 | ||
23 | /////////////////////////////////////////////////////////////////////////// | |
24 | #include "TMath.h" | |
25 | #include "Riostream.h" | |
26 | #include "TH1F.h" | |
3240a856 | 27 | #include "TH2F.h" |
21f3a443 | 28 | #include "TH3.h" |
29 | #include "TF1.h" | |
30 | #include "TTree.h" | |
31 | #include "TChain.h" | |
32 | #include "TObjString.h" | |
33 | #include "TLinearFitter.h" | |
3d7cc0b4 | 34 | #include "TGraph2D.h" |
35 | #include "TGraph.h" | |
b3453fe7 | 36 | #include "TGraphErrors.h" |
377a7d60 | 37 | #include "TMultiGraph.h" |
38 | #include "TCanvas.h" | |
39 | #include "TLatex.h" | |
40 | #include "TCut.h" | |
21f3a443 | 41 | |
42 | // | |
43 | // includes neccessary for test functions | |
44 | // | |
45 | #include "TSystem.h" | |
46 | #include "TRandom.h" | |
47 | #include "TStopwatch.h" | |
48 | #include "TTreeStream.h" | |
49 | ||
50 | #include "TStatToolkit.h" | |
51 | ||
52 | ||
53 | ClassImp(TStatToolkit) // Class implementation to enable ROOT I/O | |
54 | ||
55 | TStatToolkit::TStatToolkit() : TObject() | |
56 | { | |
57 | // | |
58 | // Default constructor | |
59 | // | |
60 | } | |
61 | /////////////////////////////////////////////////////////////////////////// | |
62 | TStatToolkit::~TStatToolkit() | |
63 | { | |
64 | // | |
65 | // Destructor | |
66 | // | |
67 | } | |
68 | ||
69 | ||
70 | //_____________________________________________________________________________ | |
71 | void TStatToolkit::EvaluateUni(Int_t nvectors, Double_t *data, Double_t &mean | |
72 | , Double_t &sigma, Int_t hh) | |
73 | { | |
74 | // | |
75 | // Robust estimator in 1D case MI version - (faster than ROOT version) | |
76 | // | |
77 | // For the univariate case | |
78 | // estimates of location and scatter are returned in mean and sigma parameters | |
79 | // the algorithm works on the same principle as in multivariate case - | |
80 | // it finds a subset of size hh with smallest sigma, and then returns mean and | |
81 | // sigma of this subset | |
82 | // | |
83 | ||
84 | if (hh==0) | |
85 | hh=(nvectors+2)/2; | |
86 | Double_t faclts[]={2.6477,2.5092,2.3826,2.2662,2.1587,2.0589,1.9660,1.879,1.7973,1.7203,1.6473}; | |
87 | Int_t *index=new Int_t[nvectors]; | |
88 | TMath::Sort(nvectors, data, index, kFALSE); | |
89 | ||
90 | Int_t nquant = TMath::Min(Int_t(Double_t(((hh*1./nvectors)-0.5)*40))+1, 11); | |
91 | Double_t factor = faclts[TMath::Max(0,nquant-1)]; | |
92 | ||
93 | Double_t sumx =0; | |
94 | Double_t sumx2 =0; | |
95 | Int_t bestindex = -1; | |
96 | Double_t bestmean = 0; | |
97 | Double_t bestsigma = (data[index[nvectors-1]]-data[index[0]]+1.); // maximal possible sigma | |
98 | bestsigma *=bestsigma; | |
99 | ||
100 | for (Int_t i=0; i<hh; i++){ | |
101 | sumx += data[index[i]]; | |
102 | sumx2 += data[index[i]]*data[index[i]]; | |
103 | } | |
104 | ||
105 | Double_t norm = 1./Double_t(hh); | |
bd7b4d18 | 106 | Double_t norm2 = (hh-1)>0 ? 1./Double_t(hh-1):1; |
21f3a443 | 107 | for (Int_t i=hh; i<nvectors; i++){ |
108 | Double_t cmean = sumx*norm; | |
109 | Double_t csigma = (sumx2 - hh*cmean*cmean)*norm2; | |
110 | if (csigma<bestsigma){ | |
111 | bestmean = cmean; | |
112 | bestsigma = csigma; | |
113 | bestindex = i-hh; | |
114 | } | |
115 | ||
116 | sumx += data[index[i]]-data[index[i-hh]]; | |
117 | sumx2 += data[index[i]]*data[index[i]]-data[index[i-hh]]*data[index[i-hh]]; | |
118 | } | |
119 | ||
120 | Double_t bstd=factor*TMath::Sqrt(TMath::Abs(bestsigma)); | |
121 | mean = bestmean; | |
122 | sigma = bstd; | |
123 | delete [] index; | |
124 | ||
125 | } | |
126 | ||
127 | ||
128 | ||
129 | void TStatToolkit::EvaluateUniExternal(Int_t nvectors, Double_t *data, Double_t &mean, Double_t &sigma, Int_t hh, Float_t externalfactor) | |
130 | { | |
131 | // Modified version of ROOT robust EvaluateUni | |
132 | // robust estimator in 1D case MI version | |
133 | // added external factor to include precision of external measurement | |
134 | // | |
135 | ||
136 | if (hh==0) | |
137 | hh=(nvectors+2)/2; | |
138 | Double_t faclts[]={2.6477,2.5092,2.3826,2.2662,2.1587,2.0589,1.9660,1.879,1.7973,1.7203,1.6473}; | |
139 | Int_t *index=new Int_t[nvectors]; | |
140 | TMath::Sort(nvectors, data, index, kFALSE); | |
141 | // | |
142 | Int_t nquant = TMath::Min(Int_t(Double_t(((hh*1./nvectors)-0.5)*40))+1, 11); | |
143 | Double_t factor = faclts[0]; | |
144 | if (nquant>0){ | |
145 | // fix proper normalization - Anja | |
146 | factor = faclts[nquant-1]; | |
147 | } | |
148 | ||
149 | // | |
150 | // | |
151 | Double_t sumx =0; | |
152 | Double_t sumx2 =0; | |
153 | Int_t bestindex = -1; | |
154 | Double_t bestmean = 0; | |
155 | Double_t bestsigma = -1; | |
156 | for (Int_t i=0; i<hh; i++){ | |
157 | sumx += data[index[i]]; | |
158 | sumx2 += data[index[i]]*data[index[i]]; | |
159 | } | |
160 | // | |
161 | Double_t kfactor = 2.*externalfactor - externalfactor*externalfactor; | |
162 | Double_t norm = 1./Double_t(hh); | |
163 | for (Int_t i=hh; i<nvectors; i++){ | |
164 | Double_t cmean = sumx*norm; | |
165 | Double_t csigma = (sumx2*norm - cmean*cmean*kfactor); | |
166 | if (csigma<bestsigma || bestsigma<0){ | |
167 | bestmean = cmean; | |
168 | bestsigma = csigma; | |
169 | bestindex = i-hh; | |
170 | } | |
171 | // | |
172 | // | |
173 | sumx += data[index[i]]-data[index[i-hh]]; | |
174 | sumx2 += data[index[i]]*data[index[i]]-data[index[i-hh]]*data[index[i-hh]]; | |
175 | } | |
176 | ||
177 | Double_t bstd=factor*TMath::Sqrt(TMath::Abs(bestsigma)); | |
178 | mean = bestmean; | |
179 | sigma = bstd; | |
180 | delete [] index; | |
181 | } | |
182 | ||
183 | ||
184 | //_____________________________________________________________________________ | |
185 | Int_t TStatToolkit::Freq(Int_t n, const Int_t *inlist | |
186 | , Int_t *outlist, Bool_t down) | |
187 | { | |
188 | // | |
189 | // Sort eleements according occurancy | |
190 | // The size of output array has is 2*n | |
191 | // | |
192 | ||
193 | Int_t * sindexS = new Int_t[n]; // temp array for sorting | |
194 | Int_t * sindexF = new Int_t[2*n]; | |
b8072cce | 195 | for (Int_t i=0;i<n;i++) sindexS[i]=0; |
196 | for (Int_t i=0;i<2*n;i++) sindexF[i]=0; | |
21f3a443 | 197 | // |
198 | TMath::Sort(n,inlist, sindexS, down); | |
199 | Int_t last = inlist[sindexS[0]]; | |
200 | Int_t val = last; | |
201 | sindexF[0] = 1; | |
202 | sindexF[0+n] = last; | |
203 | Int_t countPos = 0; | |
204 | // | |
205 | // find frequency | |
206 | for(Int_t i=1;i<n; i++){ | |
207 | val = inlist[sindexS[i]]; | |
208 | if (last == val) sindexF[countPos]++; | |
209 | else{ | |
210 | countPos++; | |
211 | sindexF[countPos+n] = val; | |
212 | sindexF[countPos]++; | |
213 | last =val; | |
214 | } | |
215 | } | |
216 | if (last==val) countPos++; | |
217 | // sort according frequency | |
218 | TMath::Sort(countPos, sindexF, sindexS, kTRUE); | |
219 | for (Int_t i=0;i<countPos;i++){ | |
220 | outlist[2*i ] = sindexF[sindexS[i]+n]; | |
221 | outlist[2*i+1] = sindexF[sindexS[i]]; | |
222 | } | |
223 | delete [] sindexS; | |
224 | delete [] sindexF; | |
225 | ||
226 | return countPos; | |
227 | ||
228 | } | |
229 | ||
230 | //___TStatToolkit__________________________________________________________________________ | |
3d7cc0b4 | 231 | void TStatToolkit::TruncatedMean(const TH1 * his, TVectorD *param, Float_t down, Float_t up, Bool_t verbose){ |
21f3a443 | 232 | // |
233 | // | |
234 | // | |
235 | Int_t nbins = his->GetNbinsX(); | |
236 | Float_t nentries = his->GetEntries(); | |
237 | Float_t sum =0; | |
238 | Float_t mean = 0; | |
239 | Float_t sigma2 = 0; | |
240 | Float_t ncumul=0; | |
241 | for (Int_t ibin=1;ibin<nbins; ibin++){ | |
242 | ncumul+= his->GetBinContent(ibin); | |
243 | Float_t fraction = Float_t(ncumul)/Float_t(nentries); | |
244 | if (fraction>down && fraction<up){ | |
245 | sum+=his->GetBinContent(ibin); | |
246 | mean+=his->GetBinCenter(ibin)*his->GetBinContent(ibin); | |
247 | sigma2+=his->GetBinCenter(ibin)*his->GetBinCenter(ibin)*his->GetBinContent(ibin); | |
248 | } | |
249 | } | |
250 | mean/=sum; | |
251 | sigma2= TMath::Sqrt(TMath::Abs(sigma2/sum-mean*mean)); | |
252 | if (param){ | |
253 | (*param)[0] = his->GetMaximum(); | |
254 | (*param)[1] = mean; | |
255 | (*param)[2] = sigma2; | |
256 | ||
257 | } | |
258 | if (verbose) printf("Mean\t%f\t Sigma2\t%f\n", mean,sigma2); | |
259 | } | |
260 | ||
32d8f622 | 261 | void TStatToolkit::LTM(TH1 * his, TVectorD *param , Float_t fraction, Bool_t verbose){ |
21f3a443 | 262 | // |
32d8f622 | 263 | // LTM : Trimmed mean on histogram - Modified version for binned data |
21f3a443 | 264 | // |
32d8f622 | 265 | // Robust statistic to estimate properties of the distribution |
266 | // See http://en.wikipedia.org/w/index.php?title=Trimmed_estimator&oldid=582847999 | |
267 | // | |
268 | // New faster version is under preparation | |
269 | // | |
270 | if (!param) return; | |
271 | (*param)[0]=0; | |
272 | (*param)[1]=0; | |
273 | (*param)[2]=0; | |
21f3a443 | 274 | Int_t nbins = his->GetNbinsX(); |
275 | Int_t nentries = (Int_t)his->GetEntries(); | |
32d8f622 | 276 | if (nentries<=0) return; |
21f3a443 | 277 | Double_t *data = new Double_t[nentries]; |
278 | Int_t npoints=0; | |
279 | for (Int_t ibin=1;ibin<nbins; ibin++){ | |
280 | Float_t entriesI = his->GetBinContent(ibin); | |
281 | Float_t xcenter= his->GetBinCenter(ibin); | |
282 | for (Int_t ic=0; ic<entriesI; ic++){ | |
283 | if (npoints<nentries){ | |
284 | data[npoints]= xcenter; | |
285 | npoints++; | |
286 | } | |
287 | } | |
288 | } | |
32d8f622 | 289 | Double_t mean, sigma; |
21f3a443 | 290 | Int_t npoints2=TMath::Min(Int_t(fraction*Float_t(npoints)),npoints-1); |
291 | npoints2=TMath::Max(Int_t(0.5*Float_t(npoints)),npoints2); | |
292 | TStatToolkit::EvaluateUni(npoints, data, mean,sigma,npoints2); | |
293 | delete [] data; | |
294 | if (verbose) printf("Mean\t%f\t Sigma2\t%f\n", mean,sigma);if (param){ | |
295 | (*param)[0] = his->GetMaximum(); | |
296 | (*param)[1] = mean; | |
297 | (*param)[2] = sigma; | |
298 | } | |
299 | } | |
300 | ||
94a43b22 | 301 | Double_t TStatToolkit::FitGaus(TH1* his, TVectorD *param, TMatrixD */*matrix*/, Float_t xmin, Float_t xmax, Bool_t verbose){ |
21f3a443 | 302 | // |
303 | // Fit histogram with gaussian function | |
304 | // | |
305 | // Prameters: | |
306 | // return value- chi2 - if negative ( not enough points) | |
307 | // his - input histogram | |
308 | // param - vector with parameters | |
309 | // xmin, xmax - range to fit - if xmin=xmax=0 - the full histogram range used | |
310 | // Fitting: | |
311 | // 1. Step - make logarithm | |
312 | // 2. Linear fit (parabola) - more robust - always converge | |
313 | // 3. In case of small statistic bins are averaged | |
314 | // | |
315 | static TLinearFitter fitter(3,"pol2"); | |
316 | TVectorD par(3); | |
317 | TVectorD sigma(3); | |
318 | TMatrixD mat(3,3); | |
319 | if (his->GetMaximum()<4) return -1; | |
320 | if (his->GetEntries()<12) return -1; | |
321 | if (his->GetRMS()<mat.GetTol()) return -1; | |
322 | Float_t maxEstimate = his->GetEntries()*his->GetBinWidth(1)/TMath::Sqrt((TMath::TwoPi()*his->GetRMS())); | |
323 | Int_t dsmooth = TMath::Nint(6./TMath::Sqrt(maxEstimate)); | |
324 | ||
325 | if (maxEstimate<1) return -1; | |
326 | Int_t nbins = his->GetNbinsX(); | |
327 | Int_t npoints=0; | |
328 | // | |
329 | ||
330 | ||
331 | if (xmin>=xmax){ | |
332 | xmin = his->GetXaxis()->GetXmin(); | |
333 | xmax = his->GetXaxis()->GetXmax(); | |
334 | } | |
335 | for (Int_t iter=0; iter<2; iter++){ | |
336 | fitter.ClearPoints(); | |
337 | npoints=0; | |
338 | for (Int_t ibin=1;ibin<nbins+1; ibin++){ | |
339 | Int_t countB=1; | |
340 | Float_t entriesI = his->GetBinContent(ibin); | |
341 | for (Int_t delta = -dsmooth; delta<=dsmooth; delta++){ | |
342 | if (ibin+delta>1 &&ibin+delta<nbins-1){ | |
343 | entriesI += his->GetBinContent(ibin+delta); | |
344 | countB++; | |
345 | } | |
346 | } | |
347 | entriesI/=countB; | |
348 | Double_t xcenter= his->GetBinCenter(ibin); | |
349 | if (xcenter<xmin || xcenter>xmax) continue; | |
350 | Double_t error=1./TMath::Sqrt(countB); | |
351 | Float_t cont=2; | |
352 | if (iter>0){ | |
353 | if (par[0]+par[1]*xcenter+par[2]*xcenter*xcenter>20) return 0; | |
354 | cont = TMath::Exp(par[0]+par[1]*xcenter+par[2]*xcenter*xcenter); | |
355 | if (cont>1.) error = 1./TMath::Sqrt(cont*Float_t(countB)); | |
356 | } | |
357 | if (entriesI>1&&cont>1){ | |
358 | fitter.AddPoint(&xcenter,TMath::Log(Float_t(entriesI)),error); | |
359 | npoints++; | |
360 | } | |
361 | } | |
362 | if (npoints>3){ | |
363 | fitter.Eval(); | |
364 | fitter.GetParameters(par); | |
365 | }else{ | |
366 | break; | |
367 | } | |
368 | } | |
369 | if (npoints<=3){ | |
370 | return -1; | |
371 | } | |
372 | fitter.GetParameters(par); | |
373 | fitter.GetCovarianceMatrix(mat); | |
374 | if (TMath::Abs(par[1])<mat.GetTol()) return -1; | |
375 | if (TMath::Abs(par[2])<mat.GetTol()) return -1; | |
376 | Double_t chi2 = fitter.GetChisquare()/Float_t(npoints); | |
377 | //fitter.GetParameters(); | |
378 | if (!param) param = new TVectorD(3); | |
cb1d20de | 379 | // if (!matrix) matrix = new TMatrixD(3,3); // Covariance matrix to be implemented |
21f3a443 | 380 | (*param)[1] = par[1]/(-2.*par[2]); |
381 | (*param)[2] = 1./TMath::Sqrt(TMath::Abs(-2.*par[2])); | |
382 | (*param)[0] = TMath::Exp(par[0]+ par[1]* (*param)[1] + par[2]*(*param)[1]*(*param)[1]); | |
383 | if (verbose){ | |
384 | par.Print(); | |
385 | mat.Print(); | |
386 | param->Print(); | |
387 | printf("Chi2=%f\n",chi2); | |
388 | TF1 * f1= new TF1("f1","[0]*exp(-(x-[1])^2/(2*[2]*[2]))",his->GetXaxis()->GetXmin(),his->GetXaxis()->GetXmax()); | |
389 | f1->SetParameter(0, (*param)[0]); | |
390 | f1->SetParameter(1, (*param)[1]); | |
391 | f1->SetParameter(2, (*param)[2]); | |
392 | f1->Draw("same"); | |
393 | } | |
394 | return chi2; | |
395 | } | |
396 | ||
cb1d20de | 397 | Double_t TStatToolkit::FitGaus(Float_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, TVectorD *param, TMatrixD */*matrix*/, Bool_t verbose){ |
21f3a443 | 398 | // |
399 | // Fit histogram with gaussian function | |
400 | // | |
401 | // Prameters: | |
402 | // nbins: size of the array and number of histogram bins | |
403 | // xMin, xMax: histogram range | |
404 | // param: paramters of the fit (0-Constant, 1-Mean, 2-Sigma) | |
405 | // matrix: covariance matrix -- not implemented yet, pass dummy matrix!!! | |
406 | // | |
407 | // Return values: | |
408 | // >0: the chi2 returned by TLinearFitter | |
409 | // -3: only three points have been used for the calculation - no fitter was used | |
410 | // -2: only two points have been used for the calculation - center of gravity was uesed for calculation | |
411 | // -1: only one point has been used for the calculation - center of gravity was uesed for calculation | |
412 | // -4: invalid result!! | |
413 | // | |
414 | // Fitting: | |
415 | // 1. Step - make logarithm | |
416 | // 2. Linear fit (parabola) - more robust - always converge | |
417 | // | |
418 | static TLinearFitter fitter(3,"pol2"); | |
419 | static TMatrixD mat(3,3); | |
420 | static Double_t kTol = mat.GetTol(); | |
421 | fitter.StoreData(kFALSE); | |
422 | fitter.ClearPoints(); | |
423 | TVectorD par(3); | |
424 | TVectorD sigma(3); | |
3d7cc0b4 | 425 | TMatrixD matA(3,3); |
21f3a443 | 426 | TMatrixD b(3,1); |
427 | Float_t rms = TMath::RMS(nBins,arr); | |
428 | Float_t max = TMath::MaxElement(nBins,arr); | |
429 | Float_t binWidth = (xMax-xMin)/(Float_t)nBins; | |
430 | ||
431 | Float_t meanCOG = 0; | |
432 | Float_t rms2COG = 0; | |
433 | Float_t sumCOG = 0; | |
434 | ||
435 | Float_t entries = 0; | |
436 | Int_t nfilled=0; | |
437 | ||
438 | for (Int_t i=0; i<nBins; i++){ | |
439 | entries+=arr[i]; | |
440 | if (arr[i]>0) nfilled++; | |
441 | } | |
442 | ||
443 | if (max<4) return -4; | |
444 | if (entries<12) return -4; | |
445 | if (rms<kTol) return -4; | |
446 | ||
447 | Int_t npoints=0; | |
448 | // | |
449 | ||
450 | // | |
451 | for (Int_t ibin=0;ibin<nBins; ibin++){ | |
452 | Float_t entriesI = arr[ibin]; | |
453 | if (entriesI>1){ | |
454 | Double_t xcenter = xMin+(ibin+0.5)*binWidth; | |
455 | ||
456 | Float_t error = 1./TMath::Sqrt(entriesI); | |
457 | Float_t val = TMath::Log(Float_t(entriesI)); | |
458 | fitter.AddPoint(&xcenter,val,error); | |
459 | if (npoints<3){ | |
3d7cc0b4 | 460 | matA(npoints,0)=1; |
461 | matA(npoints,1)=xcenter; | |
462 | matA(npoints,2)=xcenter*xcenter; | |
21f3a443 | 463 | b(npoints,0)=val; |
464 | meanCOG+=xcenter*entriesI; | |
465 | rms2COG +=xcenter*entriesI*xcenter; | |
466 | sumCOG +=entriesI; | |
467 | } | |
468 | npoints++; | |
469 | } | |
470 | } | |
471 | ||
472 | ||
473 | Double_t chi2 = 0; | |
474 | if (npoints>=3){ | |
475 | if ( npoints == 3 ){ | |
476 | //analytic calculation of the parameters for three points | |
3d7cc0b4 | 477 | matA.Invert(); |
21f3a443 | 478 | TMatrixD res(1,3); |
3d7cc0b4 | 479 | res.Mult(matA,b); |
21f3a443 | 480 | par[0]=res(0,0); |
481 | par[1]=res(0,1); | |
482 | par[2]=res(0,2); | |
483 | chi2 = -3.; | |
484 | } else { | |
485 | // use fitter for more than three points | |
486 | fitter.Eval(); | |
487 | fitter.GetParameters(par); | |
488 | fitter.GetCovarianceMatrix(mat); | |
489 | chi2 = fitter.GetChisquare()/Float_t(npoints); | |
490 | } | |
491 | if (TMath::Abs(par[1])<kTol) return -4; | |
492 | if (TMath::Abs(par[2])<kTol) return -4; | |
493 | ||
494 | if (!param) param = new TVectorD(3); | |
cb1d20de | 495 | //if (!matrix) matrix = new TMatrixD(3,3); // !!!!might be a memory leek. use dummy matrix pointer to call this function! // Covariance matrix to be implemented |
21f3a443 | 496 | |
497 | (*param)[1] = par[1]/(-2.*par[2]); | |
498 | (*param)[2] = 1./TMath::Sqrt(TMath::Abs(-2.*par[2])); | |
499 | Double_t lnparam0 = par[0]+ par[1]* (*param)[1] + par[2]*(*param)[1]*(*param)[1]; | |
500 | if ( lnparam0>307 ) return -4; | |
501 | (*param)[0] = TMath::Exp(lnparam0); | |
502 | if (verbose){ | |
503 | par.Print(); | |
504 | mat.Print(); | |
505 | param->Print(); | |
506 | printf("Chi2=%f\n",chi2); | |
507 | TF1 * f1= new TF1("f1","[0]*exp(-(x-[1])^2/(2*[2]*[2]))",xMin,xMax); | |
508 | f1->SetParameter(0, (*param)[0]); | |
509 | f1->SetParameter(1, (*param)[1]); | |
510 | f1->SetParameter(2, (*param)[2]); | |
511 | f1->Draw("same"); | |
512 | } | |
513 | return chi2; | |
514 | } | |
515 | ||
516 | if (npoints == 2){ | |
517 | //use center of gravity for 2 points | |
518 | meanCOG/=sumCOG; | |
519 | rms2COG /=sumCOG; | |
520 | (*param)[0] = max; | |
521 | (*param)[1] = meanCOG; | |
522 | (*param)[2] = TMath::Sqrt(TMath::Abs(meanCOG*meanCOG-rms2COG)); | |
523 | chi2=-2.; | |
524 | } | |
525 | if ( npoints == 1 ){ | |
526 | meanCOG/=sumCOG; | |
527 | (*param)[0] = max; | |
528 | (*param)[1] = meanCOG; | |
529 | (*param)[2] = binWidth/TMath::Sqrt(12); | |
530 | chi2=-1.; | |
531 | } | |
532 | return chi2; | |
533 | ||
534 | } | |
535 | ||
536 | ||
3d7cc0b4 | 537 | Float_t TStatToolkit::GetCOG(const Short_t *arr, Int_t nBins, Float_t xMin, Float_t xMax, Float_t *rms, Float_t *sum) |
21f3a443 | 538 | { |
539 | // | |
540 | // calculate center of gravity rms and sum for array 'arr' with nBins an a x range xMin to xMax | |
541 | // return COG; in case of failure return xMin | |
542 | // | |
543 | Float_t meanCOG = 0; | |
544 | Float_t rms2COG = 0; | |
545 | Float_t sumCOG = 0; | |
546 | Int_t npoints = 0; | |
547 | ||
548 | Float_t binWidth = (xMax-xMin)/(Float_t)nBins; | |
549 | ||
550 | for (Int_t ibin=0; ibin<nBins; ibin++){ | |
551 | Float_t entriesI = (Float_t)arr[ibin]; | |
552 | Double_t xcenter = xMin+(ibin+0.5)*binWidth; | |
553 | if ( entriesI>0 ){ | |
554 | meanCOG += xcenter*entriesI; | |
555 | rms2COG += xcenter*entriesI*xcenter; | |
556 | sumCOG += entriesI; | |
557 | npoints++; | |
558 | } | |
559 | } | |
560 | if ( sumCOG == 0 ) return xMin; | |
561 | meanCOG/=sumCOG; | |
562 | ||
563 | if ( rms ){ | |
564 | rms2COG /=sumCOG; | |
565 | (*rms) = TMath::Sqrt(TMath::Abs(meanCOG*meanCOG-rms2COG)); | |
566 | if ( npoints == 1 ) (*rms) = binWidth/TMath::Sqrt(12); | |
567 | } | |
568 | ||
569 | if ( sum ) | |
570 | (*sum) = sumCOG; | |
571 | ||
572 | return meanCOG; | |
573 | } | |
574 | ||
575 | ||
576 | ||
577 | /////////////////////////////////////////////////////////////// | |
578 | ////////////// TEST functions ///////////////////////// | |
579 | /////////////////////////////////////////////////////////////// | |
580 | ||
581 | ||
582 | ||
583 | ||
584 | ||
585 | void TStatToolkit::TestGausFit(Int_t nhistos){ | |
586 | // | |
587 | // Test performance of the parabolic - gaussian fit - compare it with | |
588 | // ROOT gauss fit | |
589 | // nhistos - number of histograms to be used for test | |
590 | // | |
32d8f622 | 591 | TTreeSRedirector *pcstream = new TTreeSRedirector("fitdebug.root","recreate"); |
21f3a443 | 592 | |
593 | Float_t *xTrue = new Float_t[nhistos]; | |
594 | Float_t *sTrue = new Float_t[nhistos]; | |
595 | TVectorD **par1 = new TVectorD*[nhistos]; | |
596 | TVectorD **par2 = new TVectorD*[nhistos]; | |
597 | TMatrixD dummy(3,3); | |
598 | ||
599 | ||
600 | TH1F **h1f = new TH1F*[nhistos]; | |
601 | TF1 *myg = new TF1("myg","gaus"); | |
602 | TF1 *fit = new TF1("fit","gaus"); | |
603 | gRandom->SetSeed(0); | |
604 | ||
605 | //init | |
606 | for (Int_t i=0;i<nhistos; i++){ | |
607 | par1[i] = new TVectorD(3); | |
608 | par2[i] = new TVectorD(3); | |
609 | h1f[i] = new TH1F(Form("h1f%d",i),Form("h1f%d",i),20,-10,10); | |
610 | xTrue[i]= gRandom->Rndm(); | |
611 | gSystem->Sleep(2); | |
612 | sTrue[i]= .75+gRandom->Rndm()*.5; | |
613 | myg->SetParameters(1,xTrue[i],sTrue[i]); | |
614 | h1f[i]->FillRandom("myg"); | |
615 | } | |
616 | ||
32d8f622 | 617 | TStopwatch s; |
21f3a443 | 618 | s.Start(); |
619 | //standard gaus fit | |
620 | for (Int_t i=0; i<nhistos; i++){ | |
621 | h1f[i]->Fit(fit,"0q"); | |
622 | (*par1[i])(0) = fit->GetParameter(0); | |
623 | (*par1[i])(1) = fit->GetParameter(1); | |
624 | (*par1[i])(2) = fit->GetParameter(2); | |
625 | } | |
626 | s.Stop(); | |
627 | printf("Gaussian fit\t"); | |
628 | s.Print(); | |
629 | ||
630 | s.Start(); | |
631 | //TStatToolkit gaus fit | |
632 | for (Int_t i=0; i<nhistos; i++){ | |
633 | TStatToolkit::FitGaus(h1f[i]->GetArray()+1,h1f[i]->GetNbinsX(),h1f[i]->GetXaxis()->GetXmin(),h1f[i]->GetXaxis()->GetXmax(),par2[i],&dummy); | |
634 | } | |
635 | ||
636 | s.Stop(); | |
637 | printf("Parabolic fit\t"); | |
638 | s.Print(); | |
639 | //write stream | |
640 | for (Int_t i=0;i<nhistos; i++){ | |
641 | Float_t xt = xTrue[i]; | |
642 | Float_t st = sTrue[i]; | |
643 | (*pcstream)<<"data" | |
644 | <<"xTrue="<<xt | |
645 | <<"sTrue="<<st | |
646 | <<"pg.="<<(par1[i]) | |
647 | <<"pa.="<<(par2[i]) | |
648 | <<"\n"; | |
649 | } | |
650 | //delete pointers | |
651 | for (Int_t i=0;i<nhistos; i++){ | |
652 | delete par1[i]; | |
653 | delete par2[i]; | |
654 | delete h1f[i]; | |
655 | } | |
656 | delete pcstream; | |
657 | delete []h1f; | |
658 | delete []xTrue; | |
659 | delete []sTrue; | |
660 | // | |
661 | delete []par1; | |
662 | delete []par2; | |
663 | ||
664 | } | |
665 | ||
666 | ||
667 | ||
668 | TGraph2D * TStatToolkit::MakeStat2D(TH3 * his, Int_t delta0, Int_t delta1, Int_t type){ | |
669 | // | |
670 | // | |
671 | // | |
672 | // delta - number of bins to integrate | |
673 | // type - 0 - mean value | |
674 | ||
675 | TAxis * xaxis = his->GetXaxis(); | |
676 | TAxis * yaxis = his->GetYaxis(); | |
677 | // TAxis * zaxis = his->GetZaxis(); | |
678 | Int_t nbinx = xaxis->GetNbins(); | |
679 | Int_t nbiny = yaxis->GetNbins(); | |
680 | char name[1000]; | |
681 | Int_t icount=0; | |
682 | TGraph2D *graph = new TGraph2D(nbinx*nbiny); | |
683 | TF1 f1("f1","gaus"); | |
684 | for (Int_t ix=0; ix<nbinx;ix++) | |
685 | for (Int_t iy=0; iy<nbiny;iy++){ | |
686 | Float_t xcenter = xaxis->GetBinCenter(ix); | |
687 | Float_t ycenter = yaxis->GetBinCenter(iy); | |
cb1d20de | 688 | snprintf(name,1000,"%s_%d_%d",his->GetName(), ix,iy); |
21f3a443 | 689 | TH1 *projection = his->ProjectionZ(name,ix-delta0,ix+delta0,iy-delta1,iy+delta1); |
690 | Float_t stat= 0; | |
691 | if (type==0) stat = projection->GetMean(); | |
692 | if (type==1) stat = projection->GetRMS(); | |
693 | if (type==2 || type==3){ | |
694 | TVectorD vec(3); | |
695 | TStatToolkit::LTM((TH1F*)projection,&vec,0.7); | |
696 | if (type==2) stat= vec[1]; | |
697 | if (type==3) stat= vec[0]; | |
698 | } | |
699 | if (type==4|| type==5){ | |
700 | projection->Fit(&f1); | |
701 | if (type==4) stat= f1.GetParameter(1); | |
702 | if (type==5) stat= f1.GetParameter(2); | |
703 | } | |
704 | //printf("%d\t%f\t%f\t%f\n", icount,xcenter, ycenter, stat); | |
705 | graph->SetPoint(icount,xcenter, ycenter, stat); | |
706 | icount++; | |
707 | } | |
708 | return graph; | |
709 | } | |
710 | ||
32d8f622 | 711 | TGraphErrors * TStatToolkit::MakeStat1D(TH2 * his, Int_t deltaBin, Double_t fraction, Int_t returnType, Int_t markerStyle, Int_t markerColor){ |
712 | // | |
713 | // function to retrieve the "mean and RMS estimate" of 2D histograms | |
714 | // | |
715 | // Robust statistic to estimate properties of the distribution | |
716 | // See http://en.wikipedia.org/wiki/Trimmed_estimator | |
717 | // | |
718 | // deltaBin - number of bins to integrate (bin+-deltaBin) | |
719 | // fraction - fraction of values for the LTM and for the gauss fit | |
720 | // returnType - | |
721 | // 0 - mean value | |
722 | // 1 - RMS | |
723 | // 2 - LTM mean | |
724 | // 3 - LTM sigma | |
725 | // 4 - Gaus fit mean - on LTM range | |
726 | // 5 - Gaus fit sigma - on LTM range | |
727 | // | |
21f3a443 | 728 | TAxis * xaxis = his->GetXaxis(); |
21f3a443 | 729 | Int_t nbinx = xaxis->GetNbins(); |
21f3a443 | 730 | char name[1000]; |
731 | Int_t icount=0; | |
32d8f622 | 732 | // |
733 | TVectorD vecX(nbinx); | |
734 | TVectorD vecXErr(nbinx); | |
735 | TVectorD vecY(nbinx); | |
736 | TVectorD vecYErr(nbinx); | |
737 | // | |
21f3a443 | 738 | TF1 f1("f1","gaus"); |
32d8f622 | 739 | TVectorD vecLTM(3); |
740 | ||
21f3a443 | 741 | for (Int_t ix=0; ix<nbinx;ix++){ |
742 | Float_t xcenter = xaxis->GetBinCenter(ix); | |
cb1d20de | 743 | snprintf(name,1000,"%s_%d",his->GetName(), ix); |
32d8f622 | 744 | TH1 *projection = his->ProjectionY(name,TMath::Max(ix-deltaBin,1),TMath::Min(ix+deltaBin,nbinx)); |
745 | Double_t stat= 0; | |
746 | Double_t err =0; | |
747 | TStatToolkit::LTM((TH1F*)projection,&vecLTM,fraction); | |
748 | // | |
749 | if (returnType==0) { | |
750 | stat = projection->GetMean(); | |
751 | err = projection->GetMeanError(); | |
21f3a443 | 752 | } |
32d8f622 | 753 | if (returnType==1) { |
754 | stat = projection->GetRMS(); | |
755 | err = projection->GetRMSError(); | |
21f3a443 | 756 | } |
32d8f622 | 757 | if (returnType==2 || returnType==3){ |
758 | if (returnType==2) stat= vecLTM[1]; | |
759 | if (returnType==3) stat= vecLTM[2]; | |
760 | } | |
761 | if (returnType==4|| returnType==5){ | |
762 | projection->Fit(&f1,"QNR","QNR"); | |
763 | if (returnType==4) { | |
764 | stat= f1.GetParameter(1); | |
765 | err=f1.GetParError(1); | |
766 | } | |
767 | if (returnType==5) { | |
768 | stat= f1.GetParameter(2); | |
769 | err=f1.GetParError(2); | |
770 | } | |
771 | } | |
772 | vecX[icount]=xcenter; | |
773 | vecY[icount]=stat; | |
774 | vecYErr[icount]=err; | |
21f3a443 | 775 | icount++; |
32d8f622 | 776 | delete projection; |
21f3a443 | 777 | } |
32d8f622 | 778 | TGraphErrors *graph = new TGraphErrors(icount,vecX.GetMatrixArray(), vecY.GetMatrixArray(),0, vecYErr.GetMatrixArray()); |
779 | graph->SetMarkerStyle(markerStyle); | |
780 | graph->SetMarkerColor(markerColor); | |
21f3a443 | 781 | return graph; |
782 | } | |
783 | ||
784 | ||
785 | ||
786 | ||
787 | ||
88b1c775 | 788 | TString* TStatToolkit::FitPlane(TTree *tree, const char* drawCommand, const char* formula, const char* cuts, Double_t & chi2, Int_t &npoints, TVectorD &fitParam, TMatrixD &covMatrix, Float_t frac, Int_t start, Int_t stop,Bool_t fix0){ |
21f3a443 | 789 | // |
790 | // fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts | |
791 | // returns chi2, fitParam and covMatrix | |
792 | // returns TString with fitted formula | |
793 | // | |
dd46129c | 794 | |
21f3a443 | 795 | TString formulaStr(formula); |
796 | TString drawStr(drawCommand); | |
797 | TString cutStr(cuts); | |
dd46129c | 798 | TString ferr("1"); |
799 | ||
800 | TString strVal(drawCommand); | |
801 | if (strVal.Contains(":")){ | |
802 | TObjArray* valTokens = strVal.Tokenize(":"); | |
803 | drawStr = valTokens->At(0)->GetName(); | |
804 | ferr = valTokens->At(1)->GetName(); | |
09d5920f | 805 | delete valTokens; |
dd46129c | 806 | } |
807 | ||
21f3a443 | 808 | |
809 | formulaStr.ReplaceAll("++", "~"); | |
810 | TObjArray* formulaTokens = formulaStr.Tokenize("~"); | |
811 | Int_t dim = formulaTokens->GetEntriesFast(); | |
812 | ||
813 | fitParam.ResizeTo(dim); | |
814 | covMatrix.ResizeTo(dim,dim); | |
815 | ||
816 | TLinearFitter* fitter = new TLinearFitter(dim+1, Form("hyp%d",dim)); | |
817 | fitter->StoreData(kTRUE); | |
818 | fitter->ClearPoints(); | |
819 | ||
820 | Int_t entries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start, start); | |
09d5920f | 821 | if (entries == -1) { |
822 | delete formulaTokens; | |
823 | return new TString("An ERROR has occured during fitting!"); | |
824 | } | |
bd7b4d18 | 825 | Double_t **values = new Double_t*[dim+1] ; |
826 | for (Int_t i=0; i<dim+1; i++) values[i]=NULL; | |
dd46129c | 827 | // |
828 | entries = tree->Draw(ferr.Data(), cutStr.Data(), "goff", stop-start, start); | |
b8072cce | 829 | if (entries == -1) { |
09d5920f | 830 | delete formulaTokens; |
b8072cce | 831 | delete []values; |
832 | return new TString("An ERROR has occured during fitting!"); | |
833 | } | |
dd46129c | 834 | Double_t *errors = new Double_t[entries]; |
835 | memcpy(errors, tree->GetV1(), entries*sizeof(Double_t)); | |
21f3a443 | 836 | |
837 | for (Int_t i = 0; i < dim + 1; i++){ | |
838 | Int_t centries = 0; | |
839 | if (i < dim) centries = tree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff", stop-start,start); | |
840 | else centries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start,start); | |
841 | ||
b8072cce | 842 | if (entries != centries) { |
843 | delete []errors; | |
844 | delete []values; | |
845 | return new TString("An ERROR has occured during fitting!"); | |
846 | } | |
21f3a443 | 847 | values[i] = new Double_t[entries]; |
848 | memcpy(values[i], tree->GetV1(), entries*sizeof(Double_t)); | |
849 | } | |
850 | ||
851 | // add points to the fitter | |
852 | for (Int_t i = 0; i < entries; i++){ | |
853 | Double_t x[1000]; | |
854 | for (Int_t j=0; j<dim;j++) x[j]=values[j][i]; | |
dd46129c | 855 | fitter->AddPoint(x, values[dim][i], errors[i]); |
21f3a443 | 856 | } |
857 | ||
858 | fitter->Eval(); | |
2c629c56 | 859 | if (frac>0.5 && frac<1){ |
860 | fitter->EvalRobust(frac); | |
88b1c775 | 861 | }else{ |
862 | if (fix0) { | |
863 | fitter->FixParameter(0,0); | |
864 | fitter->Eval(); | |
865 | } | |
2c629c56 | 866 | } |
21f3a443 | 867 | fitter->GetParameters(fitParam); |
868 | fitter->GetCovarianceMatrix(covMatrix); | |
869 | chi2 = fitter->GetChisquare(); | |
b8072cce | 870 | npoints = entries; |
21f3a443 | 871 | TString *preturnFormula = new TString(Form("( %f+",fitParam[0])), &returnFormula = *preturnFormula; |
872 | ||
873 | for (Int_t iparam = 0; iparam < dim; iparam++) { | |
874 | returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1])); | |
875 | if (iparam < dim-1) returnFormula.Append("+"); | |
876 | } | |
877 | returnFormula.Append(" )"); | |
4d61c301 | 878 | |
879 | ||
b8072cce | 880 | for (Int_t j=0; j<dim+1;j++) delete [] values[j]; |
4d61c301 | 881 | |
882 | ||
cb1d20de | 883 | delete formulaTokens; |
884 | delete fitter; | |
885 | delete[] values; | |
b8072cce | 886 | delete[] errors; |
cb1d20de | 887 | return preturnFormula; |
888 | } | |
889 | ||
890 | TString* TStatToolkit::FitPlaneConstrain(TTree *tree, const char* drawCommand, const char* formula, const char* cuts, Double_t & chi2, Int_t &npoints, TVectorD &fitParam, TMatrixD &covMatrix, Float_t frac, Int_t start, Int_t stop,Double_t constrain){ | |
891 | // | |
892 | // fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts | |
893 | // returns chi2, fitParam and covMatrix | |
894 | // returns TString with fitted formula | |
895 | // | |
896 | ||
897 | TString formulaStr(formula); | |
898 | TString drawStr(drawCommand); | |
899 | TString cutStr(cuts); | |
900 | TString ferr("1"); | |
901 | ||
902 | TString strVal(drawCommand); | |
903 | if (strVal.Contains(":")){ | |
904 | TObjArray* valTokens = strVal.Tokenize(":"); | |
905 | drawStr = valTokens->At(0)->GetName(); | |
906 | ferr = valTokens->At(1)->GetName(); | |
09d5920f | 907 | delete valTokens; |
cb1d20de | 908 | } |
909 | ||
910 | ||
911 | formulaStr.ReplaceAll("++", "~"); | |
912 | TObjArray* formulaTokens = formulaStr.Tokenize("~"); | |
913 | Int_t dim = formulaTokens->GetEntriesFast(); | |
914 | ||
915 | fitParam.ResizeTo(dim); | |
916 | covMatrix.ResizeTo(dim,dim); | |
917 | ||
918 | TLinearFitter* fitter = new TLinearFitter(dim+1, Form("hyp%d",dim)); | |
919 | fitter->StoreData(kTRUE); | |
920 | fitter->ClearPoints(); | |
921 | ||
922 | Int_t entries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start, start); | |
09d5920f | 923 | if (entries == -1) { |
924 | delete formulaTokens; | |
925 | return new TString("An ERROR has occured during fitting!"); | |
926 | } | |
cb1d20de | 927 | Double_t **values = new Double_t*[dim+1] ; |
bd7b4d18 | 928 | for (Int_t i=0; i<dim+1; i++) values[i]=NULL; |
cb1d20de | 929 | // |
930 | entries = tree->Draw(ferr.Data(), cutStr.Data(), "goff", stop-start, start); | |
b8072cce | 931 | if (entries == -1) { |
09d5920f | 932 | delete formulaTokens; |
b8072cce | 933 | delete [] values; |
934 | return new TString("An ERROR has occured during fitting!"); | |
935 | } | |
cb1d20de | 936 | Double_t *errors = new Double_t[entries]; |
937 | memcpy(errors, tree->GetV1(), entries*sizeof(Double_t)); | |
938 | ||
939 | for (Int_t i = 0; i < dim + 1; i++){ | |
940 | Int_t centries = 0; | |
941 | if (i < dim) centries = tree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff", stop-start,start); | |
942 | else centries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start,start); | |
943 | ||
b8072cce | 944 | if (entries != centries) { |
945 | delete []errors; | |
946 | delete []values; | |
09d5920f | 947 | delete formulaTokens; |
b8072cce | 948 | return new TString("An ERROR has occured during fitting!"); |
949 | } | |
cb1d20de | 950 | values[i] = new Double_t[entries]; |
951 | memcpy(values[i], tree->GetV1(), entries*sizeof(Double_t)); | |
952 | } | |
953 | ||
954 | // add points to the fitter | |
955 | for (Int_t i = 0; i < entries; i++){ | |
956 | Double_t x[1000]; | |
957 | for (Int_t j=0; j<dim;j++) x[j]=values[j][i]; | |
958 | fitter->AddPoint(x, values[dim][i], errors[i]); | |
959 | } | |
960 | if (constrain>0){ | |
961 | for (Int_t i = 0; i < dim; i++){ | |
962 | Double_t x[1000]; | |
963 | for (Int_t j=0; j<dim;j++) if (i!=j) x[j]=0; | |
964 | x[i]=1.; | |
965 | fitter->AddPoint(x, 0, constrain); | |
966 | } | |
967 | } | |
968 | ||
969 | ||
970 | fitter->Eval(); | |
971 | if (frac>0.5 && frac<1){ | |
972 | fitter->EvalRobust(frac); | |
973 | } | |
974 | fitter->GetParameters(fitParam); | |
975 | fitter->GetCovarianceMatrix(covMatrix); | |
976 | chi2 = fitter->GetChisquare(); | |
977 | npoints = entries; | |
cb1d20de | 978 | |
979 | TString *preturnFormula = new TString(Form("( %f+",fitParam[0])), &returnFormula = *preturnFormula; | |
980 | ||
981 | for (Int_t iparam = 0; iparam < dim; iparam++) { | |
982 | returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam+1])); | |
983 | if (iparam < dim-1) returnFormula.Append("+"); | |
984 | } | |
985 | returnFormula.Append(" )"); | |
986 | ||
b8072cce | 987 | for (Int_t j=0; j<dim+1;j++) delete [] values[j]; |
cb1d20de | 988 | |
989 | ||
990 | ||
991 | delete formulaTokens; | |
992 | delete fitter; | |
993 | delete[] values; | |
b8072cce | 994 | delete[] errors; |
cb1d20de | 995 | return preturnFormula; |
996 | } | |
997 | ||
998 | ||
999 | ||
1000 | TString* TStatToolkit::FitPlaneFixed(TTree *tree, const char* drawCommand, const char* formula, const char* cuts, Double_t & chi2, Int_t &npoints, TVectorD &fitParam, TMatrixD &covMatrix, Float_t frac, Int_t start, Int_t stop){ | |
1001 | // | |
1002 | // fit an arbitrary function, specified by formula into the data, specified by drawCommand and cuts | |
1003 | // returns chi2, fitParam and covMatrix | |
1004 | // returns TString with fitted formula | |
1005 | // | |
1006 | ||
1007 | TString formulaStr(formula); | |
1008 | TString drawStr(drawCommand); | |
1009 | TString cutStr(cuts); | |
1010 | TString ferr("1"); | |
1011 | ||
1012 | TString strVal(drawCommand); | |
1013 | if (strVal.Contains(":")){ | |
1014 | TObjArray* valTokens = strVal.Tokenize(":"); | |
1015 | drawStr = valTokens->At(0)->GetName(); | |
09d5920f | 1016 | ferr = valTokens->At(1)->GetName(); |
1017 | delete valTokens; | |
cb1d20de | 1018 | } |
1019 | ||
1020 | ||
1021 | formulaStr.ReplaceAll("++", "~"); | |
1022 | TObjArray* formulaTokens = formulaStr.Tokenize("~"); | |
1023 | Int_t dim = formulaTokens->GetEntriesFast(); | |
1024 | ||
1025 | fitParam.ResizeTo(dim); | |
1026 | covMatrix.ResizeTo(dim,dim); | |
1027 | TString fitString="x0"; | |
1028 | for (Int_t i=1; i<dim; i++) fitString+=Form("++x%d",i); | |
1029 | TLinearFitter* fitter = new TLinearFitter(dim, fitString.Data()); | |
1030 | fitter->StoreData(kTRUE); | |
1031 | fitter->ClearPoints(); | |
1032 | ||
1033 | Int_t entries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start, start); | |
09d5920f | 1034 | if (entries == -1) { |
1035 | delete formulaTokens; | |
1036 | return new TString("An ERROR has occured during fitting!"); | |
1037 | } | |
cb1d20de | 1038 | Double_t **values = new Double_t*[dim+1] ; |
bd7b4d18 | 1039 | for (Int_t i=0; i<dim+1; i++) values[i]=NULL; |
cb1d20de | 1040 | // |
1041 | entries = tree->Draw(ferr.Data(), cutStr.Data(), "goff", stop-start, start); | |
b8072cce | 1042 | if (entries == -1) { |
1043 | delete []values; | |
09d5920f | 1044 | delete formulaTokens; |
b8072cce | 1045 | return new TString("An ERROR has occured during fitting!"); |
1046 | } | |
cb1d20de | 1047 | Double_t *errors = new Double_t[entries]; |
1048 | memcpy(errors, tree->GetV1(), entries*sizeof(Double_t)); | |
1049 | ||
1050 | for (Int_t i = 0; i < dim + 1; i++){ | |
1051 | Int_t centries = 0; | |
1052 | if (i < dim) centries = tree->Draw(((TObjString*)formulaTokens->At(i))->GetName(), cutStr.Data(), "goff", stop-start,start); | |
1053 | else centries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff", stop-start,start); | |
1054 | ||
b8072cce | 1055 | if (entries != centries) { |
1056 | delete []errors; | |
1057 | delete []values; | |
09d5920f | 1058 | delete formulaTokens; |
b8072cce | 1059 | return new TString("An ERROR has occured during fitting!"); |
1060 | } | |
cb1d20de | 1061 | values[i] = new Double_t[entries]; |
1062 | memcpy(values[i], tree->GetV1(), entries*sizeof(Double_t)); | |
1063 | } | |
1064 | ||
1065 | // add points to the fitter | |
1066 | for (Int_t i = 0; i < entries; i++){ | |
1067 | Double_t x[1000]; | |
1068 | for (Int_t j=0; j<dim;j++) x[j]=values[j][i]; | |
1069 | fitter->AddPoint(x, values[dim][i], errors[i]); | |
1070 | } | |
1071 | ||
1072 | fitter->Eval(); | |
1073 | if (frac>0.5 && frac<1){ | |
1074 | fitter->EvalRobust(frac); | |
1075 | } | |
1076 | fitter->GetParameters(fitParam); | |
1077 | fitter->GetCovarianceMatrix(covMatrix); | |
1078 | chi2 = fitter->GetChisquare(); | |
1079 | npoints = entries; | |
cb1d20de | 1080 | |
1081 | TString *preturnFormula = new TString("("), &returnFormula = *preturnFormula; | |
1082 | ||
1083 | for (Int_t iparam = 0; iparam < dim; iparam++) { | |
1084 | returnFormula.Append(Form("%s*(%f)",((TObjString*)formulaTokens->At(iparam))->GetName(),fitParam[iparam])); | |
1085 | if (iparam < dim-1) returnFormula.Append("+"); | |
1086 | } | |
1087 | returnFormula.Append(" )"); | |
1088 | ||
1089 | ||
b8072cce | 1090 | for (Int_t j=0; j<dim+1;j++) delete [] values[j]; |
cb1d20de | 1091 | |
21f3a443 | 1092 | delete formulaTokens; |
1093 | delete fitter; | |
1094 | delete[] values; | |
b8072cce | 1095 | delete[] errors; |
21f3a443 | 1096 | return preturnFormula; |
1097 | } | |
7c9cf6e4 | 1098 | |
1099 | ||
1100 | ||
1101 | ||
1102 | ||
3d7cc0b4 | 1103 | Int_t TStatToolkit::GetFitIndex(const TString fString, const TString subString){ |
7c9cf6e4 | 1104 | // |
1105 | // fitString - ++ separated list of fits | |
1106 | // substring - ++ separated list of the requiered substrings | |
1107 | // | |
1108 | // return the last occurance of substring in fit string | |
1109 | // | |
1110 | TObjArray *arrFit = fString.Tokenize("++"); | |
1111 | TObjArray *arrSub = subString.Tokenize("++"); | |
1112 | Int_t index=-1; | |
1113 | for (Int_t i=0; i<arrFit->GetEntries(); i++){ | |
1114 | Bool_t isOK=kTRUE; | |
1115 | TString str =arrFit->At(i)->GetName(); | |
1116 | for (Int_t isub=0; isub<arrSub->GetEntries(); isub++){ | |
1117 | if (str.Contains(arrSub->At(isub)->GetName())==0) isOK=kFALSE; | |
1118 | } | |
1119 | if (isOK) index=i; | |
1120 | } | |
09d5920f | 1121 | delete arrFit; |
1122 | delete arrSub; | |
7c9cf6e4 | 1123 | return index; |
1124 | } | |
1125 | ||
1126 | ||
3d7cc0b4 | 1127 | TString TStatToolkit::FilterFit(const TString &input, const TString filter, TVectorD ¶m, TMatrixD & covar){ |
7c9cf6e4 | 1128 | // |
1129 | // Filter fit expression make sub-fit | |
1130 | // | |
1131 | TObjArray *array0= input.Tokenize("++"); | |
1132 | TObjArray *array1= filter.Tokenize("++"); | |
1133 | //TString *presult=new TString("(0"); | |
1134 | TString result="(0.0"; | |
1135 | for (Int_t i=0; i<array0->GetEntries(); i++){ | |
1136 | Bool_t isOK=kTRUE; | |
1137 | TString str(array0->At(i)->GetName()); | |
1138 | for (Int_t j=0; j<array1->GetEntries(); j++){ | |
1139 | if (str.Contains(array1->At(j)->GetName())==0) isOK=kFALSE; | |
1140 | } | |
1141 | if (isOK) { | |
1142 | result+="+"+str; | |
1143 | result+=Form("*(%f)",param[i+1]); | |
1144 | printf("%f\t%f\t%s\n",param[i+1], TMath::Sqrt(covar(i+1,i+1)),str.Data()); | |
1145 | } | |
1146 | } | |
1147 | result+="-0.)"; | |
09d5920f | 1148 | delete array0; |
1149 | delete array1; | |
7c9cf6e4 | 1150 | return result; |
1151 | } | |
1152 | ||
1153 | void TStatToolkit::Update1D(Double_t delta, Double_t sigma, Int_t s1, TMatrixD &vecXk, TMatrixD &covXk){ | |
1154 | // | |
1155 | // Update parameters and covariance - with one measurement | |
1156 | // Input: | |
1157 | // vecXk - input vector - Updated in function | |
1158 | // covXk - covariance matrix - Updated in function | |
1159 | // delta, sigma, s1 - new measurement, rms of new measurement and the index of measurement | |
1160 | const Int_t knMeas=1; | |
1161 | Int_t knElem=vecXk.GetNrows(); | |
1162 | ||
1163 | TMatrixD mat1(knElem,knElem); // update covariance matrix | |
1164 | TMatrixD matHk(1,knElem); // vector to mesurement | |
1165 | TMatrixD vecYk(knMeas,1); // Innovation or measurement residual | |
1166 | TMatrixD matHkT(knElem,knMeas); // helper matrix Hk transpose | |
1167 | TMatrixD matSk(knMeas,knMeas); // Innovation (or residual) covariance | |
1168 | TMatrixD matKk(knElem,knMeas); // Optimal Kalman gain | |
1169 | TMatrixD covXk2(knElem,knElem); // helper matrix | |
1170 | TMatrixD covXk3(knElem,knElem); // helper matrix | |
1171 | TMatrixD vecZk(1,1); | |
1172 | TMatrixD measR(1,1); | |
1173 | vecZk(0,0)=delta; | |
1174 | measR(0,0)=sigma*sigma; | |
1175 | // | |
1176 | // reset matHk | |
1177 | for (Int_t iel=0;iel<knElem;iel++) | |
1178 | for (Int_t ip=0;ip<knMeas;ip++) matHk(ip,iel)=0; | |
1179 | //mat1 | |
1180 | for (Int_t iel=0;iel<knElem;iel++) { | |
1181 | for (Int_t jel=0;jel<knElem;jel++) mat1(iel,jel)=0; | |
1182 | mat1(iel,iel)=1; | |
1183 | } | |
1184 | // | |
1185 | matHk(0, s1)=1; | |
1186 | vecYk = vecZk-matHk*vecXk; // Innovation or measurement residual | |
1187 | matHkT=matHk.T(); matHk.T(); | |
1188 | matSk = (matHk*(covXk*matHkT))+measR; // Innovation (or residual) covariance | |
1189 | matSk.Invert(); | |
1190 | matKk = (covXk*matHkT)*matSk; // Optimal Kalman gain | |
1191 | vecXk += matKk*vecYk; // updated vector | |
1192 | covXk2= (mat1-(matKk*matHk)); | |
1193 | covXk3 = covXk2*covXk; | |
1194 | covXk = covXk3; | |
1195 | Int_t nrows=covXk3.GetNrows(); | |
1196 | ||
1197 | for (Int_t irow=0; irow<nrows; irow++) | |
1198 | for (Int_t icol=0; icol<nrows; icol++){ | |
1199 | // rounding problems - make matrix again symteric | |
1200 | covXk(irow,icol)=(covXk3(irow,icol)+covXk3(icol,irow))*0.5; | |
1201 | } | |
1202 | } | |
1203 | ||
1204 | ||
1205 | ||
3d7cc0b4 | 1206 | void TStatToolkit::Constrain1D(const TString &input, const TString filter, TVectorD ¶m, TMatrixD & covar, Double_t mean, Double_t sigma){ |
7c9cf6e4 | 1207 | // |
1208 | // constrain linear fit | |
1209 | // input - string description of fit function | |
1210 | // filter - string filter to select sub fits | |
1211 | // param,covar - parameters and covariance matrix of the fit | |
1212 | // mean,sigma - new measurement uning which the fit is updated | |
1213 | // | |
ae45c94d | 1214 | |
7c9cf6e4 | 1215 | TObjArray *array0= input.Tokenize("++"); |
1216 | TObjArray *array1= filter.Tokenize("++"); | |
1217 | TMatrixD paramM(param.GetNrows(),1); | |
1218 | for (Int_t i=0; i<=array0->GetEntries(); i++){paramM(i,0)=param(i);} | |
1219 | ||
ae45c94d | 1220 | if (filter.Length()==0){ |
1221 | TStatToolkit::Update1D(mean, sigma, 0, paramM, covar);// | |
1222 | }else{ | |
1223 | for (Int_t i=0; i<array0->GetEntries(); i++){ | |
1224 | Bool_t isOK=kTRUE; | |
1225 | TString str(array0->At(i)->GetName()); | |
1226 | for (Int_t j=0; j<array1->GetEntries(); j++){ | |
1227 | if (str.Contains(array1->At(j)->GetName())==0) isOK=kFALSE; | |
1228 | } | |
1229 | if (isOK) { | |
1230 | TStatToolkit::Update1D(mean, sigma, i+1, paramM, covar);// | |
1231 | } | |
7c9cf6e4 | 1232 | } |
1233 | } | |
1234 | for (Int_t i=0; i<=array0->GetEntries(); i++){ | |
1235 | param(i)=paramM(i,0); | |
1236 | } | |
09d5920f | 1237 | delete array0; |
1238 | delete array1; | |
7c9cf6e4 | 1239 | } |
1240 | ||
ae45c94d | 1241 | TString TStatToolkit::MakeFitString(const TString &input, const TVectorD ¶m, const TMatrixD & covar, Bool_t verbose){ |
7c9cf6e4 | 1242 | // |
1243 | // | |
1244 | // | |
1245 | TObjArray *array0= input.Tokenize("++"); | |
ae45c94d | 1246 | TString result=Form("(%f",param[0]); |
1247 | printf("%f\t%f\t\n", param[0], TMath::Sqrt(covar(0,0))); | |
7c9cf6e4 | 1248 | for (Int_t i=0; i<array0->GetEntries(); i++){ |
1249 | TString str(array0->At(i)->GetName()); | |
1250 | result+="+"+str; | |
1251 | result+=Form("*(%f)",param[i+1]); | |
ae45c94d | 1252 | if (verbose) printf("%f\t%f\t%s\n", param[i+1], TMath::Sqrt(covar(i+1,i+1)),str.Data()); |
7c9cf6e4 | 1253 | } |
1254 | result+="-0.)"; | |
09d5920f | 1255 | delete array0; |
7c9cf6e4 | 1256 | return result; |
1257 | } | |
df0a2a0a | 1258 | |
1a8f4649 | 1259 | TGraphErrors * TStatToolkit::MakeGraphErrors(TTree * tree, const char * expr, const char * cut, Int_t mstyle, Int_t mcolor, Float_t msize, Float_t offset){ |
1260 | // | |
1261 | // Query a graph errors | |
1262 | // return TGraphErrors specified by expr and cut | |
1263 | // Example usage TStatToolkit::MakeGraphError(tree,"Y:X:ErrY","X>0", 25,2,0.4) | |
1264 | // tree - tree with variable | |
1265 | // expr - examp | |
1266 | const Int_t entries = tree->Draw(expr,cut,"goff"); | |
1267 | if (entries<=0) { | |
1268 | TStatToolkit t; | |
1269 | t.Error("TStatToolkit::MakeGraphError",Form("Empty or Not valid expression (%s) or cut *%s)", expr,cut)); | |
1270 | return 0; | |
1271 | } | |
1272 | if ( tree->GetV2()==0){ | |
1273 | TStatToolkit t; | |
1274 | t.Error("TStatToolkit::MakeGraphError",Form("Not valid expression (%s) ", expr)); | |
1275 | return 0; | |
1276 | } | |
1277 | TGraphErrors * graph=0; | |
1278 | if ( tree->GetV3()!=0){ | |
1279 | graph = new TGraphErrors (entries, tree->GetV2(),tree->GetV1(),0,tree->GetV3()); | |
1280 | }else{ | |
1281 | graph = new TGraphErrors (entries, tree->GetV2(),tree->GetV1(),0,0); | |
1282 | } | |
1283 | graph->SetMarkerStyle(mstyle); | |
1284 | graph->SetMarkerColor(mcolor); | |
1285 | graph->SetLineColor(mcolor); | |
1286 | if (msize>0) graph->SetMarkerSize(msize); | |
1287 | for(Int_t i=0;i<graph->GetN();i++) graph->GetX()[i]+=offset; | |
1288 | return graph; | |
1289 | ||
1290 | } | |
1291 | ||
df0a2a0a | 1292 | |
377a7d60 | 1293 | TGraph * TStatToolkit::MakeGraphSparse(TTree * tree, const char * expr, const char * cut, Int_t mstyle, Int_t mcolor, Float_t msize, Float_t offset){ |
df0a2a0a | 1294 | // |
1295 | // Make a sparse draw of the variables | |
d02b8756 | 1296 | // Format of expr : Var:Run or Var:Run:ErrorY or Var:Run:ErrorY:ErrorX |
1297 | // offset : points can slightly be shifted in x for better visibility with more graphs | |
1298 | // | |
1299 | // Written by Weilin.Yu | |
1300 | // updated & merged with QA-code by Patrick Reichelt | |
1301 | // | |
1302 | const Int_t entries = tree->Draw(expr,cut,"goff"); | |
8fb3bea1 | 1303 | if (entries<=0) { |
1304 | TStatToolkit t; | |
d02b8756 | 1305 | t.Error("TStatToolkit::MakeGraphSparse",Form("Empty or Not valid expression (%s) or cut (%s)", expr, cut)); |
8fb3bea1 | 1306 | return 0; |
1307 | } | |
df0a2a0a | 1308 | // TGraph * graph = (TGraph*)gPad->GetPrimitive("Graph"); // 2D |
d02b8756 | 1309 | |
1310 | Double_t *graphY, *graphX; | |
1311 | graphY = tree->GetV1(); | |
1312 | graphX = tree->GetV2(); | |
1313 | ||
1314 | // sort according to run number | |
eda18694 | 1315 | Int_t *index = new Int_t[entries*4]; |
d02b8756 | 1316 | TMath::Sort(entries,graphX,index,kFALSE); |
df0a2a0a | 1317 | |
d02b8756 | 1318 | // define arrays for the new graph |
1319 | Double_t *unsortedX = new Double_t[entries]; | |
1320 | Int_t *runNumber = new Int_t[entries]; | |
df0a2a0a | 1321 | Double_t count = 0.5; |
d02b8756 | 1322 | |
1323 | // evaluate arrays for the new graph according to the run-number | |
baa0041d | 1324 | Int_t icount=0; |
d02b8756 | 1325 | //first entry |
1326 | unsortedX[index[0]] = count; | |
1327 | runNumber[0] = graphX[index[0]]; | |
1328 | // loop the rest of entries | |
1329 | for(Int_t i=1;i<entries;i++) | |
1330 | { | |
1331 | if(graphX[index[i]]==graphX[index[i-1]]) | |
1332 | unsortedX[index[i]] = count; | |
1333 | else if(graphX[index[i]]!=graphX[index[i-1]]){ | |
df0a2a0a | 1334 | count++; |
baa0041d | 1335 | icount++; |
d02b8756 | 1336 | unsortedX[index[i]] = count; |
1337 | runNumber[icount]=graphX[index[i]]; | |
df0a2a0a | 1338 | } |
1339 | } | |
d02b8756 | 1340 | |
1341 | // count the number of xbins (run-wise) for the new graph | |
df0a2a0a | 1342 | const Int_t newNbins = int(count+0.5); |
1343 | Double_t *newBins = new Double_t[newNbins+1]; | |
1344 | for(Int_t i=0; i<=count+1;i++){ | |
1345 | newBins[i] = i; | |
1346 | } | |
d02b8756 | 1347 | |
1348 | // define and fill the new graph | |
b3453fe7 | 1349 | TGraph *graphNew = 0; |
d02b8756 | 1350 | if (tree->GetV3()) { |
1351 | if (tree->GetV4()) { | |
1352 | graphNew = new TGraphErrors(entries,unsortedX,graphY,tree->GetV4(),tree->GetV3()); | |
1353 | } | |
1354 | else { graphNew = new TGraphErrors(entries,unsortedX,graphY,0,tree->GetV3()); } | |
1355 | } | |
1356 | else { graphNew = new TGraphErrors(entries,unsortedX,graphY,0,0); } | |
1357 | // with "Set(...)", the x-axis is being sorted | |
df0a2a0a | 1358 | graphNew->GetXaxis()->Set(newNbins,newBins); |
d02b8756 | 1359 | |
1360 | // set the bins for the x-axis, apply shifting of points | |
df0a2a0a | 1361 | Char_t xName[50]; |
df0a2a0a | 1362 | for(Int_t i=0;i<count;i++){ |
d02b8756 | 1363 | snprintf(xName,50,"%d",runNumber[i]); |
df0a2a0a | 1364 | graphNew->GetXaxis()->SetBinLabel(i+1,xName); |
d02b8756 | 1365 | graphNew->GetX()[i]+=offset; |
df0a2a0a | 1366 | } |
d02b8756 | 1367 | |
df0a2a0a | 1368 | graphNew->GetHistogram()->SetTitle(""); |
d02b8756 | 1369 | graphNew->SetMarkerStyle(mstyle); |
b3453fe7 | 1370 | graphNew->SetMarkerColor(mcolor); |
70989f8d | 1371 | if (msize>0) graphNew->SetMarkerSize(msize); |
d02b8756 | 1372 | delete [] unsortedX; |
1373 | delete [] runNumber; | |
df0a2a0a | 1374 | delete [] index; |
1375 | delete [] newBins; | |
d02b8756 | 1376 | // |
df0a2a0a | 1377 | return graphNew; |
1378 | } | |
1379 | ||
377a7d60 | 1380 | |
1381 | ||
1382 | // | |
d02b8756 | 1383 | // functions used for the trending |
377a7d60 | 1384 | // |
1385 | ||
1386 | Int_t TStatToolkit::MakeStatAlias(TTree * tree, const char * expr, const char * cut, const char * alias) | |
1387 | { | |
1388 | // | |
1389 | // Add alias using statistical values of a given variable. | |
1390 | // (by MI, Patrick Reichelt) | |
1391 | // | |
1392 | // tree - input tree | |
1393 | // expr - variable expression | |
1394 | // cut - selection criteria | |
1395 | // Output - return number of entries used to define variable | |
1396 | // In addition mean, rms, median, and robust mean and rms (choosing fraction of data with smallest RMS) | |
1397 | // | |
d02b8756 | 1398 | /* Example usage: |
1399 | 1.) create the robust estimators for variable expr="QA.TPC.CPass1.meanTPCncl" and create a corresponding | |
1400 | aliases with the prefix alias[0]="ncl", calculated using fraction alias[1]="0.90" | |
1401 | ||
1402 | TStatToolkit::MakeStatAlias(tree,"QA.TPC.CPass1.meanTPCncl","QA.TPC.CPass1.status>0","ncl:0.9"); | |
1403 | root [4] tree->GetListOfAliases().Print() | |
1404 | OBJ: TNamed ncl_Median (130.964333+0) | |
1405 | OBJ: TNamed ncl_Mean (122.120387+0) | |
1406 | OBJ: TNamed ncl_RMS (33.509623+0) | |
1407 | OBJ: TNamed ncl_Mean90 (131.503862+0) | |
1408 | OBJ: TNamed ncl_RMS90 (3.738260+0) | |
1409 | */ | |
377a7d60 | 1410 | // |
d02b8756 | 1411 | Int_t entries = tree->Draw(expr,cut,"goff"); |
377a7d60 | 1412 | if (entries<=1){ |
1413 | printf("Expression or cut not valid:\t%s\t%s\n", expr, cut); | |
1414 | return 0; | |
1415 | } | |
1416 | // | |
1417 | TObjArray* oaAlias = TString(alias).Tokenize(":"); | |
1418 | if (oaAlias->GetEntries()<2) return 0; | |
1419 | Float_t entryFraction = atof( oaAlias->At(1)->GetName() ); | |
1420 | // | |
1421 | Double_t median = TMath::Median(entries,tree->GetV1()); | |
d02b8756 | 1422 | Double_t mean = TMath::Mean(entries,tree->GetV1()); |
377a7d60 | 1423 | Double_t rms = TMath::RMS(entries,tree->GetV1()); |
1424 | Double_t meanEF=0, rmsEF=0; | |
1425 | TStatToolkit::EvaluateUni(entries, tree->GetV1(), meanEF, rmsEF, entries*entryFraction); | |
1426 | // | |
d02b8756 | 1427 | tree->SetAlias(Form("%s_Median",oaAlias->At(0)->GetName()), Form("(%f+0)",median)); |
377a7d60 | 1428 | tree->SetAlias(Form("%s_Mean",oaAlias->At(0)->GetName()), Form("(%f+0)",mean)); |
1429 | tree->SetAlias(Form("%s_RMS",oaAlias->At(0)->GetName()), Form("(%f+0)",rms)); | |
377a7d60 | 1430 | tree->SetAlias(Form("%s_Mean%d",oaAlias->At(0)->GetName(),Int_t(entryFraction*100)), Form("(%f+0)",meanEF)); |
1431 | tree->SetAlias(Form("%s_RMS%d",oaAlias->At(0)->GetName(),Int_t(entryFraction*100)), Form("(%f+0)",rmsEF)); | |
1432 | delete oaAlias; | |
1433 | return entries; | |
1434 | } | |
1435 | ||
1436 | Int_t TStatToolkit::SetStatusAlias(TTree * tree, const char * expr, const char * cut, const char * alias) | |
1437 | { | |
1438 | // | |
1439 | // Add alias to trending tree using statistical values of a given variable. | |
1440 | // (by MI, Patrick Reichelt) | |
1441 | // | |
1442 | // format of expr : varname (e.g. meanTPCncl) | |
1443 | // format of cut : char like in TCut | |
1444 | // format of alias: alias:query:entryFraction(EF) (fraction of entries used for uniformity evaluation) | |
1445 | // e.g.: varname_Out:(abs(varname-meanEF)>6.*rmsEF):0.8 | |
d02b8756 | 1446 | // available internal variables are: 'varname, Median, Mean, MeanEF, RMS, RMSEF' |
1447 | // in the alias, 'varname' will be replaced by its content, and 'EF' by the percentage (e.g. MeanEF -> Mean80) | |
377a7d60 | 1448 | // |
1449 | /* Example usage: | |
d02b8756 | 1450 | 1.) Define robust mean (possible, but easier done with TStatToolkit::MakeStatAlias(...)) |
1451 | TStatToolkit::SetStatusAlias(tree, "meanTPCnclF", "meanTPCnclF>0", "meanTPCnclF_MeanEF:MeanEF:0.80") ; | |
377a7d60 | 1452 | root [10] tree->GetListOfAliases()->Print() |
1453 | Collection name='TList', class='TList', size=1 | |
d02b8756 | 1454 | OBJ: TNamed meanTPCnclF_Mean80 0.899308 |
377a7d60 | 1455 | 2.) create alias outlyers - 6 sigma cut |
d02b8756 | 1456 | TStatToolkit::SetStatusAlias(tree, "meanTPCnclF", "meanTPCnclF>0", "meanTPCnclF_Out:(abs(meanTPCnclF-MeanEF)>6.*RMSEF):0.8") |
1457 | meanTPCnclF_Out ==> (abs(meanTPCnclF-0.899308)>6.*0.016590) | |
377a7d60 | 1458 | 3.) the same functionality as in 2.) |
d02b8756 | 1459 | TStatToolkit::SetStatusAlias(tree, "meanTPCnclF", "meanTPCnclF>0", "varname_Out2:(abs(varname-MeanEF)>6.*RMSEF):0.8") |
1460 | meanTPCnclF_Out2 ==> (abs(meanTPCnclF-0.899308)>6.*0.016590) | |
377a7d60 | 1461 | */ |
1462 | // | |
d02b8756 | 1463 | Int_t entries = tree->Draw(expr,cut,"goff"); |
1464 | if (entries<=1){ | |
1465 | printf("Expression or cut not valid:\t%s\t%s\n", expr, cut); | |
1466 | return 0; | |
1467 | } | |
1468 | // | |
377a7d60 | 1469 | TObjArray* oaVar = TString(expr).Tokenize(":"); |
1470 | char varname[50]; | |
377a7d60 | 1471 | snprintf(varname,50,"%s", oaVar->At(0)->GetName()); |
d02b8756 | 1472 | // |
377a7d60 | 1473 | TObjArray* oaAlias = TString(alias).Tokenize(":"); |
d02b8756 | 1474 | if (oaAlias->GetEntries()<3) return 0; |
377a7d60 | 1475 | Float_t entryFraction = atof( oaAlias->At(2)->GetName() ); |
1476 | // | |
377a7d60 | 1477 | Double_t median = TMath::Median(entries,tree->GetV1()); |
d02b8756 | 1478 | Double_t mean = TMath::Mean(entries,tree->GetV1()); |
377a7d60 | 1479 | Double_t rms = TMath::RMS(entries,tree->GetV1()); |
1480 | Double_t meanEF=0, rmsEF=0; | |
1481 | TStatToolkit::EvaluateUni(entries, tree->GetV1(), meanEF, rmsEF, entries*entryFraction); | |
377a7d60 | 1482 | // |
1483 | TString sAlias( oaAlias->At(0)->GetName() ); | |
1484 | sAlias.ReplaceAll("varname",varname); | |
d02b8756 | 1485 | sAlias.ReplaceAll("MeanEF", Form("Mean%1.0f",entryFraction*100) ); |
1486 | sAlias.ReplaceAll("RMSEF", Form("RMS%1.0f",entryFraction*100) ); | |
377a7d60 | 1487 | TString sQuery( oaAlias->At(1)->GetName() ); |
1488 | sQuery.ReplaceAll("varname",varname); | |
1489 | sQuery.ReplaceAll("MeanEF", Form("%f",meanEF) ); | |
d02b8756 | 1490 | sQuery.ReplaceAll("RMSEF", Form("%f",rmsEF) ); //make sure to replace 'RMSEF' before 'RMS'... |
377a7d60 | 1491 | sQuery.ReplaceAll("Median", Form("%f",median) ); |
d02b8756 | 1492 | sQuery.ReplaceAll("Mean", Form("%f",mean) ); |
377a7d60 | 1493 | sQuery.ReplaceAll("RMS", Form("%f",rms) ); |
d02b8756 | 1494 | printf("define alias:\t%s = %s\n", sAlias.Data(), sQuery.Data()); |
377a7d60 | 1495 | // |
1496 | char query[200]; | |
1497 | char aname[200]; | |
1498 | snprintf(query,200,"%s", sQuery.Data()); | |
1499 | snprintf(aname,200,"%s", sAlias.Data()); | |
1500 | tree->SetAlias(aname, query); | |
d02b8756 | 1501 | delete oaVar; |
1502 | delete oaAlias; | |
377a7d60 | 1503 | return entries; |
1504 | } | |
1505 | ||
1506 | TMultiGraph* TStatToolkit::MakeStatusMultGr(TTree * tree, const char * expr, const char * cut, const char * alias, Int_t igr) | |
1507 | { | |
1508 | // | |
1509 | // Compute a trending multigraph that shows for which runs a variable has outliers. | |
1510 | // (by MI, Patrick Reichelt) | |
1511 | // | |
1512 | // format of expr : varname:xaxis (e.g. meanTPCncl:run) | |
1513 | // format of cut : char like in TCut | |
1514 | // format of alias: (1):(varname_Out==0):(varname_Out)[:(varname_Warning):...] | |
1515 | // in the alias, 'varname' will be replaced by its content (e.g. varname_Out -> meanTPCncl_Out) | |
d02b8756 | 1516 | // note: the aliases 'varname_Out' etc have to be defined by function TStatToolkit::SetStatusAlias(...) |
377a7d60 | 1517 | // counter igr is used to shift the multigraph in y when filling a TObjArray. |
1518 | // | |
1519 | TObjArray* oaVar = TString(expr).Tokenize(":"); | |
d02b8756 | 1520 | if (oaVar->GetEntries()<2) return 0; |
377a7d60 | 1521 | char varname[50]; |
1522 | char var_x[50]; | |
1523 | snprintf(varname,50,"%s", oaVar->At(0)->GetName()); | |
1524 | snprintf(var_x ,50,"%s", oaVar->At(1)->GetName()); | |
d02b8756 | 1525 | // |
377a7d60 | 1526 | TString sAlias(alias); |
1527 | sAlias.ReplaceAll("varname",varname); | |
1528 | TObjArray* oaAlias = TString(sAlias.Data()).Tokenize(":"); | |
d02b8756 | 1529 | if (oaAlias->GetEntries()<3) return 0; |
377a7d60 | 1530 | // |
1531 | char query[200]; | |
1532 | TMultiGraph* multGr = new TMultiGraph(); | |
1533 | Int_t marArr[6] = {24+igr%2, 20+igr%2, 20+igr%2, 20+igr%2, 22, 23}; | |
1534 | Int_t colArr[6] = {kBlack, kBlack, kRed, kOrange, kMagenta, kViolet}; | |
1535 | Double_t sizArr[6] = {1.2, 1.1, 1.0, 1.0, 1, 1}; | |
1536 | const Int_t ngr = oaAlias->GetEntriesFast(); | |
1537 | for (Int_t i=0; i<ngr; i++){ | |
1538 | if (i==2) continue; // the Fatal(Out) graph will be added in the end to be plotted on top! | |
1539 | snprintf(query,200, "%f*(%s-0.5):%s", 1.+igr, oaAlias->At(i)->GetName(), var_x); | |
d02b8756 | 1540 | multGr->Add( (TGraphErrors*) TStatToolkit::MakeGraphSparse(tree,query,cut,marArr[i],colArr[i],sizArr[i]) ); |
377a7d60 | 1541 | } |
1542 | snprintf(query,200, "%f*(%s-0.5):%s", 1.+igr, oaAlias->At(2)->GetName(), var_x); | |
d02b8756 | 1543 | multGr->Add( (TGraphErrors*) TStatToolkit::MakeGraphSparse(tree,query,cut,marArr[2],colArr[2],sizArr[2]) ); |
377a7d60 | 1544 | // |
1545 | multGr->SetName(varname); | |
1546 | multGr->SetTitle(varname); // used for y-axis labels. // details to be included! | |
d02b8756 | 1547 | delete oaVar; |
1548 | delete oaAlias; | |
377a7d60 | 1549 | return multGr; |
1550 | } | |
1551 | ||
1552 | ||
1553 | void TStatToolkit::AddStatusPad(TCanvas* c1, Float_t padratio, Float_t bottommargin) | |
1554 | { | |
1555 | // | |
1556 | // add pad to bottom of canvas for Status graphs (by Patrick Reichelt) | |
1557 | // call function "DrawStatusGraphs(...)" afterwards | |
1558 | // | |
1559 | TCanvas* c1_clone = (TCanvas*) c1->Clone("c1_clone"); | |
1560 | c1->Clear(); | |
1561 | // produce new pads | |
1562 | c1->cd(); | |
1563 | TPad* pad1 = new TPad("pad1", "pad1", 0., padratio, 1., 1.); | |
1564 | pad1->Draw(); | |
1565 | pad1->SetNumber(1); // so it can be called via "c1->cd(1);" | |
1566 | c1->cd(); | |
1567 | TPad* pad2 = new TPad("pad2", "pad2", 0., 0., 1., padratio); | |
1568 | pad2->Draw(); | |
1569 | pad2->SetNumber(2); | |
1570 | // draw original canvas into first pad | |
1571 | c1->cd(1); | |
1572 | c1_clone->DrawClonePad(); | |
1573 | pad1->SetBottomMargin(0.001); | |
1574 | pad1->SetRightMargin(0.01); | |
1575 | // set up second pad | |
1576 | c1->cd(2); | |
1577 | pad2->SetGrid(3); | |
1578 | pad2->SetTopMargin(0); | |
1579 | pad2->SetBottomMargin(bottommargin); // for the long x-axis labels (runnumbers) | |
1580 | pad2->SetRightMargin(0.01); | |
1581 | } | |
1582 | ||
1583 | ||
1584 | void TStatToolkit::DrawStatusGraphs(TObjArray* oaMultGr) | |
1585 | { | |
1586 | // | |
1587 | // draw Status graphs into active pad of canvas (by MI, Patrick Reichelt) | |
1588 | // ...into bottom pad, if called after "AddStatusPad(...)" | |
1589 | // | |
1590 | const Int_t nvars = oaMultGr->GetEntriesFast(); | |
1591 | TGraph* grAxis = (TGraph*) ((TMultiGraph*) oaMultGr->At(0))->GetListOfGraphs()->At(0); | |
1592 | grAxis->SetMaximum(0.5*nvars+0.5); | |
1593 | grAxis->SetMinimum(0); | |
1594 | grAxis->GetYaxis()->SetLabelSize(0); | |
1595 | Int_t entries = grAxis->GetN(); | |
1596 | printf("entries (via GetN()) = %d\n",entries); | |
1597 | grAxis->GetXaxis()->SetLabelSize(5.7*TMath::Min(TMath::Max(5./entries,0.01),0.03)); | |
1598 | grAxis->GetXaxis()->LabelsOption("v"); | |
1599 | grAxis->Draw("ap"); | |
1600 | // | |
1601 | // draw multigraphs & names of status variables on the y axis | |
1602 | for (Int_t i=0; i<nvars; i++){ | |
1603 | ((TMultiGraph*) oaMultGr->At(i))->Draw("p"); | |
1604 | TLatex* ylabel = new TLatex(-0.1, 0.5*i+0.5, ((TMultiGraph*) oaMultGr->At(i))->GetTitle()); | |
1605 | ylabel->SetTextAlign(32); //hor:right & vert:centered | |
1606 | ylabel->SetTextSize(0.025/gPad->GetHNDC()); | |
1607 | ylabel->Draw(); | |
1608 | } | |
1609 | } | |
376089b6 | 1610 | |
1611 | ||
1612 | void TStatToolkit::DrawHistogram(TTree * tree, const char* drawCommand, const char* cuts, const char* histoname, const char* histotitle, Int_t nsigma, Float_t fraction ) | |
1613 | { | |
1614 | // | |
1615 | // Draw histogram from TTree with robust range | |
1616 | // Only for 1D so far! | |
1617 | // | |
1618 | // Parameters: | |
1619 | // - histoname: name of histogram | |
1620 | // - histotitle: title of histgram | |
1621 | // - fraction: fraction of data to define the robust mean | |
1622 | // - nsigma: nsigma value for range | |
1623 | // | |
1624 | ||
1625 | TString drawStr(drawCommand); | |
1626 | TString cutStr(cuts); | |
1627 | Int_t dim = 1; | |
1628 | ||
376089b6 | 1629 | if(!tree) { |
1630 | cerr<<" Tree pointer is NULL!"<<endl; | |
1631 | return; | |
1632 | } | |
1633 | ||
3240a856 | 1634 | // get entries |
376089b6 | 1635 | Int_t entries = tree->Draw(drawStr.Data(), cutStr.Data(), "goff"); |
1636 | if (entries == -1) { | |
1637 | cerr<<"TTree draw returns -1"<<endl; | |
1638 | return; | |
1639 | } | |
1640 | ||
3240a856 | 1641 | // get dimension |
1642 | if(tree->GetV1()) dim = 1; | |
1643 | if(tree->GetV2()) dim = 2; | |
1644 | if(tree->GetV3()) dim = 3; | |
1645 | if(dim > 2){ | |
1646 | cerr<<"TTree has more than 2 dimensions (not yet supported)"<<endl; | |
1647 | return; | |
1648 | } | |
1649 | ||
1650 | // draw robust | |
1651 | Double_t meanX, rmsX=0; | |
1652 | Double_t meanY, rmsY=0; | |
1653 | TStatToolkit::EvaluateUni(entries, tree->GetV1(),meanX,rmsX, fraction*entries); | |
1654 | if(dim==2){ | |
1655 | TStatToolkit::EvaluateUni(entries, tree->GetV1(),meanY,rmsY, fraction*entries); | |
1656 | TStatToolkit::EvaluateUni(entries, tree->GetV2(),meanX,rmsX, fraction*entries); | |
1657 | } | |
1658 | TH1* hOut; | |
1659 | if(dim==1){ | |
1660 | hOut = new TH1F(histoname, histotitle, 200, meanX-nsigma*rmsX, meanX+nsigma*rmsX); | |
1661 | for (Int_t i=0; i<entries; i++) hOut->Fill(tree->GetV1()[i]); | |
1662 | hOut->GetXaxis()->SetTitle(tree->GetHistogram()->GetXaxis()->GetTitle()); | |
1663 | hOut->Draw(); | |
1664 | } | |
1665 | else if(dim==2){ | |
1666 | hOut = new TH2F(histoname, histotitle, 200, meanX-nsigma*rmsX, meanX+nsigma*rmsX,200, meanY-nsigma*rmsY, meanY+nsigma*rmsY); | |
1667 | for (Int_t i=0; i<entries; i++) hOut->Fill(tree->GetV2()[i],tree->GetV1()[i]); | |
1668 | hOut->GetXaxis()->SetTitle(tree->GetHistogram()->GetXaxis()->GetTitle()); | |
1669 | hOut->GetYaxis()->SetTitle(tree->GetHistogram()->GetYaxis()->GetTitle()); | |
1670 | hOut->Draw("colz"); | |
1671 | } | |
376089b6 | 1672 | |
1673 | } |