1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 // --- MUON header files ---
19 #include "AliMUONTriggerDisplay.h"
21 #include "AliMpDDLStore.h"
22 #include "AliMpVSegmentation.h"
23 #include "AliMpSegmentation.h"
24 #include "AliMpConstants.h"
26 #include "AliMpLocalBoard.h"
29 // --- AliRoot header files ---
32 // --- ROOT system ---
37 //-----------------------------------------------------------------------------
38 /// \class AliMUONTriggerDisplay
40 /// MUON base class for converting histos as a function of strip/board/slat number
41 /// into display histos showing the detection element position in the
44 /// Input histos can be given as:
45 /// - TH2 with x -> board # [1-234] and y -> strip # in board [0-15]. Option: kDisplayStrips
46 /// - TH1 with x -> board # [1-234] Option: kDisplayBoards
47 /// - TH1 with x -> slat # [0-17] Option: kDisplaySlats
52 ClassImp(AliMUONTriggerDisplay)
55 //____________________________________________________________________________
56 AliMUONTriggerDisplay::AliMUONTriggerDisplay() :
62 //__________________________________________________________________
63 AliMUONTriggerDisplay::~AliMUONTriggerDisplay()
69 //____________________________________________________________________________
70 TH2* AliMUONTriggerDisplay::GetEmptyDisplayHisto(TString displayHistoName, EDisplayType displayType,
71 Int_t cathode, Int_t chamber,
72 TString displayHistoTitle)
75 /// Return the display histogram with optimized binning
76 /// but do not fill it.
78 TH2F* displayHisto = new TH2F();
80 InitOrDisplayTriggerInfo(0x0, displayHisto, displayType,
81 cathode, chamber, displayHistoName, displayHistoTitle);
87 //____________________________________________________________________________
88 TH2* AliMUONTriggerDisplay::GetBoardNumberHisto(TString displayHistoName,
90 TString displayHistoTitle)
93 /// Return the display histogram with optimized binning
94 /// and fill it with the board number
97 const Int_t kNboards = AliMpConstants::NofLocalBoards();
98 TH1F* inputHisto = new TH1F("boardNumbers","Board Numbers",kNboards,0.5,(Float_t)kNboards + 0.5);
99 for(Int_t ibin=1; ibin<=kNboards; ibin++){
100 inputHisto->Fill(ibin,ibin);
103 TH2F* displayHisto = (TH2F*)GetEmptyDisplayHisto(displayHistoName, kDisplayBoards, 0, chamber, displayHistoTitle);
104 FillDisplayHistogram(inputHisto,displayHisto,kDisplayBoards,0,chamber,kNumbered);
108 displayHisto->SetStats(kFALSE);
113 //____________________________________________________________________________
114 TH2* AliMUONTriggerDisplay::GetDisplayHistogram(TH1* inputHisto, TString displayHistoName,
115 EDisplayType displayType, Int_t cathode,
116 Int_t chamber, TString displayHistoTitle,
117 EDisplayOption displayOpt)
120 /// Get histogram displaying the information contained in the input histogram
122 TH2* displayHisto = GetEmptyDisplayHisto(displayHistoName, displayType,
123 cathode, chamber, displayHistoTitle);
125 FillDisplayHistogram(inputHisto, displayHisto, displayType, cathode, chamber, displayOpt);
130 //____________________________________________________________________________
131 TH2* AliMUONTriggerDisplay::GetDisplayHistogram(TGraph* inputGraph, TString displayHistoName,
132 EDisplayType displayType, Int_t cathode,
133 Int_t chamber, TString displayHistoTitle,
134 EDisplayOption displayOpt)
137 /// Get histogram displaying the information contained in the input graph
139 TH2* displayHisto = GetEmptyDisplayHisto(displayHistoName, displayType,
140 cathode, chamber, displayHistoTitle);
142 FillDisplayHistogram(inputGraph, displayHisto, displayType, cathode, chamber, displayOpt);
147 //____________________________________________________________________________
148 Bool_t AliMUONTriggerDisplay::FillDisplayHistogram(TH1* inputHisto, TH2* displayHisto,
149 EDisplayType displayType, Int_t cathode,
150 Int_t chamber, EDisplayOption displayOpt)
153 /// Fill a previously initialized display histogram
154 /// with the information contained in inputHisto.
155 /// To get initialized display, please use GetEmptyDisplayHisto method
157 return InitOrDisplayTriggerInfo(inputHisto, displayHisto, displayType,
158 cathode, chamber, "", "",displayOpt);
161 //____________________________________________________________________________
162 Bool_t AliMUONTriggerDisplay::FillDisplayHistogram(TGraph* inputGraph, TH2* displayHisto,
163 EDisplayType displayType, Int_t cathode,
164 Int_t chamber, EDisplayOption displayOpt)
167 /// Fill a previously initialized display histogram
168 /// with the information contained in inputGraph.
169 /// To get initialized display, please use GetEmptyDisplayHisto method
171 return InitOrDisplayTriggerInfo(inputGraph, displayHisto, displayType,
172 cathode, chamber, "", "",displayOpt);
175 //____________________________________________________________________________
176 Bool_t AliMUONTriggerDisplay::InitOrDisplayTriggerInfo(TObject* inputObject, TH2* displayHisto,
177 EDisplayType displayType,
178 Int_t cathode, Int_t chamber,
179 TString displayHistoName, TString displayHistoTitle,
180 EDisplayOption displayOpt)
183 /// Initialize trigger information display histograms using mapping
184 /// Trigger information is displayed in a user-friendly way:
185 /// from local board and strip numbers to their position on chambers.
189 if ( ! AliMpSegmentation::Instance(kFALSE) ) {
191 if ( ! AliMpCDB::LoadDDLStore() ) {
192 AliError("Could not access mapping from OCDB !");
193 AliError("Histograms are not initialized !");
198 Int_t iCh = (chamber > AliMpConstants::NofTrackingChambers()) ? chamber - 11 : chamber;
199 Int_t iCath = cathode;
206 Float_t yOffsetLine, xOffsetLine = 0.;
208 const Float_t kResetValue=1234567.;
212 xAxisBoard.Reset(kResetValue);
214 yAxisBoard.Reset(kResetValue);
217 xAxisStrip.Reset(kResetValue);
219 yAxisStrip.Reset(kResetValue);
221 else if(!displayHisto){
222 AliWarning("Display histogram not initialized. Please initialize it first!");
226 TH1* inputHisto = dynamic_cast<TH1*>(inputObject);
228 if ( inputHisto->GetEntries() == 0 ) {
233 TGraph* inputGraph = dynamic_cast<TGraph*>(inputObject);
234 if ( inputGraph->GetN() == 0 ){
240 Float_t xWidth, yWidth, yWidthSlat=0., xWidthCol=0.;
242 Float_t x1b=0., x2b=0., y1b=0., y2b=0.;
243 Float_t xcPad, ycPad;
247 const Float_t kShiftB = 0.5;
248 const Float_t kShiftS = 0.1;
250 const Float_t kShiftX = (iCath==0) ? kShiftB : kShiftS;
251 const Float_t kShiftY = (iCath==0) ? kShiftS : kShiftB;
252 const Float_t kShiftEl = (displayType==kDisplaySlats) ? 0.01 : kShiftB;
254 Int_t iChamber = iCh + AliMpConstants::NofTrackingChambers();
255 for(Int_t iLoc = 0; iLoc < AliMpConstants::NofLocalBoards(); iLoc++) {
256 Int_t iBoard = iLoc+1;
257 Int_t detElemId = AliMpDDLStore::Instance()->GetDEfromLocalBoard(iBoard, iChamber);
259 if (!detElemId) continue;
261 AliMpLocalBoard* localBoard = AliMpDDLStore::Instance()->GetLocalBoard(iBoard, kFALSE);
264 if( !localBoard->IsNotified())
268 const AliMpVSegmentation* seg[2] = {
269 AliMpSegmentation::Instance()->GetMpSegmentation(detElemId, AliMp::GetCathodType(0)),
270 AliMpSegmentation::Instance()->GetMpSegmentation(detElemId, AliMp::GetCathodType(1))};
273 AliMpPad pad1 = seg[1]->PadByLocation(iBoard,0,kFALSE);
274 yWidthSlat = pad1.GetDimensionY();
275 AliMpPad pad0 = seg[0]->PadByLocation(iBoard,0,kFALSE);
276 xOffsetLine = TMath::Abs(pad0.GetPositionX()) + pad0.GetDimensionX();
277 xWidthCol = 2.* pad0.GetDimensionX();
280 // Get ideal global position of DetElemId center
281 slat = detElemId%100;
282 line = (4 + slat)%18;
288 yOffsetLine = (Float_t)(line - 4) * 2. * yWidthSlat;
290 for(Int_t cath=0; cath<AliMpConstants::NofCathodes(); cath++){
292 // necessary because strip info is read for each cathode
293 // board info is read only from cathode 0
296 for (Int_t iStrip = 0; iStrip < 16; ++iStrip) {
297 // get pad from electronics
300 if (cath && localBoard->GetSwitch(AliMpLocalBoard::kZeroAllYLSB)) offset = -8;
302 AliMpPad pad = seg[cath]->PadByLocation(iBoard,iStrip+offset,kFALSE);
304 if (!pad.IsValid()) continue;
306 xWidth = pad.GetDimensionX();
307 yWidth = pad.GetDimensionY();
308 xcPad = sign * (pad.GetPositionX() + xOffsetLine);
309 if(line==4) xcPad += 0.75 * sign * xWidthCol;
310 ycPad = pad.GetPositionY() + yOffsetLine;
313 x1 = xcPad - xWidth + kShiftX;
314 y1 = ycPad - yWidth + kShiftY;
315 x2 = xcPad + xWidth - kShiftX;
316 y2 = ycPad + yWidth - kShiftY;
319 AddSortedPoint(x1, xAxisStrip, kResetValue);
320 AddSortedPoint(x2, xAxisStrip, kResetValue);
322 AddSortedPoint(y1, yAxisStrip, kResetValue);
323 AddSortedPoint(y2, yAxisStrip, kResetValue);
325 else if(displayType==kDisplayStrips)
326 FillBins(inputObject, displayHisto, iBoard, iStrip, x1, x2, y1, y2, kShiftX, kShiftY, displayOpt);
330 if(iStrip+offset==0) {
331 x1b = xcPad - xWidth + kShiftEl;
332 y1b = ycPad - yWidth + kShiftEl;
334 x2b = xcPad + xWidth - kShiftEl;
335 y2b = ycPad + yWidth - kShiftEl;
339 // if iCath==0 strip info and board info are both filled -> break!
340 // if iCath==1 board info is filled at cath==0. Strip info to be filled at cath==1
342 } // loop on cathodes
346 AddSortedPoint(x1b, xAxisBoard, kResetValue);
347 AddSortedPoint(x2b, xAxisBoard, kResetValue);
349 AddSortedPoint(y1b, yAxisBoard, kResetValue);
350 AddSortedPoint(y2b, yAxisBoard, kResetValue);
352 else if(displayType==kDisplayBoards)
353 FillBins(inputObject, displayHisto, iBoard, -1, x1b, x2b, y1b, y2b, kShiftEl, kShiftEl, displayOpt);
354 else if(displayType==kDisplaySlats)
355 FillBins(inputObject, displayHisto, slat, -1, x1b, x2b, y1b, y2b, kShiftEl, kShiftEl, displayOpt);
356 } // loop on local boards
358 if ( inputObject ) return kTRUE;
360 displayHisto->Reset();
363 const Float_t kMinDiff = 0.1;
365 TArrayD* currArray[4] = {&xAxisStrip, &yAxisStrip,
366 &xAxisBoard, &yAxisBoard};
367 for(Int_t iaxis=0; iaxis<4; iaxis++){
369 while(TMath::Abs((*currArray[iaxis])[ipoint]-kResetValue)>kMinDiff) {
372 if(ipoint>currArray[iaxis]->GetSize()-2)
373 AliWarning(Form("Array size (%i) lower than the number of points!", currArray[iaxis]->GetSize()));
374 currArray[iaxis]->Set(ipoint);
379 displayHisto->SetBins(xAxisStrip.GetSize()-1, xAxisStrip.GetArray(),
380 yAxisStrip.GetSize()-1, yAxisStrip.GetArray());
384 displayHisto->SetBins(xAxisBoard.GetSize()-1, xAxisBoard.GetArray(),
385 yAxisBoard.GetSize()-1, yAxisBoard.GetArray());
389 displayHisto->SetName(displayHistoName.Data());
390 displayHisto->SetTitle(displayHistoTitle.Data());
391 displayHisto->SetXTitle("X (cm)");
392 displayHisto->SetYTitle("Y (cm)");
393 //displayHisto->SetStats(kFALSE);
399 //____________________________________________________________________________
400 Bool_t AliMUONTriggerDisplay::AddSortedPoint(Float_t currVal, TArrayD& position, const Float_t kResetValue)
403 /// Add sorted point in array according to an increasing order.
404 /// Used to build display histograms axis.
406 Int_t nEntries = position.GetSize()-1;
408 const Float_t kMinDiff = 0.1;
409 for(Int_t i=0; i<nEntries; i++){
410 if(TMath::Abs(position[i]-currVal)<kMinDiff) return kFALSE;
411 if(TMath::Abs(position[i]-kResetValue)<kMinDiff) {
412 position[i] = currVal;
415 if(currVal>position[i]) continue;
417 position[i] = currVal;
418 for(Int_t j=i+1; j<nEntries; j++){
422 if(tmp1==kResetValue) break;
430 //____________________________________________________________________________
431 void AliMUONTriggerDisplay::FillBins(TObject* inputObject, TH2* displayHisto,
432 Int_t iElement1, Int_t iElement2,
433 Float_t x1, Float_t x2, Float_t y1, Float_t y2,
434 const Float_t kShiftX, const Float_t kShiftY,
435 EDisplayOption displayOpt)
438 /// Given the bin in inputHisto, search the corresponding bins
439 /// in display histo and fill it.
441 TString className = inputObject->ClassName();
443 Float_t binContent=0;
444 if ( className.Contains("TH") ) {
445 TH1* inputHisto = dynamic_cast<TH1*>(inputObject);
446 Int_t binX = inputHisto->GetXaxis()->FindBin(iElement1);
447 if(className.Contains("TH2")) {
448 binY = inputHisto->GetYaxis()->FindBin(iElement2);
449 binContent = inputHisto->GetBinContent(binX, binY);
451 else binContent = inputHisto->GetBinContent(binX);
454 TGraph* inputGraph = dynamic_cast<TGraph*>(inputObject);
456 for ( Int_t ipt=0; ipt<inputGraph->GetN(); ipt++){
457 inputGraph->GetPoint(ipt, xpt, ypt);
458 if ( TMath::Abs(xpt - iElement1) < 0.1 ) {
466 if(displayOpt==kShowZeroes) binContent = 1e-5;
470 Int_t binX1 = displayHisto->GetXaxis()->FindBin(x1 + 0.01*kShiftX);
471 Int_t binX2 = displayHisto->GetXaxis()->FindBin(x2 - 0.01*kShiftX);
472 Int_t binY1 = displayHisto->GetYaxis()->FindBin(y1 + 0.01*kShiftY);
473 Int_t binY2 = displayHisto->GetYaxis()->FindBin(y2 - 0.01*kShiftY);
475 if(displayOpt==kNumbered) {
476 Int_t meanBin = (binX1+binX2)/2;
480 meanBin = (binY1+binY2)/2;
485 Float_t elementArea = 1.;
486 if(displayOpt == kNormalizeToArea) {
487 elementArea = (x2 - x1 + 2*kShiftX) *
488 (y2 - y1 + 2*kShiftY); // In InitOrDisplayTriggerInfo:
489 // x2 = x_c + xHalfWidth - kShiftX
490 // x1 = x_c - xHalfWidth + kShiftX
491 // so x2 - x1 + 2*kShiftX returns the element width.
495 for(Int_t ibinx=binX1; ibinx<=binX2; ibinx++){
496 for(Int_t ibiny=binY1; ibiny<=binY2; ibiny++){
497 displayHisto->SetBinContent(ibinx,ibiny,binContent/elementArea);