]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDRecParam.cxx
Added a commented out version with new digitizers.
[u/mrichter/AliRoot.git] / TRD / AliTRDRecParam.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 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 // Class containing constant reconstruction parameters                       //
21 //                                                                           //
22 // Request an instance with AliTRDRecParam::Instance()                 //
23 // Then request the needed values                                            //
24 //                                                                           //
25 ///////////////////////////////////////////////////////////////////////////////
26
27 #include "AliTRDRecParam.h"
28
29 ClassImp(AliTRDRecParam)
30
31 AliTRDRecParam* AliTRDRecParam::fgInstance = 0;
32 Bool_t AliTRDRecParam::fgTerminated = kFALSE;
33
34 //_ singleton implementation __________________________________________________
35 AliTRDRecParam* AliTRDRecParam::Instance()
36 {
37   //
38   // Singleton implementation
39   // Returns an instance of this class, it is created if neccessary
40   // 
41   
42   if (fgTerminated != kFALSE)
43     return 0;
44
45   if (fgInstance == 0)
46     fgInstance = new AliTRDRecParam();
47   
48   return fgInstance;
49 }
50
51 void AliTRDRecParam::Terminate()
52 {
53   //
54   // Singleton implementation
55   // Deletes the instance of this class and sets the terminated flag, instances cannot be requested anymore
56   // This function can be called several times.
57   //
58   
59   fgTerminated = kTRUE;
60   
61   if (fgInstance != 0)
62   {
63     delete fgInstance;
64     fgInstance = 0;
65   }
66 }
67
68 //_____________________________________________________________________________
69 AliTRDRecParam::AliTRDRecParam()
70 {
71   //
72   // constructor
73   //
74   
75   fClusMaxThresh      = 0;
76   fClusSigThresh      = 0;
77   
78   fLUT                = 0;
79   fLUTOn              = kFALSE;  
80   fLUTbin = 0;
81   
82   Init();
83 }
84
85 //_____________________________________________________________________________
86 AliTRDRecParam::~AliTRDRecParam() 
87 {
88   //
89   // destructor
90   //
91   
92   if (fLUT) {
93     delete [] fLUT;
94     fLUT    = 0;
95   }
96 }
97
98 //_____________________________________________________________________________
99 AliTRDRecParam::AliTRDRecParam(const AliTRDRecParam &p):TObject(p)
100 {
101   //
102   // copy constructor
103   //
104
105   ((AliTRDRecParam &) p).Copy(*this);
106 }
107
108
109 //_____________________________________________________________________________
110 AliTRDRecParam &AliTRDRecParam::operator=(const AliTRDRecParam &p)
111 {
112   //
113   // Assignment operator
114   //
115
116   if (this != &p) ((AliTRDRecParam &) p).Copy(*this);
117   return *this;
118 }
119
120 //_____________________________________________________________________________
121 void AliTRDRecParam::Copy(TObject &p) const
122 {
123   //
124   // Copy function
125   //
126   
127   AliTRDRecParam* target = dynamic_cast<AliTRDRecParam*> (&p);
128   if (!target)
129     return;
130   
131   target->fLUTOn              = fLUTOn;
132   target->fLUTbin             = fLUTbin;
133   if (target->fLUT)    delete [] target->fLUT;
134   target->fLUT  = new Float_t[fLUTbin];
135   for (Int_t iBin = 0; iBin < fLUTbin; iBin++) {
136     target->fLUT[iBin]  = fLUT[iBin];
137   }
138       
139   target->fClusMaxThresh      = fClusMaxThresh;
140   target->fClusSigThresh      = fClusSigThresh;
141   
142 }
143
144 //_____________________________________________________________________________
145 void AliTRDRecParam::Init() 
146 {
147   //
148   // constructor helper
149   //
150   
151   // The default parameter for the clustering
152   fClusMaxThresh = 3;
153   fClusSigThresh = 1;
154
155   // Use the LUT
156   fLUTOn         = kTRUE;
157
158   // Create the LUT
159   FillLUT();
160 }
161
162 //_____________________________________________________________________________
163 void AliTRDRecParam::FillLUT()
164 {
165   //
166   // Create the LUT
167   //
168
169   const Int_t kNlut  = 128;
170
171   fLUTbin = kNplan * kNlut;
172
173   // The lookup table from Bogdan
174   Float_t lut[kNplan][kNlut] = {  
175     {
176       0.0070, 0.0150, 0.0224, 0.0298, 0.0374, 0.0454, 0.0533, 0.0611, 
177       0.0684, 0.0755, 0.0827, 0.0900, 0.0975, 0.1049, 0.1120, 0.1187, 
178       0.1253, 0.1318, 0.1385, 0.1453, 0.1519, 0.1584, 0.1646, 0.1704, 
179       0.1762, 0.1821, 0.1879, 0.1938, 0.1996, 0.2053, 0.2108, 0.2160, 
180       0.2210, 0.2260, 0.2310, 0.2361, 0.2411, 0.2461, 0.2509, 0.2557, 
181       0.2602, 0.2646, 0.2689, 0.2732, 0.2774, 0.2816, 0.2859, 0.2901, 
182       0.2942, 0.2983, 0.3022, 0.3061, 0.3099, 0.3136, 0.3172, 0.3207, 
183       0.3242, 0.3278, 0.3312, 0.3347, 0.3382, 0.3416, 0.3450, 0.3483, 
184       0.3515, 0.3547, 0.3579, 0.3609, 0.3639, 0.3669, 0.3698, 0.3727, 
185       0.3756, 0.3785, 0.3813, 0.3842, 0.3870, 0.3898, 0.3926, 0.3952, 
186       0.3979, 0.4005, 0.4032, 0.4057, 0.4082, 0.4108, 0.4132, 0.4157, 
187       0.4181, 0.4205, 0.4228, 0.4252, 0.4275, 0.4299, 0.4322, 0.4345, 
188       0.4367, 0.4390, 0.4412, 0.4434, 0.4456, 0.4478, 0.4499, 0.4520, 
189       0.4541, 0.4562, 0.4583, 0.4603, 0.4623, 0.4643, 0.4663, 0.4683, 
190       0.4702, 0.4722, 0.4741, 0.4758, 0.4774, 0.4790, 0.4805, 0.4824, 
191       0.4844, 0.4863, 0.4883, 0.4902, 0.4921, 0.4940, 0.4959, 0.4978 
192     },
193     {
194       0.0072, 0.0156, 0.0235, 0.0313, 0.0394, 0.0478, 0.0561, 0.0642, 
195       0.0718, 0.0792, 0.0868, 0.0947, 0.1025, 0.1101, 0.1172, 0.1241, 
196       0.1309, 0.1378, 0.1449, 0.1518, 0.1586, 0.1650, 0.1710, 0.1770, 
197       0.1830, 0.1891, 0.1952, 0.2011, 0.2070, 0.2125, 0.2177, 0.2229, 
198       0.2280, 0.2332, 0.2383, 0.2435, 0.2484, 0.2533, 0.2581, 0.2627, 
199       0.2670, 0.2714, 0.2757, 0.2799, 0.2842, 0.2884, 0.2927, 0.2968, 
200       0.3008, 0.3048, 0.3086, 0.3123, 0.3159, 0.3195, 0.3231, 0.3266, 
201       0.3301, 0.3335, 0.3370, 0.3404, 0.3438, 0.3471, 0.3504, 0.3536, 
202       0.3567, 0.3598, 0.3628, 0.3657, 0.3686, 0.3715, 0.3744, 0.3772, 
203       0.3800, 0.3828, 0.3856, 0.3884, 0.3911, 0.3938, 0.3965, 0.3991, 
204       0.4016, 0.4042, 0.4067, 0.4092, 0.4116, 0.4140, 0.4164, 0.4187, 
205       0.4211, 0.4234, 0.4257, 0.4280, 0.4302, 0.4325, 0.4347, 0.4369, 
206       0.4391, 0.4413, 0.4434, 0.4456, 0.4477, 0.4497, 0.4518, 0.4538, 
207       0.4558, 0.4578, 0.4598, 0.4618, 0.4637, 0.4656, 0.4675, 0.4694, 
208       0.4713, 0.4732, 0.4750, 0.4766, 0.4781, 0.4797, 0.4813, 0.4832, 
209       0.4851, 0.4870, 0.4888, 0.4906, 0.4925, 0.4942, 0.4960, 0.4978
210     },
211     {
212       0.0075, 0.0163, 0.0246, 0.0328, 0.0415, 0.0504, 0.0592, 0.0674, 
213       0.0753, 0.0832, 0.0914, 0.0996, 0.1077, 0.1154, 0.1225, 0.1296, 
214       0.1369, 0.1442, 0.1515, 0.1585, 0.1652, 0.1714, 0.1776, 0.1839, 
215       0.1902, 0.1965, 0.2025, 0.2085, 0.2141, 0.2194, 0.2247, 0.2299, 
216       0.2352, 0.2405, 0.2457, 0.2507, 0.2557, 0.2604, 0.2649, 0.2693, 
217       0.2737, 0.2780, 0.2823, 0.2867, 0.2909, 0.2951, 0.2992, 0.3033, 
218       0.3072, 0.3110, 0.3146, 0.3182, 0.3218, 0.3253, 0.3288, 0.3323, 
219       0.3357, 0.3392, 0.3426, 0.3459, 0.3492, 0.3524, 0.3555, 0.3586, 
220       0.3616, 0.3645, 0.3674, 0.3703, 0.3731, 0.3759, 0.3787, 0.3815, 
221       0.3843, 0.3870, 0.3897, 0.3925, 0.3950, 0.3976, 0.4002, 0.4027, 
222       0.4052, 0.4076, 0.4101, 0.4124, 0.4148, 0.4171, 0.4194, 0.4217, 
223       0.4239, 0.4262, 0.4284, 0.4306, 0.4328, 0.4350, 0.4371, 0.4393, 
224       0.4414, 0.4435, 0.4455, 0.4476, 0.4496, 0.4516, 0.4536, 0.4555, 
225       0.4575, 0.4594, 0.4613, 0.4632, 0.4650, 0.4669, 0.4687, 0.4705, 
226       0.4723, 0.4741, 0.4758, 0.4773, 0.4789, 0.4804, 0.4821, 0.4839, 
227       0.4857, 0.4875, 0.4893, 0.4910, 0.4928, 0.4945, 0.4961, 0.4978
228     },
229     {
230       0.0078, 0.0171, 0.0258, 0.0345, 0.0438, 0.0532, 0.0624, 0.0708, 
231       0.0791, 0.0875, 0.0962, 0.1048, 0.1130, 0.1206, 0.1281, 0.1356, 
232       0.1432, 0.1508, 0.1582, 0.1651, 0.1716, 0.1780, 0.1845, 0.1910, 
233       0.1975, 0.2038, 0.2099, 0.2155, 0.2210, 0.2263, 0.2317, 0.2371, 
234       0.2425, 0.2477, 0.2528, 0.2578, 0.2626, 0.2671, 0.2715, 0.2759, 
235       0.2803, 0.2846, 0.2890, 0.2933, 0.2975, 0.3016, 0.3056, 0.3095, 
236       0.3132, 0.3168, 0.3204, 0.3239, 0.3274, 0.3309, 0.3344, 0.3378, 
237       0.3412, 0.3446, 0.3479, 0.3511, 0.3543, 0.3574, 0.3603, 0.3633, 
238       0.3662, 0.3690, 0.3718, 0.3747, 0.3774, 0.3802, 0.3829, 0.3857, 
239       0.3883, 0.3910, 0.3936, 0.3962, 0.3987, 0.4012, 0.4037, 0.4061, 
240       0.4085, 0.4109, 0.4132, 0.4155, 0.4177, 0.4200, 0.4222, 0.4244, 
241       0.4266, 0.4288, 0.4309, 0.4331, 0.4352, 0.4373, 0.4394, 0.4414, 
242       0.4435, 0.4455, 0.4475, 0.4494, 0.4514, 0.4533, 0.4552, 0.4571, 
243       0.4590, 0.4608, 0.4626, 0.4645, 0.4662, 0.4680, 0.4698, 0.4715, 
244       0.4733, 0.4750, 0.4766, 0.4781, 0.4796, 0.4812, 0.4829, 0.4846, 
245       0.4863, 0.4880, 0.4897, 0.4914, 0.4930, 0.4946, 0.4963, 0.4979
246     },
247     {
248       0.0081, 0.0178, 0.0270, 0.0364, 0.0463, 0.0562, 0.0656, 0.0744, 
249       0.0831, 0.0921, 0.1013, 0.1102, 0.1183, 0.1261, 0.1339, 0.1419, 
250       0.1499, 0.1576, 0.1648, 0.1715, 0.1782, 0.1849, 0.1917, 0.1984, 
251       0.2048, 0.2110, 0.2167, 0.2223, 0.2278, 0.2333, 0.2389, 0.2444, 
252       0.2497, 0.2548, 0.2598, 0.2645, 0.2691, 0.2735, 0.2780, 0.2824, 
253       0.2868, 0.2912, 0.2955, 0.2997, 0.3038, 0.3078, 0.3116, 0.3152, 
254       0.3188, 0.3224, 0.3259, 0.3294, 0.3329, 0.3364, 0.3398, 0.3432, 
255       0.3465, 0.3497, 0.3529, 0.3561, 0.3591, 0.3620, 0.3649, 0.3677, 
256       0.3705, 0.3733, 0.3761, 0.3788, 0.3816, 0.3843, 0.3869, 0.3896, 
257       0.3922, 0.3948, 0.3973, 0.3998, 0.4022, 0.4047, 0.4070, 0.4094, 
258       0.4117, 0.4139, 0.4162, 0.4184, 0.4206, 0.4227, 0.4249, 0.4270, 
259       0.4291, 0.4313, 0.4334, 0.4354, 0.4375, 0.4395, 0.4415, 0.4435, 
260       0.4455, 0.4474, 0.4493, 0.4512, 0.4531, 0.4550, 0.4568, 0.4586, 
261       0.4604, 0.4622, 0.4639, 0.4657, 0.4674, 0.4691, 0.4708, 0.4725, 
262       0.4742, 0.4758, 0.4773, 0.4788, 0.4803, 0.4819, 0.4836, 0.4852, 
263       0.4869, 0.4885, 0.4901, 0.4917, 0.4933, 0.4948, 0.4964, 0.4979
264     },
265     {
266       0.0085, 0.0189, 0.0288, 0.0389, 0.0497, 0.0603, 0.0699, 0.0792, 
267       0.0887, 0.0985, 0.1082, 0.1170, 0.1253, 0.1336, 0.1421, 0.1505, 
268       0.1587, 0.1662, 0.1733, 0.1803, 0.1874, 0.1945, 0.2014, 0.2081, 
269       0.2143, 0.2201, 0.2259, 0.2316, 0.2374, 0.2431, 0.2487, 0.2541, 
270       0.2593, 0.2642, 0.2689, 0.2735, 0.2781, 0.2826, 0.2872, 0.2917, 
271       0.2961, 0.3003, 0.3045, 0.3086, 0.3125, 0.3162, 0.3198, 0.3235, 
272       0.3270, 0.3306, 0.3342, 0.3377, 0.3411, 0.3446, 0.3479, 0.3511, 
273       0.3543, 0.3575, 0.3605, 0.3634, 0.3663, 0.3691, 0.3720, 0.3748, 
274       0.3775, 0.3803, 0.3830, 0.3857, 0.3884, 0.3911, 0.3937, 0.3962, 
275       0.3987, 0.4012, 0.4036, 0.4060, 0.4084, 0.4107, 0.4129, 0.4152, 
276       0.4174, 0.4196, 0.4218, 0.4239, 0.4261, 0.4282, 0.4303, 0.4324, 
277       0.4344, 0.4365, 0.4385, 0.4405, 0.4425, 0.4445, 0.4464, 0.4483, 
278       0.4502, 0.4521, 0.4539, 0.4558, 0.4576, 0.4593, 0.4611, 0.4629, 
279       0.4646, 0.4663, 0.4680, 0.4697, 0.4714, 0.4730, 0.4747, 0.4759, 
280       0.4769, 0.4780, 0.4790, 0.4800, 0.4811, 0.4827, 0.4843, 0.4859, 
281       0.4874, 0.4889, 0.4905, 0.4920, 0.4935, 0.4950, 0.4965, 0.4979
282     }
283   }; 
284
285   if (fLUT) delete [] fLUT;
286   fLUT = new Float_t[fLUTbin];
287
288   for (Int_t iplan = 0; iplan < kNplan; iplan++) {
289     for (Int_t ilut  = 0; ilut  < kNlut; ilut++) {
290       fLUT[iplan*kNlut+ilut] = lut[iplan][ilut];
291     }
292   }
293
294 }
295
296 //_____________________________________________________________________________
297 Double_t AliTRDRecParam::LUTposition(Int_t iplane, Double_t ampL, Double_t ampC, Double_t ampR) const
298 {
299   //
300   // Calculates the cluster position using the lookup table.
301   // Method provided by Bogdan Vulpescu.
302   //
303
304   const Int_t kNlut  = 128;
305
306   Double_t pos;
307   Double_t x = 0.0;
308   Double_t xmin;
309   Double_t xmax;
310   Double_t xwid;
311
312   Int_t    side = 0;
313   Int_t    ix;
314
315   Double_t xMin[kNplan] = { 0.006492, 0.006377, 0.006258, 0.006144, 0.006030, 0.005980 };
316   Double_t xMax[kNplan] = { 0.960351, 0.965870, 0.970445, 0.974352, 0.977667, 0.996101 };
317
318   if      (ampL > ampR) {
319     x    = (ampL - ampR) / ampC;
320     side = -1;
321   } 
322   else if (ampL < ampR) {
323     x    = (ampR - ampL) / ampC;
324     side = +1;
325   }
326
327   if (ampL != ampR) {
328
329     xmin = xMin[iplane] + 0.000005;
330     xmax = xMax[iplane] - 0.000005;
331     xwid = (xmax - xmin) / 127.0;
332
333     if      (x < xmin) {
334       pos = 0.0000;
335     } 
336     else if (x > xmax) {
337       pos = side * 0.5000;
338     } 
339     else {
340       ix  = (Int_t) ((x - xmin) / xwid);
341       pos = side * fLUT[iplane*kNlut+ix];
342     }
343        
344   } 
345   else {
346
347     pos = 0.0;
348
349   }
350
351   return pos;
352 }
353