]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/dielectron/AliDielectronCFdraw.cxx
Updates and additions: Classes for signal and spectrum extraction; saving of
[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 (!o) continue;
123     nstep+=cf->GetNStep();
124   }
125   Int_t nbins[1]={1};
126   fCfContainer=new AliCFContainer(GetName(), GetTitle(), nstep, 1, nbins);
127
128   //delete unneeded steps
129 //   for (Int_t istep=0; istep<nstep; ++istep) delete fCfContainer->GetGrid(istep);
130
131   //add step to the new container
132   Int_t istep=0;
133   for (Int_t icf=0; icf<arr->GetEntries(); ++icf){
134     AliCFContainer *cf=dynamic_cast<AliCFContainer*>(arr->At(icf));
135     if (!cf) continue;
136     for (Int_t istepCurr=0; istepCurr<cf->GetNStep(); ++istepCurr){
137       fCfContainer->SetGrid(istep, cf->GetGrid(istepCurr));
138       fCfContainer->SetStepTitle(istep,Form("%s, Pair: %s",cf->GetTitle(),cf->GetStepTitle(istepCurr)));
139       ++istep;
140     }
141   }
142   if (fEffGrid) delete fEffGrid;
143   fEffGrid=new AliCFEffGrid("eff","eff",*fCfContainer);
144 }
145
146 //________________________________________________________________
147 void AliDielectronCFdraw::SetCFContainers(const char* filename)
148 {
149   //
150   // get CF containers from file
151   //
152
153   TFile f(filename);
154   TList *l=f.GetListOfKeys();
155   Int_t entries=l->GetEntries();
156   if (entries==0) return;
157   
158   TKey *k=(TKey*)l->At(0);
159   if (!k) return;
160   TObject *o=k->ReadObj();
161   if (o->IsA()->InheritsFrom(TSeqCollection::Class())){
162     TSeqCollection *arr=static_cast<TSeqCollection*>(o);
163     SetCFContainers(arr);
164   } else if (o->IsA()==AliCFContainer::Class()){
165     fCfContainer=static_cast<AliCFContainer*>(o);
166     if (fEffGrid) delete fEffGrid;
167     fEffGrid=new AliCFEffGrid("eff","eff",*fCfContainer);
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 nim 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 void AliDielectronCFdraw::DrawEfficiency(const char* varnames, const char* nominators, Int_t denominator, const char* opt)
370 {
371   //
372   // plot efficiencies for variables. Variable may be up to 3 dim, separated by a ':'
373   // you may have several nominators, sparated by one of ',:;'
374   //
375   
376   TObjArray *arrVars=TString(varnames).Tokenize(":");
377   Int_t entries=arrVars->GetEntriesFast();
378   if (entries<1||entries>3){
379     AliError("Wrong number of variables, supported are 1 - 3 dimensions");
380     delete arrVars;
381     return;
382   }
383   
384   TIter next(arrVars);
385   TObjString *ostr=0x0;
386   Int_t ivar[3]={-1,-1,-1};
387   for (Int_t i=0; i<entries; ++i){
388     ostr=static_cast<TObjString*>(next());
389     if (ostr->GetString().IsDigit()){
390       ivar[i]=ostr->GetString().Atoi();
391     } else {
392       ivar[i]=fCfContainer->GetVar(ostr->GetName());
393     }
394   }
395
396   Int_t type=0;
397   TString optStr(opt);
398   if (optStr.Contains("2")) type=1;
399   
400   switch (entries){
401   case 1:
402     DrawEfficiency(ivar[0],nominators, denominator,opt,type);
403     break;
404   case 2:
405     DrawEfficiency(ivar[1],ivar[0], nominators, denominator,opt,type);
406     break;
407   case 3:
408     DrawEfficiency(ivar[2],ivar[1],ivar[0],nominators, denominator,opt,type);
409     break;
410   }
411   delete arrVars;
412 }
413
414 //________________________________________________________________
415 void AliDielectronCFdraw::DrawEfficiency(Int_t var, const char* nominators, Int_t denominator, const char* opt, Int_t type)
416 {
417   //
418   // Draw Efficiencies for all nominators
419   // nominators may be divided by and of ,;:
420   //
421   // if opt contains 'same' all histograms are drawn in the same pad
422   // otherwise the pad will be divided in sub pads and the histograms
423   // are drawn in each sub pad
424   //
425   
426   const Int_t ndim=1;
427   Int_t vars[ndim]={var};
428   TObjArray *arr=CollectHistosEff(ndim,vars,nominators,denominator,type);
429   TString drawOpt=opt;
430   drawOpt+="eff";
431   Draw(arr,drawOpt);
432   delete arr;
433 }
434
435 //________________________________________________________________
436 void AliDielectronCFdraw::DrawEfficiency(Int_t var0, Int_t var1, const char* nominators, Int_t denominator, const char* opt, Int_t type)
437 {
438   //
439   // Draw 2D case
440   //
441   const Int_t ndim=2;
442   Int_t vars[ndim]={var0,var1};
443   TObjArray *arr=CollectHistosEff(ndim,vars,nominators,denominator,type);
444   TString drawOpt=opt;
445   drawOpt+="eff";
446   Draw(arr,drawOpt);
447   delete arr;
448 }
449
450 //________________________________________________________________
451 void AliDielectronCFdraw::DrawEfficiency(Int_t var0, Int_t var1, Int_t var2, const char* nominators, Int_t denominator, const char* opt, Int_t type)
452 {
453   //
454   // Draw 3D case
455   //
456   const Int_t ndim=3;
457   Int_t vars[ndim]={var0,var1,var2};
458   TObjArray *arr=CollectHistosEff(ndim,vars,nominators,denominator,type);
459   TString drawOpt=opt;
460   drawOpt+="eff";
461   Draw(arr,drawOpt);
462   delete arr;
463 }
464
465 //________________________________________________________________
466 TObjArray* AliDielectronCFdraw::CollectHistosEff(Int_t dim, Int_t *vars, const char* nominators, Int_t denominator, Int_t type)
467 {
468   //
469   // Collect histos with 'dim'ension of the 'slices' separated by one of "':;'"
470   // in a TObjArray and return it
471   //
472   TObjArray *arr=TString(nominators).Tokenize(",:;");
473   TObjArray *arrHists=0x0;
474
475   if (type==0){
476     if (arr->GetEntriesFast()==0){
477     // all slices in case of 0 entries
478       arrHists=new TObjArray(fCfContainer->GetNStep());
479       fVdata.ResizeTo(arrHists->GetSize());
480       for (Int_t istep=0; istep<fCfContainer->GetNStep(); ++istep){
481         fEffGrid->CalculateEfficiency(istep,denominator);
482         TH1 *hproj=ProjectEff(dim,vars);
483         if (!hproj) continue;
484         Float_t eff=fEffGrid->GetAverage();
485         fVdata(istep)=eff;
486         hproj->SetName(Form("eff_%02d/%02d",istep,denominator));
487         hproj->SetTitle(Form("%s (%.3f)",fCfContainer->GetStepTitle(istep),eff));
488         arrHists->Add(hproj);
489       }
490     } else {
491       arrHists=new TObjArray(arr->GetEntriesFast());
492       TIter next(arr);
493       TObjString *ostr=0x0;
494       fVdata.ResizeTo(arrHists->GetSize());
495       Int_t count=0;
496       while ( (ostr=static_cast<TObjString*>(next())) ) {
497         Int_t istep=ostr->GetString().Atoi();
498         fEffGrid->CalculateEfficiency(istep,denominator);
499         TH1 *hproj=ProjectEff(dim,vars);
500         if (!hproj) continue;
501         Float_t eff=fEffGrid->GetAverage();
502         fVdata(count++)=eff;
503         hproj->SetName(Form("eff_%02d/%02d",istep,denominator));
504         hproj->SetTitle(Form("%s (%.3f)",fCfContainer->GetStepTitle(istep),eff));
505         arrHists->Add(hproj);
506       }
507     }
508   }
509
510   //second approach
511   if (type==1){
512     TH1 *hDen=Project(dim,vars,denominator);
513     Double_t entriesDen=hDen->GetEffectiveEntries();
514     if (arr->GetEntriesFast()==0){
515     // all slices in case of 0 entries
516       arrHists=new TObjArray(fCfContainer->GetNStep());
517       fVdata.ResizeTo(arrHists->GetSize());
518       for (Int_t istep=0; istep<fCfContainer->GetNStep(); ++istep){
519         TH1 *hproj=Project(dim,vars,istep);
520         if (!hproj) continue;
521         Float_t eff=0;
522         if (entriesDen>0) eff=hproj->GetEffectiveEntries()/entriesDen;
523         fVdata(istep)=eff;
524         hproj->Divide(hDen);
525         hproj->SetName(Form("eff_%02d/%02d",istep,denominator));
526         hproj->SetTitle(Form("%s (%.3f)",fCfContainer->GetStepTitle(istep),eff));
527         arrHists->Add(hproj);
528       }
529     } else {
530       arrHists=new TObjArray(arr->GetEntriesFast());
531       fVdata.ResizeTo(arrHists->GetSize());
532       TIter next(arr);
533       TObjString *ostr=0x0;
534       Int_t count=0;
535       while ( (ostr=static_cast<TObjString*>(next())) ) {
536         Int_t istep=ostr->GetString().Atoi();
537         TH1 *hproj=Project(dim,vars,istep);
538         if (!hproj) continue;
539         Float_t eff=0;
540         if (entriesDen>0) eff=hproj->GetEffectiveEntries()/entriesDen;
541         fVdata(count++)=eff;
542         hproj->Divide(hDen);
543         hproj->SetName(Form("eff_%02d/%02d",istep,denominator));
544         hproj->SetTitle(Form("%s (%.3f)",fCfContainer->GetStepTitle(istep),eff));
545         arrHists->Add(hproj);
546       }
547     }
548     delete hDen;
549   }
550   
551
552   delete arr;
553   return arrHists;
554 }
555
556 //________________________________________________________________
557 TH1* AliDielectronCFdraw::ProjectEff(Int_t ndim, Int_t *vars)
558 {
559   //
560   // Do an nim projection
561   //
562   switch (ndim){
563   case 1:
564     return fEffGrid->Project(vars[0]);
565     break;
566   case 2:
567     return fEffGrid->Project(vars[0],vars[1]);
568     break;
569   case 3:
570     return fEffGrid->Project(vars[0],vars[1],vars[2]);
571     break;
572   }
573   return 0x0;
574 }
575
576 //________________________________________________________________
577 void AliDielectronCFdraw::Draw(const TObjArray *arr, const char* opt)
578 {
579   //
580   // draw all objects in arr
581   //
582   TString optStr(opt);
583   optStr.ToLower();
584   Bool_t drawSame     = optStr.Contains("same");
585   Bool_t drawSamePlus = optStr.Contains("same+");
586   Bool_t drawEff      = optStr.Contains("eff");
587   Bool_t optLeg       = optStr.Contains("leg");
588   Bool_t optScaleMax  = optStr.Contains("max");
589   
590   if (!drawSamePlus) optStr.ReplaceAll("same","");
591   
592   optStr.ReplaceAll("+","");
593   optStr.ReplaceAll("eff","");
594   optStr.ReplaceAll("leg","");
595   optStr.ReplaceAll("max","");
596   
597   if (!gPad) new TCanvas;
598   
599   Int_t nPads = arr->GetEntriesFast();
600   if (nPads==0) return;
601   
602 //   if (nPads==1){
603 //     arr->UncheckedAt(0)->Draw(optStr.Data());
604 //     return;
605 //   }
606   
607   TCanvas *c=gPad->GetCanvas();
608   if (!gVirtualPS&&!drawSamePlus) c->Clear();
609   
610   
611   if (!drawSame){
612     //optimised division
613     Int_t nCols = (Int_t)TMath::Ceil( TMath::Sqrt(nPads) );
614     Int_t nRows = (Int_t)TMath::Ceil( (Double_t)nPads/(Double_t)nCols );
615     c->Divide(nCols,nRows);
616     for (Int_t i=0; i<nPads; ++i){
617       c->cd(i+1);
618       arr->UncheckedAt(i)->Draw(optStr.Data());
619     }
620   } else {
621     TLegend *leg=0;
622     if (drawSamePlus){
623       //find already existing legend to attach entries to it
624       TIter nextPrimitive(gPad->GetListOfPrimitives());
625       TObject *o=0x0;
626       while ((o=nextPrimitive())) if (o->IsA()==TLegend::Class()) leg=(TLegend*)o;
627     }
628     
629     if (optLeg&&!leg) leg=new TLegend(.2,.3,.99,.9);
630     Int_t addColor=0;
631     if (leg) addColor=leg->GetNRows();
632     Int_t offset=20;
633     if (nPads<7) offset=24;
634     for (Int_t i=0; i<nPads; ++i){
635       if (i==1&&!drawSamePlus) optStr+="same";
636       TH1 *hist=static_cast<TH1*>(arr->UncheckedAt(i));
637       hist->SetLineColor(i+1+addColor);
638       hist->SetLineWidth(2);
639       hist->SetMarkerColor(i+1+addColor);
640       hist->SetMarkerStyle(offset+i+addColor);
641       hist->Draw(optStr.Data());
642       
643       if (drawEff&&i==0&&!drawSamePlus) {
644         hist->GetYaxis()->SetRangeUser(0,2);
645         hist->SetYTitle("Rec. Signal / Gen. Signal");
646       }
647       
648       if (leg) leg->AddEntry(hist,hist->GetTitle(),"lp");
649     }
650     if (leg){
651       leg->SetFillColor(10);
652       leg->SetY1(.9-leg->GetNRows()*.05);
653       leg->SetY1NDC(.9-leg->GetNRows()*.05);
654       leg->SetMargin(.1);
655       if (!drawSamePlus) leg->Draw();
656     }
657     if (optScaleMax){
658       TIter nextPrimitive(gPad->GetListOfPrimitives());
659       TObject *o=0x0;
660       TH1 *hfirst=0x0;
661       Float_t max=0;
662       while ((o=nextPrimitive())) if (o->InheritsFrom(TH1::Class())){
663         TH1 *h=(TH1*)o;
664         if (!hfirst) hfirst=h;
665         if (h->GetMaximum()>max) max=h->GetMaximum();
666       }
667       max*=1.1;
668       hfirst->SetMaximum(max);
669     }
670   }
671   
672 }