1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 Revision 1.4 2000/05/08 15:53:45 cblume
19 Resolved merge conflict
21 Revision 1.3 2000/04/28 14:49:27 cblume
22 Only one declaration of iDict in MakeDigits()
24 Revision 1.1.4.1 2000/05/08 14:42:04 cblume
25 Introduced AliTRDdigitsManager
27 Revision 1.1 2000/02/28 19:00:13 cblume
32 ///////////////////////////////////////////////////////////////////////////////
34 // Creates and handles digits from TRD hits //
36 // The following effects are included: //
39 // - Gas gain including fluctuations //
40 // - Pad-response (simple Gaussian approximation) //
41 // - Electronics noise //
42 // - Electronics gain //
45 // The corresponding parameter can be adjusted via the various //
46 // Set-functions. If these parameters are not explicitly set, default //
47 // values are used (see Init-function). //
48 // To produce digits from a root-file with TRD-hits use the //
49 // slowDigitsCreate.C macro. //
51 ///////////////////////////////////////////////////////////////////////////////
58 #include "AliTRDdigitizer.h"
59 #include "AliTRDdataArrayI.h"
60 #include "AliTRDdataArrayF.h"
61 #include "AliTRDdigitsManager.h"
63 ClassImp(AliTRDdigitizer)
65 //_____________________________________________________________________________
66 AliTRDdigitizer::AliTRDdigitizer():TNamed()
69 // AliTRDdigitizer default constructor
95 //_____________________________________________________________________________
96 AliTRDdigitizer::AliTRDdigitizer(const Text_t *name, const Text_t *title)
100 // AliTRDdigitizer default constructor
114 //_____________________________________________________________________________
115 AliTRDdigitizer::~AliTRDdigitizer()
127 if (fPRF) delete fPRF;
131 //_____________________________________________________________________________
132 Int_t AliTRDdigitizer::Diffusion(Float_t driftlength, Float_t *xyz)
135 // Applies the diffusion smearing to the position of a single electron
138 Float_t driftSqrt = TMath::Sqrt(driftlength);
139 Float_t sigmaT = driftSqrt * fDiffusionT;
140 Float_t sigmaL = driftSqrt * fDiffusionL;
141 xyz[0] = gRandom->Gaus(xyz[0], sigmaL * fLorentzFactor);
142 xyz[1] = gRandom->Gaus(xyz[1], sigmaT * fLorentzFactor);
143 xyz[2] = gRandom->Gaus(xyz[2], sigmaT);
148 //_____________________________________________________________________________
149 Int_t AliTRDdigitizer::ExB(Float_t driftlength, Float_t *xyz)
152 // Applies E x B effects to the position of a single electron
156 xyz[1] = xyz[1] + fLorentzAngle * driftlength;
163 //_____________________________________________________________________________
164 void AliTRDdigitizer::Init()
167 // Initializes the digitization procedure with standard values
170 // The default parameter for the digitization
178 // Transverse and longitudinal diffusion coefficients (Xe/Isobutane)
183 // Propability for electron attachment
189 // omega * tau. (tau ~ 12 * 10^-12, B = 0.2T)
190 fLorentzAngle = 17.6 * 12.0 * 0.2 * 0.01;
192 // The pad response function
193 fPRF = new TF1("PRF","[0]*([1]+exp(-x*x/(2.0*[2])))",-2,2);
194 fPRF->SetParameter(0, 0.8872);
195 fPRF->SetParameter(1,-0.00573);
196 fPRF->SetParameter(2, 0.454 * 0.454);
200 //_____________________________________________________________________________
201 Bool_t AliTRDdigitizer::Open(const Char_t *name, Int_t nEvent)
204 // Opens a ROOT-file with TRD-hits and reads in the hit-tree
207 // Connect the AliRoot file containing Geometry, Kine, and Hits
208 fInputFile = (TFile*) gROOT->GetListOfFiles()->FindObject(name);
210 printf("AliTRDdigitizer::Open -- ");
211 printf("Open the ALIROOT-file %s.\n",name);
212 fInputFile = new TFile(name,"UPDATE");
215 printf("AliTRDdigitizer::Open -- ");
216 printf("%s is already open.\n",name);
219 gAlice = (AliRun*) fInputFile->Get("gAlice");
221 printf("AliTRDdigitizer::Open -- ");
222 printf("AliRun object found on file.\n");
225 printf("AliTRDdigitizer::Open -- ");
226 printf("Could not find AliRun object.\n");
232 // Import the Trees for the event nEvent in the file
233 Int_t nparticles = gAlice->GetEvent(fEvent);
234 if (nparticles <= 0) {
235 printf("AliTRDdigitizer::Open -- ");
236 printf("No entries in the trees for event %d.\n",fEvent);
244 //_____________________________________________________________________________
245 Bool_t AliTRDdigitizer::MakeDigits()
248 // Loops through the TRD-hits and creates the digits.
251 ///////////////////////////////////////////////////////////////
253 ///////////////////////////////////////////////////////////////
255 // Converts number of electrons to fC
256 const Float_t el2fC = 1.602E-19 * 1.0E15;
258 ///////////////////////////////////////////////////////////////
260 Int_t iRow, iCol, iTime;
264 Int_t totalSizeDigits = 0;
265 Int_t totalSizeDict0 = 0;
266 Int_t totalSizeDict1 = 0;
267 Int_t totalSizeDict2 = 0;
269 AliTRDdataArrayI *Digits;
270 AliTRDdataArrayI *Dictionary[kNDict];
273 printf("AliTRDdigitizer::MakeDigits -- ");
274 printf("No geometry defined\n");
278 // Create a digits manager
279 fDigits = new AliTRDdigitsManager();
281 // Create detector arrays to keep the signal and track numbers
282 AliTRDdataArrayF *Signal = new AliTRDdataArrayF();
283 AliTRDdataArrayI *Tracks[kNDict];
284 for (iDict = 0; iDict < kNDict; iDict++) {
285 Tracks[iDict] = new AliTRDdataArrayI();
288 // Get the pointer to the hit tree
289 TTree *HitTree = gAlice->TreeH();
291 // Get the number of entries in the hit tree
292 // (Number of primary particles creating a hit somewhere)
293 Int_t nTrack = (Int_t) HitTree->GetEntries();
295 printf("AliTRDdigitizer::MakeDigits -- ");
296 printf("Start creating digits.\n");
298 // The Lorentz factor
300 fLorentzFactor = 1.0 / (1.0 + fLorentzAngle*fLorentzAngle);
303 fLorentzFactor = 1.0;
307 Int_t chamEnd = kNcham;
308 if (fTRD->GetSensChamber() >= 0) {
309 chamBeg = fTRD->GetSensChamber();
310 chamEnd = chamBeg + 1;
313 Int_t planEnd = kNplan;
314 if (fTRD->GetSensPlane() >= 0) {
315 planBeg = fTRD->GetSensPlane();
316 planEnd = planBeg + 1;
319 Int_t sectEnd = kNsect;
320 if (fTRD->GetSensSector() >= 0) {
321 sectBeg = fTRD->GetSensSector();
322 sectEnd = sectBeg + 1;
325 Int_t count_hits = 0;
327 // Loop through all the chambers
328 for (Int_t iCham = chamBeg; iCham < chamEnd; iCham++) {
329 for (Int_t iPlan = planBeg; iPlan < planEnd; iPlan++) {
330 for (Int_t iSect = sectBeg; iSect < sectEnd; iSect++) {
334 printf("AliTRDdigitizer::MakeDigits -- ");
335 printf("Digitizing chamber %d, plane %d, sector %d.\n"
338 Int_t iDet = fGeo->GetDetector(iPlan,iCham,iSect);
339 Int_t nRowMax = fGeo->GetRowMax(iPlan,iCham,iSect);
340 Int_t nColMax = fGeo->GetColMax(iPlan);
341 Int_t nTimeMax = fGeo->GetTimeMax();
342 Float_t row0 = fGeo->GetRow0(iPlan,iCham,iSect);
343 Float_t col0 = fGeo->GetCol0(iPlan);
344 Float_t time0 = fGeo->GetTime0(iPlan);
345 Float_t rowPadSize = fGeo->GetRowPadSize();
346 Float_t colPadSize = fGeo->GetColPadSize();
347 Float_t timeBinSize = fGeo->GetTimeBinSize();
349 // Adjust the size of the detector arrays
350 Signal->Allocate(nRowMax,nColMax,nTimeMax);
351 for (iDict = 0; iDict < kNDict; iDict++) {
352 Tracks[iDict]->Allocate(nRowMax,nColMax,nTimeMax);
355 // Loop through all entries in the tree
356 for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
359 nBytes += HitTree->GetEvent(iTrack);
361 // Get the number of hits in the TRD created by this particle
362 Int_t nHit = fTRD->Hits()->GetEntriesFast();
364 // Loop through the TRD hits
365 for (Int_t iHit = 0; iHit < nHit; iHit++) {
369 AliTRDhit *Hit = (AliTRDhit *) fTRD->Hits()->UncheckedAt(iHit);
375 Int_t track = Hit->fTrack;
376 Int_t detector = Hit->fDetector;
377 Int_t plane = fGeo->GetPlane(detector);
378 Int_t sector = fGeo->GetSector(detector);
379 Int_t chamber = fGeo->GetChamber(detector);
381 if ((sector != iSect) ||
386 // Rotate the sectors on top of each other
388 fGeo->Rotate(detector,pos,rot);
390 // The hit position in pad coordinates (center pad)
391 // The pad row (z-direction)
392 Int_t rowH = (Int_t) ((rot[2] - row0) / rowPadSize);
393 // The pad column (rphi-direction)
394 Int_t colH = (Int_t) ((rot[1] - col0) / colPadSize);
396 Int_t timeH = (Int_t) ((rot[0] - time0) / timeBinSize);
398 // Array to sum up the signal in a box surrounding the
400 const Int_t timeBox = 7;
401 const Int_t colBox = 9;
402 const Int_t rowBox = 7;
403 Float_t signalSum[rowBox][colBox][timeBox];
404 for (iRow = 0; iRow < rowBox; iRow++ ) {
405 for (iCol = 0; iCol < colBox; iCol++ ) {
406 for (iTime = 0; iTime < timeBox; iTime++) {
407 signalSum[iRow][iCol][iTime] = 0;
412 // Loop over all electrons of this hit
413 Int_t nEl = (Int_t) q;
414 for (Int_t iEl = 0; iEl < nEl; iEl++) {
417 Float_t driftlength = rot[0] - time0;
418 if ((driftlength < 0) ||
419 (driftlength > kDrThick)) break;
420 Float_t driftlengthL = driftlength;
421 if (fExBOn) driftlengthL /= TMath::Sqrt(fLorentzFactor);
427 // Electron attachment
429 if (gRandom->Rndm() < (driftlengthL * fElAttachProp / 100.)) continue;
432 // Apply the diffusion smearing
434 if (!(Diffusion(driftlengthL,xyz))) continue;
437 // Apply E x B effects
439 if (!(ExB(driftlength,xyz))) continue;
442 // The electron position and the distance to the hit position
444 // The pad row (z-direction)
445 Int_t rowE = (Int_t) ((xyz[2] - row0) / rowPadSize);
446 Int_t rowD = rowH - rowE;
447 // The pad column (rphi-direction)
448 Int_t colE = (Int_t) ((xyz[1] - col0) / colPadSize);
449 Int_t colD = colH - colE;
451 Int_t timeE = (Int_t) ((xyz[0] - time0) / timeBinSize);
452 Int_t timeD = timeH - timeE;
454 // Apply the gas gain including fluctuations
455 Int_t signal = (Int_t) (-fGasGain * TMath::Log(gRandom->Rndm()));
457 // The distance of the electron to the center of the pad
458 // in units of pad width
459 Float_t dist = (xyz[1] - col0 - (colE + 0.5) * colPadSize)
462 // Sum up the signal in the different pixels
463 // and apply the pad response
464 Int_t rowIdx = rowD + (Int_t) ( rowBox / 2);
465 Int_t colIdx = colD + (Int_t) ( colBox / 2);
466 Int_t timeIdx = timeD + (Int_t) (timeBox / 2);
468 if (( rowIdx < 0) || ( rowIdx > rowBox)) {
469 printf("AliTRDdigitizer::MakeDigits -- ");
470 printf("Boundary error. rowIdx = %d (%d)\n", rowIdx, rowBox);
473 if (( colIdx < 0) || ( colIdx > colBox)) {
474 printf("AliTRDdigitizer::MakeDigits -- ");
475 printf("Boundary error. colIdx = %d (%d)\n", colIdx, colBox);
478 if ((timeIdx < 0) || (timeIdx > timeBox)) {
479 printf("AliTRDdigitizer::MakeDigits -- ");
480 printf("Boundary error. timeIdx = %d (%d)\n",timeIdx,timeBox);
483 signalSum[rowIdx][colIdx-1][timeIdx] += fPRF->Eval(dist-1.0,0,0) * signal;
484 signalSum[rowIdx][colIdx ][timeIdx] += fPRF->Eval(dist ,0,0) * signal;
485 signalSum[rowIdx][colIdx+1][timeIdx] += fPRF->Eval(dist+1.0,0,0) * signal;
489 // Add the padcluster to the detector matrix
490 for (iRow = 0; iRow < rowBox; iRow++ ) {
491 for (iCol = 0; iCol < colBox; iCol++ ) {
492 for (iTime = 0; iTime < timeBox; iTime++) {
494 Int_t rowB = rowH + iRow - (Int_t) ( rowBox / 2);
495 Int_t colB = colH + iCol - (Int_t) ( colBox / 2);
496 Int_t timeB = timeH + iTime - (Int_t) (timeBox / 2);
497 Float_t signalB = signalSum[iRow][iCol][iTime];
498 if (( rowB < 0) || ( rowB >= nRowMax)) continue;
499 if (( colB < 0) || ( colB >= nColMax)) continue;
500 if ((timeB < 0) || (timeB >= nTimeMax)) continue;
503 // Add the signal sum
504 signalB += Signal->GetData(rowB,colB,timeB);
505 Signal->SetData(rowB,colB,timeB,signalB);
506 // Store the track index in the dictionary
507 // Note: We store index+1 in order to allow the array to be compressed
508 for (iDict = 0; iDict < kNDict; iDict++) {
509 Int_t oldTrack = Tracks[iDict]->GetData(rowB,colB,timeB);
510 if (oldTrack == track+1) break;
511 if (oldTrack == -1) break;
513 Tracks[iDict]->SetData(rowB,colB,timeB,track+1);
517 if (iDict == kNDict) {
518 printf("AliTRDdigitizer::MakeDigits -- ");
519 printf("More than three tracks for one digit!\n");
531 // Add a container for the digits of this detector
532 Digits = fDigits->GetDigits(iDet);
533 // Allocate memory space for the digits buffer
534 Digits->Allocate(nRowMax,nColMax,nTimeMax);
536 // Do the same for the dictionary arrays
537 for (iDict = 0; iDict < kNDict; iDict++) {
538 Dictionary[iDict] = fDigits->GetDictionary(iDet,iDict);
539 Dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeMax);
542 // Create the digits for this chamber
543 for (iRow = 0; iRow < nRowMax; iRow++ ) {
544 for (iCol = 0; iCol < nColMax; iCol++ ) {
545 for (iTime = 0; iTime < nTimeMax; iTime++) {
547 Float_t signalAmp = Signal->GetData(iRow,iCol,iTime);
550 signalAmp = TMath::Max((Float_t) gRandom->Gaus(signalAmp,fNoise)
555 signalAmp *= fChipGain;
556 // Convert to ADC counts
557 Int_t adc = (Int_t) (signalAmp * (fADCoutRange / fADCinRange));
559 if (adc > fADCthreshold) {
563 // Store the amplitude of the digit
564 Digits->SetData(iRow,iCol,iTime,adc);
566 // Store the track index in the dictionary
567 // Note: We store index+1 in order to allow the array to be compressed
568 for (iDict = 0; iDict < kNDict; iDict++) {
569 Dictionary[iDict]->SetData(iRow,iCol,iTime
570 ,Tracks[iDict]->GetData(iRow,iCol,iTime));
579 // Compress the arrays
580 Digits->Compress(1,0);
581 for (iDict = 0; iDict < kNDict; iDict++) {
582 Dictionary[iDict]->Compress(1,0);
585 totalSizeDigits += Digits->GetSize();
586 totalSizeDict0 += Dictionary[0]->GetSize();
587 totalSizeDict1 += Dictionary[1]->GetSize();
588 totalSizeDict2 += Dictionary[2]->GetSize();
590 printf("AliTRDdigitizer::MakeDigits -- ");
591 printf("Number of digits found: %d.\n",nDigits);
595 for (iDict = 0; iDict < kNDict; iDict++) {
596 Tracks[iDict]->Reset();
603 printf("AliTRDdigitizer::MakeDigits -- ");
604 printf("Total number of analyzed hits = %d\n",count_hits);
606 printf("AliTRDdigitizer::MakeDigits -- ");
607 printf("Total digits data size = %d, %d, %d, %d\n",totalSizeDigits
616 //_____________________________________________________________________________
617 Bool_t AliTRDdigitizer::WriteDigits()
620 // Writes out the TRD-digits and the dictionaries
623 // Create the branches
624 if (!(gAlice->TreeD()->GetBranch("TRDdigits"))) {
625 if (!fDigits->MakeBranch()) return kFALSE;
628 // Store the digits and the dictionary in the tree
629 fDigits->WriteDigits();
631 // Write the new tree into the input file (use overwrite option)
633 sprintf(treeName,"TreeD%d",fEvent);
634 printf("AliTRDdigitizer::WriteDigits -- ");
635 printf("Write the digits tree %s for event %d.\n"
637 gAlice->TreeD()->Write(treeName,2);