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