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