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