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