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