]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDdigitizer.cxx
24490c6cd59badf89dadc9ec01f70f344168e70a
[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.41  2002/10/21 09:10:32  cblume
19 Fix type conversion warnings
20
21 Revision 1.40  2002/10/14 14:57:43  hristov
22 Merging the VirtualMC branch to the main development branch (HEAD)
23
24 Revision 1.33.6.3  2002/10/11 07:26:37  hristov
25 Updating VirtualMC to v3-09-02
26
27 Revision 1.39  2002/10/08 20:46:12  cblume
28 Do coupling factors before noise is applied
29
30 Revision 1.38  2002/04/30 08:30:40  cblume
31 gAlice now only read by AliRunDigitizer. Therefore it is just deleted in AliTRDmerge.C
32
33 Revision 1.37  2002/04/29 11:50:47  cblume
34 Change initialization of gAlice in the merging case
35
36 Revision 1.36  2002/04/12 12:13:23  cblume
37 Add Jiris changes
38
39 Revision 1.35  2002/03/28 14:59:07  cblume
40 Coding conventions
41
42 Revision 1.34  2002/03/25 20:00:44  cblume
43 Introduce parameter class and regions of interest for merging
44
45 Revision 1.33  2002/02/12 17:32:03  cblume
46 Rearrange the deleting of the list of sdigitsmanager
47
48 Revision 1.32  2002/02/12 16:07:21  cblume
49 Add new constructor
50
51 Revision 1.31  2002/02/11 14:27:11  cblume
52 New pad plane design, new TRF+PRF, tail cancelation, cross talk
53
54 Revision 1.30  2001/11/19 08:44:08  cblume
55 Fix bugs reported by Rene
56
57 Revision 1.29  2001/11/14 19:44:25  hristov
58 Numeric const casted (Alpha)
59
60 Revision 1.28  2001/11/14 16:35:58  cblume
61 Inherits now from AliDetector
62
63 Revision 1.27  2001/11/14 10:50:45  cblume
64 Changes in digits IO. Add merging of summable digits
65
66 Revision 1.26  2001/11/06 17:19:41  cblume
67 Add detailed geometry and simple simulator
68
69 Revision 1.25  2001/06/27 09:54:44  cblume
70 Moved fField initialization to InitDetector()
71
72 Revision 1.24  2001/05/21 16:45:47  hristov
73 Last minute changes (C.Blume)
74
75 Revision 1.23  2001/05/07 08:04:48  cblume
76 New TRF and PRF. Speedup of the code. Digits from amplification region included
77
78 Revision 1.22  2001/03/30 14:40:14  cblume
79 Update of the digitization parameter
80
81 Revision 1.21  2001/03/13 09:30:35  cblume
82 Update of digitization. Moved digit branch definition to AliTRD
83
84 Revision 1.20  2001/02/25 20:19:00  hristov
85 Minor correction: loop variable declared only once for HP, Sun
86
87 Revision 1.19  2001/02/14 18:22:26  cblume
88 Change in the geometry of the padplane
89
90 Revision 1.18  2001/01/26 19:56:57  hristov
91 Major upgrade of AliRoot code
92
93 Revision 1.17  2000/12/08 12:53:27  cblume
94 Change in Copy() function for HP-compiler
95
96 Revision 1.16  2000/12/07 12:20:46  cblume
97 Go back to array compression. Use sampled PRF to speed up digitization
98
99 Revision 1.15  2000/11/23 14:34:08  cblume
100 Fixed bug in expansion routine of arrays (initialize buffers properly)
101
102 Revision 1.14  2000/11/20 08:54:44  cblume
103 Switch off compression as default
104
105 Revision 1.13  2000/11/10 14:57:52  cblume
106 Changes in the geometry constants for the DEC compiler
107
108 Revision 1.12  2000/11/01 14:53:20  cblume
109 Merge with TRD-develop
110
111 Revision 1.1.4.9  2000/10/26 17:00:22  cblume
112 Fixed bug in CheckDetector()
113
114 Revision 1.1.4.8  2000/10/23 13:41:35  cblume
115 Added protection against Log(0) in the gas gain calulation
116
117 Revision 1.1.4.7  2000/10/17 02:27:34  cblume
118 Get rid of global constants
119
120 Revision 1.1.4.6  2000/10/16 01:16:53  cblume
121 Changed timebin 0 to be the one closest to the readout
122
123 Revision 1.1.4.5  2000/10/15 23:34:29  cblume
124 Faster version of the digitizer
125
126 Revision 1.1.4.4  2000/10/06 16:49:46  cblume
127 Made Getters const
128
129 Revision 1.1.4.3  2000/10/04 16:34:58  cblume
130 Replace include files by forward declarations
131
132 Revision 1.1.4.2  2000/09/22 14:41:10  cblume
133 Bug fix in PRF. Included time response. New structure
134
135 Revision 1.10  2000/10/05 07:27:53  cblume
136 Changes in the header-files by FCA
137
138 Revision 1.9  2000/10/02 21:28:19  fca
139 Removal of useless dependecies via forward declarations
140
141 Revision 1.8  2000/06/09 11:10:07  cblume
142 Compiler warnings and coding conventions, next round
143
144 Revision 1.7  2000/06/08 18:32:58  cblume
145 Make code compliant to coding conventions
146
147 Revision 1.6  2000/06/07 16:27:32  cblume
148 Try to remove compiler warnings on Sun and HP
149
150 Revision 1.5  2000/05/09 16:38:57  cblume
151 Removed PadResponse(). Merge problem
152
153 Revision 1.4  2000/05/08 15:53:45  cblume
154 Resolved merge conflict
155
156 Revision 1.3  2000/04/28 14:49:27  cblume
157 Only one declaration of iDict in MakeDigits()
158
159 Revision 1.1.4.1  2000/05/08 14:42:04  cblume
160 Introduced AliTRDdigitsManager
161
162 Revision 1.1  2000/02/28 19:00:13  cblume
163 Add new TRD classes
164
165 */
166
167 ///////////////////////////////////////////////////////////////////////////////
168 //                                                                           //
169 //  Creates and handles digits from TRD hits                                 //
170 //  Author: C. Blume (C.Blume@gsi.de)                                        //
171 //                                                                           //
172 //  The following effects are included:                                      //
173 //      - Diffusion                                                          //
174 //      - ExB effects                                                        //
175 //      - Gas gain including fluctuations                                    //
176 //      - Pad-response (simple Gaussian approximation)                       //
177 //      - Time-response                                                      //
178 //      - Electronics noise                                                  //
179 //      - Electronics gain                                                   //
180 //      - Digitization                                                       //
181 //      - ADC threshold                                                      //
182 //  The corresponding parameter can be adjusted via the various              //
183 //  Set-functions. If these parameters are not explicitly set, default       //
184 //  values are used (see Init-function).                                     //
185 //  As an example on how to use this class to produce digits from hits       //
186 //  have a look at the macro hits2digits.C                                   //
187 //  The production of summable digits is demonstrated in hits2sdigits.C      //
188 //  and the subsequent conversion of the s-digits into normal digits is      //
189 //  explained in sdigits2digits.C.                                           //
190 //                                                                           //
191 ///////////////////////////////////////////////////////////////////////////////
192
193 #include <stdlib.h>
194
195 #include <TMath.h>
196 #include <TVector.h>
197 #include <TRandom.h>
198 #include <TROOT.h>
199 #include <TTree.h>
200 #include <TFile.h>
201 #include <TF1.h>
202 #include <TList.h>
203 #include <TTask.h>
204
205 #include "AliRun.h"
206 #include "AliMagF.h"
207 #include "AliRunDigitizer.h"
208
209 #include "AliTRD.h"
210 #include "AliTRDhit.h"
211 #include "AliTRDdigitizer.h"
212 #include "AliTRDdataArrayI.h"
213 #include "AliTRDdataArrayF.h"
214 #include "AliTRDsegmentArray.h"
215 #include "AliTRDdigitsManager.h"
216 #include "AliTRDgeometry.h"
217 #include "AliTRDparameter.h"
218
219 ClassImp(AliTRDdigitizer)
220
221 //_____________________________________________________________________________
222 AliTRDdigitizer::AliTRDdigitizer()
223 {
224   //
225   // AliTRDdigitizer default constructor
226   //
227
228   fInputFile          = 0;
229   fDigitsManager      = 0;
230   fSDigitsManagerList = 0;
231   fSDigitsManager     = 0;
232   fTRD                = 0;
233   fGeo                = 0;
234   fPar                = 0;
235   fMasks              = 0;
236   fEvent              = 0;
237   fCompress           = kTRUE;
238   fDebug              = 0;
239   fSDigits            = kFALSE;
240   fSDigitsScale       = 0.0;
241   fMergeSignalOnly    = kFALSE;
242   fSimpleSim          = kFALSE;
243   fSimpleDet          = 0;
244  
245 }
246
247 //_____________________________________________________________________________
248 AliTRDdigitizer::AliTRDdigitizer(const Text_t *name, const Text_t *title)
249                 :AliDigitizer(name,title)
250 {
251   //
252   // AliTRDdigitizer constructor
253   //
254
255   fInputFile          = 0;
256   fDigitsManager      = 0;
257   fSDigitsManagerList = 0;
258   fSDigitsManager     = 0;
259   fTRD                = 0;
260   fGeo                = 0;
261   fPar                = 0;
262   fMasks              = 0;
263   fEvent              = 0;
264   fCompress           = kTRUE;
265   fDebug              = 0;
266   fSDigits            = kFALSE;
267   fMergeSignalOnly    = kFALSE;
268   fSimpleSim          = kFALSE;
269   fSimpleDet          = 0;
270  
271   // For the summable digits
272   fSDigitsScale       = 100.;
273
274 }
275
276 //_____________________________________________________________________________
277 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager
278                                 , const Text_t *name, const Text_t *title)
279                 :AliDigitizer(manager,name,title)
280 {
281   //
282   // AliTRDdigitizer constructor
283   //
284
285   fInputFile          = 0;
286   fDigitsManager      = 0;
287   fSDigitsManagerList = 0;
288   fSDigitsManager     = 0;
289   fTRD                = 0;
290   fGeo                = 0;
291   fPar                = 0;
292   fMasks              = 0;
293   fEvent              = 0;
294   fCompress           = kTRUE;
295   fDebug              = 0;
296   fSDigits            = kFALSE;
297   fMergeSignalOnly    = kFALSE;
298   fSimpleSim          = kFALSE;
299   fSimpleDet          = 0;
300  
301   // For the summable digits
302   fSDigitsScale       = 100.;
303
304 }
305
306 //_____________________________________________________________________________
307 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager)
308                 :AliDigitizer(manager,"AliTRDdigitizer","TRD digitizer")
309 {
310   //
311   // AliTRDdigitizer constructor
312   //
313
314   fInputFile          = 0;
315   fDigitsManager      = 0;
316   fSDigitsManagerList = 0;
317   fSDigitsManager     = 0;
318   fTRD                = 0;
319   fGeo                = 0;
320   fPar                = 0;
321   fMasks              = 0;
322   fEvent              = 0;
323   fCompress           = kTRUE;
324   fDebug              = 0;
325   fSDigits            = kFALSE;
326   fMergeSignalOnly    = kFALSE;
327   fSimpleSim          = kFALSE;
328   fSimpleDet          = 0;
329
330   // For the summable digits
331   fSDigitsScale       = 100.;
332
333 }
334
335 //_____________________________________________________________________________
336 AliTRDdigitizer::AliTRDdigitizer(const AliTRDdigitizer &d)
337 {
338   //
339   // AliTRDdigitizer copy constructor
340   //
341
342   ((AliTRDdigitizer &) d).Copy(*this);
343
344 }
345
346 //_____________________________________________________________________________
347 AliTRDdigitizer::~AliTRDdigitizer()
348 {
349   //
350   // AliTRDdigitizer destructor
351   //
352
353   if (fInputFile) {
354     fInputFile->Close();
355     delete fInputFile;
356     fInputFile = 0;
357   }
358
359   if (fDigitsManager) {
360     delete fDigitsManager;
361     fDigitsManager = 0;
362   }
363
364   if (fSDigitsManager) {
365     delete fSDigitsManager;
366     fSDigitsManager = 0;
367   }
368
369   if (fSDigitsManagerList) {
370     delete fSDigitsManagerList;
371     fSDigitsManagerList = 0;
372   }
373
374   if (fMasks) {
375     delete [] fMasks;
376     fMasks = 0;
377   }
378
379 }
380
381 //_____________________________________________________________________________
382 AliTRDdigitizer &AliTRDdigitizer::operator=(const AliTRDdigitizer &d)
383 {
384   //
385   // Assignment operator
386   //
387
388   if (this != &d) ((AliTRDdigitizer &) d).Copy(*this);
389   return *this;
390
391 }
392
393 //_____________________________________________________________________________
394 void AliTRDdigitizer::Copy(TObject &d)
395 {
396   //
397   // Copy function
398   //
399
400   ((AliTRDdigitizer &) d).fInputFile          = 0;
401   ((AliTRDdigitizer &) d).fSDigitsManagerList = 0;
402   ((AliTRDdigitizer &) d).fSDigitsManager     = 0;
403   ((AliTRDdigitizer &) d).fDigitsManager      = 0;
404   ((AliTRDdigitizer &) d).fTRD                = 0;
405   ((AliTRDdigitizer &) d).fGeo                = 0;
406   ((AliTRDdigitizer &) d).fMasks              = 0;
407   ((AliTRDdigitizer &) d).fEvent              = 0;
408   ((AliTRDdigitizer &) d).fPar                = 0;
409   ((AliTRDdigitizer &) d).fCompress           = fCompress;
410   ((AliTRDdigitizer &) d).fDebug              = fDebug  ;
411   ((AliTRDdigitizer &) d).fSDigits            = fSDigits;
412   ((AliTRDdigitizer &) d).fSDigitsScale       = fSDigitsScale;
413   ((AliTRDdigitizer &) d).fMergeSignalOnly    = fMergeSignalOnly;
414   ((AliTRDdigitizer &) d).fSimpleSim          = fSimpleSim;
415   ((AliTRDdigitizer &) d).fSimpleDet          = fSimpleDet;
416                                        
417 }
418
419 //_____________________________________________________________________________
420 void AliTRDdigitizer::Exec(Option_t* option)
421 {
422   //
423   // Executes the merging
424   //
425
426   Int_t iInput;
427
428   AliTRDdigitsManager *sdigitsManager;
429
430   TString optionString = option;
431   if (optionString.Contains("deb")) {
432     fDebug = 1;
433     if (optionString.Contains("2")) {
434       fDebug = 2;
435     }
436     printf("<AliTRDdigitizer::Exec> ");
437     printf("Called with debug option %d\n",fDebug);
438   }
439
440   // The AliRoot file is already connected by the manager
441   if (gAlice) {
442     if (fDebug > 0) {
443       printf("<AliTRDdigitizer::Exec> ");
444       printf("AliRun object found on file.\n");
445     }
446   }
447   else {
448     printf("<AliTRDdigitizer::Exec> ");
449     printf("Could not find AliRun object.\n");
450     return;
451   }
452                                                                            
453   Int_t nInput = fManager->GetNinputs();
454   fMasks = new Int_t[nInput];
455   for (iInput = 0; iInput < nInput; iInput++) {
456     fMasks[iInput] = fManager->GetMask(iInput);
457   }
458
459   // Initialization
460   InitDetector();
461
462   for (iInput = 0; iInput < nInput; iInput++) {
463
464     if (fDebug > 0) {
465       printf("<AliTRDdigitizer::Exec> ");
466       printf("Add input stream %d\n",iInput);
467     }
468
469     // check if the input tree exists
470     if (!fManager->GetInputTreeTRDS(iInput)) {
471       printf("<AliTRDdigitizer::Exec> ");
472       printf("Input stream %d does not exist\n",iInput);
473       return;
474     } 
475
476     // Read the s-digits via digits manager
477     sdigitsManager = new AliTRDdigitsManager();
478     sdigitsManager->SetDebug(fDebug);
479     sdigitsManager->SetSDigits(kTRUE);
480     sdigitsManager->ReadDigits(fManager->GetInputTreeTRDS(iInput));
481
482     // Add the s-digits to the input list 
483     AddSDigitsManager(sdigitsManager);
484
485   }
486
487   // Convert the s-digits to normal digits
488   if (fDebug > 0) {
489     printf("<AliTRDdigitizer::Exec> ");
490     printf("Do the conversion\n");
491   }
492   SDigits2Digits();
493
494   // Store the digits
495   if (fDebug > 0) {
496     printf("<AliTRDdigitizer::Exec> ");
497     printf("Write the digits\n");
498   }
499   fDigitsManager->MakeBranch(fManager->GetTreeDTRD());
500   fDigitsManager->WriteDigits();
501   if (fDebug > 0) {
502     printf("<AliTRDdigitizer::Exec> ");
503     printf("Done\n");
504   }
505
506   DeleteSDigitsManager();
507
508 }
509
510 //_____________________________________________________________________________
511 Bool_t AliTRDdigitizer::Open(const Char_t *file, Int_t nEvent)
512 {
513   //
514   // Opens a ROOT-file with TRD-hits and reads in the hit-tree
515   //
516
517   // Connect the AliRoot file containing Geometry, Kine, and Hits
518   fInputFile = (TFile *) gROOT->GetListOfFiles()->FindObject(file);
519   if (!fInputFile) {
520     if (fDebug > 0) {
521       printf("<AliTRDdigitizer::Open> ");
522       printf("Open the AliROOT-file %s.\n",file);
523     }
524     fInputFile = new TFile(file,"UPDATE");
525   }
526   else {
527     if (fDebug > 0) {
528       printf("<AliTRDdigitizer::Open> ");
529       printf("%s is already open.\n",file);
530     }
531   }
532
533   gAlice = (AliRun *) fInputFile->Get("gAlice");
534   if (gAlice) {
535     if (fDebug > 0) {
536       printf("<AliTRDdigitizer::Open> ");
537       printf("AliRun object found on file.\n");
538     }
539   }
540   else {
541     printf("<AliTRDdigitizer::Open> ");
542     printf("Could not find AliRun object.\n");
543     return kFALSE;
544   }
545
546   fEvent = nEvent;
547
548   // Import the Trees for the event nEvent in the file
549   Int_t nparticles = gAlice->GetEvent(fEvent);
550   if (nparticles <= 0) {
551     printf("<AliTRDdigitizer::Open> ");
552     printf("No entries in the trees for event %d.\n",fEvent);
553     return kFALSE;
554   }
555
556   if (InitDetector()) {
557     return MakeBranch();
558   }
559   else {
560     return kFALSE;
561   }
562
563 }
564
565 //_____________________________________________________________________________
566 Bool_t AliTRDdigitizer::InitDetector()
567 {
568   //
569   // Sets the pointer to the TRD detector and the geometry
570   //
571
572   // Get the pointer to the detector class and check for version 1
573   fTRD = (AliTRD *) gAlice->GetDetector("TRD");
574   if (!fTRD) {
575     printf("<AliTRDdigitizer::InitDetector> ");
576     printf("No TRD module found\n");
577     exit(1);
578   }
579   if (fTRD->IsVersion() != 1) {
580     printf("<AliTRDdigitizer::InitDetector> ");
581     printf("TRD must be version 1 (slow simulator).\n");
582     exit(1);
583   }
584
585   // Get the geometry
586   fGeo = fTRD->GetGeometry();
587   if (fDebug > 0) {
588     printf("<AliTRDdigitizer::InitDetector> ");
589     printf("Geometry version %d\n",fGeo->IsVersion());
590   }
591
592   // Create a digits manager
593   fDigitsManager = new AliTRDdigitsManager();
594   fDigitsManager->SetSDigits(fSDigits);
595   fDigitsManager->CreateArrays();
596   fDigitsManager->SetEvent(fEvent);
597   fDigitsManager->SetDebug(fDebug);
598
599   // The list for the input s-digits manager to be merged
600   fSDigitsManagerList = new TList();
601
602   return kTRUE;
603
604 }
605
606 //_____________________________________________________________________________
607 Bool_t AliTRDdigitizer::MakeBranch(const Char_t *file) const
608 {
609   // 
610   // Create the branches for the digits array
611   //
612
613   return fDigitsManager->MakeBranch(file);
614
615 }
616
617 //_____________________________________________________________________________
618 Bool_t AliTRDdigitizer::MakeDigits()
619 {
620   //
621   // Creates digits.
622   //
623
624   ///////////////////////////////////////////////////////////////
625   // Parameter 
626   ///////////////////////////////////////////////////////////////
627
628   // Converts number of electrons to fC
629   const Double_t kEl2fC  = 1.602E-19 * 1.0E15; 
630
631   ///////////////////////////////////////////////////////////////
632
633   // Number of pads included in the pad response
634   const Int_t kNpad  = 3;
635
636   // Number of track dictionary arrays
637   const Int_t kNDict = AliTRDdigitsManager::kNDict;
638
639   // Half the width of the amplification region
640   const Float_t kAmWidth = AliTRDgeometry::AmThick() / 2.;
641
642   Int_t   iRow, iCol, iTime, iPad;
643   Int_t   iDict  = 0;
644   Int_t   nBytes = 0;
645
646   Int_t   totalSizeDigits = 0;
647   Int_t   totalSizeDict0  = 0;
648   Int_t   totalSizeDict1  = 0;
649   Int_t   totalSizeDict2  = 0;
650
651   Int_t   timeTRDbeg = 0;
652   Int_t   timeTRDend = 1;
653
654   Float_t pos[3];
655   Float_t rot[3];
656   Float_t xyz[3];
657   Float_t padSignal[kNpad];
658   Float_t signalOld[kNpad];
659
660   AliTRDdataArrayF *signals = 0;
661   AliTRDdataArrayI *digits  = 0;
662   AliTRDdataArrayI *dictionary[kNDict];
663
664   // Create a default parameter class if none is defined
665   if (!fPar) {
666     fPar = new AliTRDparameter("TRDparameter","Standard TRD parameter");
667     if (fDebug > 0) {
668       printf("<AliTRDdigitizer::MakeDigits> ");
669       printf("Create the default parameter object\n");
670     }
671   }
672
673   // Create a container for the amplitudes
674   AliTRDsegmentArray *signalsArray 
675                      = new AliTRDsegmentArray("AliTRDdataArrayF"
676                                              ,AliTRDgeometry::Ndet());
677
678   if (fPar->TRFOn()) {
679     timeTRDbeg = ((Int_t) (-fPar->GetTRFlo() / fPar->GetTimeBinSize())) - 1;
680     timeTRDend = ((Int_t) ( fPar->GetTRFhi() / fPar->GetTimeBinSize())) - 1;
681     if (fDebug > 0) {
682       printf("<AliTRDdigitizer::MakeDigits> ");
683       printf("Sample the TRF between -%d and %d\n",timeTRDbeg,timeTRDend);
684     }
685   }
686
687   Float_t elAttachProp = fPar->GetElAttachProp() / 100.; 
688
689   if (!fGeo) {
690     printf("<AliTRDdigitizer::MakeDigits> ");
691     printf("No geometry defined\n");
692     return kFALSE;
693   }
694
695   if (fDebug > 0) {
696     printf("<AliTRDdigitizer::MakeDigits> ");
697     printf("Start creating digits.\n");
698   }
699
700   // Get the pointer to the hit tree
701   TTree *hitTree = gAlice->TreeH();
702
703   // Get the number of entries in the hit tree
704   // (Number of primary particles creating a hit somewhere)
705   Int_t nTrack = 1;
706   if (!fSimpleSim) {
707     nTrack = (Int_t) hitTree->GetEntries();
708     if (fDebug > 0) {
709       printf("<AliTRDdigitizer::MakeDigits> ");
710       printf("Found %d primary particles\n",nTrack);
711     } 
712   }
713
714   Int_t detectorOld = -1;
715   Int_t countHits   =  0; 
716
717   // Loop through all entries in the tree
718   for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
719
720     if (!fSimpleSim) {   
721       gAlice->ResetHits();
722       nBytes += hitTree->GetEvent(iTrack);
723     }
724
725     // Loop through the TRD hits
726     Int_t iHit = 0;
727     AliTRDhit *hit = (AliTRDhit *) fTRD->FirstHit(-1);
728     while (hit) {
729  
730       countHits++;
731       iHit++;
732
733               pos[0]      = hit->X();
734               pos[1]      = hit->Y();
735               pos[2]      = hit->Z();
736       Float_t q           = hit->GetCharge();
737       Int_t   track       = hit->Track();
738       Int_t   detector    = hit->GetDetector();
739       Int_t   plane       = fGeo->GetPlane(detector);
740       Int_t   sector      = fGeo->GetSector(detector);
741       Int_t   chamber     = fGeo->GetChamber(detector);
742       Int_t   nRowMax     = fPar->GetRowMax(plane,chamber,sector);
743       Int_t   nColMax     = fPar->GetColMax(plane);
744       Int_t   nTimeMax    = fPar->GetTimeMax();
745       Int_t   nTimeBefore = fPar->GetTimeBefore();
746       Int_t   nTimeAfter  = fPar->GetTimeAfter();
747       Int_t   nTimeTotal  = fPar->GetTimeTotal();
748       Float_t row0        = fPar->GetRow0(plane,chamber,sector);
749       Float_t col0        = fPar->GetCol0(plane);
750       Float_t time0       = fPar->GetTime0(plane);
751       Float_t rowPadSize  = fPar->GetRowPadSize(plane,chamber,sector);
752       Float_t colPadSize  = fPar->GetColPadSize(plane);
753       Float_t timeBinSize = fPar->GetTimeBinSize();
754       Float_t divideRow   = 1.0 / rowPadSize;
755       Float_t divideCol   = 1.0 / colPadSize;
756       Float_t divideTime  = 1.0 / timeBinSize;
757
758       if (fDebug > 1) {
759         printf("Analyze hit no. %d ",iHit);
760         printf("-----------------------------------------------------------\n");
761         hit->Dump();
762         printf("plane = %d, sector = %d, chamber = %d\n"
763               ,plane,sector,chamber);
764         printf("nRowMax = %d, nColMax = %d, nTimeMax = %d\n" 
765               ,nRowMax,nColMax,nTimeMax);
766         printf("nTimeBefore = %d, nTimeAfter = %d, nTimeTotal = %d\n"
767               ,nTimeBefore,nTimeAfter,nTimeTotal);
768         printf("row0 = %f, col0 = %f, time0 = %f\n"
769               ,row0,col0,time0);
770         printf("rowPadSize = %f, colPadSize = %f, timeBinSize = %f\n"
771                ,rowPadSize,colPadSize,timeBinSize); 
772       }
773        
774       // Don't analyze test hits and switched off detectors
775       if ((CheckDetector(plane,chamber,sector)) &&
776           (((Int_t) q) != 0)) {
777
778         if (detector != detectorOld) {
779
780           if (fDebug > 1) {
781             printf("<AliTRDdigitizer::MakeDigits> ");
782             printf("Get new container. New det = %d, Old det = %d\n"
783                   ,detector,detectorOld);
784           }
785           // Compress the old one if enabled
786           if ((fCompress) && (detectorOld > -1)) {
787             if (fDebug > 1) {
788               printf("<AliTRDdigitizer::MakeDigits> ");
789               printf("Compress the old container ...");
790             }
791             signals->Compress(1,0);
792             for (iDict = 0; iDict < kNDict; iDict++) {
793               dictionary[iDict]->Compress(1,0);
794             }
795             if (fDebug > 1) printf("done\n");
796           }
797           // Get the new container
798           signals = (AliTRDdataArrayF *) signalsArray->At(detector);
799           if (signals->GetNtime() == 0) {
800             // Allocate a new one if not yet existing
801             if (fDebug > 1) {
802               printf("<AliTRDdigitizer::MakeDigits> ");
803               printf("Allocate a new container ... ");
804             }
805             signals->Allocate(nRowMax,nColMax,nTimeTotal);
806           }
807           else if (fSimpleSim) {
808             // Clear an old one for the simple simulation
809             if (fDebug > 1) {
810               printf("<AliTRDdigitizer::MakeDigits> ");
811               printf("Clear a old container ... ");
812             }
813             signals->Clear();
814           }
815           else {
816             // Expand an existing one
817             if (fCompress) {
818               if (fDebug > 1) {
819                 printf("<AliTRDdigitizer::MakeDigits> ");
820                 printf("Expand an existing container ... ");
821               }
822               signals->Expand();
823             }
824           }
825           // The same for the dictionary
826           if (!fSimpleSim) {       
827             for (iDict = 0; iDict < kNDict; iDict++) {       
828               dictionary[iDict] = fDigitsManager->GetDictionary(detector,iDict);
829               if (dictionary[iDict]->GetNtime() == 0) {
830                 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
831               }
832               else {
833                 if (fCompress) dictionary[iDict]->Expand();
834               }
835             }
836           }      
837           if (fDebug > 1) printf("done\n");
838           detectorOld = detector;
839         }
840
841         // Rotate the sectors on top of each other
842         if (fSimpleSim) {
843           rot[0] = pos[0];
844           rot[1] = pos[1];
845           rot[2] = pos[2];
846         }
847         else {
848           fGeo->Rotate(detector,pos,rot);
849         }
850
851         // The driftlength. It is negative if the hit is in the 
852         // amplification region.
853         Float_t driftlength = time0 - rot[0];
854
855         // Take also the drift in the amplification region into account
856         // The drift length is at the moment still the same, regardless of
857         // the position relativ to the wire. This non-isochronity needs still
858         // to be implemented.
859         Float_t driftlengthL = TMath::Abs(driftlength + kAmWidth);
860         if (fPar->ExBOn()) driftlengthL /= TMath::Sqrt(fPar->GetLorentzFactor());
861
862         // Loop over all electrons of this hit
863         // TR photons produce hits with negative charge
864         Int_t nEl = ((Int_t) TMath::Abs(q));
865         for (Int_t iEl = 0; iEl < nEl; iEl++) {
866
867           xyz[0] = rot[0];
868           xyz[1] = rot[1];
869           xyz[2] = rot[2];
870
871           // Electron attachment
872           if (fPar->ElAttachOn()) {
873             if (gRandom->Rndm() < (driftlengthL * elAttachProp)) 
874               continue;
875           }
876
877           // Apply the diffusion smearing
878           if (fPar->DiffusionOn()) {
879             if (!(fPar->Diffusion(driftlengthL,xyz))) continue;
880           }
881
882           // Apply E x B effects (depends on drift direction)
883           if (fPar->ExBOn()) { 
884             if (!(fPar->ExB(driftlength+kAmWidth,xyz))) continue;   
885           }
886
887           // The electron position after diffusion and ExB in pad coordinates 
888           // The pad row (z-direction)
889           Float_t rowDist   = xyz[2] - row0;
890           Int_t   rowE      = ((Int_t) (rowDist * divideRow));
891           if ((rowE < 0) || (rowE >= nRowMax)) continue;   
892           Float_t rowOffset = ((((Float_t) rowE) + 0.5) * rowPadSize) - rowDist;
893
894           // The pad column (rphi-direction)
895           Float_t col0tilt  = fPar->Col0Tilted(col0,rowOffset,plane);
896           Float_t colDist   = xyz[1] - col0tilt;
897           Int_t   colE      = ((Int_t) (colDist * divideCol));
898           if ((colE < 0) || (colE >= nColMax)) continue;   
899           Float_t colOffset = ((((Float_t) colE) + 0.5) * colPadSize) - colDist;    
900
901           // The time bin (negative for hits in the amplification region)
902           // In the amplification region the electrons drift from both sides
903           // to the middle (anode wire plane)
904           Float_t timeDist   = time0 - xyz[0];
905           Float_t timeOffset = 0;
906           Int_t   timeE      = 0;
907           if (timeDist > 0) {
908             // The time bin
909             timeE      = ((Int_t) (timeDist * divideTime));
910             // The distance of the position to the middle of the timebin
911             timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) - timeDist;
912           }
913           else {
914             // Difference between half of the amplification gap width and
915             // the distance to the anode wire
916             Float_t anodeDist = kAmWidth - TMath::Abs(timeDist + kAmWidth);
917             // The time bin
918             timeE      = -1 * (((Int_t ) (anodeDist * divideTime)) + 1);
919             // The distance of the position to the middle of the timebin
920             timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) + anodeDist;
921           }
922  
923           // Apply the gas gain including fluctuations
924           Float_t ggRndm = 0.0;
925           do {
926             ggRndm = gRandom->Rndm();
927           } while (ggRndm <= 0);
928           Int_t signal = (Int_t) (-fPar->GetGasGain() * TMath::Log(ggRndm));
929
930           // Apply the pad response 
931           if (fPar->PRFOn()) {
932             // The distance of the electron to the center of the pad 
933             // in units of pad width
934             Float_t dist = - colOffset * divideCol;
935             if (!(fPar->PadResponse(signal,dist,plane,padSignal))) continue;
936           }
937           else {
938             padSignal[0] = 0.0;
939             padSignal[1] = signal;
940             padSignal[2] = 0.0;
941           }
942
943           // Sample the time response inside the drift region
944           // + additional time bins before and after.
945           // The sampling is done always in the middle of the time bin
946           for (Int_t iTimeBin = TMath::Max(timeE-timeTRDbeg,        -nTimeBefore) 
947                     ;iTimeBin < TMath::Min(timeE+timeTRDend,nTimeMax+nTimeAfter ) 
948                     ;iTimeBin++) {
949
950             // Apply the time response
951             Float_t timeResponse = 1.0;
952             Float_t crossTalk    = 0.0;
953             Float_t time         = (iTimeBin - timeE) * timeBinSize + timeOffset;
954             if (fPar->TRFOn()) {
955               timeResponse = fPar->TimeResponse(time);
956             }
957             if (fPar->CTOn()) {
958               crossTalk    = fPar->CrossTalk(time);
959             }
960
961             signalOld[0] = 0.0;
962             signalOld[1] = 0.0;
963             signalOld[2] = 0.0;
964
965             for (iPad = 0; iPad < kNpad; iPad++) {
966
967               Int_t colPos = colE + iPad - 1;
968               if (colPos <        0) continue;
969               if (colPos >= nColMax) break;
970
971               // Add the signals
972               // Note: The time bin number is shifted by nTimeBefore to avoid negative
973               // time bins. This has to be subtracted later.
974               Int_t iCurrentTimeBin = iTimeBin + nTimeBefore;
975               signalOld[iPad]  = signals->GetDataUnchecked(rowE,colPos,iCurrentTimeBin);
976               if( colPos != colE ) {
977                 signalOld[iPad] += padSignal[iPad] * (timeResponse + crossTalk);
978               } 
979               else {
980                 signalOld[iPad] += padSignal[iPad] * timeResponse;
981               }
982               signals->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,signalOld[iPad]);
983
984               // Store the track index in the dictionary
985               // Note: We store index+1 in order to allow the array to be compressed
986               if ((signalOld[iPad] > 0) && (!fSimpleSim)) { 
987                 for (iDict = 0; iDict < kNDict; iDict++) {
988                   Int_t oldTrack = dictionary[iDict]->GetDataUnchecked(rowE
989                                                                       ,colPos
990                                                                       ,iCurrentTimeBin);
991                   if (oldTrack == track+1) break;
992                   if (oldTrack ==       0) {
993                     dictionary[iDict]->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,track+1);
994                     break;
995                   }
996                 }
997               }
998
999             } // Loop: pads
1000
1001           } // Loop: time bins
1002
1003         } // Loop: electrons of a single hit
1004
1005       } // If: detector and test hit
1006
1007       hit = (AliTRDhit *) fTRD->NextHit();   
1008
1009     } // Loop: hits of one primary track
1010
1011   } // Loop: primary tracks
1012
1013   if (fDebug > 0) {
1014     printf("<AliTRDdigitizer::MakeDigits> ");
1015     printf("Finished analyzing %d hits\n",countHits);
1016   }
1017
1018   // The coupling factor
1019   Float_t coupling = fPar->GetPadCoupling() 
1020                    * fPar->GetTimeCoupling();
1021
1022   // The conversion factor
1023   Float_t convert  = kEl2fC
1024                    * fPar->GetChipGain();
1025
1026   // Loop through all chambers to finalize the digits
1027   Int_t iDetBeg = 0;
1028   Int_t iDetEnd = AliTRDgeometry::Ndet();
1029   if (fSimpleSim) {
1030     iDetBeg = fSimpleDet;
1031     iDetEnd = iDetBeg + 1;
1032   }
1033   for (Int_t iDet = iDetBeg; iDet < iDetEnd; iDet++) {
1034
1035     Int_t plane       = fGeo->GetPlane(iDet);
1036     Int_t sector      = fGeo->GetSector(iDet);
1037     Int_t chamber     = fGeo->GetChamber(iDet);
1038     Int_t nRowMax     = fPar->GetRowMax(plane,chamber,sector);
1039     Int_t nColMax     = fPar->GetColMax(plane);
1040     Int_t nTimeMax    = fPar->GetTimeMax();
1041     Int_t nTimeTotal  = fPar->GetTimeTotal();
1042
1043     Double_t *inADC  = new Double_t[nTimeTotal];
1044     Double_t *outADC = new Double_t[nTimeTotal];
1045
1046     if (fDebug > 0) {
1047       printf("<AliTRDdigitizer::MakeDigits> ");
1048       printf("Digitization for chamber %d\n",iDet);
1049     }
1050
1051     // Add a container for the digits of this detector
1052     digits = fDigitsManager->GetDigits(iDet);        
1053     // Allocate memory space for the digits buffer
1054     if (digits->GetNtime() == 0) {
1055       digits->Allocate(nRowMax,nColMax,nTimeTotal);
1056     }
1057     else if (fSimpleSim) {
1058       digits->Clear();
1059     }
1060  
1061     // Get the signal container
1062     signals = (AliTRDdataArrayF *) signalsArray->At(iDet);
1063     if (signals->GetNtime() == 0) {
1064       // Create missing containers
1065       signals->Allocate(nRowMax,nColMax,nTimeTotal);      
1066     }
1067     else {
1068       // Expand the container if neccessary
1069       if (fCompress) signals->Expand();
1070     }
1071     // Create the missing dictionary containers
1072     if (!fSimpleSim) {    
1073       for (iDict = 0; iDict < kNDict; iDict++) {       
1074         dictionary[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1075         if (dictionary[iDict]->GetNtime() == 0) {
1076           dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1077         }
1078       } 
1079     }
1080
1081     Int_t nDigits = 0;
1082
1083     // Don't create noise in detectors that are switched off
1084     if (CheckDetector(plane,chamber,sector)) {
1085
1086       // Create the digits for this chamber
1087       for (iRow  = 0; iRow  <  nRowMax;   iRow++ ) {
1088         for (iCol  = 0; iCol  <  nColMax;   iCol++ ) {
1089
1090           // Create summable digits
1091           if (fSDigits) {
1092
1093             for (iTime = 0; iTime < nTimeTotal; iTime++) {         
1094               Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1095               signalAmp *= fSDigitsScale;
1096               signalAmp  = TMath::Min(signalAmp,(Float_t) 1.0e9);
1097               Int_t adc  = (Int_t) signalAmp;
1098               if (adc > 0) nDigits++;
1099               digits->SetDataUnchecked(iRow,iCol,iTime,adc);
1100             }
1101
1102           }
1103           // Create normal digits
1104           else {
1105
1106             for (iTime = 0; iTime < nTimeTotal; iTime++) {         
1107               Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1108               // Pad and time coupling
1109               signalAmp *= coupling;
1110               // Add the noise, starting from minus ADC baseline in electrons
1111               Double_t baselineEl = fPar->GetADCbaseline() * (fPar->GetADCinRange()
1112                                                            / fPar->GetADCoutRange()) 
1113                                                            / convert;
1114               signalAmp  = TMath::Max((Double_t) gRandom->Gaus(signalAmp,fPar->GetNoise())
1115                                      ,-baselineEl);
1116               // Convert to mV
1117               signalAmp *= convert;
1118               // Add ADC baseline in mV
1119               signalAmp += fPar->GetADCbaseline() * (fPar->GetADCinRange()
1120                                                    / fPar->GetADCoutRange());
1121               // Convert to ADC counts. Set the overflow-bit fADCoutRange if the 
1122               // signal is larger than fADCinRange
1123               Int_t adc  = 0;
1124               if (signalAmp >= fPar->GetADCinRange()) {
1125                 adc = ((Int_t) fPar->GetADCoutRange());
1126               }
1127               else {
1128                 adc = ((Int_t) (signalAmp * (fPar->GetADCoutRange() 
1129                                            / fPar->GetADCinRange())));
1130               }
1131               inADC[iTime]  = adc;
1132               outADC[iTime] = adc;
1133             }
1134
1135             // Apply the tail cancelation via the digital filter
1136             if (fPar->TCOn()) {
1137               DeConvExp(inADC,outADC,nTimeTotal,fPar->GetTCnexp());
1138             }
1139
1140             for (iTime = 0; iTime < nTimeTotal; iTime++) {   
1141               // Store the amplitude of the digit if above threshold
1142               if (outADC[iTime] > fPar->GetADCthreshold()) {
1143                 if (fDebug > 2) {
1144                   printf("  iRow = %d, iCol = %d, iTime = %d, adc = %f\n"
1145                         ,iRow,iCol,iTime,outADC[iTime]);
1146                 }
1147                 nDigits++;
1148                 digits->SetDataUnchecked(iRow,iCol,iTime,((Int_t) outADC[iTime]));
1149               }
1150             }
1151
1152           }
1153
1154         }
1155       }
1156
1157     }
1158
1159     // Compress the arrays
1160     if (!fSimpleSim) {  
1161       digits->Compress(1,0);
1162       for (iDict = 0; iDict < kNDict; iDict++) {
1163         dictionary[iDict]->Compress(1,0);
1164       }
1165
1166       totalSizeDigits += digits->GetSize();
1167       totalSizeDict0  += dictionary[0]->GetSize();
1168       totalSizeDict1  += dictionary[1]->GetSize();
1169       totalSizeDict2  += dictionary[2]->GetSize();
1170
1171       Float_t nPixel = nRowMax * nColMax * nTimeMax;
1172       if (fDebug > 0) {
1173         printf("<AliTRDdigitizer::MakeDigits> ");
1174         printf("Found %d digits in detector %d (%3.0f).\n"
1175               ,nDigits,iDet
1176               ,100.0 * ((Float_t) nDigits) / nPixel);
1177       } 
1178
1179       if (fCompress) signals->Compress(1,0);
1180
1181     }
1182
1183     delete [] inADC;
1184     delete [] outADC;
1185
1186   }
1187
1188   if (signalsArray) {
1189     delete signalsArray;
1190     signalsArray = 0;
1191   }
1192
1193   if (fDebug > 0) {
1194     printf("<AliTRDdigitizer::MakeDigits> ");
1195     printf("Total number of analyzed hits = %d\n",countHits);
1196     if (!fSimpleSim) {    
1197       printf("<AliTRDdigitizer::MakeDigits> ");
1198       printf("Total digits data size = %d, %d, %d, %d\n",totalSizeDigits
1199                                                         ,totalSizeDict0
1200                                                         ,totalSizeDict1
1201                                                         ,totalSizeDict2);        
1202     }
1203   }
1204
1205   return kTRUE;
1206
1207 }
1208
1209 //_____________________________________________________________________________
1210 void AliTRDdigitizer::AddSDigitsManager(AliTRDdigitsManager *man)
1211 {
1212   //
1213   // Add a digits manager for s-digits to the input list.
1214   //
1215
1216   fSDigitsManagerList->Add(man);
1217
1218 }
1219
1220 //_____________________________________________________________________________
1221 void AliTRDdigitizer::DeleteSDigitsManager()
1222 {
1223   //
1224   // Removes digits manager from the input list.
1225   //
1226
1227   fSDigitsManagerList->Delete();
1228
1229 }
1230
1231 //_____________________________________________________________________________
1232 Bool_t AliTRDdigitizer::ConvertSDigits()
1233 {
1234   //
1235   // Converts s-digits to normal digits
1236   //
1237
1238   // Number of track dictionary arrays
1239   const Int_t    kNDict = AliTRDdigitsManager::kNDict;
1240
1241   // Converts number of electrons to fC
1242   const Double_t kEl2fC = 1.602E-19 * 1.0E15; 
1243
1244   Int_t iDict = 0;
1245   Int_t iRow;
1246   Int_t iCol;
1247   Int_t iTime;
1248
1249   if (!fPar) {    
1250     fPar = new AliTRDparameter("TRDparameter","Standard parameter");
1251     if (fDebug > 0) {
1252       printf("<AliTRDdigitizer::ConvertSDigits> ");
1253       printf("Create the default parameter object\n");
1254     }
1255   }
1256
1257   Double_t sDigitsScale = 1.0 / GetSDigitsScale();
1258   Double_t noise        = fPar->GetNoise();
1259   Double_t padCoupling  = fPar->GetPadCoupling();
1260   Double_t timeCoupling = fPar->GetTimeCoupling();
1261   Double_t chipGain     = fPar->GetChipGain();
1262   Double_t coupling     = padCoupling * timeCoupling;
1263   Double_t convert      = kEl2fC * chipGain;
1264   Double_t adcInRange   = fPar->GetADCinRange();
1265   Double_t adcOutRange  = fPar->GetADCoutRange();
1266   Int_t    adcThreshold = fPar->GetADCthreshold();
1267   Int_t    adcBaseline  = fPar->GetADCbaseline();   
1268
1269   AliTRDdataArrayI *digitsIn;
1270   AliTRDdataArrayI *digitsOut;
1271   AliTRDdataArrayI *dictionaryIn[kNDict];
1272   AliTRDdataArrayI *dictionaryOut[kNDict];
1273
1274   // Loop through the detectors
1275   for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1276
1277     if (fDebug > 0) {
1278       printf("<AliTRDdigitizer::ConvertSDigits> ");
1279       printf("Convert detector %d to digits.\n",iDet);
1280     }
1281
1282     Int_t plane      = fGeo->GetPlane(iDet);
1283     Int_t sector     = fGeo->GetSector(iDet);
1284     Int_t chamber    = fGeo->GetChamber(iDet);
1285     Int_t nRowMax    = fPar->GetRowMax(plane,chamber,sector);
1286     Int_t nColMax    = fPar->GetColMax(plane);
1287     Int_t nTimeTotal = fPar->GetTimeTotal();
1288
1289     Double_t *inADC  = new Double_t[nTimeTotal];
1290     Double_t *outADC = new Double_t[nTimeTotal];
1291
1292     digitsIn  = fSDigitsManager->GetDigits(iDet);
1293     digitsIn->Expand();
1294     digitsOut = fDigitsManager->GetDigits(iDet);
1295     digitsOut->Allocate(nRowMax,nColMax,nTimeTotal);
1296     for (iDict = 0; iDict < kNDict; iDict++) {
1297       dictionaryIn[iDict]  = fSDigitsManager->GetDictionary(iDet,iDict);
1298       dictionaryIn[iDict]->Expand();
1299       dictionaryOut[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1300       dictionaryOut[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1301     }
1302
1303     for (iRow  = 0; iRow  <  nRowMax;   iRow++ ) {
1304       for (iCol  = 0; iCol  <  nColMax;   iCol++ ) {
1305
1306         for (iTime = 0; iTime < nTimeTotal; iTime++) {         
1307           Double_t signal = (Double_t) digitsIn->GetDataUnchecked(iRow,iCol,iTime);
1308           signal *= sDigitsScale;
1309           // Pad and time coupling
1310           signal *= coupling;
1311           // Add the noise, starting from minus ADC baseline in electrons
1312           Double_t baselineEl = adcBaseline * (adcInRange / adcOutRange) / convert;
1313           signal  = TMath::Max((Double_t) gRandom->Gaus(signal,noise),-baselineEl);
1314           // Convert to mV
1315           signal *= convert;
1316           // add ADC baseline in mV
1317           signal += adcBaseline * (adcInRange / adcOutRange);
1318           // Convert to ADC counts. Set the overflow-bit adcOutRange if the 
1319           // signal is larger than adcInRange
1320           Int_t adc  = 0;
1321           if (signal >= adcInRange) {
1322             adc = ((Int_t) adcOutRange);
1323           }
1324           else {
1325             adc = ((Int_t) (signal * (adcOutRange / adcInRange)));
1326           }
1327           inADC[iTime]  = adc;
1328           outADC[iTime] = adc;
1329         }
1330
1331         // Apply the tail cancelation via the digital filter
1332         if (fPar->TCOn()) {
1333           DeConvExp(inADC,outADC,nTimeTotal,fPar->GetTCnexp());
1334         }
1335
1336         for (iTime = 0; iTime < nTimeTotal; iTime++) {   
1337           // Store the amplitude of the digit if above threshold
1338           if (outADC[iTime] > adcThreshold) {
1339             digitsOut->SetDataUnchecked(iRow,iCol,iTime,((Int_t) outADC[iTime]));
1340             // Copy the dictionary
1341             for (iDict = 0; iDict < kNDict; iDict++) { 
1342               Int_t track = dictionaryIn[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1343               dictionaryOut[iDict]->SetDataUnchecked(iRow,iCol,iTime,track);
1344             }
1345           }
1346         }
1347
1348       }
1349     }
1350
1351     if (fCompress) {
1352       digitsIn->Compress(1,0);
1353       digitsOut->Compress(1,0);
1354       for (iDict = 0; iDict < kNDict; iDict++) {
1355         dictionaryIn[iDict]->Compress(1,0);
1356         dictionaryOut[iDict]->Compress(1,0);
1357       }
1358     }
1359
1360     delete [] inADC;
1361     delete [] outADC;
1362
1363   }    
1364
1365   return kTRUE;
1366
1367 }
1368
1369 //_____________________________________________________________________________
1370 Bool_t AliTRDdigitizer::MergeSDigits()
1371 {
1372   //
1373   // Merges the input s-digits:
1374   //   - The amplitude of the different inputs are summed up.
1375   //   - Of the track IDs from the input dictionaries only one is
1376   //     kept for each input. This works for maximal 3 different merged inputs.
1377   //
1378
1379   // Number of track dictionary arrays
1380   const Int_t kNDict = AliTRDdigitsManager::kNDict;
1381
1382   if (!fPar) {
1383     fPar = new AliTRDparameter("TRDparameter","Standard parameter");
1384     if (fDebug > 0) {
1385       printf("<AliTRDdigitizer::MergeSDigits> ");
1386       printf("Create the default parameter object\n");
1387     }
1388   }
1389
1390   Int_t iDict = 0;
1391   Int_t jDict = 0;
1392
1393   AliTRDdataArrayI *digitsA;
1394   AliTRDdataArrayI *digitsB;
1395   AliTRDdataArrayI *dictionaryA[kNDict];
1396   AliTRDdataArrayI *dictionaryB[kNDict];
1397
1398   // Get the first s-digits
1399   fSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->First();
1400   if (!fSDigitsManager) return kFALSE;
1401
1402   // Loop through the other sets of s-digits
1403   AliTRDdigitsManager *mergeSDigitsManager;
1404   mergeSDigitsManager = (AliTRDdigitsManager *) 
1405                         fSDigitsManagerList->After(fSDigitsManager);
1406
1407   if (fDebug > 0) {
1408     if (mergeSDigitsManager) {
1409       printf("<AliTRDdigitizer::MergeSDigits> ");
1410       printf("Merge %d input files.\n",fSDigitsManagerList->GetSize());
1411     }
1412     else {
1413       printf("<AliTRDdigitizer::MergeSDigits> ");
1414       printf("Only one input file.\n");
1415     }
1416   }
1417
1418   Int_t iMerge = 0;
1419   while (mergeSDigitsManager) {
1420
1421     iMerge++;
1422
1423     // Loop through the detectors
1424     for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1425
1426       Int_t plane      = fGeo->GetPlane(iDet);
1427       Int_t sector     = fGeo->GetSector(iDet);
1428       Int_t chamber    = fGeo->GetChamber(iDet);
1429       Int_t nRowMax    = fPar->GetRowMax(plane,chamber,sector);
1430       Int_t nColMax    = fPar->GetColMax(plane);
1431       Int_t nTimeTotal = fPar->GetTimeTotal();
1432
1433       // Loop through the pixels of one detector and add the signals
1434       digitsA = fSDigitsManager->GetDigits(iDet);
1435       digitsB = mergeSDigitsManager->GetDigits(iDet);
1436       digitsA->Expand();
1437       digitsB->Expand();
1438       for (iDict = 0; iDict < kNDict; iDict++) {
1439         dictionaryA[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1440         dictionaryB[iDict] = mergeSDigitsManager->GetDictionary(iDet,iDict);
1441         dictionaryA[iDict]->Expand();
1442         dictionaryB[iDict]->Expand();
1443       }
1444
1445       // Merge only detectors that contain a signal
1446       Bool_t doMerge = kTRUE;
1447       if (fMergeSignalOnly) {
1448         if (digitsA->GetOverThreshold(0) == 0) {
1449           doMerge = kFALSE;
1450         }
1451       }
1452
1453       if (doMerge) {
1454
1455         if (fDebug > 0) {
1456           printf("<AliTRDdigitizer::MergeSDigits> ");
1457           printf("Merge detector %d of input no.%d\n",iDet,iMerge+1);
1458         }
1459
1460         for (Int_t iRow  = 0; iRow  <  nRowMax;   iRow++ ) {
1461           for (Int_t iCol  = 0; iCol  <  nColMax;   iCol++ ) {
1462             for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {         
1463
1464               // Add the amplitudes of the summable digits 
1465               Int_t ampA = digitsA->GetDataUnchecked(iRow,iCol,iTime);
1466               Int_t ampB = digitsB->GetDataUnchecked(iRow,iCol,iTime);
1467               ampA += ampB;
1468               digitsA->SetDataUnchecked(iRow,iCol,iTime,ampA);
1469
1470              // Add the mask to the track id if defined.
1471               for (iDict = 0; iDict < kNDict; iDict++) {
1472                 Int_t trackB = dictionaryB[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1473                 if ((fMasks) && (trackB > 0)) {
1474                   for (jDict = 0; jDict < kNDict; jDict++) { 
1475                     Int_t trackA = dictionaryA[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1476                     if (trackA == 0) {
1477                       trackA = trackB + fMasks[iMerge];
1478                       dictionaryA[iDict]->SetDataUnchecked(iRow,iCol,iTime,trackA);
1479                     }
1480                   }
1481                 }
1482               }
1483
1484             }
1485           }
1486         }
1487
1488       }
1489
1490       if (fCompress) {
1491         digitsA->Compress(1,0);
1492         digitsB->Compress(1,0);
1493         for (iDict = 0; iDict < kNDict; iDict++) {
1494           dictionaryA[iDict]->Compress(1,0);
1495           dictionaryB[iDict]->Compress(1,0);
1496         }
1497       }
1498
1499     }    
1500
1501     // The next set of s-digits
1502     mergeSDigitsManager = (AliTRDdigitsManager *) 
1503                           fSDigitsManagerList->After(mergeSDigitsManager);
1504
1505   }
1506
1507   return kTRUE;
1508
1509 }
1510
1511 //_____________________________________________________________________________
1512 Bool_t AliTRDdigitizer::SDigits2Digits()
1513 {
1514   //
1515   // Merges the input s-digits and converts them to normal digits
1516   //
1517
1518   if (!MergeSDigits()) return kFALSE;
1519
1520   return ConvertSDigits();
1521
1522 }
1523
1524 //_____________________________________________________________________________
1525 Bool_t AliTRDdigitizer::CheckDetector(Int_t plane, Int_t chamber, Int_t sector)
1526 {
1527   //
1528   // Checks whether a detector is enabled
1529   //
1530
1531   if (fSimpleSim) return kTRUE; 
1532
1533   if ((fTRD->GetSensChamber() >=       0) &&
1534       (fTRD->GetSensChamber() != chamber)) return kFALSE;
1535   if ((fTRD->GetSensPlane()   >=       0) &&
1536       (fTRD->GetSensPlane()   !=   plane)) return kFALSE;
1537   if ( fTRD->GetSensSector()  >=       0) {
1538     Int_t sens1 = fTRD->GetSensSector();
1539     Int_t sens2 = sens1 + fTRD->GetSensSectorRange();
1540     sens2 -= ((Int_t) (sens2 / AliTRDgeometry::Nsect())) 
1541            * AliTRDgeometry::Nsect();
1542     if (sens1 < sens2) {
1543       if ((sector < sens1) || (sector >= sens2)) return kFALSE;
1544     }
1545     else {
1546       if ((sector < sens1) && (sector >= sens2)) return kFALSE;
1547     }
1548   }
1549
1550   return kTRUE;
1551
1552 }
1553
1554 //_____________________________________________________________________________
1555 Bool_t AliTRDdigitizer::WriteDigits() const
1556 {
1557   //
1558   // Writes out the TRD-digits and the dictionaries
1559   //
1560
1561   // Store the digits and the dictionary in the tree
1562   return fDigitsManager->WriteDigits();
1563
1564 }
1565
1566 //_____________________________________________________________________________
1567 void AliTRDdigitizer::DeConvExp(Double_t *source, Double_t *target
1568                               , Int_t n, Int_t nexp) 
1569 {
1570   //
1571   // Does the deconvolution by the digital filter.
1572   //
1573   // Author:        Marcus Gutfleisch, KIP Heidelberg
1574   // Optimized for: New TRF from Venelin Angelov, simulated with CADENCE
1575   //                Pad-ground capacitance = 25 pF
1576   //                Pad-pad cross talk capacitance = 6 pF
1577   //                For 10 MHz digitization, corresponding to 20 time bins
1578   //                in the drift region
1579   //
1580
1581   Double_t rates[2];
1582   Double_t coefficients[2];
1583
1584   /* initialize (coefficient = alpha, rates = lambda) */
1585   
1586   if( nexp == 1 ) {
1587     rates[0] = 0.466998;
1588     /* no rescaling */
1589     coefficients[0] = 1.0;
1590   }
1591   if( nexp == 2 ) {
1592     rates[0] = 0.8988162;
1593     coefficients[0] = 0.11392069;
1594     rates[1] = 0.3745688;
1595     coefficients[1] = 0.8860793;
1596     /* no rescaling */
1597     Float_t sumc = coefficients[0]+coefficients[1];
1598     coefficients[0] /= sumc;
1599     coefficients[1] /= sumc;
1600   }
1601       
1602   Int_t i, k;
1603   Double_t reminder[2];
1604   Double_t correction, result;
1605
1606   /* attention: computation order is important */
1607   correction=0.0;
1608   for ( k = 0; k < nexp; k++ ) reminder[k]=0.0;
1609     
1610   for ( i = 0; i < n; i++ ) {
1611     result = ( source[i] - correction );    /* no rescaling */
1612     target[i] = result;
1613     
1614     for ( k = 0; k < nexp; k++ ) reminder[k] = rates[k] 
1615                              * ( reminder[k] + coefficients[k] * result);
1616       
1617     correction=0.0;
1618     for ( k = 0; k < nexp; k++ ) correction += reminder[k];
1619   }
1620   
1621 }
1622
1623 //_____________________________________________________________________________
1624 void AliTRDdigitizer::InitOutput(TFile *file, Int_t iEvent)
1625 {
1626   //
1627   // Initializes the output branches
1628   //
1629
1630   fEvent = iEvent;
1631   fDigitsManager->MakeTreeAndBranches(file,iEvent);
1632
1633 }