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