]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDdigitizer.cxx
Update of digits test macros
[u/mrichter/AliRoot.git] / TRD / AliTRDdigitizer.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /*
17 $Log$
18 Revision 1.1.4.9  2000/10/26 17:00:22  cblume
19 Fixed bug in CheckDetector()
20
21 Revision 1.1.4.8  2000/10/23 13:41:35  cblume
22 Added protection against Log(0) in the gas gain calulation
23
24 Revision 1.1.4.7  2000/10/17 02:27:34  cblume
25 Get rid of global constants
26
27 Revision 1.1.4.6  2000/10/16 01:16:53  cblume
28 Changed timebin 0 to be the one closest to the readout
29
30 Revision 1.1.4.5  2000/10/15 23:34:29  cblume
31 Faster version of the digitizer
32
33 Revision 1.1.4.4  2000/10/06 16:49:46  cblume
34 Made Getters const
35
36 Revision 1.1.4.3  2000/10/04 16:34:58  cblume
37 Replace include files by forward declarations
38
39 Revision 1.1.4.2  2000/09/22 14:41:10  cblume
40 Bug fix in PRF. Included time response. New structure
41
42 Revision 1.10  2000/10/05 07:27:53  cblume
43 Changes in the header-files by FCA
44
45 Revision 1.9  2000/10/02 21:28:19  fca
46 Removal of useless dependecies via forward declarations
47
48 Revision 1.8  2000/06/09 11:10:07  cblume
49 Compiler warnings and coding conventions, next round
50
51 Revision 1.7  2000/06/08 18:32:58  cblume
52 Make code compliant to coding conventions
53
54 Revision 1.6  2000/06/07 16:27:32  cblume
55 Try to remove compiler warnings on Sun and HP
56
57 Revision 1.5  2000/05/09 16:38:57  cblume
58 Removed PadResponse(). Merge problem
59
60 Revision 1.4  2000/05/08 15:53:45  cblume
61 Resolved merge conflict
62
63 Revision 1.3  2000/04/28 14:49:27  cblume
64 Only one declaration of iDict in MakeDigits()
65
66 Revision 1.1.4.1  2000/05/08 14:42:04  cblume
67 Introduced AliTRDdigitsManager
68
69 Revision 1.1  2000/02/28 19:00:13  cblume
70 Add new TRD classes
71
72 */
73
74 ///////////////////////////////////////////////////////////////////////////////
75 //                                                                           //
76 //  Creates and handles digits from TRD hits                                 //
77 //                                                                           //
78 //  The following effects are included:                                      //
79 //      - Diffusion                                                          //
80 //      - ExB effects                                                        //
81 //      - Gas gain including fluctuations                                    //
82 //      - Pad-response (simple Gaussian approximation)                       //
83 //      - Electronics noise                                                  //
84 //      - Electronics gain                                                   //
85 //      - Digitization                                                       //
86 //      - ADC threshold                                                      //
87 //  The corresponding parameter can be adjusted via the various              //
88 //  Set-functions. If these parameters are not explicitly set, default       //
89 //  values are used (see Init-function).                                     //
90 //  To produce digits from a root-file with TRD-hits use the                 //
91 //  slowDigitsCreate.C macro.                                                //
92 //                                                                           //
93 ///////////////////////////////////////////////////////////////////////////////
94
95 #include <stdlib.h>
96
97 #include <TMath.h>
98 #include <TVector.h>
99 #include <TRandom.h>
100 #include <TROOT.h>
101 #include <TTree.h>
102 #include <TFile.h>
103 #include <TF1.h>
104
105 #include "AliRun.h"
106
107 #include "AliTRD.h"
108 #include "AliTRDhit.h"
109 #include "AliTRDdigitizer.h"
110 #include "AliTRDdataArrayI.h"
111 #include "AliTRDdataArrayF.h"
112 #include "AliTRDsegmentArray.h"
113 #include "AliTRDdigitsManager.h"
114 #include "AliTRDgeometry.h"
115
116 ClassImp(AliTRDdigitizer)
117
118 //_____________________________________________________________________________
119 AliTRDdigitizer::AliTRDdigitizer():TNamed()
120 {
121   //
122   // AliTRDdigitizer default constructor
123   //
124
125   fInputFile     = NULL;
126   fDigits        = NULL;
127   fTRD           = NULL;
128   fGeo           = NULL;
129   fPRF           = NULL;
130   fTRF           = NULL;
131   fTRFint        = NULL;
132
133   fEvent         = 0;
134   fGasGain       = 0.0;
135   fNoise         = 0.0;
136   fChipGain      = 0.0;
137   fADCoutRange   = 0.0;
138   fADCinRange    = 0.0;
139   fADCthreshold  = 0;
140   fDiffusionOn   = 0;
141   fDiffusionT    = 0.0;
142   fDiffusionL    = 0.0;
143   fElAttachOn    = 0;
144   fElAttachProp  = 0.0;
145   fExBOn         = 0;
146   fOmegaTau      = 0.0;
147   fPRFOn         = 0;
148   fTRFOn         = 0;
149   fDriftVelocity = 0.0;
150
151   fTRFbin        = 0;
152   fTRFlo         = 0.0;
153   fTRFhi         = 0.0;
154   fTRFwid        = 0.0;
155   fCompress      = kFALSE;
156   fVerbose       = 1;
157
158 }
159
160 //_____________________________________________________________________________
161 AliTRDdigitizer::AliTRDdigitizer(const Text_t *name, const Text_t *title)
162                 :TNamed(name,title)
163 {
164   //
165   // AliTRDdigitizer default constructor
166   //
167
168   fInputFile     = NULL;
169   fDigits        = NULL;
170   fTRD           = NULL;
171   fGeo           = NULL;
172
173   fEvent         = 0;
174
175   fCompress      = kTRUE;
176   fVerbose       = 1;
177
178   Init();
179
180 }
181
182 //_____________________________________________________________________________
183 AliTRDdigitizer::AliTRDdigitizer(const AliTRDdigitizer &d)
184 {
185   //
186   // AliTRDdigitizer copy constructor
187   //
188
189   ((AliTRDdigitizer &) d).Copy(*this);
190
191 }
192
193 //_____________________________________________________________________________
194 AliTRDdigitizer::~AliTRDdigitizer()
195 {
196   //
197   // AliTRDdigitizer destructor
198   //
199
200   if (fInputFile) {
201     fInputFile->Close();
202     delete fInputFile;
203   }
204
205   if (fDigits) {
206     delete fDigits;
207   }
208
209   if (fPRF) delete fPRF;
210   if (fTRF) delete fTRF;
211
212 }
213
214 //_____________________________________________________________________________
215 AliTRDdigitizer &AliTRDdigitizer::operator=(const AliTRDdigitizer &d)
216 {
217   //
218   // Assignment operator
219   //
220
221   if (this != &d) ((AliTRDdigitizer &) d).Copy(*this);
222   return *this;
223
224 }
225
226 //_____________________________________________________________________________
227 void AliTRDdigitizer::Copy(TObject &d)
228 {
229   //
230   // Copy function
231   //
232
233   ((AliTRDdigitizer &) d).fInputFile     = NULL;
234   ((AliTRDdigitizer &) d).fDigits        = NULL;
235   ((AliTRDdigitizer &) d).fTRD           = NULL;
236   ((AliTRDdigitizer &) d).fGeo           = NULL;
237
238   ((AliTRDdigitizer &) d).fEvent         = 0;
239
240   ((AliTRDdigitizer &) d).fGasGain       = fGasGain;
241   ((AliTRDdigitizer &) d).fNoise         = fNoise;
242   ((AliTRDdigitizer &) d).fChipGain      = fChipGain;
243   ((AliTRDdigitizer &) d).fADCoutRange   = fADCoutRange;
244   ((AliTRDdigitizer &) d).fADCinRange    = fADCinRange;
245   ((AliTRDdigitizer &) d).fADCthreshold  = fADCthreshold;
246   ((AliTRDdigitizer &) d).fDiffusionOn   = fDiffusionOn; 
247   ((AliTRDdigitizer &) d).fDiffusionT    = fDiffusionT;
248   ((AliTRDdigitizer &) d).fDiffusionL    = fDiffusionL;
249   ((AliTRDdigitizer &) d).fElAttachOn    = fElAttachOn;
250   ((AliTRDdigitizer &) d).fElAttachProp  = fElAttachProp;
251   ((AliTRDdigitizer &) d).fExBOn         = fExBOn;
252   ((AliTRDdigitizer &) d).fOmegaTau      = fOmegaTau;
253   ((AliTRDdigitizer &) d).fLorentzFactor = fLorentzFactor;
254   ((AliTRDdigitizer &) d).fPRFOn         = fPRFOn;
255   ((AliTRDdigitizer &) d).fTRFOn         = fTRFOn;
256
257   ((AliTRDdigitizer &) d).fCompress      = fCompress;
258   ((AliTRDdigitizer &) d).fVerbose       = fVerbose;
259
260   fPRF->Copy(*((AliTRDdigitizer &) d).fPRF);
261   fTRF->Copy(*((AliTRDdigitizer &) d).fTRF);
262
263   ((AliTRDdigitizer &) d).fTRFbin        = fTRFbin;
264   ((AliTRDdigitizer &) d).fTRFlo         = fTRFlo;
265   ((AliTRDdigitizer &) d).fTRFhi         = fTRFhi;
266   ((AliTRDdigitizer &) d).fTRFwid        = fTRFwid;
267   if (((AliTRDdigitizer &) d).fTRFint) delete ((AliTRDdigitizer &) d).fTRFint;
268   ((AliTRDdigitizer &) d).fTRFint = new Float_t[fTRFbin];
269   for (Int_t iBin = 0; iBin < fTRFbin; iBin++) {
270     ((AliTRDdigitizer &) d).fTRFint[iBin] = fTRFint[iBin];
271   }                                                                             
272
273 }
274
275 //_____________________________________________________________________________
276 Int_t AliTRDdigitizer::Diffusion(Float_t driftlength, Float_t *xyz)
277 {
278   //
279   // Applies the diffusion smearing to the position of a single electron
280   //
281
282   Float_t driftSqrt = TMath::Sqrt(driftlength);
283   Float_t sigmaT = driftSqrt * fDiffusionT;
284   Float_t sigmaL = driftSqrt * fDiffusionL;
285   xyz[0] = gRandom->Gaus(xyz[0], sigmaL * fLorentzFactor);
286   xyz[1] = gRandom->Gaus(xyz[1], sigmaT * fLorentzFactor);
287   xyz[2] = gRandom->Gaus(xyz[2], sigmaT);
288
289   return 1;
290
291 }
292
293 //_____________________________________________________________________________
294 Int_t AliTRDdigitizer::ExB(Float_t driftlength, Float_t *xyz)
295 {
296   //
297   // Applies E x B effects to the position of a single electron
298   //
299
300   xyz[0] = xyz[0];
301   xyz[1] = xyz[1] + fOmegaTau * driftlength;
302   xyz[2] = xyz[2];
303
304   return 1;
305
306 }
307
308 //_____________________________________________________________________________
309 Int_t AliTRDdigitizer::PadResponse(Float_t signal, Float_t dist, Float_t *pad)
310 {
311   //
312   // Applies the pad response
313   //
314
315   if (fPRF) {
316     pad[0] = TMath::Max(fPRF->Eval(-1.0 - dist,0,0) * signal,0.0);
317     pad[1] = TMath::Max(fPRF->Eval(     - dist,0,0) * signal,0.0);
318     pad[2] = TMath::Max(fPRF->Eval( 1.0 - dist,0,0) * signal,0.0);
319     return 1;
320   }
321   else {
322     return 0;
323   }
324
325 }
326
327 //_____________________________________________________________________________
328 Float_t AliTRDdigitizer::TimeResponse(Float_t time)
329 {
330   //
331   // Applies the preamp shaper time response
332   //
333
334   Int_t iBin = ((Int_t) ((time - fTRFlo) / fTRFwid)); 
335   if ((iBin >= 0) && (iBin < fTRFbin)) {
336     return fTRFint[iBin];
337   }
338   else {
339     return 0.0;
340   }    
341
342 }
343
344 //_____________________________________________________________________________
345 void AliTRDdigitizer::Init()
346 {
347   //
348   // Initializes the digitization procedure with standard values
349   //
350
351   // The default parameter for the digitization
352   fGasGain       = 8.0E3;
353   fNoise         = 3000.;
354   fChipGain      = 10.;
355   fADCoutRange   = 255.;
356   fADCinRange    = 2000.;
357   fADCthreshold  = 1;
358
359   // Transverse and longitudinal diffusion coefficients (Xe/Isobutane)
360   fDiffusionOn   = 1;
361   fDiffusionT    = 0.060;
362   fDiffusionL    = 0.017;
363
364   // Propability for electron attachment
365   fElAttachOn    = 0;
366   fElAttachProp  = 0.0;
367
368   // E x B effects
369   fExBOn         = 0;
370   // omega * tau. (tau ~ 12 * 10^-12, B = 0.2T)
371   fOmegaTau      = 17.6 * 12.0 * 0.2 * 0.01;
372
373   // The pad response function
374   fPRFOn         = 1;
375   fPRF           = new TF1("PRF","[0]*([1]+exp(-x*x/(2.0*[2])))",-3,3);
376   fPRF->SetParameter(0, 0.8872);
377   fPRF->SetParameter(1,-0.00573);
378   fPRF->SetParameter(2, 0.454 * 0.454);
379
380   // The drift velocity (cm / mus)
381   fDriftVelocity = 1.0;
382
383   // The time response function
384   fTRFOn         = 1;
385   Float_t loTRF  = -200.0;
386   Float_t hiTRF  = 1000.0;
387   fTRF           = new TF1("TRF",TRFlandau,loTRF,hiTRF,3);
388   fTRF->SetParameter(0,  1.0 / 24.24249);
389   fTRF->SetParameter(1,  0.0);
390   fTRF->SetParameter(2, 25.0);
391   fTRFbin        = 1200;
392   fTRFlo         = loTRF * fDriftVelocity / 1000.0;
393   fTRFhi         = hiTRF * fDriftVelocity / 1000.0;
394   fTRFwid        = (fTRFhi - fTRFlo) / ((Float_t) fTRFbin);
395
396 }
397
398 //_____________________________________________________________________________
399 void AliTRDdigitizer::IntegrateTRF()
400 {
401   //
402   // Integrates the time response function over the time bin size
403   //
404
405   if (fTRFint) delete fTRFint;
406   fTRFint = new Float_t[fTRFbin];
407   Float_t hiTRF   = fTRFhi                 / fDriftVelocity * 1000.0;
408   Float_t loTRF   = fTRFlo                 / fDriftVelocity * 1000.0;
409   Float_t timeBin = fGeo->GetTimeBinSize() / fDriftVelocity * 1000.0;
410   Float_t binWidth = (hiTRF - loTRF) / ((Float_t) fTRFbin);
411   for (Int_t iBin = 0; iBin < fTRFbin; iBin++) {
412     Float_t bin = iBin * binWidth + loTRF - 0.5 * timeBin;
413     fTRFint[iBin] = fTRF->Integral(bin,bin + timeBin);
414   }
415
416 }
417
418 //_____________________________________________________________________________
419 Bool_t AliTRDdigitizer::Open(const Char_t *name, Int_t nEvent)
420 {
421   //
422   // Opens a ROOT-file with TRD-hits and reads in the hit-tree
423   //
424
425   // Connect the AliRoot file containing Geometry, Kine, and Hits
426   fInputFile = (TFile*) gROOT->GetListOfFiles()->FindObject(name);
427   if (!fInputFile) {
428     printf("AliTRDdigitizer::Open -- ");
429     printf("Open the ALIROOT-file %s.\n",name);
430     fInputFile = new TFile(name,"UPDATE");
431   }
432   else {
433     printf("AliTRDdigitizer::Open -- ");
434     printf("%s is already open.\n",name);
435   }
436
437   gAlice = (AliRun*) fInputFile->Get("gAlice");
438   if (gAlice) {
439     printf("AliTRDdigitizer::Open -- ");
440     printf("AliRun object found on file.\n");
441   }
442   else {
443     printf("AliTRDdigitizer::Open -- ");
444     printf("Could not find AliRun object.\n");
445     return kFALSE;
446   }
447
448   fEvent = nEvent;
449
450   // Import the Trees for the event nEvent in the file
451   Int_t nparticles = gAlice->GetEvent(fEvent);
452   if (nparticles <= 0) {
453     printf("AliTRDdigitizer::Open -- ");
454     printf("No entries in the trees for event %d.\n",fEvent);
455     return kFALSE;
456   }
457
458   return InitDetector();
459
460 }
461
462 //_____________________________________________________________________________
463 Bool_t AliTRDdigitizer::InitDetector()
464 {
465   //
466   // Sets the pointer to the TRD detector and the geometry
467   //
468
469   // Get the pointer to the detector class and check for version 1
470   fTRD = (AliTRD*) gAlice->GetDetector("TRD");
471   if (fTRD->IsVersion() != 1) {
472     printf("AliTRDdigitizer::Open -- ");
473     printf("TRD must be version 1 (slow simulator).\n");
474     exit(1);
475   }
476
477   // Get the geometry
478   fGeo = fTRD->GetGeometry();
479   printf("AliTRDdigitizer::Open -- ");
480   printf("Geometry version %d\n",fGeo->IsVersion());
481
482   return kTRUE;
483
484 }
485
486 //_____________________________________________________________________________
487 Bool_t AliTRDdigitizer::MakeDigits()
488 {
489   //
490   // Loops through the TRD-hits and creates the digits.
491   //
492
493   ///////////////////////////////////////////////////////////////
494   // Parameter 
495   ///////////////////////////////////////////////////////////////
496
497   // Converts number of electrons to fC
498   const Float_t kEl2fC  = 1.602E-19 * 1.0E15; 
499
500   ///////////////////////////////////////////////////////////////
501
502   // Number of pads included in the pad response
503   const Int_t kNpad  = 3;
504
505   // Number of track dictionary arrays
506   const Int_t kNDict = AliTRDdigitsManager::NDict();
507
508   Int_t   iRow, iCol, iTime;
509   Int_t   iDict;
510   Int_t   nBytes = 0;
511
512   Int_t   totalSizeDigits = 0;
513   Int_t   totalSizeDict0  = 0;
514   Int_t   totalSizeDict1  = 0;
515   Int_t   totalSizeDict2  = 0;
516
517   AliTRDdataArrayF *signals = 0;
518   AliTRDdataArrayI *digits  = 0;
519   AliTRDdataArrayI *dictionary[kNDict];
520
521   // Create a digits manager
522   fDigits = new AliTRDdigitsManager();
523
524   // Create a container for the amplitudes
525   AliTRDsegmentArray *signalsArray 
526                      = new AliTRDsegmentArray("AliTRDdataArrayF",AliTRDgeometry::Ndet());
527
528   if (!fGeo) {
529     printf("AliTRDdigitizer::MakeDigits -- ");
530     printf("No geometry defined\n");
531     return kFALSE;
532   }
533
534   printf("AliTRDdigitizer::MakeDigits -- ");
535   printf("Start creating digits.\n");
536   if (fVerbose > 0) this->Dump();
537
538   // The Lorentz factor
539   if (fExBOn) {
540     fLorentzFactor = 1.0 / (1.0 + fOmegaTau*fOmegaTau);
541   }
542   else {
543     fLorentzFactor = 1.0;
544   }
545
546   // Create the integrated TRF
547   IntegrateTRF();
548
549   // Get the pointer to the hit tree
550   TTree *HitTree = gAlice->TreeH();
551
552   // Get the number of entries in the hit tree
553   // (Number of primary particles creating a hit somewhere)
554   Int_t nTrack = (Int_t) HitTree->GetEntries();
555   if (fVerbose > 0) {
556     printf("AliTRDdigitizer::MakeDigits -- ");
557     printf("Found %d primary particles\n",nTrack);
558   } 
559
560   Int_t detectorOld = -1;
561   Int_t countHits   =  0; 
562
563   // Loop through all entries in the tree
564   for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
565
566     gAlice->ResetHits();
567     nBytes += HitTree->GetEvent(iTrack);
568
569     // Get the number of hits in the TRD created by this particle
570     Int_t nHit = fTRD->Hits()->GetEntriesFast();
571     if (fVerbose > 0) {
572       printf("AliTRDdigitizer::MakeDigits -- ");
573       printf("Found %d hits for primary particle %d\n",nHit,iTrack);
574     }
575
576     // Loop through the TRD hits  
577     for (Int_t iHit = 0; iHit < nHit; iHit++) {
578
579       countHits++;
580
581       AliTRDhit *hit = (AliTRDhit *) fTRD->Hits()->UncheckedAt(iHit);
582       Float_t pos[3];
583               pos[0]   = hit->X();
584               pos[1]   = hit->Y();
585               pos[2]   = hit->Z();
586       Float_t q        = hit->GetCharge();
587       Int_t   track    = hit->Track();
588       Int_t   detector = hit->GetDetector();
589       Int_t   plane    = fGeo->GetPlane(detector);
590       Int_t   sector   = fGeo->GetSector(detector);
591       Int_t   chamber  = fGeo->GetChamber(detector);
592
593       if (!(CheckDetector(plane,chamber,sector))) continue;
594
595       Int_t   nRowMax     = fGeo->GetRowMax(plane,chamber,sector);
596       Int_t   nColMax     = fGeo->GetColMax(plane);
597       Int_t   nTimeMax    = fGeo->GetTimeMax();
598       Float_t row0        = fGeo->GetRow0(plane,chamber,sector);
599       Float_t col0        = fGeo->GetCol0(plane);
600       Float_t time0       = fGeo->GetTime0(plane);
601       Float_t rowPadSize  = fGeo->GetRowPadSize();
602       Float_t colPadSize  = fGeo->GetColPadSize();
603       Float_t timeBinSize = fGeo->GetTimeBinSize();
604
605       if (fVerbose > 1) {
606         printf("Analyze hit no. %d ",iHit);
607         printf("-----------------------------------------------------------\n");
608         hit->Dump();
609         printf("plane = %d, sector = %d, chamber = %d\n"
610               ,plane,sector,chamber);
611         printf("nRowMax = %d, nColMax = %d, nTimeMax = %d\n" 
612               ,nRowMax,nColMax,nTimeMax);
613         printf("row0 = %f, col0 = %f, time0 = %f\n"
614               ,row0,col0,time0);
615       }
616        
617       // Get different container if the detector has changed
618       if (detector != detectorOld) {
619         if (fVerbose > 1) {
620           printf("AliTRDdigitizer::MakeDigits -- ");
621           printf("Get new container. New det = %d, Old det = %d\n"
622                 ,detector,detectorOld);
623         }
624         // Compress the old one if enabled
625         if ((fCompress) && (detectorOld > -1)) {
626           if (fVerbose > 1) {
627             printf("AliTRDdigitizer::MakeDigits -- ");
628             printf("Compress the old container ... ");
629           }
630           signals->Compress(1,0);
631           for (iDict = 0; iDict < kNDict; iDict++) {
632             dictionary[iDict]->Compress(1,0);
633           }
634           if (fVerbose > 1) printf("done\n");
635         }
636         // Get the new container
637         signals = (AliTRDdataArrayF *) signalsArray->At(detector);
638         if (signals->GetNtime() == 0) {
639           // Allocate a new one if not yet existing
640           if (fVerbose > 1) {
641             printf("AliTRDdigitizer::MakeDigits -- ");
642             printf("Allocate a new container ... ");
643           }
644           signals->Allocate(nRowMax,nColMax,nTimeMax);
645         }
646         else {
647           // Expand an existing one
648           if (fVerbose > 1) {
649             printf("AliTRDdigitizer::MakeDigits -- ");
650             printf("Expand an existing container ... ");
651           }
652           if (fCompress) signals->Expand();
653         }
654         // The same for the dictionary
655         for (iDict = 0; iDict < kNDict; iDict++) {       
656           dictionary[iDict] = fDigits->GetDictionary(detector,iDict);
657           if (dictionary[iDict]->GetNtime() == 0) {
658             dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeMax);
659           }
660           else {
661             if (fCompress) dictionary[iDict]->Expand();
662           }
663         }      
664         if (fVerbose > 1) printf("done\n");
665         detectorOld = detector;
666       }
667
668       // Rotate the sectors on top of each other
669       Float_t rot[3];
670       fGeo->Rotate(detector,pos,rot);
671
672       // The driftlength
673       Float_t driftlength = time0 - rot[0];
674       if ((driftlength < 0) || 
675           (driftlength > AliTRDgeometry::DrThick())) break;
676       Float_t driftlengthL = driftlength;
677       if (fExBOn) driftlengthL /= TMath::Sqrt(fLorentzFactor);
678
679       // The hit position in pad coordinates (center pad)
680       // The pad row (z-direction)
681       Int_t  rowH = (Int_t) ((rot[2] -  row0) /  rowPadSize);
682       // The pad column (rphi-direction)  
683       Int_t  colH = (Int_t) ((rot[1] -  col0) /  colPadSize);
684       // The time bucket
685       Int_t timeH = (Int_t) (driftlength      / timeBinSize);
686       if (fVerbose > 1) {
687         printf("rowH = %d, colH = %d, timeH = %d\n"
688               ,rowH,colH,timeH);
689       }
690
691       // Loop over all electrons of this hit
692       // TR photons produce hits with negative charge
693       Int_t nEl = ((Int_t) TMath::Abs(q));
694       for (Int_t iEl = 0; iEl < nEl; iEl++) {
695
696         Float_t xyz[3];
697         xyz[0] = rot[0];
698         xyz[1] = rot[1];
699         xyz[2] = rot[2];
700
701         // Electron attachment
702         if (fElAttachOn) {
703           if (gRandom->Rndm() < (driftlengthL * fElAttachProp / 100.)) 
704             continue;
705         }
706
707         // Apply the diffusion smearing
708         if (fDiffusionOn) {
709           if (!(Diffusion(driftlengthL,xyz))) continue;
710         }
711
712         // Apply E x B effects
713         if (fExBOn) { 
714           if (!(ExB(driftlength,xyz))) continue;   
715         }
716
717         // The electron position 
718         // The pad row (z-direction)
719         Int_t  rowE = (Int_t) ((xyz[2] -  row0) /  rowPadSize);
720         // The pad column (rphi-direction)
721         Int_t  colE = (Int_t) ((xyz[1] -  col0) /  colPadSize);
722         // The time bucket
723         Int_t timeE = (Int_t) ((time0 - xyz[0]) / timeBinSize);
724
725         if (( rowE < 0) || ( rowE >=  nRowMax)) continue;
726         if (( colE < 0) || ( colE >=  nColMax)) continue;
727         if ((timeE < 0) || (timeE >= nTimeMax)) continue;
728
729         // Apply the gas gain including fluctuations
730         Float_t ggRndm = 0.0;
731         do {
732           ggRndm = gRandom->Rndm();
733         } while (ggRndm <= 0);
734         Int_t signal = (Int_t) (-fGasGain * TMath::Log(ggRndm));
735
736         if (fVerbose > 2) {
737           printf("  electron no. %d, signal = %d\n",iEl,signal);
738           printf("  rowE = %d, colE = %d, timeE = %d\n"
739                  ,rowE,colE,timeE);
740         }
741
742         // Apply the pad response 
743         Float_t padSignal[kNpad];
744         if (fPRFOn) {
745           // The distance of the electron to the center of the pad 
746           // in units of pad width
747           Float_t dist = (xyz[1] - col0 - (colE + 0.5) * colPadSize) 
748                        / colPadSize;
749           if (!(PadResponse(signal,dist,padSignal))) continue;
750         }
751         else {
752           padSignal[0] = 0.0;
753           padSignal[1] = signal;
754           padSignal[2] = 0.0;
755         }
756
757         // The distance of the position to the beginning of the timebin
758         Float_t timeOffset = (time0 - timeE * timeBinSize) - xyz[0];
759         Int_t timeTRDbeg = 0;
760         Int_t timeTRDend = 1;
761         if (fTRFOn) {
762           timeTRDbeg =  2;
763           timeTRDend = 11;
764         }
765         for (Int_t iTimeBin = TMath::Max(timeE - timeTRDbeg,       0)
766                  ; iTimeBin < TMath::Min(timeE + timeTRDend,nTimeMax) 
767                  ; iTimeBin++) {
768
769           // Apply the time response
770           Float_t timeResponse = 1.0;
771           if (fTRFOn) {
772             Float_t time = (iTimeBin - timeE) * timeBinSize + timeOffset;
773             timeResponse = TimeResponse(time);
774           }
775
776           // Add the signals
777           Float_t signalOld[kNpad] = { 0.0, 0.0, 0.0 };
778           for (Int_t iPad = 0; iPad < kNpad; iPad++) {
779             Int_t colPos = colE + iPad - 1;
780             if (colPos <        0) continue;
781             if (colPos >= nColMax) break;
782             signalOld[iPad]  = signals->GetData(rowE,colPos,iTimeBin);
783             signalOld[iPad] += padSignal[iPad] * timeResponse;
784             signals->SetData(rowE,colPos,iTimeBin,signalOld[iPad]);
785           }
786           if (fVerbose > 3) {
787             printf("    iTimeBin = %d, timeResponse = %f\n"
788                   ,iTimeBin,timeResponse);
789             printf("    pad-signal = %f, %f, %f\n"
790                    ,signalOld[0],signalOld[1],signalOld[2]);
791           }
792
793           // Store the track index in the dictionary
794           // Note: We store index+1 in order to allow the array to be compressed
795           for (iDict = 0; iDict < kNDict; iDict++) {
796             Int_t oldTrack = dictionary[iDict]->GetData(rowE,colE,timeE);
797             if (oldTrack == track+1) break;
798             //if (oldTrack ==      -1) break;
799             if (oldTrack ==       0) {
800               dictionary[iDict]->SetData(rowE,colE,timeE,track+1);
801               if (fVerbose > 3) {
802                 printf("    track index = %d\n",track); 
803               }
804               break;
805             }
806           }
807           if ((fVerbose > 1) && (iDict == kNDict)) {
808             printf("AliTRDdigitizer::MakeDigits -- ");
809             printf("More than three tracks for one digit!\n");
810           }
811
812         }
813
814       }
815
816     }
817
818   } // All hits finished
819
820   printf("AliTRDdigitizer::MakeDigits -- ");
821   printf("Finished analyzing %d hits\n",countHits);
822
823   // Loop through all chambers to finalize the digits
824   for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
825
826     Int_t plane    = fGeo->GetPlane(iDet);
827     Int_t sector   = fGeo->GetSector(iDet);
828     Int_t chamber  = fGeo->GetChamber(iDet);
829     Int_t nRowMax  = fGeo->GetRowMax(plane,chamber,sector);
830     Int_t nColMax  = fGeo->GetColMax(plane);
831     Int_t nTimeMax = fGeo->GetTimeMax();
832
833     if (!(CheckDetector(plane,chamber,sector))) continue;
834     if (fVerbose > 0) {
835       printf("AliTRDdigitizer::MakeDigits -- ");
836       printf("Digitization for chamber %d\n",iDet);
837     }
838
839     // Add a container for the digits of this detector
840     digits = fDigits->GetDigits(iDet);        
841     // Allocate memory space for the digits buffer
842     digits->Allocate(nRowMax,nColMax,nTimeMax);
843
844     // Get the signal container
845     signals = (AliTRDdataArrayF *) signalsArray->At(iDet);
846     if (signals->GetNtime() == 0) {
847       // Create missing containers
848       signals->Allocate(nRowMax,nColMax,nTimeMax);      
849     }
850     else {
851       // Expand the container if neccessary
852       if (fCompress) signals->Expand();
853     }
854     // Create the missing dictionary containers
855     for (iDict = 0; iDict < kNDict; iDict++) {       
856       dictionary[iDict] = fDigits->GetDictionary(iDet,iDict);
857       if (dictionary[iDict]->GetNtime() == 0) {
858         dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeMax);
859       }
860     }
861
862     Int_t nDigits = 0;
863
864     // Create the digits for this chamber
865     for (iRow  = 0; iRow  <  nRowMax; iRow++ ) {
866       for (iCol  = 0; iCol  <  nColMax; iCol++ ) {
867         for (iTime = 0; iTime < nTimeMax; iTime++) {         
868
869           Float_t signalAmp = signals->GetData(iRow,iCol,iTime);
870
871           // Add the noise
872           signalAmp  = TMath::Max((Float_t) gRandom->Gaus(signalAmp,fNoise),0.0);
873           // Convert to fC
874           signalAmp *= kEl2fC;
875           // Convert to mV
876           signalAmp *= fChipGain;
877           // Convert to ADC counts. Set the overflow-bit fADCoutRange if the 
878           // signal is larger than fADCinRange
879           Int_t adc  = 0;
880           if (signalAmp >= fADCinRange) {
881             adc = ((Int_t) fADCoutRange);
882           }
883           else {
884             adc = ((Int_t) (signalAmp * (fADCoutRange / fADCinRange)));
885           }
886           if (fVerbose > 2) {
887             printf("  iRow = %d, iCol = %d, iTime = %d\n"
888                   ,iRow,iCol,iTime);
889             printf("  signal = %f, adc = %d\n",signalAmp,adc);
890           }
891
892           // Store the amplitude of the digit if above threshold
893           if (adc > fADCthreshold) {
894             nDigits++;
895             digits->SetData(iRow,iCol,iTime,adc);
896           }
897
898         }
899       }
900     }
901
902     // Compress the arrays
903     digits->Compress(1,0);
904     for (iDict = 0; iDict < kNDict; iDict++) {
905       dictionary[iDict]->Compress(1,0);
906     }
907
908     totalSizeDigits += digits->GetSize();
909     totalSizeDict0  += dictionary[0]->GetSize();
910     totalSizeDict1  += dictionary[1]->GetSize();
911     totalSizeDict2  += dictionary[2]->GetSize();
912
913     printf("AliTRDdigitizer::MakeDigits -- ");
914     printf("Found %d digits in detector %d.\n",nDigits,iDet);
915  
916     if (fCompress) signals->Compress(1,0);
917
918   }
919
920   printf("AliTRDdigitizer::MakeDigits -- ");
921   printf("Total number of analyzed hits = %d\n",countHits);
922
923   printf("AliTRDdigitizer::MakeDigits -- ");
924   printf("Total digits data size = %d, %d, %d, %d\n",totalSizeDigits
925                                                     ,totalSizeDict0
926                                                     ,totalSizeDict1
927                                                     ,totalSizeDict2);        
928
929   return kTRUE;
930
931 }
932
933 //_____________________________________________________________________________
934 Bool_t AliTRDdigitizer::CheckDetector(Int_t plane, Int_t chamber, Int_t sector)
935 {
936   //
937   // Checks whether a detector is enabled
938   //
939
940   if ((fTRD->GetSensChamber() >=       0) &&
941       (fTRD->GetSensChamber() != chamber)) return kFALSE;
942   if ((fTRD->GetSensPlane()   >=       0) &&
943       (fTRD->GetSensPlane()   !=  sector)) return kFALSE;
944   if ( fTRD->GetSensSector()  >=       0) {
945     Int_t sens1 = fTRD->GetSensSector();
946     Int_t sens2 = sens1 + fTRD->GetSensSectorRange();
947     sens2 -= ((Int_t) (sens2 / AliTRDgeometry::Nsect())) 
948            * AliTRDgeometry::Nsect();
949     if (sens1 < sens2) {
950       if ((sector < sens1) || (sector >= sens2)) return kFALSE;
951     }
952     else {
953       if ((sector < sens1) && (sector >= sens2)) return kFALSE;
954     }
955   }
956
957   return kTRUE;
958
959 }
960
961 //_____________________________________________________________________________
962 Bool_t AliTRDdigitizer::WriteDigits()
963 {
964   //
965   // Writes out the TRD-digits and the dictionaries
966   //
967
968   // Create the branches
969   if (!(gAlice->TreeD()->GetBranch("TRDdigits"))) { 
970     if (!fDigits->MakeBranch()) return kFALSE;
971   }
972
973   // Store the digits and the dictionary in the tree
974   fDigits->WriteDigits();
975
976   // Write the new tree into the input file (use overwrite option)
977   Char_t treeName[7];
978   sprintf(treeName,"TreeD%d",fEvent);
979   printf("AliTRDdigitizer::WriteDigits -- ");
980   printf("Write the digits tree %s for event %d.\n"
981         ,treeName,fEvent);
982   gAlice->TreeD()->Write(treeName,2);
983  
984   return kTRUE;
985
986 }
987
988 //_____________________________________________________________________________
989 void AliTRDdigitizer::SetPRF(TF1 *prf)
990 {
991   //
992   // Defines a new pad response function
993   //
994
995   if (fPRF) delete fPRF;
996   fPRF = prf;     
997
998 }
999
1000 //_____________________________________________________________________________
1001 void AliTRDdigitizer::SetTRF(TF1 *trf)
1002 {
1003   //
1004   // Defines a new time response function
1005   //
1006
1007   if (fTRF) delete fTRF;
1008   fTRF = trf;     
1009
1010 }
1011
1012 //_____________________________________________________________________________
1013 Double_t TRFlandau(Double_t *x, Double_t *par)
1014 {
1015
1016   Double_t xx = x[0];
1017   Double_t landau = par[0] * TMath::Landau(xx,par[1],par[2]);
1018
1019   return landau;
1020
1021 }