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.40 2002/10/14 14:57:43 hristov
19 Merging the VirtualMC branch to the main development branch (HEAD)
21 Revision 1.33.6.3 2002/10/11 07:26:37 hristov
22 Updating VirtualMC to v3-09-02
24 Revision 1.39 2002/10/08 20:46:12 cblume
25 Do coupling factors before noise is applied
27 Revision 1.38 2002/04/30 08:30:40 cblume
28 gAlice now only read by AliRunDigitizer. Therefore it is just deleted in AliTRDmerge.C
30 Revision 1.37 2002/04/29 11:50:47 cblume
31 Change initialization of gAlice in the merging case
33 Revision 1.36 2002/04/12 12:13:23 cblume
36 Revision 1.35 2002/03/28 14:59:07 cblume
39 Revision 1.34 2002/03/25 20:00:44 cblume
40 Introduce parameter class and regions of interest for merging
42 Revision 1.33 2002/02/12 17:32:03 cblume
43 Rearrange the deleting of the list of sdigitsmanager
45 Revision 1.32 2002/02/12 16:07:21 cblume
48 Revision 1.31 2002/02/11 14:27:11 cblume
49 New pad plane design, new TRF+PRF, tail cancelation, cross talk
51 Revision 1.30 2001/11/19 08:44:08 cblume
52 Fix bugs reported by Rene
54 Revision 1.29 2001/11/14 19:44:25 hristov
55 Numeric const casted (Alpha)
57 Revision 1.28 2001/11/14 16:35:58 cblume
58 Inherits now from AliDetector
60 Revision 1.27 2001/11/14 10:50:45 cblume
61 Changes in digits IO. Add merging of summable digits
63 Revision 1.26 2001/11/06 17:19:41 cblume
64 Add detailed geometry and simple simulator
66 Revision 1.25 2001/06/27 09:54:44 cblume
67 Moved fField initialization to InitDetector()
69 Revision 1.24 2001/05/21 16:45:47 hristov
70 Last minute changes (C.Blume)
72 Revision 1.23 2001/05/07 08:04:48 cblume
73 New TRF and PRF. Speedup of the code. Digits from amplification region included
75 Revision 1.22 2001/03/30 14:40:14 cblume
76 Update of the digitization parameter
78 Revision 1.21 2001/03/13 09:30:35 cblume
79 Update of digitization. Moved digit branch definition to AliTRD
81 Revision 1.20 2001/02/25 20:19:00 hristov
82 Minor correction: loop variable declared only once for HP, Sun
84 Revision 1.19 2001/02/14 18:22:26 cblume
85 Change in the geometry of the padplane
87 Revision 1.18 2001/01/26 19:56:57 hristov
88 Major upgrade of AliRoot code
90 Revision 1.17 2000/12/08 12:53:27 cblume
91 Change in Copy() function for HP-compiler
93 Revision 1.16 2000/12/07 12:20:46 cblume
94 Go back to array compression. Use sampled PRF to speed up digitization
96 Revision 1.15 2000/11/23 14:34:08 cblume
97 Fixed bug in expansion routine of arrays (initialize buffers properly)
99 Revision 1.14 2000/11/20 08:54:44 cblume
100 Switch off compression as default
102 Revision 1.13 2000/11/10 14:57:52 cblume
103 Changes in the geometry constants for the DEC compiler
105 Revision 1.12 2000/11/01 14:53:20 cblume
106 Merge with TRD-develop
108 Revision 1.1.4.9 2000/10/26 17:00:22 cblume
109 Fixed bug in CheckDetector()
111 Revision 1.1.4.8 2000/10/23 13:41:35 cblume
112 Added protection against Log(0) in the gas gain calulation
114 Revision 1.1.4.7 2000/10/17 02:27:34 cblume
115 Get rid of global constants
117 Revision 1.1.4.6 2000/10/16 01:16:53 cblume
118 Changed timebin 0 to be the one closest to the readout
120 Revision 1.1.4.5 2000/10/15 23:34:29 cblume
121 Faster version of the digitizer
123 Revision 1.1.4.4 2000/10/06 16:49:46 cblume
126 Revision 1.1.4.3 2000/10/04 16:34:58 cblume
127 Replace include files by forward declarations
129 Revision 1.1.4.2 2000/09/22 14:41:10 cblume
130 Bug fix in PRF. Included time response. New structure
132 Revision 1.10 2000/10/05 07:27:53 cblume
133 Changes in the header-files by FCA
135 Revision 1.9 2000/10/02 21:28:19 fca
136 Removal of useless dependecies via forward declarations
138 Revision 1.8 2000/06/09 11:10:07 cblume
139 Compiler warnings and coding conventions, next round
141 Revision 1.7 2000/06/08 18:32:58 cblume
142 Make code compliant to coding conventions
144 Revision 1.6 2000/06/07 16:27:32 cblume
145 Try to remove compiler warnings on Sun and HP
147 Revision 1.5 2000/05/09 16:38:57 cblume
148 Removed PadResponse(). Merge problem
150 Revision 1.4 2000/05/08 15:53:45 cblume
151 Resolved merge conflict
153 Revision 1.3 2000/04/28 14:49:27 cblume
154 Only one declaration of iDict in MakeDigits()
156 Revision 1.1.4.1 2000/05/08 14:42:04 cblume
157 Introduced AliTRDdigitsManager
159 Revision 1.1 2000/02/28 19:00:13 cblume
164 ///////////////////////////////////////////////////////////////////////////////
166 // Creates and handles digits from TRD hits //
167 // Author: C. Blume (C.Blume@gsi.de) //
169 // The following effects are included: //
172 // - Gas gain including fluctuations //
173 // - Pad-response (simple Gaussian approximation) //
174 // - Time-response //
175 // - Electronics noise //
176 // - Electronics gain //
178 // - ADC threshold //
179 // The corresponding parameter can be adjusted via the various //
180 // Set-functions. If these parameters are not explicitly set, default //
181 // values are used (see Init-function). //
182 // As an example on how to use this class to produce digits from hits //
183 // have a look at the macro hits2digits.C //
184 // The production of summable digits is demonstrated in hits2sdigits.C //
185 // and the subsequent conversion of the s-digits into normal digits is //
186 // explained in sdigits2digits.C. //
188 ///////////////////////////////////////////////////////////////////////////////
204 #include "AliRunDigitizer.h"
207 #include "AliTRDhit.h"
208 #include "AliTRDdigitizer.h"
209 #include "AliTRDdataArrayI.h"
210 #include "AliTRDdataArrayF.h"
211 #include "AliTRDsegmentArray.h"
212 #include "AliTRDdigitsManager.h"
213 #include "AliTRDgeometry.h"
214 #include "AliTRDparameter.h"
216 ClassImp(AliTRDdigitizer)
218 //_____________________________________________________________________________
219 AliTRDdigitizer::AliTRDdigitizer()
222 // AliTRDdigitizer default constructor
227 fSDigitsManagerList = 0;
238 fMergeSignalOnly = kFALSE;
244 //_____________________________________________________________________________
245 AliTRDdigitizer::AliTRDdigitizer(const Text_t *name, const Text_t *title)
246 :AliDigitizer(name,title)
249 // AliTRDdigitizer constructor
254 fSDigitsManagerList = 0;
264 fMergeSignalOnly = kFALSE;
268 // For the summable digits
269 fSDigitsScale = 100.;
273 //_____________________________________________________________________________
274 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager
275 , const Text_t *name, const Text_t *title)
276 :AliDigitizer(manager,name,title)
279 // AliTRDdigitizer constructor
284 fSDigitsManagerList = 0;
294 fMergeSignalOnly = kFALSE;
298 // For the summable digits
299 fSDigitsScale = 100.;
303 //_____________________________________________________________________________
304 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager)
305 :AliDigitizer(manager,"AliTRDdigitizer","TRD digitizer")
308 // AliTRDdigitizer constructor
313 fSDigitsManagerList = 0;
323 fMergeSignalOnly = kFALSE;
327 // For the summable digits
328 fSDigitsScale = 100.;
332 //_____________________________________________________________________________
333 AliTRDdigitizer::AliTRDdigitizer(const AliTRDdigitizer &d)
336 // AliTRDdigitizer copy constructor
339 ((AliTRDdigitizer &) d).Copy(*this);
343 //_____________________________________________________________________________
344 AliTRDdigitizer::~AliTRDdigitizer()
347 // AliTRDdigitizer destructor
356 if (fDigitsManager) {
357 delete fDigitsManager;
361 if (fSDigitsManager) {
362 delete fSDigitsManager;
366 if (fSDigitsManagerList) {
367 delete fSDigitsManagerList;
368 fSDigitsManagerList = 0;
378 //_____________________________________________________________________________
379 AliTRDdigitizer &AliTRDdigitizer::operator=(const AliTRDdigitizer &d)
382 // Assignment operator
385 if (this != &d) ((AliTRDdigitizer &) d).Copy(*this);
390 //_____________________________________________________________________________
391 void AliTRDdigitizer::Copy(TObject &d)
397 ((AliTRDdigitizer &) d).fInputFile = 0;
398 ((AliTRDdigitizer &) d).fSDigitsManagerList = 0;
399 ((AliTRDdigitizer &) d).fSDigitsManager = 0;
400 ((AliTRDdigitizer &) d).fDigitsManager = 0;
401 ((AliTRDdigitizer &) d).fTRD = 0;
402 ((AliTRDdigitizer &) d).fGeo = 0;
403 ((AliTRDdigitizer &) d).fMasks = 0;
404 ((AliTRDdigitizer &) d).fEvent = 0;
405 ((AliTRDdigitizer &) d).fPar = 0;
406 ((AliTRDdigitizer &) d).fCompress = fCompress;
407 ((AliTRDdigitizer &) d).fDebug = fDebug ;
408 ((AliTRDdigitizer &) d).fSDigits = fSDigits;
409 ((AliTRDdigitizer &) d).fSDigitsScale = fSDigitsScale;
410 ((AliTRDdigitizer &) d).fMergeSignalOnly = fMergeSignalOnly;
411 ((AliTRDdigitizer &) d).fSimpleSim = fSimpleSim;
412 ((AliTRDdigitizer &) d).fSimpleDet = fSimpleDet;
416 //_____________________________________________________________________________
417 void AliTRDdigitizer::Exec(Option_t* option)
420 // Executes the merging
425 AliTRDdigitsManager *sdigitsManager;
427 TString optionString = option;
428 if (optionString.Contains("deb")) {
430 if (optionString.Contains("2")) {
433 printf("<AliTRDdigitizer::Exec> ");
434 printf("Called with debug option %d\n",fDebug);
437 // The AliRoot file is already connected by the manager
440 printf("<AliTRDdigitizer::Exec> ");
441 printf("AliRun object found on file.\n");
445 printf("<AliTRDdigitizer::Exec> ");
446 printf("Could not find AliRun object.\n");
450 Int_t nInput = fManager->GetNinputs();
451 fMasks = new Int_t[nInput];
452 for (iInput = 0; iInput < nInput; iInput++) {
453 fMasks[iInput] = fManager->GetMask(iInput);
459 for (iInput = 0; iInput < nInput; iInput++) {
462 printf("<AliTRDdigitizer::Exec> ");
463 printf("Add input stream %d\n",iInput);
466 // check if the input tree exists
467 if (!fManager->GetInputTreeTRDS(iInput)) {
468 printf("<AliTRDdigitizer::Exec> ");
469 printf("Input stream %d does not exist\n",iInput);
473 // Read the s-digits via digits manager
474 sdigitsManager = new AliTRDdigitsManager();
475 sdigitsManager->SetDebug(fDebug);
476 sdigitsManager->SetSDigits(kTRUE);
477 sdigitsManager->ReadDigits(fManager->GetInputTreeTRDS(iInput));
479 // Add the s-digits to the input list
480 AddSDigitsManager(sdigitsManager);
484 // Convert the s-digits to normal digits
486 printf("<AliTRDdigitizer::Exec> ");
487 printf("Do the conversion\n");
493 printf("<AliTRDdigitizer::Exec> ");
494 printf("Write the digits\n");
496 fDigitsManager->MakeBranch(fManager->GetTreeDTRD());
497 fDigitsManager->WriteDigits();
499 printf("<AliTRDdigitizer::Exec> ");
503 DeleteSDigitsManager();
507 //_____________________________________________________________________________
508 Bool_t AliTRDdigitizer::Open(const Char_t *file, Int_t nEvent)
511 // Opens a ROOT-file with TRD-hits and reads in the hit-tree
514 // Connect the AliRoot file containing Geometry, Kine, and Hits
515 fInputFile = (TFile *) gROOT->GetListOfFiles()->FindObject(file);
518 printf("<AliTRDdigitizer::Open> ");
519 printf("Open the AliROOT-file %s.\n",file);
521 fInputFile = new TFile(file,"UPDATE");
525 printf("<AliTRDdigitizer::Open> ");
526 printf("%s is already open.\n",file);
530 gAlice = (AliRun *) fInputFile->Get("gAlice");
533 printf("<AliTRDdigitizer::Open> ");
534 printf("AliRun object found on file.\n");
538 printf("<AliTRDdigitizer::Open> ");
539 printf("Could not find AliRun object.\n");
545 // Import the Trees for the event nEvent in the file
546 Int_t nparticles = gAlice->GetEvent(fEvent);
547 if (nparticles <= 0) {
548 printf("<AliTRDdigitizer::Open> ");
549 printf("No entries in the trees for event %d.\n",fEvent);
553 if (InitDetector()) {
562 //_____________________________________________________________________________
563 Bool_t AliTRDdigitizer::InitDetector()
566 // Sets the pointer to the TRD detector and the geometry
569 // Get the pointer to the detector class and check for version 1
570 fTRD = (AliTRD *) gAlice->GetDetector("TRD");
572 printf("<AliTRDdigitizer::InitDetector> ");
573 printf("No TRD module found\n");
576 if (fTRD->IsVersion() != 1) {
577 printf("<AliTRDdigitizer::InitDetector> ");
578 printf("TRD must be version 1 (slow simulator).\n");
583 fGeo = fTRD->GetGeometry();
585 printf("<AliTRDdigitizer::InitDetector> ");
586 printf("Geometry version %d\n",fGeo->IsVersion());
589 // Create a digits manager
590 fDigitsManager = new AliTRDdigitsManager();
591 fDigitsManager->SetSDigits(fSDigits);
592 fDigitsManager->CreateArrays();
593 fDigitsManager->SetEvent(fEvent);
594 fDigitsManager->SetDebug(fDebug);
596 // The list for the input s-digits manager to be merged
597 fSDigitsManagerList = new TList();
603 //_____________________________________________________________________________
604 Bool_t AliTRDdigitizer::MakeBranch(const Char_t *file) const
607 // Create the branches for the digits array
610 return fDigitsManager->MakeBranch(file);
614 //_____________________________________________________________________________
615 Bool_t AliTRDdigitizer::MakeDigits()
621 ///////////////////////////////////////////////////////////////
623 ///////////////////////////////////////////////////////////////
625 // Converts number of electrons to fC
626 const Double_t kEl2fC = 1.602E-19 * 1.0E15;
628 ///////////////////////////////////////////////////////////////
630 // Number of pads included in the pad response
631 const Int_t kNpad = 3;
633 // Number of track dictionary arrays
634 const Int_t kNDict = AliTRDdigitsManager::kNDict;
636 // Half the width of the amplification region
637 const Float_t kAmWidth = AliTRDgeometry::AmThick() / 2.;
639 Int_t iRow, iCol, iTime, iPad;
643 Int_t totalSizeDigits = 0;
644 Int_t totalSizeDict0 = 0;
645 Int_t totalSizeDict1 = 0;
646 Int_t totalSizeDict2 = 0;
648 Int_t timeTRDbeg = 0;
649 Int_t timeTRDend = 1;
654 Float_t padSignal[kNpad];
655 Float_t signalOld[kNpad];
657 AliTRDdataArrayF *signals = 0;
658 AliTRDdataArrayI *digits = 0;
659 AliTRDdataArrayI *dictionary[kNDict];
661 // Create a default parameter class if none is defined
663 fPar = new AliTRDparameter("TRDparameter","Standard TRD parameter");
665 printf("<AliTRDdigitizer::MakeDigits> ");
666 printf("Create the default parameter object\n");
670 // Create a container for the amplitudes
671 AliTRDsegmentArray *signalsArray
672 = new AliTRDsegmentArray("AliTRDdataArrayF"
673 ,AliTRDgeometry::Ndet());
676 timeTRDbeg = ((Int_t) (-fPar->GetTRFlo() / fPar->GetTimeBinSize())) - 1;
677 timeTRDend = ((Int_t) ( fPar->GetTRFhi() / fPar->GetTimeBinSize())) - 1;
679 printf("<AliTRDdigitizer::MakeDigits> ");
680 printf("Sample the TRF between -%d and %d\n",timeTRDbeg,timeTRDend);
684 Float_t elAttachProp = fPar->GetElAttachProp() / 100.;
687 printf("<AliTRDdigitizer::MakeDigits> ");
688 printf("No geometry defined\n");
693 printf("<AliTRDdigitizer::MakeDigits> ");
694 printf("Start creating digits.\n");
697 // Get the pointer to the hit tree
698 TTree *hitTree = gAlice->TreeH();
700 // Get the number of entries in the hit tree
701 // (Number of primary particles creating a hit somewhere)
704 nTrack = (Int_t) hitTree->GetEntries();
706 printf("<AliTRDdigitizer::MakeDigits> ");
707 printf("Found %d primary particles\n",nTrack);
711 Int_t detectorOld = -1;
714 // Loop through all entries in the tree
715 for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
719 nBytes += hitTree->GetEvent(iTrack);
722 // Loop through the TRD hits
724 AliTRDhit *hit = (AliTRDhit *) fTRD->FirstHit(-1);
733 Float_t q = hit->GetCharge();
734 Int_t track = hit->Track();
735 Int_t detector = hit->GetDetector();
736 Int_t plane = fGeo->GetPlane(detector);
737 Int_t sector = fGeo->GetSector(detector);
738 Int_t chamber = fGeo->GetChamber(detector);
739 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
740 Int_t nColMax = fPar->GetColMax(plane);
741 Int_t nTimeMax = fPar->GetTimeMax();
742 Int_t nTimeBefore = fPar->GetTimeBefore();
743 Int_t nTimeAfter = fPar->GetTimeAfter();
744 Int_t nTimeTotal = fPar->GetTimeTotal();
745 Float_t row0 = fPar->GetRow0(plane,chamber,sector);
746 Float_t col0 = fPar->GetCol0(plane);
747 Float_t time0 = fPar->GetTime0(plane);
748 Float_t rowPadSize = fPar->GetRowPadSize(plane,chamber,sector);
749 Float_t colPadSize = fPar->GetColPadSize(plane);
750 Float_t timeBinSize = fPar->GetTimeBinSize();
751 Float_t divideRow = 1.0 / rowPadSize;
752 Float_t divideCol = 1.0 / colPadSize;
753 Float_t divideTime = 1.0 / timeBinSize;
756 printf("Analyze hit no. %d ",iHit);
757 printf("-----------------------------------------------------------\n");
759 printf("plane = %d, sector = %d, chamber = %d\n"
760 ,plane,sector,chamber);
761 printf("nRowMax = %d, nColMax = %d, nTimeMax = %d\n"
762 ,nRowMax,nColMax,nTimeMax);
763 printf("nTimeBefore = %d, nTimeAfter = %d, nTimeTotal = %d\n"
764 ,nTimeBefore,nTimeAfter,nTimeTotal);
765 printf("row0 = %f, col0 = %f, time0 = %f\n"
767 printf("rowPadSize = %f, colPadSize = %f, timeBinSize = %f\n"
768 ,rowPadSize,colPadSize,timeBinSize);
771 // Don't analyze test hits and switched off detectors
772 if ((CheckDetector(plane,chamber,sector)) &&
773 (((Int_t) q) != 0)) {
775 if (detector != detectorOld) {
778 printf("<AliTRDdigitizer::MakeDigits> ");
779 printf("Get new container. New det = %d, Old det = %d\n"
780 ,detector,detectorOld);
782 // Compress the old one if enabled
783 if ((fCompress) && (detectorOld > -1)) {
785 printf("<AliTRDdigitizer::MakeDigits> ");
786 printf("Compress the old container ...");
788 signals->Compress(1,0);
789 for (iDict = 0; iDict < kNDict; iDict++) {
790 dictionary[iDict]->Compress(1,0);
792 if (fDebug > 1) printf("done\n");
794 // Get the new container
795 signals = (AliTRDdataArrayF *) signalsArray->At(detector);
796 if (signals->GetNtime() == 0) {
797 // Allocate a new one if not yet existing
799 printf("<AliTRDdigitizer::MakeDigits> ");
800 printf("Allocate a new container ... ");
802 signals->Allocate(nRowMax,nColMax,nTimeTotal);
804 else if (fSimpleSim) {
805 // Clear an old one for the simple simulation
807 printf("<AliTRDdigitizer::MakeDigits> ");
808 printf("Clear a old container ... ");
813 // Expand an existing one
816 printf("<AliTRDdigitizer::MakeDigits> ");
817 printf("Expand an existing container ... ");
822 // The same for the dictionary
824 for (iDict = 0; iDict < kNDict; iDict++) {
825 dictionary[iDict] = fDigitsManager->GetDictionary(detector,iDict);
826 if (dictionary[iDict]->GetNtime() == 0) {
827 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
830 if (fCompress) dictionary[iDict]->Expand();
834 if (fDebug > 1) printf("done\n");
835 detectorOld = detector;
838 // Rotate the sectors on top of each other
845 fGeo->Rotate(detector,pos,rot);
848 // The driftlength. It is negative if the hit is in the
849 // amplification region.
850 Float_t driftlength = time0 - rot[0];
852 // Take also the drift in the amplification region into account
853 // The drift length is at the moment still the same, regardless of
854 // the position relativ to the wire. This non-isochronity needs still
855 // to be implemented.
856 Float_t driftlengthL = TMath::Abs(driftlength + kAmWidth);
857 if (fPar->ExBOn()) driftlengthL /= TMath::Sqrt(fPar->GetLorentzFactor());
859 // Loop over all electrons of this hit
860 // TR photons produce hits with negative charge
861 Int_t nEl = ((Int_t) TMath::Abs(q));
862 for (Int_t iEl = 0; iEl < nEl; iEl++) {
868 // Electron attachment
869 if (fPar->ElAttachOn()) {
870 if (gRandom->Rndm() < (driftlengthL * elAttachProp))
874 // Apply the diffusion smearing
875 if (fPar->DiffusionOn()) {
876 if (!(fPar->Diffusion(driftlengthL,xyz))) continue;
879 // Apply E x B effects (depends on drift direction)
881 if (!(fPar->ExB(driftlength+kAmWidth,xyz))) continue;
884 // The electron position after diffusion and ExB in pad coordinates
885 // The pad row (z-direction)
886 Float_t rowDist = xyz[2] - row0;
887 Int_t rowE = ((Int_t) (rowDist * divideRow));
888 if ((rowE < 0) || (rowE >= nRowMax)) continue;
889 Float_t rowOffset = ((((Float_t) rowE) + 0.5) * rowPadSize) - rowDist;
891 // The pad column (rphi-direction)
892 Float_t col0tilt = fPar->Col0Tilted(col0,rowOffset,plane);
893 Float_t colDist = xyz[1] - col0tilt;
894 Int_t colE = ((Int_t) (colDist * divideCol));
895 if ((colE < 0) || (colE >= nColMax)) continue;
896 Float_t colOffset = ((((Float_t) colE) + 0.5) * colPadSize) - colDist;
898 // The time bin (negative for hits in the amplification region)
899 // In the amplification region the electrons drift from both sides
900 // to the middle (anode wire plane)
901 Float_t timeDist = time0 - xyz[0];
902 Float_t timeOffset = 0;
906 timeE = ((Int_t) (timeDist * divideTime));
907 // The distance of the position to the middle of the timebin
908 timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) - timeDist;
911 // Difference between half of the amplification gap width and
912 // the distance to the anode wire
913 Float_t anodeDist = kAmWidth - TMath::Abs(timeDist + kAmWidth);
915 timeE = -1 * (((Int_t ) (anodeDist * divideTime)) + 1);
916 // The distance of the position to the middle of the timebin
917 timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) + anodeDist;
920 // Apply the gas gain including fluctuations
921 Float_t ggRndm = 0.0;
923 ggRndm = gRandom->Rndm();
924 } while (ggRndm <= 0);
925 Int_t signal = (Int_t) (-fPar->GetGasGain() * TMath::Log(ggRndm));
927 // Apply the pad response
929 // The distance of the electron to the center of the pad
930 // in units of pad width
931 Float_t dist = - colOffset * divideCol;
932 if (!(fPar->PadResponse(signal,dist,plane,padSignal))) continue;
936 padSignal[1] = signal;
940 // Sample the time response inside the drift region
941 // + additional time bins before and after.
942 // The sampling is done always in the middle of the time bin
943 for (Int_t iTimeBin = TMath::Max(timeE-timeTRDbeg, -nTimeBefore)
944 ;iTimeBin < TMath::Min(timeE+timeTRDend,nTimeMax+nTimeAfter )
947 // Apply the time response
948 Float_t timeResponse = 1.0;
949 Float_t crossTalk = 0.0;
950 Float_t time = (iTimeBin - timeE) * timeBinSize + timeOffset;
952 timeResponse = fPar->TimeResponse(time);
955 crossTalk = fPar->CrossTalk(time);
962 for (iPad = 0; iPad < kNpad; iPad++) {
964 Int_t colPos = colE + iPad - 1;
965 if (colPos < 0) continue;
966 if (colPos >= nColMax) break;
969 // Note: The time bin number is shifted by nTimeBefore to avoid negative
970 // time bins. This has to be subtracted later.
971 Int_t iCurrentTimeBin = iTimeBin + nTimeBefore;
972 signalOld[iPad] = signals->GetDataUnchecked(rowE,colPos,iCurrentTimeBin);
973 if( colPos != colE ) {
974 signalOld[iPad] += padSignal[iPad] * (timeResponse + crossTalk);
977 signalOld[iPad] += padSignal[iPad] * timeResponse;
979 signals->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,signalOld[iPad]);
981 // Store the track index in the dictionary
982 // Note: We store index+1 in order to allow the array to be compressed
983 if ((signalOld[iPad] > 0) && (!fSimpleSim)) {
984 for (iDict = 0; iDict < kNDict; iDict++) {
985 Int_t oldTrack = dictionary[iDict]->GetDataUnchecked(rowE
988 if (oldTrack == track+1) break;
990 dictionary[iDict]->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,track+1);
1000 } // Loop: electrons of a single hit
1002 } // If: detector and test hit
1004 hit = (AliTRDhit *) fTRD->NextHit();
1006 } // Loop: hits of one primary track
1008 } // Loop: primary tracks
1011 printf("<AliTRDdigitizer::MakeDigits> ");
1012 printf("Finished analyzing %d hits\n",countHits);
1015 // The coupling factor
1016 Float_t coupling = fPar->GetPadCoupling()
1017 * fPar->GetTimeCoupling();
1019 // The conversion factor
1020 Float_t convert = kEl2fC
1021 * fPar->GetChipGain();
1023 // Loop through all chambers to finalize the digits
1025 Int_t iDetEnd = AliTRDgeometry::Ndet();
1027 iDetBeg = fSimpleDet;
1028 iDetEnd = iDetBeg + 1;
1030 for (Int_t iDet = iDetBeg; iDet < iDetEnd; iDet++) {
1032 Int_t plane = fGeo->GetPlane(iDet);
1033 Int_t sector = fGeo->GetSector(iDet);
1034 Int_t chamber = fGeo->GetChamber(iDet);
1035 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1036 Int_t nColMax = fPar->GetColMax(plane);
1037 Int_t nTimeMax = fPar->GetTimeMax();
1038 Int_t nTimeTotal = fPar->GetTimeTotal();
1040 Double_t *inADC = new Double_t[nTimeTotal];
1041 Double_t *outADC = new Double_t[nTimeTotal];
1044 printf("<AliTRDdigitizer::MakeDigits> ");
1045 printf("Digitization for chamber %d\n",iDet);
1048 // Add a container for the digits of this detector
1049 digits = fDigitsManager->GetDigits(iDet);
1050 // Allocate memory space for the digits buffer
1051 if (digits->GetNtime() == 0) {
1052 digits->Allocate(nRowMax,nColMax,nTimeTotal);
1054 else if (fSimpleSim) {
1058 // Get the signal container
1059 signals = (AliTRDdataArrayF *) signalsArray->At(iDet);
1060 if (signals->GetNtime() == 0) {
1061 // Create missing containers
1062 signals->Allocate(nRowMax,nColMax,nTimeTotal);
1065 // Expand the container if neccessary
1066 if (fCompress) signals->Expand();
1068 // Create the missing dictionary containers
1070 for (iDict = 0; iDict < kNDict; iDict++) {
1071 dictionary[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1072 if (dictionary[iDict]->GetNtime() == 0) {
1073 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1080 // Don't create noise in detectors that are switched off
1081 if (CheckDetector(plane,chamber,sector)) {
1083 // Create the digits for this chamber
1084 for (iRow = 0; iRow < nRowMax; iRow++ ) {
1085 for (iCol = 0; iCol < nColMax; iCol++ ) {
1087 // Create summable digits
1090 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1091 Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1092 signalAmp *= fSDigitsScale;
1093 signalAmp = TMath::Min(signalAmp,(Float_t) 1.0e9);
1094 Int_t adc = (Int_t) signalAmp;
1095 if (adc > 0) nDigits++;
1096 digits->SetDataUnchecked(iRow,iCol,iTime,adc);
1100 // Create normal digits
1103 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1104 Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1105 // Pad and time coupling
1106 signalAmp *= coupling;
1108 signalAmp = TMath::Max((Double_t) gRandom->Gaus(signalAmp,fPar->GetNoise()),0.0);
1110 signalAmp *= convert;
1111 // Convert to ADC counts. Set the overflow-bit fADCoutRange if the
1112 // signal is larger than fADCinRange
1114 if (signalAmp >= fPar->GetADCinRange()) {
1115 adc = ((Int_t) fPar->GetADCoutRange());
1118 adc = ((Int_t) (signalAmp * (fPar->GetADCoutRange()
1119 / fPar->GetADCinRange())));
1122 outADC[iTime] = adc;
1125 // Apply the tail cancelation via the digital filter
1127 DeConvExp(inADC,outADC,nTimeTotal,fPar->GetTCnexp());
1130 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1131 // Store the amplitude of the digit if above threshold
1132 if (outADC[iTime] > fPar->GetADCthreshold()) {
1134 printf(" iRow = %d, iCol = %d, iTime = %d, adc = %f\n"
1135 ,iRow,iCol,iTime,outADC[iTime]);
1138 digits->SetDataUnchecked(iRow,iCol,iTime,((Int_t) outADC[iTime]));
1149 // Compress the arrays
1151 digits->Compress(1,0);
1152 for (iDict = 0; iDict < kNDict; iDict++) {
1153 dictionary[iDict]->Compress(1,0);
1156 totalSizeDigits += digits->GetSize();
1157 totalSizeDict0 += dictionary[0]->GetSize();
1158 totalSizeDict1 += dictionary[1]->GetSize();
1159 totalSizeDict2 += dictionary[2]->GetSize();
1161 Float_t nPixel = nRowMax * nColMax * nTimeMax;
1163 printf("<AliTRDdigitizer::MakeDigits> ");
1164 printf("Found %d digits in detector %d (%3.0f).\n"
1166 ,100.0 * ((Float_t) nDigits) / nPixel);
1169 if (fCompress) signals->Compress(1,0);
1179 delete signalsArray;
1184 printf("<AliTRDdigitizer::MakeDigits> ");
1185 printf("Total number of analyzed hits = %d\n",countHits);
1187 printf("<AliTRDdigitizer::MakeDigits> ");
1188 printf("Total digits data size = %d, %d, %d, %d\n",totalSizeDigits
1199 //_____________________________________________________________________________
1200 void AliTRDdigitizer::AddSDigitsManager(AliTRDdigitsManager *man)
1203 // Add a digits manager for s-digits to the input list.
1206 fSDigitsManagerList->Add(man);
1210 //_____________________________________________________________________________
1211 void AliTRDdigitizer::DeleteSDigitsManager()
1214 // Removes digits manager from the input list.
1217 fSDigitsManagerList->Delete();
1221 //_____________________________________________________________________________
1222 Bool_t AliTRDdigitizer::ConvertSDigits()
1225 // Converts s-digits to normal digits
1228 // Number of track dictionary arrays
1229 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1231 // Converts number of electrons to fC
1232 const Double_t kEl2fC = 1.602E-19 * 1.0E15;
1240 fPar = new AliTRDparameter("TRDparameter","Standard parameter");
1242 printf("<AliTRDdigitizer::ConvertSDigits> ");
1243 printf("Create the default parameter object\n");
1247 Double_t sDigitsScale = 1.0 / GetSDigitsScale();
1248 Double_t noise = fPar->GetNoise();
1249 Double_t padCoupling = fPar->GetPadCoupling();
1250 Double_t timeCoupling = fPar->GetTimeCoupling();
1251 Double_t chipGain = fPar->GetChipGain();
1252 Double_t coupling = padCoupling * timeCoupling;
1253 Double_t convert = kEl2fC * chipGain;
1254 Double_t adcInRange = fPar->GetADCinRange();
1255 Double_t adcOutRange = fPar->GetADCoutRange();
1256 Int_t adcThreshold = fPar->GetADCthreshold();
1258 AliTRDdataArrayI *digitsIn;
1259 AliTRDdataArrayI *digitsOut;
1260 AliTRDdataArrayI *dictionaryIn[kNDict];
1261 AliTRDdataArrayI *dictionaryOut[kNDict];
1263 // Loop through the detectors
1264 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1267 printf("<AliTRDdigitizer::ConvertSDigits> ");
1268 printf("Convert detector %d to digits.\n",iDet);
1271 Int_t plane = fGeo->GetPlane(iDet);
1272 Int_t sector = fGeo->GetSector(iDet);
1273 Int_t chamber = fGeo->GetChamber(iDet);
1274 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1275 Int_t nColMax = fPar->GetColMax(plane);
1276 Int_t nTimeTotal = fPar->GetTimeTotal();
1278 Double_t *inADC = new Double_t[nTimeTotal];
1279 Double_t *outADC = new Double_t[nTimeTotal];
1281 digitsIn = fSDigitsManager->GetDigits(iDet);
1283 digitsOut = fDigitsManager->GetDigits(iDet);
1284 digitsOut->Allocate(nRowMax,nColMax,nTimeTotal);
1285 for (iDict = 0; iDict < kNDict; iDict++) {
1286 dictionaryIn[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1287 dictionaryIn[iDict]->Expand();
1288 dictionaryOut[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1289 dictionaryOut[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1292 for (iRow = 0; iRow < nRowMax; iRow++ ) {
1293 for (iCol = 0; iCol < nColMax; iCol++ ) {
1295 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1296 Double_t signal = (Double_t) digitsIn->GetDataUnchecked(iRow,iCol,iTime);
1297 signal *= sDigitsScale;
1298 // Pad and time coupling
1301 signal = TMath::Max((Double_t) gRandom->Gaus(signal,noise),0.0);
1304 // Convert to ADC counts. Set the overflow-bit adcOutRange if the
1305 // signal is larger than adcInRange
1307 if (signal >= adcInRange) {
1308 adc = ((Int_t) adcOutRange);
1311 adc = ((Int_t) (signal * (adcOutRange / adcInRange)));
1314 outADC[iTime] = adc;
1317 // Apply the tail cancelation via the digital filter
1319 DeConvExp(inADC,outADC,nTimeTotal,fPar->GetTCnexp());
1322 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1323 // Store the amplitude of the digit if above threshold
1324 if (outADC[iTime] > adcThreshold) {
1325 digitsOut->SetDataUnchecked(iRow,iCol,iTime,((Int_t) outADC[iTime]));
1326 // Copy the dictionary
1327 for (iDict = 0; iDict < kNDict; iDict++) {
1328 Int_t track = dictionaryIn[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1329 dictionaryOut[iDict]->SetDataUnchecked(iRow,iCol,iTime,track);
1338 digitsIn->Compress(1,0);
1339 digitsOut->Compress(1,0);
1340 for (iDict = 0; iDict < kNDict; iDict++) {
1341 dictionaryIn[iDict]->Compress(1,0);
1342 dictionaryOut[iDict]->Compress(1,0);
1355 //_____________________________________________________________________________
1356 Bool_t AliTRDdigitizer::MergeSDigits()
1359 // Merges the input s-digits:
1360 // - The amplitude of the different inputs are summed up.
1361 // - Of the track IDs from the input dictionaries only one is
1362 // kept for each input. This works for maximal 3 different merged inputs.
1365 // Number of track dictionary arrays
1366 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1369 fPar = new AliTRDparameter("TRDparameter","Standard parameter");
1371 printf("<AliTRDdigitizer::MergeSDigits> ");
1372 printf("Create the default parameter object\n");
1379 AliTRDdataArrayI *digitsA;
1380 AliTRDdataArrayI *digitsB;
1381 AliTRDdataArrayI *dictionaryA[kNDict];
1382 AliTRDdataArrayI *dictionaryB[kNDict];
1384 // Get the first s-digits
1385 fSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->First();
1386 if (!fSDigitsManager) return kFALSE;
1388 // Loop through the other sets of s-digits
1389 AliTRDdigitsManager *mergeSDigitsManager;
1390 mergeSDigitsManager = (AliTRDdigitsManager *)
1391 fSDigitsManagerList->After(fSDigitsManager);
1394 if (mergeSDigitsManager) {
1395 printf("<AliTRDdigitizer::MergeSDigits> ");
1396 printf("Merge %d input files.\n",fSDigitsManagerList->GetSize());
1399 printf("<AliTRDdigitizer::MergeSDigits> ");
1400 printf("Only one input file.\n");
1405 while (mergeSDigitsManager) {
1409 // Loop through the detectors
1410 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1412 Int_t plane = fGeo->GetPlane(iDet);
1413 Int_t sector = fGeo->GetSector(iDet);
1414 Int_t chamber = fGeo->GetChamber(iDet);
1415 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1416 Int_t nColMax = fPar->GetColMax(plane);
1417 Int_t nTimeTotal = fPar->GetTimeTotal();
1419 // Loop through the pixels of one detector and add the signals
1420 digitsA = fSDigitsManager->GetDigits(iDet);
1421 digitsB = mergeSDigitsManager->GetDigits(iDet);
1424 for (iDict = 0; iDict < kNDict; iDict++) {
1425 dictionaryA[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1426 dictionaryB[iDict] = mergeSDigitsManager->GetDictionary(iDet,iDict);
1427 dictionaryA[iDict]->Expand();
1428 dictionaryB[iDict]->Expand();
1431 // Merge only detectors that contain a signal
1432 Bool_t doMerge = kTRUE;
1433 if (fMergeSignalOnly) {
1434 if (digitsA->GetOverThreshold(0) == 0) {
1442 printf("<AliTRDdigitizer::MergeSDigits> ");
1443 printf("Merge detector %d of input no.%d\n",iDet,iMerge+1);
1446 for (Int_t iRow = 0; iRow < nRowMax; iRow++ ) {
1447 for (Int_t iCol = 0; iCol < nColMax; iCol++ ) {
1448 for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {
1450 // Add the amplitudes of the summable digits
1451 Int_t ampA = digitsA->GetDataUnchecked(iRow,iCol,iTime);
1452 Int_t ampB = digitsB->GetDataUnchecked(iRow,iCol,iTime);
1454 digitsA->SetDataUnchecked(iRow,iCol,iTime,ampA);
1456 // Add the mask to the track id if defined.
1457 for (iDict = 0; iDict < kNDict; iDict++) {
1458 Int_t trackB = dictionaryB[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1459 if ((fMasks) && (trackB > 0)) {
1460 for (jDict = 0; jDict < kNDict; jDict++) {
1461 Int_t trackA = dictionaryA[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1463 trackA = trackB + fMasks[iMerge];
1464 dictionaryA[iDict]->SetDataUnchecked(iRow,iCol,iTime,trackA);
1477 digitsA->Compress(1,0);
1478 digitsB->Compress(1,0);
1479 for (iDict = 0; iDict < kNDict; iDict++) {
1480 dictionaryA[iDict]->Compress(1,0);
1481 dictionaryB[iDict]->Compress(1,0);
1487 // The next set of s-digits
1488 mergeSDigitsManager = (AliTRDdigitsManager *)
1489 fSDigitsManagerList->After(mergeSDigitsManager);
1497 //_____________________________________________________________________________
1498 Bool_t AliTRDdigitizer::SDigits2Digits()
1501 // Merges the input s-digits and converts them to normal digits
1504 if (!MergeSDigits()) return kFALSE;
1506 return ConvertSDigits();
1510 //_____________________________________________________________________________
1511 Bool_t AliTRDdigitizer::CheckDetector(Int_t plane, Int_t chamber, Int_t sector)
1514 // Checks whether a detector is enabled
1517 if (fSimpleSim) return kTRUE;
1519 if ((fTRD->GetSensChamber() >= 0) &&
1520 (fTRD->GetSensChamber() != chamber)) return kFALSE;
1521 if ((fTRD->GetSensPlane() >= 0) &&
1522 (fTRD->GetSensPlane() != plane)) return kFALSE;
1523 if ( fTRD->GetSensSector() >= 0) {
1524 Int_t sens1 = fTRD->GetSensSector();
1525 Int_t sens2 = sens1 + fTRD->GetSensSectorRange();
1526 sens2 -= ((Int_t) (sens2 / AliTRDgeometry::Nsect()))
1527 * AliTRDgeometry::Nsect();
1528 if (sens1 < sens2) {
1529 if ((sector < sens1) || (sector >= sens2)) return kFALSE;
1532 if ((sector < sens1) && (sector >= sens2)) return kFALSE;
1540 //_____________________________________________________________________________
1541 Bool_t AliTRDdigitizer::WriteDigits() const
1544 // Writes out the TRD-digits and the dictionaries
1547 // Store the digits and the dictionary in the tree
1548 return fDigitsManager->WriteDigits();
1552 //_____________________________________________________________________________
1553 void AliTRDdigitizer::DeConvExp(Double_t *source, Double_t *target
1554 , Int_t n, Int_t nexp)
1557 // Does the deconvolution by the digital filter.
1559 // Author: Marcus Gutfleisch, KIP Heidelberg
1560 // Optimized for: New TRF from Venelin Angelov, simulated with CADENCE
1561 // Pad-ground capacitance = 25 pF
1562 // Pad-pad cross talk capacitance = 6 pF
1563 // For 10 MHz digitization, corresponding to 20 time bins
1564 // in the drift region
1568 Double_t coefficients[2];
1570 /* initialize (coefficient = alpha, rates = lambda) */
1573 rates[0] = 0.466998;
1575 coefficients[0] = 1.0;
1578 rates[0] = 0.8988162;
1579 coefficients[0] = 0.11392069;
1580 rates[1] = 0.3745688;
1581 coefficients[1] = 0.8860793;
1583 Float_t sumc = coefficients[0]+coefficients[1];
1584 coefficients[0] /= sumc;
1585 coefficients[1] /= sumc;
1589 Double_t reminder[2];
1590 Double_t correction, result;
1592 /* attention: computation order is important */
1594 for ( k = 0; k < nexp; k++ ) reminder[k]=0.0;
1596 for ( i = 0; i < n; i++ ) {
1597 result = ( source[i] - correction ); /* no rescaling */
1600 for ( k = 0; k < nexp; k++ ) reminder[k] = rates[k]
1601 * ( reminder[k] + coefficients[k] * result);
1604 for ( k = 0; k < nexp; k++ ) correction += reminder[k];
1609 //_____________________________________________________________________________
1610 void AliTRDdigitizer::InitOutput(TFile *file, Int_t iEvent)
1613 // Initializes the output branches
1617 fDigitsManager->MakeTreeAndBranches(file,iEvent);