Line data Source code
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 : /*
17 : $Log$
18 : Revision 1.20.1 2007/05/19 decaro
19 : Added the following methods:
20 : GetVolumeIndices(Int_t index, Int_t *det), to get
21 : the volume indices (sector, plate, strip, padz, padx,
22 : stored respectively in det[0], det[1], det[2], det[3], det[4])
23 : from the calibration channel index;
24 : NStrip(Int_t nPlate), to get the strips number
25 : per each kind of TOF module.
26 :
27 : Revision 1.20 2007/10/08 17:52:55 decaro
28 : hole region in front of PHOS detector: update of sectors' numbers
29 :
30 : Revision 1.19 2007/10/04 14:05:09 zampolli
31 : AliTOFGeometryV5 becoming AliTOFGeometry
32 :
33 : Revision 1.18 2007/02/19 18:55:26 decaro
34 : Added getter methods for volume path (for Event Display)
35 :
36 : Revision 1.17.1 2006/12/15
37 : Added method DetToStripRF(...) to get
38 : a pad corner coordinates in its strip reference frame
39 : (A.De Caro, M.Di Stefano)
40 : Revision 1.17 2006/08/22 13:30:02 arcelli
41 : removal of effective c++ warnings (C.Zampolli)
42 :
43 : Revision 1.16 2006/04/20 22:30:50 hristov
44 : Coding conventions (Annalisa)
45 :
46 : Revision 1.15 2006/04/16 22:29:05 hristov
47 : Coding conventions (Annalisa)
48 :
49 : Revision 1.14 2006/04/05 08:35:38 hristov
50 : Coding conventions (S.Arcelli, C.Zampolli)
51 :
52 : Revision 1.13 2006/03/12 14:37:54 arcelli
53 : Changes for TOF Reconstruction using TGeo
54 :
55 : Revision 1.12 2006/02/28 10:38:00 decaro
56 : AliTOFGeometry::fAngles, AliTOFGeometry::fHeights, AliTOFGeometry::fDistances arrays: dimension definition in the right location
57 :
58 : Revision 1.11 2005/12/15 14:17:29 decaro
59 : Correction of some parameter values
60 :
61 : Revision 1.10 2005/12/15 08:55:32 decaro
62 : New TOF geometry description (V5) -G. Cara Romeo and A. De Caro
63 :
64 : Revision 1.9.1 2005/07/19 A. De Caro
65 : Created daughter-classes AliTOFGeometryV4 and AliTOFGeometryV5
66 : => moved global methods IsInsideThePad, DistanceToPad,
67 : GetPlate, GetSector, GetStrip, GetPadX, GetPadZ,
68 : GetX, GetY, GetZ, GetPadDx, GetPadDy and GetPadDz
69 : in daughter-classes
70 :
71 : Revision 1.9 2005/10/20 12:41:35 hristov
72 : Implementation of parallel tracking. It is not the default version, one can use it passing option MI from AliReconstruction to TOF (M.Ivanov)
73 :
74 : Revision 1.8 2004/11/29 08:28:01 decaro
75 : Introduction of a new TOF constant (i.e. TDC bin width)
76 :
77 : Revision 1.7 2004/11/05 07:20:08 decaro
78 : TOF library splitting and conversion of some printout messages in AliLog schema (T.Kuhr)
79 :
80 : Revision 1.6 2004/06/15 15:27:59 decaro
81 : TOF raw data: preliminary implementation and style changes
82 :
83 : Revision 1.5 2004/04/20 14:37:22 hristov
84 : Using TMath::Abs instead of fabs, arrays of variable size created/deleted correctly (HP,Sun)
85 :
86 : Revision 1.4 2004/04/13 09:42:51 decaro
87 : Track reconstruction code for TOF: updating
88 :
89 : Revision 1.3 2003/12/29 18:40:39 hristov
90 : Copy/paste error corrected
91 :
92 : Revision 1.2 2003/12/29 17:26:01 hristov
93 : Using enum to initaialize static ints in the header file, the initialization of static floats moved to the implementation file
94 :
95 : Revision 1.1 2003/12/29 15:18:03 decaro
96 : TOF geometry updating (addition of AliTOFGeometry)
97 :
98 : Revision 0.05 2004/6/11 A.De Caro
99 : Implement Global method NpadXStrip
100 : Insert four float constants (originally in AliTOF class)
101 : Revision 0.04 2004/4/05 S.Arcelli
102 : Implement Global methods IsInsideThePad
103 : DistanceToPad
104 : Revision 0.03 2003/12/14 S.Arcelli
105 : Set Phi range [-180,180]->[0,360]
106 : Revision 0.02 2003/12/10 S.Arcelli:
107 : Implement Global methods GetPos & GetDetID
108 : Revision 0.01 2003/12/04 S.Arcelli
109 : */
110 :
111 : ///////////////////////////////////////////////////////////////////////////////
112 : // //
113 : // TOF Geometry class //
114 : // //
115 : ///////////////////////////////////////////////////////////////////////////////
116 :
117 : #include "TGeoManager.h"
118 : //#include "TGeoMatrix.h"
119 : #include "TMath.h"
120 :
121 : #include "AliConst.h"
122 : #include "AliGeomManager.h"
123 : #include "AliLog.h"
124 :
125 : #include "AliTOFGeometry.h"
126 :
127 : extern TGeoManager *gGeoManager;
128 :
129 26 : ClassImp(AliTOFGeometry)
130 :
131 : const Float_t AliTOFGeometry::fgkZlenA = 370.6*2.; // length (cm) of the A module
132 : const Float_t AliTOFGeometry::fgkZlenB = 146.5; // length (cm) of the B module
133 : const Float_t AliTOFGeometry::fgkZlenC = 170.45; // length (cm) of the C module
134 : const Float_t AliTOFGeometry::fgkMaxhZtof = 370.6; // Max half z-size of TOF (cm)
135 :
136 : const Float_t AliTOFGeometry::fgkxTOF = 372.00;// Inner radius of the TOF for Reconstruction (cm)
137 : const Float_t AliTOFGeometry::fgkRmin = 371.00;// Inner radius of the TOF (cm)
138 : const Float_t AliTOFGeometry::fgkRmax = 400.05;// Outer radius of the TOF (cm)
139 :
140 : const Float_t AliTOFGeometry::fgkXPad = 2.5; // Pad size in the x direction (cm)
141 : const Float_t AliTOFGeometry::fgkZPad = 3.5; // Pad size in the z direction (cm)
142 :
143 : const Float_t AliTOFGeometry::fgkStripLength = 122.;// Strip Length (rho X phi direction) (cm)
144 :
145 : const Float_t AliTOFGeometry::fgkSigmaForTail1= 2.; //Sig1 for simulation of TDC tails
146 : const Float_t AliTOFGeometry::fgkSigmaForTail2= 0.5;//Sig2 for simulation of TDC tails
147 :
148 : const Float_t AliTOFGeometry::fgkPhiSec= 20;//sector Phi width (deg)
149 :
150 : Bool_t AliTOFGeometry::fgHoles = 1;//logical for geometry version (w/wo holes)
151 :
152 : const Float_t AliTOFGeometry::fgkTdcBin = 24.4; // time-of-flight bin width [ps]
153 : const Float_t AliTOFGeometry::fgkToTBin = 48.8; // time-over-threshold bin width [ps]
154 : const Float_t AliTOFGeometry::fgkBunchCrossingBin = fgkTdcBin * 1024; // bunch-crossing bin width [ps]
155 :
156 : const Float_t AliTOFGeometry::fgkSlewTOTMin = 10.; // min TOT for slewing correction [ns]
157 : const Float_t AliTOFGeometry::fgkSlewTOTMax = 16.; // max TOT for slewing correction [ns]
158 :
159 : const Float_t AliTOFGeometry::fgkDeadTime = 25E+03; // Single channel dead time (ps)
160 26 : const Float_t AliTOFGeometry::fgkMatchingWindow = fgkTdcBin*TMath::Power(2,13); // Matching window (ps)
161 :
162 : const Float_t AliTOFGeometry::fgkAngles[kNPlates][kMaxNstrip] = {
163 : { 43.99, 43.20, 42.40, 41.59, 40.77, 39.94, 39.11, 38.25, 37.40, 36.53,
164 : 35.65, 34.76, 33.87, 32.96, 32.05, 31.13, 30.19, 29.24, 12.33, 0.00},
165 :
166 : { 27.26, 26.28, 25.30, 24.31, 23.31, 22.31, 21.30, 20.29, 19.26, 18.24,
167 : 17.20, 16.16, 15.11, 14.05, 13.00, 11.93, 10.87, 9.80, 8.74, 0.00},
168 :
169 : { 0.00, 6.30, 5.31, 4.25, 3.19, 2.12, 1.06, 0.00, -1.06, -2.12,
170 : -3.19, -4.25, -5.31, -6.30, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00},
171 :
172 : { -8.74, -9.80, -10.87, -11.93, -13.00, -14.05, -15.11, -16.16, -17.20, -18.24,
173 : -19.26, -20.29, -21.30, -22.31, -23.31, -24.31, -25.30, -26.28, -27.26, 0.00},
174 :
175 : {-12.33, -29.24, -30.19, -31.13, -32.05, -32.96, -33.87, -34.76, -35.65, -36.53,
176 : -37.40, -38.25, -39.11, -39.94, -40.77, -41.59, -42.40, -43.20, -43.99, 0.00}
177 : };
178 :
179 : /*
180 : const Float_t AliTOFGeometry::fgkHeights[kNPlates][kMaxNstrip] = {
181 : {-8.2, -7.5, -8.2, -7.7, -8.1, -7.6, -7.7, -7.7, -7.7, -7.7,
182 : -7.5, -7.2, -7.3, -7.5, -7.6, -7.8, -8.3, -9.3, -3.1, 0.0},
183 :
184 : {-7.9, -8.1, -8.5, -9.0, -10.1, -3.9, -5.9, -7.7, -10.1, -3.6,
185 : -5.8, -8.0, -10.4, -4.4, -7.2, -10.2, -4.6, -7.4, -10.4, 0.0},
186 :
187 : {-2.5, -10.4, -5.0, -9.9, -4.8, -9.9, -4.7, -10.2, -4.7, -9.9,
188 : -4.8, -9.9, -5.0, -10.4, -2.5, 0.0, 0.0, 0.0, 0.0, 0.0},
189 :
190 : {-10.4, -7.4, -4.6, -10.2, -7.2, -4.4, -10.4, -8.0, -5.8, -3.6,
191 : -10.1, -7.7, -5.9, -3.9, -10.1, -9.0, -8.5, -8.1, -7.9, 0.0},
192 :
193 : { -3.1, -9.3, -8.3, -7.8, -7.6, -7.5, -7.3, -7.2, -7.5, -7.7,
194 : -7.7, -7.7, -7.7, -7.6, -8.1, -7.7, -8.2, -7.5, -8.2, 0.0}
195 : };
196 : */
197 : /*
198 : const Float_t AliTOFGeometry::fgkHeights[kNPlates][kMaxNstrip] = {
199 : { -8.405, -10.885, -8.405, -7.765, -8.285, -7.745, -7.865, -7.905, -7.895, -7.885,
200 : -7.705, -7.395, -7.525, -7.645, -11.285, -10.355, -8.365, -9.385, -3.255, 0.000 },
201 : { -7.905, -8.235, -8.605, -9.045, -10.205, -3.975, -5.915, -7.765, -10.205, -3.635,
202 : -5.885, -8.005, -10.505, -4.395, -7.325, -10.235, -4.655, -7.495, -10.515, 0.000 },
203 : { -2.705, -10.645, -5.165, -10.095, -4.995, -10.815, -4.835, -10.385, -4.835, -10.815,
204 : -4.995, -10.095, -5.165, -10.645, -2.705, 0.000, 0.000, 0.000, 0.000, 0.000 },
205 : { -10.515, -7.495, -4.655, -10.235, -7.325, -4.395, -10.505, -8.005, -5.885, -3.635,
206 : -10.205, -7.765, -5.915, -3.975, -10.205, -9.045, -8.605, -8.235, -7.905, 0.000 },
207 : { -3.255, -9.385, -8.365, -10.355, -11.285, -7.645, -7.525, -7.395, -7.705, -7.885,
208 : -7.895, -7.905, -7.865, -7.745, -8.285, -7.765, -8.405, -10.885, -8.405, 0.000 }
209 : };
210 : */
211 :
212 :
213 : const Float_t AliTOFGeometry::fgkHeights[kNPlates][kMaxNstrip] = {
214 : { -8.405, -7.725, -8.405, -7.765, -8.285, -7.745, -7.865, -7.905, -7.895, -7.885,
215 : -7.705, -7.395, -7.525, -7.645, -7.835, -7.965, -8.365, -9.385, -3.255, 0.000 },
216 : { -7.905, -8.235, -8.605, -9.045, -10.205, -3.975, -5.915, -7.765, -10.205, -3.635,
217 : -5.885, -8.005, -10.505, -4.395, -7.325, -10.235, -4.655, -7.495, -10.515, 0.000 },
218 : { -2.705, -10.645, -5.165, -10.095, -4.995, -10.085, -4.835, -10.385, -4.835, -10.085,
219 : -4.995, -10.095, -5.165, -10.645, -2.705, 0.000, 0.000, 0.000, 0.000, 0.000 },
220 : {-10.515, -7.495, -4.655, -10.235, -7.325, -4.395, -10.505, -8.005, -5.885, -3.635,
221 : -10.205, -7.765, -5.915, -3.975, -10.205, -9.045, -8.605, -8.235, -7.905, 0.000 },
222 : { -3.255, -9.385, -8.365, -7.965, -7.835, -7.645, -7.525, -7.395, -7.705, -7.885,
223 : -7.895, -7.905, -7.865, -7.745, -8.285, -7.765, -8.405, -7.725, -8.405, 0.000 }
224 : };
225 :
226 :
227 :
228 : const Float_t AliTOFGeometry::fgkDistances[kNPlates][kMaxNstrip] = {
229 : { 364.14, 354.88, 344.49, 335.31, 325.44, 316.51, 307.11, 297.91, 288.84, 279.89,
230 : 271.20, 262.62, 253.84, 245.20, 236.56, 228.06, 219.46, 210.63, 206.09, 0.00 },
231 : { 194.57, 186.38, 178.25, 170.13, 161.78, 156.62, 148.10, 139.72, 131.23, 125.87,
232 : 117.61, 109.44, 101.29, 95.46, 87.36, 79.37, 73.17, 65.33, 57.71, 0.00 },
233 : { 49.28, 41.35, 35.37, 27.91, 21.20, 13.94, 7.06, 0.00, -7.06, -13.94,
234 : -21.20, -27.91, -35.37, -41.35, -49.28, 0.00, 0.00, 0.00, 0.00, 0.00 },
235 : { -57.71, -65.33, -73.17, -79.37, -87.36, -95.46, -101.29, -109.44, -117.61, -125.87,
236 : -131.23, -139.72, -148.10, -156.62, -161.78, -170.13, -178.25, -186.38, -194.57, 0.00 },
237 : {-206.09, -210.63, -219.46, -228.06, -236.56, -245.20, -253.84, -262.62, -271.20, -279.89,
238 : -288.84, -297.91, -307.11, -316.51, -325.44, -335.31, -344.49, -354.88, -364.14, 0.00 }
239 : };
240 :
241 : /*
242 : const Float_t AliTOFGeometry::fgkDistances[kNPlates][kMaxNstrip] = {
243 : { 364.1, 354.9, 344.5, 335.4, 325.5, 316.6, 307.2, 298.0, 288.9, 280.0,
244 : 271.3, 262.7, 254.0, 244.8, 236.1, 227.7, 219.1, 210.3, 205.7, 0.0},
245 :
246 : { 194.2, 186.1, 177.9, 169.8, 161.5, 156.3, 147.8, 139.4, 130.9, 125.6,
247 : 117.3, 109.2, 101.1, 95.3, 87.1, 79.2, 73.0, 65.1, 57.6, 0.0},
248 :
249 : { 49.5, 41.3, 35.3, 27.8, 21.2, 13.9, 7.0, 0.0, -7.0, -13.9,
250 : -21.2, -27.8, -35.3, -41.3, -49.5, 0.0, 0.0, 0.0, 0.0, 0.0},
251 :
252 : { -57.6, -65.1, -73.0, -79.2, -87.1, -95.3, -101.1, -109.2, -117.3, -125.6,
253 : -130.9, -139.4, -147.8, -156.3, -161.5, -169.8, -177.9, -186.1, -194.2, 0.0},
254 :
255 : {-205.7, -210.3, -219.1, -227.7, -236.1, -244.8, -254.0, -262.7, -271.3, -280.0,
256 : -288.9, -298.0, -307.2, -316.6, -325.5, -335.4, -344.5, -354.9, -364.1, 0.0}
257 : };
258 : */
259 : //_____________________________________________________________________________
260 18 : AliTOFGeometry::AliTOFGeometry()
261 90 : {
262 : //
263 : // AliTOFGeometry default constructor
264 : //
265 :
266 36 : }
267 :
268 : //_____________________________________________________________________________
269 : AliTOFGeometry::~AliTOFGeometry()
270 20 : {
271 : //
272 : // AliTOFGeometry destructor
273 : //
274 20 : }
275 : //_____________________________________________________________________________
276 : void AliTOFGeometry::ImportGeometry(){
277 0 : TGeoManager::Import("geometry.root");
278 0 : }
279 : //_____________________________________________________________________________
280 : void AliTOFGeometry::GetPosPar(Int_t *det, Float_t *pos)
281 : {
282 : //
283 : // Returns space point coor (x,y,z) (cm) for Detector
284 : // Indices (iSect,iPlate,iStrip,iPadX,iPadZ)
285 : //
286 :
287 800 : pos[0]=GetX(det);
288 400 : pos[1]=GetY(det);
289 400 : pos[2]=GetZ(det);
290 :
291 400 : }
292 : //_____________________________________________________________________________
293 : void AliTOFGeometry::GetDetID( Float_t *pos, Int_t *det)
294 : {
295 : //
296 : // Returns Detector Indices (iSect,iPlate,iStrip,iPadX,iPadZ)
297 : // space point coor (x,y,z) (cm)
298 :
299 :
300 0 : det[0]=GetSector(pos);
301 0 : det[1]=GetPlate(pos);
302 0 : det[2]=GetStrip(pos);
303 0 : det[3]=GetPadZ(pos);
304 0 : det[4]=GetPadX(pos);
305 :
306 0 : }
307 : //_____________________________________________________________________________
308 :
309 : void AliTOFGeometry::DetToStripRF(Int_t nPadX, Int_t nPadZ, Float_t &x, Float_t &z) const
310 : {
311 : //
312 : // Returns the local coordinates (x, z) in strip reference frame
313 : // for the bottom corner of the pad number (nPadX, nPadZ)
314 : //
315 : /*
316 : const Float_t xCenterStrip = kNpadX * fgkXPad / 2.;
317 : const Float_t zCenterStrip = kNpadZ * fgkZPad / 2.;
318 :
319 : const Float_t xCenterPad = nPadX*fgkXPad + fgkXPad / 2.;
320 : const Float_t zCenterPad = nPadZ*fgkZPad + fgkZPad / 2.;
321 :
322 : x = xCenterPad - xCenterStrip;
323 : z = zCenterPad - zCenterStrip;
324 : */
325 :
326 :
327 0 : x = (nPadX - kNpadX*0.5) * fgkXPad;
328 0 : z = (nPadZ - kNpadZ*0.5) * fgkZPad;
329 :
330 :
331 0 : }
332 : //_____________________________________________________________________________
333 : Float_t AliTOFGeometry::DistanceToPadPar(Int_t *det, const Float_t * pos, Float_t *dist3d)
334 : {
335 : //
336 : // Returns distance of space point with coor pos (x,y,z) (cm) wrt
337 : // pad with Detector Indices idet (iSect,iPlate,iStrip,iPadX,iPadZ)
338 : //
339 :
340 : //Transform pos into Sector Frame
341 :
342 0 : Float_t x = pos[0];
343 0 : Float_t y = pos[1];
344 0 : Float_t z = pos[2];
345 :
346 0 : Float_t radius = TMath::Sqrt(x*x+y*y);
347 : //Float_t phi=TMath::ATan(y/x);
348 : //if(phi<0) phi = k2PI+phi; //2.*TMath::Pi()+phi;
349 0 : Float_t phi = TMath::Pi()+TMath::ATan2(-y,-x);
350 : // Get the local angle in the sector philoc
351 0 : Float_t angle = phi*kRaddeg-( Int_t (kRaddeg*phi/fgkPhiSec) + 0.5)*fgkPhiSec;
352 0 : Float_t xs = radius*TMath::Cos(angle/kRaddeg);
353 0 : Float_t ys = radius*TMath::Sin(angle/kRaddeg);
354 : Float_t zs = z;
355 :
356 : // Do the same for the selected pad
357 :
358 0 : Float_t g[3];
359 0 : GetPosPar(det,g);
360 :
361 0 : Float_t padRadius = TMath::Sqrt(g[0]*g[0]+g[1]*g[1]);
362 : //Float_t padPhi = TMath::ATan(g[1]/g[0]);
363 : //if(padPhi<0) padPhi = k2Pi + padPhi;
364 0 : Float_t padPhi = TMath::Pi()+TMath::ATan2(-g[1],-g[0]);
365 :
366 : // Get the local angle in the sector philoc
367 0 : Float_t padAngle = padPhi*kRaddeg-( Int_t (padPhi*kRaddeg/fgkPhiSec)+ 0.5) * fgkPhiSec;
368 0 : Float_t padxs = padRadius*TMath::Cos(padAngle/kRaddeg);
369 0 : Float_t padys = padRadius*TMath::Sin(padAngle/kRaddeg);
370 0 : Float_t padzs = g[2];
371 :
372 : //Now move to local pad coordinate frame. Translate:
373 :
374 0 : Float_t xt = xs-padxs;
375 0 : Float_t yt = ys-padys;
376 0 : Float_t zt = zs-padzs;
377 : //Now Rotate:
378 :
379 0 : Float_t alpha = GetAngles(det[1],det[2]);
380 0 : Float_t xr = xt*TMath::Cos(alpha/kRaddeg)+zt*TMath::Sin(alpha/kRaddeg);
381 : Float_t yr = yt;
382 0 : Float_t zr = -xt*TMath::Sin(alpha/kRaddeg)+zt*TMath::Cos(alpha/kRaddeg);
383 :
384 0 : Float_t dist = TMath::Sqrt(xr*xr+yr*yr+zr*zr);
385 :
386 0 : if (dist3d){
387 0 : dist3d[0] = xr;
388 0 : dist3d[1] = yr;
389 0 : dist3d[2] = zr;
390 0 : }
391 :
392 0 : return dist;
393 :
394 0 : }
395 : //_____________________________________________________________________________
396 : Bool_t AliTOFGeometry::IsInsideThePadPar(Int_t *det, const Float_t * pos)
397 : {
398 : //
399 : // Returns true if space point with coor pos (x,y,z) (cm) falls
400 : // inside pad with Detector Indices idet (iSect,iPlate,iStrip,iPadX,iPadZ)
401 : //
402 :
403 : Bool_t isInside=false;
404 :
405 : /*
406 : const Float_t khhony = 1.0 ; // heigth of HONY Layer
407 : const Float_t khpcby = 0.08 ; // heigth of PCB Layer
408 : const Float_t khrgly = 0.055 ; // heigth of RED GLASS Layer
409 : const Float_t khglfy = 0.285 ; // heigth of GLASS+FISHLINE Layer
410 : const Float_t khcpcby = 0.16 ; // heigth of PCB Central Layer
411 : //const Float_t kwcpcbz = 12.4 ; // z dimension of PCB Central Layer
412 : const Float_t khstripy = 2.*khhony+2.*khpcby+4.*khrgly+2.*khglfy+khcpcby;//3.11
413 : //const Float_t kwstripz = kwcpcbz;
414 : //const Float_t klstripx = fgkStripLength;
415 : */
416 :
417 : const Float_t kPadDepth = 0.5;//0.05;//0.11;//0.16;// // heigth of Sensitive Layer
418 :
419 : //Transform pos into Sector Frame
420 :
421 0 : Float_t x = pos[0];
422 0 : Float_t y = pos[1];
423 0 : Float_t z = pos[2];
424 :
425 0 : Float_t radius = TMath::Sqrt(x*x+y*y);
426 0 : Float_t phi = TMath::Pi()+TMath::ATan2(-y,-x);
427 :
428 : // Get the local angle in the sector philoc
429 0 : Float_t angle = phi*kRaddeg-( Int_t (kRaddeg*phi/fgkPhiSec) + 0.5) *fgkPhiSec;
430 0 : Float_t xs = radius*TMath::Cos(angle/kRaddeg);
431 0 : Float_t ys = radius*TMath::Sin(angle/kRaddeg);
432 : Float_t zs = z;
433 :
434 : // Do the same for the selected pad
435 :
436 0 : Float_t g[3];
437 0 : GetPosPar(det,g);
438 :
439 0 : Float_t padRadius = TMath::Sqrt(g[0]*g[0]+g[1]*g[1]);
440 0 : Float_t padPhi = TMath::Pi()+TMath::ATan2(-g[1],-g[0]);
441 :
442 : // Get the local angle in the sector philoc
443 0 : Float_t padAngle = padPhi*kRaddeg-( Int_t (padPhi*kRaddeg/fgkPhiSec)+ 0.5) * fgkPhiSec;
444 0 : Float_t padxs = padRadius*TMath::Cos(padAngle/kRaddeg);
445 0 : Float_t padys = padRadius*TMath::Sin(padAngle/kRaddeg);
446 0 : Float_t padzs = g[2];
447 :
448 : //Now move to local pad coordinate frame. Translate:
449 :
450 0 : Float_t xt = xs-padxs;
451 0 : Float_t yt = ys-padys;
452 0 : Float_t zt = zs-padzs;
453 :
454 : //Now Rotate:
455 :
456 0 : Float_t alpha = GetAngles(det[1],det[2]);
457 0 : Float_t xr = xt*TMath::Cos(alpha/kRaddeg)+zt*TMath::Sin(alpha/kRaddeg);
458 : Float_t yr = yt;
459 0 : Float_t zr = -xt*TMath::Sin(alpha/kRaddeg)+zt*TMath::Cos(alpha/kRaddeg);
460 :
461 0 : if(TMath::Abs(xr)<=kPadDepth*0.5 && TMath::Abs(yr)<= (fgkXPad*0.5) && TMath::Abs(zr)<= (fgkZPad*0.5))
462 0 : isInside=true;
463 0 : return isInside;
464 :
465 0 : }
466 : //_____________________________________________________________________________
467 : Bool_t AliTOFGeometry::IsInsideThePad(TGeoHMatrix *mat, const Float_t * pos, Float_t *dist3d)
468 : {
469 : //
470 : // Returns true if space point with coor pos (x,y,z) [cm] falls inside
471 : // pad identified by the matrix mat. In case dist3d!=0, dist3d vector
472 : // has been filled with the 3D distance between the impact point on
473 : // the pad and the pad centre (in the reference frame of the TOF pad
474 : // identified by the matrix mat).
475 : //
476 :
477 : const Float_t kPadDepth = 0.5; // heigth of Sensitive Layer
478 :
479 129120 : Double_t posg[3];
480 64560 : posg[0] = pos[0];
481 64560 : posg[1] = pos[1];
482 64560 : posg[2] = pos[2];
483 :
484 : // from ALICE global reference system
485 : // towards TOF pad reference system
486 64560 : Double_t posl[3] = {0., 0., 0.};
487 64560 : mat->MasterToLocal(posg,posl);
488 :
489 64560 : Float_t xr = posl[0];
490 64560 : Float_t yr = posl[1];
491 64560 : Float_t zr = posl[2];
492 :
493 : Bool_t isInside = false;
494 65293 : if (TMath::Abs(yr)<= kPadDepth*0.5 &&
495 1534 : TMath::Abs(xr)<= fgkXPad*0.5 &&
496 733 : TMath::Abs(zr)<= fgkZPad*0.5)
497 587 : isInside = true;
498 :
499 64560 : if (dist3d) {
500 : //Double_t padl[3] = {0., 0., 0.};
501 64560 : dist3d[0] = posl[0]/* - padl[0]*/;
502 64560 : dist3d[1] = posl[1]/* - padl[1]*/;
503 64560 : dist3d[2] = posl[2]/* - padl[2]*/;
504 :
505 : /*
506 : Double_t padg[3] = {0., 0., 0.};
507 : // from TOF pad local reference system
508 : // towards ALICE global reference system
509 : TGeoHMatrix inverse = mat->Inverse();
510 : inverse.MasterToLocal(padl,padg);
511 :
512 : // returns the 3d distance
513 : // between the impact point on the pad
514 : // and the pad centre (in the ALICE global reference frame)
515 : dist3d[0] = posg[0] - padg[0];
516 : dist3d[1] = posg[1] - padg[1];
517 : dist3d[2] = posg[2] - padg[2];
518 : */
519 64560 : }
520 :
521 129120 : return isInside;
522 :
523 64560 : }
524 : //_____________________________________________________________________________
525 : void AliTOFGeometry::GetVolumePath(const Int_t * ind, Char_t *path ) {
526 : //--------------------------------------------------------------------
527 : // This function returns the colume path of a given pad
528 : //--------------------------------------------------------------------
529 724 : Int_t sector = ind[0];
530 :
531 : const Int_t kSize = 100;
532 :
533 362 : Char_t string1[kSize];
534 362 : Char_t string2[kSize];
535 362 : Char_t string3[kSize];
536 :
537 : Int_t icopy=-1;
538 : icopy=sector;
539 :
540 362 : snprintf(string1,kSize,"/ALIC_1/B077_1/BSEGMO%i_1/BTOF%i_1",icopy,icopy);
541 :
542 362 : Int_t iplate=ind[1];
543 362 : Int_t istrip=ind[2];
544 365 : if( iplate==0) icopy=istrip;
545 396 : if( iplate==1) icopy=istrip+NStripC();
546 645 : if( iplate==2) icopy=istrip+NStripC()+NStripB();
547 404 : if( iplate==3) icopy=istrip+NStripC()+NStripB()+NStripA();
548 362 : if( iplate==4) icopy=istrip+NStripC()+2*NStripB()+NStripA();
549 362 : icopy++;
550 362 : snprintf(string2,kSize,"FTOA_0/FLTA_0/FSTR_%i",icopy);
551 724 : if(fgHoles && (sector==13 || sector==14 || sector==15)){
552 0 : if(iplate<2) snprintf(string2,kSize,"FTOB_0/FLTB_0/FSTR_%i",icopy);
553 0 : if(iplate>2) snprintf(string2,kSize,"FTOC_0/FLTC_0/FSTR_%i",icopy);
554 : }
555 :
556 362 : Int_t padz = ind[3]+1;
557 362 : Int_t padx = ind[4]+1;
558 362 : snprintf(string3,kSize,"FPCB_1/FSEN_1/FSEZ_%i/FPAD_%i",padz,padx);
559 362 : snprintf(path,2*kSize,"%s/%s/%s",string1,string2,string3);
560 :
561 362 : }
562 : //_____________________________________________________________________________
563 : void AliTOFGeometry::GetVolumePath(Int_t sector, Char_t *path ){
564 : //--------------------------------------------------------------------
565 : // This function returns the colume path of a given sector
566 : //--------------------------------------------------------------------
567 :
568 : const Int_t kSize = 100;
569 :
570 0 : Char_t string[kSize];
571 :
572 : Int_t icopy = sector;
573 :
574 0 : snprintf(string,kSize,"/ALIC_1/B077_1/BSEGMO%i_1/BTOF%i_1",icopy,icopy);
575 0 : snprintf(path,2*kSize,"%s",string);
576 :
577 0 : }
578 : //_____________________________________________________________________________
579 : void AliTOFGeometry::GetVolumePath(Int_t sector, Int_t plate, Int_t strip, Char_t *path ) {
580 : //--------------------------------------------------------------------
581 : // This function returns the colume path of a given strip
582 : //--------------------------------------------------------------------
583 :
584 : const Int_t kSize = 100;
585 :
586 0 : Char_t string1[kSize];
587 0 : Char_t string2[kSize];
588 0 : Char_t string3[kSize];
589 :
590 : Int_t icopy = sector;
591 :
592 0 : snprintf(string1,kSize,"/ALIC_1/B077_1/BSEGMO%i_1/BTOF%i_1",icopy,icopy);
593 :
594 0 : if(plate==0) icopy=strip;
595 0 : if(plate==1) icopy=strip+NStripC();
596 0 : if(plate==2) icopy=strip+NStripC()+NStripB();
597 0 : if(plate==3) icopy=strip+NStripC()+NStripB()+NStripA();
598 0 : if(plate==4) icopy=strip+NStripC()+2*NStripB()+NStripA();
599 0 : icopy++;
600 0 : snprintf(string2,kSize,"FTOA_0/FLTA_0/FSTR_%i",icopy);
601 0 : if(fgHoles && (sector==13 || sector==14 || sector==15)){
602 0 : if(plate<2) snprintf(string2,kSize,"FTOB_0/FLTB_0/FSTR_%i",icopy);
603 0 : if(plate>2) snprintf(string2,kSize,"FTOC_0/FLTC_0/FSTR_%i",icopy);
604 : }
605 :
606 0 : snprintf(string3,kSize,"FPCB_1/FSEN_1");
607 0 : snprintf(path,2*kSize,"%s/%s/%s",string1,string2,string3);
608 :
609 0 : }
610 : //_____________________________________________________________________________
611 : void AliTOFGeometry::GetPos(Int_t *det, Float_t *pos)
612 : {
613 : //
614 : // Returns space point coor (x,y,z) (cm) for Detector
615 : // Indices (iSect,iPlate,iStrip,iPadX,iPadZ)
616 : //
617 0 : Char_t path[200];
618 0 : GetVolumePath(det,path);
619 0 : if (!gGeoManager) {
620 0 : printf("ERROR: no TGeo\n");
621 0 : }
622 0 : gGeoManager->cd(path);
623 0 : TGeoHMatrix global;
624 0 : global = *gGeoManager->GetCurrentMatrix();
625 0 : const Double_t *tr = global.GetTranslation();
626 :
627 0 : pos[0]=tr[0];
628 0 : pos[1]=tr[1];
629 0 : pos[2]=tr[2];
630 0 : }
631 : //_____________________________________________________________________________
632 : Int_t AliTOFGeometry::GetPlate(const Float_t * pos)
633 : {
634 : //
635 : // Returns the Plate index
636 : //
637 : const Float_t kInterCentrModBorder1 = 49.5;
638 : const Float_t kInterCentrModBorder2 = 57.5;
639 : const Float_t kExterInterModBorder1 = 196.0;
640 : const Float_t kExterInterModBorder2 = 203.5;
641 :
642 : const Float_t kLengthExInModBorder = 4.7;
643 : const Float_t kLengthInCeModBorder = 7.0;
644 :
645 : //const Float_t khAlWall = 0.1;
646 : const Float_t kModuleWallThickness = 0.3;
647 : //const Float_t kHoneycombLayerThickness = 1.5;
648 :
649 : Int_t iPlate=-1;
650 :
651 0 : Float_t posLocal[3];
652 0 : for (Int_t ii=0; ii<3; ii++) posLocal[ii] = pos[ii];
653 :
654 0 : Int_t isector = GetSector(posLocal);
655 0 : if(isector == -1){
656 : //AliError("Detector Index could not be determined");
657 0 : return iPlate;
658 : }
659 :
660 : // ALICE reference frame -> B071/B074/B075 = BTO1/2/3 reference frame
661 0 : Double_t angles[6] =
662 0 : {90., 90.+(isector+0.5)*fgkPhiSec,
663 : 0., 0.,
664 : 90., (isector+0.5)*fgkPhiSec
665 : };
666 0 : Rotation(posLocal,angles);
667 :
668 0 : Float_t step[3] = {0., 0., static_cast<Float_t>((fgkRmax+fgkRmin)*0.5)};
669 0 : Translation(posLocal,step);
670 :
671 : // B071/B074/B075 = BTO1/2/3 reference frame -> FTOA = FLTA reference frame
672 0 : angles[0] = 90.;
673 0 : angles[1] = 0.;
674 0 : angles[2] = 0.;
675 0 : angles[3] = 0.;
676 0 : angles[4] = 90.;
677 0 : angles[5] =270.;
678 :
679 0 : Rotation(posLocal,angles);
680 :
681 0 : Float_t yLocal = posLocal[1];
682 0 : Float_t zLocal = posLocal[2];
683 :
684 0 : Float_t deltaRhoLoc = (fgkRmax-fgkRmin)*0.5 - kModuleWallThickness + yLocal;
685 0 : Float_t deltaZetaLoc = TMath::Abs(zLocal);
686 :
687 : Float_t deltaRHOmax = 0.;
688 :
689 0 : if (TMath::Abs(zLocal)>=kExterInterModBorder1 && TMath::Abs(zLocal)<=kExterInterModBorder2)
690 : {
691 0 : deltaRhoLoc -= kLengthExInModBorder;
692 0 : deltaZetaLoc = kExterInterModBorder2-deltaZetaLoc;
693 : deltaRHOmax = (fgkRmax - fgkRmin)*0.5 - kModuleWallThickness - 2.*kLengthExInModBorder; // old 5.35, new 4.8
694 :
695 0 : if (deltaRhoLoc > deltaZetaLoc*deltaRHOmax/(kInterCentrModBorder2-kInterCentrModBorder1)) {
696 0 : if (zLocal<0) iPlate = 0;
697 : else iPlate = 4;
698 : }
699 : else {
700 0 : if (zLocal<0) iPlate = 1;
701 : else iPlate = 3;
702 : }
703 : }
704 0 : else if (TMath::Abs(zLocal)>=kInterCentrModBorder1 && TMath::Abs(zLocal)<=kInterCentrModBorder2)
705 : {
706 0 : deltaRhoLoc -= kLengthInCeModBorder;
707 0 : deltaZetaLoc = deltaZetaLoc-kInterCentrModBorder1;
708 : deltaRHOmax = (fgkRmax - fgkRmin)*0.5 - kModuleWallThickness - 2.*kLengthInCeModBorder; // old 0.39, new 0.2
709 :
710 0 : if (deltaRhoLoc>deltaZetaLoc*deltaRHOmax/(kInterCentrModBorder2-kInterCentrModBorder1)) iPlate = 2;
711 : else {
712 0 : if (zLocal<0) iPlate = 1;
713 : else iPlate = 3;
714 : }
715 : }
716 :
717 0 : if (zLocal>-fgkZlenA*0.5 && zLocal<-kExterInterModBorder2) iPlate = 0;
718 0 : else if (zLocal>-kExterInterModBorder1 && zLocal<-kInterCentrModBorder2) iPlate = 1;
719 0 : else if (zLocal>-kInterCentrModBorder1 && zLocal< kInterCentrModBorder1) iPlate = 2;
720 0 : else if (zLocal> kInterCentrModBorder2 && zLocal< kExterInterModBorder1) iPlate = 3;
721 0 : else if (zLocal> kExterInterModBorder2 && zLocal< fgkZlenA*0.5) iPlate = 4;
722 :
723 : return iPlate;
724 :
725 0 : }
726 :
727 : //_____________________________________________________________________________
728 : Int_t AliTOFGeometry::GetSector(const Float_t * pos)
729 : {
730 : //
731 : // Returns the Sector index
732 : //
733 :
734 : Int_t iSect = -1;
735 :
736 0 : Float_t x = pos[0];
737 0 : Float_t y = pos[1];
738 0 : Float_t z = pos[2];
739 :
740 0 : Float_t rho = TMath::Sqrt(x*x + y*y);
741 :
742 0 : if (!((z>=-fgkZlenA*0.5 && z<=fgkZlenA*0.5) &&
743 0 : (rho>=(fgkRmin) && rho<=(fgkRmax)))) {
744 : //AliError("Detector Index could not be determined");
745 0 : return iSect;
746 : }
747 :
748 0 : Float_t phi = TMath::Pi() + TMath::ATan2(-y,-x);
749 :
750 0 : iSect = (Int_t) (phi*kRaddeg/fgkPhiSec);
751 :
752 : return iSect;
753 :
754 0 : }
755 : //_____________________________________________________________________________
756 : Int_t AliTOFGeometry::GetStrip(const Float_t * pos)
757 : {
758 : //
759 : // Returns the Strip index
760 : //
761 : const Float_t khhony = 1.0 ; // heigth of HONY Layer
762 : const Float_t khpcby = 0.08 ; // heigth of PCB Layer
763 : const Float_t khrgly = 0.055 ; // heigth of RED GLASS Layer
764 : const Float_t khglfy = 0.285 ; // heigth of GLASS+FISHLINE Layer
765 : const Float_t khcpcby = 0.16 ; // heigth of PCB Central Layer
766 : const Float_t kwcpcbz = 12.4 ; // z dimension of PCB Central Layer
767 : const Float_t khstripy = 2.*khhony+2.*khpcby+4.*khrgly+2.*khglfy+khcpcby;//3.11
768 : const Float_t kwstripz = kwcpcbz;
769 : const Float_t klstripx = fgkStripLength;
770 :
771 : Int_t iStrip=-1;
772 :
773 0 : Float_t posLocal[3];
774 0 : for (Int_t ii=0; ii<3; ii++) posLocal[ii] = pos[ii];
775 : // AliDebug(1,Form(" posLocal[0] = %f, posLocal[1] = %f, posLocal[2] = %f ",
776 : // posLocal[0],posLocal[1],posLocal[2]));
777 :
778 0 : Int_t isector = GetSector(posLocal);
779 0 : if(isector == -1){
780 : //AliError("Detector Index could not be determined");
781 0 : return iStrip;}
782 0 : Int_t iplate = GetPlate(posLocal);
783 0 : if(iplate == -1){
784 : //AliError("Detector Index could not be determined");
785 0 : return iStrip;}
786 :
787 : Int_t nstrips=0;
788 0 : switch (iplate) {
789 : case 0:
790 : case 4:
791 : nstrips=kNStripC;
792 0 : break;
793 : case 1:
794 : case 3:
795 : nstrips=kNStripB;
796 0 : break;
797 : case 2:
798 : nstrips=kNStripA;
799 0 : break;
800 : }
801 :
802 : // ALICE reference frame -> B071/B074/B075 = BTO1/2/3 reference frame
803 0 : Double_t angles[6] =
804 0 : {90., 90.+(isector+0.5)*fgkPhiSec,
805 : 0., 0.,
806 : 90., (isector+0.5)*fgkPhiSec
807 : };
808 0 : Rotation(posLocal,angles);
809 : // AliDebug(1,Form(" posLocal[0] = %f, posLocal[1] = %f, posLocal[2] = %f ",
810 : // posLocal[0],posLocal[1],posLocal[2]));
811 :
812 0 : Float_t step[3] = {0., 0., static_cast<Float_t>((fgkRmax+fgkRmin)*0.5)};
813 0 : Translation(posLocal,step);
814 : // AliDebug(1,Form(" posLocal[0] = %f, posLocal[1] = %f, posLocal[2] = %f ",
815 : // posLocal[0],posLocal[1],posLocal[2]));
816 :
817 : // B071/B074/B075 = BTO1/2/3 reference frame -> FTOA = FLTA reference frame
818 0 : angles[0] = 90.;
819 0 : angles[1] = 0.;
820 0 : angles[2] = 0.;
821 0 : angles[3] = 0.;
822 0 : angles[4] = 90.;
823 0 : angles[5] =270.;
824 :
825 0 : Rotation(posLocal,angles);
826 : // AliDebug(1,Form(" posLocal[0] = %f, posLocal[1] = %f, posLocal[2] = %f ",
827 : // posLocal[0],posLocal[1],posLocal[2]));
828 :
829 : // FTOA/B/C = FLTA/B/C reference frame -> FSTR reference frame
830 : Int_t totStrip=0;
831 0 : for (Int_t istrip=0; istrip<nstrips; istrip++){
832 :
833 0 : Float_t posLoc2[3]={posLocal[0],posLocal[1],posLocal[2]};
834 :
835 0 : step[0] = 0.;
836 0 : step[1] = GetHeights(iplate,istrip);
837 0 : step[2] = -GetDistances(iplate,istrip);
838 0 : Translation(posLoc2,step);
839 :
840 0 : if (GetAngles(iplate,istrip) >0.) {
841 0 : angles[0] = 90.;
842 0 : angles[1] = 0.;
843 0 : angles[2] = 90.+GetAngles(iplate,istrip);
844 0 : angles[3] = 90.;
845 0 : angles[4] = GetAngles(iplate,istrip);
846 0 : angles[5] = 90.;
847 0 : }
848 0 : else if (GetAngles(iplate,istrip)==0.) {
849 0 : angles[0] = 90.;
850 0 : angles[1] = 0.;
851 0 : angles[2] = 90.;
852 0 : angles[3] = 90.;
853 0 : angles[4] = 0;
854 0 : angles[5] = 0.;
855 0 : }
856 0 : else if (GetAngles(iplate,istrip) <0.) {
857 0 : angles[0] = 90.;
858 0 : angles[1] = 0.;
859 0 : angles[2] = 90.+GetAngles(iplate,istrip);
860 0 : angles[3] = 90.;
861 0 : angles[4] =-GetAngles(iplate,istrip);
862 0 : angles[5] = 270.;
863 0 : }
864 0 : Rotation(posLoc2,angles);
865 : // AliDebug(1,Form(" strip %2d: posLoc2[0] = %f, posLoc2[1] = %f, posLoc2[2] = %f ",
866 : // istrip, posLoc2[0],posLoc2[1],posLoc2[2]));
867 :
868 0 : if ((TMath::Abs(posLoc2[0])<=klstripx*0.5) &&
869 0 : (TMath::Abs(posLoc2[1])<=khstripy*0.5) &&
870 0 : (TMath::Abs(posLoc2[2])<=kwstripz*0.5)) {
871 : iStrip = istrip;
872 0 : totStrip++;
873 0 : for (Int_t jj=0; jj<3; jj++) posLocal[jj]=posLoc2[jj];
874 : // AliDebug(2,Form(" posLocal[0] = %f, posLocal[1] = %f, posLocal[2] = %f ",
875 : // posLocal[0],posLocal[1],posLocal[2]));
876 :
877 : // AliDebug(2,Form(" GetAngles(%1i,%2i) = %f, pos[0] = %f, pos[1] = %f, pos[2] = %f",
878 : // iplate, istrip, GetAngles(iplate,istrip), pos[0], pos[1], pos[2]));
879 0 : break;
880 : }
881 :
882 : // if (totStrip>1) AliInfo(Form("total strip number found %2i",totStrip));
883 :
884 0 : }
885 :
886 : return iStrip;
887 :
888 0 : }
889 : //_____________________________________________________________________________
890 : Int_t AliTOFGeometry::GetPadZ(const Float_t * pos)
891 : {
892 : //
893 : // Returns the Pad index along Z
894 : //
895 : //const Float_t klsensmx = kNpadX*fgkXPad; // length of Sensitive Layer
896 : //const Float_t khsensmy = 0.05;//0.11;//0.16;// heigth of Sensitive Layer
897 : //const Float_t kwsensmz = kNpadZ*fgkZPad; // width of Sensitive Layer
898 :
899 : Int_t iPadZ = -1;
900 :
901 0 : Float_t posLocal[3];
902 0 : for (Int_t ii=0; ii<3; ii++) posLocal[ii] = pos[ii];
903 :
904 0 : Int_t isector = GetSector(posLocal);
905 0 : if(isector == -1){
906 : //AliError("Detector Index could not be determined");
907 0 : return iPadZ;}
908 0 : Int_t iplate = GetPlate(posLocal);
909 0 : if(iplate == -1){
910 : //AliError("Detector Index could not be determined");
911 0 : return iPadZ;}
912 0 : Int_t istrip = GetStrip(posLocal);
913 0 : if(istrip == -1){
914 : //AliError("Detector Index could not be determined");
915 0 : return iPadZ;}
916 :
917 : // ALICE reference frame -> B071/B074/B075 = BTO1/2/3 reference frame
918 0 : Double_t angles[6] =
919 0 : {90., 90.+(isector+0.5)*fgkPhiSec,
920 : 0., 0.,
921 : 90., (isector+0.5)*fgkPhiSec
922 : };
923 0 : Rotation(posLocal,angles);
924 :
925 0 : Float_t step[3] = {0., 0., static_cast<Float_t>((fgkRmax+fgkRmin)*0.5)};
926 0 : Translation(posLocal,step);
927 :
928 : // B071/B074/B075 = BTO1/2/3 reference frame -> FTOA = FLTA reference frame
929 0 : angles[0] = 90.;
930 0 : angles[1] = 0.;
931 0 : angles[2] = 0.;
932 0 : angles[3] = 0.;
933 0 : angles[4] = 90.;
934 0 : angles[5] =270.;
935 :
936 0 : Rotation(posLocal,angles);
937 :
938 : // FTOA/B/C = FLTA/B/C reference frame -> FSTR reference frame
939 0 : step[0] = 0.;
940 0 : step[1] = GetHeights(iplate,istrip);
941 0 : step[2] = -GetDistances(iplate,istrip);
942 0 : Translation(posLocal,step);
943 :
944 0 : if (GetAngles(iplate,istrip) >0.) {
945 0 : angles[0] = 90.;
946 0 : angles[1] = 0.;
947 0 : angles[2] = 90.+GetAngles(iplate,istrip);
948 0 : angles[3] = 90.;
949 0 : angles[4] = GetAngles(iplate,istrip);
950 0 : angles[5] = 90.;
951 0 : }
952 0 : else if (GetAngles(iplate,istrip)==0.) {
953 0 : angles[0] = 90.;
954 0 : angles[1] = 0.;
955 0 : angles[2] = 90.;
956 0 : angles[3] = 90.;
957 0 : angles[4] = 0;
958 0 : angles[5] = 0.;
959 0 : }
960 0 : else if (GetAngles(iplate,istrip) <0.) {
961 0 : angles[0] = 90.;
962 0 : angles[1] = 0.;
963 0 : angles[2] = 90.+GetAngles(iplate,istrip);
964 0 : angles[3] = 90.;
965 0 : angles[4] =-GetAngles(iplate,istrip);
966 0 : angles[5] = 270.;
967 0 : }
968 0 : Rotation(posLocal,angles);
969 :
970 0 : step[0] =-0.5*kNpadX*fgkXPad;
971 0 : step[1] = 0.;
972 0 : step[2] =-0.5*kNpadZ*fgkZPad;
973 0 : Translation(posLocal,step);
974 :
975 0 : iPadZ = (Int_t)(posLocal[2]/fgkZPad);
976 0 : if (iPadZ==kNpadZ) iPadZ--;
977 0 : else if (iPadZ>kNpadZ) iPadZ=-1;
978 :
979 : return iPadZ;
980 :
981 0 : }
982 : //_____________________________________________________________________________
983 : Int_t AliTOFGeometry::GetPadX(const Float_t * pos)
984 : {
985 : //
986 : // Returns the Pad index along X
987 : //
988 : //const Float_t klsensmx = kNpadX*fgkXPad; // length of Sensitive Layer
989 : //const Float_t khsensmy = 0.05;//0.11;//0.16;// heigth of Sensitive Layer
990 : //const Float_t kwsensmz = kNpadZ*fgkZPad; // width of Sensitive Layer
991 :
992 : Int_t iPadX = -1;
993 :
994 0 : Float_t posLocal[3];
995 0 : for (Int_t ii=0; ii<3; ii++) posLocal[ii] = pos[ii];
996 :
997 0 : Int_t isector = GetSector(posLocal);
998 0 : if(isector == -1){
999 : //AliError("Detector Index could not be determined");
1000 0 : return iPadX;}
1001 0 : Int_t iplate = GetPlate(posLocal);
1002 0 : if(iplate == -1){
1003 : //AliError("Detector Index could not be determined");
1004 0 : return iPadX;}
1005 0 : Int_t istrip = GetStrip(posLocal);
1006 0 : if(istrip == -1){
1007 : //AliError("Detector Index could not be determined");
1008 0 : return iPadX;}
1009 :
1010 : // ALICE reference frame -> B071/B074/B075 = BTO1/2/3 reference frame
1011 0 : Double_t angles[6] =
1012 0 : {90., 90.+(isector+0.5)*fgkPhiSec,
1013 : 0., 0.,
1014 : 90., (isector+0.5)*fgkPhiSec
1015 : };
1016 0 : Rotation(posLocal,angles);
1017 :
1018 0 : Float_t step[3] = {0., 0., static_cast<Float_t>((fgkRmax+fgkRmin)*0.5)};
1019 0 : Translation(posLocal,step);
1020 :
1021 : // B071/B074/B075 = BTO1/2/3 reference frame -> FTOA/B/C = FLTA/B/C reference frame
1022 0 : angles[0] = 90.;
1023 0 : angles[1] = 0.;
1024 0 : angles[2] = 0.;
1025 0 : angles[3] = 0.;
1026 0 : angles[4] = 90.;
1027 0 : angles[5] =270.;
1028 :
1029 0 : Rotation(posLocal,angles);
1030 :
1031 : // FTOA/B/C = FLTA/B/C reference frame -> FSTR reference frame
1032 0 : step[0] = 0.;
1033 0 : step[1] = GetHeights(iplate,istrip);
1034 0 : step[2] = -GetDistances(iplate,istrip);
1035 0 : Translation(posLocal,step);
1036 :
1037 0 : if (GetAngles(iplate,istrip) >0.) {
1038 0 : angles[0] = 90.;
1039 0 : angles[1] = 0.;
1040 0 : angles[2] = 90.+GetAngles(iplate,istrip);
1041 0 : angles[3] = 90.;
1042 0 : angles[4] = GetAngles(iplate,istrip);
1043 0 : angles[5] = 90.;
1044 0 : }
1045 0 : else if (GetAngles(iplate,istrip)==0.) {
1046 0 : angles[0] = 90.;
1047 0 : angles[1] = 0.;
1048 0 : angles[2] = 90.;
1049 0 : angles[3] = 90.;
1050 0 : angles[4] = 0;
1051 0 : angles[5] = 0.;
1052 0 : }
1053 0 : else if (GetAngles(iplate,istrip) <0.) {
1054 0 : angles[0] = 90.;
1055 0 : angles[1] = 0.;
1056 0 : angles[2] = 90.+GetAngles(iplate,istrip);
1057 0 : angles[3] = 90.;
1058 0 : angles[4] =-GetAngles(iplate,istrip);
1059 0 : angles[5] = 270.;
1060 0 : }
1061 0 : Rotation(posLocal,angles);
1062 :
1063 0 : step[0] =-0.5*kNpadX*fgkXPad;
1064 0 : step[1] = 0.;
1065 0 : step[2] =-0.5*kNpadZ*fgkZPad;
1066 0 : Translation(posLocal,step);
1067 :
1068 0 : iPadX = (Int_t)(posLocal[0]/fgkXPad);
1069 0 : if (iPadX==kNpadX) iPadX--;
1070 0 : else if (iPadX>kNpadX) iPadX=-1;
1071 :
1072 : return iPadX;
1073 :
1074 0 : }
1075 : //_____________________________________________________________________________
1076 : Float_t AliTOFGeometry::GetX(const Int_t * det)
1077 : {
1078 : //
1079 : // Returns X coordinate (cm)
1080 : //
1081 :
1082 800 : Int_t isector = det[0];
1083 400 : Int_t iplate = det[1];
1084 400 : Int_t istrip = det[2];
1085 400 : Int_t ipadz = det[3];
1086 400 : Int_t ipadx = det[4];
1087 :
1088 : /*
1089 : // Find out distance d on the plane wrt median phi:
1090 : Float_t d = (ipadx+0.5-kNpadX*0.5)*fgkXPad;
1091 :
1092 : // The radius r in xy plane:
1093 : //Float_t r = (fgkRmin+fgkRmax)*0.5-0.01+GetHeights(iplate,istrip)+
1094 : // (ipadz-0.5)*fgkZPad*TMath::Sin(GetAngles(iplate,istrip)/kRaddeg)-0.25; ???
1095 : Float_t r = (fgkRmin+fgkRmax)*0.5-0.01+GetHeights(iplate,istrip)+
1096 : (ipadz-0.5)*fgkZPad*TMath::Sin(GetAngles(iplate,istrip)/kRaddeg);
1097 :
1098 : // local azimuthal angle in the sector philoc
1099 : Float_t philoc = TMath::ATan(d/r);
1100 : //if(philoc<0.) philoc = k2PI + philoc;
1101 :
1102 : // azimuthal angle in the global frame phi
1103 : Float_t phi = philoc*kRaddeg+(isector+0.5)*fgkPhiSec;
1104 :
1105 : Float_t xCoor = r/TMath::Cos(philoc)*TMath::Cos(phi/kRaddeg);
1106 : */
1107 :
1108 : // Pad reference frame -> FSTR reference frame
1109 400 : Float_t posLocal[3] = {0., 0., 0.};
1110 400 : Float_t step[3] = {static_cast<Float_t>(-(ipadx+0.5)*fgkXPad), 0., static_cast<Float_t>(-(ipadz+0.5)*fgkZPad)};
1111 400 : Translation(posLocal,step);
1112 :
1113 400 : step[0] = kNpadX*0.5*fgkXPad;
1114 400 : step[1] = 0.;
1115 400 : step[2] = kNpadZ*0.5*fgkZPad;
1116 : /*
1117 : Float_t posLocal[3] = {(ipadx+0.5)*fgkXPad, 0., (ipadz+0.5)*fgkZPad};
1118 : Float_t step[3]= {kNpadX*0.5*fgkXPad, 0., kNpadZ*0.5*fgkZPad};
1119 : */
1120 400 : Translation(posLocal,step);
1121 :
1122 : // FSTR reference frame -> FTOA/B/C = FLTA/B/C reference frame
1123 400 : Double_t angles[6] = {0.,0.,0.,0.,0.,0.};
1124 400 : if (GetAngles(iplate,istrip) >0.) {
1125 176 : angles[0] = 90.;
1126 176 : angles[1] = 0.;
1127 176 : angles[2] = 90.+GetAngles(iplate,istrip);
1128 176 : angles[3] = 90.;
1129 176 : angles[4] = GetAngles(iplate,istrip);
1130 176 : angles[5] = 90.;
1131 176 : }
1132 224 : else if (GetAngles(iplate,istrip)==0.) {
1133 132 : angles[0] = 90.;
1134 132 : angles[1] = 0.;
1135 132 : angles[2] = 90.;
1136 132 : angles[3] = 90.;
1137 132 : angles[4] = 0;
1138 132 : angles[5] = 0.;
1139 132 : }
1140 92 : else if (GetAngles(iplate,istrip) <0.) {
1141 92 : angles[0] = 90.;
1142 92 : angles[1] = 0.;
1143 92 : angles[2] = 90.+GetAngles(iplate,istrip);
1144 92 : angles[3] = 90.;
1145 92 : angles[4] =-GetAngles(iplate,istrip);
1146 92 : angles[5] = 270.;
1147 92 : }
1148 :
1149 400 : InverseRotation(posLocal,angles);
1150 :
1151 400 : step[0] = 0.;
1152 400 : step[1] = -GetHeights(iplate,istrip);
1153 400 : step[2] = GetDistances(iplate,istrip);
1154 400 : Translation(posLocal,step);
1155 :
1156 : // FTOA = FLTA reference frame -> B071/B074/B075 = BTO1/2/3 reference frame
1157 400 : angles[0] = 90.;
1158 400 : angles[1] = 0.;
1159 400 : angles[2] = 0.;
1160 400 : angles[3] = 0.;
1161 400 : angles[4] = 90.;
1162 400 : angles[5] =270.;
1163 :
1164 400 : InverseRotation(posLocal,angles);
1165 :
1166 : // B071/B074/B075 = BTO1/2/3 reference frame -> ALICE reference frame
1167 400 : step[0] = 0.;
1168 400 : step[1] = 0.;
1169 400 : step[2] = -((fgkRmax+fgkRmin)*0.5);
1170 400 : Translation(posLocal,step);
1171 :
1172 400 : angles[0] = 90.;
1173 400 : angles[1] = 90.+(isector+0.5)*fgkPhiSec;
1174 400 : angles[2] = 0.;
1175 400 : angles[3] = 0.;
1176 400 : angles[4] = 90.;
1177 400 : angles[5] = (isector+0.5)*fgkPhiSec;
1178 :
1179 400 : InverseRotation(posLocal,angles);
1180 :
1181 400 : Float_t xCoor = posLocal[0];
1182 :
1183 400 : return xCoor;
1184 :
1185 400 : }
1186 : //_____________________________________________________________________________
1187 : Float_t AliTOFGeometry::GetY(const Int_t * det)
1188 : {
1189 : //
1190 : // Returns Y coordinate (cm)
1191 : //
1192 :
1193 800 : Int_t isector = det[0];
1194 400 : Int_t iplate = det[1];
1195 400 : Int_t istrip = det[2];
1196 400 : Int_t ipadz = det[3];
1197 400 : Int_t ipadx = det[4];
1198 :
1199 : /*
1200 : // Find out distance d on the plane wrt median phi:
1201 : Float_t d = (ipadx+0.5-kNpadX*0.5)*fgkXPad;
1202 :
1203 : // The radius r in xy plane:
1204 : //Float_t r = (fgkRmin+fgkRmax)*0.5-0.01+GetHeights(iplate,istrip)+
1205 : // (ipadz-0.5)*fgkZPad*TMath::Sin(GetAngles(iplate,istrip)/kRaddeg)-0.25; ???
1206 : Float_t r = (fgkRmin+fgkRmax)*0.5-0.01+GetHeights(iplate,istrip)+
1207 : (ipadz-0.5)*fgkZPad*TMath::Sin(GetAngles(iplate,istrip)/kRaddeg);
1208 :
1209 : // local azimuthal angle in the sector philoc
1210 : Float_t philoc = TMath::ATan(d/r);
1211 : //if(philoc<0.) philoc = k2PI + philoc;
1212 :
1213 : // azimuthal angle in the global frame phi
1214 : Float_t phi = philoc*kRaddeg+(isector+0.5)*fgkPhiSec;
1215 :
1216 : Float_t yCoor = r/TMath::Cos(philoc)*TMath::Sin(phi/kRaddeg);
1217 : */
1218 :
1219 : // Pad reference frame -> FSTR reference frame
1220 400 : Float_t posLocal[3] = {0., 0., 0.};
1221 400 : Float_t step[3] = {static_cast<Float_t>(-(ipadx+0.5)*fgkXPad), 0., static_cast<Float_t>(-(ipadz+0.5)*fgkZPad)};
1222 400 : Translation(posLocal,step);
1223 :
1224 400 : step[0] = kNpadX*0.5*fgkXPad;
1225 400 : step[1] = 0.;
1226 400 : step[2] = kNpadZ*0.5*fgkZPad;
1227 : /*
1228 : Float_t posLocal[3] = {(ipadx+0.5)*fgkXPad, 0., (ipadz+0.5)*fgkZPad};
1229 : Float_t step[3]= {kNpadX*0.5*fgkXPad, 0., kNpadZ*0.5*fgkZPad};
1230 : */
1231 400 : Translation(posLocal,step);
1232 :
1233 : // FSTR reference frame -> FTOA/B/C = FLTA/B/C reference frame
1234 :
1235 400 : Double_t angles[6] = {0.,0.,0.,0.,0.,0.};
1236 400 : if (GetAngles(iplate,istrip) >0.) {
1237 176 : angles[0] = 90.;
1238 176 : angles[1] = 0.;
1239 176 : angles[2] = 90.+GetAngles(iplate,istrip);
1240 176 : angles[3] = 90.;
1241 176 : angles[4] = GetAngles(iplate,istrip);
1242 176 : angles[5] = 90.;
1243 176 : }
1244 224 : else if (GetAngles(iplate,istrip)==0.) {
1245 132 : angles[0] = 90.;
1246 132 : angles[1] = 0.;
1247 132 : angles[2] = 90.;
1248 132 : angles[3] = 90.;
1249 132 : angles[4] = 0;
1250 132 : angles[5] = 0.;
1251 132 : }
1252 92 : else if (GetAngles(iplate,istrip) <0.) {
1253 92 : angles[0] = 90.;
1254 92 : angles[1] = 0.;
1255 92 : angles[2] = 90.+GetAngles(iplate,istrip);
1256 92 : angles[3] = 90.;
1257 92 : angles[4] =-GetAngles(iplate,istrip);
1258 92 : angles[5] = 270.;
1259 92 : }
1260 :
1261 400 : InverseRotation(posLocal,angles);
1262 :
1263 400 : step[0] = 0.;
1264 400 : step[1] = -GetHeights(iplate,istrip);
1265 400 : step[2] = GetDistances(iplate,istrip);
1266 400 : Translation(posLocal,step);
1267 :
1268 : // FTOA = FLTA reference frame -> B071/B074/B075 = BTO1/2/3 reference frame
1269 400 : angles[0] = 90.;
1270 400 : angles[1] = 0.;
1271 400 : angles[2] = 0.;
1272 400 : angles[3] = 0.;
1273 400 : angles[4] = 90.;
1274 400 : angles[5] =270.;
1275 :
1276 400 : InverseRotation(posLocal,angles);
1277 :
1278 : // B071/B074/B075 = BTO1/2/3 reference frame -> ALICE reference frame
1279 400 : step[0] = 0.;
1280 400 : step[1] = 0.;
1281 400 : step[2] = -((fgkRmax+fgkRmin)*0.5);
1282 400 : Translation(posLocal,step);
1283 :
1284 400 : angles[0] = 90.;
1285 400 : angles[1] = 90.+(isector+0.5)*fgkPhiSec;
1286 400 : angles[2] = 0.;
1287 400 : angles[3] = 0.;
1288 400 : angles[4] = 90.;
1289 400 : angles[5] = (isector+0.5)*fgkPhiSec;
1290 :
1291 400 : InverseRotation(posLocal,angles);
1292 :
1293 400 : Float_t yCoor = posLocal[1];
1294 :
1295 400 : return yCoor;
1296 :
1297 400 : }
1298 :
1299 : //_____________________________________________________________________________
1300 : Float_t AliTOFGeometry::GetZ(const Int_t * det)
1301 : {
1302 : //
1303 : // Returns Z coordinate (cm)
1304 : //
1305 :
1306 800 : Int_t isector = det[0];
1307 400 : Int_t iplate = det[1];
1308 400 : Int_t istrip = det[2];
1309 400 : Int_t ipadz = det[3];
1310 400 : Int_t ipadx = det[4];
1311 :
1312 : /*
1313 : Float_t zCoor = GetDistances(iplate,istrip) +
1314 : (0.5-ipadz) * fgkZPad * TMath::Cos(GetAngles(iplate,istrip)*kDegrad);
1315 : */
1316 :
1317 : // Pad reference frame -> FSTR reference frame
1318 400 : Float_t posLocal[3] = {0., 0., 0.};
1319 400 : Float_t step[3] = {static_cast<Float_t>(-(ipadx+0.5)*fgkXPad), 0., static_cast<Float_t>(-(ipadz+0.5)*fgkZPad)};
1320 400 : Translation(posLocal,step);
1321 :
1322 400 : step[0] = kNpadX*0.5*fgkXPad;
1323 400 : step[1] = 0.;
1324 400 : step[2] = kNpadZ*0.5*fgkZPad;
1325 : /*
1326 : Float_t posLocal[3] = {(ipadx+0.5)*fgkXPad, 0., (ipadz+0.5)*fgkZPad};
1327 : Float_t step[3]= {kNpadX*0.5*fgkXPad, 0., kNpadZ*0.5*fgkZPad};
1328 : */
1329 400 : Translation(posLocal,step);
1330 :
1331 : // FSTR reference frame -> FTOA/B/C = FLTA/B/C reference frame
1332 400 : Double_t angles[6] = {0.,0.,0.,0.,0.,0.};
1333 400 : if (GetAngles(iplate,istrip) >0.) {
1334 176 : angles[0] = 90.;
1335 176 : angles[1] = 0.;
1336 176 : angles[2] = 90.+GetAngles(iplate,istrip);
1337 176 : angles[3] = 90.;
1338 176 : angles[4] = GetAngles(iplate,istrip);
1339 176 : angles[5] = 90.;
1340 176 : }
1341 224 : else if (GetAngles(iplate,istrip)==0.) {
1342 132 : angles[0] = 90.;
1343 132 : angles[1] = 0.;
1344 132 : angles[2] = 90.;
1345 132 : angles[3] = 90.;
1346 132 : angles[4] = 0;
1347 132 : angles[5] = 0.;
1348 132 : }
1349 92 : else if (GetAngles(iplate,istrip) <0.) {
1350 92 : angles[0] = 90.;
1351 92 : angles[1] = 0.;
1352 92 : angles[2] = 90.+GetAngles(iplate,istrip);
1353 92 : angles[3] = 90.;
1354 92 : angles[4] =-GetAngles(iplate,istrip);
1355 92 : angles[5] = 270.;
1356 92 : }
1357 :
1358 400 : InverseRotation(posLocal,angles);
1359 :
1360 400 : step[0] = 0.;
1361 400 : step[1] = -GetHeights(iplate,istrip);
1362 400 : step[2] = GetDistances(iplate,istrip);
1363 400 : Translation(posLocal,step);
1364 :
1365 : // FTOA = FLTA reference frame -> B071/B074/B075 = BTO1/2/3 reference frame
1366 400 : angles[0] = 90.;
1367 400 : angles[1] = 0.;
1368 400 : angles[2] = 0.;
1369 400 : angles[3] = 0.;
1370 400 : angles[4] = 90.;
1371 400 : angles[5] =270.;
1372 :
1373 400 : InverseRotation(posLocal,angles);
1374 :
1375 : // B071/B074/B075 = BTO1/2/3 reference frame -> ALICE reference frame
1376 400 : step[0] = 0.;
1377 400 : step[1] = 0.;
1378 400 : step[2] = -((fgkRmax+fgkRmin)*0.5);
1379 400 : Translation(posLocal,step);
1380 :
1381 400 : angles[0] = 90.;
1382 400 : angles[1] = 90.+(isector+0.5)*fgkPhiSec;
1383 400 : angles[2] = 0.;
1384 400 : angles[3] = 0.;
1385 400 : angles[4] = 90.;
1386 400 : angles[5] = (isector+0.5)*fgkPhiSec;
1387 :
1388 400 : InverseRotation(posLocal,angles);
1389 :
1390 400 : Float_t zCoor = posLocal[2];
1391 :
1392 400 : return zCoor;
1393 :
1394 400 : }
1395 : //_____________________________________________________________________________
1396 :
1397 : void AliTOFGeometry::DetToSectorRF(Int_t vol[5], Double_t coord[4][3])
1398 : {
1399 : //
1400 : // Returns the local coordinates (x, y, z) in sector reference frame
1401 : // for the 4 corners of each sector pad (vol[1], vol[2], vol[3], vol[4])
1402 : //
1403 :
1404 0 : if (!gGeoManager) printf("ERROR: no TGeo\n");
1405 :
1406 : // ALICE -> TOF Sector
1407 0 : Char_t path1[200];
1408 0 : GetVolumePath(vol[0],path1);
1409 0 : gGeoManager->cd(path1);
1410 0 : TGeoHMatrix aliceToSector;
1411 0 : aliceToSector = *gGeoManager->GetCurrentMatrix();
1412 :
1413 : // TOF Sector -> ALICE
1414 : //TGeoHMatrix sectorToALICE = aliceToSector.Inverse();
1415 :
1416 : // ALICE -> TOF Pad
1417 0 : Char_t path2[200];
1418 0 : GetVolumePath(vol,path2);
1419 0 : gGeoManager->cd(path2);
1420 0 : TGeoHMatrix aliceToPad;
1421 0 : aliceToPad = *gGeoManager->GetCurrentMatrix();
1422 :
1423 : // TOF Pad -> ALICE
1424 0 : TGeoHMatrix padToALICE = aliceToPad.Inverse();
1425 :
1426 : // TOF Pad -> TOF Sector
1427 0 : TGeoHMatrix padToSector = padToALICE*aliceToSector;
1428 :
1429 : // TOF Sector -> TOF Pad
1430 : //TGeoHMatrix sectorToPad = sectorToALICE*aliceToPad;
1431 :
1432 : // coordinates of the pad bottom corner
1433 0 : Double_t **cornerPad = new Double_t*[4];
1434 0 : for (Int_t ii=0; ii<4; ii++) cornerPad[ii] = new Double_t[3];
1435 :
1436 0 : cornerPad[0][0] = -fgkXPad/2.;
1437 0 : cornerPad[0][1] = 0.;
1438 0 : cornerPad[0][2] = -fgkZPad/2.;
1439 :
1440 0 : cornerPad[1][0] = fgkXPad/2.;
1441 0 : cornerPad[1][1] = 0.;
1442 0 : cornerPad[1][2] = -fgkZPad/2.;
1443 :
1444 0 : cornerPad[2][0] = fgkXPad/2.;
1445 0 : cornerPad[2][1] = 0.;
1446 0 : cornerPad[2][2] = fgkZPad/2.;
1447 :
1448 0 : cornerPad[3][0] = -fgkXPad/2.;
1449 0 : cornerPad[3][1] = 0.;
1450 0 : cornerPad[3][2] = fgkZPad/2.;
1451 :
1452 0 : for(Int_t aa=0; aa<4; aa++) for(Int_t bb=0; bb<3; bb++) coord[aa][bb]=0.;
1453 :
1454 0 : for (Int_t jj=0; jj<4; jj++) padToSector.MasterToLocal(&cornerPad[jj][0], &coord[jj][0]);
1455 :
1456 0 : delete [] cornerPad;
1457 :
1458 : //sectorToPad.LocalToMaster(cornerPad, coord);
1459 :
1460 0 : }
1461 : //_____________________________________________________________________________
1462 : Float_t AliTOFGeometry::GetPadDx(const Float_t * pos)
1463 : {
1464 : //
1465 : // Returns the x coordinate in the Pad reference frame
1466 : //
1467 :
1468 : Float_t xpad = -2.;
1469 :
1470 0 : Float_t posLocal[3];
1471 0 : for (Int_t ii=0; ii<3; ii++) posLocal[ii] = pos[ii];
1472 :
1473 0 : Int_t isector = GetSector(posLocal);
1474 0 : if(isector == -1){
1475 : //AliError("Detector Index could not be determined");
1476 0 : return xpad;}
1477 0 : Int_t iplate = GetPlate(posLocal);
1478 0 : if(iplate == -1){
1479 : //AliError("Detector Index could not be determined");
1480 0 : return xpad;}
1481 0 : Int_t istrip = GetStrip(posLocal);
1482 0 : if(istrip == -1){
1483 : //AliError("Detector Index could not be determined");
1484 0 : return xpad;}
1485 0 : Int_t ipadz = GetPadZ(posLocal);
1486 0 : if(ipadz == -1){
1487 : //AliError("Detector Index could not be determined");
1488 0 : return xpad;}
1489 0 : Int_t ipadx = GetPadX(posLocal);
1490 0 : if(ipadx == -1){
1491 : //AliError("Detector Index could not be determined");
1492 0 : return xpad;}
1493 :
1494 : // ALICE reference frame -> B071/B074/B075 = BTO1/2/3 reference frame
1495 0 : Double_t angles[6] =
1496 0 : {90., 90.+(isector+0.5)*fgkPhiSec,
1497 : 0., 0.,
1498 : 90., (isector+0.5)*fgkPhiSec
1499 : };
1500 0 : Rotation(posLocal,angles);
1501 :
1502 0 : Float_t step[3] = {0., 0., static_cast<Float_t>((fgkRmax+fgkRmin)*0.5)};
1503 0 : Translation(posLocal,step);
1504 :
1505 : // B071/B074/B075 = BTO1/2/3 reference frame -> FTOA/B/C = FLTA/B/C reference frame
1506 0 : angles[0] = 90.;
1507 0 : angles[1] = 0.;
1508 0 : angles[2] = 0.;
1509 0 : angles[3] = 0.;
1510 0 : angles[4] = 90.;
1511 0 : angles[5] =270.;
1512 :
1513 0 : Rotation(posLocal,angles);
1514 :
1515 : // FTOA/B/C = FLTA/B/C reference frame -> FSTR reference frame
1516 0 : step[0] = 0.;
1517 0 : step[1] = GetHeights(iplate,istrip);
1518 0 : step[2] = -GetDistances(iplate,istrip);
1519 0 : Translation(posLocal,step);
1520 :
1521 0 : if (GetAngles(iplate,istrip) >0.) {
1522 0 : angles[0] = 90.;
1523 0 : angles[1] = 0.;
1524 0 : angles[2] = 90.+GetAngles(iplate,istrip);
1525 0 : angles[3] = 90.;
1526 0 : angles[4] = GetAngles(iplate,istrip);
1527 0 : angles[5] = 90.;
1528 0 : }
1529 0 : else if (GetAngles(iplate,istrip)==0.) {
1530 0 : angles[0] = 90.;
1531 0 : angles[1] = 0.;
1532 0 : angles[2] = 90.;
1533 0 : angles[3] = 90.;
1534 0 : angles[4] = 0;
1535 0 : angles[5] = 0.;
1536 0 : }
1537 0 : else if (GetAngles(iplate,istrip) <0.) {
1538 0 : angles[0] = 90.;
1539 0 : angles[1] = 0.;
1540 0 : angles[2] = 90.+GetAngles(iplate,istrip);
1541 0 : angles[3] = 90.;
1542 0 : angles[4] =-GetAngles(iplate,istrip);
1543 0 : angles[5] = 270.;
1544 0 : }
1545 0 : Rotation(posLocal,angles);
1546 :
1547 0 : step[0] =-0.5*kNpadX*fgkXPad;
1548 0 : step[1] = 0.;
1549 0 : step[2] =-0.5*kNpadZ*fgkZPad;
1550 0 : Translation(posLocal,step);
1551 :
1552 0 : step[0] = (ipadx+0.5)*fgkXPad;
1553 0 : step[1] = 0.;
1554 0 : step[2] = (ipadz+0.5)*fgkZPad;
1555 0 : Translation(posLocal,step);
1556 :
1557 0 : xpad=posLocal[0];
1558 :
1559 : return xpad;
1560 :
1561 0 : }
1562 : //_____________________________________________________________________________
1563 : Float_t AliTOFGeometry::GetPadDy(const Float_t * pos)
1564 : {
1565 : //
1566 : // Returns the y coordinate in the Pad reference frame
1567 : //
1568 :
1569 : Float_t ypad = -2.;
1570 :
1571 0 : Float_t posLocal[3];
1572 0 : for (Int_t ii=0; ii<3; ii++) posLocal[ii] = pos[ii];
1573 :
1574 0 : Int_t isector = GetSector(posLocal);
1575 0 : if(isector == -1){
1576 : //AliError("Detector Index could not be determined");
1577 0 : return ypad;}
1578 0 : Int_t iplate = GetPlate(posLocal);
1579 0 : if(iplate == -1){
1580 : //AliError("Detector Index could not be determined");
1581 0 : return ypad;}
1582 0 : Int_t istrip = GetStrip(posLocal);
1583 0 : if(istrip == -1){
1584 : //AliError("Detector Index could not be determined");
1585 0 : return ypad;}
1586 0 : Int_t ipadz = GetPadZ(posLocal);
1587 0 : if(ipadz == -1){
1588 : //AliError("Detector Index could not be determined");
1589 0 : return ypad;}
1590 0 : Int_t ipadx = GetPadX(posLocal);
1591 0 : if(ipadx == -1){
1592 : //AliError("Detector Index could not be determined");
1593 0 : return ypad;}
1594 :
1595 : // ALICE reference frame -> B071/B074/B075 = BTO1/2/3 reference frame
1596 0 : Double_t angles[6] =
1597 0 : {90., 90.+(isector+0.5)*fgkPhiSec,
1598 : 0., 0.,
1599 : 90., (isector+0.5)*fgkPhiSec
1600 : };
1601 0 : Rotation(posLocal,angles);
1602 :
1603 0 : Float_t step[3] = {0., 0., static_cast<Float_t>((fgkRmax+fgkRmin)*0.5)};
1604 0 : Translation(posLocal,step);
1605 :
1606 : // B071/B074/B075 = BTO1/2/3 reference frame -> FTOA/B/C = FLTA/B/C reference frame
1607 0 : angles[0] = 90.;
1608 0 : angles[1] = 0.;
1609 0 : angles[2] = 0.;
1610 0 : angles[3] = 0.;
1611 0 : angles[4] = 90.;
1612 0 : angles[5] =270.;
1613 :
1614 0 : Rotation(posLocal,angles);
1615 :
1616 : // FTOA/B/C = FLTA/B/C reference frame -> FSTR reference frame
1617 0 : step[0] = 0.;
1618 0 : step[1] = GetHeights(iplate,istrip);
1619 0 : step[2] = -GetDistances(iplate,istrip);
1620 0 : Translation(posLocal,step);
1621 :
1622 0 : if (GetAngles(iplate,istrip) >0.) {
1623 0 : angles[0] = 90.;
1624 0 : angles[1] = 0.;
1625 0 : angles[2] = 90.+GetAngles(iplate,istrip);
1626 0 : angles[3] = 90.;
1627 0 : angles[4] = GetAngles(iplate,istrip);
1628 0 : angles[5] = 90.;
1629 0 : }
1630 0 : else if (GetAngles(iplate,istrip)==0.) {
1631 0 : angles[0] = 90.;
1632 0 : angles[1] = 0.;
1633 0 : angles[2] = 90.;
1634 0 : angles[3] = 90.;
1635 0 : angles[4] = 0;
1636 0 : angles[5] = 0.;
1637 0 : }
1638 0 : else if (GetAngles(iplate,istrip) <0.) {
1639 0 : angles[0] = 90.;
1640 0 : angles[1] = 0.;
1641 0 : angles[2] = 90.+GetAngles(iplate,istrip);
1642 0 : angles[3] = 90.;
1643 0 : angles[4] =-GetAngles(iplate,istrip);
1644 0 : angles[5] = 270.;
1645 0 : }
1646 0 : Rotation(posLocal,angles);
1647 :
1648 0 : step[0] =-0.5*kNpadX*fgkXPad;
1649 0 : step[1] = 0.;
1650 0 : step[2] =-0.5*kNpadZ*fgkZPad;
1651 0 : Translation(posLocal,step);
1652 :
1653 0 : step[0] = (ipadx+0.5)*fgkXPad;
1654 0 : step[1] = 0.;
1655 0 : step[2] = (ipadz+0.5)*fgkZPad;
1656 0 : Translation(posLocal,step);
1657 :
1658 0 : ypad=posLocal[1];
1659 :
1660 : return ypad;
1661 :
1662 0 : }
1663 : //_____________________________________________________________________________
1664 : Float_t AliTOFGeometry::GetPadDz(const Float_t * pos)
1665 : {
1666 : //
1667 : // Returns the z coordinate in the Pad reference frame
1668 : //
1669 :
1670 : Float_t zpad = -2.;
1671 :
1672 0 : Float_t posLocal[3];
1673 0 : for (Int_t ii=0; ii<3; ii++) posLocal[ii] = pos[ii];
1674 :
1675 0 : Int_t isector = GetSector(posLocal);
1676 0 : if(isector == -1){
1677 : //AliError("Detector Index could not be determined");
1678 0 : return zpad;}
1679 0 : Int_t iplate = GetPlate(posLocal);
1680 0 : if(iplate == -1){
1681 : //AliError("Detector Index could not be determined");
1682 0 : return zpad;}
1683 0 : Int_t istrip = GetStrip(posLocal);
1684 0 : if(istrip == -1){
1685 : //AliError("Detector Index could not be determined");
1686 0 : return zpad;}
1687 0 : Int_t ipadz = GetPadZ(posLocal);
1688 0 : if(ipadz == -1){
1689 : //AliError("Detector Index could not be determined");
1690 0 : return zpad;}
1691 0 : Int_t ipadx = GetPadX(posLocal);
1692 0 : if(ipadx == -1){
1693 : //AliError("Detector Index could not be determined");
1694 0 : return zpad;}
1695 :
1696 : // ALICE reference frame -> B071/B074/B075 = BTO1/2/3 reference frame
1697 0 : Double_t angles[6] =
1698 0 : {90., 90.+(isector+0.5)*fgkPhiSec,
1699 : 0., 0.,
1700 : 90., (isector+0.5)*fgkPhiSec
1701 : };
1702 0 : Rotation(posLocal,angles);
1703 :
1704 0 : Float_t step[3] = {0., 0., static_cast<Float_t>((fgkRmax+fgkRmin)*0.5)};
1705 0 : Translation(posLocal,step);
1706 :
1707 : // B071/B074/B075 = BTO1/2/3 reference frame -> FTOA/B/C = FLTA/B/C reference frame
1708 0 : angles[0] = 90.;
1709 0 : angles[1] = 0.;
1710 0 : angles[2] = 0.;
1711 0 : angles[3] = 0.;
1712 0 : angles[4] = 90.;
1713 0 : angles[5] =270.;
1714 :
1715 0 : Rotation(posLocal,angles);
1716 :
1717 : // FTOA/B/C = FLTA/B/C reference frame -> FSTR reference frame
1718 0 : step[0] = 0.;
1719 0 : step[1] = GetHeights(iplate,istrip);
1720 0 : step[2] = -GetDistances(iplate,istrip);
1721 0 : Translation(posLocal,step);
1722 :
1723 0 : if (GetAngles(iplate,istrip) >0.) {
1724 0 : angles[0] = 90.;
1725 0 : angles[1] = 0.;
1726 0 : angles[2] = 90.+GetAngles(iplate,istrip);
1727 0 : angles[3] = 90.;
1728 0 : angles[4] = GetAngles(iplate,istrip);
1729 0 : angles[5] = 90.;
1730 0 : }
1731 0 : else if (GetAngles(iplate,istrip)==0.) {
1732 0 : angles[0] = 90.;
1733 0 : angles[1] = 0.;
1734 0 : angles[2] = 90.;
1735 0 : angles[3] = 90.;
1736 0 : angles[4] = 0;
1737 0 : angles[5] = 0.;
1738 0 : }
1739 0 : else if (GetAngles(iplate,istrip) <0.) {
1740 0 : angles[0] = 90.;
1741 0 : angles[1] = 0.;
1742 0 : angles[2] = 90.+GetAngles(iplate,istrip);
1743 0 : angles[3] = 90.;
1744 0 : angles[4] =-GetAngles(iplate,istrip);
1745 0 : angles[5] = 270.;
1746 0 : }
1747 0 : Rotation(posLocal,angles);
1748 :
1749 0 : step[0] =-0.5*kNpadX*fgkXPad;
1750 0 : step[1] = 0.;
1751 0 : step[2] =-0.5*kNpadZ*fgkZPad;
1752 0 : Translation(posLocal,step);
1753 :
1754 0 : step[0] = (ipadx+0.5)*fgkXPad;
1755 0 : step[1] = 0.;
1756 0 : step[2] = (ipadz+0.5)*fgkZPad;
1757 0 : Translation(posLocal,step);
1758 :
1759 0 : zpad=posLocal[2];
1760 :
1761 : return zpad;
1762 :
1763 0 : }
1764 : //_____________________________________________________________________________
1765 :
1766 : void AliTOFGeometry::Translation(Float_t *xyz, Float_t translationVector[3])
1767 : {
1768 : //
1769 : // Return the vector xyz translated by translationVector vector
1770 : //
1771 :
1772 : Int_t ii=0;
1773 :
1774 43200 : for (ii=0; ii<3; ii++)
1775 14400 : xyz[ii] -= translationVector[ii];
1776 :
1777 : return;
1778 :
1779 4800 : }
1780 : //_____________________________________________________________________________
1781 :
1782 : void AliTOFGeometry::Rotation(Float_t *xyz, Double_t rotationAngles[6])
1783 : {
1784 : //
1785 : // Return the vector xyz rotated according to the rotationAngles angles
1786 : //
1787 :
1788 : Int_t ii=0;
1789 : /*
1790 : TRotMatrix *matrix = new TRotMatrix("matrix","matrix", angles[0], angles[1],
1791 : angles[2], angles[3],
1792 : angles[4], angles[5]);
1793 : */
1794 :
1795 0 : for (ii=0; ii<6; ii++) rotationAngles[ii]*=kDegrad;
1796 :
1797 0 : Float_t xyzDummy[3] = {0., 0., 0.};
1798 :
1799 0 : for (ii=0; ii<3; ii++) {
1800 0 : xyzDummy[ii] =
1801 0 : xyz[0]*TMath::Sin(rotationAngles[2*ii])*TMath::Cos(rotationAngles[2*ii+1]) +
1802 0 : xyz[1]*TMath::Sin(rotationAngles[2*ii])*TMath::Sin(rotationAngles[2*ii+1]) +
1803 0 : xyz[2]*TMath::Cos(rotationAngles[2*ii]);
1804 : }
1805 :
1806 0 : for (ii=0; ii<3; ii++) xyz[ii]=xyzDummy[ii];
1807 :
1808 : return;
1809 :
1810 0 : }
1811 : //_____________________________________________________________________________
1812 : void AliTOFGeometry::InverseRotation(Float_t *xyz, Double_t rotationAngles[6])
1813 : {
1814 : //
1815 : // Rotates the vector xyz acordint to the rotationAngles
1816 : //
1817 :
1818 : Int_t ii=0;
1819 :
1820 54000 : for (ii=0; ii<6; ii++) rotationAngles[ii]*=kDegrad;
1821 :
1822 3600 : Float_t xyzDummy[3] = {0., 0., 0.};
1823 :
1824 3600 : xyzDummy[0] =
1825 10800 : xyz[0]*TMath::Sin(rotationAngles[0])*TMath::Cos(rotationAngles[1]) +
1826 7200 : xyz[1]*TMath::Sin(rotationAngles[2])*TMath::Cos(rotationAngles[3]) +
1827 3600 : xyz[2]*TMath::Sin(rotationAngles[4])*TMath::Cos(rotationAngles[5]);
1828 :
1829 3600 : xyzDummy[1] =
1830 10800 : xyz[0]*TMath::Sin(rotationAngles[0])*TMath::Sin(rotationAngles[1]) +
1831 7200 : xyz[1]*TMath::Sin(rotationAngles[2])*TMath::Sin(rotationAngles[3]) +
1832 3600 : xyz[2]*TMath::Sin(rotationAngles[4])*TMath::Sin(rotationAngles[5]);
1833 :
1834 3600 : xyzDummy[2] =
1835 10800 : xyz[0]*TMath::Cos(rotationAngles[0]) +
1836 7200 : xyz[1]*TMath::Cos(rotationAngles[2]) +
1837 3600 : xyz[2]*TMath::Cos(rotationAngles[4]);
1838 :
1839 28800 : for (ii=0; ii<3; ii++) xyz[ii]=xyzDummy[ii];
1840 :
1841 : return;
1842 :
1843 3600 : }
1844 : //_____________________________________________________________________________
1845 :
1846 : Int_t AliTOFGeometry::GetIndex(const Int_t * detId)
1847 : {
1848 : //Retrieve calibration channel index
1849 1150 : Int_t isector = detId[0];
1850 575 : if (isector >= kNSectors){
1851 0 : printf("Wrong sector number in TOF (%d) !\n",isector);
1852 0 : return -1;
1853 : }
1854 575 : Int_t iplate = detId[1];
1855 575 : if (iplate >= kNPlates){
1856 0 : printf("Wrong plate number in TOF (%d) !\n",iplate);
1857 0 : return -1;
1858 : }
1859 575 : Int_t istrip = detId[2];
1860 575 : Int_t stripOffset = GetStripNumberPerSM(iplate,istrip);
1861 575 : if (stripOffset==-1) {
1862 0 : printf("Wrong strip number per SM in TOF (%d) !\n",stripOffset);
1863 0 : return -1;
1864 : }
1865 :
1866 575 : Int_t ipadz = detId[3];
1867 575 : Int_t ipadx = detId[4];
1868 :
1869 1150 : Int_t idet = ((2*(kNStripC+kNStripB)+kNStripA)*kNpadZ*kNpadX)*isector +
1870 1150 : (stripOffset*kNpadZ*kNpadX)+
1871 1150 : (kNpadX)*ipadz+
1872 : ipadx;
1873 : return idet;
1874 575 : }
1875 : //_____________________________________________________________________________
1876 :
1877 : void AliTOFGeometry::GetVolumeIndices(Int_t index, Int_t *detId)
1878 : {
1879 : //
1880 : // Retrieve volume indices from the calibration channel index
1881 : //
1882 :
1883 0 : detId[0] = index/NpadXStrip()/NStripXSector();
1884 :
1885 : Int_t dummyStripPerModule =
1886 0 : ( index - ( NStripXSector()*NpadXStrip()*detId[0]) ) / NpadXStrip();
1887 0 : if (dummyStripPerModule<kNStripC) {
1888 0 : detId[1] = 0;
1889 0 : detId[2] = dummyStripPerModule;
1890 0 : }
1891 0 : else if (dummyStripPerModule>=kNStripC && dummyStripPerModule<kNStripC+kNStripB) {
1892 0 : detId[1] = 1;
1893 0 : detId[2] = dummyStripPerModule-kNStripC;
1894 0 : }
1895 0 : else if (dummyStripPerModule>=kNStripC+kNStripB && dummyStripPerModule<kNStripC+kNStripB+kNStripA) {
1896 0 : detId[1] = 2;
1897 0 : detId[2] = dummyStripPerModule-kNStripC-kNStripB;
1898 0 : }
1899 0 : else if (dummyStripPerModule>=kNStripC+kNStripB+kNStripA && dummyStripPerModule<kNStripC+kNStripB+kNStripA+kNStripB) {
1900 0 : detId[1] = 3;
1901 0 : detId[2] = dummyStripPerModule-kNStripC-kNStripB-kNStripA;
1902 0 : }
1903 0 : else if (dummyStripPerModule>=kNStripC+kNStripB+kNStripA+kNStripB && dummyStripPerModule<NStripXSector()) {
1904 0 : detId[1] = 4;
1905 0 : detId[2] = dummyStripPerModule-kNStripC-kNStripB-kNStripA-kNStripB;
1906 0 : }
1907 :
1908 0 : Int_t padPerStrip = ( index - ( NStripXSector()*NpadXStrip()*detId[0]) ) - dummyStripPerModule*NpadXStrip();
1909 :
1910 0 : detId[3] = padPerStrip / kNpadX; // padZ
1911 0 : detId[4] = padPerStrip - detId[3]*kNpadX; // padX
1912 :
1913 0 : }
1914 : //_____________________________________________________________________________
1915 :
1916 : Int_t AliTOFGeometry::NStrip(Int_t nPlate)
1917 : {
1918 : //
1919 : // Returns the strips number for the plate number 'nPlate'
1920 : //
1921 :
1922 : Int_t nStrips = kNStripC;
1923 :
1924 0 : switch(nPlate) {
1925 : case 2:
1926 : nStrips = kNStripA;
1927 0 : break;
1928 : case 1:
1929 : case 3:
1930 : nStrips = kNStripB;
1931 0 : break;
1932 : case 0:
1933 : case 4:
1934 : default:
1935 : nStrips = kNStripC;
1936 0 : break;
1937 : }
1938 :
1939 0 : return nStrips;
1940 :
1941 : }
1942 : //-------------------------------------------------------------------------
1943 :
1944 : UShort_t AliTOFGeometry::GetAliSensVolIndex(Int_t isector, Int_t iplate, Int_t istrip)
1945 : {
1946 : //
1947 : // Get the index of the TOF alignable volume in the AliGeomManager order.
1948 : //
1949 :
1950 0 : Int_t index = GetStripNumber(isector, iplate, istrip);
1951 :
1952 0 : UShort_t volIndex = AliGeomManager::LayerToVolUID(AliGeomManager::kTOF,index);
1953 :
1954 0 : return volIndex;
1955 :
1956 : }
1957 : //-------------------------------------------------------------------------
1958 :
1959 : Int_t AliTOFGeometry::GetStripNumber(Int_t isector, Int_t iplate, Int_t istrip)
1960 : {
1961 : //
1962 : // Get the serial number of the TOF strip number istrip [0,14/18],
1963 : // in the module number iplate [0,4],
1964 : // in the TOF SM number isector [0,17].
1965 : // This number will range in [0,1637].
1966 : //
1967 :
1968 0 : Bool_t check = (isector >= kNSectors);
1969 :
1970 0 : if (check)
1971 0 : printf("E-AliTOFGeometry::GetStripNumber: Wrong sector number in TOF (%d)!\n",isector);
1972 :
1973 : Int_t index = -1;
1974 0 : Int_t stripInSM = GetStripNumberPerSM(iplate, istrip);
1975 0 : if (!check && stripInSM!=-1)
1976 0 : index = (2*(kNStripC+kNStripB)+kNStripA)*isector + stripInSM;
1977 :
1978 0 : return index;
1979 :
1980 : }
1981 : //-------------------------------------------------------------------------
1982 :
1983 : void AliTOFGeometry::GetStripAndModule(Int_t iStripPerSM, Int_t &iplate, Int_t &istrip)
1984 : {
1985 : //
1986 : // Convert the serial number of the TOF strip number iStripPerSM [0,90]
1987 : // in module number iplate [0,4] and strip number istrip [0,14/18].
1988 : //
1989 :
1990 0 : if (iStripPerSM<0 || iStripPerSM>=kNStripC+kNStripB+kNStripA+kNStripB+kNStripC) {
1991 0 : iplate = -1;
1992 0 : istrip = -1;
1993 0 : }
1994 0 : else if (iStripPerSM<kNStripC) {
1995 0 : iplate = 0;
1996 0 : istrip = iStripPerSM;
1997 0 : }
1998 0 : else if (iStripPerSM>=kNStripC && iStripPerSM<kNStripC+kNStripB) {
1999 0 : iplate = 1;
2000 0 : istrip = iStripPerSM-kNStripC;
2001 0 : }
2002 0 : else if (iStripPerSM>=kNStripC+kNStripB && iStripPerSM<kNStripC+kNStripB+kNStripA) {
2003 0 : iplate = 2;
2004 0 : istrip = iStripPerSM-kNStripC-kNStripB;
2005 0 : }
2006 0 : else if (iStripPerSM>=kNStripC+kNStripB+kNStripA && iStripPerSM<kNStripC+kNStripB+kNStripA+kNStripB) {
2007 0 : iplate = 3;
2008 0 : istrip = iStripPerSM-kNStripC-kNStripB-kNStripA;
2009 0 : }
2010 0 : else if (iStripPerSM>=kNStripC+kNStripB+kNStripA+kNStripB && iStripPerSM<kNStripC+kNStripB+kNStripA+kNStripB+kNStripC) {
2011 0 : iplate = 4;
2012 0 : istrip = iStripPerSM-kNStripC-kNStripB-kNStripA-kNStripB;
2013 0 : }
2014 :
2015 :
2016 0 : }
2017 : //-------------------------------------------------------------------------
2018 :
2019 : Int_t AliTOFGeometry::GetStripNumberPerSM(Int_t iplate, Int_t istrip)
2020 : {
2021 : //
2022 : // Get the serial number of the TOF strip number istrip [0,14/18],
2023 : // in the module number iplate [0,4].
2024 : // This number will range in [0,90].
2025 : //
2026 :
2027 : Int_t index = -1;
2028 :
2029 575 : Bool_t check = (
2030 1150 : (iplate<0 || iplate>=kNPlates)
2031 : ||
2032 : (
2033 978 : (iplate==2 && (istrip<0 || istrip>=kNStripA))
2034 : ||
2035 747 : (iplate!=2 && (istrip<0 || istrip>=kNStripC))
2036 : )
2037 : );
2038 :
2039 575 : if (iplate<0 || iplate>=kNPlates)
2040 0 : printf("E-AliTOFGeometry::GetStripNumberPerSM: Wrong plate number in TOF (%1d)!\n",iplate);
2041 :
2042 : if (
2043 978 : (iplate==2 && (istrip<0 || istrip>=kNStripA))
2044 : ||
2045 747 : (iplate!=2 && (istrip<0 || istrip>=kNStripC))
2046 : )
2047 0 : printf("E-AliTOFGeometry::GetStripNumberPerSM: Wrong strip number in TOF "
2048 : "(strip=%2d in the plate=%1d)!\n",istrip,iplate);
2049 :
2050 : Int_t stripOffset = 0;
2051 1150 : switch (iplate) {
2052 : case 0:
2053 : stripOffset = 0;
2054 27 : break;
2055 : case 1:
2056 : stripOffset = kNStripC;
2057 79 : break;
2058 : case 2:
2059 : stripOffset = kNStripC+kNStripB;
2060 403 : break;
2061 : case 3:
2062 : stripOffset = kNStripC+kNStripB+kNStripA;
2063 66 : break;
2064 : case 4:
2065 : stripOffset = kNStripC+kNStripB+kNStripA+kNStripB;
2066 0 : break;
2067 : };
2068 :
2069 1150 : if (!check) index = stripOffset + istrip;
2070 :
2071 575 : return index;
2072 :
2073 : }
2074 : //-------------------------------------------------------------------------
2075 :
2076 : void AliTOFGeometry::PadRF2TrackingRF(Float_t *ctrackPos, Float_t *differenceT)
2077 : {
2078 : //
2079 : // To convert the 3D distance ctrackPos, referred to the ALICE RF,
2080 : // into the 3D distance differenceT, referred to the tracking RF
2081 : // in case ctrakPos belongs to a TOF sensitive volume.
2082 : //
2083 :
2084 0 : for (Int_t ii=0; ii<3; ii++) differenceT[ii] = 999.;
2085 :
2086 0 : AliDebug(1,Form(" track position in ALICE global Ref. frame -> %f, %f, %f",
2087 : ctrackPos[0],ctrackPos[1],ctrackPos[2]));
2088 :
2089 0 : Int_t detId[5] = {-1,-1,-1,-1,-1};
2090 :
2091 0 : detId[0] = GetSector(ctrackPos);
2092 0 : if (detId[0]==-1) {
2093 0 : AliWarning(Form("This point does not belong to any TOF sector"));
2094 0 : return;
2095 : }
2096 :
2097 0 : detId[1] = GetPlate(ctrackPos);
2098 0 : if (detId[1]==-1) {
2099 0 : AliWarning(Form("This point does not belong to any TOF module"));
2100 0 : return;
2101 : }
2102 :
2103 0 : detId[2] = GetStrip(ctrackPos);
2104 0 : if (detId[2]==-1) {
2105 0 : AliWarning(Form("This point does not belong to any TOF strip"));
2106 0 : return;
2107 : }
2108 :
2109 0 : detId[3] = GetPadZ(ctrackPos);
2110 0 : if (detId[3]==-1) {
2111 0 : AliWarning(Form("This point does not belong to any TOF pad-row"));
2112 0 : return;
2113 : }
2114 :
2115 0 : detId[4] = GetPadX(ctrackPos);
2116 0 : if (detId[4]==-1) {
2117 0 : AliWarning(Form("This point does not belong to any TOF pad"));
2118 0 : return;
2119 : }
2120 :
2121 :
2122 : UShort_t alignableStripIndex =
2123 0 : GetAliSensVolIndex(detId[0],detId[1],detId[2]);
2124 0 : AliDebug(1,Form(" sector = %2d, plate = %1d, strip = %2d (padZ = %1d, padX = %2d) "
2125 : "---> stripIndex = %4d",
2126 : detId[0], detId[1], detId[2], detId[3], detId[4], alignableStripIndex));
2127 :
2128 : // pad centre coordinates in the strip ref. frame
2129 0 : Double_t padCentreL[3] = {(detId[4]-AliTOFGeometry::NpadX()/2)*AliTOFGeometry::XPad()
2130 0 : +AliTOFGeometry::XPad()/2.,
2131 : 0.,
2132 0 : (detId[3]-AliTOFGeometry::NpadZ()/2)*AliTOFGeometry::XPad()
2133 0 : +AliTOFGeometry::XPad()/2.};
2134 : // pad centre coordinates in the strip tracking frame
2135 0 : Double_t padCentreT[3] = {0., 0., 0.};
2136 0 : TGeoHMatrix l2t = *AliGeomManager::GetTracking2LocalMatrix(alignableStripIndex);
2137 0 : l2t.MasterToLocal(padCentreL,padCentreT);
2138 :
2139 :
2140 0 : Char_t path[200];
2141 : // pad centre coordinates in its ref. frame
2142 0 : Double_t padCentreL2[3] = {0., 0., 0.};
2143 : // pad centre coordinates in the ALICE global ref. frame
2144 0 : Double_t padCentreG[3] = {0., 0., 0.};
2145 0 : GetVolumePath(detId,path);
2146 0 : gGeoManager->cd(path);
2147 0 : TGeoHMatrix g2l = *gGeoManager->GetCurrentMatrix();
2148 0 : TGeoHMatrix l2g = g2l.Inverse();
2149 0 : l2g.MasterToLocal(padCentreL2,padCentreG);
2150 :
2151 :
2152 0 : Char_t path2[200];
2153 : // strip centre coordinates in its ref. frame
2154 0 : Double_t stripCentreL[3] = {0., 0., 0.};
2155 : // strip centre coordinates in the ALICE global ref. frame
2156 0 : Double_t stripCentreG[3] = {0., 0., 0.};
2157 0 : GetVolumePath(detId[0],detId[1],detId[2],path2);
2158 0 : gGeoManager->cd(path2);
2159 0 : TGeoHMatrix g2lb = *gGeoManager->GetCurrentMatrix();
2160 0 : TGeoHMatrix l2gb = g2lb.Inverse();
2161 0 : l2gb.MasterToLocal(stripCentreL,stripCentreG);
2162 :
2163 0 : TGeoHMatrix g2t = 0;
2164 0 : AliGeomManager::GetTrackingMatrix(alignableStripIndex, g2t);
2165 :
2166 : // track position in the ALICE global ref. frame
2167 0 : Double_t posG[3];
2168 0 : for (Int_t ii=0; ii<3; ii++) posG[ii] = (Double_t)ctrackPos[ii];
2169 :
2170 : // strip centre coordinates in the tracking ref. frame
2171 0 : Double_t stripCentreT[3] = {0., 0., 0.};
2172 : // track position in the tracking ref. frame
2173 0 : Double_t posT[3] = {0., 0., 0.};
2174 0 : g2t.MasterToLocal(posG,posT);
2175 0 : g2t.MasterToLocal(stripCentreG,stripCentreT);
2176 :
2177 0 : for (Int_t ii=0; ii<3; ii++)
2178 0 : AliDebug(1,Form(" track position in ALICE global and tracking RFs -> posG[%d] = %f --- posT[%d] = %f",
2179 : ii, posG[ii], ii, posT[ii]));
2180 0 : for (Int_t ii=0; ii<3; ii++)
2181 0 : AliDebug(1,Form(" pad centre coordinates in its, the ALICE global and tracking RFs -> "
2182 : "padCentreL[%d] = %f --- padCentreG[%d] = %f --- padCentreT[%d] = %f",
2183 : ii, padCentreL[ii],
2184 : ii, padCentreG[ii],
2185 : ii, padCentreT[ii]));
2186 0 : for (Int_t ii=0; ii<3; ii++)
2187 0 : AliDebug(1,Form(" strip centre coordinates in its, the ALICE global and tracking RFs -> "
2188 : "stripCentreL[%d] = %f --- stripCentreG[%d] = %f --- stripCentreT[%d] = %f",
2189 : ii, stripCentreL[ii],
2190 : ii, stripCentreG[ii],
2191 : ii, stripCentreT[ii]));
2192 0 : for (Int_t ii=0; ii<3; ii++)
2193 0 : AliDebug(1,Form(" difference between the track position and the pad centre in the tracking RF "
2194 : "-> posT[%d]-padCentreT[%d] = %f",
2195 : ii,ii,
2196 : posT[ii]-padCentreT[ii]));
2197 :
2198 0 : for (Int_t ii=0; ii<3; ii++) differenceT[ii] = (Float_t)(posT[ii]-padCentreT[ii]);
2199 :
2200 0 : }
2201 : //-------------------------------------------------------------------------
2202 :
2203 : Int_t AliTOFGeometry::GetTOFsupermodule(Int_t index)
2204 : {
2205 : // Return the TOF supermodule where TOF channel index is located
2206 :
2207 0 : if (index<0 || index>=NPadXSector()*NSectors()) return -1;
2208 0 : else return index/NpadXStrip()/NStripXSector();
2209 :
2210 0 : }
|