]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/dielectron/AliDielectronCFdraw.cxx
framework update; new classes for track rotation (for background), cuts grouping...
[u/mrichter/AliRoot.git] / PWG3 / dielectron / AliDielectronCFdraw.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 //       Dielectron Correction framework draw helper                     //
18 //                                                                       //
19 /*
20
21
22
23
24
25
26
27
28
29 */
30 //                                                                       //
31 ///////////////////////////////////////////////////////////////////////////
32
33 #include <TSeqCollection.h>
34 #include <TObjArray.h>
35 #include <TKey.h>
36 #include <TList.h>
37 #include <TClass.h>
38 #include <TObject.h>
39 #include <TVirtualPS.h>
40 #include <TFile.h>
41 #include <TString.h>
42 #include <TObjString.h>
43 #include <TVectorD.h>
44 #include <TMath.h>
45 #include <TH1.h>
46 #include <TH2.h>
47 #include <TH3.h>
48 #include <TPad.h>
49 #include <TCanvas.h>
50 #include <TLegend.h>
51 #include <AliCFEffGrid.h>
52
53 #include <AliLog.h>
54
55 #include "AliDielectronCFdraw.h"
56
57 ClassImp(AliDielectronCFdraw)
58
59 AliDielectronCFdraw::AliDielectronCFdraw() :
60   TNamed(),
61   fCfContainer(0x0),
62   fEffGrid(0x0),
63   fVdata(0)
64 {
65   //
66   // Ctor
67   //
68 }
69
70 //________________________________________________________________
71 AliDielectronCFdraw::AliDielectronCFdraw(const char*name, const char* title) :
72   TNamed(name,title),
73   fCfContainer(0x0),
74   fEffGrid(0x0),
75   fVdata(0)
76 {
77   //
78   // Named Ctor
79   //
80   
81 }
82
83 //________________________________________________________________
84 AliDielectronCFdraw::AliDielectronCFdraw(AliCFContainer *cont) :
85   TNamed(cont->GetName(), cont->GetTitle()),
86   fCfContainer(cont),
87   fEffGrid(new AliCFEffGrid("eff","eff",*cont)),
88   fVdata(0)
89 {
90   //
91   // directly set the CF container
92   //
93
94 }
95
96 //________________________________________________________________
97 AliDielectronCFdraw::AliDielectronCFdraw(const char* filename) :
98   TNamed(),
99   fCfContainer(0x0),
100   fEffGrid(0x0),
101   fVdata(0)
102 {
103   //
104   // get CF container(s) from file 'filename'
105   //
106   SetCFContainers(filename);
107 }
108
109 //________________________________________________________________
110 void AliDielectronCFdraw::SetCFContainers(const TSeqCollection *arr)
111 {
112   //
113   // Merge CF Container out of several containers
114   //
115
116   TIter next(arr);
117   TObject *o=0x0;
118
119   Int_t nstep=0;
120   while ( (o=next()) ){
121     AliCFContainer *cf=dynamic_cast<AliCFContainer*>(o);
122     if (!cf) continue;
123     nstep+=cf->GetNStep();
124   }
125   if (nstep==0) return;
126   Int_t nbins[1]={1};
127   fCfContainer=new AliCFContainer(GetName(), GetTitle(), nstep, 1, nbins);
128
129   //delete unneeded steps
130 //   for (Int_t istep=0; istep<nstep; ++istep) delete fCfContainer->GetGrid(istep);
131
132   //add step to the new container
133   Int_t istep=0;
134   for (Int_t icf=0; icf<arr->GetEntries(); ++icf){
135     AliCFContainer *cf=dynamic_cast<AliCFContainer*>(arr->At(icf));
136     if (!cf) continue;
137     for (Int_t istepCurr=0; istepCurr<cf->GetNStep(); ++istepCurr){
138       fCfContainer->SetGrid(istep, cf->GetGrid(istepCurr));
139       fCfContainer->SetStepTitle(istep,Form("%s, Pair: %s",cf->GetTitle(),cf->GetStepTitle(istepCurr)));
140       ++istep;
141     }
142   }
143   if (fEffGrid) delete fEffGrid;
144   fEffGrid=new AliCFEffGrid("eff","eff",*fCfContainer);
145 }
146
147 //________________________________________________________________
148 void AliDielectronCFdraw::SetCFContainers(const char* filename)
149 {
150   //
151   // get CF containers from file
152   //
153
154   TFile f(filename);
155   TList *l=f.GetListOfKeys();
156   TIter nextKey(l);
157   TKey *k=0x0;
158   while ( (k=static_cast<TKey*>(nextKey())) ){
159     TObject *o=k->ReadObj();
160     if (o->IsA()->InheritsFrom(TSeqCollection::Class())){
161       TSeqCollection *arr=static_cast<TSeqCollection*>(o);
162       SetCFContainers(arr);
163     } else if (o->IsA()==AliCFContainer::Class()){
164       fCfContainer=static_cast<AliCFContainer*>(o);
165       if (fEffGrid) delete fEffGrid;
166       fEffGrid=new AliCFEffGrid("eff","eff",*fCfContainer);
167     }
168   }
169 }
170
171 //________________________________________________________________
172 void AliDielectronCFdraw::SetRangeUser(Int_t ivar, Double_t min, Double_t max, const char* slices)
173 {
174   //
175   // Set range of cut steps defined in slices
176   // Steps may be separated by one the the characteres ,;:
177   //
178   TObjArray *arr=TString(slices).Tokenize(",:;");
179
180   if (arr->GetEntriesFast()==0){
181     // all slices in case of 0 entries
182     for (Int_t istep=0; istep<fCfContainer->GetNStep(); ++istep){
183       fCfContainer->SetRangeUser(ivar,min,max,istep);
184     }
185   } else {
186     TIter next(arr);
187     TObjString *ostr=0x0;
188     while ( (ostr=static_cast<TObjString*>(next())) ) {
189       Int_t istep=ostr->GetString().Atoi();
190       fCfContainer->SetRangeUser(ivar,min,max,istep);
191     }
192   }
193   delete arr;
194 }
195
196 //________________________________________________________________
197 void AliDielectronCFdraw::UnsetRangeUser(Int_t ivar, const char* slices)
198 {
199   //
200   // Unset range of cut steps defined in slices
201   // Steps may be separated by one the the characteres ,;:
202   //
203   TObjArray *arr=TString(slices).Tokenize(",:;");
204   
205   if (arr->GetEntriesFast()==0){
206     // all slices in case of 0 entries
207     for (Int_t istep=0; istep<fCfContainer->GetNStep(); ++istep){
208       fCfContainer->GetAxis(ivar,istep)->SetRange(0,0);
209     }
210   } else {
211     TIter next(arr);
212     TObjString *ostr=0x0;
213     while ( (ostr=static_cast<TObjString*>(next())) ) {
214       Int_t istep=ostr->GetString().Atoi();
215       fCfContainer->GetAxis(ivar,istep)->SetRange(0,0);
216     }
217   }
218   delete arr;
219 }
220
221 //________________________________________________________________
222 void AliDielectronCFdraw::Draw(const Option_t* varnames, const char* opt, const char* slices)
223 {
224   //
225   // Draw 'variables' of 'slices'
226   // for multidimensional draw variables may be separated by a ':'
227   // slice numbers may be separated by any of ,:;
228   //
229   // variables may be called by either their name or number
230   //
231
232   TObjArray *arrVars=TString(varnames).Tokenize(":");
233   Int_t entries=arrVars->GetEntriesFast();
234   if (entries<1||entries>3){
235     AliError("Wrong number of variables, supported are 1 - 3 dimensions");
236     delete arrVars;
237     return;
238   }
239   
240   TIter next(arrVars);
241   TObjString *ostr=0x0;
242   Int_t ivar[3]={-1,-1,-1};
243   for (Int_t i=0; i<entries; ++i){
244     ostr=static_cast<TObjString*>(next());
245     if (ostr->GetString().IsDigit()){
246       ivar[i]=ostr->GetString().Atoi();
247     } else {
248       ivar[i]=fCfContainer->GetVar(ostr->GetName());
249     }
250   }
251
252   switch (entries){
253   case 1:
254     Draw(ivar[0],opt,slices);
255     break;
256   case 2:
257     Draw(ivar[1],ivar[0],opt,slices);
258     break;
259   case 3:
260     Draw(ivar[2],ivar[1],ivar[0],opt,slices);
261     break;
262   }
263   delete arrVars;
264 }
265
266 //________________________________________________________________
267 void AliDielectronCFdraw::Draw(Int_t var, const char* opt, const char* slices)
268 {
269   //
270   // Draw variable var for all slices
271   // slices may be divided by and of ,;:
272   //
273   // if opt contains 'same' all histograms are drawn in the same pad
274   // otherwise the pad will be divided in sub pads and the histograms
275   // are drawn in each sub pad
276   //
277
278   const Int_t ndim=1;
279   Int_t vars[ndim]={var};
280   TObjArray *arr=CollectHistosProj(ndim,vars,slices);
281   Draw(arr,opt);
282   delete arr; 
283 }
284
285 //________________________________________________________________
286 void AliDielectronCFdraw::Draw(Int_t var0, Int_t var1, const char* opt, const char* slices)
287 {
288   //
289   // Draw 2D case
290   //
291   const Int_t ndim=2;
292   Int_t vars[ndim]={var0,var1};
293   TObjArray *arr=CollectHistosProj(ndim,vars,slices);
294   Draw(arr,opt);
295   delete arr;
296 }
297
298 //________________________________________________________________
299 void AliDielectronCFdraw::Draw(Int_t var0, Int_t var1, Int_t var2, const char* opt, const char* slices)
300 {
301   //
302   // Draw 3D case
303   //
304   const Int_t ndim=3;
305   Int_t vars[ndim]={var0,var1,var2};
306   TObjArray *arr=CollectHistosProj(ndim,vars,slices);
307   Draw(arr,opt);
308   delete arr;
309 }
310
311 //________________________________________________________________
312 TObjArray* AliDielectronCFdraw::CollectHistosProj(Int_t dim, Int_t *vars, const char* slices)
313 {
314   //
315   // Collect histos with 'dim'ension of the 'slices' separated by one of "':;'"
316   // in a TObjArray and return it
317   //
318   TObjArray *arr=TString(slices).Tokenize(",:;");
319   TObjArray *arrHists=0x0;
320   if (arr->GetEntriesFast()==0){
321     // all slices in case of 0 entries
322     arrHists=new TObjArray(fCfContainer->GetNStep());
323     for (Int_t istep=0; istep<fCfContainer->GetNStep(); ++istep){
324       TH1 *hproj=Project(dim,vars,istep);
325       if (!hproj) continue;
326       hproj->SetName(Form("proj_%02d",istep));
327       hproj->SetTitle(fCfContainer->GetStepTitle(istep));
328       arrHists->Add(hproj);
329     }
330   } else {
331     arrHists=new TObjArray(arr->GetEntriesFast());
332     TIter next(arr);
333     TObjString *ostr=0x0;
334     while ( (ostr=static_cast<TObjString*>(next())) ) {
335       Int_t istep=ostr->GetString().Atoi();
336       TH1 *hproj=Project(dim,vars,istep);
337       if (!hproj) continue;
338       hproj->SetName(Form("proj_%02d",istep));
339       hproj->SetTitle(fCfContainer->GetStepTitle(istep));
340       arrHists->Add(hproj);
341     }
342   }
343   delete arr;
344
345   return arrHists;
346 }
347
348 //________________________________________________________________
349 TH1* AliDielectronCFdraw::Project(Int_t ndim, Int_t *vars, Int_t slice)
350 {
351   //
352   // Do an ndim projection
353   //
354   switch (ndim){
355   case 1:
356     return fCfContainer->Project(vars[0],slice);
357     break;
358   case 2:
359     return fCfContainer->Project(vars[0],vars[1],slice);
360     break;
361   case 3:
362     return fCfContainer->Project(vars[0],vars[1],vars[2],slice);
363     break;
364   }
365   return 0x0;
366 }
367
368 //________________________________________________________________
369 TH1* AliDielectronCFdraw::Project(const Option_t* var, Int_t slice)
370 {
371   //
372   // translate variable names and do projection
373   //
374   TObjArray *arrVars=TString(var).Tokenize(":");
375   Int_t entries=arrVars->GetEntriesFast();
376   if (entries<1||entries>3){
377     AliError("Wrong number of variables, supported are 1 - 3 dimensions");
378     delete arrVars;
379     return 0x0;
380   }
381   
382   TIter next(arrVars);
383   TObjString *ostr=0x0;
384   Int_t ivar[3]={-1,-1,-1};
385   for (Int_t i=0; i<entries; ++i){
386     ostr=static_cast<TObjString*>(next());
387     if (ostr->GetString().IsDigit()){
388       ivar[i]=ostr->GetString().Atoi();
389     } else {
390       ivar[i]=fCfContainer->GetVar(ostr->GetName());
391     }
392   }
393   delete arrVars;
394   return Project(entries,ivar,slice);
395 }
396
397 //________________________________________________________________
398 void AliDielectronCFdraw::DrawEfficiency(const char* varnames, const char* nominators, Int_t denominator, const char* opt)
399 {
400   //
401   // plot efficiencies for variables. Variable may be up to 3 dim, separated by a ':'
402   // you may have several nominators, sparated by one of ',:;'
403   //
404   
405   TObjArray *arrVars=TString(varnames).Tokenize(":");
406   Int_t entries=arrVars->GetEntriesFast();
407   if (entries<1||entries>3){
408     AliError("Wrong number of variables, supported are 1 - 3 dimensions");
409     delete arrVars;
410     return;
411   }
412   
413   TIter next(arrVars);
414   TObjString *ostr=0x0;
415   Int_t ivar[3]={-1,-1,-1};
416   for (Int_t i=0; i<entries; ++i){
417     ostr=static_cast<TObjString*>(next());
418     if (ostr->GetString().IsDigit()){
419       ivar[i]=ostr->GetString().Atoi();
420     } else {
421       ivar[i]=fCfContainer->GetVar(ostr->GetName());
422     }
423   }
424
425   Int_t type=0;
426   TString optStr(opt);
427   if (optStr.Contains("2")) type=1;
428   
429   switch (entries){
430   case 1:
431     DrawEfficiency(ivar[0],nominators, denominator,opt,type);
432     break;
433   case 2:
434     DrawEfficiency(ivar[1],ivar[0], nominators, denominator,opt,type);
435     break;
436   case 3:
437     DrawEfficiency(ivar[2],ivar[1],ivar[0],nominators, denominator,opt,type);
438     break;
439   }
440   delete arrVars;
441 }
442
443 //________________________________________________________________
444 void AliDielectronCFdraw::DrawEfficiency(Int_t var, const char* nominators, Int_t denominator, const char* opt, Int_t type)
445 {
446   //
447   // Draw Efficiencies for all nominators
448   // nominators may be divided by and of ,;:
449   //
450   // if opt contains 'same' all histograms are drawn in the same pad
451   // otherwise the pad will be divided in sub pads and the histograms
452   // are drawn in each sub pad
453   //
454   
455   const Int_t ndim=1;
456   Int_t vars[ndim]={var};
457   TObjArray *arr=CollectHistosEff(ndim,vars,nominators,denominator,type);
458   TString drawOpt=opt;
459   drawOpt+="eff";
460   Draw(arr,drawOpt);
461   delete arr;
462 }
463
464 //________________________________________________________________
465 void AliDielectronCFdraw::DrawEfficiency(Int_t var0, Int_t var1, const char* nominators, Int_t denominator, const char* opt, Int_t type)
466 {
467   //
468   // Draw 2D case
469   //
470   const Int_t ndim=2;
471   Int_t vars[ndim]={var0,var1};
472   TObjArray *arr=CollectHistosEff(ndim,vars,nominators,denominator,type);
473   TString drawOpt=opt;
474   drawOpt+="eff";
475   Draw(arr,drawOpt);
476   delete arr;
477 }
478
479 //________________________________________________________________
480 void AliDielectronCFdraw::DrawEfficiency(Int_t var0, Int_t var1, Int_t var2, const char* nominators, Int_t denominator, const char* opt, Int_t type)
481 {
482   //
483   // Draw 3D case
484   //
485   const Int_t ndim=3;
486   Int_t vars[ndim]={var0,var1,var2};
487   TObjArray *arr=CollectHistosEff(ndim,vars,nominators,denominator,type);
488   TString drawOpt=opt;
489   drawOpt+="eff";
490   Draw(arr,drawOpt);
491   delete arr;
492 }
493
494 //________________________________________________________________
495 TObjArray* AliDielectronCFdraw::CollectHistosEff(Int_t dim, Int_t *vars, const char* nominators, Int_t denominator, Int_t type)
496 {
497   //
498   // Collect histos with 'dim'ension of the 'slices' separated by one of "':;'"
499   // in a TObjArray and return it
500   //
501   TObjArray *arr=TString(nominators).Tokenize(",:;");
502   TObjArray *arrHists=0x0;
503
504   if (type==0){
505     if (arr->GetEntriesFast()==0){
506     // all slices in case of 0 entries
507       arrHists=new TObjArray(fCfContainer->GetNStep());
508       fVdata.ResizeTo(arrHists->GetSize());
509       for (Int_t istep=0; istep<fCfContainer->GetNStep(); ++istep){
510         fEffGrid->CalculateEfficiency(istep,denominator);
511         TH1 *hproj=ProjectEff(dim,vars);
512         if (!hproj) continue;
513         Float_t eff=fEffGrid->GetAverage();
514         fVdata(istep)=eff;
515         hproj->SetName(Form("eff_%02d/%02d",istep,denominator));
516         hproj->SetTitle(Form("%s (%.3f)",fCfContainer->GetStepTitle(istep),eff));
517         arrHists->Add(hproj);
518       }
519     } else {
520       arrHists=new TObjArray(arr->GetEntriesFast());
521       TIter next(arr);
522       TObjString *ostr=0x0;
523       fVdata.ResizeTo(arrHists->GetSize());
524       Int_t count=0;
525       while ( (ostr=static_cast<TObjString*>(next())) ) {
526         Int_t istep=ostr->GetString().Atoi();
527         fEffGrid->CalculateEfficiency(istep,denominator);
528         TH1 *hproj=ProjectEff(dim,vars);
529         if (!hproj) continue;
530         Float_t eff=fEffGrid->GetAverage();
531         fVdata(count++)=eff;
532         hproj->SetName(Form("eff_%02d/%02d",istep,denominator));
533         hproj->SetTitle(Form("%s (%.3f)",fCfContainer->GetStepTitle(istep),eff));
534         arrHists->Add(hproj);
535       }
536     }
537   }
538
539   //second approach
540   if (type==1){
541     TH1 *hDen=Project(dim,vars,denominator);
542     Double_t entriesDen=hDen->GetEffectiveEntries();
543     if (arr->GetEntriesFast()==0){
544     // all slices in case of 0 entries
545       arrHists=new TObjArray(fCfContainer->GetNStep());
546       fVdata.ResizeTo(arrHists->GetSize());
547       for (Int_t istep=0; istep<fCfContainer->GetNStep(); ++istep){
548         TH1 *hproj=Project(dim,vars,istep);
549         if (!hproj) continue;
550         Float_t eff=0;
551         if (entriesDen>0) eff=hproj->GetEffectiveEntries()/entriesDen;
552         fVdata(istep)=eff;
553         hproj->Divide(hproj,hDen,1,1,"B");
554         hproj->SetName(Form("eff_%02d/%02d",istep,denominator));
555         hproj->SetTitle(Form("%s (%.3f)",fCfContainer->GetStepTitle(istep),eff));
556         arrHists->Add(hproj);
557       }
558     } else {
559       arrHists=new TObjArray(arr->GetEntriesFast());
560       fVdata.ResizeTo(arrHists->GetSize());
561       TIter next(arr);
562       TObjString *ostr=0x0;
563       Int_t count=0;
564       while ( (ostr=static_cast<TObjString*>(next())) ) {
565         Int_t istep=ostr->GetString().Atoi();
566         TH1 *hproj=Project(dim,vars,istep);
567         if (!hproj) continue;
568         Float_t eff=0;
569         if (entriesDen>0) eff=hproj->GetEffectiveEntries()/entriesDen;
570         fVdata(count++)=eff;
571         hproj->Divide(hproj,hDen,1,1,"B");
572         hproj->SetName(Form("eff_%02d/%02d",istep,denominator));
573         hproj->SetTitle(Form("%s (%.3f)",fCfContainer->GetStepTitle(istep),eff));
574         arrHists->Add(hproj);
575       }
576     }
577     delete hDen;
578   }
579   
580
581   delete arr;
582   return arrHists;
583 }
584
585 //________________________________________________________________
586 TH1* AliDielectronCFdraw::ProjectEff(Int_t ndim, Int_t *vars)
587 {
588   //
589   // Do an nim projection
590   //
591   switch (ndim){
592   case 1:
593     return fEffGrid->Project(vars[0]);
594     break;
595   case 2:
596     return fEffGrid->Project(vars[0],vars[1]);
597     break;
598   case 3:
599     return fEffGrid->Project(vars[0],vars[1],vars[2]);
600     break;
601   }
602   return 0x0;
603 }
604
605 //________________________________________________________________
606 void AliDielectronCFdraw::Draw(const TObjArray *arr, const char* opt)
607 {
608   //
609   // draw all objects in arr
610   //
611   TString optStr(opt);
612   optStr.ToLower();
613   Bool_t drawSame     = optStr.Contains("same");
614   Bool_t drawSamePlus = optStr.Contains("same+");
615   Bool_t drawEff      = optStr.Contains("eff");
616   Bool_t optLeg       = optStr.Contains("leg");
617   Bool_t optScaleMax  = optStr.Contains("max");
618   
619   if (!drawSamePlus) optStr.ReplaceAll("same","");
620   
621   optStr.ReplaceAll("+","");
622   optStr.ReplaceAll("eff","");
623   optStr.ReplaceAll("leg","");
624   optStr.ReplaceAll("max","");
625   
626   if (!gPad) new TCanvas;
627   
628   Int_t nPads = arr->GetEntriesFast();
629   if (nPads==0) return;
630   
631 //   if (nPads==1){
632 //     arr->UncheckedAt(0)->Draw(optStr.Data());
633 //     return;
634 //   }
635   
636   TCanvas *c=gPad->GetCanvas();
637   if (!gVirtualPS&&!drawSamePlus&&!drawSame&&nPads>1) c->Clear();
638   
639   if (!drawSame&&nPads>1){
640     //optimised division
641     Int_t nCols = (Int_t)TMath::Ceil( TMath::Sqrt(nPads) );
642     Int_t nRows = (Int_t)TMath::Ceil( (Double_t)nPads/(Double_t)nCols );
643     c->Divide(nCols,nRows);
644     for (Int_t i=0; i<nPads; ++i){
645       c->cd(i+1);
646       arr->UncheckedAt(i)->Draw(optStr.Data());
647     }
648   } else {
649     TLegend *leg=0;
650     if (drawSamePlus){
651       //find already existing legend to attach entries to it
652       TIter nextPrimitive(gPad->GetListOfPrimitives());
653       TObject *o=0x0;
654       while ((o=nextPrimitive())) if (o->IsA()==TLegend::Class()) leg=(TLegend*)o;
655     }
656     
657     if (optLeg&&!leg) leg=new TLegend(.2,.3,.99,.9);
658     Int_t addColor=0;
659     if (leg) addColor=leg->GetNRows();
660     Int_t offset=20;
661     if (nPads<7) offset=24;
662     for (Int_t i=0; i<nPads; ++i){
663       if (i==1&&!drawSamePlus) optStr+="same";
664       TH1 *hist=static_cast<TH1*>(arr->UncheckedAt(i));
665       hist->SetLineColor(i+1+addColor);
666       hist->SetLineWidth(2);
667       hist->SetMarkerColor(i+1+addColor);
668       hist->SetMarkerStyle(offset+i+addColor);
669       hist->Draw(optStr.Data());
670       
671       if (drawEff&&i==0&&!drawSamePlus) {
672         hist->GetYaxis()->SetRangeUser(0,2);
673         hist->SetYTitle("Rec. Signal / Gen. Signal");
674       }
675       
676       if (leg) leg->AddEntry(hist,hist->GetTitle(),"lp");
677     }
678     if (leg){
679       leg->SetFillColor(10);
680       leg->SetY1(.9-leg->GetNRows()*.05);
681       leg->SetY1NDC(.9-leg->GetNRows()*.05);
682       leg->SetMargin(.1);
683       if (!drawSamePlus) leg->Draw();
684     }
685     if (optScaleMax){
686       TIter nextPrimitive(gPad->GetListOfPrimitives());
687       TObject *o=0x0;
688       TH1 *hfirst=0x0;
689       Float_t max=0;
690       while ((o=nextPrimitive())) if (o->InheritsFrom(TH1::Class())){
691         TH1 *h=(TH1*)o;
692         if (!hfirst) hfirst=h;
693         if (h->GetMaximum()>max) max=h->GetMaximum();
694       }
695       max*=1.1;
696       hfirst->SetMaximum(max);
697     }
698   }
699   
700 }