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