]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALHistoUtilities.cxx
removed static from variable declaration in AliEMCALGeometry::GetPhiCenterOfSM
[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   // Oct 15, 2007
148   if(h==0) return;
149   if(name  && strlen(name))  h->SetName(Form("%s%s",h->GetName(),name));
150   if(title && strlen(title)) h->SetTitle(Form("%s%s",h->GetTitle(),title));
151 }
152
153 void AliEMCALHistoUtilities::AddToNameAndTitleToList(TList *l, const char *name, const char *title)
154 {
155   // Oct 15, 2007
156   if(l==0) return;
157   if(name || title) {
158     for(int i=0; i<l->GetSize(); i++) {
159       TObject *o = l->At(i);
160       if(o->InheritsFrom("TH1")) {
161         TH1 *h = (TH1*)o;
162         AddToNameAndTitle(h, name, title);
163       }
164     }
165   }
166 }
167
168 void AliEMCALHistoUtilities::ResetListOfHists(TList *l)
169 {
170   // Oct 15, 2007
171   if(l==0) return;
172   for(int i=0; i<l->GetSize(); i++) {
173     TH1F* h = (TH1F*)l->At(i);
174     h->Reset(); 
175   }
176 }
177
178 TLatex *AliEMCALHistoUtilities::Lat(const char *text, Float_t x,Float_t y, Int_t align, Float_t tsize, short tcolor)
179
180   // Oct 15, 2007
181   double y1=y;
182   TLatex *latex = new TLatex;
183   latex->SetTextAlign(align);
184   latex->SetTextSize(tsize);
185   latex->SetTextColor(tcolor);
186   latex->DrawLatex(x, y1, text);
187   printf("<I> AliEMCALHistoUtilities::lat() -> %s gPad->GetLogy() %i\n", 
188   text, gPad->GetLogy());
189   return latex;
190 }
191
192 TGraph *AliEMCALHistoUtilities::DrawGraph(Int_t n, Double_t *x, Double_t *y, Int_t markerColor, 
193 Int_t markerStyle, const char* opt, const char* tit, const char* xTit,const char* yTit, Int_t ifun,  
194 const char *optFit, const char *fun)
195 {
196   /* Drawing options 
197   chopt='L' :  A simple polyline between every points is drawn
198   chopt='F' :  A fill area is drawn ('CF' draw a smooth fill area)
199   chopt='A' :  Axis are drawn around the graph
200   chopt='C' :  A smooth Curve is drawn
201   chopt='*' :  A Star is plotted at each point
202   chopt='P' :  Idem with the current marker
203   chopt='B' :  A Bar chart is drawn at each point
204   chopt='1' :  ylow=rwymin
205   chopt='X+' : The X-axis is drawn on the top side of the plot.
206   chopt='Y+' : The Y-axis is drawn on the right side of the plot.
207
208     Fitting options
209    The list of fit options is given in parameter option.
210       option = "W"  Set all errors to 1
211              = "U" Use a User specified fitting algorithm (via SetFCN)
212              = "Q" Quiet mode (minimum printing)
213              = "V" Verbose mode (default is between Q and V)
214              = "B" Use this option when you want to fix one or more parameters
215                    and the fitting function is like "gaus","expo","poln","landau".
216              = "R" Use the Range specified in the function range
217              = "N" Do not store the graphics function, do not draw
218              = "0" Do not plot the result of the fit. By default the fitted function
219                    is drawn unless the option"N" above is specified.
220              = "+" Add this new fitted function to the list of fitted functions
221                    (by default, any previous function is deleted)
222              = "C" In case of linear fitting, not calculate the chisquare
223                     (saves time)
224              = "F" If fitting a polN, switch to minuit fitter
225              = "ROB" In case of linear fitting, compute the LTS regression
226                      coefficients (robust(resistant) regression), using
227                      the default fraction of good points
228                "ROB=0.x" - compute the LTS regression coefficients, using
229                            0.x as a fraction of good points
230
231   */
232   printf("AliEMCALHistoUtilities::drawGraph started \n");
233   printf("Drawing opt |%s| : Fitting opt |%s|\n", opt, optFit);
234
235     TGraph *gr = new TGraph(n, x, y);
236     gr->SetMarkerColor(markerColor);
237     gr->SetLineColor(markerColor);
238     gr->SetMarkerStyle(markerStyle);
239     gr->SetTitle(tit);
240     gr->Draw(opt);
241
242     TString ctmp(opt);
243     if(ctmp.Contains("A")) {
244        gr->GetHistogram()->SetXTitle(xTit);
245        gr->GetHistogram()->SetYTitle(yTit);
246     }
247     ctmp = optFit; 
248     if(ifun>=0) {
249       TString sf("pol"); sf += ifun;
250       gr->Fit(sf.Data(),optFit);
251       printf("\n ** Fit by Polynomial of degree %i : %s **\n", ifun, sf.Data());
252     } else if(!ctmp.Contains("no",TString::kIgnoreCase)){
253       printf("\n ** ifun %i : %s **\n", ifun, fun);
254       gr->Fit(fun, optFit);
255     }
256
257     return gr;
258 }
259
260 TGraphErrors *AliEMCALHistoUtilities::DrawGraphErrors(const Int_t n,Double_t *x,Double_t *y,Double_t *ex, 
261 Double_t *ey, Int_t markerColor,  Int_t markerStyle, const char* opt, const char* tit, 
262 const char* xTit,char* yTit, Int_t ifun, const char *optFit, const char *fun)
263 {
264   // Oct 15, 2007
265   printf("AliEMCALHistoUtilities::drawGraphErrors started \n");
266   printf("Drawing opt |%s| : ifun %i: Fitting opt |%s|, fun |%s|\n", 
267          opt, ifun, optFit, fun);
268
269   TGraphErrors *gr = new TGraphErrors(n, x,y,ex,ey);
270   gr->SetMarkerColor(markerColor);
271   gr->SetLineColor(markerColor);
272   gr->SetMarkerStyle(markerStyle);
273   if(tit&&strlen(tit)>0) gr->SetTitle(tit);
274
275   TString ctmp(opt);
276   if(ctmp.Contains("A")) {
277      gr->GetHistogram()->SetXTitle(xTit);
278      gr->GetHistogram()->SetYTitle(yTit);
279   }
280   if(ifun>0) {
281     if(ifun != 999) {
282       TString sf("pol"); sf += ifun;
283       gr->Fit(sf.Data(),optFit);
284       printf("\n ** Fit by Polynomial of degree %i : %s **\n", ifun, sf.Data());
285     } else {
286       gr->Fit(fun, optFit);
287       printf("\n ** Fit by %s **\n", fun);
288     }
289   } else {
290     if(strlen(optFit)) {
291       printf("\n ** ifun %i : %s **\n", ifun, fun);
292       gr->Fit(fun, optFit);
293     }
294   }
295
296   gr->Draw(opt);
297
298   return gr;
299 }
300
301 TF1* AliEMCALHistoUtilities::GetResolutionFunction(const char *opt, TString &latexName)
302 {
303   // Oct 15, 2007
304   TString sopt(opt);
305   sopt.ToUpper();
306   TF1 *fres=0;
307   if      (sopt.Contains("FRES1")) {
308     fres = new TF1("fres","[0]+[1]/sqrt(x)", 0.0, 101.);
309     latexName = "#frac{#sigma_{E}}{E} = A+#frac{B}{#sqrt{E}}";
310   } else if(sopt.Contains("FRES2")) { 
311     fres = new TF1("fres","sqrt([0]*[0]+[1]*[1]/x)", 0.0, 101.);
312     latexName = "#sqrt{A^{2}+#frac{B^{2}}{E}}";
313   }
314   if(fres) {
315     fres->SetParName(0,"A");
316     fres->SetParName(1,"B");
317
318     fres->SetParameter(0, 2.0);
319     fres->SetParameter(1, 6.6);
320     fres->SetLineWidth(2);
321     fres->SetLineColor(kRed);
322   }
323   return fres;
324 }
325
326 void AliEMCALHistoUtilities::InitChain(TChain *chain, const char* nameListOfFiles, Int_t nFileMax)
327 {
328   // Read name of files from text file with nameListOfFiles and added to the chain
329   if(chain==0 || nameListOfFiles==0) return;
330  
331   ifstream in;
332   in.open(nameListOfFiles);
333   if (!in) {
334     cout << "Input file "<<nameListOfFiles<<" cannot be opened.\n";
335     return;
336   }
337
338   Int_t nfiles = 0; // number of files in chain
339   char nf[200];     // name of current file
340   while (!in.getline(nf,200).eof()) {
341     if(!in.good()) break;
342     chain->Add(nf);
343     nfiles++;
344     cout<<nfiles<<" "<<nf<<endl;
345     if(nFileMax && nfiles>=nFileMax) break;
346   }
347   cout << " \n ********** <I> Accepted file "<< nfiles << "********* \n"<<endl;
348   //  chainEsd->Print();
349   //  chain->Lookup();
350 }
351
352 AliRunLoader* AliEMCALHistoUtilities::InitKinematics(const Int_t nev, const char* galiceName)
353 {
354   // Oct 15, 2007
355   static AliRunLoader *rl = 0;
356   if(rl == 0 || nev%1000==0) {
357     if(rl)  {
358       rl->UnloadgAlice();
359       delete rl;
360     }
361     rl = AliRunLoader::Open(galiceName,AliConfig::GetDefaultEventFolderName(),"read");
362     rl->LoadgAlice(); // obligatory
363   }
364   if(rl) {
365     rl->GetEvent(nev%1000);
366     rl->LoadKinematics();
367     /* Get what you need after that
368       rl->LoadHits();
369       AliStack *stack=rl->Stack();
370       AliESDCaloCluster * clus = esd->GetCaloCluster(n);
371       Int_t label = clus->GetLabel(); // what is this 
372       TParticle *primary=stack->Particle(label); 
373     */
374   }
375   return rl;
376 }
377
378 Double_t AliEMCALHistoUtilities::GetMomentum(const char* nameListOfFiles)
379 {
380   // Get momentum value from string  - like /....100GEV/.. 
381   Double_t p=-1; // negative if undefined 
382   TString st(nameListOfFiles);
383   if(st.Length()>=5) {
384     st.ToUpper();
385     Ssiz_t ind1 = st.Index("GEV"), ind2=0;
386     if(ind1>0) {
387       for(Int_t i=2; i<=10; i++) {
388         ind2 = st.Index("/",ind1-i);
389         if(ind2>0 && ind2<ind1) break;
390       }
391       TString mom  = st(ind2+1, ind1-ind2-1);
392       if(mom.IsFloat()) p = mom.Atof();
393       printf(" dir |%s| : mom |%s| : p %f : ind2,1 %i->%i\n", st.Data(), mom.Data(), p, ind2, ind1);
394     }
395   }
396   return p;
397 }
398
399 int AliEMCALHistoUtilities::ParseString(const TString &topt, TObjArray &Opt)
400
401   // Moved from AliEMCALGeometry
402   // Feb 06, 2006
403   Ssiz_t begin, index, end, end2;
404   begin = index = end = end2 = 0;
405   TRegexp separator("[^ ;,\\t\\s/]+");
406   while ( (begin < topt.Length()) && (index != kNPOS) ) {
407     // loop over given options
408     index = topt.Index(separator,&end,begin);
409     if (index >= 0 && end >= 1) {
410       TString substring(topt(index,end));
411       Opt.Add(new TObjString(substring.Data()));
412     }
413     begin += end+1;
414   }
415   return Opt.GetEntries();
416 }
417
418 // Analysis utilites
419 Bool_t AliEMCALHistoUtilities::GetLorentzVectorFromESDCluster(TLorentzVector &v, const AliESDCaloCluster* cl)
420 {
421   // May 8, 2007
422   static Double_t e=0.0;
423   static Float_t pos[3];
424   static TVector3 gpos;
425   if(cl==0) return kFALSE;
426   
427   e = Double_t(cl->E());
428   if(e<=0.0) {
429     printf(" negative cluster energy : %f \n", e);
430     return kFALSE;
431   }
432   cl->GetPosition(pos);
433   gpos.SetXYZ(Double_t(pos[0]), Double_t(pos[1]), Double_t(pos[2]));
434   gpos.SetMag(e);
435   v.SetVectM(gpos, 0.0);
436
437   return kTRUE;
438 }
439
440 Bool_t AliEMCALHistoUtilities::GetLorentzVectorFromRecPoint(TLorentzVector &v, const AliEMCALRecPoint  *rp)
441 {
442   // Jun 20, 2007
443   static Double_t e=0.0;
444   static TVector3 gpos;
445   if(rp==0) return kFALSE;
446   
447   e = Double_t(rp->GetPointEnergy());
448   if(e<=0.0) {
449     printf(" negative rec.point energy : %f \n", e);
450     return kFALSE;
451   }
452   rp->GetGlobalPosition(gpos);
453   gpos.SetMag(e);
454   v.SetVectM(gpos, 0.0);
455
456   return kTRUE;
457 }
458 //
459 //// Drawing:
460 //
461 void AliEMCALHistoUtilities::DrawHist(TH1* hid,int lineWidth,int lineColor,const char* opt, int lineStyle)
462 {
463   // Oct 15, 2007
464   TString sopt;
465   sopt = opt;
466   if(!hid) return;
467   printf(" Hist. %s : option |%s| \n", hid->GetName(), opt);
468   if(hid && hid->GetEntries()>=1.) {
469     if(lineWidth) hid->SetLineWidth(lineWidth);
470     if(lineColor) hid->SetLineColor(lineColor);
471     if(lineStyle) hid->SetLineStyle(lineStyle);
472     if(sopt.Contains("stat",TString::kIgnoreCase)) hid->SetStats(kTRUE);
473     hid->Draw(opt);
474   } else {
475     if   (strcmp(opt,"empty")==0) hid->Draw();
476     else printf(" has fewer entries %i or hid is zero\n", (int) hid->GetEntries());
477   }
478 }
479
480 //
481 //// Fitting:
482 //
483 TF1* AliEMCALHistoUtilities::Gausi(const char *addName,double xmi,double xma,double N,double mean,double sig,double width)
484
485   // Fit by gaus where first parameter is the number of events under ga
486   // Oct 15, 2007
487   TString name("gi");
488   name += addName;
489   TF1 *f = new TF1(name.Data(), Gi, xmi, xma, 4); 
490   f->SetParNames("INTEGRAL","MEAN","SIGMA","WIDTH");
491
492   f->SetParameter(0,N);
493   f->SetParameter(1,mean);
494   f->SetParameter(2,sig);
495
496   f->FixParameter(3,width); // width of histogramm bin
497   return f;
498 }
499
500 TF1* AliEMCALHistoUtilities::Gausi(const char *addName, double xmi, double xma, TH1 *h) 
501 {
502   // You can change parameters after that if you don't like this assumption
503   if(h) return Gausi(addName, xmi, xma, h->Integral(), h->GetMean(),h->GetRMS(), h->GetBinWidth(1));
504   else  return 0; 
505 }
506
507 TF1* AliEMCALHistoUtilities::GausiPol2(const char *addName,double xmi,double xma, TF1 *g, TF1* bg)
508
509   // Fit by gaus where first parameter is the number of events under ga
510   TString name("giPol2");
511   name += addName;
512   TF1 *f = new TF1(name.Data(), GiPol2, xmi, xma, 7); 
513   f->SetParNames("INTEGRAL","MEAN","SIGMA","WIDTH","a0","a1","a2");
514
515   if(g) {
516     for(int i=0; i<4; i++) f->SetParameter(i, g->GetParameter(i));
517     f->FixParameter(3,g->GetParameter(3));
518   }
519
520   if(bg) {
521     for(int i=4; i<7; i++) f->SetParameter(i, bg->GetParameter(i+4));
522   }
523   f->SetLineColor(kRed);
524   return f;
525 }
526
527 Double_t AliEMCALHistoUtilities::Gi(Double_t *x, Double_t *par)
528
529   // First parameter is integral (number events under gaus)
530   // Forth parameter (par[3]) is width of histogram bin 
531   // gaus(0) is a substitute for [0]*exp(-0.5*((x-[1])/[2])**2)
532
533   static Double_t c = TMath::Sqrt(2.*TMath::Pi()), y=0.0, f=0.0; // sqrt(2.*pi)
534
535   y  = (x[0]-par[1])/par[2];
536   f  = par[0]*par[3]/(par[2]*c) * TMath::Exp(-0.5*y*y);
537
538   return f;
539 }
540
541 Double_t AliEMCALHistoUtilities::GiPol2(Double_t *x, Double_t *par)
542
543   // First parameter is integral (number events under gaus)
544   // Forth parameter (par[3]) is width of histogram bin 
545   // gaus(0) is a substitute for [0]*exp(-0.5*((x-[1])/[2])**2)
546   // + pol2 -> 7 parameters
547   static Double_t c = TMath::Sqrt(2.*TMath::Pi()), y=0.0, f=0.0; // sqrt(2.*pi)
548
549   y  = (x[0]-par[1])/par[2];
550   f  = par[0]*par[3]/(par[2]*c) * TMath::Exp(-0.5*y*y);
551
552   f += par[4] + par[5]*x[0] + par[6]*x[0]*x[0];
553
554   return f;
555 }
556
557 // Calibration stuff
558 Double_t AliEMCALHistoUtilities::GetCorrectionCoefficientForGamma1(const Double_t eRec)
559 {
560   // Correction to rec.energy - Jul 15, 2007
561   // E(gamma) = corr * E(eRec);
562   /* corr = p0*(eRec-7.5)+p1*(eRec-7.5)*(eRec-7.5)+p2; 0.0<eRec<10.0
563    1  p0           6.07157e-05   1.15179e-04  -0.00000e+00   1.20997e-03
564    2  p1           1.50019e-04   3.13566e-05  -0.00000e+00   1.22531e-02
565    3  p2           9.99019e-01   4.08251e-04  -0.00000e+00   1.40606e-03
566
567      corr = p3 + p4*eRec + p5*eRec*eRec
568    1  p3           9.97135e-01   5.31970e-04   1.37962e-09   1.68120e-08
569    2  p4           3.15740e-04   2.53371e-05   1.11475e-11   1.74192e-04
570    3  p5          -1.35383e-06   2.19495e-07  -5.82864e-13   4.52277e-02
571   */
572   static Double_t p0=6.07157e-05, p1=1.50019e-04, p2=9.99019e-01;
573   static Double_t p3=9.97135e-01, p4=3.15740e-04, p5=-1.35383e-06;
574   static Double_t corr=1.0;
575   if(eRec>=0.0 && eRec <=10.0) {
576     corr = p0*(eRec-7.5) + p1*(eRec-7.5)*(eRec-7.5) + p2;
577   } else if(eRec>10.0 && eRec <=105.0) {
578     corr = p3 + p4*eRec + p5*eRec*eRec;
579   }
580   //printf(" eRec %f | corr %f \n", eRec, corr);
581   return corr;
582 }
583
584 Double_t AliEMCALHistoUtilities::GetCorrectedEnergyForGamma1(const Double_t eRec)
585 {
586   return GetCorrectionCoefficientForGamma1(eRec) * eRec;
587 }