]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONTriggerEfficiencyCells.cxx
Changes for removal of AliMpManuList (Laurent)
[u/mrichter/AliRoot.git] / MUON / AliMUONTriggerEfficiencyCells.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 #include "AliMUONTriggerEfficiencyCells.h"
17 #include "AliMpConstants.h"
18 #include "AliMpDEManager.h"
19
20 #include "AliLog.h"
21
22 #include "TRandom.h"
23 #include "Riostream.h"
24 #include "TSystem.h"
25 #include "TFile.h"
26 #include "TH1F.h"
27 #include "TMath.h"
28
29 #include "TH2F.h"
30 #include "TStyle.h"
31 #include "TPaveLabel.h"
32 #include "TCanvas.h"
33
34 #include <fstream>
35 #include <cassert>
36
37 //-----------------------------------------------------------------------------
38 /// \class AliMUONTriggerEfficiencyCells
39 /// A class to store and give access to the trigger chamber efficiency.
40 ///
41 /// Efficiency is stored per cathode on local boards, or, alternatively,
42 /// on "cells" of a given size.
43 ///
44 /// The main method of this class is IsTriggered().
45 ///
46 /// $ALICE_ROOT/MUON/data/TriggerChamberefficiencyCells.dat contains efficiency 
47 /// for each chamber (i.e. DetElement). 
48 ///
49 /// In the case of local boards, efficiency is stored from left to right
50 /// per increasing board number (from 1 to 234)
51 ///
52 /// Otherwise, he efficiency cells goes from right to left and 
53 /// from bottom to top of the chamber, namely, the efficiencies tabulated in the 
54 /// file refers to the following reference frame:
55 ///
56 /// <pre>
57 /// x
58 /// <----------------------------------|
59 ///                                    |
60 ///    ---------------------------     |
61 ///   | 0.97 | 0.97 | 0.97 | 0.97 |    |
62 ///    ---------------------------     |
63 ///   | 0.97 | 0.97 | 0.97 | 0.97 |    |
64 ///    ---------------------------     |
65 ///   | 0.97 | 0.97 | 0.97 | 0.97 |    |
66 ///    ---------------------------     |
67 ///                                    |
68 ///                                   \/ y
69 /// </pre>
70 ///
71 ///  In both cases, the file can be edited in order to change efficiency
72 ///  in a chosen local board/region of the chamber.
73 ///
74 ///
75 /// But please note that this object is also available from the CDB 
76 /// (generated using the MUONCDB.C macro)
77 ///
78 /// \author Diego Stocco; INFN Torino
79 //-----------------------------------------------------------------------------
80
81 /// \cond CLASSIMP
82 ClassImp(AliMUONTriggerEfficiencyCells)
83 /// \endcond
84
85 //__________________________________________________________________________
86 AliMUONTriggerEfficiencyCells::AliMUONTriggerEfficiencyCells()
87 :
88 TObject(),
89 fCountHistoList(0x0),
90 fNoCountHistoList(0x0)
91 {
92 ///  Default constructor.
93   CheckConstants();
94   Reset();
95   InitHistos();
96 }
97
98 //__________________________________________________________________________
99 AliMUONTriggerEfficiencyCells::AliMUONTriggerEfficiencyCells(const char* filename)
100 :
101 TObject(),
102 fCountHistoList(0x0),
103 fNoCountHistoList(0x0)
104 {
105 ///  Constructor using an ASCII file.
106   CheckConstants();
107   Reset();
108   ReadFile(filename);
109 }
110
111 AliMUONTriggerEfficiencyCells::AliMUONTriggerEfficiencyCells(TList *countHistoList,
112                                                              TList *noCountHistoList)
113 :
114 TObject(),
115 fCountHistoList(countHistoList),
116 fNoCountHistoList(noCountHistoList)
117 {
118 ///  Constructor using an ASCII file.
119   CheckConstants();
120   Reset();
121   InitHistos();
122   FillHistosFromList();
123 }
124
125 //_____________________________________________________________________________
126 AliMUONTriggerEfficiencyCells::AliMUONTriggerEfficiencyCells(const AliMUONTriggerEfficiencyCells& other)
127 :
128 TObject(other),
129 fCountHistoList(other.fCountHistoList),
130 fNoCountHistoList(other.fNoCountHistoList)
131 {
132 /// Copy constructor
133
134   for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
135     for(Int_t slat=0; slat<fgkNslats; slat++){
136       fCellContent[chCath][slat] = other.fCellContent[chCath][slat];
137     }
138     fCellSize[chCath] = other.fCellSize[chCath];
139     fCellNumber[chCath] = other.fCellNumber[chCath];
140     fBoardEfficiency[chCath] = other.fBoardEfficiency[chCath];
141     fSlatEfficiency[chCath] = other.fSlatEfficiency[chCath];
142   }
143 }
144
145 //_____________________________________________________________________________
146 AliMUONTriggerEfficiencyCells& AliMUONTriggerEfficiencyCells::operator=(const AliMUONTriggerEfficiencyCells& other)
147 {
148   /// Asignment operator
149   // check assignement to self
150   if (this == &other)
151     return *this;
152
153   for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
154     for(Int_t slat=0; slat<fgkNslats; slat++){
155       fCellContent[chCath][slat] = other.fCellContent[chCath][slat];
156     }
157     fCellSize[chCath] = other.fCellSize[chCath];
158     fCellNumber[chCath] = other.fCellNumber[chCath];
159     fBoardEfficiency[chCath] = other.fBoardEfficiency[chCath];
160     fSlatEfficiency[chCath] = other.fSlatEfficiency[chCath];
161   }
162
163   fCountHistoList = other.fCountHistoList;
164   fNoCountHistoList = other.fNoCountHistoList;
165
166   return *this;
167 }
168
169 //__________________________________________________________________________
170 AliMUONTriggerEfficiencyCells::~AliMUONTriggerEfficiencyCells()
171 {
172 ///  Destructor.
173 }
174
175
176 //__________________________________________________________________________
177 void AliMUONTriggerEfficiencyCells::GetCellEfficiency(Int_t detElemId, Float_t x, Float_t y, Float_t &eff1, Float_t &eff2) const
178 {
179 ///  Get the efficiencies of the 2 cathodes at a given location (x,y)
180
181   Int_t chamber = FindChamberIndex(detElemId);
182   Int_t slat = FindSlatIndex(detElemId);
183   TArrayI cell = CellByCoord(detElemId,x,y);
184   eff1 = 0.0;
185   eff2 = 0.0;
186   if(cell.At(0)>=0 && cell.At(1)>=0)
187   {
188     eff1 = fCellContent[chamber][slat][cell.At(0)][cell.At(1)];
189     eff2 = fCellContent[fgkNchambers+chamber][slat][cell.At(0)][cell.At(1)];
190   }
191 }
192
193
194 //__________________________________________________________________________
195 void AliMUONTriggerEfficiencyCells::GetCellEfficiency(Int_t detElemId, Int_t localBoard, Float_t &eff1, Float_t &eff2) const
196 {
197 ///  Get the efficiencies of the 2 cathodes at a given local board
198
199   Int_t chamber = FindChamberIndex(detElemId);
200   Int_t bin = fBoardEfficiency[chamber]->FindBin(localBoard);
201   eff1 = fBoardEfficiency[chamber]->GetBinContent(bin);
202   eff2 = fBoardEfficiency[fgkNchambers+chamber]->GetBinContent(bin);
203 }
204
205
206 //__________________________________________________________________________
207 void 
208 AliMUONTriggerEfficiencyCells::IsTriggered(Int_t detElemId, Float_t x, Float_t y, Bool_t &trig1, Bool_t &trig2) const
209 {
210 ///  Whether or not a given location (x,y) has a chance to trig, on each cathode.
211
212   Float_t eff1 = 0.0;
213   Float_t eff2 = 0.0;
214   GetCellEfficiency(detElemId, x, y, eff1, eff2);
215   trig1 = kTRUE; 
216   trig2 = kTRUE;
217   if(gRandom->Rndm()>eff1)trig1 = kFALSE;
218   if(gRandom->Rndm()>eff2)trig2 = kFALSE;
219 }
220
221
222 //__________________________________________________________________________
223 void 
224 AliMUONTriggerEfficiencyCells::IsTriggered(Int_t detElemId, Int_t localBoard, Bool_t &trig1, Bool_t &trig2) const
225 {
226 ///  Whether or not a given local board has a chance to trig, on each cathode.
227
228   Float_t eff1 = 0.0;
229   Float_t eff2 = 0.0;
230   GetCellEfficiency(detElemId, localBoard, eff1, eff2);
231   trig1 = kTRUE; 
232   trig2 = kTRUE;
233   if(gRandom->Rndm()>eff1)trig1 = kFALSE;
234   if(gRandom->Rndm()>eff2)trig2 = kFALSE;
235 }
236
237 //__________________________________________________________________________
238 TArrayI AliMUONTriggerEfficiencyCells::CellByCoord(Int_t detElemId, Float_t x, Float_t y) const
239 {
240 ///  Get the efficiencies at a given location.
241
242   Int_t chamber = FindChamberIndex(detElemId);
243   Int_t slat = FindSlatIndex(detElemId);
244   Int_t cell[fgkNcathodes]={-1,-1};
245   Float_t maxX = fCellSize[chamber][slat]*((Float_t)fCellNumber[chamber][slat]);
246   Float_t maxY = fCellSize[fgkNchambers+chamber][slat]*((Float_t)fCellNumber[fgkNchambers+chamber][slat]);
247   if(x>=0 & x<maxX & y>=0 & y<maxY)
248   {
249     cell[0] = (Int_t)(x/fCellSize[chamber][slat]);
250     cell[1] = (Int_t)(y/fCellSize[fgkNchambers+chamber][slat]);
251   }
252   return TArrayI(fgkNcathodes,cell);
253 }
254
255 //__________________________________________________________________________
256 void AliMUONTriggerEfficiencyCells::ReadFile(const char* filename)
257 {
258 ///  Reads a file containing the efficiency map.
259
260   TString fileName = gSystem->ExpandPathName(filename);
261   if(fileName.EndsWith(".root")){
262       ReadHistoBoards(fileName.Data());
263       return;
264   }
265
266   InitHistos();
267   ifstream file(fileName.Data());
268   char dat[50];
269   if (file.good()){
270       file >> dat;
271       if(!strcmp(dat,"localBoards"))ReadFileBoards(file);
272       else ReadFileXY(file);
273       file.close();
274   } else {
275       AliWarning(Form("Can't read file %s",fileName.Data()));
276   }
277 }
278
279
280 //__________________________________________________________________________
281 void AliMUONTriggerEfficiencyCells::ReadFileXY(ifstream &file)
282 {
283 ///  Structure of file (.dat) containing geometrical efficency
284     Int_t datInt=0, detEl=0, chamber=0, rpc=0, chCath=0;
285     Float_t datFloat=0.0;
286     char dat[50];
287
288     while (file >> dat) {
289         file >> detEl;
290         chamber = FindChamberIndex(detEl);
291         rpc = FindSlatIndex(detEl);
292         file >> dat;
293         for(Int_t i=0; i<fgkNcathodes; i++){
294             chCath = fgkNchambers*i + chamber;
295             file >> datInt;
296             fCellNumber[chCath][rpc] = datInt;
297             file >> dat;
298         }
299         for(Int_t i=0; i<fgkNcathodes; i++){
300             chCath = fgkNchambers*i + chamber;
301             file >> datFloat;
302             fCellSize[chCath][rpc] = datFloat;
303             if(i==0)file >> dat;
304         }
305         for(Int_t cath=0; cath<fgkNcathodes; cath++){
306             chCath = fgkNchambers*cath + chamber;
307             file >> dat;
308             file >> datInt;
309             for(Int_t iy=0; iy<fCellNumber[fgkNchambers+chamber][rpc]; iy++){
310                 for(Int_t ix=0; ix<fCellNumber[chamber][rpc]; ix++){
311                     file >> datFloat;
312                     fCellContent[chCath][rpc][ix][iy] = datFloat;
313                 }
314             }
315         }
316     }
317 }
318
319 //__________________________________________________________________________
320 void AliMUONTriggerEfficiencyCells::ReadFileBoards(ifstream &file)
321 {
322 ///  Structure of file (.dat) containing local board efficency
323   Int_t datInt=0, detEl=0, chamber=0, chCath=0, bin=0;
324     Float_t datFloat=0.0;
325     char dat[50];
326
327     while (file >> dat) {
328             file >> detEl;
329             chamber = FindChamberIndex(detEl);
330             //rpc = FindSlatIndex(detEl);
331             for(Int_t cath=0; cath<fgkNcathodes; cath++){
332                 chCath = fgkNchambers*cath + chamber;
333                 file >> dat;
334                 file >> datInt;
335                 for(Int_t board=1; board<=AliMpConstants::NofLocalBoards(); board++){
336                     file >> datFloat;
337                     bin = fBoardEfficiency[chCath]->FindBin(board);
338                     fBoardEfficiency[chCath]->SetBinContent(bin, datFloat);
339                 }
340             }
341     }
342 }
343
344
345 //__________________________________________________________________________
346 void AliMUONTriggerEfficiencyCells::ReadHistoBoards(const char *filename)
347 {
348 ///  Structure of file (.root) containing local board efficency
349     TFile *file = new TFile(filename, "read");
350     if(!file) {
351         AliWarning(Form("Can't read file %s",filename));
352         return;
353     }
354     char histoName[30];
355     char *cathCode[fgkNcathodes] = {"bendPlane", "nonBendPlane"};
356
357     for(Int_t ch=0; ch<fgkNchambers; ch++){
358         for(Int_t cath=0; cath<fgkNcathodes; cath++){
359             sprintf(histoName, "%sBoardEffChamber%i", cathCode[cath], 11+ch);
360             if(!(TH1F *)file->Get(histoName)) {
361                 AliWarning(Form("Can't find histo %s in file %s",histoName, filename));
362                 continue;
363             }
364             Int_t chCath = fgkNchambers*cath + ch;
365             fBoardEfficiency[chCath] = (TH1F *)file->Get(histoName);
366         }
367     }
368 }
369
370
371 //_____________________________________________________________________________
372 void AliMUONTriggerEfficiencyCells::CheckConstants() const
373 {
374 /// Check consistence of redefined constants 
375
376   assert(fgkNcathodes == AliMpConstants::NofCathodes());    
377   assert(fgkNchambers == AliMpConstants::NofTriggerChambers());    
378   assert(fgkNplanes == AliMpConstants::NofTriggerChambers() * fgkNcathodes);    
379 }
380
381
382 //__________________________________________________________________________
383 Int_t AliMUONTriggerEfficiencyCells::FindChamberIndex(Int_t detElemId) const
384 {
385 ///  From detElemId to chamber number
386   Int_t iChamber = AliMpDEManager::GetChamberId(detElemId);
387   Int_t chamber = iChamber-AliMpConstants::NofTrackingChambers();
388   return chamber;
389 }
390
391
392 //__________________________________________________________________________
393 Int_t AliMUONTriggerEfficiencyCells::FindSlatIndex(Int_t detElemId) const
394 {
395 ///  From detElemId to slat index.
396   Int_t slat = detElemId%100;
397   return slat;
398 }
399
400
401 //__________________________________________________________________________
402 TVector2 AliMUONTriggerEfficiencyCells::ChangeReferenceFrame(Float_t x, Float_t y, Float_t x0, Float_t y0)
403 {
404 /// (x0,y0) position of the local reference frame (center of the chamber)
405
406     Float_t x1 = x0-x;//reflection of axis
407     Float_t y1 = y+y0;
408     return TVector2(x1,y1);
409 }
410
411 //__________________________________________________________________________
412 void
413 AliMUONTriggerEfficiencyCells::Reset()
414 {
415 ///  Sets our internal array contents to zero.
416
417     for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
418         fCellSize[chCath].Set(fgkNslats);
419         fCellNumber[chCath].Set(fgkNslats);
420         for(Int_t slat=0; slat<fgkNslats; slat++){
421             fCellContent[chCath][slat].ResizeTo(fgkNcells,fgkNcells);
422         }
423     }
424
425     for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
426         fCellSize[chCath].Reset();
427         fCellNumber[chCath].Reset();
428         for(Int_t slat=0; slat<fgkNslats; slat++){
429             fCellContent[chCath][slat].Zero();
430         }
431     }
432
433     for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
434       fBoardEfficiency[chCath] = 0x0;
435       fSlatEfficiency[chCath] = 0x0;
436     }
437     
438 }
439
440
441 //__________________________________________________________________________
442 void
443 AliMUONTriggerEfficiencyCells::InitHistos()
444 {
445 ///  Sets our internal array contents to zero.
446
447   const Int_t kNumOfBoards = AliMpConstants::NofLocalBoards();
448   Int_t chCath=0;
449   Char_t histoName[40];
450
451   Char_t *cathCode[fgkNcathodes] = {"bendPlane", "nonBendPlane"};
452
453   for(Int_t ch=0; ch<fgkNchambers; ch++){
454     for(Int_t cath=0; cath<fgkNcathodes; cath++){
455       chCath = fgkNchambers*cath + ch;
456       sprintf(histoName, "%sBoardEffChamber%i", cathCode[cath], 11+ch);
457       fBoardEfficiency[chCath] = new TH1F(histoName, histoName, kNumOfBoards, 1-0.5, kNumOfBoards+1.-0.5);
458       sprintf(histoName, "%sSlatEffChamber%i", cathCode[cath], 11+ch);
459       fSlatEfficiency[chCath] = new TH1F(histoName, histoName, fgkNslats, 0-0.5, fgkNslats-0.5);
460     }
461   }
462 }
463
464
465 //__________________________________________________________________________
466 void
467 AliMUONTriggerEfficiencyCells::FillHistosFromList()
468 {
469 ///  Sets our internal array contents to zero.
470
471   Int_t nHistoBins=0;
472   TH1F *histoNum = 0x0, *histoDen=0x0, *currHisto = 0x0;
473   TString slatName = "Slat", boardName = "Board", histoName;
474   Int_t iHistoBoard = -1, iHistoSlat = -1;
475   Float_t efficiency, efficiencyError;
476
477   Int_t nentries = fCountHistoList->GetEntries();
478
479   for(Int_t iEntry=0; iEntry<nentries; iEntry++){
480     histoNum = (TH1F*)fCountHistoList->At(iEntry);
481     histoDen = (TH1F*)fNoCountHistoList->At(iEntry);
482
483     if(!histoNum) {
484       AliWarning("Histogram not found in fCountHistoList. Skip to next");
485       continue;
486     }
487     if(!histoDen) {
488       AliWarning("Histogram not found in fNoCountHistoList. Skip to next");
489       continue;
490     }
491
492     histoName = histoNum->GetName();
493     nHistoBins = histoNum->GetNbinsX();
494
495     if(histoName.Contains(boardName)){
496       iHistoBoard++;
497       currHisto = fBoardEfficiency[iHistoBoard];
498     }
499     else if(histoName.Contains(slatName)){
500       iHistoSlat++;
501       currHisto = fSlatEfficiency[iHistoSlat];
502     }
503     else continue;
504
505     for(Int_t iBin=1; iBin<=nHistoBins; iBin++){
506       CalculateEfficiency((Int_t)histoNum->GetBinContent(iBin), (Int_t)histoNum->GetBinContent(iBin) + (Int_t)histoDen->GetBinContent(iBin), efficiency, efficiencyError, kFALSE);
507
508       currHisto->SetBinContent(iBin, efficiency);
509       currHisto->SetBinError(iBin, efficiencyError);
510     }
511   }
512 }
513
514
515
516 //_____________________________________________________________________________
517 void AliMUONTriggerEfficiencyCells::CalculateEfficiency(Int_t trigger44, Int_t trigger34,
518                                                         Float_t &efficiency, Float_t &error,
519                                                         Bool_t failuresAsInput)
520 {
521     //
522     /// Returns the efficiency.
523     //
524
525     efficiency=-9.;
526     error=0.;
527     if(trigger34>0){
528         efficiency=(Double_t)trigger44/((Double_t)trigger34);
529         if(failuresAsInput)efficiency=1.-(Double_t)trigger44/((Double_t)trigger34);
530     }
531     Double_t q = TMath::Abs(1-efficiency);
532     if(efficiency<0)error=0.0;
533     else error = TMath::Sqrt(efficiency*q/((Double_t)trigger34));
534 }
535
536 //_____________________________________________________________________________
537 void AliMUONTriggerEfficiencyCells::DisplayEfficiency(Bool_t perSlat)
538 {
539   //
540   /// Display calculated efficiency.
541   //
542
543   const Int_t kNumOfBoards = AliMpConstants::NofLocalBoards();
544
545   Bool_t isInitSlat = kFALSE, isInitBoard = kFALSE;
546   for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
547     if(fBoardEfficiency[chCath]->GetEntries()>0) isInitBoard = kTRUE;
548     if(fSlatEfficiency[chCath]->GetEntries()>0) isInitSlat = kTRUE;
549   }
550
551   if(!isInitBoard){
552     printf("Trigger efficiency not initialized per board.\nDisplay not yet implemented.\n");
553     return;
554   }
555   if(!isInitSlat && perSlat){
556     printf("Trigger efficiency not initialized for slat.\nPlease try option kFALSE.\n");
557     return;
558   }
559
560   Int_t side, col, line, nbx, slat;
561   Float_t xCenter, yCenter, zCenter, xWidth, yWidth;
562   Float_t x1Label, x2Label, y1Label, y2Label;
563   Int_t x1, y1, x2, y2, board=0;
564   Char_t name[8], text[200];
565
566   gStyle->SetPalette(1);
567
568   Char_t *cathCode[fgkNcathodes] = {"bendPlane", "nonBendPlane"};
569
570   Float_t boardsX = 257.00;  // cm
571   Float_t boardsY = 307.00;  // cm
572
573   TH2F *histo[fgkNplanes];
574   TPaveLabel *boardLabel[fgkNplanes][234];
575   assert(kNumOfBoards==234);    
576   TArrayI boardsPerColumn[9];
577   for(Int_t iLine=0; iLine<9; iLine++){
578     boardsPerColumn[iLine].Set(7);
579     boardsPerColumn[iLine].Reset();
580   }
581
582   Char_t histoName[40], histoTitle[90], labelTxt[5];
583
584   Float_t efficiency, efficiencyError;
585
586   for(Int_t cath=0; cath<fgkNcathodes; cath++){
587     for(Int_t ch=0; ch<fgkNchambers; ch++){
588       Int_t chCath = fgkNchambers*cath + ch;
589       sprintf(histoName, "%sChamber%i", cathCode[cath], 11+ch);
590       sprintf(histoTitle, "Chamber %i: efficiency %s", 11+ch, cathCode[cath]);
591       histo[chCath] = new TH2F(histoName, histoTitle, (Int_t)boardsX, -boardsX, boardsX, (Int_t)boardsY, -boardsY, boardsY);
592       histo[chCath]->SetXTitle("X (cm)");
593       histo[chCath]->SetYTitle("Y (cm)");
594     }
595   }
596
597   TString mapspath = gSystem->Getenv("ALICE_ROOT");
598   mapspath.Append("/MUON/data");
599
600   sprintf(text,"%s/guimapp11.txt",mapspath.Data());
601
602   FileStat_t fs;
603   if(gSystem->GetPathInfo(text,fs)){
604     AliWarning(Form("Map file %s not found. Nothing done",text));
605     return;
606   }
607
608   FILE *fmap = fopen(text,"r");
609
610   for (Int_t ib = 0; ib < kNumOfBoards; ib++) {
611     fscanf(fmap,"%d   %d   %d   %d   %f   %f   %f   %f   %f   %s   \n",&side,&col,&line,&nbx,&xCenter,&yCenter,&xWidth,&yWidth,&zCenter,&name[0]);
612     if(side==0) continue;
613     boardsPerColumn[line-1][col-1]++;
614   }
615
616   rewind(fmap);
617
618   for (Int_t ib = 0; ib < kNumOfBoards; ib++) {
619     fscanf(fmap,"%d   %d   %d   %d   %f   %f   %f   %f   %f   %s   \n",&side,&col,&line,&nbx,&xCenter,&yCenter,&xWidth,&yWidth,&zCenter,&name[0]);
620
621     board=0;
622     for(Int_t iCol=1; iCol<=col; iCol++){
623       Int_t lastLine = 9;
624       if(iCol==col) lastLine = line-1;
625       for(Int_t iLine=1; iLine<=lastLine; iLine++){
626         board += boardsPerColumn[iLine-1][iCol-1];
627       }
628     }
629     if(side==0) board += kNumOfBoards/2;
630     board += nbx - 1;
631
632     slat = (line+13)%fgkNslats;
633     if(side==0) slat = 14-line;
634
635     for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
636       x1 = histo[chCath]->GetXaxis()->FindBin(xCenter-xWidth/2.)+1;
637       y1 = histo[chCath]->GetYaxis()->FindBin(yCenter-yWidth/2.)+1;
638       x2 = histo[chCath]->GetXaxis()->FindBin(xCenter+xWidth/2.)-1;
639       y2 = histo[chCath]->GetYaxis()->FindBin(yCenter+yWidth/2.)-1;
640
641       x1Label = xCenter-xWidth/2.;
642       y1Label = yCenter-yWidth/2.;
643       x2Label = xCenter+xWidth/2.;
644       y2Label =  yCenter+yWidth/2.;
645       sprintf(labelTxt,"%3d", board+1);
646       if(perSlat){
647         x1Label = 140.;
648         x2Label = x1Label + 40.;
649         y1Label = -285. + ((Float_t)(line - 1)) * 68;
650         y2Label = y1Label + 34.;
651         if(side==0){
652           x1Label = -x2Label;
653           x2Label = x1Label + 40.;
654         }
655         sprintf(labelTxt,"%2d", slat);
656       }
657     
658       boardLabel[chCath][board] = new TPaveLabel(x1Label, y1Label, x2Label, y2Label, labelTxt);
659       boardLabel[chCath][board]->SetFillStyle(0);
660       boardLabel[chCath][board]->SetBorderSize(0);
661
662       Int_t histoBin = board+1;
663       if(perSlat) histoBin = slat+1;
664
665       if(!perSlat){
666         efficiency = fBoardEfficiency[chCath]->GetBinContent(histoBin);
667         efficiencyError = fBoardEfficiency[chCath]->GetBinError(histoBin);
668       }
669       else {
670         efficiency = fSlatEfficiency[chCath]->GetBinContent(histoBin);
671         efficiencyError = fSlatEfficiency[chCath]->GetBinError(histoBin);
672       }
673
674       for(Int_t binX=x1; binX<=x2; binX++){
675         for(Int_t binY=y1; binY<=y2; binY++){
676           histo[chCath]->SetBinContent(binX, binY, efficiency);
677           histo[chCath]->SetBinError(binX, binY, efficiencyError);
678         }
679       }
680     }
681   }
682
683   fclose(fmap);
684
685   TCanvas *can[8];
686   for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
687     sprintf(histoName, "%sCan", histo[chCath]->GetName());
688     sprintf(histoTitle, "%s", histo[chCath]->GetTitle());
689     can[chCath] = new TCanvas(histoName, histoTitle, 100+10*chCath, 10*chCath, 700, 700);
690     can[chCath]->SetRightMargin(0.14);
691     can[chCath]->SetLeftMargin(0.12);
692     histo[chCath]->GetZaxis()->SetRangeUser(0.,1.);
693     histo[chCath]->GetYaxis()->SetTitleOffset(1.4);
694     histo[chCath]->SetStats(kFALSE);
695     histo[chCath]->Draw("COLZ");
696     for (Int_t board = 0; board < kNumOfBoards; board++) {
697       boardLabel[chCath][board]->Draw("same");
698     }
699   }
700 }
701
702
703 //__________________________________________________________________________
704 Bool_t AliMUONTriggerEfficiencyCells::SumRunEfficiency(const AliMUONTriggerEfficiencyCells &other)
705 {
706 ///  Sums results from different runs and gives the efficiency
707   if(!fCountHistoList || !fNoCountHistoList) {
708     AliWarning("Histograms for efficiency calculations not implemented in object");
709     return kFALSE;
710   }
711   if(!other.fCountHistoList || !other.fNoCountHistoList) {
712     AliWarning("Histograms for efficiency calculations not implemented in object passed as argument");
713     return kFALSE;
714   }
715
716   Int_t nentries = fCountHistoList->GetEntries();
717   TH1F *currNum = 0x0, *currDen = 0x0, *otherNum = 0x0, *otherDen = 0x0;
718
719   for(Int_t iEntry=0; iEntry<nentries; iEntry++){
720     currNum = (TH1F*)fCountHistoList->At(iEntry);
721     currDen = (TH1F*)fNoCountHistoList->At(iEntry);
722     otherNum = (TH1F*)other.fCountHistoList->At(iEntry);
723     otherDen = (TH1F*)other.fNoCountHistoList->At(iEntry);
724     currNum->Add(otherNum);
725     currDen->Add(otherDen);
726   }
727
728   FillHistosFromList();
729   return kTRUE;
730 }