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