]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALHistoUtilities.cxx
cout replaced by AliDebug (Plamen)
[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 #include <TArrayF.h>
50
51 #include "AliESDCaloCluster.h"
52 #include "AliEMCALRecPoint.h"
53 #include "AliRunLoader.h"
54 #include "AliHeader.h"
55 #include "AliGenEventHeader.h"
56 #include "AliGenPythiaEventHeader.h"
57
58 using namespace std;
59
60 ClassImp(AliEMCALHistoUtilities)
61
62 AliEMCALHistoUtilities::AliEMCALHistoUtilities(const char *name, const char *tit) : TNamed(name,tit)
63 {
64   // constructor
65 }
66
67 AliEMCALHistoUtilities::~AliEMCALHistoUtilities()
68 {
69         //destructor
70 }  
71
72 TList* AliEMCALHistoUtilities::MoveHistsToList(const char* name, Bool_t putToBrowser)
73 {
74   // Move HIST to list
75   gROOT->cd();
76   TIter nextHist(gDirectory->GetList());
77   TList *list = new TList;
78   list->SetName(name);
79   TObject *objHist;
80   while((objHist=nextHist())){
81     if (!objHist->InheritsFrom("TH1")) continue;
82     ((TH1*)objHist)->SetDirectory(0); // Remove from gROOT
83     list->Add(objHist);
84     gDirectory->GetList()->Remove(objHist);
85   }
86   if(putToBrowser) gROOT->GetListOfBrowsables()->Add((TObject*)list);
87   return list;
88 }
89
90 void AliEMCALHistoUtilities::FillH1(TList *l, Int_t ind, Double_t x, Double_t w, Double_t error)
91 {
92   //fill 1d histogram
93   static TH1* hid=0;
94   static int  bin=0;
95   if(l == 0) return;
96
97   if(ind>=0 && ind < l->GetSize()){
98     hid = (TH1*)l->At(ind);
99     if(error <= 0.0) { // standard way
100       hid->Fill(x,w);
101     } else{
102       bin = hid->FindBin(x);
103       hid->SetBinContent(bin,w);
104       hid->SetBinError(bin,error);
105     }
106   }
107 }
108
109 void AliEMCALHistoUtilities::FillH2(TList *l, Int_t ind, Double_t x, Double_t y, Double_t w)
110 {
111   //fill 2d histogram
112   static TH2* hid=0;
113   if(l == 0) return;
114   if(ind>=0 && ind < l->GetSize()){
115     hid = (TH2*)l->At(ind);
116     hid->Fill(x,y,w);
117   }
118 }
119
120 int AliEMCALHistoUtilities::SaveListOfHists(TList *mylist,const char* name,Bool_t kSingleKey,const char* opt)
121 {
122   //write histograms to file
123   printf(" Name of out file |%s|\n", name); 
124   int save = 0;
125   if(mylist && mylist->GetSize() && strlen(name)){
126     TString nf(name); 
127     if(nf.Contains(".root") == kFALSE) nf += ".root";
128     TFile file(nf.Data(),opt);
129     TIter nextHist(mylist);
130     TObject* objHist=0;
131     int nh=0;
132     if(kSingleKey) {
133        file.cd();
134        mylist->Write(mylist->GetName(),TObject::kSingleKey);
135        mylist->ls();
136        save = 1;
137     } else {
138       while((objHist=nextHist())) { // loop over list 
139         if(objHist->InheritsFrom("TH1")) {
140           TH1* hid = (TH1*)objHist;
141           file.cd();
142           hid->Write();
143           nh++;
144           printf("Save hist. %s \n",hid ->GetName());
145         }
146       }
147       printf("%i hists. save to file -> %s\n", nh,file.GetName());
148       if(nh>0) save = 1;
149     }
150     file.Close();
151   } else {
152     printf("AliEMCALHistoUtilities::SaveListOfHists : N O  S A V I N G \n");
153   }
154   return save;
155 }
156
157 void AliEMCALHistoUtilities::AddToNameAndTitle(TH1 *h, const char *name, const char *title)
158 {
159   // Oct 15, 2007
160   if(h==0) return;
161   if(name  && strlen(name))  h->SetName(Form("%s%s",h->GetName(),name));
162   if(title && strlen(title)) h->SetTitle(Form("%s%s",h->GetTitle(),title));
163 }
164
165 void AliEMCALHistoUtilities::AddToNameAndTitleToList(TList *l, const char *name, const char *title)
166 {
167   // Oct 15, 2007
168   if(l==0) return;
169   if(name || title) {
170     for(int i=0; i<l->GetSize(); i++) {
171       TObject *o = l->At(i);
172       if(o->InheritsFrom("TH1")) {
173         TH1 *h = (TH1*)o;
174         AddToNameAndTitle(h, name, title);
175       }
176     }
177   }
178 }
179
180 void AliEMCALHistoUtilities::ResetListOfHists(TList *l)
181 {
182   // Oct 15, 2007
183   if(l==0) return;
184   for(int i=0; i<l->GetSize(); i++) {
185     TH1F* h = (TH1F*)l->At(i);
186     h->Reset(); 
187   }
188 }
189
190 TLatex *AliEMCALHistoUtilities::Lat(const char *text, Float_t x,Float_t y, Int_t align, Float_t tsize, short tcolor)
191
192   // Oct 15, 2007
193   double y1=y;
194   TLatex *latex = new TLatex;
195   latex->SetTextAlign(align);
196   latex->SetTextSize(tsize);
197   latex->SetTextColor(tcolor);
198   latex->DrawLatex(x, y1, text);
199   printf("<I> AliEMCALHistoUtilities::lat() -> %s gPad->GetLogy() %i\n", 
200   text, gPad->GetLogy());
201   return latex;
202 }
203
204 TGraph *AliEMCALHistoUtilities::DrawGraph(Int_t n, Double_t *x, Double_t *y, Int_t markerColor, 
205 Int_t markerStyle, const char* opt, const char* tit, const char* xTit,const char* yTit, Int_t ifun,  
206 const char *optFit, const char *fun)
207 {
208   /* Drawing options 
209   chopt='L' :  A simple polyline between every points is drawn
210   chopt='F' :  A fill area is drawn ('CF' draw a smooth fill area)
211   chopt='A' :  Axis are drawn around the graph
212   chopt='C' :  A smooth Curve is drawn
213   chopt='*' :  A Star is plotted at each point
214   chopt='P' :  Idem with the current marker
215   chopt='B' :  A Bar chart is drawn at each point
216   chopt='1' :  ylow=rwymin
217   chopt='X+' : The X-axis is drawn on the top side of the plot.
218   chopt='Y+' : The Y-axis is drawn on the right side of the plot.
219
220     Fitting options
221    The list of fit options is given in parameter option.
222       option = "W"  Set all errors to 1
223              = "U" Use a User specified fitting algorithm (via SetFCN)
224              = "Q" Quiet mode (minimum printing)
225              = "V" Verbose mode (default is between Q and V)
226              = "B" Use this option when you want to fix one or more parameters
227                    and the fitting function is like "gaus","expo","poln","landau".
228              = "R" Use the Range specified in the function range
229              = "N" Do not store the graphics function, do not draw
230              = "0" Do not plot the result of the fit. By default the fitted function
231                    is drawn unless the option"N" above is specified.
232              = "+" Add this new fitted function to the list of fitted functions
233                    (by default, any previous function is deleted)
234              = "C" In case of linear fitting, not calculate the chisquare
235                     (saves time)
236              = "F" If fitting a polN, switch to minuit fitter
237              = "ROB" In case of linear fitting, compute the LTS regression
238                      coefficients (robust(resistant) regression), using
239                      the default fraction of good points
240                "ROB=0.x" - compute the LTS regression coefficients, using
241                            0.x as a fraction of good points
242
243   */
244   printf("AliEMCALHistoUtilities::drawGraph started \n");
245   printf("Drawing opt |%s| : Fitting opt |%s|\n", opt, optFit);
246
247     TGraph *gr = new TGraph(n, x, y);
248     gr->SetMarkerColor(markerColor);
249     gr->SetLineColor(markerColor);
250     gr->SetMarkerStyle(markerStyle);
251     gr->SetTitle(tit);
252     gr->Draw(opt);
253
254     TString ctmp(opt);
255     if(ctmp.Contains("A")) {
256        gr->GetHistogram()->SetXTitle(xTit);
257        gr->GetHistogram()->SetYTitle(yTit);
258     }
259     ctmp = optFit; 
260     if(ifun>=0) {
261       TString sf("pol"); sf += ifun;
262       gr->Fit(sf.Data(),optFit);
263       printf("\n ** Fit by Polynomial of degree %i : %s **\n", ifun, sf.Data());
264     } else if(!ctmp.Contains("no",TString::kIgnoreCase)){
265       printf("\n ** ifun %i : %s **\n", ifun, fun);
266       gr->Fit(fun, optFit);
267     }
268
269     return gr;
270 }
271
272 TGraphErrors *AliEMCALHistoUtilities::DrawGraphErrors(const Int_t n,Double_t *x,Double_t *y,Double_t *ex, 
273 Double_t *ey, Int_t markerColor,  Int_t markerStyle, const char* opt, const char* tit, 
274 const char* xTit,const char* yTit, Int_t ifun, const char *optFit, const char *fun)
275 {
276   // Oct 15, 2007
277   printf("AliEMCALHistoUtilities::drawGraphErrors started \n");
278   printf("Drawing opt |%s| : ifun %i: Fitting opt |%s|, fun |%s|\n", 
279          opt, ifun, optFit, fun);
280
281   TGraphErrors *gr = new TGraphErrors(n, x,y,ex,ey);
282   gr->SetMarkerColor(markerColor);
283   gr->SetLineColor(markerColor);
284   gr->SetMarkerStyle(markerStyle);
285   if(tit&&strlen(tit)>0) gr->SetTitle(tit);
286
287   TString ctmp(opt);
288   if(ctmp.Contains("A")) {
289      gr->GetHistogram()->SetXTitle(xTit);
290      gr->GetHistogram()->SetYTitle(yTit);
291   }
292   if(ifun>0) {
293     if(ifun != 999) {
294       TString sf("pol"); sf += ifun;
295       gr->Fit(sf.Data(),optFit);
296       printf("\n ** Fit by Polynomial of degree %i : %s **\n", ifun, sf.Data());
297     } else {
298       gr->Fit(fun, optFit);
299       printf("\n ** Fit by %s **\n", fun);
300     }
301   } else {
302     if(strlen(optFit)) {
303       printf("\n ** ifun %i : %s **\n", ifun, fun);
304       gr->Fit(fun, optFit);
305     }
306   }
307
308   gr->Draw(opt);
309
310   return gr;
311 }
312
313 TF1* AliEMCALHistoUtilities::GetResolutionFunction(const char *opt, TString &latexName)
314 {
315   // Oct 15, 2007
316   TString sopt(opt);
317   sopt.ToUpper();
318   TF1 *fres=0;
319   if      (sopt.Contains("FRES1")) {
320     fres = new TF1("fres","[0]+[1]/sqrt(x)", 0.0, 101.);
321     latexName = "#frac{#sigma_{E}}{E} = A+#frac{B}{#sqrt{E}}";
322   } else if(sopt.Contains("FRES2")) { 
323     fres = new TF1("fres","sqrt([0]*[0]+[1]*[1]/x)", 0.0, 101.);
324     latexName = "#sqrt{A^{2}+#frac{B^{2}}{E}}";
325   }
326   if(fres) {
327     fres->SetParName(0,"A");
328     fres->SetParName(1,"B");
329
330     fres->SetParameter(0, 2.0);
331     fres->SetParameter(1, 6.6);
332     fres->SetLineWidth(2);
333     fres->SetLineColor(kRed);
334   }
335   return fres;
336 }
337
338 void AliEMCALHistoUtilities::InitChain(TChain *chain, const char* nameListOfFiles, Int_t nFileMax)
339 {
340   // Read name of files from text file with nameListOfFiles and added to the chain
341   if(chain==0 || nameListOfFiles==0) return;
342  
343   ifstream fin;
344   fin.open(nameListOfFiles);
345   if (!fin.is_open()) {
346     cout << "Input file "<<nameListOfFiles<<" cannot be opened.\n";
347     return;
348   }
349
350   Int_t nfiles = 0; // number of files in chain
351   char nf[200];     // name of current file
352   while (!fin.getline(nf,200).eof()) {
353     if(!fin.good()) break;
354     chain->Add(nf);
355     nfiles++;
356     cout<<nfiles<<" "<<nf<<endl;
357     if(nFileMax && nfiles>=nFileMax) break;
358   }
359   fin.close();
360   //
361   cout << " \n ********** <I> Accepted file "<< nfiles << "********* \n"<<endl;
362   //  chainEsd->Print();
363   //  chain->Lookup();
364 }
365
366 AliRunLoader* AliEMCALHistoUtilities::InitKinematics(const Int_t nev, const char* galiceName)
367 {
368   // Oct 15, 2007
369   // nev == 0 - new file
370   static AliRunLoader *rl = 0;
371
372   if((rl == 0 || nev==0) && galiceName) {
373     //printf("<I> AliEMCALHistoUtilities::InitKinematics() : nev %i : rl %p : %s (IN)\n", 
374         // nev, rl, galiceName);  
375     if(rl)  {
376       rl->UnloadgAlice();
377       delete rl;
378     }
379     rl = AliRunLoader::Open(galiceName,AliConfig::GetDefaultEventFolderName(),"read");
380     rl->LoadgAlice(); // obligatory
381     //printf("<I> AliEMCALHistoUtilities::InitKinematics() : nev %i : rl %p : %s (OUT)\n", 
382          //nev, rl, galiceName);  
383   }
384   if(rl) {
385     rl->GetEvent(nev);
386     rl->LoadKinematics();
387     /* Get what you need after that
388       rl->LoadHits();
389       AliStack *stack=rl->Stack();
390       AliESDCaloCluster * clus = esd->GetCaloCluster(n);
391       Int_t label = clus->GetLabel(); // what is this 
392       TParticle *primary=stack->Particle(label); 
393     */
394   }
395   return rl;
396 }
397
398 AliRunLoader* AliEMCALHistoUtilities::GetRunLoader(const Int_t nev, const Char_t* galiceName,
399                                        const Char_t* eventFolderName, AliRunLoader* rlOld)
400 {
401   // Nov 26, 2007
402   // nev == 0 - new file
403   AliRunLoader *rl = 0;
404
405   if(nev==0 && galiceName) {
406     printf("<I> AliEMCALHistoUtilities::GetLoader() : nev %i : %s (IN)\n", 
407            nev, galiceName);  
408     if(rlOld) {
409       rlOld->UnloadgAlice();
410       delete rlOld;
411     }
412     TString folderName(eventFolderName);
413     if(folderName.Length() == 0) folderName = AliConfig::GetDefaultEventFolderName();
414     rl = AliRunLoader::Open(galiceName, folderName.Data(),"read");
415     rl->LoadgAlice(); // obligatory
416     printf("<I> AliEMCALHistoUtilities::GetLoader() : nev %i : %s : %s (OUT)\n", 
417            nev, galiceName, folderName.Data());  
418   } else {
419     rl = rlOld;
420   }
421   if(rl) rl->LoadgAlice(); // obligatory
422   // if(rl) rl->GetEvent(nev);
423   return rl;
424 }
425
426 Double_t AliEMCALHistoUtilities::GetMomentum(const char* nameListOfFiles)
427 {
428   // Get momentum value from string  - like /....100GEV/.. 
429   Double_t p=-1; // negative if undefined 
430   TString st(nameListOfFiles);
431   if(st.Length()>=5) {
432     st.ToUpper();
433     Ssiz_t ind1 = st.Index("GEV"), ind2=0;
434     if(ind1>0) {
435       for(Int_t i=2; i<=10; i++) {
436         ind2 = st.Index("/",ind1-i);
437         if(ind2>0 && ind2<ind1) break;
438       }
439       TString mom  = st(ind2+1, ind1-ind2-1);
440       if(mom.IsFloat()) p = mom.Atof();
441       printf(" dir |%s| : mom |%s| : p %f : ind2,1 %i->%i\n", st.Data(), mom.Data(), p, ind2, ind1);
442     }
443   }
444   return p;
445 }
446
447 int AliEMCALHistoUtilities::ParseString(const TString &topt, TObjArray &Opt)
448
449   // Moved from AliEMCALGeometry
450   // Feb 06, 2006
451   Ssiz_t begin, index, end, end2;
452   begin = index = end = end2 = 0;
453   TRegexp separator("[^ ;,\\t\\s/]+");
454   while ( (begin < topt.Length()) && (index != kNPOS) ) {
455     // loop over given options
456     index = topt.Index(separator,&end,begin);
457     if (index >= 0 && end >= 1) {
458       TString substring(topt(index,end));
459       Opt.Add(new TObjString(substring.Data()));
460     }
461     begin += end+1;
462   }
463   return Opt.GetEntries();
464 }
465
466 // Analysis utilites
467 Bool_t AliEMCALHistoUtilities::GetLorentzVectorFromESDCluster(TLorentzVector &v, const AliESDCaloCluster* cl)
468 {
469   // May 8, 2007
470   static Double_t e=0.0;
471   static Float_t pos[3];
472   static TVector3 gpos;
473   if(cl==0) return kFALSE;
474   
475   e = Double_t(cl->E());
476   if(e<=0.0) {
477     printf(" negative cluster energy : %f \n", e);
478     return kFALSE;
479   }
480   cl->GetPosition(pos);
481   gpos.SetXYZ(Double_t(pos[0]), Double_t(pos[1]), Double_t(pos[2]));
482   gpos.SetMag(e);
483   v.SetVectM(gpos, 0.0);
484
485   return kTRUE;
486 }
487
488 Bool_t AliEMCALHistoUtilities::GetLorentzVectorFromRecPoint(TLorentzVector &v, const AliEMCALRecPoint  *rp)
489 {
490   // Jun 20, 2007
491   static Double_t e=0.0;
492   static TVector3 gpos;
493   if(rp==0) return kFALSE;
494   
495   e = Double_t(rp->GetPointEnergy());
496   if(e<=0.0) {
497     printf(" negative rec.point energy : %f \n", e);
498     return kFALSE;
499   }
500   rp->GetGlobalPosition(gpos);
501   gpos.SetMag(e);
502   v.SetVectM(gpos, 0.0);
503
504   return kTRUE;
505 }
506 //
507 //// Drawing:
508 //
509 void AliEMCALHistoUtilities::DrawHist(TH1* hid,int lineWidth,int lineColor,const char* opt, int lineStyle)
510 {
511   // Oct 15, 2007
512   TString sopt;
513   sopt = opt;
514   if(!hid) return;
515   printf(" Hist. %s : option |%s| \n", hid->GetName(), opt);
516   if(hid && hid->GetEntries()>=1.) {
517     if(lineWidth) hid->SetLineWidth(lineWidth);
518     if(lineColor) hid->SetLineColor(lineColor);
519     if(lineStyle) hid->SetLineStyle(lineStyle);
520     if(sopt.Contains("stat",TString::kIgnoreCase)) hid->SetStats(kTRUE);
521     hid->Draw(opt);
522   } else {
523     if   (strcmp(opt,"empty")==0) hid->Draw();
524     else printf(" has fewer entries %i or hid is zero\n", (int) hid->GetEntries());
525   }
526 }
527
528 //
529 //// Fitting:
530 //
531 TF1* AliEMCALHistoUtilities::Gausi(const char *addName,double xmi,double xma,double N,double mean,double sig,double width)
532
533   // Fit by gaus where first parameter is the number of events under ga
534   // Oct 15, 2007
535   TString name("gi");
536   name += addName;
537   TF1 *f = new TF1(name.Data(), Gi, xmi, xma, 4); 
538   f->SetParNames("INTEGRAL","MEAN","SIGMA","WIDTH");
539
540   f->SetParameter(0,N);
541   f->SetParameter(1,mean);
542   f->SetParameter(2,sig);
543
544   f->FixParameter(3,width); // width of histogramm bin
545   return f;
546 }
547
548 TF1* AliEMCALHistoUtilities::Gausi(const char *addName, double xmi, double xma, TH1 *h) 
549 {
550   // You can change parameters after that if you don't like this assumption
551   if(h) return Gausi(addName, xmi, xma, h->Integral(), h->GetMean(),h->GetRMS(), h->GetBinWidth(1));
552   else  return 0; 
553 }
554
555 TF1* AliEMCALHistoUtilities::GausiPol2(const char *addName,double xmi,double xma, TF1 *g, TF1* bg)
556
557   // Fit by gaus where first parameter is the number of events under ga
558   TString name("giPol2");
559   name += addName;
560   TF1 *f = new TF1(name.Data(), GiPol2, xmi, xma, 7); 
561   f->SetParNames("INTEGRAL","MEAN","SIGMA","WIDTH","a0","a1","a2");
562
563   if(g) {
564     for(int i=0; i<4; i++) f->SetParameter(i, g->GetParameter(i));
565     f->FixParameter(3,g->GetParameter(3));
566   }
567
568   if(bg) {
569     for(int i=4; i<7; i++) f->SetParameter(i, bg->GetParameter(i+4));
570   }
571   f->SetLineColor(kRed);
572   return f;
573 }
574
575 Double_t AliEMCALHistoUtilities::Gi(Double_t *x, Double_t *par)
576
577   // First parameter is integral (number events under gaus)
578   // Forth parameter (par[3]) is width of histogram bin 
579   // gaus(0) is a substitute for [0]*exp(-0.5*((x-[1])/[2])**2)
580
581   static Double_t c = TMath::Sqrt(2.*TMath::Pi()), y=0.0, f=0.0; // sqrt(2.*pi)
582
583   y  = (x[0]-par[1])/par[2];
584   f  = par[0]*par[3]/(par[2]*c) * TMath::Exp(-0.5*y*y);
585
586   return f;
587 }
588
589 Double_t AliEMCALHistoUtilities::GiPol2(Double_t *x, Double_t *par)
590
591   // First parameter is integral (number events under gaus)
592   // Forth parameter (par[3]) is width of histogram bin 
593   // gaus(0) is a substitute for [0]*exp(-0.5*((x-[1])/[2])**2)
594   // + pol2 -> 7 parameters
595   static Double_t c = TMath::Sqrt(2.*TMath::Pi()), y=0.0, f=0.0; // sqrt(2.*pi)
596
597   y  = (x[0]-par[1])/par[2];
598   f  = par[0]*par[3]/(par[2]*c) * TMath::Exp(-0.5*y*y);
599
600   f += par[4] + par[5]*x[0] + par[6]*x[0]*x[0];
601
602   return f;
603 }
604
605 // Calibration stuff
606 Double_t AliEMCALHistoUtilities::GetCorrectionCoefficientForGamma1(const Double_t eRec)
607 {
608   // Correction to rec.energy - Jul 15, 2007
609   // E(gamma) = corr * E(eRec);
610   /* corr = p0*(eRec-7.5)+p1*(eRec-7.5)*(eRec-7.5)+p2; 0.0<eRec<10.0
611    1  p0           6.07157e-05   1.15179e-04  -0.00000e+00   1.20997e-03
612    2  p1           1.50019e-04   3.13566e-05  -0.00000e+00   1.22531e-02
613    3  p2           9.99019e-01   4.08251e-04  -0.00000e+00   1.40606e-03
614
615      corr = p3 + p4*eRec + p5*eRec*eRec
616    1  p3           9.97135e-01   5.31970e-04   1.37962e-09   1.68120e-08
617    2  p4           3.15740e-04   2.53371e-05   1.11475e-11   1.74192e-04
618    3  p5          -1.35383e-06   2.19495e-07  -5.82864e-13   4.52277e-02
619   */
620   static Double_t p0=6.07157e-05, p1=1.50019e-04, p2=9.99019e-01;
621   static Double_t p3=9.97135e-01, p4=3.15740e-04, p5=-1.35383e-06;
622   static Double_t corr=1.0;
623   if(eRec>=0.0 && eRec <=10.0) {
624     corr = p0*(eRec-7.5) + p1*(eRec-7.5)*(eRec-7.5) + p2;
625   } else if(eRec>10.0 && eRec <=105.0) {
626     corr = p3 + p4*eRec + p5*eRec*eRec;
627   }
628   //printf(" eRec %f | corr %f \n", eRec, corr);
629   return corr;
630 }
631
632 Double_t AliEMCALHistoUtilities::GetCorrectedEnergyForGamma1(const Double_t eRec)
633 {
634   return GetCorrectionCoefficientForGamma1(eRec) * eRec;
635 }
636
637 // Trigger 
638 TList* AliEMCALHistoUtilities::GetTriggersListOfHists(const Int_t scale, const Int_t nTrig, const Bool_t toBrowser)
639 {
640   // Oct 22, 2007 - trigger technical assurance
641   gROOT->cd();
642   TH1::AddDirectory(1);
643
644   new TH1F("00_hXpos2x2", "X coord. of max Amp 2x2",100, -500., +500.);
645   new TH1F("01_hYpos2x2", "Y coord. of max Amp 2x2",100, -500., +500.);
646   new TH1F("02_hZpos2x2", "Z coord. of max Amp 2x2",100, -500., +500.);
647   new TH1F("03_hXposnxn", "X coord. of max Amp NXN",100, -500., +500.);
648   new TH1F("04_hYposnxn", "Y coord. of max Amp NXN",100, -500., +500.);
649   new TH1F("05_hZposnxn", "Z coord. of max Amp NXN",100, -500., +500.);
650   // May 7, 2008 - jet trigger
651   new TH1F("06_hJetTriggerPhi", "%phi of COG of jet trigger patch", 110, 80., 190.);
652   new TH1F("07_hJetTriggerEta", "%eta of COG of jet trigger patch", 70, -0.7, +0.7);
653   // 
654   new TH1F("08_hMaxAmp2x2", "max Amp 2x2", 1000, 0.0, pow(2.,14.));
655   new TH1F("09_hAmpOutOf2x2", "Amp out of patch 2x2", 1000, 0.0, pow(2.,14.));
656   new TH1F("10_hMaxAmpnxn", "max Amp NXN", 1000, 0.0, pow(2.,14.));
657   new TH1F("11_hAmpOutOfnxn", "Amp out of patch nxn", 1000, 0.0, pow(2.,14.));
658   // May 7, 2008 - jet trigger
659   for(Int_t i=0; i<nTrig; i++) {
660     new TH1F(Form("%2.2i_hJetPatchAmp%2.2i",i+12,i), Form("jet patch amplitude : jet trig %i",i), 
661     1000, 0.0, pow(2.,14.));
662   }
663   // For checking
664   Double_t maxEdigit=100., maxN=1000., maxPC = 200.;
665   if(scale==1) {
666     maxN  *= 10.;
667     maxPC *= 10.;
668   }
669   Int_t ind = 12+nTrig; 
670   new TH1F(Form("%2.2i_hDigitsAmp",ind++), "amplitude of digits (PC)  ", 1001, -0.5, 1000.5);
671   new TH1F(Form("%2.2i_hDigitsE",ind++), " energy of digits (PC)", 1000, 0.0, maxEdigit);
672   new TH1F(Form("%2.2i_hNDigitsInPCs",ind++), " number of digits in PC's", 1000, 0.0, maxN);
673   new TH1F(Form("%2.2i_hEinInPCs",ind++), " energy in PC's", 200, 0.0, maxPC);
674
675   return MoveHistsToList("TriggerLiOfHists", toBrowser);
676 }
677
678 void AliEMCALHistoUtilities::FillTriggersListOfHists(TList *l, TArrayF *triggerPosition, TArrayF *triggerAmplitudes)
679 {
680   // Oct 22, 2007 - trigger technical assurance
681   if(l==0) {
682     printf("<E> FillTriggersListOfHists() : list of hists undefined. \n");
683     return;
684   }
685   for(int i=0; i<triggerPosition->GetSize(); i++) {
686     FillH1(l, i, double(triggerPosition->At(i)));
687   }
688
689   for(int i=0; i<triggerAmplitudes->GetSize(); i++) {
690     FillH1(l, triggerPosition->GetSize() + i, double(triggerAmplitudes->At(i)) );
691   }
692 }
693
694 // Jet(s) kinematics 
695 TList* AliEMCALHistoUtilities::GetJetsListOfHists(Int_t njet, Bool_t toBrowser)
696 {
697   // Oct 30, 2007
698   gROOT->cd();
699   TH1::AddDirectory(1);
700   new TH1F("00_nJet", "number of jets in Pythia",5, -0.5, +4.5);
701   int ic=1;
702   for(Int_t ij=0; ij<njet; ij++)
703   {
704     new TH1F(Form("%2.2i_EtaOfJet%i",ic++,ij), Form("#eta of jet (%i)",ij),80, -2., 2.);
705     new TH1F(Form("%2.2i_ThetaOfJet%i",ic++,ij), Form("#theta(degree) of jet (%i)",ij),
706     180, 0.0, 180.);
707     new TH1F(Form("%2.2i_PhiOfJet%i",ic++,ij), Form("#phi(degree) of jet (%i)",ij),
708     180, 0.0, 360.);
709     new TH1F(Form("%2.2i_PtOfJet%i",ic++,ij), Form(" p_{T} of jet (%i)",ij),
710     200, 0.0, 200.);
711     new TH1F(Form("%2.2i_EOfJet%i",ic++,ij), Form(" E of jet (%i)",ij),
712     200, 0.0, 200.);
713   }
714
715   return MoveHistsToList("JetLiOfHists", toBrowser);
716 }
717
718 void  AliEMCALHistoUtilities::FillJetKineListOfHists(TList *l, AliRunLoader* rl, TLorentzVector &goodJet)
719 {
720   // Oct 30, 2007; Nov 07;
721   // goodJet - output; only one jet with EMCAL acceptance
722
723   // Bad case
724   goodJet.SetPxPyPzE(0., 0., 0., 0.); 
725
726   if(l==0 || rl==0) return;
727   //try to get trigger jet info
728   AliHeader* alih = rl->GetHeader();
729   if (alih == 0) return;
730
731   AliGenEventHeader * genh = alih->GenEventHeader();
732   if (genh == 0) return;
733
734   AliGenPythiaEventHeader* genhpy = dynamic_cast<AliGenPythiaEventHeader *>(genh);
735   if(genhpy==0) return; // Not Pythia
736
737   Int_t nj = genhpy->NTriggerJets();
738   FillH1(l, 0, double(nj));
739
740   TLorentzVector* jets[4];
741   nj = nj>4?4:nj;
742
743   int ic=1;
744   for (Int_t i=0; i< nj; i++) {
745      Float_t p[4];
746      genhpy->TriggerJet(i,p);
747      jets[i] = new TLorentzVector(p);
748      FillH1(l, ic++, jets[i]->Eta());
749      FillH1(l, ic++, jets[i]->Theta()*TMath::RadToDeg());
750      FillH1(l, ic++, TVector2::Phi_0_2pi(jets[i]->Phi())*TMath::RadToDeg());
751      FillH1(l, ic++, jets[i]->Pt());
752      FillH1(l, ic++, jets[i]->E());
753
754      //printf(" %i pj %f %f %f %f : ic %i\n", i, p[0], p[1], p[2], p[3], ic);
755      if(ic >= l->GetSize()) break;
756   }
757   if(nj==1) {
758     Double_t eta=jets[0]->Eta(), phi=TVector2::Phi_0_2pi(jets[0]->Phi())*TMath::RadToDeg();
759     if(TMath::Abs(eta)<0.5 && (phi>90.&&phi<180.)) {
760       goodJet = (*jets[0]);
761     }
762   }
763 }