]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDdigitizer.cxx
Add detailed geometry and simple simulator
[u/mrichter/AliRoot.git] / TRD / AliTRDdigitizer.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.25  2001/06/27 09:54:44  cblume
19 Moved fField initialization to InitDetector()
20
21 Revision 1.24  2001/05/21 16:45:47  hristov
22 Last minute changes (C.Blume)
23
24 Revision 1.23  2001/05/07 08:04:48  cblume
25 New TRF and PRF. Speedup of the code. Digits from amplification region included
26
27 Revision 1.22  2001/03/30 14:40:14  cblume
28 Update of the digitization parameter
29
30 Revision 1.21  2001/03/13 09:30:35  cblume
31 Update of digitization. Moved digit branch definition to AliTRD
32
33 Revision 1.20  2001/02/25 20:19:00  hristov
34 Minor correction: loop variable declared only once for HP, Sun
35
36 Revision 1.19  2001/02/14 18:22:26  cblume
37 Change in the geometry of the padplane
38
39 Revision 1.18  2001/01/26 19:56:57  hristov
40 Major upgrade of AliRoot code
41
42 Revision 1.17  2000/12/08 12:53:27  cblume
43 Change in Copy() function for HP-compiler
44
45 Revision 1.16  2000/12/07 12:20:46  cblume
46 Go back to array compression. Use sampled PRF to speed up digitization
47
48 Revision 1.15  2000/11/23 14:34:08  cblume
49 Fixed bug in expansion routine of arrays (initialize buffers properly)
50
51 Revision 1.14  2000/11/20 08:54:44  cblume
52 Switch off compression as default
53
54 Revision 1.13  2000/11/10 14:57:52  cblume
55 Changes in the geometry constants for the DEC compiler
56
57 Revision 1.12  2000/11/01 14:53:20  cblume
58 Merge with TRD-develop
59
60 Revision 1.1.4.9  2000/10/26 17:00:22  cblume
61 Fixed bug in CheckDetector()
62
63 Revision 1.1.4.8  2000/10/23 13:41:35  cblume
64 Added protection against Log(0) in the gas gain calulation
65
66 Revision 1.1.4.7  2000/10/17 02:27:34  cblume
67 Get rid of global constants
68
69 Revision 1.1.4.6  2000/10/16 01:16:53  cblume
70 Changed timebin 0 to be the one closest to the readout
71
72 Revision 1.1.4.5  2000/10/15 23:34:29  cblume
73 Faster version of the digitizer
74
75 Revision 1.1.4.4  2000/10/06 16:49:46  cblume
76 Made Getters const
77
78 Revision 1.1.4.3  2000/10/04 16:34:58  cblume
79 Replace include files by forward declarations
80
81 Revision 1.1.4.2  2000/09/22 14:41:10  cblume
82 Bug fix in PRF. Included time response. New structure
83
84 Revision 1.10  2000/10/05 07:27:53  cblume
85 Changes in the header-files by FCA
86
87 Revision 1.9  2000/10/02 21:28:19  fca
88 Removal of useless dependecies via forward declarations
89
90 Revision 1.8  2000/06/09 11:10:07  cblume
91 Compiler warnings and coding conventions, next round
92
93 Revision 1.7  2000/06/08 18:32:58  cblume
94 Make code compliant to coding conventions
95
96 Revision 1.6  2000/06/07 16:27:32  cblume
97 Try to remove compiler warnings on Sun and HP
98
99 Revision 1.5  2000/05/09 16:38:57  cblume
100 Removed PadResponse(). Merge problem
101
102 Revision 1.4  2000/05/08 15:53:45  cblume
103 Resolved merge conflict
104
105 Revision 1.3  2000/04/28 14:49:27  cblume
106 Only one declaration of iDict in MakeDigits()
107
108 Revision 1.1.4.1  2000/05/08 14:42:04  cblume
109 Introduced AliTRDdigitsManager
110
111 Revision 1.1  2000/02/28 19:00:13  cblume
112 Add new TRD classes
113
114 */
115
116 ///////////////////////////////////////////////////////////////////////////////
117 //                                                                           //
118 //  Creates and handles digits from TRD hits                                 //
119 //                                                                           //
120 //  The following effects are included:                                      //
121 //      - Diffusion                                                          //
122 //      - ExB effects                                                        //
123 //      - Gas gain including fluctuations                                    //
124 //      - Pad-response (simple Gaussian approximation)                       //
125 //      - Electronics noise                                                  //
126 //      - Electronics gain                                                   //
127 //      - Digitization                                                       //
128 //      - ADC threshold                                                      //
129 //  The corresponding parameter can be adjusted via the various              //
130 //  Set-functions. If these parameters are not explicitly set, default       //
131 //  values are used (see Init-function).                                     //
132 //  To produce digits from a root-file with TRD-hits use the                 //
133 //  slowDigitsCreate.C macro.                                                //
134 //                                                                           //
135 ///////////////////////////////////////////////////////////////////////////////
136
137 #include <stdlib.h>
138
139 #include <TMath.h>
140 #include <TVector.h>
141 #include <TRandom.h>
142 #include <TROOT.h>
143 #include <TTree.h>
144 #include <TFile.h>
145 #include <TF1.h>
146
147 #include "AliRun.h"
148 #include "AliMagF.h"
149
150 #include "AliTRD.h"
151 #include "AliTRDhit.h"
152 #include "AliTRDdigitizer.h"
153 #include "AliTRDdataArrayI.h"
154 #include "AliTRDdataArrayF.h"
155 #include "AliTRDsegmentArray.h"
156 #include "AliTRDdigitsManager.h"
157 #include "AliTRDgeometry.h"
158
159 ClassImp(AliTRDdigitizer)
160
161 //_____________________________________________________________________________
162 AliTRDdigitizer::AliTRDdigitizer():TNamed()
163 {
164   //
165   // AliTRDdigitizer default constructor
166   //
167
168   fInputFile      = NULL;
169   fDigits         = NULL;
170   fTRD            = NULL;
171   fGeo            = NULL;
172   fPRFsmp         = NULL;
173   fTRFsmp         = NULL;
174
175   fEvent          = 0;
176   fGasGain        = 0.0;
177   fNoise          = 0.0;
178   fChipGain       = 0.0;
179   fSinRange       = 0.0;
180   fSoutRange      = 0.0;
181   fADCoutRange    = 0.0;
182   fADCinRange     = 0.0;
183   fADCthreshold   = 0;
184   fDiffusionOn    = 0;
185   fDiffusionT     = 0.0;
186   fDiffusionL     = 0.0;
187   fElAttachOn     = 0;
188   fElAttachProp   = 0.0;
189   fExBOn          = 0;
190   fOmegaTau       = 0.0;
191   fPRFOn          = 0;
192   fTRFOn          = 0;
193   fDriftVelocity  = 0.0;
194   fPadCoupling    = 0.0;
195   fTimeCoupling   = 0.0;
196   fTimeBinWidth   = 0.0;
197   fField          = 0.0;
198
199   fPRFbin         = 0;
200   fPRFlo          = 0.0;
201   fPRFhi          = 0.0;
202   fPRFwid         = 0.0;
203   fPRFpad         = 0;
204   fTRFbin         = 0;
205   fTRFlo          = 0.0;
206   fTRFhi          = 0.0;
207   fTRFwid         = 0.0;
208
209   fCompress       = kTRUE;
210   fVerbose        = 1;
211   fSDigits        = kFALSE;
212
213 }
214
215 //_____________________________________________________________________________
216 AliTRDdigitizer::AliTRDdigitizer(const Text_t *name, const Text_t *title)
217                 :TNamed(name,title)
218 {
219   //
220   // AliTRDdigitizer default constructor
221   //
222
223   fInputFile     = NULL;
224   fDigits        = NULL;
225   fTRD           = NULL;
226   fGeo           = NULL;
227   fPRFsmp        = NULL;
228   fTRFsmp        = NULL;
229
230   fEvent         = 0;
231
232   fCompress      = kTRUE;
233   fVerbose       = 1;
234   fSDigits       = kFALSE;
235
236   Init();
237
238 }
239
240 //_____________________________________________________________________________
241 AliTRDdigitizer::AliTRDdigitizer(const AliTRDdigitizer &d)
242 {
243   //
244   // AliTRDdigitizer copy constructor
245   //
246
247   ((AliTRDdigitizer &) d).Copy(*this);
248
249 }
250
251 //_____________________________________________________________________________
252 AliTRDdigitizer::~AliTRDdigitizer()
253 {
254   //
255   // AliTRDdigitizer destructor
256   //
257
258   if (fInputFile) {
259     fInputFile->Close();
260     delete fInputFile;
261   }
262
263   if (fDigits) {
264     delete fDigits;
265   }
266
267 }
268
269 //_____________________________________________________________________________
270 AliTRDdigitizer &AliTRDdigitizer::operator=(const AliTRDdigitizer &d)
271 {
272   //
273   // Assignment operator
274   //
275
276   if (this != &d) ((AliTRDdigitizer &) d).Copy(*this);
277   return *this;
278
279 }
280
281 //_____________________________________________________________________________
282 void AliTRDdigitizer::Copy(TObject &d)
283 {
284   //
285   // Copy function
286   //
287
288   Int_t iBin;
289
290   ((AliTRDdigitizer &) d).fInputFile      = NULL;
291   ((AliTRDdigitizer &) d).fDigits         = NULL;
292   ((AliTRDdigitizer &) d).fTRD            = NULL;
293   ((AliTRDdigitizer &) d).fGeo            = NULL;
294
295   ((AliTRDdigitizer &) d).fEvent          = 0;
296
297   ((AliTRDdigitizer &) d).fGasGain        = fGasGain;
298   ((AliTRDdigitizer &) d).fNoise          = fNoise;
299   ((AliTRDdigitizer &) d).fChipGain       = fChipGain;
300   ((AliTRDdigitizer &) d).fSoutRange      = fSoutRange;
301   ((AliTRDdigitizer &) d).fSinRange       = fSinRange;
302   ((AliTRDdigitizer &) d).fADCoutRange    = fADCoutRange;
303   ((AliTRDdigitizer &) d).fADCinRange     = fADCinRange;
304   ((AliTRDdigitizer &) d).fADCthreshold   = fADCthreshold;
305   ((AliTRDdigitizer &) d).fDiffusionOn    = fDiffusionOn; 
306   ((AliTRDdigitizer &) d).fDiffusionT     = fDiffusionT;
307   ((AliTRDdigitizer &) d).fDiffusionL     = fDiffusionL;
308   ((AliTRDdigitizer &) d).fElAttachOn     = fElAttachOn;
309   ((AliTRDdigitizer &) d).fElAttachProp   = fElAttachProp;
310   ((AliTRDdigitizer &) d).fExBOn          = fExBOn;
311   ((AliTRDdigitizer &) d).fOmegaTau       = fOmegaTau;
312   ((AliTRDdigitizer &) d).fLorentzFactor  = fLorentzFactor;
313   ((AliTRDdigitizer &) d).fDriftVelocity  = fDriftVelocity;
314   ((AliTRDdigitizer &) d).fPadCoupling    = fPadCoupling;
315   ((AliTRDdigitizer &) d).fTimeCoupling   = fTimeCoupling;
316   ((AliTRDdigitizer &) d).fTimeBinWidth   = fTimeBinWidth;
317   ((AliTRDdigitizer &) d).fField          = fField;
318   ((AliTRDdigitizer &) d).fPRFOn          = fPRFOn;
319   ((AliTRDdigitizer &) d).fTRFOn          = fTRFOn;
320
321   ((AliTRDdigitizer &) d).fCompress       = fCompress;
322   ((AliTRDdigitizer &) d).fVerbose        = fVerbose;
323   ((AliTRDdigitizer &) d).fSDigits        = fSDigits;
324
325   ((AliTRDdigitizer &) d).fPRFbin         = fPRFbin;
326   ((AliTRDdigitizer &) d).fPRFlo          = fPRFlo;
327   ((AliTRDdigitizer &) d).fPRFhi          = fPRFhi;
328   ((AliTRDdigitizer &) d).fPRFwid         = fPRFwid;
329   ((AliTRDdigitizer &) d).fPRFpad         = fPRFpad;
330   if (((AliTRDdigitizer &) d).fPRFsmp) delete ((AliTRDdigitizer &) d).fPRFsmp;
331   ((AliTRDdigitizer &) d).fPRFsmp = new Float_t[fPRFbin];
332   for (iBin = 0; iBin < fPRFbin; iBin++) {
333     ((AliTRDdigitizer &) d).fPRFsmp[iBin] = fPRFsmp[iBin];
334   }                                                                             
335   ((AliTRDdigitizer &) d).fTRFbin         = fTRFbin;
336   ((AliTRDdigitizer &) d).fTRFlo          = fTRFlo;
337   ((AliTRDdigitizer &) d).fTRFhi          = fTRFhi;
338   ((AliTRDdigitizer &) d).fTRFwid         = fTRFwid;
339   if (((AliTRDdigitizer &) d).fTRFsmp) delete ((AliTRDdigitizer &) d).fTRFsmp;
340   ((AliTRDdigitizer &) d).fTRFsmp = new Float_t[fTRFbin];
341   for (iBin = 0; iBin < fTRFbin; iBin++) {
342     ((AliTRDdigitizer &) d).fTRFsmp[iBin] = fTRFsmp[iBin];
343   }                                      
344                                        
345 }
346
347 //_____________________________________________________________________________
348 Int_t AliTRDdigitizer::Diffusion(Float_t driftlength, Float_t *xyz)
349 {
350   //
351   // Applies the diffusion smearing to the position of a single electron
352   //
353
354   Float_t driftSqrt = TMath::Sqrt(driftlength);
355   Float_t sigmaT = driftSqrt * fDiffusionT;
356   Float_t sigmaL = driftSqrt * fDiffusionL;
357   xyz[0] = gRandom->Gaus(xyz[0], sigmaL * fLorentzFactor);
358   xyz[1] = gRandom->Gaus(xyz[1], sigmaT * fLorentzFactor);
359   xyz[2] = gRandom->Gaus(xyz[2], sigmaT);
360
361   return 1;
362
363 }
364
365 //_____________________________________________________________________________
366 Int_t AliTRDdigitizer::ExB(Float_t driftlength, Float_t *xyz)
367 {
368   //
369   // Applies E x B effects to the position of a single electron
370   //
371
372   xyz[0] = xyz[0];
373   xyz[1] = xyz[1] + fOmegaTau * driftlength;
374   xyz[2] = xyz[2];
375
376   return 1;
377
378 }
379
380 //_____________________________________________________________________________
381 Int_t AliTRDdigitizer::PadResponse(Float_t signal, Float_t dist, Float_t *pad)
382 {
383   //
384   // Applies the pad response
385   //
386
387   Int_t iBin =  ((Int_t) (( - dist - fPRFlo) / fPRFwid));
388
389   Int_t iBin0 = iBin - fPRFpad;
390   Int_t iBin1 = iBin;
391   Int_t iBin2 = iBin + fPRFpad;
392
393   if ((iBin0 >= 0) && (iBin2 < fPRFbin)) {
394
395     pad[0] = signal * fPRFsmp[iBin0];
396     pad[1] = signal * fPRFsmp[iBin1];
397     pad[2] = signal * fPRFsmp[iBin2];
398
399     return 1;
400
401   }
402   else {
403
404     return 0;
405
406   }
407
408 }
409
410 //_____________________________________________________________________________
411 Float_t AliTRDdigitizer::TimeResponse(Float_t time)
412 {
413   //
414   // Applies the preamp shaper time response
415   //
416
417   Int_t iBin = ((Int_t) ((time - fTRFlo) / fTRFwid)); 
418   if ((iBin >= 0) && (iBin < fTRFbin)) {
419     return fTRFsmp[iBin];
420   }
421   else {
422     return 0.0;
423   }    
424
425 }
426
427 //_____________________________________________________________________________
428 void AliTRDdigitizer::Init()
429 {
430   //
431   // Initializes the digitization procedure with standard values
432   //
433
434   // The default parameter for the digitization
435   fGasGain        = 2800.;
436   fChipGain       = 6.1;
437   fNoise          = 1000.;
438   fADCoutRange    = 1023.;          // 10-bit ADC
439   fADCinRange     = 1000.;          // 1V input range
440   fADCthreshold   = 1;
441
442   // For the summable digits
443   fSinRange       = 1000000.;
444   fSoutRange      = 1000000.;
445
446   // The drift velocity (cm / mus)
447   fDriftVelocity  = 1.5;
448
449   // Diffusion on
450   fDiffusionOn    = 1;
451
452   // E x B effects
453   fExBOn          = 0;
454
455   // Propability for electron attachment
456   fElAttachOn     = 0;
457   fElAttachProp   = 0.0;
458
459   // The pad response function
460   fPRFOn          =  1;
461
462   // The time response function
463   fTRFOn          =  1;
464
465   // The pad coupling factor (same number as for the TPC)
466   fPadCoupling    = 0.5;
467
468   // The time coupling factor (same number as for the TPC)
469   fTimeCoupling   = 0.4;
470
471 }
472
473 //_____________________________________________________________________________
474 void AliTRDdigitizer::ReInit()
475 {
476   //
477   // Reinitializes the digitization procedure after a change in the parameter
478   //
479
480   if (!fGeo) {
481     printf("AliTRDdigitizer::ReInit -- ");
482     printf("No geometry defined. Run InitDetector() first\n");
483     exit(1);
484   }
485
486   // Calculate the time bin width in ns
487   fTimeBinWidth   = fGeo->GetTimeBinSize() / fDriftVelocity * 1000.0;
488
489   // The range and the binwidth for the sampled TRF 
490   fTRFbin = 100;
491   // Start 0.2 mus before the signal
492   fTRFlo  = -0.2 * fDriftVelocity;
493   // End the maximum driftlength after the signal 
494   fTRFhi  = AliTRDgeometry::DrThick() 
495           + fGeo->GetTimeAfter() * fGeo->GetTimeBinSize();
496   fTRFwid = (fTRFhi - fTRFlo) / ((Float_t) fTRFbin);
497
498   // Transverse and longitudinal diffusion coefficients (Xe/CO2)
499   fDiffusionT     = GetDiffusionT(fDriftVelocity,fField);
500   fDiffusionL     = GetDiffusionL(fDriftVelocity,fField);
501
502   // omega * tau.= tan(Lorentz-angle)
503   fOmegaTau       = GetOmegaTau(fDriftVelocity,fField);
504
505   // The Lorentz factor
506   if (fExBOn) {
507     fLorentzFactor = 1.0 / (1.0 + fOmegaTau*fOmegaTau);
508   }
509   else {
510     fLorentzFactor = 1.0;
511   }
512
513 }
514
515 //_____________________________________________________________________________
516 void AliTRDdigitizer::SampleTRF()
517 {
518   //
519   // Samples the time response function
520   // It is defined according to Vasiles simulation of the preamp shaper
521   // output and includes the effect of the ion tail (based on Tariqs 
522   // Garfield simulation) and a shaping time of 125 ns FWHM
523   //
524
525   Int_t   ipos1;
526   Int_t   ipos2;
527   Float_t diff;
528
529   const Int_t kNpasa = 200;
530   Float_t time[kNpasa]   = {  -0.280000, -0.270000, -0.260000, -0.250000
531                             , -0.240000, -0.230000, -0.220000, -0.210000
532                             , -0.200000, -0.190000, -0.180000, -0.170000
533                             , -0.160000, -0.150000, -0.140000, -0.130000
534                             , -0.120000, -0.110000, -0.100000, -0.090000
535                             , -0.080000, -0.070000, -0.060000, -0.050000
536                             , -0.040000, -0.030000, -0.020000, -0.010000
537                             , -0.000000,  0.010000,  0.020000,  0.030000
538                             ,  0.040000,  0.050000,  0.060000,  0.070000
539                             ,  0.080000,  0.090000,  0.100000,  0.110000
540                             ,  0.120000,  0.130000,  0.140000,  0.150000
541                             ,  0.160000,  0.170000,  0.180000,  0.190000
542                             ,  0.200000,  0.210000,  0.220000,  0.230000
543                             ,  0.240000,  0.250000,  0.260000,  0.270000
544                             ,  0.280000,  0.290000,  0.300000,  0.310000
545                             ,  0.320000,  0.330000,  0.340000,  0.350000
546                             ,  0.360000,  0.370000,  0.380000,  0.390000
547                             ,  0.400000,  0.410000,  0.420000,  0.430000
548                             ,  0.440000,  0.450000,  0.460000,  0.470000
549                             ,  0.480000,  0.490000,  0.500000,  0.510000
550                             ,  0.520000,  0.530000,  0.540000,  0.550000
551                             ,  0.560000,  0.570000,  0.580000,  0.590000
552                             ,  0.600000,  0.610000,  0.620000,  0.630000
553                             ,  0.640000,  0.650000,  0.660000,  0.670000
554                             ,  0.680000,  0.690000,  0.700000,  0.710000
555                             ,  0.720000,  0.730000,  0.740000,  0.750000
556                             ,  0.760000,  0.770000,  0.780000,  0.790000
557                             ,  0.800000,  0.810000,  0.820000,  0.830000
558                             ,  0.840000,  0.850000,  0.860000,  0.870000
559                             ,  0.880000,  0.890000,  0.900000,  0.910000
560                             ,  0.920000,  0.930000,  0.940000,  0.950000
561                             ,  0.960000,  0.970000,  0.980000,  0.990000
562                             ,  1.000000,  1.010000,  1.020000,  1.030000
563                             ,  1.040000,  1.050000,  1.060000,  1.070000
564                             ,  1.080000,  1.090000,  1.100000,  1.110000
565                             ,  1.120000,  1.130000,  1.140000,  1.150000
566                             ,  1.160000,  1.170000,  1.180000,  1.190000
567                             ,  1.200000,  1.210000,  1.220000,  1.230000
568                             ,  1.240000,  1.250000,  1.260000,  1.270000
569                             ,  1.280000,  1.290000,  1.300000,  1.310000
570                             ,  1.320000,  1.330000,  1.340000,  1.350000
571                             ,  1.360000,  1.370000,  1.380000,  1.390000
572                             ,  1.400000,  1.410000,  1.420000,  1.430000
573                             ,  1.440000,  1.450000,  1.460000,  1.470000
574                             ,  1.480000,  1.490000,  1.500000,  1.510000
575                             ,  1.520000,  1.530000,  1.540000,  1.550000
576                             ,  1.560000,  1.570000,  1.580000,  1.590000
577                             ,  1.600000,  1.610000,  1.620000,  1.630000
578                             ,  1.640000,  1.650000,  1.660000,  1.670000
579                             ,  1.680000,  1.690000,  1.700000,  1.710000 };
580
581   Float_t signal[kNpasa] = {   0.000000,  0.000000,  0.000000,  0.000000 
582                             ,  0.000000,  0.000000,  0.000000,  0.000000
583                             ,  0.000000,  0.000000,  0.000000,  0.000000
584                             ,  0.000000,  0.000000,  0.000000,  0.000098
585                             ,  0.003071,  0.020056,  0.066053,  0.148346
586                             ,  0.263120,  0.398496,  0.540226,  0.674436
587                             ,  0.790977,  0.883083,  0.947744,  0.985714
588                             ,  0.999248,  0.992105,  0.967669,  0.930827
589                             ,  0.884586,  0.833083,  0.778571,  0.723684
590                             ,  0.669173,  0.617293,  0.567669,  0.521805
591                             ,  0.479699,  0.440977,  0.405639,  0.373985
592                             ,  0.345526,  0.320038,  0.297256,  0.276917
593                             ,  0.258797,  0.242632,  0.228195,  0.215301
594                             ,  0.203759,  0.193383,  0.184023,  0.175564
595                             ,  0.167895,  0.160940,  0.154549,  0.148722
596                             ,  0.143308,  0.138346,  0.133722,  0.129398
597                             ,  0.125376,  0.121617,  0.118045,  0.114699
598                             ,  0.111541,  0.108571,  0.105714,  0.103008
599                             ,  0.100414,  0.097970,  0.095602,  0.093346
600                             ,  0.091165,  0.089060,  0.087068,  0.085150
601                             ,  0.083308,  0.081541,  0.079812,  0.078158
602                             ,  0.076541,  0.075000,  0.073496,  0.072068
603                             ,  0.070677,  0.069286,  0.068008,  0.066729
604                             ,  0.065489,  0.064286,  0.063120,  0.061992
605                             ,  0.060902,  0.059850,  0.058797,  0.057820
606                             ,  0.056842,  0.055902,  0.054962,  0.054060
607                             ,  0.053158,  0.052293,  0.051466,  0.050639
608                             ,  0.049850,  0.049060,  0.048308,  0.047556
609                             ,  0.046842,  0.046128,  0.045451,  0.044774
610                             ,  0.044098,  0.043459,  0.042820,  0.042218
611                             ,  0.041617,  0.041015,  0.040451,  0.039887
612                             ,  0.039323,  0.038797,  0.038271,  0.037744
613                             ,  0.037237,  0.036744,  0.036259,  0.035786
614                             ,  0.035323,  0.034872,  0.034429,  0.033996
615                             ,  0.033575,  0.033162,  0.032756,  0.032361
616                             ,  0.031974,  0.031594,  0.031222,  0.030857
617                             ,  0.030496,  0.030143,  0.029793,  0.029451
618                             ,  0.029109,  0.028774,  0.028444,  0.028113
619                             ,  0.027793,  0.027477,  0.027165,  0.026861
620                             ,  0.026564,  0.026271,  0.025981,  0.025699
621                             ,  0.025421,  0.025147,  0.024880,  0.024613
622                             ,  0.024353,  0.024094,  0.023842,  0.023590
623                             ,  0.023346,  0.023102,  0.022865,  0.022628
624                             ,  0.022398,  0.022173,  0.021951,  0.021733
625                             ,  0.021519,  0.021308,  0.021098,  0.020891
626                             ,  0.020688,  0.020485,  0.020286,  0.020090
627                             ,  0.019895,  0.019707,  0.019519,  0.019335
628                             ,  0.019150,  0.018974,  0.018797,  0.018624
629                             ,  0.018451,  0.018282,  0.018113,  0.017947
630                             ,  0.017782,  0.017617,  0.017455,  0.017297 };
631
632   //for (Int_t ipasa = 0; ipasa < kNpasa; ipasa++) {
633   //  time[ipasa] += 0.13; 
634   //  time[ipasa] *= 0.5;
635   //}
636
637   if (fTRFsmp) delete fTRFsmp;
638   fTRFsmp = new Float_t[fTRFbin];
639
640   Float_t loTRF    = TMath::Max(fTRFlo / fDriftVelocity,time[0]);
641   Float_t hiTRF    = TMath::Min(fTRFhi / fDriftVelocity,time[kNpasa-1]);
642   Float_t binWidth = (hiTRF - loTRF) / ((Float_t) fTRFbin);
643
644   // Take the linear interpolation
645   for (Int_t iBin = 0; iBin < fTRFbin; iBin++) {
646
647     Float_t bin = (((Float_t) iBin) + 0.5) * binWidth + loTRF;
648     ipos1 = ipos2 = 0;
649     diff  = 0;
650     do {
651       diff = bin - time[ipos2++];
652     } while (diff > 0);
653     ipos2--;
654     if (ipos2 > kNpasa) ipos2 = kNpasa - 1;
655     ipos1 = ipos2 - 1;
656
657     fTRFsmp[iBin] = signal[ipos2] 
658                   + diff * (signal[ipos2] - signal[ipos1]) 
659                          / (  time[ipos2] -   time[ipos1]);
660
661   }
662
663 }
664
665 //_____________________________________________________________________________
666 void AliTRDdigitizer::SamplePRF()
667 {
668   //
669   // Samples the pad response function
670   //
671
672   const Int_t kPRFbin = 61;
673   Float_t prf[kPRFbin] = { 0.002340, 0.003380, 0.004900, 0.007080, 0.010220
674                          , 0.014740, 0.021160, 0.030230, 0.042800, 0.059830
675                          , 0.082030, 0.109700, 0.142550, 0.179840, 0.220610
676                          , 0.263980, 0.309180, 0.355610, 0.402790, 0.450350
677                          , 0.497930, 0.545190, 0.591740, 0.637100, 0.680610
678                          , 0.721430, 0.758400, 0.790090, 0.814720, 0.830480
679                          , 0.835930, 0.830480, 0.814710, 0.790070, 0.758390
680                          , 0.721410, 0.680590, 0.637080, 0.591730, 0.545180
681                          , 0.497920, 0.450340, 0.402790, 0.355610, 0.309190
682                          , 0.263990, 0.220630, 0.179850, 0.142570, 0.109720
683                          , 0.082040, 0.059830, 0.042820, 0.030230, 0.021170
684                          , 0.014740, 0.010230, 0.007080, 0.004900, 0.003380
685                          , 0.002340 };
686
687   fPRFlo  = -1.5;
688   fPRFhi  =  1.5;
689   fPRFbin = kPRFbin;
690   fPRFwid = (fPRFhi - fPRFlo) / ((Float_t) fPRFbin);
691   fPRFpad = ((Int_t) (1.0 / fPRFwid));
692
693   if (fPRFsmp) delete fPRFsmp;
694   fPRFsmp = new Float_t[fPRFbin];
695   for (Int_t iBin = 0; iBin < fPRFbin; iBin++) {
696     fPRFsmp[iBin] = prf[iBin];
697   }
698
699 }
700
701 //_____________________________________________________________________________
702 Bool_t AliTRDdigitizer::Open(const Char_t *name, Int_t nEvent)
703 {
704   //
705   // Opens a ROOT-file with TRD-hits and reads in the hit-tree
706   //
707
708   // Connect the AliRoot file containing Geometry, Kine, and Hits
709   fInputFile = (TFile*) gROOT->GetListOfFiles()->FindObject(name);
710   if (!fInputFile) {
711     printf("AliTRDdigitizer::Open -- ");
712     printf("Open the ALIROOT-file %s.\n",name);
713     fInputFile = new TFile(name,"UPDATE");
714   }
715   else {
716     printf("AliTRDdigitizer::Open -- ");
717     printf("%s is already open.\n",name);
718   }
719
720   gAlice = (AliRun*) fInputFile->Get("gAlice");
721   if (gAlice) {
722     printf("AliTRDdigitizer::Open -- ");
723     printf("AliRun object found on file.\n");
724   }
725   else {
726     printf("AliTRDdigitizer::Open -- ");
727     printf("Could not find AliRun object.\n");
728     return kFALSE;
729   }
730
731   fEvent = nEvent;
732
733   // Import the Trees for the event nEvent in the file
734   Int_t nparticles = gAlice->GetEvent(fEvent);
735   if (nparticles <= 0) {
736     printf("AliTRDdigitizer::Open -- ");
737     printf("No entries in the trees for event %d.\n",fEvent);
738     return kFALSE;
739   }
740
741   return InitDetector();
742
743 }
744
745 //_____________________________________________________________________________
746 Bool_t AliTRDdigitizer::InitDetector()
747 {
748   //
749   // Sets the pointer to the TRD detector and the geometry
750   //
751
752   // Get the pointer to the detector class and check for version 1
753   fTRD = (AliTRD*) gAlice->GetDetector("TRD");
754   if (fTRD->IsVersion() != 1) {
755     printf("AliTRDdigitizer::InitDetector -- ");
756     printf("TRD must be version 1 (slow simulator).\n");
757     exit(1);
758   }
759
760   // Get the geometry
761   fGeo = fTRD->GetGeometry();
762   printf("AliTRDdigitizer::InitDetector -- ");
763   printf("Geometry version %d\n",fGeo->IsVersion());
764
765   // The magnetic field strength in Tesla
766   fField = 0.2 * gAlice->Field()->Factor();
767
768   ReInit();
769
770   return kTRUE;
771
772 }
773
774 //_____________________________________________________________________________
775 Bool_t AliTRDdigitizer::SumSDigits()
776 {
777   //
778   // Sums up the summable digits and creates final digits
779   // Not yet implemented
780   //
781
782   return kFALSE;
783
784 }
785
786 //_____________________________________________________________________________
787 Bool_t AliTRDdigitizer::MakeDigits()
788 {
789   //
790   // Creates digits.
791   //
792
793   ///////////////////////////////////////////////////////////////
794   // Parameter 
795   ///////////////////////////////////////////////////////////////
796
797   // Converts number of electrons to fC
798   const Double_t kEl2fC  = 1.602E-19 * 1.0E15; 
799
800   ///////////////////////////////////////////////////////////////
801
802   // Number of pads included in the pad response
803   const Int_t kNpad  = 3;
804
805   // Number of track dictionary arrays
806   const Int_t kNDict = AliTRDdigitsManager::kNDict;
807
808   // Half the width of the amplification region
809   const Float_t kAmWidth = AliTRDgeometry::AmThick() / 2.;
810
811   Int_t   iRow, iCol, iTime, iPad;
812   Int_t   iDict  = 0;
813   Int_t   nBytes = 0;
814
815   Int_t   totalSizeDigits = 0;
816   Int_t   totalSizeDict0  = 0;
817   Int_t   totalSizeDict1  = 0;
818   Int_t   totalSizeDict2  = 0;
819
820   Int_t   timeTRDbeg = 0;
821   Int_t   timeTRDend = 1;
822
823   Float_t pos[3];
824   Float_t rot[3];
825   Float_t xyz[3];
826   Float_t padSignal[kNpad];
827   Float_t signalOld[kNpad];
828
829   AliTRDdataArrayF *signals = 0;
830   AliTRDdataArrayI *digits  = 0;
831   AliTRDdataArrayI *dictionary[kNDict];
832
833   // Create a digits manager
834   fDigits = new AliTRDdigitsManager();
835
836   // Create a container for the amplitudes
837   AliTRDsegmentArray *signalsArray 
838                      = new AliTRDsegmentArray("AliTRDdataArrayF",AliTRDgeometry::Ndet());
839
840   if (fTRFOn) {
841     timeTRDbeg = ((Int_t) (-fTRFlo / fGeo->GetTimeBinSize())) - 1;
842     timeTRDend = ((Int_t) ( fTRFhi / fGeo->GetTimeBinSize())) - 1;
843     printf("AliTRDdigitizer::MakeDigits -- ");
844     printf("Sample the TRF between -%d and %d\n",timeTRDbeg,timeTRDend);
845   }
846
847   Float_t elAttachProp = fElAttachProp / 100.; 
848
849   // Create the sampled PRF
850   SamplePRF();
851
852   // Create the sampled TRF
853   SampleTRF();
854
855   if (!fGeo) {
856     printf("AliTRDdigitizer::MakeDigits -- ");
857     printf("No geometry defined\n");
858     return kFALSE;
859   }
860
861   printf("AliTRDdigitizer::MakeDigits -- ");
862   printf("Start creating digits.\n");
863   if (fVerbose > 0) this->Dump();
864
865   // Get the pointer to the hit tree
866   TTree *HitTree = gAlice->TreeH();
867
868   // Get the number of entries in the hit tree
869   // (Number of primary particles creating a hit somewhere)
870   Int_t nTrack = (Int_t) HitTree->GetEntries();
871   if (fVerbose > 0) {
872     printf("AliTRDdigitizer::MakeDigits -- ");
873     printf("Found %d primary particles\n",nTrack);
874   } 
875
876   Int_t detectorOld = -1;
877   Int_t countHits   =  0; 
878
879   // Loop through all entries in the tree
880   for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
881
882     gAlice->ResetHits();
883     nBytes += HitTree->GetEvent(iTrack);
884
885     // Get the number of hits in the TRD created by this particle
886     Int_t nHit = fTRD->Hits()->GetEntriesFast();
887     if (fVerbose > 0) {
888       printf("AliTRDdigitizer::MakeDigits -- ");
889       printf("Found %d hits for primary particle %d\n",nHit,iTrack);
890     }
891
892     // Loop through the TRD hits  
893     for (Int_t iHit = 0; iHit < nHit; iHit++) {
894
895       countHits++;
896
897       AliTRDhit *hit = (AliTRDhit *) fTRD->Hits()->UncheckedAt(iHit);
898               pos[0]   = hit->X();
899               pos[1]   = hit->Y();
900               pos[2]   = hit->Z();
901       Float_t q        = hit->GetCharge();
902       Int_t   track    = hit->Track();
903       Int_t   detector = hit->GetDetector();
904       Int_t   plane    = fGeo->GetPlane(detector);
905       Int_t   sector   = fGeo->GetSector(detector);
906       Int_t   chamber  = fGeo->GetChamber(detector);
907
908       if (!(CheckDetector(plane,chamber,sector))) continue;
909
910       Int_t   nRowMax     = fGeo->GetRowMax(plane,chamber,sector);
911       Int_t   nColMax     = fGeo->GetColMax(plane);
912       Int_t   nTimeMax    = fGeo->GetTimeMax();
913       Int_t   nTimeBefore = fGeo->GetTimeBefore();
914       Int_t   nTimeAfter  = fGeo->GetTimeAfter();
915       Int_t   nTimeTotal  = fGeo->GetTimeTotal();
916       Float_t row0        = fGeo->GetRow0(plane,chamber,sector);
917       Float_t col0        = fGeo->GetCol0(plane);
918       Float_t time0       = fGeo->GetTime0(plane);
919       Float_t rowPadSize  = fGeo->GetRowPadSize(plane,chamber,sector);
920       Float_t colPadSize  = fGeo->GetColPadSize(plane);
921       Float_t timeBinSize = fGeo->GetTimeBinSize();
922       Float_t divideRow   = 1.0 / rowPadSize;
923       Float_t divideCol   = 1.0 / colPadSize;
924       Float_t divideTime  = 1.0 / timeBinSize;
925
926       if (fVerbose > 1) {
927         printf("Analyze hit no. %d ",iHit);
928         printf("-----------------------------------------------------------\n");
929         hit->Dump();
930         printf("plane = %d, sector = %d, chamber = %d\n"
931               ,plane,sector,chamber);
932         printf("nRowMax = %d, nColMax = %d, nTimeMax = %d\n" 
933               ,nRowMax,nColMax,nTimeMax);
934         printf("nTimeBefore = %d, nTimeAfter = %d, nTimeTotal = %d\n"
935               ,nTimeBefore,nTimeAfter,nTimeTotal);
936         printf("row0 = %f, col0 = %f, time0 = %f\n"
937               ,row0,col0,time0);
938         printf("rowPadSize = %f, colPadSize = %f, timeBinSize = %f\n"
939                ,rowPadSize,colPadSize,timeBinSize); 
940       }
941        
942       // Don't analyze test hits
943       if (hit->FromTest()) continue;
944
945       if (detector != detectorOld) {
946
947         if (fVerbose > 1) {
948           printf("AliTRDdigitizer::MakeDigits -- ");
949           printf("Get new container. New det = %d, Old det = %d\n"
950                 ,detector,detectorOld);
951         }
952         // Compress the old one if enabled
953         if ((fCompress) && (detectorOld > -1)) {
954           if (fVerbose > 1) {
955             printf("AliTRDdigitizer::MakeDigits -- ");
956             printf("Compress the old container ...");
957           }
958           signals->Compress(1,0);
959           for (iDict = 0; iDict < kNDict; iDict++) {
960             dictionary[iDict]->Compress(1,0);
961           }
962           if (fVerbose > 1) printf("done\n");
963         }
964         // Get the new container
965         signals = (AliTRDdataArrayF *) signalsArray->At(detector);
966         if (signals->GetNtime() == 0) {
967           // Allocate a new one if not yet existing
968           if (fVerbose > 1) {
969             printf("AliTRDdigitizer::MakeDigits -- ");
970             printf("Allocate a new container ... ");
971           }
972           signals->Allocate(nRowMax,nColMax,nTimeTotal);
973         }
974         else {
975           // Expand an existing one
976           if (fCompress) {
977             if (fVerbose > 1) {
978               printf("AliTRDdigitizer::MakeDigits -- ");
979               printf("Expand an existing container ... ");
980             }
981             signals->Expand();
982           }
983         }
984         // The same for the dictionary
985         for (iDict = 0; iDict < kNDict; iDict++) {       
986           dictionary[iDict] = fDigits->GetDictionary(detector,iDict);
987           if (dictionary[iDict]->GetNtime() == 0) {
988             dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
989           }
990           else {
991             if (fCompress) dictionary[iDict]->Expand();
992           }
993         }      
994         if (fVerbose > 1) printf("done\n");
995         detectorOld = detector;
996       }
997
998       // Rotate the sectors on top of each other
999       fGeo->Rotate(detector,pos,rot);
1000
1001       // The driftlength. It is negative if the hit is in the 
1002       // amplification region.
1003       Float_t driftlength = time0 - rot[0];
1004
1005       // Take also the drift in the amplification region into account
1006       // The drift length is at the moment still the same, regardless of
1007       // the position relativ to the wire. This non-isochronity needs still
1008       // to be implemented.
1009       Float_t driftlengthL = TMath::Abs(driftlength + kAmWidth);
1010       if (fExBOn) driftlengthL /= TMath::Sqrt(fLorentzFactor);
1011
1012       // Loop over all electrons of this hit
1013       // TR photons produce hits with negative charge
1014       Int_t nEl = ((Int_t) TMath::Abs(q));
1015       for (Int_t iEl = 0; iEl < nEl; iEl++) {
1016
1017         xyz[0] = rot[0];
1018         xyz[1] = rot[1];
1019         xyz[2] = rot[2];
1020
1021         // Electron attachment
1022         if (fElAttachOn) {
1023           if (gRandom->Rndm() < (driftlengthL * elAttachProp)) 
1024             continue;
1025         }
1026
1027         // Apply the diffusion smearing
1028         if (fDiffusionOn) {
1029           if (!(Diffusion(driftlengthL,xyz))) continue;
1030         }
1031
1032         // Apply E x B effects (depends on drift direction)
1033         if (fExBOn) { 
1034           if (!(ExB(driftlength+kAmWidth,xyz))) continue;   
1035         }
1036
1037         // The electron position after diffusion and ExB in pad coordinates 
1038         // The pad row (z-direction)
1039         Int_t  rowE = ((Int_t) ((xyz[2] -  row0) * divideRow));
1040         if ((rowE < 0) || (rowE >= nRowMax)) continue;   
1041
1042         // The pad column (rphi-direction)
1043         Int_t  colE = ((Int_t) ((xyz[1] -  col0) * divideCol));
1044         if ((colE < 0) || (colE >= nColMax)) continue;   
1045   
1046         // The time bin (negative for hits in the amplification region)
1047         // In the amplification region the electrons drift from both sides
1048         // to the middle (anode wire plane)
1049         Float_t timeDist   = time0 - xyz[0];
1050         Float_t timeOffset = 0;
1051         Int_t   timeE      = 0;
1052         if (timeDist > 0) {
1053           // The time bin
1054           timeE      = ((Int_t) (timeDist * divideTime));
1055           // The distance of the position to the middle of the timebin
1056           timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) - timeDist;
1057         }
1058         else {
1059           // Difference between half of the amplification gap width and
1060           // the distance to the anode wire
1061           Float_t anodeDist = kAmWidth - TMath::Abs(timeDist + kAmWidth);
1062           // The time bin
1063           timeE      = -1 * (((Int_t ) (anodeDist * divideTime)) + 1);
1064           // The distance of the position to the middle of the timebin
1065           timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) + anodeDist;
1066         }
1067  
1068         // Apply the gas gain including fluctuations
1069         Float_t ggRndm = 0.0;
1070         do {
1071           ggRndm = gRandom->Rndm();
1072         } while (ggRndm <= 0);
1073         Int_t signal = (Int_t) (-fGasGain * TMath::Log(ggRndm));
1074
1075         // Apply the pad response 
1076         if (fPRFOn) {
1077           // The distance of the electron to the center of the pad 
1078           // in units of pad width
1079           Float_t dist = (xyz[1] - col0 - (colE + 0.5) * colPadSize) 
1080                        * divideCol;
1081           if (!(PadResponse(signal,dist,padSignal))) continue;
1082         }
1083         else {
1084           padSignal[0] = 0.0;
1085           padSignal[1] = signal;
1086           padSignal[2] = 0.0;
1087         }
1088
1089         // Sample the time response inside the drift region
1090         // + additional time bins before and after.
1091         // The sampling is done always in the middle of the time bin
1092         for (Int_t iTimeBin = TMath::Max(timeE-timeTRDbeg,        -nTimeBefore) 
1093                   ;iTimeBin < TMath::Min(timeE+timeTRDend,nTimeMax+nTimeAfter ) 
1094                   ;iTimeBin++) {
1095
1096           // Apply the time response
1097           Float_t timeResponse = 1.0;
1098           if (fTRFOn) {
1099             Float_t time = (iTimeBin - timeE) * timeBinSize + timeOffset;
1100             timeResponse = TimeResponse(time);
1101           }
1102
1103           signalOld[0] = 0.0;
1104           signalOld[1] = 0.0;
1105           signalOld[2] = 0.0;
1106
1107           for (iPad = 0; iPad < kNpad; iPad++) {
1108
1109             Int_t colPos = colE + iPad - 1;
1110             if (colPos <        0) continue;
1111             if (colPos >= nColMax) break;
1112
1113             // Add the signals
1114             // Note: The time bin number is shifted by nTimeBefore to avoid negative
1115             // time bins. This has to be subtracted later.
1116             Int_t iCurrentTimeBin = iTimeBin + nTimeBefore;
1117             signalOld[iPad]  = signals->GetDataUnchecked(rowE,colPos,iCurrentTimeBin);
1118             signalOld[iPad] += padSignal[iPad] * timeResponse;
1119             signals->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,signalOld[iPad]);
1120
1121             // Store the track index in the dictionary
1122             // Note: We store index+1 in order to allow the array to be compressed
1123             if (signalOld[iPad] > 0) {
1124               for (iDict = 0; iDict < kNDict; iDict++) {
1125                 Int_t oldTrack = dictionary[iDict]->GetDataUnchecked(rowE
1126                                                                     ,colPos
1127                                                                     ,iCurrentTimeBin);
1128                 if (oldTrack == track+1) break;
1129                 if (oldTrack ==       0) {
1130                   dictionary[iDict]->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,track+1);
1131                   break;
1132                 }
1133               }
1134             }
1135
1136           }
1137
1138         }
1139
1140       }
1141
1142     }
1143
1144   } // All hits finished
1145
1146   printf("AliTRDdigitizer::MakeDigits -- ");
1147   printf("Finished analyzing %d hits\n",countHits);
1148
1149   // The total conversion factor
1150   Float_t convert = kEl2fC * fPadCoupling * fTimeCoupling * fChipGain;
1151
1152   // Loop through all chambers to finalize the digits
1153   for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1154
1155     Int_t plane       = fGeo->GetPlane(iDet);
1156     Int_t sector      = fGeo->GetSector(iDet);
1157     Int_t chamber     = fGeo->GetChamber(iDet);
1158     Int_t nRowMax     = fGeo->GetRowMax(plane,chamber,sector);
1159     Int_t nColMax     = fGeo->GetColMax(plane);
1160     Int_t nTimeMax    = fGeo->GetTimeMax();
1161     Int_t nTimeTotal  = fGeo->GetTimeTotal();
1162
1163     if (fVerbose > 0) {
1164       printf("AliTRDdigitizer::MakeDigits -- ");
1165       printf("Digitization for chamber %d\n",iDet);
1166     }
1167
1168     // Add a container for the digits of this detector
1169     digits = fDigits->GetDigits(iDet);        
1170     // Allocate memory space for the digits buffer
1171     digits->Allocate(nRowMax,nColMax,nTimeTotal);
1172
1173     // Get the signal container
1174     signals = (AliTRDdataArrayF *) signalsArray->At(iDet);
1175     if (signals->GetNtime() == 0) {
1176       // Create missing containers
1177       signals->Allocate(nRowMax,nColMax,nTimeTotal);      
1178     }
1179     else {
1180       // Expand the container if neccessary
1181       if (fCompress) signals->Expand();
1182     }
1183     // Create the missing dictionary containers
1184     for (iDict = 0; iDict < kNDict; iDict++) {       
1185       dictionary[iDict] = fDigits->GetDictionary(iDet,iDict);
1186       if (dictionary[iDict]->GetNtime() == 0) {
1187         dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1188       }
1189     }
1190
1191     Int_t nDigits = 0;
1192
1193     // Don't create noise in detectors that are switched off
1194     if (CheckDetector(plane,chamber,sector)) {
1195
1196       // Create the digits for this chamber
1197       for (iRow  = 0; iRow  <  nRowMax;   iRow++ ) {
1198         for (iCol  = 0; iCol  <  nColMax;   iCol++ ) {
1199           for (iTime = 0; iTime < nTimeTotal; iTime++) {         
1200
1201             // Create summable digits
1202             if (fSDigits) {
1203
1204               Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1205               Int_t adc  = 0;
1206               if (signalAmp >= fSinRange) {
1207                 adc = ((Int_t) fSoutRange);
1208               }
1209               else {
1210                 adc = ((Int_t) (signalAmp * (fSoutRange / fSinRange)));
1211               }
1212               nDigits++;
1213               digits->SetDataUnchecked(iRow,iCol,iTime,adc);
1214
1215             }
1216             // Create normal digits
1217             else {
1218
1219               Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1220
1221               // Add the noise
1222               signalAmp  = TMath::Max((Double_t) gRandom->Gaus(signalAmp,fNoise),0.0);
1223               // Convert to mV
1224               signalAmp *= convert;
1225               // Convert to ADC counts. Set the overflow-bit fADCoutRange if the 
1226               // signal is larger than fADCinRange
1227               Int_t adc  = 0;
1228               if (signalAmp >= fADCinRange) {
1229                 adc = ((Int_t) fADCoutRange);
1230               }
1231               else {
1232                 adc = ((Int_t) (signalAmp * (fADCoutRange / fADCinRange)));
1233               }
1234
1235               // Store the amplitude of the digit if above threshold
1236               if (adc > fADCthreshold) {
1237                 if (fVerbose > 2) {
1238                   printf("  iRow = %d, iCol = %d, iTime = %d\n"
1239                         ,iRow,iCol,iTime);
1240                   printf("  signal = %f, adc = %d\n",signalAmp,adc);
1241                 }
1242                 nDigits++;
1243                 digits->SetDataUnchecked(iRow,iCol,iTime,adc);
1244               }
1245
1246             }
1247
1248           }
1249         }
1250       }
1251
1252     }
1253
1254     // Compress the arrays
1255     digits->Compress(1,0);
1256     for (iDict = 0; iDict < kNDict; iDict++) {
1257       dictionary[iDict]->Compress(1,0);
1258     }
1259
1260     totalSizeDigits += digits->GetSize();
1261     totalSizeDict0  += dictionary[0]->GetSize();
1262     totalSizeDict1  += dictionary[1]->GetSize();
1263     totalSizeDict2  += dictionary[2]->GetSize();
1264
1265     Float_t nPixel = nRowMax * nColMax * nTimeMax;
1266     printf("AliTRDdigitizer::MakeDigits -- ");
1267     printf("Found %d digits in detector %d (%3.0f).\n"
1268           ,nDigits,iDet
1269           ,100.0 * ((Float_t) nDigits) / nPixel);
1270  
1271     if (fCompress) signals->Compress(1,0);
1272
1273   }
1274
1275   printf("AliTRDdigitizer::MakeDigits -- ");
1276   printf("Total number of analyzed hits = %d\n",countHits);
1277
1278   printf("AliTRDdigitizer::MakeDigits -- ");
1279   printf("Total digits data size = %d, %d, %d, %d\n",totalSizeDigits
1280                                                     ,totalSizeDict0
1281                                                     ,totalSizeDict1
1282                                                     ,totalSizeDict2);        
1283
1284   return kTRUE;
1285
1286 }
1287
1288 //_____________________________________________________________________________
1289 Bool_t AliTRDdigitizer::Merge(TTree *trees, Int_t *mask, Int_t nin, Int_t event)
1290 {
1291   //
1292   // Merges the summable digits of different events
1293   //
1294
1295   return kTRUE;
1296
1297 }
1298
1299 //_____________________________________________________________________________
1300 Bool_t AliTRDdigitizer::CheckDetector(Int_t plane, Int_t chamber, Int_t sector)
1301 {
1302   //
1303   // Checks whether a detector is enabled
1304   //
1305
1306   if ((fTRD->GetSensChamber() >=       0) &&
1307       (fTRD->GetSensChamber() != chamber)) return kFALSE;
1308   if ((fTRD->GetSensPlane()   >=       0) &&
1309       (fTRD->GetSensPlane()   !=   plane)) return kFALSE;
1310   if ( fTRD->GetSensSector()  >=       0) {
1311     Int_t sens1 = fTRD->GetSensSector();
1312     Int_t sens2 = sens1 + fTRD->GetSensSectorRange();
1313     sens2 -= ((Int_t) (sens2 / AliTRDgeometry::Nsect())) 
1314            * AliTRDgeometry::Nsect();
1315     if (sens1 < sens2) {
1316       if ((sector < sens1) || (sector >= sens2)) return kFALSE;
1317     }
1318     else {
1319       if ((sector < sens1) && (sector >= sens2)) return kFALSE;
1320     }
1321   }
1322
1323   return kTRUE;
1324
1325 }
1326
1327 //_____________________________________________________________________________
1328 Bool_t AliTRDdigitizer::WriteDigits()
1329 {
1330   //
1331   // Writes out the TRD-digits and the dictionaries
1332   //
1333
1334   // Create the branches
1335   if (!(gAlice->TreeD()->GetBranch("TRDdigits"))) { 
1336     return kFALSE;
1337   }
1338
1339   // Store the digits and the dictionary in the tree
1340   fDigits->WriteDigits();
1341
1342   // Write the new tree into the input file (use overwrite option)
1343   Char_t treeName[15];
1344   sprintf(treeName,"TreeD%d",fEvent);
1345   printf("AliTRDdigitizer::WriteDigits -- ");
1346   printf("Write the digits tree %s for event %d.\n"
1347         ,treeName,fEvent);
1348   gAlice->TreeD()->Write(treeName,TObject::kOverwrite);
1349  
1350   return kTRUE;
1351
1352 }
1353
1354 //_____________________________________________________________________________
1355 Float_t AliTRDdigitizer::GetDiffusionL(Float_t vd, Float_t b)
1356 {
1357   //
1358   // Returns the longitudinal diffusion coefficient for a given drift 
1359   // velocity <vd> and a B-field <b> for Xe/CO2 (15%).
1360   // The values are according to a GARFIELD simulation.
1361   //
1362
1363   const Int_t kNb = 5;
1364   Float_t p0[kNb] = {  0.007440,  0.007493,  0.007513,  0.007672,  0.007831 };
1365   Float_t p1[kNb] = {  0.019252,  0.018912,  0.018636,  0.018012,  0.017343 };
1366   Float_t p2[kNb] = { -0.005042, -0.004926, -0.004867, -0.004650, -0.004424 };
1367   Float_t p3[kNb] = {  0.000195,  0.000189,  0.000195,  0.000182,  0.000169 };
1368
1369   Int_t ib = ((Int_t) (10 * (b - 0.15)));
1370   ib       = TMath::Max(  0,ib);
1371   ib       = TMath::Min(kNb,ib);
1372
1373   Float_t diff = p0[ib] 
1374                + p1[ib] * vd
1375                + p2[ib] * vd*vd
1376                + p3[ib] * vd*vd*vd;
1377
1378   return diff;
1379
1380 }
1381
1382 //_____________________________________________________________________________
1383 Float_t AliTRDdigitizer::GetDiffusionT(Float_t vd, Float_t b)
1384 {
1385   //
1386   // Returns the transverse diffusion coefficient for a given drift 
1387   // velocity <vd> and a B-field <b> for Xe/CO2 (15%).
1388   // The values are according to a GARFIELD simulation.
1389   //
1390
1391   const Int_t kNb = 5;
1392   Float_t p0[kNb] = {  0.009550,  0.009599,  0.009674,  0.009757,  0.009850 };
1393   Float_t p1[kNb] = {  0.006667,  0.006539,  0.006359,  0.006153,  0.005925 };
1394   Float_t p2[kNb] = { -0.000853, -0.000798, -0.000721, -0.000635, -0.000541 };
1395   Float_t p3[kNb] = {  0.000131,  0.000122,  0.000111,  0.000098,  0.000085 };
1396
1397   Int_t ib = ((Int_t) (10 * (b - 0.15)));
1398   ib       = TMath::Max(  0,ib);
1399   ib       = TMath::Min(kNb,ib);
1400
1401   Float_t diff = p0[ib] 
1402                + p1[ib] * vd
1403                + p2[ib] * vd*vd
1404                + p3[ib] * vd*vd*vd;
1405
1406   return diff;
1407
1408 }
1409
1410 //_____________________________________________________________________________
1411 Float_t AliTRDdigitizer::GetOmegaTau(Float_t vd, Float_t b)
1412 {
1413   //
1414   // Returns omega*tau (tan(Lorentz-angle)) for a given drift velocity <vd> 
1415   // and a B-field <b> for Xe/CO2 (15%).
1416   // The values are according to a GARFIELD simulation.
1417   //
1418
1419   const Int_t kNb = 5;
1420   Float_t p0[kNb] = {  0.004810,  0.007412,  0.010252,  0.013409,  0.016888 };
1421   Float_t p1[kNb] = {  0.054875,  0.081534,  0.107333,  0.131983,  0.155455 };
1422   Float_t p2[kNb] = { -0.008682, -0.012896, -0.016987, -0.020880, -0.024623 };
1423   Float_t p3[kNb] = {  0.000155,  0.000238,  0.000330,  0.000428,  0.000541 };
1424
1425   Int_t ib = ((Int_t) (10 * (b - 0.15)));
1426   ib       = TMath::Max(  0,ib);
1427   ib       = TMath::Min(kNb,ib);
1428
1429   Float_t alphaL = p0[ib] 
1430                  + p1[ib] * vd
1431                  + p2[ib] * vd*vd
1432                  + p3[ib] * vd*vd*vd;
1433
1434   return TMath::Tan(alphaL);
1435
1436 }
1437