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