]>
Commit | Line | Data |
---|---|---|
4c0d7bc7 | 1 | // ---------------------------------------------------------------------- |
2 | // AliBWTools | |
3 | // | |
4 | // This class provides some tools which can be useful in the analsis | |
5 | // of spectra, to fit or transform histograms. See the comments of the | |
6 | // individual methods for details | |
7 | // | |
8 | // Author: M. Floris (CERN) | |
9 | // ---------------------------------------------------------------------- | |
10 | ||
11 | #include "AliBWTools.h" | |
12 | #include "TH1D.h" | |
13 | #include "TF1.h" | |
14 | #include "TH1.h" | |
15 | #include "TMath.h" | |
16 | #include "TGraphErrors.h" | |
17 | #include "TVirtualFitter.h" | |
18 | #include "TMinuit.h" | |
19 | #include "AliLog.h" | |
20 | #include <iostream> | |
21 | ||
22 | using namespace std; | |
23 | ||
4c0d7bc7 | 24 | |
25 | ClassImp(AliBWTools) | |
26 | ||
27 | AliBWTools::AliBWTools() { | |
28 | // ctor | |
29 | } | |
30 | ||
31 | AliBWTools::~AliBWTools(){ | |
32 | // dtor | |
33 | } | |
34 | ||
3beb22a2 | 35 | TH1 * AliBWTools::GetdNdmtFromdNdpt(const TH1 * hpt, Double_t mass) { |
4c0d7bc7 | 36 | // convert the x axis from pt to mt. Assumes you have 1/pt dNdpt in the histo you start with |
37 | ||
38 | Int_t nbins = hpt->GetNbinsX(); | |
39 | Float_t * xbins = new Float_t[nbins+1]; | |
40 | for(Int_t ibins = 0; ibins <= nbins; ibins++){ | |
41 | xbins[ibins] = TMath::Sqrt(hpt->GetBinLowEdge(ibins+1)*hpt->GetBinLowEdge(ibins+1) + | |
0e019ec6 | 42 | mass *mass) - mass; |
77233d47 | 43 | // // xbins[ibins] = TMath::Sqrt(hpt->GetBinLowEdge(ibins+1)*hpt->GetBinLowEdge(ibins+1) + |
44 | // // mass *mass); | |
4c0d7bc7 | 45 | // cout << ibins << " "<< xbins[ibins] << endl; |
46 | ||
47 | } | |
48 | ||
49 | TH1D * hmt = new TH1D(TString(hpt->GetName())+"_mt", | |
50 | TString(hpt->GetName())+"_mt", | |
51 | nbins, xbins); | |
52 | for(Int_t ibins = 1; ibins <= nbins; ibins++){ | |
53 | hmt->SetBinContent(ibins, hpt->GetBinContent(ibins)); | |
54 | hmt->SetBinError(ibins, hpt->GetBinError(ibins)); | |
55 | ||
56 | } | |
57 | ||
58 | hmt->SetXTitle("m_{t} - m_{0} (GeV/c^{2})"); | |
59 | hmt->SetYTitle("1/m_{t} dN/dm_{t} (a.u.)"); | |
0e019ec6 | 60 | hmt->SetMarkerStyle(hpt->GetMarkerStyle()); |
61 | hmt->SetMarkerColor(hpt->GetMarkerColor()); | |
62 | hmt->SetLineColor(hpt->GetLineColor()); | |
63 | ||
4c0d7bc7 | 64 | return hmt; |
65 | ||
66 | } | |
67 | ||
9ca56c22 | 68 | TH1 * AliBWTools::GetdNdptFromdNdmt(const TH1 * hmt, Double_t mass) { |
69 | // convert the x axis from mt to pt. Assumes you have 1/mt dNdmt in the histo you start with | |
70 | ||
71 | Int_t nbins = hmt->GetNbinsX(); | |
72 | Float_t * xbins = new Float_t[nbins+1]; | |
73 | for(Int_t ibins = 0; ibins <= nbins; ibins++){ | |
74 | xbins[ibins] = TMath::Sqrt((hmt->GetBinLowEdge(ibins+1)+mass)*(hmt->GetBinLowEdge(ibins+1)+mass) - | |
75 | mass *mass); | |
76 | xbins[ibins] = Float_t(TMath::Nint(xbins[ibins]*100))/100; | |
77 | // // xbins[ibins] = TMath::Sqrt(hmt->GetBinLowEdge(ibins+1)*hmt->GetBinLowEdge(ibins+1) + | |
78 | // // mass *mass); | |
79 | cout << ibins << " "<< xbins[ibins] << endl; | |
80 | ||
81 | } | |
82 | ||
83 | TH1D * hptL = new TH1D(TString(hmt->GetName())+"_pt", | |
84 | TString(hmt->GetName())+"_pt", | |
85 | nbins, xbins); | |
86 | for(Int_t ibins = 1; ibins <= nbins; ibins++){ | |
87 | hptL->SetBinContent(ibins, hmt->GetBinContent(ibins)); | |
88 | hptL->SetBinError(ibins, hmt->GetBinError(ibins)); | |
89 | ||
90 | } | |
91 | ||
92 | hptL->SetXTitle("p_{t} (GeV/c)"); | |
93 | hptL->SetYTitle("1/p_{t} dN/dp_{t} (a.u.)"); | |
94 | hptL->SetMarkerStyle(hmt->GetMarkerStyle()); | |
95 | hptL->SetMarkerColor(hmt->GetMarkerColor()); | |
96 | hptL->SetLineColor(hmt->GetLineColor()); | |
97 | ||
98 | return hptL; | |
99 | ||
100 | } | |
101 | ||
102 | ||
3beb22a2 | 103 | TH1 * AliBWTools::GetdNdPtFromOneOverPt(const TH1 * h1Pt) { |
4c0d7bc7 | 104 | |
105 | // convert an histo from 1/pt dNdpt to dNdpt, using the pt at the center of the bin | |
106 | ||
107 | ||
108 | TH1 * hPt = (TH1 *) h1Pt->Clone((TString(h1Pt->GetName()) + "_inv").Data()); | |
109 | hPt->Reset(); | |
110 | ||
111 | Int_t nbinx = hPt->GetNbinsX(); | |
112 | ||
113 | for(Int_t ibinx = 1; ibinx <= nbinx; ibinx++){ | |
114 | ||
115 | Double_t cont = h1Pt->GetBinContent(ibinx); | |
116 | Double_t err = h1Pt->GetBinError(ibinx); | |
117 | ||
118 | Double_t pt = h1Pt->GetBinCenter(ibinx); | |
119 | ||
c3c47c65 | 120 | if(pt > 0) { |
4c0d7bc7 | 121 | cont *= pt; |
122 | err *= pt; | |
123 | } else { | |
124 | cont = 0; | |
125 | err = 0; | |
126 | } | |
127 | ||
128 | hPt->SetBinContent(ibinx, cont); | |
129 | hPt->SetBinError(ibinx, err); | |
130 | ||
131 | } | |
132 | ||
133 | hPt->SetXTitle("p_{t} (GeV)"); | |
134 | hPt->SetYTitle("dN/dp_{t} (GeV^{-2})"); | |
135 | ||
136 | return hPt; | |
137 | ||
138 | } | |
139 | ||
140 | ||
141 | ||
142 | ||
3beb22a2 | 143 | TH1 * AliBWTools::GetOneOverPtdNdPt(const TH1 * hPt) { |
4c0d7bc7 | 144 | |
145 | // convert an histo from dNdpt to 1/pt dNdpt, using the pt at the center of the bin | |
146 | ||
147 | TH1 * h1Pt = (TH1 *) hPt->Clone((TString(hPt->GetName()) + "_inv").Data()); | |
148 | h1Pt->Reset(); | |
149 | ||
150 | Int_t nbinx = h1Pt->GetNbinsX(); | |
151 | ||
152 | for(Int_t ibinx = 1; ibinx <= nbinx; ibinx++){ | |
153 | ||
154 | Double_t cont = hPt->GetBinContent(ibinx); | |
155 | Double_t err = hPt->GetBinError(ibinx); | |
156 | ||
157 | Double_t pt = hPt->GetBinCenter(ibinx); | |
158 | ||
c3c47c65 | 159 | if(pt > 0) { |
4c0d7bc7 | 160 | cont /= pt; |
161 | err /= pt; | |
162 | } else { | |
163 | cont = 0; | |
164 | err = 0; | |
165 | } | |
166 | ||
167 | h1Pt->SetBinContent(ibinx, cont); | |
168 | h1Pt->SetBinError(ibinx, err); | |
169 | ||
170 | } | |
171 | ||
172 | h1Pt->SetXTitle("p_{t} (GeV)"); | |
173 | h1Pt->SetYTitle("1/p_{t} dN/dp_{t} (GeV^{-2})"); | |
174 | ||
175 | return h1Pt; | |
176 | ||
177 | } | |
178 | ||
179 | ||
3beb22a2 | 180 | TGraphErrors * AliBWTools::GetGraphFromHisto(const TH1F * h, Bool_t binWidth) { |
4c0d7bc7 | 181 | // Convert a histo to a graph |
182 | // if binWidth is true ex is set to the bin width of the histos, otherwise it is set to zero | |
183 | Int_t nbin = h->GetNbinsX(); | |
184 | ||
185 | TGraphErrors * g = new TGraphErrors(); | |
186 | Int_t ipoint = 0; | |
187 | for(Int_t ibin = 1; ibin <= nbin; ibin++){ | |
188 | Double_t xerr = binWidth ? h->GetBinWidth(ibin)/2 : 0; | |
3beb22a2 | 189 | if (h->GetBinContent(ibin)) { |
4c0d7bc7 | 190 | g->SetPoint (ipoint, h->GetBinCenter(ibin), h->GetBinContent(ibin)); |
191 | g->SetPointError(ipoint, xerr, h->GetBinError(ibin)); | |
192 | ipoint++; | |
193 | } | |
194 | } | |
195 | ||
196 | g->SetMarkerStyle(h->GetMarkerStyle()); | |
197 | g->SetMarkerColor(h->GetMarkerColor()); | |
198 | g->SetLineColor(h->GetLineColor()); | |
2210b002 | 199 | g->SetLineStyle(h->GetLineStyle()); |
200 | g->SetLineWidth(h->GetLineWidth()); | |
4c0d7bc7 | 201 | |
202 | g->SetTitle(h->GetTitle()); | |
203 | g->SetName(TString("g_")+h->GetName()); | |
204 | ||
205 | return g; | |
206 | ||
207 | } | |
208 | ||
3beb22a2 | 209 | TH1F * AliBWTools::GetHistoFromGraph(const TGraphErrors * g, const TH1F* hTemplate) { |
4c0d7bc7 | 210 | |
211 | // convert a graph to histo with the binning of hTemplate. | |
212 | // Warning: the template should be chosen | |
213 | // properly: if you have two graph points in the same histo bin this | |
214 | // won't work! | |
215 | ||
216 | TH1F * h = (TH1F*) hTemplate->Clone(TString("h_")+g->GetName()); | |
217 | h->Reset(); | |
218 | Int_t npoint = g->GetN(); | |
219 | // g->Print(); | |
220 | for(Int_t ipoint = 0; ipoint < npoint; ipoint++){ | |
221 | Float_t x = g->GetX() [ipoint]; | |
222 | Float_t y = g->GetY() [ipoint]; | |
223 | Float_t ey = g->GetEY()[ipoint]; | |
224 | Int_t bin = h->FindBin(x); | |
225 | // cout << "bin: "<< bin << " -> " << x << ", "<< y <<", " << ey << endl; | |
226 | ||
227 | h->SetBinContent(bin,y); | |
228 | h->SetBinError (bin,ey); | |
229 | } | |
230 | ||
231 | h->SetMarkerStyle(g->GetMarkerStyle()); | |
232 | h->SetMarkerColor(g->GetMarkerColor()); | |
233 | h->SetLineColor (g->GetLineColor()); | |
234 | ||
235 | ||
236 | return h; | |
237 | } | |
238 | ||
3beb22a2 | 239 | TGraphErrors * AliBWTools::ConcatenateGraphs(const TGraphErrors * g1,const TGraphErrors * g2){ |
4c0d7bc7 | 240 | |
241 | // concatenates two graphs | |
242 | ||
243 | Int_t npoint1=g1->GetN(); | |
244 | Int_t npoint2=g2->GetN(); | |
245 | ||
246 | TGraphErrors * gClone = (TGraphErrors*) g1->Clone(); | |
247 | ||
248 | // for(Int_t ipoint = 0; ipoint < npoint1; ipoint++){ | |
249 | // gClone->SetPointError(ipoint,0,g1->GetEY()[ipoint]); | |
250 | ||
251 | // } | |
252 | for(Int_t ipoint = 0; ipoint < npoint2; ipoint++){ | |
253 | gClone->SetPoint(ipoint+npoint1,g2->GetX()[ipoint],g2->GetY()[ipoint]); | |
254 | gClone->SetPointError(ipoint+npoint1,g2->GetEX()[ipoint],g2->GetEY()[ipoint]); | |
255 | // gClone->SetPointError(ipoint+npoint1,0,g2->GetEY()[ipoint]); | |
256 | } | |
257 | ||
258 | gClone->GetHistogram()->GetXaxis()->SetTimeDisplay(1); | |
259 | gClone->SetTitle(TString(gClone->GetTitle())+" + "+g2->GetTitle()); | |
260 | gClone->SetName(TString(gClone->GetName())+"_"+g2->GetName()); | |
261 | ||
262 | return gClone; | |
263 | } | |
264 | ||
265 | ||
3beb22a2 | 266 | TH1F * AliBWTools::Combine3HistosWithErrors(const TH1 * h1, const TH1 * h2, const TH1* h3, |
3725266c | 267 | TH1 * he1, TH1 * he2, TH1 * he3, |
3beb22a2 | 268 | const TH1* htemplate, Int_t statFrom, |
3725266c | 269 | Float_t renorm1, Float_t renorm2, Float_t renorm3, |
270 | TH1 ** hSyst, Bool_t errorFromBinContent) { | |
ee08b77d | 271 | |
272 | // Combines 3 histos (h1,h2,h3), weighting by the errors provided in | |
273 | // he1,he2,he3, supposed to be the independent systematic errors. | |
274 | // he1,he2,he3 are also assumed to have the same binning as h1,h2,h3 | |
275 | // The combined histo must fit the template provided (no check is performed on this) | |
276 | // The histogram are supposed to come from the same (nearly) sample | |
277 | // of tracks. statFrom tells the stat error of which of the 3 is | |
278 | // going to be assigned to the combined | |
279 | // Optionally, it is possible to rescale any of the histograms. | |
3725266c | 280 | // if hSyst is give, the histo is filled with combined syst error vs pt |
281 | // if errorFromBinContent is true, the weights are taken from the he* content rather than errors | |
ee08b77d | 282 | TH1F * hcomb = (TH1F*) htemplate->Clone(TString("hComb_")+h1->GetName()+"_"+h2->GetName()+h3->GetName()); |
283 | ||
284 | // TODO: I should have used an array for h*local... | |
285 | ||
286 | // Clone histos locally to rescale them | |
287 | TH1F * h1local = (TH1F*) h1->Clone(); | |
288 | TH1F * h2local = (TH1F*) h2->Clone(); | |
289 | TH1F * h3local = (TH1F*) h3->Clone(); | |
290 | h1local->Scale(renorm1); | |
291 | h2local->Scale(renorm2); | |
292 | h3local->Scale(renorm3); | |
293 | ||
3beb22a2 | 294 | const TH1 * hStatError = 0; |
ee08b77d | 295 | if (statFrom == 0) hStatError = h1; |
296 | else if (statFrom == 1) hStatError = h2; | |
297 | else if (statFrom == 2) hStatError = h3; | |
37e1d042 | 298 | else { |
299 | Printf("AliBWTools::Combine3HistosWithErrors: wrong value of the statFrom parameter"); | |
300 | return NULL; | |
301 | } | |
ee08b77d | 302 | Printf("AliBWTools::Combine3HistosWithErrors: improve error on combined"); |
303 | // Loop over all bins and take weighted mean of all points | |
304 | Int_t nBinComb = hcomb->GetNbinsX(); | |
305 | for(Int_t ibin = 1; ibin <= nBinComb; ibin++){ | |
306 | Int_t ibin1 = h1local->FindBin(hcomb->GetBinCenter(ibin)); | |
307 | Int_t ibin2 = h2local->FindBin(hcomb->GetBinCenter(ibin)); | |
308 | Int_t ibin3 = h3local->FindBin(hcomb->GetBinCenter(ibin)); | |
309 | Int_t ibinError = -1; // bin used to get stat error | |
310 | ||
311 | if (statFrom == 0) ibinError = ibin1; | |
312 | else if (statFrom == 1) ibinError = ibin2; | |
313 | else if (statFrom == 2) ibinError = ibin3; | |
314 | else Printf("AliBWTools::Combine3HistosWithErrors: wrong value of the statFrom parameter"); | |
315 | ||
316 | ||
317 | Double_t y[3]; | |
318 | Double_t ye[3]; | |
319 | y[0] = h1local->GetBinContent(ibin1); | |
320 | y[1] = h2local->GetBinContent(ibin2); | |
321 | y[2] = h3local->GetBinContent(ibin3); | |
3725266c | 322 | if (errorFromBinContent) { |
323 | ye[0] = he1->GetBinContent(he1->FindBin(hcomb->GetBinCenter(ibin))); | |
324 | ye[1] = he2->GetBinContent(he2->FindBin(hcomb->GetBinCenter(ibin))); | |
325 | ye[2] = he3->GetBinContent(he3->FindBin(hcomb->GetBinCenter(ibin))); | |
326 | } else { | |
327 | ye[0] = he1->GetBinError(ibin1); | |
328 | ye[1] = he2->GetBinError(ibin2); | |
329 | ye[2] = he3->GetBinError(ibin3); | |
330 | } | |
ee08b77d | 331 | // Set error to 0 if content is 0 (means it was not filled) |
3beb22a2 | 332 | if(!h1local->GetBinContent(ibin1)) ye[0] = 0; |
333 | if(!h2local->GetBinContent(ibin2)) ye[1] = 0; | |
334 | if(!h3local->GetBinContent(ibin3)) ye[2] = 0; | |
ee08b77d | 335 | |
336 | // Compute weighted mean | |
3725266c | 337 | // cout << "Bins: "<< ibin1 << " " << ibin2 << " " << ibin3 << endl; |
ee08b77d | 338 | Double_t mean, err; |
339 | WeightedMean(3,y,ye,mean,err); | |
340 | ||
341 | ||
342 | // Fill combined | |
ee08b77d | 343 | hcomb->SetBinContent(ibin,mean); |
344 | Double_t statError = 0; | |
3beb22a2 | 345 | if (hStatError->GetBinContent(ibinError)) { |
ee08b77d | 346 | // cout << "def" << endl; |
347 | statError = hStatError->GetBinError(ibinError)/hStatError->GetBinContent(ibinError)*hcomb->GetBinContent(ibin); | |
348 | } | |
3beb22a2 | 349 | else if (h1local->GetBinContent(ibin1)) { |
ee08b77d | 350 | // cout << "1" << endl; |
351 | statError = h1local->GetBinError(ibin1)/h1local->GetBinContent(ibin1)*hcomb->GetBinContent(ibin); | |
352 | } | |
3beb22a2 | 353 | else if (h2local->GetBinContent(ibin2)) { |
ee08b77d | 354 | // cout << "2" << endl; |
355 | statError = h2local->GetBinError(ibin2)/h2local->GetBinContent(ibin2)*hcomb->GetBinContent(ibin); | |
356 | } | |
3beb22a2 | 357 | else if (h3local->GetBinContent(ibin3)) { |
ee08b77d | 358 | // cout << "3" << endl; |
359 | statError = h3local->GetBinError(ibin3)/h3local->GetBinContent(ibin3)*hcomb->GetBinContent(ibin); | |
360 | } | |
361 | hcomb->SetBinError (ibin,statError); | |
3725266c | 362 | if(hSyst) (*hSyst)->SetBinContent(ibin,err); |
ee08b77d | 363 | // cout << "BIN " << ibin << " " << mean << " " << statError << endl; |
364 | ||
365 | } | |
366 | ||
367 | hcomb->SetMarkerStyle(hStatError->GetMarkerStyle()); | |
368 | hcomb->SetMarkerColor(hStatError->GetMarkerColor()); | |
369 | hcomb->SetLineColor (hStatError->GetLineColor()); | |
370 | ||
371 | hcomb->SetXTitle(hStatError->GetXaxis()->GetTitle()); | |
372 | hcomb->SetYTitle(hStatError->GetYaxis()->GetTitle()); | |
373 | ||
374 | delete h1local; | |
375 | delete h2local; | |
376 | delete h3local; | |
377 | ||
378 | return hcomb; | |
379 | } | |
380 | ||
381 | ||
ddecf727 | 382 | void AliBWTools::GetMeanDataAndExtrapolation(const TH1 * hData, TF1 * fExtrapolation, Double_t &mean, Double_t &error, Float_t min, Float_t max){ |
383 | // Computes the mean of the combined data and extrapolation in a | |
384 | // given range, use data where they are available and the function | |
385 | // where data are not available | |
386 | // ERROR is from DATA ONLY returned in this version! | |
387 | // | |
388 | Printf("AliBWTools::GetMeanDataAndExtrapolation: WARNING from data only"); | |
389 | Float_t minData = GetLowestNotEmptyBinEdge (hData); | |
58eeda84 | 390 | Int_t minDataBin = GetLowestNotEmptyBin (hData); |
ddecf727 | 391 | Float_t maxData = GetHighestNotEmptyBinEdge(hData); |
58eeda84 | 392 | Int_t maxDataBin = GetHighestNotEmptyBin (hData); |
ddecf727 | 393 | Double_t integral = 0; |
394 | mean = 0; | |
395 | error = 0; | |
396 | if (min < minData) { | |
397 | // Compute integral exploiting root function to calculate moments, "unnormalizing" them | |
398 | mean += fExtrapolation->Mean(min,minData)*fExtrapolation->Integral(min,minData); | |
399 | integral += fExtrapolation->Integral(min,minData); | |
400 | cout << " Low "<< mean << " " << integral << endl; | |
401 | ||
402 | } | |
403 | ||
404 | if (max > maxData) { | |
405 | // Compute integral exploiting root function to calculate moments, "unnormalizing" them | |
406 | mean += fExtrapolation->Mean(maxData,max)*fExtrapolation->Integral(maxData,max); | |
407 | integral += fExtrapolation->Integral(maxData,max); | |
408 | cout << " Hi "<< mean << " " << integral << endl; | |
409 | } | |
410 | Float_t err2 = 0; | |
411 | ||
412 | for(Int_t ibin = minDataBin; ibin <= maxDataBin; ibin++){ | |
413 | if(hData->GetBinCenter(ibin) < min) continue; | |
414 | if(hData->GetBinCenter(ibin) > max) continue; | |
415 | mean = mean + (hData->GetBinCenter(ibin) * hData->GetBinWidth(ibin)* hData->GetBinContent(ibin)); | |
416 | err2 = err2 + TMath::Power(hData->GetBinError(ibin) * hData->GetBinCenter(ibin) * hData->GetBinWidth(ibin),2); | |
417 | integral = integral + hData->GetBinContent(ibin) * hData->GetBinWidth(ibin); | |
418 | } | |
419 | cout << " Data "<< mean << " " << integral << endl; | |
420 | ||
421 | mean = mean/integral; | |
ffcaf7be | 422 | error = TMath::Sqrt(err2)/integral; |
ddecf727 | 423 | |
4c0d7bc7 | 424 | |
ddecf727 | 425 | } |
426 | ||
427 | TH1F * AliBWTools::CombineHistos(const TH1 * h1, TH1 * h2, const TH1* htemplate, Float_t renorm1){ | |
4c0d7bc7 | 428 | // Combine two histos. This assumes the histos have the same binning |
429 | // in the overlapping region. It computes the arithmetic mean in the | |
430 | // overlapping region and assigns as an error the relative error | |
431 | // h2. TO BE IMPROVED | |
432 | ||
433 | TH1F * hcomb = (TH1F*) htemplate->Clone(TString(h1->GetName())+"_"+h2->GetName()); | |
434 | ||
435 | TH1F * h1local = (TH1F*) h1->Clone(); | |
436 | h1local->Scale(renorm1); | |
437 | ||
438 | Int_t nBinComb = hcomb->GetNbinsX(); | |
439 | for(Int_t ibin = 1; ibin <= nBinComb; ibin++){ | |
440 | Int_t ibin1 = h1local->FindBin(hcomb->GetBinCenter(ibin)); | |
441 | Int_t ibin2 = h2->FindBin(hcomb->GetBinCenter(ibin)); | |
442 | ||
3beb22a2 | 443 | if (!h1local->GetBinContent(ibin1) && !h2->GetBinContent(ibin2) ) { |
4c0d7bc7 | 444 | // None has data: go to next bin |
445 | hcomb->SetBinContent(ibin,0); | |
446 | hcomb->SetBinError (ibin,0); | |
3beb22a2 | 447 | } else if(h1local->GetBinContent(ibin1) && !h2->GetBinContent(ibin2)) { |
4c0d7bc7 | 448 | // take data from h1local: |
449 | hcomb->SetBinContent(ibin,h1local->GetBinContent(ibin1)); | |
450 | hcomb->SetBinError (ibin,h1local->GetBinError(ibin1)); | |
3beb22a2 | 451 | } else if(!h1local->GetBinContent(ibin1) && h2->GetBinContent(ibin2)) { |
4c0d7bc7 | 452 | // take data from h2: |
453 | hcomb->SetBinContent(ibin,h2->GetBinContent(ibin2)); | |
454 | hcomb->SetBinError (ibin,h2->GetBinError(ibin2)); | |
455 | } else { | |
456 | hcomb->SetBinContent(ibin,(h1local->GetBinContent(ibin1) +h2->GetBinContent(ibin2))/2); | |
457 | // hcomb->SetBinError (ibin,h1local->GetBinError(ibin1)/h1local->GetBinContent(ibin1)*hcomb->GetBinContent(ibin)); | |
458 | hcomb->SetBinError (ibin,h2->GetBinError(ibin2)/h2->GetBinContent(ibin2)*hcomb->GetBinContent(ibin)); | |
459 | } | |
460 | ||
461 | ||
462 | } | |
463 | ||
464 | ||
465 | hcomb->SetMarkerStyle(h1local->GetMarkerStyle()); | |
466 | hcomb->SetMarkerColor(h1local->GetMarkerColor()); | |
467 | hcomb->SetLineColor (h1local->GetLineColor()); | |
468 | ||
469 | hcomb->SetXTitle(h1local->GetXaxis()->GetTitle()); | |
470 | hcomb->SetYTitle(h1local->GetYaxis()->GetTitle()); | |
471 | delete h1local; | |
472 | return hcomb; | |
473 | ||
474 | } | |
475 | ||
476 | ||
3beb22a2 | 477 | void AliBWTools::GetFromHistoGraphDifferentX(const TH1F * h, TF1 * f, TGraphErrors ** gBarycentre, TGraphErrors ** gXlw) { |
4c0d7bc7 | 478 | |
479 | // Computes the baycentre in each bin and the xlw as defined in NIMA | |
480 | // 355 - 541 using f. Returs 2 graphs with the same y content of h | |
481 | // but with a different x (barycentre and xlw) | |
482 | ||
483 | Int_t nbin = h->GetNbinsX(); | |
484 | ||
485 | (*gBarycentre) = new TGraphErrors(); | |
486 | (*gXlw) = new TGraphErrors(); | |
487 | ||
488 | Int_t ipoint = 0; | |
489 | for(Int_t ibin = 1; ibin <= nbin; ibin++){ | |
490 | Float_t min = h->GetBinLowEdge(ibin); | |
491 | Float_t max = h->GetBinLowEdge(ibin+1); | |
492 | Double_t xerr = 0; | |
493 | Double_t xbar = f->Mean(min,max); | |
c3c47c65 | 494 | // compute xLW |
4c0d7bc7 | 495 | Double_t temp = 1./(max-min) * f->Integral(min,max); |
496 | Double_t epsilon = 0.000000001; | |
497 | Double_t increment = 0.0000000001; | |
c3c47c65 | 498 | Double_t xLW = min; |
4c0d7bc7 | 499 | |
c3c47c65 | 500 | while ((f->Eval(xLW)- temp) > epsilon) { |
501 | xLW += increment; | |
502 | if(xLW > max) { | |
503 | Printf("Cannot find xLW"); | |
4c0d7bc7 | 504 | break; |
505 | } | |
506 | } | |
507 | ||
508 | if (h->GetBinContent(ibin)!=0) { | |
509 | (*gBarycentre)->SetPoint (ipoint, xbar, h->GetBinContent(ibin)); | |
510 | (*gBarycentre)->SetPointError(ipoint, xerr, h->GetBinError(ibin)); | |
c3c47c65 | 511 | (*gXlw) ->SetPoint (ipoint, xLW, h->GetBinContent(ibin)); |
4c0d7bc7 | 512 | (*gXlw) ->SetPointError(ipoint, xerr, h->GetBinError(ibin)); |
513 | ipoint++; | |
514 | } | |
515 | } | |
516 | ||
517 | (*gBarycentre)->SetMarkerStyle(h->GetMarkerStyle()); | |
518 | (*gBarycentre)->SetMarkerColor(h->GetMarkerColor()); | |
519 | (*gBarycentre)->SetLineColor(h->GetLineColor()); | |
520 | ||
521 | (*gBarycentre)->SetTitle(h->GetTitle()); | |
522 | (*gBarycentre)->SetName(TString("g_")+h->GetName()); | |
523 | ||
524 | (*gXlw)->SetMarkerStyle(h->GetMarkerStyle()); | |
525 | (*gXlw)->SetMarkerColor(h->GetMarkerColor()); | |
526 | (*gXlw)->SetLineColor(h->GetLineColor()); | |
527 | (*gXlw)->SetTitle(h->GetTitle()); | |
528 | (*gXlw)->SetName(TString("g_")+h->GetName()); | |
529 | ||
530 | ||
531 | } | |
532 | ||
533 | ||
ddecf727 | 534 | Float_t AliBWTools::GetMean(TH1F * h, Float_t min, Float_t max, Float_t * error) { |
4c0d7bc7 | 535 | |
3472a634 | 536 | // Get the mean of histo in a range; root is not reliable in sub |
537 | // ranges with variable binning. | |
4c0d7bc7 | 538 | Int_t minBin = h->FindBin(min); |
ddecf727 | 539 | Int_t maxBin = h->FindBin(max-0.00001); |
4c0d7bc7 | 540 | |
541 | Float_t mean = 0 ; | |
542 | Float_t integral = 0; | |
ddecf727 | 543 | Float_t err2 = 0; |
544 | for(Int_t ibin = minBin; ibin <= maxBin; ibin++){ | |
4c0d7bc7 | 545 | mean = mean + (h->GetBinCenter(ibin) * h->GetBinWidth(ibin)* h->GetBinContent(ibin)); |
ddecf727 | 546 | err2 = err2 + TMath::Power(h->GetBinError(ibin) * h->GetBinCenter(ibin) * h->GetBinWidth(ibin),2); |
4c0d7bc7 | 547 | integral = integral + h->GetBinContent(ibin) * h->GetBinWidth(ibin); |
548 | } | |
549 | ||
ddecf727 | 550 | Float_t value = mean/integral; |
551 | if (error) (*error) = TMath::Sqrt(err2); | |
552 | return value; | |
4c0d7bc7 | 553 | |
554 | ||
555 | } | |
556 | ||
3472a634 | 557 | void AliBWTools::GetMean(TF1 * func, Float_t &mean, Float_t &error, Float_t min, Float_t max, Int_t normPar) { |
4c0d7bc7 | 558 | |
3472a634 | 559 | // Get the mean of function in a range; If normPar is >= 0, it means |
560 | // that the function is defined such that par[normPar] is its | |
561 | // integral. In this case the error on meanpt can be calculated | |
562 | // correctly. Otherwise, the function is normalized in get moment, | |
563 | // but the error is not computed correctly. | |
4c0d7bc7 | 564 | |
3472a634 | 565 | return GetMoment("fMean", TString("x*")+func->GetExpFormula(), func, mean, error, min, max, normPar); |
4c0d7bc7 | 566 | |
567 | } | |
568 | ||
3472a634 | 569 | void AliBWTools::GetMeanSquare(TF1 * func, Float_t &mean, Float_t &error, Float_t min, Float_t max, Int_t normPar) { |
4c0d7bc7 | 570 | |
3472a634 | 571 | // Get the mean2 of function in a range; If normPar is >= 0, it means |
572 | // that the function is defined such that par[normPar] is its | |
573 | // integral. In this case the error on meanpt can be calculated | |
574 | // correctly. Otherwise, the function is normalized in get moment, | |
575 | // but the error is not computed correctly. | |
4c0d7bc7 | 576 | |
3472a634 | 577 | return GetMoment("fMean2", TString("x*x*")+func->GetExpFormula(), func, mean, error, min, max, normPar); |
4c0d7bc7 | 578 | |
579 | ||
580 | } | |
581 | ||
3472a634 | 582 | void AliBWTools::GetMoment(TString name, TString var, TF1 * func, Float_t &mean, Float_t &error, Float_t min, Float_t max, Int_t normPar) { |
4c0d7bc7 | 583 | |
584 | // returns the integral of a function derived from func by prefixing some operation. | |
3472a634 | 585 | // the derived function MUST have the same parameter in the same order |
4c0d7bc7 | 586 | // Used as a base method for mean and mean2 |
3472a634 | 587 | // If normPar is >= 0, it means that the function is defined such |
588 | // that par[normPar] is its integral. In this case the error on | |
589 | // meanpt can be calculated correctly. Otherwise, the function is | |
590 | // normalized using its numerical integral, but the error is not computed | |
591 | // correctly. | |
592 | ||
593 | // TODO: | |
594 | // - improve to propagate error also in the case you need the | |
595 | // numerical integrals (will be rather slow) | |
596 | // - this version assumes that func is defined using a | |
597 | // TFormula. Generalize to the case of a C++ function | |
598 | ||
599 | if (normPar<0) Printf("AliBWTools::GetMoment: Warning: If normalization is required, the error may bot be correct"); | |
600 | if (!strcmp(func->GetExpFormula(),"")) Printf("AliBWTools::GetMoment: Warning: Empty formula in the base function"); | |
4c0d7bc7 | 601 | Int_t npar = func->GetNpar(); |
602 | ||
3472a634 | 603 | // Definition changes according to the value of normPar |
604 | TF1 * f = normPar < 0 ? | |
605 | new TF1(name, var) : // not normalized | |
606 | new TF1(name, var+Form("/[%d]",normPar)); // normalized with par normPar | |
607 | ||
608 | // integr is used to normalize if no parameter is provided | |
609 | Double_t integr = normPar < 0 ? func->Integral(min,max) : 1; | |
4c0d7bc7 | 610 | |
611 | // The parameter of the function used to compute the mean should be | |
612 | // the same as the parent function: fixed if needed and they should | |
613 | // also have the same errors. | |
614 | ||
615 | // cout << "npar :" << npar << endl; | |
616 | ||
617 | for(Int_t ipar = 0; ipar < npar; ipar++){ | |
618 | Double_t parmin, parmax; | |
619 | Double_t value = func->GetParameter(ipar); | |
3472a634 | 620 | f->SetParameter(ipar,value); |
4c0d7bc7 | 621 | func->GetParLimits(ipar, parmin, parmax); |
622 | if ( parmin == parmax ) { | |
3beb22a2 | 623 | // if ( parmin || (parmin == 1 && !value) ) { // not sure we I check parmin == 1 here. |
624 | if ( parmin || (TMath::Abs(parmin-1)<0.000001 && !value) ) { // not sure we I check parmin == 1 here. Changed like this because of coding conventions. Does it still work? FIXME | |
4c0d7bc7 | 625 | f->FixParameter(ipar,func->GetParameter(ipar)); |
626 | // cout << "Fixing " << ipar << "("<<value<<","<<parmin<<","<<parmax<<")"<<endl; | |
627 | } | |
628 | else { | |
629 | f->SetParError (ipar,func->GetParError(ipar) ); | |
630 | // cout << "Setting Err" << ipar << "("<<func->GetParError(ipar)<<")"<<endl; | |
631 | } | |
632 | } | |
633 | else { | |
634 | f->SetParError (ipar,func->GetParError(ipar) ); | |
635 | // cout << "Setting Err" << ipar << "("<<func->GetParError(ipar)<<")"<<endl; | |
636 | } | |
637 | } | |
3472a634 | 638 | |
639 | // f->Print(); | |
640 | // cout << "----" << endl; | |
641 | // func->Print(); | |
642 | ||
643 | mean = normPar < 0 ? f->Integral (min,max)/integr : f->Integral (min,max); | |
644 | error = normPar < 0 ? f->IntegralError(min,max)/integr : f->IntegralError(min,max); | |
4c0d7bc7 | 645 | // cout << "Mean " << mean <<"+-"<< error<< endl; |
646 | // cout << "Integral Error " << error << endl; | |
647 | ||
648 | } | |
649 | ||
4c0d7bc7 | 650 | |
651 | ||
652 | Bool_t AliBWTools::Fit (TH1 * h1, TF1* func, Float_t min, Float_t max) { | |
653 | ||
654 | // Fits h1 with func, preforms several checks on the quality of the | |
655 | // fit and tries to improve it. If the fit is not good enough, it | |
656 | // returs false. | |
657 | ||
658 | Double_t amin; Double_t edm; Double_t errdef; Int_t nvpar; Int_t nparx; | |
659 | TVirtualFitter *fitter; | |
660 | cout << "--- Fitting : " << h1->GetName() << " ["<< h1->GetTitle() <<"] ---"<< endl; | |
661 | ||
662 | h1-> Fit(func,"IME0","",min,max); | |
663 | Int_t fitResult = h1-> Fit(func,"IME0","",min,max); | |
664 | // h1-> Fit(func,"0","",min,max); | |
665 | // Int_t fitResult = h1-> Fit(func,"0IE","",min,max); | |
666 | ||
667 | ||
668 | // From TH1: | |
669 | // The fitStatus is 0 if the fit is OK (i.e no error occurred). The | |
670 | // value of the fit status code is negative in case of an error not | |
671 | // connected with the minimization procedure, for example when a wrong | |
672 | // function is used. Otherwise the return value is the one returned | |
673 | // from the minimization procedure. When TMinuit (default case) or | |
674 | // Minuit2 are used as minimizer the status returned is : fitStatus = | |
675 | // migradResult + 10*minosResult + 100*hesseResult + | |
676 | // 1000*improveResult. TMinuit will return 0 (for migrad, minos, | |
677 | // hesse or improve) in case of success and 4 in case of error (see | |
678 | // the documentation of TMinuit::mnexcm). So for example, for an error | |
679 | // only in Minos but not in Migrad a fitStatus of 40 will be returned. | |
680 | // Minuit2 will return also 0 in case of success and different values | |
681 | // in migrad minos or hesse depending on the error. See in this case | |
682 | // the documentation of Minuit2Minimizer::Minimize for the | |
683 | // migradResult, Minuit2Minimizer::GetMinosError for the minosResult | |
684 | // and Minuit2Minimizer::Hesse for the hesseResult. If other | |
685 | // minimizers are used see their specific documentation for the status | |
686 | // code returned. For example in the case of Fumili, for the status | |
687 | // returned see TFumili::Minimize. | |
688 | ||
689 | ||
690 | if( gMinuit->fLimset ) { | |
691 | Printf("ERROR: AliBWTools: Parameters at limits"); | |
692 | return kFALSE; | |
693 | } | |
694 | ||
695 | ||
696 | /// | |
697 | fitter = TVirtualFitter::GetFitter(); | |
698 | Int_t fitStat = fitter->GetStats(amin, edm, errdef, nvpar, nparx); | |
699 | ||
700 | if( ( (fitStat < 3 && gMinuit->fCstatu != "UNCHANGED ")|| (edm > 1e6) || (fitResult !=0 && fitResult < 4000) ) && | |
701 | TString(gMinuit->fCstatu) != "SUCCESSFUL" && | |
702 | TString(gMinuit->fCstatu) != "CONVERGED " ) { | |
703 | if(fitStat < 3 && gMinuit->fCstatu != "UNCHANGED ") { | |
704 | Printf("WARNING: AliBWTools: Cannot properly compute errors"); | |
705 | if (fitStat == 0) Printf(" not calculated at all"); | |
706 | if (fitStat == 1) Printf(" approximation only, not accurate"); | |
707 | if (fitStat == 2) Printf(" full matrix, but forced positive-definite"); | |
708 | } | |
709 | ||
710 | if (edm > 1e6) { | |
711 | ||
712 | Printf("WARNING: AliBWTools: Huge EDM (%f): Fit probably failed!", edm); | |
713 | } | |
714 | if (fitResult != 0) { | |
715 | Printf("WARNING: AliBWTools: Fit Result (%d)", fitResult); | |
716 | } | |
717 | ||
718 | Printf ("AliBWTools: Trying Again with Strategy = 2"); | |
719 | ||
720 | gMinuit->Command("SET STRATEGY 2"); // more effort | |
721 | fitResult = 0; | |
722 | fitResult = h1-> Fit(func,"0","",min,max); | |
723 | fitResult = h1-> Fit(func,"IME0","",min,max); | |
724 | fitResult = h1-> Fit(func,"IME0","",min,max); | |
725 | ||
726 | fitter = TVirtualFitter::GetFitter(); | |
727 | ||
728 | fitStat = fitter->GetStats(amin, edm, errdef, nvpar, nparx); | |
729 | ||
730 | if(fitStat < 3 && gMinuit->fCstatu != "UNCHANGED ") { | |
731 | Printf("ERROR: AliBWTools: Cannot properly compute errors"); | |
732 | if (fitStat == 0) Printf(" not calculated at all"); | |
733 | if (fitStat == 1) Printf(" approximation only, not accurate"); | |
734 | if (fitStat == 2) Printf(" full matrix, but forced positive-definite"); | |
735 | cout << "[" <<gMinuit->fCstatu<<"]" << endl; | |
736 | return kFALSE; | |
737 | } | |
738 | ||
739 | if (edm > 1e6) { | |
740 | Printf("ERROR: AliBWTools: Huge EDM (%f): Fit probably failed!", edm); | |
741 | ||
742 | return kFALSE; | |
743 | } | |
744 | if (fitResult != 0) { | |
745 | Printf("ERROR: AliBWTools: Fit Result (%d)", fitResult); | |
746 | ||
747 | return kFALSE; | |
748 | } | |
749 | ||
750 | gMinuit->Command("SET STRATEGY 1"); // back to normal value | |
751 | ||
752 | } | |
753 | ||
754 | cout << "---- FIT OK ----" << endl; | |
755 | ||
756 | return kTRUE; | |
757 | ||
758 | } | |
759 | ||
3beb22a2 | 760 | Int_t AliBWTools::GetLowestNotEmptyBin(const TH1*h) { |
4c0d7bc7 | 761 | |
762 | // Return the index of the lowest non empty bin in the histo h | |
763 | ||
764 | Int_t nbin = h->GetNbinsX(); | |
765 | for(Int_t ibin = 1; ibin <= nbin; ibin++){ | |
766 | if(h->GetBinContent(ibin)>0) return ibin; | |
767 | } | |
768 | ||
769 | return -1; | |
770 | ||
771 | } | |
772 | ||
3beb22a2 | 773 | Int_t AliBWTools::GetHighestNotEmptyBin(const TH1*h) { |
4c0d7bc7 | 774 | |
775 | // Return the index of the highest non empty bin in the histo h | |
776 | ||
777 | Int_t nbin = h->GetNbinsX(); | |
778 | for(Int_t ibin = nbin; ibin > 0; ibin--){ | |
779 | if(h->GetBinContent(ibin)>0) return ibin; | |
780 | } | |
781 | ||
782 | return -1; | |
783 | ||
784 | } | |
785 | ||
3beb22a2 | 786 | void AliBWTools::GetResiduals(const TGraphErrors * gdata, const TF1 * func, TH1F ** hres, TGraphErrors ** gres) { |
4c0d7bc7 | 787 | |
788 | // Returns a graph of residuals vs point and the res/err distribution | |
789 | ||
790 | Int_t npoint = gdata->GetN(); | |
791 | ||
792 | (*gres) =new TGraphErrors(); | |
793 | (*hres) = new TH1F(TString("hres_")+gdata->GetName()+"-"+func->GetName(), | |
794 | TString("hres_")+gdata->GetName()+"-"+func->GetName(), | |
795 | 20,-5,5); | |
796 | ||
797 | ||
798 | for(Int_t ipoint = 0; ipoint < npoint; ipoint++){ | |
799 | Float_t x = gdata->GetX()[ipoint]; | |
800 | Float_t res = (gdata->GetY()[ipoint] - func->Eval(x))/func->Eval(x); | |
801 | Float_t err = gdata->GetEY()[ipoint]/func->Eval(x); | |
802 | (*hres)->Fill(res/err); | |
803 | (*gres)->SetPoint(ipoint, x, res/err); | |
804 | // (*gres)->SetPointError(ipoint, 0, err); | |
805 | ||
806 | } | |
807 | ||
808 | (*gres)->SetMarkerStyle(gdata->GetMarkerStyle()); | |
809 | (*gres)->SetMarkerColor(gdata->GetMarkerColor()); | |
810 | (*gres)->SetLineColor (gdata->GetLineColor()); | |
811 | (*gres)->GetHistogram()->GetYaxis()->SetTitle("(data-function)/function"); | |
812 | (*hres)->SetMarkerStyle(gdata->GetMarkerStyle()); | |
813 | (*hres)->SetMarkerColor(gdata->GetMarkerColor()); | |
814 | (*hres)->SetLineColor (gdata->GetLineColor()); | |
815 | ||
816 | ||
817 | ||
818 | } | |
819 | ||
3beb22a2 | 820 | void AliBWTools::GetResiduals(const TH1F* hdata, const TF1 * func, TH1F ** hres, TH1F ** hresVsBin) { |
4c0d7bc7 | 821 | |
822 | // Returns an histo of residuals bin by bin and the res/err distribution | |
823 | ||
824 | if (!func) { | |
825 | Printf("AliBWTools::GetResiduals: No function provided"); | |
826 | return; | |
827 | } | |
828 | if (!hdata) { | |
829 | Printf("AliBWTools::GetResiduals: No data provided"); | |
830 | return; | |
831 | } | |
832 | ||
833 | (*hresVsBin) = (TH1F*) hdata->Clone(TString("hres_")+hdata->GetName()); | |
834 | (*hresVsBin)->Reset(); | |
835 | (*hres) = new TH1F(TString("hres_")+hdata->GetName()+"-"+func->GetName(), | |
836 | TString("hres_")+hdata->GetName()+"-"+func->GetName(), | |
837 | 20,-5,5); | |
838 | ||
839 | Int_t nbin = hdata->GetNbinsX(); | |
840 | for(Int_t ibin = 1; ibin <= nbin; ibin++){ | |
3beb22a2 | 841 | if(!hdata->GetBinContent(ibin)) continue; |
4c0d7bc7 | 842 | Float_t res = (hdata->GetBinContent(ibin) - func->Eval(hdata->GetBinCenter(ibin)) ) / |
843 | func->Eval(hdata->GetBinCenter(ibin)); | |
844 | Float_t err = hdata->GetBinError (ibin) / func->Eval(hdata->GetBinCenter(ibin)); | |
845 | (*hresVsBin)->SetBinContent(ibin,res); | |
846 | (*hresVsBin)->SetBinError (ibin,err); | |
847 | (*hres)->Fill(res/err); | |
848 | ||
849 | } | |
850 | ||
851 | (*hresVsBin)->SetMarkerStyle(hdata->GetMarkerStyle()); | |
852 | (*hresVsBin)->SetMarkerColor(hdata->GetMarkerColor()); | |
853 | (*hresVsBin)->SetLineColor (hdata->GetLineColor() ); | |
854 | (*hresVsBin)->GetYaxis()->SetTitle("(data-function)/function"); | |
855 | (*hres)->SetMarkerStyle(hdata->GetMarkerStyle()); | |
856 | (*hres)->SetMarkerColor(hdata->GetMarkerColor()); | |
857 | (*hres)->SetLineColor (hdata->GetLineColor() ); | |
858 | ||
859 | } | |
860 | ||
3beb22a2 | 861 | void AliBWTools::GetYield(TH1* h, TF1 * f, Double_t &yield, Double_t &yieldError, Float_t min, Float_t max, |
4c0d7bc7 | 862 | Double_t *partialYields, Double_t *partialYieldsErrors){ |
863 | ||
864 | // Returns the yield extracted from the data in the histo where | |
865 | // there are points and from the fit to extrapolate, in the given | |
866 | // range. | |
867 | ||
868 | // Partial yields are also returned if the corresponding pointers are non null | |
869 | ||
870 | Int_t bin1 = h->FindBin(min); | |
871 | Int_t bin2 = h->FindBin(max); | |
872 | Float_t bin1Edge = GetLowestNotEmptyBinEdge (h); | |
873 | Float_t bin2Edge = GetHighestNotEmptyBinEdge(h); | |
874 | ||
875 | Double_t integralFromHistoError ; | |
57b5143c | 876 | Double_t integralFromHisto = DoIntegral(h,bin1,bin2,-1,-1,-1,-1,integralFromHistoError,"width",1); |
4c0d7bc7 | 877 | |
878 | Double_t integralBelow = min < bin1Edge ? f->Integral(min,bin1Edge) : 0; | |
879 | Double_t integralBelowError = min < bin1Edge ? f->IntegralError(min,bin1Edge) : 0; | |
880 | Double_t integralAbove = max > bin2Edge ? f->Integral(bin2Edge,max) : 0; | |
169d33f6 | 881 | Double_t integralAboveError = max > bin2Edge ? f->IntegralError(bin2Edge,max) : 0; |
4c0d7bc7 | 882 | |
3beb22a2 | 883 | // cout << "GetYield INFO" << endl; |
884 | // cout << " " << bin1Edge << " " << bin2Edge << endl; | |
885 | // cout << " " << integralFromHisto << " " << integralBelow << " " << integralAbove << endl; | |
886 | // cout << " " << integralFromHistoError << " " << integralBelowError << " " << integralAboveError << endl; | |
4c0d7bc7 | 887 | |
888 | if(partialYields) { | |
889 | partialYields[0] = integralFromHisto; | |
890 | partialYields[1] = integralBelow; | |
891 | partialYields[2] = integralAbove; | |
892 | } | |
893 | if(partialYieldsErrors) { | |
894 | partialYieldsErrors[0] = integralFromHistoError; | |
895 | partialYieldsErrors[1] = integralBelowError; | |
896 | partialYieldsErrors[2] = integralAboveError; | |
897 | } | |
898 | yield = integralFromHisto+integralBelow+integralAbove; | |
899 | yieldError = TMath::Sqrt(integralFromHistoError*integralFromHistoError+ | |
900 | integralBelowError*integralBelowError+ | |
901 | integralAboveError*integralAboveError); | |
902 | ||
903 | } | |
904 | ||
3beb22a2 | 905 | TGraphErrors * AliBWTools::DivideGraphByFunc(const TGraphErrors * g, const TF1 * f, Bool_t invert){ |
4c0d7bc7 | 906 | |
907 | // Divides g/f. If invert == true => f/g | |
908 | ||
909 | TGraphErrors * gRatio = new TGraphErrors(); | |
910 | Int_t npoint = g->GetN(); | |
911 | for(Int_t ipoint = 0; ipoint < npoint; ipoint++){ | |
912 | Double_t x = g->GetX()[ipoint]; | |
913 | Double_t ratio = invert ? f->Eval(x)/g->GetY()[ipoint] :g->GetY()[ipoint]/f->Eval(x); | |
914 | gRatio->SetPoint (ipoint, x, ratio); | |
915 | gRatio->SetPointError(ipoint, 0, g->GetEY()[ipoint]/f->Eval(x)); | |
916 | // cout << x << " " << g->GetY()[ipoint] << " " << f->Eval(x) << endl; | |
917 | ||
918 | } | |
919 | gRatio->SetMarkerStyle(20); | |
920 | //gRatio->Print(); | |
921 | return gRatio; | |
922 | ||
923 | } | |
924 | ||
3beb22a2 | 925 | TGraphErrors * AliBWTools::DivideGraphByHisto(const TGraphErrors * g, TH1 * h, Bool_t invert){ |
4c0d7bc7 | 926 | |
927 | // Divides g/h. If invert == true => h/g | |
928 | ||
169d33f6 | 929 | Bool_t skipError = kFALSE; |
930 | if(!strcmp(g->ClassName(),"TGraph")) skipError = kTRUE; | |
931 | if(!strcmp(g->ClassName(),"TGraphAsymmErrors")) skipError = kTRUE; | |
932 | if(skipError) { | |
933 | Printf("WARNING: Skipping graph errors"); | |
934 | } | |
4c0d7bc7 | 935 | TGraphErrors * gRatio = new TGraphErrors(); |
936 | Int_t npoint = g->GetN(); | |
937 | for(Int_t ipoint = 0; ipoint < npoint; ipoint++){ | |
938 | Double_t xj = g->GetX()[ipoint]; | |
939 | Double_t yj = g->GetY()[ipoint]; | |
169d33f6 | 940 | Double_t yje = skipError ? 0 : g->GetEY()[ipoint]; |
4c0d7bc7 | 941 | |
942 | Int_t binData = h->FindBin(xj); | |
943 | Double_t yd = h->GetBinContent(binData); | |
944 | Double_t yde = h->GetBinError(binData); | |
945 | Double_t xd = h->GetBinCenter(binData); | |
946 | ||
4c0d7bc7 | 947 | |
948 | ||
3beb22a2 | 949 | if (!yd) continue; |
4c0d7bc7 | 950 | |
951 | if (TMath::Abs(xd-xj)/TMath::Abs(xd) > 0.01){ | |
952 | Printf( "WARNING: bin center (%f) and x graph (%f) are more than 1 %% away, skipping",xd,xj ); | |
953 | continue; | |
954 | ||
955 | } | |
956 | ||
957 | Double_t ratio = invert ? yd/yj : yj/yd; | |
958 | ||
959 | gRatio->SetPoint(ipoint, xj, ratio); | |
960 | gRatio->SetPointError(ipoint, 0, TMath::Sqrt(yde*yde/yd/yd + yje*yje/yj/yj)*ratio); | |
961 | // gRatio->SetPointError(ipoint, 0, yje/yj * ratio); | |
962 | } | |
963 | ||
964 | return gRatio; | |
965 | ||
966 | ||
967 | } | |
968 | ||
969 | TH1F * AliBWTools::DivideHistoByFunc(TH1F * h, TF1 * f, Bool_t invert){ | |
970 | ||
971 | // Divides h/f. If invert == true => f/g | |
972 | // Performs the integral of f on the bin range to perform the ratio | |
973 | // Returns a histo with the same binnig as h | |
974 | ||
975 | // Prepare histo for ratio | |
976 | TH1F * hRatio = (TH1F*) h->Clone(TString("hRatio_")+h->GetName()+"_"+f->GetName()); | |
977 | hRatio->Reset(); | |
978 | // Set y title | |
979 | if(!invert) hRatio->SetYTitle(TString(h->GetName())+"/"+f->GetName()); | |
980 | else hRatio->SetYTitle(TString(f->GetName())+"/"+h->GetName()); | |
981 | ||
982 | // Loop over all bins | |
983 | Int_t nbin = hRatio->GetNbinsX(); | |
984 | ||
985 | for(Int_t ibin = 1; ibin <= nbin; ibin++){ | |
986 | Double_t yhisto = h->GetBinContent(ibin); | |
987 | Double_t yerror = h->GetBinError(ibin); | |
988 | Double_t xmin = h->GetBinLowEdge(ibin); | |
989 | Double_t xmax = h->GetBinLowEdge(ibin+1); | |
990 | Double_t yfunc = f->Integral(xmin,xmax)/(xmax-xmin); | |
991 | Double_t ratio = invert ? yfunc/yhisto : yhisto/yfunc ; | |
992 | Double_t error = yerror/yfunc ; | |
993 | hRatio->SetBinContent(ibin,ratio); | |
994 | hRatio->SetBinError (ibin,error); | |
995 | } | |
996 | ||
997 | return hRatio; | |
998 | ||
999 | } | |
ee08b77d | 1000 | |
1001 | void AliBWTools::WeightedMean(Int_t npoints, const Double_t *x, const Double_t *xerr, Double_t &mean, Double_t &meanerr){ | |
1002 | ||
3beb22a2 | 1003 | // Performs the weighted mean of npoints numbers in x with errors in xerr |
ee08b77d | 1004 | |
1005 | mean = 0; | |
1006 | meanerr = 0; | |
1007 | ||
1008 | Double_t sumweight = 0; | |
1009 | ||
1010 | for (Int_t ipoint = 0; ipoint < npoints; ipoint++){ | |
1011 | ||
1012 | Double_t xerr2 = xerr[ipoint]*xerr[ipoint]; | |
1013 | if(xerr2>0){ | |
1014 | // cout << "xe2 " << xerr2 << endl; | |
1015 | Double_t weight = 1. / xerr2; | |
1016 | sumweight += weight; | |
1017 | mean += weight * x[ipoint]; | |
3725266c | 1018 | }// else cout << " Skipping " << ipoint << endl; |
ee08b77d | 1019 | |
1020 | } | |
1021 | ||
1022 | ||
1023 | if(sumweight){ | |
1024 | mean /= sumweight; | |
1025 | meanerr = TMath::Sqrt(1./ sumweight); | |
1026 | } | |
1027 | else { | |
3725266c | 1028 | // cout << " No sumweight" << endl; |
ee08b77d | 1029 | mean = 0; |
1030 | meanerr = 0; | |
1031 | } | |
1032 | ||
1033 | ||
1034 | } | |
3725266c | 1035 | |
57b5143c | 1036 | TH1 * AliBWTools::GetRelativeError(TH1 * h){ |
1037 | // Returns an histogram with the same binning as h, filled with the relative error bin by bin | |
1038 | TH1 * hrel = (TH1*) h->Clone(TString(h->GetName())+"_rel"); | |
1039 | hrel->Reset(); | |
1040 | Int_t nbin = hrel->GetNbinsX(); | |
1041 | for(Int_t ibin = 1; ibin <= nbin; ibin++){ | |
1042 | hrel->SetBinContent(ibin,h->GetBinError(ibin)/h->GetBinContent(ibin)); | |
1043 | hrel->SetBinError(ibin,0); | |
1044 | } | |
1045 | ||
1046 | return hrel; | |
1047 | } | |
1048 | ||
1049 | ||
058137a6 | 1050 | void AliBWTools::GetValueAndError(TH1 * hdest, const TH1 * hvalue, const TH1 * herror, Bool_t isPercentError) { |
3725266c | 1051 | |
1052 | // Put into source, bin-by-bin, the values from hvalue and the | |
1053 | // errors from content from herror. | |
1054 | // Used mainly to combine histos of systemati errors with their spectra | |
1055 | // Set isPercentError to kTRUE if the error is given in % | |
1056 | ||
1057 | if(hdest == NULL){ | |
1058 | Printf("AliBWTools::GetValueAndError Errror: hdest is null"); | |
1059 | return; | |
1060 | } | |
1061 | ||
1062 | ||
1063 | Int_t nbin = hdest->GetNbinsX(); | |
1064 | Int_t nBinSourceVal = hvalue->GetNbinsX(); | |
1065 | Int_t nBinSourceErr = herror->GetNbinsX(); | |
1066 | ||
1067 | for(Int_t iBinDest = 1; iBinDest <= nbin; iBinDest++){ | |
1068 | Float_t lowPtDest=hdest->GetBinLowEdge(iBinDest); | |
1069 | Float_t binWidDest=hdest->GetBinWidth(iBinDest); | |
1070 | // Loop over Source bins and find overlapping bins to Dest | |
1071 | // First value then error | |
1072 | // Value | |
1073 | Bool_t foundValue = kFALSE; | |
1074 | for(Int_t iBinSourceVal=1; iBinSourceVal<=nBinSourceVal; iBinSourceVal++){ | |
1075 | Float_t lowPtSource= hvalue->GetBinLowEdge(iBinSourceVal) ; | |
1076 | Float_t binWidSource= hvalue->GetBinWidth(iBinSourceVal); | |
1077 | if(TMath::Abs(lowPtDest-lowPtSource)<0.001 && TMath::Abs(binWidSource-binWidDest)<0.001){ | |
1078 | Double_t content = hvalue->GetBinContent(iBinSourceVal); | |
1079 | hdest->SetBinContent(iBinDest, content); | |
1080 | foundValue = kTRUE; | |
1081 | break; | |
1082 | } | |
1083 | } | |
1084 | // if (!foundValue){ | |
1085 | // Printf("AliBWTools::GetValueAndError: Error: cannot find matching value source bin for destination %d",iBinDest); | |
1086 | // } | |
1087 | ||
1088 | // Error | |
1089 | Bool_t foundError = kFALSE; | |
1090 | for(Int_t iBinSourceErr=1; iBinSourceErr<=nBinSourceErr; iBinSourceErr++){ | |
1091 | Float_t lowPtSource= herror->GetBinLowEdge(iBinSourceErr) ; | |
1092 | Float_t binWidSource= herror->GetBinWidth(iBinSourceErr); | |
1093 | if(TMath::Abs(lowPtDest-lowPtSource)<0.001 && TMath::Abs(binWidSource-binWidDest)<0.001){ | |
1094 | Double_t error = herror->GetBinContent(iBinSourceErr); | |
1095 | // cout << "-> " << iBinDest << " " << error << " " << hdest->GetBinContent(iBinDest) << endl; | |
1096 | ||
1097 | hdest->SetBinError(iBinDest, isPercentError ? error * hdest->GetBinContent(iBinDest) : error); | |
1098 | foundError=kTRUE; | |
1099 | break; | |
1100 | } | |
1101 | } | |
1102 | // if (!foundError ){ | |
1103 | // Printf("AliBWTools::GetValueAndError: Error: cannot find matching error source bin for destination %d",iBinDest); | |
1104 | // } | |
1105 | } | |
1106 | ||
1107 | ||
1108 | } | |
1109 | ||
058137a6 | 1110 | void AliBWTools::AddHisto(TH1 * hdest, const TH1* hsource, Bool_t getMirrorBins) { |
3725266c | 1111 | |
1112 | // Adds hsource to hdest bin by bin, even if they have a different | |
1113 | // binning If getMirrorBins is true, it takes the negative bins | |
1114 | // (Needed because sometimes the TPC uses the positive axis for | |
1115 | // negative particles and the possitive axis for positive | |
1116 | // particles). | |
1117 | ||
1118 | ||
1119 | if (hdest == NULL) { | |
1120 | Printf("Error: hdest is NULL\n"); | |
37e1d042 | 1121 | return; |
3725266c | 1122 | } |
1123 | if (hsource == NULL) { | |
1124 | Printf("Error: hsource is NULL\n"); | |
37e1d042 | 1125 | return; |
3725266c | 1126 | } |
1127 | ||
1128 | Int_t nBinSource = hsource->GetNbinsX(); | |
1129 | Int_t nBinDest = hdest->GetNbinsX(); | |
1130 | ||
1131 | // Loop over destination bins, | |
1132 | for(Int_t iBinDest=1; iBinDest<=nBinDest; iBinDest++){ | |
1133 | Float_t lowPtDest=hdest->GetBinLowEdge(iBinDest); | |
1134 | Float_t binWidDest=hdest->GetBinWidth(iBinDest); | |
1135 | // Loop over Source bins and find overlapping bins to Dest | |
1136 | Bool_t found = kFALSE; | |
1137 | for(Int_t iBinSource=1; iBinSource<=nBinSource; iBinSource++){ | |
1138 | Float_t lowPtSource= getMirrorBins ? -hsource->GetBinLowEdge(iBinSource)+hsource->GetBinWidth(iBinSource) : hsource->GetBinLowEdge(iBinSource) ; | |
1139 | Float_t binWidSource= hsource->GetBinWidth(iBinSource) ; | |
1140 | if(TMath::Abs(lowPtDest-lowPtSource)<0.001 && TMath::Abs(binWidSource-binWidDest)<0.001){ | |
1141 | Float_t dest=hdest->GetBinContent(iBinDest); | |
1142 | Float_t source=hsource->GetBinContent(iBinSource); | |
1143 | Float_t edest=hdest->GetBinError(iBinDest); | |
1144 | Float_t esource=hsource->GetBinError(iBinSource); | |
1145 | Double_t cont=dest+source; | |
1146 | Double_t econt=TMath::Sqrt(edest*edest+esource*esource); | |
1147 | hdest->SetBinContent(iBinDest,cont); | |
1148 | hdest->SetBinError (iBinDest,econt); | |
1149 | found = kTRUE; | |
1150 | ||
1151 | break; | |
1152 | } | |
1153 | } | |
77233d47 | 1154 | // if (!found){ |
1155 | // Printf("Error: cannot find matching source bin for destination %d",iBinDest); | |
1156 | // } | |
3725266c | 1157 | } |
1158 | ||
1159 | ||
1160 | } | |
1161 | ||
058137a6 | 1162 | void AliBWTools::GetHistoCombinedErrors(TH1 * hdest, const TH1 * h1) { |
1163 | ||
1164 | // Combine the errors of hdest with the errors of h1, summing in | |
1165 | // quadrature. Results are put in hdest. Histograms are assumed to | |
1166 | // have the same binning | |
1167 | ||
1168 | Int_t nbin = hdest->GetNbinsX(); | |
1169 | for(Int_t ibin = 0; ibin < nbin; ibin++){ | |
1170 | Double_t e1 = hdest->GetBinError(ibin); | |
1171 | Double_t e2 = h1->GetBinError(ibin); | |
1172 | hdest->SetBinError(ibin, TMath::Sqrt(e1*e1+e2*e2)); | |
1173 | } | |
1174 | ||
1175 | ||
1176 | } | |
d8e6677d | 1177 | |
1178 | TH1F * AliBWTools::DivideHistosDifferentBins(const TH1F* h1, const TH1F* h2) { | |
1179 | // Divides 2 histos even if they have a different binning. Finds | |
1180 | // overlapping bins and divides them | |
1181 | ||
1182 | // 1. clone histo | |
1183 | TH1F * hRatio = new TH1F(*h1); | |
1184 | Int_t nBinsH1=h1->GetNbinsX(); | |
1185 | Int_t nBinsH2=h2->GetNbinsX(); | |
1186 | // Loop over H1 bins, | |
1187 | for(Int_t iBin=1; iBin<=nBinsH1; iBin++){ | |
1188 | hRatio->SetBinContent(iBin,0.); | |
1189 | hRatio->SetBinContent(iBin,0.); | |
1190 | Float_t lowPtH1=h1->GetBinLowEdge(iBin); | |
1191 | Float_t binWidH1=h1->GetBinWidth(iBin); | |
1192 | // Loop over H2 bins and find overlapping bins to H1 | |
1193 | for(Int_t jBin=1; jBin<=nBinsH2; jBin++){ | |
1194 | Float_t lowPtH2=h2->GetBinLowEdge(jBin); | |
1195 | Float_t binWidH2=h2->GetBinWidth(jBin); | |
1196 | if(TMath::Abs(lowPtH1-lowPtH2)<0.001 && TMath::Abs(binWidH2-binWidH1)<0.001){ | |
1197 | Float_t numer=h1->GetBinContent(iBin); | |
1198 | Float_t denom=h2->GetBinContent(jBin); | |
1199 | Float_t enumer=h1->GetBinError(iBin); | |
1200 | Float_t edenom=h2->GetBinError(jBin); | |
1201 | Double_t ratio=0.; | |
1202 | Double_t eratio=0.; | |
1203 | if(numer>0. && denom>0.){ | |
1204 | ratio=numer/denom; | |
1205 | eratio=ratio*TMath::Sqrt((enumer/numer)*(enumer/numer)+(edenom/denom)*(edenom/denom)); | |
1206 | } | |
1207 | hRatio->SetBinContent(iBin,ratio); | |
1208 | hRatio->SetBinError(iBin,eratio); | |
1209 | break; | |
1210 | } | |
1211 | } | |
1212 | } | |
1213 | return hRatio; | |
1214 | } | |
57b5143c | 1215 | |
1216 | Double_t AliBWTools::DoIntegral(TH1* h, Int_t binx1, Int_t binx2, Int_t biny1, Int_t biny2, Int_t binz1, Int_t binz2, Double_t & error , | |
1217 | Option_t *option, Bool_t doError) | |
1218 | { | |
1219 | // function to compute integral and optionally the error between the limits | |
1220 | // specified by the bin number values working for all histograms (1D, 2D and 3D) | |
1221 | // copied from TH! to fix a bug still present in 5-27-06b | |
1222 | Int_t nbinsx = h->GetNbinsX(); | |
1223 | if (binx1 < 0) binx1 = 0; | |
1224 | if (binx2 > nbinsx+1 || binx2 < binx1) binx2 = nbinsx+1; | |
1225 | if (h->GetDimension() > 1) { | |
1226 | Int_t nbinsy = h->GetNbinsY(); | |
1227 | if (biny1 < 0) biny1 = 0; | |
1228 | if (biny2 > nbinsy+1 || biny2 < biny1) biny2 = nbinsy+1; | |
1229 | } else { | |
1230 | biny1 = 0; biny2 = 0; | |
1231 | } | |
1232 | if (h->GetDimension() > 2) { | |
1233 | Int_t nbinsz = h->GetNbinsZ(); | |
1234 | if (binz1 < 0) binz1 = 0; | |
1235 | if (binz2 > nbinsz+1 || binz2 < binz1) binz2 = nbinsz+1; | |
1236 | } else { | |
1237 | binz1 = 0; binz2 = 0; | |
1238 | } | |
1239 | ||
1240 | // - Loop on bins in specified range | |
1241 | TString opt = option; | |
1242 | opt.ToLower(); | |
1243 | Bool_t width = kFALSE; | |
1244 | if (opt.Contains("width")) width = kTRUE; | |
1245 | ||
1246 | ||
1247 | Double_t dx = 1.; | |
1248 | Double_t dy = 1.; | |
1249 | Double_t dz = 1.; | |
1250 | Double_t integral = 0; | |
1251 | Double_t igerr2 = 0; | |
1252 | for (Int_t binx = binx1; binx <= binx2; ++binx) { | |
1253 | if (width) dx = h->GetXaxis()->GetBinWidth(binx); | |
1254 | for (Int_t biny = biny1; biny <= biny2; ++biny) { | |
1255 | if (width) dy = h->GetYaxis()->GetBinWidth(biny); | |
1256 | for (Int_t binz = binz1; binz <= binz2; ++binz) { | |
1257 | if (width) dz = h->GetZaxis()->GetBinWidth(binz); | |
1258 | Int_t bin = h->GetBin(binx, biny, binz); | |
1259 | if (width) integral += h->GetBinContent(bin)*dx*dy*dz; | |
1260 | else integral += h->GetBinContent(bin); | |
1261 | if (doError) { | |
1262 | if (width) igerr2 += h->GetBinError(bin)*h->GetBinError(bin)*dx*dy*dz*dx*dy*dz; | |
1263 | else igerr2 += h->GetBinError(bin)*h->GetBinError(bin); | |
1264 | } | |
1265 | // cout << h->GetBinContent(bin) << " " << h->GetBinError(bin) << " " << dx*dy*dz << " " << integral << " +- " << igerr2 << endl; | |
1266 | ||
1267 | } | |
1268 | } | |
1269 | } | |
1270 | ||
1271 | if (doError) error = TMath::Sqrt(igerr2); | |
1272 | return integral; | |
1273 | } |