]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDdigitizer.cxx
Store raw pad signals in AliTRDcluster (suggested by Xianguo)
[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 //                                                                        //
22 //  Authors: C. Blume (blume@ikf.uni-frankfurt.de)                        //
23 //           C. Lippmann                                                  //
24 //           B. Vulpescu                                                  //
25 //                                                                        //
26 //  The following effects are included:                                   //
27 //      - Diffusion                                                       //
28 //      - ExB effects                                                     //
29 //      - Gas gain including fluctuations                                 //
30 //      - Pad-response (simple Gaussian approximation)                    //
31 //      - Time-response                                                   //
32 //      - Electronics noise                                               //
33 //      - Electronics gain                                                //
34 //      - Digitization                                                    //
35 //      - Zero suppression                                                //
36 //                                                                        //
37 ////////////////////////////////////////////////////////////////////////////
38
39 #include <TGeoManager.h>
40 #include <TList.h>
41 #include <TMath.h>
42 #include <TRandom.h>
43 #include <TTree.h>
44
45 #include "AliRun.h"
46 #include "AliMC.h"
47 #include "AliRunLoader.h"
48 #include "AliLoader.h"
49 #include "AliConfig.h"
50 #include "AliDigitizationInput.h"
51 #include "AliRunLoader.h"
52 #include "AliLoader.h"
53 #include "AliLog.h"
54
55 #include "AliTRD.h"
56 #include "AliTRDhit.h"
57 #include "AliTRDdigitizer.h"
58 #include "AliTRDarrayDictionary.h"
59 #include "AliTRDarrayADC.h"
60 #include "AliTRDarraySignal.h"
61 #include "AliTRDdigitsManager.h"
62 #include "AliTRDgeometry.h"
63 #include "AliTRDpadPlane.h"
64 #include "AliTRDcalibDB.h"
65 #include "AliTRDSimParam.h"
66 #include "AliTRDCommonParam.h"
67 #include "AliTRDfeeParam.h"
68 #include "AliTRDmcmSim.h"
69 #include "AliTRDdigitsParam.h"
70
71 #include "Cal/AliTRDCalROC.h"
72 #include "Cal/AliTRDCalDet.h"
73 #include "Cal/AliTRDCalOnlineGainTableROC.h"
74
75 ClassImp(AliTRDdigitizer)
76
77 //_____________________________________________________________________________
78 AliTRDdigitizer::AliTRDdigitizer()
79   :AliDigitizer()
80   ,fRunLoader(0)
81   ,fDigitsManager(0)
82   ,fSDigitsManager(0)
83   ,fSDigitsManagerList(0)
84   ,fTRD(0)
85   ,fGeo(0)
86   ,fMcmSim(new AliTRDmcmSim)
87   ,fEvent(0)
88   ,fMasks(0)
89   ,fCompress(kTRUE)
90   ,fSDigits(kFALSE)
91   ,fMergeSignalOnly(kFALSE)
92 {
93   //
94   // AliTRDdigitizer default constructor
95   //
96   
97 }
98
99 //_____________________________________________________________________________
100 AliTRDdigitizer::AliTRDdigitizer(const Text_t *name, const Text_t *title)              
101   :AliDigitizer(name,title)
102   ,fRunLoader(0)
103   ,fDigitsManager(0)
104   ,fSDigitsManager(0)
105   ,fSDigitsManagerList(0)
106   ,fTRD(0)
107   ,fGeo(0)
108   ,fMcmSim(new AliTRDmcmSim)
109   ,fEvent(0)
110   ,fMasks(0)
111   ,fCompress(kTRUE)
112   ,fSDigits(kFALSE)
113   ,fMergeSignalOnly(kFALSE)
114 {
115   //
116   // AliTRDdigitizer constructor
117   //
118
119 }
120
121 //_____________________________________________________________________________
122 AliTRDdigitizer::AliTRDdigitizer(AliDigitizationInput* digInput
123                                , const Text_t *name, const Text_t *title)
124   :AliDigitizer(digInput,name,title)
125   ,fRunLoader(0)
126   ,fDigitsManager(0)
127   ,fSDigitsManager(0)
128   ,fSDigitsManagerList(0)
129   ,fTRD(0)
130   ,fGeo(0)
131   ,fMcmSim(new AliTRDmcmSim)
132   ,fEvent(0)
133   ,fMasks(0)
134   ,fCompress(kTRUE)
135   ,fSDigits(kFALSE)
136   ,fMergeSignalOnly(kFALSE)
137 {
138   //
139   // AliTRDdigitizer constructor
140   //
141
142 }
143
144 //_____________________________________________________________________________
145 AliTRDdigitizer::AliTRDdigitizer(AliDigitizationInput* digInput)
146   :AliDigitizer(digInput,"AliTRDdigitizer","TRD digitizer")
147   ,fRunLoader(0)
148   ,fDigitsManager(0)
149   ,fSDigitsManager(0)
150   ,fSDigitsManagerList(0)
151   ,fTRD(0)
152   ,fGeo(0)
153   ,fMcmSim(new AliTRDmcmSim)
154   ,fEvent(0)
155   ,fMasks(0)
156   ,fCompress(kTRUE)
157   ,fSDigits(kFALSE)
158   ,fMergeSignalOnly(kFALSE)
159 {
160   //
161   // AliTRDdigitizer constructor
162   //
163
164 }
165
166 //_____________________________________________________________________________
167 AliTRDdigitizer::AliTRDdigitizer(const AliTRDdigitizer &d)
168   :AliDigitizer(d)
169   ,fRunLoader(0)
170   ,fDigitsManager(0)
171   ,fSDigitsManager(0)
172   ,fSDigitsManagerList(0)
173   ,fTRD(0)
174   ,fGeo(0)
175   ,fMcmSim(new AliTRDmcmSim)
176   ,fEvent(0)
177   ,fMasks(0)
178   ,fCompress(d.fCompress)
179   ,fSDigits(d.fSDigits)
180   ,fMergeSignalOnly(d.fMergeSignalOnly)
181 {
182   //
183   // AliTRDdigitizer copy constructor
184   //
185
186 }
187
188 //_____________________________________________________________________________
189 AliTRDdigitizer::~AliTRDdigitizer()
190 {
191   //
192   // AliTRDdigitizer destructor
193   //
194
195   delete fDigitsManager;
196   fDigitsManager      = 0;
197
198   // s-digitsmanager will be deleted via list
199   fSDigitsManager     = 0;
200   if (fSDigitsManagerList) {
201     fSDigitsManagerList->Delete();
202     delete fSDigitsManagerList;
203   }
204   fSDigitsManagerList = 0;
205
206   delete [] fMasks;
207   fMasks       = 0;
208
209   delete fMcmSim;
210   fMcmSim = 0;
211
212   delete fGeo;
213   fGeo = 0;
214
215 }
216
217 //_____________________________________________________________________________
218 AliTRDdigitizer &AliTRDdigitizer::operator=(const AliTRDdigitizer &d)
219 {
220   //
221   // Assignment operator
222   //
223
224   if (this != &d) {
225     ((AliTRDdigitizer &) d).Copy(*this);
226   }
227
228   return *this;
229
230 }
231
232 //_____________________________________________________________________________
233 void AliTRDdigitizer::Copy(TObject &d) const
234 {
235   //
236   // Copy function
237   //
238
239   ((AliTRDdigitizer &) d).fRunLoader          = 0;
240   ((AliTRDdigitizer &) d).fDigitsManager      = 0;
241   ((AliTRDdigitizer &) d).fSDigitsManager     = 0;
242   ((AliTRDdigitizer &) d).fSDigitsManagerList = 0;
243   ((AliTRDdigitizer &) d).fTRD                = 0;
244   ((AliTRDdigitizer &) d).fGeo                = 0;
245   ((AliTRDdigitizer &) d).fEvent              = 0;
246   ((AliTRDdigitizer &) d).fMasks              = 0;
247   ((AliTRDdigitizer &) d).fCompress           = fCompress;
248   ((AliTRDdigitizer &) d).fSDigits            = fSDigits;
249   ((AliTRDdigitizer &) d).fMergeSignalOnly    = fMergeSignalOnly;
250
251 }
252
253 //_____________________________________________________________________________
254 void AliTRDdigitizer::Digitize(const Option_t* option)
255 {
256   //
257   // Executes the merging
258   //
259
260   Int_t iInput;
261
262   AliTRDdigitsManager *sdigitsManager;
263
264   TString optionString = option;
265   if (optionString.Contains("deb")) {
266     AliLog::SetClassDebugLevel("AliTRDdigitizer",1);
267     AliInfo("Called with debug option");
268   }
269
270   // The AliRoot file is already connected by the manager
271   AliRunLoader *inrl = 0x0;
272   
273   if (gAlice) {
274     AliDebug(1,"AliRun object found on file.");
275   }
276   else {
277     inrl = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(0));
278     inrl->LoadgAlice();
279     gAlice = inrl->GetAliRun();
280     if (!gAlice) {
281       AliError("Could not find AliRun object.");
282       return;
283     }
284   }
285                                                                            
286   Int_t nInput = fDigInput->GetNinputs();
287   fMasks       = new Int_t[nInput];
288   for (iInput = 0; iInput < nInput; iInput++) {
289     fMasks[iInput] = fDigInput->GetMask(iInput);
290   }
291
292   //
293   // Initialization
294   //
295
296   AliRunLoader *orl = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());
297
298   if (InitDetector()) {
299
300     AliLoader *ogime = orl->GetLoader("TRDLoader");
301
302     TTree *tree = 0;
303     if (fSDigits) { 
304       // If we produce SDigits
305       tree = ogime->TreeS();
306       if (!tree) {
307         ogime->MakeTree("S");
308         tree = ogime->TreeS();
309       }
310     }
311     else {
312       // If we produce Digits
313       tree = ogime->TreeD();
314       if (!tree) {
315         ogime->MakeTree("D");
316         tree = ogime->TreeD();
317       }
318     }
319
320     MakeBranch(tree);
321
322   }
323  
324   for (iInput = 0; iInput < nInput; iInput++) {
325
326     AliDebug(1,Form("Add input stream %d",iInput));
327
328     // Check if the input tree exists
329     inrl = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(iInput));
330     AliLoader *gime = inrl->GetLoader("TRDLoader");
331
332     TTree *treees = gime->TreeS();
333     if (treees == 0x0) {
334       if (gime->LoadSDigits()) {
335         AliError(Form("Error Occured while loading S. Digits for input %d.",iInput));
336         return;
337       }
338       treees = gime->TreeS();
339     }
340     
341     if (treees == 0x0) {
342       AliError(Form("Input stream %d does not exist",iInput));
343       return;
344     } 
345
346     // Read the s-digits via digits manager
347     sdigitsManager = new AliTRDdigitsManager();
348     sdigitsManager->SetSDigits(kTRUE);
349     
350     AliRunLoader *rl = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(iInput));
351     AliLoader *gimme = rl->GetLoader("TRDLoader");
352     if (!gimme->TreeS()) 
353       {
354         gimme->LoadSDigits();
355       }
356
357     sdigitsManager->ReadDigits(gimme->TreeS());
358    
359     // Add the s-digits to the input list 
360     AddSDigitsManager(sdigitsManager);
361
362   }
363
364   // Convert the s-digits to normal digits
365   AliDebug(1,"Do the conversion");
366   SDigits2Digits();
367
368   // Store the digits
369   AliDebug(1,"Write the digits");
370   fDigitsManager->WriteDigits();
371
372   // Write parameters
373   orl->CdGAFile();
374
375   // Clean up
376   DeleteSDigitsManager();
377
378   AliDebug(1,"Done");
379
380 }
381
382 //_____________________________________________________________________________
383 Bool_t AliTRDdigitizer::Open(const Char_t *file, Int_t nEvent)
384 {
385   //
386   // Opens a ROOT-file with TRD-hits and reads in the hit-tree
387   //
388   // Connect the AliRoot file containing Geometry, Kine, and Hits
389   //  
390
391   TString evfoldname = AliConfig::GetDefaultEventFolderName();
392
393   fRunLoader = AliRunLoader::GetRunLoader(evfoldname);
394   if (!fRunLoader) {
395     fRunLoader = AliRunLoader::Open(file,evfoldname,"UPDATE");
396   }  
397   if (!fRunLoader) {
398     AliError(Form("Can not open session for file %s.",file));
399     return kFALSE;
400   }
401    
402   if (!fRunLoader->GetAliRun()) {
403     fRunLoader->LoadgAlice();
404   }
405   gAlice = fRunLoader->GetAliRun();
406   
407   if (gAlice) {
408     AliDebug(1,"AliRun object found on file.");
409   }
410   else {
411     AliError("Could not find AliRun object.");
412     return kFALSE;
413   }
414
415   fEvent = nEvent;
416
417   AliLoader *loader = fRunLoader->GetLoader("TRDLoader");
418   if (!loader) {
419     AliError("Can not get TRD loader from Run Loader");
420     return kFALSE;
421   }
422   
423   if (InitDetector()) {
424     TTree *tree = 0;
425     if (fSDigits) { 
426       // If we produce SDigits
427       tree = loader->TreeS();
428       if (!tree) {
429         loader->MakeTree("S");
430         tree = loader->TreeS();
431       }
432     }
433     else {
434       // If we produce Digits
435       tree = loader->TreeD();
436       if (!tree) {
437         loader->MakeTree("D");
438         tree = loader->TreeD();
439       }
440     }
441     return MakeBranch(tree);
442   }
443   else {
444     return kFALSE;
445   }
446
447 }
448
449 //_____________________________________________________________________________
450 Bool_t AliTRDdigitizer::Open(AliRunLoader * const runLoader, Int_t nEvent)
451 {
452   //
453   // Opens a ROOT-file with TRD-hits and reads in the hit-tree
454   //
455   // Connect the AliRoot file containing Geometry, Kine, and Hits
456   //  
457
458   fRunLoader = runLoader;
459   if (!fRunLoader) {
460     AliError("RunLoader does not exist");
461     return kFALSE;
462   }
463    
464   if (!fRunLoader->GetAliRun()) {
465     fRunLoader->LoadgAlice();
466   }
467   gAlice = fRunLoader->GetAliRun();
468   
469   if (gAlice) {
470     AliDebug(1,"AliRun object found on file.");
471   }
472   else {
473     AliError("Could not find AliRun object.");
474     return kFALSE;
475   }
476
477   fEvent = nEvent;
478
479   AliLoader *loader = fRunLoader->GetLoader("TRDLoader");
480   if (!loader) {
481     AliError("Can not get TRD loader from Run Loader");
482     return kFALSE;
483   }
484   
485   if (InitDetector()) {
486     TTree *tree = 0;
487     if (fSDigits) { 
488       // If we produce SDigits
489       tree = loader->TreeS();
490       if (!tree) {
491         loader->MakeTree("S");
492         tree = loader->TreeS();
493       }
494     }
495     else {
496       // If we produce Digits
497       tree = loader->TreeD();
498       if (!tree) {
499         loader->MakeTree("D");
500         tree = loader->TreeD();
501       }
502     }
503     return MakeBranch(tree);
504   }
505   else {
506     return kFALSE;
507   }
508
509 }
510
511 //_____________________________________________________________________________
512 Bool_t AliTRDdigitizer::InitDetector()
513 {
514   //
515   // Sets the pointer to the TRD detector and the geometry
516   //
517
518   // Get the pointer to the detector class and check for version 1
519   fTRD = (AliTRD *) gAlice->GetDetector("TRD");
520   if (!fTRD) {
521     AliFatal("No TRD module found");
522     exit(1);
523   }
524   if (fTRD->IsVersion() != 1) {
525     AliFatal("TRD must be version 1 (slow simulator)");
526     exit(1);
527   }
528
529   // Get the geometry
530   fGeo = new AliTRDgeometry();
531
532   // Create a digits manager
533   if (fDigitsManager) {
534     delete fDigitsManager;
535   }
536   fDigitsManager = new AliTRDdigitsManager();
537   fDigitsManager->SetSDigits(fSDigits);
538   fDigitsManager->CreateArrays();
539   fDigitsManager->SetEvent(fEvent);
540
541   // The list for the input s-digits manager to be merged
542   if (fSDigitsManagerList) {
543     fSDigitsManagerList->Delete();
544   } 
545   else {
546     fSDigitsManagerList = new TList();
547   }
548
549   return kTRUE;
550
551 }
552
553 //_____________________________________________________________________________
554 Bool_t AliTRDdigitizer::MakeBranch(TTree *tree) const
555 {
556   // 
557   // Create the branches for the digits array
558   //
559
560   return fDigitsManager->MakeBranch(tree);
561
562 }
563
564 //_____________________________________________________________________________
565 void AliTRDdigitizer::AddSDigitsManager(AliTRDdigitsManager *man)
566 {
567   //
568   // Add a digits manager for s-digits to the input list.
569   //
570
571   fSDigitsManagerList->Add(man);
572
573 }
574
575 //_____________________________________________________________________________
576 void AliTRDdigitizer::DeleteSDigitsManager()
577 {
578   //
579   // Removes digits manager from the input list.
580   //
581
582   fSDigitsManagerList->Delete();
583
584 }
585
586 //_____________________________________________________________________________
587 Bool_t AliTRDdigitizer::MakeDigits()
588 {
589   //
590   // Creates digits.
591   //
592
593   AliDebug(1,"Start creating digits");
594
595   if (!fGeo) {
596     AliError("No geometry defined");
597     return kFALSE;
598   }
599
600   AliTRDcalibDB *calibration = AliTRDcalibDB::Instance();
601   if (!calibration) {
602     AliFatal("Could not get calibration object");
603     return kFALSE;
604   }
605
606   const Int_t kNdet  = AliTRDgeometry::Ndet();
607
608   Float_t **hits = new Float_t*[kNdet];
609   Int_t    *nhit = new Int_t[kNdet];
610   memset(nhit,0,kNdet*sizeof(Int_t));
611
612   AliTRDarraySignal *signals = 0x0;
613
614   // Check the number of time bins from simParam against OCDB,
615   // if OCDB value is not supposed to be used.
616   // As default, the value from OCDB is taken
617   if (AliTRDSimParam::Instance()->GetNTBoverwriteOCDB()) {
618     if (calibration->GetNumberOfTimeBinsDCS() != AliTRDSimParam::Instance()->GetNTimeBins()) {
619       AliWarning(Form("Number of time bins is different to OCDB value [SIM=%d, OCDB=%d]"
620                      ,AliTRDSimParam::Instance()->GetNTimeBins()
621                      ,calibration->GetNumberOfTimeBinsDCS()));
622     }
623     // Save the values for the raw data headers
624     fDigitsManager->GetDigitsParam()->SetNTimeBinsAll(AliTRDSimParam::Instance()->GetNTimeBins());
625   }
626   else {
627     // Save the values for the raw data headers
628     fDigitsManager->GetDigitsParam()->SetNTimeBinsAll(calibration->GetNumberOfTimeBinsDCS());
629   }
630
631   // Save the values for the raw data headers
632   fDigitsManager->GetDigitsParam()->SetADCbaselineAll(AliTRDSimParam::Instance()->GetADCbaseline());
633  
634   // Sort all hits according to detector number
635   if (!SortHits(hits,nhit)) {
636     AliError("Sorting hits failed");
637     delete [] hits;
638     delete [] nhit;
639     return kFALSE;
640   }
641
642   // Loop through all detectors
643   for (Int_t det = 0; det < kNdet; det++) {
644
645     // Detectors that are switched off, not installed, etc.
646     if ((!calibration->IsChamberNoData(det))    &&
647         ( fGeo->ChamberInGeometry(det))         &&
648         (nhit[det] > 0)) {
649
650       signals = new AliTRDarraySignal();
651           
652       // Convert the hits of the current detector to detector signals
653       if (!ConvertHits(det,hits[det],nhit[det],signals)) {
654         AliError(Form("Conversion of hits failed for detector=%d",det));
655         delete [] hits;
656         delete [] nhit;
657         delete signals;
658         signals = 0x0;
659         return kFALSE;
660       }
661
662       // Convert the detector signals to digits or s-digits
663       if (!ConvertSignals(det,signals)) {
664         AliError(Form("Conversion of signals failed for detector=%d",det));
665         delete [] hits;
666         delete [] nhit;
667         delete signals;
668         signals = 0x0;
669         return kFALSE;
670       }
671
672       // Delete the signals array
673       delete signals;
674       signals = 0x0;
675
676     } // if: detector status
677
678     delete [] hits[det];
679
680   } // for: detector
681
682   if (!fSDigits) {
683     if (AliDataLoader *trklLoader 
684           = AliRunLoader::Instance()->GetLoader("TRDLoader")->GetDataLoader("tracklets")) {
685       if (trklLoader->Tree())
686         trklLoader->WriteData("OVERWRITE");
687     }
688   }
689     
690   delete [] hits;
691   delete [] nhit;
692
693   return kTRUE;
694
695 }
696
697 //_____________________________________________________________________________
698 Bool_t AliTRDdigitizer::SortHits(Float_t **hits, Int_t *nhit)
699 {
700   //
701   // Read all the hits and sorts them according to detector number
702   // in the output array <hits>.
703   //
704
705   AliDebug(1,"Start sorting hits");
706
707   const Int_t kNdet = AliTRDgeometry::Ndet();
708   // Size of the hit vector
709   const Int_t kNhit = 6;
710
711   Float_t *xyz      = 0;
712   Int_t    nhitTrk  = 0;
713
714   Int_t   *lhit     = new Int_t[kNdet];
715   memset(lhit,0,kNdet*sizeof(Int_t));
716
717   for (Int_t det = 0; det < kNdet; det++) {
718     hits[det] = 0x0;
719   }
720
721   AliLoader *gimme = fRunLoader->GetLoader("TRDLoader");
722   if (!gimme->TreeH()) {
723     gimme->LoadHits();
724   }
725   TTree *hitTree = gimme->TreeH();
726   if (hitTree == 0x0) {
727     AliError("Can not get TreeH");
728     delete [] lhit;
729     return kFALSE;
730   }
731   fTRD->SetTreeAddress();
732
733   // Get the number of entries in the hit tree
734   // (Number of primary particles creating a hit somewhere)
735   Int_t nTrk = (Int_t) hitTree->GetEntries();
736   AliDebug(1,Form("Found %d tracks",nTrk));
737
738   // Loop through all the tracks in the tree
739   for (Int_t iTrk = 0; iTrk < nTrk; iTrk++) {
740
741     gAlice->GetMCApp()->ResetHits();
742     hitTree->GetEvent(iTrk);
743
744     if (!fTRD->Hits()) {
745       AliError(Form("No hits array for track = %d",iTrk));
746       continue;
747     }
748
749     // Number of hits for this track
750     nhitTrk = fTRD->Hits()->GetEntriesFast();
751
752     Int_t hitCnt = 0;
753     // Loop through the TRD hits
754     AliTRDhit *hit = (AliTRDhit *) fTRD->FirstHit(-1);
755     while (hit) {
756
757       hitCnt++;
758
759       // Don't analyze test hits
760       if (((Int_t) hit->GetCharge()) != 0) {
761
762         Int_t   trk  = hit->Track();
763         Int_t   det  = hit->GetDetector();
764         Int_t   q    = hit->GetCharge();
765         Float_t x    = hit->X();
766         Float_t y    = hit->Y();
767         Float_t z    = hit->Z();
768         Float_t time = hit->GetTime();
769
770         if (nhit[det] == lhit[det]) {
771           // Inititialization of new detector
772           xyz        = new Float_t[kNhit*(nhitTrk+lhit[det])];
773           if (hits[det]) {
774             memcpy(xyz,hits[det],sizeof(Float_t)*kNhit*lhit[det]);
775             delete [] hits[det];
776           }
777           lhit[det] += nhitTrk;
778           hits[det]  = xyz;
779         }
780         else {
781           xyz        = hits[det];
782         }
783         xyz[nhit[det]*kNhit+0] = x;
784         xyz[nhit[det]*kNhit+1] = y;
785         xyz[nhit[det]*kNhit+2] = z;
786         xyz[nhit[det]*kNhit+3] = q;
787         xyz[nhit[det]*kNhit+4] = trk;
788         xyz[nhit[det]*kNhit+5] = time;
789         nhit[det]++;
790
791       } // if: charge != 0
792
793       hit = (AliTRDhit *) fTRD->NextHit();   
794
795     } // for: hits of one track
796
797   } // for: tracks
798
799   delete [] lhit;
800
801   return kTRUE;
802
803 }
804
805 //_____________________________________________________________________________
806 Bool_t AliTRDdigitizer::ConvertHits(Int_t det
807                                   , const Float_t * const hits
808                                   , Int_t nhit
809                                   , AliTRDarraySignal *signals)
810 {
811   //
812   // Converts the detectorwise sorted hits to detector signals
813   //
814
815   AliDebug(1,Form("Start converting hits for detector=%d (nhits=%d)",det,nhit));
816
817   // Number of pads included in the pad response
818   const Int_t kNpad      = 3;
819   // Number of track dictionary arrays
820   const Int_t kNdict     = AliTRDdigitsManager::kNDict;
821   // Size of the hit vector
822   const Int_t kNhit      = 6;
823
824   // Width of the amplification region
825   const Float_t kAmWidth = AliTRDgeometry::AmThick();
826   // Width of the drift region
827   const Float_t kDrWidth = AliTRDgeometry::DrThick();
828   // Drift + amplification region 
829   const Float_t kDrMin   =          - 0.5 * kAmWidth;
830   const Float_t kDrMax   = kDrWidth + 0.5 * kAmWidth;
831   
832   Int_t    iPad          = 0;
833   Int_t    dict          = 0;
834   Int_t    timeBinTRFend = 1;
835
836   Double_t pos[3];
837   Double_t loc[3];
838   Double_t padSignal[kNpad];
839   Double_t signalOld[kNpad];
840
841   AliTRDarrayDictionary *dictionary[kNdict];
842
843   AliTRDSimParam    *simParam    = AliTRDSimParam::Instance();
844   AliTRDCommonParam *commonParam = AliTRDCommonParam::Instance();
845   AliTRDcalibDB     *calibration = AliTRDcalibDB::Instance();
846
847   if (!commonParam) {
848     AliFatal("Could not get common parameterss");
849     return kFALSE;
850   }
851   if (!simParam) {
852     AliFatal("Could not get simulation parameters");
853     return kFALSE;
854   }  
855   if (!calibration) {
856     AliFatal("Could not get calibration object");  
857     return kFALSE;
858   }
859
860   // Get the detector wise calibration objects
861   AliTRDCalROC       *calVdriftROC      = 0;
862   Float_t             calVdriftDetValue = 0.0;
863   const AliTRDCalDet *calVdriftDet      = calibration->GetVdriftDet();  
864   AliTRDCalROC       *calT0ROC          = 0;
865   Float_t             calT0DetValue     = 0.0;
866   const AliTRDCalDet *calT0Det          = calibration->GetT0Det();  
867   Double_t            calExBDetValue    = 0.0;
868   const AliTRDCalDet *calExBDet         = calibration->GetExBDet();
869
870   if (simParam->TRFOn()) {
871     timeBinTRFend = ((Int_t) (simParam->GetTRFhi() 
872                   * commonParam->GetSamplingFrequency())) - 1;
873   }
874
875   Int_t   nTimeTotal   = fDigitsManager->GetDigitsParam()->GetNTimeBins(det);
876   Float_t samplingRate = commonParam->GetSamplingFrequency();
877   Float_t elAttachProp = simParam->GetElAttachProp() / 100.0; 
878
879   AliTRDpadPlane *padPlane = fGeo->GetPadPlane(det);
880   Int_t   layer   = fGeo->GetLayer(det);         //update
881   Float_t row0    = padPlane->GetRow0ROC();
882   Int_t   nRowMax = padPlane->GetNrows();
883   Int_t   nColMax = padPlane->GetNcols();
884
885   // Create a new array for the signals
886   signals->Allocate(nRowMax,nColMax,nTimeTotal);
887
888   // Create a new array for the dictionary
889   for (dict = 0; dict < kNdict; dict++) {       
890     dictionary[dict] = (AliTRDarrayDictionary *) fDigitsManager->GetDictionary(det,dict);
891     dictionary[dict]->Allocate(nRowMax,nColMax,nTimeTotal);
892   }      
893
894   // Loop through the hits in this detector
895   for (Int_t hit = 0; hit < nhit; hit++) {
896
897     pos[0]          = hits[hit*kNhit+0];
898     pos[1]          = hits[hit*kNhit+1];
899     pos[2]          = hits[hit*kNhit+2];
900     Float_t q       = hits[hit*kNhit+3];
901     Float_t hittime = hits[hit*kNhit+5];
902     Int_t   track   = ((Int_t) hits[hit*kNhit+4]);
903
904     Int_t   inDrift = 1;
905
906     // Find the current volume with the geo manager
907     gGeoManager->SetCurrentPoint(pos);
908     gGeoManager->FindNode();      
909     if (strstr(gGeoManager->GetPath(),"/UK")) {
910       inDrift = 0;
911     }
912
913     // Get the calibration objects
914     calVdriftROC      = calibration->GetVdriftROC(det);
915     calVdriftDetValue = calVdriftDet->GetValue(det);
916     calT0ROC          = calibration->GetT0ROC(det);
917     calT0DetValue     = calT0Det->GetValue(det);
918     calExBDetValue    = calExBDet->GetValue(det);
919
920     // Go to the local coordinate system:
921     // loc[0] - col  direction in amplification or driftvolume
922     // loc[1] - row  direction in amplification or driftvolume
923     // loc[2] - time direction in amplification or driftvolume
924     gGeoManager->MasterToLocal(pos,loc);
925     if (inDrift) {
926       // Relative to middle of amplification region
927       loc[2] = loc[2] - kDrWidth/2.0 - kAmWidth/2.0;
928     } 
929
930     // The driftlength [cm] (w/o diffusion yet !).
931     // It is negative if the hit is between pad plane and anode wires.
932     Double_t driftlength = -1.0 * loc[2];
933
934     // Stupid patch to take care of TR photons that are absorbed
935     // outside the chamber volume. A real fix would actually need
936     // a more clever implementation of the TR hit generation
937     if (q < 0.0) {
938       if ((loc[1] < padPlane->GetRowEndROC()) ||
939           (loc[1] > padPlane->GetRow0ROC())) {
940         continue;
941       }
942       if ((driftlength < kDrMin) ||
943           (driftlength > kDrMax)) {
944         continue;
945       }
946     }
947
948     // Get row and col of unsmeared electron to retrieve drift velocity
949     // The pad row (z-direction)
950     Int_t    rowE         = padPlane->GetPadRowNumberROC(loc[1]);
951     if (rowE < 0) {
952       continue;
953     }
954     Double_t rowOffset    = padPlane->GetPadRowOffsetROC(rowE,loc[1]);
955     // The pad column (rphi-direction)
956     Double_t offsetTilt   = padPlane->GetTiltOffset(rowOffset);
957     Int_t    colE         = padPlane->GetPadColNumber(loc[0]+offsetTilt);
958     if (colE < 0) {
959       continue;   
960     }
961     Double_t colOffset    = 0.0;
962
963     // Normalized drift length
964     Float_t  driftvelocity  = calVdriftDetValue * calVdriftROC->GetValue(colE,rowE);
965     Double_t absdriftlength = TMath::Abs(driftlength);
966     if (commonParam->ExBOn()) {
967       absdriftlength /= TMath::Sqrt(1.0 / (1.0 + calExBDetValue*calExBDetValue));
968     }
969
970     // Loop over all electrons of this hit
971     // TR photons produce hits with negative charge
972     Int_t nEl = ((Int_t) TMath::Abs(q));
973     for (Int_t iEl = 0; iEl < nEl; iEl++) {
974
975       // Now the real local coordinate system of the ROC
976       // column direction: locC
977       // row direction:    locR 
978       // time direction:   locT
979       // locR and locC are identical to the coordinates of the corresponding
980       // volumina of the drift or amplification region.
981       // locT is defined relative to the wire plane (i.e. middle of amplification
982       // region), meaning locT = 0, and is negative for hits coming from the
983       // drift region. 
984       Double_t locC = loc[0];
985       Double_t locR = loc[1];
986       Double_t locT = loc[2];
987
988       // Electron attachment
989       if (simParam->ElAttachOn()) {
990         if (gRandom->Rndm() < (absdriftlength * elAttachProp)) {
991           continue;
992         }
993       }
994           
995       // Apply the diffusion smearing
996       if (simParam->DiffusionOn()) {
997         if (!(Diffusion(driftvelocity,absdriftlength,calExBDetValue,locR,locC,locT))) {
998           continue;
999         }
1000       }
1001
1002       // Apply E x B effects (depends on drift direction)
1003       if (commonParam->ExBOn()) {
1004         locC = locC + calExBDetValue * driftlength;
1005       }
1006
1007       // The electron position after diffusion and ExB in pad coordinates.
1008       // The pad row (z-direction)
1009       rowE       = padPlane->GetPadRowNumberROC(locR);
1010       if (rowE < 0) continue;
1011       rowOffset  = padPlane->GetPadRowOffsetROC(rowE,locR);
1012
1013       // The pad column (rphi-direction)
1014       offsetTilt = padPlane->GetTiltOffset(rowOffset);
1015       colE       = padPlane->GetPadColNumber(locC+offsetTilt);
1016       if (colE < 0) continue;         
1017       colOffset  = padPlane->GetPadColOffset(colE,locC+offsetTilt);
1018           
1019       // Also re-retrieve drift velocity because col and row may have changed
1020       driftvelocity = calVdriftDetValue * calVdriftROC->GetValue(colE,rowE);
1021       Float_t t0    = calT0DetValue     + calT0ROC->GetValue(colE,rowE);
1022
1023       // Convert the position to drift time [mus], using either constant drift velocity or
1024       // time structure of drift cells (non-isochronity, GARFIELD calculation).
1025       // Also add absolute time of hits to take pile-up events into account properly
1026       Double_t drifttime;
1027       if (simParam->TimeStructOn()) {
1028         // Get z-position with respect to anode wire
1029         Double_t zz  =  row0 - locR + padPlane->GetAnodeWireOffset();
1030         zz -= ((Int_t)(2 * zz)) / 2.0;
1031         if (zz > 0.25) {
1032           zz  = 0.5 - zz;
1033         }
1034         // Use drift time map (GARFIELD)
1035         drifttime = commonParam->TimeStruct(driftvelocity,0.5*kAmWidth-1.0*locT,zz)
1036                   + hittime;
1037       } 
1038       else {
1039         // Use constant drift velocity
1040         drifttime = TMath::Abs(locT) / driftvelocity
1041                   + hittime;
1042       }
1043
1044       // Apply the gas gain including fluctuations
1045       Double_t ggRndm = 0.0;
1046       do {
1047         ggRndm = gRandom->Rndm();
1048       } while (ggRndm <= 0);
1049       Double_t signal = -(simParam->GetGasGain()) * TMath::Log(ggRndm);
1050
1051       // Apply the pad response 
1052       if (simParam->PRFOn()) {
1053         // The distance of the electron to the center of the pad 
1054         // in units of pad width
1055         Double_t dist = (colOffset - 0.5*padPlane->GetColSize(colE))
1056                       / padPlane->GetColSize(colE);
1057         // This is a fixed parametrization, i.e. not dependent on
1058         // calibration values !
1059         if (!(calibration->PadResponse(signal,dist,layer,padSignal))) continue;
1060       }
1061       else {
1062         padSignal[0] = 0.0;
1063         padSignal[1] = signal;
1064         padSignal[2] = 0.0;
1065       }
1066
1067       // The time bin (always positive), with t0 distortion
1068       Double_t timeBinIdeal = drifttime * samplingRate + t0;
1069       // Protection 
1070       if (TMath::Abs(timeBinIdeal) > 2*nTimeTotal) {
1071         timeBinIdeal = 2 * nTimeTotal;
1072       }
1073       Int_t    timeBinTruncated = ((Int_t) timeBinIdeal);
1074       // The distance of the position to the middle of the timebin
1075       Double_t timeOffset       = ((Float_t) timeBinTruncated 
1076                                 + 0.5 - timeBinIdeal) / samplingRate;
1077           
1078       // Sample the time response inside the drift region
1079       // + additional time bins before and after.
1080       // The sampling is done always in the middle of the time bin
1081       for (Int_t iTimeBin = TMath::Max(timeBinTruncated,0)
1082           ;iTimeBin < TMath::Min(timeBinTruncated+timeBinTRFend,nTimeTotal)
1083           ;iTimeBin++) {
1084
1085         // Apply the time response
1086         Double_t timeResponse = 1.0;
1087         Double_t crossTalk    = 0.0;
1088         Double_t time         = (iTimeBin - timeBinTruncated) / samplingRate + timeOffset;
1089
1090         if (simParam->TRFOn()) {
1091           timeResponse = simParam->TimeResponse(time);
1092         }
1093         if (simParam->CTOn()) {
1094           crossTalk    = simParam->CrossTalk(time);
1095         }
1096
1097         signalOld[0] = 0.0;
1098         signalOld[1] = 0.0;
1099         signalOld[2] = 0.0;
1100
1101         for (iPad = 0; iPad < kNpad; iPad++) {
1102
1103           Int_t colPos = colE + iPad - 1;
1104           if (colPos <        0) continue;
1105           if (colPos >= nColMax) break;
1106
1107           // Add the signals
1108           signalOld[iPad] = signals->GetData(rowE,colPos,iTimeBin);
1109
1110           if (colPos != colE) {
1111             // Cross talk added to non-central pads
1112             signalOld[iPad] += padSignal[iPad] 
1113                              * (timeResponse + crossTalk);
1114           } 
1115           else {
1116             // W/o cross talk at central pad
1117             signalOld[iPad] += padSignal[iPad] 
1118                              * timeResponse;
1119           }
1120
1121           signals->SetData(rowE,colPos,iTimeBin,signalOld[iPad]);
1122
1123           // Store the track index in the dictionary
1124           // Note: We store index+1 in order to allow the array to be compressed
1125           // Note2: Taking out the +1 in track
1126           if (signalOld[iPad] > 0.0) { 
1127             for (dict = 0; dict < kNdict; dict++) {
1128               Int_t oldTrack = dictionary[dict]->GetData(rowE,colPos,iTimeBin); 
1129               if (oldTrack == track) break;
1130               if (oldTrack ==  -1 ) {
1131                 dictionary[dict]->SetData(rowE,colPos,iTimeBin,track);
1132                 break;
1133               }
1134             }
1135           }
1136
1137         } // Loop: pads
1138
1139       } // Loop: time bins
1140
1141     } // Loop: electrons of a single hit
1142
1143   } // Loop: hits
1144
1145   AliDebug(2,Form("Finished analyzing %d hits",nhit));
1146
1147   return kTRUE;
1148
1149 }
1150
1151 //_____________________________________________________________________________
1152 Bool_t AliTRDdigitizer::ConvertSignals(Int_t det, AliTRDarraySignal *signals)
1153 {
1154   //
1155   // Convert signals to digits
1156   //
1157
1158   AliDebug(1,Form("Start converting the signals for detector %d",det));
1159
1160   if (fSDigits) {
1161     // Convert the signal array to s-digits
1162     if (!Signal2SDigits(det,signals)) {
1163       return kFALSE;
1164     }
1165   }
1166   else {
1167     // Convert the signal array to digits
1168     if (!Signal2ADC(det,signals)) {
1169       return kFALSE;
1170     }
1171     // Run digital processing for digits
1172     RunDigitalProcessing(det);
1173   }
1174
1175   // Compress the arrays
1176   CompressOutputArrays(det);   
1177
1178   return kTRUE;
1179
1180 }
1181
1182 //_____________________________________________________________________________
1183 Bool_t AliTRDdigitizer::Signal2ADC(Int_t det, AliTRDarraySignal *signals)
1184 {
1185   //
1186   // Converts the sampled electron signals to ADC values for a given chamber
1187   //
1188
1189   AliDebug(1,Form("Start converting signals to ADC values for detector=%d",det));
1190
1191   AliTRDcalibDB     *calibration = AliTRDcalibDB::Instance();
1192   if (!calibration) {
1193     AliFatal("Could not get calibration object");
1194     return kFALSE;
1195   }
1196
1197   AliTRDSimParam    *simParam    = AliTRDSimParam::Instance();
1198   if (!simParam) {
1199     AliFatal("Could not get simulation parameters");
1200     return kFALSE;
1201   }  
1202
1203   // Converts number of electrons to fC
1204   const Double_t kEl2fC = 1.602e-19 * 1.0e15; 
1205
1206   // Coupling factor
1207   Double_t coupling     = simParam->GetPadCoupling() 
1208                         * simParam->GetTimeCoupling();
1209   // Electronics conversion factor
1210   Double_t convert      = kEl2fC 
1211                         * simParam->GetChipGain();
1212   // ADC conversion factor
1213   Double_t adcConvert   = simParam->GetADCoutRange()
1214                         / simParam->GetADCinRange();
1215   // The electronics baseline in mV
1216   Double_t baseline     = simParam->GetADCbaseline() 
1217                         / adcConvert;
1218   // The electronics baseline in electrons
1219   Double_t baselineEl   = baseline
1220                         / convert;
1221
1222   Int_t row  = 0;
1223   Int_t col  = 0;
1224   Int_t time = 0;
1225
1226   Int_t nRowMax    = fGeo->GetPadPlane(det)->GetNrows();
1227   Int_t nColMax    = fGeo->GetPadPlane(det)->GetNcols();
1228   Int_t nTimeTotal = fDigitsManager->GetDigitsParam()->GetNTimeBins(det);
1229   if (fSDigitsManager->GetDigitsParam()->GetNTimeBins(det)) {
1230     nTimeTotal = fSDigitsManager->GetDigitsParam()->GetNTimeBins(det);
1231   }
1232   else {
1233     AliFatal("Could not get number of time bins");
1234     return kFALSE;
1235   }
1236
1237   // The gain factor calibration objects
1238   const AliTRDCalDet *calGainFactorDet      = calibration->GetGainFactorDet();  
1239   AliTRDCalROC       *calGainFactorROC      = 0x0;
1240   Float_t             calGainFactorDetValue = 0.0;
1241
1242   AliTRDarrayADC     *digits = 0x0;
1243
1244   if (!signals) {
1245     AliError(Form("Signals array for detector %d does not exist\n",det));
1246     return kFALSE;
1247   }
1248   if (signals->HasData()) {
1249     // Expand the container if neccessary
1250     signals->Expand();   
1251   }
1252   else {
1253     // Create missing containers
1254     signals->Allocate(nRowMax,nColMax,nTimeTotal);      
1255   }
1256
1257   // Get the container for the digits of this detector
1258   if (fDigitsManager->HasSDigits()) {
1259     AliError("Digits manager has s-digits");
1260     return kFALSE;
1261   }
1262
1263   digits = (AliTRDarrayADC *) fDigitsManager->GetDigits(det);
1264   // Allocate memory space for the digits buffer
1265   if (!digits->HasData()) {
1266     digits->Allocate(nRowMax,nColMax,nTimeTotal);
1267   }
1268
1269   // Get the calibration objects
1270   calGainFactorROC      = calibration->GetGainFactorROC(det);
1271   calGainFactorDetValue = calGainFactorDet->GetValue(det);
1272
1273   // Get the online gain factors
1274   //AliTRDCalOnlineGainTableROC *onlineGainFactorROC 
1275   //  = calibration->GetOnlineGainTableROC(det);
1276
1277   // Create the digits for this chamber
1278   for (row  = 0; row  <  nRowMax; row++ ) {
1279     for (col  = 0; col  <  nColMax; col++ ) {
1280
1281       // Check whether pad is masked
1282       // Bridged pads are not considered yet!!!
1283       if (calibration->IsPadMasked(det,col,row) || 
1284           calibration->IsPadNotConnected(det,col,row)) {
1285         continue;
1286       }
1287
1288       // The gain factors
1289       Float_t padgain = calGainFactorDetValue 
1290                       * calGainFactorROC->GetValue(col,row);
1291       if (padgain <= 0) {
1292         AliError(Form("Not a valid gain %f, %d %d %d",padgain,det,col,row));
1293       }
1294
1295       for (time = 0; time < nTimeTotal; time++) {
1296
1297         // Get the signal amplitude
1298         Float_t signalAmp = signals->GetData(row,col,time);
1299         // Pad and time coupling
1300         signalAmp *= coupling;
1301         // Gain factors
1302         signalAmp *= padgain;
1303
1304         // Add the noise, starting from minus ADC baseline in electrons
1305         signalAmp  = TMath::Max((Double_t) gRandom->Gaus(signalAmp,simParam->GetNoise())
1306                                ,-baselineEl);
1307
1308         // Convert to mV
1309         signalAmp *= convert;
1310         // Add ADC baseline in mV
1311         signalAmp += baseline;
1312
1313         // Convert to ADC counts. Set the overflow-bit fADCoutRange if the
1314         // signal is larger than fADCinRange
1315         Short_t adc  = 0;
1316         if (signalAmp >= simParam->GetADCinRange()) {
1317           adc = ((Short_t) simParam->GetADCoutRange());
1318         }
1319         else {
1320           adc = TMath::Nint(signalAmp * adcConvert);
1321         }
1322
1323         // Saving all digits
1324         digits->SetData(row,col,time,adc);
1325
1326       } // for: time
1327
1328     } // for: col
1329   } // for: row
1330
1331   return kTRUE;
1332
1333 }
1334
1335 //_____________________________________________________________________________
1336 Bool_t AliTRDdigitizer::Signal2SDigits(Int_t det, AliTRDarraySignal *signals)
1337 {
1338   //
1339   // Converts the sampled electron signals to s-digits
1340   //
1341
1342   AliDebug(1,Form("Start converting signals to s-digits for detector=%d",det));
1343
1344   AliTRDcalibDB     *calibration = AliTRDcalibDB::Instance();
1345   if (!calibration) {
1346     AliFatal("Could not get calibration object");
1347     return kFALSE;
1348   }
1349
1350   Int_t row  = 0;
1351   Int_t col  = 0;
1352   Int_t time = 0;
1353
1354   Int_t nRowMax    = fGeo->GetPadPlane(det)->GetNrows();
1355   Int_t nColMax    = fGeo->GetPadPlane(det)->GetNcols();
1356   Int_t nTimeTotal = fDigitsManager->GetDigitsParam()->GetNTimeBins(det);
1357
1358   // Get the container for the digits of this detector
1359   if (!fDigitsManager->HasSDigits()) {
1360     AliError("Digits manager has no s-digits");
1361     return kFALSE;
1362   }
1363
1364   AliTRDarraySignal *digits = (AliTRDarraySignal *) fDigitsManager->GetSDigits(det);
1365   // Allocate memory space for the digits buffer
1366   if (!digits->HasData()) {
1367     digits->Allocate(nRowMax,nColMax,nTimeTotal);
1368   }
1369
1370   // Create the sdigits for this chamber
1371   for (row  = 0; row  <  nRowMax; row++ ) {
1372     for (col  = 0; col  <  nColMax; col++ ) {
1373       for (time = 0; time < nTimeTotal; time++) {         
1374         digits->SetData(row,col,time,signals->GetData(row,col,time));
1375       } // for: time
1376     } // for: col
1377   } // for: row
1378   
1379   return kTRUE;
1380
1381 }
1382
1383 //_____________________________________________________________________________
1384 Bool_t AliTRDdigitizer::Digits2SDigits(AliTRDdigitsManager * const manDig
1385                                      , AliTRDdigitsManager * const manSDig)
1386 {
1387   //
1388   // Converts digits into s-digits. Needed for embedding into real data.
1389   //
1390
1391   AliDebug(1,"Start converting digits to s-digits");
1392
1393   if (!fGeo) {
1394     fGeo = new AliTRDgeometry();
1395   }
1396
1397   AliTRDcalibDB     *calibration = AliTRDcalibDB::Instance();
1398   if (!calibration) {
1399     AliFatal("Could not get calibration object");
1400     return kFALSE;
1401   }
1402
1403   AliTRDSimParam    *simParam    = AliTRDSimParam::Instance();
1404   if (!simParam) {
1405     AliFatal("Could not get simulation parameters");
1406     return kFALSE;
1407   }  
1408
1409   // Converts number of electrons to fC
1410   const Double_t kEl2fC = 1.602e-19 * 1.0e15; 
1411
1412   // Coupling factor
1413   Double_t coupling     = simParam->GetPadCoupling() 
1414                         * simParam->GetTimeCoupling();
1415   // Electronics conversion factor
1416   Double_t convert      = kEl2fC 
1417                         * simParam->GetChipGain();
1418   // ADC conversion factor
1419   Double_t adcConvert   = simParam->GetADCoutRange()
1420                         / simParam->GetADCinRange();
1421   // The electronics baseline in mV
1422   Double_t baseline     = simParam->GetADCbaseline() 
1423                         / adcConvert;
1424   // The electronics baseline in electrons
1425   //Double_t baselineEl   = baseline
1426   //                      / convert;
1427
1428   // The gainfactor calibration objects
1429   //const AliTRDCalDet *calGainFactorDet      = calibration->GetGainFactorDet();  
1430   //AliTRDCalROC       *calGainFactorROC      = 0;
1431   //Float_t             calGainFactorDetValue = 0.0;
1432
1433   Int_t row  = 0;
1434   Int_t col  = 0;
1435   Int_t time = 0;
1436
1437   for (Int_t det = 0; det < AliTRDgeometry::Ndet(); det++) {
1438
1439     Int_t nRowMax    = fGeo->GetPadPlane(det)->GetNrows();
1440     Int_t nColMax    = fGeo->GetPadPlane(det)->GetNcols();
1441     Int_t nTimeTotal = manDig->GetDigitsParam()->GetNTimeBins(det);
1442
1443     // Get the calibration objects
1444     //calGainFactorROC      = calibration->GetGainFactorROC(det);
1445     //calGainFactorDetValue = calGainFactorDet->GetValue(det);
1446
1447     // Get the digits
1448     AliTRDarrayADC *digits = (AliTRDarrayADC *) manDig->GetDigits(det);
1449
1450     if (!manSDig->HasSDigits()) {
1451       AliError("SDigits manager has no s-digits");
1452       return kFALSE;
1453     }
1454     // Get the s-digits
1455     AliTRDarraySignal     *sdigits = (AliTRDarraySignal *)     manSDig->GetSDigits(det);
1456     AliTRDarrayDictionary *tracks0 = (AliTRDarrayDictionary *) manSDig->GetDictionary(det,0);
1457     AliTRDarrayDictionary *tracks1 = (AliTRDarrayDictionary *) manSDig->GetDictionary(det,1);
1458     AliTRDarrayDictionary *tracks2 = (AliTRDarrayDictionary *) manSDig->GetDictionary(det,2);
1459     // Allocate memory space for the digits buffer
1460     sdigits->Allocate(nRowMax,nColMax,nTimeTotal);
1461     tracks0->Allocate(nRowMax,nColMax,nTimeTotal);
1462     tracks1->Allocate(nRowMax,nColMax,nTimeTotal);
1463     tracks2->Allocate(nRowMax,nColMax,nTimeTotal);
1464
1465     // Keep the digits param
1466     manSDig->GetDigitsParam()->SetNTimeBinsAll(manDig->GetDigitsParam()->GetNTimeBins(0));
1467     manSDig->GetDigitsParam()->SetADCbaselineAll(manDig->GetDigitsParam()->GetADCbaseline(0));
1468
1469     if (digits->HasData()) {
1470
1471       digits->Expand();
1472
1473       // Create the sdigits for this chamber
1474       for (row  = 0; row  <  nRowMax; row++ ) {
1475         for (col  = 0; col  <  nColMax; col++ ) {
1476
1477           // The gain factors
1478           //Float_t padgain = calGainFactorDetValue 
1479           //                * calGainFactorROC->GetValue(col,row);
1480
1481           for (time = 0; time < nTimeTotal; time++) {
1482
1483             Short_t  adcVal = digits->GetData(row,col,time);
1484             Double_t signal = (Double_t) adcVal;
1485             // ADC -> signal in mV
1486             signal /= adcConvert;
1487             // Subtract baseline in mV
1488             signal -= baseline;
1489             // Signal in mV -> signal in #electrons
1490             signal /= convert;
1491             // Gain factor
1492             //signal /= padgain; // Not needed for real data
1493             // Pad and time coupling
1494             signal /= coupling;
1495
1496             sdigits->SetData(row,col,time,signal);
1497             tracks0->SetData(row,col,time,0);
1498             tracks1->SetData(row,col,time,0);
1499             tracks2->SetData(row,col,time,0);
1500
1501           } // for: time
1502
1503         } // for: col
1504       } // for: row
1505   
1506     } // if: has data
1507
1508     sdigits->Compress(0);
1509     tracks0->Compress();
1510     tracks1->Compress();
1511     tracks2->Compress();
1512
1513     // No compress just remove
1514     manDig->RemoveDigits(det);
1515     manDig->RemoveDictionaries(det);      
1516
1517   } // for: det
1518
1519   return kTRUE;
1520
1521 }
1522
1523 //_____________________________________________________________________________
1524 Bool_t AliTRDdigitizer::SDigits2Digits()
1525 {
1526   //
1527   // Merges the input s-digits and converts them to normal digits
1528   //
1529
1530   if (!MergeSDigits()) {
1531     return kFALSE;
1532   }
1533
1534   return ConvertSDigits();
1535
1536 }
1537
1538 //_____________________________________________________________________________
1539 Bool_t AliTRDdigitizer::MergeSDigits()
1540 {
1541   //
1542   // Merges the input s-digits:
1543   //   - The amplitude of the different inputs are summed up.
1544   //   - Of the track IDs from the input dictionaries only one is
1545   //     kept for each input. This works for maximal 3 different merged inputs.
1546   //
1547
1548   // Number of track dictionary arrays
1549   const Int_t kNDict = AliTRDdigitsManager::kNDict;
1550
1551   AliTRDSimParam    *simParam    = AliTRDSimParam::Instance();
1552   if (!simParam) {
1553     AliFatal("Could not get simulation parameters");
1554     return kFALSE;
1555   }
1556   
1557   AliTRDcalibDB     *calibration = AliTRDcalibDB::Instance();
1558   if (!calibration) {
1559     AliFatal("Could not get calibration object");
1560     return kFALSE;
1561   }
1562   
1563   Int_t iDict = 0;
1564   Int_t jDict = 0;
1565
1566   AliTRDarraySignal     *digitsA;
1567   AliTRDarraySignal     *digitsB;
1568   AliTRDarrayDictionary *dictionaryA[kNDict];
1569   AliTRDarrayDictionary *dictionaryB[kNDict];
1570
1571   AliTRDdigitsManager   *mergeSDigitsManager = 0x0;
1572   // Get the first s-digits
1573   fSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->First();
1574   if (!fSDigitsManager) { 
1575     AliError("No SDigits manager");
1576     return kFALSE;
1577   }
1578
1579   // Loop through the other sets of s-digits
1580   mergeSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->After(fSDigitsManager);
1581
1582   if (mergeSDigitsManager) {
1583     AliDebug(1,Form("Merge %d input files.",fSDigitsManagerList->GetSize()));
1584   }
1585   else {
1586     AliDebug(1,"Only one input file.");
1587   }
1588   
1589   Int_t iMerge = 0;
1590
1591   while (mergeSDigitsManager) {
1592
1593     iMerge++;
1594
1595     // Loop through the detectors
1596     for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1597
1598       Int_t nTimeTotal = fSDigitsManager->GetDigitsParam()->GetNTimeBins(iDet);
1599       if (mergeSDigitsManager->GetDigitsParam()->GetNTimeBins(iDet) != nTimeTotal) {
1600         AliError(Form("Mismatch in the number of time bins [%d,%d] in detector %d"
1601                      ,nTimeTotal
1602                      ,mergeSDigitsManager->GetDigitsParam()->GetNTimeBins(iDet)
1603                      ,iDet));
1604         return kFALSE;
1605       }
1606
1607       Int_t nRowMax = fGeo->GetPadPlane(iDet)->GetNrows();
1608       Int_t nColMax = fGeo->GetPadPlane(iDet)->GetNcols();
1609           
1610       // Loop through the pixels of one detector and add the signals
1611       digitsA = (AliTRDarraySignal *) fSDigitsManager->GetSDigits(iDet);    
1612       digitsB = (AliTRDarraySignal *) mergeSDigitsManager->GetSDigits(iDet); 
1613       digitsA->Expand();  
1614       if (!digitsA->HasData()) continue;
1615       digitsB->Expand();    
1616       if (!digitsB->HasData()) continue;
1617           
1618       for (iDict = 0; iDict < kNDict; iDict++) {
1619         dictionaryA[iDict] = (AliTRDarrayDictionary *) fSDigitsManager->GetDictionary(iDet,iDict);
1620         dictionaryB[iDict] = (AliTRDarrayDictionary *) mergeSDigitsManager->GetDictionary(iDet,iDict);
1621         dictionaryA[iDict]->Expand();  
1622         dictionaryB[iDict]->Expand();
1623       }
1624
1625       // Merge only detectors that contain a signal
1626       Bool_t doMerge = kTRUE;
1627       if (fMergeSignalOnly) {
1628         if (digitsA->GetOverThreshold(0) == 0) {                             
1629           doMerge = kFALSE;
1630         }
1631       }
1632           
1633       if (doMerge) {
1634               
1635         AliDebug(1,Form("Merge detector %d of input no.%d",iDet,iMerge+1));
1636               
1637         for (Int_t iRow  = 0; iRow  <  nRowMax;   iRow++ ) {
1638           for (Int_t iCol  = 0; iCol  <  nColMax;   iCol++ ) {
1639             for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {
1640                 
1641               // Add the amplitudes of the summable digits 
1642               Float_t ampA = digitsA->GetData(iRow,iCol,iTime);
1643               Float_t ampB = digitsB->GetData(iRow,iCol,iTime);
1644               ampA += ampB;
1645               digitsA->SetData(iRow,iCol,iTime,ampA);
1646
1647               // Add the mask to the track id if defined.
1648               for (iDict = 0; iDict < kNDict; iDict++) {
1649                 Int_t trackB = dictionaryB[iDict]->GetData(iRow,iCol,iTime);
1650                 if ((fMasks) && (trackB > 0))  {
1651                   for (jDict = 0; jDict < kNDict; jDict++) { 
1652                     Int_t trackA = dictionaryA[iDict]->GetData(iRow,iCol,iTime); 
1653                     if (trackA == 0) {
1654                       trackA = trackB + fMasks[iMerge];
1655                       dictionaryA[iDict]->SetData(iRow,iCol,iTime,trackA);  
1656                     } // if:  track A == 0
1657                   } // for: jDict
1658                 } // if:  fMasks and trackB > 0
1659               } // for: iDict
1660
1661             } // for: iTime
1662           } // for: iCol
1663         } // for: iRow
1664
1665       } // if:  doMerge
1666
1667       mergeSDigitsManager->RemoveDigits(iDet);
1668       mergeSDigitsManager->RemoveDictionaries(iDet);
1669   
1670       if (fCompress) {
1671         digitsA->Compress(0); 
1672         for (iDict = 0; iDict < kNDict; iDict++) {                                     
1673           dictionaryA[iDict]->Compress();
1674         }
1675       }
1676       
1677     } // for: detectors    
1678       
1679     // The next set of s-digits
1680     mergeSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->After(mergeSDigitsManager);
1681     
1682   } // while: mergeDigitsManagers
1683   
1684   return kTRUE;
1685
1686 }
1687
1688 //_____________________________________________________________________________
1689 Bool_t AliTRDdigitizer::ConvertSDigits()
1690 {
1691   //
1692   // Converts s-digits to normal digits
1693   //
1694
1695   AliTRDarraySignal *digitsIn = 0x0;
1696
1697   if (!fSDigitsManager->HasSDigits()) {
1698     AliError("No s-digits in digits manager");
1699     return kFALSE;
1700   }
1701
1702   // Loop through the detectors
1703   for (Int_t det = 0; det < AliTRDgeometry::Ndet(); det++) {
1704
1705     // Get the merged s-digits (signals)
1706     digitsIn = (AliTRDarraySignal *) fSDigitsManager->GetSDigits(det);
1707     if (!digitsIn->HasData()) {
1708       AliDebug(2,Form("No digits for det=%d",det));
1709       continue;
1710     }
1711     
1712     // Convert the merged sdigits to digits
1713     if (!Signal2ADC(det,digitsIn)) {
1714       continue;
1715     }
1716
1717     // Copy the dictionary information to the output array
1718     if (!CopyDictionary(det)) {
1719       continue;
1720     }
1721
1722     // Delete 
1723     fSDigitsManager->RemoveDigits(det);
1724     fSDigitsManager->RemoveDictionaries(det);
1725
1726     // Run digital processing
1727     RunDigitalProcessing(det);
1728
1729     // Compress the arrays
1730     CompressOutputArrays(det);
1731
1732   } // for: detector numbers
1733
1734   if (AliDataLoader *trklLoader = AliRunLoader::Instance()->GetLoader("TRDLoader")->GetDataLoader("tracklets")) {
1735     if (trklLoader->Tree())
1736       trklLoader->WriteData("OVERWRITE");
1737   }
1738
1739   // Save the values for the raw data headers
1740   if (AliTRDSimParam::Instance()->GetNTBoverwriteOCDB()) {
1741     fDigitsManager->GetDigitsParam()->SetNTimeBinsAll(AliTRDSimParam::Instance()->GetNTimeBins());
1742   }
1743   else {
1744     fDigitsManager->GetDigitsParam()->SetNTimeBinsAll(AliTRDcalibDB::Instance()->GetNumberOfTimeBinsDCS());
1745   }
1746   fDigitsManager->GetDigitsParam()->SetADCbaselineAll(AliTRDSimParam::Instance()->GetADCbaseline());
1747
1748   return kTRUE;
1749
1750 }
1751
1752 //_____________________________________________________________________________
1753 Bool_t AliTRDdigitizer::CopyDictionary(Int_t det)
1754 {
1755   //
1756   // Copies the dictionary information from the s-digits arrays
1757   // to the output arrays
1758   //
1759
1760   AliTRDcalibDB     *calibration = AliTRDcalibDB::Instance();
1761   if (!calibration) {
1762     AliFatal("Could not get calibration object");
1763     return kFALSE;
1764   }
1765
1766   AliDebug(1,Form("Start copying dictionaries for detector=%d",det));
1767
1768   const Int_t kNDict = AliTRDdigitsManager::kNDict;
1769   AliTRDarrayDictionary *dictionaryIn[kNDict];
1770   AliTRDarrayDictionary *dictionaryOut[kNDict];
1771
1772   Int_t nRowMax    = fGeo->GetPadPlane(det)->GetNrows();
1773   Int_t nColMax    = fGeo->GetPadPlane(det)->GetNcols();
1774   Int_t nTimeTotal = fSDigitsManager->GetDigitsParam()->GetNTimeBins(det);
1775
1776   Int_t row  = 0;
1777   Int_t col  = 0;
1778   Int_t time = 0;
1779   Int_t dict = 0;
1780
1781   for (dict = 0; dict < kNDict; dict++) {
1782
1783     dictionaryIn[dict]  = (AliTRDarrayDictionary *) fSDigitsManager->GetDictionary(det,dict);
1784     dictionaryIn[dict]->Expand();
1785     dictionaryOut[dict] = (AliTRDarrayDictionary *) fDigitsManager->GetDictionary(det,dict);
1786     dictionaryOut[dict]->Allocate(nRowMax,nColMax,nTimeTotal);
1787
1788     for (row = 0; row < nRowMax; row++) {
1789       for (col = 0; col < nColMax; col++) {
1790         for (time = 0; time < nTimeTotal; time++) {
1791           Int_t track = dictionaryIn[dict]->GetData(row,col,time);
1792           dictionaryOut[dict]->SetData(row,col,time,track);
1793         } // for: time
1794       } // for: col
1795     } // for: row
1796     
1797   } // for: dictionaries
1798   
1799   return kTRUE;
1800
1801 }
1802
1803 //_____________________________________________________________________________
1804 void AliTRDdigitizer::CompressOutputArrays(Int_t det)
1805 {
1806   //
1807   // Compress the output arrays
1808   //
1809
1810   const Int_t kNDict = AliTRDdigitsManager::kNDict;
1811   AliTRDarrayDictionary *dictionary = 0x0;
1812
1813   if (fCompress) {
1814
1815     if (!fSDigits) {
1816       AliTRDarrayADC *digits = 0x0;  
1817       digits = (AliTRDarrayADC *) fDigitsManager->GetDigits(det);
1818       digits->Compress();
1819     }
1820
1821     if (fSDigits) {
1822       AliTRDarraySignal *digits = 0x0; 
1823       digits = (AliTRDarraySignal *) fDigitsManager->GetSDigits(det);
1824       digits->Compress(0);
1825     }
1826
1827     for (Int_t dict = 0; dict < kNDict; dict++) {
1828       dictionary = (AliTRDarrayDictionary *) fDigitsManager->GetDictionary(det,dict);
1829       dictionary->Compress();
1830     }
1831
1832   }
1833
1834 }
1835
1836 //_____________________________________________________________________________
1837 Bool_t AliTRDdigitizer::WriteDigits() const
1838 {
1839   //
1840   // Writes out the TRD-digits and the dictionaries
1841   //
1842
1843   // Write parameters
1844   fRunLoader->CdGAFile();
1845
1846   // Store the digits and the dictionary in the tree
1847   return fDigitsManager->WriteDigits();
1848
1849 }
1850
1851 //_____________________________________________________________________________
1852 void AliTRDdigitizer::InitOutput(Int_t iEvent)
1853 {
1854   //
1855   // Initializes the output branches
1856   //
1857
1858   fEvent = iEvent;
1859    
1860   if (!fRunLoader) {
1861     AliError("Run Loader is NULL");
1862     return;  
1863   }
1864
1865   AliLoader *loader = fRunLoader->GetLoader("TRDLoader");
1866   if (!loader) {
1867     AliError("Can not get TRD loader from Run Loader");
1868     return;
1869   }
1870
1871   TTree *tree = 0;
1872   
1873   if (fSDigits) { 
1874     // If we produce SDigits
1875     tree = loader->TreeS();
1876     if (!tree) {
1877       loader->MakeTree("S");
1878       tree = loader->TreeS();
1879     }
1880   }
1881   else {
1882     // If we produce Digits
1883     tree = loader->TreeD();
1884     if (!tree) {
1885       loader->MakeTree("D");
1886       tree = loader->TreeD();
1887     }
1888   }
1889   fDigitsManager->SetEvent(iEvent);
1890   fDigitsManager->MakeBranch(tree);
1891
1892 }
1893   
1894 //_____________________________________________________________________________
1895 Int_t AliTRDdigitizer::Diffusion(Float_t vdrift, Double_t absdriftlength
1896                                , Double_t exbvalue
1897                                , Double_t &lRow, Double_t &lCol, Double_t &lTime)
1898 {
1899   //
1900   // Applies the diffusion smearing to the position of a single electron.
1901   // Depends on absolute drift length.
1902   //
1903   
1904   Float_t diffL = 0.0;
1905   Float_t diffT = 0.0;
1906
1907   if (AliTRDCommonParam::Instance()->GetDiffCoeff(diffL,diffT,vdrift)) {
1908
1909     Float_t driftSqrt = TMath::Sqrt(absdriftlength);
1910     Float_t sigmaT    = driftSqrt * diffT;
1911     Float_t sigmaL    = driftSqrt * diffL;
1912     lRow  = gRandom->Gaus(lRow ,sigmaT);
1913     if (AliTRDCommonParam::Instance()->ExBOn()) {
1914       lCol  = gRandom->Gaus(lCol ,sigmaT * 1.0 / (1.0 + exbvalue*exbvalue));
1915       lTime = gRandom->Gaus(lTime,sigmaL * 1.0 / (1.0 + exbvalue*exbvalue));
1916     }
1917     else {
1918       lCol  = gRandom->Gaus(lCol ,sigmaT);
1919       lTime = gRandom->Gaus(lTime,sigmaL);
1920     }
1921
1922     return 1;
1923
1924   }
1925   else {
1926
1927     return 0;
1928
1929   }
1930
1931 }
1932   
1933 //_____________________________________________________________________________
1934 void AliTRDdigitizer::RunDigitalProcessing(Int_t det)
1935 {
1936   //
1937   // Run the digital processing in the TRAP
1938   //
1939
1940   AliTRDfeeParam *feeParam = AliTRDfeeParam::Instance();
1941
1942   AliTRDarrayADC *digits = fDigitsManager->GetDigits(det);
1943   if (!digits)
1944     return;
1945
1946   //Call the methods in the mcm class using the temporary array as input  
1947   // process the data in the same order as in hardware
1948   for (Int_t side = 0; side <= 1; side++) {
1949     for(Int_t rob = side; rob < digits->GetNrow() / 2; rob += 2) {
1950       for(Int_t mcm = 0; mcm < 16; mcm++) {
1951         fMcmSim->Init(det, rob, mcm); 
1952         fMcmSim->SetDataByPad(digits, fDigitsManager);
1953         fMcmSim->Filter();
1954         if (feeParam->GetTracklet()) {
1955           fMcmSim->Tracklet();
1956           fMcmSim->StoreTracklets();
1957         }
1958         fMcmSim->ZSMapping();
1959         fMcmSim->WriteData(digits);
1960       }
1961     }
1962   }
1963 }
1964