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