]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/AliMUONTriggerEfficiencyCells.cxx
No more misaligned_geometry
[u/mrichter/AliRoot.git] / MUON / AliMUONTriggerEfficiencyCells.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 #include "AliMUONTriggerEfficiencyCells.h"
17 #include "AliMpConstants.h"
18 #include "AliMpDEManager.h"
19
20 #include "AliLog.h"
21
22 #include "TRandom.h"
23 #include "Riostream.h"
24 #include "TSystem.h"
25 #include "TFile.h"
26 #include "TH1F.h"
27 #include "TMatrixF.h"
28
29 #include <fstream>
30 #include <cassert>
31
32 /// \class AliMUONTriggerEfficiencyCells
33 /// A class to store and give access to the trigger chamber efficiency.
34 ///
35 /// Efficiency is stored per cathode on local boards, or, alternatively,
36 /// on "cells" of a given size.
37 ///
38 /// The main method of this class is IsTriggered().
39 ///
40 /// $ALICE_ROOT/MUON/data/TriggerChamberefficiencyCells.dat contains efficiency 
41 /// for each chamber (i.e. DetElement). 
42 ///
43 /// In the case of local boards, efficiency is stored from left to right
44 /// per increasing board number (from 1 to 234)
45 ///
46 /// Otherwise, he efficiency cells goes from right to left and 
47 /// from bottom to top of the chamber, namely, the efficiencies tabulated in the 
48 /// file refers to the following reference frame:
49 ///
50 /// <pre>
51 /// x
52 /// <----------------------------------|
53 ///                                    |
54 ///    ---------------------------     |
55 ///   | 0.97 | 0.97 | 0.97 | 0.97 |    |
56 ///    ---------------------------     |
57 ///   | 0.97 | 0.97 | 0.97 | 0.97 |    |
58 ///    ---------------------------     |
59 ///   | 0.97 | 0.97 | 0.97 | 0.97 |    |
60 ///    ---------------------------     |
61 ///                                    |
62 ///                                   \/ y
63 /// </pre>
64 ///
65 ///  In both cases, the file can be edited in order to change efficiency
66 ///  in a chosen local board/region of the chamber.
67 ///
68 ///
69 /// But please note that this object is also available from the CDB 
70 /// (generated using the MUONCDB.C macro)
71 ///
72 /// \author Diego Stocco; INFN Torino
73
74 /// \cond CLASSIMP
75 ClassImp(AliMUONTriggerEfficiencyCells)
76 /// \endcond
77
78 //__________________________________________________________________________
79 AliMUONTriggerEfficiencyCells::AliMUONTriggerEfficiencyCells()
80 :
81 TObject()
82 {
83 ///  Default constructor.
84   CheckConstants();
85   Reset();
86 }
87
88 //__________________________________________________________________________
89 AliMUONTriggerEfficiencyCells::AliMUONTriggerEfficiencyCells(const char* filename)
90 :
91 TObject()
92 {
93 ///  Constructor using an ASCII file.
94   CheckConstants();
95   Reset();
96   ReadFile(filename);
97 }
98
99
100 //__________________________________________________________________________
101 AliMUONTriggerEfficiencyCells::~AliMUONTriggerEfficiencyCells()
102 {
103 ///  Destructor. Does nothing ;-)
104 }
105
106
107 //__________________________________________________________________________
108 void AliMUONTriggerEfficiencyCells::GetCellEfficiency(Int_t detElemId, Float_t x, Float_t y, Float_t &eff1, Float_t &eff2) const
109 {
110 ///  Get the efficiencies of the 2 cathodes at a given location (x,y)
111
112   Int_t chamber = FindChamberIndex(detElemId);
113   Int_t slat = FindSlatIndex(detElemId);
114   TArrayI cell = CellByCoord(detElemId,x,y);
115   eff1 = 0.0;
116   eff2 = 0.0;
117   if(cell.At(0)>=0 && cell.At(1)>=0)
118   {
119     eff1 = fCellContent[chamber][slat][cell.At(0)][cell.At(1)];
120     eff2 = fCellContent[fgkNchambers+chamber][slat][cell.At(0)][cell.At(1)];
121   }
122 }
123
124
125 //__________________________________________________________________________
126 void AliMUONTriggerEfficiencyCells::GetCellEfficiency(Int_t detElemId, Int_t localBoard, Float_t &eff1, Float_t &eff2) const
127 {
128 ///  Get the efficiencies of the 2 cathodes at a given local board
129
130   Int_t chamber = FindChamberIndex(detElemId);
131   eff1 = fBoardContent[chamber][localBoard-1];
132   eff2 = fBoardContent[fgkNchambers+chamber][localBoard-1];
133 }
134
135
136 //__________________________________________________________________________
137 void 
138 AliMUONTriggerEfficiencyCells::IsTriggered(Int_t detElemId, Float_t x, Float_t y, Bool_t &trig1, Bool_t &trig2) const
139 {
140 ///  Whether or not a given location (x,y) has a chance to trig, on each cathode.
141
142   Float_t eff1 = 0.0;
143   Float_t eff2 = 0.0;
144   GetCellEfficiency(detElemId, x, y, eff1, eff2);
145   trig1 = kTRUE; 
146   trig2 = kTRUE;
147   if(gRandom->Rndm()>eff1)trig1 = kFALSE;
148   if(gRandom->Rndm()>eff2)trig2 = kFALSE;
149 }
150
151
152 //__________________________________________________________________________
153 void 
154 AliMUONTriggerEfficiencyCells::IsTriggered(Int_t detElemId, Int_t localBoard, Bool_t &trig1, Bool_t &trig2) const
155 {
156 ///  Whether or not a given local board has a chance to trig, on each cathode.
157
158   Float_t eff1 = 0.0;
159   Float_t eff2 = 0.0;
160   GetCellEfficiency(detElemId, localBoard, eff1, eff2);
161   trig1 = kTRUE; 
162   trig2 = kTRUE;
163   if(gRandom->Rndm()>eff1)trig1 = kFALSE;
164   if(gRandom->Rndm()>eff2)trig2 = kFALSE;
165 }
166
167 //__________________________________________________________________________
168 TArrayI AliMUONTriggerEfficiencyCells::CellByCoord(Int_t detElemId, Float_t x, Float_t y) const
169 {
170 ///  Get the efficiencies at a given location.
171
172   Int_t chamber = FindChamberIndex(detElemId);
173   Int_t slat = FindSlatIndex(detElemId);
174   Int_t cell[fgkNcathodes]={-1,-1};
175   Float_t maxX = fCellSize[chamber][slat]*((Float_t)fCellNumber[chamber][slat]);
176   Float_t maxY = fCellSize[fgkNchambers+chamber][slat]*((Float_t)fCellNumber[fgkNchambers+chamber][slat]);
177   if(x>=0 & x<maxX & y>=0 & y<maxY)
178   {
179     cell[0] = (Int_t)(x/fCellSize[chamber][slat]);
180     cell[1] = (Int_t)(y/fCellSize[fgkNchambers+chamber][slat]);
181   }
182   return TArrayI(fgkNcathodes,cell);
183 }
184
185 //__________________________________________________________________________
186 void AliMUONTriggerEfficiencyCells::ReadFile(const char* filename)
187 {
188 ///  Reads a file containing the efficiency map.
189
190   TString fileName = gSystem->ExpandPathName(filename);
191   if(fileName.EndsWith(".root")){
192       ReadHistoBoards(fileName.Data());
193       return;
194   }
195
196   ifstream file(fileName.Data());
197   char dat[50];
198   if (file.good()){
199       file >> dat;
200       if(!strcmp(dat,"localBoards"))ReadFileBoards(file);
201       else ReadFileXY(file);
202       file.close();
203   } else {
204       AliWarning(Form("Can't read file %s",fileName.Data()));
205   }
206 }
207
208
209 //__________________________________________________________________________
210 void AliMUONTriggerEfficiencyCells::ReadFileXY(ifstream &file)
211 {
212 ///  Structure of file (.dat) containing geometrical efficency
213     Int_t datInt=0, detEl=0, chamber=0, rpc=0, chCath=0;
214     Float_t datFloat=0.0;
215     char dat[50];
216
217     while (file >> dat) {
218         file >> detEl;
219         chamber = FindChamberIndex(detEl);
220         rpc = FindSlatIndex(detEl);
221         file >> dat;
222         for(Int_t i=0; i<fgkNcathodes; i++){
223             chCath = fgkNchambers*i + chamber;
224             file >> datInt;
225             fCellNumber[chCath][rpc] = datInt;
226             file >> dat;
227         }
228         for(Int_t i=0; i<fgkNcathodes; i++){
229             chCath = fgkNchambers*i + chamber;
230             file >> datFloat;
231             fCellSize[chCath][rpc] = datFloat;
232             if(i==0)file >> dat;
233         }
234         for(Int_t cath=0; cath<fgkNcathodes; cath++){
235             chCath = fgkNchambers*cath + chamber;
236             file >> dat;
237             file >> datInt;
238             for(Int_t iy=0; iy<fCellNumber[fgkNchambers+chamber][rpc]; iy++){
239                 for(Int_t ix=0; ix<fCellNumber[chamber][rpc]; ix++){
240                     file >> datFloat;
241                     fCellContent[chCath][rpc][ix][iy] = datFloat;
242                 }
243             }
244         }
245     }
246 }
247
248 //__________________________________________________________________________
249 void AliMUONTriggerEfficiencyCells::ReadFileBoards(ifstream &file)
250 {
251 ///  Structure of file (.dat) containing local board efficency
252     Int_t datInt=0, detEl=0, chamber=0, chCath=0;
253     Float_t datFloat=0.0;
254     char dat[50];
255
256     while (file >> dat) {
257             file >> detEl;
258             chamber = FindChamberIndex(detEl);
259             //rpc = FindSlatIndex(detEl);
260             for(Int_t cath=0; cath<fgkNcathodes; cath++){
261                 chCath = fgkNchambers*cath + chamber;
262                 file >> dat;
263                 file >> datInt;
264                 for(Int_t board=0; board<AliMpConstants::NofLocalBoards(); board++){
265                     file >> datFloat;
266                     fBoardContent[chCath][board] = datFloat;
267                 }
268             }
269     }
270 }
271
272
273 //__________________________________________________________________________
274 void AliMUONTriggerEfficiencyCells::ReadHistoBoards(const char *filename)
275 {
276 ///  Structure of file (.root) containing local board efficency
277     TFile *file = new TFile(filename, "read");
278     if(!file) {
279         AliWarning(Form("Can't read file %s",filename));
280         return;
281     }
282     char histoName[30];
283     char *cathCode[fgkNcathodes] = {"bendPlane", "nonBendPlane"};
284     TH1F *histo = 0x0;
285     for(Int_t ch=0; ch<fgkNchambers; ch++){
286         for(Int_t cath=0; cath<fgkNcathodes; cath++){
287             sprintf(histoName, "%sBoardEffChamber%i", cathCode[cath], 11+ch);
288             histo = (TH1F *)file->Get(histoName);
289             if(!(TH1F *)file->Get(histoName)) {
290                 AliWarning(Form("Can't find histo %s in file %s",histoName, filename));
291                 continue;
292             }
293             Int_t chCath = fgkNchambers*cath + ch;
294             for(Int_t board=0; board<AliMpConstants::NofLocalBoards(); board++){
295                 Int_t bin = histo->FindBin(board+1);
296                 fBoardContent[chCath][board] = histo->GetBinContent(bin);
297             }
298         }
299     }
300 }
301
302
303 //_____________________________________________________________________________
304 void AliMUONTriggerEfficiencyCells::CheckConstants() const
305 {
306 /// Check consistence of redefined constants 
307
308   assert(fgkNcathodes == AliMpConstants::NofCathodes());    
309   assert(fgkNchambers == AliMpConstants::NofTriggerChambers());    
310   assert(fgkNplanes == AliMpConstants::NofTriggerChambers() * fgkNcathodes);    
311 }
312
313
314 //__________________________________________________________________________
315 Int_t AliMUONTriggerEfficiencyCells::FindChamberIndex(Int_t detElemId) const
316 {
317 ///  From detElemId to chamber number
318   Int_t iChamber = AliMpDEManager::GetChamberId(detElemId);
319   Int_t chamber = iChamber-AliMpConstants::NofTrackingChambers();
320   return chamber;
321 }
322
323
324 //__________________________________________________________________________
325 Int_t AliMUONTriggerEfficiencyCells::FindSlatIndex(Int_t detElemId) const
326 {
327 ///  From detElemId to slat index.
328   Int_t slat = detElemId%100;
329   return slat;
330 }
331
332
333 //__________________________________________________________________________
334 TVector2 AliMUONTriggerEfficiencyCells::ChangeReferenceFrame(Float_t x, Float_t y, Float_t x0, Float_t y0)
335 {
336 /// (x0,y0) position of the local reference frame (center of the chamber)
337
338     Float_t x1 = x0-x;//reflection of axis
339     Float_t y1 = y+y0;
340     return TVector2(x1,y1);
341 }
342
343 //__________________________________________________________________________
344 void
345 AliMUONTriggerEfficiencyCells::Reset()
346 {
347 ///  Sets our internal array contents to zero.
348
349     for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
350         fBoardContent[chCath].Set(AliMpConstants::NofLocalBoards());
351         fCellSize[chCath].Set(fgkNslats);
352         fCellNumber[chCath].Set(fgkNslats);
353         for(Int_t slat=0; slat<fgkNslats; slat++){
354             fCellContent[chCath][slat].ResizeTo(fgkNcells,fgkNcells);
355         }
356     }
357
358     for(Int_t chCath=0; chCath<fgkNplanes; chCath++){
359         fBoardContent[chCath].Reset();
360         fCellSize[chCath].Reset();
361         fCellNumber[chCath].Reset();
362         for(Int_t slat=0; slat<fgkNslats; slat++){
363             fCellContent[chCath][slat].Zero();
364         }
365     }
366 }
367
368