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