]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONTriggerDisplay.cxx
Go from pointer to ifstream to ifstream.
[u/mrichter/AliRoot.git] / MUON / AliMUONTriggerDisplay.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 "AliMUONTriggerDisplay.h"
20
21 #include "AliMpDDLStore.h"
22 #include "AliMpVSegmentation.h"
23 #include "AliMpSegmentation.h"
24 #include "AliMpConstants.h"
25 #include "AliMpPad.h"
26 #include "AliMpLocalBoard.h"
27 #include "AliMpCDB.h"
28
29 // --- AliRoot header files ---
30 #include "AliLog.h"
31
32 // --- ROOT system ---
33 #include <TH1.h> 
34 #include <TH2.h>
35 #include <TGraph.h>
36
37 //-----------------------------------------------------------------------------
38 /// \class AliMUONTriggerDisplay
39 ///
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 
42 /// trigger chamber.
43 ///
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
48 ///
49 /// \author D. Stocco
50
51 /// \cond CLASSIMP
52 ClassImp(AliMUONTriggerDisplay)
53 /// \endcond
54            
55 //____________________________________________________________________________ 
56 AliMUONTriggerDisplay::AliMUONTriggerDisplay() : 
57 TObject()
58 {
59   /// ctor
60 }
61
62 //__________________________________________________________________
63 AliMUONTriggerDisplay::~AliMUONTriggerDisplay()
64 {
65   /// dtor
66 }
67
68
69 //____________________________________________________________________________ 
70 TH2* AliMUONTriggerDisplay::GetEmptyDisplayHisto(TString displayHistoName, EDisplayType displayType,
71                                                  Int_t cathode, Int_t chamber,
72                                                  TString displayHistoTitle)
73 {
74   //
75   /// Return the display histogram with optimized binning
76   /// but do not fill it.
77   //
78   TH2F* displayHisto = new TH2F();
79   
80   InitOrDisplayTriggerInfo(0x0, displayHisto, displayType,
81                            cathode, chamber, displayHistoName, displayHistoTitle);
82   
83   return displayHisto;
84 }
85
86
87 //____________________________________________________________________________ 
88 TH2* AliMUONTriggerDisplay::GetBoardNumberHisto(TString displayHistoName,
89                                                 Int_t chamber,
90                                                 TString displayHistoTitle)
91 {
92   //
93   /// Return the display histogram with optimized binning
94   /// and fill it with the board number
95   //
96
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);
101   }
102
103   TH2F* displayHisto = (TH2F*)GetEmptyDisplayHisto(displayHistoName, kDisplayBoards, 0, chamber, displayHistoTitle);
104   FillDisplayHistogram(inputHisto,displayHisto,kDisplayBoards,0,chamber,kNumbered);
105   
106   delete inputHisto;
107
108   displayHisto->SetStats(kFALSE);
109   return displayHisto;
110 }
111
112
113 //____________________________________________________________________________ 
114 TH2* AliMUONTriggerDisplay::GetDisplayHistogram(TH1* inputHisto, TString displayHistoName,
115                                                 EDisplayType displayType, Int_t cathode,
116                                                 Int_t chamber, TString displayHistoTitle,
117                                                 EDisplayOption displayOpt)
118 {
119   //
120   /// Get histogram displaying the information contained in the input histogram
121   //
122   TH2* displayHisto = GetEmptyDisplayHisto(displayHistoName, displayType,
123                                            cathode, chamber, displayHistoTitle);
124
125   FillDisplayHistogram(inputHisto, displayHisto, displayType, cathode, chamber, displayOpt);
126
127   return displayHisto;
128 }
129
130 //____________________________________________________________________________ 
131 TH2* AliMUONTriggerDisplay::GetDisplayHistogram(TGraph* inputGraph, TString displayHistoName,
132                                                 EDisplayType displayType, Int_t cathode,
133                                                 Int_t chamber, TString displayHistoTitle,
134                                                 EDisplayOption displayOpt)
135 {
136   //
137   /// Get histogram displaying the information contained in the input graph
138   //
139   TH2* displayHisto = GetEmptyDisplayHisto(displayHistoName, displayType,
140                                            cathode, chamber, displayHistoTitle);
141
142   FillDisplayHistogram(inputGraph, displayHisto, displayType, cathode, chamber, displayOpt);
143
144   return displayHisto;
145 }
146
147 //____________________________________________________________________________ 
148 Bool_t AliMUONTriggerDisplay::FillDisplayHistogram(TH1* inputHisto, TH2* displayHisto,
149                                                    EDisplayType displayType, Int_t cathode,
150                                                    Int_t chamber, EDisplayOption displayOpt)
151 {
152   //
153   /// Fill a previously initialized display histogram 
154   /// with the information contained in inputHisto.
155   /// To get initialized display, please use GetEmptyDisplayHisto method
156   //
157   return InitOrDisplayTriggerInfo(inputHisto, displayHisto, displayType,
158                                   cathode, chamber, "", "",displayOpt);
159 }
160
161 //____________________________________________________________________________ 
162 Bool_t AliMUONTriggerDisplay::FillDisplayHistogram(TGraph* inputGraph, TH2* displayHisto,
163                                                    EDisplayType displayType, Int_t cathode,
164                                                    Int_t chamber, EDisplayOption displayOpt)
165 {
166   //
167   /// Fill a previously initialized display histogram 
168   /// with the information contained in inputGraph.
169   /// To get initialized display, please use GetEmptyDisplayHisto method
170   //
171   return InitOrDisplayTriggerInfo(inputGraph, displayHisto, displayType,
172                                   cathode, chamber, "", "",displayOpt);
173 }
174
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)
181 {
182   //
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.
186   //
187
188   // Load mapping
189   if ( ! AliMpSegmentation::Instance(kFALSE) ) {
190     /// Load mapping
191     if ( ! AliMpCDB::LoadDDLStore() ) {
192       AliError("Could not access mapping from OCDB !");
193       AliError("Histograms are not initialized !");
194       return kFALSE;
195     }
196   }
197   
198   Int_t iCh = (chamber > AliMpConstants::NofTrackingChambers()) ? chamber - 11 : chamber;
199   Int_t iCath = cathode;
200   
201   TArrayD xAxisStrip;
202   TArrayD yAxisStrip;
203   TArrayD xAxisBoard;
204   TArrayD yAxisBoard;
205   
206   Float_t yOffsetLine, xOffsetLine = 0.;
207   
208   const Float_t kResetValue=1234567.;
209   
210   if(!inputObject){
211     xAxisBoard.Set(55);
212     xAxisBoard.Reset(kResetValue);
213     yAxisBoard.Set(50);
214     yAxisBoard.Reset(kResetValue);
215     
216     xAxisStrip.Set(420);
217     xAxisStrip.Reset(kResetValue);
218     yAxisStrip.Set(710);
219     yAxisStrip.Reset(kResetValue);
220   }
221   else if(!displayHisto){
222     AliWarning("Display histogram not initialized. Please initialize it first!");
223     return kFALSE;
224   }
225   else {
226     TH1* inputHisto = dynamic_cast<TH1*>(inputObject);
227     if ( inputHisto ) {
228       if ( inputHisto->GetEntries() == 0 ) {
229         return kTRUE;
230       }
231     }
232     else {
233       TGraph* inputGraph = dynamic_cast<TGraph*>(inputObject);
234       if ( inputGraph ) {
235         if ( inputGraph->GetN() == 0 ){
236           return kTRUE;
237         }
238       }
239       else {
240         AliWarning("The object should inherit from TH1 or TGraph!");
241         return kFALSE;
242       }
243     }
244   }
245   
246   Float_t xWidth, yWidth, yWidthSlat=0., xWidthCol=0.;
247   Float_t x1,x2,y1,y2;
248   Float_t x1b=0., x2b=0., y1b=0., y2b=0.;
249   Float_t xcPad, ycPad;
250   Int_t line=0, slat;
251   Float_t sign = 1.;
252   
253   const Float_t kShiftB = 0.5;
254   const Float_t kShiftS = 0.1;
255   
256   const Float_t kShiftX = (iCath==0) ? kShiftB : kShiftS;
257   const Float_t kShiftY = (iCath==0) ? kShiftS : kShiftB;
258   const Float_t kShiftEl = (displayType==kDisplaySlats) ? 0.01 : kShiftB;
259   
260   Int_t iChamber = iCh + AliMpConstants::NofTrackingChambers();
261   for(Int_t iLoc = 0; iLoc < AliMpConstants::NofLocalBoards(); iLoc++) {
262     Int_t iBoard = iLoc+1;
263     Int_t detElemId = AliMpDDLStore::Instance()->GetDEfromLocalBoard(iBoard, iChamber);
264     
265     if (!detElemId) continue;
266     
267     AliMpLocalBoard* localBoard = AliMpDDLStore::Instance()->GetLocalBoard(iBoard, kFALSE);
268     
269     // skip copy cards
270     if( !localBoard->IsNotified())
271       continue;
272     
273     // get segmentation
274     const AliMpVSegmentation* seg[2] = {
275       AliMpSegmentation::Instance()->GetMpSegmentation(detElemId, AliMp::GetCathodType(0)),
276       AliMpSegmentation::Instance()->GetMpSegmentation(detElemId, AliMp::GetCathodType(1))};
277     
278     if(iLoc==0){
279       AliMpPad pad1 = seg[1]->PadByLocation(iBoard,0,kFALSE);
280       yWidthSlat = pad1.GetDimensionY();
281       AliMpPad pad0 = seg[0]->PadByLocation(iBoard,0,kFALSE);
282       xOffsetLine = TMath::Abs(pad0.GetPositionX()) + pad0.GetDimensionX();
283       xWidthCol = 2.* pad0.GetDimensionX();
284     }
285     
286     // Get ideal global position of DetElemId center
287     slat = detElemId%100;
288     line = (4 + slat)%18;
289     sign = 1.;
290     if(line>8) {
291       line = 17 - line;
292       sign = -1.;
293     }
294     yOffsetLine = (Float_t)(line - 4) * 2. * yWidthSlat;
295           
296     Int_t nLocations = 1;
297     
298     for(Int_t cath=0; cath<AliMpConstants::NofCathodes(); cath++){
299       // Loop on cathodes:
300       // necessary because strip info is read for each cathode
301       // board info is read only from cathode 0
302       
303       // loop over strips
304       for (Int_t ibitxy = 0; ibitxy < 16; ++ibitxy) {
305         // get pad from electronics
306         
307         Int_t offset = 0;
308         if (cath && localBoard->GetSwitch(AliMpLocalBoard::kZeroAllYLSB)) offset = -8;
309         
310         AliMpPad pad = seg[cath]->PadByLocation(iBoard,ibitxy+offset,kFALSE);
311         
312         if (!pad.IsValid()) continue;
313                 
314         xWidth = pad.GetDimensionX();
315         yWidth = pad.GetDimensionY();
316         xcPad = sign * (pad.GetPositionX() + xOffsetLine);
317         if(line==4) xcPad += 0.75 * sign * xWidthCol;
318         ycPad = pad.GetPositionY() + yOffsetLine;
319         nLocations = pad.GetNofLocations();
320         Int_t iStrip = pad.GetLocalBoardChannel(0);
321         
322         if(cath==0){
323           if(iStrip==0) {
324             x1b = xcPad - xWidth + kShiftEl;
325             y1b = ycPad - yWidth + kShiftEl;
326           }
327           x2b = xcPad + xWidth - kShiftEl;
328           y2b = ycPad + yWidth - kShiftEl;
329         }
330         
331         // For non-bending plane fill only the first board covered by the strip
332         // i.e. avoid filling many times the same information
333         if ( cath == 1 && pad.GetLocalBoardId(0) != iBoard ) continue;
334
335         
336         if(cath==iCath){
337           x1 = xcPad - xWidth + kShiftX;
338           y1 = ycPad - yWidth + kShiftY;
339           x2 = xcPad + xWidth - kShiftX;
340           y2 = ycPad + yWidth - kShiftY;
341           
342           if(!inputObject){
343             AddSortedPoint(x1, xAxisStrip, kResetValue);
344             AddSortedPoint(x2, xAxisStrip, kResetValue);
345             
346             AddSortedPoint(y1, yAxisStrip, kResetValue);
347             AddSortedPoint(y2, yAxisStrip, kResetValue);
348           }
349           else if(displayType==kDisplayStrips)
350             FillBins(inputObject, displayHisto, iBoard, iStrip, x1, x2, y1, y2, kShiftX, kShiftY, displayOpt);
351         }
352                 
353       } // loop on strips
354       
355       // if iCath==0 strip info and board info are both filled -> break!
356       // if iCath==1 board info is filled at cath==0. Strip info to be filled at cath==1
357       if(iCath==0) break;
358     } // loop on cathodes
359     
360     if(!inputObject){
361       // Per board
362       AddSortedPoint(x1b, xAxisBoard, kResetValue);
363       AddSortedPoint(x2b, xAxisBoard, kResetValue);
364       
365       AddSortedPoint(y1b, yAxisBoard, kResetValue);
366       AddSortedPoint(y2b, yAxisBoard, kResetValue);
367     }
368     else if(displayType==kDisplayBoards)
369       FillBins(inputObject, displayHisto, iBoard, -nLocations, x1b, x2b, y1b, y2b, kShiftEl, kShiftEl, displayOpt);
370     else if(displayType==kDisplaySlats)
371       FillBins(inputObject, displayHisto, slat, -1, x1b, x2b, y1b, y2b, kShiftEl, kShiftEl, displayOpt);
372   } // loop on local boards
373   
374   if ( inputObject ) return kTRUE;
375   
376   displayHisto->Reset();
377   
378   // Book histos
379   const Float_t kMinDiff = 0.1;
380   
381   TArrayD* currArray[4] = {&xAxisStrip, &yAxisStrip,
382     &xAxisBoard, &yAxisBoard};
383   for(Int_t iaxis=0; iaxis<4; iaxis++){
384     Int_t ipoint=0;
385     while(TMath::Abs((*currArray[iaxis])[ipoint]-kResetValue)>kMinDiff) {
386       ipoint++;
387     }
388     if(ipoint>currArray[iaxis]->GetSize()-2)
389       AliWarning(Form("Array size (%i) lower than the number of points!", currArray[iaxis]->GetSize()));
390     currArray[iaxis]->Set(ipoint);
391   }
392   
393   switch(displayType){
394     case kDisplayStrips:
395       displayHisto->SetBins(xAxisStrip.GetSize()-1, xAxisStrip.GetArray(),
396                             yAxisStrip.GetSize()-1, yAxisStrip.GetArray());
397       break;
398     case kDisplayBoards:
399     case kDisplaySlats:
400       displayHisto->SetBins(xAxisBoard.GetSize()-1, xAxisBoard.GetArray(),
401                             yAxisBoard.GetSize()-1, yAxisBoard.GetArray());
402       break;
403   }
404   
405   displayHisto->SetName(displayHistoName.Data());
406   displayHisto->SetTitle(displayHistoTitle.Data());
407   displayHisto->SetXTitle("X (cm)");
408   displayHisto->SetYTitle("Y (cm)");
409   //displayHisto->SetStats(kFALSE);
410   
411   return kTRUE;
412 }
413
414
415 //____________________________________________________________________________ 
416 Bool_t AliMUONTriggerDisplay::AddSortedPoint(Float_t currVal, TArrayD& position, const Float_t kResetValue)
417 {
418   //
419   /// Add sorted point in array according to an increasing order.
420   /// Used to build display histograms axis.
421   //
422   Int_t nEntries = position.GetSize()-1;
423   Float_t tmp1, tmp2;
424   const Float_t kMinDiff = 0.1;
425   for(Int_t i=0; i<nEntries; i++){
426     if(TMath::Abs(position[i]-currVal)<kMinDiff) return kFALSE;
427     if(TMath::Abs(position[i]-kResetValue)<kMinDiff) {
428       position[i] = currVal;
429       return kTRUE;
430     }
431     if(currVal>position[i]) continue;
432     tmp1 = position[i];
433     position[i] = currVal;
434     for(Int_t j=i+1; j<nEntries; j++){
435       tmp2 = position[j];
436       position[j] = tmp1;
437       tmp1 = tmp2;
438       if(tmp1==kResetValue) break;
439     }
440     return kTRUE;
441   }
442   return kFALSE;
443 }
444
445
446 //____________________________________________________________________________ 
447 void AliMUONTriggerDisplay::FillBins(TObject* inputObject, TH2* displayHisto,
448                                      Int_t iElement1, Int_t iElement2,
449                                      Float_t x1, Float_t x2, Float_t y1, Float_t y2,
450                                      const Float_t kShiftX, const Float_t kShiftY,
451                                      EDisplayOption displayOpt)
452 {
453   //
454   /// Given the bin in inputHisto, search the corresponding bins
455   /// in display histo and fill it.
456   //
457   Int_t binY=0;
458   Float_t binContent=0;
459   TH1* inputHisto = dynamic_cast<TH1*>(inputObject);
460   if ( inputHisto ) {
461     Int_t binX = inputHisto->GetXaxis()->FindBin(iElement1);
462     if ( inputObject->IsA()->InheritsFrom("TH2") ) {
463       binY = inputHisto->GetYaxis()->FindBin(iElement2);
464       binContent = inputHisto->GetBinContent(binX, binY);
465     }
466     else binContent = inputHisto->GetBinContent(binX);
467   }
468   else {
469     TGraph* inputGraph = dynamic_cast<TGraph*>(inputObject);
470     if ( inputGraph ) {
471       Double_t xpt, ypt;
472       for ( Int_t ipt=0; ipt<inputGraph->GetN(); ipt++){
473         inputGraph->GetPoint(ipt, xpt, ypt);
474         if ( TMath::Abs(xpt - iElement1) < 0.1 ) {
475           binContent = ypt;
476           break;
477         }
478       }
479     }
480     else return;
481   }  
482
483   if(binContent==0) {
484     if(displayOpt==kShowZeroes) binContent = 1e-5;
485     else return;
486   }
487
488   Int_t binX1 = displayHisto->GetXaxis()->FindBin(x1 + 0.01*kShiftX);
489   Int_t binX2 = displayHisto->GetXaxis()->FindBin(x2 - 0.01*kShiftX);
490   Int_t binY1 = displayHisto->GetYaxis()->FindBin(y1 + 0.01*kShiftY);
491   Int_t binY2 = displayHisto->GetYaxis()->FindBin(y2 - 0.01*kShiftY);
492
493   if(displayOpt==kNumbered) {
494     Int_t meanBin = (binX1+binX2)/2;
495     binX1 = meanBin;
496     binX2 = meanBin;
497     
498     meanBin = (binY1+binY2)/2;
499     binY1 = meanBin;
500     binY2 = meanBin;
501   }
502
503   Float_t elementArea = 1.;
504   if(displayOpt == kNormalizeToArea) {
505     elementArea = (x2 - x1 + 2*kShiftX) * 
506       (y2 - y1 + 2*kShiftY); // In InitOrDisplayTriggerInfo: 
507                              // x2 = x_c + xHalfWidth - kShiftX
508                              // x1 = x_c - xHalfWidth + kShiftX
509                              // so x2 - x1 + 2*kShiftX returns the element width.
510           
511           // If iElement2 is less than 0, then its meaning is the 
512           // number of boards covered by one strip
513           // the area has therefore to be multiplied accordingly
514           // This fixes the problem when filling the trigger rate per boards in non-bending plane
515           // which is overestimated since the segmentation is always given by the bending-plane
516           if ( iElement2 < 0 ) elementArea *= -(Double_t)iElement2;
517
518   }
519
520   for(Int_t ibinx=binX1; ibinx<=binX2; ibinx++){
521     for(Int_t ibiny=binY1; ibiny<=binY2; ibiny++){
522       displayHisto->SetBinContent(ibinx,ibiny,binContent/elementArea);
523     }
524   }
525 }
526