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