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.42 2002/12/13 13:31:44 cblume
19 Add ADC offset and change default digitazation parameters
21 Revision 1.41 2002/10/21 09:10:32 cblume
22 Fix type conversion warnings
24 Revision 1.40 2002/10/14 14:57:43 hristov
25 Merging the VirtualMC branch to the main development branch (HEAD)
27 Revision 1.33.6.3 2002/10/11 07:26:37 hristov
28 Updating VirtualMC to v3-09-02
30 Revision 1.39 2002/10/08 20:46:12 cblume
31 Do coupling factors before noise is applied
33 Revision 1.38 2002/04/30 08:30:40 cblume
34 gAlice now only read by AliRunDigitizer. Therefore it is just deleted in AliTRDmerge.C
36 Revision 1.37 2002/04/29 11:50:47 cblume
37 Change initialization of gAlice in the merging case
39 Revision 1.36 2002/04/12 12:13:23 cblume
42 Revision 1.35 2002/03/28 14:59:07 cblume
45 Revision 1.34 2002/03/25 20:00:44 cblume
46 Introduce parameter class and regions of interest for merging
48 Revision 1.33 2002/02/12 17:32:03 cblume
49 Rearrange the deleting of the list of sdigitsmanager
51 Revision 1.32 2002/02/12 16:07:21 cblume
54 Revision 1.31 2002/02/11 14:27:11 cblume
55 New pad plane design, new TRF+PRF, tail cancelation, cross talk
57 Revision 1.30 2001/11/19 08:44:08 cblume
58 Fix bugs reported by Rene
60 Revision 1.29 2001/11/14 19:44:25 hristov
61 Numeric const casted (Alpha)
63 Revision 1.28 2001/11/14 16:35:58 cblume
64 Inherits now from AliDetector
66 Revision 1.27 2001/11/14 10:50:45 cblume
67 Changes in digits IO. Add merging of summable digits
69 Revision 1.26 2001/11/06 17:19:41 cblume
70 Add detailed geometry and simple simulator
72 Revision 1.25 2001/06/27 09:54:44 cblume
73 Moved fField initialization to InitDetector()
75 Revision 1.24 2001/05/21 16:45:47 hristov
76 Last minute changes (C.Blume)
78 Revision 1.23 2001/05/07 08:04:48 cblume
79 New TRF and PRF. Speedup of the code. Digits from amplification region included
81 Revision 1.22 2001/03/30 14:40:14 cblume
82 Update of the digitization parameter
84 Revision 1.21 2001/03/13 09:30:35 cblume
85 Update of digitization. Moved digit branch definition to AliTRD
87 Revision 1.20 2001/02/25 20:19:00 hristov
88 Minor correction: loop variable declared only once for HP, Sun
90 Revision 1.19 2001/02/14 18:22:26 cblume
91 Change in the geometry of the padplane
93 Revision 1.18 2001/01/26 19:56:57 hristov
94 Major upgrade of AliRoot code
96 Revision 1.17 2000/12/08 12:53:27 cblume
97 Change in Copy() function for HP-compiler
99 Revision 1.16 2000/12/07 12:20:46 cblume
100 Go back to array compression. Use sampled PRF to speed up digitization
102 Revision 1.15 2000/11/23 14:34:08 cblume
103 Fixed bug in expansion routine of arrays (initialize buffers properly)
105 Revision 1.14 2000/11/20 08:54:44 cblume
106 Switch off compression as default
108 Revision 1.13 2000/11/10 14:57:52 cblume
109 Changes in the geometry constants for the DEC compiler
111 Revision 1.12 2000/11/01 14:53:20 cblume
112 Merge with TRD-develop
114 Revision 1.1.4.9 2000/10/26 17:00:22 cblume
115 Fixed bug in CheckDetector()
117 Revision 1.1.4.8 2000/10/23 13:41:35 cblume
118 Added protection against Log(0) in the gas gain calulation
120 Revision 1.1.4.7 2000/10/17 02:27:34 cblume
121 Get rid of global constants
123 Revision 1.1.4.6 2000/10/16 01:16:53 cblume
124 Changed timebin 0 to be the one closest to the readout
126 Revision 1.1.4.5 2000/10/15 23:34:29 cblume
127 Faster version of the digitizer
129 Revision 1.1.4.4 2000/10/06 16:49:46 cblume
132 Revision 1.1.4.3 2000/10/04 16:34:58 cblume
133 Replace include files by forward declarations
135 Revision 1.1.4.2 2000/09/22 14:41:10 cblume
136 Bug fix in PRF. Included time response. New structure
138 Revision 1.10 2000/10/05 07:27:53 cblume
139 Changes in the header-files by FCA
141 Revision 1.9 2000/10/02 21:28:19 fca
142 Removal of useless dependecies via forward declarations
144 Revision 1.8 2000/06/09 11:10:07 cblume
145 Compiler warnings and coding conventions, next round
147 Revision 1.7 2000/06/08 18:32:58 cblume
148 Make code compliant to coding conventions
150 Revision 1.6 2000/06/07 16:27:32 cblume
151 Try to remove compiler warnings on Sun and HP
153 Revision 1.5 2000/05/09 16:38:57 cblume
154 Removed PadResponse(). Merge problem
156 Revision 1.4 2000/05/08 15:53:45 cblume
157 Resolved merge conflict
159 Revision 1.3 2000/04/28 14:49:27 cblume
160 Only one declaration of iDict in MakeDigits()
162 Revision 1.1.4.1 2000/05/08 14:42:04 cblume
163 Introduced AliTRDdigitsManager
165 Revision 1.1 2000/02/28 19:00:13 cblume
170 ///////////////////////////////////////////////////////////////////////////////
172 // Creates and handles digits from TRD hits //
173 // Author: C. Blume (C.Blume@gsi.de) //
175 // The following effects are included: //
178 // - Gas gain including fluctuations //
179 // - Pad-response (simple Gaussian approximation) //
180 // - Time-response //
181 // - Electronics noise //
182 // - Electronics gain //
184 // - ADC threshold //
185 // The corresponding parameter can be adjusted via the various //
186 // Set-functions. If these parameters are not explicitly set, default //
187 // values are used (see Init-function). //
188 // As an example on how to use this class to produce digits from hits //
189 // have a look at the macro hits2digits.C //
190 // The production of summable digits is demonstrated in hits2sdigits.C //
191 // and the subsequent conversion of the s-digits into normal digits is //
192 // explained in sdigits2digits.C. //
194 ///////////////////////////////////////////////////////////////////////////////
210 #include "AliRunDigitizer.h"
213 #include "AliTRDhit.h"
214 #include "AliTRDdigitizer.h"
215 #include "AliTRDdataArrayI.h"
216 #include "AliTRDdataArrayF.h"
217 #include "AliTRDsegmentArray.h"
218 #include "AliTRDdigitsManager.h"
219 #include "AliTRDgeometry.h"
220 #include "AliTRDparameter.h"
222 ClassImp(AliTRDdigitizer)
224 //_____________________________________________________________________________
225 AliTRDdigitizer::AliTRDdigitizer()
228 // AliTRDdigitizer default constructor
233 fSDigitsManagerList = 0;
244 fMergeSignalOnly = kFALSE;
250 //_____________________________________________________________________________
251 AliTRDdigitizer::AliTRDdigitizer(const Text_t *name, const Text_t *title)
252 :AliDigitizer(name,title)
255 // AliTRDdigitizer constructor
260 fSDigitsManagerList = 0;
270 fMergeSignalOnly = kFALSE;
274 // For the summable digits
275 fSDigitsScale = 100.;
279 //_____________________________________________________________________________
280 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager
281 , const Text_t *name, const Text_t *title)
282 :AliDigitizer(manager,name,title)
285 // AliTRDdigitizer constructor
290 fSDigitsManagerList = 0;
300 fMergeSignalOnly = kFALSE;
304 // For the summable digits
305 fSDigitsScale = 100.;
309 //_____________________________________________________________________________
310 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager)
311 :AliDigitizer(manager,"AliTRDdigitizer","TRD digitizer")
314 // AliTRDdigitizer constructor
319 fSDigitsManagerList = 0;
329 fMergeSignalOnly = kFALSE;
333 // For the summable digits
334 fSDigitsScale = 100.;
338 //_____________________________________________________________________________
339 AliTRDdigitizer::AliTRDdigitizer(const AliTRDdigitizer &d)
342 // AliTRDdigitizer copy constructor
345 ((AliTRDdigitizer &) d).Copy(*this);
349 //_____________________________________________________________________________
350 AliTRDdigitizer::~AliTRDdigitizer()
353 // AliTRDdigitizer destructor
362 if (fDigitsManager) {
363 delete fDigitsManager;
367 if (fSDigitsManager) {
368 delete fSDigitsManager;
372 if (fSDigitsManagerList) {
373 delete fSDigitsManagerList;
374 fSDigitsManagerList = 0;
384 //_____________________________________________________________________________
385 AliTRDdigitizer &AliTRDdigitizer::operator=(const AliTRDdigitizer &d)
388 // Assignment operator
391 if (this != &d) ((AliTRDdigitizer &) d).Copy(*this);
396 //_____________________________________________________________________________
397 void AliTRDdigitizer::Copy(TObject &d)
403 ((AliTRDdigitizer &) d).fInputFile = 0;
404 ((AliTRDdigitizer &) d).fSDigitsManagerList = 0;
405 ((AliTRDdigitizer &) d).fSDigitsManager = 0;
406 ((AliTRDdigitizer &) d).fDigitsManager = 0;
407 ((AliTRDdigitizer &) d).fTRD = 0;
408 ((AliTRDdigitizer &) d).fGeo = 0;
409 ((AliTRDdigitizer &) d).fMasks = 0;
410 ((AliTRDdigitizer &) d).fEvent = 0;
411 ((AliTRDdigitizer &) d).fPar = 0;
412 ((AliTRDdigitizer &) d).fCompress = fCompress;
413 ((AliTRDdigitizer &) d).fDebug = fDebug ;
414 ((AliTRDdigitizer &) d).fSDigits = fSDigits;
415 ((AliTRDdigitizer &) d).fSDigitsScale = fSDigitsScale;
416 ((AliTRDdigitizer &) d).fMergeSignalOnly = fMergeSignalOnly;
417 ((AliTRDdigitizer &) d).fSimpleSim = fSimpleSim;
418 ((AliTRDdigitizer &) d).fSimpleDet = fSimpleDet;
422 //_____________________________________________________________________________
423 void AliTRDdigitizer::Exec(Option_t* option)
426 // Executes the merging
431 AliTRDdigitsManager *sdigitsManager;
433 TString optionString = option;
434 if (optionString.Contains("deb")) {
436 if (optionString.Contains("2")) {
439 printf("<AliTRDdigitizer::Exec> ");
440 printf("Called with debug option %d\n",fDebug);
443 // The AliRoot file is already connected by the manager
446 printf("<AliTRDdigitizer::Exec> ");
447 printf("AliRun object found on file.\n");
451 printf("<AliTRDdigitizer::Exec> ");
452 printf("Could not find AliRun object.\n");
456 Int_t nInput = fManager->GetNinputs();
457 fMasks = new Int_t[nInput];
458 for (iInput = 0; iInput < nInput; iInput++) {
459 fMasks[iInput] = fManager->GetMask(iInput);
465 for (iInput = 0; iInput < nInput; iInput++) {
468 printf("<AliTRDdigitizer::Exec> ");
469 printf("Add input stream %d\n",iInput);
472 // check if the input tree exists
473 if (!fManager->GetInputTreeTRDS(iInput)) {
474 printf("<AliTRDdigitizer::Exec> ");
475 printf("Input stream %d does not exist\n",iInput);
479 // Read the s-digits via digits manager
480 sdigitsManager = new AliTRDdigitsManager();
481 sdigitsManager->SetDebug(fDebug);
482 sdigitsManager->SetSDigits(kTRUE);
483 sdigitsManager->ReadDigits(fManager->GetInputTreeTRDS(iInput));
485 // Add the s-digits to the input list
486 AddSDigitsManager(sdigitsManager);
490 // Convert the s-digits to normal digits
492 printf("<AliTRDdigitizer::Exec> ");
493 printf("Do the conversion\n");
499 printf("<AliTRDdigitizer::Exec> ");
500 printf("Write the digits\n");
502 fDigitsManager->MakeBranch(fManager->GetTreeDTRD());
503 fDigitsManager->WriteDigits();
505 printf("<AliTRDdigitizer::Exec> ");
509 DeleteSDigitsManager();
513 //_____________________________________________________________________________
514 Bool_t AliTRDdigitizer::Open(const Char_t *file, Int_t nEvent)
517 // Opens a ROOT-file with TRD-hits and reads in the hit-tree
520 // Connect the AliRoot file containing Geometry, Kine, and Hits
521 fInputFile = (TFile *) gROOT->GetListOfFiles()->FindObject(file);
524 printf("<AliTRDdigitizer::Open> ");
525 printf("Open the AliROOT-file %s.\n",file);
527 fInputFile = new TFile(file,"UPDATE");
531 printf("<AliTRDdigitizer::Open> ");
532 printf("%s is already open.\n",file);
536 gAlice = (AliRun *) fInputFile->Get("gAlice");
539 printf("<AliTRDdigitizer::Open> ");
540 printf("AliRun object found on file.\n");
544 printf("<AliTRDdigitizer::Open> ");
545 printf("Could not find AliRun object.\n");
551 // Import the Trees for the event nEvent in the file
552 Int_t nparticles = gAlice->GetEvent(fEvent);
553 if (nparticles <= 0) {
554 printf("<AliTRDdigitizer::Open> ");
555 printf("No entries in the trees for event %d.\n",fEvent);
559 if (InitDetector()) {
568 //_____________________________________________________________________________
569 Bool_t AliTRDdigitizer::InitDetector()
572 // Sets the pointer to the TRD detector and the geometry
575 // Get the pointer to the detector class and check for version 1
576 fTRD = (AliTRD *) gAlice->GetDetector("TRD");
578 printf("<AliTRDdigitizer::InitDetector> ");
579 printf("No TRD module found\n");
582 if (fTRD->IsVersion() != 1) {
583 printf("<AliTRDdigitizer::InitDetector> ");
584 printf("TRD must be version 1 (slow simulator).\n");
589 fGeo = fTRD->GetGeometry();
591 printf("<AliTRDdigitizer::InitDetector> ");
592 printf("Geometry version %d\n",fGeo->IsVersion());
595 // Create a digits manager
596 fDigitsManager = new AliTRDdigitsManager();
597 fDigitsManager->SetSDigits(fSDigits);
598 fDigitsManager->CreateArrays();
599 fDigitsManager->SetEvent(fEvent);
600 fDigitsManager->SetDebug(fDebug);
602 // The list for the input s-digits manager to be merged
603 fSDigitsManagerList = new TList();
609 //_____________________________________________________________________________
610 Bool_t AliTRDdigitizer::MakeBranch(const Char_t *file) const
613 // Create the branches for the digits array
616 return fDigitsManager->MakeBranch(file);
620 //_____________________________________________________________________________
621 Bool_t AliTRDdigitizer::MakeDigits()
627 ///////////////////////////////////////////////////////////////
629 ///////////////////////////////////////////////////////////////
631 // Converts number of electrons to fC
632 const Double_t kEl2fC = 1.602E-19 * 1.0E15;
634 ///////////////////////////////////////////////////////////////
636 // Number of pads included in the pad response
637 const Int_t kNpad = 3;
639 // Number of track dictionary arrays
640 const Int_t kNDict = AliTRDdigitsManager::kNDict;
642 // Half the width of the amplification region
643 const Float_t kAmWidth = AliTRDgeometry::AmThick() / 2.;
645 Int_t iRow, iCol, iTime, iPad;
649 Int_t totalSizeDigits = 0;
650 Int_t totalSizeDict0 = 0;
651 Int_t totalSizeDict1 = 0;
652 Int_t totalSizeDict2 = 0;
654 Int_t timeTRDbeg = 0;
655 Int_t timeTRDend = 1;
660 Float_t padSignal[kNpad];
661 Float_t signalOld[kNpad];
663 AliTRDdataArrayF *signals = 0;
664 AliTRDdataArrayI *digits = 0;
665 AliTRDdataArrayI *dictionary[kNDict];
667 // Reset the digits arrays
668 fDigitsManager->ResetArrays();
670 // Create a default parameter class if none is defined
672 fPar = new AliTRDparameter("TRDparameter","Standard TRD parameter");
674 printf("<AliTRDdigitizer::MakeDigits> ");
675 printf("Create the default parameter object\n");
679 // Create a container for the amplitudes
680 AliTRDsegmentArray *signalsArray
681 = new AliTRDsegmentArray("AliTRDdataArrayF"
682 ,AliTRDgeometry::Ndet());
685 timeTRDbeg = ((Int_t) (-fPar->GetTRFlo() / fPar->GetTimeBinSize())) - 1;
686 timeTRDend = ((Int_t) ( fPar->GetTRFhi() / fPar->GetTimeBinSize())) - 1;
688 printf("<AliTRDdigitizer::MakeDigits> ");
689 printf("Sample the TRF between -%d and %d\n",timeTRDbeg,timeTRDend);
693 Float_t elAttachProp = fPar->GetElAttachProp() / 100.;
696 printf("<AliTRDdigitizer::MakeDigits> ");
697 printf("No geometry defined\n");
702 printf("<AliTRDdigitizer::MakeDigits> ");
703 printf("Start creating digits.\n");
706 // Get the pointer to the hit tree
707 TTree *hitTree = gAlice->TreeH();
709 // Get the number of entries in the hit tree
710 // (Number of primary particles creating a hit somewhere)
713 nTrack = (Int_t) hitTree->GetEntries();
715 printf("<AliTRDdigitizer::MakeDigits> ");
716 printf("Found %d primary particles\n",nTrack);
720 Int_t detectorOld = -1;
723 // Loop through all entries in the tree
724 for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
728 nBytes += hitTree->GetEvent(iTrack);
731 // Loop through the TRD hits
733 AliTRDhit *hit = (AliTRDhit *) fTRD->FirstHit(-1);
742 Float_t q = hit->GetCharge();
743 Int_t track = hit->Track();
744 Int_t detector = hit->GetDetector();
745 Int_t plane = fGeo->GetPlane(detector);
746 Int_t sector = fGeo->GetSector(detector);
747 Int_t chamber = fGeo->GetChamber(detector);
748 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
749 Int_t nColMax = fPar->GetColMax(plane);
750 Int_t nTimeMax = fPar->GetTimeMax();
751 Int_t nTimeBefore = fPar->GetTimeBefore();
752 Int_t nTimeAfter = fPar->GetTimeAfter();
753 Int_t nTimeTotal = fPar->GetTimeTotal();
754 Float_t row0 = fPar->GetRow0(plane,chamber,sector);
755 Float_t col0 = fPar->GetCol0(plane);
756 Float_t time0 = fPar->GetTime0(plane);
757 Float_t rowPadSize = fPar->GetRowPadSize(plane,chamber,sector);
758 Float_t colPadSize = fPar->GetColPadSize(plane);
759 Float_t timeBinSize = fPar->GetTimeBinSize();
760 Float_t divideRow = 1.0 / rowPadSize;
761 Float_t divideCol = 1.0 / colPadSize;
762 Float_t divideTime = 1.0 / timeBinSize;
765 printf("Analyze hit no. %d ",iHit);
766 printf("-----------------------------------------------------------\n");
768 printf("plane = %d, sector = %d, chamber = %d\n"
769 ,plane,sector,chamber);
770 printf("nRowMax = %d, nColMax = %d, nTimeMax = %d\n"
771 ,nRowMax,nColMax,nTimeMax);
772 printf("nTimeBefore = %d, nTimeAfter = %d, nTimeTotal = %d\n"
773 ,nTimeBefore,nTimeAfter,nTimeTotal);
774 printf("row0 = %f, col0 = %f, time0 = %f\n"
776 printf("rowPadSize = %f, colPadSize = %f, timeBinSize = %f\n"
777 ,rowPadSize,colPadSize,timeBinSize);
780 // Don't analyze test hits and switched off detectors
781 if ((CheckDetector(plane,chamber,sector)) &&
782 (((Int_t) q) != 0)) {
784 if (detector != detectorOld) {
787 printf("<AliTRDdigitizer::MakeDigits> ");
788 printf("Get new container. New det = %d, Old det = %d\n"
789 ,detector,detectorOld);
791 // Compress the old one if enabled
792 if ((fCompress) && (detectorOld > -1)) {
794 printf("<AliTRDdigitizer::MakeDigits> ");
795 printf("Compress the old container ...");
797 signals->Compress(1,0);
798 for (iDict = 0; iDict < kNDict; iDict++) {
799 dictionary[iDict]->Compress(1,0);
801 if (fDebug > 1) printf("done\n");
803 // Get the new container
804 signals = (AliTRDdataArrayF *) signalsArray->At(detector);
805 if (signals->GetNtime() == 0) {
806 // Allocate a new one if not yet existing
808 printf("<AliTRDdigitizer::MakeDigits> ");
809 printf("Allocate a new container ... ");
811 signals->Allocate(nRowMax,nColMax,nTimeTotal);
813 else if (fSimpleSim) {
814 // Clear an old one for the simple simulation
816 printf("<AliTRDdigitizer::MakeDigits> ");
817 printf("Clear a old container ... ");
822 // Expand an existing one
825 printf("<AliTRDdigitizer::MakeDigits> ");
826 printf("Expand an existing container ... ");
831 // The same for the dictionary
833 for (iDict = 0; iDict < kNDict; iDict++) {
834 dictionary[iDict] = fDigitsManager->GetDictionary(detector,iDict);
835 if (dictionary[iDict]->GetNtime() == 0) {
836 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
839 if (fCompress) dictionary[iDict]->Expand();
843 if (fDebug > 1) printf("done\n");
844 detectorOld = detector;
847 // Rotate the sectors on top of each other
854 fGeo->Rotate(detector,pos,rot);
857 // The driftlength. It is negative if the hit is in the
858 // amplification region.
859 Float_t driftlength = time0 - rot[0];
861 // Take also the drift in the amplification region into account
862 // The drift length is at the moment still the same, regardless of
863 // the position relativ to the wire. This non-isochronity needs still
864 // to be implemented.
865 Float_t driftlengthL = TMath::Abs(driftlength + kAmWidth);
866 if (fPar->ExBOn()) driftlengthL /= TMath::Sqrt(fPar->GetLorentzFactor());
868 // Loop over all electrons of this hit
869 // TR photons produce hits with negative charge
870 Int_t nEl = ((Int_t) TMath::Abs(q));
871 for (Int_t iEl = 0; iEl < nEl; iEl++) {
877 // Electron attachment
878 if (fPar->ElAttachOn()) {
879 if (gRandom->Rndm() < (driftlengthL * elAttachProp))
883 // Apply the diffusion smearing
884 if (fPar->DiffusionOn()) {
885 if (!(fPar->Diffusion(driftlengthL,xyz))) continue;
888 // Apply E x B effects (depends on drift direction)
890 if (!(fPar->ExB(driftlength+kAmWidth,xyz))) continue;
893 // The electron position after diffusion and ExB in pad coordinates
894 // The pad row (z-direction)
895 Float_t rowDist = xyz[2] - row0;
896 Int_t rowE = ((Int_t) (rowDist * divideRow));
897 if ((rowE < 0) || (rowE >= nRowMax)) continue;
898 Float_t rowOffset = ((((Float_t) rowE) + 0.5) * rowPadSize) - rowDist;
900 // The pad column (rphi-direction)
901 Float_t col0tilt = fPar->Col0Tilted(col0,rowOffset,plane);
902 Float_t colDist = xyz[1] - col0tilt;
903 Int_t colE = ((Int_t) (colDist * divideCol));
904 if ((colE < 0) || (colE >= nColMax)) continue;
905 Float_t colOffset = ((((Float_t) colE) + 0.5) * colPadSize) - colDist;
907 // The time bin (negative for hits in the amplification region)
908 // In the amplification region the electrons drift from both sides
909 // to the middle (anode wire plane)
910 Float_t timeDist = time0 - xyz[0];
911 Float_t timeOffset = 0;
915 timeE = ((Int_t) (timeDist * divideTime));
916 // The distance of the position to the middle of the timebin
917 timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) - timeDist;
920 // Difference between half of the amplification gap width and
921 // the distance to the anode wire
922 Float_t anodeDist = kAmWidth - TMath::Abs(timeDist + kAmWidth);
924 timeE = -1 * (((Int_t ) (anodeDist * divideTime)) + 1);
925 // The distance of the position to the middle of the timebin
926 timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) + anodeDist;
929 // Apply the gas gain including fluctuations
930 Float_t ggRndm = 0.0;
932 ggRndm = gRandom->Rndm();
933 } while (ggRndm <= 0);
934 Int_t signal = (Int_t) (-fPar->GetGasGain() * TMath::Log(ggRndm));
936 // Apply the pad response
938 // The distance of the electron to the center of the pad
939 // in units of pad width
940 Float_t dist = - colOffset * divideCol;
941 if (!(fPar->PadResponse(signal,dist,plane,padSignal))) continue;
945 padSignal[1] = signal;
949 // Sample the time response inside the drift region
950 // + additional time bins before and after.
951 // The sampling is done always in the middle of the time bin
952 for (Int_t iTimeBin = TMath::Max(timeE-timeTRDbeg, -nTimeBefore)
953 ;iTimeBin < TMath::Min(timeE+timeTRDend,nTimeMax+nTimeAfter )
956 // Apply the time response
957 Float_t timeResponse = 1.0;
958 Float_t crossTalk = 0.0;
959 Float_t time = (iTimeBin - timeE) * timeBinSize + timeOffset;
961 timeResponse = fPar->TimeResponse(time);
964 crossTalk = fPar->CrossTalk(time);
971 for (iPad = 0; iPad < kNpad; iPad++) {
973 Int_t colPos = colE + iPad - 1;
974 if (colPos < 0) continue;
975 if (colPos >= nColMax) break;
978 // Note: The time bin number is shifted by nTimeBefore to avoid negative
979 // time bins. This has to be subtracted later.
980 Int_t iCurrentTimeBin = iTimeBin + nTimeBefore;
981 signalOld[iPad] = signals->GetDataUnchecked(rowE,colPos,iCurrentTimeBin);
982 if( colPos != colE ) {
983 signalOld[iPad] += padSignal[iPad] * (timeResponse + crossTalk);
986 signalOld[iPad] += padSignal[iPad] * timeResponse;
988 signals->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,signalOld[iPad]);
990 // Store the track index in the dictionary
991 // Note: We store index+1 in order to allow the array to be compressed
992 if ((signalOld[iPad] > 0) && (!fSimpleSim)) {
993 for (iDict = 0; iDict < kNDict; iDict++) {
994 Int_t oldTrack = dictionary[iDict]->GetDataUnchecked(rowE
997 if (oldTrack == track+1) break;
999 dictionary[iDict]->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,track+1);
1007 } // Loop: time bins
1009 } // Loop: electrons of a single hit
1011 } // If: detector and test hit
1013 hit = (AliTRDhit *) fTRD->NextHit();
1015 } // Loop: hits of one primary track
1017 } // Loop: primary tracks
1020 printf("<AliTRDdigitizer::MakeDigits> ");
1021 printf("Finished analyzing %d hits\n",countHits);
1024 // The coupling factor
1025 Float_t coupling = fPar->GetPadCoupling()
1026 * fPar->GetTimeCoupling();
1028 // The conversion factor
1029 Float_t convert = kEl2fC
1030 * fPar->GetChipGain();
1032 // Loop through all chambers to finalize the digits
1034 Int_t iDetEnd = AliTRDgeometry::Ndet();
1036 iDetBeg = fSimpleDet;
1037 iDetEnd = iDetBeg + 1;
1039 for (Int_t iDet = iDetBeg; iDet < iDetEnd; iDet++) {
1041 Int_t plane = fGeo->GetPlane(iDet);
1042 Int_t sector = fGeo->GetSector(iDet);
1043 Int_t chamber = fGeo->GetChamber(iDet);
1044 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1045 Int_t nColMax = fPar->GetColMax(plane);
1046 Int_t nTimeMax = fPar->GetTimeMax();
1047 Int_t nTimeTotal = fPar->GetTimeTotal();
1049 Double_t *inADC = new Double_t[nTimeTotal];
1050 Double_t *outADC = new Double_t[nTimeTotal];
1053 printf("<AliTRDdigitizer::MakeDigits> ");
1054 printf("Digitization for chamber %d\n",iDet);
1057 // Add a container for the digits of this detector
1058 digits = fDigitsManager->GetDigits(iDet);
1059 // Allocate memory space for the digits buffer
1060 if (digits->GetNtime() == 0) {
1061 digits->Allocate(nRowMax,nColMax,nTimeTotal);
1063 else if (fSimpleSim) {
1067 // Get the signal container
1068 signals = (AliTRDdataArrayF *) signalsArray->At(iDet);
1069 if (signals->GetNtime() == 0) {
1070 // Create missing containers
1071 signals->Allocate(nRowMax,nColMax,nTimeTotal);
1074 // Expand the container if neccessary
1075 if (fCompress) signals->Expand();
1077 // Create the missing dictionary containers
1079 for (iDict = 0; iDict < kNDict; iDict++) {
1080 dictionary[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1081 if (dictionary[iDict]->GetNtime() == 0) {
1082 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1089 // Don't create noise in detectors that are switched off
1090 if (CheckDetector(plane,chamber,sector)) {
1092 // Create the digits for this chamber
1093 for (iRow = 0; iRow < nRowMax; iRow++ ) {
1094 for (iCol = 0; iCol < nColMax; iCol++ ) {
1096 // Create summable digits
1099 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1100 Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1101 signalAmp *= fSDigitsScale;
1102 signalAmp = TMath::Min(signalAmp,(Float_t) 1.0e9);
1103 Int_t adc = (Int_t) signalAmp;
1104 if (adc > 0) nDigits++;
1105 digits->SetDataUnchecked(iRow,iCol,iTime,adc);
1109 // Create normal digits
1112 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1113 Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1114 // Pad and time coupling
1115 signalAmp *= coupling;
1116 // Add the noise, starting from minus ADC baseline in electrons
1117 Double_t baselineEl = fPar->GetADCbaseline() * (fPar->GetADCinRange()
1118 / fPar->GetADCoutRange())
1120 signalAmp = TMath::Max((Double_t) gRandom->Gaus(signalAmp,fPar->GetNoise())
1123 signalAmp *= convert;
1124 // Add ADC baseline in mV
1125 signalAmp += fPar->GetADCbaseline() * (fPar->GetADCinRange()
1126 / fPar->GetADCoutRange());
1127 // Convert to ADC counts. Set the overflow-bit fADCoutRange if the
1128 // signal is larger than fADCinRange
1130 if (signalAmp >= fPar->GetADCinRange()) {
1131 adc = ((Int_t) fPar->GetADCoutRange());
1134 adc = ((Int_t) (signalAmp * (fPar->GetADCoutRange()
1135 / fPar->GetADCinRange())));
1138 outADC[iTime] = adc;
1141 // Apply the tail cancelation via the digital filter
1143 DeConvExp(inADC,outADC,nTimeTotal,fPar->GetTCnexp());
1146 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1147 // Store the amplitude of the digit if above threshold
1148 if (outADC[iTime] > fPar->GetADCthreshold()) {
1150 printf(" iRow = %d, iCol = %d, iTime = %d, adc = %f\n"
1151 ,iRow,iCol,iTime,outADC[iTime]);
1154 digits->SetDataUnchecked(iRow,iCol,iTime,((Int_t) outADC[iTime]));
1165 // Compress the arrays
1168 digits->Compress(1,0);
1169 for (iDict = 0; iDict < kNDict; iDict++) {
1170 dictionary[iDict]->Compress(1,0);
1173 totalSizeDigits += digits->GetSize();
1174 totalSizeDict0 += dictionary[0]->GetSize();
1175 totalSizeDict1 += dictionary[1]->GetSize();
1176 totalSizeDict2 += dictionary[2]->GetSize();
1178 Float_t nPixel = nRowMax * nColMax * nTimeMax;
1180 printf("<AliTRDdigitizer::MakeDigits> ");
1181 printf("Found %d digits in detector %d (%3.0f).\n"
1183 ,100.0 * ((Float_t) nDigits) / nPixel);
1186 if (fCompress) signals->Compress(1,0);
1196 delete signalsArray;
1201 printf("<AliTRDdigitizer::MakeDigits> ");
1202 printf("Total number of analyzed hits = %d\n",countHits);
1204 printf("<AliTRDdigitizer::MakeDigits> ");
1205 printf("Total digits data size = %d, %d, %d, %d\n",totalSizeDigits
1216 //_____________________________________________________________________________
1217 void AliTRDdigitizer::AddSDigitsManager(AliTRDdigitsManager *man)
1220 // Add a digits manager for s-digits to the input list.
1223 fSDigitsManagerList->Add(man);
1227 //_____________________________________________________________________________
1228 void AliTRDdigitizer::DeleteSDigitsManager()
1231 // Removes digits manager from the input list.
1234 fSDigitsManagerList->Delete();
1238 //_____________________________________________________________________________
1239 Bool_t AliTRDdigitizer::ConvertSDigits()
1242 // Converts s-digits to normal digits
1245 // Number of track dictionary arrays
1246 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1248 // Converts number of electrons to fC
1249 const Double_t kEl2fC = 1.602E-19 * 1.0E15;
1257 fPar = new AliTRDparameter("TRDparameter","Standard parameter");
1259 printf("<AliTRDdigitizer::ConvertSDigits> ");
1260 printf("Create the default parameter object\n");
1264 Double_t sDigitsScale = 1.0 / GetSDigitsScale();
1265 Double_t noise = fPar->GetNoise();
1266 Double_t padCoupling = fPar->GetPadCoupling();
1267 Double_t timeCoupling = fPar->GetTimeCoupling();
1268 Double_t chipGain = fPar->GetChipGain();
1269 Double_t coupling = padCoupling * timeCoupling;
1270 Double_t convert = kEl2fC * chipGain;
1271 Double_t adcInRange = fPar->GetADCinRange();
1272 Double_t adcOutRange = fPar->GetADCoutRange();
1273 Int_t adcThreshold = fPar->GetADCthreshold();
1274 Int_t adcBaseline = fPar->GetADCbaseline();
1276 AliTRDdataArrayI *digitsIn;
1277 AliTRDdataArrayI *digitsOut;
1278 AliTRDdataArrayI *dictionaryIn[kNDict];
1279 AliTRDdataArrayI *dictionaryOut[kNDict];
1281 // Loop through the detectors
1282 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1285 printf("<AliTRDdigitizer::ConvertSDigits> ");
1286 printf("Convert detector %d to digits.\n",iDet);
1289 Int_t plane = fGeo->GetPlane(iDet);
1290 Int_t sector = fGeo->GetSector(iDet);
1291 Int_t chamber = fGeo->GetChamber(iDet);
1292 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1293 Int_t nColMax = fPar->GetColMax(plane);
1294 Int_t nTimeTotal = fPar->GetTimeTotal();
1296 Double_t *inADC = new Double_t[nTimeTotal];
1297 Double_t *outADC = new Double_t[nTimeTotal];
1299 digitsIn = fSDigitsManager->GetDigits(iDet);
1301 digitsOut = fDigitsManager->GetDigits(iDet);
1302 digitsOut->Allocate(nRowMax,nColMax,nTimeTotal);
1303 for (iDict = 0; iDict < kNDict; iDict++) {
1304 dictionaryIn[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1305 dictionaryIn[iDict]->Expand();
1306 dictionaryOut[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1307 dictionaryOut[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1310 for (iRow = 0; iRow < nRowMax; iRow++ ) {
1311 for (iCol = 0; iCol < nColMax; iCol++ ) {
1313 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1314 Double_t signal = (Double_t) digitsIn->GetDataUnchecked(iRow,iCol,iTime);
1315 signal *= sDigitsScale;
1316 // Pad and time coupling
1318 // Add the noise, starting from minus ADC baseline in electrons
1319 Double_t baselineEl = adcBaseline * (adcInRange / adcOutRange) / convert;
1320 signal = TMath::Max((Double_t) gRandom->Gaus(signal,noise),-baselineEl);
1323 // add ADC baseline in mV
1324 signal += adcBaseline * (adcInRange / adcOutRange);
1325 // Convert to ADC counts. Set the overflow-bit adcOutRange if the
1326 // signal is larger than adcInRange
1328 if (signal >= adcInRange) {
1329 adc = ((Int_t) adcOutRange);
1332 adc = ((Int_t) (signal * (adcOutRange / adcInRange)));
1335 outADC[iTime] = adc;
1338 // Apply the tail cancelation via the digital filter
1340 DeConvExp(inADC,outADC,nTimeTotal,fPar->GetTCnexp());
1343 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1344 // Store the amplitude of the digit if above threshold
1345 if (outADC[iTime] > adcThreshold) {
1346 digitsOut->SetDataUnchecked(iRow,iCol,iTime,((Int_t) outADC[iTime]));
1347 // Copy the dictionary
1348 for (iDict = 0; iDict < kNDict; iDict++) {
1349 Int_t track = dictionaryIn[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1350 dictionaryOut[iDict]->SetDataUnchecked(iRow,iCol,iTime,track);
1359 digitsIn->Compress(1,0);
1360 digitsOut->Compress(1,0);
1361 for (iDict = 0; iDict < kNDict; iDict++) {
1362 dictionaryIn[iDict]->Compress(1,0);
1363 dictionaryOut[iDict]->Compress(1,0);
1376 //_____________________________________________________________________________
1377 Bool_t AliTRDdigitizer::MergeSDigits()
1380 // Merges the input s-digits:
1381 // - The amplitude of the different inputs are summed up.
1382 // - Of the track IDs from the input dictionaries only one is
1383 // kept for each input. This works for maximal 3 different merged inputs.
1386 // Number of track dictionary arrays
1387 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1390 fPar = new AliTRDparameter("TRDparameter","Standard parameter");
1392 printf("<AliTRDdigitizer::MergeSDigits> ");
1393 printf("Create the default parameter object\n");
1400 AliTRDdataArrayI *digitsA;
1401 AliTRDdataArrayI *digitsB;
1402 AliTRDdataArrayI *dictionaryA[kNDict];
1403 AliTRDdataArrayI *dictionaryB[kNDict];
1405 // Get the first s-digits
1406 fSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->First();
1407 if (!fSDigitsManager) return kFALSE;
1409 // Loop through the other sets of s-digits
1410 AliTRDdigitsManager *mergeSDigitsManager;
1411 mergeSDigitsManager = (AliTRDdigitsManager *)
1412 fSDigitsManagerList->After(fSDigitsManager);
1415 if (mergeSDigitsManager) {
1416 printf("<AliTRDdigitizer::MergeSDigits> ");
1417 printf("Merge %d input files.\n",fSDigitsManagerList->GetSize());
1420 printf("<AliTRDdigitizer::MergeSDigits> ");
1421 printf("Only one input file.\n");
1426 while (mergeSDigitsManager) {
1430 // Loop through the detectors
1431 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1433 Int_t plane = fGeo->GetPlane(iDet);
1434 Int_t sector = fGeo->GetSector(iDet);
1435 Int_t chamber = fGeo->GetChamber(iDet);
1436 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1437 Int_t nColMax = fPar->GetColMax(plane);
1438 Int_t nTimeTotal = fPar->GetTimeTotal();
1440 // Loop through the pixels of one detector and add the signals
1441 digitsA = fSDigitsManager->GetDigits(iDet);
1442 digitsB = mergeSDigitsManager->GetDigits(iDet);
1445 for (iDict = 0; iDict < kNDict; iDict++) {
1446 dictionaryA[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1447 dictionaryB[iDict] = mergeSDigitsManager->GetDictionary(iDet,iDict);
1448 dictionaryA[iDict]->Expand();
1449 dictionaryB[iDict]->Expand();
1452 // Merge only detectors that contain a signal
1453 Bool_t doMerge = kTRUE;
1454 if (fMergeSignalOnly) {
1455 if (digitsA->GetOverThreshold(0) == 0) {
1463 printf("<AliTRDdigitizer::MergeSDigits> ");
1464 printf("Merge detector %d of input no.%d\n",iDet,iMerge+1);
1467 for (Int_t iRow = 0; iRow < nRowMax; iRow++ ) {
1468 for (Int_t iCol = 0; iCol < nColMax; iCol++ ) {
1469 for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {
1471 // Add the amplitudes of the summable digits
1472 Int_t ampA = digitsA->GetDataUnchecked(iRow,iCol,iTime);
1473 Int_t ampB = digitsB->GetDataUnchecked(iRow,iCol,iTime);
1475 digitsA->SetDataUnchecked(iRow,iCol,iTime,ampA);
1477 // Add the mask to the track id if defined.
1478 for (iDict = 0; iDict < kNDict; iDict++) {
1479 Int_t trackB = dictionaryB[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1480 if ((fMasks) && (trackB > 0)) {
1481 for (jDict = 0; jDict < kNDict; jDict++) {
1482 Int_t trackA = dictionaryA[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1484 trackA = trackB + fMasks[iMerge];
1485 dictionaryA[iDict]->SetDataUnchecked(iRow,iCol,iTime,trackA);
1498 digitsA->Compress(1,0);
1499 digitsB->Compress(1,0);
1500 for (iDict = 0; iDict < kNDict; iDict++) {
1501 dictionaryA[iDict]->Compress(1,0);
1502 dictionaryB[iDict]->Compress(1,0);
1508 // The next set of s-digits
1509 mergeSDigitsManager = (AliTRDdigitsManager *)
1510 fSDigitsManagerList->After(mergeSDigitsManager);
1518 //_____________________________________________________________________________
1519 Bool_t AliTRDdigitizer::SDigits2Digits()
1522 // Merges the input s-digits and converts them to normal digits
1525 if (!MergeSDigits()) return kFALSE;
1527 return ConvertSDigits();
1531 //_____________________________________________________________________________
1532 Bool_t AliTRDdigitizer::CheckDetector(Int_t plane, Int_t chamber, Int_t sector)
1535 // Checks whether a detector is enabled
1538 if (fSimpleSim) return kTRUE;
1540 if ((fTRD->GetSensChamber() >= 0) &&
1541 (fTRD->GetSensChamber() != chamber)) return kFALSE;
1542 if ((fTRD->GetSensPlane() >= 0) &&
1543 (fTRD->GetSensPlane() != plane)) return kFALSE;
1544 if ( fTRD->GetSensSector() >= 0) {
1545 Int_t sens1 = fTRD->GetSensSector();
1546 Int_t sens2 = sens1 + fTRD->GetSensSectorRange();
1547 sens2 -= ((Int_t) (sens2 / AliTRDgeometry::Nsect()))
1548 * AliTRDgeometry::Nsect();
1549 if (sens1 < sens2) {
1550 if ((sector < sens1) || (sector >= sens2)) return kFALSE;
1553 if ((sector < sens1) && (sector >= sens2)) return kFALSE;
1561 //_____________________________________________________________________________
1562 Bool_t AliTRDdigitizer::WriteDigits() const
1565 // Writes out the TRD-digits and the dictionaries
1568 // Store the digits and the dictionary in the tree
1569 return fDigitsManager->WriteDigits();
1573 //_____________________________________________________________________________
1574 void AliTRDdigitizer::DeConvExp(Double_t *source, Double_t *target
1575 , Int_t n, Int_t nexp)
1578 // Does the deconvolution by the digital filter.
1580 // Author: Marcus Gutfleisch, KIP Heidelberg
1581 // Optimized for: New TRF from Venelin Angelov, simulated with CADENCE
1582 // Pad-ground capacitance = 25 pF
1583 // Pad-pad cross talk capacitance = 6 pF
1584 // For 10 MHz digitization, corresponding to 20 time bins
1585 // in the drift region
1589 Double_t coefficients[2];
1591 /* initialize (coefficient = alpha, rates = lambda) */
1594 rates[0] = 0.466998;
1596 coefficients[0] = 1.0;
1599 rates[0] = 0.8988162;
1600 coefficients[0] = 0.11392069;
1601 rates[1] = 0.3745688;
1602 coefficients[1] = 0.8860793;
1604 Float_t sumc = coefficients[0]+coefficients[1];
1605 coefficients[0] /= sumc;
1606 coefficients[1] /= sumc;
1610 Double_t reminder[2];
1611 Double_t correction, result;
1613 /* attention: computation order is important */
1615 for ( k = 0; k < nexp; k++ ) reminder[k]=0.0;
1617 for ( i = 0; i < n; i++ ) {
1618 result = ( source[i] - correction ); /* no rescaling */
1621 for ( k = 0; k < nexp; k++ ) reminder[k] = rates[k]
1622 * ( reminder[k] + coefficients[k] * result);
1625 for ( k = 0; k < nexp; k++ ) correction += reminder[k];
1630 //_____________________________________________________________________________
1631 void AliTRDdigitizer::InitOutput(TFile *file, Int_t iEvent)
1634 // Initializes the output branches
1638 fDigitsManager->MakeTreeAndBranches(file,iEvent);