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