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