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