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