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