]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONTriggerQADataMakerRec.cxx
In AliMUONTriggerQADataMakerRec:
[u/mrichter/AliRoot.git] / MUON / AliMUONTriggerQADataMakerRec.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: AliMUONTriggerQADataMakerRec.cxx 35760 2009-10-21 21:45:42Z ivana $
17
18 // --- MUON header files ---
19 #include "AliMUONTriggerQADataMakerRec.h"
20
21 //-----------------------------------------------------------------------------
22 /// \class AliMUONTriggerQADataMakerRec
23 ///
24 /// MUON class for quality assurance data (histo) maker
25 ///
26 /// \author C. Finck, D. Stocco, L. Aphecetche, A. Blanc
27
28 /// \cond CLASSIMP
29 ClassImp(AliMUONTriggerQADataMakerRec)
30 /// \endcond
31            
32 #include "AliCodeTimer.h"
33 #include "AliMUONConstants.h"
34 #include "AliMpConstants.h"
35 #include "AliMUONTriggerDisplay.h"
36 #include "TH2.h"
37 #include "TH1F.h"
38 #include "TString.h"
39 #include "AliRecoParam.h"
40 #include "AliMUONDigitStoreV2R.h"
41 #include "AliMUONTriggerStoreV1.h"
42 #include "AliMpCDB.h"
43 #include "AliMUONRawStreamTriggerHP.h"
44 #include "AliMpDDLStore.h"
45 #include "AliMpTriggerCrate.h"
46 #include "AliMpLocalBoard.h"
47 #include "AliQAv1.h"
48 #include "AliRawReader.h"
49 #include "AliMUONDigitMaker.h"
50 #include "AliMUONLocalTrigger.h"
51 #include "AliMUONRecoParam.h"
52 #include "AliMUONTriggerElectronics.h"
53 #include "AliMUONCalibrationData.h"
54 #include "AliDCSValue.h"
55 #include "AliMpDCSNamer.h"
56 #include "AliMpDEManager.h"
57 #include "AliMpDEIterator.h"
58 #include "AliCDBManager.h"
59 #include "TTree.h"
60 #include "AliMUONGlobalTriggerBoard.h"
61 #include "AliMUONGlobalTrigger.h"
62 #include "AliMUONGlobalCrateConfig.h"
63 #include "AliMUONQAIndices.h"
64
65 namespace
66 {
67   Double_t ProtectedSqrt(Double_t x)
68   {
69     return ( x > 0.0 ? TMath::Sqrt(x) : 0.0 );
70   }
71 }
72 //____________________________________________________________________________ 
73 AliMUONTriggerQADataMakerRec::AliMUONTriggerQADataMakerRec(AliQADataMakerRec* master) : 
74 AliMUONVQADataMakerRec(master),
75 fNumberOf34Dec(0x0),
76 fNumberOf44Dec(0x0),
77 fDigitMaker(new AliMUONDigitMaker(kFALSE)),
78 fCalibrationData(new AliMUONCalibrationData(AliCDBManager::Instance()->GetRun())),
79 fTriggerProcessor(new AliMUONTriggerElectronics(fCalibrationData)),
80 fDigitStore(0x0)
81 {
82     /// ctor
83 }
84
85
86 //__________________________________________________________________
87 AliMUONTriggerQADataMakerRec::~AliMUONTriggerQADataMakerRec()
88 {
89     /// dtor
90   delete fDigitMaker;
91   delete fDigitStore;
92   delete fTriggerProcessor;
93   delete fCalibrationData;
94   delete fNumberOf34Dec;
95   delete fNumberOf44Dec;
96 }
97
98 //____________________________________________________________________________ 
99 void AliMUONTriggerQADataMakerRec::EndOfDetectorCycleESDs(Int_t /*specie*/, TObjArray** /*list*/)
100 {
101   /// Normalize ESD histograms
102 }
103   
104 //____________________________________________________________________________ 
105 void AliMUONTriggerQADataMakerRec::EndOfDetectorCycleRecPoints(Int_t /*specie*/, TObjArray** /*list*/)
106 {
107   /// Normalize RecPoints histograms
108   
109 }
110
111
112 //____________________________________________________________________________ 
113 void AliMUONTriggerQADataMakerRec::EndOfDetectorCycleRaws(Int_t /*specie*/, TObjArray** /*list*/)
114 {
115   /// create Raws histograms in Raws subdir
116
117   DisplayTriggerInfo();
118
119   // Normalize RawData histos
120   Float_t nbevent = GetRawsData(AliMUONQAIndices::kTriggerRawNAnalyzedEvents)->GetBinContent(1);
121   Int_t histoRawsIndex[] = {
122     AliMUONQAIndices::kTriggerErrorSummary,
123     AliMUONQAIndices::kTriggerCalibSummary,
124     AliMUONQAIndices::kTriggerReadOutErrors,
125     AliMUONQAIndices::kTriggerGlobalOutput
126   };
127   Int_t histoRawsScaledIndex[] = {
128     AliMUONQAIndices::kTriggerErrorSummaryNorm,
129     AliMUONQAIndices::kTriggerCalibSummaryNorm,
130     AliMUONQAIndices::kTriggerReadOutErrorsNorm,
131     AliMUONQAIndices::kTriggerGlobalOutputNorm
132   };
133   const Int_t kNrawsHistos = sizeof(histoRawsIndex)/sizeof(histoRawsIndex[0]);
134   Float_t scaleFactor[kNrawsHistos] = {100., 100., 100., 1.};
135   for(Int_t ihisto=0; ihisto<kNrawsHistos; ihisto++){
136     TH1* inputHisto = GetRawsData(histoRawsIndex[ihisto]);
137     TH1* scaledHisto = GetRawsData(histoRawsScaledIndex[ihisto]);
138     if ( scaledHisto && inputHisto &&  nbevent > 0 ) {
139       scaledHisto->Reset();
140       scaledHisto->Add(inputHisto);
141       scaledHisto->Scale(scaleFactor[ihisto]/nbevent);
142     }
143   } // loop on histos
144
145   TH1* hYCopy = GetRawsData(AliMUONQAIndices::kTriggerErrorLocalYCopy); //number of YCopy error per board
146   TH1* hYCopyTests = GetRawsData(AliMUONQAIndices::kTriggerErrorLocalYCopyTest); //contains the number of YCopy test per board
147   TH1* hYCopyNorm = GetRawsData(AliMUONQAIndices::kTriggerErrorLocalYCopyNorm); 
148   hYCopyNorm->Reset();
149   hYCopyNorm->Divide(hYCopy, hYCopyTests, 100., 1.);
150      
151   Float_t mean = hYCopyNorm->Integral();
152       
153   TH1* hSummary = GetRawsData(AliMUONQAIndices::kTriggerErrorSummary);
154   hSummary->SetBinContent(AliMUONQAIndices::kAlgoLocalYCopy+1,mean/192.); //put the mean of the % of YCopy error in the kTriggerError's corresponding bin
155
156   ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerLocalRatio4434))->Divide(fNumberOf44Dec,fNumberOf34Dec);
157
158   FillRatio4434Histos(1);
159
160   //reset bins temporary used to store informations
161   ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerRatio4434AllEvents))->SetBinContent(0,0); 
162   Int_t nbins =  ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerRatio4434AllEvents))->GetNbinsX();
163   ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerRatio4434AllEvents))->SetBinContent(nbins+1,0);
164
165   ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerLocalRatio4434))->SetMaximum(1.1);
166   ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerRatio4434AllEvents))->SetMaximum(1.1);
167   ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerRatio4434SinceLastUpdate))->SetMaximum(1.1);
168 }
169
170 //____________________________________________________________________________ 
171 void AliMUONTriggerQADataMakerRec::InitRaws()
172 {
173     /// create Raws histograms in Raws subdir
174         
175   AliCodeTimerAuto("",0);
176   
177   const Bool_t expert   = kTRUE ; 
178   const Bool_t saveCorr = kTRUE ; 
179   const Bool_t image    = kTRUE ; 
180  
181   TString boardName = "Local board Id";
182
183   Int_t nbLocalBoard = AliMUONConstants::NTriggerCircuit();
184
185   TH1F* histo1D = 0x0;
186   TH2F* histo2D = 0x0;
187
188   AliMUONTriggerDisplay triggerDisplay;
189
190   TString histoName, histoTitle;
191   if ( GetRecoParam()->GetEventSpecie() == AliRecoParam::kCalib ) {
192     histo1D = new TH1F("hTriggerScalersTime", "Acquisition time from trigger scalers", 1, 0.5, 1.5);
193     histo1D->GetXaxis()->SetBinLabel(1, "One-bin histogram: bin is filled at each scaler event.");
194     histo1D->GetYaxis()->SetTitle("Cumulated scaler time (s)");
195     Add2RawsList(histo1D, AliMUONQAIndices::kTriggerScalersTime, expert, !image, !saveCorr);
196
197     for(Int_t iCath=0; iCath<AliMpConstants::NofCathodes(); iCath++){
198       TString cathName = ( iCath==0 ) ? "BendPlane" : "NonBendPlane";
199       for(Int_t iChamber=0; iChamber<AliMpConstants::NofTriggerChambers(); iChamber++){
200         histoName = Form("hTriggerScalers%sChamber%i", cathName.Data(), 11+iChamber);
201         histoTitle = Form("Chamber %i - %s: trigger scaler counts", 11+iChamber, cathName.Data());
202         histo2D = new TH2F(histoName.Data(), histoTitle.Data(),
203                            nbLocalBoard, 0.5, (Float_t)nbLocalBoard + 0.5,
204                            16, -0.5, 15.5);
205         histo2D->GetXaxis()->SetTitle(boardName.Data());
206         histo2D->GetYaxis()->SetTitle("Strip"); 
207         histo2D->SetOption("COLZ");     
208         Add2RawsList(histo2D, AliMUONQAIndices::kTriggerScalers + AliMpConstants::NofTriggerChambers()*iCath + iChamber, expert, !image, !saveCorr);
209       } // loop on chambers
210     } // loop on cathodes
211         
212     for(Int_t iCath=0; iCath<AliMpConstants::NofCathodes(); iCath++){
213       TString cathName = ( iCath==0 ) ? "BendPlane" : "NonBendPlane";
214       for(Int_t iChamber=0; iChamber<AliMpConstants::NofTriggerChambers(); iChamber++){
215         histoName = Form("hTriggerScalersDisplay%sChamber%i", cathName.Data(), 11+iChamber);
216         histoTitle = Form("Chamber %i - %s: Hit rate from scalers (Hz/cm^{2})", 11+iChamber, cathName.Data());
217         histo2D = (TH2F*)triggerDisplay.GetEmptyDisplayHisto(histoName, AliMUONTriggerDisplay::kDisplayStrips, 
218                                                              iCath, iChamber, histoTitle);
219         histo2D->SetOption("COLZ");
220         Add2RawsList(histo2D, AliMUONQAIndices::kTriggerScalersDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber, expert, !image, !saveCorr);
221       } // loop on chambers
222     } // loop on cathodes    
223
224     TString axisLabel[AliMUONQAIndices::kNtrigCalibSummaryBins] = {"#splitline{Dead}{Channels}", "#splitline{Dead}{Local Boards}", "#splitline{Dead}{Regional Boards}", "#splitline{Dead}{Global Board}", "#splitline{Noisy}{Strips}"};
225
226     TH1F* histoCalib = new TH1F("hTriggerCalibSummaryAll", "MTR calibration summary counts", AliMUONQAIndices::kNtrigCalibSummaryBins, -0.5, (Float_t)AliMUONQAIndices::kNtrigCalibSummaryBins - 0.5);
227     for (Int_t ibin=1; ibin<=AliMUONQAIndices::kNtrigCalibSummaryBins; ibin++){
228       histoCalib->GetXaxis()->SetBinLabel(ibin, axisLabel[ibin-1].Data());
229     }
230     histoCalib->SetFillColor(kBlue);
231     histoCalib->GetYaxis()->SetTitle("Counts");
232     // Copy of previous histo for scaling purposes
233     TH1F* histoCalibNorm = (TH1F*)histoCalib->Clone("hTriggerCalibSummary");
234     histoCalibNorm->SetTitle("MTR calibration summary");
235     histoCalibNorm->SetOption("bar2");
236     histoCalibNorm->GetYaxis()->SetTitle("Percentage per event (%)");
237     // Adding both histos after cloning to avoid problems with the expert bit
238     Add2RawsList(histoCalib,     AliMUONQAIndices::kTriggerCalibSummary,      expert, !image, !saveCorr);
239     Add2RawsList(histoCalibNorm, AliMUONQAIndices::kTriggerCalibSummaryNorm, !expert,  image, !saveCorr);
240   } // Calibration reco param
241         
242   const char *globalXaxisName[6] = {"US HPt", "US LPt", "LS HPt", "LS LPt", "SGL HPt", "SGL LPt"};
243   const char *allLevelXaxisName[AliMUONQAIndices::kNtrigAlgoErrorBins] = {"Local algo X", "Local algo Y", "Local LUT","Local Y Copy" , "Local2Regional", "Regional", "Regional2Global", "GlobalFromInGlobal", "GlobalFromInLocal", "GlobalFromOutLocal"};
244   const char *readoutErrNames[AliMUONQAIndices::kNtrigStructErrorBins]={"Local","Regional","Global","DARC"};
245
246   TString errorAxisTitle = "Number of errors";
247
248   histo1D = new TH1F("hTriggerErrorLocalXPos", "ErrorLocalXPos",nbLocalBoard,0.5,(Float_t)nbLocalBoard+0.5);
249   histo1D->GetXaxis()->SetTitle(boardName.Data());
250   histo1D->GetYaxis()->SetTitle(errorAxisTitle.Data());
251   Add2RawsList(histo1D, AliMUONQAIndices::kTriggerErrorLocalXPos, expert, !image, !saveCorr);
252
253   histo1D = new TH1F("hTriggerErrorLocalYPos", "ErrorLocalYPos",nbLocalBoard,0.5,(Float_t)nbLocalBoard+0.5);
254   histo1D->GetXaxis()->SetTitle(boardName.Data());
255   histo1D->GetYaxis()->SetTitle(errorAxisTitle.Data());
256   Add2RawsList(histo1D, AliMUONQAIndices::kTriggerErrorLocalYPos, expert, !image, !saveCorr);
257
258   histo1D = new TH1F("hTriggerErrorLocalDev", "ErrorLocalDev",nbLocalBoard,0.5,(Float_t)nbLocalBoard+0.5);
259   histo1D->GetXaxis()->SetTitle(boardName.Data());
260   histo1D->GetYaxis()->SetTitle(errorAxisTitle.Data());
261   Add2RawsList(histo1D, AliMUONQAIndices::kTriggerErrorLocalDev, expert, !image, !saveCorr);
262
263   histo1D = new TH1F("hTriggerErrorLocalTriggerDec", "ErrorLocalTriggerDec",nbLocalBoard,0.5,(Float_t)nbLocalBoard+0.5);
264   histo1D->GetXaxis()->SetTitle(boardName.Data());
265   histo1D->GetYaxis()->SetTitle(errorAxisTitle.Data());
266   Add2RawsList(histo1D, AliMUONQAIndices::kTriggerErrorLocalTriggerDec, expert, !image, !saveCorr);
267
268   histo1D = new TH1F("hTriggerErrorLocalLPtLSB", "ErrorLocalLPtLSB",nbLocalBoard,0.5,(Float_t)nbLocalBoard+0.5);
269   histo1D->GetXaxis()->SetTitle(boardName.Data());
270   histo1D->GetYaxis()->SetTitle(errorAxisTitle.Data());
271   Add2RawsList(histo1D, AliMUONQAIndices::kTriggerErrorLocalLPtLSB, expert, !image, !saveCorr);
272
273   histo1D = new TH1F("hTriggerErrorLocalLPtMSB", "ErrorLocalLPtMSB",nbLocalBoard,0.5,(Float_t)nbLocalBoard+0.5);
274   histo1D->GetXaxis()->SetTitle(boardName.Data());
275   histo1D->GetYaxis()->SetTitle(errorAxisTitle.Data());
276   Add2RawsList(histo1D, AliMUONQAIndices::kTriggerErrorLocalLPtMSB, expert, !image, !saveCorr);
277
278   histo1D = new TH1F("hTriggerErrorLocalHPtLSB", "ErrorLocalHPtLSB",nbLocalBoard,0.5,(Float_t)nbLocalBoard+0.5);
279   histo1D->GetXaxis()->SetTitle(boardName.Data());
280   histo1D->GetYaxis()->SetTitle(errorAxisTitle.Data());
281   Add2RawsList(histo1D, AliMUONQAIndices::kTriggerErrorLocalHPtLSB, expert, !image, !saveCorr);
282
283   histo1D = new TH1F("hTriggerErrorLocalHPtMSB", "ErrorLocalHPtMSB",nbLocalBoard,0.5,(Float_t)nbLocalBoard+0.5);
284   histo1D->GetXaxis()->SetTitle(boardName.Data());
285   histo1D->GetYaxis()->SetTitle(errorAxisTitle.Data());
286   Add2RawsList(histo1D, AliMUONQAIndices::kTriggerErrorLocalHPtMSB, expert, !image, !saveCorr);
287
288   histo1D = new TH1F("hTriggerErrorLocalTrigY", "ErrorLocalTrigY",nbLocalBoard,0.5,(Float_t)nbLocalBoard+0.5);
289   histo1D->GetXaxis()->SetTitle(boardName.Data());
290   histo1D->GetYaxis()->SetTitle(errorAxisTitle.Data());
291   Add2RawsList(histo1D, AliMUONQAIndices::kTriggerErrorLocalTrigY, expert, !image, !saveCorr);
292
293   histo1D = new TH1F("hTriggerRatio4434Local", "Ratio4434Local",nbLocalBoard,0.5,(Float_t)nbLocalBoard+0.5);
294   histo1D->GetXaxis()->SetTitle(boardName.Data());
295   histo1D->GetYaxis()->SetTitle("ratio 44/34");
296   Add2RawsList(histo1D, AliMUONQAIndices::kTriggerLocalRatio4434, expert, !image, !saveCorr);                                               
297   histo1D = new TH1F("hTriggerRatio4434AllEvents", "Ratio4434AllEvents",1,0,1);
298   histo1D->GetXaxis()->SetTitle("Event number");
299   histo1D->GetYaxis()->SetTitle("ratio 44/34");
300   histo1D->SetLineColor(4);                           
301   Add2RawsList(histo1D, AliMUONQAIndices::kTriggerRatio4434AllEvents, expert, !image, !saveCorr);                                               
302   histo1D = new TH1F("hTriggerRatio4434SinceLastUpdate", "Ratio4434SinceLastUpdate",1,0,1);
303   histo1D->GetXaxis()->SetTitle("Event number");
304   histo1D->GetYaxis()->SetTitle("ratio 44/34");                           
305   Add2RawsList(histo1D, AliMUONQAIndices::kTriggerRatio4434SinceLastUpdate, expert, !image, !saveCorr);
306
307   histo1D = new TH1F("hTriggerErrorLocal2RegionalLPtLSB", "ErrorLocal2RegionalLPtLSB",nbLocalBoard,0.5,(Float_t)nbLocalBoard+0.5);
308   histo1D->GetXaxis()->SetTitle(boardName.Data());
309   histo1D->GetYaxis()->SetTitle(errorAxisTitle.Data());
310   Add2RawsList(histo1D, AliMUONQAIndices::kTriggerErrorLocal2RegionalLPtLSB, expert, !image, !saveCorr);
311
312   histo1D = new TH1F("hTriggerErrorLocal2RegionalLPtMSB", "ErrorLocal2RegionalLPtMSB",nbLocalBoard,0.5,(Float_t)nbLocalBoard+0.5);
313   histo1D->GetXaxis()->SetTitle(boardName.Data());
314   histo1D->GetYaxis()->SetTitle(errorAxisTitle.Data());
315   Add2RawsList(histo1D, AliMUONQAIndices::kTriggerErrorLocal2RegionalLPtMSB, expert, !image, !saveCorr);
316
317   histo1D = new TH1F("hTriggerErrorLocal2RegionalHPtLSB", "ErrorLocal2RegionalHPtLSB",nbLocalBoard,0.5,(Float_t)nbLocalBoard+0.5);
318   histo1D->GetXaxis()->SetTitle(boardName.Data());
319   histo1D->GetYaxis()->SetTitle(errorAxisTitle.Data());
320   Add2RawsList(histo1D, AliMUONQAIndices::kTriggerErrorLocal2RegionalHPtLSB, expert, !image, !saveCorr);
321
322   histo1D = new TH1F("hTriggerErrorLocal2RegionalHPtMSB", "ErrorLocal2RegionalHPtMSB",nbLocalBoard,0.5,(Float_t)nbLocalBoard+0.5);
323   histo1D->GetXaxis()->SetTitle(boardName.Data());
324   histo1D->GetYaxis()->SetTitle(errorAxisTitle.Data());
325   Add2RawsList(histo1D, AliMUONQAIndices::kTriggerErrorLocal2RegionalHPtMSB, expert, !image, !saveCorr);
326
327   histo1D = new TH1F("hTriggerErrorOutGlobalFromInGlobal", "ErrorOutGlobalFromInGlobal",6,-0.5,6-0.5);
328   histo1D->GetYaxis()->SetTitle(errorAxisTitle.Data());
329   for (int ibin=0;ibin<6;ibin++){
330     histo1D->GetXaxis()->SetBinLabel(ibin+1,globalXaxisName[ibin]);
331   }
332   Add2RawsList(histo1D, AliMUONQAIndices::kTriggerErrorOutGlobalFromInGlobal, expert, !image, !saveCorr);
333
334   histo1D = new TH1F("hTriggerErrorOutGlobalFromInLocal", "ErrorOutGlobalFromInLocal",6,-0.5,6-0.5);
335   histo1D->GetYaxis()->SetTitle(errorAxisTitle.Data());
336   for (int ibin=0;ibin<6;ibin++){
337     histo1D->GetXaxis()->SetBinLabel(ibin+1,globalXaxisName[ibin]);
338   }
339   Add2RawsList(histo1D, AliMUONQAIndices::kTriggerErrorOutGlobalFromInLocal, expert, !image, !saveCorr);
340
341   TH1F* histoAlgoErr = new TH1F("hTriggerAlgoNumOfErrors", "Trigger Algorithm total errors",AliMUONQAIndices::kNtrigAlgoErrorBins,-0.5,(Float_t)AliMUONQAIndices::kNtrigAlgoErrorBins-0.5);
342   histoAlgoErr->GetYaxis()->SetTitle("Number of events with errors");
343   for (int ibin=0;ibin<AliMUONQAIndices::kNtrigAlgoErrorBins;ibin++){
344     histoAlgoErr->GetXaxis()->SetBinLabel(ibin+1,allLevelXaxisName[ibin]);
345   }
346   histoAlgoErr->SetFillColor(kBlue);
347   // Copy of previous histo for scaling purposes
348   TH1F* histoAlgoErrNorm = (TH1F*)histoAlgoErr->Clone("hTriggerAlgoErrors");
349   histoAlgoErrNorm->SetOption("bar2");
350   histoAlgoErrNorm->SetTitle("Trigger algorithm errors");
351   histoAlgoErrNorm->GetYaxis()->SetTitle("% of events with errors");
352   // Adding both histos after cloning to avoid problems with the expert bit
353   Add2RawsList(histoAlgoErr,     AliMUONQAIndices::kTriggerErrorSummary,      expert, !image, !saveCorr);
354   Add2RawsList(histoAlgoErrNorm, AliMUONQAIndices::kTriggerErrorSummaryNorm, !expert,  image, !saveCorr);  
355
356   histo1D = new TH1F("hTriggerTriggeredBoards", "Triggered boards", nbLocalBoard, 0.5, (Float_t)nbLocalBoard + 0.5);
357   Add2RawsList(histo1D, AliMUONQAIndices::kTriggeredBoards, expert, !image, !saveCorr);
358
359   histo2D = (TH2F*)triggerDisplay.GetEmptyDisplayHisto("hTriggerFiredBoardsDisplay", AliMUONTriggerDisplay::kDisplayBoards,
360                                                        0, 0, "Local board triggers / event");
361   histo2D->SetOption("COLZ");
362   Add2RawsList(histo2D, AliMUONQAIndices::kTriggerBoardsDisplay, expert, !image, !saveCorr);
363
364   TH1F* histoYCopyErr = new TH1F("hTriggerErrorLocalYCopy", "Number of YCopy errors",nbLocalBoard,0.5,(Float_t)nbLocalBoard+0.5);
365   histoYCopyErr->GetXaxis()->SetTitle(boardName.Data());
366   histoYCopyErr->GetYaxis()->SetTitle(errorAxisTitle.Data());
367   // Copy of previous histo for scaling purposes
368   TH1F* histoYCopyErrTest = (TH1F*)histoYCopyErr->Clone("hTriggerErrorLocalYCopyTest");
369   histoYCopyErrTest->SetTitle("Number of YCopy tested");
370   // Copy of previous histo for scaling purposes
371   TH1F* histoYCopyErrNorm = (TH1F*)histoYCopyErr->Clone("hTriggerErrorLocalYCopyNorm");
372   histoYCopyErrNorm->SetTitle("% of YCopy errors");
373   // Adding both histos after cloning to avoid problems with the expert bit
374   Add2RawsList(histoYCopyErr,     AliMUONQAIndices::kTriggerErrorLocalYCopy,     expert, !image, !saveCorr);
375   Add2RawsList(histoYCopyErrTest, AliMUONQAIndices::kTriggerErrorLocalYCopyTest, expert, !image, !saveCorr);
376   Add2RawsList(histoYCopyErrNorm, AliMUONQAIndices::kTriggerErrorLocalYCopyNorm, expert, !image, !saveCorr);
377
378   TH1F* histoROerr = new TH1F("hTriggerReadoutNumOfErrors","Trigger Read-Out total errors", AliMUONQAIndices::kNtrigStructErrorBins, -0.5, (Float_t)AliMUONQAIndices::kNtrigStructErrorBins-0.5);
379   histoROerr->GetYaxis()->SetTitle("Fraction of errors");
380   histoROerr->SetFillColor(kBlue);
381   for (int ibin=0;ibin<AliMUONQAIndices::kNtrigStructErrorBins;ibin++){
382     histoROerr->GetXaxis()->SetBinLabel(ibin+1,readoutErrNames[ibin]);
383   }
384   // Copy of previous histo for scaling purposes
385   TH1F* histoROerrNorm = (TH1F*)histoROerr->Clone("hTriggerReadoutErrors");
386   histoROerrNorm->SetTitle("Trigger Read-Out errors");
387   histoROerrNorm->SetOption("bar2");
388   histoROerrNorm->GetYaxis()->SetTitle("% of errors per event");
389   // Adding both histos after cloning to avoid problems with the expert bit
390   Add2RawsList(histoROerr,     AliMUONQAIndices::kTriggerReadOutErrors,      expert, !image, !saveCorr);
391   Add2RawsList(histoROerrNorm, AliMUONQAIndices::kTriggerReadOutErrorsNorm, !expert,  image, !saveCorr);
392
393   TH1F* histoGlobalMult = new TH1F("hTriggerGlobalOutMultiplicity","Trigger global outputs multiplicity", 6, -0.5, 6.-0.5);
394   histoGlobalMult->GetYaxis()->SetTitle("Number of triggers"); 
395   histoGlobalMult->GetXaxis()->SetTitle("Global output");
396   for (int ibin=0;ibin<6;ibin++){
397     histoGlobalMult->GetXaxis()->SetBinLabel(ibin+1,globalXaxisName[ibin]);
398   }        
399   histoGlobalMult->SetFillColor(kBlue);
400   // Copy of previous histo for scaling purposes
401   TH1F* histoGlobalMultNorm = (TH1F*)histoGlobalMult->Clone("hTriggerGlobalOutMultiplicityPerEvt");
402   histoGlobalMultNorm->SetTitle("Trigger global outputs multiplicity per event");
403   histoGlobalMultNorm->SetOption("bar2");
404   histoGlobalMultNorm->SetBarWidth(0.5);
405   histoGlobalMultNorm->SetBarOffset(0.25);
406   histoGlobalMultNorm->GetYaxis()->SetTitle("Triggers per event");
407   // Adding both histos after cloning to avoid problems with the expert bit
408   Add2RawsList(histoGlobalMult,     AliMUONQAIndices::kTriggerGlobalOutput,     expert, !image, !saveCorr);
409   Add2RawsList(histoGlobalMultNorm, AliMUONQAIndices::kTriggerGlobalOutputNorm, expert, !image, !saveCorr);
410
411   histo1D = new TH1F("hTriggerRawNAnalyzedEvents", "Number of analyzed events per specie", 1, 0.5, 1.5);
412   Int_t esindex = AliRecoParam::AConvert(CurrentEventSpecie());
413   histo1D->GetXaxis()->SetBinLabel(1, AliRecoParam::GetEventSpecieName(esindex));
414   histo1D->GetYaxis()->SetTitle("Number of analyzed events");
415   Add2RawsList(histo1D, AliMUONQAIndices::kTriggerRawNAnalyzedEvents, expert, !image, !saveCorr);
416
417   fNumberOf34Dec = new TH1F("hTriggerNumberOf34Dec", "hNumberOf34Dec",nbLocalBoard,0.5,(Float_t)nbLocalBoard+0.5);
418   fNumberOf34Dec->SetDirectory(0); // Detach histo from file to avoid double delete
419   fNumberOf44Dec = new TH1F("hTriggerNumberOf44Dec", "hNumberOf44Dec",nbLocalBoard,0.5,(Float_t)nbLocalBoard+0.5);
420   fNumberOf44Dec->SetDirectory(0); // Detach histo from file to avoid double delete
421 }
422
423 //__________________________________________________________________
424 void AliMUONTriggerQADataMakerRec::InitDigits() 
425 {
426   /// Initialized Digits spectra 
427   const Bool_t expert   = kTRUE ; 
428   const Bool_t image    = kTRUE ; 
429   
430   TH1I* h0 = new TH1I("hTriggerDigitsDetElem", "Detection element distribution in Digits;Detection element Id;Counts",  400, 1100, 1500); 
431   Add2DigitsList(h0, 0, !expert, image);
432
433
434 //____________________________________________________________________________ 
435 void AliMUONTriggerQADataMakerRec::InitRecPoints()
436 {
437         /// create Reconstructed Points histograms in RecPoints subdir for the
438         /// MUON Trigger subsystem.
439
440   const Bool_t expert   = kTRUE ; 
441   const Bool_t image    = kTRUE ; 
442
443   TH1F* histo1D = 0x0;
444
445   histo1D = new TH1F("hTriggerNAnalyzedEvents", "Number of analyzed events per specie", 1, 0.5, 1.5);
446   Int_t esindex = AliRecoParam::AConvert(CurrentEventSpecie());
447   histo1D->GetXaxis()->SetBinLabel(1, AliRecoParam::GetEventSpecieName(esindex));
448   histo1D->GetYaxis()->SetTitle("Number of analyzed events");
449   Add2RecPointsList(histo1D, AliMUONQAIndices::kTriggerNAnalyzedEvents, expert, !image);
450
451   histo1D = new TH1F("hTriggerTrippedChambers", "Trigger RPCs in trip", 418, 1100-0.5, 1417+0.5);
452   histo1D->GetXaxis()->SetTitle("DetElemId");
453   histo1D->GetYaxis()->SetTitle("# of trips");
454   histo1D->SetFillColor(kRed);
455   histo1D->SetLineColor(kRed);
456   Add2RecPointsList(histo1D, AliMUONQAIndices::kTriggerRPCtrips, !expert, image);
457
458   FillTriggerDCSHistos();       
459 }
460
461
462 //____________________________________________________________________________ 
463 void AliMUONTriggerQADataMakerRec::InitESDs()
464 {
465   /// Empty implementation
466 }
467
468 //____________________________________________________________________________
469 void AliMUONTriggerQADataMakerRec::MakeRaws(AliRawReader* rawReader)
470 {
471         /// make QA for rawdata trigger
472
473     AliCodeTimerAuto("",0);
474         
475     // Init Local/Regional/Global decision with fake values
476
477     UInt_t globalInput[4];
478     for (Int_t bit=0; bit<4; bit++){
479         globalInput[bit]=0;
480     }
481
482     //for (Int_t reg=0;reg<16;reg++){
483     //fTriggerOutputRegionalData[reg]=0;
484     //for (Int_t bit=0;bit<4;bit++){
485     //fTriggerInputGlobalDataLPt[reg][bit]=0;
486     //fTriggerInputGlobalDataHPt[reg][bit]=0;
487     //}
488     //}
489
490     AliMUONDigitStoreV2R digitStore;
491
492     AliMUONDigitStoreV2R digitStoreAll;
493     TArrayS xyPatternAll[2];
494     for(Int_t icath=0; icath<AliMpConstants::NofCathodes(); icath++){
495       xyPatternAll[icath].Set(AliMpConstants::NofTriggerChambers());
496       xyPatternAll[icath].Reset(1);
497     }
498     
499     AliMUONTriggerStoreV1 recoTriggerStore;
500
501     AliMUONTriggerStoreV1 inputTriggerStore;
502
503     AliMUONGlobalTrigger inputGlobalTrigger;
504
505     UShort_t maxNcounts = 0xFFFF;
506     
507     // Get trigger Local, Regional, Global in/outputs and scalers
508
509     Int_t loCircuit=0;
510     AliMpCDB::LoadDDLStore();
511
512     const AliMUONRawStreamTriggerHP::AliHeader*          darcHeader  = 0x0;
513     const AliMUONRawStreamTriggerHP::AliRegionalHeader*  regHeader   = 0x0;
514     const AliMUONRawStreamTriggerHP::AliLocalStruct*     localStruct = 0x0;
515
516     Int_t nDeadLocal = 0, nDeadRegional = 0, nDeadGlobal = 0, nNoisyStrips = 0;
517
518     // When a crate is not present, the loop on boards is not performed
519     // This should allow to correctly count the local boards
520     Int_t countNotifiedBoards = 0, countAllBoards = 0;
521
522     Bool_t containTriggerData = kFALSE;
523     AliMUONRawStreamTriggerHP rawStreamTrig(rawReader);
524     while (rawStreamTrig.NextDDL()) 
525     {
526       containTriggerData = kTRUE;
527
528       Bool_t scalerEvent =  rawReader->GetDataHeader()->GetL1TriggerMessage() & 0x1;
529
530       Bool_t fillScalerHistos = ( scalerEvent && 
531                                   ( GetRecoParam()->GetEventSpecie() == AliRecoParam::kCalib ) );
532
533       if ( scalerEvent != fillScalerHistos ) {
534         Int_t esindex = AliRecoParam::AConvert(CurrentEventSpecie());
535         AliWarning(Form("Scaler event found but event specie is %s. Scaler histos will not be filled", AliRecoParam::GetEventSpecieName(esindex)));
536       }
537
538       darcHeader = rawStreamTrig.GetHeaders();
539
540       if (darcHeader->GetGlobalFlag()){
541         if ( fillScalerHistos ) {
542           UInt_t nOfClocks = darcHeader->GetGlobalClock();
543           Double_t nOfSeconds = ((Double_t) nOfClocks) / 40e6; // 1 clock each 25 ns
544           ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerScalersTime))->Fill(1., nOfSeconds);
545         }
546
547         //Get Global datas
548         inputGlobalTrigger.SetFromGlobalResponse(darcHeader->GetGlobalOutput());
549         Bool_t resp[6] = {inputGlobalTrigger.PairUnlikeHpt(), inputGlobalTrigger.PairUnlikeLpt(),
550                           inputGlobalTrigger.PairLikeHpt(), inputGlobalTrigger.PairLikeLpt(),
551                           inputGlobalTrigger.SingleHpt(), inputGlobalTrigger.SingleLpt()}; 
552         for (Int_t bit=0; bit<6; bit++){
553           if ( ! resp[bit] ){
554             if ( fillScalerHistos )
555               nDeadGlobal++;
556           }
557           else
558             ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerGlobalOutput))->Fill(bit);
559         }
560
561         //for (Int_t Bit=0; Bit<32; Bit++){
562         //fTriggerInputGlobalDataLPt[Bit/4][Bit%4]=((darcHeader->GetGlobalInput(0)>>Bit)&1);
563         //fTriggerInputGlobalDataLPt[Bit/4+8][Bit%4]=((darcHeader->GetGlobalInput(1)>>Bit)&1);
564         //fTriggerInputGlobalDataHPt[Bit/4][Bit%4]=((darcHeader->GetGlobalInput(2)>>Bit)&1);
565         //fTriggerInputGlobalDataHPt[Bit/4+8][Bit%4]=((darcHeader->GetGlobalInput(3)>>Bit)&1);
566         //}
567
568         for (Int_t i=0; i<4; i++){
569           globalInput[i]=darcHeader->GetGlobalInput(i);
570         }
571       }
572
573       Int_t nReg = rawStreamTrig.GetRegionalHeaderCount();
574
575       for(Int_t iReg = 0; iReg < nReg ;iReg++)
576       {   //reg loop
577
578         //Int_t regId=rawStreamTrig.GetDDL()*8+iReg;
579
580         // crate info  
581           AliMpTriggerCrate* crate = AliMpDDLStore::Instance()->GetTriggerCrate(rawStreamTrig.GetDDL(), iReg);
582
583           regHeader =  rawStreamTrig.GetRegionalHeader(iReg);
584
585           //Get regional outputs -> not checked, hardware read-out doesn't work
586           //fTriggerOutputRegionalData[regId]=Int_t(regHeader->GetOutput());
587           // if ( ! fTriggerOutputRegionalData[regId] )
588           // nDeadRegional++;
589         Int_t nBoardsInReg = 0; // Not necessary when regional output will work
590
591         // loop over local structures
592         Int_t nLocal = regHeader->GetLocalStructCount();
593
594         for(Int_t iLocal = 0; iLocal < nLocal; iLocal++) 
595         {
596             
597             localStruct = regHeader->GetLocalStruct(iLocal);
598
599           // if card exist
600           if (!localStruct) continue;
601           
602           loCircuit = crate->GetLocalBoardId(localStruct->GetId());
603
604           if ( !loCircuit ) continue; // empty slot
605
606           AliMpLocalBoard* localBoard = AliMpDDLStore::Instance()->GetLocalBoard(loCircuit, false);
607
608           nBoardsInReg++; // Not necessary when regional output will work
609           countAllBoards++;
610
611           // skip copy cards
612           if( !localBoard->IsNotified()) 
613             continue;
614
615           AliMUONLocalTrigger inputLocalTrigger;
616           inputLocalTrigger.SetLocalStruct(loCircuit, *localStruct);
617           inputTriggerStore.Add(inputLocalTrigger);
618
619           countNotifiedBoards++;  
620
621           TArrayS xyPattern[2];   
622           localStruct->GetXPattern(xyPattern[0]);
623           localStruct->GetYPattern(xyPattern[1]);
624           fDigitMaker->TriggerDigits(loCircuit, xyPattern, digitStore);
625           if ( fillScalerHistos ) // Compute total number of strips
626             fDigitMaker->TriggerDigits(loCircuit, xyPatternAll, digitStoreAll);
627
628           //Get electronic Decisions from data
629
630           //Get regional inputs -> not checked, hardware read-out doesn't work
631           //fTriggerInputRegionalDataLPt[0][loCircuit]=Int_t(((regHeader->GetInput(0))>>(2*iLocal))&1);
632           //fTriggerInputRegionalDataLPt[1][loCircuit]=Int_t(((regHeader->GetInput(1))>>((2*iLocal)+1))&1);
633
634           //Get local in/outputs
635           if (Int_t(localStruct->GetDec())!=0){
636               ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggeredBoards))->Fill(loCircuit);
637           }
638           else if ( fillScalerHistos ){
639             nDeadLocal++;
640           }
641
642           // loop over strips
643           if ( fillScalerHistos ) {
644             Int_t cathode = localStruct->GetComptXY()%2;
645             for (Int_t ibitxy = 0; ibitxy < 16; ++ibitxy) {
646               if (ibitxy==0){
647                 AliDebug(AliQAv1::GetQADebugLevel(),"Filling trigger scalers");
648               }
649
650               UShort_t scalerVal[4] = {
651                 localStruct->GetXY1(ibitxy),
652                 localStruct->GetXY2(ibitxy),
653                 localStruct->GetXY3(ibitxy),
654                 localStruct->GetXY4(ibitxy)
655               };
656
657               for(Int_t ich=0; ich<AliMpConstants::NofTriggerChambers(); ich++){
658                 if ( scalerVal[ich] > 0 )
659                   ((TH2F*)GetRawsData(AliMUONQAIndices::kTriggerScalers + AliMpConstants::NofTriggerChambers()*cathode + ich))
660                     ->Fill(loCircuit, ibitxy, 2*(Float_t)scalerVal[ich]);
661
662                 if ( scalerVal[ich] >= maxNcounts )
663                   nNoisyStrips++;
664               } // loop on chamber
665             } // loop on strips
666           } // scaler event
667         } // iLocal
668         if ( nBoardsInReg == 0 )
669           nDeadRegional++; // Not necessary when regional output will work
670       } // iReg
671
672       Float_t readoutErrors[AliMUONQAIndices::kNtrigStructErrorBins] = {
673         ((Float_t)rawStreamTrig.GetLocalEoWErrors())/((Float_t)countAllBoards),
674         ((Float_t)rawStreamTrig.GetRegEoWErrors())/16.,
675         ((Float_t)rawStreamTrig.GetGlobalEoWErrors())/6.,
676         ((Float_t)rawStreamTrig.GetDarcEoWErrors())/2.
677       };
678     
679       for (Int_t ibin=0; ibin<AliMUONQAIndices::kNtrigStructErrorBins; ibin++){
680         if ( readoutErrors[ibin] > 0 )
681           ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerReadOutErrors))->Fill(ibin, readoutErrors[ibin]);
682       }
683     } // NextDDL
684
685     if ( ! containTriggerData ) return;
686     GetRawsData(AliMUONQAIndices::kTriggerRawNAnalyzedEvents)->Fill(1.);
687
688     nDeadLocal += AliMUONConstants::NTriggerCircuit() - countNotifiedBoards;
689     Int_t nStripsTot = digitStoreAll.GetSize();
690     if ( nStripsTot > 0 ) { // The value is != 0 only for scaler events
691       Float_t fraction[AliMUONQAIndices::kNtrigCalibSummaryBins] = {
692         ((Float_t)(nStripsTot - digitStore.GetSize())) / ((Float_t)nStripsTot),
693         //(Float_t)nDeadLocal / ((Float_t)countNotifiedBoards),
694         (Float_t)nDeadLocal / ((Float_t)AliMUONConstants::NTriggerCircuit()),
695         (Float_t)nDeadRegional / 16.,
696         (Float_t)nDeadGlobal / 6., // Number of bits of global response
697         (Float_t)nNoisyStrips / ((Float_t)nStripsTot),
698       };
699
700       for(Int_t ibin = 0; ibin < AliMUONQAIndices::kNtrigCalibSummaryBins; ibin++){
701         if ( fraction[ibin] > 0. )
702           ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerCalibSummary))->Fill(ibin, fraction[ibin]);
703       }
704     }
705
706   fTriggerProcessor->Digits2Trigger(digitStore,recoTriggerStore);
707
708   AliMUONGlobalTrigger* recoGlobalTriggerFromLocal;
709   recoGlobalTriggerFromLocal = recoTriggerStore.Global();
710
711   //Reconstruct Global decision from Global inputs
712   UChar_t recoResp = RawTriggerInGlobal2OutGlobal(globalInput);
713   AliMUONGlobalTrigger recoGlobalTriggerFromGlobal;
714   recoGlobalTriggerFromGlobal.SetFromGlobalResponse(recoResp);
715
716   // Compare data and reconstructed decisions and fill histos
717   RawTriggerMatchOutLocal(inputTriggerStore, recoTriggerStore);
718   //Fill ratio 44/34 histos
719   FillRatio4434Histos(fgkUpdateRatio4434);
720   //RawTriggerMatchOutLocalInRegional(); // Not tested, hardware read-out doesn't work
721   RawTriggerMatchOutGlobal(inputGlobalTrigger, recoGlobalTriggerFromGlobal, 'G');
722   // Global, reconstruction from Local inputs: compare data and reconstructed decisions and fill histos
723   RawTriggerMatchOutGlobal(inputGlobalTrigger, *recoGlobalTriggerFromLocal, 'L');
724   // Global, reconstruction from Global inputs: compare data and reconstructed decisions and fill histos
725 }
726
727 //__________________________________________________________________
728 void AliMUONTriggerQADataMakerRec::MakeDigits(TTree* digitsTree)         
729 {
730   /// makes data from Digits
731
732   // Do nothing in case of calibration event
733   if ( GetRecoParam()->GetEventSpecie() == AliRecoParam::kCalib ) return;
734
735   if (!fDigitStore)
736     fDigitStore = AliMUONVDigitStore::Create(*digitsTree);
737
738   fDigitStore->Clear();
739   fDigitStore->Connect(*digitsTree, false);
740   digitsTree->GetEvent(0);
741   
742   TIter next(fDigitStore->CreateIterator());
743   
744   AliMUONVDigit* dig = 0x0;
745   
746   while ( ( dig = static_cast<AliMUONVDigit*>(next()) ) )
747     {
748     GetDigitsData(0)->Fill(dig->DetElemId());
749     }
750 }
751
752 //____________________________________________________________________________
753 void AliMUONTriggerQADataMakerRec::MakeRecPoints(TTree* /*clustersTree*/)
754 {
755   /// Fill histogram with total number of analyzed events for normalization purposes
756
757   // Do nothing in case of calibration event
758   if ( GetRecoParam()->GetEventSpecie() == AliRecoParam::kCalib ) return;
759         
760   GetRecPointsData(AliMUONQAIndices::kTriggerNAnalyzedEvents)->Fill(1.);
761 }
762
763 //____________________________________________________________________________
764 void AliMUONTriggerQADataMakerRec::MakeESDs(AliESDEvent* /*esd*/)
765 {  
766   /// Empty implementation
767 }
768
769
770 //____________________________________________________________________________ 
771 void AliMUONTriggerQADataMakerRec::DisplayTriggerInfo()
772 {
773   //
774   /// Display trigger information in a user-friendly way:
775   /// from local board and strip numbers to their position on chambers
776   //
777
778   AliMUONTriggerDisplay triggerDisplay;
779   
780   TH2F* histoStrips=0x0;
781   TH2F* histoDisplayStrips=0x0;
782   if ( GetRawsData(AliMUONQAIndices::kTriggerScalers) ) {
783     AliMUONTriggerDisplay::EDisplayOption displayOption = AliMUONTriggerDisplay::kNormalizeToArea;
784     for (Int_t iCath = 0; iCath < AliMpConstants::NofCathodes(); iCath++)
785       {    
786         for (Int_t iChamber = 0; iChamber < AliMpConstants::NofTriggerChambers(); iChamber++)
787           {
788             histoStrips = (TH2F*)GetRawsData(AliMUONQAIndices::kTriggerScalers + AliMpConstants::NofTriggerChambers()*iCath + iChamber);
789
790             if(histoStrips->GetEntries()==0) continue; // No events found => No need to display
791
792             histoDisplayStrips = (TH2F*)GetRawsData(AliMUONQAIndices::kTriggerScalersDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber);
793
794             triggerDisplay.FillDisplayHistogram(histoStrips, histoDisplayStrips,
795                                                 AliMUONTriggerDisplay::kDisplayStrips, iCath, iChamber, displayOption);
796
797             Float_t scaleValue = ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerScalersTime))->GetBinContent(1);
798             if(scaleValue>0.) histoDisplayStrips->Scale(1./scaleValue);
799           } // iChamber
800       } // iCath
801   }
802
803   if ( GetRawsData(AliMUONQAIndices::kTriggeredBoards) ){
804     TH1F* histoBoards = (TH1F*)GetRawsData(AliMUONQAIndices::kTriggeredBoards);
805     TH2F* histoDisplayBoards = (TH2F*)GetRawsData(AliMUONQAIndices::kTriggerBoardsDisplay);
806     triggerDisplay.FillDisplayHistogram(histoBoards, histoDisplayBoards, AliMUONTriggerDisplay::kDisplayBoards, 0, 0);
807     Float_t scaleValue = GetRawsData(AliMUONQAIndices::kTriggerRawNAnalyzedEvents)->GetBinContent(1);
808     if(scaleValue>0.) histoDisplayBoards->Scale(1./scaleValue);
809   }
810 }
811
812
813 //_____________________________________________________________________________
814 Bool_t 
815 AliMUONTriggerQADataMakerRec::FillTriggerDCSHistos()
816 {
817   /// Get HV and currents values for one trigger chamber
818   
819   AliCodeTimerAuto("",0);
820
821   TMap* triggerDcsMap = fCalibrationData->TriggerDCS();
822
823   if ( !triggerDcsMap ) 
824   {
825     AliError("Cannot fill DCS histos, as triggerDcsMap is NULL");
826     return kFALSE;
827   }
828
829   const Double_t kMaxDelay = 3.;
830   const Double_t kMaxVariation = 25.; // Volts
831   const Int_t kDefaultNpoints = 200;
832   Double_t scaleFactor = 1./1000.;
833
834   Bool_t error = kFALSE;
835   Bool_t expert   = kTRUE;
836   Bool_t image    = kTRUE;
837
838   AliMpDEIterator deIt;
839   
840   AliMpDCSNamer triggerDcsNamer("TRIGGER");
841
842   TH2F* currHisto = 0x0;
843   Int_t histoIndex = 0;
844   TString histoName, histoTitle;
845
846   TArrayD axisSlat(18+1);
847   for(Int_t islat=0; islat<=18; islat++){
848     axisSlat[islat] = -0.5 + (Float_t)islat;
849   }
850
851   TArrayD axisTimeAux[4], axisTime[4], axisTimeDE(kDefaultNpoints);
852   TArrayI index[4], npoints(4);
853
854   // Build axis of times
855   npoints.Reset();
856   for(Int_t ich=0; ich<4; ich++){
857     axisTimeAux[ich].Set(kDefaultNpoints);
858     axisTimeAux[ich].Reset(-1.);
859   }
860
861   deIt.First();
862   while ( !deIt.IsDone() )
863   {
864     Int_t detElemId = deIt.CurrentDEId();
865     TObjArray* values = GetDCSValues(AliMpDCSNamer::kDCSHV, detElemId, triggerDcsMap, triggerDcsNamer);
866
867     if ( values ) {
868
869       AliDebug(AliQAv1::GetQADebugLevel(), Form("DetElemId %i", detElemId));
870
871       axisTimeDE.Reset(-1.);
872
873       TIter next(values);
874       AliDCSValue* val = 0x0;
875       Double_t previousVal = -999.;
876       Int_t npointsde = 0;
877       while ( ( val = static_cast<AliDCSValue*>(next()) ) )
878       {
879         if ( npointsde + 1 > kDefaultNpoints ) {
880           axisTimeDE.Set(npointsde + 1);
881         }
882
883         Double_t currVal = val->GetFloat();
884         Double_t currTime = (Double_t)val->GetTimeStamp();
885         if (npointsde > 0 ){
886           if ( TMath::Abs( currVal - previousVal ) < kMaxVariation && 
887                TMath::Abs( currTime - axisTimeDE[npointsde-1] ) < 40 ) continue;
888         }
889
890         axisTimeDE[npointsde] = currTime;
891         previousVal = currVal;
892         npointsde++;
893       } // loop on values
894
895       //      AliDebug(AliQAv1::GetQADebugLevel(), Form("Adding DE point %2i  (%2i)  %.2f  (%i)\n", previousBin, npointsde, axisTimeDE[previousBin], nTimesPerBin));
896
897       Int_t iChamber = AliMpDEManager::GetChamberId(detElemId);
898       Int_t ich = iChamber - AliMpConstants::NofTrackingChambers();
899
900       for(Int_t ipde=0; ipde<npointsde; ipde++){
901
902         if ( npoints[ich] + 1 > kDefaultNpoints ) {
903           axisTimeAux[ich].Set(npoints[ich] + 1);
904         }
905
906         for(Int_t ipoint = 0; ipoint < axisTimeAux[ich].GetSize(); ipoint++){
907           if (axisTimeAux[ich][ipoint] < 0.) {
908             axisTimeAux[ich][ipoint] = axisTimeDE[ipde];
909             npoints[ich]++;
910             AliDebug(AliQAv1::GetQADebugLevel(), Form("Adding point %2i  %.0f\n", ipoint, axisTimeAux[ich][ipoint]));
911             break;
912           }
913           if ( TMath::Abs( axisTimeDE[ipde] - axisTimeAux[ich][ipoint]) < kMaxDelay ) {
914             axisTimeAux[ich][ipoint] = TMath::Min(axisTimeAux[ich][ipoint], axisTimeDE[ipde]);
915             break;
916           }
917         } // loop on points
918       } // loop on reorganized values
919
920     } // if ( values ) 
921     deIt.Next();
922   } // loop on DetElemId
923
924   for(Int_t ich=0; ich<4; ich++){
925     axisTimeAux[ich].Set(npoints[ich]);
926     index[ich].Set(npoints[ich]);
927     TMath::Sort(npoints[ich], axisTimeAux[ich].GetArray(), index[ich].GetArray(), kFALSE);
928
929     axisTime[ich].Set(npoints[ich]+1);
930     for(Int_t ipoint = 0; ipoint < axisTimeAux[ich].GetSize(); ipoint++){
931       axisTime[ich][ipoint] = axisTimeAux[ich][index[ich][ipoint]];
932     }
933     Double_t minStartEndWidth = 0.1 * (axisTime[ich][npoints[ich]-1] - axisTime[ich][0]);
934     axisTime[ich][npoints[ich]] = axisTime[ich][npoints[ich]-1] + minStartEndWidth;
935     if ( npoints[ich] >= 1)
936       if ( axisTime[ich][1] - axisTime[ich][0] < minStartEndWidth )
937         axisTime[ich][0] = axisTime[ich][1] - minStartEndWidth;
938   }
939
940
941   // Loop again on detection elements: create and fill histos
942   deIt.First();
943   while ( !deIt.IsDone() )
944   {
945     Int_t detElemId = deIt.CurrentDEId();
946     TObjArray* values = GetDCSValues(AliMpDCSNamer::kDCSHV, detElemId, triggerDcsMap, triggerDcsNamer);
947
948     if ( values ) {
949       Int_t iChamber = AliMpDEManager::GetChamberId(detElemId);
950       Int_t ich = iChamber - AliMpConstants::NofTrackingChambers();
951
952       histoIndex = AliMUONQAIndices::kTriggerRPChv + ich;
953       histoName = Form("hTriggerRPCHVChamber%i", 11+ich);
954       histoTitle = Form("Chamber %i: RPC HV (kV)", 11+ich);
955
956       currHisto = (TH2F*)GetRecPointsData(histoIndex);
957
958       if(!currHisto){
959         currHisto  = new TH2F(histoName.Data(), histoTitle.Data(),
960                               npoints[ich], axisTime[ich].GetArray(),
961                               18, axisSlat.GetArray());
962         currHisto->GetXaxis()->SetTitle("Time");
963         currHisto->GetXaxis()->SetTimeDisplay(1);
964         //currHisto->GetXaxis()->SetTimeFormat("%d%b%y %H:%M:%S");
965         currHisto->GetXaxis()->SetLabelSize(0.03);
966         currHisto->GetYaxis()->SetTitle("RPC");
967         currHisto->SetOption("TEXT45COLZ");
968         Add2RecPointsList(currHisto, histoIndex, expert, !image);
969       }
970
971       Int_t slat = detElemId%100;
972       Int_t slatBin = currHisto->GetYaxis()->FindBin(slat);
973
974       TIter next(values);
975       AliDCSValue* val = 0x0;
976       Double_t sumValuesPerBin = 0.;
977       Int_t nValuesPerBin = 0;
978       Int_t previousBin = -1;
979       Double_t previousTime = -1., previousVal = -999., sumVal = 0., sumTime = 0.;
980       Bool_t isTrip = kFALSE;
981       Int_t nPointsForSlope = 0;
982       while ( ( val = static_cast<AliDCSValue*>(next()) ) )
983       {
984         Double_t currTime = (Double_t)val->GetTimeStamp();
985         Int_t currentBin = currHisto->GetXaxis()->FindBin(currTime+0.5);
986         Double_t currVal = val->GetFloat();
987         Double_t deltaVal = currVal - previousVal;
988         Bool_t isRepeated = kFALSE;
989         if ( previousTime > 0 ){
990           isRepeated = ( TMath::Abs( currVal - previousVal ) < kMaxVariation && 
991                          TMath::Abs( currTime - previousTime ) < 40 );
992
993           // Check for trips
994           sumTime += currTime - previousTime;
995           sumVal += deltaVal;
996           nPointsForSlope++;
997
998           if ( sumTime > 0. && nPointsForSlope >= 3 ){
999             Double_t slope = sumVal / sumTime;
1000             if ( slope < -10. ) // going down of more than 10V/s
1001               isTrip = kTRUE;
1002           }
1003
1004           if ( deltaVal * sumVal < 0. ) {
1005             sumTime = 0.;
1006             sumVal = 0.;
1007             nPointsForSlope = 0;
1008           }
1009         }
1010
1011         if ( ! isRepeated ) {
1012           if ( currentBin != previousBin ) {
1013             if ( previousBin >= 0 ) {
1014               currHisto->SetBinContent(previousBin, slatBin, scaleFactor*sumValuesPerBin/((Double_t)nValuesPerBin));
1015               sumValuesPerBin = 0.;
1016               nValuesPerBin = 0;
1017             }
1018             previousBin = currentBin;
1019           }
1020         }
1021           
1022         sumValuesPerBin += currVal;
1023         nValuesPerBin++;
1024         previousTime = currTime;
1025         previousVal = currVal;
1026       } // loop on values
1027       currHisto->SetBinContent(previousBin, slatBin, scaleFactor*sumValuesPerBin/((Double_t)nValuesPerBin)); // Fill last value
1028       if ( isTrip ) ((TH1F*)GetRecPointsData(AliMUONQAIndices::kTriggerRPCtrips))->Fill(detElemId);
1029     } // if ( values ) 
1030     deIt.Next();
1031   } // loop on detElem
1032   return error;
1033 }
1034
1035
1036 //____________________________________________________________________________ 
1037 TObjArray* 
1038 AliMUONTriggerQADataMakerRec::GetDCSValues(Int_t iMeas, Int_t detElemId,
1039                                            TMap* triggerDcsMap, AliMpDCSNamer& triggerDcsNamer)
1040 {
1041   //
1042   /// Get values of DCS data points from the map
1043   //
1044
1045   if ( AliMpDEManager::GetStationType(detElemId) != AliMp::kStationTrigger) return 0x0;
1046
1047   TString currAlias = triggerDcsNamer.DCSChannelName(detElemId, 0, iMeas);
1048
1049   TPair* triggerDcsPair = static_cast<TPair*>(triggerDcsMap->FindObject(currAlias.Data()));
1050
1051   if (!triggerDcsPair)
1052   {
1053     AliError(Form("Did not find expected alias (%s) for DE %d\n",
1054                   currAlias.Data(),detElemId));
1055     return 0x0;
1056   }
1057
1058   TObjArray* values = static_cast<TObjArray*>(triggerDcsPair->Value());
1059   if (!values)
1060   {
1061     AliError(Form("Could not get values for alias %s\n",currAlias.Data()));
1062     return 0x0;
1063   }
1064
1065   return values;
1066 }
1067
1068
1069 //____________________________________________________________________________ 
1070 UChar_t AliMUONTriggerQADataMakerRec::RawTriggerInGlobal2OutGlobal(UInt_t globalInput[4])
1071 {
1072   //
1073   /// Reconstruct Global Trigger decision using Global Inputs
1074   //
1075
1076     AliCodeTimerAuto("",0);
1077
1078     AliMUONGlobalCrateConfig* globalConfig = fCalibrationData->GlobalTriggerCrateConfig();
1079
1080     AliMUONGlobalTriggerBoard globalTriggerBoard;
1081     globalTriggerBoard.Reset();
1082     for (Int_t i = 0; i < 4; i++) {
1083         globalTriggerBoard.Mask(i,globalConfig->GetGlobalMask(i));
1084     }
1085
1086     globalTriggerBoard.RecomputeRegional(globalInput);
1087     globalTriggerBoard.Response();
1088     return globalTriggerBoard.GetResponse();
1089
1090 }
1091
1092 //____________________________________________________________________________ 
1093 void AliMUONTriggerQADataMakerRec::RawTriggerMatchOutLocal(AliMUONVTriggerStore& inputTriggerStore,
1094                                                            AliMUONVTriggerStore& recoTriggerStore)
1095 {
1096   //
1097   /// Match data and reconstructed Local Trigger decision
1098   //
1099
1100   AliCodeTimerAuto("",0);
1101
1102   Bool_t skipBoard[234];
1103   memset(skipBoard,0,AliMUONConstants::NTriggerCircuit()*sizeof(Bool_t));
1104
1105   Bool_t errorInYCopy = kFALSE;
1106
1107   // First search for YCopy errors.
1108   Int_t loCircuit = -1;
1109   TIter next(recoTriggerStore.CreateLocalIterator());
1110   AliMUONLocalTrigger* recoLocalTrigger, *inputLocalTrigger;
1111   while ( ( recoLocalTrigger = static_cast<AliMUONLocalTrigger*>(next()) ) )
1112   {  
1113     loCircuit = recoLocalTrigger->LoCircuit();
1114     Int_t iboard = loCircuit - 1;
1115
1116     ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerErrorLocalYCopyTest))->Fill(loCircuit);
1117   
1118     inputLocalTrigger = inputTriggerStore.FindLocal(loCircuit);
1119
1120     Int_t recoTrigPattern[4]  = {recoLocalTrigger->GetY1Pattern(), recoLocalTrigger->GetY2Pattern(), recoLocalTrigger->GetY3Pattern(), recoLocalTrigger->GetY4Pattern()};
1121     Int_t inputTrigPattern[4] = {inputLocalTrigger->GetY1Pattern(), inputLocalTrigger->GetY2Pattern(), inputLocalTrigger->GetY3Pattern(), inputLocalTrigger->GetY4Pattern()};
1122
1123     AliMpLocalBoard* localBoardMp = AliMpDDLStore::Instance()->GetLocalBoard(loCircuit); // get local board object for switch value
1124
1125     Bool_t errorInCopyBoard = kFALSE;
1126     for(Int_t ich=0; ich<4; ich++){
1127       if ( recoTrigPattern[ich] != inputTrigPattern[ich] ){
1128         skipBoard[iboard] = kTRUE;
1129         if ( ich >=2 ){
1130           if ( localBoardMp->GetSwitch(AliMpLocalBoard::kOR0) )
1131             skipBoard[iboard+1] = kTRUE;
1132           if ( localBoardMp->GetSwitch(AliMpLocalBoard::kOR1) )
1133             skipBoard[iboard-1] = kTRUE;
1134         }
1135         errorInCopyBoard = kTRUE;
1136         errorInYCopy = kTRUE;
1137       }
1138     } // loop on chambers
1139     if ( errorInCopyBoard )
1140       ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerErrorLocalYCopy))->Fill(loCircuit);    
1141   } // loop on local boards
1142
1143   if (errorInYCopy)
1144     ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerErrorSummary))->Fill(AliMUONQAIndices::kAlgoLocalYCopy);
1145   
1146   Bool_t errorInXPosDev = kFALSE;
1147   Bool_t errorInYPosTrigY = kFALSE;
1148   Bool_t errorInLUT = kFALSE;
1149
1150   next.Reset();
1151   Bool_t respBendPlane, respNonBendPlane;
1152   while ( ( recoLocalTrigger = static_cast<AliMUONLocalTrigger*>(next()) ) )
1153   {  
1154     loCircuit = recoLocalTrigger->LoCircuit();
1155     Int_t iboard = loCircuit - 1;
1156   
1157     // Fill ratio 44/34 histos (if not scaler event)
1158     if ( GetRecoParam()->GetEventSpecie() != AliRecoParam::kCalib ) {
1159       Bool_t is34 = ( recoLocalTrigger->GetLoDecision() != 0 );
1160       Bool_t is44 = fTriggerProcessor->ModifiedLocalResponse(loCircuit, respBendPlane, respNonBendPlane, kTRUE);
1161       if ( is34 ) fNumberOf34Dec->Fill(loCircuit);
1162       if ( is44 ) fNumberOf44Dec->Fill(loCircuit);
1163
1164       if ( is44 && ! is34 )
1165         AliWarning("Event satisfies the 4/4 conditions but not the 3/4");
1166     }
1167     
1168     inputLocalTrigger = inputTriggerStore.FindLocal(loCircuit);
1169
1170     if ( recoLocalTrigger->LoStripX() != inputLocalTrigger->LoStripX() ) {
1171       ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerErrorLocalXPos))->Fill(loCircuit);
1172       errorInXPosDev = kTRUE;
1173     }
1174
1175     if ( recoLocalTrigger->GetDeviation() != inputLocalTrigger->GetDeviation() ) {
1176       ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerErrorLocalDev))->Fill(loCircuit);
1177       errorInXPosDev = kTRUE;
1178     }
1179
1180     // Skip following checks in case we previously found YCopy error and YPos or trigY errors
1181     if ( (!skipBoard[iboard]) || ( (recoLocalTrigger->LoStripY() == inputLocalTrigger->LoStripY()) && (recoLocalTrigger->LoTrigY() == inputLocalTrigger->LoTrigY())) ) {
1182         
1183         if ( recoLocalTrigger->GetLoDecision() != inputLocalTrigger->GetLoDecision() ) {
1184             ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerErrorLocalTriggerDec))->Fill(loCircuit);
1185         }
1186         
1187         // Test Hpt and LPT
1188         Int_t recoLut[2]  = { recoLocalTrigger->LoLpt(),  recoLocalTrigger->LoHpt() };
1189         Int_t inputLut[2] = {inputLocalTrigger->LoLpt(), inputLocalTrigger->LoHpt() };
1190         Int_t currIndex[2][2] = {{AliMUONQAIndices::kTriggerErrorLocalLPtLSB, AliMUONQAIndices::kTriggerErrorLocalLPtMSB},
1191                                  {AliMUONQAIndices::kTriggerErrorLocalHPtMSB, AliMUONQAIndices::kTriggerErrorLocalHPtMSB}};
1192         for (Int_t ilut=0; ilut<2; ilut++){
1193             Int_t bitDiff = recoLut[ilut]^inputLut[ilut];
1194             if ( bitDiff == 0 ) continue;
1195             for (Int_t ibit=0; ibit<2; ibit++){
1196                 Bool_t isBitDifferent = (bitDiff>>ibit)&1;
1197                 if ( isBitDifferent ){
1198                     ((TH1F*)GetRawsData(currIndex[ilut][ibit]))->Fill(loCircuit);
1199                     errorInLUT = kTRUE;
1200                 }
1201             }
1202         }
1203     }
1204     
1205  
1206     // Skip following checks in case we previously found YCopy errors
1207     if ( skipBoard[iboard] ) continue;
1208
1209     if ( recoLocalTrigger->LoStripY() != inputLocalTrigger->LoStripY() ) {
1210       ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerErrorLocalYPos))->Fill(loCircuit);
1211       errorInYPosTrigY = kTRUE;
1212     }
1213
1214     if ( recoLocalTrigger->LoTrigY() != inputLocalTrigger->LoTrigY()  ) {
1215       ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerErrorLocalTrigY))->Fill(loCircuit); 
1216       errorInYPosTrigY = kTRUE;
1217     }
1218   } // loop on local boards
1219
1220   if (errorInXPosDev)
1221     ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerErrorSummary))->Fill(AliMUONQAIndices::kAlgoLocalX);
1222
1223   if (errorInLUT)
1224     ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerErrorSummary))->Fill(AliMUONQAIndices::kAlgoLocalLUT);
1225
1226   if (errorInYPosTrigY)
1227     ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerErrorSummary))->Fill(AliMUONQAIndices::kAlgoLocalY);
1228
1229 }
1230 /*
1231 //____________________________________________________________________________ 
1232 void AliMUONTriggerQADataMakerRec::RawTriggerMatchOutLocalInRegional()
1233 {
1234   //
1235   /// Match Local outputs and Regional inputs
1236   /// Not tested, hardware read-out doesn't work
1237   //
1238
1239     for (int localId=1;localId<235;localId++){
1240         if(fTriggerOutputLocalDataLPtDec[0][localId]!=fTriggerInputRegionalDataLPt[0][localId]){
1241             ((TH1F*)GetRawsData(kTriggerErrorLocal2RegionalLPtLSB))->Fill(localId);
1242         }
1243         if(fTriggerOutputLocalDataLPtDec[1][localId]!=fTriggerInputRegionalDataLPt[1][localId]){
1244             ((TH1F*)GetRawsData(kTriggerErrorLocal2RegionalLPtMSB))->Fill(localId);
1245         }
1246         if(fTriggerOutputLocalDataHPtDec[0][localId]!=fTriggerInputRegionalDataHPt[0][localId]){
1247             ((TH1F*)GetRawsData(kTriggerErrorLocal2RegionalHPtLSB))->Fill(localId);
1248         }
1249         if(fTriggerOutputLocalDataHPtDec[1][localId]!=fTriggerInputRegionalDataHPt[1][localId]){
1250             ((TH1F*)GetRawsData(kTriggerErrorLocal2RegionalHPtMSB))->Fill(localId);
1251         }
1252     }
1253 }
1254 */
1255
1256
1257 //____________________________________________________________________________ 
1258 void AliMUONTriggerQADataMakerRec::RawTriggerMatchOutGlobal(AliMUONGlobalTrigger& inputGlobalTrigger, 
1259                                                                         AliMUONGlobalTrigger& recoGlobalTrigger, 
1260                                                                         Char_t histo)
1261 {
1262   //
1263   /// Match data and reconstructed Global Trigger decision for a reconstruction from Global inputs.
1264   /// histo='G': fill FromGlobalInput histo='F': fill from Local input;
1265   //
1266
1267   if ( recoGlobalTrigger.GetGlobalResponse() == inputGlobalTrigger.GetGlobalResponse() )
1268     return;
1269   Int_t histoToFill;
1270   Int_t binToFill;
1271   
1272   if (histo=='G'){
1273       histoToFill=AliMUONQAIndices::kTriggerErrorOutGlobalFromInGlobal;
1274       binToFill=AliMUONQAIndices::kAlgoGlobalFromGlobal;
1275   }else{
1276       if (histo=='L'){
1277           histoToFill=AliMUONQAIndices::kTriggerErrorOutGlobalFromInLocal;
1278           binToFill=AliMUONQAIndices::kAlgoGlobalFromLocal;
1279       }else{
1280           AliWarning(Form("Global histos not filled, 3rd argument must be 'G' or 'F'"));
1281           return;
1282       } 
1283   }
1284
1285   ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerErrorSummary))->Fill(binToFill);
1286
1287   Bool_t inputResp[6] = {inputGlobalTrigger.PairUnlikeHpt(), inputGlobalTrigger.PairUnlikeLpt(),
1288                          inputGlobalTrigger.PairLikeHpt(), inputGlobalTrigger.PairLikeLpt(),
1289                          inputGlobalTrigger.SingleHpt(), inputGlobalTrigger.SingleLpt()};
1290
1291   Bool_t recoResp[6] = {recoGlobalTrigger.PairUnlikeHpt(), recoGlobalTrigger.PairUnlikeLpt(),
1292                         recoGlobalTrigger.PairLikeHpt(), recoGlobalTrigger.PairLikeLpt(),
1293                         recoGlobalTrigger.SingleHpt(), recoGlobalTrigger.SingleLpt()};
1294
1295   for (int bit=0;bit<6;bit++){
1296     if ( recoResp[bit] != inputResp[bit] )
1297       ((TH1F*)GetRawsData(histoToFill))->Fill(bit);
1298   }
1299 }
1300
1301 //____________________________________________________________________________ 
1302 void AliMUONTriggerQADataMakerRec::FillRatio4434Histos(Int_t evtInterval)
1303 {
1304   /// Fill ratio 44/34 histos
1305
1306   Int_t numEvent = Int_t(((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerRawNAnalyzedEvents))->GetBinContent(1));
1307
1308   // Fill every fgkUpdateRatio4434 events
1309   if (numEvent % evtInterval != 0)
1310     return;
1311
1312   Float_t totalNumberOf44 = fNumberOf44Dec->GetSumOfWeights();
1313   Float_t totalNumberOf34 = fNumberOf34Dec->GetSumOfWeights();
1314
1315   if ( totalNumberOf34 == 0 )
1316     return;
1317
1318   TH1F* histoAllEvents = (TH1F*)GetRawsData(AliMUONQAIndices::kTriggerRatio4434AllEvents);
1319   Int_t nbins =  histoAllEvents->GetNbinsX();
1320   Float_t maxBin = histoAllEvents->GetXaxis()->GetBinLowEdge(nbins+1);
1321
1322   if ( numEvent - maxBin < 1) return;
1323
1324   // Use the underflow and overflow to store the number of 34 and 44
1325   // in previous event
1326   Float_t previousNumOf34 = histoAllEvents->GetBinContent(0);
1327   Float_t previousNumOf44 = histoAllEvents->GetBinContent(nbins+1);
1328
1329   Float_t numOf34Update = totalNumberOf34 - previousNumOf34;
1330   Float_t numOf44Update = totalNumberOf44 - previousNumOf44;
1331
1332   if ( numOf34Update == 0 && numOf44Update == 0 )
1333     return;
1334
1335   Int_t newNbins = ( (Int_t)maxBin % fgkUpdateRatio4434 ) ? nbins : nbins+1;
1336   TString cloneName;
1337
1338   Int_t hIndex[2] = {AliMUONQAIndices::kTriggerRatio4434AllEvents, AliMUONQAIndices::kTriggerRatio4434SinceLastUpdate};
1339   
1340   for (Int_t ihisto=0; ihisto<2; ihisto++){
1341     TH1F* currHisto = (TH1F*)GetRawsData(hIndex[ihisto]);
1342     cloneName = Form("%sClone", currHisto->GetName());
1343     TArrayD newAxis(newNbins+1);
1344     for (Int_t ibin=0; ibin<newNbins; ibin++){
1345       newAxis[ibin] = currHisto->GetXaxis()->GetBinLowEdge(ibin+1);
1346     }
1347     newAxis[newNbins] = numEvent;
1348     TH1F* copyHisto = (TH1F*)currHisto->Clone(cloneName.Data());
1349     //currHisto->SetBins(newNbins, 0., fgkUpdateRatio4434*newNbins);
1350     currHisto->SetBins(newNbins, newAxis.GetArray());
1351     for (Int_t ibin=1; ibin<newNbins; ibin++){
1352       currHisto->SetBinContent(ibin, copyHisto->GetBinContent(ibin));
1353       currHisto->SetBinError(ibin, copyHisto->GetBinError(ibin));
1354     }
1355     delete copyHisto;
1356   }
1357
1358   Float_t ratio4434 = totalNumberOf44/totalNumberOf34;
1359   Float_t errorRatio4434 = ProtectedSqrt(totalNumberOf44*(1-ratio4434))/totalNumberOf34;
1360     
1361   histoAllEvents->SetBinContent(newNbins,ratio4434);
1362   histoAllEvents->SetBinError(newNbins,errorRatio4434);
1363
1364   Float_t ratio4434Update = 0.;
1365   Float_t errorRatio4434Update = 0.;
1366
1367   if(numOf34Update!=0){
1368     ratio4434Update = numOf44Update/numOf34Update;
1369     if ( numOf44Update > numOf34Update )
1370       AliWarning("Number of 4/4 is higher than number of 3/4");
1371     errorRatio4434Update = ProtectedSqrt(numOf44Update*(1-ratio4434Update))/numOf34Update;
1372   }
1373   
1374   ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerRatio4434SinceLastUpdate))->SetBinContent(newNbins,ratio4434Update);
1375   ((TH1F*)GetRawsData(AliMUONQAIndices::kTriggerRatio4434SinceLastUpdate))->SetBinError(newNbins,errorRatio4434Update);
1376
1377   histoAllEvents->SetBinContent(0,totalNumberOf34);
1378   histoAllEvents->SetBinContent(newNbins+1,totalNumberOf44);
1379
1380 }
1381