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