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