]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONQADataMakerRec.cxx
- Bug fix: in the creation/reading of the trigger DCS values, the HV and currents...
[u/mrichter/AliRoot.git] / MUON / AliMUONQADataMakerRec.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$
17
18 // --- MUON header files ---
19 #include "AliMUONQADataMakerRec.h"
20
21 #include "AliMUON2DMap.h"
22 #include "AliMUONCluster.h"  
23 #include "AliMUONConstants.h"  
24 #include "AliMUONDDLTrigger.h"
25 #include "AliMUONDarcHeader.h"
26 #include "AliMUONDigitMaker.h"
27 #include "AliMUONLocalStruct.h"
28 #include "AliMUONLocalTrigger.h"
29 #include "AliMUONRawStreamTracker.h"
30 #include "AliMUONRawStreamTrigger.h"
31 #include "AliMUONRegHeader.h"
32 #include "AliMUONTrackerDataMaker.h"
33 #include "AliMUONTriggerDisplay.h"
34 #include "AliMUONVCluster.h"
35 #include "AliMUONVClusterStore.h"
36 #include "AliMUONVDigit.h"
37 #include "AliMUONVDigitStore.h"
38 #include "AliMUONVTrackerData.h"
39 #include "AliMUONVTriggerStore.h"
40 #include "AliMUONTrack.h"
41 #include "AliMUONTrackParam.h"
42 #include "AliMUONESDInterface.h"
43 #include "AliMUONCalibrationData.h"
44 #include "AliMpBusPatch.h"
45 #include "AliMpCDB.h"
46 #include "AliMpConstants.h"
47 #include "AliMpDDLStore.h"
48 #include "AliMpDEIterator.h"
49 #include "AliMpDEManager.h"
50 #include "AliMpLocalBoard.h"
51 #include "AliMpStationType.h"
52 #include "AliMpTriggerCrate.h"
53 #include "AliMpDCSNamer.h"
54 #include "AliRawEventHeaderBase.h"
55
56 // --- AliRoot header files ---
57 #include "AliCDBManager.h"
58 #include "AliCDBStorage.h"
59 #include "AliESDEvent.h"
60 #include "AliESDMuonTrack.h"
61 #include "AliESDMuonCluster.h"
62 #include "AliLog.h"
63 #include "AliRawReader.h"
64 #include "AliQAChecker.h"
65 #include "AliCodeTimer.h"
66 #include "AliDCSValue.h"
67 #include "AliMUONVDigit.h"
68
69 // --- ROOT system ---
70 #include <TClonesArray.h>
71 #include <TFile.h> 
72 #include <TH1F.h> 
73 #include <TH1I.h> 
74 #include <TH2F.h>
75 #include <Riostream.h>
76
77 //-----------------------------------------------------------------------------
78 /// \class AliMUONQADataMakerRec
79 ///
80 /// MUON base class for quality assurance data (histo) maker
81 ///
82 /// \author C. Finck, D. Stocco, L. Aphecetche
83
84 /// \cond CLASSIMP
85 ClassImp(AliMUONQADataMakerRec)
86 /// \endcond
87            
88 //____________________________________________________________________________ 
89 AliMUONQADataMakerRec::AliMUONQADataMakerRec() : 
90 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kMUON), "MUON Quality Assurance Data Maker"),
91 fDigitStore(0x0),
92 fTriggerStore(0x0),
93 fDigitMaker(0x0),
94 fClusterStore(0x0),
95 fTrackerDataMaker(0x0)
96 {
97     /// ctor
98         
99   AliDebug(AliQAv1::GetQADebugLevel(),"");
100
101         Ctor();
102 }
103
104 //____________________________________________________________________________ 
105 void
106 AliMUONQADataMakerRec::Ctor()
107 {
108         /// Init some members
109         fDigitStore = AliMUONVDigitStore::Create("AliMUONDigitStoreV1");
110         fDigitMaker = new AliMUONDigitMaker(kTRUE);
111 }
112
113 //____________________________________________________________________________ 
114 AliMUONQADataMakerRec::AliMUONQADataMakerRec(const AliMUONQADataMakerRec& qadm) :
115 AliQADataMakerRec(qadm),
116 fDigitStore(0x0),
117 fTriggerStore(0x0),
118 fDigitMaker(0x0),
119 fClusterStore(0x0),
120 fTrackerDataMaker(0x0)
121 {
122     ///copy ctor 
123
124   AliDebug(AliQAv1::GetQADebugLevel(),"");
125
126
127     SetName((const char*)qadm.GetName()) ; 
128     SetTitle((const char*)qadm.GetTitle()); 
129
130         // Do not copy the digit store and digit maker, but create its own ones
131         
132         Ctor();
133         
134 }
135
136 //__________________________________________________________________
137 AliMUONQADataMakerRec& AliMUONQADataMakerRec::operator = (const AliMUONQADataMakerRec& qadm )
138 {
139   /// Assignment operator
140
141   AliDebug(AliQAv1::GetQADebugLevel(),"");
142
143   // check assignment to self
144   if (this == &qadm) return *this;
145
146   this->~AliMUONQADataMakerRec();
147   new(this) AliMUONQADataMakerRec(qadm);
148   return *this;
149 }
150
151 //__________________________________________________________________
152 AliMUONQADataMakerRec::~AliMUONQADataMakerRec()
153 {
154     /// dtor
155   
156   AliDebug(AliQAv1::GetQADebugLevel(),"");
157
158   AliCodeTimerAuto("");
159   
160   delete fDigitStore;
161   delete fTriggerStore;
162   delete fDigitMaker;
163   delete fClusterStore;
164   delete fTrackerDataMaker;
165 }
166
167 //____________________________________________________________________________ 
168 void AliMUONQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray** list)
169 {
170   ///Detector specific actions at end of cycle
171   
172   AliCodeTimerAuto("");
173
174   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) 
175   {
176     if (! IsValidEventSpecie(specie, list)  ) 
177       continue ;
178     SetEventSpecie(AliRecoParam::ConvertIndex(specie)) ; 
179     if ( task == AliQAv1::kRAWS && fTrackerDataMaker ) 
180       {
181
182         if ( !GetRawsData(kTrackerBusPatchOccupancy) ) continue;
183         
184         TIter next(list[specie]);
185         TObject* o;
186         Bool_t alreadyThere(kFALSE);
187         while ( ( o = next() ) && !alreadyThere )
188           {
189             TString classname(o->ClassName());
190             if ( classname.Contains("TrackerData") ) alreadyThere = kTRUE;
191           }
192         if (!alreadyThere && fTrackerDataMaker) 
193           {
194           AliDebug(AliQAv1::GetQADebugLevel(), "Adding fTrackerDataMaker to the list of qa objects");
195             list[specie]->AddAt(fTrackerDataMaker->Data(),(Int_t)kTrackerData);
196           }
197           if ( fTrackerDataMaker ) 
198             {
199               TH1* hbp = GetRawsData(kTrackerBusPatchOccupancy);
200               hbp->Reset();
201               TIter nextBP(AliMpDDLStore::Instance()->CreateBusPatchIterator());
202               AliMpBusPatch* bp(0x0);
203               AliMUONVTrackerData* data = fTrackerDataMaker->Data();
204               Int_t occDim = 2;
205       
206               while ( ( bp = static_cast<AliMpBusPatch*>(nextBP())) )
207                 {
208                   Int_t busPatchId = bp->GetId();
209                   Int_t bin = hbp->FindBin(busPatchId);
210                   hbp->SetBinContent(bin,data->BusPatch(busPatchId,occDim));
211                 }
212             }
213       }
214     
215     // Normalize RecPoints histos
216     if ( task == AliQAv1::kRECPOINTS ) {
217       
218       if (!GetRecPointsData(kTrackerClusterChargePerChMean)) continue;
219
220       TH1* hTrackerClusterChargePerChMean = GetRecPointsData(kTrackerClusterChargePerChMean);
221       TH1* hTrackerClusterChargePerChSigma = GetRecPointsData(kTrackerClusterChargePerChSigma);
222       TH1* hTrackerClusterMultiplicityPerChMean = GetRecPointsData(kTrackerClusterMultiplicityPerChMean);
223       TH1* hTrackerClusterMultiplicityPerChSigma = GetRecPointsData(kTrackerClusterMultiplicityPerChSigma);
224       TH1* hTrackerClusterChargePerDEMean = GetRecPointsData(kTrackerClusterChargePerDEMean);
225       TH1* hTrackerClusterMultiplicityPerDEMean = GetRecPointsData(kTrackerClusterMultiplicityPerDEMean);
226       
227       // loop over chambers
228       for (Int_t iCh = 0; iCh < AliMUONConstants::NTrackingCh(); iCh++) {
229         
230         TH1* hTrackerClusterChargePerChamber = GetRecPointsData(kTrackerClusterChargePerChamber+iCh);
231         Double_t sigmaCharge = hTrackerClusterChargePerChamber->GetRMS();
232         hTrackerClusterChargePerChMean->SetBinContent(iCh+1, hTrackerClusterChargePerChamber->GetMean());
233         hTrackerClusterChargePerChMean->SetBinError(iCh+1, hTrackerClusterChargePerChamber->GetMeanError());
234         hTrackerClusterChargePerChSigma->SetBinContent(iCh+1, sigmaCharge);
235         hTrackerClusterChargePerChSigma->SetBinError(iCh+1, hTrackerClusterChargePerChamber->GetRMSError());
236         
237         TH1* hTrackerClusterMultiplicityPerChamber = GetRecPointsData(kTrackerClusterMultiplicityPerChamber+iCh);
238         Double_t sigmaSize = hTrackerClusterMultiplicityPerChamber->GetRMS();
239         hTrackerClusterMultiplicityPerChMean->SetBinContent(iCh+1, hTrackerClusterMultiplicityPerChamber->GetMean());
240         hTrackerClusterMultiplicityPerChMean->SetBinError(iCh+1, hTrackerClusterMultiplicityPerChamber->GetMeanError());
241         hTrackerClusterMultiplicityPerChSigma->SetBinContent(iCh+1, sigmaSize);
242         hTrackerClusterMultiplicityPerChSigma->SetBinError(iCh+1, hTrackerClusterMultiplicityPerChamber->GetRMSError());
243         
244         // loop over DE into chamber iCh
245         AliMpDEIterator it;
246         it.First(iCh);
247         while ( !it.IsDone()) {
248           
249           Int_t iDE = it.CurrentDEId();
250           
251           TH1* hTrackerClusterChargePerDE = GetRecPointsData(kTrackerClusterChargePerDE+iDE);
252           hTrackerClusterChargePerDEMean->SetBinContent(iDE+1, hTrackerClusterChargePerDE->GetMean());
253           Double_t nClusters = hTrackerClusterChargePerDE->GetEntries();
254           if (nClusters > 1) hTrackerClusterChargePerDEMean->SetBinError(iDE+1, sigmaCharge/TMath::Sqrt(nClusters));
255           else hTrackerClusterChargePerDEMean->SetBinError(iDE+1, hTrackerClusterChargePerChamber->GetXaxis()->GetXmax());
256           
257           TH1* hTrackerClusterMultiplicityPerDE = GetRecPointsData(kTrackerClusterMultiplicityPerDE+iDE);
258           hTrackerClusterMultiplicityPerDEMean->SetBinContent(iDE+1, hTrackerClusterMultiplicityPerDE->GetMean());
259           nClusters = hTrackerClusterMultiplicityPerDE->GetEntries();
260           if (nClusters > 1) hTrackerClusterMultiplicityPerDEMean->SetBinError(iDE+1, sigmaSize/TMath::Sqrt(nClusters));
261           else hTrackerClusterMultiplicityPerDEMean->SetBinError(iDE+1, hTrackerClusterMultiplicityPerChamber->GetXaxis()->GetXmax());
262           
263           it.Next();
264         }
265         
266       }
267       
268     }
269     
270     // Normalize ESD histos
271     if ( task == AliQAv1::kESDS ) {
272       
273       if (!GetESDsData(kESDnClustersPerTrack)) continue;
274       
275       Double_t nTracks = GetESDsData(kESDnClustersPerTrack)->GetEntries();
276       if (nTracks <= 0) continue;
277       
278       TH1* hESDnClustersPerCh = GetESDsData(kESDnClustersPerCh);
279       TH1* hESDnClustersPerDE = GetESDsData(kESDnClustersPerDE);
280       TH1* hESDClusterChargePerChMean = GetESDsData(kESDClusterChargePerChMean);
281       TH1* hESDClusterChargePerChSigma = GetESDsData(kESDClusterChargePerChSigma);
282       TH1* hESDClusterSizePerChMean = GetESDsData(kESDClusterSizePerChMean);
283       TH1* hESDClusterSizePerChSigma = GetESDsData(kESDClusterSizePerChSigma);
284       TH1* hESDResidualXPerChMean = GetESDsData(kESDResidualXPerChMean);
285       TH1* hESDResidualXPerChSigma = GetESDsData(kESDResidualXPerChSigma);
286       TH1* hESDResidualYPerChMean = GetESDsData(kESDResidualYPerChMean);
287       TH1* hESDResidualYPerChSigma = GetESDsData(kESDResidualYPerChSigma);
288       TH1* hESDLocalChi2XPerChMean = GetESDsData(kESDLocalChi2XPerChMean);
289       TH1* hESDLocalChi2YPerChMean = GetESDsData(kESDLocalChi2YPerChMean);
290       TH1* hESDLocalChi2PerChMean = GetESDsData(kESDLocalChi2PerChMean);
291       TH1* hESDClusterChargePerDE = GetESDsData(kESDClusterChargePerDE);
292       TH1* hESDClusterSizePerDE = GetESDsData(kESDClusterSizePerDE);
293       TH1* hESDResidualXPerDEMean = GetESDsData(kESDResidualXPerDEMean);
294       TH1* hESDResidualXPerDESigma = GetESDsData(kESDResidualXPerDESigma);
295       TH1* hESDResidualYPerDEMean = GetESDsData(kESDResidualYPerDEMean);
296       TH1* hESDResidualYPerDESigma = GetESDsData(kESDResidualYPerDESigma);
297       TH1* hESDLocalChi2XPerDEMean = GetESDsData(kESDLocalChi2XPerDEMean);
298       TH1* hESDLocalChi2YPerDEMean = GetESDsData(kESDLocalChi2YPerDEMean);
299       TH1* hESDLocalChi2PerDEMean = GetESDsData(kESDLocalChi2PerDEMean);
300       TH1* hESDnTotClustersPerCh = GetESDsData(kESDnTotClustersPerCh);
301       TH1* hESDnTotClustersPerDE = GetESDsData(kESDnTotClustersPerDE);
302       TH1* hESDnTotFullClustersPerDE = GetESDsData(kESDnTotFullClustersPerDE);
303       TH1* hESDSumClusterChargePerDE = GetESDsData(kESDSumClusterChargePerDE);
304       TH1* hESDSumClusterSizePerDE = GetESDsData(kESDSumClusterSizePerDE);
305       TH1* hESDSumResidualXPerDE = GetESDsData(kESDSumResidualXPerDE);
306       TH1* hESDSumResidualYPerDE = GetESDsData(kESDSumResidualYPerDE);
307       TH1* hESDSumResidualX2PerDE = GetESDsData(kESDSumResidualX2PerDE);
308       TH1* hESDSumResidualY2PerDE = GetESDsData(kESDSumResidualY2PerDE);
309       TH1* hESDSumLocalChi2XPerDE = GetESDsData(kESDSumLocalChi2XPerDE);
310       TH1* hESDSumLocalChi2YPerDE = GetESDsData(kESDSumLocalChi2YPerDE);
311       TH1* hESDSumLocalChi2PerDE = GetESDsData(kESDSumLocalChi2PerDE);
312       
313       hESDnClustersPerCh->Reset();
314       hESDnClustersPerDE->Reset();
315       hESDnClustersPerCh->Add(hESDnTotClustersPerCh, 1./nTracks);
316       hESDnClustersPerDE->Add(hESDnTotClustersPerDE, 1./nTracks);
317       
318       // loop over chambers
319       for (Int_t iCh = 0; iCh < AliMUONConstants::NTrackingCh(); iCh++) {
320         
321         TH1* hESDClusterChargeInCh = GetESDsData(kESDClusterChargeInCh+iCh);
322         Double_t sigmaCharge = hESDClusterChargeInCh->GetRMS();
323         hESDClusterChargePerChMean->SetBinContent(iCh+1, hESDClusterChargeInCh->GetMean());
324         hESDClusterChargePerChMean->SetBinError(iCh+1, hESDClusterChargeInCh->GetMeanError());
325         hESDClusterChargePerChSigma->SetBinContent(iCh+1, sigmaCharge);
326         hESDClusterChargePerChSigma->SetBinError(iCh+1, hESDClusterChargeInCh->GetRMSError());
327         
328         TH1* hESDClusterSizeInCh = GetESDsData(kESDClusterSizeInCh+iCh);
329         Double_t sigmaSize = hESDClusterSizeInCh->GetRMS();
330         hESDClusterSizePerChMean->SetBinContent(iCh+1, hESDClusterSizeInCh->GetMean());
331         hESDClusterSizePerChMean->SetBinError(iCh+1, hESDClusterSizeInCh->GetMeanError());
332         hESDClusterSizePerChSigma->SetBinContent(iCh+1, sigmaSize);
333         hESDClusterSizePerChSigma->SetBinError(iCh+1, hESDClusterSizeInCh->GetRMSError());
334         
335         TH1* hESDResidualXInCh = GetESDsData(kESDResidualXInCh+iCh);
336         Double_t sigmaResidualX = hESDResidualXInCh->GetRMS();
337         hESDResidualXPerChMean->SetBinContent(iCh+1, hESDResidualXInCh->GetMean());
338         hESDResidualXPerChMean->SetBinError(iCh+1, hESDResidualXInCh->GetMeanError());
339         hESDResidualXPerChSigma->SetBinContent(iCh+1, sigmaResidualX);
340         hESDResidualXPerChSigma->SetBinError(iCh+1, hESDResidualXInCh->GetRMSError());
341         
342         TH1* hESDResidualYInCh = GetESDsData(kESDResidualYInCh+iCh);
343         Double_t sigmaResidualY = hESDResidualYInCh->GetRMS();
344         hESDResidualYPerChMean->SetBinContent(iCh+1, hESDResidualYInCh->GetMean());
345         hESDResidualYPerChMean->SetBinError(iCh+1, hESDResidualYInCh->GetMeanError());
346         hESDResidualYPerChSigma->SetBinContent(iCh+1, sigmaResidualY);
347         hESDResidualYPerChSigma->SetBinError(iCh+1, hESDResidualYInCh->GetRMSError());
348         
349         TH1* hESDLocalChi2XInCh = GetESDsData(kESDLocalChi2XInCh+iCh);
350         Double_t sigmaLocalChi2X = hESDLocalChi2XInCh->GetRMS();
351         hESDLocalChi2XPerChMean->SetBinContent(iCh+1, hESDLocalChi2XInCh->GetMean());
352         hESDLocalChi2XPerChMean->SetBinError(iCh+1, hESDLocalChi2XInCh->GetMeanError());
353         
354         TH1* hESDLocalChi2YInCh = GetESDsData(kESDLocalChi2YInCh+iCh);
355         Double_t sigmaLocalChi2Y = hESDLocalChi2YInCh->GetRMS();
356         hESDLocalChi2YPerChMean->SetBinContent(iCh+1, hESDLocalChi2YInCh->GetMean());
357         hESDLocalChi2YPerChMean->SetBinError(iCh+1, hESDLocalChi2YInCh->GetMeanError());
358         
359         TH1* hESDLocalChi2InCh = GetESDsData(kESDLocalChi2InCh+iCh);
360         Double_t sigmaLocalChi2 = hESDLocalChi2InCh->GetRMS();
361         hESDLocalChi2PerChMean->SetBinContent(iCh+1, hESDLocalChi2InCh->GetMean());
362         hESDLocalChi2PerChMean->SetBinError(iCh+1, hESDLocalChi2InCh->GetMeanError());
363         
364         // loop over DE into chamber iCh
365         AliMpDEIterator it;
366         it.First(iCh);
367         while ( !it.IsDone()) {
368           
369           Int_t iDE = it.CurrentDEId();
370           
371           Double_t nClusters = hESDnTotClustersPerDE->GetBinContent(iDE+1);
372           if (nClusters > 1) {
373             
374             hESDClusterChargePerDE->SetBinContent(iDE+1, hESDSumClusterChargePerDE->GetBinContent(iDE+1)/nClusters);
375             hESDClusterChargePerDE->SetBinError(iDE+1, sigmaCharge/TMath::Sqrt(nClusters));
376             
377             Double_t meanResX = hESDSumResidualXPerDE->GetBinContent(iDE+1)/nClusters;
378             hESDResidualXPerDEMean->SetBinContent(iDE+1, meanResX);
379             hESDResidualXPerDEMean->SetBinError(iDE+1, sigmaResidualX/TMath::Sqrt(nClusters));
380             hESDResidualXPerDESigma->SetBinContent(iDE+1, TMath::Sqrt(hESDSumResidualX2PerDE->GetBinContent(iDE+1)/nClusters - meanResX*meanResX));
381             hESDResidualXPerDESigma->SetBinError(iDE+1, sigmaResidualX/TMath::Sqrt(2.*nClusters));
382             
383             Double_t meanResY = hESDSumResidualYPerDE->GetBinContent(iDE+1)/nClusters;
384             hESDResidualYPerDEMean->SetBinContent(iDE+1, meanResY);
385             hESDResidualYPerDEMean->SetBinError(iDE+1, sigmaResidualY/TMath::Sqrt(nClusters));
386             hESDResidualYPerDESigma->SetBinContent(iDE+1, TMath::Sqrt(hESDSumResidualY2PerDE->GetBinContent(iDE+1)/nClusters - meanResY*meanResY));
387             hESDResidualYPerDESigma->SetBinError(iDE+1, sigmaResidualY/TMath::Sqrt(2.*nClusters));
388             
389             hESDLocalChi2XPerDEMean->SetBinContent(iDE+1, hESDSumLocalChi2XPerDE->GetBinContent(iDE+1)/nClusters);
390             hESDLocalChi2XPerDEMean->SetBinError(iDE+1, sigmaLocalChi2X/TMath::Sqrt(nClusters));
391             
392             hESDLocalChi2YPerDEMean->SetBinContent(iDE+1, hESDSumLocalChi2YPerDE->GetBinContent(iDE+1)/nClusters);
393             hESDLocalChi2YPerDEMean->SetBinError(iDE+1, sigmaLocalChi2Y/TMath::Sqrt(nClusters));
394             
395             hESDLocalChi2PerDEMean->SetBinContent(iDE+1, hESDSumLocalChi2PerDE->GetBinContent(iDE+1)/nClusters);
396             hESDLocalChi2PerDEMean->SetBinError(iDE+1, sigmaLocalChi2/TMath::Sqrt(nClusters));
397             
398           } else {
399             
400             hESDClusterChargePerDE->SetBinContent(iDE+1, hESDSumClusterChargePerDE->GetBinContent(iDE+1));
401             hESDClusterChargePerDE->SetBinError(iDE+1, hESDClusterChargeInCh->GetXaxis()->GetXmax());
402             
403             hESDResidualXPerDEMean->SetBinContent(iDE+1, hESDSumResidualXPerDE->GetBinContent(iDE+1));
404             hESDResidualXPerDEMean->SetBinError(iDE+1, hESDResidualXInCh->GetXaxis()->GetXmax());
405             hESDResidualXPerDESigma->SetBinContent(iDE+1, 0.);
406             hESDResidualXPerDESigma->SetBinError(iDE+1, hESDResidualXInCh->GetXaxis()->GetXmax());
407             
408             hESDResidualYPerDEMean->SetBinContent(iDE+1, hESDSumResidualYPerDE->GetBinContent(iDE+1));
409             hESDResidualYPerDEMean->SetBinError(iDE+1, hESDResidualYInCh->GetXaxis()->GetXmax());
410             hESDResidualYPerDESigma->SetBinContent(iDE+1, 0.);
411             hESDResidualYPerDESigma->SetBinError(iDE+1, hESDResidualYInCh->GetXaxis()->GetXmax());
412             
413             hESDLocalChi2XPerDEMean->SetBinContent(iDE+1, hESDSumLocalChi2XPerDE->GetBinContent(iDE+1));
414             hESDLocalChi2XPerDEMean->SetBinError(iDE+1, hESDLocalChi2XInCh->GetXaxis()->GetXmax());
415             
416             hESDLocalChi2YPerDEMean->SetBinContent(iDE+1, hESDSumLocalChi2YPerDE->GetBinContent(iDE+1));
417             hESDLocalChi2YPerDEMean->SetBinError(iDE+1, hESDLocalChi2YInCh->GetXaxis()->GetXmax());
418             
419             hESDLocalChi2PerDEMean->SetBinContent(iDE+1, hESDSumLocalChi2PerDE->GetBinContent(iDE+1));
420             hESDLocalChi2PerDEMean->SetBinError(iDE+1, hESDLocalChi2InCh->GetXaxis()->GetXmax());
421             
422           }
423           
424           Double_t nFullClusters = hESDnTotFullClustersPerDE->GetBinContent(iDE+1);
425           if (nFullClusters > 1) {
426             
427             hESDClusterSizePerDE->SetBinContent(iDE+1, hESDSumClusterSizePerDE->GetBinContent(iDE+1)/nFullClusters);
428             hESDClusterSizePerDE->SetBinError(iDE+1, sigmaSize/TMath::Sqrt(nFullClusters));
429             
430           } else {
431             
432             hESDClusterSizePerDE->SetBinContent(iDE+1, hESDSumClusterSizePerDE->GetBinContent(iDE+1));
433             hESDClusterSizePerDE->SetBinError(iDE+1, hESDClusterSizeInCh->GetXaxis()->GetXmax());
434             
435           }
436           
437           it.Next();
438         }
439         
440       }
441       
442     }
443     
444     // Display trigger histos in a more user friendly way
445     DisplayTriggerInfo(task);
446
447   } // loop on specie
448   
449   // do the QA checking
450   AliQAChecker::Instance()->Run(AliQAv1::kMUON, task, list) ;
451 }
452
453 //____________________________________________________________________________ 
454 void AliMUONQADataMakerRec::InitRaws()
455 {
456     /// create Raws histograms in Raws subdir
457         
458   AliCodeTimerAuto("");
459   
460   const Bool_t expert   = kTRUE ; 
461   const Bool_t saveCorr = kTRUE ; 
462   const Bool_t image    = kTRUE ; 
463  
464   TString boardName = "Local board Id";
465   
466   TString histoName, histoTitle;
467   for(Int_t iCath=0; iCath<AliMpConstants::NofCathodes(); iCath++){
468     TString cathName = ( iCath==0 ) ? "BendPlane" : "NonBendPlane";
469     for(Int_t iChamber=0; iChamber<AliMpConstants::NofTriggerChambers(); iChamber++){
470       histoName = Form("hTriggerScalers%sChamber%i", cathName.Data(), 11+iChamber);
471       histoTitle = Form("Chamber %i - %s: trigger scaler counts", 11+iChamber, cathName.Data());
472       TH2F* h3 = new TH2F(histoName.Data(), histoTitle.Data(),
473                           234, 0.5, 234.5,
474                           16, -0.5, 15.5);
475       h3->GetXaxis()->SetTitle(boardName.Data());
476       h3->GetYaxis()->SetTitle("Strip");        
477       Add2RawsList(h3, kTriggerScalers + AliMpConstants::NofTriggerChambers()*iCath + iChamber, expert, !image, !saveCorr);
478     }
479   }
480         
481         AliMUONTriggerDisplay triggerDisplay;
482         for(Int_t iCath=0; iCath<AliMpConstants::NofCathodes(); iCath++){
483                 TString cathName = ( iCath==0 ) ? "BendPlane" : "NonBendPlane";
484                 for(Int_t iChamber=0; iChamber<AliMpConstants::NofTriggerChambers(); iChamber++){
485                         histoName = Form("hScalersDisplay%sChamber%i", cathName.Data(), 11+iChamber);
486                         histoTitle = Form("Chamber %i - %s: Hit rate from scalers (Hz/cm^{2})", 11+iChamber, cathName.Data());
487                         TH2F* h5 = (TH2F*)triggerDisplay.GetEmptyDisplayHisto(histoName, AliMUONTriggerDisplay::kDisplayStrips, 
488                                                                               iCath, iChamber, histoTitle);
489                         Add2RawsList(h5, kTriggerScalersDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber, !expert, image, !saveCorr);
490                 }
491         }
492
493         FillTriggerDCSHistos();
494
495         TH1F* h10 = new TH1F("hTriggerScalersTime", "Acquisition time from trigger scalers", 1, 0.5, 1.5);
496         h10->GetXaxis()->SetBinLabel(1, "One-bin histogram: bin is filled at each scaler event.");
497         h10->GetYaxis()->SetTitle("Cumulated scaler time (s)");
498         Add2RawsList(h10, kTriggerScalersTime, !expert, !image, !saveCorr);
499         
500   Int_t nbp(0);
501   TIter next(AliMpDDLStore::Instance()->CreateBusPatchIterator());
502   while (next())
503   {
504     ++nbp;
505   }
506   
507   TH1* hbp = new TH1F("hTrackerBusPatchOccupancy","Occupancy of bus patches",
508                       nbp,-0.5,nbp-0.5);
509   
510   Add2RawsList(hbp,kTrackerBusPatchOccupancy, !expert, image, !saveCorr);
511
512   const Bool_t histogram(kFALSE);
513
514   fTrackerDataMaker = new AliMUONTrackerDataMaker(GetMUONRecoParam(),
515                                                   AliCDBManager::Instance()->GetRun(),
516                                                   0x0,
517                                                   "",
518                                                   "NOGAIN",
519                                                   histogram,
520                                                   0.0,0.0);
521                 
522   fTrackerDataMaker->Data()->DisableChannelLevel(); // to save up disk space, we only store starting at the manu level
523         
524   fTrackerDataMaker->SetRunning(kTRUE);
525 }
526
527 //__________________________________________________________________
528 void AliMUONQADataMakerRec::InitDigits() 
529 {
530   /// Initialized Digits spectra 
531   const Bool_t expert   = kTRUE ; 
532   const Bool_t image    = kTRUE ; 
533   
534   TH1I* h0 = new TH1I("hDigitsDetElem", "Detection element distribution in Digits;Detection element Id;Counts",  1400, 100, 1500); 
535   Add2DigitsList(h0, 0, !expert, image);
536   
537   TH1I* h1 = new TH1I("hDigitsADC", "ADC distribution in Digits;ACD value;Counts", 4096, 0, 4095); 
538   Add2DigitsList(h1, 1, !expert, image);    
539
540
541 //____________________________________________________________________________ 
542 void AliMUONQADataMakerRec::InitRecPoints()
543 {
544         /// create Reconstructed Points histograms in RecPoints subdir
545   
546   AliCodeTimerAuto("");
547   
548         InitRecPointsTrigger();
549         InitRecPointsTracker();
550 }
551
552 //____________________________________________________________________________ 
553 void AliMUONQADataMakerRec::InitRecPointsTracker()
554 {
555   /// create Reconstructed Points histograms in RecPoints subdir for the
556   /// MUON tracker subsystem.
557   const Bool_t expert   = kTRUE ; 
558   const Bool_t image    = kTRUE ; 
559   
560   AliCodeTimerAuto("");
561   
562   TH1I *h1I;
563   TH1F *h1F;
564   TH2F *h2F;
565   
566   // histograms per chamber
567   Int_t nCh = AliMpConstants::NofTrackingChambers();
568   for ( Int_t i = 0; i < nCh; ++i ) 
569   {
570     h1I = new TH1I(Form("hTrackerClusterMultiplicityForChamber%d",i+1), Form("cluster size distribution in chamber %d;size (n_{pads};Counts)",i+1), 100,0,100);
571     Add2RecPointsList(h1I,kTrackerClusterMultiplicityPerChamber+i, expert, !image);
572     
573     h1I = new TH1I(Form("hTrackerClusterChargeForChamber%d",i+1), Form("cluster charge distribution in chamber %d;charge (fC);Counts",i+1), 100,0,1000);
574     Add2RecPointsList(h1I,kTrackerClusterChargePerChamber+i, expert, !image);
575     
576     Float_t rMax = AliMUONConstants::Rmax(i/2);
577     h2F = new TH2F(Form("hTrackerClusterHitMapForChamber%d",i+1), Form("cluster position distribution in chamber %d;X (cm);Y (cm)",i+1), 100, -rMax, rMax, 100, -rMax, rMax);
578     Add2RecPointsList(h2F, kTrackerClusterHitMapPerChamber+i, expert, !image);
579   }
580   
581   // summary histograms per chamber
582   h1I = new TH1I("hTrackerNumberOfClustersPerChamber", "Number of clusters per chamber;chamber ID;n_{clusters}", nCh,-0.5,nCh-0.5);
583   Add2RecPointsList(h1I,kTrackerNumberOfClustersPerChamber, !expert, image);
584   
585   h1F = new TH1F("hTrackerClusterMultiplicityPerChMean", "cluster mean size per chamber;chamber ID;<size> (n_{pads})", nCh,-0.5,nCh-0.5);
586   h1F->SetOption("P");
587   h1F->SetMarkerStyle(kFullDotMedium);
588   h1F->SetMarkerColor(kRed);
589   Add2RecPointsList(h1F, kTrackerClusterMultiplicityPerChMean, !expert, image);
590   
591   h1F = new TH1F("hTrackerClusterMultiplicityPerChSigma", "cluster size dispersion per chamber;chamber ID;#sigma_{size} (n_{pads})", nCh,-0.5,nCh-0.5);
592   h1F->SetOption("P");
593   h1F->SetMarkerStyle(kFullDotMedium);
594   h1F->SetMarkerColor(kRed);
595   Add2RecPointsList(h1F, kTrackerClusterMultiplicityPerChSigma, !expert, image);
596   
597   h1F = new TH1F("hTrackerClusterChargePerChMean", "cluster mean charge per chamber;chamber ID;<charge> (fC)", nCh,-0.5,nCh-0.5);
598   h1F->SetOption("P");
599   h1F->SetMarkerStyle(kFullDotMedium);
600   h1F->SetMarkerColor(kRed);
601   Add2RecPointsList(h1F, kTrackerClusterChargePerChMean, !expert, image);
602   
603   h1F = new TH1F("hTrackerClusterChargePerChSigma", "cluster charge dispersion per chamber;chamber ID;#sigma_{charge} (fC)", nCh,-0.5,nCh-0.5);
604   h1F->SetOption("P");
605   h1F->SetMarkerStyle(kFullDotMedium);
606   h1F->SetMarkerColor(kRed);
607   Add2RecPointsList(h1F, kTrackerClusterChargePerChSigma, !expert, image);
608   
609   // histograms per DE
610   Int_t ndes(0);
611   AliMpDEIterator it;
612   it.First();
613   while ( !it.IsDone())
614   {
615     Int_t detElemId = it.CurrentDEId();
616     
617     if ( AliMpDEManager::GetStationType(detElemId) != AliMp::kStationTrigger )
618     {
619       ndes = TMath::Max(ndes,detElemId);
620       
621       h1I = new TH1I(Form("hTrackerClusterMultiplicityForDE%04d",detElemId), Form("cluster size distribution in detection element %d;size (n_{pads})",detElemId), 100,0,100);
622       Add2RecPointsList(h1I,kTrackerClusterMultiplicityPerDE+detElemId, expert, !image);
623       
624       h1I = new TH1I(Form("hTrackerClusterChargeForDE%04d",detElemId), Form("cluster charge distribution in detection element %d;charge (fC)",detElemId), 100,0,1000);
625       Add2RecPointsList(h1I,kTrackerClusterChargePerDE+detElemId, expert, !image);
626     }
627     
628     it.Next();
629   }
630   
631   // summary histograms per DE
632   h1I = new TH1I("hTrackerNumberOfClustersPerDE", "Number of clusters per detection element;DetElem ID;n_{clusters}", ndes+1,-0.5,ndes+0.5);
633   Add2RecPointsList(h1I, kTrackerNumberOfClustersPerDE, !expert, image);
634   
635   h1F = new TH1F("hTrackerClusterMultiplicityPerDEMean", "cluster mean size per DE;DetElem ID;<size> (n_{pads})", ndes+1,-0.5,ndes+0.5);
636   h1F->SetOption("P");
637   h1F->SetMarkerStyle(kFullDotMedium);
638   h1F->SetMarkerColor(kRed);
639   Add2RecPointsList(h1F, kTrackerClusterMultiplicityPerDEMean, !expert, image);
640   
641   h1F = new TH1F("hTrackerClusterChargePerDEMean", "cluster mean charge per DE;DetElem ID;<charge> (fC)", ndes+1,-0.5,ndes+0.5);
642   h1F->SetOption("P");
643   h1F->SetMarkerStyle(kFullDotMedium);
644   h1F->SetMarkerColor(kRed);
645   Add2RecPointsList(h1F, kTrackerClusterChargePerDEMean, !expert, image);
646 }
647
648 //____________________________________________________________________________ 
649 void AliMUONQADataMakerRec::InitRecPointsTrigger()
650 {
651         /// create Reconstructed Points histograms in RecPoints subdir for the
652         /// MUON Trigger subsystem.
653         
654   const Bool_t expert   = kTRUE ; 
655   const Bool_t image    = kTRUE ; 
656
657   TString boardName = "Local board Id";
658
659   TString histoName, histoTitle;
660   for(Int_t iCath=0; iCath<AliMpConstants::NofCathodes(); iCath++){
661     TString cathName = ( iCath==0 ) ? "BendPlane" : "NonBendPlane";
662     for(Int_t iChamber=0; iChamber<AliMpConstants::NofTriggerChambers(); iChamber++){
663       histoName = Form("hTriggerDigits%sChamber%i", cathName.Data(), 11+iChamber);
664       histoTitle = Form("Chamber %i - %s: counts per strip", 11+iChamber, cathName.Data());
665       TH2F* h0 = new TH2F(histoName.Data(), histoTitle.Data(),
666                           234, 0.5, 234.5,
667                           16, -0.5, 15.5);
668       h0->GetXaxis()->SetTitle(boardName.Data());
669       h0->GetYaxis()->SetTitle("Strip");
670       Add2RecPointsList(h0, kTriggerDigits + AliMpConstants::NofTriggerChambers()*iCath + iChamber, expert, !image);
671     }
672   }
673
674     TH1F* h2 = new TH1F("hTriggeredBoards", "Triggered boards", 234, 0.5, 234.5);
675     Add2RecPointsList(h2, kTriggeredBoards, expert, !image);
676
677     AliMUONTriggerDisplay triggerDisplay;
678     for(Int_t iCath=0; iCath<AliMpConstants::NofCathodes(); iCath++){
679       TString cathName = ( iCath==0 ) ? "BendPlane" : "NonBendPlane";
680       for(Int_t iChamber=0; iChamber<AliMpConstants::NofTriggerChambers(); iChamber++){
681         histoName = Form("hTriggerDigitsDisplay%sChamber%i", cathName.Data(), 11+iChamber);
682         histoTitle = Form("Chamber %i - %s: Occupancy per strip (counts/event)", 11+iChamber, cathName.Data());
683         TH2F* h3 = (TH2F*)triggerDisplay.GetEmptyDisplayHisto(histoName, AliMUONTriggerDisplay::kDisplayStrips, 
684                                                               iCath, iChamber, histoTitle);
685         Add2RecPointsList(h3, kTriggerDigitsDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber, !expert, image);
686       }
687     }
688
689     TH2F* h4 = (TH2F*)triggerDisplay.GetEmptyDisplayHisto("hFiredBoardsDisplay", AliMUONTriggerDisplay::kDisplayBoards,
690                                                           0, 0, "Local board triggers / event");
691     Add2RecPointsList(h4, kTriggerBoardsDisplay, !expert, image);
692
693     TH1F* h5 = new TH1F("hNAnalyzedEvents", "Number of analyzed events per specie", 1, 0.5, 1.5);
694     Int_t esindex = AliRecoParam::AConvert(fEventSpecie);
695     h5->GetXaxis()->SetBinLabel(1, AliRecoParam::GetEventSpecieName(esindex));
696     h5->GetYaxis()->SetTitle("Number of analyzed events");
697     Add2RecPointsList(h5, kNAnalyzedEvents, !expert, image);
698 }
699
700
701 //____________________________________________________________________________ 
702 void AliMUONQADataMakerRec::InitESDs()
703 {
704   ///create ESDs histograms in ESDs subdir
705   
706   const Bool_t expert   = kTRUE ; 
707   const Bool_t image    = kTRUE ; 
708   
709   Int_t nCh = AliMUONConstants::NTrackingCh();
710   Int_t nDE = 1100;
711   
712   // track info
713   TH1F* hESDnTracks = new TH1F("hESDnTracks", "number of tracks;n_{tracks}", 20, 0., 20.);
714   Add2ESDsList(hESDnTracks, kESDnTracks, !expert, image);
715
716   TH1F* hESDMatchTrig = new TH1F("hESDMatchTrig", "number of tracks matched with trigger;n_{tracks}", 20, 0., 20.);
717   Add2ESDsList(hESDMatchTrig, kESDMatchTrig, !expert, image);
718   
719   TH1F* hESDMomentum = new TH1F("hESDMomentum", "momentum distribution;p (GeV/c)", 300, 0., 300);
720   Add2ESDsList(hESDMomentum, kESDMomentum, !expert, image);
721
722   TH1F* hESDPt = new TH1F("hESDPt", "transverse momentum distribution;p_{t} (GeV/c)", 200, 0., 50);
723   Add2ESDsList(hESDPt, kESDPt, !expert, image);
724
725   TH1F* hESDRapidity = new TH1F("hESDRapidity", "rapidity distribution;rapidity", 200, -4.5, -2.);
726   Add2ESDsList(hESDRapidity, kESDRapidity, !expert, image);
727
728   TH1F* hESDChi2 = new TH1F("hESDChi2", "normalized #chi^{2} distribution;#chi^{2} / ndf", 500, 0., 50.);
729   Add2ESDsList(hESDChi2, kESDChi2, !expert, image);
730   
731   TH1F* hESDProbChi2 = new TH1F("hESDProbChi2", "distribution of probability of #chi^{2};prob(#chi^{2})", 100, 0., 1.);
732   Add2ESDsList(hESDProbChi2, kESDProbChi2, !expert, image);
733   
734   TH1F* hESDThetaX = new TH1F("hESDThetaX", "#theta_{X} distribution;#theta_{X} (degree)", 360, -180., 180);
735   Add2ESDsList(hESDThetaX, kESDThetaX, !expert, image);
736   
737   TH1F* hESDThetaY = new TH1F("hESDThetaY", "#theta_{Y} distribution;#theta_{Y} (degree)", 360, -180., 180);
738   Add2ESDsList(hESDThetaY, kESDThetaY, !expert, image);
739   
740   // cluster info
741   for (Int_t i = 0; i < nCh; i++) {
742     Float_t rMax = AliMUONConstants::Rmax(i/2);
743     TH2F* hESDClusterHitMap = new TH2F(Form("hESDClusterHitMap%d",i+1), Form("cluster position distribution in chamber %d;X (cm);Y (cm)",i+1),
744                                        100, -rMax, rMax, 100, -rMax, rMax);
745     Add2ESDsList(hESDClusterHitMap, kESDClusterHitMap+i, expert, !image);
746   }
747   
748   TH1F* hESDnClustersPerTrack = new TH1F("hESDnClustersPerTrack", "number of associated clusters per track;n_{clusters}", 20, 0., 20.);
749   Add2ESDsList(hESDnClustersPerTrack, kESDnClustersPerTrack, !expert, image);
750   
751   TH1F* hESDnClustersPerCh = new TH1F("hESDnClustersPerCh", "averaged number of clusters per chamber per track;chamber ID;<n_{clusters}>", nCh, -0.5, nCh-0.5);
752   hESDnClustersPerCh->SetFillColor(kRed);
753   Add2ESDsList(hESDnClustersPerCh, kESDnClustersPerCh, !expert, image);
754   
755   TH1F* hESDnClustersPerDE = new TH1F("hESDnClustersPerDE", "averaged number of clusters per DE per track;DetElem ID;<n_{clusters}>", nDE+1, -0.5, nDE+0.5);
756   hESDnClustersPerDE->SetFillColor(kRed);
757   Add2ESDsList(hESDnClustersPerDE, kESDnClustersPerDE, !expert, image);
758   
759   for (Int_t i = 0; i < nCh; i++) {
760     TH1F* hESDClusterChargeInCh = new TH1F(Form("hESDClusterChargeInCh%d",i+1), Form("cluster charge distribution in chamber %d;charge (fC)",i+1), 100, 0., 1000.);
761     Add2ESDsList(hESDClusterChargeInCh, kESDClusterChargeInCh+i, expert, !image);
762   }
763   
764   TH1F* hESDClusterChargePerChMean = new TH1F("hESDClusterChargePerChMean", "cluster mean charge per chamber;chamber ID;<charge> (fC)", nCh, -0.5, nCh-0.5);
765   hESDClusterChargePerChMean->SetOption("P");
766   hESDClusterChargePerChMean->SetMarkerStyle(kFullDotMedium);
767   hESDClusterChargePerChMean->SetMarkerColor(kRed);
768   Add2ESDsList(hESDClusterChargePerChMean, kESDClusterChargePerChMean, !expert, image);
769   
770   TH1F* hESDClusterChargePerChSigma = new TH1F("hESDClusterChargePerChSigma", "cluster charge dispersion per chamber;chamber ID;#sigma_{charge} (fC)", nCh, -0.5, nCh-0.5);
771   hESDClusterChargePerChSigma->SetOption("P");
772   hESDClusterChargePerChSigma->SetMarkerStyle(kFullDotMedium);
773   hESDClusterChargePerChSigma->SetMarkerColor(kRed);
774   Add2ESDsList(hESDClusterChargePerChSigma, kESDClusterChargePerChSigma, !expert, image);
775   
776   TH1F* hESDClusterChargePerDE = new TH1F("hESDClusterChargePerDE", "cluster mean charge per DE;DetElem ID;<charge> (fC)", nDE+1, -0.5, nDE+0.5);
777   hESDClusterChargePerDE->SetOption("P");
778   hESDClusterChargePerDE->SetMarkerStyle(kFullDotMedium);
779   hESDClusterChargePerDE->SetMarkerColor(kRed);
780   Add2ESDsList(hESDClusterChargePerDE, kESDClusterChargePerDE, !expert, image);
781   
782   for (Int_t i = 0; i < nCh; i++) {
783     TH1F* hESDClusterSizeInCh = new TH1F(Form("hESDClusterSizeInCh%d",i+1), Form("cluster size distribution in chamber %d;size (n_{pads})",i+1), 200, 0., 200.);
784     Add2ESDsList(hESDClusterSizeInCh, kESDClusterSizeInCh+i, expert, !image);
785   }
786   
787   TH1F* hESDClusterSizePerChMean = new TH1F("hESDClusterSizePerChMean", "cluster mean size per chamber;chamber ID;<size> (n_{pads})", nCh, -0.5, nCh-0.5);
788   hESDClusterSizePerChMean->SetOption("P");
789   hESDClusterSizePerChMean->SetMarkerStyle(kFullDotMedium);
790   hESDClusterSizePerChMean->SetMarkerColor(kRed);
791   Add2ESDsList(hESDClusterSizePerChMean, kESDClusterSizePerChMean, !expert, image);
792   
793   TH1F* hESDClusterSizePerChSigma = new TH1F("hESDClusterSizePerChSigma", "cluster size dispersion per chamber;chamber ID;#sigma_{size} (n_{pads})", nCh, -0.5, nCh-0.5);
794   hESDClusterSizePerChSigma->SetOption("P");
795   hESDClusterSizePerChSigma->SetMarkerStyle(kFullDotMedium);
796   hESDClusterSizePerChSigma->SetMarkerColor(kRed);
797   Add2ESDsList(hESDClusterSizePerChSigma, kESDClusterSizePerChSigma, !expert, image);
798   
799   TH1F* hESDClusterSizePerDE = new TH1F("hESDClusterSizePerDE", "cluster mean size per DE;DetElem ID;<size> (n_{pads})", nDE+1, -0.5, nDE+0.5);
800   hESDClusterSizePerDE->SetOption("P");
801   hESDClusterSizePerDE->SetMarkerStyle(kFullDotMedium);
802   hESDClusterSizePerDE->SetMarkerColor(kRed);
803   Add2ESDsList(hESDClusterSizePerDE, kESDClusterSizePerDE, !expert, image);
804   
805   // cluster - track info
806   for (Int_t i = 0; i < nCh; i++) {
807     TH1F* hESDResidualXInCh = new TH1F(Form("hESDResidualXInCh%d",i+1), Form("cluster-track residual-X distribution in chamber %d;#Delta_{X} (cm)",i+1), 1000, -5., 5.);
808     Add2ESDsList(hESDResidualXInCh, kESDResidualXInCh+i, expert, !image);
809     
810     TH1F* hESDResidualYInCh = new TH1F(Form("hESDResidualYInCh%d",i+1), Form("cluster-track residual-Y distribution in chamber %d;#Delta_{Y} (cm)",i+1), 1000, -1., 1.);
811     Add2ESDsList(hESDResidualYInCh, kESDResidualYInCh+i, expert, !image);
812     
813     TH1F* hESDLocalChi2XInCh = new TH1F(Form("hESDLocalChi2XInCh%d",i+1), Form("local chi2-X distribution in chamber %d;local #chi^{2}_{X}",i+1), 1000, 0., 25);
814     Add2ESDsList(hESDLocalChi2XInCh, kESDLocalChi2XInCh+i, expert, !image);
815     
816     TH1F* hESDLocalChi2YInCh = new TH1F(Form("hESDLocalChi2YInCh%d",i+1), Form("local chi2-Y distribution in chamber %d;local #chi^{2}_{Y}",i+1), 1000, 0., 25);
817     Add2ESDsList(hESDLocalChi2YInCh, kESDLocalChi2YInCh+i, expert, !image);
818     
819     TH1F* hESDLocalChi2InCh = new TH1F(Form("hESDLocalChi2InCh%d",i+1), Form("local chi2 (~0.5*(#chi^{2}_{X}+#chi^{2}_{Y})) distribution in chamber %d;local #chi^{2}",i+1), 1000, 0., 25);
820     Add2ESDsList(hESDLocalChi2InCh, kESDLocalChi2InCh+i, expert, !image);
821   }
822   
823   TH1F* hESDResidualXPerChMean = new TH1F("hESDResidualXPerChMean", "cluster-track residual-X per Ch: mean;chamber ID;<#Delta_{X}> (cm)", nCh, -0.5, nCh-0.5);
824   hESDResidualXPerChMean->SetOption("P");
825   hESDResidualXPerChMean->SetMarkerStyle(kFullDotMedium);
826   hESDResidualXPerChMean->SetMarkerColor(kRed);
827   Add2ESDsList(hESDResidualXPerChMean, kESDResidualXPerChMean, !expert, image);
828   
829   TH1F* hESDResidualYPerChMean = new TH1F("hESDResidualYPerChMean", "cluster-track residual-Y per Ch: mean;chamber ID;<#Delta_{Y}> (cm)", nCh, -0.5, nCh-0.5);
830   hESDResidualYPerChMean->SetOption("P");
831   hESDResidualYPerChMean->SetMarkerStyle(kFullDotMedium);
832   hESDResidualYPerChMean->SetMarkerColor(kRed);
833   Add2ESDsList(hESDResidualYPerChMean, kESDResidualYPerChMean, !expert, image);
834   
835   TH1F* hESDResidualXPerChSigma = new TH1F("hESDResidualXPerChSigma", "cluster-track residual-X per Ch: sigma;chamber ID;#sigma_{X} (cm)", nCh, -0.5, nCh-0.5);
836   hESDResidualXPerChSigma->SetOption("P");
837   hESDResidualXPerChSigma->SetMarkerStyle(kFullDotMedium);
838   hESDResidualXPerChSigma->SetMarkerColor(kRed);
839   Add2ESDsList(hESDResidualXPerChSigma, kESDResidualXPerChSigma, !expert, image);
840   
841   TH1F* hESDResidualYPerChSigma = new TH1F("hESDResidualYPerChSigma", "cluster-track residual-Y per Ch: sigma;chamber ID;#sigma_{Y} (cm)", nCh, -0.5, nCh-0.5);
842   hESDResidualYPerChSigma->SetOption("P");
843   hESDResidualYPerChSigma->SetMarkerStyle(kFullDotMedium);
844   hESDResidualYPerChSigma->SetMarkerColor(kRed);
845   Add2ESDsList(hESDResidualYPerChSigma, kESDResidualYPerChSigma, !expert, image);
846   
847   TH1F* hESDLocalChi2XPerChMean = new TH1F("hESDLocalChi2XPerCh", "local chi2-X per Ch: mean;chamber ID;<local #chi^{2}_{X}>", nCh, -0.5, nCh-0.5);
848   hESDLocalChi2XPerChMean->SetOption("P");
849   hESDLocalChi2XPerChMean->SetMarkerStyle(kFullDotMedium);
850   hESDLocalChi2XPerChMean->SetMarkerColor(kRed);
851   Add2ESDsList(hESDLocalChi2XPerChMean, kESDLocalChi2XPerChMean, !expert, image);
852   
853   TH1F* hESDLocalChi2YPerChMean = new TH1F("hESDLocalChi2YPerCh", "local chi2-Y per Ch: mean;chamber ID;<local #chi^{2}_{Y}>", nCh, -0.5, nCh-0.5);
854   hESDLocalChi2YPerChMean->SetOption("P");
855   hESDLocalChi2YPerChMean->SetMarkerStyle(kFullDotMedium);
856   hESDLocalChi2YPerChMean->SetMarkerColor(kRed);
857   Add2ESDsList(hESDLocalChi2YPerChMean, kESDLocalChi2YPerChMean, !expert, image);
858   
859   TH1F* hESDLocalChi2PerChMean = new TH1F("hESDLocalChi2PerCh", "local chi2 (~0.5*(#chi^{2}_{X}+#chi^{2}_{Y})) per Ch: mean;chamber ID;<local #chi^{2}>", nCh, -0.5, nCh-0.5);
860   hESDLocalChi2PerChMean->SetOption("P");
861   hESDLocalChi2PerChMean->SetMarkerStyle(kFullDotMedium);
862   hESDLocalChi2PerChMean->SetMarkerColor(kRed);
863   Add2ESDsList(hESDLocalChi2PerChMean, kESDLocalChi2PerChMean, !expert, image);
864   
865   TH1F* hESDResidualXPerDEMean = new TH1F("hESDResidualXPerDEMean", "cluster-track residual-X per DE: mean;DetElem ID;<#Delta_{X}> (cm)", nDE+1, -0.5, nDE+0.5);
866   hESDResidualXPerDEMean->SetOption("P");
867   hESDResidualXPerDEMean->SetMarkerStyle(kFullDotMedium);
868   hESDResidualXPerDEMean->SetMarkerColor(kRed);
869   Add2ESDsList(hESDResidualXPerDEMean, kESDResidualXPerDEMean, !expert, image);
870   
871   TH1F* hESDResidualYPerDEMean = new TH1F("hESDResidualYPerDEMean", "cluster-track residual-Y per DE: mean;DetElem ID;<#Delta_{Y}> (cm)", nDE+1, -0.5, nDE+0.5);
872   hESDResidualYPerDEMean->SetOption("P");
873   hESDResidualYPerDEMean->SetMarkerStyle(kFullDotMedium);
874   hESDResidualYPerDEMean->SetMarkerColor(kRed);
875   Add2ESDsList(hESDResidualYPerDEMean, kESDResidualYPerDEMean, !expert, image);
876   
877   TH1F* hESDResidualXPerDESigma = new TH1F("hESDResidualXPerDESigma", "cluster-track residual-X per DE: sigma;DetElem ID;#sigma_{X} (cm)", nDE+1, -0.5, nDE+0.5);
878   hESDResidualXPerDESigma->SetOption("P");
879   hESDResidualXPerDESigma->SetMarkerStyle(kFullDotMedium);
880   hESDResidualXPerDESigma->SetMarkerColor(kRed);
881   Add2ESDsList(hESDResidualXPerDESigma, kESDResidualXPerDESigma, !expert, image);
882   
883   TH1F* hESDResidualYPerDESigma = new TH1F("hESDResidualYPerDESigma", "cluster-track residual-Y per DE: sigma;DetElem ID;#sigma_{Y} (cm)", nDE+1, -0.5, nDE+0.5);
884   hESDResidualYPerDESigma->SetOption("P");
885   hESDResidualYPerDESigma->SetMarkerStyle(kFullDotMedium);
886   hESDResidualYPerDESigma->SetMarkerColor(kRed);
887   Add2ESDsList(hESDResidualYPerDESigma, kESDResidualYPerDESigma, !expert, image);
888   
889   TH1F* hESDLocalChi2XPerDEMean = new TH1F("hESDLocalChi2XPerDE", "local chi2-X per DE: mean;DetElem ID;<local #chi^{2}_{X}>", nDE+1, -0.5, nDE+0.5);
890   hESDLocalChi2XPerDEMean->SetOption("P");
891   hESDLocalChi2XPerDEMean->SetMarkerStyle(kFullDotMedium);
892   hESDLocalChi2XPerDEMean->SetMarkerColor(kRed);
893   Add2ESDsList(hESDLocalChi2XPerDEMean, kESDLocalChi2XPerDEMean, !expert, image);
894   
895   TH1F* hESDLocalChi2YPerDEMean = new TH1F("hESDLocalChi2YPerDE", "local chi2-Y per DE: mean;DetElem ID;<local #chi^{2}_{Y}>", nDE+1, -0.5, nDE+0.5);
896   hESDLocalChi2YPerDEMean->SetOption("P");
897   hESDLocalChi2YPerDEMean->SetMarkerStyle(kFullDotMedium);
898   hESDLocalChi2YPerDEMean->SetMarkerColor(kRed);
899   Add2ESDsList(hESDLocalChi2YPerDEMean, kESDLocalChi2YPerDEMean, !expert, image);
900   
901   TH1F* hESDLocalChi2PerDEMean = new TH1F("hESDLocalChi2PerDE", "local chi2 (~0.5*(#chi^{2}_{X}+#chi^{2}_{Y})) per DE: mean;DetElem ID;<local #chi^{2}>", nDE+1, -0.5, nDE+0.5);
902   hESDLocalChi2PerDEMean->SetOption("P");
903   hESDLocalChi2PerDEMean->SetMarkerStyle(kFullDotMedium);
904   hESDLocalChi2PerDEMean->SetMarkerColor(kRed);
905   Add2ESDsList(hESDLocalChi2PerDEMean, kESDLocalChi2PerDEMean, !expert, image);
906   
907   // intermediate histograms
908   TH1F* hESDnTotClustersPerCh = new TH1F("hESDnTotClustersPerCh", "total number of associated clusters per chamber;chamber ID;#Sigma(n_{clusters})", nCh, -0.5, nCh-0.5);
909   Add2ESDsList(hESDnTotClustersPerCh, kESDnTotClustersPerCh, expert, !image);
910   TH1F* hESDnTotClustersPerDE = new TH1F("hESDnTotClustersPerDE", "total number of associated clusters per DE;DetElem ID;#Sigma(n_{clusters})", nDE+1, -0.5, nDE+0.5);
911   Add2ESDsList(hESDnTotClustersPerDE, kESDnTotClustersPerDE, expert, !image);
912   TH1F* hESDnTotFullClustersPerDE = new TH1F("hESDnTotFullClustersPerDE", "total number of associated clusters containing pad info per DE;DetElem ID;#Sigma(n_{full clusters})", nDE+1, -0.5, nDE+0.5);
913   Add2ESDsList(hESDnTotFullClustersPerDE, kESDnTotFullClustersPerDE, expert, !image);
914   TH1F* hESDSumClusterChargePerDE = new TH1F("hESDSumClusterChargePerDE", "sum of cluster charge per DE;DetElem ID;#Sigma(charge) (fC)", nDE+1, -0.5, nDE+0.5);
915   Add2ESDsList(hESDSumClusterChargePerDE, kESDSumClusterChargePerDE, expert, !image);
916   TH1F* hESDSumClusterSizePerDE = new TH1F("hESDSumClusterSizePerDE", "sum of cluster size per DE;DetElem ID;#Sigma(size) (n_{pads})", nDE+1, -0.5, nDE+0.5);
917   Add2ESDsList(hESDSumClusterSizePerDE, kESDSumClusterSizePerDE, expert, !image);
918   TH1F* hESDSumResidualXPerDE = new TH1F("hESDSumResidualXPerDE", "sum of cluster-track residual-X per DE;DetElem ID;#Sigma(#Delta_{X}) (cm)", nDE+1, -0.5, nDE+0.5);
919   Add2ESDsList(hESDSumResidualXPerDE, kESDSumResidualXPerDE, expert, !image);
920   TH1F* hESDSumResidualYPerDE = new TH1F("hESDSumResidualYPerDE", "sum of cluster-track residual-Y per DE;DetElem ID;#Sigma(#Delta_{Y}) (cm)", nDE+1, -0.5, nDE+0.5);
921   Add2ESDsList(hESDSumResidualYPerDE, kESDSumResidualYPerDE, expert, !image);
922   TH1F* hESDSumResidualX2PerDE = new TH1F("hESDSumResidualX2PerDE", "sum of cluster-track residual-X**2 per DE;DetElem ID;#Sigma(#Delta_{X}^{2}) (cm^{2})", nDE+1, -0.5, nDE+0.5);
923   Add2ESDsList(hESDSumResidualX2PerDE, kESDSumResidualX2PerDE, expert, !image);
924   TH1F* hESDSumResidualY2PerDE = new TH1F("hESDSumResidualY2PerDE", "sum of cluster-track residual-Y**2 per DE;DetElem ID;#Sigma(#Delta_{Y}^{2}) (cm^{2})", nDE+1, -0.5, nDE+0.5);
925   Add2ESDsList(hESDSumResidualY2PerDE, kESDSumResidualY2PerDE, expert, !image);
926   TH1F* hESDSumLocalChi2XPerDE = new TH1F("hESDSumLocalChi2XPerDE", "sum of local chi2-X per DE;DetElem ID;#Sigma(local #chi^{2}_{X})", nDE+1, -0.5, nDE+0.5);
927   Add2ESDsList(hESDSumLocalChi2XPerDE, kESDSumLocalChi2XPerDE, expert, !image);
928   TH1F* hESDSumLocalChi2YPerDE = new TH1F("hESDSumLocalChi2YPerDE", "sum of local chi2-Y per DE;DetElem ID;#Sigma(local #chi^{2}_{Y})", nDE+1, -0.5, nDE+0.5);
929   Add2ESDsList(hESDSumLocalChi2YPerDE, kESDSumLocalChi2YPerDE, expert, !image);
930   TH1F* hESDSumLocalChi2PerDE = new TH1F("hESDSumLocalChi2PerDE", "sum of local chi2 (~0.5*(#chi^{2}_{X}+#chi^{2}_{Y})) per DE;DetElem ID;#Sigma(local #chi^{2})", nDE+1, -0.5, nDE+0.5);
931   Add2ESDsList(hESDSumLocalChi2PerDE, kESDSumLocalChi2PerDE, expert, !image);
932 }
933
934 //____________________________________________________________________________
935 void AliMUONQADataMakerRec::MakeRaws(AliRawReader* rawReader)
936 {
937     /// make QA for rawdata
938
939     // Check id histograms already created for this Event Specie
940
941   AliDebug(AliQAv1::GetQADebugLevel(),
942            Form("RAW event type %s", AliRawEventHeaderBase::GetTypeName(rawReader->GetType())));
943   
944   if ( rawReader->GetType() == AliRawEventHeaderBase::kPhysicsEvent ) 
945   {
946     rawReader->Reset();
947     MakeRawsTracker(rawReader);
948   }
949   
950   if ( rawReader->GetType() == AliRawEventHeaderBase::kPhysicsEvent ||
951        rawReader->GetType() == AliRawEventHeaderBase::kCalibrationEvent ) 
952   {
953     rawReader->Reset();    
954     MakeRawsTrigger(rawReader);
955   }
956 }
957
958 //____________________________________________________________________________
959 void AliMUONQADataMakerRec::MakeRawsTracker(AliRawReader* rawReader)
960 {
961         /// make QA for rawdata tracker
962         
963   /// forces init
964   GetRawsData(kTrackerBusPatchOccupancy);
965   
966         ((AliMUONTrackerDataMaker*)fTrackerDataMaker)->SetRawReader(rawReader);
967         
968         fTrackerDataMaker->ProcessEvent();
969 }
970
971 //____________________________________________________________________________
972 void AliMUONQADataMakerRec::MakeRawsTrigger(AliRawReader* rawReader)
973 {
974         /// make QA for rawdata trigger
975         
976     // Get trigger scalers
977
978     Int_t loCircuit=0;
979     AliMpCDB::LoadDDLStore();
980
981     AliMUONRawStreamTrigger rawStreamTrig(rawReader);
982     while (rawStreamTrig.NextDDL()) 
983     {
984       // If not a scaler event, do nothing
985       Bool_t scalerEvent =  rawReader->GetDataHeader()->GetL1TriggerMessage() & 0x1;
986       if ( !scalerEvent ) continue;
987
988       AliDebug(AliQAv1::GetQADebugLevel(),"Filling trigger scalers");
989
990       AliMUONDDLTrigger* ddlTrigger = rawStreamTrig.GetDDLTrigger();
991       AliMUONDarcHeader* darcHeader = ddlTrigger->GetDarcHeader();
992
993       if (darcHeader->GetGlobalFlag()){
994         UInt_t nOfClocks = darcHeader->GetGlobalClock();
995         Double_t nOfSeconds = ((Double_t) nOfClocks) / 40e6; // 1 clock each 25 ns
996         ((TH1F*)GetRawsData(kTriggerScalersTime))->Fill(1., nOfSeconds);
997       }
998
999       Int_t nReg = darcHeader->GetRegHeaderEntries();
1000     
1001       for(Int_t iReg = 0; iReg < nReg ;iReg++)
1002       {   //reg loop
1003
1004         // crate info  
1005         AliMpTriggerCrate* crate = AliMpDDLStore::Instance()->
1006           GetTriggerCrate(rawStreamTrig.GetDDL(), iReg);
1007
1008         AliMUONRegHeader* regHeader =  darcHeader->GetRegHeaderEntry(iReg);
1009
1010         // loop over local structures
1011         Int_t nLocal = regHeader->GetLocalEntries();
1012         for(Int_t iLocal = 0; iLocal < nLocal; iLocal++) 
1013         {
1014           AliMUONLocalStruct* localStruct = regHeader->GetLocalEntry(iLocal);
1015         
1016           // if card exist
1017           if (!localStruct) continue;
1018           
1019           loCircuit = crate->GetLocalBoardId(localStruct->GetId());
1020
1021           if ( !loCircuit ) continue; // empty slot
1022
1023           AliMpLocalBoard* localBoard = AliMpDDLStore::Instance()->GetLocalBoard(loCircuit, false);
1024
1025           // skip copy cards
1026           if( !localBoard->IsNotified()) 
1027             continue;
1028
1029           Int_t cathode = localStruct->GetComptXY()%2;
1030
1031           // loop over strips
1032           for (Int_t ibitxy = 0; ibitxy < 16; ++ibitxy) {
1033             if(localStruct->GetXY1(ibitxy) > 0)
1034               ((TH2F*)GetRawsData(kTriggerScalers + AliMpConstants::NofTriggerChambers()*cathode + 0))
1035                 ->Fill(loCircuit, ibitxy, 2*localStruct->GetXY1(ibitxy));
1036             if(localStruct->GetXY2(ibitxy) > 0)
1037               ((TH2F*)GetRawsData(kTriggerScalers + AliMpConstants::NofTriggerChambers()*cathode + 1))
1038                 ->Fill(loCircuit, ibitxy, 2*localStruct->GetXY2(ibitxy));
1039             if(localStruct->GetXY3(ibitxy) > 0)
1040               ((TH2F*)GetRawsData(kTriggerScalers + AliMpConstants::NofTriggerChambers()*cathode + 2))
1041                 ->Fill(loCircuit, ibitxy, 2*localStruct->GetXY3(ibitxy));
1042             if(localStruct->GetXY4(ibitxy) > 0)
1043               ((TH2F*)GetRawsData(kTriggerScalers + AliMpConstants::NofTriggerChambers()*cathode + 3))
1044                 ->Fill(loCircuit, ibitxy, 2*localStruct->GetXY4(ibitxy));
1045           } // loop on strips
1046         } // iLocal
1047       } // iReg
1048     } // NextDDL
1049 }
1050
1051 //__________________________________________________________________
1052 void AliMUONQADataMakerRec::MakeDigits(TTree* digitsTree)         
1053 {
1054   /// makes data from Digits
1055
1056   // Do nothing in case of calibration event
1057   if ( GetRecoParam()->GetEventSpecie() == AliRecoParam::kCalib ) return;
1058
1059   if (!fDigitStore)
1060     fDigitStore = AliMUONVDigitStore::Create(*digitsTree);
1061   fDigitStore->Connect(*digitsTree, false);
1062   digitsTree->GetEvent(0);
1063   
1064   TIter next(fDigitStore->CreateIterator());
1065   
1066   AliMUONVDigit* dig = 0x0;
1067   
1068   while ( ( dig = static_cast<AliMUONVDigit*>(next()) ) )
1069     {
1070     GetDigitsData(0)->Fill(dig->DetElemId());
1071     GetDigitsData(1)->Fill(dig->ADC());
1072     }
1073 }
1074
1075 //____________________________________________________________________________
1076 void AliMUONQADataMakerRec::MakeRecPoints(TTree* clustersTree)
1077 {
1078         /// Fill histograms from treeR
1079
1080   // Do nothing in case of calibration event
1081   if ( GetRecoParam()->GetEventSpecie() == AliRecoParam::kCalib ) return;
1082         
1083   GetRecPointsData(kNAnalyzedEvents)->Fill(1.);
1084
1085   MakeRecPointsTracker(clustersTree);
1086         MakeRecPointsTrigger(clustersTree);
1087 }
1088
1089 //____________________________________________________________________________
1090 void AliMUONQADataMakerRec::MakeRecPointsTracker(TTree* clustersTree)
1091 {
1092         /// Fill histograms related to tracker clusters 
1093         
1094         // First things first : do we have clusters in the TreeR ?
1095         // In "normal" production mode, it should be perfectly normal
1096         // *not* to have them.
1097         // But if for some reason we de-activated the combined tracking,
1098         // then we have clusters in TreeR, so let's take that opportunity
1099         // to QA them...
1100         
1101         if (!fClusterStore)
1102         {
1103                 AliCodeTimerAuto("ClusterStore creation");
1104                 fClusterStore = AliMUONVClusterStore::Create(*clustersTree);
1105                 if (!fClusterStore) 
1106                 {
1107                         return;
1108                 }
1109         }
1110         
1111         AliCodeTimerAuto("");
1112         
1113         fClusterStore->Connect(*clustersTree,kFALSE);
1114         clustersTree->GetEvent(0);
1115
1116         TIter next(fClusterStore->CreateIterator());
1117         AliMUONVCluster* cluster;
1118         
1119         while ( ( cluster = static_cast<AliMUONVCluster*>(next()) ) )
1120         {
1121                 Int_t detElemId = cluster->GetDetElemId();
1122                 Int_t chamberId = AliMpDEManager::GetChamberId(detElemId);
1123                 
1124                 GetRecPointsData(kTrackerNumberOfClustersPerDE)->Fill(detElemId);
1125                 GetRecPointsData(kTrackerClusterChargePerDE+detElemId)->Fill(cluster->GetCharge());
1126                 GetRecPointsData(kTrackerClusterMultiplicityPerDE+detElemId)->Fill(cluster->GetNDigits());
1127
1128                 GetRecPointsData(kTrackerNumberOfClustersPerChamber)->Fill(chamberId);
1129                 GetRecPointsData(kTrackerClusterChargePerChamber+chamberId)->Fill(cluster->GetCharge());
1130                 GetRecPointsData(kTrackerClusterMultiplicityPerChamber+chamberId)->Fill(cluster->GetNDigits());
1131                 GetRecPointsData(kTrackerClusterHitMapPerChamber+chamberId)->Fill(cluster->GetX(),cluster->GetY());
1132                 
1133         }
1134         
1135         fClusterStore->Clear();
1136 }
1137
1138 //____________________________________________________________________________
1139 void AliMUONQADataMakerRec::MakeRecPointsTrigger(TTree* clustersTree)
1140 {
1141         /// makes data from trigger response
1142       
1143     // Fired pads info
1144     fDigitStore->Clear();
1145
1146     if (!fTriggerStore) fTriggerStore = AliMUONVTriggerStore::Create(*clustersTree);
1147     fTriggerStore->Clear();
1148     fTriggerStore->Connect(*clustersTree, false);
1149     clustersTree->GetEvent(0);
1150
1151     AliMUONLocalTrigger* locTrg;
1152     TIter nextLoc(fTriggerStore->CreateLocalIterator());
1153
1154     while ( ( locTrg = static_cast<AliMUONLocalTrigger*>(nextLoc()) ) ) 
1155     {
1156       if (locTrg->IsNull()) continue;
1157    
1158       TArrayS xyPattern[2];
1159       locTrg->GetXPattern(xyPattern[0]);
1160       locTrg->GetYPattern(xyPattern[1]);
1161
1162       Int_t nBoard = locTrg->LoCircuit();
1163
1164       Bool_t xTrig=locTrg->IsTrigX();
1165       Bool_t yTrig=locTrg->IsTrigY();
1166     
1167       if (xTrig && yTrig)
1168         ((TH1F*)GetRecPointsData(kTriggeredBoards))->Fill(nBoard);
1169
1170       fDigitMaker->TriggerDigits(nBoard, xyPattern, *fDigitStore);
1171     }
1172
1173     TIter nextDigit(fDigitStore->CreateIterator());
1174     AliMUONVDigit* mDigit;
1175     while ( ( mDigit = static_cast<AliMUONVDigit*>(nextDigit()) ) )
1176     {
1177       Int_t detElemId = mDigit->DetElemId();
1178       Int_t ch = detElemId/100;
1179       Int_t localBoard = mDigit->ManuId();
1180       Int_t channel = mDigit->ManuChannel();
1181       Int_t cathode = mDigit->Cathode();
1182       Int_t iChamber = ch - 11;
1183       
1184       ((TH2F*)GetRecPointsData(kTriggerDigits + AliMpConstants::NofTriggerChambers()*cathode + iChamber))
1185         ->Fill(localBoard, channel);
1186     }
1187 }
1188
1189 //____________________________________________________________________________
1190 void AliMUONQADataMakerRec::MakeESDs(AliESDEvent* esd)
1191 {
1192   /// make QA data from ESDs
1193
1194   // Do nothing in case of calibration event
1195   if ( GetRecoParam()->GetEventSpecie() == AliRecoParam::kCalib ) return;
1196   
1197   // load ESD event in the interface
1198   AliMUONESDInterface esdInterface;
1199   if (GetMUONRecoParam()) AliMUONESDInterface::ResetTracker(GetMUONRecoParam());
1200   else AliError("Unable to get recoParam: use default ones for residual calculation");
1201   esdInterface.LoadEvent(*esd);
1202   
1203   GetESDsData(kESDnTracks)->Fill(esdInterface.GetNTracks());
1204   
1205   Int_t nTrackMatchTrig = 0;
1206   
1207   // loop over tracks
1208   Int_t nTracks = (Int_t) esd->GetNumberOfMuonTracks(); 
1209   for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack) {
1210     
1211     // get the ESD track and skip "ghosts"
1212     AliESDMuonTrack* esdTrack = esd->GetMuonTrack(iTrack);
1213     if (!esdTrack->ContainTrackerData()) continue;
1214     
1215     // get corresponding MUON track
1216     AliMUONTrack* track = esdInterface.FindTrack(esdTrack->GetUniqueID());
1217     
1218     if (esdTrack->ContainTriggerData()) nTrackMatchTrig++;
1219     
1220     GetESDsData(kESDMomentum)->Fill(esdTrack->P());
1221     GetESDsData(kESDPt)->Fill(esdTrack->Pt());
1222     GetESDsData(kESDRapidity)->Fill(esdTrack->Y());
1223     GetESDsData(kESDChi2)->Fill(track->GetNormalizedChi2());
1224     GetESDsData(kESDProbChi2)->Fill(TMath::Prob(track->GetGlobalChi2(),track->GetNDF()));
1225     GetESDsData(kESDThetaX)->Fill(esdTrack->GetThetaXUncorrected() / TMath::Pi() * 180.);
1226     GetESDsData(kESDThetaY)->Fill(esdTrack->GetThetaYUncorrected() / TMath::Pi() * 180.);
1227     GetESDsData(kESDnClustersPerTrack)->Fill(track->GetNClusters());
1228     
1229     // loop over clusters
1230     AliMUONTrackParam* trackParam = static_cast<AliMUONTrackParam*>(track->GetTrackParamAtCluster()->First());
1231     while (trackParam) {
1232       
1233       AliMUONVCluster* cluster = trackParam->GetClusterPtr();
1234       Int_t chId = cluster->GetChamberId();
1235       Int_t deID = cluster->GetDetElemId();
1236       Double_t residualX = cluster->GetX() - trackParam->GetNonBendingCoor();
1237       Double_t residualY = cluster->GetY() - trackParam->GetBendingCoor();
1238       Double_t sigmaResidualX2 = cluster->GetErrX2() - trackParam->GetCovariances()(0,0);
1239       Double_t sigmaResidualY2 = cluster->GetErrY2() - trackParam->GetCovariances()(2,2);
1240       Double_t localChi2X = (sigmaResidualX2 > 0.) ? residualX*residualX/sigmaResidualX2 : 0.;
1241       Double_t localChi2Y = (sigmaResidualY2 > 0.) ? residualY*residualY/sigmaResidualY2 : 0.;
1242       Double_t localChi2 = 0.5 * trackParam->GetLocalChi2();
1243       
1244       GetESDsData(kESDClusterHitMap+chId)->Fill(cluster->GetX(), cluster->GetY());
1245       
1246       GetESDsData(kESDnTotClustersPerCh)->Fill(chId);
1247       GetESDsData(kESDnTotClustersPerDE)->Fill(deID);
1248       
1249       GetESDsData(kESDClusterChargeInCh+chId)->Fill(cluster->GetCharge());
1250       GetESDsData(kESDSumClusterChargePerDE)->Fill(deID, cluster->GetCharge());
1251       
1252       if (cluster->GetNDigits() > 0) { // discard clusters with pad not stored in ESD
1253         GetESDsData(kESDnTotFullClustersPerDE)->Fill(deID);
1254         GetESDsData(kESDClusterSizeInCh+chId)->Fill(cluster->GetNDigits());
1255         GetESDsData(kESDSumClusterSizePerDE)->Fill(deID, cluster->GetNDigits());
1256       }
1257       
1258       GetESDsData(kESDResidualXInCh+chId)->Fill(residualX);
1259       GetESDsData(kESDResidualYInCh+chId)->Fill(residualY);
1260       GetESDsData(kESDSumResidualXPerDE)->Fill(deID, residualX);
1261       GetESDsData(kESDSumResidualYPerDE)->Fill(deID, residualY);
1262       GetESDsData(kESDSumResidualX2PerDE)->Fill(deID, residualX*residualX);
1263       GetESDsData(kESDSumResidualY2PerDE)->Fill(deID, residualY*residualY);
1264       
1265       GetESDsData(kESDLocalChi2XInCh+chId)->Fill(localChi2X);
1266       GetESDsData(kESDLocalChi2YInCh+chId)->Fill(localChi2Y);
1267       GetESDsData(kESDLocalChi2InCh+chId)->Fill(localChi2);
1268       GetESDsData(kESDSumLocalChi2XPerDE)->Fill(deID, localChi2X);
1269       GetESDsData(kESDSumLocalChi2YPerDE)->Fill(deID, localChi2Y);
1270       GetESDsData(kESDSumLocalChi2PerDE)->Fill(deID, localChi2);
1271       
1272       trackParam = static_cast<AliMUONTrackParam*>(track->GetTrackParamAtCluster()->After(trackParam));
1273     }
1274     
1275   }
1276   
1277   GetESDsData(kESDMatchTrig)->Fill(nTrackMatchTrig);
1278   
1279 }
1280
1281 //____________________________________________________________________________ 
1282 void AliMUONQADataMakerRec::StartOfDetectorCycle()
1283 {
1284     /// Detector specific actions at start of cycle
1285   
1286 }
1287
1288 //____________________________________________________________________________ 
1289 void AliMUONQADataMakerRec::DisplayTriggerInfo(AliQAv1::TASKINDEX_t task)
1290 {
1291   //
1292   /// Display trigger information in a user-friendly way:
1293   /// from local board and strip numbers to their position on chambers
1294   //
1295   
1296
1297   if(task!=AliQAv1::kRECPOINTS && task!=AliQAv1::kRAWS) return;
1298
1299   // check we get histograms, otherwise return right now
1300   if ( task == AliQAv1::kRECPOINTS )
1301   {
1302     if ( !GetRecPointsData(kTriggerDigits) ) return;
1303   }
1304   
1305   if ( task == AliQAv1::kRAWS ) 
1306   {
1307     if ( !GetRawsData(kTriggerScalers) ) return;
1308   }
1309   
1310   AliMUONTriggerDisplay triggerDisplay;
1311   
1312   TH2F* histoStrips=0x0;
1313   TH2F* histoDisplayStrips=0x0;
1314   AliMUONTriggerDisplay::EDisplayOption displayOption = AliMUONTriggerDisplay::kDefaultDisplay;
1315   if(task == AliQAv1::kRAWS) displayOption = AliMUONTriggerDisplay::kNormalizeToArea;
1316
1317   for (Int_t iCath = 0; iCath < AliMpConstants::NofCathodes(); iCath++)
1318   {    
1319     for (Int_t iChamber = 0; iChamber < AliMpConstants::NofTriggerChambers(); iChamber++)
1320     {
1321       if(task==AliQAv1::kRECPOINTS){
1322         histoStrips = (TH2F*)GetRecPointsData(kTriggerDigits + AliMpConstants::NofTriggerChambers()*iCath + iChamber);
1323       }
1324       else if(task==AliQAv1::kRAWS){
1325         histoStrips = (TH2F*)GetRawsData(kTriggerScalers + AliMpConstants::NofTriggerChambers()*iCath + iChamber);
1326       }
1327
1328       if(histoStrips->GetEntries()==0) continue; // No events found => No need to display
1329
1330       if(task==AliQAv1::kRECPOINTS){
1331         histoDisplayStrips = (TH2F*)GetRecPointsData(kTriggerDigitsDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber);
1332       }
1333       else if(task==AliQAv1::kRAWS){
1334         histoDisplayStrips = (TH2F*)GetRawsData(kTriggerScalersDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber);
1335       }
1336
1337       triggerDisplay.FillDisplayHistogram(histoStrips, histoDisplayStrips,
1338                                           AliMUONTriggerDisplay::kDisplayStrips, iCath, iChamber, displayOption);
1339
1340       Float_t scaleValue = 0.;
1341       if(task==AliQAv1::kRAWS) {
1342         scaleValue = ((TH1F*)GetRawsData(kTriggerScalersTime))->GetBinContent(1);
1343       }
1344       else if ( task == AliQAv1::kRECPOINTS ) {
1345         scaleValue = GetRecPointsData(kNAnalyzedEvents)->GetBinContent(1);
1346       }
1347       if(scaleValue>0.) histoDisplayStrips->Scale(1./scaleValue);
1348     } // iChamber
1349   } // iCath
1350
1351   if(task==AliQAv1::kRECPOINTS){
1352     TH1F* histoBoards = (TH1F*)GetRecPointsData(kTriggeredBoards);
1353     TH2F* histoDisplayBoards = (TH2F*)GetRecPointsData(kTriggerBoardsDisplay);
1354     triggerDisplay.FillDisplayHistogram(histoBoards, histoDisplayBoards, AliMUONTriggerDisplay::kDisplayBoards, 0, 0);
1355     Float_t scaleValue = GetRecPointsData(kNAnalyzedEvents)->GetBinContent(1);
1356     if(scaleValue>0.) histoDisplayBoards->Scale(1./scaleValue);
1357   }
1358 }
1359
1360
1361 //_____________________________________________________________________________
1362 Bool_t 
1363 AliMUONQADataMakerRec::FillTriggerDCSHistos()
1364 {
1365   /// Get HV and currents values for one trigger chamber
1366   
1367   AliCodeTimerAuto("");
1368
1369   AliMUONCalibrationData calibrationData(AliCDBManager::Instance()->GetRun());
1370
1371   TMap* triggerDcsMap = calibrationData.TriggerDCS();
1372
1373   if ( !triggerDcsMap ) 
1374   {
1375     AliError("Cannot fill DCS histos, as triggerDcsMap is NULL");
1376     return kFALSE;
1377   }
1378
1379   AliMpDEIterator deIt;
1380   
1381   deIt.First();
1382     
1383   AliMpDCSNamer triggerDcsNamer("TRIGGER");
1384
1385   TH2* currHisto = 0x0;
1386   Int_t histoIndex = 0;
1387   TString histoName, histoTitle;
1388
1389   Bool_t error = kFALSE;
1390   Bool_t expert   = kTRUE;
1391   Bool_t saveCorr = kTRUE;
1392   Bool_t image    = kTRUE;
1393   
1394   while ( !deIt.IsDone() )
1395   {
1396     Int_t detElemId = deIt.CurrentDEId();
1397     
1398     if ( AliMpDEManager::GetStationType(detElemId) == AliMp::kStationTrigger) {
1399
1400       Int_t iChamber = AliMpDEManager::GetChamberId(detElemId);
1401       Int_t slat = detElemId%100;
1402
1403       for(Int_t iMeas=0; iMeas<AliMpDCSNamer::kNDCSMeas; iMeas++){
1404         TString currAlias = triggerDcsNamer.DCSChannelName(detElemId, 0, iMeas);
1405
1406         AliDebug(AliQAv1::GetQADebugLevel(), Form("\nDetElemId %i   dcsAlias %s", detElemId, currAlias.Data()));
1407
1408         TPair* triggerDcsPair = static_cast<TPair*>(triggerDcsMap->FindObject(currAlias.Data()));
1409
1410         if (!triggerDcsPair)
1411         {
1412           AliError(Form("Did not find expected alias (%s) for DE %d",
1413                         currAlias.Data(),detElemId));  
1414           error = kTRUE;
1415         }
1416         else
1417         {
1418           TObjArray* values = static_cast<TObjArray*>(triggerDcsPair->Value());
1419           if (!values)
1420           {
1421             AliError(Form("Could not get values for alias %s",currAlias.Data()));
1422             error = kTRUE;
1423           }
1424           else
1425           {
1426             TIter next(values);
1427             AliDCSValue* val = 0x0;
1428
1429             Int_t ich = iChamber - AliMpConstants::NofTrackingChambers();
1430
1431             switch(iMeas){
1432             case AliMpDCSNamer::kDCSI:
1433               histoIndex = kTriggerRPCi + ich;
1434               histoName = Form("hRPCIChamber%i", 11+ich);
1435               histoTitle = Form("Chamber %i: RPC Currents (#muA)", 11+ich);
1436               break;
1437             case AliMpDCSNamer::kDCSHV:
1438               histoIndex = kTriggerRPChv + ich;
1439               histoName = Form("hRPCHVChamber%i", 11+ich);
1440               histoTitle = Form("Chamber %i: RPC HV (V)", 11+ich);
1441               break;
1442             }
1443
1444             currHisto = (TH2F*) GetRawsData(histoIndex);
1445
1446             if(!currHisto){
1447               Int_t npoints = values->GetEntries();
1448               TArrayF axisSlat(18+1), axisTime(npoints+1);
1449               for(Int_t islat=0; islat<=18; islat++){
1450                 axisSlat[islat] = -0.5 + (Float_t)islat;
1451               }
1452               for(Int_t ientry=0; ientry<npoints; ientry++){
1453                 val = static_cast<AliDCSValue*>(values->At(ientry));
1454                 axisTime[ientry] = val->GetTimeStamp();
1455               }
1456               axisTime[npoints] = val->GetTimeStamp() + 100;
1457
1458               currHisto  = new TH2F(histoName.Data(), histoTitle.Data(),
1459                                     npoints, axisTime.GetArray(),
1460                                     18, axisSlat.GetArray());
1461               currHisto->GetXaxis()->SetTitle("Time");
1462               currHisto->GetXaxis()->SetTimeDisplay(1);
1463               //currHisto->GetXaxis()->SetTimeFormat("%d%b%y %H:%M:%S");
1464               currHisto->GetYaxis()->SetTitle("RPC");
1465               Add2RawsList(currHisto, histoIndex, !expert, image, !saveCorr);
1466             }
1467
1468             while ( ( val = static_cast<AliDCSValue*>(next()) ) )
1469             {
1470               Float_t hvi = val->GetFloat();
1471
1472               AliDebug(AliQAv1::GetQADebugLevel(), Form("Value %f", hvi));
1473  
1474               currHisto->Fill(1.0001 * val->GetTimeStamp(), slat, hvi);
1475             } // loop on values
1476           } // if (!values)
1477         } // if (!triggerDcsPair)
1478       } // loop on measured types (HV and currents)
1479     } // if (stationType == kStationTrigger)
1480
1481     deIt.Next();
1482   }
1483   return error;
1484 }
1485
1486 //____________________________________________________________________________ 
1487 AliMUONVTrackerData* AliMUONQADataMakerRec::GetTrackerData() const
1488
1489 /// Return tracker data
1490   
1491   return fTrackerDataMaker->Data(); 
1492   
1493 }