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