]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONQADataMakerRec.cxx
Adding check for mapping; try to load it if it is not yet available
[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
17 // --- ROOT system ---
18 #include <TClonesArray.h>
19 #include <TFile.h> 
20 #include <TH1F.h> 
21 #include <TH1I.h> 
22 #include <TH2F.h>
23 #include <TH3F.h> 
24 #include <TLorentzVector.h>
25
26 // --- AliRoot header files ---
27 #include "AliESDEvent.h"
28 #include "AliMUONConstants.h"  
29 #include "AliLog.h"
30 #include "AliRawReader.h"
31 #include "AliQAChecker.h"
32 #include "AliMpBusPatch.h"
33 #include "AliMUONCluster.h"  
34 #include "AliMUONRawStreamTracker.h"
35 #include "AliMUONRawStreamTrigger.h"
36
37 #include "AliMUONVClusterStore.h"
38 #include "AliMUONVCluster.h"
39 #include "AliESDMuonTrack.h"
40 #include "AliESDMuonCluster.h"
41
42 #include "AliMUONDigitMaker.h"
43 #include "AliMUONVDigitStore.h"
44 #include "AliMUONVTriggerStore.h"
45 #include "AliMUONVDigit.h"
46 #include "AliMUONLocalTrigger.h"
47
48 #include "AliMUONDDLTrigger.h"
49 #include "AliMUONRegHeader.h"
50 #include "AliMUONDarcHeader.h"
51 #include "AliMpTriggerCrate.h"
52 #include "AliMpDDLStore.h"
53 #include "AliMUONLocalStruct.h"
54 #include "AliMpLocalBoard.h"
55 #include "AliMpCDB.h"
56
57 #include "AliMpPad.h"
58 #include "AliMUONGeometryTransformer.h"
59 #include "AliMpVSegmentation.h"
60 #include "AliMpSegmentation.h"
61 #include "AliMpConstants.h"
62
63 #include "AliMUONQADataMakerRec.h"
64
65 //-----------------------------------------------------------------------------
66 /// \class AliMUONQADataMakerRec
67 ///
68 /// MUON base class for quality assurance data (histo) maker
69 ///
70 /// \author C. Finck
71
72 /// \cond CLASSIMP
73 ClassImp(AliMUONQADataMakerRec)
74 /// \endcond
75            
76 //____________________________________________________________________________ 
77 AliMUONQADataMakerRec::AliMUONQADataMakerRec() : 
78     AliQADataMakerRec(AliQA::GetDetName(AliQA::kMUON), "MUON Quality Assurance Data Maker"),
79     fDigitStore(0x0),
80     fTriggerStore(0x0),
81     fDigitMaker(0x0)
82 {
83     /// ctor
84   fDigitStore = AliMUONVDigitStore::Create("AliMUONDigitStoreV1");
85   fDigitMaker = new AliMUONDigitMaker(kTRUE,kFALSE);
86 }
87
88 //____________________________________________________________________________ 
89 AliMUONQADataMakerRec::AliMUONQADataMakerRec(const AliMUONQADataMakerRec& qadm) :
90     AliQADataMakerRec(qadm),
91     fDigitStore(0x0),
92     fTriggerStore(0x0),
93     fDigitMaker(0x0)
94 {
95     ///copy ctor 
96     SetName((const char*)qadm.GetName()) ; 
97     SetTitle((const char*)qadm.GetTitle()); 
98
99     // Do not copy the digit store and digit maker, but create its own ones
100     fDigitStore = AliMUONVDigitStore::Create("AliMUONDigitStoreV1");
101     fDigitMaker = new AliMUONDigitMaker(kTRUE,kFALSE);
102 }
103
104 //__________________________________________________________________
105 AliMUONQADataMakerRec& AliMUONQADataMakerRec::operator = (const AliMUONQADataMakerRec& qadm )
106 {
107   /// Assignment operator
108
109   // check assignment to self
110   if (this == &qadm) return *this;
111
112   this->~AliMUONQADataMakerRec();
113   new(this) AliMUONQADataMakerRec(qadm);
114   return *this;
115 }
116
117 //__________________________________________________________________
118 AliMUONQADataMakerRec::~AliMUONQADataMakerRec()
119 {
120     /// dtor
121     delete fDigitStore;
122     delete fTriggerStore;
123     delete fDigitMaker;
124 }
125
126 //____________________________________________________________________________ 
127 void AliMUONQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray* list)
128 {
129     ///Detector specific actions at end of cycle
130
131     // Display trigger histos in a more user friendly way
132     DisplayTriggerInfo(task);
133
134     // do the QA checking
135     AliQAChecker::Instance()->Run(AliQA::kMUON, task, list) ;
136 }
137
138 //____________________________________________________________________________ 
139 void AliMUONQADataMakerRec::InitRaws()
140 {
141     /// create Raws histograms in Raws subdir
142     TH1I* h0 = new TH1I("hRawBusPatch", "buspatch distribution",  1932, 1, 1932); 
143     Add2RawsList(h0, kRawBusPatch);
144
145     TH1I* h1 = new TH1I("hRawCharge", "Charge distribution in rawdata", 4096, 0, 4095); 
146     Add2RawsList(h1, kRawCharge);
147                 
148     for (Int_t iDDL = 0; iDDL < 20; ++iDDL) 
149     {
150       TH1F* h2 = new TH1F(Form("%s%d", "hRawBusPatchDDL", iDDL), Form("%s %d","RAW Buspatch distribution for DDL", iDDL), 50, 0, 49); 
151       Add2RawsList(h2, kRawBuspatchDDL+iDDL);
152     }
153
154     TH3F* h3 = new TH3F("hTriggerScalersBendPlane", "Trigger scalers in bending plane",
155                         4, 10.5, 14.5,
156                         234, 0.5, 234.5,
157                         16, -0.5, 15.5);
158     h3->GetXaxis()->SetTitle("Chamber");
159     h3->GetYaxis()->SetTitle("Board");
160     h3->GetZaxis()->SetTitle("Strip");
161     Add2RawsList(h3, kTriggerScalersBP);
162
163     TH3F* h4 = new TH3F("hTriggerScalersNonBendPlane", "Trigger scalers in non-bending plane",
164                         4, 10.5, 14.5,
165                         234, 0.5, 234.5,
166                         16, -0.5, 15.5);
167     h4->GetXaxis()->SetTitle("Chamber");
168     h4->GetYaxis()->SetTitle("Board");
169     h4->GetZaxis()->SetTitle("Strip");
170     Add2RawsList(h4, kTriggerScalersNBP);
171
172     InitDisplayHistos(AliQA::kRAWS);
173 }
174
175 //____________________________________________________________________________ 
176 void AliMUONQADataMakerRec::InitRecPoints()
177 {
178     /// create Reconstructed Points histograms in RecPoints subdir
179     TH3F* h0 = new TH3F("hTriggerDigitsBendPlane", "Trigger digits in bending plane",
180                         4, 10.5, 14.5,
181                         234, 0.5, 234.5,
182                         16, -0.5, 15.5);
183     h0->GetXaxis()->SetTitle("Chamber");
184     h0->GetYaxis()->SetTitle("Board");
185     h0->GetZaxis()->SetTitle("Strip");
186     Add2RecPointsList(h0, kTriggerDigitsBendPlane);
187
188     TH3F* h1 = new TH3F("hTriggerDigitsNonBendPlane", "Trigger digits in non-bending plane",
189                         4, 10.5, 14.5,
190                         234, 0.5, 234.5,
191                         16, -0.5, 15.5);
192     h1->GetXaxis()->SetTitle("Chamber");
193     h1->GetYaxis()->SetTitle("Board");
194     h1->GetZaxis()->SetTitle("Strip");
195     Add2RecPointsList(h1, kTriggerDigitsNonBendPlane);
196
197     TH1F* h2 = new TH1F("hTriggeredBoards", "Triggered boards", 234, 0.5, 234.5);
198     Add2RecPointsList(h2, kTriggeredBoards);
199
200     InitDisplayHistos(AliQA::kRECPOINTS);
201 }
202
203
204 //____________________________________________________________________________ 
205 void AliMUONQADataMakerRec::InitESDs()
206 {
207     ///create ESDs histograms in ESDs subdir
208   TH1F* h0 = new TH1F("hESDnTracks", "ESDs track number distribution", 30, 0., 30.);  
209   Add2ESDsList(h0, kESDnTracks);
210
211   TH1F* h1 = new TH1F("hESDMomentum", "ESDs P distribution", 300, 0., 300) ; 
212   Add2ESDsList(h1, kESDMomentum);
213
214   TH1F* h2 = new TH1F("hESDPt", "ESDs Pt distribution", 200, 0., 50) ; 
215   Add2ESDsList(h2, kESDPt);
216
217   TH1F* h3 = new TH1F("hESDRapidity", "ESDs rapidity distribution", 200, -4.5,-2.) ; 
218   Add2ESDsList(h3, kESDRapidity);
219
220   for (Int_t i = 0; i < AliMUONConstants::NTrackingCh(); ++i) 
221   {
222     TH2F* h4 = new TH2F(Form("%s%d", "hESDClusterHitMap", i), 
223                      Form("%s %d", "ESD Clusters hit distribution for chamber", i),
224                      100, -1*AliMUONConstants::Rmax(i/2), AliMUONConstants::Rmax(i/2),
225                      100, -1*AliMUONConstants::Rmax(i/2), AliMUONConstants::Rmax(i/2)); 
226     Add2ESDsList(h4, kESDClusterHitMap+i);
227   }
228 }
229
230 //____________________________________________________________________________
231 void AliMUONQADataMakerRec::MakeRaws(AliRawReader* rawReader)
232 {
233     /// make QA for rawdata
234     Int_t busPatchId;
235     UShort_t manuId;  
236     UChar_t channelId;
237     UShort_t charge;
238
239     rawReader->Reset();
240     AliMUONRawStreamTracker rawStreamTrack(rawReader);
241     rawStreamTrack.First();
242     while( rawStreamTrack.Next(busPatchId, manuId, channelId, charge) ) {
243   
244       GetRawsData(kRawBusPatch)->Fill(busPatchId);
245       GetRawsData(kRawCharge)->Fill(charge);
246       Int_t iDDL = rawStreamTrack.GetCurentDDL();
247       GetRawsData(kRawBuspatchDDL + iDDL)->Fill(AliMpBusPatch::GetLocalBusID(busPatchId, iDDL));
248                 
249                   
250     } // Next digit
251
252
253     // Get trigger scalers
254
255     Int_t loCircuit=0;
256     AliMpCDB::LoadDDLStore();
257
258     AliMUONRawStreamTrigger rawStreamTrig(rawReader);
259     while (rawStreamTrig.NextDDL()) 
260     {
261       // If not a scaler event, do nothing
262       Bool_t scalerEvent =  rawReader->GetDataHeader()->GetL1TriggerMessage() & 0x1;
263       if(!scalerEvent) break;
264
265       AliMUONDDLTrigger* ddlTrigger = rawStreamTrig.GetDDLTrigger();
266       AliMUONDarcHeader* darcHeader = ddlTrigger->GetDarcHeader();
267
268       Int_t nReg = darcHeader->GetRegHeaderEntries();
269     
270       for(Int_t iReg = 0; iReg < nReg ;iReg++)
271       {   //reg loop
272
273         // crate info  
274         AliMpTriggerCrate* crate = AliMpDDLStore::Instance()->
275           GetTriggerCrate(rawStreamTrig.GetDDL(), iReg);
276
277         AliMUONRegHeader* regHeader =  darcHeader->GetRegHeaderEntry(iReg);
278
279         // loop over local structures
280         Int_t nLocal = regHeader->GetLocalEntries();
281         for(Int_t iLocal = 0; iLocal < nLocal; iLocal++) 
282         {
283           AliMUONLocalStruct* localStruct = regHeader->GetLocalEntry(iLocal);
284         
285           // if card exist
286           if (!localStruct) continue;
287           
288           loCircuit = crate->GetLocalBoardId(localStruct->GetId());
289
290           if ( !loCircuit ) continue; // empty slot
291
292           AliMpLocalBoard* localBoard = AliMpDDLStore::Instance()->GetLocalBoard(loCircuit, false);
293
294           // skip copy cards
295           if( !localBoard->IsNotified()) 
296             continue;
297
298           Int_t cathode = localStruct->GetComptXY()%2;
299
300           ERaw hindex = (cathode==0) ? kTriggerScalersBP : kTriggerScalersNBP;
301
302           // loop over strips
303           for (Int_t ibitxy = 0; ibitxy < 16; ++ibitxy) {
304             if(localStruct->GetXY1(ibitxy) > 0)
305               ((TH3F*)GetRawsData(hindex))->Fill(11+0, loCircuit, ibitxy, 2*localStruct->GetXY1(ibitxy));
306             if(localStruct->GetXY2(ibitxy) > 0)
307               ((TH3F*)GetRawsData(hindex))->Fill(11+1, loCircuit, ibitxy, 2*localStruct->GetXY2(ibitxy));
308             if(localStruct->GetXY3(ibitxy) > 0)
309               ((TH3F*)GetRawsData(hindex))->Fill(11+2, loCircuit, ibitxy, 2*localStruct->GetXY3(ibitxy));
310             if(localStruct->GetXY4(ibitxy) > 0)
311               ((TH3F*)GetRawsData(hindex))->Fill(11+3, loCircuit, ibitxy, 2*localStruct->GetXY4(ibitxy));
312           } // loop on strips
313         } // iLocal
314       } // iReg
315     } // NextDDL
316 }
317
318 //____________________________________________________________________________
319 void AliMUONQADataMakerRec::MakeRecPoints(TTree* clustersTree)
320 {
321   
322     /// makes data from trigger response
323       
324     // Fired pads info
325     fDigitStore->Clear();
326
327     if (!fTriggerStore) fTriggerStore = AliMUONVTriggerStore::Create(*clustersTree);
328     fTriggerStore->Clear();
329     fTriggerStore->Connect(*clustersTree, false);
330     clustersTree->GetEvent(0);
331
332     AliMUONLocalTrigger* locTrg;
333     TIter nextLoc(fTriggerStore->CreateLocalIterator());
334
335     while ( ( locTrg = static_cast<AliMUONLocalTrigger*>(nextLoc()) ) ) 
336     {
337       if (locTrg->IsNull()) continue;
338    
339       TArrayS xyPattern[2];
340       locTrg->GetXPattern(xyPattern[0]);
341       locTrg->GetYPattern(xyPattern[1]);
342
343       Int_t nBoard = locTrg->LoCircuit();
344
345       Bool_t xTrig=locTrg->IsTrigX();
346       Bool_t yTrig=locTrg->IsTrigY();
347     
348       if (xTrig && yTrig)
349         ((TH1F*)GetRecPointsData(kTriggeredBoards))->Fill(nBoard);
350
351       fDigitMaker->TriggerDigits(nBoard, xyPattern, *fDigitStore);
352     }
353
354     TIter nextDigit(fDigitStore->CreateIterator());
355     AliMUONVDigit* mDigit;
356     while ( ( mDigit = static_cast<AliMUONVDigit*>(nextDigit()) ) )
357     {
358       Int_t detElemId = mDigit->DetElemId();
359       Int_t ch = detElemId/100;
360       Int_t localBoard = mDigit->ManuId();
361       Int_t channel = mDigit->ManuChannel();
362       Int_t cathode = mDigit->Cathode();
363       ERecPoints hindex 
364         = ( cathode == 0 ) ? kTriggerDigitsBendPlane : kTriggerDigitsNonBendPlane;
365       
366       ((TH3F*)GetRecPointsData(hindex))->Fill(ch, localBoard, channel);
367     }
368 }
369
370 //____________________________________________________________________________
371 void AliMUONQADataMakerRec::MakeESDs(AliESDEvent* esd)
372 {
373     /// make QA data from ESDs
374
375     TLorentzVector v1;
376
377     Int_t nTracks = (Int_t)esd->GetNumberOfMuonTracks() ; 
378     GetESDsData(0)->Fill(nTracks);
379
380     for (Int_t iTrack = 0; iTrack < nTracks; ++iTrack) {
381
382       AliESDMuonTrack* muonTrack = esd->GetMuonTrack(iTrack);
383       muonTrack->LorentzP(v1);
384
385       GetESDsData(1)->Fill(v1.P());
386       GetESDsData(2)->Fill(v1.Pt());
387       GetESDsData(3)->Fill(v1.Rapidity());
388
389       TClonesArray clusters =  muonTrack->GetClusters();
390
391       for (Int_t iCluster = 0; iCluster <clusters.GetEntriesFast(); ++iCluster) {
392         AliESDMuonCluster* cluster = (AliESDMuonCluster*)clusters.At(iCluster);
393         GetESDsData(kESDClusterHitMap+cluster->GetChamberId())
394           ->Fill(cluster->GetX(), cluster->GetY());
395       }
396     }
397 }
398
399 //____________________________________________________________________________ 
400 void AliMUONQADataMakerRec::StartOfDetectorCycle()
401 {
402     /// Detector specific actions at start of cycle
403   
404 }
405
406 //____________________________________________________________________________ 
407 void AliMUONQADataMakerRec::DisplayTriggerInfo(AliQA::TASKINDEX_t task)
408 {
409   //
410   /// Display trigger information in a user-friendly way:
411   /// from local board and strip numbers to their position on chambers
412   //
413   if(task!=AliQA::kRECPOINTS && task!=AliQA::kRAWS) return;
414
415   AliMpCDB::LoadDDLStore();
416
417   Int_t detElemId;
418   Float_t xg1, yg1, zg1, xg2, yg2, zg2, binContent;
419   Float_t xLocal1=0., yLocal1=0., xLocal2=0., yLocal2=0., xWidth=0., yWidth=0.;
420   Float_t x1Label, y1Label, x2Label, y2Label;
421   Int_t x1Int, x2Int, y1Int, y2Int;
422   Int_t binCh, binBoard, binStrip;
423
424   AliMUONGeometryTransformer transform;
425   transform.LoadGeometryData();
426
427   TH3F* histoStrips=0x0;
428   TH2F* histoDisplayStrips=0x0;
429   TH2F* histoDisplayBoards=0x0;
430   TH1F* histoBoards=0x0;
431
432   const Float_t kShift = 0.15;
433
434   for (Int_t iCath = 0; iCath < 2; ++iCath)
435   {
436
437     if(task==AliQA::kRECPOINTS){
438       ERecPoints hindex 
439         = ( iCath == 0 ) ? kTriggerDigitsBendPlane : kTriggerDigitsNonBendPlane;
440       histoStrips = (TH3F*)GetRecPointsData(hindex);
441       histoBoards = (TH1F*)GetRecPointsData(kTriggeredBoards);
442     }
443     else if(task==AliQA::kRAWS){
444       ERaw hindex 
445         = ( iCath == 0 ) ? kTriggerScalersBP : kTriggerScalersNBP;
446       histoStrips = (TH3F*)GetRawsData(hindex);
447       if(histoStrips->GetEntries()==0) return; // No scalers found
448     }
449     
450     // Load mapping
451     if ( ! AliMpSegmentation::Instance(kFALSE) ) {
452       /// Load mapping
453       if ( ! AliMpCDB::LoadDDLStore() ) {
454         AliFatal("Could not access mapping from OCDB !");
455       }
456     }  
457
458     for (Int_t iChamber = 0; iChamber < 4; ++iChamber)
459     {
460       Int_t iCh = iChamber + AliMpConstants::NofTrackingChambers();
461
462       if(task==AliQA::kRECPOINTS){
463         histoDisplayStrips = (TH2F*)GetRecPointsData(kTriggerDigitsDisplay + 4*iCath + iChamber);
464         histoDisplayBoards = (TH2F*)GetRecPointsData(kTriggerBoardsDisplay);
465       }
466       else if(task==AliQA::kRAWS){
467         histoDisplayStrips = (TH2F*)GetRawsData(kTriggerScalersDisplay + 4*iCath + iChamber);
468       }
469
470       for(Int_t iBoard=1; iBoard<=234; iBoard++){
471         // get detElemId
472         detElemId = AliMpDDLStore::Instance()->GetDEfromLocalBoard(iBoard, iCh);
473
474         if (!detElemId) continue;
475
476         AliMpLocalBoard* localBoard = AliMpDDLStore::Instance()->GetLocalBoard(iBoard, false);
477
478         // skip copy cards
479         if( !localBoard->IsNotified()) 
480           continue;
481
482         const AliMpVSegmentation* seg 
483           = AliMpSegmentation::Instance()
484           ->GetMpSegmentation(detElemId, AliMp::GetCathodType(iCath));  
485         
486         // loop over strips
487         for (Int_t ibitxy = 0; ibitxy < 16; ++ibitxy) 
488         {
489           AliMpPad pad = seg->PadByLocation(AliMpIntPair(iBoard,ibitxy),kFALSE);
490                         
491           if (!pad.IsValid())
492             continue;
493
494           if(iCath==0){ // Geometry info from bending plane only
495             if(ibitxy==0) {
496               xLocal1 = pad.Position().X();
497               yLocal1 = pad.Position().Y();
498               xWidth = pad.Dimensions().X();
499               yWidth = pad.Dimensions().Y();
500             }
501             xLocal2 = pad.Position().X();
502             yLocal2 = pad.Position().Y();
503           }
504
505           Float_t dimX = pad.Dimensions().X();
506           Float_t dimY = pad.Dimensions().Y();
507
508           Float_t stripX = pad.Position().X();
509           Float_t stripY = pad.Position().Y();
510                    
511           transform.Local2Global(detElemId, stripX, stripY, 0, xg1, yg1, zg1);
512
513           x1Int = histoDisplayStrips->GetXaxis()->FindBin(xg1 - dimX + kShift);
514           y1Int = histoDisplayStrips->GetYaxis()->FindBin(yg1 - dimY + kShift);
515           x2Int = histoDisplayStrips->GetXaxis()->FindBin(xg1 + dimX - kShift);
516           y2Int = histoDisplayStrips->GetYaxis()->FindBin(yg1 + dimY - kShift);
517
518           binCh = histoStrips->GetXaxis()->FindBin(iCh+1);
519           binBoard = histoStrips->GetYaxis()->FindBin(iBoard);
520           binStrip = histoStrips->GetZaxis()->FindBin(ibitxy);
521           binContent = histoStrips->GetBinContent(binCh, binBoard, binStrip);
522
523           if(binContent==0) continue;
524
525           for(Int_t binX=x1Int; binX<=x2Int; binX++){
526             for(Int_t binY=y1Int; binY<=y2Int; binY++){       
527               histoDisplayStrips->SetBinContent(binX, binY, binContent);
528             }
529           }
530         }// ibitxy
531
532         if(task != AliQA::kRECPOINTS) continue;
533         if(iCath>0 || iChamber>0) continue; // Geometry info from bending plane only
534         transform.Local2Global(detElemId, xLocal1, yLocal1, 0, xg1, yg1, zg1);
535         transform.Local2Global(detElemId, xLocal2, yLocal2, 0, xg2, yg2, zg2);
536
537         // Per board
538         x1Label = TMath::Min(xg1,xg2) - xWidth + kShift;
539         y1Label = TMath::Min(yg1,yg2) - yWidth + kShift;
540         x2Label = TMath::Max(xg1,xg2) + xWidth - kShift;
541         y2Label = TMath::Max(yg1,yg2) + yWidth - kShift;
542
543         x1Int = histoDisplayBoards->GetXaxis()->FindBin(x1Label);
544         y1Int = histoDisplayBoards->GetYaxis()->FindBin(y1Label);
545         x2Int = histoDisplayBoards->GetXaxis()->FindBin(x2Label);
546         y2Int = histoDisplayBoards->GetYaxis()->FindBin(y2Label);
547
548         binBoard = histoBoards->GetXaxis()->FindBin(iBoard);
549         binContent = histoBoards->GetBinContent(binBoard);
550
551         if(binContent==0) continue;
552
553         for(Int_t binX=x1Int; binX<=x2Int; binX++){
554           for(Int_t binY=y1Int; binY<=y2Int; binY++){
555             histoDisplayBoards->SetBinContent(binX, binY, binContent);
556           }
557         }
558       } // iBoard
559     } // iChamber
560   } // iCath
561 }
562
563
564 //____________________________________________________________________________ 
565 void AliMUONQADataMakerRec::InitDisplayHistos(AliQA::TASKINDEX_t task)
566 {
567   //
568   /// Initialize trigger information display histograms,
569   /// using mapping and geometry
570   //
571   AliMpCDB::LoadDDLStore();
572
573   const Int_t kNumOfBoards = AliMpConstants::NofLocalBoards();
574
575   AliMUONGeometryTransformer transform;
576   transform.LoadGeometryData();
577
578   TString cathCode[2] = {"BendPlane", "NonBendPlane"};
579
580   TArrayF xLocal1[4], yLocal1[4], xLocal2[4], yLocal2[4], xWidth[4], yWidth[4];
581
582   TArrayF xAxisStrip[8];
583   TArrayF yAxisStrip[8];
584   TArrayF xAxisBoard[8];
585   TArrayF yAxisBoard[8];
586
587   TH2F* histoStrips=0x0;
588   TH2F* histoBoards=0x0;
589
590   const Float_t kResetValue=1234567.;
591    
592   for(Int_t ch=0; ch<4; ch++){
593     xLocal1[ch].Set(kNumOfBoards);
594     yLocal1[ch].Set(kNumOfBoards);
595     xLocal2[ch].Set(kNumOfBoards);
596     yLocal2[ch].Set(kNumOfBoards);
597     xWidth[ch].Set(kNumOfBoards);
598     yWidth[ch].Set(kNumOfBoards);
599   }
600
601   for(Int_t cath=0; cath<2; cath++){
602     for(Int_t ch=0; ch<4; ch++){
603       Int_t chCath = 4*cath + ch;
604       xAxisBoard[chCath].Set(60);
605       xAxisBoard[chCath].Reset(kResetValue);
606       yAxisBoard[chCath].Set(60);
607       yAxisBoard[chCath].Reset(kResetValue);
608
609       xAxisStrip[chCath].Set(700);
610       xAxisStrip[chCath].Reset(kResetValue);
611       yAxisStrip[chCath].Set(700);
612       yAxisStrip[chCath].Reset(kResetValue);
613     }
614   }
615    
616   Float_t xg1, yg1, zg1, xg2, yg2, zg2;
617
618   TString histoName, histoTitle;
619
620   const Float_t kShift = 0.;
621
622   // Load mapping
623   if ( ! AliMpSegmentation::Instance(kFALSE) ) {
624     /// Load mapping
625     if ( ! AliMpCDB::LoadDDLStore() ) {
626       AliFatal("Could not access mapping from OCDB !");
627     }
628   }  
629
630   for(Int_t iCath=0; iCath<2; iCath++){
631     for (Int_t iChamber = 0; iChamber < 4; ++iChamber) {
632       Int_t iCh = iChamber + AliMpConstants::NofTrackingChambers();
633       for(Int_t iLoc = 0; iLoc < 234; iLoc++) {  
634         Int_t iBoard = iLoc+1;
635         Int_t detElemId = AliMpDDLStore::Instance()->GetDEfromLocalBoard(iBoard, iCh);
636
637         if (!detElemId) continue;
638
639         AliMpLocalBoard* localBoard = AliMpDDLStore::Instance()->GetLocalBoard(iBoard, kFALSE);
640
641         // skip copy cards
642         if( !localBoard->IsNotified()) 
643           continue;
644
645         // get segmentation
646         const AliMpVSegmentation* seg = AliMpSegmentation::Instance()
647           ->GetMpSegmentation(detElemId,
648                               AliMp::GetCathodType(iCath));
649
650         Int_t chCath = 4*iCath + iChamber;
651         // loop over strips
652         for (Int_t ibitxy = 0; ibitxy < 16; ++ibitxy) {
653           // get pad from electronics
654           AliMpPad pad = seg->PadByLocation(AliMpIntPair(iBoard,ibitxy),kFALSE);
655
656           if (!pad.IsValid()) continue;
657
658           if(iCath==0){
659             if(ibitxy==0) {
660               xLocal1[iChamber][iLoc] = pad.Position().X();
661               yLocal1[iChamber][iLoc] = pad.Position().Y();
662               xWidth[iChamber][iLoc] = pad.Dimensions().X();
663               yWidth[iChamber][iLoc] = pad.Dimensions().Y();
664             }
665             xLocal2[iChamber][iLoc] = pad.Position().X();
666             yLocal2[iChamber][iLoc] = pad.Position().Y();
667           }
668
669           // Check fired pads
670           Float_t dimX = pad.Dimensions().X();
671           Float_t dimY = pad.Dimensions().Y();
672      
673           Float_t stripX = pad.Position().X();
674           Float_t stripY = pad.Position().Y();
675
676           transform.Local2Global(detElemId, stripX, stripY, 0, xg1, yg1, zg1);
677
678           AddSortedPoint(xg1 - dimX + kShift, xAxisStrip[chCath], kResetValue);
679           AddSortedPoint(xg1 + dimX - kShift, xAxisStrip[chCath], kResetValue);
680
681           AddSortedPoint(yg1 - dimY + kShift, yAxisStrip[chCath], kResetValue);
682           AddSortedPoint(yg1 + dimY - kShift, yAxisStrip[chCath], kResetValue);
683         } // loop on strips  
684
685         transform.Local2Global(detElemId, xLocal1[iChamber][iLoc], yLocal1[iChamber][iLoc], 0, xg1, yg1, zg1);
686         transform.Local2Global(detElemId, xLocal2[iChamber][iLoc], yLocal2[iChamber][iLoc], 0, xg2, yg2, zg2);
687
688         // Per board
689         AddSortedPoint(TMath::Min(xg1,xg2) - xWidth[iChamber][iLoc] + kShift, xAxisBoard[chCath], kResetValue);
690         AddSortedPoint(TMath::Max(xg1,xg2) + xWidth[iChamber][iLoc] - kShift, xAxisBoard[chCath], kResetValue);
691
692         AddSortedPoint(TMath::Min(yg1,yg2) - yWidth[iChamber][iLoc] + kShift, yAxisBoard[chCath], kResetValue);
693         AddSortedPoint(TMath::Max(yg1,yg2) + yWidth[iChamber][iLoc] - kShift, yAxisBoard[chCath], kResetValue);
694       } // loop on local boards 
695     } // loop on chambers
696   } // loop on cathodes
697
698   const Float_t kMinDiff = 0.1;
699
700   // Book histos
701   for(Int_t iCath=0; iCath<2; iCath++){
702     for(Int_t iChamber=0; iChamber<4; iChamber++){
703       Int_t chCath = 4*iCath + iChamber;
704       Int_t ipoint=0;
705       while(TMath::Abs(xAxisStrip[chCath][ipoint]-kResetValue)>kMinDiff) { ipoint++; }
706       xAxisStrip[chCath].Set(ipoint);
707
708       ipoint = 0;
709       while(TMath::Abs(yAxisStrip[chCath][ipoint]-kResetValue)>kMinDiff) { ipoint++; }
710       yAxisStrip[chCath].Set(ipoint);
711
712       ipoint = 0;
713       while(TMath::Abs(xAxisBoard[chCath][ipoint]-kResetValue)>kMinDiff) { ipoint++; }
714       xAxisBoard[chCath].Set(ipoint);
715
716       ipoint = 0;
717       while(TMath::Abs(yAxisBoard[chCath][ipoint]-kResetValue)>kMinDiff) { ipoint++; }
718       yAxisBoard[chCath].Set(ipoint);
719
720       if(task==AliQA::kRECPOINTS){
721         histoName = Form("hTriggerDigits%sChamber%i", cathCode[iCath].Data(), 11+iChamber);
722         histoTitle = Form("Chamber %i: Fired pads %s", 11+iChamber, cathCode[iCath].Data());
723       }
724       else if(task==AliQA::kRAWS){
725         histoName = Form("hScalers%sChamber%i", cathCode[iCath].Data(), 11+iChamber);
726         histoTitle = Form("Chamber %i: Scalers %s", 11+iChamber, cathCode[iCath].Data());
727       }
728       histoStrips = new TH2F(histoName.Data(), histoTitle.Data(), 
729                              xAxisStrip[chCath].GetSize()-1, xAxisStrip[chCath].GetArray(),
730                              yAxisStrip[chCath].GetSize()-1, yAxisStrip[chCath].GetArray());
731       histoStrips->SetXTitle("X (cm)");
732       histoStrips->SetYTitle("Y (cm)");
733
734       if(task==AliQA::kRECPOINTS){
735         Add2RecPointsList(histoStrips, kTriggerDigitsDisplay + 4*iCath + iChamber);
736       }
737       else if(task==AliQA::kRAWS){
738         Add2RawsList(histoStrips, kTriggerScalersDisplay + 4*iCath + iChamber);
739       }
740
741       if(task != AliQA::kRECPOINTS) continue;
742       if(iCath>0 || iChamber>0) continue;
743
744       histoName = "hFiredBoardsDisplay";
745       histoTitle = "Fired boards";
746       histoBoards = new TH2F(histoName.Data(), histoTitle.Data(), 
747                              xAxisBoard[chCath].GetSize()-1, xAxisBoard[chCath].GetArray(),
748                              yAxisBoard[chCath].GetSize()-1, yAxisBoard[chCath].GetArray());
749       histoBoards->SetXTitle("X (cm)");
750       histoBoards->SetYTitle("Y (cm)");
751
752       Add2RecPointsList(histoBoards, kTriggerBoardsDisplay + 4*iCath + iChamber);
753     } // loop on chamber
754   } // loop on cathode
755 }
756
757
758 //____________________________________________________________________________ 
759 Bool_t AliMUONQADataMakerRec::AddSortedPoint(Float_t currVal, TArrayF& position, const Float_t kResetValue)
760 {
761   //
762   /// Add sorted point in array according to an increasing order.
763   /// Used to build display histograms axis.
764   //
765   Int_t nEntries = position.GetSize();
766   Float_t tmp1, tmp2;
767   const Float_t kMinDiff = 0.1;
768   for(Int_t i=0; i<nEntries; i++){
769     if(TMath::Abs(position[i]-currVal)<kMinDiff) return kFALSE;
770     if(TMath::Abs(position[i]-kResetValue)<kMinDiff) {
771       position[i] = currVal;
772       return kTRUE;
773     }
774     if(currVal>position[i]) continue;
775     tmp1 = position[i];
776     position[i] = currVal;
777     for(Int_t j=i+1; j<nEntries; j++){
778       tmp2 = position[j];
779       position[j] = tmp1;
780       tmp1 = tmp2;
781       if(tmp1==kResetValue) break;
782     }
783     return kTRUE;
784   }
785   return kFALSE;
786 }
787