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