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