]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDdigitizer.cxx
Add new constructor
[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.31  2002/02/11 14:27:11  cblume
19 New pad plane design, new TRF+PRF, tail cancelation, cross talk
20
21 Revision 1.30  2001/11/19 08:44:08  cblume
22 Fix bugs reported by Rene
23
24 Revision 1.29  2001/11/14 19:44:25  hristov
25 Numeric const casted (Alpha)
26
27 Revision 1.28  2001/11/14 16:35:58  cblume
28 Inherits now from AliDetector
29
30 Revision 1.27  2001/11/14 10:50:45  cblume
31 Changes in digits IO. Add merging of summable digits
32
33 Revision 1.26  2001/11/06 17:19:41  cblume
34 Add detailed geometry and simple simulator
35
36 Revision 1.25  2001/06/27 09:54:44  cblume
37 Moved fField initialization to InitDetector()
38
39 Revision 1.24  2001/05/21 16:45:47  hristov
40 Last minute changes (C.Blume)
41
42 Revision 1.23  2001/05/07 08:04:48  cblume
43 New TRF and PRF. Speedup of the code. Digits from amplification region included
44
45 Revision 1.22  2001/03/30 14:40:14  cblume
46 Update of the digitization parameter
47
48 Revision 1.21  2001/03/13 09:30:35  cblume
49 Update of digitization. Moved digit branch definition to AliTRD
50
51 Revision 1.20  2001/02/25 20:19:00  hristov
52 Minor correction: loop variable declared only once for HP, Sun
53
54 Revision 1.19  2001/02/14 18:22:26  cblume
55 Change in the geometry of the padplane
56
57 Revision 1.18  2001/01/26 19:56:57  hristov
58 Major upgrade of AliRoot code
59
60 Revision 1.17  2000/12/08 12:53:27  cblume
61 Change in Copy() function for HP-compiler
62
63 Revision 1.16  2000/12/07 12:20:46  cblume
64 Go back to array compression. Use sampled PRF to speed up digitization
65
66 Revision 1.15  2000/11/23 14:34:08  cblume
67 Fixed bug in expansion routine of arrays (initialize buffers properly)
68
69 Revision 1.14  2000/11/20 08:54:44  cblume
70 Switch off compression as default
71
72 Revision 1.13  2000/11/10 14:57:52  cblume
73 Changes in the geometry constants for the DEC compiler
74
75 Revision 1.12  2000/11/01 14:53:20  cblume
76 Merge with TRD-develop
77
78 Revision 1.1.4.9  2000/10/26 17:00:22  cblume
79 Fixed bug in CheckDetector()
80
81 Revision 1.1.4.8  2000/10/23 13:41:35  cblume
82 Added protection against Log(0) in the gas gain calulation
83
84 Revision 1.1.4.7  2000/10/17 02:27:34  cblume
85 Get rid of global constants
86
87 Revision 1.1.4.6  2000/10/16 01:16:53  cblume
88 Changed timebin 0 to be the one closest to the readout
89
90 Revision 1.1.4.5  2000/10/15 23:34:29  cblume
91 Faster version of the digitizer
92
93 Revision 1.1.4.4  2000/10/06 16:49:46  cblume
94 Made Getters const
95
96 Revision 1.1.4.3  2000/10/04 16:34:58  cblume
97 Replace include files by forward declarations
98
99 Revision 1.1.4.2  2000/09/22 14:41:10  cblume
100 Bug fix in PRF. Included time response. New structure
101
102 Revision 1.10  2000/10/05 07:27:53  cblume
103 Changes in the header-files by FCA
104
105 Revision 1.9  2000/10/02 21:28:19  fca
106 Removal of useless dependecies via forward declarations
107
108 Revision 1.8  2000/06/09 11:10:07  cblume
109 Compiler warnings and coding conventions, next round
110
111 Revision 1.7  2000/06/08 18:32:58  cblume
112 Make code compliant to coding conventions
113
114 Revision 1.6  2000/06/07 16:27:32  cblume
115 Try to remove compiler warnings on Sun and HP
116
117 Revision 1.5  2000/05/09 16:38:57  cblume
118 Removed PadResponse(). Merge problem
119
120 Revision 1.4  2000/05/08 15:53:45  cblume
121 Resolved merge conflict
122
123 Revision 1.3  2000/04/28 14:49:27  cblume
124 Only one declaration of iDict in MakeDigits()
125
126 Revision 1.1.4.1  2000/05/08 14:42:04  cblume
127 Introduced AliTRDdigitsManager
128
129 Revision 1.1  2000/02/28 19:00:13  cblume
130 Add new TRD classes
131
132 */
133
134 ///////////////////////////////////////////////////////////////////////////////
135 //                                                                           //
136 //  Creates and handles digits from TRD hits                                 //
137 //  Author: C. Blume (C.Blume@gsi.de)                                        //
138 //                                                                           //
139 //  The following effects are included:                                      //
140 //      - Diffusion                                                          //
141 //      - ExB effects                                                        //
142 //      - Gas gain including fluctuations                                    //
143 //      - Pad-response (simple Gaussian approximation)                       //
144 //      - Time-response                                                      //
145 //      - Electronics noise                                                  //
146 //      - Electronics gain                                                   //
147 //      - Digitization                                                       //
148 //      - ADC threshold                                                      //
149 //  The corresponding parameter can be adjusted via the various              //
150 //  Set-functions. If these parameters are not explicitly set, default       //
151 //  values are used (see Init-function).                                     //
152 //  As an example on how to use this class to produce digits from hits       //
153 //  have a look at the macro hits2digits.C                                   //
154 //  The production of summable digits is demonstrated in hits2sdigits.C      //
155 //  and the subsequent conversion of the s-digits into normal digits is      //
156 //  explained in sdigits2digits.C.                                           //
157 //                                                                           //
158 ///////////////////////////////////////////////////////////////////////////////
159
160 #include <stdlib.h>
161
162 #include <TMath.h>
163 #include <TVector.h>
164 #include <TRandom.h>
165 #include <TROOT.h>
166 #include <TTree.h>
167 #include <TFile.h>
168 #include <TF1.h>
169 #include <TList.h>
170 #include <TTask.h>
171
172 #include "AliRun.h"
173 #include "AliMagF.h"
174 #include "AliRunDigitizer.h"
175
176 #include "AliTRD.h"
177 #include "AliTRDhit.h"
178 #include "AliTRDdigitizer.h"
179 #include "AliTRDdataArrayI.h"
180 #include "AliTRDdataArrayF.h"
181 #include "AliTRDsegmentArray.h"
182 #include "AliTRDdigitsManager.h"
183 #include "AliTRDgeometry.h"
184
185 ClassImp(AliTRDdigitizer)
186
187 //_____________________________________________________________________________
188 AliTRDdigitizer::AliTRDdigitizer()
189 {
190   //
191   // AliTRDdigitizer default constructor
192   //
193
194   fInputFile          = NULL;
195   fDigitsManager      = NULL;
196   fSDigitsManagerList = NULL;
197   fSDigitsManager     = NULL;
198   fTRD                = NULL;
199   fGeo                = NULL;
200   fPRFsmp             = NULL;
201   fTRFsmp             = NULL;
202   fCTsmp              = NULL;
203
204   fMasks              = 0;
205
206   fEvent              = 0;
207   fGasGain            = 0.0;
208   fNoise              = 0.0;
209   fChipGain           = 0.0;
210   fADCoutRange        = 0.0;
211   fADCinRange         = 0.0;
212   fADCthreshold       = 0;
213   fDiffusionOn        = 0;
214   fDiffusionT         = 0.0;
215   fDiffusionL         = 0.0;
216   fElAttachOn         = 0;
217   fElAttachProp       = 0.0;
218   fExBOn              = 0;
219   fOmegaTau           = 0.0;
220   fPRFOn              = 0;
221   fTRFOn              = 0;
222   fCTOn               = 0;
223   fTCOn               = 0;
224   fDriftVelocity      = 0.0;
225   fPadCoupling        = 0.0;
226   fTimeCoupling       = 0.0;
227   fTimeBinWidth       = 0.0;
228   fField              = 0.0;
229   fTiltingAngle       = 0.0;
230
231   fPRFbin             = 0;
232   fPRFlo              = 0.0;
233   fPRFhi              = 0.0;
234   fPRFwid             = 0.0;
235   fPRFpad             = 0;
236   fTRFbin             = 0;
237   fTRFlo              = 0.0;
238   fTRFhi              = 0.0;
239   fTRFwid             = 0.0;
240   fTCnexp             = 0;
241
242   fCompress           = kTRUE;
243   fDebug              = 0;
244   fSDigits            = kFALSE;
245   fSDigitsScale       = 0.0;
246
247 }
248
249 //_____________________________________________________________________________
250 AliTRDdigitizer::AliTRDdigitizer(const Text_t *name, const Text_t *title)
251                 :AliDigitizer(name,title)
252 {
253   //
254   // AliTRDdigitizer constructor
255   //
256
257   fInputFile          = NULL;
258
259   fDigitsManager      = NULL;
260   fSDigitsManager     = NULL;
261   fSDigitsManagerList = NULL;
262
263   fTRD                = NULL;
264   fGeo                = NULL;
265   fPRFsmp             = NULL;
266   fTRFsmp             = NULL;
267   fCTsmp              = NULL;
268
269   fMasks              = 0;
270
271   fEvent              = 0;
272
273   fCompress           = kTRUE;
274   fDebug              = 0;
275   fSDigits            = kFALSE;
276
277   Init();
278
279 }
280
281 //_____________________________________________________________________________
282 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager
283                                 , const Text_t *name, const Text_t *title)
284                 :AliDigitizer(manager,name,title)
285 {
286   //
287   // AliTRDdigitizer constructor
288   //
289
290   fInputFile          = NULL;
291
292   fDigitsManager      = NULL;
293   fSDigitsManager     = NULL;
294   fSDigitsManagerList = NULL;
295
296   fTRD                = NULL;
297   fGeo                = NULL;
298   fPRFsmp             = NULL;
299   fTRFsmp             = NULL;
300   fCTsmp              = NULL;
301
302   fMasks              = 0;
303
304   fEvent              = 0;
305
306   fCompress           = kTRUE;
307   fDebug              = 0;
308   fSDigits            = kFALSE;
309
310   Init();
311
312 }
313
314 //_____________________________________________________________________________
315 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager)
316                 :AliDigitizer(manager,"AliTRDdigitizer","TRD digitizer")
317 {
318   //
319   // AliTRDdigitizer constructor
320   //
321
322   fInputFile          = NULL;
323
324   fDigitsManager      = NULL;
325   fSDigitsManager     = NULL;
326   fSDigitsManagerList = NULL;
327
328   fTRD                = NULL;
329   fGeo                = NULL;
330   fPRFsmp             = NULL;
331   fTRFsmp             = NULL;
332   fCTsmp              = NULL;
333
334   fMasks              = 0;
335
336   fEvent              = 0;
337
338   fCompress           = kTRUE;
339   fDebug              = 0;
340   fSDigits            = kFALSE;
341
342   Init();
343
344 }
345
346 //_____________________________________________________________________________
347 AliTRDdigitizer::AliTRDdigitizer(const AliTRDdigitizer &d)
348 {
349   //
350   // AliTRDdigitizer copy constructor
351   //
352
353   ((AliTRDdigitizer &) d).Copy(*this);
354
355 }
356
357 //_____________________________________________________________________________
358 AliTRDdigitizer::~AliTRDdigitizer()
359 {
360   //
361   // AliTRDdigitizer destructor
362   //
363
364   if (fInputFile) {
365     fInputFile->Close();
366     delete fInputFile;
367     fInputFile = NULL;
368   }
369
370   if (fDigitsManager) {
371     delete fDigitsManager;
372     fDigitsManager = NULL;
373   }
374
375   if (fSDigitsManager) {
376     delete fSDigitsManager;
377     fSDigitsManager = NULL;
378   }
379
380   if (fSDigitsManagerList) {
381     fSDigitsManagerList->Delete();
382     delete fSDigitsManagerList;
383     fSDigitsManagerList = NULL;
384   }
385
386   if (fTRFsmp) {
387     delete [] fTRFsmp;
388     fTRFsmp = NULL;
389   }
390
391   if (fPRFsmp) {
392     delete [] fPRFsmp;
393     fPRFsmp = NULL;
394   }
395
396   if (fCTsmp) {
397     delete [] fCTsmp;
398     fCTsmp  = NULL;
399   }
400
401   if (fMasks) {
402     delete [] fMasks;
403     fMasks = 0;
404   }
405
406 }
407
408 //_____________________________________________________________________________
409 AliTRDdigitizer &AliTRDdigitizer::operator=(const AliTRDdigitizer &d)
410 {
411   //
412   // Assignment operator
413   //
414
415   if (this != &d) ((AliTRDdigitizer &) d).Copy(*this);
416   return *this;
417
418 }
419
420 //_____________________________________________________________________________
421 void AliTRDdigitizer::Copy(TObject &d)
422 {
423   //
424   // Copy function
425   //
426
427   Int_t iBin;
428
429   ((AliTRDdigitizer &) d).fInputFile          = NULL;
430   ((AliTRDdigitizer &) d).fSDigitsManagerList = NULL;
431   ((AliTRDdigitizer &) d).fSDigitsManager     = NULL;
432   ((AliTRDdigitizer &) d).fDigitsManager      = NULL;
433   ((AliTRDdigitizer &) d).fTRD                = NULL;
434   ((AliTRDdigitizer &) d).fGeo                = NULL;
435
436   ((AliTRDdigitizer &) d).fMasks              = 0;
437
438   ((AliTRDdigitizer &) d).fEvent              = 0;
439
440   ((AliTRDdigitizer &) d).fGasGain            = fGasGain;
441   ((AliTRDdigitizer &) d).fNoise              = fNoise;
442   ((AliTRDdigitizer &) d).fChipGain           = fChipGain;
443   ((AliTRDdigitizer &) d).fADCoutRange        = fADCoutRange;
444   ((AliTRDdigitizer &) d).fADCinRange         = fADCinRange;
445   ((AliTRDdigitizer &) d).fADCthreshold       = fADCthreshold;
446   ((AliTRDdigitizer &) d).fDiffusionOn        = fDiffusionOn; 
447   ((AliTRDdigitizer &) d).fDiffusionT         = fDiffusionT;
448   ((AliTRDdigitizer &) d).fDiffusionL         = fDiffusionL;
449   ((AliTRDdigitizer &) d).fElAttachOn         = fElAttachOn;
450   ((AliTRDdigitizer &) d).fElAttachProp       = fElAttachProp;
451   ((AliTRDdigitizer &) d).fExBOn              = fExBOn;
452   ((AliTRDdigitizer &) d).fOmegaTau           = fOmegaTau;
453   ((AliTRDdigitizer &) d).fLorentzFactor      = fLorentzFactor;
454   ((AliTRDdigitizer &) d).fDriftVelocity      = fDriftVelocity;
455   ((AliTRDdigitizer &) d).fPadCoupling        = fPadCoupling;
456   ((AliTRDdigitizer &) d).fTimeCoupling       = fTimeCoupling;
457   ((AliTRDdigitizer &) d).fTimeBinWidth       = fTimeBinWidth;
458   ((AliTRDdigitizer &) d).fField              = fField;
459   ((AliTRDdigitizer &) d).fPRFOn              = fPRFOn;
460   ((AliTRDdigitizer &) d).fTRFOn              = fTRFOn;
461   ((AliTRDdigitizer &) d).fCTOn               = fCTOn;
462   ((AliTRDdigitizer &) d).fTCOn               = fTCOn;
463   ((AliTRDdigitizer &) d).fTiltingAngle       = fTiltingAngle;
464
465   ((AliTRDdigitizer &) d).fCompress           = fCompress;
466   ((AliTRDdigitizer &) d).fDebug              = fDebug  ;
467   ((AliTRDdigitizer &) d).fSDigits            = fSDigits;
468   ((AliTRDdigitizer &) d).fSDigitsScale       = fSDigitsScale;
469
470   ((AliTRDdigitizer &) d).fPRFbin             = fPRFbin;
471   ((AliTRDdigitizer &) d).fPRFlo              = fPRFlo;
472   ((AliTRDdigitizer &) d).fPRFhi              = fPRFhi;
473   ((AliTRDdigitizer &) d).fPRFwid             = fPRFwid;
474   ((AliTRDdigitizer &) d).fPRFpad             = fPRFpad;
475   if (((AliTRDdigitizer &) d).fPRFsmp) delete [] ((AliTRDdigitizer &) d).fPRFsmp;
476   ((AliTRDdigitizer &) d).fPRFsmp = new Float_t[fPRFbin];
477   for (iBin = 0; iBin < fPRFbin; iBin++) {
478     ((AliTRDdigitizer &) d).fPRFsmp[iBin] = fPRFsmp[iBin];
479   }                                                                             
480   ((AliTRDdigitizer &) d).fTRFbin             = fTRFbin;
481   ((AliTRDdigitizer &) d).fTRFlo              = fTRFlo;
482   ((AliTRDdigitizer &) d).fTRFhi              = fTRFhi;
483   ((AliTRDdigitizer &) d).fTRFwid             = fTRFwid;
484   if (((AliTRDdigitizer &) d).fTRFsmp) delete [] ((AliTRDdigitizer &) d).fTRFsmp;
485   ((AliTRDdigitizer &) d).fTRFsmp = new Float_t[fTRFbin];
486   for (iBin = 0; iBin < fTRFbin; iBin++) {
487     ((AliTRDdigitizer &) d).fTRFsmp[iBin] = fTRFsmp[iBin];
488   }                                      
489
490   if (((AliTRDdigitizer &) d).fCTsmp)  delete [] ((AliTRDdigitizer &) d).fCTsmp;
491   ((AliTRDdigitizer &) d).fCTsmp  = new Float_t[fTRFbin];
492   for (iBin = 0; iBin < fTRFbin; iBin++) {
493     ((AliTRDdigitizer &) d).fCTsmp[iBin]  = fCTsmp[iBin];
494   }                                      
495
496   ((AliTRDdigitizer &) d).fTCnexp             = fTCnexp;
497                                       
498 }
499
500 //_____________________________________________________________________________
501 Float_t AliTRDdigitizer::CrossTalk(Float_t time)
502 {
503   //
504   // Applies the pad-pad capacitive cross talk
505   //
506
507   Int_t iBin = ((Int_t) ((time - fTRFlo) / fTRFwid)); 
508   if ((iBin >= 0) && (iBin < fTRFbin)) {
509     return fCTsmp[iBin];
510   }
511   else {
512     return 0.0;
513   }    
514
515 }
516
517 //_____________________________________________________________________________
518 Int_t AliTRDdigitizer::Diffusion(Float_t driftlength, Float_t *xyz)
519 {
520   //
521   // Applies the diffusion smearing to the position of a single electron
522   //
523
524   Float_t driftSqrt = TMath::Sqrt(driftlength);
525   Float_t sigmaT = driftSqrt * fDiffusionT;
526   Float_t sigmaL = driftSqrt * fDiffusionL;
527   xyz[0] = gRandom->Gaus(xyz[0], sigmaL * fLorentzFactor);
528   xyz[1] = gRandom->Gaus(xyz[1], sigmaT * fLorentzFactor);
529   xyz[2] = gRandom->Gaus(xyz[2], sigmaT);
530
531   return 1;
532
533 }
534
535 //_____________________________________________________________________________
536 Int_t AliTRDdigitizer::ExB(Float_t driftlength, Float_t *xyz)
537 {
538   //
539   // Applies E x B effects to the position of a single electron
540   //
541
542   xyz[0] = xyz[0];
543   xyz[1] = xyz[1] + fOmegaTau * driftlength;
544   xyz[2] = xyz[2];
545
546   return 1;
547
548 }
549
550 //_____________________________________________________________________________
551 Int_t AliTRDdigitizer::PadResponse(Float_t signal, Float_t dist
552                                  , Int_t plane, Float_t *pad)
553 {
554   //
555   // Applies the pad response
556   //
557
558   const Int_t kNplan = AliTRDgeometry::kNplan;
559
560   Int_t iBin  = ((Int_t) (( - dist - fPRFlo) / fPRFwid));
561   Int_t iOff  = plane * fPRFbin;
562
563   Int_t iBin0 = iBin - fPRFpad + iOff;
564   Int_t iBin1 = iBin           + iOff;
565   Int_t iBin2 = iBin + fPRFpad + iOff;
566
567   if ((iBin0 >= 0) && (iBin2 < (fPRFbin*kNplan))) {
568
569     pad[0] = signal * fPRFsmp[iBin0];
570     pad[1] = signal * fPRFsmp[iBin1];
571     pad[2] = signal * fPRFsmp[iBin2];
572
573     return 1;
574
575   }
576   else {
577
578     return 0;
579
580   }
581
582 }
583
584 //_____________________________________________________________________________
585 Float_t AliTRDdigitizer::TimeResponse(Float_t time)
586 {
587   //
588   // Applies the preamp shaper time response
589   //
590
591   Int_t iBin = ((Int_t) ((time - fTRFlo) / fTRFwid)); 
592   if ((iBin >= 0) && (iBin < fTRFbin)) {
593     return fTRFsmp[iBin];
594   }
595   else {
596     return 0.0;
597   }    
598
599 }
600
601 //_____________________________________________________________________________
602 Float_t AliTRDdigitizer::Col0Tilted(Float_t col0, Float_t rowOffset
603                                   , Int_t plane)
604 {
605   //
606   // Calculates col0 for tilted pads
607   //
608
609   Float_t diff = fTiltingAngle * rowOffset;
610   return (col0 + TMath::Power(-1.0,plane) * diff);
611
612 }
613
614 //_____________________________________________________________________________
615 void AliTRDdigitizer::Exec(Option_t* option)
616 {
617   //
618   // Executes the merging
619   //
620
621   Int_t iInput;
622
623   AliTRDdigitsManager *sdigitsManager;
624
625   TString optionString = option;
626   if (optionString.Contains("deb")) {
627     fDebug = 1;
628     if (optionString.Contains("2")) {
629       fDebug = 2;
630     }
631     printf("<AliTRDdigitizer::Exec> ");
632     printf("Called with debug option %d\n",fDebug);
633   }
634
635   Int_t nInput = fManager->GetNinputs();
636   fMasks = new Int_t[nInput];
637   for (iInput = 0; iInput < nInput; iInput++) {
638     fMasks[iInput] = fManager->GetMask(iInput);
639   }
640   
641   // Set the event number
642   fEvent = gAlice->GetEvNumber();
643
644   // Initialization
645   InitDetector();
646
647   for (iInput = 0; iInput < nInput; iInput++) {
648
649     if (fDebug > 0) {
650       printf("<AliTRDdigitizer::Exec> ");
651       printf("Add input stream %d\n",iInput);
652     }
653
654     // Read the s-digits via digits manager
655     sdigitsManager = new AliTRDdigitsManager();
656     sdigitsManager->SetDebug(fDebug);
657     sdigitsManager->SetSDigits(kTRUE);
658     sdigitsManager->ReadDigits(fManager->GetInputTreeTRDS(iInput));
659
660     // Add the s-digits to the input list 
661     AddSDigitsManager(sdigitsManager);
662
663   }
664
665   // Convert the s-digits to normal digits
666   if (fDebug > 0) {
667     printf("<AliTRDdigitizer::Exec> ");
668     printf("Do the conversion\n");
669   }
670   SDigits2Digits();
671
672   // Store the digits
673   if (fDebug > 0) {
674     printf("<AliTRDdigitizer::Exec> ");
675     printf("Write the digits\n");
676   }
677   fDigitsManager->MakeBranch(fManager->GetTreeDTRD());
678   fDigitsManager->WriteDigits();
679   if (fDebug > 0) {
680     printf("<AliTRDdigitizer::Exec> ");
681     printf("Done\n");
682   }
683
684 }
685
686 //_____________________________________________________________________________
687 Bool_t AliTRDdigitizer::Init()
688 {
689   //
690   // Initializes the digitization procedure with standard values
691   //
692
693   // The default parameter for the digitization
694   fGasGain        = 2800.;
695   fChipGain       = 6.1;
696   fNoise          = 1000.;
697   fADCoutRange    = 1023.;          // 10-bit ADC
698   fADCinRange     = 1000.;          // 1V input range
699   fADCthreshold   = 1;
700
701   // For the summable digits
702   fSDigitsScale   = 100.;
703
704   // The drift velocity (cm / mus)
705   fDriftVelocity  = 1.5;
706
707   // Diffusion on
708   fDiffusionOn    = 1;
709
710   // E x B effects
711   fExBOn          = 0;
712
713   // Propability for electron attachment
714   fElAttachOn     = 0;
715   fElAttachProp   = 0.0;
716
717   // The pad response function
718   fPRFOn          = 1;
719
720   // The time response function
721   fTRFOn          = 1;
722
723   // The cross talk
724   fCTOn           = 0;
725
726   // The tail cancelation
727   fTCOn           = 1;
728   
729   // The number of exponentials
730   fTCnexp         = 2;
731
732   // The pad coupling factor (same number as for the TPC)
733   fPadCoupling    = 0.5;
734
735   // The time coupling factor (same number as for the TPC)
736   fTimeCoupling   = 0.4;
737
738   // The tilting angle for the readout pads
739   SetTiltingAngle(5.0);
740
741   return kTRUE;
742
743 }
744
745 //_____________________________________________________________________________
746 Bool_t AliTRDdigitizer::ReInit()
747 {
748   //
749   // Reinitializes the digitization procedure after a change in the parameter
750   //
751
752   if (!fGeo) {
753     printf("AliTRDdigitizer::ReInit -- ");
754     printf("No geometry defined. Run InitDetector() first\n");
755     return kFALSE;
756   }
757
758   // Calculate the time bin width in ns
759   fTimeBinWidth   = fGeo->GetTimeBinSize() / fDriftVelocity * 1000.0;
760
761   // The range and the binwidth for the sampled TRF 
762   fTRFbin = 100;
763   // Start 0.2 mus before the signal
764   fTRFlo  = -0.2 * fDriftVelocity;
765   // End the maximum driftlength after the signal 
766   fTRFhi  = AliTRDgeometry::DrThick() 
767           + fGeo->GetTimeAfter() * fGeo->GetTimeBinSize();
768   fTRFwid = (fTRFhi - fTRFlo) / ((Float_t) fTRFbin);
769
770   // Transverse and longitudinal diffusion coefficients (Xe/CO2)
771   fDiffusionT     = GetDiffusionT(fDriftVelocity,fField);
772   fDiffusionL     = GetDiffusionL(fDriftVelocity,fField);
773
774   // omega * tau.= tan(Lorentz-angle)
775   fOmegaTau       = GetOmegaTau(fDriftVelocity,fField);
776
777   // The Lorentz factor
778   if (fExBOn) {
779     fLorentzFactor = 1.0 / (1.0 + fOmegaTau*fOmegaTau);
780   }
781   else {
782     fLorentzFactor = 1.0;
783   }
784
785   return kTRUE;
786
787 }
788
789 //_____________________________________________________________________________
790 void AliTRDdigitizer::SampleTRF()
791 {
792   //
793   // Samples the time response function
794   //
795   // New TRF from Venelin Angelov, simulated with CADENCE
796   // Pad-ground capacitance = 25 pF
797   // Pad-pad cross talk capacitance = 6 pF   
798   //
799
800   Int_t   ipos1;
801   Int_t   ipos2;
802   Float_t diff;
803
804   const Int_t kNpasa     = 252;
805
806   Float_t time[kNpasa]   = { -0.220000, -0.210000, -0.200000, -0.190000 
807                            , -0.180000, -0.170000, -0.160000, -0.150000 
808                            , -0.140000, -0.130000, -0.120000, -0.110000 
809                            , -0.100000, -0.090000, -0.080000, -0.070000 
810                            , -0.060000, -0.050000, -0.040000, -0.030000 
811                            , -0.020000, -0.010000, -0.000000,  0.010000 
812                            ,  0.020000,  0.030000,  0.040000,  0.050000 
813                            ,  0.060000,  0.070000,  0.080000,  0.090000 
814                            ,  0.100000,  0.110000,  0.120000,  0.130000 
815                            ,  0.140000,  0.150000,  0.160000,  0.170000 
816                            ,  0.180000,  0.190000,  0.200000,  0.210000 
817                            ,  0.220000,  0.230000,  0.240000,  0.250000 
818                            ,  0.260000,  0.270000,  0.280000,  0.290000 
819                            ,  0.300000,  0.310000,  0.320000,  0.330000 
820                            ,  0.340000,  0.350000,  0.360000,  0.370000 
821                            ,  0.380000,  0.390000,  0.400000,  0.410000 
822                            ,  0.420000,  0.430000,  0.440000,  0.450000 
823                            ,  0.460000,  0.470000,  0.480000,  0.490000 
824                            ,  0.500000,  0.510000,  0.520000,  0.530000 
825                            ,  0.540000,  0.550000,  0.560000,  0.570000 
826                            ,  0.580000,  0.590000,  0.600000,  0.610000 
827                            ,  0.620000,  0.630000,  0.640000,  0.650000 
828                            ,  0.660000,  0.670000,  0.680000,  0.690000 
829                            ,  0.700000,  0.710000,  0.720000,  0.730000 
830                            ,  0.740000,  0.750000,  0.760000,  0.770000 
831                            ,  0.780000,  0.790000,  0.800000,  0.810000 
832                            ,  0.820000,  0.830000,  0.840000,  0.850000 
833                            ,  0.860000,  0.870000,  0.880000,  0.890000 
834                            ,  0.900000,  0.910000,  0.920000,  0.930000 
835                            ,  0.940000,  0.950000,  0.960000,  0.970000 
836                            ,  0.980000,  0.990000,  1.000000,  1.010000 
837                            ,  1.020000,  1.030000,  1.040000,  1.050000 
838                            ,  1.060000,  1.070000,  1.080000,  1.090000 
839                            ,  1.100000,  1.110000,  1.120000,  1.130000 
840                            ,  1.140000,  1.150000,  1.160000,  1.170000 
841                            ,  1.180000,  1.190000,  1.200000,  1.210000 
842                            ,  1.220000,  1.230000,  1.240000,  1.250000 
843                            ,  1.260000,  1.270000,  1.280000,  1.290000 
844                            ,  1.300000,  1.310000,  1.320000,  1.330000 
845                            ,  1.340000,  1.350000,  1.360000,  1.370000 
846                            ,  1.380000,  1.390000,  1.400000,  1.410000 
847                            ,  1.420000,  1.430000,  1.440000,  1.450000 
848                            ,  1.460000,  1.470000,  1.480000,  1.490000 
849                            ,  1.500000,  1.510000,  1.520000,  1.530000 
850                            ,  1.540000,  1.550000,  1.560000,  1.570000 
851                            ,  1.580000,  1.590000,  1.600000,  1.610000 
852                            ,  1.620000,  1.630000,  1.640000,  1.650000 
853                            ,  1.660000,  1.670000,  1.680000,  1.690000 
854                            ,  1.700000,  1.710000,  1.720000,  1.730000 
855                            ,  1.740000,  1.750000,  1.760000,  1.770000 
856                            ,  1.780000,  1.790000,  1.800000,  1.810000 
857                            ,  1.820000,  1.830000,  1.840000,  1.850000 
858                            ,  1.860000,  1.870000,  1.880000,  1.890000 
859                            ,  1.900000,  1.910000,  1.920000,  1.930000 
860                            ,  1.940000,  1.950000,  1.960000,  1.970000 
861                            ,  1.980000,  1.990000,  2.000000,  2.010000 
862                            ,  2.020000,  2.030000,  2.040000,  2.050000 
863                            ,  2.060000,  2.070000,  2.080000,  2.090000 
864                            ,  2.100000,  2.110000,  2.120000,  2.130000 
865                            ,  2.140000,  2.150000,  2.160000,  2.170000 
866                            ,  2.180000,  2.190000,  2.200000,  2.210000 
867                            ,  2.220000,  2.230000,  2.240000,  2.250000 
868                            ,  2.260000,  2.270000,  2.280000,  2.290000 };
869
870   Float_t signal[kNpasa] = {  0.000000,  0.000000,  0.000000,  0.000000 
871                            ,  0.000000,  0.000000,  0.000000,  0.000396 
872                            ,  0.005096,  0.022877,  0.061891,  0.126614 
873                            ,  0.215798,  0.324406,  0.444507,  0.566817 
874                            ,  0.683465,  0.787089,  0.873159,  0.937146 
875                            ,  0.979049,  0.999434,  1.000000,  0.983579 
876                            ,  0.954134,  0.913364,  0.866365,  0.813703 
877                            ,  0.759910,  0.706116,  0.653454,  0.603624 
878                            ,  0.556625,  0.514156,  0.475085,  0.439977 
879                            ,  0.408834,  0.380578,  0.355549,  0.333352 
880                            ,  0.313647,  0.296093,  0.280351,  0.266195 
881                            ,  0.253397,  0.241789,  0.231257,  0.221574 
882                            ,  0.212627,  0.204417,  0.196772,  0.189581 
883                            ,  0.182956,  0.176784,  0.171008,  0.165515 
884                            ,  0.160419,  0.155606,  0.151076,  0.146716 
885                            ,  0.142639,  0.138845,  0.135221,  0.131767 
886                            ,  0.128482,  0.125368,  0.122424,  0.119592 
887                            ,  0.116931,  0.114326,  0.111891,  0.109513 
888                            ,  0.107248,  0.105096,  0.103058,  0.101019 
889                            ,  0.099151,  0.097282,  0.095527,  0.093715 
890                            ,  0.092129,  0.090544,  0.088958,  0.087429 
891                            ,  0.086014,  0.084598,  0.083239,  0.081880 
892                            ,  0.080634,  0.079388,  0.078143,  0.077010 
893                            ,  0.075878,  0.074745,  0.073669,  0.072593 
894                            ,  0.071574,  0.070612,  0.069649,  0.068686 
895                            ,  0.067780,  0.066874,  0.066025,  0.065176 
896                            ,  0.064326,  0.063533,  0.062684,  0.061948 
897                            ,  0.061212,  0.060419,  0.059740,  0.059003 
898                            ,  0.058324,  0.057644,  0.057022,  0.056342 
899                            ,  0.055663,  0.055096,  0.054473,  0.053851 
900                            ,  0.053284,  0.052718,  0.052152,  0.051585 
901                            ,  0.051019,  0.050566,  0.050000,  0.049490 
902                            ,  0.048981,  0.048528,  0.048018,  0.047508 
903                            ,  0.047055,  0.046602,  0.046149,  0.045696 
904                            ,  0.045300,  0.044904,  0.044451,  0.044054 
905                            ,  0.043658,  0.043205,  0.042865,  0.042469 
906                            ,  0.042072,  0.041733,  0.041336,  0.040997 
907                            ,  0.040657,  0.040260,  0.039921,  0.039581 
908                            ,  0.039241,  0.038958,  0.038618,  0.038335 
909                            ,  0.037995,  0.037656,  0.037373,  0.037089 
910                            ,  0.036806,  0.036467,  0.036183,  0.035900 
911                            ,  0.035617,  0.035334,  0.035108,  0.034824 
912                            ,  0.034541,  0.034315,  0.034032,  0.033805 
913                            ,  0.033522,  0.033296,  0.033069,  0.032786 
914                            ,  0.032559,  0.032333,  0.032106,  0.031880 
915                            ,  0.031653,  0.031427,  0.031200,  0.030974 
916                            ,  0.030804,  0.030578,  0.030351,  0.030125 
917                            ,  0.029955,  0.029785,  0.029558,  0.029332 
918                            ,  0.029162,  0.028992,  0.028766,  0.028596 
919                            ,  0.028426,  0.028199,  0.028086,  0.027860 
920                            ,  0.027746,  0.027633,  0.027463,  0.027293 
921                            ,  0.027180,  0.027067,  0.026954,  0.026954 
922                            ,  0.026840,  0.026727,  0.026727,  0.026614 
923                            ,  0.026614,  0.026614,  0.026557,  0.026501 
924                            ,  0.026501,  0.026501,  0.026501,  0.026501 
925                            ,  0.026501,  0.026501,  0.026501,  0.026387 
926                            ,  0.026387,  0.026387,  0.026387,  0.026387 
927                            ,  0.026387,  0.026387,  0.026387,  0.026387 
928                            ,  0.026387,  0.026387,  0.026387,  0.026387 
929                            ,  0.026387,  0.026274,  0.026274,  0.026274 
930                            ,  0.026274,  0.026274,  0.026274,  0.026274 
931                            ,  0.026274,  0.026274,  0.026274,  0.026274 
932                            ,  0.026274,  0.026274,  0.026274,  0.026161 };
933
934   Float_t xtalk[kNpasa]  = {  0.000000,  0.000000,  0.000000,  0.000000 
935                            ,  0.000000,  0.000000,  0.000000,  0.000113 
936                            ,  0.000793,  0.003058,  0.007305,  0.013194 
937                            ,  0.019706,  0.025821,  0.030634,  0.033465 
938                            ,  0.034145,  0.032729,  0.029615,  0.025198 
939                            ,  0.019989,  0.014496,  0.009003,  0.003964 
940                            , -0.000510, -0.004190, -0.007191, -0.009400 
941                            , -0.010872, -0.011835, -0.012288, -0.012288 
942                            , -0.012005, -0.011495, -0.010872, -0.010136 
943                            , -0.009343, -0.008607, -0.007871, -0.007191 
944                            , -0.006512, -0.005946, -0.005379, -0.004926 
945                            , -0.004473, -0.004077, -0.003737, -0.003398 
946                            , -0.003114, -0.002831, -0.002605, -0.002378 
947                            , -0.002208, -0.002039, -0.001869, -0.001699 
948                            , -0.001585, -0.001472, -0.001359, -0.001246 
949                            , -0.001132, -0.001019, -0.001019, -0.000906 
950                            , -0.000906, -0.000793, -0.000793, -0.000680 
951                            , -0.000680, -0.000680, -0.000566, -0.000566 
952                            , -0.000566, -0.000566, -0.000453, -0.000453 
953                            , -0.000453, -0.000453, -0.000453, -0.000453 
954                            , -0.000340, -0.000340, -0.000340, -0.000340 
955                            , -0.000340, -0.000340, -0.000340, -0.000340 
956                            , -0.000340, -0.000340, -0.000340, -0.000340 
957                            , -0.000340, -0.000227, -0.000227, -0.000227 
958                            , -0.000227, -0.000227, -0.000227, -0.000227 
959                            , -0.000227, -0.000227, -0.000227, -0.000227 
960                            , -0.000227, -0.000227, -0.000227, -0.000227 
961                            , -0.000227, -0.000227, -0.000227, -0.000227 
962                            , -0.000227, -0.000227, -0.000227, -0.000227 
963                            , -0.000227, -0.000227, -0.000227, -0.000227 
964                            , -0.000227, -0.000227, -0.000227, -0.000227 
965                            , -0.000227, -0.000227, -0.000227, -0.000227 
966                            , -0.000227, -0.000227, -0.000227, -0.000113 
967                            , -0.000113, -0.000113, -0.000113, -0.000113 
968                            , -0.000113, -0.000113, -0.000113, -0.000113 
969                            , -0.000113, -0.000113, -0.000113, -0.000113 
970                            , -0.000113, -0.000113, -0.000113, -0.000113 
971                            , -0.000113, -0.000113, -0.000113, -0.000113 
972                            , -0.000113, -0.000113, -0.000113, -0.000113 
973                            , -0.000113, -0.000113, -0.000113, -0.000113 
974                            , -0.000113, -0.000113, -0.000113, -0.000113 
975                            , -0.000113, -0.000113, -0.000113, -0.000113 
976                            , -0.000113, -0.000113, -0.000113, -0.000113 
977                            , -0.000113, -0.000113, -0.000113, -0.000113 
978                            , -0.000113, -0.000113, -0.000113, -0.000113 
979                            , -0.000113, -0.000113, -0.000113, -0.000113 
980                            , -0.000113, -0.000113, -0.000113, -0.000113 
981                            , -0.000113, -0.000113, -0.000113, -0.000113 
982                            , -0.000113, -0.000113, -0.000113, -0.000113 
983                            , -0.000113, -0.000113, -0.000113, -0.000113 
984                            , -0.000113, -0.000113, -0.000113, -0.000113 
985                            , -0.000113, -0.000113, -0.000113, -0.000113 
986                            , -0.000113, -0.000113, -0.000113,  0.000000 
987                            ,  0.000000,  0.000000,  0.000000,  0.000000 
988                            ,  0.000000,  0.000000,  0.000000,  0.000000 
989                            ,  0.000000,  0.000000,  0.000000,  0.000000 
990                            ,  0.000000,  0.000000,  0.000000,  0.000000 
991                            ,  0.000000,  0.000000,  0.000000,  0.000000 
992                            ,  0.000000,  0.000000,  0.000000,  0.000000 
993                            ,  0.000000,  0.000000,  0.000000,  0.000000 
994                            ,  0.000000,  0.000000,  0.000000,  0.000000 
995                            ,  0.000000,  0.000000,  0.000000,  0.000000 
996                            ,  0.000000,  0.000000,  0.000000,  0.000000 };
997
998   // increase CrossTalk to measurements
999   for (Int_t ipasa = 0; ipasa < kNpasa; ipasa++) {
1000     xtalk[ipasa] *= 1.75;
1001   }
1002
1003   if (fTRFsmp) delete [] fTRFsmp;
1004   fTRFsmp = new Float_t[fTRFbin];
1005   if (fCTsmp)  delete [] fCTsmp;
1006   fCTsmp  = new Float_t[fTRFbin];
1007
1008   Float_t loTRF    = TMath::Max(fTRFlo / fDriftVelocity,time[0]);
1009   Float_t hiTRF    = TMath::Min(fTRFhi / fDriftVelocity,time[kNpasa-1]);
1010   Float_t binWidth = (hiTRF - loTRF) / ((Float_t) fTRFbin);
1011
1012   // Take the linear interpolation
1013   for (Int_t iBin = 0; iBin < fTRFbin; iBin++) {
1014
1015     Float_t bin = (((Float_t) iBin) + 0.5) * binWidth + loTRF;
1016     ipos1 = ipos2 = 0;
1017     diff  = 0;
1018     do {
1019       diff = bin - time[ipos2++];
1020     } while (diff > 0);
1021     ipos2--;
1022     if (ipos2 >= kNpasa) ipos2 = kNpasa - 1;
1023     ipos1 = ipos2 - 1;
1024
1025     fTRFsmp[iBin] = signal[ipos2] 
1026                   + diff * (signal[ipos2] - signal[ipos1]) 
1027                          / (  time[ipos2] -   time[ipos1]);
1028
1029     fCTsmp[iBin]  = xtalk[ipos2] 
1030                   + diff * (xtalk[ipos2]  -  xtalk[ipos1]) 
1031                          / (  time[ipos2] -   time[ipos1]);
1032
1033   }
1034
1035 }
1036
1037 //_____________________________________________________________________________
1038 void AliTRDdigitizer::SamplePRF()
1039 {
1040   //
1041   // Samples the pad response function
1042   //
1043
1044   const Int_t kNplan  = AliTRDgeometry::kNplan;
1045   const Int_t kPRFbin = 61;
1046
1047   Float_t prf[kNplan][kPRFbin] = { { 0.018570, 0.022270, 0.026710, 0.032010
1048                                    , 0.038350, 0.045920, 0.054930, 0.065650
1049                                    , 0.078370, 0.093420, 0.111150, 0.131940
1050                                    , 0.156160, 0.184160, 0.216220, 0.252470
1051                                    , 0.292860, 0.337030, 0.384330, 0.433750
1052                                    , 0.484010, 0.533630, 0.581150, 0.625200
1053                                    , 0.664710, 0.698860, 0.727130, 0.749230
1054                                    , 0.765050, 0.774540, 0.777700, 0.774540
1055                                    , 0.765050, 0.749230, 0.727130, 0.698860
1056                                    , 0.664710, 0.625200, 0.581150, 0.533630
1057                                    , 0.484010, 0.433750, 0.384330, 0.337030
1058                                    , 0.292860, 0.252470, 0.216220, 0.184160
1059                                    , 0.156160, 0.131940, 0.111150, 0.093420
1060                                    , 0.078370, 0.065650, 0.054930, 0.045920
1061                                    , 0.038350, 0.032010, 0.026710, 0.022270
1062                                    , 0.018570                               }
1063                                  , { 0.015730, 0.019040, 0.023030, 0.027840
1064                                    , 0.033650, 0.040650, 0.049060, 0.059160
1065                                    , 0.071260, 0.085710, 0.102910, 0.123270
1066                                    , 0.147240, 0.175220, 0.207590, 0.244540
1067                                    , 0.286090, 0.331920, 0.381350, 0.433290
1068                                    , 0.486290, 0.538710, 0.588870, 0.635280
1069                                    , 0.676760, 0.712460, 0.741890, 0.764810
1070                                    , 0.781150, 0.790930, 0.794180, 0.790930
1071                                    , 0.781150, 0.764810, 0.741890, 0.712460
1072                                    , 0.676760, 0.635280, 0.588870, 0.538710
1073                                    , 0.486290, 0.433290, 0.381350, 0.331920
1074                                    , 0.286090, 0.244540, 0.207590, 0.175220
1075                                    , 0.147240, 0.123270, 0.102910, 0.085710
1076                                    , 0.071260, 0.059160, 0.049060, 0.040650
1077                                    , 0.033650, 0.027840, 0.023030, 0.019040
1078                                    , 0.015730                               }
1079                                  , { 0.013330, 0.016270, 0.019850, 0.024210
1080                                    , 0.029510, 0.035960, 0.043790, 0.053280
1081                                    , 0.064740, 0.078580, 0.095190, 0.115070
1082                                    , 0.138700, 0.166570, 0.199120, 0.236660
1083                                    , 0.279260, 0.326660, 0.378140, 0.432540
1084                                    , 0.488260, 0.543440, 0.596200, 0.644900
1085                                    , 0.688240, 0.725380, 0.755840, 0.779470
1086                                    , 0.796260, 0.806280, 0.809610, 0.806280
1087                                    , 0.796260, 0.779470, 0.755840, 0.725380
1088                                    , 0.688240, 0.644900, 0.596200, 0.543440
1089                                    , 0.488260, 0.432540, 0.378140, 0.326660
1090                                    , 0.279260, 0.236660, 0.199120, 0.166570
1091                                    , 0.138700, 0.115070, 0.095190, 0.078580
1092                                    , 0.064740, 0.053280, 0.043790, 0.035960
1093                                    , 0.029510, 0.024210, 0.019850, 0.016270
1094                                    , 0.013330                               }
1095                                  , { 0.011280, 0.013890, 0.017090, 0.021030
1096                                    , 0.025870, 0.031800, 0.039060, 0.047940
1097                                    , 0.058790, 0.071980, 0.087990, 0.107330
1098                                    , 0.130550, 0.158220, 0.190850, 0.228870
1099                                    , 0.272410, 0.321270, 0.374740, 0.431560
1100                                    , 0.489960, 0.547870, 0.603180, 0.654080
1101                                    , 0.699190, 0.737640, 0.769030, 0.793260
1102                                    , 0.810410, 0.820620, 0.824010, 0.820620
1103                                    , 0.810410, 0.793260, 0.769030, 0.737640
1104                                    , 0.699190, 0.654080, 0.603180, 0.547870
1105                                    , 0.489960, 0.431560, 0.374740, 0.321270
1106                                    , 0.272410, 0.228870, 0.190850, 0.158220
1107                                    , 0.130550, 0.107330, 0.087990, 0.071980
1108                                    , 0.058790, 0.047940, 0.039060, 0.031800
1109                                    , 0.025870, 0.021030, 0.017090, 0.013890
1110                                    , 0.011280                               }
1111                                  , { 0.009550, 0.011860, 0.014720, 0.018270
1112                                    , 0.022660, 0.028100, 0.034820, 0.043120
1113                                    , 0.053340, 0.065900, 0.081280, 0.100040
1114                                    , 0.122800, 0.150180, 0.182800, 0.221170
1115                                    , 0.265550, 0.315790, 0.371180, 0.430370
1116                                    , 0.491430, 0.552030, 0.609840, 0.662860
1117                                    , 0.709630, 0.749290, 0.781490, 0.806220
1118                                    , 0.823650, 0.834000, 0.837430, 0.834000
1119                                    , 0.823650, 0.806220, 0.781490, 0.749290
1120                                    , 0.709630, 0.662860, 0.609840, 0.552030
1121                                    , 0.491430, 0.430370, 0.371180, 0.315790
1122                                    , 0.265550, 0.221170, 0.182800, 0.150180
1123                                    , 0.122800, 0.100040, 0.081280, 0.065900
1124                                    , 0.053340, 0.043120, 0.034820, 0.028100
1125                                    , 0.022660, 0.018270, 0.014720, 0.011860
1126                                    , 0.009550                               }
1127                                  , { 0.008080, 0.010120, 0.012670, 0.015860
1128                                    , 0.019840, 0.024820, 0.031030, 0.038760
1129                                    , 0.048370, 0.060300, 0.075040, 0.093200
1130                                    , 0.115430, 0.142450, 0.174980, 0.213610
1131                                    , 0.258720, 0.310250, 0.367480, 0.429010
1132                                    , 0.492690, 0.555950, 0.616210, 0.671280
1133                                    , 0.719600, 0.760350, 0.793250, 0.818380
1134                                    , 0.836020, 0.846460, 0.849920, 0.846460
1135                                    , 0.836020, 0.818380, 0.793250, 0.760350
1136                                    , 0.719600, 0.671280, 0.616210, 0.555950
1137                                    , 0.492690, 0.429010, 0.367480, 0.310250
1138                                    , 0.258720, 0.213610, 0.174980, 0.142450
1139                                    , 0.115430, 0.093200, 0.075040, 0.060300
1140                                    , 0.048370, 0.038760, 0.031030, 0.024820
1141                                    , 0.019840, 0.015860, 0.012670, 0.010120
1142                                    , 0.008080                               } };
1143
1144   fPRFlo  = -1.5;
1145   fPRFhi  =  1.5;
1146   fPRFbin = kPRFbin;
1147   fPRFwid = (fPRFhi - fPRFlo) / ((Float_t) fPRFbin);
1148   fPRFpad = ((Int_t) (1.0 / fPRFwid));
1149
1150   if (fPRFsmp) delete [] fPRFsmp;
1151   fPRFsmp = new Float_t[kNplan*fPRFbin];
1152   for (Int_t iPla = 0; iPla < kNplan; iPla++) {
1153     for (Int_t iBin = 0; iBin < fPRFbin; iBin++) {
1154       fPRFsmp[iPla*kPRFbin+iBin] = prf[iPla][iBin];
1155     }
1156   } 
1157
1158 }
1159
1160 //_____________________________________________________________________________
1161 Bool_t AliTRDdigitizer::Open(const Char_t *file, Int_t nEvent)
1162 {
1163   //
1164   // Opens a ROOT-file with TRD-hits and reads in the hit-tree
1165   //
1166
1167   // Connect the AliRoot file containing Geometry, Kine, and Hits
1168   fInputFile = (TFile*) gROOT->GetListOfFiles()->FindObject(file);
1169   if (!fInputFile) {
1170     if (fDebug > 0) {
1171       printf("<AliTRDdigitizer::Open> ");
1172       printf("Open the AliROOT-file %s.\n",file);
1173     }
1174     fInputFile = new TFile(file,"UPDATE");
1175   }
1176   else {
1177     if (fDebug > 0) {
1178       printf("<AliTRDdigitizer::Open> ");
1179       printf("%s is already open.\n",file);
1180     }
1181   }
1182
1183   gAlice = (AliRun*) fInputFile->Get("gAlice");
1184   if (gAlice) {
1185     if (fDebug > 0) {
1186       printf("<AliTRDdigitizer::Open> ");
1187       printf("AliRun object found on file.\n");
1188     }
1189   }
1190   else {
1191     printf("<AliTRDdigitizer::Open> ");
1192     printf("Could not find AliRun object.\n");
1193     return kFALSE;
1194   }
1195
1196   fEvent = nEvent;
1197
1198   // Import the Trees for the event nEvent in the file
1199   Int_t nparticles = gAlice->GetEvent(fEvent);
1200   if (nparticles <= 0) {
1201     printf("<AliTRDdigitizer::Open> ");
1202     printf("No entries in the trees for event %d.\n",fEvent);
1203     return kFALSE;
1204   }
1205
1206   if (InitDetector()) {
1207     return MakeBranch();
1208   }
1209   else {
1210     return kFALSE;
1211   }
1212
1213 }
1214
1215 //_____________________________________________________________________________
1216 Bool_t AliTRDdigitizer::InitDetector()
1217 {
1218   //
1219   // Sets the pointer to the TRD detector and the geometry
1220   //
1221
1222   // Get the pointer to the detector class and check for version 1
1223   fTRD = (AliTRD*) gAlice->GetDetector("TRD");
1224   if (fTRD->IsVersion() != 1) {
1225     printf("<AliTRDdigitizer::InitDetector> ");
1226     printf("TRD must be version 1 (slow simulator).\n");
1227     exit(1);
1228   }
1229
1230   // Get the geometry
1231   fGeo = fTRD->GetGeometry();
1232   if (fDebug > 0) {
1233     printf("<AliTRDdigitizer::InitDetector> ");
1234     printf("Geometry version %d\n",fGeo->IsVersion());
1235   }
1236
1237   // The magnetic field strength in Tesla
1238   fField = 0.2 * gAlice->Field()->Factor();
1239
1240   // Create a digits manager
1241   fDigitsManager = new AliTRDdigitsManager();
1242   fDigitsManager->SetSDigits(fSDigits);
1243   fDigitsManager->CreateArrays();
1244   fDigitsManager->SetEvent(fEvent);
1245   fDigitsManager->SetDebug(fDebug);
1246
1247   // The list for the input s-digits manager to be merged
1248   fSDigitsManagerList = new TList();
1249
1250   return ReInit();
1251
1252 }
1253
1254 //_____________________________________________________________________________
1255 Bool_t AliTRDdigitizer::MakeBranch(const Char_t *file)
1256 {
1257   // 
1258   // Create the branches for the digits array
1259   //
1260
1261   return fDigitsManager->MakeBranch(file);
1262
1263 }
1264
1265 //_____________________________________________________________________________
1266 Bool_t AliTRDdigitizer::MakeDigits()
1267 {
1268   //
1269   // Creates digits.
1270   //
1271
1272   ///////////////////////////////////////////////////////////////
1273   // Parameter 
1274   ///////////////////////////////////////////////////////////////
1275
1276   // Converts number of electrons to fC
1277   const Double_t kEl2fC  = 1.602E-19 * 1.0E15; 
1278
1279   ///////////////////////////////////////////////////////////////
1280
1281   // Number of pads included in the pad response
1282   const Int_t kNpad  = 3;
1283
1284   // Number of track dictionary arrays
1285   const Int_t kNDict = AliTRDdigitsManager::kNDict;
1286
1287   // Half the width of the amplification region
1288   const Float_t kAmWidth = AliTRDgeometry::AmThick() / 2.;
1289
1290   Int_t   iRow, iCol, iTime, iPad;
1291   Int_t   iDict  = 0;
1292   Int_t   nBytes = 0;
1293
1294   Int_t   totalSizeDigits = 0;
1295   Int_t   totalSizeDict0  = 0;
1296   Int_t   totalSizeDict1  = 0;
1297   Int_t   totalSizeDict2  = 0;
1298
1299   Int_t   timeTRDbeg = 0;
1300   Int_t   timeTRDend = 1;
1301
1302   Float_t pos[3];
1303   Float_t rot[3];
1304   Float_t xyz[3];
1305   Float_t padSignal[kNpad];
1306   Float_t signalOld[kNpad];
1307
1308   AliTRDdataArrayF *signals = 0;
1309   AliTRDdataArrayI *digits  = 0;
1310   AliTRDdataArrayI *dictionary[kNDict];
1311
1312   // Create a container for the amplitudes
1313   AliTRDsegmentArray *signalsArray 
1314                      = new AliTRDsegmentArray("AliTRDdataArrayF"
1315                                              ,AliTRDgeometry::Ndet());
1316
1317   if (fTRFOn) {
1318     timeTRDbeg = ((Int_t) (-fTRFlo / fGeo->GetTimeBinSize())) - 1;
1319     timeTRDend = ((Int_t) ( fTRFhi / fGeo->GetTimeBinSize())) - 1;
1320     if (fDebug > 0) {
1321       printf("<AliTRDdigitizer::MakeDigits> ");
1322       printf("Sample the TRF between -%d and %d\n",timeTRDbeg,timeTRDend);
1323     }
1324   }
1325
1326   Float_t elAttachProp = fElAttachProp / 100.; 
1327
1328   // Create the sampled PRF
1329   SamplePRF();
1330
1331   // Create the sampled TRF
1332   SampleTRF();
1333
1334   if (!fGeo) {
1335     printf("<AliTRDdigitizer::MakeDigits> ");
1336     printf("No geometry defined\n");
1337     return kFALSE;
1338   }
1339
1340   if (fDebug > 0) {
1341     printf("<AliTRDdigitizer::MakeDigits> ");
1342     printf("Start creating digits.\n");
1343   }
1344
1345   // Get the pointer to the hit tree
1346   TTree *HitTree = gAlice->TreeH();
1347
1348   // Get the number of entries in the hit tree
1349   // (Number of primary particles creating a hit somewhere)
1350   Int_t nTrack = (Int_t) HitTree->GetEntries();
1351   if (fDebug > 0) {
1352     printf("<AliTRDdigitizer::MakeDigits> ");
1353     printf("Found %d primary particles\n",nTrack);
1354   } 
1355
1356   Int_t detectorOld = -1;
1357   Int_t countHits   =  0; 
1358
1359   // Loop through all entries in the tree
1360   for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
1361
1362     gAlice->ResetHits();
1363     nBytes += HitTree->GetEvent(iTrack);
1364
1365     // Loop through the TRD hits
1366     Int_t iHit = 0;
1367     AliTRDhit *hit = (AliTRDhit *) fTRD->FirstHit(-1);
1368     while (hit) {
1369  
1370       countHits++;
1371       iHit++;
1372
1373               pos[0]      = hit->X();
1374               pos[1]      = hit->Y();
1375               pos[2]      = hit->Z();
1376       Float_t q           = hit->GetCharge();
1377       Int_t   track       = hit->Track();
1378       Int_t   detector    = hit->GetDetector();
1379       Int_t   plane       = fGeo->GetPlane(detector);
1380       Int_t   sector      = fGeo->GetSector(detector);
1381       Int_t   chamber     = fGeo->GetChamber(detector);
1382       Int_t   nRowMax     = fGeo->GetRowMax(plane,chamber,sector);
1383       Int_t   nColMax     = fGeo->GetColMax(plane);
1384       Int_t   nTimeMax    = fGeo->GetTimeMax();
1385       Int_t   nTimeBefore = fGeo->GetTimeBefore();
1386       Int_t   nTimeAfter  = fGeo->GetTimeAfter();
1387       Int_t   nTimeTotal  = fGeo->GetTimeTotal();
1388       Float_t row0        = fGeo->GetRow0(plane,chamber,sector);
1389       Float_t col0        = fGeo->GetCol0(plane);
1390       Float_t time0       = fGeo->GetTime0(plane);
1391       Float_t rowPadSize  = fGeo->GetRowPadSize(plane,chamber,sector);
1392       Float_t colPadSize  = fGeo->GetColPadSize(plane);
1393       Float_t timeBinSize = fGeo->GetTimeBinSize();
1394       Float_t divideRow   = 1.0 / rowPadSize;
1395       Float_t divideCol   = 1.0 / colPadSize;
1396       Float_t divideTime  = 1.0 / timeBinSize;
1397
1398       if (fDebug > 1) {
1399         printf("Analyze hit no. %d ",iHit);
1400         printf("-----------------------------------------------------------\n");
1401         hit->Dump();
1402         printf("plane = %d, sector = %d, chamber = %d\n"
1403               ,plane,sector,chamber);
1404         printf("nRowMax = %d, nColMax = %d, nTimeMax = %d\n" 
1405               ,nRowMax,nColMax,nTimeMax);
1406         printf("nTimeBefore = %d, nTimeAfter = %d, nTimeTotal = %d\n"
1407               ,nTimeBefore,nTimeAfter,nTimeTotal);
1408         printf("row0 = %f, col0 = %f, time0 = %f\n"
1409               ,row0,col0,time0);
1410         printf("rowPadSize = %f, colPadSize = %f, timeBinSize = %f\n"
1411                ,rowPadSize,colPadSize,timeBinSize); 
1412       }
1413        
1414       // Don't analyze test hits and switched off detectors
1415       if ((CheckDetector(plane,chamber,sector)) &&
1416           (((Int_t) q) != 0)) {
1417
1418         if (detector != detectorOld) {
1419
1420           if (fDebug > 1) {
1421             printf("<AliTRDdigitizer::MakeDigits> ");
1422             printf("Get new container. New det = %d, Old det = %d\n"
1423                   ,detector,detectorOld);
1424           }
1425           // Compress the old one if enabled
1426           if ((fCompress) && (detectorOld > -1)) {
1427             if (fDebug > 1) {
1428               printf("<AliTRDdigitizer::MakeDigits> ");
1429               printf("Compress the old container ...");
1430             }
1431             signals->Compress(1,0);
1432             for (iDict = 0; iDict < kNDict; iDict++) {
1433               dictionary[iDict]->Compress(1,0);
1434             }
1435             if (fDebug > 1) printf("done\n");
1436           }
1437           // Get the new container
1438           signals = (AliTRDdataArrayF *) signalsArray->At(detector);
1439           if (signals->GetNtime() == 0) {
1440             // Allocate a new one if not yet existing
1441             if (fDebug > 1) {
1442               printf("<AliTRDdigitizer::MakeDigits> ");
1443               printf("Allocate a new container ... ");
1444             }
1445             signals->Allocate(nRowMax,nColMax,nTimeTotal);
1446           }
1447           else {
1448             // Expand an existing one
1449             if (fCompress) {
1450               if (fDebug > 1) {
1451                 printf("<AliTRDdigitizer::MakeDigits> ");
1452                 printf("Expand an existing container ... ");
1453               }
1454               signals->Expand();
1455             }
1456           }
1457           // The same for the dictionary
1458           for (iDict = 0; iDict < kNDict; iDict++) {       
1459             dictionary[iDict] = fDigitsManager->GetDictionary(detector,iDict);
1460             if (dictionary[iDict]->GetNtime() == 0) {
1461               dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1462             }
1463             else {
1464               if (fCompress) dictionary[iDict]->Expand();
1465             }
1466           }      
1467           if (fDebug > 1) printf("done\n");
1468           detectorOld = detector;
1469         }
1470
1471         // Rotate the sectors on top of each other
1472         fGeo->Rotate(detector,pos,rot);
1473
1474         // The driftlength. It is negative if the hit is in the 
1475         // amplification region.
1476         Float_t driftlength = time0 - rot[0];
1477
1478         // Take also the drift in the amplification region into account
1479         // The drift length is at the moment still the same, regardless of
1480         // the position relativ to the wire. This non-isochronity needs still
1481         // to be implemented.
1482         Float_t driftlengthL = TMath::Abs(driftlength + kAmWidth);
1483         if (fExBOn) driftlengthL /= TMath::Sqrt(fLorentzFactor);
1484
1485         // Loop over all electrons of this hit
1486         // TR photons produce hits with negative charge
1487         Int_t nEl = ((Int_t) TMath::Abs(q));
1488         for (Int_t iEl = 0; iEl < nEl; iEl++) {
1489
1490           xyz[0] = rot[0];
1491           xyz[1] = rot[1];
1492           xyz[2] = rot[2];
1493
1494           // Electron attachment
1495           if (fElAttachOn) {
1496             if (gRandom->Rndm() < (driftlengthL * elAttachProp)) 
1497               continue;
1498           }
1499
1500           // Apply the diffusion smearing
1501           if (fDiffusionOn) {
1502             if (!(Diffusion(driftlengthL,xyz))) continue;
1503           }
1504
1505           // Apply E x B effects (depends on drift direction)
1506           if (fExBOn) { 
1507             if (!(ExB(driftlength+kAmWidth,xyz))) continue;   
1508           }
1509
1510           // The electron position after diffusion and ExB in pad coordinates 
1511           // The pad row (z-direction)
1512           Float_t rowDist   = xyz[2] - row0;
1513           Int_t   rowE      = ((Int_t) (rowDist * divideRow));
1514           if ((rowE < 0) || (rowE >= nRowMax)) continue;   
1515           Float_t rowOffset = ((((Float_t) rowE) + 0.5) * rowPadSize) - rowDist;
1516
1517           // The pad column (rphi-direction)
1518           Float_t col0tilt  =  Col0Tilted(col0,rowOffset,plane);
1519           Float_t colDist   = xyz[1] - col0tilt;
1520           Int_t   colE      = ((Int_t) (colDist * divideCol));
1521           if ((colE < 0) || (colE >= nColMax)) continue;   
1522           Float_t colOffset = ((((Float_t) colE) + 0.5) * colPadSize) - colDist;    
1523
1524           // The time bin (negative for hits in the amplification region)
1525           // In the amplification region the electrons drift from both sides
1526           // to the middle (anode wire plane)
1527           Float_t timeDist   = time0 - xyz[0];
1528           Float_t timeOffset = 0;
1529           Int_t   timeE      = 0;
1530           if (timeDist > 0) {
1531             // The time bin
1532             timeE      = ((Int_t) (timeDist * divideTime));
1533             // The distance of the position to the middle of the timebin
1534             timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) - timeDist;
1535           }
1536           else {
1537             // Difference between half of the amplification gap width and
1538             // the distance to the anode wire
1539             Float_t anodeDist = kAmWidth - TMath::Abs(timeDist + kAmWidth);
1540             // The time bin
1541             timeE      = -1 * (((Int_t ) (anodeDist * divideTime)) + 1);
1542             // The distance of the position to the middle of the timebin
1543             timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) + anodeDist;
1544           }
1545  
1546           // Apply the gas gain including fluctuations
1547           Float_t ggRndm = 0.0;
1548           do {
1549             ggRndm = gRandom->Rndm();
1550           } while (ggRndm <= 0);
1551           Int_t signal = (Int_t) (-fGasGain * TMath::Log(ggRndm));
1552
1553           // Apply the pad response 
1554           if (fPRFOn) {
1555             // The distance of the electron to the center of the pad 
1556             // in units of pad width
1557             Float_t dist = - colOffset * divideCol;
1558             if (!(PadResponse(signal,dist,plane,padSignal))) continue;
1559           }
1560           else {
1561             padSignal[0] = 0.0;
1562             padSignal[1] = signal;
1563             padSignal[2] = 0.0;
1564           }
1565
1566           // Sample the time response inside the drift region
1567           // + additional time bins before and after.
1568           // The sampling is done always in the middle of the time bin
1569           for (Int_t iTimeBin = TMath::Max(timeE-timeTRDbeg,        -nTimeBefore) 
1570                     ;iTimeBin < TMath::Min(timeE+timeTRDend,nTimeMax+nTimeAfter ) 
1571                     ;iTimeBin++) {
1572
1573             // Apply the time response
1574             Float_t timeResponse = 1.0;
1575             Float_t crossTalk    = 0.0;
1576             Float_t time         = (iTimeBin - timeE) * timeBinSize + timeOffset;
1577             if (fTRFOn) {
1578               timeResponse = TimeResponse(time);
1579             }
1580             if (fCTOn) {
1581               crossTalk    = CrossTalk(time);
1582             }
1583
1584             signalOld[0] = 0.0;
1585             signalOld[1] = 0.0;
1586             signalOld[2] = 0.0;
1587
1588             for (iPad = 0; iPad < kNpad; iPad++) {
1589
1590               Int_t colPos = colE + iPad - 1;
1591               if (colPos <        0) continue;
1592               if (colPos >= nColMax) break;
1593
1594               // Add the signals
1595               // Note: The time bin number is shifted by nTimeBefore to avoid negative
1596               // time bins. This has to be subtracted later.
1597               Int_t iCurrentTimeBin = iTimeBin + nTimeBefore;
1598               signalOld[iPad]  = signals->GetDataUnchecked(rowE,colPos,iCurrentTimeBin);
1599               if( colPos != colE ) {
1600                 signalOld[iPad] += padSignal[iPad] * (timeResponse + crossTalk);
1601               } 
1602               else {
1603                 signalOld[iPad] += padSignal[iPad] * timeResponse;
1604               }
1605               signals->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,signalOld[iPad]);
1606
1607               // Store the track index in the dictionary
1608               // Note: We store index+1 in order to allow the array to be compressed
1609               if (signalOld[iPad] > 0) {
1610                 for (iDict = 0; iDict < kNDict; iDict++) {
1611                   Int_t oldTrack = dictionary[iDict]->GetDataUnchecked(rowE
1612                                                                       ,colPos
1613                                                                       ,iCurrentTimeBin);
1614                   if (oldTrack == track+1) break;
1615                   if (oldTrack ==       0) {
1616                     dictionary[iDict]->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,track+1);
1617                     break;
1618                   }
1619                 }
1620               }
1621
1622             } // Loop: pads
1623
1624           } // Loop: time bins
1625
1626         } // Loop: electrons of a single hit
1627
1628       } // If: detector and test hit
1629
1630       hit = (AliTRDhit *) fTRD->NextHit();   
1631
1632     } // Loop: hits of one primary track
1633
1634   } // Loop: primary tracks
1635
1636   if (fDebug > 0) {
1637     printf("<AliTRDdigitizer::MakeDigits> ");
1638     printf("Finished analyzing %d hits\n",countHits);
1639   }
1640
1641   // The total conversion factor
1642   Float_t convert = kEl2fC * fPadCoupling * fTimeCoupling * fChipGain;
1643
1644   // Loop through all chambers to finalize the digits
1645   for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1646
1647     Int_t plane       = fGeo->GetPlane(iDet);
1648     Int_t sector      = fGeo->GetSector(iDet);
1649     Int_t chamber     = fGeo->GetChamber(iDet);
1650     Int_t nRowMax     = fGeo->GetRowMax(plane,chamber,sector);
1651     Int_t nColMax     = fGeo->GetColMax(plane);
1652     Int_t nTimeMax    = fGeo->GetTimeMax();
1653     Int_t nTimeTotal  = fGeo->GetTimeTotal();
1654
1655     Double_t *inADC  = new Double_t[nTimeTotal];
1656     Double_t *outADC = new Double_t[nTimeTotal];
1657
1658     if (fDebug > 0) {
1659       printf("<AliTRDdigitizer::MakeDigits> ");
1660       printf("Digitization for chamber %d\n",iDet);
1661     }
1662
1663     // Add a container for the digits of this detector
1664     digits = fDigitsManager->GetDigits(iDet);        
1665     // Allocate memory space for the digits buffer
1666     digits->Allocate(nRowMax,nColMax,nTimeTotal);
1667
1668     // Get the signal container
1669     signals = (AliTRDdataArrayF *) signalsArray->At(iDet);
1670     if (signals->GetNtime() == 0) {
1671       // Create missing containers
1672       signals->Allocate(nRowMax,nColMax,nTimeTotal);      
1673     }
1674     else {
1675       // Expand the container if neccessary
1676       if (fCompress) signals->Expand();
1677     }
1678     // Create the missing dictionary containers
1679     for (iDict = 0; iDict < kNDict; iDict++) {       
1680       dictionary[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1681       if (dictionary[iDict]->GetNtime() == 0) {
1682         dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1683       }
1684     }
1685
1686     Int_t nDigits = 0;
1687
1688     // Don't create noise in detectors that are switched off
1689     if (CheckDetector(plane,chamber,sector)) {
1690
1691       // Create the digits for this chamber
1692       for (iRow  = 0; iRow  <  nRowMax;   iRow++ ) {
1693         for (iCol  = 0; iCol  <  nColMax;   iCol++ ) {
1694
1695           // Create summable digits
1696           if (fSDigits) {
1697
1698             for (iTime = 0; iTime < nTimeTotal; iTime++) {         
1699               Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1700               signalAmp *= fSDigitsScale;
1701               signalAmp  = TMath::Min(signalAmp,(Float_t) 1.0e9);
1702               Int_t adc  = (Int_t) signalAmp;
1703               nDigits++;
1704               digits->SetDataUnchecked(iRow,iCol,iTime,adc);
1705             }
1706
1707           }
1708           // Create normal digits
1709           else {
1710
1711             for (iTime = 0; iTime < nTimeTotal; iTime++) {         
1712               Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1713               // Add the noise
1714               signalAmp  = TMath::Max((Double_t) gRandom->Gaus(signalAmp,fNoise),0.0);
1715               // Convert to mV
1716               signalAmp *= convert;
1717               // Convert to ADC counts. Set the overflow-bit fADCoutRange if the 
1718               // signal is larger than fADCinRange
1719               Int_t adc  = 0;
1720               if (signalAmp >= fADCinRange) {
1721                 adc = ((Int_t) fADCoutRange);
1722               }
1723               else {
1724                 adc = ((Int_t) (signalAmp * (fADCoutRange / fADCinRange)));
1725               }
1726               inADC[iTime]  = adc;
1727               outADC[iTime] = adc;
1728             }
1729
1730             // Apply the tail cancelation via the digital filter
1731             if (fTCOn) {
1732               DeConvExp(inADC,outADC,nTimeTotal,fTCnexp);
1733             }
1734
1735             for (iTime = 0; iTime < nTimeTotal; iTime++) {   
1736               // Store the amplitude of the digit if above threshold
1737               if (outADC[iTime] > fADCthreshold) {
1738                 if (fDebug > 2) {
1739                   printf("  iRow = %d, iCol = %d, iTime = %d, adc = %f\n"
1740                         ,iRow,iCol,iTime,outADC[iTime]);
1741                 }
1742                 nDigits++;
1743                 digits->SetDataUnchecked(iRow,iCol,iTime,outADC[iTime]);
1744               }
1745             }
1746
1747           }
1748
1749         }
1750       }
1751
1752     }
1753
1754     // Compress the arrays
1755     digits->Compress(1,0);
1756     for (iDict = 0; iDict < kNDict; iDict++) {
1757       dictionary[iDict]->Compress(1,0);
1758     }
1759
1760     totalSizeDigits += digits->GetSize();
1761     totalSizeDict0  += dictionary[0]->GetSize();
1762     totalSizeDict1  += dictionary[1]->GetSize();
1763     totalSizeDict2  += dictionary[2]->GetSize();
1764
1765     Float_t nPixel = nRowMax * nColMax * nTimeMax;
1766     if (fDebug > 0) {
1767       printf("<AliTRDdigitizer::MakeDigits> ");
1768       printf("Found %d digits in detector %d (%3.0f).\n"
1769             ,nDigits,iDet
1770             ,100.0 * ((Float_t) nDigits) / nPixel);
1771     } 
1772
1773     if (fCompress) signals->Compress(1,0);
1774
1775     delete [] inADC;
1776     delete [] outADC;
1777
1778   }
1779
1780   if (fDebug > 0) {
1781     printf("<AliTRDdigitizer::MakeDigits> ");
1782     printf("Total number of analyzed hits = %d\n",countHits);
1783     printf("<AliTRDdigitizer::MakeDigits> ");
1784     printf("Total digits data size = %d, %d, %d, %d\n",totalSizeDigits
1785                                                       ,totalSizeDict0
1786                                                       ,totalSizeDict1
1787                                                       ,totalSizeDict2);        
1788   }
1789
1790   return kTRUE;
1791
1792 }
1793
1794 //_____________________________________________________________________________
1795 void AliTRDdigitizer::AddSDigitsManager(AliTRDdigitsManager *man)
1796 {
1797   //
1798   // Add a digits manager for s-digits to the input list.
1799   //
1800
1801   fSDigitsManagerList->Add(man);
1802
1803 }
1804
1805 //_____________________________________________________________________________
1806 Bool_t AliTRDdigitizer::ConvertSDigits()
1807 {
1808   //
1809   // Converts s-digits to normal digits
1810   //
1811
1812   // Number of track dictionary arrays
1813   const Int_t    kNDict = AliTRDdigitsManager::kNDict;
1814
1815   // Converts number of electrons to fC
1816   const Double_t kEl2fC = 1.602E-19 * 1.0E15; 
1817
1818   Int_t iDict = 0;
1819   Int_t iRow;
1820   Int_t iCol;
1821   Int_t iTime;
1822
1823   if (fDebug > 0) {
1824     this->Dump();
1825   }
1826
1827   Double_t sDigitsScale = 1.0 / GetSDigitsScale();
1828   Double_t noise        = GetNoise();
1829   Double_t padCoupling  = GetPadCoupling();
1830   Double_t timeCoupling = GetTimeCoupling();
1831   Double_t chipGain     = GetChipGain();
1832   Double_t convert      = kEl2fC * padCoupling * timeCoupling * chipGain;;
1833   Double_t adcInRange   = GetADCinRange();
1834   Double_t adcOutRange  = GetADCoutRange();
1835   Int_t    adcThreshold = GetADCthreshold();
1836
1837   AliTRDdataArrayI *digitsIn;
1838   AliTRDdataArrayI *digitsOut;
1839   AliTRDdataArrayI *dictionaryIn[kNDict];
1840   AliTRDdataArrayI *dictionaryOut[kNDict];
1841
1842   // Loop through the detectors
1843   for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1844
1845     if (fDebug > 0) {
1846       printf("<AliTRDdigitizer::ConvertSDigits> ");
1847       printf("Convert detector %d to digits.\n",iDet);
1848     }
1849
1850     Int_t plane      = fGeo->GetPlane(iDet);
1851     Int_t sector     = fGeo->GetSector(iDet);
1852     Int_t chamber    = fGeo->GetChamber(iDet);
1853     Int_t nRowMax    = fGeo->GetRowMax(plane,chamber,sector);
1854     Int_t nColMax    = fGeo->GetColMax(plane);
1855     Int_t nTimeTotal = fGeo->GetTimeTotal();
1856
1857     Double_t *inADC  = new Double_t[nTimeTotal];
1858     Double_t *outADC = new Double_t[nTimeTotal];
1859
1860     digitsIn  = fSDigitsManager->GetDigits(iDet);
1861     digitsIn->Expand();
1862     digitsOut = fDigitsManager->GetDigits(iDet);
1863     digitsOut->Allocate(nRowMax,nColMax,nTimeTotal);
1864     for (iDict = 0; iDict < kNDict; iDict++) {
1865       dictionaryIn[iDict]  = fSDigitsManager->GetDictionary(iDet,iDict);
1866       dictionaryIn[iDict]->Expand();
1867       dictionaryOut[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1868       dictionaryOut[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1869     }
1870
1871     for (iRow  = 0; iRow  <  nRowMax;   iRow++ ) {
1872       for (iCol  = 0; iCol  <  nColMax;   iCol++ ) {
1873
1874         for (iTime = 0; iTime < nTimeTotal; iTime++) {         
1875           Double_t signal = (Double_t) digitsIn->GetDataUnchecked(iRow,iCol,iTime);
1876           signal *= sDigitsScale;
1877           // Add the noise
1878           signal  = TMath::Max((Double_t) gRandom->Gaus(signal,noise),0.0);
1879           // Convert to mV
1880           signal *= convert;
1881           // Convert to ADC counts. Set the overflow-bit adcOutRange if the 
1882           // signal is larger than adcInRange
1883           Int_t adc  = 0;
1884           if (signal >= adcInRange) {
1885             adc = ((Int_t) adcOutRange);
1886           }
1887           else {
1888             adc = ((Int_t) (signal * (adcOutRange / adcInRange)));
1889           }
1890           inADC[iTime]  = adc;
1891           outADC[iTime] = adc;
1892         }
1893
1894         // Apply the tail cancelation via the digital filter
1895         if (fTCOn) {
1896           DeConvExp(inADC,outADC,nTimeTotal,fTCnexp);
1897         }
1898
1899         for (iTime = 0; iTime < nTimeTotal; iTime++) {   
1900           // Store the amplitude of the digit if above threshold
1901           if (outADC[iTime] > adcThreshold) {
1902             digitsOut->SetDataUnchecked(iRow,iCol,iTime,outADC[iTime]);
1903             // Copy the dictionary
1904             for (iDict = 0; iDict < kNDict; iDict++) { 
1905               Int_t track = dictionaryIn[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1906               dictionaryOut[iDict]->SetDataUnchecked(iRow,iCol,iTime,track);
1907             }
1908           }
1909         }
1910
1911       }
1912     }
1913
1914     if (fCompress) {
1915       digitsIn->Compress(1,0);
1916       digitsOut->Compress(1,0);
1917       for (iDict = 0; iDict < kNDict; iDict++) {
1918         dictionaryIn[iDict]->Compress(1,0);
1919         dictionaryOut[iDict]->Compress(1,0);
1920       }
1921     }
1922
1923     delete [] inADC;
1924     delete [] outADC;
1925
1926   }    
1927
1928   return kTRUE;
1929
1930 }
1931
1932 //_____________________________________________________________________________
1933 Bool_t AliTRDdigitizer::MergeSDigits()
1934 {
1935   //
1936   // Merges the input s-digits:
1937   //   - The amplitude of the different inputs are summed up.
1938   //   - Of the track IDs from the input dictionaries only one is
1939   //     kept for each input. This works for maximal 3 different merged inputs.
1940   //
1941
1942   // Number of track dictionary arrays
1943   const Int_t kNDict = AliTRDdigitsManager::kNDict;
1944
1945   Int_t iDict = 0;
1946   Int_t jDict = 0;
1947
1948   AliTRDdataArrayI *digitsA;
1949   AliTRDdataArrayI *digitsB;
1950   AliTRDdataArrayI *dictionaryA[kNDict];
1951   AliTRDdataArrayI *dictionaryB[kNDict];
1952
1953   // Get the first s-digits
1954   fSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->First();
1955   if (!fSDigitsManager) return kFALSE;
1956
1957   // Loop through the other sets of s-digits
1958   AliTRDdigitsManager *mergeSDigitsManager;
1959   mergeSDigitsManager = (AliTRDdigitsManager *) 
1960                         fSDigitsManagerList->After(fSDigitsManager);
1961
1962   if (fDebug > 0) {
1963     if (mergeSDigitsManager) {
1964       printf("<AliTRDdigitizer::MergeSDigits> ");
1965       printf("Merge %d input files.\n",fSDigitsManagerList->GetSize());
1966     }
1967     else {
1968       printf("<AliTRDdigitizer::MergeSDigits> ");
1969       printf("Only one input file.\n");
1970     }
1971   }
1972
1973   Int_t iMerge = 0;
1974   while (mergeSDigitsManager) {
1975
1976     iMerge++;
1977
1978     // Loop through the detectors
1979     for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1980
1981       Int_t plane      = fGeo->GetPlane(iDet);
1982       Int_t sector     = fGeo->GetSector(iDet);
1983       Int_t chamber    = fGeo->GetChamber(iDet);
1984       Int_t nRowMax    = fGeo->GetRowMax(plane,chamber,sector);
1985       Int_t nColMax    = fGeo->GetColMax(plane);
1986       Int_t nTimeTotal = fGeo->GetTimeTotal();
1987
1988       // Loop through the pixels of one detector and add the signals
1989       digitsA = fSDigitsManager->GetDigits(iDet);
1990       digitsB = mergeSDigitsManager->GetDigits(iDet);
1991       digitsA->Expand();
1992       digitsB->Expand();
1993       for (iDict = 0; iDict < kNDict; iDict++) {
1994         dictionaryA[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1995         dictionaryB[iDict] = mergeSDigitsManager->GetDictionary(iDet,iDict);
1996         dictionaryA[iDict]->Expand();
1997         dictionaryB[iDict]->Expand();
1998       }
1999
2000       if (fDebug > 0) {
2001         printf("<AliTRDdigitizer::MergeSDigits> ");
2002         printf("Merge detector %d of input no.%d\n",iDet,iMerge+1);
2003       }
2004
2005       for (Int_t iRow  = 0; iRow  <  nRowMax;   iRow++ ) {
2006         for (Int_t iCol  = 0; iCol  <  nColMax;   iCol++ ) {
2007           for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {         
2008
2009             // Add the amplitudes of the summable digits 
2010             Int_t ampA = digitsA->GetDataUnchecked(iRow,iCol,iTime);
2011             Int_t ampB = digitsB->GetDataUnchecked(iRow,iCol,iTime);
2012             ampA += ampB;
2013             digitsA->SetDataUnchecked(iRow,iCol,iTime,ampA);
2014
2015             // Add the mask to the track id if defined.
2016             for (iDict = 0; iDict < kNDict; iDict++) {
2017               Int_t trackB = dictionaryB[iDict]->GetDataUnchecked(iRow,iCol,iTime);
2018               if ((fMasks) && (trackB > 0)) {
2019                 for (jDict = 0; jDict < kNDict; jDict++) { 
2020                   Int_t trackA = dictionaryA[iDict]->GetDataUnchecked(iRow,iCol,iTime);
2021                   if (trackA == 0) {
2022                     trackA = trackB + fMasks[iMerge];
2023                     dictionaryA[iDict]->SetDataUnchecked(iRow,iCol,iTime,trackA);
2024                   }
2025                 }
2026               }
2027             }
2028
2029           }
2030         }
2031       }
2032
2033       if (fCompress) {
2034         digitsA->Compress(1,0);
2035         digitsB->Compress(1,0);
2036         for (iDict = 0; iDict < kNDict; iDict++) {
2037           dictionaryA[iDict]->Compress(1,0);
2038           dictionaryB[iDict]->Compress(1,0);
2039         }
2040       }
2041
2042     }    
2043
2044     // The next set of s-digits
2045     mergeSDigitsManager = (AliTRDdigitsManager *) 
2046                           fSDigitsManagerList->After(mergeSDigitsManager);
2047
2048   }
2049
2050   return kTRUE;
2051
2052 }
2053
2054 //_____________________________________________________________________________
2055 Bool_t AliTRDdigitizer::SDigits2Digits()
2056 {
2057   //
2058   // Merges the input s-digits and converts them to normal digits
2059   //
2060
2061   if (!MergeSDigits()) return kFALSE;
2062
2063   return ConvertSDigits();
2064
2065 }
2066
2067 //_____________________________________________________________________________
2068 Bool_t AliTRDdigitizer::CheckDetector(Int_t plane, Int_t chamber, Int_t sector)
2069 {
2070   //
2071   // Checks whether a detector is enabled
2072   //
2073
2074   if ((fTRD->GetSensChamber() >=       0) &&
2075       (fTRD->GetSensChamber() != chamber)) return kFALSE;
2076   if ((fTRD->GetSensPlane()   >=       0) &&
2077       (fTRD->GetSensPlane()   !=   plane)) return kFALSE;
2078   if ( fTRD->GetSensSector()  >=       0) {
2079     Int_t sens1 = fTRD->GetSensSector();
2080     Int_t sens2 = sens1 + fTRD->GetSensSectorRange();
2081     sens2 -= ((Int_t) (sens2 / AliTRDgeometry::Nsect())) 
2082            * AliTRDgeometry::Nsect();
2083     if (sens1 < sens2) {
2084       if ((sector < sens1) || (sector >= sens2)) return kFALSE;
2085     }
2086     else {
2087       if ((sector < sens1) && (sector >= sens2)) return kFALSE;
2088     }
2089   }
2090
2091   return kTRUE;
2092
2093 }
2094
2095 //_____________________________________________________________________________
2096 Bool_t AliTRDdigitizer::WriteDigits()
2097 {
2098   //
2099   // Writes out the TRD-digits and the dictionaries
2100   //
2101
2102   // Store the digits and the dictionary in the tree
2103   return fDigitsManager->WriteDigits();
2104
2105 }
2106
2107 //_____________________________________________________________________________
2108 void AliTRDdigitizer::SetTiltingAngle(Float_t v)
2109 {
2110   //
2111   // Set the tilting angle for the readout pads
2112   //
2113
2114   fTiltingAngle = TMath::Tan(TMath::Pi()/180.0 * v);
2115
2116 }
2117
2118 //_____________________________________________________________________________
2119 Float_t AliTRDdigitizer::GetTiltingAngle() const
2120 {
2121   //
2122   // Get the tilting angle for the readout pads
2123   //
2124
2125   return TMath::ATan(180.0/TMath::Pi() * fTiltingAngle);
2126
2127 }
2128
2129 //_____________________________________________________________________________
2130 Float_t AliTRDdigitizer::GetDiffusionL(Float_t vd, Float_t b)
2131 {
2132   //
2133   // Returns the longitudinal diffusion coefficient for a given drift 
2134   // velocity <vd> and a B-field <b> for Xe/CO2 (15%).
2135   // The values are according to a GARFIELD simulation.
2136   //
2137
2138   const Int_t kNb = 5;
2139   Float_t p0[kNb] = {  0.007440,  0.007493,  0.007513,  0.007672,  0.007831 };
2140   Float_t p1[kNb] = {  0.019252,  0.018912,  0.018636,  0.018012,  0.017343 };
2141   Float_t p2[kNb] = { -0.005042, -0.004926, -0.004867, -0.004650, -0.004424 };
2142   Float_t p3[kNb] = {  0.000195,  0.000189,  0.000195,  0.000182,  0.000169 };
2143
2144   Int_t ib = ((Int_t) (10 * (b - 0.15)));
2145   ib       = TMath::Max(  0,ib);
2146   ib       = TMath::Min(kNb,ib);
2147
2148   Float_t diff = p0[ib] 
2149                + p1[ib] * vd
2150                + p2[ib] * vd*vd
2151                + p3[ib] * vd*vd*vd;
2152
2153   return diff;
2154
2155 }
2156
2157 //_____________________________________________________________________________
2158 Float_t AliTRDdigitizer::GetDiffusionT(Float_t vd, Float_t b)
2159 {
2160   //
2161   // Returns the transverse diffusion coefficient for a given drift 
2162   // velocity <vd> and a B-field <b> for Xe/CO2 (15%).
2163   // The values are according to a GARFIELD simulation.
2164   //
2165
2166   const Int_t kNb = 5;
2167   Float_t p0[kNb] = {  0.009550,  0.009599,  0.009674,  0.009757,  0.009850 };
2168   Float_t p1[kNb] = {  0.006667,  0.006539,  0.006359,  0.006153,  0.005925 };
2169   Float_t p2[kNb] = { -0.000853, -0.000798, -0.000721, -0.000635, -0.000541 };
2170   Float_t p3[kNb] = {  0.000131,  0.000122,  0.000111,  0.000098,  0.000085 };
2171
2172   Int_t ib = ((Int_t) (10 * (b - 0.15)));
2173   ib       = TMath::Max(  0,ib);
2174   ib       = TMath::Min(kNb,ib);
2175
2176   Float_t diff = p0[ib] 
2177                + p1[ib] * vd
2178                + p2[ib] * vd*vd
2179                + p3[ib] * vd*vd*vd;
2180
2181   return diff;
2182
2183 }
2184
2185 //_____________________________________________________________________________
2186 Float_t AliTRDdigitizer::GetOmegaTau(Float_t vd, Float_t b)
2187 {
2188   //
2189   // Returns omega*tau (tan(Lorentz-angle)) for a given drift velocity <vd> 
2190   // and a B-field <b> for Xe/CO2 (15%).
2191   // The values are according to a GARFIELD simulation.
2192   //
2193
2194   const Int_t kNb = 5;
2195   Float_t p0[kNb] = {  0.004810,  0.007412,  0.010252,  0.013409,  0.016888 };
2196   Float_t p1[kNb] = {  0.054875,  0.081534,  0.107333,  0.131983,  0.155455 };
2197   Float_t p2[kNb] = { -0.008682, -0.012896, -0.016987, -0.020880, -0.024623 };
2198   Float_t p3[kNb] = {  0.000155,  0.000238,  0.000330,  0.000428,  0.000541 };
2199
2200   Int_t ib = ((Int_t) (10 * (b - 0.15)));
2201   ib       = TMath::Max(  0,ib);
2202   ib       = TMath::Min(kNb,ib);
2203
2204   Float_t alphaL = p0[ib] 
2205                  + p1[ib] * vd
2206                  + p2[ib] * vd*vd
2207                  + p3[ib] * vd*vd*vd;
2208
2209   return TMath::Tan(alphaL);
2210
2211 }
2212
2213 //_____________________________________________________________________________
2214 void AliTRDdigitizer::DeConvExp(Double_t *source, Double_t *target
2215                               , Int_t n, Int_t nexp) 
2216 {
2217   //
2218   // Does the deconvolution by the digital filter.
2219   //
2220   // Author:        Marcus Gutfleisch, KIP Heidelberg
2221   // Optimized for: New TRF from Venelin Angelov, simulated with CADENCE
2222   //                Pad-ground capacitance = 25 pF
2223   //                Pad-pad cross talk capacitance = 6 pF
2224   //                For 10 MHz digitization, corresponding to 20 time bins
2225   //                in the drift region
2226   //
2227
2228   Double_t rates[2];
2229   Double_t coefficients[2];
2230
2231   /* initialize (coefficient = alpha, rates = lambda) */
2232   
2233   if( nexp == 1 ) {
2234     rates[0] = 0.466998;
2235     /* no rescaling */
2236     coefficients[0] = 1.0;
2237   }
2238   if( nexp == 2 ) {
2239     rates[0] = 0.8988162;
2240     coefficients[0] = 0.11392069;
2241     rates[1] = 0.3745688;
2242     coefficients[1] = 0.8860793;
2243     /* no rescaling */
2244     Float_t sumc = coefficients[0]+coefficients[1];
2245     coefficients[0] /= sumc;
2246     coefficients[1] /= sumc;
2247   }
2248       
2249   Int_t i, k;
2250   Double_t reminder[2];
2251   Double_t correction, result;
2252
2253   /* attention: computation order is important */
2254   correction=0.0;
2255   for ( k = 0; k < nexp; k++ ) reminder[k]=0.0;
2256     
2257   for ( i = 0; i < n; i++ ) {
2258     result = ( source[i] - correction );    /* no rescaling */
2259     target[i] = result;
2260     
2261     for ( k = 0; k < nexp; k++ ) reminder[k] = rates[k] 
2262                              * ( reminder[k] + coefficients[k] * result);
2263       
2264     correction=0.0;
2265     for ( k = 0; k < nexp; k++ ) correction += reminder[k];
2266   }
2267   
2268 }
2269