]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONTriggerEfficiencyCells.cxx
In AliMUONCDB:
[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 // $Id$
17
18 #include "AliMUONTriggerEfficiencyCells.h"
19 #include "AliMpConstants.h"
20
21 // Classes for display
22 #include "AliMUONTriggerDisplay.h"
23 #include "AliCDBManager.h"
24 #include "AliMpCDB.h"
25 #include "AliMpDDLStore.h"
26 #include "AliMpLocalBoard.h"
27 #include "AliMpPad.h"
28 #include "AliMpVSegmentation.h"
29 #include "AliMpSegmentation.h"
30
31 #include "AliLog.h"
32
33 #include "TRandom.h"
34 #include "Riostream.h"
35 #include "TSystem.h"
36 #include "TFile.h"
37 #include "TH1F.h"
38 #include "TMath.h"
39
40 #include "TH2F.h"
41 #include "TH3F.h"
42 #include "TF1.h"
43 #include "TStyle.h"
44 #include "TPaveLabel.h"
45 #include "TCanvas.h"
46
47 #include <fstream>
48 #include <cassert>
49
50 //-----------------------------------------------------------------------------
51 /// \class AliMUONTriggerEfficiencyCells
52 /// A class to store and give access to the trigger chamber efficiency.
53 ///
54 /// Efficiency is stored per cathode on local boards
55 ///
56 /// The main method of this class is IsTriggered().
57 ///
58 /// $ALICE_ROOT/MUON/data/efficiencyCells.dat contains efficiency 
59 /// for each chamber (i.e. DetElement). 
60 ///
61 /// In the case of local boards, efficiency is stored from left to right
62 /// per increasing board number (from 1 to 234)
63 ///
64 /// The file can be edited in order to change efficiency
65 /// in a chosen local board/region of the chamber.
66 ///
67 ///
68 /// But please note that this object is also available from the CDB 
69 ///
70 /// \author Diego Stocco; INFN Torino
71 //-----------------------------------------------------------------------------
72
73 /// \cond CLASSIMP
74 ClassImp(AliMUONTriggerEfficiencyCells)
75 /// \endcond
76
77 //__________________________________________________________________________
78 AliMUONTriggerEfficiencyCells::AliMUONTriggerEfficiencyCells()
79 :
80 TObject(),
81 fCountHistoList(0x0),
82 fNoCountHistoList(0x0),
83 fFiredStrips(0x0),
84 fDisplayHistoList(0x0),
85 fBoardLabelList(0x0),
86 fFiredFitHistoList(0x0),
87 fFiredDisplayHistoList(0x0)
88 {
89 ///  Default constructor.
90   CheckConstants();
91   Reset();
92   InitHistos();
93 }
94
95 //__________________________________________________________________________
96 AliMUONTriggerEfficiencyCells::AliMUONTriggerEfficiencyCells(const Char_t* filename)
97 :
98 TObject(),
99 fCountHistoList(0x0),
100 fNoCountHistoList(0x0),
101 fFiredStrips(0x0),
102 fDisplayHistoList(0x0),
103 fBoardLabelList(0x0),
104 fFiredFitHistoList(0x0),
105 fFiredDisplayHistoList(0x0)
106 {
107 ///  Constructor using an ASCII file.
108   CheckConstants();
109   Reset();
110   InitHistos();
111   ReadFile(filename);
112 }
113
114
115 //__________________________________________________________________________
116 AliMUONTriggerEfficiencyCells::AliMUONTriggerEfficiencyCells(TList *countHistoList,
117                                                              TList *noCountHistoList)
118 :
119 TObject(),
120 fCountHistoList(countHistoList),
121 fNoCountHistoList(noCountHistoList),
122 fFiredStrips(0x0),
123 fDisplayHistoList(0x0),
124 fBoardLabelList(0x0),
125 fFiredFitHistoList(0x0),
126 fFiredDisplayHistoList(0x0)
127 {
128 ///  Constructor using a list of histograms with counts.
129   CheckConstants();
130   Reset();
131   InitHistos();
132   FillHistosFromList();
133 }
134
135 //_____________________________________________________________________________
136 AliMUONTriggerEfficiencyCells::AliMUONTriggerEfficiencyCells(const AliMUONTriggerEfficiencyCells& other)
137 :
138 TObject(other),
139 fCountHistoList(other.fCountHistoList),
140 fNoCountHistoList(other.fNoCountHistoList),
141 fFiredStrips(other.fFiredStrips),
142 fDisplayHistoList(other.fDisplayHistoList),
143 fBoardLabelList(other.fBoardLabelList),
144 fFiredFitHistoList(other.fFiredFitHistoList),
145 fFiredDisplayHistoList(other.fFiredDisplayHistoList)
146 {
147 /// Copy constructor
148
149   for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
150     fBoardEfficiency[chCath] = other.fBoardEfficiency[chCath];
151     fSlatEfficiency[chCath] = other.fSlatEfficiency[chCath];
152   }
153 }
154
155 //_____________________________________________________________________________
156 AliMUONTriggerEfficiencyCells& AliMUONTriggerEfficiencyCells::operator=(const AliMUONTriggerEfficiencyCells& other)
157 {
158   /// Asignment operator
159   // check assignement to self
160   if (this == &other)
161     return *this;
162
163   for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
164     fBoardEfficiency[chCath] = other.fBoardEfficiency[chCath];
165     fSlatEfficiency[chCath] = other.fSlatEfficiency[chCath];
166   }
167
168   fCountHistoList = other.fCountHistoList;
169   fNoCountHistoList = other.fNoCountHistoList;
170   fFiredStrips = other.fFiredStrips;
171
172   fDisplayHistoList = other.fDisplayHistoList;
173   fBoardLabelList = other.fBoardLabelList;
174   fFiredFitHistoList = other.fFiredFitHistoList;
175   fFiredDisplayHistoList = other.fFiredDisplayHistoList;
176     
177   return *this;
178 }
179
180 //__________________________________________________________________________
181 AliMUONTriggerEfficiencyCells::~AliMUONTriggerEfficiencyCells()
182 {
183 ///  Destructor.
184
185   delete [] fBoardEfficiency;
186   delete [] fSlatEfficiency;
187   Reset();
188 }
189
190
191 //__________________________________________________________________________
192 void AliMUONTriggerEfficiencyCells::GetCellEfficiency(Int_t detElemId, Int_t localBoard, Float_t &eff1, Float_t &eff2) const
193 {
194 ///  Get the efficiencies of the 2 cathodes at a given local board
195
196   Int_t chamber = FindChamberIndex(detElemId);
197   Int_t bin = fBoardEfficiency[chamber]->FindBin(localBoard);
198   eff1 = fBoardEfficiency[chamber]->GetBinContent(bin);
199   eff2 = fBoardEfficiency[fgkNchambers+chamber]->GetBinContent(bin);
200 }
201
202
203 //__________________________________________________________________________
204 void 
205 AliMUONTriggerEfficiencyCells::IsTriggered(Int_t detElemId, Int_t localBoard, Bool_t &trig1, Bool_t &trig2) const
206 {
207 ///  Whether or not a given local board has a chance to trig, on each cathode.
208
209   Float_t eff1 = 0.0;
210   Float_t eff2 = 0.0;
211   GetCellEfficiency(detElemId, localBoard, eff1, eff2);
212   trig1 = kTRUE; 
213   trig2 = kTRUE;
214   if(gRandom->Rndm()>eff1)trig1 = kFALSE;
215   if(gRandom->Rndm()>eff2)trig2 = kFALSE;
216 }
217
218
219 //__________________________________________________________________________
220 void AliMUONTriggerEfficiencyCells::ReadFile(const Char_t* filename)
221 {
222 ///  Reads a file containing the efficiency map.
223
224   TString fileName = gSystem->ExpandPathName(filename);
225   if(fileName.EndsWith(".root")){
226       ReadHistoBoards(fileName.Data());
227       return;
228   }
229
230   ifstream file(fileName.Data());
231   Char_t dat[50];
232   if (file.good()){
233       file >> dat;
234       if(!strcmp(dat,"localBoards"))ReadFileBoards(file);
235       else AliWarning("File .dat in wrong format");
236       file.close();
237   } else {
238     AliError(Form("Can't read file %s",fileName.Data()));
239   }
240 }
241
242
243 //__________________________________________________________________________
244 void AliMUONTriggerEfficiencyCells::ReadFileBoards(ifstream &file)
245 {
246 ///  Structure of file (.dat) containing local board efficency
247   Int_t datInt=0, detEl=0, chamber=0, chCath=0, bin=0;
248     Float_t datFloat=0.0;
249     Char_t dat[50];
250
251     while (file >> dat) {
252             file >> detEl;
253             chamber = FindChamberIndex(detEl);
254             for(Int_t cath=0; cath<fgkNcathodes; cath++){
255                 chCath = GetPlane(chamber, cath);
256                 file >> dat;
257                 file >> datInt;
258                 for(Int_t board=1; board<=AliMpConstants::NofLocalBoards(); board++){
259                     file >> datFloat;
260                     bin = fBoardEfficiency[chCath]->FindBin(board);
261                     fBoardEfficiency[chCath]->SetBinContent(bin, datFloat);
262                 }
263             }
264     }
265 }
266
267
268 //__________________________________________________________________________
269 void AliMUONTriggerEfficiencyCells::ReadHistoBoards(const Char_t *filename)
270 {
271 ///  Structure of file (.root) containing local board efficency
272     TFile *file = new TFile(filename, "read");
273     if(!file || !file->IsOpen()) {
274       AliError(Form("Can't read file %s",filename));
275       return;
276     }
277
278     TString histoName;
279     TString cathCode[fgkNcathodes] = {"bendPlane", "nonBendPlane"};
280     enum {kAllChEff, kChNonEff, kNumOfHistoTypes};
281     TString histoTypeName[2] = {"CountInCh", "NonCountInCh"};
282
283     if(!fCountHistoList) fCountHistoList = new TList();
284     else fCountHistoList->Delete();
285     if(!fNoCountHistoList) fNoCountHistoList = new TList();
286     else fNoCountHistoList->Delete();
287
288     TList *currList[2] = {fCountHistoList, fNoCountHistoList};
289
290     TH1F *histo = 0x0;
291     
292     for(Int_t cath=0; cath<fgkNcathodes; cath++){
293       for(Int_t hType=0; hType<kNumOfHistoTypes; hType++){
294         histoName = Form("%sChamber%s", cathCode[cath].Data(), histoTypeName[hType].Data());
295         histo = (TH1F*)file->Get(histoName.Data());
296         currList[hType]->Add(histo);
297       }
298     }
299
300     for(Int_t cath=0; cath<fgkNcathodes; cath++){
301       for(Int_t ch=0; ch<fgkNchambers; ch++){
302         for(Int_t hType=0; hType<kNumOfHistoTypes; hType++){
303           histoName = Form("%sSlat%s%i", cathCode[cath].Data(), histoTypeName[hType].Data(), 11+ch);
304           histo = (TH1F*)file->Get(histoName.Data());
305           currList[hType]->Add(histo);
306         }
307       }
308     }
309
310     for(Int_t cath=0; cath<fgkNcathodes; cath++){
311       for(Int_t ch=0; ch<fgkNchambers; ch++){
312         for(Int_t hType=0; hType<kNumOfHistoTypes; hType++){
313           histoName = Form("%sBoard%s%i", cathCode[cath].Data(), histoTypeName[hType].Data(), 11+ch);
314           histo = (TH1F*)file->Get(histoName.Data());
315           currList[hType]->Add(histo);
316         }
317       }
318     }
319
320     FillHistosFromList();
321 }
322
323
324 //_____________________________________________________________________________
325 void AliMUONTriggerEfficiencyCells::CheckConstants() const
326 {
327 /// Check consistence of redefined constants 
328
329   assert(fgkNcathodes == AliMpConstants::NofCathodes());    
330   assert(fgkNchambers == AliMpConstants::NofTriggerChambers());    
331   assert(fgkNplanes == AliMpConstants::NofTriggerChambers() * fgkNcathodes);    
332 }
333
334
335 //__________________________________________________________________________
336 Int_t AliMUONTriggerEfficiencyCells::FindChamberIndex(Int_t detElemId) const
337 {
338 ///  From detElemId to chamber number
339
340   // Int_t iChamber = AliMpDEManager::GetChamberId(detElemId);
341   Int_t iChamber = detElemId/100 - 1;
342   Int_t chamber = iChamber-AliMpConstants::NofTrackingChambers();
343   return chamber;
344 }
345
346
347 //__________________________________________________________________________
348 void
349 AliMUONTriggerEfficiencyCells::Reset(Bool_t resetAll)
350 {
351 ///  Sets our internal array contents to zero.
352
353   if ( resetAll ){
354     for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
355       if ( fBoardEfficiency[chCath] ) {
356         fBoardEfficiency[chCath] = 0x0;
357       }
358       if ( fSlatEfficiency[chCath] ) {
359         fSlatEfficiency[chCath] = 0x0;
360       }
361     }
362
363     if ( fCountHistoList )
364       delete fCountHistoList;
365     if ( fNoCountHistoList )
366       delete fNoCountHistoList;
367   }
368
369   if ( fFiredStrips )
370     delete fFiredStrips;
371   if ( fDisplayHistoList )
372     delete fDisplayHistoList;
373   if ( fBoardLabelList )
374     delete fBoardLabelList;
375   if ( fFiredFitHistoList )
376     delete fFiredFitHistoList;
377   if ( fFiredDisplayHistoList )
378     delete fFiredDisplayHistoList;
379 }
380
381
382 //__________________________________________________________________________
383 void
384 AliMUONTriggerEfficiencyCells::InitHistos()
385 {
386 ///  Sets our internal array contents to zero.
387
388   const Int_t kNumOfBoards = AliMpConstants::NofLocalBoards();
389   const Int_t kNslats = 18;
390   Int_t chCath=0;
391   TString histoName;
392
393   TString cathCode[fgkNcathodes] = {"bendPlane", "nonBendPlane"};
394
395   for(Int_t ch=0; ch<fgkNchambers; ch++){
396     for(Int_t cath=0; cath<fgkNcathodes; cath++){
397       chCath = GetPlane(ch, cath);
398       histoName = Form("%sBoardEffChamber%i", cathCode[cath].Data(), 11+ch);
399       fBoardEfficiency[chCath] = new TH1F(histoName.Data(), histoName.Data(), kNumOfBoards, 1-0.5, kNumOfBoards+1.-0.5);
400       histoName = Form("%sSlatEffChamber%i", cathCode[cath].Data(), 11+ch);
401       fSlatEfficiency[chCath] = new TH1F(histoName.Data(), histoName.Data(), kNslats, 0-0.5, kNslats-0.5);
402     }
403   }
404 }
405
406
407 //__________________________________________________________________________
408 void
409 AliMUONTriggerEfficiencyCells::FillHistosFromList(Bool_t useMeanValues)
410 {
411 ///  Fills internal histos from list.
412
413   Int_t nHistoBins=0;
414   TH1F *histoNum = 0x0, *histoDen=0x0, *currHisto = 0x0;
415   TString slatName = "SLAT", boardName = "BOARD", histoName;
416   Float_t efficiency, efficiencyError;
417   TString chName;
418
419   Int_t nentries = fCountHistoList->GetEntries();
420
421   if ( useMeanValues )
422     AliInfo("Boards filled with the average efficiency of the RPC");
423
424   for(Int_t iEntry=0; iEntry<nentries; iEntry++){
425     histoNum = (TH1F*)fCountHistoList->At(iEntry);
426     histoDen = (TH1F*)fNoCountHistoList->At(iEntry);
427
428     if(!histoNum) {
429       AliWarning("Histogram not found in fCountHistoList. Skip to next");
430       continue;
431     }
432     if(!histoDen) {
433       AliWarning("Histogram not found in fNoCountHistoList. Skip to next");
434       continue;
435     }
436
437     histoName = histoNum->GetName();
438     histoName.ToUpper();
439     nHistoBins = histoNum->GetNbinsX();
440
441
442     Int_t currCh = 0;
443     for(Int_t ich = 11; ich<=14; ich++){
444       chName = Form("%i", ich);
445       if ( histoName.Contains(chName.Data())){
446         currCh = ich-11;
447         break;
448       }
449     }
450     Int_t currCath = ( histoName.Contains("NONBEND") ) ? 1 : 0;
451
452     Int_t chCath = GetPlane(currCh, currCath);
453
454     if(histoName.Contains(boardName)){
455       if ( useMeanValues ) continue;
456       currHisto = fBoardEfficiency[chCath];
457     }
458     else if(histoName.Contains(slatName)){
459       currHisto = fSlatEfficiency[chCath];
460     }
461     else
462       continue;
463
464     for(Int_t iBin=1; iBin<=nHistoBins; iBin++){
465       CalculateEfficiency((Int_t)histoNum->GetBinContent(iBin), (Int_t)histoNum->GetBinContent(iBin) + (Int_t)histoDen->GetBinContent(iBin), efficiency, efficiencyError, kFALSE);
466
467       currHisto->SetBinContent(iBin, efficiency);
468       currHisto->SetBinError(iBin, efficiencyError);
469     }
470
471     if ( useMeanValues ){
472       currHisto = fBoardEfficiency[chCath];
473       Int_t nBinsBoards = currHisto->GetNbinsX();
474       Int_t currChamber = currCh + AliMpConstants::NofTrackingChambers();
475       for ( Int_t iBinBoard = 1; iBinBoard<= nBinsBoards; iBinBoard++){
476         Int_t detElemId = AliMpDDLStore::Instance()->GetDEfromLocalBoard(iBinBoard, currChamber);
477         Int_t iBin = histoNum->FindBin(detElemId%100);
478
479         CalculateEfficiency((Int_t)histoNum->GetBinContent(iBin), (Int_t)histoNum->GetBinContent(iBin) + (Int_t)histoDen->GetBinContent(iBin), efficiency, efficiencyError, kFALSE);
480
481         currHisto->SetBinContent(iBinBoard, efficiency);
482         currHisto->SetBinError(iBinBoard, efficiencyError);
483       }
484     }
485   }
486 }
487
488
489 //_____________________________________________________________________________
490 void AliMUONTriggerEfficiencyCells::CalculateEfficiency(Int_t trigger44, Int_t trigger34,
491                                                         Float_t &efficiency, Float_t &error,
492                                                         Bool_t failuresAsInput)
493 {
494     //
495     /// Returns the efficiency.
496     //
497
498     efficiency=-9.;
499     error=0.;
500     if(trigger34>0){
501         efficiency=(Double_t)trigger44/((Double_t)trigger34);
502         if(failuresAsInput)efficiency=1.-(Double_t)trigger44/((Double_t)trigger34);
503     }
504     Double_t q = TMath::Abs(1-efficiency);
505     if(efficiency<0)error=0.0;
506     else error = TMath::Sqrt(efficiency*q/((Double_t)trigger34));
507 }
508
509
510 //_____________________________________________________________________________
511 void AliMUONTriggerEfficiencyCells::CheckFiredStrips()
512 {
513   //
514   /// Check for fired strips participating to efficiency
515   /// calculation (when available).
516   /// Strips inside a local board should be quite homogeneously hit
517   /// If not, this could be a problem of electronics (i.e. ADULT board off).
518   //
519
520   if(!fFiredStrips) {
521     AliWarning("List of fired pads not present. Check not performable.");
522     return;
523   }
524
525   GetListsForCheck();
526
527   TString histoName;
528
529   // Check fired pads (when available)
530   if(fFiredFitHistoList){
531     TH1F *histo1D = 0x0;
532     TF1 *fitFunc = 0x0;
533     TCanvas *histoFiredCan[20];
534     Int_t nEntries = fFiredFitHistoList->GetEntries();
535     for(Int_t iEntry=0; iEntry<nEntries; iEntry++){
536       histo1D = (TH1F*)fFiredFitHistoList->At(iEntry);
537       printf("Problems found in %s\n", histo1D->GetTitle());
538     }
539     Int_t nPrintCan = nEntries;
540     if(nPrintCan>20) {
541       AliWarning("Too many boards with problems: only 20 will be shown");
542       nPrintCan = 20;
543     }
544     for(Int_t iCan=0; iCan<nPrintCan; iCan++){
545       histo1D = (TH1F*)fFiredFitHistoList->At(iCan);
546       histoName = histo1D->GetName();
547       histoName.Append("Can");
548       histoFiredCan[iCan] = new TCanvas(histoName.Data(), histoName.Data(), 100+10*iCan, 10*iCan, 700, 700);
549       histoFiredCan[iCan]->SetRightMargin(0.14);
550       histoFiredCan[iCan]->SetLeftMargin(0.12);
551       histo1D->Draw("E");
552       fitFunc = histo1D->GetFunction("pol0");
553       fitFunc->SetLineColor(2);
554       fitFunc->Draw("same");
555     }
556     if(nEntries==0){
557       printf("\nAll local boards seem ok!!\n\n");
558     }
559   }
560 }
561
562
563 //_____________________________________________________________________________
564 void AliMUONTriggerEfficiencyCells::DisplayEfficiency(Bool_t perSlat)
565 {
566   //
567   /// Display calculated efficiency.
568   //
569
570   Bool_t isInitSlat = kFALSE, isInitBoard = kFALSE;
571   for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
572     if(fBoardEfficiency[chCath]->GetEntries()>0) isInitBoard = kTRUE;
573     if(fSlatEfficiency[chCath]->GetEntries()>0) isInitSlat = kTRUE;
574   }
575
576   if(!isInitBoard){
577     AliWarning("Trigger efficiency not initialized per board.\nDisplay not yet implemented.");
578     return;
579   }
580   if(!isInitSlat && perSlat){
581     AliWarning("Trigger efficiency not initialized for slat.\nPlease try option kFALSE.");
582     return;
583   }
584   if ( !AliCDBManager::Instance()->GetDefaultStorage() ){
585     AliWarning("Please set default CDB storage (needed for mapping).");
586     return;
587   }
588   if ( AliCDBManager::Instance()->GetRun() < 0 ){
589     AliWarning("Please set CDB run number (needed for mapping).");
590     return;
591   }
592   
593   GetListsForCheck();
594
595   //const Int_t kNumOfBoards = AliMpConstants::NofLocalBoards();
596
597   TH2F *histo = 0x0;
598   TString histoName, histoTitle;
599
600   // Plot fired strips (when available)
601   if(fFiredDisplayHistoList){
602     TCanvas *displayFiredCan[fgkNplanes];
603     for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
604       histo = (TH2F*)fFiredDisplayHistoList->At(chCath);
605       histoName = Form("%sCan", histo->GetName());
606       histoTitle = Form("%s", histo->GetTitle());
607       displayFiredCan[chCath] = new TCanvas(histoName.Data(), histoTitle.Data(), 100+10*chCath, 10*chCath, 700, 700);
608       displayFiredCan[chCath]->SetRightMargin(0.14);
609       displayFiredCan[chCath]->SetLeftMargin(0.12);
610       histo->GetYaxis()->SetTitleOffset(1.4);
611       histo->SetStats(kFALSE);
612       histo->Draw("COLZ");
613     }
614   }
615
616   // Plot efficiency
617   if(fDisplayHistoList){
618     TCanvas *can[fgkNplanes];
619     for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
620       Int_t currChCath = chCath;
621       if(perSlat==kTRUE) currChCath += fgkNplanes;
622       histo = (TH2F*)fDisplayHistoList->At(currChCath);
623       histoName = Form("%sCan", histo->GetName());
624       histoTitle = Form("%s", histo->GetTitle());
625       can[chCath] = new TCanvas(histoName.Data(), histoTitle.Data(), 100+10*chCath, 10*chCath, 700, 700);
626       can[chCath]->SetRightMargin(0.14);
627       can[chCath]->SetLeftMargin(0.12);
628       histo->GetZaxis()->SetRangeUser(0.,1.);
629       histo->GetYaxis()->SetTitleOffset(1.4);
630       histo->SetStats(kFALSE);
631       histo->Draw("COLZ");
632       if(perSlat==kTRUE) continue;
633       histo = (TH2F*)fBoardLabelList->At(currChCath);
634       histo->Draw("textsame");
635     }
636   }
637 }
638
639
640 //__________________________________________________________________________
641 Bool_t AliMUONTriggerEfficiencyCells::GetListsForCheck()
642 {
643   //
644   /// Getting histograms for efficiency, 
645   /// map of fired strips entering efficiency calculations,
646   /// fits for checking switched-off elements in chambers.
647   //
648   const Float_t kChi2RedMax = 1.5;
649   const Float_t kDummyFired = 1e-5;
650
651   Reset(kFALSE);
652
653   fDisplayHistoList = new TList();
654   fDisplayHistoList->SetOwner();
655   fBoardLabelList = new TList();
656   fBoardLabelList->SetOwner();
657   if ( fFiredStrips ){
658     fFiredFitHistoList = new TList();
659     fFiredFitHistoList->SetOwner();
660     fFiredDisplayHistoList = new TList();
661     fFiredDisplayHistoList->SetOwner();
662   }
663
664   TH3F* padFired = 0x0;
665   TH2F* displayHisto = 0x0;
666   TH1F* histoFired[fgkNplanes][234];
667   TF1* fitFunc = 0x0;
668   Bool_t isStripOffInBoard[fgkNplanes][234];
669   Int_t bin=0;
670   TString histoName, histoTitle;
671   TString cathName[fgkNcathodes] = {"bendPlane", "nonBendPlane"};
672
673   AliMUONTriggerDisplay triggerDisplay;
674   // Book histos
675   for(Int_t iCath=0; iCath<fgkNcathodes; iCath++){
676     if(fFiredStrips) padFired = (TH3F*)fFiredStrips->At(iCath);
677     for(Int_t iCh=0; iCh<fgkNchambers; iCh++){
678       Int_t chCath = GetPlane(iCh, iCath);
679       Int_t currCh = 11 + iCh;
680       histoName = Form("%sChamber%i", cathName[iCath].Data(), currCh);
681       histoTitle = Form("Chamber %i: efficiency %s", currCh, cathName[iCath].Data());
682       displayHisto = 
683         (TH2F*)triggerDisplay.GetDisplayHistogram(fBoardEfficiency[chCath], histoName,
684                                                   AliMUONTriggerDisplay::kDisplayBoards,
685                                                   iCath,currCh,histoTitle,
686                                                   AliMUONTriggerDisplay::kShowZeroes);
687       fDisplayHistoList->Add(displayHisto);
688
689       histoName = Form("labels%sChamber%i", cathName[iCath].Data(), currCh);
690       displayHisto = 
691         (TH2F*)triggerDisplay.GetBoardNumberHisto(histoName,currCh);
692       fBoardLabelList->Add(displayHisto);
693
694       if(!fFiredStrips) continue;
695       histoName = Form("firedPads%sChamber%i", cathName[iCath].Data(), currCh);
696       histoTitle = Form("Chamber %i: Fired pads %s", currCh, cathName[iCath].Data());
697       bin = padFired->GetXaxis()->FindBin(currCh);
698       padFired->GetXaxis()->SetRange(bin,bin);
699       displayHisto = 
700         (TH2F*)triggerDisplay.GetDisplayHistogram(padFired->Project3D("zy"),histoName,
701                                                   AliMUONTriggerDisplay::kDisplayStrips,
702                                                   iCath,currCh,histoTitle);
703       fFiredDisplayHistoList->Add(displayHisto);
704       
705       for(Int_t ib=0; ib<AliMpConstants::NofLocalBoards(); ib++){
706         histoName = Form("%sChamber%iBoard%i", cathName[iCath].Data(), currCh, ib+1);
707         histoTitle = Form("Chamber %i: fired pads %s board = %i", currCh, cathName[iCath].Data(), ib+1);
708         histoFired[chCath][ib] = new TH1F(histoName.Data(), histoTitle.Data(), 16, -0.5, 15.5);
709         histoFired[chCath][ib]->SetXTitle("board");
710         isStripOffInBoard[chCath][ib] = kFALSE;
711       } // loop on board
712     } // loop on chamber
713   } // loop on cathode
714
715   for(Int_t iCath=0; iCath<fgkNcathodes; iCath++){
716     for(Int_t iCh=0; iCh<fgkNchambers; iCh++){
717       Int_t chCath = GetPlane(iCh, iCath);
718       Int_t currCh = 11+iCh;
719       histoName = Form("%sChamber%iSlatEff", cathName[iCath].Data(), currCh);
720       histoTitle = Form("Chamber %i: efficiency %s per slat", currCh, cathName[iCath].Data());
721       displayHisto = 
722         (TH2F*)triggerDisplay.GetDisplayHistogram(fSlatEfficiency[chCath], histoName,
723                                                   AliMUONTriggerDisplay::kDisplaySlats,
724                                                   iCath,currCh,histoTitle);
725       fDisplayHistoList->Add(displayHisto);
726     }
727   }
728
729   if(!fFiredStrips) return kTRUE;
730
731   AliMpCDB::LoadDDLStore();
732
733   // Check fired pads (when available)
734   for(Int_t iLoc = 0; iLoc < AliMpConstants::NofLocalBoards(); iLoc++) {  
735     Int_t iBoard = iLoc+1;
736     for(Int_t iCh=0; iCh<AliMpConstants::NofChambers(); iCh++){
737       Int_t iChamber = iCh + AliMpConstants::NofTrackingChambers();
738
739       Int_t detElemId = AliMpDDLStore::Instance()->GetDEfromLocalBoard(iBoard, iChamber);
740
741       if (!detElemId) continue;
742
743       AliMpLocalBoard* localBoard = AliMpDDLStore::Instance()->GetLocalBoard(iBoard, kFALSE);
744
745       // skip copy cards
746       if( !localBoard->IsNotified()) 
747         continue;
748       
749       for(Int_t iCath=0; iCath<AliMpConstants::NofCathodes(); iCath++){
750         Int_t chCath = GetPlane(iCh, iCath);
751         // loop over strips
752         const AliMpVSegmentation* seg = 
753           AliMpSegmentation::Instance()->GetMpSegmentation(detElemId, AliMp::GetCathodType(iCath));
754         Int_t nStrips=0;
755         for (Int_t iStrip = 0; iStrip < 16; ++iStrip) {
756           AliMpPad pad = seg->PadByLocation(iBoard,iStrip,kFALSE);
757           if (!pad.IsValid()) continue;
758           nStrips++;
759           padFired = (TH3F*)fFiredStrips->At(iCath);
760           Float_t nFired = padFired->GetBinContent(11+iCh, iBoard, iStrip);
761           if(nFired==0.) nFired = kDummyFired;
762           histoFired[chCath][iLoc]->Fill(iStrip, nFired);
763         }
764
765         histoFired[chCath][iLoc]->Fit("pol0","Q0R","",0., (Float_t)nStrips-1.);
766         fitFunc = histoFired[chCath][iLoc]->GetFunction("pol0");
767         Float_t chi2 = fitFunc->GetChisquare();
768         Float_t ndf = (Float_t)fitFunc->GetNDF();
769         Float_t reducedChi2 = chi2/ndf;
770         if(reducedChi2>kChi2RedMax) {
771           isStripOffInBoard[chCath][iLoc] = kTRUE;
772           fFiredFitHistoList->Add(histoFired[chCath][iLoc]);
773         }
774       } // loop on cathodes
775     } // loop on chambers
776   } // loop on local boards 
777
778   return kTRUE;
779 }
780
781 //__________________________________________________________________________
782 Bool_t AliMUONTriggerEfficiencyCells::SumRunEfficiency(const AliMUONTriggerEfficiencyCells &other)
783 {
784 ///  Sums results from different runs and gives the efficiency
785   if(!fCountHistoList || !fNoCountHistoList) {
786     AliWarning("Histograms for efficiency calculations not implemented in object");
787     return kFALSE;
788   }
789   if(!other.fCountHistoList || !other.fNoCountHistoList) {
790     AliWarning("Histograms for efficiency calculations not implemented in object passed as argument");
791     return kFALSE;
792   }
793
794   Int_t nentries = fCountHistoList->GetEntries();
795   TH1F *currNum = 0x0, *currDen = 0x0, *otherNum = 0x0, *otherDen = 0x0;
796
797   for(Int_t iEntry=0; iEntry<nentries; iEntry++){
798     currNum = (TH1F*)fCountHistoList->At(iEntry);
799     currDen = (TH1F*)fNoCountHistoList->At(iEntry);
800     otherNum = (TH1F*)other.fCountHistoList->At(iEntry);
801     otherDen = (TH1F*)other.fNoCountHistoList->At(iEntry);
802     currNum->Add(otherNum);
803     currDen->Add(otherDen);
804   }
805
806   FillHistosFromList();
807
808   if(!fFiredStrips) {
809     AliWarning("Histograms for fired region check not implemented in object");
810     return kFALSE;
811   }
812   if(!other.fFiredStrips) {
813     AliWarning("Histograms for fired region check not implemented in object passed as argument");
814     return kFALSE;
815   }
816   
817   TH3F *currFired = 0x0, *otherFired = 0x0;
818   nentries = fFiredStrips->GetEntries();
819   for(Int_t iEntry=0; iEntry<nentries; iEntry++){
820     currFired  = (TH3F*)fFiredStrips->At(iEntry);
821     otherFired = (TH3F*)fFiredStrips->At(iEntry);
822     currFired->Add(otherFired);
823   }
824     
825   return kTRUE;
826 }
827
828
829 //__________________________________________________________________________
830 Bool_t AliMUONTriggerEfficiencyCells::LowStatisticsSettings(Bool_t useMeanValues)
831 {
832   //
833   /// In case of low statistics, fill the local board efficiency with
834   /// the average value of the RPC
835   //
836
837   if ( ! fCountHistoList ){
838     AliWarning("Cannot find histograms for RPC efficiency calculation. Nothing done!");
839     return kFALSE;
840   }
841   if ( !AliCDBManager::Instance()->GetDefaultStorage() ){
842     AliWarning("Please set default CDB storage (needed for mapping).");
843     return kFALSE;
844   }
845   if ( AliCDBManager::Instance()->GetRun() < 0 ){
846     AliWarning("Please set CDB run number (needed for mapping).");
847     return kFALSE;
848   }
849
850   AliMpCDB::LoadDDLStore();
851   
852   for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
853     fBoardEfficiency[chCath]->Reset();
854     fSlatEfficiency[chCath]->Reset();
855   }
856
857   FillHistosFromList(useMeanValues);
858
859   return kTRUE;
860 }