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.33 2002/02/12 17:32:03 cblume
19 Rearrange the deleting of the list of sdigitsmanager
21 Revision 1.32 2002/02/12 16:07:21 cblume
24 Revision 1.31 2002/02/11 14:27:11 cblume
25 New pad plane design, new TRF+PRF, tail cancelation, cross talk
27 Revision 1.30 2001/11/19 08:44:08 cblume
28 Fix bugs reported by Rene
30 Revision 1.29 2001/11/14 19:44:25 hristov
31 Numeric const casted (Alpha)
33 Revision 1.28 2001/11/14 16:35:58 cblume
34 Inherits now from AliDetector
36 Revision 1.27 2001/11/14 10:50:45 cblume
37 Changes in digits IO. Add merging of summable digits
39 Revision 1.26 2001/11/06 17:19:41 cblume
40 Add detailed geometry and simple simulator
42 Revision 1.25 2001/06/27 09:54:44 cblume
43 Moved fField initialization to InitDetector()
45 Revision 1.24 2001/05/21 16:45:47 hristov
46 Last minute changes (C.Blume)
48 Revision 1.23 2001/05/07 08:04:48 cblume
49 New TRF and PRF. Speedup of the code. Digits from amplification region included
51 Revision 1.22 2001/03/30 14:40:14 cblume
52 Update of the digitization parameter
54 Revision 1.21 2001/03/13 09:30:35 cblume
55 Update of digitization. Moved digit branch definition to AliTRD
57 Revision 1.20 2001/02/25 20:19:00 hristov
58 Minor correction: loop variable declared only once for HP, Sun
60 Revision 1.19 2001/02/14 18:22:26 cblume
61 Change in the geometry of the padplane
63 Revision 1.18 2001/01/26 19:56:57 hristov
64 Major upgrade of AliRoot code
66 Revision 1.17 2000/12/08 12:53:27 cblume
67 Change in Copy() function for HP-compiler
69 Revision 1.16 2000/12/07 12:20:46 cblume
70 Go back to array compression. Use sampled PRF to speed up digitization
72 Revision 1.15 2000/11/23 14:34:08 cblume
73 Fixed bug in expansion routine of arrays (initialize buffers properly)
75 Revision 1.14 2000/11/20 08:54:44 cblume
76 Switch off compression as default
78 Revision 1.13 2000/11/10 14:57:52 cblume
79 Changes in the geometry constants for the DEC compiler
81 Revision 1.12 2000/11/01 14:53:20 cblume
82 Merge with TRD-develop
84 Revision 1.1.4.9 2000/10/26 17:00:22 cblume
85 Fixed bug in CheckDetector()
87 Revision 1.1.4.8 2000/10/23 13:41:35 cblume
88 Added protection against Log(0) in the gas gain calulation
90 Revision 1.1.4.7 2000/10/17 02:27:34 cblume
91 Get rid of global constants
93 Revision 1.1.4.6 2000/10/16 01:16:53 cblume
94 Changed timebin 0 to be the one closest to the readout
96 Revision 1.1.4.5 2000/10/15 23:34:29 cblume
97 Faster version of the digitizer
99 Revision 1.1.4.4 2000/10/06 16:49:46 cblume
102 Revision 1.1.4.3 2000/10/04 16:34:58 cblume
103 Replace include files by forward declarations
105 Revision 1.1.4.2 2000/09/22 14:41:10 cblume
106 Bug fix in PRF. Included time response. New structure
108 Revision 1.10 2000/10/05 07:27:53 cblume
109 Changes in the header-files by FCA
111 Revision 1.9 2000/10/02 21:28:19 fca
112 Removal of useless dependecies via forward declarations
114 Revision 1.8 2000/06/09 11:10:07 cblume
115 Compiler warnings and coding conventions, next round
117 Revision 1.7 2000/06/08 18:32:58 cblume
118 Make code compliant to coding conventions
120 Revision 1.6 2000/06/07 16:27:32 cblume
121 Try to remove compiler warnings on Sun and HP
123 Revision 1.5 2000/05/09 16:38:57 cblume
124 Removed PadResponse(). Merge problem
126 Revision 1.4 2000/05/08 15:53:45 cblume
127 Resolved merge conflict
129 Revision 1.3 2000/04/28 14:49:27 cblume
130 Only one declaration of iDict in MakeDigits()
132 Revision 1.1.4.1 2000/05/08 14:42:04 cblume
133 Introduced AliTRDdigitsManager
135 Revision 1.1 2000/02/28 19:00:13 cblume
140 ///////////////////////////////////////////////////////////////////////////////
142 // Creates and handles digits from TRD hits //
143 // Author: C. Blume (C.Blume@gsi.de) //
145 // The following effects are included: //
148 // - Gas gain including fluctuations //
149 // - Pad-response (simple Gaussian approximation) //
150 // - Time-response //
151 // - Electronics noise //
152 // - Electronics gain //
154 // - ADC threshold //
155 // The corresponding parameter can be adjusted via the various //
156 // Set-functions. If these parameters are not explicitly set, default //
157 // values are used (see Init-function). //
158 // As an example on how to use this class to produce digits from hits //
159 // have a look at the macro hits2digits.C //
160 // The production of summable digits is demonstrated in hits2sdigits.C //
161 // and the subsequent conversion of the s-digits into normal digits is //
162 // explained in sdigits2digits.C. //
164 ///////////////////////////////////////////////////////////////////////////////
180 #include "AliRunDigitizer.h"
183 #include "AliTRDhit.h"
184 #include "AliTRDdigitizer.h"
185 #include "AliTRDdataArrayI.h"
186 #include "AliTRDdataArrayF.h"
187 #include "AliTRDsegmentArray.h"
188 #include "AliTRDdigitsManager.h"
189 #include "AliTRDgeometry.h"
190 #include "AliTRDparameter.h"
192 ClassImp(AliTRDdigitizer)
194 //_____________________________________________________________________________
195 AliTRDdigitizer::AliTRDdigitizer()
198 // AliTRDdigitizer default constructor
203 fSDigitsManagerList = 0;
214 fMergeSignalOnly = kFALSE;
218 //_____________________________________________________________________________
219 AliTRDdigitizer::AliTRDdigitizer(const Text_t *name, const Text_t *title)
220 :AliDigitizer(name,title)
223 // AliTRDdigitizer constructor
228 fSDigitsManagerList = 0;
238 fMergeSignalOnly = kFALSE;
240 // For the summable digits
241 fSDigitsScale = 100.;
245 //_____________________________________________________________________________
246 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager
247 , const Text_t *name, const Text_t *title)
248 :AliDigitizer(manager,name,title)
251 // AliTRDdigitizer constructor
256 fSDigitsManagerList = 0;
266 fMergeSignalOnly = kFALSE;
268 // For the summable digits
269 fSDigitsScale = 100.;
273 //_____________________________________________________________________________
274 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager)
275 :AliDigitizer(manager,"AliTRDdigitizer","TRD digitizer")
278 // AliTRDdigitizer constructor
283 fSDigitsManagerList = 0;
293 fMergeSignalOnly = kFALSE;
295 // For the summable digits
296 fSDigitsScale = 100.;
300 //_____________________________________________________________________________
301 AliTRDdigitizer::AliTRDdigitizer(const AliTRDdigitizer &d)
304 // AliTRDdigitizer copy constructor
307 ((AliTRDdigitizer &) d).Copy(*this);
311 //_____________________________________________________________________________
312 AliTRDdigitizer::~AliTRDdigitizer()
315 // AliTRDdigitizer destructor
324 if (fDigitsManager) {
325 delete fDigitsManager;
329 if (fSDigitsManager) {
330 delete fSDigitsManager;
334 if (fSDigitsManagerList) {
335 delete fSDigitsManagerList;
336 fSDigitsManagerList = 0;
346 //_____________________________________________________________________________
347 AliTRDdigitizer &AliTRDdigitizer::operator=(const AliTRDdigitizer &d)
350 // Assignment operator
353 if (this != &d) ((AliTRDdigitizer &) d).Copy(*this);
358 //_____________________________________________________________________________
359 void AliTRDdigitizer::Copy(TObject &d)
365 ((AliTRDdigitizer &) d).fInputFile = 0;
366 ((AliTRDdigitizer &) d).fSDigitsManagerList = 0;
367 ((AliTRDdigitizer &) d).fSDigitsManager = 0;
368 ((AliTRDdigitizer &) d).fDigitsManager = 0;
369 ((AliTRDdigitizer &) d).fTRD = 0;
370 ((AliTRDdigitizer &) d).fGeo = 0;
371 ((AliTRDdigitizer &) d).fMasks = 0;
372 ((AliTRDdigitizer &) d).fEvent = 0;
373 ((AliTRDdigitizer &) d).fPar = 0;
374 ((AliTRDdigitizer &) d).fCompress = fCompress;
375 ((AliTRDdigitizer &) d).fDebug = fDebug ;
376 ((AliTRDdigitizer &) d).fSDigits = fSDigits;
377 ((AliTRDdigitizer &) d).fSDigitsScale = fSDigitsScale;
378 ((AliTRDdigitizer &) d).fMergeSignalOnly = fMergeSignalOnly;
382 //_____________________________________________________________________________
383 void AliTRDdigitizer::Exec(Option_t* option)
386 // Executes the merging
391 AliTRDdigitsManager *sdigitsManager;
393 TString optionString = option;
394 if (optionString.Contains("deb")) {
396 if (optionString.Contains("2")) {
399 printf("<AliTRDdigitizer::Exec> ");
400 printf("Called with debug option %d\n",fDebug);
403 // Connect the AliRoot file containing Geometry, Kine, and Hits
404 fInputFile = (TFile *) fManager->GetInputTreeTRDS(0)->GetCurrentFile();
407 printf("<AliTRDdigitizer::Exec> ");
408 printf("Cannot open the input file %s.\n",fInputFile->GetName());
416 gAlice = (AliRun *) fInputFile->Get("gAlice");
419 printf("<AliTRDdigitizer::Exec> ");
420 printf("AliRun object found on file.\n");
424 printf("<AliTRDdigitizer::Exec> ");
425 printf("Could not find AliRun object.\n");
428 Int_t nInput = fManager->GetNinputs();
429 fMasks = new Int_t[nInput];
430 for (iInput = 0; iInput < nInput; iInput++) {
431 fMasks[iInput] = fManager->GetMask(iInput);
437 for (iInput = 0; iInput < nInput; iInput++) {
440 printf("<AliTRDdigitizer::Exec> ");
441 printf("Add input stream %d\n",iInput);
444 // Read the s-digits via digits manager
445 sdigitsManager = new AliTRDdigitsManager();
446 sdigitsManager->SetDebug(fDebug);
447 sdigitsManager->SetSDigits(kTRUE);
448 sdigitsManager->ReadDigits(fManager->GetInputTreeTRDS(iInput));
450 // Add the s-digits to the input list
451 AddSDigitsManager(sdigitsManager);
455 // Convert the s-digits to normal digits
457 printf("<AliTRDdigitizer::Exec> ");
458 printf("Do the conversion\n");
464 printf("<AliTRDdigitizer::Exec> ");
465 printf("Write the digits\n");
467 fDigitsManager->MakeBranch(fManager->GetTreeDTRD());
468 fDigitsManager->WriteDigits();
470 printf("<AliTRDdigitizer::Exec> ");
474 DeleteSDigitsManager();
478 //_____________________________________________________________________________
479 Bool_t AliTRDdigitizer::Open(const Char_t *file, Int_t nEvent)
482 // Opens a ROOT-file with TRD-hits and reads in the hit-tree
485 // Connect the AliRoot file containing Geometry, Kine, and Hits
486 fInputFile = (TFile *) gROOT->GetListOfFiles()->FindObject(file);
489 printf("<AliTRDdigitizer::Open> ");
490 printf("Open the AliROOT-file %s.\n",file);
492 fInputFile = new TFile(file,"UPDATE");
496 printf("<AliTRDdigitizer::Open> ");
497 printf("%s is already open.\n",file);
501 gAlice = (AliRun *) fInputFile->Get("gAlice");
504 printf("<AliTRDdigitizer::Open> ");
505 printf("AliRun object found on file.\n");
509 printf("<AliTRDdigitizer::Open> ");
510 printf("Could not find AliRun object.\n");
516 // Import the Trees for the event nEvent in the file
517 Int_t nparticles = gAlice->GetEvent(fEvent);
518 if (nparticles <= 0) {
519 printf("<AliTRDdigitizer::Open> ");
520 printf("No entries in the trees for event %d.\n",fEvent);
524 if (InitDetector()) {
533 //_____________________________________________________________________________
534 Bool_t AliTRDdigitizer::InitDetector()
537 // Sets the pointer to the TRD detector and the geometry
540 // Get the pointer to the detector class and check for version 1
541 fTRD = (AliTRD *) gAlice->GetDetector("TRD");
542 if (fTRD->IsVersion() != 1) {
543 printf("<AliTRDdigitizer::InitDetector> ");
544 printf("TRD must be version 1 (slow simulator).\n");
549 fGeo = fTRD->GetGeometry();
551 printf("<AliTRDdigitizer::InitDetector> ");
552 printf("Geometry version %d\n",fGeo->IsVersion());
555 // Create a digits manager
556 fDigitsManager = new AliTRDdigitsManager();
557 fDigitsManager->SetSDigits(fSDigits);
558 fDigitsManager->CreateArrays();
559 fDigitsManager->SetEvent(fEvent);
560 fDigitsManager->SetDebug(fDebug);
562 // The list for the input s-digits manager to be merged
563 fSDigitsManagerList = new TList();
569 //_____________________________________________________________________________
570 Bool_t AliTRDdigitizer::MakeBranch(const Char_t *file)
573 // Create the branches for the digits array
576 return fDigitsManager->MakeBranch(file);
580 //_____________________________________________________________________________
581 Bool_t AliTRDdigitizer::MakeDigits()
587 ///////////////////////////////////////////////////////////////
589 ///////////////////////////////////////////////////////////////
591 // Converts number of electrons to fC
592 const Double_t kEl2fC = 1.602E-19 * 1.0E15;
594 ///////////////////////////////////////////////////////////////
596 // Number of pads included in the pad response
597 const Int_t kNpad = 3;
599 // Number of track dictionary arrays
600 const Int_t kNDict = AliTRDdigitsManager::kNDict;
602 // Half the width of the amplification region
603 const Float_t kAmWidth = AliTRDgeometry::AmThick() / 2.;
605 Int_t iRow, iCol, iTime, iPad;
609 Int_t totalSizeDigits = 0;
610 Int_t totalSizeDict0 = 0;
611 Int_t totalSizeDict1 = 0;
612 Int_t totalSizeDict2 = 0;
614 Int_t timeTRDbeg = 0;
615 Int_t timeTRDend = 1;
620 Float_t padSignal[kNpad];
621 Float_t signalOld[kNpad];
623 AliTRDdataArrayF *signals = 0;
624 AliTRDdataArrayI *digits = 0;
625 AliTRDdataArrayI *dictionary[kNDict];
627 // Create a default parameter class if none is defined
629 fPar = new AliTRDparameter("TRDparameter","Standard TRD parameter");
631 printf("<AliTRDdigitizer::MakeDigits> ");
632 printf("Create the default parameter object\n");
636 // Create a container for the amplitudes
637 AliTRDsegmentArray *signalsArray
638 = new AliTRDsegmentArray("AliTRDdataArrayF"
639 ,AliTRDgeometry::Ndet());
642 timeTRDbeg = ((Int_t) (-fPar->GetTRFlo() / fPar->GetTimeBinSize())) - 1;
643 timeTRDend = ((Int_t) ( fPar->GetTRFhi() / fPar->GetTimeBinSize())) - 1;
645 printf("<AliTRDdigitizer::MakeDigits> ");
646 printf("Sample the TRF between -%d and %d\n",timeTRDbeg,timeTRDend);
650 Float_t elAttachProp = fPar->GetElAttachProp() / 100.;
653 printf("<AliTRDdigitizer::MakeDigits> ");
654 printf("No geometry defined\n");
659 printf("<AliTRDdigitizer::MakeDigits> ");
660 printf("Start creating digits.\n");
663 // Get the pointer to the hit tree
664 TTree *HitTree = gAlice->TreeH();
666 // Get the number of entries in the hit tree
667 // (Number of primary particles creating a hit somewhere)
668 Int_t nTrack = (Int_t) HitTree->GetEntries();
670 printf("<AliTRDdigitizer::MakeDigits> ");
671 printf("Found %d primary particles\n",nTrack);
674 Int_t detectorOld = -1;
677 // Loop through all entries in the tree
678 for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
681 nBytes += HitTree->GetEvent(iTrack);
683 // Loop through the TRD hits
685 AliTRDhit *hit = (AliTRDhit *) fTRD->FirstHit(-1);
694 Float_t q = hit->GetCharge();
695 Int_t track = hit->Track();
696 Int_t detector = hit->GetDetector();
697 Int_t plane = fGeo->GetPlane(detector);
698 Int_t sector = fGeo->GetSector(detector);
699 Int_t chamber = fGeo->GetChamber(detector);
700 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
701 Int_t nColMax = fPar->GetColMax(plane);
702 Int_t nTimeMax = fPar->GetTimeMax();
703 Int_t nTimeBefore = fPar->GetTimeBefore();
704 Int_t nTimeAfter = fPar->GetTimeAfter();
705 Int_t nTimeTotal = fPar->GetTimeTotal();
706 Float_t row0 = fPar->GetRow0(plane,chamber,sector);
707 Float_t col0 = fPar->GetCol0(plane);
708 Float_t time0 = fPar->GetTime0(plane);
709 Float_t rowPadSize = fPar->GetRowPadSize(plane,chamber,sector);
710 Float_t colPadSize = fPar->GetColPadSize(plane);
711 Float_t timeBinSize = fPar->GetTimeBinSize();
712 Float_t divideRow = 1.0 / rowPadSize;
713 Float_t divideCol = 1.0 / colPadSize;
714 Float_t divideTime = 1.0 / timeBinSize;
717 printf("Analyze hit no. %d ",iHit);
718 printf("-----------------------------------------------------------\n");
720 printf("plane = %d, sector = %d, chamber = %d\n"
721 ,plane,sector,chamber);
722 printf("nRowMax = %d, nColMax = %d, nTimeMax = %d\n"
723 ,nRowMax,nColMax,nTimeMax);
724 printf("nTimeBefore = %d, nTimeAfter = %d, nTimeTotal = %d\n"
725 ,nTimeBefore,nTimeAfter,nTimeTotal);
726 printf("row0 = %f, col0 = %f, time0 = %f\n"
728 printf("rowPadSize = %f, colPadSize = %f, timeBinSize = %f\n"
729 ,rowPadSize,colPadSize,timeBinSize);
732 // Don't analyze test hits and switched off detectors
733 if ((CheckDetector(plane,chamber,sector)) &&
734 (((Int_t) q) != 0)) {
736 if (detector != detectorOld) {
739 printf("<AliTRDdigitizer::MakeDigits> ");
740 printf("Get new container. New det = %d, Old det = %d\n"
741 ,detector,detectorOld);
743 // Compress the old one if enabled
744 if ((fCompress) && (detectorOld > -1)) {
746 printf("<AliTRDdigitizer::MakeDigits> ");
747 printf("Compress the old container ...");
749 signals->Compress(1,0);
750 for (iDict = 0; iDict < kNDict; iDict++) {
751 dictionary[iDict]->Compress(1,0);
753 if (fDebug > 1) printf("done\n");
755 // Get the new container
756 signals = (AliTRDdataArrayF *) signalsArray->At(detector);
757 if (signals->GetNtime() == 0) {
758 // Allocate a new one if not yet existing
760 printf("<AliTRDdigitizer::MakeDigits> ");
761 printf("Allocate a new container ... ");
763 signals->Allocate(nRowMax,nColMax,nTimeTotal);
766 // Expand an existing one
769 printf("<AliTRDdigitizer::MakeDigits> ");
770 printf("Expand an existing container ... ");
775 // The same for the dictionary
776 for (iDict = 0; iDict < kNDict; iDict++) {
777 dictionary[iDict] = fDigitsManager->GetDictionary(detector,iDict);
778 if (dictionary[iDict]->GetNtime() == 0) {
779 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
782 if (fCompress) dictionary[iDict]->Expand();
785 if (fDebug > 1) printf("done\n");
786 detectorOld = detector;
789 // Rotate the sectors on top of each other
790 fGeo->Rotate(detector,pos,rot);
792 // The driftlength. It is negative if the hit is in the
793 // amplification region.
794 Float_t driftlength = time0 - rot[0];
796 // Take also the drift in the amplification region into account
797 // The drift length is at the moment still the same, regardless of
798 // the position relativ to the wire. This non-isochronity needs still
799 // to be implemented.
800 Float_t driftlengthL = TMath::Abs(driftlength + kAmWidth);
801 if (fPar->ExBOn()) driftlengthL /= TMath::Sqrt(fPar->GetLorentzFactor());
803 // Loop over all electrons of this hit
804 // TR photons produce hits with negative charge
805 Int_t nEl = ((Int_t) TMath::Abs(q));
806 for (Int_t iEl = 0; iEl < nEl; iEl++) {
812 // Electron attachment
813 if (fPar->ElAttachOn()) {
814 if (gRandom->Rndm() < (driftlengthL * elAttachProp))
818 // Apply the diffusion smearing
819 if (fPar->DiffusionOn()) {
820 if (!(fPar->Diffusion(driftlengthL,xyz))) continue;
823 // Apply E x B effects (depends on drift direction)
825 if (!(fPar->ExB(driftlength+kAmWidth,xyz))) continue;
828 // The electron position after diffusion and ExB in pad coordinates
829 // The pad row (z-direction)
830 Float_t rowDist = xyz[2] - row0;
831 Int_t rowE = ((Int_t) (rowDist * divideRow));
832 if ((rowE < 0) || (rowE >= nRowMax)) continue;
833 Float_t rowOffset = ((((Float_t) rowE) + 0.5) * rowPadSize) - rowDist;
835 // The pad column (rphi-direction)
836 Float_t col0tilt = fPar->Col0Tilted(col0,rowOffset,plane);
837 Float_t colDist = xyz[1] - col0tilt;
838 Int_t colE = ((Int_t) (colDist * divideCol));
839 if ((colE < 0) || (colE >= nColMax)) continue;
840 Float_t colOffset = ((((Float_t) colE) + 0.5) * colPadSize) - colDist;
842 // The time bin (negative for hits in the amplification region)
843 // In the amplification region the electrons drift from both sides
844 // to the middle (anode wire plane)
845 Float_t timeDist = time0 - xyz[0];
846 Float_t timeOffset = 0;
850 timeE = ((Int_t) (timeDist * divideTime));
851 // The distance of the position to the middle of the timebin
852 timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) - timeDist;
855 // Difference between half of the amplification gap width and
856 // the distance to the anode wire
857 Float_t anodeDist = kAmWidth - TMath::Abs(timeDist + kAmWidth);
859 timeE = -1 * (((Int_t ) (anodeDist * divideTime)) + 1);
860 // The distance of the position to the middle of the timebin
861 timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) + anodeDist;
864 // Apply the gas gain including fluctuations
865 Float_t ggRndm = 0.0;
867 ggRndm = gRandom->Rndm();
868 } while (ggRndm <= 0);
869 Int_t signal = (Int_t) (-fPar->GetGasGain() * TMath::Log(ggRndm));
871 // Apply the pad response
873 // The distance of the electron to the center of the pad
874 // in units of pad width
875 Float_t dist = - colOffset * divideCol;
876 if (!(fPar->PadResponse(signal,dist,plane,padSignal))) continue;
880 padSignal[1] = signal;
884 // Sample the time response inside the drift region
885 // + additional time bins before and after.
886 // The sampling is done always in the middle of the time bin
887 for (Int_t iTimeBin = TMath::Max(timeE-timeTRDbeg, -nTimeBefore)
888 ;iTimeBin < TMath::Min(timeE+timeTRDend,nTimeMax+nTimeAfter )
891 // Apply the time response
892 Float_t timeResponse = 1.0;
893 Float_t crossTalk = 0.0;
894 Float_t time = (iTimeBin - timeE) * timeBinSize + timeOffset;
896 timeResponse = fPar->TimeResponse(time);
899 crossTalk = fPar->CrossTalk(time);
906 for (iPad = 0; iPad < kNpad; iPad++) {
908 Int_t colPos = colE + iPad - 1;
909 if (colPos < 0) continue;
910 if (colPos >= nColMax) break;
913 // Note: The time bin number is shifted by nTimeBefore to avoid negative
914 // time bins. This has to be subtracted later.
915 Int_t iCurrentTimeBin = iTimeBin + nTimeBefore;
916 signalOld[iPad] = signals->GetDataUnchecked(rowE,colPos,iCurrentTimeBin);
917 if( colPos != colE ) {
918 signalOld[iPad] += padSignal[iPad] * (timeResponse + crossTalk);
921 signalOld[iPad] += padSignal[iPad] * timeResponse;
923 signals->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,signalOld[iPad]);
925 // Store the track index in the dictionary
926 // Note: We store index+1 in order to allow the array to be compressed
927 if (signalOld[iPad] > 0) {
928 for (iDict = 0; iDict < kNDict; iDict++) {
929 Int_t oldTrack = dictionary[iDict]->GetDataUnchecked(rowE
932 if (oldTrack == track+1) break;
934 dictionary[iDict]->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,track+1);
944 } // Loop: electrons of a single hit
946 } // If: detector and test hit
948 hit = (AliTRDhit *) fTRD->NextHit();
950 } // Loop: hits of one primary track
952 } // Loop: primary tracks
955 printf("<AliTRDdigitizer::MakeDigits> ");
956 printf("Finished analyzing %d hits\n",countHits);
959 // The total conversion factor
960 Float_t convert = kEl2fC * fPar->GetPadCoupling()
961 * fPar->GetTimeCoupling()
962 * fPar->GetChipGain();
964 // Loop through all chambers to finalize the digits
965 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
967 Int_t plane = fGeo->GetPlane(iDet);
968 Int_t sector = fGeo->GetSector(iDet);
969 Int_t chamber = fGeo->GetChamber(iDet);
970 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
971 Int_t nColMax = fPar->GetColMax(plane);
972 Int_t nTimeMax = fPar->GetTimeMax();
973 Int_t nTimeTotal = fPar->GetTimeTotal();
975 Double_t *inADC = new Double_t[nTimeTotal];
976 Double_t *outADC = new Double_t[nTimeTotal];
979 printf("<AliTRDdigitizer::MakeDigits> ");
980 printf("Digitization for chamber %d\n",iDet);
983 // Add a container for the digits of this detector
984 digits = fDigitsManager->GetDigits(iDet);
985 // Allocate memory space for the digits buffer
986 digits->Allocate(nRowMax,nColMax,nTimeTotal);
988 // Get the signal container
989 signals = (AliTRDdataArrayF *) signalsArray->At(iDet);
990 if (signals->GetNtime() == 0) {
991 // Create missing containers
992 signals->Allocate(nRowMax,nColMax,nTimeTotal);
995 // Expand the container if neccessary
996 if (fCompress) signals->Expand();
998 // Create the missing dictionary containers
999 for (iDict = 0; iDict < kNDict; iDict++) {
1000 dictionary[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1001 if (dictionary[iDict]->GetNtime() == 0) {
1002 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1008 // Don't create noise in detectors that are switched off
1009 if (CheckDetector(plane,chamber,sector)) {
1011 // Create the digits for this chamber
1012 for (iRow = 0; iRow < nRowMax; iRow++ ) {
1013 for (iCol = 0; iCol < nColMax; iCol++ ) {
1015 // Create summable digits
1018 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1019 Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1020 signalAmp *= fSDigitsScale;
1021 signalAmp = TMath::Min(signalAmp,(Float_t) 1.0e9);
1022 Int_t adc = (Int_t) signalAmp;
1023 if (adc > 0) nDigits++;
1024 digits->SetDataUnchecked(iRow,iCol,iTime,adc);
1028 // Create normal digits
1031 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1032 Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1034 signalAmp = TMath::Max((Double_t) gRandom->Gaus(signalAmp,fPar->GetNoise()),0.0);
1036 signalAmp *= convert;
1037 // Convert to ADC counts. Set the overflow-bit fADCoutRange if the
1038 // signal is larger than fADCinRange
1040 if (signalAmp >= fPar->GetADCinRange()) {
1041 adc = ((Int_t) fPar->GetADCoutRange());
1044 adc = ((Int_t) (signalAmp * (fPar->GetADCoutRange()
1045 / fPar->GetADCinRange())));
1048 outADC[iTime] = adc;
1051 // Apply the tail cancelation via the digital filter
1053 DeConvExp(inADC,outADC,nTimeTotal,fPar->GetTCnexp());
1056 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1057 // Store the amplitude of the digit if above threshold
1058 if (outADC[iTime] > fPar->GetADCthreshold()) {
1060 printf(" iRow = %d, iCol = %d, iTime = %d, adc = %f\n"
1061 ,iRow,iCol,iTime,outADC[iTime]);
1064 digits->SetDataUnchecked(iRow,iCol,iTime,outADC[iTime]);
1075 // Compress the arrays
1076 digits->Compress(1,0);
1077 for (iDict = 0; iDict < kNDict; iDict++) {
1078 dictionary[iDict]->Compress(1,0);
1081 totalSizeDigits += digits->GetSize();
1082 totalSizeDict0 += dictionary[0]->GetSize();
1083 totalSizeDict1 += dictionary[1]->GetSize();
1084 totalSizeDict2 += dictionary[2]->GetSize();
1086 Float_t nPixel = nRowMax * nColMax * nTimeMax;
1088 printf("<AliTRDdigitizer::MakeDigits> ");
1089 printf("Found %d digits in detector %d (%3.0f).\n"
1091 ,100.0 * ((Float_t) nDigits) / nPixel);
1094 if (fCompress) signals->Compress(1,0);
1102 printf("<AliTRDdigitizer::MakeDigits> ");
1103 printf("Total number of analyzed hits = %d\n",countHits);
1104 printf("<AliTRDdigitizer::MakeDigits> ");
1105 printf("Total digits data size = %d, %d, %d, %d\n",totalSizeDigits
1115 //_____________________________________________________________________________
1116 void AliTRDdigitizer::AddSDigitsManager(AliTRDdigitsManager *man)
1119 // Add a digits manager for s-digits to the input list.
1122 fSDigitsManagerList->Add(man);
1126 //_____________________________________________________________________________
1127 void AliTRDdigitizer::DeleteSDigitsManager()
1130 // Removes digits manager from the input list.
1133 fSDigitsManagerList->Delete();
1137 //_____________________________________________________________________________
1138 Bool_t AliTRDdigitizer::ConvertSDigits()
1141 // Converts s-digits to normal digits
1144 // Number of track dictionary arrays
1145 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1147 // Converts number of electrons to fC
1148 const Double_t kEl2fC = 1.602E-19 * 1.0E15;
1156 fPar = new AliTRDparameter("TRDparameter","Standard parameter");
1158 printf("<AliTRDdigitizer::ConvertSDigits> ");
1159 printf("Create the default parameter object\n");
1163 Double_t sDigitsScale = 1.0 / GetSDigitsScale();
1164 Double_t noise = fPar->GetNoise();
1165 Double_t padCoupling = fPar->GetPadCoupling();
1166 Double_t timeCoupling = fPar->GetTimeCoupling();
1167 Double_t chipGain = fPar->GetChipGain();
1168 Double_t convert = kEl2fC * padCoupling * timeCoupling * chipGain;;
1169 Double_t adcInRange = fPar->GetADCinRange();
1170 Double_t adcOutRange = fPar->GetADCoutRange();
1171 Int_t adcThreshold = fPar->GetADCthreshold();
1173 AliTRDdataArrayI *digitsIn;
1174 AliTRDdataArrayI *digitsOut;
1175 AliTRDdataArrayI *dictionaryIn[kNDict];
1176 AliTRDdataArrayI *dictionaryOut[kNDict];
1178 // Loop through the detectors
1179 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1182 printf("<AliTRDdigitizer::ConvertSDigits> ");
1183 printf("Convert detector %d to digits.\n",iDet);
1186 Int_t plane = fGeo->GetPlane(iDet);
1187 Int_t sector = fGeo->GetSector(iDet);
1188 Int_t chamber = fGeo->GetChamber(iDet);
1189 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1190 Int_t nColMax = fPar->GetColMax(plane);
1191 Int_t nTimeTotal = fPar->GetTimeTotal();
1193 Double_t *inADC = new Double_t[nTimeTotal];
1194 Double_t *outADC = new Double_t[nTimeTotal];
1196 digitsIn = fSDigitsManager->GetDigits(iDet);
1198 digitsOut = fDigitsManager->GetDigits(iDet);
1199 digitsOut->Allocate(nRowMax,nColMax,nTimeTotal);
1200 for (iDict = 0; iDict < kNDict; iDict++) {
1201 dictionaryIn[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1202 dictionaryIn[iDict]->Expand();
1203 dictionaryOut[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1204 dictionaryOut[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1207 for (iRow = 0; iRow < nRowMax; iRow++ ) {
1208 for (iCol = 0; iCol < nColMax; iCol++ ) {
1210 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1211 Double_t signal = (Double_t) digitsIn->GetDataUnchecked(iRow,iCol,iTime);
1212 signal *= sDigitsScale;
1214 signal = TMath::Max((Double_t) gRandom->Gaus(signal,noise),0.0);
1217 // Convert to ADC counts. Set the overflow-bit adcOutRange if the
1218 // signal is larger than adcInRange
1220 if (signal >= adcInRange) {
1221 adc = ((Int_t) adcOutRange);
1224 adc = ((Int_t) (signal * (adcOutRange / adcInRange)));
1227 outADC[iTime] = adc;
1230 // Apply the tail cancelation via the digital filter
1232 DeConvExp(inADC,outADC,nTimeTotal,fPar->GetTCnexp());
1235 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1236 // Store the amplitude of the digit if above threshold
1237 if (outADC[iTime] > adcThreshold) {
1238 digitsOut->SetDataUnchecked(iRow,iCol,iTime,outADC[iTime]);
1239 // Copy the dictionary
1240 for (iDict = 0; iDict < kNDict; iDict++) {
1241 Int_t track = dictionaryIn[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1242 dictionaryOut[iDict]->SetDataUnchecked(iRow,iCol,iTime,track);
1251 digitsIn->Compress(1,0);
1252 digitsOut->Compress(1,0);
1253 for (iDict = 0; iDict < kNDict; iDict++) {
1254 dictionaryIn[iDict]->Compress(1,0);
1255 dictionaryOut[iDict]->Compress(1,0);
1268 //_____________________________________________________________________________
1269 Bool_t AliTRDdigitizer::MergeSDigits()
1272 // Merges the input s-digits:
1273 // - The amplitude of the different inputs are summed up.
1274 // - Of the track IDs from the input dictionaries only one is
1275 // kept for each input. This works for maximal 3 different merged inputs.
1278 // Number of track dictionary arrays
1279 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1282 fPar = new AliTRDparameter("TRDparameter","Standard parameter");
1284 printf("<AliTRDdigitizer::MergeSDigits> ");
1285 printf("Create the default parameter object\n");
1292 AliTRDdataArrayI *digitsA;
1293 AliTRDdataArrayI *digitsB;
1294 AliTRDdataArrayI *dictionaryA[kNDict];
1295 AliTRDdataArrayI *dictionaryB[kNDict];
1297 // Get the first s-digits
1298 fSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->First();
1299 if (!fSDigitsManager) return kFALSE;
1301 // Loop through the other sets of s-digits
1302 AliTRDdigitsManager *mergeSDigitsManager;
1303 mergeSDigitsManager = (AliTRDdigitsManager *)
1304 fSDigitsManagerList->After(fSDigitsManager);
1307 if (mergeSDigitsManager) {
1308 printf("<AliTRDdigitizer::MergeSDigits> ");
1309 printf("Merge %d input files.\n",fSDigitsManagerList->GetSize());
1312 printf("<AliTRDdigitizer::MergeSDigits> ");
1313 printf("Only one input file.\n");
1318 while (mergeSDigitsManager) {
1322 // Loop through the detectors
1323 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1325 Int_t plane = fGeo->GetPlane(iDet);
1326 Int_t sector = fGeo->GetSector(iDet);
1327 Int_t chamber = fGeo->GetChamber(iDet);
1328 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1329 Int_t nColMax = fPar->GetColMax(plane);
1330 Int_t nTimeTotal = fPar->GetTimeTotal();
1332 // Loop through the pixels of one detector and add the signals
1333 digitsA = fSDigitsManager->GetDigits(iDet);
1334 digitsB = mergeSDigitsManager->GetDigits(iDet);
1337 for (iDict = 0; iDict < kNDict; iDict++) {
1338 dictionaryA[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1339 dictionaryB[iDict] = mergeSDigitsManager->GetDictionary(iDet,iDict);
1340 dictionaryA[iDict]->Expand();
1341 dictionaryB[iDict]->Expand();
1344 // Merge only detectors that contain a signal
1345 Bool_t doMerge = kTRUE;
1346 if (fMergeSignalOnly) {
1347 if (digitsA->GetOverThreshold(0) == 0) {
1355 printf("<AliTRDdigitizer::MergeSDigits> ");
1356 printf("Merge detector %d of input no.%d\n",iDet,iMerge+1);
1359 for (Int_t iRow = 0; iRow < nRowMax; iRow++ ) {
1360 for (Int_t iCol = 0; iCol < nColMax; iCol++ ) {
1361 for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {
1363 // Add the amplitudes of the summable digits
1364 Int_t ampA = digitsA->GetDataUnchecked(iRow,iCol,iTime);
1365 Int_t ampB = digitsB->GetDataUnchecked(iRow,iCol,iTime);
1367 digitsA->SetDataUnchecked(iRow,iCol,iTime,ampA);
1369 // Add the mask to the track id if defined.
1370 for (iDict = 0; iDict < kNDict; iDict++) {
1371 Int_t trackB = dictionaryB[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1372 if ((fMasks) && (trackB > 0)) {
1373 for (jDict = 0; jDict < kNDict; jDict++) {
1374 Int_t trackA = dictionaryA[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1376 trackA = trackB + fMasks[iMerge];
1377 dictionaryA[iDict]->SetDataUnchecked(iRow,iCol,iTime,trackA);
1390 digitsA->Compress(1,0);
1391 digitsB->Compress(1,0);
1392 for (iDict = 0; iDict < kNDict; iDict++) {
1393 dictionaryA[iDict]->Compress(1,0);
1394 dictionaryB[iDict]->Compress(1,0);
1400 // The next set of s-digits
1401 mergeSDigitsManager = (AliTRDdigitsManager *)
1402 fSDigitsManagerList->After(mergeSDigitsManager);
1410 //_____________________________________________________________________________
1411 Bool_t AliTRDdigitizer::SDigits2Digits()
1414 // Merges the input s-digits and converts them to normal digits
1417 if (!MergeSDigits()) return kFALSE;
1419 return ConvertSDigits();
1423 //_____________________________________________________________________________
1424 Bool_t AliTRDdigitizer::CheckDetector(Int_t plane, Int_t chamber, Int_t sector)
1427 // Checks whether a detector is enabled
1430 if ((fTRD->GetSensChamber() >= 0) &&
1431 (fTRD->GetSensChamber() != chamber)) return kFALSE;
1432 if ((fTRD->GetSensPlane() >= 0) &&
1433 (fTRD->GetSensPlane() != plane)) return kFALSE;
1434 if ( fTRD->GetSensSector() >= 0) {
1435 Int_t sens1 = fTRD->GetSensSector();
1436 Int_t sens2 = sens1 + fTRD->GetSensSectorRange();
1437 sens2 -= ((Int_t) (sens2 / AliTRDgeometry::Nsect()))
1438 * AliTRDgeometry::Nsect();
1439 if (sens1 < sens2) {
1440 if ((sector < sens1) || (sector >= sens2)) return kFALSE;
1443 if ((sector < sens1) && (sector >= sens2)) return kFALSE;
1451 //_____________________________________________________________________________
1452 Bool_t AliTRDdigitizer::WriteDigits()
1455 // Writes out the TRD-digits and the dictionaries
1458 // Store the digits and the dictionary in the tree
1459 return fDigitsManager->WriteDigits();
1463 //_____________________________________________________________________________
1464 void AliTRDdigitizer::DeConvExp(Double_t *source, Double_t *target
1465 , Int_t n, Int_t nexp)
1468 // Does the deconvolution by the digital filter.
1470 // Author: Marcus Gutfleisch, KIP Heidelberg
1471 // Optimized for: New TRF from Venelin Angelov, simulated with CADENCE
1472 // Pad-ground capacitance = 25 pF
1473 // Pad-pad cross talk capacitance = 6 pF
1474 // For 10 MHz digitization, corresponding to 20 time bins
1475 // in the drift region
1479 Double_t coefficients[2];
1481 /* initialize (coefficient = alpha, rates = lambda) */
1484 rates[0] = 0.466998;
1486 coefficients[0] = 1.0;
1489 rates[0] = 0.8988162;
1490 coefficients[0] = 0.11392069;
1491 rates[1] = 0.3745688;
1492 coefficients[1] = 0.8860793;
1494 Float_t sumc = coefficients[0]+coefficients[1];
1495 coefficients[0] /= sumc;
1496 coefficients[1] /= sumc;
1500 Double_t reminder[2];
1501 Double_t correction, result;
1503 /* attention: computation order is important */
1505 for ( k = 0; k < nexp; k++ ) reminder[k]=0.0;
1507 for ( i = 0; i < n; i++ ) {
1508 result = ( source[i] - correction ); /* no rescaling */
1511 for ( k = 0; k < nexp; k++ ) reminder[k] = rates[k]
1512 * ( reminder[k] + coefficients[k] * result);
1515 for ( k = 0; k < nexp; k++ ) correction += reminder[k];