]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDparameter.cxx
Macro to copy gAlice, TreeE, TreeK to a file
[u/mrichter/AliRoot.git] / TRD / AliTRDparameter.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 /*
17 $Log$
18 Revision 1.3  2002/03/28 14:59:07  cblume
19 Coding conventions
20
21 Revision 1.2  2002/03/28 10:00:36  hristov
22 Some additional initialisation
23
24 Revision 1.1  2002/03/25 20:01:18  cblume
25 Introduce parameter class
26
27 */
28
29 ///////////////////////////////////////////////////////////////////////////////
30 //                                                                           //
31 //  TRD parameter class                                                      //
32 //                                                                           //
33 ///////////////////////////////////////////////////////////////////////////////
34
35 #include "AliRun.h"
36 #include "AliMagF.h"
37
38 #include "AliTRD.h"
39 #include "AliTRDparameter.h"
40 #include "AliTRDgeometryFull.h"
41
42 ClassImp(AliTRDparameter)
43
44 //_____________________________________________________________________________
45 AliTRDparameter::AliTRDparameter():TNamed()
46 {
47   //
48   // AliTRDparameter default constructor
49   //
50
51   fGeo                = 0;
52   fPRFsmp             = 0;
53   fTRFsmp             = 0;
54   fCTsmp              = 0;
55   fGasGain            = 0.0;
56   fNoise              = 0.0;
57   fChipGain           = 0.0;
58   fADCoutRange        = 0.0;
59   fADCinRange         = 0.0;
60   fADCthreshold       = 0;
61   fDiffusionOn        = 0;
62   fDiffusionT         = 0.0;
63   fDiffusionL         = 0.0;
64   fElAttachOn         = 0;
65   fElAttachProp       = 0.0;
66   fExBOn              = 0;
67   fOmegaTau           = 0.0;
68   fPRFOn              = 0;
69   fTRFOn              = 0;
70   fCTOn               = 0;
71   fTCOn               = 0;
72   fDriftVelocity      = 0.0;
73   fPadCoupling        = 0.0;
74   fTimeCoupling       = 0.0;
75   fTimeBinWidth       = 0.0;
76   fField              = 0.0;
77   fTiltingAngle       = 0.0;
78   fPRFbin             = 0;
79   fPRFlo              = 0.0;
80   fPRFhi              = 0.0;
81   fPRFwid             = 0.0;
82   fPRFpad             = 0;
83   fTRFbin             = 0;
84   fTRFlo              = 0.0;
85   fTRFhi              = 0.0;
86   fTRFwid             = 0.0;
87   fTCnexp             = 0;
88
89   fLUTOn              = 0;  
90   fLUT                = 0;
91   fClusMaxThresh      = 0;
92   fClusSigThresh      = 0;
93
94 }
95
96 //_____________________________________________________________________________
97 AliTRDparameter::AliTRDparameter(const Text_t *name, const Text_t *title)
98                 :TNamed(name,title)
99 {
100   //
101   // AliTRDparameter constructor
102   //
103
104   fGeo                = new AliTRDgeometryFull();
105   fPRFsmp             = 0;
106   fTRFsmp             = 0;
107   fCTsmp              = 0;
108   fGasGain            = 0.0;
109   fNoise              = 0.0;
110   fChipGain           = 0.0;
111   fADCoutRange        = 0.0;
112   fADCinRange         = 0.0;
113   fADCthreshold       = 0;
114   fDiffusionOn        = 0;
115   fDiffusionT         = 0.0;
116   fDiffusionL         = 0.0;
117   fElAttachOn         = 0;
118   fElAttachProp       = 0.0;
119   fExBOn              = 0;
120   fOmegaTau           = 0.0;
121   fPRFOn              = 0;
122   fTRFOn              = 0;
123   fCTOn               = 0;
124   fTCOn               = 0;
125   fDriftVelocity      = 0.0;
126   fPadCoupling        = 0.0;
127   fTimeCoupling       = 0.0;
128   fTimeBinWidth       = 0.0;
129   fField              = 0.0;
130   fTiltingAngle       = 0.0;
131   fPRFbin             = 0;
132   fPRFlo              = 0.0;
133   fPRFhi              = 0.0;
134   fPRFwid             = 0.0;
135   fPRFpad             = 0;
136   fTRFbin             = 0;
137   fTRFlo              = 0.0;
138   fTRFhi              = 0.0;
139   fTRFwid             = 0.0;
140   fTCnexp             = 0;
141
142   fLUTOn              = 0;  
143   fLUT                = 0;
144   fClusMaxThresh      = 0;
145   fClusSigThresh      = 0;
146
147   Init();
148
149 }
150
151
152 //_____________________________________________________________________________
153 AliTRDparameter::AliTRDparameter(const AliTRDparameter &p)
154 {
155   //
156   // AliTRDparameter copy constructor
157   //
158
159   ((AliTRDparameter &) p).Copy(*this);
160
161 }
162
163 ///_____________________________________________________________________________
164 AliTRDparameter::~AliTRDparameter()
165 {
166   //
167   // AliTRDparameter destructor
168   //
169
170   if (fTRFsmp) {
171     delete [] fTRFsmp;
172     fTRFsmp = 0;
173   }
174
175   if (fPRFsmp) {
176     delete [] fPRFsmp;
177     fPRFsmp = 0;
178   }
179
180   if (fCTsmp) {
181     delete [] fCTsmp;
182     fCTsmp  = 0;
183   }
184
185   if (fLUT) {
186     delete [] fLUT;
187     fLUT    = 0;
188   }
189
190   if (fGeo) {
191     delete fGeo;
192     fGeo    = 0;
193   }
194
195 }
196
197 //_____________________________________________________________________________
198 AliTRDparameter &AliTRDparameter::operator=(const AliTRDparameter &p)
199 {
200   //
201   // Assignment operator
202   //
203
204   if (this != &p) ((AliTRDparameter &) p).Copy(*this);
205   return *this;
206
207 }
208
209 //_____________________________________________________________________________
210 void AliTRDparameter::Copy(TObject &p)
211 {
212   //
213   // Copy function
214   //
215
216   Int_t iBin;
217
218   ((AliTRDparameter &) p).fGasGain            = fGasGain;
219   ((AliTRDparameter &) p).fNoise              = fNoise;
220   ((AliTRDparameter &) p).fChipGain           = fChipGain;
221   ((AliTRDparameter &) p).fADCoutRange        = fADCoutRange;
222   ((AliTRDparameter &) p).fADCinRange         = fADCinRange;
223   ((AliTRDparameter &) p).fADCthreshold       = fADCthreshold;
224   ((AliTRDparameter &) p).fDiffusionOn        = fDiffusionOn; 
225   ((AliTRDparameter &) p).fDiffusionT         = fDiffusionT;
226   ((AliTRDparameter &) p).fDiffusionL         = fDiffusionL;
227   ((AliTRDparameter &) p).fElAttachOn         = fElAttachOn;
228   ((AliTRDparameter &) p).fElAttachProp       = fElAttachProp;
229   ((AliTRDparameter &) p).fExBOn              = fExBOn;
230   ((AliTRDparameter &) p).fOmegaTau           = fOmegaTau;
231   ((AliTRDparameter &) p).fLorentzFactor      = fLorentzFactor;
232   ((AliTRDparameter &) p).fDriftVelocity      = fDriftVelocity;
233   ((AliTRDparameter &) p).fPadCoupling        = fPadCoupling;
234   ((AliTRDparameter &) p).fTimeCoupling       = fTimeCoupling;
235   ((AliTRDparameter &) p).fTimeBinWidth       = fTimeBinWidth;
236   ((AliTRDparameter &) p).fField              = fField;
237   ((AliTRDparameter &) p).fPRFOn              = fPRFOn;
238   ((AliTRDparameter &) p).fTRFOn              = fTRFOn;
239   ((AliTRDparameter &) p).fCTOn               = fCTOn;
240   ((AliTRDparameter &) p).fTCOn               = fTCOn;
241   ((AliTRDparameter &) p).fTiltingAngle       = fTiltingAngle;
242   ((AliTRDparameter &) p).fPRFbin             = fPRFbin;
243   ((AliTRDparameter &) p).fPRFlo              = fPRFlo;
244   ((AliTRDparameter &) p).fPRFhi              = fPRFhi;
245   ((AliTRDparameter &) p).fPRFwid             = fPRFwid;
246   ((AliTRDparameter &) p).fPRFpad             = fPRFpad;
247   if (((AliTRDparameter &) p).fPRFsmp) delete [] ((AliTRDparameter &) p).fPRFsmp;
248   ((AliTRDparameter &) p).fPRFsmp = new Float_t[fPRFbin];
249   for (iBin = 0; iBin < fPRFbin; iBin++) {
250     ((AliTRDparameter &) p).fPRFsmp[iBin] = fPRFsmp[iBin];
251   }                                                                             
252   ((AliTRDparameter &) p).fTRFbin             = fTRFbin;
253   ((AliTRDparameter &) p).fTRFlo              = fTRFlo;
254   ((AliTRDparameter &) p).fTRFhi              = fTRFhi;
255   ((AliTRDparameter &) p).fTRFwid             = fTRFwid;
256   if (((AliTRDparameter &) p).fTRFsmp) delete [] ((AliTRDparameter &) p).fTRFsmp;
257   ((AliTRDparameter &) p).fTRFsmp = new Float_t[fTRFbin];
258   for (iBin = 0; iBin < fTRFbin; iBin++) {
259     ((AliTRDparameter &) p).fTRFsmp[iBin] = fTRFsmp[iBin];
260   }                                      
261   if (((AliTRDparameter &) p).fCTsmp)  delete [] ((AliTRDparameter &) p).fCTsmp;
262   ((AliTRDparameter &) p).fCTsmp  = new Float_t[fTRFbin];
263   for (iBin = 0; iBin < fTRFbin; iBin++) {
264     ((AliTRDparameter &) p).fCTsmp[iBin]  = fCTsmp[iBin];
265   }                                      
266   ((AliTRDparameter &) p).fTCnexp             = fTCnexp;
267
268   ((AliTRDparameter &) p).fLUTOn              = fLUTOn;
269   ((AliTRDparameter &) p).fLUTbin             = fLUTbin;
270   if (((AliTRDparameter &) p).fLUT)    delete [] ((AliTRDparameter &) p).fLUT;
271   ((AliTRDparameter &) p).fLUT  = new Float_t[fLUTbin];
272   for (iBin = 0; iBin < fLUTbin; iBin++) {
273     ((AliTRDparameter &) p).fLUT[iBin]  = fLUT[iBin];
274   }                                      
275   ((AliTRDparameter &) p).fClusMaxThresh      = fClusMaxThresh;
276   ((AliTRDparameter &) p).fClusSigThresh      = fClusSigThresh;
277
278 }
279
280 //_____________________________________________________________________________
281 void AliTRDparameter::Init()
282 {
283   //
284   // Initializes the parameter
285   //
286   // The maximum number of pads
287   // and the position of pad 0,0,0
288   //
289   // chambers seen from the top:
290   //     +----------------------------+
291   //     |                            |
292   //     |                            |      ^
293   //     |                            |  rphi|
294   //     |                            |      |
295   //     |0                           |      |
296   //     +----------------------------+      +------>
297   //                                             z
298   // chambers seen from the side:            ^
299   //     +----------------------------+ drift|
300   //     |0                           |      |
301   //     |                            |      |
302   //     +----------------------------+      +------>
303   //                                             z
304   //
305   // IMPORTANT: time bin 0 is now the first one in the drift region
306   // closest to the readout !!!
307   //
308
309   //
310   // ----------------------------------------------------------------------------
311   // The pad definition
312   // ----------------------------------------------------------------------------
313   //
314
315   // The pad size in column direction (rphi-direction)
316   SetColPadSize(0,0.65);
317   SetColPadSize(1,0.68);
318   SetColPadSize(2,0.71);
319   SetColPadSize(3,0.74);
320   SetColPadSize(4,0.77);
321   SetColPadSize(5,0.80);
322
323   // The pad row (z-direction)
324   SetNRowPad();
325
326   // The number of time bins. Default is 100 ns timbin size
327   SetNTimeBin(15);
328
329   // Additional time bins before and after the drift region.
330   // Default is to only sample the drift region
331   SetExpandTimeBin(0,0);
332
333   //
334   // ----------------------------------------------------------------------------
335   // The digitization parameter
336   // ----------------------------------------------------------------------------
337   //
338
339   // The default parameter for the digitization
340   fGasGain        = 4500.;
341   fChipGain       = 12.4;
342   fNoise          = 1000.;
343   fADCoutRange    = 1023.;          // 10-bit ADC
344   fADCinRange     = 1000.;          // 1V input range
345   fADCthreshold   = 1;
346
347   // The drift velocity (cm / mus)
348   fDriftVelocity  = 1.5;
349
350   // Diffusion on
351   fDiffusionOn    = 1;
352
353   // E x B effects
354   fExBOn          = 0;
355
356   // Propability for electron attachment
357   fElAttachOn     = 0;
358   fElAttachProp   = 0.0;
359
360   // The pad response function
361   fPRFOn          = 1;
362
363   // The time response function
364   fTRFOn          = 1;
365
366   // The cross talk
367   fCTOn           = 0;
368
369   // The tail cancelation
370   fTCOn           = 1;
371   
372   // The number of exponentials
373   fTCnexp         = 2;
374
375   // The pad coupling factor (same number as for the TPC)
376   fPadCoupling    = 0.5;
377
378   // The time coupling factor (same number as for the TPC)
379   fTimeCoupling   = 0.4;
380
381   // The tilting angle for the readout pads
382   SetTiltingAngle(5.0);
383
384   // The magnetic field strength in Tesla
385   fField           = 0.2 * gAlice->Field()->Factor();
386
387   //
388   // ----------------------------------------------------------------------------
389   // The clusterization parameter
390   // ----------------------------------------------------------------------------
391   //
392
393   // The default parameter for the clustering
394   fClusMaxThresh = 3;
395   fClusSigThresh = 1;
396
397   // Use the LUT
398   fLUTOn         = 1;  
399
400   ReInit();
401
402 }
403
404 //_____________________________________________________________________________
405 void AliTRDparameter::ReInit()
406 {
407   //
408   // Reinitializes the parameter class after a change
409   //
410
411   // Calculate the time bin width in ns
412   fTimeBinWidth   = fTimeBinSize / fDriftVelocity * 1000.0;
413
414   // The range and the binwidth for the sampled TRF 
415   fTRFbin = 100;
416   // Start 0.2 mus before the signal
417   fTRFlo  = -0.2 * fDriftVelocity;
418   // End the maximum driftlength after the signal 
419   fTRFhi  = AliTRDgeometry::DrThick() 
420           + fTimeAfter * fTimeBinSize;
421   fTRFwid = (fTRFhi - fTRFlo) / ((Float_t) fTRFbin);
422
423   // Transverse and longitudinal diffusion coefficients (Xe/CO2)
424   fDiffusionT     = GetDiffusionT(fDriftVelocity,fField);
425   fDiffusionL     = GetDiffusionL(fDriftVelocity,fField);
426
427   // omega * tau.= tan(Lorentz-angle)
428   fOmegaTau       = GetOmegaTau(fDriftVelocity,fField);
429
430   // The Lorentz factor
431   if (fExBOn) {
432     fLorentzFactor = 1.0 / (1.0 + fOmegaTau*fOmegaTau);
433   }
434   else {
435     fLorentzFactor = 1.0;
436   }
437
438   // Create the sampled PRF
439   SamplePRF();
440
441   // Create the sampled TRF
442   SampleTRF();
443
444   // Create the LUT
445   FillLUT();
446
447 }
448
449 //_____________________________________________________________________________
450 void AliTRDparameter::SetNRowPad(const Int_t p, const Int_t c, const Int_t npad)
451 {
452   //
453   // Redefines the number of pads in raw direction for
454   // a given plane and chamber number
455   //
456
457   for (Int_t isect = 0; isect < AliTRDgeometry::Nsect(); isect++) {
458
459     fRowMax[p][c][isect] = npad;
460
461     fRowPadSize[p][c][isect] = (fGeo->GetChamberLength(p,c) 
462                                 - 2.* AliTRDgeometry::RpadW())
463                              / ((Float_t) npad);
464
465   }
466
467 }
468
469 //_____________________________________________________________________________
470 void AliTRDparameter::SetNRowPad()
471 {
472   //
473   // Defines the number of pads in row direction
474   //
475
476   Int_t isect;
477   Int_t icham;
478   Int_t iplan;
479
480   Int_t rowMax[kNplan][kNcham] = { { 16, 16, 12, 16, 16 }
481                                  , { 16, 16, 12, 16, 16 }
482                                  , { 16, 16, 12, 16, 16 }
483                                  , { 16, 16, 12, 16, 16 }
484                                  , { 14, 16, 12, 16, 14 }
485                                  , { 13, 16, 12, 16, 13 } };
486
487   Float_t rpadW = AliTRDgeometry::RpadW();
488
489   for (isect = 0; isect < kNsect; isect++) {
490     for (icham = 0; icham < kNcham; icham++) {
491       for (iplan = 0; iplan < kNplan; iplan++) {
492
493         fRowMax[iplan][icham][isect]     = rowMax[iplan][icham];
494
495         fRowPadSize[iplan][icham][isect] = (fGeo->GetChamberLength(iplan,icham) 
496                                             - 2.*rpadW)
497                                          / ((Float_t) rowMax[iplan][icham]);
498
499         Float_t row0 = rpadW - fGeo->GetChamberLength(iplan,0)
500                              - fGeo->GetChamberLength(iplan,1)
501                              - fGeo->GetChamberLength(iplan,2) / 2.;
502         for (Int_t ic = 0; ic < icham; ic++) {
503           row0 += fGeo->GetChamberLength(iplan,icham);
504         }
505         fRow0[iplan][icham][isect]       = row0;
506
507       }
508     }
509   }
510
511 }
512
513 //_____________________________________________________________________________
514 void AliTRDparameter::SetColPadSize(const Int_t p, const Float_t s)
515 {
516   //
517   // Redefines the pad size in column direction
518   //
519
520   Float_t cpadW  = AliTRDgeometry::CpadW();
521
522   fColPadSize[p] = s;
523   fCol0[p]       = - fGeo->GetChamberWidth(p)/2. + cpadW;
524   fColMax[p]     = ((Int_t) ((fGeo->GetChamberWidth(p) - 2.*cpadW) / s));
525
526 }
527
528 //_____________________________________________________________________________
529 void AliTRDparameter::SetNTimeBin(const Int_t nbin)
530 {
531   //
532   // Redefines the number of time bins in the drift region.
533   // The time bin width is defined by the length of the
534   // drift region divided by <nbin>.
535   //
536
537   fTimeMax     = nbin;
538   fTimeBinSize = AliTRDgeometry::DrThick() / ((Float_t) fTimeMax);
539   for (Int_t iplan = 0; iplan < AliTRDgeometry::Nplan(); iplan++) {
540     fTime0[iplan] = AliTRDgeometry::Rmin()
541                   + AliTRDgeometry::CraHght()
542                   + AliTRDgeometry::CdrHght()
543                   + iplan * (AliTRDgeometry::Cheight() 
544                            + AliTRDgeometry::Cspace());
545   }
546
547 }
548
549 //_____________________________________________________________________________
550 Float_t AliTRDparameter::CrossTalk(Float_t time) const
551 {
552   //
553   // Applies the pad-pad capacitive cross talk
554   //
555
556   Int_t iBin = ((Int_t) ((time - fTRFlo) / fTRFwid)); 
557   if ((iBin >= 0) && (iBin < fTRFbin)) {
558     return fCTsmp[iBin];
559   }
560   else {
561     return 0.0;
562   }    
563
564 }
565
566 //_____________________________________________________________________________
567 Int_t AliTRDparameter::Diffusion(Float_t driftlength, Float_t *xyz)
568 {
569   //
570   // Applies the diffusion smearing to the position of a single electron
571   //
572
573   Float_t driftSqrt = TMath::Sqrt(driftlength);
574   Float_t sigmaT = driftSqrt * fDiffusionT;
575   Float_t sigmaL = driftSqrt * fDiffusionL;
576   xyz[0] = gRandom->Gaus(xyz[0], sigmaL * fLorentzFactor);
577   xyz[1] = gRandom->Gaus(xyz[1], sigmaT * fLorentzFactor);
578   xyz[2] = gRandom->Gaus(xyz[2], sigmaT);
579
580   return 1;
581
582 }
583
584 //_____________________________________________________________________________
585 Int_t AliTRDparameter::ExB(Float_t driftlength, Float_t *xyz) const
586 {
587   //
588   // Applies E x B effects to the position of a single electron
589   //
590
591   xyz[0] = xyz[0];
592   xyz[1] = xyz[1] + fOmegaTau * driftlength;
593   xyz[2] = xyz[2];
594
595   return 1;
596
597 }
598
599 //_____________________________________________________________________________
600 Int_t AliTRDparameter::PadResponse(Float_t signal, Float_t dist
601                                  , Int_t plane, Float_t *pad) const
602 {
603   //
604   // Applies the pad response
605   //
606
607   const Int_t kNplan = AliTRDgeometry::kNplan;
608
609   Int_t iBin  = ((Int_t) (( - dist - fPRFlo) / fPRFwid));
610   Int_t iOff  = plane * fPRFbin;
611
612   Int_t iBin0 = iBin - fPRFpad + iOff;
613   Int_t iBin1 = iBin           + iOff;
614   Int_t iBin2 = iBin + fPRFpad + iOff;
615
616   pad[0] = 0.0;
617   pad[1] = 0.0;
618   pad[2] = 0.0;
619   if ((iBin1 >= 0) && (iBin1 < (fPRFbin*kNplan))) {
620
621     if (iBin0 >= 0) {
622       pad[0] = signal * fPRFsmp[iBin0];
623     }
624     pad[1] = signal * fPRFsmp[iBin1];
625     if (iBin2 < (fPRFbin*kNplan)) {
626       pad[2] = signal * fPRFsmp[iBin2];
627     }
628
629     return 1;
630
631   }
632   else {
633
634     return 0;
635
636   }
637
638 }
639
640 //_____________________________________________________________________________
641 Float_t AliTRDparameter::TimeResponse(Float_t time) const
642 {
643   //
644   // Applies the preamp shaper time response
645   //
646
647   Int_t iBin = ((Int_t) ((time - fTRFlo) / fTRFwid)); 
648   if ((iBin >= 0) && (iBin < fTRFbin)) {
649     return fTRFsmp[iBin];
650   }
651   else {
652     return 0.0;
653   }    
654
655 }
656
657 //_____________________________________________________________________________
658 Float_t AliTRDparameter::Col0Tilted(Float_t col0, Float_t rowOffset
659                                   , Int_t plane)
660 {
661   //
662   // Calculates col0 for tilted pads
663   //
664
665   Float_t diff = fTiltingAngle * rowOffset;
666   return (col0 + TMath::Power(-1.0,plane) * diff);
667
668 }
669
670 //_____________________________________________________________________________
671 void AliTRDparameter::SampleTRF()
672 {
673   //
674   // Samples the time response function
675   //
676   // New TRF from Venelin Angelov, simulated with CADENCE
677   // Pad-ground capacitance = 25 pF
678   // Pad-pad cross talk capacitance = 6 pF   
679   //
680
681   Int_t   ipos1;
682   Int_t   ipos2;
683   Float_t diff;
684
685   const Int_t kNpasa     = 252;
686
687   Float_t time[kNpasa]   = { -0.220000, -0.210000, -0.200000, -0.190000 
688                            , -0.180000, -0.170000, -0.160000, -0.150000 
689                            , -0.140000, -0.130000, -0.120000, -0.110000 
690                            , -0.100000, -0.090000, -0.080000, -0.070000 
691                            , -0.060000, -0.050000, -0.040000, -0.030000 
692                            , -0.020000, -0.010000, -0.000000,  0.010000 
693                            ,  0.020000,  0.030000,  0.040000,  0.050000 
694                            ,  0.060000,  0.070000,  0.080000,  0.090000 
695                            ,  0.100000,  0.110000,  0.120000,  0.130000 
696                            ,  0.140000,  0.150000,  0.160000,  0.170000 
697                            ,  0.180000,  0.190000,  0.200000,  0.210000 
698                            ,  0.220000,  0.230000,  0.240000,  0.250000 
699                            ,  0.260000,  0.270000,  0.280000,  0.290000 
700                            ,  0.300000,  0.310000,  0.320000,  0.330000 
701                            ,  0.340000,  0.350000,  0.360000,  0.370000 
702                            ,  0.380000,  0.390000,  0.400000,  0.410000 
703                            ,  0.420000,  0.430000,  0.440000,  0.450000 
704                            ,  0.460000,  0.470000,  0.480000,  0.490000 
705                            ,  0.500000,  0.510000,  0.520000,  0.530000 
706                            ,  0.540000,  0.550000,  0.560000,  0.570000 
707                            ,  0.580000,  0.590000,  0.600000,  0.610000 
708                            ,  0.620000,  0.630000,  0.640000,  0.650000 
709                            ,  0.660000,  0.670000,  0.680000,  0.690000 
710                            ,  0.700000,  0.710000,  0.720000,  0.730000 
711                            ,  0.740000,  0.750000,  0.760000,  0.770000 
712                            ,  0.780000,  0.790000,  0.800000,  0.810000 
713                            ,  0.820000,  0.830000,  0.840000,  0.850000 
714                            ,  0.860000,  0.870000,  0.880000,  0.890000 
715                            ,  0.900000,  0.910000,  0.920000,  0.930000 
716                            ,  0.940000,  0.950000,  0.960000,  0.970000 
717                            ,  0.980000,  0.990000,  1.000000,  1.010000 
718                            ,  1.020000,  1.030000,  1.040000,  1.050000 
719                            ,  1.060000,  1.070000,  1.080000,  1.090000 
720                            ,  1.100000,  1.110000,  1.120000,  1.130000 
721                            ,  1.140000,  1.150000,  1.160000,  1.170000 
722                            ,  1.180000,  1.190000,  1.200000,  1.210000 
723                            ,  1.220000,  1.230000,  1.240000,  1.250000 
724                            ,  1.260000,  1.270000,  1.280000,  1.290000 
725                            ,  1.300000,  1.310000,  1.320000,  1.330000 
726                            ,  1.340000,  1.350000,  1.360000,  1.370000 
727                            ,  1.380000,  1.390000,  1.400000,  1.410000 
728                            ,  1.420000,  1.430000,  1.440000,  1.450000 
729                            ,  1.460000,  1.470000,  1.480000,  1.490000 
730                            ,  1.500000,  1.510000,  1.520000,  1.530000 
731                            ,  1.540000,  1.550000,  1.560000,  1.570000 
732                            ,  1.580000,  1.590000,  1.600000,  1.610000 
733                            ,  1.620000,  1.630000,  1.640000,  1.650000 
734                            ,  1.660000,  1.670000,  1.680000,  1.690000 
735                            ,  1.700000,  1.710000,  1.720000,  1.730000 
736                            ,  1.740000,  1.750000,  1.760000,  1.770000 
737                            ,  1.780000,  1.790000,  1.800000,  1.810000 
738                            ,  1.820000,  1.830000,  1.840000,  1.850000 
739                            ,  1.860000,  1.870000,  1.880000,  1.890000 
740                            ,  1.900000,  1.910000,  1.920000,  1.930000 
741                            ,  1.940000,  1.950000,  1.960000,  1.970000 
742                            ,  1.980000,  1.990000,  2.000000,  2.010000 
743                            ,  2.020000,  2.030000,  2.040000,  2.050000 
744                            ,  2.060000,  2.070000,  2.080000,  2.090000 
745                            ,  2.100000,  2.110000,  2.120000,  2.130000 
746                            ,  2.140000,  2.150000,  2.160000,  2.170000 
747                            ,  2.180000,  2.190000,  2.200000,  2.210000 
748                            ,  2.220000,  2.230000,  2.240000,  2.250000 
749                            ,  2.260000,  2.270000,  2.280000,  2.290000 };
750
751   Float_t signal[kNpasa] = {  0.000000,  0.000000,  0.000000,  0.000000 
752                            ,  0.000000,  0.000000,  0.000000,  0.000396 
753                            ,  0.005096,  0.022877,  0.061891,  0.126614 
754                            ,  0.215798,  0.324406,  0.444507,  0.566817 
755                            ,  0.683465,  0.787089,  0.873159,  0.937146 
756                            ,  0.979049,  0.999434,  1.000000,  0.983579 
757                            ,  0.954134,  0.913364,  0.866365,  0.813703 
758                            ,  0.759910,  0.706116,  0.653454,  0.603624 
759                            ,  0.556625,  0.514156,  0.475085,  0.439977 
760                            ,  0.408834,  0.380578,  0.355549,  0.333352 
761                            ,  0.313647,  0.296093,  0.280351,  0.266195 
762                            ,  0.253397,  0.241789,  0.231257,  0.221574 
763                            ,  0.212627,  0.204417,  0.196772,  0.189581 
764                            ,  0.182956,  0.176784,  0.171008,  0.165515 
765                            ,  0.160419,  0.155606,  0.151076,  0.146716 
766                            ,  0.142639,  0.138845,  0.135221,  0.131767 
767                            ,  0.128482,  0.125368,  0.122424,  0.119592 
768                            ,  0.116931,  0.114326,  0.111891,  0.109513 
769                            ,  0.107248,  0.105096,  0.103058,  0.101019 
770                            ,  0.099151,  0.097282,  0.095527,  0.093715 
771                            ,  0.092129,  0.090544,  0.088958,  0.087429 
772                            ,  0.086014,  0.084598,  0.083239,  0.081880 
773                            ,  0.080634,  0.079388,  0.078143,  0.077010 
774                            ,  0.075878,  0.074745,  0.073669,  0.072593 
775                            ,  0.071574,  0.070612,  0.069649,  0.068686 
776                            ,  0.067780,  0.066874,  0.066025,  0.065176 
777                            ,  0.064326,  0.063533,  0.062684,  0.061948 
778                            ,  0.061212,  0.060419,  0.059740,  0.059003 
779                            ,  0.058324,  0.057644,  0.057022,  0.056342 
780                            ,  0.055663,  0.055096,  0.054473,  0.053851 
781                            ,  0.053284,  0.052718,  0.052152,  0.051585 
782                            ,  0.051019,  0.050566,  0.050000,  0.049490 
783                            ,  0.048981,  0.048528,  0.048018,  0.047508 
784                            ,  0.047055,  0.046602,  0.046149,  0.045696 
785                            ,  0.045300,  0.044904,  0.044451,  0.044054 
786                            ,  0.043658,  0.043205,  0.042865,  0.042469 
787                            ,  0.042072,  0.041733,  0.041336,  0.040997 
788                            ,  0.040657,  0.040260,  0.039921,  0.039581 
789                            ,  0.039241,  0.038958,  0.038618,  0.038335 
790                            ,  0.037995,  0.037656,  0.037373,  0.037089 
791                            ,  0.036806,  0.036467,  0.036183,  0.035900 
792                            ,  0.035617,  0.035334,  0.035108,  0.034824 
793                            ,  0.034541,  0.034315,  0.034032,  0.033805 
794                            ,  0.033522,  0.033296,  0.033069,  0.032786 
795                            ,  0.032559,  0.032333,  0.032106,  0.031880 
796                            ,  0.031653,  0.031427,  0.031200,  0.030974 
797                            ,  0.030804,  0.030578,  0.030351,  0.030125 
798                            ,  0.029955,  0.029785,  0.029558,  0.029332 
799                            ,  0.029162,  0.028992,  0.028766,  0.028596 
800                            ,  0.028426,  0.028199,  0.028086,  0.027860 
801                            ,  0.027746,  0.027633,  0.027463,  0.027293 
802                            ,  0.027180,  0.027067,  0.026954,  0.026954 
803                            ,  0.026840,  0.026727,  0.026727,  0.026614 
804                            ,  0.026614,  0.026614,  0.026557,  0.026501 
805                            ,  0.026501,  0.026501,  0.026501,  0.026501 
806                            ,  0.026501,  0.026501,  0.026501,  0.026387 
807                            ,  0.026387,  0.026387,  0.026387,  0.026387 
808                            ,  0.026387,  0.026387,  0.026387,  0.026387 
809                            ,  0.026387,  0.026387,  0.026387,  0.026387 
810                            ,  0.026387,  0.026274,  0.026274,  0.026274 
811                            ,  0.026274,  0.026274,  0.026274,  0.026274 
812                            ,  0.026274,  0.026274,  0.026274,  0.026274 
813                            ,  0.026274,  0.026274,  0.026274,  0.026161 };
814
815   Float_t xtalk[kNpasa]  = {  0.000000,  0.000000,  0.000000,  0.000000 
816                            ,  0.000000,  0.000000,  0.000000,  0.000113 
817                            ,  0.000793,  0.003058,  0.007305,  0.013194 
818                            ,  0.019706,  0.025821,  0.030634,  0.033465 
819                            ,  0.034145,  0.032729,  0.029615,  0.025198 
820                            ,  0.019989,  0.014496,  0.009003,  0.003964 
821                            , -0.000510, -0.004190, -0.007191, -0.009400 
822                            , -0.010872, -0.011835, -0.012288, -0.012288 
823                            , -0.012005, -0.011495, -0.010872, -0.010136 
824                            , -0.009343, -0.008607, -0.007871, -0.007191 
825                            , -0.006512, -0.005946, -0.005379, -0.004926 
826                            , -0.004473, -0.004077, -0.003737, -0.003398 
827                            , -0.003114, -0.002831, -0.002605, -0.002378 
828                            , -0.002208, -0.002039, -0.001869, -0.001699 
829                            , -0.001585, -0.001472, -0.001359, -0.001246 
830                            , -0.001132, -0.001019, -0.001019, -0.000906 
831                            , -0.000906, -0.000793, -0.000793, -0.000680 
832                            , -0.000680, -0.000680, -0.000566, -0.000566 
833                            , -0.000566, -0.000566, -0.000453, -0.000453 
834                            , -0.000453, -0.000453, -0.000453, -0.000453 
835                            , -0.000340, -0.000340, -0.000340, -0.000340 
836                            , -0.000340, -0.000340, -0.000340, -0.000340 
837                            , -0.000340, -0.000340, -0.000340, -0.000340 
838                            , -0.000340, -0.000227, -0.000227, -0.000227 
839                            , -0.000227, -0.000227, -0.000227, -0.000227 
840                            , -0.000227, -0.000227, -0.000227, -0.000227 
841                            , -0.000227, -0.000227, -0.000227, -0.000227 
842                            , -0.000227, -0.000227, -0.000227, -0.000227 
843                            , -0.000227, -0.000227, -0.000227, -0.000227 
844                            , -0.000227, -0.000227, -0.000227, -0.000227 
845                            , -0.000227, -0.000227, -0.000227, -0.000227 
846                            , -0.000227, -0.000227, -0.000227, -0.000227 
847                            , -0.000227, -0.000227, -0.000227, -0.000113 
848                            , -0.000113, -0.000113, -0.000113, -0.000113 
849                            , -0.000113, -0.000113, -0.000113, -0.000113 
850                            , -0.000113, -0.000113, -0.000113, -0.000113 
851                            , -0.000113, -0.000113, -0.000113, -0.000113 
852                            , -0.000113, -0.000113, -0.000113, -0.000113 
853                            , -0.000113, -0.000113, -0.000113, -0.000113 
854                            , -0.000113, -0.000113, -0.000113, -0.000113 
855                            , -0.000113, -0.000113, -0.000113, -0.000113 
856                            , -0.000113, -0.000113, -0.000113, -0.000113 
857                            , -0.000113, -0.000113, -0.000113, -0.000113 
858                            , -0.000113, -0.000113, -0.000113, -0.000113 
859                            , -0.000113, -0.000113, -0.000113, -0.000113 
860                            , -0.000113, -0.000113, -0.000113, -0.000113 
861                            , -0.000113, -0.000113, -0.000113, -0.000113 
862                            , -0.000113, -0.000113, -0.000113, -0.000113 
863                            , -0.000113, -0.000113, -0.000113, -0.000113 
864                            , -0.000113, -0.000113, -0.000113, -0.000113 
865                            , -0.000113, -0.000113, -0.000113, -0.000113 
866                            , -0.000113, -0.000113, -0.000113, -0.000113 
867                            , -0.000113, -0.000113, -0.000113,  0.000000 
868                            ,  0.000000,  0.000000,  0.000000,  0.000000 
869                            ,  0.000000,  0.000000,  0.000000,  0.000000 
870                            ,  0.000000,  0.000000,  0.000000,  0.000000 
871                            ,  0.000000,  0.000000,  0.000000,  0.000000 
872                            ,  0.000000,  0.000000,  0.000000,  0.000000 
873                            ,  0.000000,  0.000000,  0.000000,  0.000000 
874                            ,  0.000000,  0.000000,  0.000000,  0.000000 
875                            ,  0.000000,  0.000000,  0.000000,  0.000000 
876                            ,  0.000000,  0.000000,  0.000000,  0.000000 
877                            ,  0.000000,  0.000000,  0.000000,  0.000000 };
878
879   // increase CrossTalk to measurements
880   for (Int_t ipasa = 0; ipasa < kNpasa; ipasa++) {
881     xtalk[ipasa] *= 1.75;
882   }
883
884   if (fTRFsmp) delete [] fTRFsmp;
885   fTRFsmp = new Float_t[fTRFbin];
886   if (fCTsmp)  delete [] fCTsmp;
887   fCTsmp  = new Float_t[fTRFbin];
888
889   Float_t loTRF    = TMath::Max(fTRFlo / fDriftVelocity,time[0]);
890   Float_t hiTRF    = TMath::Min(fTRFhi / fDriftVelocity,time[kNpasa-1]);
891   Float_t binWidth = (hiTRF - loTRF) / ((Float_t) fTRFbin);
892
893   // Take the linear interpolation
894   for (Int_t iBin = 0; iBin < fTRFbin; iBin++) {
895
896     Float_t bin = (((Float_t) iBin) + 0.5) * binWidth + loTRF;
897     ipos1 = ipos2 = 0;
898     diff  = 0;
899     do {
900       diff = bin - time[ipos2++];
901     } while (diff > 0);
902     ipos2--;
903     if (ipos2 >= kNpasa) ipos2 = kNpasa - 1;
904     ipos1 = ipos2 - 1;
905
906     fTRFsmp[iBin] = signal[ipos2] 
907                   + diff * (signal[ipos2] - signal[ipos1]) 
908                          / (  time[ipos2] -   time[ipos1]);
909
910     fCTsmp[iBin]  = xtalk[ipos2] 
911                   + diff * (xtalk[ipos2]  -  xtalk[ipos1]) 
912                          / (  time[ipos2] -   time[ipos1]);
913
914   }
915
916 }
917
918 //_____________________________________________________________________________
919 void AliTRDparameter::SamplePRF()
920 {
921   //
922   // Samples the pad response function
923   //
924
925   const Int_t kNplan  = AliTRDgeometry::kNplan;
926   const Int_t kPRFbin = 61;
927
928   Float_t prf[kNplan][kPRFbin] = { { 0.018570, 0.022270, 0.026710, 0.032010
929                                    , 0.038350, 0.045920, 0.054930, 0.065650
930                                    , 0.078370, 0.093420, 0.111150, 0.131940
931                                    , 0.156160, 0.184160, 0.216220, 0.252470
932                                    , 0.292860, 0.337030, 0.384330, 0.433750
933                                    , 0.484010, 0.533630, 0.581150, 0.625200
934                                    , 0.664710, 0.698860, 0.727130, 0.749230
935                                    , 0.765050, 0.774540, 0.777700, 0.774540
936                                    , 0.765050, 0.749230, 0.727130, 0.698860
937                                    , 0.664710, 0.625200, 0.581150, 0.533630
938                                    , 0.484010, 0.433750, 0.384330, 0.337030
939                                    , 0.292860, 0.252470, 0.216220, 0.184160
940                                    , 0.156160, 0.131940, 0.111150, 0.093420
941                                    , 0.078370, 0.065650, 0.054930, 0.045920
942                                    , 0.038350, 0.032010, 0.026710, 0.022270
943                                    , 0.018570                               }
944                                  , { 0.015730, 0.019040, 0.023030, 0.027840
945                                    , 0.033650, 0.040650, 0.049060, 0.059160
946                                    , 0.071260, 0.085710, 0.102910, 0.123270
947                                    , 0.147240, 0.175220, 0.207590, 0.244540
948                                    , 0.286090, 0.331920, 0.381350, 0.433290
949                                    , 0.486290, 0.538710, 0.588870, 0.635280
950                                    , 0.676760, 0.712460, 0.741890, 0.764810
951                                    , 0.781150, 0.790930, 0.794180, 0.790930
952                                    , 0.781150, 0.764810, 0.741890, 0.712460
953                                    , 0.676760, 0.635280, 0.588870, 0.538710
954                                    , 0.486290, 0.433290, 0.381350, 0.331920
955                                    , 0.286090, 0.244540, 0.207590, 0.175220
956                                    , 0.147240, 0.123270, 0.102910, 0.085710
957                                    , 0.071260, 0.059160, 0.049060, 0.040650
958                                    , 0.033650, 0.027840, 0.023030, 0.019040
959                                    , 0.015730                               }
960                                  , { 0.013330, 0.016270, 0.019850, 0.024210
961                                    , 0.029510, 0.035960, 0.043790, 0.053280
962                                    , 0.064740, 0.078580, 0.095190, 0.115070
963                                    , 0.138700, 0.166570, 0.199120, 0.236660
964                                    , 0.279260, 0.326660, 0.378140, 0.432540
965                                    , 0.488260, 0.543440, 0.596200, 0.644900
966                                    , 0.688240, 0.725380, 0.755840, 0.779470
967                                    , 0.796260, 0.806280, 0.809610, 0.806280
968                                    , 0.796260, 0.779470, 0.755840, 0.725380
969                                    , 0.688240, 0.644900, 0.596200, 0.543440
970                                    , 0.488260, 0.432540, 0.378140, 0.326660
971                                    , 0.279260, 0.236660, 0.199120, 0.166570
972                                    , 0.138700, 0.115070, 0.095190, 0.078580
973                                    , 0.064740, 0.053280, 0.043790, 0.035960
974                                    , 0.029510, 0.024210, 0.019850, 0.016270
975                                    , 0.013330                               }
976                                  , { 0.011280, 0.013890, 0.017090, 0.021030
977                                    , 0.025870, 0.031800, 0.039060, 0.047940
978                                    , 0.058790, 0.071980, 0.087990, 0.107330
979                                    , 0.130550, 0.158220, 0.190850, 0.228870
980                                    , 0.272410, 0.321270, 0.374740, 0.431560
981                                    , 0.489960, 0.547870, 0.603180, 0.654080
982                                    , 0.699190, 0.737640, 0.769030, 0.793260
983                                    , 0.810410, 0.820620, 0.824010, 0.820620
984                                    , 0.810410, 0.793260, 0.769030, 0.737640
985                                    , 0.699190, 0.654080, 0.603180, 0.547870
986                                    , 0.489960, 0.431560, 0.374740, 0.321270
987                                    , 0.272410, 0.228870, 0.190850, 0.158220
988                                    , 0.130550, 0.107330, 0.087990, 0.071980
989                                    , 0.058790, 0.047940, 0.039060, 0.031800
990                                    , 0.025870, 0.021030, 0.017090, 0.013890
991                                    , 0.011280                               }
992                                  , { 0.009550, 0.011860, 0.014720, 0.018270
993                                    , 0.022660, 0.028100, 0.034820, 0.043120
994                                    , 0.053340, 0.065900, 0.081280, 0.100040
995                                    , 0.122800, 0.150180, 0.182800, 0.221170
996                                    , 0.265550, 0.315790, 0.371180, 0.430370
997                                    , 0.491430, 0.552030, 0.609840, 0.662860
998                                    , 0.709630, 0.749290, 0.781490, 0.806220
999                                    , 0.823650, 0.834000, 0.837430, 0.834000
1000                                    , 0.823650, 0.806220, 0.781490, 0.749290
1001                                    , 0.709630, 0.662860, 0.609840, 0.552030
1002                                    , 0.491430, 0.430370, 0.371180, 0.315790
1003                                    , 0.265550, 0.221170, 0.182800, 0.150180
1004                                    , 0.122800, 0.100040, 0.081280, 0.065900
1005                                    , 0.053340, 0.043120, 0.034820, 0.028100
1006                                    , 0.022660, 0.018270, 0.014720, 0.011860
1007                                    , 0.009550                               }
1008                                  , { 0.008080, 0.010120, 0.012670, 0.015860
1009                                    , 0.019840, 0.024820, 0.031030, 0.038760
1010                                    , 0.048370, 0.060300, 0.075040, 0.093200
1011                                    , 0.115430, 0.142450, 0.174980, 0.213610
1012                                    , 0.258720, 0.310250, 0.367480, 0.429010
1013                                    , 0.492690, 0.555950, 0.616210, 0.671280
1014                                    , 0.719600, 0.760350, 0.793250, 0.818380
1015                                    , 0.836020, 0.846460, 0.849920, 0.846460
1016                                    , 0.836020, 0.818380, 0.793250, 0.760350
1017                                    , 0.719600, 0.671280, 0.616210, 0.555950
1018                                    , 0.492690, 0.429010, 0.367480, 0.310250
1019                                    , 0.258720, 0.213610, 0.174980, 0.142450
1020                                    , 0.115430, 0.093200, 0.075040, 0.060300
1021                                    , 0.048370, 0.038760, 0.031030, 0.024820
1022                                    , 0.019840, 0.015860, 0.012670, 0.010120
1023                                    , 0.008080                               } };
1024
1025   // More sampling precision with linear interpolation
1026   fPRFlo  = -1.5;
1027   fPRFhi  =  1.5;
1028   Float_t pad[kPRFbin];
1029   Int_t   sPRFbin = kPRFbin;  
1030   Float_t sPRFwid = (fPRFhi - fPRFlo) / ((Float_t) sPRFbin);
1031   for (Int_t iPad = 0; iPad < sPRFbin; iPad++) {
1032     pad[iPad] = ((Float_t) iPad + 0.5) * sPRFwid + fPRFlo;
1033   }
1034   fPRFbin = 500;  
1035   fPRFwid = (fPRFhi - fPRFlo) / ((Float_t) fPRFbin);
1036   fPRFpad = ((Int_t) (1.0 / fPRFwid));
1037
1038   if (fPRFsmp) delete [] fPRFsmp;
1039   fPRFsmp = new Float_t[kNplan*fPRFbin];
1040
1041   Int_t   ipos1;
1042   Int_t   ipos2;
1043   Float_t diff;
1044
1045   for (Int_t iPla = 0; iPla < kNplan; iPla++) {
1046
1047     for (Int_t iBin = 0; iBin < fPRFbin; iBin++) {
1048
1049       Float_t bin = (((Float_t) iBin) + 0.5) * fPRFwid + fPRFlo;
1050       ipos1 = ipos2 = 0;
1051       diff  = 0;
1052       do {
1053         diff = bin - pad[ipos2++];
1054       } while ((diff > 0) && (ipos2 < kPRFbin));
1055       if      (ipos2 == kPRFbin) {
1056         fPRFsmp[iPla*fPRFbin+iBin] = prf[iPla][ipos2-1];
1057       }
1058       else if (ipos2 == 1) {
1059         fPRFsmp[iPla*fPRFbin+iBin] = prf[iPla][ipos2-1];
1060       }
1061       else {
1062         ipos2--;
1063         if (ipos2 >= kPRFbin) ipos2 = kPRFbin - 1;
1064         ipos1 = ipos2 - 1;
1065         fPRFsmp[iPla*fPRFbin+iBin] = prf[iPla][ipos2] 
1066                                    + diff * (prf[iPla][ipos2] - prf[iPla][ipos1]) 
1067                                           / sPRFwid;
1068       }
1069
1070     }
1071   } 
1072
1073 }
1074
1075 //_____________________________________________________________________________
1076 void AliTRDparameter::FillLUT()
1077 {
1078   //
1079   // Create the LUT
1080   //
1081
1082   const Int_t kNplan = AliTRDgeometry::kNplan;
1083   const Int_t kNlut  = 128;
1084
1085   fLUTbin = kNplan * kNlut;
1086
1087   // The lookup table from Bogdan
1088   Float_t lut[kNplan][kNlut] = {  
1089     {
1090       0.0070, 0.0150, 0.0224, 0.0298, 0.0374, 0.0454, 0.0533, 0.0611, 
1091       0.0684, 0.0755, 0.0827, 0.0900, 0.0975, 0.1049, 0.1120, 0.1187, 
1092       0.1253, 0.1318, 0.1385, 0.1453, 0.1519, 0.1584, 0.1646, 0.1704, 
1093       0.1762, 0.1821, 0.1879, 0.1938, 0.1996, 0.2053, 0.2108, 0.2160, 
1094       0.2210, 0.2260, 0.2310, 0.2361, 0.2411, 0.2461, 0.2509, 0.2557, 
1095       0.2602, 0.2646, 0.2689, 0.2732, 0.2774, 0.2816, 0.2859, 0.2901, 
1096       0.2942, 0.2983, 0.3022, 0.3061, 0.3099, 0.3136, 0.3172, 0.3207, 
1097       0.3242, 0.3278, 0.3312, 0.3347, 0.3382, 0.3416, 0.3450, 0.3483, 
1098       0.3515, 0.3547, 0.3579, 0.3609, 0.3639, 0.3669, 0.3698, 0.3727, 
1099       0.3756, 0.3785, 0.3813, 0.3842, 0.3870, 0.3898, 0.3926, 0.3952, 
1100       0.3979, 0.4005, 0.4032, 0.4057, 0.4082, 0.4108, 0.4132, 0.4157, 
1101       0.4181, 0.4205, 0.4228, 0.4252, 0.4275, 0.4299, 0.4322, 0.4345, 
1102       0.4367, 0.4390, 0.4412, 0.4434, 0.4456, 0.4478, 0.4499, 0.4520, 
1103       0.4541, 0.4562, 0.4583, 0.4603, 0.4623, 0.4643, 0.4663, 0.4683, 
1104       0.4702, 0.4722, 0.4741, 0.4758, 0.4774, 0.4790, 0.4805, 0.4824, 
1105       0.4844, 0.4863, 0.4883, 0.4902, 0.4921, 0.4940, 0.4959, 0.4978 
1106     },
1107     {
1108       0.0072, 0.0156, 0.0235, 0.0313, 0.0394, 0.0478, 0.0561, 0.0642, 
1109       0.0718, 0.0792, 0.0868, 0.0947, 0.1025, 0.1101, 0.1172, 0.1241, 
1110       0.1309, 0.1378, 0.1449, 0.1518, 0.1586, 0.1650, 0.1710, 0.1770, 
1111       0.1830, 0.1891, 0.1952, 0.2011, 0.2070, 0.2125, 0.2177, 0.2229, 
1112       0.2280, 0.2332, 0.2383, 0.2435, 0.2484, 0.2533, 0.2581, 0.2627, 
1113       0.2670, 0.2714, 0.2757, 0.2799, 0.2842, 0.2884, 0.2927, 0.2968, 
1114       0.3008, 0.3048, 0.3086, 0.3123, 0.3159, 0.3195, 0.3231, 0.3266, 
1115       0.3301, 0.3335, 0.3370, 0.3404, 0.3438, 0.3471, 0.3504, 0.3536, 
1116       0.3567, 0.3598, 0.3628, 0.3657, 0.3686, 0.3715, 0.3744, 0.3772, 
1117       0.3800, 0.3828, 0.3856, 0.3884, 0.3911, 0.3938, 0.3965, 0.3991, 
1118       0.4016, 0.4042, 0.4067, 0.4092, 0.4116, 0.4140, 0.4164, 0.4187, 
1119       0.4211, 0.4234, 0.4257, 0.4280, 0.4302, 0.4325, 0.4347, 0.4369, 
1120       0.4391, 0.4413, 0.4434, 0.4456, 0.4477, 0.4497, 0.4518, 0.4538, 
1121       0.4558, 0.4578, 0.4598, 0.4618, 0.4637, 0.4656, 0.4675, 0.4694, 
1122       0.4713, 0.4732, 0.4750, 0.4766, 0.4781, 0.4797, 0.4813, 0.4832, 
1123       0.4851, 0.4870, 0.4888, 0.4906, 0.4925, 0.4942, 0.4960, 0.4978
1124     },
1125     {
1126       0.0075, 0.0163, 0.0246, 0.0328, 0.0415, 0.0504, 0.0592, 0.0674, 
1127       0.0753, 0.0832, 0.0914, 0.0996, 0.1077, 0.1154, 0.1225, 0.1296, 
1128       0.1369, 0.1442, 0.1515, 0.1585, 0.1652, 0.1714, 0.1776, 0.1839, 
1129       0.1902, 0.1965, 0.2025, 0.2085, 0.2141, 0.2194, 0.2247, 0.2299, 
1130       0.2352, 0.2405, 0.2457, 0.2507, 0.2557, 0.2604, 0.2649, 0.2693, 
1131       0.2737, 0.2780, 0.2823, 0.2867, 0.2909, 0.2951, 0.2992, 0.3033, 
1132       0.3072, 0.3110, 0.3146, 0.3182, 0.3218, 0.3253, 0.3288, 0.3323, 
1133       0.3357, 0.3392, 0.3426, 0.3459, 0.3492, 0.3524, 0.3555, 0.3586, 
1134       0.3616, 0.3645, 0.3674, 0.3703, 0.3731, 0.3759, 0.3787, 0.3815, 
1135       0.3843, 0.3870, 0.3897, 0.3925, 0.3950, 0.3976, 0.4002, 0.4027, 
1136       0.4052, 0.4076, 0.4101, 0.4124, 0.4148, 0.4171, 0.4194, 0.4217, 
1137       0.4239, 0.4262, 0.4284, 0.4306, 0.4328, 0.4350, 0.4371, 0.4393, 
1138       0.4414, 0.4435, 0.4455, 0.4476, 0.4496, 0.4516, 0.4536, 0.4555, 
1139       0.4575, 0.4594, 0.4613, 0.4632, 0.4650, 0.4669, 0.4687, 0.4705, 
1140       0.4723, 0.4741, 0.4758, 0.4773, 0.4789, 0.4804, 0.4821, 0.4839, 
1141       0.4857, 0.4875, 0.4893, 0.4910, 0.4928, 0.4945, 0.4961, 0.4978
1142     },
1143     {
1144       0.0078, 0.0171, 0.0258, 0.0345, 0.0438, 0.0532, 0.0624, 0.0708, 
1145       0.0791, 0.0875, 0.0962, 0.1048, 0.1130, 0.1206, 0.1281, 0.1356, 
1146       0.1432, 0.1508, 0.1582, 0.1651, 0.1716, 0.1780, 0.1845, 0.1910, 
1147       0.1975, 0.2038, 0.2099, 0.2155, 0.2210, 0.2263, 0.2317, 0.2371, 
1148       0.2425, 0.2477, 0.2528, 0.2578, 0.2626, 0.2671, 0.2715, 0.2759, 
1149       0.2803, 0.2846, 0.2890, 0.2933, 0.2975, 0.3016, 0.3056, 0.3095, 
1150       0.3132, 0.3168, 0.3204, 0.3239, 0.3274, 0.3309, 0.3344, 0.3378, 
1151       0.3412, 0.3446, 0.3479, 0.3511, 0.3543, 0.3574, 0.3603, 0.3633, 
1152       0.3662, 0.3690, 0.3718, 0.3747, 0.3774, 0.3802, 0.3829, 0.3857, 
1153       0.3883, 0.3910, 0.3936, 0.3962, 0.3987, 0.4012, 0.4037, 0.4061, 
1154       0.4085, 0.4109, 0.4132, 0.4155, 0.4177, 0.4200, 0.4222, 0.4244, 
1155       0.4266, 0.4288, 0.4309, 0.4331, 0.4352, 0.4373, 0.4394, 0.4414, 
1156       0.4435, 0.4455, 0.4475, 0.4494, 0.4514, 0.4533, 0.4552, 0.4571, 
1157       0.4590, 0.4608, 0.4626, 0.4645, 0.4662, 0.4680, 0.4698, 0.4715, 
1158       0.4733, 0.4750, 0.4766, 0.4781, 0.4796, 0.4812, 0.4829, 0.4846, 
1159       0.4863, 0.4880, 0.4897, 0.4914, 0.4930, 0.4946, 0.4963, 0.4979
1160     },
1161     {
1162       0.0081, 0.0178, 0.0270, 0.0364, 0.0463, 0.0562, 0.0656, 0.0744, 
1163       0.0831, 0.0921, 0.1013, 0.1102, 0.1183, 0.1261, 0.1339, 0.1419, 
1164       0.1499, 0.1576, 0.1648, 0.1715, 0.1782, 0.1849, 0.1917, 0.1984, 
1165       0.2048, 0.2110, 0.2167, 0.2223, 0.2278, 0.2333, 0.2389, 0.2444, 
1166       0.2497, 0.2548, 0.2598, 0.2645, 0.2691, 0.2735, 0.2780, 0.2824, 
1167       0.2868, 0.2912, 0.2955, 0.2997, 0.3038, 0.3078, 0.3116, 0.3152, 
1168       0.3188, 0.3224, 0.3259, 0.3294, 0.3329, 0.3364, 0.3398, 0.3432, 
1169       0.3465, 0.3497, 0.3529, 0.3561, 0.3591, 0.3620, 0.3649, 0.3677, 
1170       0.3705, 0.3733, 0.3761, 0.3788, 0.3816, 0.3843, 0.3869, 0.3896, 
1171       0.3922, 0.3948, 0.3973, 0.3998, 0.4022, 0.4047, 0.4070, 0.4094, 
1172       0.4117, 0.4139, 0.4162, 0.4184, 0.4206, 0.4227, 0.4249, 0.4270, 
1173       0.4291, 0.4313, 0.4334, 0.4354, 0.4375, 0.4395, 0.4415, 0.4435, 
1174       0.4455, 0.4474, 0.4493, 0.4512, 0.4531, 0.4550, 0.4568, 0.4586, 
1175       0.4604, 0.4622, 0.4639, 0.4657, 0.4674, 0.4691, 0.4708, 0.4725, 
1176       0.4742, 0.4758, 0.4773, 0.4788, 0.4803, 0.4819, 0.4836, 0.4852, 
1177       0.4869, 0.4885, 0.4901, 0.4917, 0.4933, 0.4948, 0.4964, 0.4979
1178     },
1179     {
1180       0.0085, 0.0189, 0.0288, 0.0389, 0.0497, 0.0603, 0.0699, 0.0792, 
1181       0.0887, 0.0985, 0.1082, 0.1170, 0.1253, 0.1336, 0.1421, 0.1505, 
1182       0.1587, 0.1662, 0.1733, 0.1803, 0.1874, 0.1945, 0.2014, 0.2081, 
1183       0.2143, 0.2201, 0.2259, 0.2316, 0.2374, 0.2431, 0.2487, 0.2541, 
1184       0.2593, 0.2642, 0.2689, 0.2735, 0.2781, 0.2826, 0.2872, 0.2917, 
1185       0.2961, 0.3003, 0.3045, 0.3086, 0.3125, 0.3162, 0.3198, 0.3235, 
1186       0.3270, 0.3306, 0.3342, 0.3377, 0.3411, 0.3446, 0.3479, 0.3511, 
1187       0.3543, 0.3575, 0.3605, 0.3634, 0.3663, 0.3691, 0.3720, 0.3748, 
1188       0.3775, 0.3803, 0.3830, 0.3857, 0.3884, 0.3911, 0.3937, 0.3962, 
1189       0.3987, 0.4012, 0.4036, 0.4060, 0.4084, 0.4107, 0.4129, 0.4152, 
1190       0.4174, 0.4196, 0.4218, 0.4239, 0.4261, 0.4282, 0.4303, 0.4324, 
1191       0.4344, 0.4365, 0.4385, 0.4405, 0.4425, 0.4445, 0.4464, 0.4483, 
1192       0.4502, 0.4521, 0.4539, 0.4558, 0.4576, 0.4593, 0.4611, 0.4629, 
1193       0.4646, 0.4663, 0.4680, 0.4697, 0.4714, 0.4730, 0.4747, 0.4759, 
1194       0.4769, 0.4780, 0.4790, 0.4800, 0.4811, 0.4827, 0.4843, 0.4859, 
1195       0.4874, 0.4889, 0.4905, 0.4920, 0.4935, 0.4950, 0.4965, 0.4979
1196     }
1197   }; 
1198
1199   if (fLUT) delete [] fLUT;
1200   fLUT = new Float_t[fLUTbin];
1201
1202   for (Int_t iplan = 0; iplan < kNplan; iplan++) {
1203     for (Int_t ilut  = 0; ilut  < kNlut; ilut++) {
1204       fLUT[iplan*kNplan+ilut] = lut[iplan][ilut];
1205     }
1206   }
1207
1208 }
1209
1210 //_____________________________________________________________________________
1211 void AliTRDparameter::SetTiltingAngle(Float_t v)
1212 {
1213   //
1214   // Set the tilting angle for the readout pads
1215   //
1216
1217   fTiltingAngle = TMath::Tan(TMath::Pi()/180.0 * v);
1218
1219 }
1220
1221 //_____________________________________________________________________________
1222 Float_t AliTRDparameter::GetTiltingAngle() const
1223 {
1224   //
1225   // Get the tilting angle for the readout pads
1226   //
1227
1228   return TMath::ATan(180.0/TMath::Pi() * fTiltingAngle);
1229
1230 }
1231
1232 //_____________________________________________________________________________
1233 Float_t AliTRDparameter::GetDiffusionL(Float_t vd, Float_t b)
1234 {
1235   //
1236   // Returns the longitudinal diffusion coefficient for a given drift 
1237   // velocity <vd> and a B-field <b> for Xe/CO2 (15%).
1238   // The values are according to a GARFIELD simulation.
1239   //
1240
1241   const Int_t kNb = 5;
1242   Float_t p0[kNb] = {  0.007440,  0.007493,  0.007513,  0.007672,  0.007831 };
1243   Float_t p1[kNb] = {  0.019252,  0.018912,  0.018636,  0.018012,  0.017343 };
1244   Float_t p2[kNb] = { -0.005042, -0.004926, -0.004867, -0.004650, -0.004424 };
1245   Float_t p3[kNb] = {  0.000195,  0.000189,  0.000195,  0.000182,  0.000169 };
1246
1247   Int_t ib = ((Int_t) (10 * (b - 0.15)));
1248   ib       = TMath::Max(  0,ib);
1249   ib       = TMath::Min(kNb,ib);
1250
1251   Float_t diff = p0[ib] 
1252                + p1[ib] * vd
1253                + p2[ib] * vd*vd
1254                + p3[ib] * vd*vd*vd;
1255
1256   return diff;
1257
1258 }
1259
1260 //_____________________________________________________________________________
1261 Float_t AliTRDparameter::GetDiffusionT(Float_t vd, Float_t b)
1262 {
1263   //
1264   // Returns the transverse diffusion coefficient for a given drift 
1265   // velocity <vd> and a B-field <b> for Xe/CO2 (15%).
1266   // The values are according to a GARFIELD simulation.
1267   //
1268
1269   const Int_t kNb = 5;
1270   Float_t p0[kNb] = {  0.009550,  0.009599,  0.009674,  0.009757,  0.009850 };
1271   Float_t p1[kNb] = {  0.006667,  0.006539,  0.006359,  0.006153,  0.005925 };
1272   Float_t p2[kNb] = { -0.000853, -0.000798, -0.000721, -0.000635, -0.000541 };
1273   Float_t p3[kNb] = {  0.000131,  0.000122,  0.000111,  0.000098,  0.000085 };
1274
1275   Int_t ib = ((Int_t) (10 * (b - 0.15)));
1276   ib       = TMath::Max(  0,ib);
1277   ib       = TMath::Min(kNb,ib);
1278
1279   Float_t diff = p0[ib] 
1280                + p1[ib] * vd
1281                + p2[ib] * vd*vd
1282                + p3[ib] * vd*vd*vd;
1283
1284   return diff;
1285
1286 }
1287
1288 //_____________________________________________________________________________
1289 Float_t AliTRDparameter::GetOmegaTau(Float_t vd, Float_t b)
1290 {
1291   //
1292   // Returns omega*tau (tan(Lorentz-angle)) for a given drift velocity <vd> 
1293   // and a B-field <b> for Xe/CO2 (15%).
1294   // The values are according to a GARFIELD simulation.
1295   //
1296
1297   const Int_t kNb = 5;
1298   Float_t p0[kNb] = {  0.004810,  0.007412,  0.010252,  0.013409,  0.016888 };
1299   Float_t p1[kNb] = {  0.054875,  0.081534,  0.107333,  0.131983,  0.155455 };
1300   Float_t p2[kNb] = { -0.008682, -0.012896, -0.016987, -0.020880, -0.024623 };
1301   Float_t p3[kNb] = {  0.000155,  0.000238,  0.000330,  0.000428,  0.000541 };
1302
1303   Int_t ib = ((Int_t) (10 * (b - 0.15)));
1304   ib       = TMath::Max(  0,ib);
1305   ib       = TMath::Min(kNb,ib);
1306
1307   Float_t alphaL = p0[ib] 
1308                  + p1[ib] * vd
1309                  + p2[ib] * vd*vd
1310                  + p3[ib] * vd*vd*vd;
1311
1312   return TMath::Tan(alphaL);
1313
1314 }
1315
1316 //_____________________________________________________________________________
1317 Double_t AliTRDparameter::LUTposition(Int_t iplane, Double_t ampL 
1318                                                   , Double_t ampC
1319                                                   , Double_t ampR) const
1320 {
1321   //
1322   // Calculates the cluster position using the lookup table.
1323   // Method provided by Bogdan Vulpescu.
1324   //
1325
1326   const Int_t kNplan = AliTRDgeometry::kNplan;
1327
1328   Double_t pos;
1329   Double_t x = 0.0;
1330   Double_t xmin;
1331   Double_t xmax;
1332   Double_t xwid;
1333
1334   Int_t    side = 0;
1335   Int_t    ix;
1336
1337   Double_t xMin[kNplan] = { 0.006492, 0.006377, 0.006258, 0.006144, 0.006030, 0.005980 };
1338   Double_t xMax[kNplan] = { 0.960351, 0.965870, 0.970445, 0.974352, 0.977667, 0.996101 };
1339
1340   if      (ampL > ampR) {
1341     x    = (ampL - ampR) / ampC;
1342     side = -1;
1343   } 
1344   else if (ampL < ampR) {
1345     x    = (ampR - ampL) / ampC;
1346     side = +1;
1347   }
1348
1349   if (ampL != ampR) {
1350
1351     xmin = xMin[iplane] + 0.000005;
1352     xmax = xMax[iplane] - 0.000005;
1353     xwid = (xmax - xmin) / 127.0;
1354
1355     if      (x < xmin) {
1356       pos = 0.0000;
1357     } 
1358     else if (x > xmax) {
1359       pos = side * 0.5000;
1360     } 
1361     else {
1362       ix  = (Int_t) ((x - xmin) / xwid);
1363       pos = side * fLUT[iplane*kNplan+ix];
1364     }
1365        
1366   } 
1367   else {
1368
1369     pos = 0.0;
1370
1371   }
1372
1373   return pos;
1374
1375 }