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.36 2002/04/12 12:13:23 cblume
21 Revision 1.35 2002/03/28 14:59:07 cblume
24 Revision 1.34 2002/03/25 20:00:44 cblume
25 Introduce parameter class and regions of interest for merging
27 Revision 1.33 2002/02/12 17:32:03 cblume
28 Rearrange the deleting of the list of sdigitsmanager
30 Revision 1.32 2002/02/12 16:07:21 cblume
33 Revision 1.31 2002/02/11 14:27:11 cblume
34 New pad plane design, new TRF+PRF, tail cancelation, cross talk
36 Revision 1.30 2001/11/19 08:44:08 cblume
37 Fix bugs reported by Rene
39 Revision 1.29 2001/11/14 19:44:25 hristov
40 Numeric const casted (Alpha)
42 Revision 1.28 2001/11/14 16:35:58 cblume
43 Inherits now from AliDetector
45 Revision 1.27 2001/11/14 10:50:45 cblume
46 Changes in digits IO. Add merging of summable digits
48 Revision 1.26 2001/11/06 17:19:41 cblume
49 Add detailed geometry and simple simulator
51 Revision 1.25 2001/06/27 09:54:44 cblume
52 Moved fField initialization to InitDetector()
54 Revision 1.24 2001/05/21 16:45:47 hristov
55 Last minute changes (C.Blume)
57 Revision 1.23 2001/05/07 08:04:48 cblume
58 New TRF and PRF. Speedup of the code. Digits from amplification region included
60 Revision 1.22 2001/03/30 14:40:14 cblume
61 Update of the digitization parameter
63 Revision 1.21 2001/03/13 09:30:35 cblume
64 Update of digitization. Moved digit branch definition to AliTRD
66 Revision 1.20 2001/02/25 20:19:00 hristov
67 Minor correction: loop variable declared only once for HP, Sun
69 Revision 1.19 2001/02/14 18:22:26 cblume
70 Change in the geometry of the padplane
72 Revision 1.18 2001/01/26 19:56:57 hristov
73 Major upgrade of AliRoot code
75 Revision 1.17 2000/12/08 12:53:27 cblume
76 Change in Copy() function for HP-compiler
78 Revision 1.16 2000/12/07 12:20:46 cblume
79 Go back to array compression. Use sampled PRF to speed up digitization
81 Revision 1.15 2000/11/23 14:34:08 cblume
82 Fixed bug in expansion routine of arrays (initialize buffers properly)
84 Revision 1.14 2000/11/20 08:54:44 cblume
85 Switch off compression as default
87 Revision 1.13 2000/11/10 14:57:52 cblume
88 Changes in the geometry constants for the DEC compiler
90 Revision 1.12 2000/11/01 14:53:20 cblume
91 Merge with TRD-develop
93 Revision 1.1.4.9 2000/10/26 17:00:22 cblume
94 Fixed bug in CheckDetector()
96 Revision 1.1.4.8 2000/10/23 13:41:35 cblume
97 Added protection against Log(0) in the gas gain calulation
99 Revision 1.1.4.7 2000/10/17 02:27:34 cblume
100 Get rid of global constants
102 Revision 1.1.4.6 2000/10/16 01:16:53 cblume
103 Changed timebin 0 to be the one closest to the readout
105 Revision 1.1.4.5 2000/10/15 23:34:29 cblume
106 Faster version of the digitizer
108 Revision 1.1.4.4 2000/10/06 16:49:46 cblume
111 Revision 1.1.4.3 2000/10/04 16:34:58 cblume
112 Replace include files by forward declarations
114 Revision 1.1.4.2 2000/09/22 14:41:10 cblume
115 Bug fix in PRF. Included time response. New structure
117 Revision 1.10 2000/10/05 07:27:53 cblume
118 Changes in the header-files by FCA
120 Revision 1.9 2000/10/02 21:28:19 fca
121 Removal of useless dependecies via forward declarations
123 Revision 1.8 2000/06/09 11:10:07 cblume
124 Compiler warnings and coding conventions, next round
126 Revision 1.7 2000/06/08 18:32:58 cblume
127 Make code compliant to coding conventions
129 Revision 1.6 2000/06/07 16:27:32 cblume
130 Try to remove compiler warnings on Sun and HP
132 Revision 1.5 2000/05/09 16:38:57 cblume
133 Removed PadResponse(). Merge problem
135 Revision 1.4 2000/05/08 15:53:45 cblume
136 Resolved merge conflict
138 Revision 1.3 2000/04/28 14:49:27 cblume
139 Only one declaration of iDict in MakeDigits()
141 Revision 1.1.4.1 2000/05/08 14:42:04 cblume
142 Introduced AliTRDdigitsManager
144 Revision 1.1 2000/02/28 19:00:13 cblume
149 ///////////////////////////////////////////////////////////////////////////////
151 // Creates and handles digits from TRD hits //
152 // Author: C. Blume (C.Blume@gsi.de) //
154 // The following effects are included: //
157 // - Gas gain including fluctuations //
158 // - Pad-response (simple Gaussian approximation) //
159 // - Time-response //
160 // - Electronics noise //
161 // - Electronics gain //
163 // - ADC threshold //
164 // The corresponding parameter can be adjusted via the various //
165 // Set-functions. If these parameters are not explicitly set, default //
166 // values are used (see Init-function). //
167 // As an example on how to use this class to produce digits from hits //
168 // have a look at the macro hits2digits.C //
169 // The production of summable digits is demonstrated in hits2sdigits.C //
170 // and the subsequent conversion of the s-digits into normal digits is //
171 // explained in sdigits2digits.C. //
173 ///////////////////////////////////////////////////////////////////////////////
189 #include "AliRunDigitizer.h"
192 #include "AliTRDhit.h"
193 #include "AliTRDdigitizer.h"
194 #include "AliTRDdataArrayI.h"
195 #include "AliTRDdataArrayF.h"
196 #include "AliTRDsegmentArray.h"
197 #include "AliTRDdigitsManager.h"
198 #include "AliTRDgeometry.h"
199 #include "AliTRDparameter.h"
201 ClassImp(AliTRDdigitizer)
203 //_____________________________________________________________________________
204 AliTRDdigitizer::AliTRDdigitizer()
207 // AliTRDdigitizer default constructor
212 fSDigitsManagerList = 0;
223 fMergeSignalOnly = kFALSE;
227 //_____________________________________________________________________________
228 AliTRDdigitizer::AliTRDdigitizer(const Text_t *name, const Text_t *title)
229 :AliDigitizer(name,title)
232 // AliTRDdigitizer constructor
237 fSDigitsManagerList = 0;
247 fMergeSignalOnly = kFALSE;
249 // For the summable digits
250 fSDigitsScale = 100.;
254 //_____________________________________________________________________________
255 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager
256 , const Text_t *name, const Text_t *title)
257 :AliDigitizer(manager,name,title)
260 // AliTRDdigitizer constructor
265 fSDigitsManagerList = 0;
275 fMergeSignalOnly = kFALSE;
277 // For the summable digits
278 fSDigitsScale = 100.;
282 //_____________________________________________________________________________
283 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager)
284 :AliDigitizer(manager,"AliTRDdigitizer","TRD digitizer")
287 // AliTRDdigitizer constructor
292 fSDigitsManagerList = 0;
302 fMergeSignalOnly = kFALSE;
304 // For the summable digits
305 fSDigitsScale = 100.;
309 //_____________________________________________________________________________
310 AliTRDdigitizer::AliTRDdigitizer(const AliTRDdigitizer &d)
313 // AliTRDdigitizer copy constructor
316 ((AliTRDdigitizer &) d).Copy(*this);
320 //_____________________________________________________________________________
321 AliTRDdigitizer::~AliTRDdigitizer()
324 // AliTRDdigitizer destructor
333 if (fDigitsManager) {
334 delete fDigitsManager;
338 if (fSDigitsManager) {
339 delete fSDigitsManager;
343 if (fSDigitsManagerList) {
344 delete fSDigitsManagerList;
345 fSDigitsManagerList = 0;
355 //_____________________________________________________________________________
356 AliTRDdigitizer &AliTRDdigitizer::operator=(const AliTRDdigitizer &d)
359 // Assignment operator
362 if (this != &d) ((AliTRDdigitizer &) d).Copy(*this);
367 //_____________________________________________________________________________
368 void AliTRDdigitizer::Copy(TObject &d)
374 ((AliTRDdigitizer &) d).fInputFile = 0;
375 ((AliTRDdigitizer &) d).fSDigitsManagerList = 0;
376 ((AliTRDdigitizer &) d).fSDigitsManager = 0;
377 ((AliTRDdigitizer &) d).fDigitsManager = 0;
378 ((AliTRDdigitizer &) d).fTRD = 0;
379 ((AliTRDdigitizer &) d).fGeo = 0;
380 ((AliTRDdigitizer &) d).fMasks = 0;
381 ((AliTRDdigitizer &) d).fEvent = 0;
382 ((AliTRDdigitizer &) d).fPar = 0;
383 ((AliTRDdigitizer &) d).fCompress = fCompress;
384 ((AliTRDdigitizer &) d).fDebug = fDebug ;
385 ((AliTRDdigitizer &) d).fSDigits = fSDigits;
386 ((AliTRDdigitizer &) d).fSDigitsScale = fSDigitsScale;
387 ((AliTRDdigitizer &) d).fMergeSignalOnly = fMergeSignalOnly;
391 //_____________________________________________________________________________
392 void AliTRDdigitizer::Exec(Option_t* option)
395 // Executes the merging
400 AliTRDdigitsManager *sdigitsManager;
402 TString optionString = option;
403 if (optionString.Contains("deb")) {
405 if (optionString.Contains("2")) {
408 printf("<AliTRDdigitizer::Exec> ");
409 printf("Called with debug option %d\n",fDebug);
412 // The AliRoot file is already connected by the manager
415 printf("<AliTRDdigitizer::Exec> ");
416 printf("AliRun object found on file.\n");
420 printf("<AliTRDdigitizer::Exec> ");
421 printf("Could not find AliRun object.\n");
425 Int_t nInput = fManager->GetNinputs();
426 fMasks = new Int_t[nInput];
427 for (iInput = 0; iInput < nInput; iInput++) {
428 fMasks[iInput] = fManager->GetMask(iInput);
434 for (iInput = 0; iInput < nInput; iInput++) {
437 printf("<AliTRDdigitizer::Exec> ");
438 printf("Add input stream %d\n",iInput);
441 // Read the s-digits via digits manager
442 sdigitsManager = new AliTRDdigitsManager();
443 sdigitsManager->SetDebug(fDebug);
444 sdigitsManager->SetSDigits(kTRUE);
445 sdigitsManager->ReadDigits(fManager->GetInputTreeTRDS(iInput));
447 // Add the s-digits to the input list
448 AddSDigitsManager(sdigitsManager);
452 // Convert the s-digits to normal digits
454 printf("<AliTRDdigitizer::Exec> ");
455 printf("Do the conversion\n");
461 printf("<AliTRDdigitizer::Exec> ");
462 printf("Write the digits\n");
464 fDigitsManager->MakeBranch(fManager->GetTreeDTRD());
465 fDigitsManager->WriteDigits();
467 printf("<AliTRDdigitizer::Exec> ");
471 DeleteSDigitsManager();
475 //_____________________________________________________________________________
476 Bool_t AliTRDdigitizer::Open(const Char_t *file, Int_t nEvent)
479 // Opens a ROOT-file with TRD-hits and reads in the hit-tree
482 // Connect the AliRoot file containing Geometry, Kine, and Hits
483 fInputFile = (TFile *) gROOT->GetListOfFiles()->FindObject(file);
486 printf("<AliTRDdigitizer::Open> ");
487 printf("Open the AliROOT-file %s.\n",file);
489 fInputFile = new TFile(file,"UPDATE");
493 printf("<AliTRDdigitizer::Open> ");
494 printf("%s is already open.\n",file);
498 gAlice = (AliRun *) fInputFile->Get("gAlice");
501 printf("<AliTRDdigitizer::Open> ");
502 printf("AliRun object found on file.\n");
506 printf("<AliTRDdigitizer::Open> ");
507 printf("Could not find AliRun object.\n");
513 // Import the Trees for the event nEvent in the file
514 Int_t nparticles = gAlice->GetEvent(fEvent);
515 if (nparticles <= 0) {
516 printf("<AliTRDdigitizer::Open> ");
517 printf("No entries in the trees for event %d.\n",fEvent);
521 if (InitDetector()) {
530 //_____________________________________________________________________________
531 Bool_t AliTRDdigitizer::InitDetector()
534 // Sets the pointer to the TRD detector and the geometry
537 // Get the pointer to the detector class and check for version 1
538 fTRD = (AliTRD *) gAlice->GetDetector("TRD");
540 printf("<AliTRDdigitizer::InitDetector> ");
541 printf("No TRD module found\n");
544 if (fTRD->IsVersion() != 1) {
545 printf("<AliTRDdigitizer::InitDetector> ");
546 printf("TRD must be version 1 (slow simulator).\n");
551 fGeo = fTRD->GetGeometry();
553 printf("<AliTRDdigitizer::InitDetector> ");
554 printf("Geometry version %d\n",fGeo->IsVersion());
557 // Create a digits manager
558 fDigitsManager = new AliTRDdigitsManager();
559 fDigitsManager->SetSDigits(fSDigits);
560 fDigitsManager->CreateArrays();
561 fDigitsManager->SetEvent(fEvent);
562 fDigitsManager->SetDebug(fDebug);
564 // The list for the input s-digits manager to be merged
565 fSDigitsManagerList = new TList();
571 //_____________________________________________________________________________
572 Bool_t AliTRDdigitizer::MakeBranch(const Char_t *file) const
575 // Create the branches for the digits array
578 return fDigitsManager->MakeBranch(file);
582 //_____________________________________________________________________________
583 Bool_t AliTRDdigitizer::MakeDigits()
589 ///////////////////////////////////////////////////////////////
591 ///////////////////////////////////////////////////////////////
593 // Converts number of electrons to fC
594 const Double_t kEl2fC = 1.602E-19 * 1.0E15;
596 ///////////////////////////////////////////////////////////////
598 // Number of pads included in the pad response
599 const Int_t kNpad = 3;
601 // Number of track dictionary arrays
602 const Int_t kNDict = AliTRDdigitsManager::kNDict;
604 // Half the width of the amplification region
605 const Float_t kAmWidth = AliTRDgeometry::AmThick() / 2.;
607 Int_t iRow, iCol, iTime, iPad;
611 Int_t totalSizeDigits = 0;
612 Int_t totalSizeDict0 = 0;
613 Int_t totalSizeDict1 = 0;
614 Int_t totalSizeDict2 = 0;
616 Int_t timeTRDbeg = 0;
617 Int_t timeTRDend = 1;
622 Float_t padSignal[kNpad];
623 Float_t signalOld[kNpad];
625 AliTRDdataArrayF *signals = 0;
626 AliTRDdataArrayI *digits = 0;
627 AliTRDdataArrayI *dictionary[kNDict];
629 // Create a default parameter class if none is defined
631 fPar = new AliTRDparameter("TRDparameter","Standard TRD parameter");
633 printf("<AliTRDdigitizer::MakeDigits> ");
634 printf("Create the default parameter object\n");
638 // Create a container for the amplitudes
639 AliTRDsegmentArray *signalsArray
640 = new AliTRDsegmentArray("AliTRDdataArrayF"
641 ,AliTRDgeometry::Ndet());
644 timeTRDbeg = ((Int_t) (-fPar->GetTRFlo() / fPar->GetTimeBinSize())) - 1;
645 timeTRDend = ((Int_t) ( fPar->GetTRFhi() / fPar->GetTimeBinSize())) - 1;
647 printf("<AliTRDdigitizer::MakeDigits> ");
648 printf("Sample the TRF between -%d and %d\n",timeTRDbeg,timeTRDend);
652 Float_t elAttachProp = fPar->GetElAttachProp() / 100.;
655 printf("<AliTRDdigitizer::MakeDigits> ");
656 printf("No geometry defined\n");
661 printf("<AliTRDdigitizer::MakeDigits> ");
662 printf("Start creating digits.\n");
665 // Get the pointer to the hit tree
666 TTree *hitTree = gAlice->TreeH();
668 // Get the number of entries in the hit tree
669 // (Number of primary particles creating a hit somewhere)
670 Int_t nTrack = (Int_t) hitTree->GetEntries();
672 printf("<AliTRDdigitizer::MakeDigits> ");
673 printf("Found %d primary particles\n",nTrack);
676 Int_t detectorOld = -1;
679 // Loop through all entries in the tree
680 for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
683 nBytes += hitTree->GetEvent(iTrack);
685 // Loop through the TRD hits
687 AliTRDhit *hit = (AliTRDhit *) fTRD->FirstHit(-1);
696 Float_t q = hit->GetCharge();
697 Int_t track = hit->Track();
698 Int_t detector = hit->GetDetector();
699 Int_t plane = fGeo->GetPlane(detector);
700 Int_t sector = fGeo->GetSector(detector);
701 Int_t chamber = fGeo->GetChamber(detector);
702 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
703 Int_t nColMax = fPar->GetColMax(plane);
704 Int_t nTimeMax = fPar->GetTimeMax();
705 Int_t nTimeBefore = fPar->GetTimeBefore();
706 Int_t nTimeAfter = fPar->GetTimeAfter();
707 Int_t nTimeTotal = fPar->GetTimeTotal();
708 Float_t row0 = fPar->GetRow0(plane,chamber,sector);
709 Float_t col0 = fPar->GetCol0(plane);
710 Float_t time0 = fPar->GetTime0(plane);
711 Float_t rowPadSize = fPar->GetRowPadSize(plane,chamber,sector);
712 Float_t colPadSize = fPar->GetColPadSize(plane);
713 Float_t timeBinSize = fPar->GetTimeBinSize();
714 Float_t divideRow = 1.0 / rowPadSize;
715 Float_t divideCol = 1.0 / colPadSize;
716 Float_t divideTime = 1.0 / timeBinSize;
719 printf("Analyze hit no. %d ",iHit);
720 printf("-----------------------------------------------------------\n");
722 printf("plane = %d, sector = %d, chamber = %d\n"
723 ,plane,sector,chamber);
724 printf("nRowMax = %d, nColMax = %d, nTimeMax = %d\n"
725 ,nRowMax,nColMax,nTimeMax);
726 printf("nTimeBefore = %d, nTimeAfter = %d, nTimeTotal = %d\n"
727 ,nTimeBefore,nTimeAfter,nTimeTotal);
728 printf("row0 = %f, col0 = %f, time0 = %f\n"
730 printf("rowPadSize = %f, colPadSize = %f, timeBinSize = %f\n"
731 ,rowPadSize,colPadSize,timeBinSize);
734 // Don't analyze test hits and switched off detectors
735 if ((CheckDetector(plane,chamber,sector)) &&
736 (((Int_t) q) != 0)) {
738 if (detector != detectorOld) {
741 printf("<AliTRDdigitizer::MakeDigits> ");
742 printf("Get new container. New det = %d, Old det = %d\n"
743 ,detector,detectorOld);
745 // Compress the old one if enabled
746 if ((fCompress) && (detectorOld > -1)) {
748 printf("<AliTRDdigitizer::MakeDigits> ");
749 printf("Compress the old container ...");
751 signals->Compress(1,0);
752 for (iDict = 0; iDict < kNDict; iDict++) {
753 dictionary[iDict]->Compress(1,0);
755 if (fDebug > 1) printf("done\n");
757 // Get the new container
758 signals = (AliTRDdataArrayF *) signalsArray->At(detector);
759 if (signals->GetNtime() == 0) {
760 // Allocate a new one if not yet existing
762 printf("<AliTRDdigitizer::MakeDigits> ");
763 printf("Allocate a new container ... ");
765 signals->Allocate(nRowMax,nColMax,nTimeTotal);
768 // Expand an existing one
771 printf("<AliTRDdigitizer::MakeDigits> ");
772 printf("Expand an existing container ... ");
777 // The same for the dictionary
778 for (iDict = 0; iDict < kNDict; iDict++) {
779 dictionary[iDict] = fDigitsManager->GetDictionary(detector,iDict);
780 if (dictionary[iDict]->GetNtime() == 0) {
781 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
784 if (fCompress) dictionary[iDict]->Expand();
787 if (fDebug > 1) printf("done\n");
788 detectorOld = detector;
791 // Rotate the sectors on top of each other
792 fGeo->Rotate(detector,pos,rot);
794 // The driftlength. It is negative if the hit is in the
795 // amplification region.
796 Float_t driftlength = time0 - rot[0];
798 // Take also the drift in the amplification region into account
799 // The drift length is at the moment still the same, regardless of
800 // the position relativ to the wire. This non-isochronity needs still
801 // to be implemented.
802 Float_t driftlengthL = TMath::Abs(driftlength + kAmWidth);
803 if (fPar->ExBOn()) driftlengthL /= TMath::Sqrt(fPar->GetLorentzFactor());
805 // Loop over all electrons of this hit
806 // TR photons produce hits with negative charge
807 Int_t nEl = ((Int_t) TMath::Abs(q));
808 for (Int_t iEl = 0; iEl < nEl; iEl++) {
814 // Electron attachment
815 if (fPar->ElAttachOn()) {
816 if (gRandom->Rndm() < (driftlengthL * elAttachProp))
820 // Apply the diffusion smearing
821 if (fPar->DiffusionOn()) {
822 if (!(fPar->Diffusion(driftlengthL,xyz))) continue;
825 // Apply E x B effects (depends on drift direction)
827 if (!(fPar->ExB(driftlength+kAmWidth,xyz))) continue;
830 // The electron position after diffusion and ExB in pad coordinates
831 // The pad row (z-direction)
832 Float_t rowDist = xyz[2] - row0;
833 Int_t rowE = ((Int_t) (rowDist * divideRow));
834 if ((rowE < 0) || (rowE >= nRowMax)) continue;
835 Float_t rowOffset = ((((Float_t) rowE) + 0.5) * rowPadSize) - rowDist;
837 // The pad column (rphi-direction)
838 Float_t col0tilt = fPar->Col0Tilted(col0,rowOffset,plane);
839 Float_t colDist = xyz[1] - col0tilt;
840 Int_t colE = ((Int_t) (colDist * divideCol));
841 if ((colE < 0) || (colE >= nColMax)) continue;
842 Float_t colOffset = ((((Float_t) colE) + 0.5) * colPadSize) - colDist;
844 // The time bin (negative for hits in the amplification region)
845 // In the amplification region the electrons drift from both sides
846 // to the middle (anode wire plane)
847 Float_t timeDist = time0 - xyz[0];
848 Float_t timeOffset = 0;
852 timeE = ((Int_t) (timeDist * divideTime));
853 // The distance of the position to the middle of the timebin
854 timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) - timeDist;
857 // Difference between half of the amplification gap width and
858 // the distance to the anode wire
859 Float_t anodeDist = kAmWidth - TMath::Abs(timeDist + kAmWidth);
861 timeE = -1 * (((Int_t ) (anodeDist * divideTime)) + 1);
862 // The distance of the position to the middle of the timebin
863 timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) + anodeDist;
866 // Apply the gas gain including fluctuations
867 Float_t ggRndm = 0.0;
869 ggRndm = gRandom->Rndm();
870 } while (ggRndm <= 0);
871 Int_t signal = (Int_t) (-fPar->GetGasGain() * TMath::Log(ggRndm));
873 // Apply the pad response
875 // The distance of the electron to the center of the pad
876 // in units of pad width
877 Float_t dist = - colOffset * divideCol;
878 if (!(fPar->PadResponse(signal,dist,plane,padSignal))) continue;
882 padSignal[1] = signal;
886 // Sample the time response inside the drift region
887 // + additional time bins before and after.
888 // The sampling is done always in the middle of the time bin
889 for (Int_t iTimeBin = TMath::Max(timeE-timeTRDbeg, -nTimeBefore)
890 ;iTimeBin < TMath::Min(timeE+timeTRDend,nTimeMax+nTimeAfter )
893 // Apply the time response
894 Float_t timeResponse = 1.0;
895 Float_t crossTalk = 0.0;
896 Float_t time = (iTimeBin - timeE) * timeBinSize + timeOffset;
898 timeResponse = fPar->TimeResponse(time);
901 crossTalk = fPar->CrossTalk(time);
908 for (iPad = 0; iPad < kNpad; iPad++) {
910 Int_t colPos = colE + iPad - 1;
911 if (colPos < 0) continue;
912 if (colPos >= nColMax) break;
915 // Note: The time bin number is shifted by nTimeBefore to avoid negative
916 // time bins. This has to be subtracted later.
917 Int_t iCurrentTimeBin = iTimeBin + nTimeBefore;
918 signalOld[iPad] = signals->GetDataUnchecked(rowE,colPos,iCurrentTimeBin);
919 if( colPos != colE ) {
920 signalOld[iPad] += padSignal[iPad] * (timeResponse + crossTalk);
923 signalOld[iPad] += padSignal[iPad] * timeResponse;
925 signals->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,signalOld[iPad]);
927 // Store the track index in the dictionary
928 // Note: We store index+1 in order to allow the array to be compressed
929 if (signalOld[iPad] > 0) {
930 for (iDict = 0; iDict < kNDict; iDict++) {
931 Int_t oldTrack = dictionary[iDict]->GetDataUnchecked(rowE
934 if (oldTrack == track+1) break;
936 dictionary[iDict]->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,track+1);
946 } // Loop: electrons of a single hit
948 } // If: detector and test hit
950 hit = (AliTRDhit *) fTRD->NextHit();
952 } // Loop: hits of one primary track
954 } // Loop: primary tracks
957 printf("<AliTRDdigitizer::MakeDigits> ");
958 printf("Finished analyzing %d hits\n",countHits);
961 // The total conversion factor
962 Float_t convert = kEl2fC * fPar->GetPadCoupling()
963 * fPar->GetTimeCoupling()
964 * fPar->GetChipGain();
966 // Loop through all chambers to finalize the digits
967 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
969 Int_t plane = fGeo->GetPlane(iDet);
970 Int_t sector = fGeo->GetSector(iDet);
971 Int_t chamber = fGeo->GetChamber(iDet);
972 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
973 Int_t nColMax = fPar->GetColMax(plane);
974 Int_t nTimeMax = fPar->GetTimeMax();
975 Int_t nTimeTotal = fPar->GetTimeTotal();
977 Double_t *inADC = new Double_t[nTimeTotal];
978 Double_t *outADC = new Double_t[nTimeTotal];
981 printf("<AliTRDdigitizer::MakeDigits> ");
982 printf("Digitization for chamber %d\n",iDet);
985 // Add a container for the digits of this detector
986 digits = fDigitsManager->GetDigits(iDet);
987 // Allocate memory space for the digits buffer
988 digits->Allocate(nRowMax,nColMax,nTimeTotal);
990 // Get the signal container
991 signals = (AliTRDdataArrayF *) signalsArray->At(iDet);
992 if (signals->GetNtime() == 0) {
993 // Create missing containers
994 signals->Allocate(nRowMax,nColMax,nTimeTotal);
997 // Expand the container if neccessary
998 if (fCompress) signals->Expand();
1000 // Create the missing dictionary containers
1001 for (iDict = 0; iDict < kNDict; iDict++) {
1002 dictionary[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1003 if (dictionary[iDict]->GetNtime() == 0) {
1004 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1010 // Don't create noise in detectors that are switched off
1011 if (CheckDetector(plane,chamber,sector)) {
1013 // Create the digits for this chamber
1014 for (iRow = 0; iRow < nRowMax; iRow++ ) {
1015 for (iCol = 0; iCol < nColMax; iCol++ ) {
1017 // Create summable digits
1020 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1021 Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1022 signalAmp *= fSDigitsScale;
1023 signalAmp = TMath::Min(signalAmp,(Float_t) 1.0e9);
1024 Int_t adc = (Int_t) signalAmp;
1025 if (adc > 0) nDigits++;
1026 digits->SetDataUnchecked(iRow,iCol,iTime,adc);
1030 // Create normal digits
1033 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1034 Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1036 signalAmp = TMath::Max((Double_t) gRandom->Gaus(signalAmp,fPar->GetNoise()),0.0);
1038 signalAmp *= convert;
1039 // Convert to ADC counts. Set the overflow-bit fADCoutRange if the
1040 // signal is larger than fADCinRange
1042 if (signalAmp >= fPar->GetADCinRange()) {
1043 adc = ((Int_t) fPar->GetADCoutRange());
1046 adc = ((Int_t) (signalAmp * (fPar->GetADCoutRange()
1047 / fPar->GetADCinRange())));
1050 outADC[iTime] = adc;
1053 // Apply the tail cancelation via the digital filter
1055 DeConvExp(inADC,outADC,nTimeTotal,fPar->GetTCnexp());
1058 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1059 // Store the amplitude of the digit if above threshold
1060 if (outADC[iTime] > fPar->GetADCthreshold()) {
1062 printf(" iRow = %d, iCol = %d, iTime = %d, adc = %f\n"
1063 ,iRow,iCol,iTime,outADC[iTime]);
1066 digits->SetDataUnchecked(iRow,iCol,iTime,outADC[iTime]);
1077 // Compress the arrays
1078 digits->Compress(1,0);
1079 for (iDict = 0; iDict < kNDict; iDict++) {
1080 dictionary[iDict]->Compress(1,0);
1083 totalSizeDigits += digits->GetSize();
1084 totalSizeDict0 += dictionary[0]->GetSize();
1085 totalSizeDict1 += dictionary[1]->GetSize();
1086 totalSizeDict2 += dictionary[2]->GetSize();
1088 Float_t nPixel = nRowMax * nColMax * nTimeMax;
1090 printf("<AliTRDdigitizer::MakeDigits> ");
1091 printf("Found %d digits in detector %d (%3.0f).\n"
1093 ,100.0 * ((Float_t) nDigits) / nPixel);
1096 if (fCompress) signals->Compress(1,0);
1104 printf("<AliTRDdigitizer::MakeDigits> ");
1105 printf("Total number of analyzed hits = %d\n",countHits);
1106 printf("<AliTRDdigitizer::MakeDigits> ");
1107 printf("Total digits data size = %d, %d, %d, %d\n",totalSizeDigits
1117 //_____________________________________________________________________________
1118 void AliTRDdigitizer::AddSDigitsManager(AliTRDdigitsManager *man)
1121 // Add a digits manager for s-digits to the input list.
1124 fSDigitsManagerList->Add(man);
1128 //_____________________________________________________________________________
1129 void AliTRDdigitizer::DeleteSDigitsManager()
1132 // Removes digits manager from the input list.
1135 fSDigitsManagerList->Delete();
1139 //_____________________________________________________________________________
1140 Bool_t AliTRDdigitizer::ConvertSDigits()
1143 // Converts s-digits to normal digits
1146 // Number of track dictionary arrays
1147 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1149 // Converts number of electrons to fC
1150 const Double_t kEl2fC = 1.602E-19 * 1.0E15;
1158 fPar = new AliTRDparameter("TRDparameter","Standard parameter");
1160 printf("<AliTRDdigitizer::ConvertSDigits> ");
1161 printf("Create the default parameter object\n");
1165 Double_t sDigitsScale = 1.0 / GetSDigitsScale();
1166 Double_t noise = fPar->GetNoise();
1167 Double_t padCoupling = fPar->GetPadCoupling();
1168 Double_t timeCoupling = fPar->GetTimeCoupling();
1169 Double_t chipGain = fPar->GetChipGain();
1170 Double_t convert = kEl2fC * padCoupling * timeCoupling * chipGain;;
1171 Double_t adcInRange = fPar->GetADCinRange();
1172 Double_t adcOutRange = fPar->GetADCoutRange();
1173 Int_t adcThreshold = fPar->GetADCthreshold();
1175 AliTRDdataArrayI *digitsIn;
1176 AliTRDdataArrayI *digitsOut;
1177 AliTRDdataArrayI *dictionaryIn[kNDict];
1178 AliTRDdataArrayI *dictionaryOut[kNDict];
1180 // Loop through the detectors
1181 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1184 printf("<AliTRDdigitizer::ConvertSDigits> ");
1185 printf("Convert detector %d to digits.\n",iDet);
1188 Int_t plane = fGeo->GetPlane(iDet);
1189 Int_t sector = fGeo->GetSector(iDet);
1190 Int_t chamber = fGeo->GetChamber(iDet);
1191 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1192 Int_t nColMax = fPar->GetColMax(plane);
1193 Int_t nTimeTotal = fPar->GetTimeTotal();
1195 Double_t *inADC = new Double_t[nTimeTotal];
1196 Double_t *outADC = new Double_t[nTimeTotal];
1198 digitsIn = fSDigitsManager->GetDigits(iDet);
1200 digitsOut = fDigitsManager->GetDigits(iDet);
1201 digitsOut->Allocate(nRowMax,nColMax,nTimeTotal);
1202 for (iDict = 0; iDict < kNDict; iDict++) {
1203 dictionaryIn[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1204 dictionaryIn[iDict]->Expand();
1205 dictionaryOut[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1206 dictionaryOut[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1209 for (iRow = 0; iRow < nRowMax; iRow++ ) {
1210 for (iCol = 0; iCol < nColMax; iCol++ ) {
1212 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1213 Double_t signal = (Double_t) digitsIn->GetDataUnchecked(iRow,iCol,iTime);
1214 signal *= sDigitsScale;
1216 signal = TMath::Max((Double_t) gRandom->Gaus(signal,noise),0.0);
1219 // Convert to ADC counts. Set the overflow-bit adcOutRange if the
1220 // signal is larger than adcInRange
1222 if (signal >= adcInRange) {
1223 adc = ((Int_t) adcOutRange);
1226 adc = ((Int_t) (signal * (adcOutRange / adcInRange)));
1229 outADC[iTime] = adc;
1232 // Apply the tail cancelation via the digital filter
1234 DeConvExp(inADC,outADC,nTimeTotal,fPar->GetTCnexp());
1237 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1238 // Store the amplitude of the digit if above threshold
1239 if (outADC[iTime] > adcThreshold) {
1240 digitsOut->SetDataUnchecked(iRow,iCol,iTime,outADC[iTime]);
1241 // Copy the dictionary
1242 for (iDict = 0; iDict < kNDict; iDict++) {
1243 Int_t track = dictionaryIn[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1244 dictionaryOut[iDict]->SetDataUnchecked(iRow,iCol,iTime,track);
1253 digitsIn->Compress(1,0);
1254 digitsOut->Compress(1,0);
1255 for (iDict = 0; iDict < kNDict; iDict++) {
1256 dictionaryIn[iDict]->Compress(1,0);
1257 dictionaryOut[iDict]->Compress(1,0);
1270 //_____________________________________________________________________________
1271 Bool_t AliTRDdigitizer::MergeSDigits()
1274 // Merges the input s-digits:
1275 // - The amplitude of the different inputs are summed up.
1276 // - Of the track IDs from the input dictionaries only one is
1277 // kept for each input. This works for maximal 3 different merged inputs.
1280 // Number of track dictionary arrays
1281 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1284 fPar = new AliTRDparameter("TRDparameter","Standard parameter");
1286 printf("<AliTRDdigitizer::MergeSDigits> ");
1287 printf("Create the default parameter object\n");
1294 AliTRDdataArrayI *digitsA;
1295 AliTRDdataArrayI *digitsB;
1296 AliTRDdataArrayI *dictionaryA[kNDict];
1297 AliTRDdataArrayI *dictionaryB[kNDict];
1299 // Get the first s-digits
1300 fSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->First();
1301 if (!fSDigitsManager) return kFALSE;
1303 // Loop through the other sets of s-digits
1304 AliTRDdigitsManager *mergeSDigitsManager;
1305 mergeSDigitsManager = (AliTRDdigitsManager *)
1306 fSDigitsManagerList->After(fSDigitsManager);
1309 if (mergeSDigitsManager) {
1310 printf("<AliTRDdigitizer::MergeSDigits> ");
1311 printf("Merge %d input files.\n",fSDigitsManagerList->GetSize());
1314 printf("<AliTRDdigitizer::MergeSDigits> ");
1315 printf("Only one input file.\n");
1320 while (mergeSDigitsManager) {
1324 // Loop through the detectors
1325 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1327 Int_t plane = fGeo->GetPlane(iDet);
1328 Int_t sector = fGeo->GetSector(iDet);
1329 Int_t chamber = fGeo->GetChamber(iDet);
1330 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1331 Int_t nColMax = fPar->GetColMax(plane);
1332 Int_t nTimeTotal = fPar->GetTimeTotal();
1334 // Loop through the pixels of one detector and add the signals
1335 digitsA = fSDigitsManager->GetDigits(iDet);
1336 digitsB = mergeSDigitsManager->GetDigits(iDet);
1339 for (iDict = 0; iDict < kNDict; iDict++) {
1340 dictionaryA[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1341 dictionaryB[iDict] = mergeSDigitsManager->GetDictionary(iDet,iDict);
1342 dictionaryA[iDict]->Expand();
1343 dictionaryB[iDict]->Expand();
1346 // Merge only detectors that contain a signal
1347 Bool_t doMerge = kTRUE;
1348 if (fMergeSignalOnly) {
1349 if (digitsA->GetOverThreshold(0) == 0) {
1357 printf("<AliTRDdigitizer::MergeSDigits> ");
1358 printf("Merge detector %d of input no.%d\n",iDet,iMerge+1);
1361 for (Int_t iRow = 0; iRow < nRowMax; iRow++ ) {
1362 for (Int_t iCol = 0; iCol < nColMax; iCol++ ) {
1363 for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {
1365 // Add the amplitudes of the summable digits
1366 Int_t ampA = digitsA->GetDataUnchecked(iRow,iCol,iTime);
1367 Int_t ampB = digitsB->GetDataUnchecked(iRow,iCol,iTime);
1369 digitsA->SetDataUnchecked(iRow,iCol,iTime,ampA);
1371 // Add the mask to the track id if defined.
1372 for (iDict = 0; iDict < kNDict; iDict++) {
1373 Int_t trackB = dictionaryB[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1374 if ((fMasks) && (trackB > 0)) {
1375 for (jDict = 0; jDict < kNDict; jDict++) {
1376 Int_t trackA = dictionaryA[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1378 trackA = trackB + fMasks[iMerge];
1379 dictionaryA[iDict]->SetDataUnchecked(iRow,iCol,iTime,trackA);
1392 digitsA->Compress(1,0);
1393 digitsB->Compress(1,0);
1394 for (iDict = 0; iDict < kNDict; iDict++) {
1395 dictionaryA[iDict]->Compress(1,0);
1396 dictionaryB[iDict]->Compress(1,0);
1402 // The next set of s-digits
1403 mergeSDigitsManager = (AliTRDdigitsManager *)
1404 fSDigitsManagerList->After(mergeSDigitsManager);
1412 //_____________________________________________________________________________
1413 Bool_t AliTRDdigitizer::SDigits2Digits()
1416 // Merges the input s-digits and converts them to normal digits
1419 if (!MergeSDigits()) return kFALSE;
1421 return ConvertSDigits();
1425 //_____________________________________________________________________________
1426 Bool_t AliTRDdigitizer::CheckDetector(Int_t plane, Int_t chamber, Int_t sector)
1429 // Checks whether a detector is enabled
1432 if ((fTRD->GetSensChamber() >= 0) &&
1433 (fTRD->GetSensChamber() != chamber)) return kFALSE;
1434 if ((fTRD->GetSensPlane() >= 0) &&
1435 (fTRD->GetSensPlane() != plane)) return kFALSE;
1436 if ( fTRD->GetSensSector() >= 0) {
1437 Int_t sens1 = fTRD->GetSensSector();
1438 Int_t sens2 = sens1 + fTRD->GetSensSectorRange();
1439 sens2 -= ((Int_t) (sens2 / AliTRDgeometry::Nsect()))
1440 * AliTRDgeometry::Nsect();
1441 if (sens1 < sens2) {
1442 if ((sector < sens1) || (sector >= sens2)) return kFALSE;
1445 if ((sector < sens1) && (sector >= sens2)) return kFALSE;
1453 //_____________________________________________________________________________
1454 Bool_t AliTRDdigitizer::WriteDigits() const
1457 // Writes out the TRD-digits and the dictionaries
1460 // Store the digits and the dictionary in the tree
1461 return fDigitsManager->WriteDigits();
1465 //_____________________________________________________________________________
1466 void AliTRDdigitizer::DeConvExp(Double_t *source, Double_t *target
1467 , Int_t n, Int_t nexp)
1470 // Does the deconvolution by the digital filter.
1472 // Author: Marcus Gutfleisch, KIP Heidelberg
1473 // Optimized for: New TRF from Venelin Angelov, simulated with CADENCE
1474 // Pad-ground capacitance = 25 pF
1475 // Pad-pad cross talk capacitance = 6 pF
1476 // For 10 MHz digitization, corresponding to 20 time bins
1477 // in the drift region
1481 Double_t coefficients[2];
1483 /* initialize (coefficient = alpha, rates = lambda) */
1486 rates[0] = 0.466998;
1488 coefficients[0] = 1.0;
1491 rates[0] = 0.8988162;
1492 coefficients[0] = 0.11392069;
1493 rates[1] = 0.3745688;
1494 coefficients[1] = 0.8860793;
1496 Float_t sumc = coefficients[0]+coefficients[1];
1497 coefficients[0] /= sumc;
1498 coefficients[1] /= sumc;
1502 Double_t reminder[2];
1503 Double_t correction, result;
1505 /* attention: computation order is important */
1507 for ( k = 0; k < nexp; k++ ) reminder[k]=0.0;
1509 for ( i = 0; i < n; i++ ) {
1510 result = ( source[i] - correction ); /* no rescaling */
1513 for ( k = 0; k < nexp; k++ ) reminder[k] = rates[k]
1514 * ( reminder[k] + coefficients[k] * result);
1517 for ( k = 0; k < nexp; k++ ) correction += reminder[k];
1522 //_____________________________________________________________________________
1523 void AliTRDdigitizer::InitOutput(TFile *file, Int_t iEvent)
1526 // Initializes the output branches
1530 fDigitsManager->MakeTreeAndBranches(file,iEvent);