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