]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDdigitizer.cxx
Replaced Fill3() by Fill()
[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  2000/02/28 19:00:13  cblume
19 Add new TRD classes
20
21 */
22
23 ///////////////////////////////////////////////////////////////////////////////
24 //                                                                           //
25 //  Creates and handles digits from TRD hits                                 //
26 //                                                                           //
27 //  The following effects are included:                                      //
28 //      - Diffusion                                                          //
29 //      - ExB effects                                                        //
30 //      - Gas gain including fluctuations                                    //
31 //      - Pad-response (simple Gaussian approximation)                       //
32 //      - Electronics noise                                                  //
33 //      - Electronics gain                                                   //
34 //      - Digitization                                                       //
35 //      - ADC threshold                                                      //
36 //  The corresponding parameter can be adjusted via the various              //
37 //  Set-functions. If these parameters are not explicitly set, default       //
38 //  values are used (see Init-function).                                     //
39 //  To produce digits from a root-file with TRD-hits use the                 //
40 //  slowDigitsCreate.C macro.                                                //
41 //                                                                           //
42 ///////////////////////////////////////////////////////////////////////////////
43
44 #include <TMath.h>
45 #include <TVector.h>
46 #include <TRandom.h>
47
48 #include "AliTRD.h"
49 #include "AliTRDdigitizer.h"
50 #include "AliTRDmatrix.h"
51
52 ClassImp(AliTRDdigitizer)
53
54 //_____________________________________________________________________________
55 AliTRDdigitizer::AliTRDdigitizer():TNamed()
56 {
57   //
58   // AliTRDdigitizer default constructor
59   //
60
61   fInputFile = NULL;
62   fEvent     = 0;
63
64 }
65
66 //_____________________________________________________________________________
67 AliTRDdigitizer::AliTRDdigitizer(const Text_t *name, const Text_t *title)
68                 :TNamed(name,title)
69 {
70   //
71   // AliTRDdigitizer default constructor
72   //
73
74   fInputFile   = NULL;
75   fEvent       = 0;
76
77   fDigitsArray   = new AliTRDsegmentArray(kNsect*kNplan*kNcham);
78   for (Int_t iDict = 0; iDict < kNDict; iDict++) {
79     fDictionary[iDict] = new AliTRDsegmentArray(kNsect*kNplan*kNcham);
80   }
81
82   Init();
83
84 }
85
86 //_____________________________________________________________________________
87 AliTRDdigitizer::~AliTRDdigitizer()
88 {
89
90   if (fInputFile) {
91     fInputFile->Close();
92     delete fInputFile;
93   }
94
95   if (fDigitsArray) {
96     fDigitsArray->Delete();
97     delete fDigitsArray;
98   }
99
100   for (Int_t iDict = 0; iDict < kNDict; iDict++) {
101     fDictionary[iDict]->Delete();
102     delete fDictionary[iDict];
103   }
104
105 }
106
107 //_____________________________________________________________________________
108 Int_t AliTRDdigitizer::Diffusion(Float_t driftlength, Float_t *xyz)
109 {
110   //
111   // Applies the diffusion smearing to the position of a single electron
112   //
113
114   Float_t driftSqrt = TMath::Sqrt(driftlength);
115   Float_t sigmaT = driftSqrt * fDiffusionT;
116   Float_t sigmaL = driftSqrt * fDiffusionL;
117   xyz[0] = gRandom->Gaus(xyz[0], sigmaL * fLorentzFactor);
118   xyz[1] = gRandom->Gaus(xyz[1], sigmaT * fLorentzFactor);
119   xyz[2] = gRandom->Gaus(xyz[2], sigmaT);
120   return 1;
121
122 }
123
124 //_____________________________________________________________________________
125 Int_t AliTRDdigitizer::ExB(Float_t driftlength, Float_t *xyz)
126 {
127   //
128   // Applies E x B effects to the position of a single electron
129   //
130
131   xyz[0] = xyz[0];
132   xyz[1] = xyz[1] + fLorentzAngle * driftlength;
133   xyz[2] = xyz[2];
134
135   return 1;
136
137 }
138
139 //_____________________________________________________________________________
140 void AliTRDdigitizer::Init()
141 {
142   //
143   // Initializes the digitization procedure with standard values
144   //
145
146   // The default parameter for the digitization
147   fGasGain       = 2.0E3;
148   fNoise         = 3000.;
149   fChipGain      = 10.;
150   fADCoutRange   = 255.;
151   fADCinRange    = 2000.;
152   fADCthreshold  = 1;
153
154   // Transverse and longitudinal diffusion coefficients (Xe/Isobutane)
155   fDiffusionOn   = 1;
156   fDiffusionT    = 0.060;
157   fDiffusionL    = 0.017;
158
159   // Propability for electron attachment
160   fElAttachOn    = 0;
161   fElAttachProp  = 0.0;
162
163   // E x B effects
164   fExBOn         = 0;
165   // omega * tau. (tau ~ 12 * 10^-12, B = 0.2T)
166   fLorentzAngle  = 17.6 * 12.0 * 0.2 * 0.01;
167
168 }
169
170 //_____________________________________________________________________________
171 Bool_t AliTRDdigitizer::Open(const Char_t *name, Int_t nEvent)
172 {
173   //
174   // Opens a ROOT-file with TRD-hits and reads in the hit-tree
175   //
176
177   // Connect the AliRoot file containing Geometry, Kine, and Hits
178   fInputFile = (TFile*) gROOT->GetListOfFiles()->FindObject(name);
179   if (!fInputFile) {
180     printf("AliTRDdigitizer::Open -- ");
181     printf("Open the ALIROOT-file %s.\n",name);
182     fInputFile = new TFile(name,"UPDATE");
183   }
184   else {
185     printf("AliTRDdigitizer::Open -- ");
186     printf("%s is already open.\n",name);
187   }
188
189   // Get AliRun object from file or create it if not on file
190   //if (!gAlice) {
191     gAlice = (AliRun*) fInputFile->Get("gAlice");
192     if (gAlice) {
193       printf("AliTRDdigitizer::Open -- ");
194       printf("AliRun object found on file.\n");
195     }
196     else {
197       printf("AliTRDdigitizer::Open -- ");
198       printf("Could not find AliRun object.\n");
199       return kFALSE;
200     }
201   //}
202
203   fEvent = nEvent;
204
205   // Import the Trees for the event nEvent in the file
206   Int_t nparticles = gAlice->GetEvent(fEvent);
207   if (nparticles <= 0) {
208     printf("AliTRDdigitizer::Open -- ");
209     printf("No entries in the trees for event %d.\n",fEvent);
210     return kFALSE;
211   }
212
213   return kTRUE;
214
215 }
216
217 //_____________________________________________________________________________
218 Float_t AliTRDdigitizer::PadResponse(Float_t x)
219 {
220   //
221   // The pad response for the chevron pads. 
222   // We use a simple Gaussian approximation which should be good
223   // enough for our purpose.
224   //
225
226   // The parameters for the response function
227   const Float_t aa  =  0.8872;
228   const Float_t bb  = -0.00573;
229   const Float_t cc  =  0.454;
230   const Float_t cc2 =  cc*cc;
231
232   Float_t pr = aa * (bb + TMath::Exp(-x*x / (2. * cc2)));
233
234   return (pr);
235
236 }
237
238 //_____________________________________________________________________________
239 Bool_t AliTRDdigitizer::MakeDigits()
240 {
241   //
242   // Loops through the TRD-hits and creates the digits.
243   //
244
245   // Get the pointer to the detector class and check for version 1
246   AliTRD *TRD = (AliTRD*) gAlice->GetDetector("TRD");
247   if (TRD->IsVersion() != 1) {
248     printf("AliTRDdigitizer::MakeDigits -- ");
249     printf("TRD must be version 1 (slow simulator).\n");
250     return kFALSE; 
251   }
252
253   // Get the geometry
254   AliTRDgeometry *Geo = TRD->GetGeometry();
255   printf("AliTRDdigitizer::MakeDigits -- ");
256   printf("Geometry version %d\n",Geo->IsVersion());
257
258   printf("AliTRDdigitizer::MakeDigits -- ");
259   printf("Start creating digits.\n");
260
261   ///////////////////////////////////////////////////////////////
262   // Parameter 
263   ///////////////////////////////////////////////////////////////
264
265   // Converts number of electrons to fC
266   const Float_t el2fC  = 1.602E-19 * 1.0E15; 
267
268   ///////////////////////////////////////////////////////////////
269
270   Int_t   iRow, iCol, iTime;
271   Int_t   nBytes = 0;
272
273   Int_t   totalSizeDigits = 0;
274   Int_t   totalSizeDict0  = 0;
275   Int_t   totalSizeDict1  = 0;
276   Int_t   totalSizeDict2  = 0;
277
278   AliTRDhit       *Hit;
279   AliTRDdataArray *Digits;
280   AliTRDdataArray *Dictionary[kNDict];
281
282   // Get the pointer to the hit tree
283   TTree *HitTree    = gAlice->TreeH();
284
285   // The Lorentz factor
286   if (fExBOn) {
287     fLorentzFactor = 1.0 / (1.0 + fLorentzAngle*fLorentzAngle);
288   }
289   else {
290     fLorentzFactor = 1.0;
291   }
292
293   // Get the number of entries in the hit tree
294   // (Number of primary particles creating a hit somewhere)
295   Int_t nTrack = (Int_t) HitTree->GetEntries();
296
297   Int_t chamBeg = 0;
298   Int_t chamEnd = kNcham;
299   if (TRD->GetSensChamber() >= 0) {
300     chamBeg = TRD->GetSensChamber();
301     chamEnd = chamEnd + 1;
302   }
303   Int_t planBeg = 0;
304   Int_t planEnd = kNplan;
305   if (TRD->GetSensPlane()   >= 0) {
306     planBeg = TRD->GetSensPlane();
307     planEnd = planBeg + 1;
308   }
309   Int_t sectBeg = 0;
310   Int_t sectEnd = kNsect;
311   if (TRD->GetSensSector()  >= 0) {
312     sectBeg = TRD->GetSensSector();
313     sectEnd = sectBeg + 1;
314   }
315
316   // Loop through all the chambers
317   for (Int_t iCham = chamBeg; iCham < chamEnd; iCham++) {
318     for (Int_t iPlan = planBeg; iPlan < planEnd; iPlan++) {
319       for (Int_t iSect = sectBeg; iSect < sectEnd; iSect++) {
320
321         Int_t nDigits = 0;
322
323         Int_t iDet = Geo->GetDetector(iPlan,iCham,iSect);
324
325         printf("AliTRDdigitizer::MakeDigits -- ");
326         printf("Digitizing chamber %d, plane %d, sector %d.\n"
327               ,iCham,iPlan,iSect);
328
329         Int_t   nRowMax     = Geo->GetRowMax(iPlan,iCham,iSect);
330         Int_t   nColMax     = Geo->GetColMax(iPlan);
331         Int_t   nTimeMax    = Geo->GetTimeMax();
332         Float_t row0        = Geo->GetRow0(iPlan,iCham,iSect);
333         Float_t col0        = Geo->GetCol0(iPlan);
334         Float_t time0       = Geo->GetTime0(iPlan);
335         Float_t rowPadSize  = Geo->GetRowPadSize();
336         Float_t colPadSize  = Geo->GetColPadSize();
337         Float_t timeBinSize = Geo->GetTimeBinSize();
338
339         // Create a detector matrix to keep the signal and track numbers
340         AliTRDmatrix *Matrix = new AliTRDmatrix(nRowMax,nColMax,nTimeMax
341                                                ,iSect,iCham,iPlan);
342
343         // Loop through all entries in the tree
344         for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
345
346           gAlice->ResetHits();
347           nBytes += HitTree->GetEvent(iTrack);
348
349           // Get the number of hits in the TRD created by this particle
350           Int_t nHit = TRD->Hits()->GetEntriesFast();
351
352           // Loop through the TRD hits  
353           for (Int_t iHit = 0; iHit < nHit; iHit++) {
354
355             if (!(Hit = (AliTRDhit *) TRD->Hits()->UncheckedAt(iHit))) 
356               continue;
357
358             Float_t pos[3];
359                     pos[0]   = Hit->fX;
360                     pos[1]   = Hit->fY;
361                     pos[2]   = Hit->fZ;
362             Float_t q        = Hit->fQ;
363             Int_t   track    = Hit->fTrack;
364             Int_t   detector = Hit->fDetector;
365             Int_t   plane    = Geo->GetPlane(detector);
366             Int_t   sector   = Geo->GetSector(detector);
367             Int_t   chamber  = Geo->GetChamber(detector);
368
369             if ((sector  != iSect) ||
370                 (plane   != iPlan) ||
371                 (chamber != iCham)) 
372               continue;
373
374             // Rotate the sectors on top of each other
375             Float_t rot[3];
376             Geo->Rotate(detector,pos,rot);
377
378             // The hit position in pad coordinates (center pad)
379             // The pad row (z-direction)
380             Int_t  rowH = (Int_t) ((rot[2] -  row0) /  rowPadSize);
381             // The pad column (rphi-direction)  
382             Int_t  colH = (Int_t) ((rot[1] -  col0) /  colPadSize);
383             // The time bucket
384             Int_t timeH = (Int_t) ((rot[0] - time0) / timeBinSize);
385
386             // Array to sum up the signal in a box surrounding the
387             // hit postition
388             const Int_t timeBox = 7;
389             const Int_t  colBox = 9;
390             const Int_t  rowBox = 7;
391             Float_t signalSum[rowBox][colBox][timeBox];
392             for (iRow  = 0;  iRow <  rowBox; iRow++ ) {
393               for (iCol  = 0;  iCol <  colBox; iCol++ ) {
394                 for (iTime = 0; iTime < timeBox; iTime++) {
395                   signalSum[iRow][iCol][iTime] = 0;
396                 }
397               }
398             }
399
400             // Loop over all electrons of this hit
401             Int_t nEl = (Int_t) q;
402             for (Int_t iEl = 0; iEl < nEl; iEl++) {
403
404               // The driftlength
405               Float_t driftlength = rot[0] - time0;
406               if ((driftlength <        0) || 
407                   (driftlength > kDrThick)) break;
408               Float_t driftlengthL = driftlength;
409               if (fExBOn) driftlengthL /= TMath::Sqrt(fLorentzFactor);
410               Float_t xyz[3];
411               xyz[0] = rot[0];
412               xyz[1] = rot[1];
413               xyz[2] = rot[2];
414
415               // Electron attachment
416               if (fElAttachOn) {
417                 if (gRandom->Rndm() < (driftlengthL * fElAttachProp / 100.)) continue;
418               }
419
420               // Apply the diffusion smearing
421               if (fDiffusionOn) {
422                 if (!(Diffusion(driftlengthL,xyz))) continue;
423               }
424
425               // Apply E x B effects
426               if (fExBOn) { 
427                 if (!(ExB(driftlength,xyz))) continue;   
428               }
429
430               // The electron position and the distance to the hit position
431               // in pad units
432               // The pad row (z-direction)
433               Int_t  rowE = (Int_t) ((xyz[2] -  row0) /  rowPadSize);
434               Int_t  rowD =  rowH -  rowE;
435               // The pad column (rphi-direction)
436               Int_t  colE = (Int_t) ((xyz[1] -  col0) /  colPadSize);
437               Int_t  colD =  colH -  colE;
438               // The time bucket
439               Int_t timeE = (Int_t) ((xyz[0] - time0) / timeBinSize);
440               Int_t timeD = timeH - timeE;
441
442               // Apply the gas gain including fluctuations
443               Int_t signal = (Int_t) (-fGasGain * TMath::Log(gRandom->Rndm()));
444
445               // The distance of the electron to the center of the pad 
446               // in units of pad width
447               Float_t dist = (xyz[1] - col0 - (colE + 0.5) * colPadSize) 
448                            / colPadSize;
449
450               // Sum up the signal in the different pixels
451               // and apply the pad response
452               Int_t  rowIdx =  rowD + (Int_t) ( rowBox / 2);
453               Int_t  colIdx =  colD + (Int_t) ( colBox / 2);
454               Int_t timeIdx = timeD + (Int_t) (timeBox / 2);
455
456               if (( rowIdx < 0) || ( rowIdx >  rowBox)) {
457                 printf("AliTRDdigitizer::MakeDigits -- ");
458                 printf("Boundary error. rowIdx = %d (%d)\n", rowIdx, rowBox);
459                 continue;
460               }
461               if (( colIdx < 0) || ( colIdx >  colBox)) {
462                 printf("AliTRDdigitizer::MakeDigits -- ");
463                 printf("Boundary error. colIdx = %d (%d)\n", colIdx, colBox);
464                 continue;
465               }
466               if ((timeIdx < 0) || (timeIdx > timeBox)) {
467                 printf("AliTRDdigitizer::MakeDigits -- ");
468                 printf("Boundary error. timeIdx = %d (%d)\n",timeIdx,timeBox);
469                 continue;
470               }
471               signalSum[rowIdx][colIdx-1][timeIdx] += PadResponse(dist-1.) * signal;
472               signalSum[rowIdx][colIdx  ][timeIdx] += PadResponse(dist   ) * signal;
473               signalSum[rowIdx][colIdx+1][timeIdx] += PadResponse(dist+1.) * signal;
474
475             }
476
477             // Add the padcluster to the detector matrix
478             for (iRow  = 0;  iRow <  rowBox; iRow++ ) {
479               for (iCol  = 0;  iCol <  colBox; iCol++ ) {
480                 for (iTime = 0; iTime < timeBox; iTime++) {
481
482                   Int_t  rowB =  rowH + iRow  - (Int_t) ( rowBox / 2); 
483                   Int_t  colB =  colH + iCol  - (Int_t) ( colBox / 2);
484                   Int_t timeB = timeH + iTime - (Int_t) (timeBox / 2);
485
486                   Float_t signalB = signalSum[iRow][iCol][iTime];
487                   if (signalB > 0.0) {
488                     Matrix->AddSignal(rowB,colB,timeB,signalB);
489                     if (!(Matrix->AddTrack(rowB,colB,timeB,track))) { 
490                       printf("AliTRDdigitizer::MakeDigits -- ");
491                       printf("More than three tracks in a pixel!\n");
492                     }
493                   }
494
495                 }
496               }
497             }
498
499           }
500
501         }
502
503         // Add a container for the digits of this detector
504         Digits = (AliTRDdataArray *) fDigitsArray->At(iDet);        
505         // Allocate memory space for the digits buffer
506         Digits->Allocate(nRowMax,nColMax,nTimeMax);
507
508         for (Int_t iDict = 0; iDict < kNDict; iDict++) {
509           Dictionary[iDict] = (AliTRDdataArray *) fDictionary[iDict]->At(iDet);
510           Dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeMax);
511         }
512
513         // Create the hits for this chamber
514         for (iRow  = 0; iRow  <  nRowMax; iRow++ ) {
515           for (iCol  = 0; iCol  <  nColMax; iCol++ ) {
516             for (iTime = 0; iTime < nTimeMax; iTime++) {         
517
518               Float_t signalAmp = Matrix->GetSignal(iRow,iCol,iTime);
519
520               // Add the noise
521               signalAmp  = TMath::Max((Float_t) gRandom->Gaus(signalAmp,fNoise)
522                                      ,(Float_t) 0.0);
523               // Convert to fC
524               signalAmp *= el2fC;
525               // Convert to mV
526               signalAmp *= fChipGain;
527               // Convert to ADC counts
528               Int_t adc  = (Int_t) (signalAmp * (fADCoutRange / fADCinRange));
529
530               // Store the amplitude of the digit
531               Digits->SetData(iRow,iCol,iTime,adc);
532
533               // Store the track index in the dictionary
534               // Note: We store index+1 in order to allow the array to be compressed
535               for (Int_t iDict = 0; iDict < kNDict; iDict++) {
536                 Dictionary[iDict]->SetData(iRow,iCol,iTime
537                                           ,Matrix->GetTrack(iRow,iCol,iTime,iDict)+1);
538               }
539
540               if (adc > fADCthreshold) nDigits++;
541
542             }
543           }
544         }
545
546         // Compress the arrays
547         Digits->Compress(1,fADCthreshold);
548         for (Int_t iDict = 0; iDict < kNDict; iDict++) {
549           Dictionary[iDict]->Compress(1,0);
550         }
551
552         totalSizeDigits += Digits->GetSize();
553         totalSizeDict0  += Dictionary[0]->GetSize();
554         totalSizeDict1  += Dictionary[1]->GetSize();
555         totalSizeDict2  += Dictionary[2]->GetSize();
556
557         printf("AliTRDdigitizer::MakeDigits -- ");
558         printf("Number of digits found: %d.\n",nDigits);
559
560         // Clean up
561         if (Matrix) delete Matrix;
562
563       }
564     }
565   }
566
567   printf("AliTRDdigitizer::MakeDigits -- ");
568   printf("Total digits data size = %d, %d, %d, %d\n",totalSizeDigits
569                                                     ,totalSizeDict0
570                                                     ,totalSizeDict1
571                                                     ,totalSizeDict2);        
572
573   return kTRUE;
574
575 }
576
577 //_____________________________________________________________________________
578 Bool_t AliTRDdigitizer::MakeBranch()
579 {
580   //
581   // Creates the branches for the digits and the dictionary
582   //
583
584   Int_t buffersize = 64000;
585
586   Bool_t status = kTRUE;
587
588   if (gAlice->TreeD()) {
589
590     // Make the branch for the digits
591     if (fDigitsArray) {
592       const AliTRDdataArray *Digits = 
593            (AliTRDdataArray *) fDigitsArray->At(0);
594       if (Digits) {
595         gAlice->TreeD()->Branch("TRDdigits",Digits->IsA()->GetName()
596                                            ,&Digits,buffersize,1);
597         printf("AliTRDdigitizer::MakeBranch -- ");
598         printf("Making branch TRDdigits\n");
599       }
600       else {
601         status = kFALSE;
602       }
603     }
604     else {
605       status = kFALSE;
606     }
607
608     // Make the branches for the dictionaries
609     for (Int_t iDict = 0; iDict < kNDict; iDict++) {
610
611       Char_t branchname[15];
612       sprintf(branchname,"TRDdictionary%d",iDict);
613       if (fDictionary[iDict]) {
614         const AliTRDdataArray *Dictionary = 
615              (AliTRDdataArray *) fDictionary[iDict]->At(0);
616         if (Dictionary) {
617           gAlice->TreeD()->Branch(branchname,Dictionary->IsA()->GetName()
618                                             ,&Dictionary,buffersize,1);
619           printf("AliTRDdigitizer::MakeBranch -- ");
620           printf("Making branch %s\n",branchname);
621         }
622         else {
623           status = kFALSE;
624         }
625       }
626       else {
627         status = kFALSE;
628       }
629     }
630
631   }
632   else {
633     status = kFALSE;
634   }
635
636   return status;
637
638 }
639
640 //_____________________________________________________________________________
641 Bool_t AliTRDdigitizer::WriteDigits()
642 {
643   //
644   // Writes out the TRD-digits and the dictionaries
645   //
646
647   // Create the branches
648   if (!(gAlice->TreeD()->GetBranch("TRDdigits"))) { 
649     if (!MakeBranch()) return kFALSE;
650   }
651
652   // Store the contents of the segment array in the tree
653   if (!fDigitsArray->StoreArray("TRDdigits")) {
654     printf("AliTRDdigitizer::WriteDigits -- ");
655     printf("Error while storing digits in branch TRDdigits\n");
656     return kFALSE;
657   }
658   for (Int_t iDict = 0; iDict < kNDict; iDict++) {
659     Char_t branchname[15];
660     sprintf(branchname,"TRDdictionary%d",iDict);
661     if (!fDictionary[iDict]->StoreArray(branchname)) {
662       printf("AliTRDdigitizer::WriteDigits -- ");
663       printf("Error while storing dictionary in branch %s\n",branchname);
664       return kFALSE;
665     }
666   }
667
668   // Write the new tree into the input file (use overwrite option)
669   Char_t treeName[7];
670   sprintf(treeName,"TreeD%d",fEvent);
671   printf("AliTRDdigitizer::WriteDigits -- ");
672   printf("Write the digits tree %s for event %d.\n"
673         ,treeName,fEvent);
674   gAlice->TreeD()->Write(treeName,2);
675  
676   return kTRUE;
677
678 }
679
680 ClassImp(AliTRDdigit)
681
682 //_____________________________________________________________________________
683 AliTRDdigit::AliTRDdigit(Int_t *digits):AliDigitNew()
684 {
685   //
686   // Create a TRD digit
687   //
688
689   // Store the volume hierarchy
690   fDetector  = digits[0];
691
692   // Store the row, pad, and time bucket number
693   fRow       = digits[1];
694   fCol       = digits[2];
695   fTime      = digits[3];
696
697   // Store the signal amplitude
698   fAmplitude = digits[4];
699
700 }