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