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