]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/dielectron/AliDielectronHistos.cxx
small fix to remove TF1Helper which sneaked in again, inadvertently...
[u/mrichter/AliRoot.git] / PWG3 / dielectron / AliDielectronHistos.cxx
1 /*************************************************************************
2 * Copyright(c) 1998-2009, 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 //
17 // Generic Histogram container with support for groups and filling of groups by passing
18 // a vector of data
19 //
20 // Authors: 
21 //   Jens Wiechula <Jens.Wiechula@cern.ch> 
22 // 
23
24 #include <TH1.h>
25 #include <TH1F.h>
26 #include <TH2.h>
27 #include <TH3.h>
28 #include <TCollection.h>
29 #include <THashList.h>
30 #include <TString.h>
31 #include <TObjString.h>
32 #include <TObjArray.h>
33 #include <TFile.h>
34 #include <TError.h>
35 #include <TCanvas.h>
36 #include <TMath.h>
37 #include <TROOT.h>
38 #include <TLegend.h>
39 #include <TKey.h>
40 #include <TAxis.h>
41 #include <TVirtualPS.h>
42 // #include <TVectorD.h>
43
44 #include "AliDielectronHistos.h"
45
46
47 ClassImp(AliDielectronHistos)
48
49
50 AliDielectronHistos::AliDielectronHistos() :
51 //   TCollection(),
52   TNamed("AliDielectronHistos","Dielectron Histogram Container"),
53   fHistoList(),
54   fList(0x0),
55   fReservedWords(new TString)
56 {
57   //
58   // Default constructor
59   //
60   fHistoList.SetOwner(kTRUE);
61   fHistoList.SetName("Dielectron_Histos");
62 }
63
64 //_____________________________________________________________________________
65 AliDielectronHistos::AliDielectronHistos(const char* name, const char* title) :
66 //   TCollection(),
67   TNamed(name, title),
68   fHistoList(),
69   fList(0x0),
70   fReservedWords(new TString)
71 {
72   //
73   // TNamed constructor
74   //
75   fHistoList.SetOwner(kTRUE);
76   fHistoList.SetName(name);
77 }
78
79 //_____________________________________________________________________________
80 AliDielectronHistos::~AliDielectronHistos()
81 {
82   //
83   // Destructor
84   //
85   fHistoList.Clear();
86   if (fList) fList->Clear();
87   delete fReservedWords;
88 }
89
90 //_____________________________________________________________________________
91 void AliDielectronHistos::UserHistogram(const char* histClass,const char *name, const char* title,
92                                         Int_t nbinsX, Double_t xmin, Double_t xmax,
93                                         UInt_t valTypeX, Bool_t logBinX)
94 {
95   //
96   // Default histogram creation 1D case
97   //
98   if (!IsHistogramOk(histClass,name)) return;
99
100   Double_t *binLimX=0x0;
101   
102   if (logBinX) {
103     binLimX=MakeLogBinning(nbinsX, xmin, xmax);
104   } else {
105     binLimX=MakeLinBinning(nbinsX, xmin, xmax);
106   }
107
108   TH1* hist=new TH1F(name,title,nbinsX,binLimX);
109
110   delete [] binLimX;
111   
112   Bool_t isReserved=fReservedWords->Contains(histClass);
113   if (isReserved)
114     UserHistogramReservedWords(histClass, hist, valTypeX);
115   else
116     UserHistogram(histClass, hist, valTypeX);
117 }
118
119 //_____________________________________________________________________________
120 void AliDielectronHistos::UserHistogram(const char* histClass,const char *name, const char* title,
121                                         Int_t nbinsX, Double_t xmin, Double_t xmax,
122                                         Int_t nbinsY, Double_t ymin, Double_t ymax,
123                                         UInt_t valTypeX, UInt_t valTypeY,
124                                         Bool_t logBinX, Bool_t logBinY)
125 {
126   //
127   // Default histogram creation 2D case
128   //
129   if (!IsHistogramOk(histClass,name)) return;
130
131   Double_t *binLimX=0x0;
132   Double_t *binLimY=0x0;
133   
134   if (logBinX) {
135     binLimX=MakeLogBinning(nbinsX, xmin, xmax);
136   } else {
137     binLimX=MakeLinBinning(nbinsX, xmin, xmax);
138   }
139   if (logBinY) {
140     binLimY=MakeLogBinning(nbinsY, ymin, ymax);
141   } else {
142     binLimY=MakeLinBinning(nbinsY, ymin, ymax);
143   }
144   
145   TH1* hist=new TH2F(name,title,nbinsX,binLimX,nbinsY,binLimY);
146   
147   delete [] binLimX;
148   delete [] binLimY;
149   
150   
151   Bool_t isReserved=fReservedWords->Contains(histClass);
152   if (isReserved)
153     UserHistogramReservedWords(histClass, hist, valTypeX+100*valTypeY);
154   else
155     UserHistogram(histClass, hist, valTypeX+100*valTypeY);
156 }
157
158
159 //_____________________________________________________________________________
160 void AliDielectronHistos::UserHistogram(const char* histClass,const char *name, const char* title,
161                                         Int_t nbinsX, Double_t xmin, Double_t xmax,
162                                         Int_t nbinsY, Double_t ymin, Double_t ymax,
163                                         Int_t nbinsZ, Double_t zmin, Double_t zmax,
164                                         UInt_t valTypeX, UInt_t valTypeY, UInt_t valTypeZ,
165                                         Bool_t logBinX, Bool_t logBinY, Bool_t logBinZ)
166 {
167   //
168   // Default histogram creation 3D case
169   //
170   if (!IsHistogramOk(histClass,name)) return;
171
172   Double_t *binLimX=0x0;
173   Double_t *binLimY=0x0;
174   Double_t *binLimZ=0x0;
175   
176   if (logBinX) {
177     binLimX=MakeLogBinning(nbinsX, xmin, xmax);
178   } else {
179     binLimX=MakeLinBinning(nbinsX, xmin, xmax);
180   }
181   
182   if (logBinY) {
183     binLimY=MakeLogBinning(nbinsY, ymin, ymax);
184   } else {
185     binLimY=MakeLinBinning(nbinsY, ymin, ymax);
186   }
187   
188   if (logBinZ) {
189     binLimZ=MakeLogBinning(nbinsZ, zmin, zmax);
190   } else {
191     binLimZ=MakeLinBinning(nbinsZ, zmin, zmax);
192   }
193   
194   TH1* hist=new TH3F(name,title,nbinsX,binLimX,nbinsY,binLimY,nbinsZ,binLimZ);
195   
196   delete [] binLimX;
197   delete [] binLimY;
198   delete [] binLimZ;
199   
200   Bool_t isReserved=fReservedWords->Contains(histClass);
201   if (isReserved)
202     UserHistogramReservedWords(histClass, hist, valTypeX+100*valTypeY+10000*valTypeZ);
203   else
204     UserHistogram(histClass, hist, valTypeX+100*valTypeY+10000*valTypeZ);
205 }
206
207 //_____________________________________________________________________________
208 void AliDielectronHistos::UserHistogram(const char* histClass, TH1* hist, UInt_t valTypes)
209 {
210   //
211   // Add any type of user histogram
212   //
213
214   //special case for the calss Pair. where histograms will be created for all pair classes
215   Bool_t isReserved=fReservedWords->Contains(histClass);
216   if (isReserved) {
217     UserHistogramReservedWords(histClass, hist, valTypes);
218     return;
219   }
220   
221   if (!IsHistogramOk(histClass,hist->GetName())) return;
222   
223   THashList *classTable=(THashList*)fHistoList.FindObject(histClass);
224   hist->SetDirectory(0);
225   hist->SetUniqueID(valTypes);
226   classTable->Add(hist);
227 }
228
229 //_____________________________________________________________________________
230 void AliDielectronHistos::AddClass(const char* histClass)
231 {
232   //
233   // Add a class of histograms
234   // Several classes can be added by separating them by a ';' e.g. 'class1;class2;class3'
235   //
236   TString hists(histClass);
237   TObjArray *arr=hists.Tokenize(";");
238   TIter next(arr);
239   TObject *o=0;
240   while ( (o=next()) ){
241     if (fHistoList.FindObject(o->GetName())){
242       Warning("AddClass","Cannot create class '%s' it already exists.",histClass);
243       continue;
244     }
245     if (fReservedWords->Contains(o->GetName())){
246       Error("AddClass","Pair is a reserved word, please use another name");
247       continue;
248     }
249     THashList *table=new THashList;
250     table->SetOwner(kTRUE);
251     table->SetName(o->GetName());
252     fHistoList.Add(table);
253   }
254   delete arr;
255 }
256
257 //_____________________________________________________________________________
258 void AliDielectronHistos::Fill(const char* histClass, const char* name, Double_t xval)
259 {
260   //
261   // Fill function 1D case
262   //
263   THashList *classTable=(THashList*)fHistoList.FindObject(histClass);
264   TH1* hist=0;
265   if (!classTable || !(hist=(TH1*)classTable->FindObject(name)) ){
266     Warning("Fill","Cannot fill histogram. Either class '%s' or histogram '%s' not existing.",histClass,name);
267     return;
268   }
269   hist->Fill(xval);
270 }
271
272 //_____________________________________________________________________________
273 void AliDielectronHistos::Fill(const char* histClass, const char* name, Double_t xval, Double_t yval)
274 {
275   //
276   // Fill function 2D case
277   //
278   THashList *classTable=(THashList*)fHistoList.FindObject(histClass);
279   TH2* hist=0;
280   if (!classTable || !(hist=(TH2*)classTable->FindObject(name)) ){
281     Warning("UserHistogram","Cannot fill histogram. Either class '%s' or histogram '%s' not existing.",histClass,name);
282     return;
283   }
284   hist->Fill(xval,yval);
285 }
286
287 //_____________________________________________________________________________
288 void AliDielectronHistos::Fill(const char* histClass, const char* name, Double_t xval, Double_t yval, Double_t zval)
289 {
290   //
291   // Fill function 3D case
292   //
293   THashList *classTable=(THashList*)fHistoList.FindObject(histClass);
294   TH3* hist=0;
295   if (!classTable || !(hist=(TH3*)classTable->FindObject(name)) ){
296     Warning("UserHistogram","Cannot fill histogram. Either class '%s' or histogram '%s' not existing.",histClass,name);
297     return;
298   }
299   hist->Fill(xval,yval,zval);
300 }
301
302 //_____________________________________________________________________________
303 void AliDielectronHistos::FillClass(const char* histClass, Int_t nValues, const Double_t *values)
304 {
305   //
306   // Fill class 'histClass' (by name)
307   //
308   
309   THashList *classTable=(THashList*)fHistoList.FindObject(histClass);
310   if (!classTable){
311     Warning("FillClass","Cannot fill class '%s' its not defined.",histClass);
312     return;
313   }
314   
315   TIter nextHist(classTable);
316   TH1 *hist=0;
317   while ( (hist=(TH1*)nextHist()) ){
318     UInt_t valueTypes=hist->GetUniqueID();
319     if (valueTypes==(UInt_t)kNoAutoFill) continue;
320     UInt_t value1=valueTypes%100;        //last two digits
321     UInt_t value2=valueTypes/100%100;    //second last two digits
322     UInt_t value3=valueTypes/10000%100;  //third last two digits
323     if (value1>=(UInt_t)nValues||value2>=(UInt_t)nValues||value3>=(UInt_t)nValues) {
324       Warning("FillClass","One of the values is out of range. Not filling histogram '%s/%s'.", histClass, hist->GetName());
325       continue;
326     }
327     switch (hist->GetDimension()){
328     case 1:
329       hist->Fill(values[value1]);
330       break;
331     case 2:
332       ((TH2*)hist)->Fill(values[value1],values[value2]);
333       break;
334     case 3:
335       ((TH3*)hist)->Fill(values[value1],values[value2],values[value3]);
336       break;
337     }
338   }
339 }
340
341 //_____________________________________________________________________________
342 // void AliDielectronHistos::FillClass(const char* histClass, const TVectorD &vals)
343 // {
344 //   //
345 //   //
346 //   //
347 //   FillClass(histClass, vals.GetNrows(), vals.GetMatrixArray());
348 // }
349
350 //_____________________________________________________________________________
351 void AliDielectronHistos::UserHistogramReservedWords(const char* histClass, TH1 *hist, UInt_t valTypes)
352 {
353   //
354   // Creation of histogram for all pair types
355   //
356   TString title(hist->GetTitle());
357   // Same Event Like Sign
358   TIter nextClass(&fHistoList);
359   THashList *l=0;
360   while ( (l=static_cast<THashList*>(nextClass())) ){
361     TString name(l->GetName());
362     if (name.Contains(histClass)){
363       TH1 *h=static_cast<TH1*>(hist->Clone());
364       h->SetDirectory(0);
365       h->SetTitle(Form("%s %s",title.Data(),l->GetName()));
366       UserHistogram(l->GetName(),h,valTypes);
367     }
368   }
369   delete hist;
370 }
371
372 //_____________________________________________________________________________
373 void AliDielectronHistos::DumpToFile(const char* file)
374 {
375   //
376   // Dump the histogram list to a newly created root file
377   //
378   TFile f(file,"recreate");
379   fHistoList.Write(fHistoList.GetName(),TObject::kSingleKey);
380   f.Close();
381 }
382
383 //_____________________________________________________________________________
384 TH1* AliDielectronHistos::GetHistogram(const char* histClass, const char* name) const
385 {
386   //
387   // return histogram 'name' in 'histClass'
388   //
389   THashList *classTable=(THashList*)fHistoList.FindObject(histClass);
390   if (!classTable) return 0x0;
391   return (TH1*)classTable->FindObject(name);
392 }
393
394 //_____________________________________________________________________________
395 TH1* AliDielectronHistos::GetHistogram(const char* cutClass, const char* histClass, const char* name) const
396 {
397   //
398   // return histogram from list of list of histograms
399   // this function is thought for retrieving histograms if a list of AliDielectronHistos is set
400   //
401   
402   if (!fList) return 0x0;
403   THashList *h=dynamic_cast<THashList*>(fList->FindObject(cutClass));
404   if (!h)return 0x0;
405   THashList *classTable=dynamic_cast<THashList*>(h->FindObject(histClass));
406   if (!classTable) return 0x0;
407   return (TH1*)classTable->FindObject(name);
408 }
409
410 //_____________________________________________________________________________
411 void AliDielectronHistos::Draw(const Option_t* option)
412 {
413   //
414   // Draw histograms
415   //
416
417   TString drawStr(option);
418   TObjArray *arr=drawStr.Tokenize(";");
419   arr->SetOwner();
420   TIter nextOpt(arr);
421
422   TString drawClasses;
423   TObjString *ostr=0x0;
424
425   TString currentOpt;
426   TString testOpt;
427   while ( (ostr=(TObjString*)nextOpt()) ){
428     currentOpt=ostr->GetString();
429     currentOpt.Remove(TString::kBoth,'\t');
430     currentOpt.Remove(TString::kBoth,' ');
431
432     testOpt="classes=";
433     if ( currentOpt.Contains(testOpt.Data()) ){
434       drawClasses=currentOpt(testOpt.Length(),currentOpt.Length());
435     }
436   }
437
438   delete arr;
439   drawStr.ToLower();
440   //optionsfList
441 //   Bool_t same=drawOpt.Contains("same"); //FIXME not yet implemented
442
443   TCanvas *c=0x0;
444   if (gVirtualPS) {
445     if (!gPad){
446       Error("Draw","When writing to a file you have to create a canvas before opening the file!!!");
447       return;
448     }
449     c=gPad->GetCanvas();
450     c->cd();
451 //     c=new TCanvas;
452   }
453   
454   TIter nextClass(&fHistoList);
455   THashList *classTable=0;
456 //   Bool_t first=kTRUE;
457   while ( (classTable=(THashList*)nextClass()) ){
458     //test classes option
459     if (!drawClasses.IsNull() && !drawClasses.Contains(classTable->GetName())) continue;
460     //optimised division
461     Int_t nPads = classTable->GetEntries();
462     Int_t nCols = (Int_t)TMath::Ceil( TMath::Sqrt(nPads) );
463     Int_t nRows = (Int_t)TMath::Ceil( (Double_t)nPads/(Double_t)nCols );
464
465     //create canvas
466     if (!gVirtualPS){
467       TString canvasName;
468       canvasName.Form("c%s_%s",GetName(),classTable->GetName());
469       c=(TCanvas*)gROOT->FindObject(canvasName.Data());
470       if (!c) c=new TCanvas(canvasName.Data(),Form("%s: %s",GetName(),classTable->GetName()));
471       c->Clear();
472     } else {
473 //       if (first){
474 //         first=kFALSE;
475 //         if (nPads>1) gVirtualPS->NewPage();
476 //       } else {
477         if (nPads>1) c->Clear();
478 //       }
479     }
480     if (nCols>1||nRows>1) c->Divide(nCols,nRows);
481     
482     //loop over histograms and draw them
483     TIter nextHist(classTable);
484     Int_t iPad=0;
485     TH1 *h=0;
486     while ( (h=(TH1*)nextHist()) ){
487       TString drawOpt;
488       if ( (h->InheritsFrom(TH2::Class())) ) drawOpt="colz";
489       if (nCols>1||nRows>1) c->cd(++iPad);
490       if ( TMath::Abs(h->GetXaxis()->GetBinWidth(1)-h->GetXaxis()->GetBinWidth(2))>1e-10 ) gPad->SetLogx();
491       if ( TMath::Abs(h->GetYaxis()->GetBinWidth(1)-h->GetYaxis()->GetBinWidth(2))>1e-10 ) gPad->SetLogy();
492       if ( TMath::Abs(h->GetZaxis()->GetBinWidth(1)-h->GetZaxis()->GetBinWidth(2))>1e-10 ) gPad->SetLogz();
493       TString histOpt=h->GetOption();
494       histOpt.ToLower();
495       if (histOpt.Contains("logx")) gPad->SetLogx();
496       if (histOpt.Contains("logy")) gPad->SetLogy();
497       if (histOpt.Contains("logz")) gPad->SetLogz();
498       histOpt.ReplaceAll("logx","");
499       histOpt.ReplaceAll("logy","");
500       histOpt.ReplaceAll("logz","");
501       h->Draw(drawOpt.Data());
502     }
503     if (gVirtualPS) {
504       c->Update();
505     }
506     
507   }
508 //   if (gVirtualPS) delete c;
509 }
510
511 //_____________________________________________________________________________
512 void AliDielectronHistos::Print(const Option_t* option) const
513 {
514   //
515   // Print classes and histograms
516   //
517   TString optString(option);
518
519   if (optString.IsNull()) PrintStructure();
520
521
522
523 }
524
525 //_____________________________________________________________________________
526 void AliDielectronHistos::PrintStructure() const
527 {
528   //
529   // Print classes and histograms in the class to stdout
530   //
531   if (!fList){
532     TIter nextClass(&fHistoList);
533     THashList *classTable=0;
534     while ( (classTable=(THashList*)nextClass()) ){
535       TIter nextHist(classTable);
536       TObject *o=0;
537       printf("+ %s\n",classTable->GetName());
538       while ( (o=nextHist()) )
539         printf("| ->%s\n",o->GetName());
540     }
541   } else {
542     TIter nextCutClass(fList);
543     THashList *cutClass=0x0;
544     while ( (cutClass=(THashList*)nextCutClass()) ) {
545       printf("+ %s\n",cutClass->GetName());
546       TIter nextClass(cutClass);
547       THashList *classTable=0;
548       while ( (classTable=(THashList*)nextClass()) ){
549         TIter nextHist(classTable);
550         TObject *o=0;
551         printf("|  + %s\n",classTable->GetName());
552         while ( (o=nextHist()) )
553           printf("|  | ->%s\n",o->GetName());
554       }
555       
556     }
557   }
558 }
559
560 //_____________________________________________________________________________
561 void AliDielectronHistos::SetHistogramList(THashList &list, Bool_t setOwner/*=kTRUE*/)
562 {
563   //
564   // set histogram classes and histograms to this instance. It will take onwnership!
565   //
566   ResetHistogramList();
567   TString name(GetName());
568   if (name == "AliDielectronHistos") SetName(list.GetName());
569   TIter next(&list);
570   TObject *o;
571   while ( (o=next()) ){
572     fHistoList.Add(o);
573   }
574   if (setOwner){
575     list.SetOwner(kFALSE);
576     fHistoList.SetOwner(kTRUE);
577   }
578 }
579
580 //_____________________________________________________________________________
581 Bool_t AliDielectronHistos::SetCutClass(const char* cutClass)
582 {
583   //
584   // Assign histogram list according to cutClass
585   //
586
587   if (!fList) return kFALSE;
588   ResetHistogramList();
589   THashList *h=dynamic_cast<THashList*>(fList->FindObject(cutClass));
590   if (!h) {
591     Warning("SetCutClass","cutClass '%s' not found", cutClass);
592     return kFALSE;
593   }
594   SetHistogramList(*h,kFALSE);
595   return kTRUE;
596 }
597
598 //_____________________________________________________________________________
599 Bool_t AliDielectronHistos::IsHistogramOk(const char* histClass, const char* name)
600 {
601   //
602   // check whether the histogram class exists and the histogram itself does not exist yet
603   //
604   Bool_t isReserved=fReservedWords->Contains(histClass);
605   if (!fHistoList.FindObject(histClass)&&!isReserved){
606     Warning("IsHistogramOk","Cannot create histogram. Class '%s' not defined. Please create it using AddClass before.",histClass);
607     return kFALSE;
608   }
609   if (GetHistogram(histClass,name)){
610     Warning("IsHistogramOk","Cannot create histogram '%s' in class '%s': It already exists!",name,histClass);
611     return kFALSE;
612   }
613   return kTRUE;
614 }
615
616 // //_____________________________________________________________________________
617 // TIterator* AliDielectronHistos::MakeIterator(Bool_t dir) const
618 // {
619 //   //
620 //   //
621 //   //
622 //   return new TListIter(&fHistoList, dir);
623 // }
624
625 //_____________________________________________________________________________
626 void AliDielectronHistos::ReadFromFile(const char* file)
627 {
628   //
629   // Read histos from file
630   //
631   TFile f(file);
632   TIter nextKey(f.GetListOfKeys());
633   TKey *key=0;
634   while ( (key=(TKey*)nextKey()) ){
635     TObject *o=f.Get(key->GetName());
636     THashList *list=dynamic_cast<THashList*>(o);
637     if (!list) continue;
638     SetHistogramList(*list);
639     break;
640   }
641   f.Close();
642 }
643
644 //_____________________________________________________________________________
645 void AliDielectronHistos::DrawSame(const char* histName, const Option_t *opt)
646 {
647   //
648   // Draw all histograms with the same name into one canvas
649   // if option contains 'leg' a legend will be created with the class name as caption
650   // if option contains 'can' a new canvas is created
651   //
652
653   TString optString(opt);
654   optString.ToLower();
655   Bool_t optLeg=optString.Contains("leg");
656   Bool_t optCan=optString.Contains("can");
657
658   TLegend *leg=0;
659   TCanvas *c=0;
660   if (optCan){
661     c=(TCanvas*)gROOT->FindObject(Form("c%s",histName));
662     if (!c) c=new TCanvas(Form("c%s",histName),Form("All '%s' histograms",histName));
663     c->Clear();
664     c->cd();
665   }
666
667   if (optLeg) leg=new TLegend(.8,.3,.99,.9);
668   
669   Int_t i=0;
670   TIter next(&fHistoList);
671   THashList *classTable=0;
672   Double_t max=-1e10;
673   TH1 *hFirst=0x0;
674   while ( (classTable=(THashList*)next()) ){
675     if ( TH1 *h=(TH1*)classTable->FindObject(histName) ){
676       if (i==0) hFirst=h;
677       h->SetLineColor(i+1);
678       h->SetMarkerColor(i+1);
679       h->Draw(i>0?"same":"");
680       if (leg) leg->AddEntry(h,classTable->GetName(),"lp");
681       ++i;
682       max=TMath::Max(max,h->GetMaximum());
683     }
684   }
685   if (leg){
686     leg->SetFillColor(10);
687     leg->SetY1(.9-i*.05);
688     leg->Draw();
689   }
690   if (hFirst&&(hFirst->GetYaxis()->GetXmax()<max)){
691     hFirst->SetMaximum(max);
692   }
693 }
694
695 //_____________________________________________________________________________
696 void AliDielectronHistos::SetReservedWords(const char* words)
697 {
698   //
699   // set reserved words
700   //
701   
702   (*fReservedWords)=words;
703 }
704
705 //_____________________________________________________________________________
706 Double_t* AliDielectronHistos::MakeLogBinning(Int_t nbinsX, Double_t xmin, Double_t xmax) const
707 {
708   //
709   // Make logarithmic binning
710   // the user has to delete the array afterwards!!!
711   //
712
713   //check limits
714   if (xmin<1e-20 || xmax<1e-20){
715     Error("MakeLogBinning","For Log binning xmin and xmax must be > 1e-20. Using linear binning instead!");
716     return MakeLinBinning(nbinsX, xmin, xmax);
717   }
718   Double_t *binLim=new Double_t[nbinsX+1];
719   Double_t first=xmin;
720   Double_t last=xmax;
721   Double_t expMax=TMath::Log(last/first);
722   for (Int_t i=0; i<nbinsX+1; ++i){
723     binLim[i]=first*TMath::Exp(expMax/nbinsX*(Double_t)i);
724   }
725   return binLim;
726 }
727
728 //_____________________________________________________________________________
729 Double_t* AliDielectronHistos::MakeLinBinning(Int_t nbinsX, Double_t xmin, Double_t xmax) const
730 {
731   //
732   // Make logarithmic binning
733   // the user has to delete the array afterwards!!!
734   //
735   Double_t *binLim=new Double_t[nbinsX+1];
736   Double_t first=xmin;
737   Double_t last=xmax;
738   Double_t binWidth=(last-first)/nbinsX;
739   for (Int_t i=0; i<nbinsX+1; ++i){
740     binLim[i]=first+binWidth*(Double_t)i;
741   }
742   return binLim;
743 }
744