]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONQADataMakerRec.cxx
Introducing event specie in QA (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 "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 "AliMpBusPatch.h"
44 #include "AliMpCDB.h"
45 #include "AliMpConstants.h"
46 #include "AliMpDDLStore.h"
47 #include "AliMpDEIterator.h"
48 #include "AliMpDEManager.h"
49 #include "AliMpLocalBoard.h"
50 #include "AliMpStationType.h"
51 #include "AliMpTriggerCrate.h"
52 #include "AliRawEventHeaderBase.h"
53
54 // --- AliRoot header files ---
55 #include "AliCDBManager.h"
56 #include "AliCDBStorage.h"
57 #include "AliESDEvent.h"
58 #include "AliESDMuonTrack.h"
59 #include "AliESDMuonCluster.h"
60 #include "AliLog.h"
61 #include "AliRawReader.h"
62 #include "AliQAChecker.h"
63 #include "AliCodeTimer.h"
64
65 // --- ROOT system ---
66 #include <TClonesArray.h>
67 #include <TFile.h> 
68 #include <TH1F.h> 
69 #include <TH1I.h> 
70 #include <TH2F.h>
71 #include <TH3F.h> 
72 #include <Riostream.h>
73
74 //-----------------------------------------------------------------------------
75 /// \class AliMUONQADataMakerRec
76 ///
77 /// MUON base class for quality assurance data (histo) maker
78 ///
79 /// \author C. Finck, D. Stocco, L. Aphecetche
80
81 /// \cond CLASSIMP
82 ClassImp(AliMUONQADataMakerRec)
83 /// \endcond
84            
85 //____________________________________________________________________________ 
86 AliMUONQADataMakerRec::AliMUONQADataMakerRec() : 
87 AliQADataMakerRec(AliQA::GetDetName(AliQA::kMUON), "MUON Quality Assurance Data Maker"),
88 fIsInitRaws(kFALSE),
89 fIsInitRecPointsTracker(kFALSE),
90 fIsInitRecPointsTrigger(kFALSE),
91 fIsInitESDs(kFALSE),
92 fDigitStore(0x0),
93 fTriggerStore(0x0),
94 fDigitMaker(0x0),
95 fClusterStore(0x0),
96 fTrackerDataMaker(0x0)
97 {
98     /// ctor
99         
100   AliDebug(1,"");
101
102         Ctor();
103 }
104
105 //____________________________________________________________________________ 
106 void
107 AliMUONQADataMakerRec::Ctor()
108 {
109         /// Init some members
110         fDigitStore = AliMUONVDigitStore::Create("AliMUONDigitStoreV1");
111         fDigitMaker = new AliMUONDigitMaker(kTRUE);
112 }
113
114 //____________________________________________________________________________ 
115 AliMUONQADataMakerRec::AliMUONQADataMakerRec(const AliMUONQADataMakerRec& qadm) :
116 AliQADataMakerRec(qadm),
117 fIsInitRaws(kFALSE),
118 fIsInitRecPointsTracker(kFALSE),
119 fIsInitRecPointsTrigger(kFALSE),
120 fIsInitESDs(kFALSE),
121 fDigitStore(0x0),
122 fTriggerStore(0x0),
123 fDigitMaker(0x0),
124 fClusterStore(0x0),
125 fTrackerDataMaker(0x0)
126 {
127     ///copy ctor 
128
129     AliDebug(1,"");
130
131
132     SetName((const char*)qadm.GetName()) ; 
133     SetTitle((const char*)qadm.GetTitle()); 
134
135         // Do not copy the digit store and digit maker, but create its own ones
136         
137         Ctor();
138         
139 }
140
141 //__________________________________________________________________
142 AliMUONQADataMakerRec& AliMUONQADataMakerRec::operator = (const AliMUONQADataMakerRec& qadm )
143 {
144   /// Assignment operator
145
146   AliDebug(1,"");
147
148   // check assignment to self
149   if (this == &qadm) return *this;
150
151   this->~AliMUONQADataMakerRec();
152   new(this) AliMUONQADataMakerRec(qadm);
153   return *this;
154 }
155
156 //__________________________________________________________________
157 AliMUONQADataMakerRec::~AliMUONQADataMakerRec()
158 {
159     /// dtor
160   
161   AliDebug(1,"");
162
163   AliCodeTimerAuto("");
164   
165   delete fDigitStore;
166   delete fTriggerStore;
167   delete fDigitMaker;
168   delete fClusterStore;
169   delete fTrackerDataMaker;
170 }
171
172 //____________________________________________________________________________ 
173 void AliMUONQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray** list)
174 {
175   ///Detector specific actions at end of cycle
176   
177   AliCodeTimerAuto("");
178   
179   // Display trigger histos in a more user friendly way
180   DisplayTriggerInfo(task);
181   
182   for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
183     SetEventSpecie(specie) ; 
184     if ( task == AliQA::kRAWS && fTrackerDataMaker ) 
185       {
186         TIter next(list[specie]);
187         TObject* o;
188         Bool_t alreadyThere(kFALSE);
189         while ( ( o = next() ) && !alreadyThere )
190           {
191             TString classname(o->ClassName());
192             if ( classname.Contains("TrackerData") ) alreadyThere = kTRUE;
193           }
194         if (!alreadyThere && fTrackerDataMaker) 
195           {
196             AliInfo("Adding fTrackerDataMaker to the list of qa objects");
197             list[specie]->AddAt(fTrackerDataMaker->Data(),(Int_t)kTrackerData);
198           }
199           if ( fTrackerDataMaker ) 
200             {
201               TH1* hbp = GetRawsData(kTrackerBusPatchOccupancy);
202               hbp->Reset();
203               TIter nextBP(AliMpDDLStore::Instance()->CreateBusPatchIterator());
204               AliMpBusPatch* bp(0x0);
205               AliMUONVTrackerData* data = fTrackerDataMaker->Data();
206               Int_t occDim = 2;
207       
208               while ( ( bp = static_cast<AliMpBusPatch*>(nextBP())) )
209                 {
210                   Int_t busPatchId = bp->GetId();
211                   Int_t bin = hbp->FindBin(busPatchId);
212                   hbp->SetBinContent(bin,data->BusPatch(busPatchId,occDim));
213                 }
214             }
215       }
216   
217     if ( task == AliQA::kESDS ) {
218       // Normalize ESD histos
219       TH1* h;
220       Int_t bin;
221       AliMpDEIterator it;
222       it.First();
223       while ( !it.IsDone()) {
224     
225         Int_t detElemId = it.CurrentDEId();
226     
227         if ( AliMpDEManager::GetStationType(detElemId) != AliMp::kStationTrigger ) {
228       
229           h = GetESDsData(kESDnClustersPerDE);
230           Double_t nClusters = h->GetBinContent(h->GetXaxis()->FindFixBin((Double_t)detElemId));
231       
232         if (nClusters > 0) {
233         
234           h = GetESDsData(kESDClusterChargePerDE);
235           bin = h->GetXaxis()->FindFixBin((Double_t)detElemId);
236           h->SetBinContent(bin, h->GetBinContent(bin)/nClusters);
237         
238           h = GetESDsData(kESDClusterMultPerDE);
239           bin = h->GetXaxis()->FindFixBin((Double_t)detElemId);
240           h->SetBinContent(bin, h->GetBinContent(bin)/nClusters);
241         
242           h = GetESDsData(kESDResidualXPerDEMean);
243           bin = h->GetXaxis()->FindFixBin((Double_t)detElemId);
244           Double_t meanResX = h->GetBinContent(bin)/nClusters;
245           h->SetBinContent(bin, meanResX);
246         
247           h = GetESDsData(kESDResidualYPerDEMean);
248           bin = h->GetXaxis()->FindFixBin((Double_t)detElemId);
249           Double_t meanResY = h->GetBinContent(bin)/nClusters;
250           h->SetBinContent(bin, meanResY);
251         
252           h = GetESDsData(kESDResidualXPerDESigma);
253           bin = h->GetXaxis()->FindFixBin((Double_t)detElemId);
254           if (nClusters > 1) h->SetBinContent(bin, TMath::Sqrt(h->GetBinContent(bin)/nClusters - meanResX*meanResX));
255           else h->SetBinContent(bin, 0.);
256         
257           h = GetESDsData(kESDResidualYPerDESigma);
258           bin = h->GetXaxis()->FindFixBin((Double_t)detElemId);
259           if (nClusters > 1) h->SetBinContent(bin, TMath::Sqrt(h->GetBinContent(bin)/nClusters - meanResY*meanResY));
260           else h->SetBinContent(bin, 0.);
261         }
262       
263       }
264     
265       it.Next();
266     }
267   
268     Double_t nTracks = GetESDsData(kESDnClustersPerTrack)->GetEntries();
269     if (nTracks > 0) {
270       GetESDsData(kESDnClustersPerCh)->Scale(1./nTracks);
271       GetESDsData(kESDnClustersPerDE)->Scale(1./nTracks);
272     }
273   }
274   }
275   // do the QA checking
276   AliQAChecker::Instance()->Run(AliQA::kMUON, task, list) ;
277 }
278
279 //____________________________________________________________________________ 
280 void AliMUONQADataMakerRec::InitRaws()
281 {
282     /// create Raws histograms in Raws subdir
283         
284   AliCodeTimerAuto("");
285   
286   Bool_t forExpert(kTRUE);
287   
288         if ( ! AliCDBManager::Instance()->GetDefaultStorage() )
289         {
290                 AliError("CDB default storage not set. Cannot work.");
291                 fIsInitRaws=kFALSE;
292         }
293         
294         TH3F* h3 = new TH3F("hTriggerScalersBendPlane", "Trigger scalers in bending plane",
295                                                                                         4, 10.5, 14.5,
296                                                                                         234, 0.5, 234.5,
297                                                                                         16, -0.5, 15.5);
298         h3->GetXaxis()->SetTitle("Chamber");
299         h3->GetYaxis()->SetTitle("Board");
300         h3->GetZaxis()->SetTitle("Strip");
301         Add2RawsList(h3, kTriggerScalersBP,forExpert);
302         
303         TH3F* h4 = new TH3F("hTriggerScalersNonBendPlane", "Trigger scalers in non-bending plane",
304                                                                                         4, 10.5, 14.5,
305                                                                                         234, 0.5, 234.5,
306                                                                                         16, -0.5, 15.5);
307         h4->GetXaxis()->SetTitle("Chamber");
308         h4->GetYaxis()->SetTitle("Board");
309         h4->GetZaxis()->SetTitle("Strip");
310         Add2RawsList(h4, kTriggerScalersNBP,forExpert);
311         
312         AliMUONTriggerDisplay triggerDisplay;
313         TString histoName, histoTitle;
314         for(Int_t iCath=0; iCath<AliMpConstants::NofCathodes(); iCath++){
315                 TString cathName = ( iCath==0 ) ? "BendPlane" : "NonBendPlane";
316                 for(Int_t iChamber=0; iChamber<AliMpConstants::NofTriggerChambers(); iChamber++){
317                         histoName = Form("hScalers%sChamber%i", cathName.Data(), 11+iChamber);
318                         histoTitle = Form("Chamber %i: Scalers %s", 11+iChamber, cathName.Data());
319                         TH2F* h5 = (TH2F*)triggerDisplay.GetEmptyDisplayHisto(histoName, AliMUONTriggerDisplay::kDisplayStrips, 
320                                                                               iCath, iChamber, histoTitle);
321                         Add2RawsList(h5, kTriggerScalersDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber,forExpert);
322                 }
323         }
324         
325   Int_t nbp(0);
326   TIter next(AliMpDDLStore::Instance()->CreateBusPatchIterator());
327   while (next())
328   {
329     ++nbp;
330   }
331   
332   TH1* hbp = new TH1F("hTrackerBusPatchOccupancy","Occupancy of bus patches",
333                       nbp,-0.5,nbp-0.5);
334   
335   Add2RawsList(hbp,kTrackerBusPatchOccupancy,!forExpert);
336   
337         fIsInitRaws = kTRUE;
338 }
339
340 //____________________________________________________________________________ 
341 void AliMUONQADataMakerRec::InitRecPoints()
342 {
343         /// create Reconstructed Points histograms in RecPoints subdir
344   
345   AliCodeTimerAuto("");
346   
347         InitRecPointsTrigger();
348         InitRecPointsTracker();
349 }
350
351 //____________________________________________________________________________ 
352 void AliMUONQADataMakerRec::InitRecPointsTracker()
353 {
354         /// create Reconstructed Points histograms in RecPoints subdir for the
355         /// MUON tracker subsystem.
356
357   AliCodeTimerAuto("");
358   
359   Bool_t forExpert(kTRUE);
360   
361         AliMpDEIterator it;
362         
363         it.First();
364         
365         Int_t ndes(0);
366         
367         while ( !it.IsDone())
368         {
369                 Int_t detElemId = it.CurrentDEId();
370                 
371                 it.Next();
372
373                 if ( AliMpDEManager::GetStationType(detElemId) != AliMp::kStationTrigger )
374                 {
375                         ndes = TMath::Max(ndes,detElemId);
376
377                         TH1* h = new TH1I(Form("hTrackerClusterMultiplicityForDE%04d",detElemId),
378                                           Form("Multiplicity of the clusters in detection element %d",detElemId),
379                                           100,0,100);
380                         
381                         h->GetXaxis()->SetTitle("Detection Element Id");
382                         
383                         Add2RecPointsList(h,kTrackerClusterMultiplicityPerDE+detElemId,forExpert);
384                         
385                         h =  new TH1I(Form("hTrackerClusterChargeForDE%04d",detElemId),
386                                       Form("Charge of the clusters in detection element %d",detElemId),
387                                       100,0,1000);
388
389                         h->GetXaxis()->SetTitle("Detection Element Id");
390
391                         Add2RecPointsList(h,kTrackerClusterChargePerDE+detElemId,forExpert);
392
393                 }
394
395         }
396
397         TH1* h = new TH1I("hTrackerNumberOfClustersPerDE","Number of clusters per detection element",
398                                                                                 ndes, 0.5, ndes + 0.5);
399
400         h->GetXaxis()->SetTitle("Detection Element Id");
401
402         Add2RecPointsList(h, kTrackerNumberOfClustersPerDE,!forExpert);
403
404         for ( Int_t i = 0; i < AliMpConstants::NofTrackingChambers(); ++i ) 
405         {
406                 TH1* h1 = new TH1I("hTrackerNumberOfClustersPerChamber","Number of clusters per chamber",
407                                    AliMpConstants::NofTrackingChambers(),-0.5,AliMpConstants::NofTrackingChambers()-0.5);
408                 Add2RecPointsList(h1,kTrackerNumberOfClustersPerChamber,forExpert);
409                 h1 = new TH1I(Form("hTrackerClusterMultiplicityForChamber%d",i),
410                               Form("Cluster multiplicity for chamber %d",i),
411                               100,0,100);
412                 Add2RecPointsList(h1,kTrackerClusterMultiplicityPerChamber+i,forExpert);
413                 h1 = new TH1I(Form("hTrackerClusterChargeForChamber%d",i),
414                               Form("Cluster charge for chamber %d",i),
415                               100,0,1000);
416                 Add2RecPointsList(h1,kTrackerClusterChargePerChamber+i,forExpert);
417         }
418         
419         fIsInitRecPointsTracker=kTRUE;
420 }
421
422 //____________________________________________________________________________ 
423 void AliMUONQADataMakerRec::InitRecPointsTrigger()
424 {
425         /// create Reconstructed Points histograms in RecPoints subdir for the
426         /// MUON Trigger subsystem.
427         
428   Bool_t forExpert(kTRUE);
429   
430     TH3F* h0 = new TH3F("hTriggerDigitsBendPlane", "Trigger digits in bending plane",
431                         4, 10.5, 14.5,
432                         234, 0.5, 234.5,
433                         16, -0.5, 15.5);
434     h0->GetXaxis()->SetTitle("Chamber");
435     h0->GetYaxis()->SetTitle("Board");
436     h0->GetZaxis()->SetTitle("Strip");
437     Add2RecPointsList(h0, kTriggerDigitsBendPlane,forExpert);
438
439     TH3F* h1 = new TH3F("hTriggerDigitsNonBendPlane", "Trigger digits in non-bending plane",
440                         4, 10.5, 14.5,
441                         234, 0.5, 234.5,
442                         16, -0.5, 15.5);
443     h1->GetXaxis()->SetTitle("Chamber");
444     h1->GetYaxis()->SetTitle("Board");
445     h1->GetZaxis()->SetTitle("Strip");
446     Add2RecPointsList(h1, kTriggerDigitsNonBendPlane,forExpert);
447
448     TH1F* h2 = new TH1F("hTriggeredBoards", "Triggered boards", 234, 0.5, 234.5);
449     Add2RecPointsList(h2, kTriggeredBoards,forExpert);
450
451     AliMUONTriggerDisplay triggerDisplay;
452     TString histoName, histoTitle;
453     for(Int_t iCath=0; iCath<AliMpConstants::NofCathodes(); iCath++){
454       TString cathName = ( iCath==0 ) ? "BendPlane" : "NonBendPlane";
455       for(Int_t iChamber=0; iChamber<AliMpConstants::NofTriggerChambers(); iChamber++){
456         histoName = Form("hTriggerDigits%sChamber%i", cathName.Data(), 11+iChamber);
457         histoTitle = Form("Chamber %i: Fired pads %s", 11+iChamber, cathName.Data());
458         TH2F* h3 = (TH2F*)triggerDisplay.GetEmptyDisplayHisto(histoName, AliMUONTriggerDisplay::kDisplayStrips, 
459                                                               iCath, iChamber, histoTitle);
460         Add2RecPointsList(h3, kTriggerDigitsDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber,forExpert);
461       }
462     }
463
464     TH2F* h4 = (TH2F*)triggerDisplay.GetEmptyDisplayHisto("hFiredBoardsDisplay", AliMUONTriggerDisplay::kDisplayBoards,
465                                                           0, 0, "Fired boards");
466     Add2RecPointsList(h4, kTriggerBoardsDisplay,forExpert);
467         
468         fIsInitRecPointsTrigger = kTRUE;
469 }
470
471
472 //____________________________________________________________________________ 
473 void AliMUONQADataMakerRec::InitESDs()
474 {
475   ///create ESDs histograms in ESDs subdir
476   
477   Bool_t forExpert(kTRUE);
478   
479   Int_t nCh = AliMUONConstants::NTrackingCh();
480   Int_t nDE = 1100;
481   
482   // track info
483   TH1F* hESDnTracks = new TH1F("hESDnTracks", "number of tracks", 20, 0., 20.);
484   Add2ESDsList(hESDnTracks, kESDnTracks,!forExpert);
485
486   TH1F* hESDMatchTrig = new TH1F("hESDMatchTrig", "number of tracks matched with trigger", 20, 0., 20.);
487   Add2ESDsList(hESDMatchTrig, kESDMatchTrig,!forExpert);
488   
489   TH1F* hESDMomentum = new TH1F("hESDMomentum", "P distribution", 300, 0., 300);
490   Add2ESDsList(hESDMomentum, kESDMomentum,forExpert);
491
492   TH1F* hESDPt = new TH1F("hESDPt", "Pt distribution", 200, 0., 50);
493   Add2ESDsList(hESDPt, kESDPt,forExpert);
494
495   TH1F* hESDRapidity = new TH1F("hESDRapidity", "rapidity distribution", 200, -4.5, -2.);
496   Add2ESDsList(hESDRapidity, kESDRapidity,forExpert);
497
498   TH1F* hESDChi2 = new TH1F("hESDChi2", "normalized chi2 distribution", 500, 0., 50.);
499   Add2ESDsList(hESDChi2, kESDChi2,forExpert);
500   
501   // cluster info
502   for (Int_t i = 0; i < nCh; i++) {
503     Float_t rMax = AliMUONConstants::Rmax(i/2);
504     TH2F* hESDClusterHitMap = new TH2F(Form("hESDClusterHitMap%d",i+1), Form("cluster position distribution in chamber %d",i+1),
505                                        100, -rMax, rMax, 100, -rMax, rMax);
506     Add2ESDsList(hESDClusterHitMap, kESDClusterHitMap+i,forExpert);
507   }
508   
509   TH1F* hESDnClustersPerTrack = new TH1F("hESDnClustersPerTrack", "number of clusters per track", 20, 0., 20.);
510   Add2ESDsList(hESDnClustersPerTrack, kESDnClustersPerTrack,!forExpert);
511   
512   TH1F* hESDnClustersPerCh = new TH1F("hESDnClustersPerCh", "number of clusters per chamber per track;chamber ID", nCh, 0, nCh);
513   hESDnClustersPerCh->SetFillColor(kRed);
514   Add2ESDsList(hESDnClustersPerCh, kESDnClustersPerCh,forExpert);
515   
516   TH1F* hESDnClustersPerDE = new TH1F("hESDnClustersPerDE", "number of clusters per DE per track;DetElem ID", nDE, 0, nDE);
517   hESDnClustersPerDE->SetFillColor(kRed);
518   Add2ESDsList(hESDnClustersPerDE, kESDnClustersPerDE,forExpert);
519   
520   TH1F* hESDClusterCharge = new TH1F("hESDClusterCharge", "cluster charge distribution", 500, 0., 5000.);
521   Add2ESDsList(hESDClusterCharge, kESDClusterCharge,forExpert);
522   
523   for (Int_t i = 0; i < nCh; i++) {
524     TH1F* hESDClusterChargeInCh = new TH1F(Form("hESDClusterChargeInCh%d",i+1), Form("cluster charge distribution in chamber %d",i+1), 500, 0., 5000.);
525     Add2ESDsList(hESDClusterChargeInCh, kESDClusterChargeInCh+i,forExpert);
526   }
527   
528   TH1F* hESDClusterChargePerDE = new TH1F("hESDClusterChargePerDE", "cluster mean charge per DE;DetElem ID", nDE, 0, nDE);
529   hESDClusterChargePerDE->SetOption("P");
530   hESDClusterChargePerDE->SetMarkerStyle(kFullDotMedium);
531   hESDClusterChargePerDE->SetMarkerColor(kRed);
532   Add2ESDsList(hESDClusterChargePerDE, kESDClusterChargePerDE,forExpert);
533   
534   TH1F* hESDClusterMult = new TH1F("hESDClusterMult", "cluster multiplicity distribution", 200, 0., 200.);
535   Add2ESDsList(hESDClusterMult, kESDClusterMult,forExpert);
536   
537   for (Int_t i = 0; i < nCh; i++) {
538     TH1F* hESDClusterMultInCh = new TH1F(Form("hESDClusterMultInCh%d",i+1), Form("cluster multiplicity distribution in chamber %d",i+1), 200, 0., 200.);
539     Add2ESDsList(hESDClusterMultInCh, kESDClusterMultInCh+i,forExpert);
540   }
541   
542   TH1F* hESDClusterMultPerDE = new TH1F("hESDClusterMultPerDE", "cluster mean multiplicity per DE;DetElem ID", nDE, 0, nDE);
543   hESDClusterMultPerDE->SetOption("P");
544   hESDClusterMultPerDE->SetMarkerStyle(kFullDotMedium);
545   hESDClusterMultPerDE->SetMarkerColor(kRed);
546   Add2ESDsList(hESDClusterMultPerDE, kESDClusterMultPerDE,forExpert);
547   
548   // cluster - track info
549   TH1F* hESDResidualX = new TH1F("hESDResidualX", "cluster-track residual-X distribution", 1000, -5., 5.);
550   Add2ESDsList(hESDResidualX, kESDResidualX,forExpert);
551   
552   TH1F* hESDResidualY = new TH1F("hESDResidualY", "cluster-track residual-Y distribution", 1000, -1., 1.);
553   Add2ESDsList(hESDResidualY, kESDResidualY,forExpert);
554   
555   for (Int_t i = 0; i < nCh; i++) {
556     TH1F* hESDResidualXInCh = new TH1F(Form("hESDResidualXInCh%d",i+1), Form("cluster-track residual-X distribution in chamber %d",i+1), 1000, -5., 5.);
557     Add2ESDsList(hESDResidualXInCh, kESDResidualXInCh+i,forExpert);
558     
559     TH1F* hESDResidualYInCh = new TH1F(Form("hESDResidualYInCh%d",i+1), Form("cluster-track residual-Y distribution in chamber %d",i+1), 1000, -1., 1.);
560     Add2ESDsList(hESDResidualYInCh, kESDResidualYInCh+i,forExpert);
561   }
562   
563   TH1F* hESDResidualXPerDEMean = new TH1F("hESDResidualXPerDEMean", "cluster-track residual-X per DE: mean;DetElem ID", nDE, 0, nDE);
564   hESDResidualXPerDEMean->SetOption("P");
565   hESDResidualXPerDEMean->SetMarkerStyle(kFullDotMedium);
566   hESDResidualXPerDEMean->SetMarkerColor(kRed);
567   Add2ESDsList(hESDResidualXPerDEMean, kESDResidualXPerDEMean,forExpert);
568   
569   TH1F* hESDResidualYPerDEMean = new TH1F("hESDResidualYPerDEMean", "cluster-track residual-Y per DE: mean;DetElem ID", nDE, 0, nDE);
570   hESDResidualYPerDEMean->SetOption("P");
571   hESDResidualYPerDEMean->SetMarkerStyle(kFullDotMedium);
572   hESDResidualYPerDEMean->SetMarkerColor(kRed);
573   Add2ESDsList(hESDResidualYPerDEMean, kESDResidualYPerDEMean,forExpert);
574   
575   TH1F* hESDResidualXPerDESigma = new TH1F("hESDResidualXPerDESigma", "cluster-track residual-X per DE: sigma;DetElem ID", nDE, 0, nDE);
576   hESDResidualXPerDESigma->SetOption("P");
577   hESDResidualXPerDESigma->SetMarkerStyle(kFullDotMedium);
578   hESDResidualXPerDESigma->SetMarkerColor(kRed);
579   Add2ESDsList(hESDResidualXPerDESigma, kESDResidualXPerDESigma,forExpert);
580   
581   TH1F* hESDResidualYPerDESigma = new TH1F("hESDResidualYPerDESigma", "cluster-track residual-Y per DE: sigma;DetElem ID", nDE, 0, nDE);
582   hESDResidualYPerDESigma->SetOption("P");
583   hESDResidualYPerDESigma->SetMarkerStyle(kFullDotMedium);
584   hESDResidualYPerDESigma->SetMarkerColor(kRed);
585   Add2ESDsList(hESDResidualYPerDESigma, kESDResidualYPerDESigma,forExpert);
586   
587   fIsInitESDs =  kTRUE;
588 }
589
590 //____________________________________________________________________________
591 void AliMUONQADataMakerRec::MakeRaws(AliRawReader* rawReader)
592 {
593     /// make QA for rawdata
594
595     if ( ! fIsInitRaws ) {
596       AliWarningStream() 
597         << "Skipping function due to a failure in Init" << endl;
598       return;
599     }    
600
601   if ( rawReader->GetType() == AliRawEventHeaderBase::kPhysicsEvent ) 
602   {
603     rawReader->Reset();
604     MakeRawsTracker(rawReader);
605   }
606   
607   rawReader->Reset();    
608   MakeRawsTrigger(rawReader);
609 }
610
611 //____________________________________________________________________________
612 void AliMUONQADataMakerRec::MakeRawsTracker(AliRawReader* rawReader)
613 {
614         /// make QA for rawdata tracker
615   
616         if (!fTrackerDataMaker) 
617         {
618                 const Bool_t histogram(kFALSE);
619                 const Bool_t fastDecoder(kTRUE);
620     
621 //    fTrackerDataMaker = new AliMUONTrackerRawDataMaker(rawReader,histogram,fastDecoder,takeRawReaderOwnership);
622
623                 fTrackerDataMaker = new AliMUONTrackerCalibratedDataMaker(GetRecoParam(),
624                                                                           AliCDBManager::Instance()->GetRun(),
625                                                                           rawReader,
626                                                                           AliCDBManager::Instance()->GetDefaultStorage()->GetURI(),
627                                                                           "NOGAIN",
628                                                                           histogram,
629                                                                           0.0,0.0,
630                                                                           fastDecoder);
631                 
632                 fTrackerDataMaker->Data()->DisableChannelLevel(); // to save up disk space, we only store starting at the manu level
633                 
634                 fTrackerDataMaker->SetRunning(kTRUE);
635         }
636         
637         ((AliMUONTrackerCalibratedDataMaker*)fTrackerDataMaker)->SetRawReader(rawReader);
638         
639         fTrackerDataMaker->ProcessEvent();
640 }
641
642 //____________________________________________________________________________
643 void AliMUONQADataMakerRec::MakeRawsTrigger(AliRawReader* rawReader)
644 {
645         /// make QA for rawdata trigger
646         
647     // Get trigger scalers
648
649     Int_t loCircuit=0;
650     AliMpCDB::LoadDDLStore();
651
652     AliMUONRawStreamTrigger rawStreamTrig(rawReader);
653     while (rawStreamTrig.NextDDL()) 
654     {
655       // If not a scaler event, do nothing
656       Bool_t scalerEvent =  rawReader->GetDataHeader()->GetL1TriggerMessage() & 0x1;
657       if(!scalerEvent) break;
658
659       AliMUONDDLTrigger* ddlTrigger = rawStreamTrig.GetDDLTrigger();
660       AliMUONDarcHeader* darcHeader = ddlTrigger->GetDarcHeader();
661
662       Int_t nReg = darcHeader->GetRegHeaderEntries();
663     
664       for(Int_t iReg = 0; iReg < nReg ;iReg++)
665       {   //reg loop
666
667         // crate info  
668         AliMpTriggerCrate* crate = AliMpDDLStore::Instance()->
669           GetTriggerCrate(rawStreamTrig.GetDDL(), iReg);
670
671         AliMUONRegHeader* regHeader =  darcHeader->GetRegHeaderEntry(iReg);
672
673         // loop over local structures
674         Int_t nLocal = regHeader->GetLocalEntries();
675         for(Int_t iLocal = 0; iLocal < nLocal; iLocal++) 
676         {
677           AliMUONLocalStruct* localStruct = regHeader->GetLocalEntry(iLocal);
678         
679           // if card exist
680           if (!localStruct) continue;
681           
682           loCircuit = crate->GetLocalBoardId(localStruct->GetId());
683
684           if ( !loCircuit ) continue; // empty slot
685
686           AliMpLocalBoard* localBoard = AliMpDDLStore::Instance()->GetLocalBoard(loCircuit, false);
687
688           // skip copy cards
689           if( !localBoard->IsNotified()) 
690             continue;
691
692           Int_t cathode = localStruct->GetComptXY()%2;
693
694           ERaw hindex = (cathode==0) ? kTriggerScalersBP : kTriggerScalersNBP;
695
696           // loop over strips
697           for (Int_t ibitxy = 0; ibitxy < 16; ++ibitxy) {
698             if(localStruct->GetXY1(ibitxy) > 0)
699               ((TH3F*)GetRawsData(hindex))->Fill(11+0, loCircuit, ibitxy, 2*localStruct->GetXY1(ibitxy));
700             if(localStruct->GetXY2(ibitxy) > 0)
701               ((TH3F*)GetRawsData(hindex))->Fill(11+1, loCircuit, ibitxy, 2*localStruct->GetXY2(ibitxy));
702             if(localStruct->GetXY3(ibitxy) > 0)
703               ((TH3F*)GetRawsData(hindex))->Fill(11+2, loCircuit, ibitxy, 2*localStruct->GetXY3(ibitxy));
704             if(localStruct->GetXY4(ibitxy) > 0)
705               ((TH3F*)GetRawsData(hindex))->Fill(11+3, loCircuit, ibitxy, 2*localStruct->GetXY4(ibitxy));
706           } // loop on strips
707         } // iLocal
708       } // iReg
709     } // NextDDL
710 }
711
712 //____________________________________________________________________________
713 void AliMUONQADataMakerRec::MakeRecPoints(TTree* clustersTree)
714 {
715         /// Fill histograms from treeR
716         
717         if (fIsInitRecPointsTracker) MakeRecPointsTracker(clustersTree);
718         if (fIsInitRecPointsTrigger) MakeRecPointsTrigger(clustersTree);
719 }
720
721 //____________________________________________________________________________
722 void AliMUONQADataMakerRec::MakeRecPointsTracker(TTree* clustersTree)
723 {
724         /// Fill histograms related to tracker clusters 
725         
726         // First things first : do we have clusters in the TreeR ?
727         // In "normal" production mode, it should be perfectly normal
728         // *not* to have them.
729         // But if for some reason we de-activated the combined tracking,
730         // then we have clusters in TreeR, so let's take that opportunity
731         // to QA them...
732         
733         if (!fClusterStore)
734         {
735                 AliCodeTimerAuto("ClusterStore creation");
736                 fClusterStore = AliMUONVClusterStore::Create(*clustersTree);
737                 if (!fClusterStore) 
738                 {
739                         fIsInitRecPointsTracker = kFALSE;
740                         return;
741                 }
742         }
743         
744         AliCodeTimerAuto("");
745         
746         fClusterStore->Connect(*clustersTree,kFALSE);
747         clustersTree->GetEvent(0);
748
749         TIter next(fClusterStore->CreateIterator());
750         AliMUONVCluster* cluster;
751         
752         while ( ( cluster = static_cast<AliMUONVCluster*>(next()) ) )
753         {
754                 Int_t detElemId = cluster->GetDetElemId();
755                 Int_t chamberId = AliMpDEManager::GetChamberId(detElemId);
756                 
757                 GetRecPointsData(kTrackerNumberOfClustersPerDE)->Fill(detElemId);
758                 GetRecPointsData(kTrackerClusterChargePerDE+detElemId)->Fill(cluster->GetCharge());
759                 GetRecPointsData(kTrackerClusterMultiplicityPerDE+detElemId)->Fill(cluster->GetNDigits());
760
761                 GetRecPointsData(kTrackerNumberOfClustersPerChamber)->Fill(chamberId);
762                 GetRecPointsData(kTrackerClusterChargePerChamber+chamberId)->Fill(cluster->GetCharge());
763                 GetRecPointsData(kTrackerClusterMultiplicityPerChamber+chamberId)->Fill(cluster->GetNDigits());
764                 
765         }
766         
767         fClusterStore->Clear();
768 }
769
770 //____________________________________________________________________________
771 void AliMUONQADataMakerRec::MakeRecPointsTrigger(TTree* clustersTree)
772 {
773         /// makes data from trigger response
774       
775     // Fired pads info
776     fDigitStore->Clear();
777
778     if (!fTriggerStore) fTriggerStore = AliMUONVTriggerStore::Create(*clustersTree);
779     fTriggerStore->Clear();
780     fTriggerStore->Connect(*clustersTree, false);
781     clustersTree->GetEvent(0);
782
783     AliMUONLocalTrigger* locTrg;
784     TIter nextLoc(fTriggerStore->CreateLocalIterator());
785
786     while ( ( locTrg = static_cast<AliMUONLocalTrigger*>(nextLoc()) ) ) 
787     {
788       if (locTrg->IsNull()) continue;
789    
790       TArrayS xyPattern[2];
791       locTrg->GetXPattern(xyPattern[0]);
792       locTrg->GetYPattern(xyPattern[1]);
793
794       Int_t nBoard = locTrg->LoCircuit();
795
796       Bool_t xTrig=locTrg->IsTrigX();
797       Bool_t yTrig=locTrg->IsTrigY();
798     
799       if (xTrig && yTrig)
800         ((TH1F*)GetRecPointsData(kTriggeredBoards))->Fill(nBoard);
801
802       fDigitMaker->TriggerDigits(nBoard, xyPattern, *fDigitStore);
803     }
804
805     TIter nextDigit(fDigitStore->CreateIterator());
806     AliMUONVDigit* mDigit;
807     while ( ( mDigit = static_cast<AliMUONVDigit*>(nextDigit()) ) )
808     {
809       Int_t detElemId = mDigit->DetElemId();
810       Int_t ch = detElemId/100;
811       Int_t localBoard = mDigit->ManuId();
812       Int_t channel = mDigit->ManuChannel();
813       Int_t cathode = mDigit->Cathode();
814       ERecPoints hindex 
815         = ( cathode == 0 ) ? kTriggerDigitsBendPlane : kTriggerDigitsNonBendPlane;
816       
817       ((TH3F*)GetRecPointsData(hindex))->Fill(ch, localBoard, channel);
818     }
819 }
820
821 //____________________________________________________________________________
822 void AliMUONQADataMakerRec::MakeESDs(AliESDEvent* esd)
823 {
824   /// make QA data from ESDs
825   
826   if ( ! fIsInitESDs ) {
827     AliWarningStream() 
828     << "Skipping function due to a failure in Init" << endl;
829     return;
830   }    
831   
832   // load ESD event in the interface
833   AliMUONESDInterface esdInterface;
834   if (GetRecoParam()) AliMUONESDInterface::ResetTracker(GetRecoParam());
835   else AliError("Unable to get recoParam: use default ones for residual calculation");
836   esdInterface.LoadEvent(*esd);
837   
838   GetESDsData(kESDnTracks)->Fill(esdInterface.GetNTracks());
839   
840   Int_t nTrackMatchTrig = 0;
841   
842   // loop over tracks
843   Int_t nTracks = (Int_t) esd->GetNumberOfMuonTracks(); 
844   for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack) {
845     
846     // get the ESD track and skip "ghosts"
847     AliESDMuonTrack* esdTrack = esd->GetMuonTrack(iTrack);
848     if (!esdTrack->ContainTrackerData()) continue;
849     
850     // get corresponding MUON track
851     AliMUONTrack* track = esdInterface.FindTrack(esdTrack->GetUniqueID());
852     
853     if (esdTrack->ContainTriggerData()) nTrackMatchTrig++;
854     
855     GetESDsData(kESDMomentum)->Fill(esdTrack->P());
856     GetESDsData(kESDPt)->Fill(esdTrack->Pt());
857     GetESDsData(kESDRapidity)->Fill(esdTrack->Y());
858     GetESDsData(kESDChi2)->Fill(track->GetNormalizedChi2());
859     GetESDsData(kESDnClustersPerTrack)->Fill(track->GetNClusters());
860     
861     // loop over clusters
862     AliMUONTrackParam* trackParam = static_cast<AliMUONTrackParam*>(track->GetTrackParamAtCluster()->First());
863     while (trackParam) {
864       
865       AliMUONVCluster* cluster = trackParam->GetClusterPtr();
866       Int_t chId = cluster->GetChamberId();
867       Int_t deID = cluster->GetDetElemId();
868       Double_t residualX = cluster->GetX() - trackParam->GetNonBendingCoor();
869       Double_t residualY = cluster->GetY() - trackParam->GetBendingCoor();
870       
871       GetESDsData(kESDClusterHitMap+chId)->Fill(cluster->GetX(), cluster->GetY());
872       
873       GetESDsData(kESDnClustersPerCh)->Fill(chId);
874       GetESDsData(kESDnClustersPerDE)->Fill(deID);
875       
876       GetESDsData(kESDClusterCharge)->Fill(cluster->GetCharge());
877       GetESDsData(kESDClusterChargeInCh+chId)->Fill(cluster->GetCharge());
878       GetESDsData(kESDClusterChargePerDE)->Fill(deID, cluster->GetCharge());
879       
880       if (cluster->GetNDigits() > 0) { // discard clusters with pad not stored in ESD
881         GetESDsData(kESDClusterMult)->Fill(cluster->GetNDigits());
882         GetESDsData(kESDClusterMultInCh+chId)->Fill(cluster->GetNDigits());
883         GetESDsData(kESDClusterMultPerDE)->Fill(deID, cluster->GetNDigits());
884       }
885       
886       GetESDsData(kESDResidualX)->Fill(residualX);
887       GetESDsData(kESDResidualY)->Fill(residualY);
888       GetESDsData(kESDResidualXInCh+chId)->Fill(residualX);
889       GetESDsData(kESDResidualYInCh+chId)->Fill(residualY);
890       GetESDsData(kESDResidualXPerDEMean)->Fill(deID, residualX);
891       GetESDsData(kESDResidualYPerDEMean)->Fill(deID, residualY);
892       GetESDsData(kESDResidualXPerDESigma)->Fill(deID, residualX*residualX);
893       GetESDsData(kESDResidualYPerDESigma)->Fill(deID, residualY*residualY);
894       
895       trackParam = static_cast<AliMUONTrackParam*>(track->GetTrackParamAtCluster()->After(trackParam));
896     }
897     
898   }
899   
900   GetESDsData(kESDMatchTrig)->Fill(nTrackMatchTrig);
901   
902 }
903
904 //____________________________________________________________________________ 
905 void AliMUONQADataMakerRec::StartOfDetectorCycle()
906 {
907     /// Detector specific actions at start of cycle
908   
909 }
910
911 //____________________________________________________________________________ 
912 void AliMUONQADataMakerRec::DisplayTriggerInfo(AliQA::TASKINDEX_t task)
913 {
914   //
915   /// Display trigger information in a user-friendly way:
916   /// from local board and strip numbers to their position on chambers
917   //
918   if(task!=AliQA::kRECPOINTS && task!=AliQA::kRAWS) return;
919
920   AliMUONTriggerDisplay triggerDisplay;
921   
922   TH3F* histoStrips=0x0;
923   TH2F* histoDisplayStrips=0x0;
924
925   for (Int_t iCath = 0; iCath < AliMpConstants::NofCathodes(); iCath++)
926   {
927     if(task==AliQA::kRECPOINTS){
928       ERecPoints hindex 
929         = ( iCath == 0 ) ? kTriggerDigitsBendPlane : kTriggerDigitsNonBendPlane;
930       histoStrips = (TH3F*)GetRecPointsData(hindex);
931     }
932     else if(task==AliQA::kRAWS){
933       ERaw hindex 
934         = ( iCath == 0 ) ? kTriggerScalersBP : kTriggerScalersNBP;
935       histoStrips = (TH3F*)GetRawsData(hindex);
936       if(histoStrips->GetEntries()==0) return; // No scalers found
937     }
938     
939     for (Int_t iChamber = 0; iChamber < AliMpConstants::NofTriggerChambers(); iChamber++)
940     {
941       if(task==AliQA::kRECPOINTS){
942         histoDisplayStrips = (TH2F*)GetRecPointsData(kTriggerDigitsDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber);
943       }
944       else if(task==AliQA::kRAWS){
945         histoDisplayStrips = (TH2F*)GetRawsData(kTriggerScalersDisplay + AliMpConstants::NofTriggerChambers()*iCath + iChamber);
946       }
947       Int_t bin = histoStrips->GetXaxis()->FindBin(11+iChamber);
948       histoStrips->GetXaxis()->SetRange(bin,bin);
949       TH2F* inputHisto = (TH2F*)histoStrips->Project3D("zy");
950       triggerDisplay.FillDisplayHistogram(inputHisto, histoDisplayStrips, AliMUONTriggerDisplay::kDisplayStrips, iCath, iChamber);
951     } // iChamber
952   } // iCath
953
954   if(task!=AliQA::kRECPOINTS) return;
955   TH1F* histoBoards = (TH1F*)GetRecPointsData(kTriggeredBoards);
956   TH2F* histoDisplayBoards = (TH2F*)GetRecPointsData(kTriggerBoardsDisplay);
957   triggerDisplay.FillDisplayHistogram(histoBoards, histoDisplayBoards, AliMUONTriggerDisplay::kDisplayBoards, 0, 0);
958 }