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