]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGPP/MUON/lite/AliAnalysisTaskTrigChEff.cxx
1) Adding class AliAnalysisMuonUtility which contains static methods allowing to...
[u/mrichter/AliRoot.git] / PWGPP / MUON / lite / AliAnalysisTaskTrigChEff.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2007, 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 #include "AliAnalysisTaskTrigChEff.h"
19
20 // ROOT includes
21 #include "TH1.h"
22 #include "TH2.h"
23 #include "TCanvas.h"
24 #include "TLegend.h"
25 #include "TMath.h"
26 #include "TObjString.h"
27 #include "TObjArray.h"
28 #include "TStyle.h"
29 #include "TROOT.h"
30 #include "TGraphAsymmErrors.h"
31 #include "TList.h"
32 #include "TFile.h"
33
34 // STEER includes
35 #include "AliVParticle.h"
36 #include "AliAODTrack.h"
37 #include "AliESDMuonTrack.h"
38
39 // ANALYSIS includes
40 #include "AliAnalysisManager.h"
41 #include "AliAnalysisDataSlot.h"
42 #include "AliAnalysisDataContainer.h"
43
44 // PWG3 includes
45 #include "AliVAnalysisMuon.h"
46 #include "AliMuonEventCuts.h"
47 #include "AliMuonTrackCuts.h"
48 #include "AliMergeableCollection.h"
49 #include "AliAnalysisMuonUtility.h"
50
51
52 /// \cond CLASSIMP
53 ClassImp(AliAnalysisTaskTrigChEff) // Class implementation in ROOT context
54 /// \endcond
55
56
57 //________________________________________________________________________
58 AliAnalysisTaskTrigChEff::AliAnalysisTaskTrigChEff() :
59   AliVAnalysisMuon(),
60   fTrackSelKeys(0x0),
61   fCountTypeKeys(0x0),
62   fHistoTypeKeys(0x0),
63   fEffMethodKeys(0x0),
64   fUseGhosts(kFALSE),
65   fList(0x0)
66 {
67   /// Default ctor.
68 }
69
70 //________________________________________________________________________
71 AliAnalysisTaskTrigChEff::AliAnalysisTaskTrigChEff(const char *name, const AliMuonTrackCuts& cuts) :
72   AliVAnalysisMuon(name, cuts),
73   fTrackSelKeys(0x0),
74   fCountTypeKeys(0x0),
75   fHistoTypeKeys(0x0),
76   fEffMethodKeys(0x0),
77   fUseGhosts(kFALSE),
78   fList(0x0)
79 {
80   //
81   /// Constructor.
82   //
83
84   InitLocalKeys();
85   
86   DefineOutput(2, TList::Class());
87 }
88
89
90 //________________________________________________________________________
91 AliAnalysisTaskTrigChEff::~AliAnalysisTaskTrigChEff()
92 {
93   //
94   /// Destructor
95   //
96   delete fTrackSelKeys;
97   delete fCountTypeKeys;
98   delete fHistoTypeKeys;
99   delete fEffMethodKeys;
100   if ( ! AliAnalysisManager::GetAnalysisManager() || ! AliAnalysisManager::GetAnalysisManager()->IsProofMode() ) {
101     delete fList;
102   }
103 }
104
105 //___________________________________________________________________________
106 TList* AliAnalysisTaskTrigChEff::GetEffHistoList(TString physSel, TString trigClassNames, TString centrality, TString trackSelection)
107 {
108   /// Get the list of efficiency objects by merging the 
109   // results from the histogram collection
110   
111   TList* outList = new TList();
112   outList->SetOwner();
113   FillEffHistoList(physSel, trigClassNames, centrality, trackSelection, outList);
114   return outList;
115 }
116
117 //___________________________________________________________________________
118 Bool_t AliAnalysisTaskTrigChEff::FillEffHistoList(TString physSel, TString trigClassNames, TString centrality, TString trackSelection, TList* outList)
119 {
120   /// Fill the list of objects for the efficiency calculation
121   /// merging the splitted output of the fHistogramCollection
122   /// The obtained list can be converted in the efficiency map used in simulations
123   /// in a backward compatible way
124   
125   if ( ! fMergeableCollection ) return kFALSE;
126   if ( ! outList ) outList = fList;
127   TString histoName = "";
128   TString histoPattern = "";
129   TH1* histo = 0x0;
130   Bool_t isOk = kTRUE;
131   for ( Int_t icount=0; icount<kNcounts; ++icount ) {
132     histoName = GetHistoName(kHchamberEff, icount, -1, -1, -1);
133     histoPattern = Form("%s%s", histoName.Data(), trackSelection.Data());
134     histo = (TH1*)GetSum(physSel, trigClassNames, centrality, histoPattern);
135     if ( histo ) {
136       histo->SetName(histoName.Data());
137       histo->SetTitle(histoName.Data());
138       histo->SetDirectory(0);
139       outList->Add(histo);
140     }
141     else isOk = kFALSE;
142   }
143   for ( Int_t icount=0; icount<kNcounts; ++icount ) {
144     for ( Int_t ich=0; ich<4; ++ich ) {
145       histoName = GetHistoName(kHslatEff, icount, ich, -1, -1);
146       histoPattern = Form("%s%s", histoName.Data(), trackSelection.Data());
147       histo = (TH1*)GetSum(physSel, trigClassNames, centrality, histoPattern);
148       if ( histo ) {
149         histo->SetName(histoName.Data());
150         histo->SetTitle(histoName.Data());
151         histo->SetDirectory(0);
152         outList->Add(histo);
153       }
154       else isOk = kFALSE;
155     }
156   }
157   for ( Int_t icount=0; icount<kNcounts; ++icount ) {
158     for ( Int_t ich=0; ich<4; ++ich ) {
159       histoName = GetHistoName(kHboardEff, icount, ich, -1, -1);
160       histoPattern = Form("%s%s", histoName.Data(), trackSelection.Data());
161       histo = (TH1*)GetSum(physSel, trigClassNames, centrality, histoPattern);
162       if ( histo ) {
163         histo->SetName(histoName.Data());
164         histo->SetTitle(histoName.Data());
165         histo->SetDirectory(0);
166         outList->Add(histo);
167       }
168       else isOk = kFALSE;
169     }
170   }
171   
172   histoName = GetHistoName(kHcheckBoard, -1, -1, -1, -1);
173   histoPattern = Form("%s%s", histoName.Data(), trackSelection.Data());
174   histo = (TH1*)GetSum(physSel, trigClassNames, centrality, histoPattern);
175   if ( histo ) {
176     histo->SetName(histoName.Data());
177     histo->SetTitle(histoName.Data());
178     histo->SetDirectory(0);
179     outList->Add(histo);
180   }
181   else isOk = kFALSE;
182   
183   return isOk;
184 }
185
186 //___________________________________________________________________________
187 void AliAnalysisTaskTrigChEff::FinishTaskOutput()
188 {
189   //
190   /// Merge Apt, Lpt and Hpt
191   /// Fill the final efficiency object (for backward compatibility)
192   //
193
194   TString histoName = "";
195   for ( Int_t isel=0; isel<kNselections; ++isel ) {
196     for ( Int_t itrig=0; itrig<GetAllSelectedTrigClasses()->GetEntries(); ++itrig ) {
197       for ( Int_t icent=1; icent<=GetCentralityClasses()->GetNbins(); ++icent ) {
198         for ( Int_t imethod=0; imethod<kNeffMethods; ++imethod ) {
199           for ( Int_t itype=0; itype<kNhistoTypes; ++itype ) {
200             for ( Int_t icount=-1; icount<kNcounts; ++icount ) {
201               for ( Int_t ich=-1; ich<4; ++ich ) {
202                 for ( Int_t imatch=kMatchApt; imatch<kMatchHpt; ++imatch ) {
203                   TH1* histo = 0x0;
204                   for ( Int_t jmatch=imatch+1; jmatch<=kMatchHpt; ++jmatch ) {
205                     histoName = GetHistoName(itype, icount, ich, jmatch, imethod);
206                     TH1* histoAdd = (TH1*)fMergeableCollection->GetObject(Form("/%s/%s/%s/",fPhysSelKeys->At(isel)->GetName(), GetAllSelectedTrigClasses()->At(itrig)->GetName(), GetCentralityClasses()->GetBinLabel(icent)), histoName);
207                     if ( ! histoAdd ) continue;
208                     histoName = GetHistoName(itype, icount, ich, imatch, imethod);
209                     if ( ! histo ) histo = (TH1*)GetMergeableObject(fPhysSelKeys->At(isel)->GetName(), GetAllSelectedTrigClasses()->At(itrig)->GetName(), GetCentralityClasses()->GetBinLabel(icent), histoName);
210                     AliDebug(2,Form("Adding %s (%g) to %s (%g)", histoAdd->GetName(), histoAdd->Integral(), histo->GetName(), histo->Integral()));
211                     histo->Add(histoAdd);
212                   } // loop on higher pt matching
213                 } // loop on match trigger
214               } // loop on chamber
215             } // loop on count type
216           } // loop on histo type
217         } // loop on eff method
218       } // loop on centrality
219     } // loop on trigger classes
220   } // loop on physics selection
221
222   TString physSel = fPhysSelKeys->At(kPhysSelPass)->GetName();
223   TString trigClass = "ANY";
224   TString centrality = "all";
225   TString histoPattern = Form("%s%s",fTrackSelKeys->At(kMatchApt)->GetName(),fEffMethodKeys->At(kEffFromTrack)->GetName());
226   
227   FillEffHistoList(physSel, trigClass, centrality, histoPattern, fList);
228
229   AliVAnalysisMuon::FinishTaskOutput();
230 }
231
232
233 //___________________________________________________________________________
234 void AliAnalysisTaskTrigChEff::InitLocalKeys()
235 {
236   //
237   /// Initialyze objects
238   //
239   
240   TString trackSelNames = "MatchNoTrackSelCut MatchApt MatchLpt MatchHpt"; 
241   fTrackSelKeys = trackSelNames.Tokenize(" ");
242   
243   TString countTypeNames = "bendPlaneCount nonBendPlaneCount bothPlanesCount allTracksCount";
244   fCountTypeKeys = countTypeNames.Tokenize(" ");
245   
246   TString histoTypeKeys = "Chamber Slat Board checkRejectedBoard";
247   fHistoTypeKeys = histoTypeKeys.Tokenize(" ");
248   
249   TString effMethodKeys = "FromTrk FromTrg";
250   fEffMethodKeys = effMethodKeys.Tokenize(" ");
251 }
252
253 //___________________________________________________________________________
254 void AliAnalysisTaskTrigChEff::MyUserCreateOutputObjects()
255 {
256   //
257   /// Create prototype objects for mergeable collection
258   //
259   
260   const Char_t* yAxisTitle = "counts";
261   
262   Int_t nChamberBins = 4;
263   Float_t chamberLow = 11.-0.5, chamberHigh = (Float_t)(nChamberBins)+11.-0.5;
264   const Char_t* chamberName = "chamber";
265   
266   Int_t nSlatBins = 18;
267   Float_t slatLow = 0.-0.5, slatHigh = (Float_t)nSlatBins-0.5;
268   const Char_t* slatName = "slat";
269   
270   Int_t nBoardBins = 234;
271   Float_t boardLow = 1.-0.5, boardHigh = (Float_t)nBoardBins+1.-0.5;
272   const Char_t* boardName = "board";
273   
274   TString baseName, histoName, histoTitle;
275   
276   TH1* histo;
277   TH2F* histo2D;
278   
279   for ( Int_t imethod=0; imethod<kNeffMethods; ++imethod ) {
280     for ( Int_t itrackSel = 0; itrackSel<kNtrackSel; ++itrackSel ) {
281       for ( Int_t icount=0; icount<kNcounts; ++icount ) {
282         histoName = GetHistoName(kHchamberEff, icount, -1, itrackSel, imethod);
283         histo = new TH1F(histoName, histoName,
284                          nChamberBins, chamberLow, chamberHigh);
285         histo->GetXaxis()->SetTitle(chamberName);
286         histo->GetYaxis()->SetTitle(yAxisTitle);
287         AddObjectToCollection(histo);
288       } // loop on counts
289       
290       for ( Int_t icount=0; icount<kNcounts; ++icount ) {
291         for ( Int_t ich=0; ich<4; ++ich ) {
292           histoName = GetHistoName(kHslatEff, icount, ich, itrackSel, imethod);
293           histo = new TH1F(histoName, histoName,
294                            nSlatBins, slatLow, slatHigh);
295           histo->GetXaxis()->SetTitle(slatName);
296           histo->GetYaxis()->SetTitle(yAxisTitle);
297           AddObjectToCollection(histo);
298         } // loop on chamber
299       } // loop on counts
300       
301       for ( Int_t icount=0; icount<kNcounts; ++icount ) {
302         for ( Int_t ich=0; ich<4; ++ich ) {
303           histoName = GetHistoName(kHboardEff, icount, ich, itrackSel, imethod);
304           histo = new TH1F(histoName, histoName,
305                            nBoardBins, boardLow, boardHigh);
306           histo->GetXaxis()->SetTitle(boardName);
307           histo->GetYaxis()->SetTitle(yAxisTitle);
308           AddObjectToCollection(histo);
309         } // loop on chamber
310       } // loop on counts
311       
312       histoName = GetHistoName(kHcheckBoard, -1, -1, itrackSel, imethod);
313       histo2D = new TH2F(histoName.Data(), "Rejected tracks motivation", 
314                          5, 20.5, 25.5, nBoardBins, boardLow, boardHigh);
315       histo2D->GetXaxis()->SetBinLabel(1,"Many pads");
316       histo2D->GetXaxis()->SetBinLabel(2,"Few pads");
317       histo2D->GetXaxis()->SetBinLabel(3,"Outside geom");
318       histo2D->GetXaxis()->SetBinLabel(4,"Tracker track");
319       histo2D->GetXaxis()->SetBinLabel(5,"Masked board");
320       histo2D->GetYaxis()->SetTitle(boardName);
321       AddObjectToCollection(histo2D);
322     } // loop on track selection
323   } // loop on eff method
324
325   fMuonTrackCuts->Print();
326   
327   fList = new TList();
328   fList->SetOwner();
329   PostData(2, fList);
330
331 }
332
333 //___________________________________________________________________________
334 TString AliAnalysisTaskTrigChEff::GetHistoName(Int_t itype, Int_t icount, Int_t ichamber, Int_t itrackSel, Int_t imethod)
335 {
336   /// Get histogram index
337   TString histoName = "";
338   if ( itype < kHcheckBoard && icount >= 0 ) histoName += ((TObjString*)fCountTypeKeys->At(icount))->GetString();
339   histoName += ((TObjString*)fHistoTypeKeys->At(itype))->GetString();
340   if ( ichamber >= 0 ) histoName += Form("Ch%i", 11+ichamber);
341   if ( itrackSel >= 0 ) histoName += ((TObjString*)fTrackSelKeys->At(itrackSel))->GetString();
342   if ( imethod >= 0 ) histoName += ((TObjString*)fEffMethodKeys->At(imethod))->GetString();
343   return histoName;
344 }
345
346 //________________________________________________________________________
347 void AliAnalysisTaskTrigChEff::ProcessEvent(TString physSel, const TObjArray& selectTrigClasses, TString centrality)
348 {
349   //
350   /// Fill histogram
351   //
352
353   Int_t slat = 0, board = 0;
354   UInt_t pattern = 0;
355   TString histoName = "";
356
357   TArrayI othersEfficient(4);
358
359   AliVParticle* track = 0x0;
360
361   Int_t itrackSel = -1;
362
363   Int_t nTracks = AliAnalysisMuonUtility::GetNTracks(InputEvent());
364   for ( Int_t itrack = 0; itrack < nTracks; ++itrack ) {
365     track = AliAnalysisMuonUtility::GetTrack(itrack,InputEvent());
366
367     Bool_t matchTracker = AliAnalysisMuonUtility::IsMuonTrack(track);
368     if ( ! matchTracker && ! fUseGhosts ) continue;
369
370     Int_t matchTrig = AliAnalysisMuonUtility::GetMatchTrigger(track);
371     itrackSel = matchTrig;
372     UInt_t selection = fMuonTrackCuts->GetSelectionMask(track);
373     Bool_t isSelected = ( ( selection & fMuonTrackCuts->GetFilterMask() ) == fMuonTrackCuts->GetFilterMask() );
374     if ( matchTrig == 1 && ( ( selection & AliMuonTrackCuts::kMuMatchApt ) == 0 ) ) isSelected = kFALSE;
375     if ( matchTrig == 2 && ( ( selection & AliMuonTrackCuts::kMuMatchLpt ) == 0 ) ) isSelected = kFALSE;
376     if ( matchTrig == 3 && ( ( selection & AliMuonTrackCuts::kMuMatchHpt ) == 0 ) ) isSelected = kFALSE;
377     if ( ! isSelected ) itrackSel = kNoSelCutApt;
378
379     for ( Int_t imethod=0; imethod<kNeffMethods; ++imethod ) {
380       if ( imethod == kEffFromTrack ) {
381         if ( track->P() < 10. ) continue;
382         pattern = AliAnalysisMuonUtility::GetMUONTrigHitsMapTrk(track);
383         board = AliESDMuonTrack::GetCrossedBoard(pattern);
384       }
385       else {
386         pattern = AliAnalysisMuonUtility::GetMUONTrigHitsMapTrg(track);
387         board = ( AliAnalysisMuonUtility::IsAODEvent(InputEvent()) ) ? AliESDMuonTrack::GetCrossedBoard(pattern) : ((AliESDMuonTrack*)track)->LoCircuit();
388       }
389       
390       Int_t effFlag = AliESDMuonTrack::GetEffFlag(pattern);
391       
392       if ( effFlag < AliESDMuonTrack::kChEff ) {
393         for ( Int_t itrig=0; itrig<selectTrigClasses.GetEntries(); ++itrig ) {
394           TString trigClassName = ((TObjString*)selectTrigClasses.At(itrig))->GetString();
395           
396           histoName = GetHistoName(kHcheckBoard, -1, -1, itrackSel, imethod);
397           ((TH2F*)GetMergeableObject(physSel, trigClassName, centrality, histoName))->Fill(AliESDMuonTrack::GetSlatOrInfo(pattern), board);
398         }
399         continue; // Track not good for efficiency calculation
400       }
401       
402       othersEfficient.Reset(1);
403       for ( Int_t cath=0; cath<2; ++cath ) {
404         for ( Int_t ich=0; ich<4; ++ich ) {
405           if ( ! AliESDMuonTrack::IsChamberHit(pattern, cath, ich) ) {
406             for ( Int_t jch=0; jch<4; jch++ ) {
407               if ( jch != ich ) othersEfficient[jch] = 0;
408             } // loop on other chambers
409             break;
410           } // if chamber not efficient
411         } // loop on chambers
412       } // loop on cathodes
413       
414       Bool_t rejectTrack = kTRUE;
415       for ( Int_t ich=0; ich<4; ++ich ) {
416         if ( othersEfficient[ich] > 0 ) {
417           rejectTrack = kFALSE;
418           break;
419         }
420       }
421       
422       if ( rejectTrack ) continue;
423       
424       slat = AliESDMuonTrack::GetSlatOrInfo(pattern);
425       
426       for ( Int_t ich=0; ich<4; ++ich ) {
427         if ( ! othersEfficient[ich] )
428           continue; // Reject track if the info of the chamber under study 
429         // is necessary to create the track itself
430         
431         Int_t iChamber = 11 + ich;
432         
433         Bool_t hitsBend = AliESDMuonTrack::IsChamberHit(pattern, 0, ich);
434         Bool_t hitsNonBend = AliESDMuonTrack::IsChamberHit(pattern, 1, ich);
435         
436         Bool_t fillHisto[kNcounts] = {
437           hitsBend,
438           hitsNonBend,
439           ( hitsBend && hitsNonBend ),
440           kTRUE
441         };
442         
443         for (Int_t icount=0; icount<kNcounts; ++icount){
444           if ( ! fillHisto[icount] ) continue;
445           for ( Int_t itrig=0; itrig<selectTrigClasses.GetEntries(); ++itrig ) {
446             TString trigClassName = ((TObjString*)selectTrigClasses.At(itrig))->GetString();
447             
448             histoName = GetHistoName(kHchamberEff, icount, -1, itrackSel, imethod);
449             ((TH1*)GetMergeableObject(physSel, trigClassName, centrality, histoName))->Fill(iChamber);
450             
451             if ( effFlag < AliESDMuonTrack::kSlatEff ) continue; // Track crossed different slats
452             histoName = GetHistoName(kHslatEff, icount, ich, itrackSel, imethod);
453             ((TH1*)GetMergeableObject(physSel, trigClassName, centrality, histoName))->Fill(slat);
454             
455             if ( effFlag < AliESDMuonTrack::kBoardEff ) continue; // Track crossed different boards
456             histoName = GetHistoName(kHboardEff, icount, ich, itrackSel, imethod);
457             ((TH1*)GetMergeableObject(physSel, trigClassName, centrality, histoName))->Fill(board);
458           } // loop on trigger classes
459         } // loop on chambers
460       } // loop on count types
461     } // loop on tracks
462   } // loop on eff methods
463   
464   PostData(2, fList);
465 }
466
467
468 //________________________________________________________________________
469 void AliAnalysisTaskTrigChEff::Terminate(Option_t *)
470 {
471   //
472   /// Draw some histogram at the end.
473   //
474
475   AliVAnalysisMuon::Terminate("");
476   
477   if ( ! fMergeableCollection ) return;
478   
479
480   Int_t xshift = 100;
481   Int_t yshift = 20;
482   Int_t igroup1 = -1;
483   Int_t igroup2 = 0;
484
485   TObjArray* physSel =  ((TObjString*)fTerminateOptions->At(0))->GetString().Tokenize(" ");
486   physSel->SetOwner();
487   TObjArray* trigClasses = ((TObjString*)fTerminateOptions->At(1))->GetString().Tokenize(" ");
488   trigClasses->SetOwner();
489   TObjArray* centrality = ((TObjString*)fTerminateOptions->At(2))->GetString().Tokenize(" ");
490   centrality->SetOwner();
491   TString furtherOpt = ((TObjString*)fTerminateOptions->At(3))->GetString();
492
493   TString currName = "";
494   TObjArray* optArr = furtherOpt.Tokenize(" ");
495   TObjArray trackSel, methodSel;
496   trackSel.SetOwner();
497   methodSel.SetOwner();
498   TString outFileOpt = "";
499   for ( Int_t iopt=0; iopt<optArr->GetEntries(); iopt++ ) {
500     currName = optArr->At(iopt)->GetName();
501     if ( currName.Contains(".root") ) outFileOpt = currName;
502     else if ( currName.Contains("Match") ) trackSel.Add(new TObjString(currName));
503     else if ( currName.Contains("From") ) methodSel.Add(new TObjString(currName));
504   }
505   delete optArr;
506   
507   if ( trackSel.GetEntries() == 0 ) trackSel.Add(new TObjString(fTrackSelKeys->At(kMatchApt)->GetName()));
508   if ( methodSel.GetEntries() == 0 ) methodSel.Add(new TObjString(fEffMethodKeys->At(kEffFromTrack)->GetName()));
509
510   furtherOpt.ToUpper();
511   
512   Int_t chosenType = ( furtherOpt.Contains("BOARD") ) ? kHboardEff : kHslatEff;
513
514   igroup1++;
515   igroup2 = 0;
516   TString histoName = "", yAxisTitle = "";
517   
518   TH1 *num = 0x0;
519   TH1 *den = 0x0;
520   TGraphAsymmErrors* effGraph = 0x0;
521   
522   ////////////////
523   // Show tests //
524   ////////////////
525   
526   for ( Int_t icount=0; icount<kNcounts-1; ++icount ) {
527     currName = Form("%s%s_can", fHistoTypeKeys->At(chosenType)->GetName(), fCountTypeKeys->At(icount)->GetName());
528     TCanvas* can = new TCanvas(currName.Data(), currName.Data(), igroup1*xshift,igroup2*yshift,600,600);
529     can->Divide(2,2);
530     TLegend* leg = new TLegend(0.6, 0.6, 0.9, 0.9);
531     leg->SetBorderSize(1);
532     for ( Int_t ich=0; ich<4; ++ich ) {
533       TGraph* refGraph = 0x0;
534       can->cd(ich+1);
535       gPad->SetRightMargin(0.03);
536       Int_t icolor = 1;
537       Int_t istyle = 0;
538       TString drawOpt = "AP";
539       for ( Int_t isel=0; isel<physSel->GetEntries(); ++isel ) {
540         for ( Int_t itrig=0; itrig<trigClasses->GetEntries(); ++itrig ) {
541           for ( Int_t icent=0; icent<centrality->GetEntries(); ++icent ) {
542             for ( Int_t imethodSel=0; imethodSel<methodSel.GetEntries(); ++imethodSel ) {
543               for ( Int_t itrackSel=0; itrackSel<trackSel.GetEntries(); ++itrackSel ) {
544                 histoName = GetHistoName(chosenType, kAllTracks, ich, -1, -1); // partial name
545                 histoName += Form("%s%s",trackSel.At(itrackSel)->GetName(),methodSel.At(imethodSel)->GetName());
546                 den = (TH1*)GetSum(physSel->At(isel)->GetName(), trigClasses->At(itrig)->GetName(), centrality->At(icent)->GetName(), histoName.Data());
547                 if ( ! den ) continue;
548                 histoName = GetHistoName(chosenType, icount, ich, -1, -1); // partial name
549                 histoName += Form("%s%s",trackSel.At(itrackSel)->GetName(),methodSel.At(imethodSel)->GetName());
550                 num = (TH1*)GetSum(physSel->At(isel)->GetName(), trigClasses->At(itrig)->GetName(), centrality->At(icent)->GetName(), histoName.Data());
551                 if ( ! num ) continue;
552                 effGraph = new TGraphAsymmErrors(num, den);
553                 currName = Form("%s_%s_%s_%s_%s", physSel->At(isel)->GetName(), trigClasses->At(itrig)->GetName(), centrality->At(icent)->GetName(), trackSel.At(itrackSel)->GetName(), methodSel.At(imethodSel)->GetName());
554                 effGraph->SetTitle(currName.Data());
555
556                 Double_t ymin = 0.;
557                 Double_t ymax = 1.1;
558                 yAxisTitle = "Efficiency";
559
560                 if ( furtherOpt.Contains("DIFF") || furtherOpt.Contains("PULL") ) {
561                   if ( ! refGraph ) {
562                     refGraph = effGraph;
563                     continue;
564                   }
565                   Double_t currX, currY, baseX, baseY, errYlow = 0., errYhigh = 0., newY = 0.;
566                   Double_t refVal = 1., errY = 0.; 
567                   for ( Int_t ipoint=0; ipoint<effGraph->GetN(); ipoint++ ) {
568                     refGraph->GetPoint(ipoint, baseX, baseY);
569                     effGraph->GetPoint(ipoint, currX, currY);
570                     Double_t errX = effGraph->GetErrorXlow(ipoint);
571                     if ( furtherOpt.Contains("DIFF") ) {
572                       refVal = ( baseY > 0. ) ? baseY : 1.;
573                       newY = ( currY - baseY ) / refVal;
574                       Double_t errYlow1 = effGraph->GetErrorYlow(ipoint);
575                       Double_t errYlow2 = refGraph->GetErrorYlow(ipoint);
576                       Double_t errYhigh1 = effGraph->GetErrorYhigh(ipoint);
577                       Double_t errYhigh2 = refGraph->GetErrorYhigh(ipoint);
578                       errYlow = TMath::Sqrt(errYlow1*errYlow1 + errYlow2*errYlow2) / refVal;
579                       errYhigh = TMath::Sqrt(errYhigh1*errYhigh1 + errYhigh2*errYhigh2) / refVal;
580                       //yAxisTitle = Form("(%s - %s) / %s", effGraph->GetTitle(), refGraph->GetTitle(), refGraph->GetTitle());
581                       yAxisTitle = "(eff - ref ) / ref";
582                       effGraph->SetTitle(Form("Rel. diff. w.r.t. %s", refGraph->GetTitle()));
583                       ymin = -0.1;
584                       ymax = 0.1;
585                     }
586                     else if ( furtherOpt.Contains("PULL") ) {
587                       errY = 0.5 * ( effGraph->GetErrorYlow(ipoint) + effGraph->GetErrorYhigh(ipoint));
588                       newY = ( errY > 0. ) ? ( currY - baseY ) / errY : 0.;
589                       errYlow = 1.;
590                       errYhigh = 1.;
591                       //yAxisTitle = Form("( %s - %s ) / err", effGraph->GetTitle(), refGraph->GetTitle());
592                       yAxisTitle = "(eff - ref ) / err";
593                       effGraph->SetTitle(Form("Pull w.r.t. %s", refGraph->GetTitle()));
594                       ymin = -4.;
595                       ymax = 4.;
596                     }
597                     effGraph->SetPoint(ipoint, currX, newY);
598                     effGraph->SetPointError(ipoint, errX, errX, errYlow, errYhigh);
599                   } // loop on points
600                 }
601                 effGraph->GetYaxis()->SetRangeUser(ymin, ymax);
602                 effGraph->GetYaxis()->SetTitle(yAxisTitle.Data());
603                 effGraph->SetLineColor(icolor);
604                 effGraph->SetMarkerColor(icolor);
605                 effGraph->SetMarkerStyle(20 + istyle);
606                 //effGraph->SetMarkerSize(0.3);
607                 icolor++;
608                 if ( icolor == 5 || icolor == 10 ) icolor++;
609                 istyle++;
610                 effGraph->Draw(drawOpt.Data());
611                 drawOpt = "P";
612                 if ( ich < 3 ) continue;
613                 leg->AddEntry(effGraph, currName.Data(), "lp");
614               } // loop on match trigger
615             } // loop on eff method
616           } // loop on centrality
617         } // loop on trigger classes
618       } // loop on physics selection
619     } // loop on chamber
620     leg->Draw("same");
621   } // loop on count type
622   //} // loop on detection element type
623
624   delete physSel;
625   delete trigClasses;
626   delete centrality;
627   
628    
629   fList = dynamic_cast<TList*>(GetOutputData(2));
630   if ( fList ) {
631   
632     ///////////////////////////
633     // Show final efficiency //
634     ///////////////////////////
635     TString baseName[3] = {"Chamber", "RPC", "Board"};
636     Int_t baseIndex[3] = {kHchamberEff, kHslatEff, kHboardEff};
637     TString effName[kNcounts-1] = {"BendPlane", "NonBendPlane", "BothPlanes"};
638     for ( Int_t itype=0; itype<3; itype++ ) {
639       for ( Int_t icount=0; icount<kNcounts-1; icount++ ){
640         TString canName = Form("efficiencyPer%s_%s",baseName[itype].Data(),effName[icount].Data());
641         TCanvas* can = new TCanvas(canName.Data(),canName.Data(),10*(1+kNcounts*itype+icount),10*(1+kNcounts*itype+icount),310,310);
642         can->SetFillColor(10); can->SetHighLightColor(10);
643         can->SetLeftMargin(0.15); can->SetBottomMargin(0.15);  
644         if ( itype > 0 )
645           can->Divide(2,2);
646         
647         for ( Int_t ich=-1; ich<4; ich++ ) {
648           histoName = GetHistoName(baseIndex[itype], icount, ich, -1, -1);
649           num = (TH1*)fList->FindObject(histoName.Data());
650           histoName = GetHistoName(baseIndex[itype], kNcounts-1, ich, -1, -1);
651           den = (TH1*)fList->FindObject(histoName.Data());
652           if ( ! num || ! den ) continue;
653           effGraph = new TGraphAsymmErrors(num, den);
654           effGraph->GetYaxis()->SetRangeUser(0., 1.1);
655           effGraph->GetYaxis()->SetTitle("Efficiency");
656           effGraph->GetXaxis()->SetTitle(baseName[itype].Data());
657           can->cd(ich+1);
658           effGraph->Draw("AP");
659           if ( itype == 0 ) break;
660         } // loop on chamber
661       } // loop on count types
662     } // loop on histo
663   }
664
665   
666   if ( ! outFileOpt.IsNull() ) {
667     TObjArray* outFileOptList = outFileOpt.Tokenize("?");
668     AliInfo(Form("Creating file %s", outFileOptList->At(0)->GetName()));
669     TList* effList = GetEffHistoList(outFileOptList->At(1)->GetName(), outFileOptList->At(2)->GetName(), outFileOptList->At(3)->GetName(), outFileOptList->At(4)->GetName());
670     effList->SetName(GetOutputSlot(2)->GetContainer()->GetName());
671     TFile* outFile = TFile::Open(outFileOptList->At(0)->GetName(), "RECREATE");
672     effList->Write(effList->GetName(),TObject::kSingleKey);
673     outFile->Close();
674     delete effList;
675     delete outFileOptList;
676   }
677 }