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