]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONQADataMakerRec.cxx
update QA of RecPoints with histograms taken from QA of ESD, namely:
[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 "AliMUONTrackerCalibratedDataMaker.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(AliQA::GetDetName(AliQA::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(AliQA::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 == AliQA::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 == AliQA::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 == AliQA::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(AliQA::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   const Bool_t fastDecoder(kTRUE);
532
533   fTrackerDataMaker = new AliMUONTrackerCalibratedDataMaker(GetRecoParam(),
534                                                             AliCDBManager::Instance()->GetRun(),
535                                                             0x0,
536                                                             AliCDBManager::Instance()->GetDefaultStorage()->GetURI(),
537                                                             "NOGAIN",
538                                                             histogram,
539                                                             0.0,0.0,
540                                                             fastDecoder);
541                 
542   fTrackerDataMaker->Data()->DisableChannelLevel(); // to save up disk space, we only store starting at the manu level
543         
544   fTrackerDataMaker->SetRunning(kTRUE);
545   
546         fIsInitRaws = kTRUE;
547 }
548
549 //____________________________________________________________________________ 
550 void AliMUONQADataMakerRec::InitRecPoints()
551 {
552         /// create Reconstructed Points histograms in RecPoints subdir
553   
554   AliCodeTimerAuto("");
555   
556         InitRecPointsTrigger();
557         InitRecPointsTracker();
558 }
559
560 //____________________________________________________________________________ 
561 void AliMUONQADataMakerRec::InitRecPointsTracker()
562 {
563   /// create Reconstructed Points histograms in RecPoints subdir for the
564   /// MUON tracker subsystem.
565   
566   AliCodeTimerAuto("");
567   
568   Bool_t forExpert(kTRUE);
569   TH1I *h1I;
570   TH1F *h1F;
571   TH2F *h2F;
572   
573   // histograms per chamber
574   Int_t nCh = AliMpConstants::NofTrackingChambers();
575   for ( Int_t i = 0; i < nCh; ++i ) 
576   {
577     h1I = new TH1I(Form("hTrackerClusterMultiplicityForChamber%d",i+1), Form("cluster size distribution in chamber %d;size (n_{pads})",i+1), 100,0,100);
578     Add2RecPointsList(h1I,kTrackerClusterMultiplicityPerChamber+i,forExpert);
579     
580     h1I = new TH1I(Form("hTrackerClusterChargeForChamber%d",i+1), Form("cluster charge distribution in chamber %d;charge (ADC counts)",i+1), 500,0,5000);
581     Add2RecPointsList(h1I,kTrackerClusterChargePerChamber+i,forExpert);
582     
583     Float_t rMax = AliMUONConstants::Rmax(i/2);
584     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);
585     Add2RecPointsList(h2F, kTrackerClusterHitMapPerChamber+i,forExpert);
586   }
587   
588   // summary histograms per chamber
589   h1I = new TH1I("hTrackerNumberOfClustersPerChamber", "Number of clusters per chamber;chamber ID;n_{clusters}", nCh,-0.5,nCh-0.5);
590   Add2RecPointsList(h1I,kTrackerNumberOfClustersPerChamber,forExpert);
591   
592   h1F = new TH1F("hTrackerClusterMultiplicityPerChMean", "cluster mean size per chamber;chamber ID;<size> (n_{pads})", nCh,-0.5,nCh-0.5);
593   h1F->SetOption("P");
594   h1F->SetMarkerStyle(kFullDotMedium);
595   h1F->SetMarkerColor(kRed);
596   Add2RecPointsList(h1F, kTrackerClusterMultiplicityPerChMean,forExpert);
597   
598   h1F = new TH1F("hTrackerClusterMultiplicityPerChSigma", "cluster size dispersion per chamber;chamber ID;#sigma_{size} (n_{pads})", nCh,-0.5,nCh-0.5);
599   h1F->SetOption("P");
600   h1F->SetMarkerStyle(kFullDotMedium);
601   h1F->SetMarkerColor(kRed);
602   Add2RecPointsList(h1F, kTrackerClusterMultiplicityPerChSigma,forExpert);
603   
604   h1F = new TH1F("hTrackerClusterChargePerChMean", "cluster mean charge per chamber;chamber ID;<charge> (ADC counts)", nCh,-0.5,nCh-0.5);
605   h1F->SetOption("P");
606   h1F->SetMarkerStyle(kFullDotMedium);
607   h1F->SetMarkerColor(kRed);
608   Add2RecPointsList(h1F, kTrackerClusterChargePerChMean,forExpert);
609   
610   h1F = new TH1F("hTrackerClusterChargePerChSigma", "cluster charge dispersion per chamber;chamber ID;#sigma_{charge} (ADC counts)", nCh,-0.5,nCh-0.5);
611   h1F->SetOption("P");
612   h1F->SetMarkerStyle(kFullDotMedium);
613   h1F->SetMarkerColor(kRed);
614   Add2RecPointsList(h1F, kTrackerClusterChargePerChSigma,forExpert);
615   
616   // histograms per DE
617   Int_t ndes(0);
618   AliMpDEIterator it;
619   it.First();
620   while ( !it.IsDone())
621   {
622     Int_t detElemId = it.CurrentDEId();
623     
624     if ( AliMpDEManager::GetStationType(detElemId) != AliMp::kStationTrigger )
625     {
626       ndes = TMath::Max(ndes,detElemId);
627       
628       h1I = new TH1I(Form("hTrackerClusterMultiplicityForDE%04d",detElemId), Form("cluster size distribution in detection element %d;size (n_{pads})",detElemId), 100,0,100);
629       Add2RecPointsList(h1I,kTrackerClusterMultiplicityPerDE+detElemId,forExpert);
630       
631       h1I = new TH1I(Form("hTrackerClusterChargeForDE%04d",detElemId), Form("cluster charge distribution in detection element %d;charge (ADC counts)",detElemId), 500,0,5000);
632       Add2RecPointsList(h1I,kTrackerClusterChargePerDE+detElemId,forExpert);
633     }
634     
635     it.Next();
636   }
637   
638   // summary histograms per DE
639   h1I = new TH1I("hTrackerNumberOfClustersPerDE", "Number of clusters per detection element;DetElem ID;n_{clusters}", ndes+1,-0.5,ndes+0.5);
640   Add2RecPointsList(h1I, kTrackerNumberOfClustersPerDE,!forExpert);
641   
642   h1F = new TH1F("hTrackerClusterMultiplicityPerDEMean", "cluster mean size per DE;DetElem ID;<size> (n_{pads})", ndes+1,-0.5,ndes+0.5);
643   h1F->SetOption("P");
644   h1F->SetMarkerStyle(kFullDotMedium);
645   h1F->SetMarkerColor(kRed);
646   Add2RecPointsList(h1F, kTrackerClusterMultiplicityPerDEMean,forExpert);
647   
648   h1F = new TH1F("hTrackerClusterChargePerDEMean", "cluster mean charge per DE;DetElem ID;<charge> (ADC counts)", ndes+1,-0.5,ndes+0.5);
649   h1F->SetOption("P");
650   h1F->SetMarkerStyle(kFullDotMedium);
651   h1F->SetMarkerColor(kRed);
652   Add2RecPointsList(h1F, kTrackerClusterChargePerDEMean,forExpert);
653   
654   fIsInitRecPointsTracker=kTRUE;
655 }
656
657 //____________________________________________________________________________ 
658 void AliMUONQADataMakerRec::InitRecPointsTrigger()
659 {
660         /// create Reconstructed Points histograms in RecPoints subdir for the
661         /// MUON Trigger subsystem.
662         
663   Bool_t forExpert(kTRUE);
664   
665     TH3F* h0 = new TH3F("hTriggerDigitsBendPlane", "Trigger digits in bending plane",
666                         4, 10.5, 14.5,
667                         234, 0.5, 234.5,
668                         16, -0.5, 15.5);
669     h0->GetXaxis()->SetTitle("Chamber");
670     h0->GetYaxis()->SetTitle("Board");
671     h0->GetZaxis()->SetTitle("Strip");
672     Add2RecPointsList(h0, kTriggerDigitsBendPlane,forExpert);
673
674     TH3F* h1 = new TH3F("hTriggerDigitsNonBendPlane", "Trigger digits in non-bending plane",
675                         4, 10.5, 14.5,
676                         234, 0.5, 234.5,
677                         16, -0.5, 15.5);
678     h1->GetXaxis()->SetTitle("Chamber");
679     h1->GetYaxis()->SetTitle("Board");
680     h1->GetZaxis()->SetTitle("Strip");
681     Add2RecPointsList(h1, kTriggerDigitsNonBendPlane,forExpert);
682
683     TH1F* h2 = new TH1F("hTriggeredBoards", "Triggered boards", 234, 0.5, 234.5);
684     Add2RecPointsList(h2, kTriggeredBoards,forExpert);
685
686     AliMUONTriggerDisplay triggerDisplay;
687     TString histoName, histoTitle;
688     for(Int_t iCath=0; iCath<AliMpConstants::NofCathodes(); iCath++){
689       TString cathName = ( iCath==0 ) ? "BendPlane" : "NonBendPlane";
690       for(Int_t iChamber=0; iChamber<AliMpConstants::NofTriggerChambers(); iChamber++){
691         histoName = Form("hTriggerDigits%sChamber%i", cathName.Data(), 11+iChamber);
692         histoTitle = Form("Chamber %i: Fired pads %s", 11+iChamber, cathName.Data());
693         TH2F* h3 = (TH2F*)triggerDisplay.GetEmptyDisplayHisto(histoName, AliMUONTriggerDisplay::kDisplayStrips, 
694                                                               iCath, iChamber, histoTitle);
695         Add2RecPointsList(h3, kTriggerDigitsDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber,!forExpert);
696       }
697     }
698
699     TH2F* h4 = (TH2F*)triggerDisplay.GetEmptyDisplayHisto("hFiredBoardsDisplay", AliMUONTriggerDisplay::kDisplayBoards,
700                                                           0, 0, "Fired boards");
701     Add2RecPointsList(h4, kTriggerBoardsDisplay,!forExpert);
702         
703         fIsInitRecPointsTrigger = kTRUE;
704 }
705
706
707 //____________________________________________________________________________ 
708 void AliMUONQADataMakerRec::InitESDs()
709 {
710   ///create ESDs histograms in ESDs subdir
711   
712   Bool_t forExpert(kTRUE);
713   
714   Int_t nCh = AliMUONConstants::NTrackingCh();
715   Int_t nDE = 1100;
716   
717   // track info
718   TH1F* hESDnTracks = new TH1F("hESDnTracks", "number of tracks;n_{tracks}", 20, 0., 20.);
719   Add2ESDsList(hESDnTracks, kESDnTracks,!forExpert);
720
721   TH1F* hESDMatchTrig = new TH1F("hESDMatchTrig", "number of tracks matched with trigger;n_{tracks}", 20, 0., 20.);
722   Add2ESDsList(hESDMatchTrig, kESDMatchTrig,!forExpert);
723   
724   TH1F* hESDMomentum = new TH1F("hESDMomentum", "momentum distribution;p (GeV/c)", 300, 0., 300);
725   Add2ESDsList(hESDMomentum, kESDMomentum,forExpert);
726
727   TH1F* hESDPt = new TH1F("hESDPt", "transverse momentum distribution;p_{t} (GeV/c)", 200, 0., 50);
728   Add2ESDsList(hESDPt, kESDPt,forExpert);
729
730   TH1F* hESDRapidity = new TH1F("hESDRapidity", "rapidity distribution;rapidity", 200, -4.5, -2.);
731   Add2ESDsList(hESDRapidity, kESDRapidity,forExpert);
732
733   TH1F* hESDChi2 = new TH1F("hESDChi2", "normalized #chi^{2} distribution;#chi^{2} / ndf", 500, 0., 50.);
734   Add2ESDsList(hESDChi2, kESDChi2,forExpert);
735   
736   TH1F* hESDProbChi2 = new TH1F("hESDProbChi2", "distribution of probability of #chi^{2};prob(#chi^{2})", 100, 0., 1.);
737   Add2ESDsList(hESDProbChi2, kESDProbChi2,forExpert);
738   
739   // cluster info
740   for (Int_t i = 0; i < nCh; i++) {
741     Float_t rMax = AliMUONConstants::Rmax(i/2);
742     TH2F* hESDClusterHitMap = new TH2F(Form("hESDClusterHitMap%d",i+1), Form("cluster position distribution in chamber %d;X (cm);Y (cm)",i+1),
743                                        100, -rMax, rMax, 100, -rMax, rMax);
744     Add2ESDsList(hESDClusterHitMap, kESDClusterHitMap+i,forExpert);
745   }
746   
747   TH1F* hESDnClustersPerTrack = new TH1F("hESDnClustersPerTrack", "number of associated clusters per track;n_{clusters}", 20, 0., 20.);
748   Add2ESDsList(hESDnClustersPerTrack, kESDnClustersPerTrack,!forExpert);
749   
750   TH1F* hESDnClustersPerCh = new TH1F("hESDnClustersPerCh", "averaged number of clusters per chamber per track;chamber ID;<n_{clusters}>", nCh, -0.5, nCh-0.5);
751   hESDnClustersPerCh->SetFillColor(kRed);
752   Add2ESDsList(hESDnClustersPerCh, kESDnClustersPerCh,forExpert);
753   
754   TH1F* hESDnClustersPerDE = new TH1F("hESDnClustersPerDE", "averaged number of clusters per DE per track;DetElem ID;<n_{clusters}>", nDE+1, -0.5, nDE+0.5);
755   hESDnClustersPerDE->SetFillColor(kRed);
756   Add2ESDsList(hESDnClustersPerDE, kESDnClustersPerDE,forExpert);
757   
758   for (Int_t i = 0; i < nCh; i++) {
759     TH1F* hESDClusterChargeInCh = new TH1F(Form("hESDClusterChargeInCh%d",i+1), Form("cluster charge distribution in chamber %d;charge (ADC counts)",i+1), 500, 0., 5000.);
760     Add2ESDsList(hESDClusterChargeInCh, kESDClusterChargeInCh+i,forExpert);
761   }
762   
763   TH1F* hESDClusterChargePerChMean = new TH1F("hESDClusterChargePerChMean", "cluster mean charge per chamber;chamber ID;<charge> (ADC counts)", nCh, -0.5, nCh-0.5);
764   hESDClusterChargePerChMean->SetOption("P");
765   hESDClusterChargePerChMean->SetMarkerStyle(kFullDotMedium);
766   hESDClusterChargePerChMean->SetMarkerColor(kRed);
767   Add2ESDsList(hESDClusterChargePerChMean, kESDClusterChargePerChMean,forExpert);
768   
769   TH1F* hESDClusterChargePerChSigma = new TH1F("hESDClusterChargePerChSigma", "cluster charge dispersion per chamber;chamber ID;#sigma_{charge} (ADC counts)", nCh, -0.5, nCh-0.5);
770   hESDClusterChargePerChSigma->SetOption("P");
771   hESDClusterChargePerChSigma->SetMarkerStyle(kFullDotMedium);
772   hESDClusterChargePerChSigma->SetMarkerColor(kRed);
773   Add2ESDsList(hESDClusterChargePerChSigma, kESDClusterChargePerChSigma,forExpert);
774   
775   TH1F* hESDClusterChargePerDE = new TH1F("hESDClusterChargePerDE", "cluster mean charge per DE;DetElem ID;<charge> (ADC counts)", nDE+1, -0.5, nDE+0.5);
776   hESDClusterChargePerDE->SetOption("P");
777   hESDClusterChargePerDE->SetMarkerStyle(kFullDotMedium);
778   hESDClusterChargePerDE->SetMarkerColor(kRed);
779   Add2ESDsList(hESDClusterChargePerDE, kESDClusterChargePerDE,forExpert);
780   
781   for (Int_t i = 0; i < nCh; i++) {
782     TH1F* hESDClusterSizeInCh = new TH1F(Form("hESDClusterSizeInCh%d",i+1), Form("cluster size distribution in chamber %d;size (n_{pads})",i+1), 200, 0., 200.);
783     Add2ESDsList(hESDClusterSizeInCh, kESDClusterSizeInCh+i,forExpert);
784   }
785   
786   TH1F* hESDClusterSizePerChMean = new TH1F("hESDClusterSizePerChMean", "cluster mean size per chamber;chamber ID;<size> (n_{pads})", nCh, -0.5, nCh-0.5);
787   hESDClusterSizePerChMean->SetOption("P");
788   hESDClusterSizePerChMean->SetMarkerStyle(kFullDotMedium);
789   hESDClusterSizePerChMean->SetMarkerColor(kRed);
790   Add2ESDsList(hESDClusterSizePerChMean, kESDClusterSizePerChMean,forExpert);
791   
792   TH1F* hESDClusterSizePerChSigma = new TH1F("hESDClusterSizePerChSigma", "cluster size dispersion per chamber;chamber ID;#sigma_{size} (n_{pads})", nCh, -0.5, nCh-0.5);
793   hESDClusterSizePerChSigma->SetOption("P");
794   hESDClusterSizePerChSigma->SetMarkerStyle(kFullDotMedium);
795   hESDClusterSizePerChSigma->SetMarkerColor(kRed);
796   Add2ESDsList(hESDClusterSizePerChSigma, kESDClusterSizePerChSigma,forExpert);
797   
798   TH1F* hESDClusterSizePerDE = new TH1F("hESDClusterSizePerDE", "cluster mean size per DE;DetElem ID;<size> (n_{pads})", nDE+1, -0.5, nDE+0.5);
799   hESDClusterSizePerDE->SetOption("P");
800   hESDClusterSizePerDE->SetMarkerStyle(kFullDotMedium);
801   hESDClusterSizePerDE->SetMarkerColor(kRed);
802   Add2ESDsList(hESDClusterSizePerDE, kESDClusterSizePerDE,forExpert);
803   
804   // cluster - track info
805   for (Int_t i = 0; i < nCh; i++) {
806     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.);
807     Add2ESDsList(hESDResidualXInCh, kESDResidualXInCh+i,forExpert);
808     
809     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.);
810     Add2ESDsList(hESDResidualYInCh, kESDResidualYInCh+i,forExpert);
811     
812     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);
813     Add2ESDsList(hESDLocalChi2XInCh, kESDLocalChi2XInCh+i,forExpert);
814     
815     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);
816     Add2ESDsList(hESDLocalChi2YInCh, kESDLocalChi2YInCh+i,forExpert);
817   }
818   
819   TH1F* hESDResidualXPerChMean = new TH1F("hESDResidualXPerChMean", "cluster-track residual-X per Ch: mean;chamber ID;<#Delta_{X}> (cm)", nCh, -0.5, nCh-0.5);
820   hESDResidualXPerChMean->SetOption("P");
821   hESDResidualXPerChMean->SetMarkerStyle(kFullDotMedium);
822   hESDResidualXPerChMean->SetMarkerColor(kRed);
823   Add2ESDsList(hESDResidualXPerChMean, kESDResidualXPerChMean,forExpert);
824   
825   TH1F* hESDResidualYPerChMean = new TH1F("hESDResidualYPerChMean", "cluster-track residual-Y per Ch: mean;chamber ID;<#Delta_{Y}> (cm)", nCh, -0.5, nCh-0.5);
826   hESDResidualYPerChMean->SetOption("P");
827   hESDResidualYPerChMean->SetMarkerStyle(kFullDotMedium);
828   hESDResidualYPerChMean->SetMarkerColor(kRed);
829   Add2ESDsList(hESDResidualYPerChMean, kESDResidualYPerChMean,forExpert);
830   
831   TH1F* hESDResidualXPerChSigma = new TH1F("hESDResidualXPerChSigma", "cluster-track residual-X per Ch: sigma;chamber ID;#sigma_{X} (cm)", nCh, -0.5, nCh-0.5);
832   hESDResidualXPerChSigma->SetOption("P");
833   hESDResidualXPerChSigma->SetMarkerStyle(kFullDotMedium);
834   hESDResidualXPerChSigma->SetMarkerColor(kRed);
835   Add2ESDsList(hESDResidualXPerChSigma, kESDResidualXPerChSigma,forExpert);
836   
837   TH1F* hESDResidualYPerChSigma = new TH1F("hESDResidualYPerChSigma", "cluster-track residual-Y per Ch: sigma;chamber ID;#sigma_{Y} (cm)", nCh, -0.5, nCh-0.5);
838   hESDResidualYPerChSigma->SetOption("P");
839   hESDResidualYPerChSigma->SetMarkerStyle(kFullDotMedium);
840   hESDResidualYPerChSigma->SetMarkerColor(kRed);
841   Add2ESDsList(hESDResidualYPerChSigma, kESDResidualYPerChSigma,forExpert);
842   
843   TH1F* hESDLocalChi2XPerChMean = new TH1F("hESDLocalChi2XPerChMean", "local chi2-X per Ch: mean;chamber ID;<local #chi^{2}_{X}>", nCh, -0.5, nCh-0.5);
844   hESDLocalChi2XPerChMean->SetOption("P");
845   hESDLocalChi2XPerChMean->SetMarkerStyle(kFullDotMedium);
846   hESDLocalChi2XPerChMean->SetMarkerColor(kRed);
847   Add2ESDsList(hESDLocalChi2XPerChMean, kESDLocalChi2XPerChMean,forExpert);
848   
849   TH1F* hESDLocalChi2YPerChMean = new TH1F("hESDLocalChi2YPerChMean", "local chi2-Y per Ch: mean;chamber ID;<local #chi^{2}_{Y}>", nCh, -0.5, nCh-0.5);
850   hESDLocalChi2YPerChMean->SetOption("P");
851   hESDLocalChi2YPerChMean->SetMarkerStyle(kFullDotMedium);
852   hESDLocalChi2YPerChMean->SetMarkerColor(kRed);
853   Add2ESDsList(hESDLocalChi2YPerChMean, kESDLocalChi2YPerChMean,forExpert);
854   
855   TH1F* hESDResidualXPerDEMean = new TH1F("hESDResidualXPerDEMean", "cluster-track residual-X per DE: mean;DetElem ID;<#Delta_{X}> (cm)", nDE+1, -0.5, nDE+0.5);
856   hESDResidualXPerDEMean->SetOption("P");
857   hESDResidualXPerDEMean->SetMarkerStyle(kFullDotMedium);
858   hESDResidualXPerDEMean->SetMarkerColor(kRed);
859   Add2ESDsList(hESDResidualXPerDEMean, kESDResidualXPerDEMean,forExpert);
860   
861   TH1F* hESDResidualYPerDEMean = new TH1F("hESDResidualYPerDEMean", "cluster-track residual-Y per DE: mean;DetElem ID;<#Delta_{Y}> (cm)", nDE+1, -0.5, nDE+0.5);
862   hESDResidualYPerDEMean->SetOption("P");
863   hESDResidualYPerDEMean->SetMarkerStyle(kFullDotMedium);
864   hESDResidualYPerDEMean->SetMarkerColor(kRed);
865   Add2ESDsList(hESDResidualYPerDEMean, kESDResidualYPerDEMean,forExpert);
866   
867   TH1F* hESDResidualXPerDESigma = new TH1F("hESDResidualXPerDESigma", "cluster-track residual-X per DE: sigma;DetElem ID;#sigma_{X} (cm)", nDE+1, -0.5, nDE+0.5);
868   hESDResidualXPerDESigma->SetOption("P");
869   hESDResidualXPerDESigma->SetMarkerStyle(kFullDotMedium);
870   hESDResidualXPerDESigma->SetMarkerColor(kRed);
871   Add2ESDsList(hESDResidualXPerDESigma, kESDResidualXPerDESigma,forExpert);
872   
873   TH1F* hESDResidualYPerDESigma = new TH1F("hESDResidualYPerDESigma", "cluster-track residual-Y per DE: sigma;DetElem ID;#sigma_{Y} (cm)", nDE+1, -0.5, nDE+0.5);
874   hESDResidualYPerDESigma->SetOption("P");
875   hESDResidualYPerDESigma->SetMarkerStyle(kFullDotMedium);
876   hESDResidualYPerDESigma->SetMarkerColor(kRed);
877   Add2ESDsList(hESDResidualYPerDESigma, kESDResidualYPerDESigma,forExpert);
878   
879   TH1F* hESDLocalChi2XPerDEMean = new TH1F("hESDLocalChi2XPerDEMean", "local chi2-X per DE: mean;DetElem ID;<local #chi^{2}_{X}>", nDE+1, -0.5, nDE+0.5);
880   hESDLocalChi2XPerDEMean->SetOption("P");
881   hESDLocalChi2XPerDEMean->SetMarkerStyle(kFullDotMedium);
882   hESDLocalChi2XPerDEMean->SetMarkerColor(kRed);
883   Add2ESDsList(hESDLocalChi2XPerDEMean, kESDLocalChi2XPerDEMean,forExpert);
884   
885   TH1F* hESDLocalChi2YPerDEMean = new TH1F("hESDLocalChi2YPerDEMean", "local chi2-Y per DE: mean;DetElem ID;<local #chi^{2}_{Y}>", nDE+1, -0.5, nDE+0.5);
886   hESDLocalChi2YPerDEMean->SetOption("P");
887   hESDLocalChi2YPerDEMean->SetMarkerStyle(kFullDotMedium);
888   hESDLocalChi2YPerDEMean->SetMarkerColor(kRed);
889   Add2ESDsList(hESDLocalChi2YPerDEMean, kESDLocalChi2YPerDEMean,forExpert);
890   
891   // intermediate histograms
892   TH1F* hESDnTotClustersPerCh = new TH1F("hESDnTotClustersPerCh", "total number of associated clusters per chamber;chamber ID;#Sigma(n_{clusters})", nCh, -0.5, nCh-0.5);
893   Add2ESDsList(hESDnTotClustersPerCh, kESDnTotClustersPerCh, forExpert);
894   TH1F* hESDnTotClustersPerDE = new TH1F("hESDnTotClustersPerDE", "total number of associated clusters per DE;DetElem ID;#Sigma(n_{clusters})", nDE+1, -0.5, nDE+0.5);
895   Add2ESDsList(hESDnTotClustersPerDE, kESDnTotClustersPerDE, forExpert);
896   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);
897   Add2ESDsList(hESDnTotFullClustersPerDE, kESDnTotFullClustersPerDE, forExpert);
898   TH1F* hESDSumClusterChargePerDE = new TH1F("hESDSumClusterChargePerDE", "sum of cluster charge per DE;DetElem ID;#Sigma(charge) (ADC counts)", nDE+1, -0.5, nDE+0.5);
899   Add2ESDsList(hESDSumClusterChargePerDE, kESDSumClusterChargePerDE, forExpert);
900   TH1F* hESDSumClusterSizePerDE = new TH1F("hESDSumClusterSizePerDE", "sum of cluster size per DE;DetElem ID;#Sigma(size) (n_{pads})", nDE+1, -0.5, nDE+0.5);
901   Add2ESDsList(hESDSumClusterSizePerDE, kESDSumClusterSizePerDE, forExpert);
902   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);
903   Add2ESDsList(hESDSumResidualXPerDE, kESDSumResidualXPerDE, forExpert);
904   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);
905   Add2ESDsList(hESDSumResidualYPerDE, kESDSumResidualYPerDE, forExpert);
906   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);
907   Add2ESDsList(hESDSumResidualX2PerDE, kESDSumResidualX2PerDE, forExpert);
908   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);
909   Add2ESDsList(hESDSumResidualY2PerDE, kESDSumResidualY2PerDE, forExpert);
910   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);
911   Add2ESDsList(hESDSumLocalChi2XPerDE, kESDSumLocalChi2XPerDE, forExpert);
912   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);
913   Add2ESDsList(hESDSumLocalChi2YPerDE, kESDSumLocalChi2YPerDE, forExpert);
914   
915   fIsInitESDs =  kTRUE;
916 }
917
918 //____________________________________________________________________________
919 void AliMUONQADataMakerRec::MakeRaws(AliRawReader* rawReader)
920 {
921     /// make QA for rawdata
922
923     if ( ! fIsInitRaws ) {
924       AliWarningStream() 
925         << "Skipping function due to a failure in Init" << endl;
926       return;
927     }    
928
929   if ( rawReader->GetType() == AliRawEventHeaderBase::kPhysicsEvent ) 
930   {
931     rawReader->Reset();
932     MakeRawsTracker(rawReader);
933   }
934   
935   rawReader->Reset();    
936   MakeRawsTrigger(rawReader);
937 }
938
939 //____________________________________________________________________________
940 void AliMUONQADataMakerRec::MakeRawsTracker(AliRawReader* rawReader)
941 {
942         /// make QA for rawdata tracker
943         
944         ((AliMUONTrackerCalibratedDataMaker*)fTrackerDataMaker)->SetRawReader(rawReader);
945         
946         fTrackerDataMaker->ProcessEvent();
947 }
948
949 //____________________________________________________________________________
950 void AliMUONQADataMakerRec::MakeRawsTrigger(AliRawReader* rawReader)
951 {
952         /// make QA for rawdata trigger
953         
954     // Get trigger scalers
955
956     Int_t loCircuit=0;
957     AliMpCDB::LoadDDLStore();
958
959     AliMUONRawStreamTrigger rawStreamTrig(rawReader);
960     while (rawStreamTrig.NextDDL()) 
961     {
962       // If not a scaler event, do nothing
963       Bool_t scalerEvent =  rawReader->GetDataHeader()->GetL1TriggerMessage() & 0x1;
964       if(!scalerEvent) break;
965
966       AliMUONDDLTrigger* ddlTrigger = rawStreamTrig.GetDDLTrigger();
967       AliMUONDarcHeader* darcHeader = ddlTrigger->GetDarcHeader();
968
969       if (darcHeader->GetGlobalFlag()){
970         UInt_t nOfClocks = darcHeader->GetGlobalClock();
971         Double_t nOfSeconds = ((Double_t) nOfClocks) / 40e6; // 1 clock each 25 ns
972         ((TH1F*)GetRawsData(kTriggerScalersTime))->Fill(1., nOfSeconds);
973       }
974
975       Int_t nReg = darcHeader->GetRegHeaderEntries();
976     
977       for(Int_t iReg = 0; iReg < nReg ;iReg++)
978       {   //reg loop
979
980         // crate info  
981         AliMpTriggerCrate* crate = AliMpDDLStore::Instance()->
982           GetTriggerCrate(rawStreamTrig.GetDDL(), iReg);
983
984         AliMUONRegHeader* regHeader =  darcHeader->GetRegHeaderEntry(iReg);
985
986         // loop over local structures
987         Int_t nLocal = regHeader->GetLocalEntries();
988         for(Int_t iLocal = 0; iLocal < nLocal; iLocal++) 
989         {
990           AliMUONLocalStruct* localStruct = regHeader->GetLocalEntry(iLocal);
991         
992           // if card exist
993           if (!localStruct) continue;
994           
995           loCircuit = crate->GetLocalBoardId(localStruct->GetId());
996
997           if ( !loCircuit ) continue; // empty slot
998
999           AliMpLocalBoard* localBoard = AliMpDDLStore::Instance()->GetLocalBoard(loCircuit, false);
1000
1001           // skip copy cards
1002           if( !localBoard->IsNotified()) 
1003             continue;
1004
1005           Int_t cathode = localStruct->GetComptXY()%2;
1006
1007           ERaw hindex = (cathode==0) ? kTriggerScalersBP : kTriggerScalersNBP;
1008
1009           // loop over strips
1010           for (Int_t ibitxy = 0; ibitxy < 16; ++ibitxy) {
1011             if(localStruct->GetXY1(ibitxy) > 0)
1012               ((TH3F*)GetRawsData(hindex))->Fill(11+0, loCircuit, ibitxy, 2*localStruct->GetXY1(ibitxy));
1013             if(localStruct->GetXY2(ibitxy) > 0)
1014               ((TH3F*)GetRawsData(hindex))->Fill(11+1, loCircuit, ibitxy, 2*localStruct->GetXY2(ibitxy));
1015             if(localStruct->GetXY3(ibitxy) > 0)
1016               ((TH3F*)GetRawsData(hindex))->Fill(11+2, loCircuit, ibitxy, 2*localStruct->GetXY3(ibitxy));
1017             if(localStruct->GetXY4(ibitxy) > 0)
1018               ((TH3F*)GetRawsData(hindex))->Fill(11+3, loCircuit, ibitxy, 2*localStruct->GetXY4(ibitxy));
1019           } // loop on strips
1020         } // iLocal
1021       } // iReg
1022     } // NextDDL
1023 }
1024
1025 //____________________________________________________________________________
1026 void AliMUONQADataMakerRec::MakeRecPoints(TTree* clustersTree)
1027 {
1028         /// Fill histograms from treeR
1029         
1030         if (fIsInitRecPointsTracker) MakeRecPointsTracker(clustersTree);
1031         if (fIsInitRecPointsTrigger) MakeRecPointsTrigger(clustersTree);
1032 }
1033
1034 //____________________________________________________________________________
1035 void AliMUONQADataMakerRec::MakeRecPointsTracker(TTree* clustersTree)
1036 {
1037         /// Fill histograms related to tracker clusters 
1038         
1039         // First things first : do we have clusters in the TreeR ?
1040         // In "normal" production mode, it should be perfectly normal
1041         // *not* to have them.
1042         // But if for some reason we de-activated the combined tracking,
1043         // then we have clusters in TreeR, so let's take that opportunity
1044         // to QA them...
1045         
1046         if (!fClusterStore)
1047         {
1048                 AliCodeTimerAuto("ClusterStore creation");
1049                 fClusterStore = AliMUONVClusterStore::Create(*clustersTree);
1050                 if (!fClusterStore) 
1051                 {
1052                         fIsInitRecPointsTracker = kFALSE;
1053                         return;
1054                 }
1055         }
1056         
1057         AliCodeTimerAuto("");
1058         
1059         fClusterStore->Connect(*clustersTree,kFALSE);
1060         clustersTree->GetEvent(0);
1061
1062         TIter next(fClusterStore->CreateIterator());
1063         AliMUONVCluster* cluster;
1064         
1065         while ( ( cluster = static_cast<AliMUONVCluster*>(next()) ) )
1066         {
1067                 Int_t detElemId = cluster->GetDetElemId();
1068                 Int_t chamberId = AliMpDEManager::GetChamberId(detElemId);
1069                 
1070                 GetRecPointsData(kTrackerNumberOfClustersPerDE)->Fill(detElemId);
1071                 GetRecPointsData(kTrackerClusterChargePerDE+detElemId)->Fill(cluster->GetCharge());
1072                 GetRecPointsData(kTrackerClusterMultiplicityPerDE+detElemId)->Fill(cluster->GetNDigits());
1073
1074                 GetRecPointsData(kTrackerNumberOfClustersPerChamber)->Fill(chamberId);
1075                 GetRecPointsData(kTrackerClusterChargePerChamber+chamberId)->Fill(cluster->GetCharge());
1076                 GetRecPointsData(kTrackerClusterMultiplicityPerChamber+chamberId)->Fill(cluster->GetNDigits());
1077                 GetRecPointsData(kTrackerClusterHitMapPerChamber+chamberId)->Fill(cluster->GetX(),cluster->GetY());
1078                 
1079         }
1080         
1081         fClusterStore->Clear();
1082 }
1083
1084 //____________________________________________________________________________
1085 void AliMUONQADataMakerRec::MakeRecPointsTrigger(TTree* clustersTree)
1086 {
1087         /// makes data from trigger response
1088       
1089     // Fired pads info
1090     fDigitStore->Clear();
1091
1092     if (!fTriggerStore) fTriggerStore = AliMUONVTriggerStore::Create(*clustersTree);
1093     fTriggerStore->Clear();
1094     fTriggerStore->Connect(*clustersTree, false);
1095     clustersTree->GetEvent(0);
1096
1097     AliMUONLocalTrigger* locTrg;
1098     TIter nextLoc(fTriggerStore->CreateLocalIterator());
1099
1100     while ( ( locTrg = static_cast<AliMUONLocalTrigger*>(nextLoc()) ) ) 
1101     {
1102       if (locTrg->IsNull()) continue;
1103    
1104       TArrayS xyPattern[2];
1105       locTrg->GetXPattern(xyPattern[0]);
1106       locTrg->GetYPattern(xyPattern[1]);
1107
1108       Int_t nBoard = locTrg->LoCircuit();
1109
1110       Bool_t xTrig=locTrg->IsTrigX();
1111       Bool_t yTrig=locTrg->IsTrigY();
1112     
1113       if (xTrig && yTrig)
1114         ((TH1F*)GetRecPointsData(kTriggeredBoards))->Fill(nBoard);
1115
1116       fDigitMaker->TriggerDigits(nBoard, xyPattern, *fDigitStore);
1117     }
1118
1119     TIter nextDigit(fDigitStore->CreateIterator());
1120     AliMUONVDigit* mDigit;
1121     while ( ( mDigit = static_cast<AliMUONVDigit*>(nextDigit()) ) )
1122     {
1123       Int_t detElemId = mDigit->DetElemId();
1124       Int_t ch = detElemId/100;
1125       Int_t localBoard = mDigit->ManuId();
1126       Int_t channel = mDigit->ManuChannel();
1127       Int_t cathode = mDigit->Cathode();
1128       ERecPoints hindex 
1129         = ( cathode == 0 ) ? kTriggerDigitsBendPlane : kTriggerDigitsNonBendPlane;
1130       
1131       ((TH3F*)GetRecPointsData(hindex))->Fill(ch, localBoard, channel);
1132     }
1133 }
1134
1135 //____________________________________________________________________________
1136 void AliMUONQADataMakerRec::MakeESDs(AliESDEvent* esd)
1137 {
1138   /// make QA data from ESDs
1139   
1140   if ( ! fIsInitESDs ) {
1141     AliWarningStream() 
1142     << "Skipping function due to a failure in Init" << endl;
1143     return;
1144   }    
1145   
1146   // load ESD event in the interface
1147   AliMUONESDInterface esdInterface;
1148   if (GetRecoParam()) AliMUONESDInterface::ResetTracker(GetRecoParam());
1149   else AliError("Unable to get recoParam: use default ones for residual calculation");
1150   esdInterface.LoadEvent(*esd);
1151   
1152   GetESDsData(kESDnTracks)->Fill(esdInterface.GetNTracks());
1153   
1154   Int_t nTrackMatchTrig = 0;
1155   
1156   // loop over tracks
1157   Int_t nTracks = (Int_t) esd->GetNumberOfMuonTracks(); 
1158   for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack) {
1159     
1160     // get the ESD track and skip "ghosts"
1161     AliESDMuonTrack* esdTrack = esd->GetMuonTrack(iTrack);
1162     if (!esdTrack->ContainTrackerData()) continue;
1163     
1164     // get corresponding MUON track
1165     AliMUONTrack* track = esdInterface.FindTrack(esdTrack->GetUniqueID());
1166     
1167     if (esdTrack->ContainTriggerData()) nTrackMatchTrig++;
1168     
1169     GetESDsData(kESDMomentum)->Fill(esdTrack->P());
1170     GetESDsData(kESDPt)->Fill(esdTrack->Pt());
1171     GetESDsData(kESDRapidity)->Fill(esdTrack->Y());
1172     GetESDsData(kESDChi2)->Fill(track->GetNormalizedChi2());
1173     GetESDsData(kESDProbChi2)->Fill(TMath::Prob(track->GetGlobalChi2(),track->GetNDF()));
1174     GetESDsData(kESDnClustersPerTrack)->Fill(track->GetNClusters());
1175     
1176     // loop over clusters
1177     AliMUONTrackParam* trackParam = static_cast<AliMUONTrackParam*>(track->GetTrackParamAtCluster()->First());
1178     while (trackParam) {
1179       
1180       AliMUONVCluster* cluster = trackParam->GetClusterPtr();
1181       Int_t chId = cluster->GetChamberId();
1182       Int_t deID = cluster->GetDetElemId();
1183       Double_t residualX = cluster->GetX() - trackParam->GetNonBendingCoor();
1184       Double_t residualY = cluster->GetY() - trackParam->GetBendingCoor();
1185       Double_t localChi2X = (cluster->GetErrX2() > 0.) ? residualX*residualX/cluster->GetErrX2() : 0.;
1186       Double_t localChi2Y = (cluster->GetErrY2() > 0.) ? residualY*residualY/cluster->GetErrY2() : 0.;
1187       
1188       GetESDsData(kESDClusterHitMap+chId)->Fill(cluster->GetX(), cluster->GetY());
1189       
1190       GetESDsData(kESDnTotClustersPerCh)->Fill(chId);
1191       GetESDsData(kESDnTotClustersPerDE)->Fill(deID);
1192       
1193       GetESDsData(kESDClusterChargeInCh+chId)->Fill(cluster->GetCharge());
1194       GetESDsData(kESDSumClusterChargePerDE)->Fill(deID, cluster->GetCharge());
1195       
1196       if (cluster->GetNDigits() > 0) { // discard clusters with pad not stored in ESD
1197         GetESDsData(kESDnTotFullClustersPerDE)->Fill(deID);
1198         GetESDsData(kESDClusterSizeInCh+chId)->Fill(cluster->GetNDigits());
1199         GetESDsData(kESDSumClusterSizePerDE)->Fill(deID, cluster->GetNDigits());
1200       }
1201       
1202       GetESDsData(kESDResidualXInCh+chId)->Fill(residualX);
1203       GetESDsData(kESDResidualYInCh+chId)->Fill(residualY);
1204       GetESDsData(kESDSumResidualXPerDE)->Fill(deID, residualX);
1205       GetESDsData(kESDSumResidualYPerDE)->Fill(deID, residualY);
1206       GetESDsData(kESDSumResidualX2PerDE)->Fill(deID, residualX*residualX);
1207       GetESDsData(kESDSumResidualY2PerDE)->Fill(deID, residualY*residualY);
1208       
1209       GetESDsData(kESDLocalChi2XInCh+chId)->Fill(localChi2X);
1210       GetESDsData(kESDLocalChi2YInCh+chId)->Fill(localChi2Y);
1211       GetESDsData(kESDSumLocalChi2XPerDE)->Fill(deID, localChi2X);
1212       GetESDsData(kESDSumLocalChi2YPerDE)->Fill(deID, localChi2Y);
1213       
1214       trackParam = static_cast<AliMUONTrackParam*>(track->GetTrackParamAtCluster()->After(trackParam));
1215     }
1216     
1217   }
1218   
1219   GetESDsData(kESDMatchTrig)->Fill(nTrackMatchTrig);
1220   
1221 }
1222
1223 //____________________________________________________________________________ 
1224 void AliMUONQADataMakerRec::StartOfDetectorCycle()
1225 {
1226     /// Detector specific actions at start of cycle
1227   
1228 }
1229
1230 //____________________________________________________________________________ 
1231 void AliMUONQADataMakerRec::DisplayTriggerInfo(AliQA::TASKINDEX_t task)
1232 {
1233   //
1234   /// Display trigger information in a user-friendly way:
1235   /// from local board and strip numbers to their position on chambers
1236   //
1237
1238   if(task!=AliQA::kRECPOINTS && task!=AliQA::kRAWS) return;
1239
1240   AliMUONTriggerDisplay triggerDisplay;
1241   
1242   TH3F* histoStrips=0x0;
1243   TH2F* histoDisplayStrips=0x0;
1244
1245   for (Int_t iCath = 0; iCath < AliMpConstants::NofCathodes(); iCath++)
1246   {
1247     if(task==AliQA::kRECPOINTS){
1248       ERecPoints hindex 
1249         = ( iCath == 0 ) ? kTriggerDigitsBendPlane : kTriggerDigitsNonBendPlane;
1250       histoStrips = (TH3F*)GetRecPointsData(hindex);
1251     }
1252     else if(task==AliQA::kRAWS){
1253       ERaw hindex 
1254         = ( iCath == 0 ) ? kTriggerScalersBP : kTriggerScalersNBP;
1255       histoStrips = (TH3F*)GetRawsData(hindex);
1256       if(histoStrips->GetEntries()==0) continue; // No scalers found
1257     }
1258     
1259     for (Int_t iChamber = 0; iChamber < AliMpConstants::NofTriggerChambers(); iChamber++)
1260     {
1261       if(task==AliQA::kRECPOINTS){
1262         histoDisplayStrips = (TH2F*)GetRecPointsData(kTriggerDigitsDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber);
1263       }
1264       else if(task==AliQA::kRAWS){
1265         histoDisplayStrips = (TH2F*)GetRawsData(kTriggerScalersDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber);
1266       }
1267       Int_t bin = histoStrips->GetXaxis()->FindBin(11+iChamber);
1268       histoStrips->GetXaxis()->SetRange(bin,bin);
1269       TH2F* inputHisto = (TH2F*)histoStrips->Project3D("zy");
1270       triggerDisplay.FillDisplayHistogram(inputHisto, histoDisplayStrips, AliMUONTriggerDisplay::kDisplayStrips, iCath, iChamber);
1271       if(task==AliQA::kRAWS) {
1272         Float_t acquisitionTime = ((TH1F*)GetRawsData(kTriggerScalersTime))->GetBinContent(1);
1273         histoDisplayStrips->Scale(1./acquisitionTime);
1274       }
1275     } // iChamber
1276   } // iCath
1277
1278   if(task==AliQA::kRAWS){    
1279     TH2F* histoI  = (TH2F*) GetRawsData(kTriggerRPCi);
1280     TH2F* histoHV = (TH2F*) GetRawsData(kTriggerRPChv);
1281     FillTriggerDCSHistos();
1282     for (Int_t iChamber = 0; iChamber < AliMpConstants::NofTriggerChambers(); iChamber++) {
1283       Int_t bin = histoI->GetXaxis()->FindBin(11+iChamber);
1284       TH2F* histoDisplayI = (TH2F*)GetRawsData(kTriggerIDisplay + iChamber);
1285       triggerDisplay.FillDisplayHistogram(histoI->ProjectionY("_px", bin, bin), histoDisplayI, AliMUONTriggerDisplay::kDisplaySlats, 0, iChamber);
1286       TH2F* histoDisplayHV = (TH2F*)GetRawsData(kTriggerHVDisplay + iChamber);
1287       bin = histoHV->GetXaxis()->FindBin(11+iChamber);
1288       triggerDisplay.FillDisplayHistogram(histoHV->ProjectionY("_px", bin, bin), histoDisplayHV, AliMUONTriggerDisplay::kDisplaySlats, 0, iChamber);
1289     }
1290   }
1291
1292   if(task==AliQA::kRECPOINTS){
1293     TH1F* histoBoards = (TH1F*)GetRecPointsData(kTriggeredBoards);
1294     TH2F* histoDisplayBoards = (TH2F*)GetRecPointsData(kTriggerBoardsDisplay);
1295     triggerDisplay.FillDisplayHistogram(histoBoards, histoDisplayBoards, AliMUONTriggerDisplay::kDisplayBoards, 0, 0);
1296   }
1297 }
1298
1299
1300 //_____________________________________________________________________________
1301 Bool_t 
1302 AliMUONQADataMakerRec::FillTriggerDCSHistos()
1303 {
1304   /// Get HV and currents values for one trigger chamber
1305   
1306   AliCodeTimerAuto("");
1307
1308   AliMUONCalibrationData calibrationData(AliCDBManager::Instance()->GetRun());
1309
1310   TMap* triggerDcsMap = calibrationData.TriggerDCS();
1311
1312   if ( !triggerDcsMap ) 
1313   {
1314     AliError("Cannot fill DCS histos, as triggerDcsMap is NULL");
1315     return kFALSE;
1316   }
1317
1318   AliMpDEIterator deIt;
1319   
1320   deIt.First();
1321     
1322   AliMpDCSNamer triggerDcsNamer("TRIGGER");
1323
1324   TH2* currHisto = 0x0;
1325
1326   Bool_t error = kFALSE;
1327   
1328   while ( !deIt.IsDone() )
1329   {
1330     Int_t detElemId = deIt.CurrentDEId();
1331     
1332     if ( AliMpDEManager::GetStationType(detElemId) == AliMp::kStationTrigger) {
1333
1334       Int_t iChamber = AliMpDEManager::GetChamberId(detElemId);
1335       Int_t slat = detElemId%100;
1336
1337       for(Int_t iMeas=0; iMeas<AliMpDCSNamer::kNDCSMeas; iMeas++){
1338         TString currAlias = triggerDcsNamer.DCSChannelName(detElemId, 0, iMeas);
1339
1340         AliDebug(2, Form("\nDetElemId %i   dcsAlias %s", detElemId, currAlias.Data()));
1341
1342         TPair* triggerDcsPair = static_cast<TPair*>(triggerDcsMap->FindObject(currAlias.Data()));
1343
1344         if (!triggerDcsPair)
1345         {
1346           AliError(Form("Did not find expected alias (%s) for DE %d",
1347                         currAlias.Data(),detElemId));  
1348           error = kTRUE;
1349         }
1350         else
1351         {
1352           TObjArray* values = static_cast<TObjArray*>(triggerDcsPair->Value());
1353           if (!values)
1354           {
1355             AliError(Form("Could not get values for alias %s",currAlias.Data()));
1356             error = kTRUE;
1357           }
1358           else
1359           {
1360             TIter next(values);
1361             AliDCSValue* val;
1362
1363             while ( ( val = static_cast<AliDCSValue*>(next()) ) )
1364             {
1365               Float_t hvi = val->GetFloat();
1366
1367               AliDebug(2, Form("Value %f", hvi));
1368
1369               switch(iMeas){
1370               case AliMpDCSNamer::kDCSI:
1371                 currHisto = (TH2F*) GetRawsData(kTriggerRPCi);
1372                 break;
1373               case AliMpDCSNamer::kDCSHV:
1374                 currHisto = (TH2F*) GetRawsData(kTriggerRPChv);
1375                 break;
1376               } 
1377               Int_t binX = currHisto->GetXaxis()->FindBin(iChamber+1);
1378               Int_t binY = currHisto->GetYaxis()->FindBin(slat);
1379               currHisto->SetBinContent(binX, binY, hvi);
1380             } // loop on values
1381           } // if (!values)
1382         } // if (!triggerDcsPair)
1383       } // loop on measured types (HV and currents)
1384     } // if (stationType == kStationTrigger)
1385
1386     deIt.Next();
1387   }
1388   return error;
1389 }
1390
1391 //____________________________________________________________________________ 
1392 AliMUONVTrackerData* AliMUONQADataMakerRec::GetTrackerData() const
1393
1394   
1395   return fTrackerDataMaker->Data(); 
1396   
1397 }