]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALHistoUtilities.cxx
added pi0 calibration, linearity, shower profile
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALHistoUtilities.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2002, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //_________________________________________________________________________
19 // This is a set of histogram
20 // utilities for the EMCAL
21 // to make some common
22 // functions easier
23 //
24 //*-- Authors: J.L. Klay (LLNL) & Aleksei Pavlinov (WSU) 
25
26 #include "AliEMCALHistoUtilities.h"
27
28 #include <iostream>
29 #include <iomanip>
30 #include <fstream>
31
32 #include <TROOT.h>
33 #include <TPad.h>
34 #include <TFile.h>
35 #include <TH1.h>
36 #include <TH2.h>
37 #include <TF1.h>
38 #include <TGraph.h>
39 #include <TGraphErrors.h>
40 #include <TLatex.h>
41 #include <TChain.h>
42 #include <TList.h>
43 #include <TObjArray.h>
44 #include <TObjString.h>
45 #include <TRegexp.h>
46 #include <TString.h>
47 #include <TLorentzVector.h>
48 #include <Gtypes.h> // color, line style and so on
49
50 #include "AliESDCaloCluster.h"
51 #include "AliEMCALRecPoint.h"
52 #include "AliRunLoader.h"
53
54 using namespace std;
55
56 ClassImp(AliEMCALHistoUtilities)
57
58 AliEMCALHistoUtilities::AliEMCALHistoUtilities(const char *name, const char *tit) : TNamed(name,tit)
59 {
60   // constructor
61 }
62
63 AliEMCALHistoUtilities::~AliEMCALHistoUtilities()
64 {
65         //destructor
66 }  
67
68 TList* AliEMCALHistoUtilities::MoveHistsToList(const char* name, Bool_t putToBrowser)
69 {
70   // Move HIST to list
71   gROOT->cd();
72   TIter nextHist(gDirectory->GetList());
73   TList *list = new TList;
74   list->SetName(name);
75   TObject *objHist;
76   while((objHist=nextHist())){
77     if (!objHist->InheritsFrom("TH1")) continue;
78     ((TH1*)objHist)->SetDirectory(0); // Remove from gROOT
79     list->Add(objHist);
80     gDirectory->GetList()->Remove(objHist);
81   }
82   if(putToBrowser) gROOT->GetListOfBrowsables()->Add((TObject*)list);
83   return list;
84 }
85
86 void AliEMCALHistoUtilities::FillH1(TList *l, Int_t ind, Double_t x, Double_t w)
87 {
88   //fill 1d histogram
89   static TH1* hid=0;
90   if(l == 0) return;
91   if(ind>=0 && ind < l->GetSize()){
92     hid = (TH1*)l->At(ind);
93     hid->Fill(x,w);
94   }
95 }
96
97 void AliEMCALHistoUtilities::FillH2(TList *l, Int_t ind, Double_t x, Double_t y, Double_t w)
98 {
99   //fill 2d histogram
100   static TH2* hid=0;
101   if(l == 0) return;
102   if(ind>=0 && ind < l->GetSize()){
103     hid = (TH2*)l->At(ind);
104     hid->Fill(x,y,w);
105   }
106 }
107
108 int AliEMCALHistoUtilities::SaveListOfHists(TList *mylist,const char* name,Bool_t kSingleKey,const char* opt)
109 {
110   //write histograms to file
111   printf(" Name of out file |%s|\n", name); 
112   int save = 0;
113   if(mylist && mylist->GetSize() && strlen(name)){
114     TString nf(name); 
115     if(nf.Contains(".root") == kFALSE) nf += ".root";
116     TFile file(nf.Data(),opt);
117     TIter nextHist(mylist);
118     TObject* objHist=0;
119     int nh=0;
120     if(kSingleKey) {
121        file.cd();
122        mylist->Write(mylist->GetName(),TObject::kSingleKey);
123        mylist->ls();
124        save = 1;
125     } else {
126       while((objHist=nextHist())) { // loop over list 
127         if(objHist->InheritsFrom("TH1")) {
128           TH1* hid = (TH1*)objHist;
129           file.cd();
130           hid->Write();
131           nh++;
132           printf("Save hist. %s \n",hid ->GetName());
133         }
134       }
135       printf("%i hists. save to file -> %s\n", nh,file.GetName());
136       if(nh>0) save = 1;
137     }
138     file.Close();
139   } else {
140     printf("AliEMCALHistoUtilities::SaveListOfHists : N O  S A V I N G \n");
141   }
142   return save;
143 }
144
145 void AliEMCALHistoUtilities::AddToNameAndTitle(TH1 *h, const char *name, const char *title)
146 {
147   if(h==0) return;
148   if(name  && strlen(name))  h->SetName(Form("%s%s",h->GetName(),name));
149   if(title && strlen(title)) h->SetTitle(Form("%s%s",h->GetTitle(),title));
150 }
151
152 void AliEMCALHistoUtilities::AddToNameAndTitleToList(TList *l, const char *name, const char *title)
153 {
154   if(l==0) return;
155   if(name || title) {
156     for(int i=0; i<l->GetSize(); i++) {
157       TObject *o = l->At(i);
158       if(o->InheritsFrom("TH1")) {
159         TH1 *h = (TH1*)o;
160         AddToNameAndTitle(h, name, title);
161       }
162     }
163   }
164 }
165
166 void AliEMCALHistoUtilities::ResetListOfHists(TList *l)
167 {
168   if(l==0) return;
169
170   for(int i=0; i<l->GetSize(); i++) {
171     TH1F* h = (TH1F*)l->At(i);
172     h->Reset(); 
173   }
174 }
175
176 TLatex *AliEMCALHistoUtilities::Lat(const char *text, Float_t x,Float_t y, Int_t align, Float_t tsize, short tcolor)
177
178   double y1=y;
179   TLatex *latex = new TLatex;
180   latex->SetTextAlign(align);
181   latex->SetTextSize(tsize);
182   latex->SetTextColor(tcolor);
183   latex->DrawLatex(x, y1, text);
184   printf("<I> AliEMCALHistoUtilities::lat() -> %s gPad->GetLogy() %i\n", 
185   text, gPad->GetLogy());
186   return latex;
187 }
188
189 TGraph *AliEMCALHistoUtilities::DrawGraph(Int_t n, Double_t *x, Double_t *y, Int_t markerColor, 
190 Int_t markerStyle, const char* opt, const char* tit, const char* xTit,const char* yTit, Int_t ifun,  
191 const char *optFit, const char *fun)
192 {
193   /* Drawing options 
194   chopt='L' :  A simple polyline between every points is drawn
195   chopt='F' :  A fill area is drawn ('CF' draw a smooth fill area)
196   chopt='A' :  Axis are drawn around the graph
197   chopt='C' :  A smooth Curve is drawn
198   chopt='*' :  A Star is plotted at each point
199   chopt='P' :  Idem with the current marker
200   chopt='B' :  A Bar chart is drawn at each point
201   chopt='1' :  ylow=rwymin
202   chopt='X+' : The X-axis is drawn on the top side of the plot.
203   chopt='Y+' : The Y-axis is drawn on the right side of the plot.
204
205     Fitting options
206    The list of fit options is given in parameter option.
207       option = "W"  Set all errors to 1
208              = "U" Use a User specified fitting algorithm (via SetFCN)
209              = "Q" Quiet mode (minimum printing)
210              = "V" Verbose mode (default is between Q and V)
211              = "B" Use this option when you want to fix one or more parameters
212                    and the fitting function is like "gaus","expo","poln","landau".
213              = "R" Use the Range specified in the function range
214              = "N" Do not store the graphics function, do not draw
215              = "0" Do not plot the result of the fit. By default the fitted function
216                    is drawn unless the option"N" above is specified.
217              = "+" Add this new fitted function to the list of fitted functions
218                    (by default, any previous function is deleted)
219              = "C" In case of linear fitting, not calculate the chisquare
220                     (saves time)
221              = "F" If fitting a polN, switch to minuit fitter
222              = "ROB" In case of linear fitting, compute the LTS regression
223                      coefficients (robust(resistant) regression), using
224                      the default fraction of good points
225                "ROB=0.x" - compute the LTS regression coefficients, using
226                            0.x as a fraction of good points
227
228   */
229   printf("AliEMCALHistoUtilities::drawGraph started \n");
230   printf("Drawing opt |%s| : Fitting opt |%s|\n", opt, optFit);
231
232     TGraph *gr = new TGraph(n, x, y);
233     gr->SetMarkerColor(markerColor);
234     gr->SetLineColor(markerColor);
235     gr->SetMarkerStyle(markerStyle);
236     gr->SetTitle(tit);
237     gr->Draw(opt);
238
239     TString ctmp(opt);
240     if(ctmp.Contains("A")) {
241        gr->GetHistogram()->SetXTitle(xTit);
242        gr->GetHistogram()->SetYTitle(yTit);
243     }
244     ctmp = optFit; 
245     if(ifun>=0) {
246       TString sf("pol"); sf += ifun;
247       gr->Fit(sf.Data(),optFit);
248       printf("\n ** Fit by Polynomial of degree %i : %s **\n", ifun, sf.Data());
249     } else if(!ctmp.Contains("no",TString::kIgnoreCase)){
250       printf("\n ** ifun %i : %s **\n", ifun, fun);
251       gr->Fit(fun, optFit);
252     }
253
254     return gr;
255 }
256
257 TGraphErrors *AliEMCALHistoUtilities::DrawGraphErrors(const Int_t n,Double_t *x,Double_t *y,Double_t *ex, 
258 Double_t *ey, Int_t markerColor,  Int_t markerStyle, const char* opt, const char* tit, 
259 const char* xTit,char* yTit, Int_t ifun, const char *optFit, const char *fun)
260 {
261   printf("AliEMCALHistoUtilities::drawGraphErrors started \n");
262   printf("Drawing opt |%s| : ifun %i: Fitting opt |%s|, fun |%s|\n", 
263          opt, ifun, optFit, fun);
264
265   TGraphErrors *gr = new TGraphErrors(n, x,y,ex,ey);
266   gr->SetMarkerColor(markerColor);
267   gr->SetLineColor(markerColor);
268   gr->SetMarkerStyle(markerStyle);
269   if(tit&&strlen(tit)>0) gr->SetTitle(tit);
270
271   TString ctmp(opt);
272   if(ctmp.Contains("A")) {
273      gr->GetHistogram()->SetXTitle(xTit);
274      gr->GetHistogram()->SetYTitle(yTit);
275   }
276   if(ifun>0) {
277     if(ifun != 999) {
278       TString sf("pol"); sf += ifun;
279       gr->Fit(sf.Data(),optFit);
280       printf("\n ** Fit by Polynomial of degree %i : %s **\n", ifun, sf.Data());
281     } else {
282       gr->Fit(fun, optFit);
283       printf("\n ** Fit by %s **\n", fun);
284     }
285   } else {
286     if(strlen(optFit)) {
287       printf("\n ** ifun %i : %s **\n", ifun, fun);
288       gr->Fit(fun, optFit);
289     }
290   }
291
292   gr->Draw(opt);
293
294   return gr;
295 }
296
297 TF1* AliEMCALHistoUtilities::GetResolutionFunction(const char *opt, TString &latexName)
298 {
299   TString OPT(opt);
300   OPT.ToUpper();
301   TF1 *fres=0;
302   if      (OPT.Contains("FRES1")) {
303     fres = new TF1("fres","[0]+[1]/sqrt(x)", 0.0, 101.);
304     latexName = "#frac{#sigma_{E}}{E} = A+#frac{B}{#sqrt{E}}";
305   } else if(OPT.Contains("FRES2")) { 
306     fres = new TF1("fres","sqrt([0]*[0]+[1]*[1]/x)", 0.0, 101.);
307     latexName = "#sqrt{A^{2}+#frac{B^{2}}{E}}";
308   }
309   if(fres) {
310     fres->SetParName(0,"A");
311     fres->SetParName(1,"B");
312
313     fres->SetParameter(0, 2.0);
314     fres->SetParameter(1, 6.6);
315     fres->SetLineWidth(2);
316     fres->SetLineColor(kRed);
317   }
318   return fres;
319 }
320
321 void AliEMCALHistoUtilities::InitChain(TChain *chain, const char* nameListOfFiles, Int_t nFileMax)
322 {
323   // Read name of files from text file with nameListOfFiles and added to the chain
324   if(chain==0 || nameListOfFiles==0) return;
325  
326   ifstream in;
327   in.open(nameListOfFiles);
328   if (!in) {
329     cout << "Input file "<<nameListOfFiles<<" cannot be opened.\n";
330     return;
331   }
332
333   Int_t nfiles = 0; // number of files in chain
334   char nf[200];     // name of current file
335   while (!in.getline(nf,200).eof()) {
336     if(!in.good()) break;
337     chain->Add(nf);
338     nfiles++;
339     cout<<nfiles<<" "<<nf<<endl;
340     if(nFileMax && nfiles>=nFileMax) break;
341   }
342   cout << " \n ********** <I> Accepted file "<< nfiles << "********* \n"<<endl;
343   //  chainEsd->Print();
344   //  chain->Lookup();
345 }
346
347 AliRunLoader* AliEMCALHistoUtilities::InitKinematics(const Int_t nev, const char* galiceName)
348 {
349   static AliRunLoader *RL = 0;
350   if(RL == 0 || nev%1000==0) {
351     if(RL)  {
352       RL->UnloadgAlice();
353       delete RL;
354     }
355     RL = AliRunLoader::Open(galiceName,AliConfig::GetDefaultEventFolderName(),"read");
356     RL->LoadgAlice(); // obligatory
357   }
358   if(RL) {
359     RL->GetEvent(nev%1000);
360     RL->LoadKinematics();
361     /* Get what you need after that
362       RL->LoadHits();
363       AliStack *stack=RL->Stack();
364       AliESDCaloCluster * clus = esd->GetCaloCluster(n);
365       Int_t label = clus->GetLabel(); // what is this 
366       TParticle *primary=stack->Particle(label); 
367     */
368   }
369   return RL;
370 }
371
372 Double_t AliEMCALHistoUtilities::GetMomentum(const char* nameListOfFiles)
373 {
374   // Get momentum value from string  - like /....100GEV/.. 
375   Double_t p=-1; // negative if undefined 
376   TString st(nameListOfFiles);
377   if(st.Length()>=5) {
378     st.ToUpper();
379     Ssiz_t ind1 = st.Index("GEV"), ind2=0;
380     if(ind1>0) {
381       for(Int_t i=2; i<=10; i++) {
382         ind2 = st.Index("/",ind1-i);
383         if(ind2>0 && ind2<ind1) break;
384       }
385       TString mom  = st(ind2+1, ind1-ind2-1);
386       if(mom.IsFloat()) p = mom.Atof();
387       printf(" dir |%s| : mom |%s| : p %f : ind2,1 %i->%i\n", st.Data(), mom.Data(), p, ind2, ind1);
388     }
389   }
390   return p;
391 }
392
393 int AliEMCALHistoUtilities::ParseString(const TString &topt, TObjArray &Opt)
394
395   // Moved from AliEMCALGeometry
396   // Feb 06, 2006
397   Ssiz_t begin, index, end, end2;
398   begin = index = end = end2 = 0;
399   TRegexp separator("[^ ;,\\t\\s/]+");
400   while ( (begin < topt.Length()) && (index != kNPOS) ) {
401     // loop over given options
402     index = topt.Index(separator,&end,begin);
403     if (index >= 0 && end >= 1) {
404       TString substring(topt(index,end));
405       Opt.Add(new TObjString(substring.Data()));
406     }
407     begin += end+1;
408   }
409   return Opt.GetEntries();
410 }
411
412 // Analysis utilites
413 Bool_t AliEMCALHistoUtilities::GetLorentzVectorFromESDCluster(TLorentzVector &v, const AliESDCaloCluster* cl)
414 {
415   // May 8, 2007
416   static Double_t e=0.0;
417   static Float_t pos[3];
418   static TVector3 gpos;
419   if(cl==0) return kFALSE;
420   
421   e = Double_t(cl->E());
422   if(e<=0.0) {
423     printf(" negative cluster energy : %f \n", e);
424     return kFALSE;
425   }
426   cl->GetPosition(pos);
427   gpos.SetXYZ(Double_t(pos[0]), Double_t(pos[1]), Double_t(pos[2]));
428   gpos.SetMag(e);
429   v.SetVectM(gpos, 0.0);
430
431   return kTRUE;
432 }
433
434 Bool_t AliEMCALHistoUtilities::GetLorentzVectorFromRecPoint(TLorentzVector &v, const AliEMCALRecPoint  *rp)
435 {
436   // Jun 20, 2007
437   static Double_t e=0.0;
438   static TVector3 gpos;
439   if(rp==0) return kFALSE;
440   
441   e = Double_t(rp->GetPointEnergy());
442   if(e<=0.0) {
443     printf(" negative rec.point energy : %f \n", e);
444     return kFALSE;
445   }
446   rp->GetGlobalPosition(gpos);
447   gpos.SetMag(e);
448   v.SetVectM(gpos, 0.0);
449
450   return kTRUE;
451 }
452 //
453 //// Drawing:
454 //
455 void AliEMCALHistoUtilities::DrawHist(TH1* hid,int lineWidth,int lineColor,const char* opt, int lineStyle)
456 {
457   TString OPT;
458   OPT = opt;
459   if(!hid) return;
460   printf(" Hist. %s : option |%s| \n", hid->GetName(), opt);
461   if(hid && hid->GetEntries()>=1.) {
462     if(lineWidth) hid->SetLineWidth(lineWidth);
463     if(lineColor) hid->SetLineColor(lineColor);
464     if(lineStyle) hid->SetLineStyle(lineStyle);
465     if(OPT.Contains("stat",TString::kIgnoreCase)) hid->SetStats(kTRUE);
466     hid->Draw(opt);
467   } else {
468     if   (strcmp(opt,"empty")==0) hid->Draw();
469     else printf(" has fewer entries %i or hid is zero\n", (int) hid->GetEntries());
470   }
471 }
472
473 //
474 //// Fitting:
475 //
476 TF1* AliEMCALHistoUtilities::Gausi(const char *addName,double xmi,double xma,double N,double mean,double sig,double width)
477 { // Fit by gaus where first parameter is the number of events under ga
478   TString name("gi");
479   name += addName;
480   TF1 *F = new TF1(name.Data(), Gi, xmi, xma, 4); 
481   F->SetParNames("INTEGRAL","MEAN","SIGMA","WIDTH");
482
483   F->SetParameter(0,N);
484   F->SetParameter(1,mean);
485   F->SetParameter(2,sig);
486
487   F->FixParameter(3,width); // width of histogramm bin
488   return F;
489 }
490
491 TF1* AliEMCALHistoUtilities::Gausi(const char *addName, double xmi, double xma, TH1 *h) 
492 {
493   // You can change parameters after that if you don't like this assumption
494   if(h) return Gausi(addName, xmi, xma, h->Integral(), h->GetMean(),h->GetRMS(), h->GetBinWidth(1));
495   else  return 0; 
496 }
497
498 TF1* AliEMCALHistoUtilities::GausiPol2(const char *addName,double xmi,double xma, TF1 *g, TF1* bg)
499 { // Fit by gaus where first parameter is the number of events under ga
500   TString name("giPol2");
501   name += addName;
502   TF1 *F = new TF1(name.Data(), GiPol2, xmi, xma, 7); 
503   F->SetParNames("INTEGRAL","MEAN","SIGMA","WIDTH","a0","a1","a2");
504
505   if(g) {
506     for(int i=0; i<4; i++) F->SetParameter(i, g->GetParameter(i));
507     F->FixParameter(3,g->GetParameter(3));
508   }
509
510   if(bg) {
511     for(int i=4; i<7; i++) F->SetParameter(i, bg->GetParameter(i+4));
512   }
513   F->SetLineColor(kRed);
514   return F;
515 }
516
517 Double_t AliEMCALHistoUtilities::Gi(Double_t *x, Double_t *par)
518
519   // First parameter is integral (number events under gaus)
520   // Forth parameter (par[3]) is width of histogram bin 
521   // gaus(0) is a substitute for [0]*exp(-0.5*((x-[1])/[2])**2)
522
523   static Double_t C = TMath::Sqrt(2.*TMath::Pi()), y=0.0, f=0.0; // sqrt(2.*pi)
524
525   y  = (x[0]-par[1])/par[2];
526   f  = par[0]*par[3]/(par[2]*C) * TMath::Exp(-0.5*y*y);
527
528   return f;
529 }
530
531 Double_t AliEMCALHistoUtilities::GiPol2(Double_t *x, Double_t *par)
532
533   // First parameter is integral (number events under gaus)
534   // Forth parameter (par[3]) is width of histogram bin 
535   // gaus(0) is a substitute for [0]*exp(-0.5*((x-[1])/[2])**2)
536   // + pol2 -> 7 parameters
537   static Double_t C = TMath::Sqrt(2.*TMath::Pi()), y=0.0, f=0.0; // sqrt(2.*pi)
538
539   y  = (x[0]-par[1])/par[2];
540   f  = par[0]*par[3]/(par[2]*C) * TMath::Exp(-0.5*y*y);
541
542   f += par[4] + par[5]*x[0] + par[6]*x[0]*x[0];
543
544   return f;
545 }
546
547 // Calibration stuff
548 Double_t AliEMCALHistoUtilities::GetCorrectionCoefficientForGamma_1(const Double_t eRec)
549 {
550   // Correction to rec.energy - Jul 15, 2007
551   // E(gamma) = corr * E(eRec);
552   /* corr = p0*(eRec-7.5)+p1*(eRec-7.5)*(eRec-7.5)+p2; 0.0<eRec<10.0
553    1  p0           6.07157e-05   1.15179e-04  -0.00000e+00   1.20997e-03
554    2  p1           1.50019e-04   3.13566e-05  -0.00000e+00   1.22531e-02
555    3  p2           9.99019e-01   4.08251e-04  -0.00000e+00   1.40606e-03
556
557      corr = p3 + p4*eRec + p5*eRec*eRec
558    1  p3           9.97135e-01   5.31970e-04   1.37962e-09   1.68120e-08
559    2  p4           3.15740e-04   2.53371e-05   1.11475e-11   1.74192e-04
560    3  p5          -1.35383e-06   2.19495e-07  -5.82864e-13   4.52277e-02
561   */
562   static Double_t p0=6.07157e-05, p1=1.50019e-04, p2=9.99019e-01;
563   static Double_t p3=9.97135e-01, p4=3.15740e-04, p5=-1.35383e-06;
564   static Double_t corr=1.0;
565   if(eRec>=0.0 && eRec <=10.0) {
566     corr = p0*(eRec-7.5) + p1*(eRec-7.5)*(eRec-7.5) + p2;
567   } else if(eRec>10.0 && eRec <=105.0) {
568     corr = p3 + p4*eRec + p5*eRec*eRec;
569   }
570   //printf(" eRec %f | corr %f \n", eRec, corr);
571   return corr;
572 }
573
574 Double_t AliEMCALHistoUtilities::GetCorrectedEnergyForGamma_1(const Double_t eRec)
575 {
576   return GetCorrectionCoefficientForGamma_1(eRec) * eRec;
577 }