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.35 2002/03/28 14:59:07 cblume
21 Revision 1.34 2002/03/25 20:00:44 cblume
22 Introduce parameter class and regions of interest for merging
24 Revision 1.33 2002/02/12 17:32:03 cblume
25 Rearrange the deleting of the list of sdigitsmanager
27 Revision 1.32 2002/02/12 16:07:21 cblume
30 Revision 1.31 2002/02/11 14:27:11 cblume
31 New pad plane design, new TRF+PRF, tail cancelation, cross talk
33 Revision 1.30 2001/11/19 08:44:08 cblume
34 Fix bugs reported by Rene
36 Revision 1.29 2001/11/14 19:44:25 hristov
37 Numeric const casted (Alpha)
39 Revision 1.28 2001/11/14 16:35:58 cblume
40 Inherits now from AliDetector
42 Revision 1.27 2001/11/14 10:50:45 cblume
43 Changes in digits IO. Add merging of summable digits
45 Revision 1.26 2001/11/06 17:19:41 cblume
46 Add detailed geometry and simple simulator
48 Revision 1.25 2001/06/27 09:54:44 cblume
49 Moved fField initialization to InitDetector()
51 Revision 1.24 2001/05/21 16:45:47 hristov
52 Last minute changes (C.Blume)
54 Revision 1.23 2001/05/07 08:04:48 cblume
55 New TRF and PRF. Speedup of the code. Digits from amplification region included
57 Revision 1.22 2001/03/30 14:40:14 cblume
58 Update of the digitization parameter
60 Revision 1.21 2001/03/13 09:30:35 cblume
61 Update of digitization. Moved digit branch definition to AliTRD
63 Revision 1.20 2001/02/25 20:19:00 hristov
64 Minor correction: loop variable declared only once for HP, Sun
66 Revision 1.19 2001/02/14 18:22:26 cblume
67 Change in the geometry of the padplane
69 Revision 1.18 2001/01/26 19:56:57 hristov
70 Major upgrade of AliRoot code
72 Revision 1.17 2000/12/08 12:53:27 cblume
73 Change in Copy() function for HP-compiler
75 Revision 1.16 2000/12/07 12:20:46 cblume
76 Go back to array compression. Use sampled PRF to speed up digitization
78 Revision 1.15 2000/11/23 14:34:08 cblume
79 Fixed bug in expansion routine of arrays (initialize buffers properly)
81 Revision 1.14 2000/11/20 08:54:44 cblume
82 Switch off compression as default
84 Revision 1.13 2000/11/10 14:57:52 cblume
85 Changes in the geometry constants for the DEC compiler
87 Revision 1.12 2000/11/01 14:53:20 cblume
88 Merge with TRD-develop
90 Revision 1.1.4.9 2000/10/26 17:00:22 cblume
91 Fixed bug in CheckDetector()
93 Revision 1.1.4.8 2000/10/23 13:41:35 cblume
94 Added protection against Log(0) in the gas gain calulation
96 Revision 1.1.4.7 2000/10/17 02:27:34 cblume
97 Get rid of global constants
99 Revision 1.1.4.6 2000/10/16 01:16:53 cblume
100 Changed timebin 0 to be the one closest to the readout
102 Revision 1.1.4.5 2000/10/15 23:34:29 cblume
103 Faster version of the digitizer
105 Revision 1.1.4.4 2000/10/06 16:49:46 cblume
108 Revision 1.1.4.3 2000/10/04 16:34:58 cblume
109 Replace include files by forward declarations
111 Revision 1.1.4.2 2000/09/22 14:41:10 cblume
112 Bug fix in PRF. Included time response. New structure
114 Revision 1.10 2000/10/05 07:27:53 cblume
115 Changes in the header-files by FCA
117 Revision 1.9 2000/10/02 21:28:19 fca
118 Removal of useless dependecies via forward declarations
120 Revision 1.8 2000/06/09 11:10:07 cblume
121 Compiler warnings and coding conventions, next round
123 Revision 1.7 2000/06/08 18:32:58 cblume
124 Make code compliant to coding conventions
126 Revision 1.6 2000/06/07 16:27:32 cblume
127 Try to remove compiler warnings on Sun and HP
129 Revision 1.5 2000/05/09 16:38:57 cblume
130 Removed PadResponse(). Merge problem
132 Revision 1.4 2000/05/08 15:53:45 cblume
133 Resolved merge conflict
135 Revision 1.3 2000/04/28 14:49:27 cblume
136 Only one declaration of iDict in MakeDigits()
138 Revision 1.1.4.1 2000/05/08 14:42:04 cblume
139 Introduced AliTRDdigitsManager
141 Revision 1.1 2000/02/28 19:00:13 cblume
146 ///////////////////////////////////////////////////////////////////////////////
148 // Creates and handles digits from TRD hits //
149 // Author: C. Blume (C.Blume@gsi.de) //
151 // The following effects are included: //
154 // - Gas gain including fluctuations //
155 // - Pad-response (simple Gaussian approximation) //
156 // - Time-response //
157 // - Electronics noise //
158 // - Electronics gain //
160 // - ADC threshold //
161 // The corresponding parameter can be adjusted via the various //
162 // Set-functions. If these parameters are not explicitly set, default //
163 // values are used (see Init-function). //
164 // As an example on how to use this class to produce digits from hits //
165 // have a look at the macro hits2digits.C //
166 // The production of summable digits is demonstrated in hits2sdigits.C //
167 // and the subsequent conversion of the s-digits into normal digits is //
168 // explained in sdigits2digits.C. //
170 ///////////////////////////////////////////////////////////////////////////////
186 #include "AliRunDigitizer.h"
189 #include "AliTRDhit.h"
190 #include "AliTRDdigitizer.h"
191 #include "AliTRDdataArrayI.h"
192 #include "AliTRDdataArrayF.h"
193 #include "AliTRDsegmentArray.h"
194 #include "AliTRDdigitsManager.h"
195 #include "AliTRDgeometry.h"
196 #include "AliTRDparameter.h"
198 ClassImp(AliTRDdigitizer)
200 //_____________________________________________________________________________
201 AliTRDdigitizer::AliTRDdigitizer()
204 // AliTRDdigitizer default constructor
209 fSDigitsManagerList = 0;
220 fMergeSignalOnly = kFALSE;
224 //_____________________________________________________________________________
225 AliTRDdigitizer::AliTRDdigitizer(const Text_t *name, const Text_t *title)
226 :AliDigitizer(name,title)
229 // AliTRDdigitizer constructor
234 fSDigitsManagerList = 0;
244 fMergeSignalOnly = kFALSE;
246 // For the summable digits
247 fSDigitsScale = 100.;
251 //_____________________________________________________________________________
252 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager
253 , const Text_t *name, const Text_t *title)
254 :AliDigitizer(manager,name,title)
257 // AliTRDdigitizer constructor
262 fSDigitsManagerList = 0;
272 fMergeSignalOnly = kFALSE;
274 // For the summable digits
275 fSDigitsScale = 100.;
279 //_____________________________________________________________________________
280 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager)
281 :AliDigitizer(manager,"AliTRDdigitizer","TRD digitizer")
284 // AliTRDdigitizer constructor
289 fSDigitsManagerList = 0;
299 fMergeSignalOnly = kFALSE;
301 // For the summable digits
302 fSDigitsScale = 100.;
306 //_____________________________________________________________________________
307 AliTRDdigitizer::AliTRDdigitizer(const AliTRDdigitizer &d)
310 // AliTRDdigitizer copy constructor
313 ((AliTRDdigitizer &) d).Copy(*this);
317 //_____________________________________________________________________________
318 AliTRDdigitizer::~AliTRDdigitizer()
321 // AliTRDdigitizer destructor
330 if (fDigitsManager) {
331 delete fDigitsManager;
335 if (fSDigitsManager) {
336 delete fSDigitsManager;
340 if (fSDigitsManagerList) {
341 delete fSDigitsManagerList;
342 fSDigitsManagerList = 0;
352 //_____________________________________________________________________________
353 AliTRDdigitizer &AliTRDdigitizer::operator=(const AliTRDdigitizer &d)
356 // Assignment operator
359 if (this != &d) ((AliTRDdigitizer &) d).Copy(*this);
364 //_____________________________________________________________________________
365 void AliTRDdigitizer::Copy(TObject &d)
371 ((AliTRDdigitizer &) d).fInputFile = 0;
372 ((AliTRDdigitizer &) d).fSDigitsManagerList = 0;
373 ((AliTRDdigitizer &) d).fSDigitsManager = 0;
374 ((AliTRDdigitizer &) d).fDigitsManager = 0;
375 ((AliTRDdigitizer &) d).fTRD = 0;
376 ((AliTRDdigitizer &) d).fGeo = 0;
377 ((AliTRDdigitizer &) d).fMasks = 0;
378 ((AliTRDdigitizer &) d).fEvent = 0;
379 ((AliTRDdigitizer &) d).fPar = 0;
380 ((AliTRDdigitizer &) d).fCompress = fCompress;
381 ((AliTRDdigitizer &) d).fDebug = fDebug ;
382 ((AliTRDdigitizer &) d).fSDigits = fSDigits;
383 ((AliTRDdigitizer &) d).fSDigitsScale = fSDigitsScale;
384 ((AliTRDdigitizer &) d).fMergeSignalOnly = fMergeSignalOnly;
388 //_____________________________________________________________________________
389 void AliTRDdigitizer::Exec(Option_t* option)
392 // Executes the merging
397 AliTRDdigitsManager *sdigitsManager;
399 TString optionString = option;
400 if (optionString.Contains("deb")) {
402 if (optionString.Contains("2")) {
405 printf("<AliTRDdigitizer::Exec> ");
406 printf("Called with debug option %d\n",fDebug);
409 // Connect the AliRoot file containing Geometry, Kine, and Hits
410 fInputFile = (TFile *) fManager->GetInputTreeTRDS(0)->GetCurrentFile();
413 printf("<AliTRDdigitizer::Exec> ");
414 printf("Cannot open the input file %s.\n",fInputFile->GetName());
422 gAlice = (AliRun *) fInputFile->Get("gAlice");
425 printf("<AliTRDdigitizer::Exec> ");
426 printf("AliRun object found on file.\n");
430 printf("<AliTRDdigitizer::Exec> ");
431 printf("Could not find AliRun object.\n");
434 Int_t nInput = fManager->GetNinputs();
435 fMasks = new Int_t[nInput];
436 for (iInput = 0; iInput < nInput; iInput++) {
437 fMasks[iInput] = fManager->GetMask(iInput);
443 for (iInput = 0; iInput < nInput; iInput++) {
446 printf("<AliTRDdigitizer::Exec> ");
447 printf("Add input stream %d\n",iInput);
450 // Read the s-digits via digits manager
451 sdigitsManager = new AliTRDdigitsManager();
452 sdigitsManager->SetDebug(fDebug);
453 sdigitsManager->SetSDigits(kTRUE);
454 sdigitsManager->ReadDigits(fManager->GetInputTreeTRDS(iInput));
456 // Add the s-digits to the input list
457 AddSDigitsManager(sdigitsManager);
461 // Convert the s-digits to normal digits
463 printf("<AliTRDdigitizer::Exec> ");
464 printf("Do the conversion\n");
470 printf("<AliTRDdigitizer::Exec> ");
471 printf("Write the digits\n");
473 fDigitsManager->MakeBranch(fManager->GetTreeDTRD());
474 fDigitsManager->WriteDigits();
476 printf("<AliTRDdigitizer::Exec> ");
480 DeleteSDigitsManager();
484 //_____________________________________________________________________________
485 Bool_t AliTRDdigitizer::Open(const Char_t *file, Int_t nEvent)
488 // Opens a ROOT-file with TRD-hits and reads in the hit-tree
491 // Connect the AliRoot file containing Geometry, Kine, and Hits
492 fInputFile = (TFile *) gROOT->GetListOfFiles()->FindObject(file);
495 printf("<AliTRDdigitizer::Open> ");
496 printf("Open the AliROOT-file %s.\n",file);
498 fInputFile = new TFile(file,"UPDATE");
502 printf("<AliTRDdigitizer::Open> ");
503 printf("%s is already open.\n",file);
507 gAlice = (AliRun *) fInputFile->Get("gAlice");
510 printf("<AliTRDdigitizer::Open> ");
511 printf("AliRun object found on file.\n");
515 printf("<AliTRDdigitizer::Open> ");
516 printf("Could not find AliRun object.\n");
522 // Import the Trees for the event nEvent in the file
523 Int_t nparticles = gAlice->GetEvent(fEvent);
524 if (nparticles <= 0) {
525 printf("<AliTRDdigitizer::Open> ");
526 printf("No entries in the trees for event %d.\n",fEvent);
530 if (InitDetector()) {
539 //_____________________________________________________________________________
540 Bool_t AliTRDdigitizer::InitDetector()
543 // Sets the pointer to the TRD detector and the geometry
546 // Get the pointer to the detector class and check for version 1
547 fTRD = (AliTRD *) gAlice->GetDetector("TRD");
548 if (fTRD->IsVersion() != 1) {
549 printf("<AliTRDdigitizer::InitDetector> ");
550 printf("TRD must be version 1 (slow simulator).\n");
555 fGeo = fTRD->GetGeometry();
557 printf("<AliTRDdigitizer::InitDetector> ");
558 printf("Geometry version %d\n",fGeo->IsVersion());
561 // Create a digits manager
562 fDigitsManager = new AliTRDdigitsManager();
563 fDigitsManager->SetSDigits(fSDigits);
564 fDigitsManager->CreateArrays();
565 fDigitsManager->SetEvent(fEvent);
566 fDigitsManager->SetDebug(fDebug);
568 // The list for the input s-digits manager to be merged
569 fSDigitsManagerList = new TList();
575 //_____________________________________________________________________________
576 Bool_t AliTRDdigitizer::MakeBranch(const Char_t *file) const
579 // Create the branches for the digits array
582 return fDigitsManager->MakeBranch(file);
586 //_____________________________________________________________________________
587 Bool_t AliTRDdigitizer::MakeDigits()
593 ///////////////////////////////////////////////////////////////
595 ///////////////////////////////////////////////////////////////
597 // Converts number of electrons to fC
598 const Double_t kEl2fC = 1.602E-19 * 1.0E15;
600 ///////////////////////////////////////////////////////////////
602 // Number of pads included in the pad response
603 const Int_t kNpad = 3;
605 // Number of track dictionary arrays
606 const Int_t kNDict = AliTRDdigitsManager::kNDict;
608 // Half the width of the amplification region
609 const Float_t kAmWidth = AliTRDgeometry::AmThick() / 2.;
611 Int_t iRow, iCol, iTime, iPad;
615 Int_t totalSizeDigits = 0;
616 Int_t totalSizeDict0 = 0;
617 Int_t totalSizeDict1 = 0;
618 Int_t totalSizeDict2 = 0;
620 Int_t timeTRDbeg = 0;
621 Int_t timeTRDend = 1;
626 Float_t padSignal[kNpad];
627 Float_t signalOld[kNpad];
629 AliTRDdataArrayF *signals = 0;
630 AliTRDdataArrayI *digits = 0;
631 AliTRDdataArrayI *dictionary[kNDict];
633 // Create a default parameter class if none is defined
635 fPar = new AliTRDparameter("TRDparameter","Standard TRD parameter");
637 printf("<AliTRDdigitizer::MakeDigits> ");
638 printf("Create the default parameter object\n");
642 // Create a container for the amplitudes
643 AliTRDsegmentArray *signalsArray
644 = new AliTRDsegmentArray("AliTRDdataArrayF"
645 ,AliTRDgeometry::Ndet());
648 timeTRDbeg = ((Int_t) (-fPar->GetTRFlo() / fPar->GetTimeBinSize())) - 1;
649 timeTRDend = ((Int_t) ( fPar->GetTRFhi() / fPar->GetTimeBinSize())) - 1;
651 printf("<AliTRDdigitizer::MakeDigits> ");
652 printf("Sample the TRF between -%d and %d\n",timeTRDbeg,timeTRDend);
656 Float_t elAttachProp = fPar->GetElAttachProp() / 100.;
659 printf("<AliTRDdigitizer::MakeDigits> ");
660 printf("No geometry defined\n");
665 printf("<AliTRDdigitizer::MakeDigits> ");
666 printf("Start creating digits.\n");
669 // Get the pointer to the hit tree
670 TTree *hitTree = gAlice->TreeH();
672 // Get the number of entries in the hit tree
673 // (Number of primary particles creating a hit somewhere)
674 Int_t nTrack = (Int_t) hitTree->GetEntries();
676 printf("<AliTRDdigitizer::MakeDigits> ");
677 printf("Found %d primary particles\n",nTrack);
680 Int_t detectorOld = -1;
683 // Loop through all entries in the tree
684 for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
687 nBytes += hitTree->GetEvent(iTrack);
689 // Loop through the TRD hits
691 AliTRDhit *hit = (AliTRDhit *) fTRD->FirstHit(-1);
700 Float_t q = hit->GetCharge();
701 Int_t track = hit->Track();
702 Int_t detector = hit->GetDetector();
703 Int_t plane = fGeo->GetPlane(detector);
704 Int_t sector = fGeo->GetSector(detector);
705 Int_t chamber = fGeo->GetChamber(detector);
706 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
707 Int_t nColMax = fPar->GetColMax(plane);
708 Int_t nTimeMax = fPar->GetTimeMax();
709 Int_t nTimeBefore = fPar->GetTimeBefore();
710 Int_t nTimeAfter = fPar->GetTimeAfter();
711 Int_t nTimeTotal = fPar->GetTimeTotal();
712 Float_t row0 = fPar->GetRow0(plane,chamber,sector);
713 Float_t col0 = fPar->GetCol0(plane);
714 Float_t time0 = fPar->GetTime0(plane);
715 Float_t rowPadSize = fPar->GetRowPadSize(plane,chamber,sector);
716 Float_t colPadSize = fPar->GetColPadSize(plane);
717 Float_t timeBinSize = fPar->GetTimeBinSize();
718 Float_t divideRow = 1.0 / rowPadSize;
719 Float_t divideCol = 1.0 / colPadSize;
720 Float_t divideTime = 1.0 / timeBinSize;
723 printf("Analyze hit no. %d ",iHit);
724 printf("-----------------------------------------------------------\n");
726 printf("plane = %d, sector = %d, chamber = %d\n"
727 ,plane,sector,chamber);
728 printf("nRowMax = %d, nColMax = %d, nTimeMax = %d\n"
729 ,nRowMax,nColMax,nTimeMax);
730 printf("nTimeBefore = %d, nTimeAfter = %d, nTimeTotal = %d\n"
731 ,nTimeBefore,nTimeAfter,nTimeTotal);
732 printf("row0 = %f, col0 = %f, time0 = %f\n"
734 printf("rowPadSize = %f, colPadSize = %f, timeBinSize = %f\n"
735 ,rowPadSize,colPadSize,timeBinSize);
738 // Don't analyze test hits and switched off detectors
739 if ((CheckDetector(plane,chamber,sector)) &&
740 (((Int_t) q) != 0)) {
742 if (detector != detectorOld) {
745 printf("<AliTRDdigitizer::MakeDigits> ");
746 printf("Get new container. New det = %d, Old det = %d\n"
747 ,detector,detectorOld);
749 // Compress the old one if enabled
750 if ((fCompress) && (detectorOld > -1)) {
752 printf("<AliTRDdigitizer::MakeDigits> ");
753 printf("Compress the old container ...");
755 signals->Compress(1,0);
756 for (iDict = 0; iDict < kNDict; iDict++) {
757 dictionary[iDict]->Compress(1,0);
759 if (fDebug > 1) printf("done\n");
761 // Get the new container
762 signals = (AliTRDdataArrayF *) signalsArray->At(detector);
763 if (signals->GetNtime() == 0) {
764 // Allocate a new one if not yet existing
766 printf("<AliTRDdigitizer::MakeDigits> ");
767 printf("Allocate a new container ... ");
769 signals->Allocate(nRowMax,nColMax,nTimeTotal);
772 // Expand an existing one
775 printf("<AliTRDdigitizer::MakeDigits> ");
776 printf("Expand an existing container ... ");
781 // The same for the dictionary
782 for (iDict = 0; iDict < kNDict; iDict++) {
783 dictionary[iDict] = fDigitsManager->GetDictionary(detector,iDict);
784 if (dictionary[iDict]->GetNtime() == 0) {
785 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
788 if (fCompress) dictionary[iDict]->Expand();
791 if (fDebug > 1) printf("done\n");
792 detectorOld = detector;
795 // Rotate the sectors on top of each other
796 fGeo->Rotate(detector,pos,rot);
798 // The driftlength. It is negative if the hit is in the
799 // amplification region.
800 Float_t driftlength = time0 - rot[0];
802 // Take also the drift in the amplification region into account
803 // The drift length is at the moment still the same, regardless of
804 // the position relativ to the wire. This non-isochronity needs still
805 // to be implemented.
806 Float_t driftlengthL = TMath::Abs(driftlength + kAmWidth);
807 if (fPar->ExBOn()) driftlengthL /= TMath::Sqrt(fPar->GetLorentzFactor());
809 // Loop over all electrons of this hit
810 // TR photons produce hits with negative charge
811 Int_t nEl = ((Int_t) TMath::Abs(q));
812 for (Int_t iEl = 0; iEl < nEl; iEl++) {
818 // Electron attachment
819 if (fPar->ElAttachOn()) {
820 if (gRandom->Rndm() < (driftlengthL * elAttachProp))
824 // Apply the diffusion smearing
825 if (fPar->DiffusionOn()) {
826 if (!(fPar->Diffusion(driftlengthL,xyz))) continue;
829 // Apply E x B effects (depends on drift direction)
831 if (!(fPar->ExB(driftlength+kAmWidth,xyz))) continue;
834 // The electron position after diffusion and ExB in pad coordinates
835 // The pad row (z-direction)
836 Float_t rowDist = xyz[2] - row0;
837 Int_t rowE = ((Int_t) (rowDist * divideRow));
838 if ((rowE < 0) || (rowE >= nRowMax)) continue;
839 Float_t rowOffset = ((((Float_t) rowE) + 0.5) * rowPadSize) - rowDist;
841 // The pad column (rphi-direction)
842 Float_t col0tilt = fPar->Col0Tilted(col0,rowOffset,plane);
843 Float_t colDist = xyz[1] - col0tilt;
844 Int_t colE = ((Int_t) (colDist * divideCol));
845 if ((colE < 0) || (colE >= nColMax)) continue;
846 Float_t colOffset = ((((Float_t) colE) + 0.5) * colPadSize) - colDist;
848 // The time bin (negative for hits in the amplification region)
849 // In the amplification region the electrons drift from both sides
850 // to the middle (anode wire plane)
851 Float_t timeDist = time0 - xyz[0];
852 Float_t timeOffset = 0;
856 timeE = ((Int_t) (timeDist * divideTime));
857 // The distance of the position to the middle of the timebin
858 timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) - timeDist;
861 // Difference between half of the amplification gap width and
862 // the distance to the anode wire
863 Float_t anodeDist = kAmWidth - TMath::Abs(timeDist + kAmWidth);
865 timeE = -1 * (((Int_t ) (anodeDist * divideTime)) + 1);
866 // The distance of the position to the middle of the timebin
867 timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) + anodeDist;
870 // Apply the gas gain including fluctuations
871 Float_t ggRndm = 0.0;
873 ggRndm = gRandom->Rndm();
874 } while (ggRndm <= 0);
875 Int_t signal = (Int_t) (-fPar->GetGasGain() * TMath::Log(ggRndm));
877 // Apply the pad response
879 // The distance of the electron to the center of the pad
880 // in units of pad width
881 Float_t dist = - colOffset * divideCol;
882 if (!(fPar->PadResponse(signal,dist,plane,padSignal))) continue;
886 padSignal[1] = signal;
890 // Sample the time response inside the drift region
891 // + additional time bins before and after.
892 // The sampling is done always in the middle of the time bin
893 for (Int_t iTimeBin = TMath::Max(timeE-timeTRDbeg, -nTimeBefore)
894 ;iTimeBin < TMath::Min(timeE+timeTRDend,nTimeMax+nTimeAfter )
897 // Apply the time response
898 Float_t timeResponse = 1.0;
899 Float_t crossTalk = 0.0;
900 Float_t time = (iTimeBin - timeE) * timeBinSize + timeOffset;
902 timeResponse = fPar->TimeResponse(time);
905 crossTalk = fPar->CrossTalk(time);
912 for (iPad = 0; iPad < kNpad; iPad++) {
914 Int_t colPos = colE + iPad - 1;
915 if (colPos < 0) continue;
916 if (colPos >= nColMax) break;
919 // Note: The time bin number is shifted by nTimeBefore to avoid negative
920 // time bins. This has to be subtracted later.
921 Int_t iCurrentTimeBin = iTimeBin + nTimeBefore;
922 signalOld[iPad] = signals->GetDataUnchecked(rowE,colPos,iCurrentTimeBin);
923 if( colPos != colE ) {
924 signalOld[iPad] += padSignal[iPad] * (timeResponse + crossTalk);
927 signalOld[iPad] += padSignal[iPad] * timeResponse;
929 signals->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,signalOld[iPad]);
931 // Store the track index in the dictionary
932 // Note: We store index+1 in order to allow the array to be compressed
933 if (signalOld[iPad] > 0) {
934 for (iDict = 0; iDict < kNDict; iDict++) {
935 Int_t oldTrack = dictionary[iDict]->GetDataUnchecked(rowE
938 if (oldTrack == track+1) break;
940 dictionary[iDict]->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,track+1);
950 } // Loop: electrons of a single hit
952 } // If: detector and test hit
954 hit = (AliTRDhit *) fTRD->NextHit();
956 } // Loop: hits of one primary track
958 } // Loop: primary tracks
961 printf("<AliTRDdigitizer::MakeDigits> ");
962 printf("Finished analyzing %d hits\n",countHits);
965 // The total conversion factor
966 Float_t convert = kEl2fC * fPar->GetPadCoupling()
967 * fPar->GetTimeCoupling()
968 * fPar->GetChipGain();
970 // Loop through all chambers to finalize the digits
971 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
973 Int_t plane = fGeo->GetPlane(iDet);
974 Int_t sector = fGeo->GetSector(iDet);
975 Int_t chamber = fGeo->GetChamber(iDet);
976 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
977 Int_t nColMax = fPar->GetColMax(plane);
978 Int_t nTimeMax = fPar->GetTimeMax();
979 Int_t nTimeTotal = fPar->GetTimeTotal();
981 Double_t *inADC = new Double_t[nTimeTotal];
982 Double_t *outADC = new Double_t[nTimeTotal];
985 printf("<AliTRDdigitizer::MakeDigits> ");
986 printf("Digitization for chamber %d\n",iDet);
989 // Add a container for the digits of this detector
990 digits = fDigitsManager->GetDigits(iDet);
991 // Allocate memory space for the digits buffer
992 digits->Allocate(nRowMax,nColMax,nTimeTotal);
994 // Get the signal container
995 signals = (AliTRDdataArrayF *) signalsArray->At(iDet);
996 if (signals->GetNtime() == 0) {
997 // Create missing containers
998 signals->Allocate(nRowMax,nColMax,nTimeTotal);
1001 // Expand the container if neccessary
1002 if (fCompress) signals->Expand();
1004 // Create the missing dictionary containers
1005 for (iDict = 0; iDict < kNDict; iDict++) {
1006 dictionary[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1007 if (dictionary[iDict]->GetNtime() == 0) {
1008 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1014 // Don't create noise in detectors that are switched off
1015 if (CheckDetector(plane,chamber,sector)) {
1017 // Create the digits for this chamber
1018 for (iRow = 0; iRow < nRowMax; iRow++ ) {
1019 for (iCol = 0; iCol < nColMax; iCol++ ) {
1021 // Create summable digits
1024 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1025 Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1026 signalAmp *= fSDigitsScale;
1027 signalAmp = TMath::Min(signalAmp,(Float_t) 1.0e9);
1028 Int_t adc = (Int_t) signalAmp;
1029 if (adc > 0) nDigits++;
1030 digits->SetDataUnchecked(iRow,iCol,iTime,adc);
1034 // Create normal digits
1037 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1038 Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1040 signalAmp = TMath::Max((Double_t) gRandom->Gaus(signalAmp,fPar->GetNoise()),0.0);
1042 signalAmp *= convert;
1043 // Convert to ADC counts. Set the overflow-bit fADCoutRange if the
1044 // signal is larger than fADCinRange
1046 if (signalAmp >= fPar->GetADCinRange()) {
1047 adc = ((Int_t) fPar->GetADCoutRange());
1050 adc = ((Int_t) (signalAmp * (fPar->GetADCoutRange()
1051 / fPar->GetADCinRange())));
1054 outADC[iTime] = adc;
1057 // Apply the tail cancelation via the digital filter
1059 DeConvExp(inADC,outADC,nTimeTotal,fPar->GetTCnexp());
1062 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1063 // Store the amplitude of the digit if above threshold
1064 if (outADC[iTime] > fPar->GetADCthreshold()) {
1066 printf(" iRow = %d, iCol = %d, iTime = %d, adc = %f\n"
1067 ,iRow,iCol,iTime,outADC[iTime]);
1070 digits->SetDataUnchecked(iRow,iCol,iTime,outADC[iTime]);
1081 // Compress the arrays
1082 digits->Compress(1,0);
1083 for (iDict = 0; iDict < kNDict; iDict++) {
1084 dictionary[iDict]->Compress(1,0);
1087 totalSizeDigits += digits->GetSize();
1088 totalSizeDict0 += dictionary[0]->GetSize();
1089 totalSizeDict1 += dictionary[1]->GetSize();
1090 totalSizeDict2 += dictionary[2]->GetSize();
1092 Float_t nPixel = nRowMax * nColMax * nTimeMax;
1094 printf("<AliTRDdigitizer::MakeDigits> ");
1095 printf("Found %d digits in detector %d (%3.0f).\n"
1097 ,100.0 * ((Float_t) nDigits) / nPixel);
1100 if (fCompress) signals->Compress(1,0);
1108 printf("<AliTRDdigitizer::MakeDigits> ");
1109 printf("Total number of analyzed hits = %d\n",countHits);
1110 printf("<AliTRDdigitizer::MakeDigits> ");
1111 printf("Total digits data size = %d, %d, %d, %d\n",totalSizeDigits
1121 //_____________________________________________________________________________
1122 void AliTRDdigitizer::AddSDigitsManager(AliTRDdigitsManager *man)
1125 // Add a digits manager for s-digits to the input list.
1128 fSDigitsManagerList->Add(man);
1132 //_____________________________________________________________________________
1133 void AliTRDdigitizer::DeleteSDigitsManager()
1136 // Removes digits manager from the input list.
1139 fSDigitsManagerList->Delete();
1143 //_____________________________________________________________________________
1144 Bool_t AliTRDdigitizer::ConvertSDigits()
1147 // Converts s-digits to normal digits
1150 // Number of track dictionary arrays
1151 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1153 // Converts number of electrons to fC
1154 const Double_t kEl2fC = 1.602E-19 * 1.0E15;
1162 fPar = new AliTRDparameter("TRDparameter","Standard parameter");
1164 printf("<AliTRDdigitizer::ConvertSDigits> ");
1165 printf("Create the default parameter object\n");
1169 Double_t sDigitsScale = 1.0 / GetSDigitsScale();
1170 Double_t noise = fPar->GetNoise();
1171 Double_t padCoupling = fPar->GetPadCoupling();
1172 Double_t timeCoupling = fPar->GetTimeCoupling();
1173 Double_t chipGain = fPar->GetChipGain();
1174 Double_t convert = kEl2fC * padCoupling * timeCoupling * chipGain;;
1175 Double_t adcInRange = fPar->GetADCinRange();
1176 Double_t adcOutRange = fPar->GetADCoutRange();
1177 Int_t adcThreshold = fPar->GetADCthreshold();
1179 AliTRDdataArrayI *digitsIn;
1180 AliTRDdataArrayI *digitsOut;
1181 AliTRDdataArrayI *dictionaryIn[kNDict];
1182 AliTRDdataArrayI *dictionaryOut[kNDict];
1184 // Loop through the detectors
1185 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1188 printf("<AliTRDdigitizer::ConvertSDigits> ");
1189 printf("Convert detector %d to digits.\n",iDet);
1192 Int_t plane = fGeo->GetPlane(iDet);
1193 Int_t sector = fGeo->GetSector(iDet);
1194 Int_t chamber = fGeo->GetChamber(iDet);
1195 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1196 Int_t nColMax = fPar->GetColMax(plane);
1197 Int_t nTimeTotal = fPar->GetTimeTotal();
1199 Double_t *inADC = new Double_t[nTimeTotal];
1200 Double_t *outADC = new Double_t[nTimeTotal];
1202 digitsIn = fSDigitsManager->GetDigits(iDet);
1204 digitsOut = fDigitsManager->GetDigits(iDet);
1205 digitsOut->Allocate(nRowMax,nColMax,nTimeTotal);
1206 for (iDict = 0; iDict < kNDict; iDict++) {
1207 dictionaryIn[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1208 dictionaryIn[iDict]->Expand();
1209 dictionaryOut[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1210 dictionaryOut[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1213 for (iRow = 0; iRow < nRowMax; iRow++ ) {
1214 for (iCol = 0; iCol < nColMax; iCol++ ) {
1216 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1217 Double_t signal = (Double_t) digitsIn->GetDataUnchecked(iRow,iCol,iTime);
1218 signal *= sDigitsScale;
1220 signal = TMath::Max((Double_t) gRandom->Gaus(signal,noise),0.0);
1223 // Convert to ADC counts. Set the overflow-bit adcOutRange if the
1224 // signal is larger than adcInRange
1226 if (signal >= adcInRange) {
1227 adc = ((Int_t) adcOutRange);
1230 adc = ((Int_t) (signal * (adcOutRange / adcInRange)));
1233 outADC[iTime] = adc;
1236 // Apply the tail cancelation via the digital filter
1238 DeConvExp(inADC,outADC,nTimeTotal,fPar->GetTCnexp());
1241 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1242 // Store the amplitude of the digit if above threshold
1243 if (outADC[iTime] > adcThreshold) {
1244 digitsOut->SetDataUnchecked(iRow,iCol,iTime,outADC[iTime]);
1245 // Copy the dictionary
1246 for (iDict = 0; iDict < kNDict; iDict++) {
1247 Int_t track = dictionaryIn[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1248 dictionaryOut[iDict]->SetDataUnchecked(iRow,iCol,iTime,track);
1257 digitsIn->Compress(1,0);
1258 digitsOut->Compress(1,0);
1259 for (iDict = 0; iDict < kNDict; iDict++) {
1260 dictionaryIn[iDict]->Compress(1,0);
1261 dictionaryOut[iDict]->Compress(1,0);
1274 //_____________________________________________________________________________
1275 Bool_t AliTRDdigitizer::MergeSDigits()
1278 // Merges the input s-digits:
1279 // - The amplitude of the different inputs are summed up.
1280 // - Of the track IDs from the input dictionaries only one is
1281 // kept for each input. This works for maximal 3 different merged inputs.
1284 // Number of track dictionary arrays
1285 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1288 fPar = new AliTRDparameter("TRDparameter","Standard parameter");
1290 printf("<AliTRDdigitizer::MergeSDigits> ");
1291 printf("Create the default parameter object\n");
1298 AliTRDdataArrayI *digitsA;
1299 AliTRDdataArrayI *digitsB;
1300 AliTRDdataArrayI *dictionaryA[kNDict];
1301 AliTRDdataArrayI *dictionaryB[kNDict];
1303 // Get the first s-digits
1304 fSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->First();
1305 if (!fSDigitsManager) return kFALSE;
1307 // Loop through the other sets of s-digits
1308 AliTRDdigitsManager *mergeSDigitsManager;
1309 mergeSDigitsManager = (AliTRDdigitsManager *)
1310 fSDigitsManagerList->After(fSDigitsManager);
1313 if (mergeSDigitsManager) {
1314 printf("<AliTRDdigitizer::MergeSDigits> ");
1315 printf("Merge %d input files.\n",fSDigitsManagerList->GetSize());
1318 printf("<AliTRDdigitizer::MergeSDigits> ");
1319 printf("Only one input file.\n");
1324 while (mergeSDigitsManager) {
1328 // Loop through the detectors
1329 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1331 Int_t plane = fGeo->GetPlane(iDet);
1332 Int_t sector = fGeo->GetSector(iDet);
1333 Int_t chamber = fGeo->GetChamber(iDet);
1334 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1335 Int_t nColMax = fPar->GetColMax(plane);
1336 Int_t nTimeTotal = fPar->GetTimeTotal();
1338 // Loop through the pixels of one detector and add the signals
1339 digitsA = fSDigitsManager->GetDigits(iDet);
1340 digitsB = mergeSDigitsManager->GetDigits(iDet);
1343 for (iDict = 0; iDict < kNDict; iDict++) {
1344 dictionaryA[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1345 dictionaryB[iDict] = mergeSDigitsManager->GetDictionary(iDet,iDict);
1346 dictionaryA[iDict]->Expand();
1347 dictionaryB[iDict]->Expand();
1350 // Merge only detectors that contain a signal
1351 Bool_t doMerge = kTRUE;
1352 if (fMergeSignalOnly) {
1353 if (digitsA->GetOverThreshold(0) == 0) {
1361 printf("<AliTRDdigitizer::MergeSDigits> ");
1362 printf("Merge detector %d of input no.%d\n",iDet,iMerge+1);
1365 for (Int_t iRow = 0; iRow < nRowMax; iRow++ ) {
1366 for (Int_t iCol = 0; iCol < nColMax; iCol++ ) {
1367 for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {
1369 // Add the amplitudes of the summable digits
1370 Int_t ampA = digitsA->GetDataUnchecked(iRow,iCol,iTime);
1371 Int_t ampB = digitsB->GetDataUnchecked(iRow,iCol,iTime);
1373 digitsA->SetDataUnchecked(iRow,iCol,iTime,ampA);
1375 // Add the mask to the track id if defined.
1376 for (iDict = 0; iDict < kNDict; iDict++) {
1377 Int_t trackB = dictionaryB[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1378 if ((fMasks) && (trackB > 0)) {
1379 for (jDict = 0; jDict < kNDict; jDict++) {
1380 Int_t trackA = dictionaryA[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1382 trackA = trackB + fMasks[iMerge];
1383 dictionaryA[iDict]->SetDataUnchecked(iRow,iCol,iTime,trackA);
1396 digitsA->Compress(1,0);
1397 digitsB->Compress(1,0);
1398 for (iDict = 0; iDict < kNDict; iDict++) {
1399 dictionaryA[iDict]->Compress(1,0);
1400 dictionaryB[iDict]->Compress(1,0);
1406 // The next set of s-digits
1407 mergeSDigitsManager = (AliTRDdigitsManager *)
1408 fSDigitsManagerList->After(mergeSDigitsManager);
1416 //_____________________________________________________________________________
1417 Bool_t AliTRDdigitizer::SDigits2Digits()
1420 // Merges the input s-digits and converts them to normal digits
1423 if (!MergeSDigits()) return kFALSE;
1425 return ConvertSDigits();
1429 //_____________________________________________________________________________
1430 Bool_t AliTRDdigitizer::CheckDetector(Int_t plane, Int_t chamber, Int_t sector)
1433 // Checks whether a detector is enabled
1436 if ((fTRD->GetSensChamber() >= 0) &&
1437 (fTRD->GetSensChamber() != chamber)) return kFALSE;
1438 if ((fTRD->GetSensPlane() >= 0) &&
1439 (fTRD->GetSensPlane() != plane)) return kFALSE;
1440 if ( fTRD->GetSensSector() >= 0) {
1441 Int_t sens1 = fTRD->GetSensSector();
1442 Int_t sens2 = sens1 + fTRD->GetSensSectorRange();
1443 sens2 -= ((Int_t) (sens2 / AliTRDgeometry::Nsect()))
1444 * AliTRDgeometry::Nsect();
1445 if (sens1 < sens2) {
1446 if ((sector < sens1) || (sector >= sens2)) return kFALSE;
1449 if ((sector < sens1) && (sector >= sens2)) return kFALSE;
1457 //_____________________________________________________________________________
1458 Bool_t AliTRDdigitizer::WriteDigits() const
1461 // Writes out the TRD-digits and the dictionaries
1464 // Store the digits and the dictionary in the tree
1465 return fDigitsManager->WriteDigits();
1469 //_____________________________________________________________________________
1470 void AliTRDdigitizer::DeConvExp(Double_t *source, Double_t *target
1471 , Int_t n, Int_t nexp)
1474 // Does the deconvolution by the digital filter.
1476 // Author: Marcus Gutfleisch, KIP Heidelberg
1477 // Optimized for: New TRF from Venelin Angelov, simulated with CADENCE
1478 // Pad-ground capacitance = 25 pF
1479 // Pad-pad cross talk capacitance = 6 pF
1480 // For 10 MHz digitization, corresponding to 20 time bins
1481 // in the drift region
1485 Double_t coefficients[2];
1487 /* initialize (coefficient = alpha, rates = lambda) */
1490 rates[0] = 0.466998;
1492 coefficients[0] = 1.0;
1495 rates[0] = 0.8988162;
1496 coefficients[0] = 0.11392069;
1497 rates[1] = 0.3745688;
1498 coefficients[1] = 0.8860793;
1500 Float_t sumc = coefficients[0]+coefficients[1];
1501 coefficients[0] /= sumc;
1502 coefficients[1] /= sumc;
1506 Double_t reminder[2];
1507 Double_t correction, result;
1509 /* attention: computation order is important */
1511 for ( k = 0; k < nexp; k++ ) reminder[k]=0.0;
1513 for ( i = 0; i < n; i++ ) {
1514 result = ( source[i] - correction ); /* no rescaling */
1517 for ( k = 0; k < nexp; k++ ) reminder[k] = rates[k]
1518 * ( reminder[k] + coefficients[k] * result);
1521 for ( k = 0; k < nexp; k++ ) correction += reminder[k];
1526 //_____________________________________________________________________________
1527 void AliTRDdigitizer::InitOutput(TFile *file, Int_t iEvent)
1530 // Initializes the output branches
1534 fDigitsManager->MakeTreeAndBranches(file,iEvent);