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.38 2002/04/30 08:30:40 cblume
19 gAlice now only read by AliRunDigitizer. Therefore it is just deleted in AliTRDmerge.C
21 Revision 1.37 2002/04/29 11:50:47 cblume
22 Change initialization of gAlice in the merging case
24 Revision 1.36 2002/04/12 12:13:23 cblume
27 Revision 1.35 2002/03/28 14:59:07 cblume
30 Revision 1.34 2002/03/25 20:00:44 cblume
31 Introduce parameter class and regions of interest for merging
33 Revision 1.33 2002/02/12 17:32:03 cblume
34 Rearrange the deleting of the list of sdigitsmanager
36 Revision 1.32 2002/02/12 16:07:21 cblume
39 Revision 1.31 2002/02/11 14:27:11 cblume
40 New pad plane design, new TRF+PRF, tail cancelation, cross talk
42 Revision 1.30 2001/11/19 08:44:08 cblume
43 Fix bugs reported by Rene
45 Revision 1.29 2001/11/14 19:44:25 hristov
46 Numeric const casted (Alpha)
48 Revision 1.28 2001/11/14 16:35:58 cblume
49 Inherits now from AliDetector
51 Revision 1.27 2001/11/14 10:50:45 cblume
52 Changes in digits IO. Add merging of summable digits
54 Revision 1.26 2001/11/06 17:19:41 cblume
55 Add detailed geometry and simple simulator
57 Revision 1.25 2001/06/27 09:54:44 cblume
58 Moved fField initialization to InitDetector()
60 Revision 1.24 2001/05/21 16:45:47 hristov
61 Last minute changes (C.Blume)
63 Revision 1.23 2001/05/07 08:04:48 cblume
64 New TRF and PRF. Speedup of the code. Digits from amplification region included
66 Revision 1.22 2001/03/30 14:40:14 cblume
67 Update of the digitization parameter
69 Revision 1.21 2001/03/13 09:30:35 cblume
70 Update of digitization. Moved digit branch definition to AliTRD
72 Revision 1.20 2001/02/25 20:19:00 hristov
73 Minor correction: loop variable declared only once for HP, Sun
75 Revision 1.19 2001/02/14 18:22:26 cblume
76 Change in the geometry of the padplane
78 Revision 1.18 2001/01/26 19:56:57 hristov
79 Major upgrade of AliRoot code
81 Revision 1.17 2000/12/08 12:53:27 cblume
82 Change in Copy() function for HP-compiler
84 Revision 1.16 2000/12/07 12:20:46 cblume
85 Go back to array compression. Use sampled PRF to speed up digitization
87 Revision 1.15 2000/11/23 14:34:08 cblume
88 Fixed bug in expansion routine of arrays (initialize buffers properly)
90 Revision 1.14 2000/11/20 08:54:44 cblume
91 Switch off compression as default
93 Revision 1.13 2000/11/10 14:57:52 cblume
94 Changes in the geometry constants for the DEC compiler
96 Revision 1.12 2000/11/01 14:53:20 cblume
97 Merge with TRD-develop
99 Revision 1.1.4.9 2000/10/26 17:00:22 cblume
100 Fixed bug in CheckDetector()
102 Revision 1.1.4.8 2000/10/23 13:41:35 cblume
103 Added protection against Log(0) in the gas gain calulation
105 Revision 1.1.4.7 2000/10/17 02:27:34 cblume
106 Get rid of global constants
108 Revision 1.1.4.6 2000/10/16 01:16:53 cblume
109 Changed timebin 0 to be the one closest to the readout
111 Revision 1.1.4.5 2000/10/15 23:34:29 cblume
112 Faster version of the digitizer
114 Revision 1.1.4.4 2000/10/06 16:49:46 cblume
117 Revision 1.1.4.3 2000/10/04 16:34:58 cblume
118 Replace include files by forward declarations
120 Revision 1.1.4.2 2000/09/22 14:41:10 cblume
121 Bug fix in PRF. Included time response. New structure
123 Revision 1.10 2000/10/05 07:27:53 cblume
124 Changes in the header-files by FCA
126 Revision 1.9 2000/10/02 21:28:19 fca
127 Removal of useless dependecies via forward declarations
129 Revision 1.8 2000/06/09 11:10:07 cblume
130 Compiler warnings and coding conventions, next round
132 Revision 1.7 2000/06/08 18:32:58 cblume
133 Make code compliant to coding conventions
135 Revision 1.6 2000/06/07 16:27:32 cblume
136 Try to remove compiler warnings on Sun and HP
138 Revision 1.5 2000/05/09 16:38:57 cblume
139 Removed PadResponse(). Merge problem
141 Revision 1.4 2000/05/08 15:53:45 cblume
142 Resolved merge conflict
144 Revision 1.3 2000/04/28 14:49:27 cblume
145 Only one declaration of iDict in MakeDigits()
147 Revision 1.1.4.1 2000/05/08 14:42:04 cblume
148 Introduced AliTRDdigitsManager
150 Revision 1.1 2000/02/28 19:00:13 cblume
155 ///////////////////////////////////////////////////////////////////////////////
157 // Creates and handles digits from TRD hits //
158 // Author: C. Blume (C.Blume@gsi.de) //
160 // The following effects are included: //
163 // - Gas gain including fluctuations //
164 // - Pad-response (simple Gaussian approximation) //
165 // - Time-response //
166 // - Electronics noise //
167 // - Electronics gain //
169 // - ADC threshold //
170 // The corresponding parameter can be adjusted via the various //
171 // Set-functions. If these parameters are not explicitly set, default //
172 // values are used (see Init-function). //
173 // As an example on how to use this class to produce digits from hits //
174 // have a look at the macro hits2digits.C //
175 // The production of summable digits is demonstrated in hits2sdigits.C //
176 // and the subsequent conversion of the s-digits into normal digits is //
177 // explained in sdigits2digits.C. //
179 ///////////////////////////////////////////////////////////////////////////////
195 #include "AliRunDigitizer.h"
198 #include "AliTRDhit.h"
199 #include "AliTRDdigitizer.h"
200 #include "AliTRDdataArrayI.h"
201 #include "AliTRDdataArrayF.h"
202 #include "AliTRDsegmentArray.h"
203 #include "AliTRDdigitsManager.h"
204 #include "AliTRDgeometry.h"
205 #include "AliTRDparameter.h"
207 ClassImp(AliTRDdigitizer)
209 //_____________________________________________________________________________
210 AliTRDdigitizer::AliTRDdigitizer()
213 // AliTRDdigitizer default constructor
218 fSDigitsManagerList = 0;
229 fMergeSignalOnly = kFALSE;
235 //_____________________________________________________________________________
236 AliTRDdigitizer::AliTRDdigitizer(const Text_t *name, const Text_t *title)
237 :AliDigitizer(name,title)
240 // AliTRDdigitizer constructor
245 fSDigitsManagerList = 0;
255 fMergeSignalOnly = kFALSE;
259 // For the summable digits
260 fSDigitsScale = 100.;
264 //_____________________________________________________________________________
265 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager
266 , const Text_t *name, const Text_t *title)
267 :AliDigitizer(manager,name,title)
270 // AliTRDdigitizer constructor
275 fSDigitsManagerList = 0;
285 fMergeSignalOnly = kFALSE;
289 // For the summable digits
290 fSDigitsScale = 100.;
294 //_____________________________________________________________________________
295 AliTRDdigitizer::AliTRDdigitizer(AliRunDigitizer *manager)
296 :AliDigitizer(manager,"AliTRDdigitizer","TRD digitizer")
299 // AliTRDdigitizer constructor
304 fSDigitsManagerList = 0;
314 fMergeSignalOnly = kFALSE;
318 // For the summable digits
319 fSDigitsScale = 100.;
323 //_____________________________________________________________________________
324 AliTRDdigitizer::AliTRDdigitizer(const AliTRDdigitizer &d)
327 // AliTRDdigitizer copy constructor
330 ((AliTRDdigitizer &) d).Copy(*this);
334 //_____________________________________________________________________________
335 AliTRDdigitizer::~AliTRDdigitizer()
338 // AliTRDdigitizer destructor
347 if (fDigitsManager) {
348 delete fDigitsManager;
352 if (fSDigitsManager) {
353 delete fSDigitsManager;
357 if (fSDigitsManagerList) {
358 delete fSDigitsManagerList;
359 fSDigitsManagerList = 0;
369 //_____________________________________________________________________________
370 AliTRDdigitizer &AliTRDdigitizer::operator=(const AliTRDdigitizer &d)
373 // Assignment operator
376 if (this != &d) ((AliTRDdigitizer &) d).Copy(*this);
381 //_____________________________________________________________________________
382 void AliTRDdigitizer::Copy(TObject &d)
388 ((AliTRDdigitizer &) d).fInputFile = 0;
389 ((AliTRDdigitizer &) d).fSDigitsManagerList = 0;
390 ((AliTRDdigitizer &) d).fSDigitsManager = 0;
391 ((AliTRDdigitizer &) d).fDigitsManager = 0;
392 ((AliTRDdigitizer &) d).fTRD = 0;
393 ((AliTRDdigitizer &) d).fGeo = 0;
394 ((AliTRDdigitizer &) d).fMasks = 0;
395 ((AliTRDdigitizer &) d).fEvent = 0;
396 ((AliTRDdigitizer &) d).fPar = 0;
397 ((AliTRDdigitizer &) d).fCompress = fCompress;
398 ((AliTRDdigitizer &) d).fDebug = fDebug ;
399 ((AliTRDdigitizer &) d).fSDigits = fSDigits;
400 ((AliTRDdigitizer &) d).fSDigitsScale = fSDigitsScale;
401 ((AliTRDdigitizer &) d).fMergeSignalOnly = fMergeSignalOnly;
402 ((AliTRDdigitizer &) d).fSimpleSim = fSimpleSim;
403 ((AliTRDdigitizer &) d).fSimpleDet = fSimpleDet;
407 //_____________________________________________________________________________
408 void AliTRDdigitizer::Exec(Option_t* option)
411 // Executes the merging
416 AliTRDdigitsManager *sdigitsManager;
418 TString optionString = option;
419 if (optionString.Contains("deb")) {
421 if (optionString.Contains("2")) {
424 printf("<AliTRDdigitizer::Exec> ");
425 printf("Called with debug option %d\n",fDebug);
428 // The AliRoot file is already connected by the manager
431 printf("<AliTRDdigitizer::Exec> ");
432 printf("AliRun object found on file.\n");
436 printf("<AliTRDdigitizer::Exec> ");
437 printf("Could not find AliRun object.\n");
441 Int_t nInput = fManager->GetNinputs();
442 fMasks = new Int_t[nInput];
443 for (iInput = 0; iInput < nInput; iInput++) {
444 fMasks[iInput] = fManager->GetMask(iInput);
450 for (iInput = 0; iInput < nInput; iInput++) {
453 printf("<AliTRDdigitizer::Exec> ");
454 printf("Add input stream %d\n",iInput);
457 // check if the input tree exists
458 if (!fManager->GetInputTreeTRDS(iInput)) {
459 printf("<AliTRDdigitizer::Exec> ");
460 printf("Input stream %d does not exist\n",iInput);
464 // Read the s-digits via digits manager
465 sdigitsManager = new AliTRDdigitsManager();
466 sdigitsManager->SetDebug(fDebug);
467 sdigitsManager->SetSDigits(kTRUE);
468 sdigitsManager->ReadDigits(fManager->GetInputTreeTRDS(iInput));
470 // Add the s-digits to the input list
471 AddSDigitsManager(sdigitsManager);
475 // Convert the s-digits to normal digits
477 printf("<AliTRDdigitizer::Exec> ");
478 printf("Do the conversion\n");
484 printf("<AliTRDdigitizer::Exec> ");
485 printf("Write the digits\n");
487 fDigitsManager->MakeBranch(fManager->GetTreeDTRD());
488 fDigitsManager->WriteDigits();
490 printf("<AliTRDdigitizer::Exec> ");
494 DeleteSDigitsManager();
498 //_____________________________________________________________________________
499 Bool_t AliTRDdigitizer::Open(const Char_t *file, Int_t nEvent)
502 // Opens a ROOT-file with TRD-hits and reads in the hit-tree
505 // Connect the AliRoot file containing Geometry, Kine, and Hits
506 fInputFile = (TFile *) gROOT->GetListOfFiles()->FindObject(file);
509 printf("<AliTRDdigitizer::Open> ");
510 printf("Open the AliROOT-file %s.\n",file);
512 fInputFile = new TFile(file,"UPDATE");
516 printf("<AliTRDdigitizer::Open> ");
517 printf("%s is already open.\n",file);
521 gAlice = (AliRun *) fInputFile->Get("gAlice");
524 printf("<AliTRDdigitizer::Open> ");
525 printf("AliRun object found on file.\n");
529 printf("<AliTRDdigitizer::Open> ");
530 printf("Could not find AliRun object.\n");
536 // Import the Trees for the event nEvent in the file
537 Int_t nparticles = gAlice->GetEvent(fEvent);
538 if (nparticles <= 0) {
539 printf("<AliTRDdigitizer::Open> ");
540 printf("No entries in the trees for event %d.\n",fEvent);
544 if (InitDetector()) {
553 //_____________________________________________________________________________
554 Bool_t AliTRDdigitizer::InitDetector()
557 // Sets the pointer to the TRD detector and the geometry
560 // Get the pointer to the detector class and check for version 1
561 fTRD = (AliTRD *) gAlice->GetDetector("TRD");
563 printf("<AliTRDdigitizer::InitDetector> ");
564 printf("No TRD module found\n");
567 if (fTRD->IsVersion() != 1) {
568 printf("<AliTRDdigitizer::InitDetector> ");
569 printf("TRD must be version 1 (slow simulator).\n");
574 fGeo = fTRD->GetGeometry();
576 printf("<AliTRDdigitizer::InitDetector> ");
577 printf("Geometry version %d\n",fGeo->IsVersion());
580 // Create a digits manager
581 fDigitsManager = new AliTRDdigitsManager();
582 fDigitsManager->SetSDigits(fSDigits);
583 fDigitsManager->CreateArrays();
584 fDigitsManager->SetEvent(fEvent);
585 fDigitsManager->SetDebug(fDebug);
587 // The list for the input s-digits manager to be merged
588 fSDigitsManagerList = new TList();
594 //_____________________________________________________________________________
595 Bool_t AliTRDdigitizer::MakeBranch(const Char_t *file) const
598 // Create the branches for the digits array
601 return fDigitsManager->MakeBranch(file);
605 //_____________________________________________________________________________
606 Bool_t AliTRDdigitizer::MakeDigits()
612 ///////////////////////////////////////////////////////////////
614 ///////////////////////////////////////////////////////////////
616 // Converts number of electrons to fC
617 const Double_t kEl2fC = 1.602E-19 * 1.0E15;
619 ///////////////////////////////////////////////////////////////
621 // Number of pads included in the pad response
622 const Int_t kNpad = 3;
624 // Number of track dictionary arrays
625 const Int_t kNDict = AliTRDdigitsManager::kNDict;
627 // Half the width of the amplification region
628 const Float_t kAmWidth = AliTRDgeometry::AmThick() / 2.;
630 Int_t iRow, iCol, iTime, iPad;
634 Int_t totalSizeDigits = 0;
635 Int_t totalSizeDict0 = 0;
636 Int_t totalSizeDict1 = 0;
637 Int_t totalSizeDict2 = 0;
639 Int_t timeTRDbeg = 0;
640 Int_t timeTRDend = 1;
645 Float_t padSignal[kNpad];
646 Float_t signalOld[kNpad];
648 AliTRDdataArrayF *signals = 0;
649 AliTRDdataArrayI *digits = 0;
650 AliTRDdataArrayI *dictionary[kNDict];
652 // Create a default parameter class if none is defined
654 fPar = new AliTRDparameter("TRDparameter","Standard TRD parameter");
656 printf("<AliTRDdigitizer::MakeDigits> ");
657 printf("Create the default parameter object\n");
661 // Create a container for the amplitudes
662 AliTRDsegmentArray *signalsArray
663 = new AliTRDsegmentArray("AliTRDdataArrayF"
664 ,AliTRDgeometry::Ndet());
667 timeTRDbeg = ((Int_t) (-fPar->GetTRFlo() / fPar->GetTimeBinSize())) - 1;
668 timeTRDend = ((Int_t) ( fPar->GetTRFhi() / fPar->GetTimeBinSize())) - 1;
670 printf("<AliTRDdigitizer::MakeDigits> ");
671 printf("Sample the TRF between -%d and %d\n",timeTRDbeg,timeTRDend);
675 Float_t elAttachProp = fPar->GetElAttachProp() / 100.;
678 printf("<AliTRDdigitizer::MakeDigits> ");
679 printf("No geometry defined\n");
684 printf("<AliTRDdigitizer::MakeDigits> ");
685 printf("Start creating digits.\n");
688 // Get the pointer to the hit tree
689 TTree *hitTree = gAlice->TreeH();
691 // Get the number of entries in the hit tree
692 // (Number of primary particles creating a hit somewhere)
695 nTrack = (Int_t) hitTree->GetEntries();
697 printf("<AliTRDdigitizer::MakeDigits> ");
698 printf("Found %d primary particles\n",nTrack);
702 Int_t detectorOld = -1;
705 // Loop through all entries in the tree
706 for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
710 nBytes += hitTree->GetEvent(iTrack);
713 // Loop through the TRD hits
715 AliTRDhit *hit = (AliTRDhit *) fTRD->FirstHit(-1);
724 Float_t q = hit->GetCharge();
725 Int_t track = hit->Track();
726 Int_t detector = hit->GetDetector();
727 Int_t plane = fGeo->GetPlane(detector);
728 Int_t sector = fGeo->GetSector(detector);
729 Int_t chamber = fGeo->GetChamber(detector);
730 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
731 Int_t nColMax = fPar->GetColMax(plane);
732 Int_t nTimeMax = fPar->GetTimeMax();
733 Int_t nTimeBefore = fPar->GetTimeBefore();
734 Int_t nTimeAfter = fPar->GetTimeAfter();
735 Int_t nTimeTotal = fPar->GetTimeTotal();
736 Float_t row0 = fPar->GetRow0(plane,chamber,sector);
737 Float_t col0 = fPar->GetCol0(plane);
738 Float_t time0 = fPar->GetTime0(plane);
739 Float_t rowPadSize = fPar->GetRowPadSize(plane,chamber,sector);
740 Float_t colPadSize = fPar->GetColPadSize(plane);
741 Float_t timeBinSize = fPar->GetTimeBinSize();
742 Float_t divideRow = 1.0 / rowPadSize;
743 Float_t divideCol = 1.0 / colPadSize;
744 Float_t divideTime = 1.0 / timeBinSize;
747 printf("Analyze hit no. %d ",iHit);
748 printf("-----------------------------------------------------------\n");
750 printf("plane = %d, sector = %d, chamber = %d\n"
751 ,plane,sector,chamber);
752 printf("nRowMax = %d, nColMax = %d, nTimeMax = %d\n"
753 ,nRowMax,nColMax,nTimeMax);
754 printf("nTimeBefore = %d, nTimeAfter = %d, nTimeTotal = %d\n"
755 ,nTimeBefore,nTimeAfter,nTimeTotal);
756 printf("row0 = %f, col0 = %f, time0 = %f\n"
758 printf("rowPadSize = %f, colPadSize = %f, timeBinSize = %f\n"
759 ,rowPadSize,colPadSize,timeBinSize);
762 // Don't analyze test hits and switched off detectors
763 if ((CheckDetector(plane,chamber,sector)) &&
764 (((Int_t) q) != 0)) {
766 if (detector != detectorOld) {
769 printf("<AliTRDdigitizer::MakeDigits> ");
770 printf("Get new container. New det = %d, Old det = %d\n"
771 ,detector,detectorOld);
773 // Compress the old one if enabled
774 if ((fCompress) && (detectorOld > -1)) {
776 printf("<AliTRDdigitizer::MakeDigits> ");
777 printf("Compress the old container ...");
779 signals->Compress(1,0);
780 for (iDict = 0; iDict < kNDict; iDict++) {
781 dictionary[iDict]->Compress(1,0);
783 if (fDebug > 1) printf("done\n");
785 // Get the new container
786 signals = (AliTRDdataArrayF *) signalsArray->At(detector);
787 if (signals->GetNtime() == 0) {
788 // Allocate a new one if not yet existing
790 printf("<AliTRDdigitizer::MakeDigits> ");
791 printf("Allocate a new container ... ");
793 signals->Allocate(nRowMax,nColMax,nTimeTotal);
795 else if (fSimpleSim) {
796 // Clear an old one for the simple simulation
798 printf("<AliTRDdigitizer::MakeDigits> ");
799 printf("Clear a old container ... ");
804 // Expand an existing one
807 printf("<AliTRDdigitizer::MakeDigits> ");
808 printf("Expand an existing container ... ");
813 // The same for the dictionary
815 for (iDict = 0; iDict < kNDict; iDict++) {
816 dictionary[iDict] = fDigitsManager->GetDictionary(detector,iDict);
817 if (dictionary[iDict]->GetNtime() == 0) {
818 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
821 if (fCompress) dictionary[iDict]->Expand();
825 if (fDebug > 1) printf("done\n");
826 detectorOld = detector;
829 // Rotate the sectors on top of each other
836 fGeo->Rotate(detector,pos,rot);
839 // The driftlength. It is negative if the hit is in the
840 // amplification region.
841 Float_t driftlength = time0 - rot[0];
843 // Take also the drift in the amplification region into account
844 // The drift length is at the moment still the same, regardless of
845 // the position relativ to the wire. This non-isochronity needs still
846 // to be implemented.
847 Float_t driftlengthL = TMath::Abs(driftlength + kAmWidth);
848 if (fPar->ExBOn()) driftlengthL /= TMath::Sqrt(fPar->GetLorentzFactor());
850 // Loop over all electrons of this hit
851 // TR photons produce hits with negative charge
852 Int_t nEl = ((Int_t) TMath::Abs(q));
853 for (Int_t iEl = 0; iEl < nEl; iEl++) {
859 // Electron attachment
860 if (fPar->ElAttachOn()) {
861 if (gRandom->Rndm() < (driftlengthL * elAttachProp))
865 // Apply the diffusion smearing
866 if (fPar->DiffusionOn()) {
867 if (!(fPar->Diffusion(driftlengthL,xyz))) continue;
870 // Apply E x B effects (depends on drift direction)
872 if (!(fPar->ExB(driftlength+kAmWidth,xyz))) continue;
875 // The electron position after diffusion and ExB in pad coordinates
876 // The pad row (z-direction)
877 Float_t rowDist = xyz[2] - row0;
878 Int_t rowE = ((Int_t) (rowDist * divideRow));
879 if ((rowE < 0) || (rowE >= nRowMax)) continue;
880 Float_t rowOffset = ((((Float_t) rowE) + 0.5) * rowPadSize) - rowDist;
882 // The pad column (rphi-direction)
883 Float_t col0tilt = fPar->Col0Tilted(col0,rowOffset,plane);
884 Float_t colDist = xyz[1] - col0tilt;
885 Int_t colE = ((Int_t) (colDist * divideCol));
886 if ((colE < 0) || (colE >= nColMax)) continue;
887 Float_t colOffset = ((((Float_t) colE) + 0.5) * colPadSize) - colDist;
889 // The time bin (negative for hits in the amplification region)
890 // In the amplification region the electrons drift from both sides
891 // to the middle (anode wire plane)
892 Float_t timeDist = time0 - xyz[0];
893 Float_t timeOffset = 0;
897 timeE = ((Int_t) (timeDist * divideTime));
898 // The distance of the position to the middle of the timebin
899 timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) - timeDist;
902 // Difference between half of the amplification gap width and
903 // the distance to the anode wire
904 Float_t anodeDist = kAmWidth - TMath::Abs(timeDist + kAmWidth);
906 timeE = -1 * (((Int_t ) (anodeDist * divideTime)) + 1);
907 // The distance of the position to the middle of the timebin
908 timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) + anodeDist;
911 // Apply the gas gain including fluctuations
912 Float_t ggRndm = 0.0;
914 ggRndm = gRandom->Rndm();
915 } while (ggRndm <= 0);
916 Int_t signal = (Int_t) (-fPar->GetGasGain() * TMath::Log(ggRndm));
918 // Apply the pad response
920 // The distance of the electron to the center of the pad
921 // in units of pad width
922 Float_t dist = - colOffset * divideCol;
923 if (!(fPar->PadResponse(signal,dist,plane,padSignal))) continue;
927 padSignal[1] = signal;
931 // Sample the time response inside the drift region
932 // + additional time bins before and after.
933 // The sampling is done always in the middle of the time bin
934 for (Int_t iTimeBin = TMath::Max(timeE-timeTRDbeg, -nTimeBefore)
935 ;iTimeBin < TMath::Min(timeE+timeTRDend,nTimeMax+nTimeAfter )
938 // Apply the time response
939 Float_t timeResponse = 1.0;
940 Float_t crossTalk = 0.0;
941 Float_t time = (iTimeBin - timeE) * timeBinSize + timeOffset;
943 timeResponse = fPar->TimeResponse(time);
946 crossTalk = fPar->CrossTalk(time);
953 for (iPad = 0; iPad < kNpad; iPad++) {
955 Int_t colPos = colE + iPad - 1;
956 if (colPos < 0) continue;
957 if (colPos >= nColMax) break;
960 // Note: The time bin number is shifted by nTimeBefore to avoid negative
961 // time bins. This has to be subtracted later.
962 Int_t iCurrentTimeBin = iTimeBin + nTimeBefore;
963 signalOld[iPad] = signals->GetDataUnchecked(rowE,colPos,iCurrentTimeBin);
964 if( colPos != colE ) {
965 signalOld[iPad] += padSignal[iPad] * (timeResponse + crossTalk);
968 signalOld[iPad] += padSignal[iPad] * timeResponse;
970 signals->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,signalOld[iPad]);
972 // Store the track index in the dictionary
973 // Note: We store index+1 in order to allow the array to be compressed
974 if ((signalOld[iPad] > 0) && (!fSimpleSim)) {
975 for (iDict = 0; iDict < kNDict; iDict++) {
976 Int_t oldTrack = dictionary[iDict]->GetDataUnchecked(rowE
979 if (oldTrack == track+1) break;
981 dictionary[iDict]->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,track+1);
991 } // Loop: electrons of a single hit
993 } // If: detector and test hit
995 hit = (AliTRDhit *) fTRD->NextHit();
997 } // Loop: hits of one primary track
999 } // Loop: primary tracks
1002 printf("<AliTRDdigitizer::MakeDigits> ");
1003 printf("Finished analyzing %d hits\n",countHits);
1006 // The coupling factor
1007 Float_t coupling = fPar->GetPadCoupling()
1008 * fPar->GetTimeCoupling();
1010 // The conversion factor
1011 Float_t convert = kEl2fC
1012 * fPar->GetChipGain();
1014 // Loop through all chambers to finalize the digits
1016 Int_t iDetEnd = AliTRDgeometry::Ndet();
1018 iDetBeg = fSimpleDet;
1019 iDetEnd = iDetBeg + 1;
1021 for (Int_t iDet = iDetBeg; iDet < iDetEnd; iDet++) {
1023 Int_t plane = fGeo->GetPlane(iDet);
1024 Int_t sector = fGeo->GetSector(iDet);
1025 Int_t chamber = fGeo->GetChamber(iDet);
1026 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1027 Int_t nColMax = fPar->GetColMax(plane);
1028 Int_t nTimeMax = fPar->GetTimeMax();
1029 Int_t nTimeTotal = fPar->GetTimeTotal();
1031 Double_t *inADC = new Double_t[nTimeTotal];
1032 Double_t *outADC = new Double_t[nTimeTotal];
1035 printf("<AliTRDdigitizer::MakeDigits> ");
1036 printf("Digitization for chamber %d\n",iDet);
1039 // Add a container for the digits of this detector
1040 digits = fDigitsManager->GetDigits(iDet);
1041 // Allocate memory space for the digits buffer
1042 if (digits->GetNtime() == 0) {
1043 digits->Allocate(nRowMax,nColMax,nTimeTotal);
1045 else if (fSimpleSim) {
1049 // Get the signal container
1050 signals = (AliTRDdataArrayF *) signalsArray->At(iDet);
1051 if (signals->GetNtime() == 0) {
1052 // Create missing containers
1053 signals->Allocate(nRowMax,nColMax,nTimeTotal);
1056 // Expand the container if neccessary
1057 if (fCompress) signals->Expand();
1059 // Create the missing dictionary containers
1061 for (iDict = 0; iDict < kNDict; iDict++) {
1062 dictionary[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1063 if (dictionary[iDict]->GetNtime() == 0) {
1064 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1071 // Don't create noise in detectors that are switched off
1072 if (CheckDetector(plane,chamber,sector)) {
1074 // Create the digits for this chamber
1075 for (iRow = 0; iRow < nRowMax; iRow++ ) {
1076 for (iCol = 0; iCol < nColMax; iCol++ ) {
1078 // Create summable digits
1081 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1082 Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1083 signalAmp *= fSDigitsScale;
1084 signalAmp = TMath::Min(signalAmp,(Float_t) 1.0e9);
1085 Int_t adc = (Int_t) signalAmp;
1086 if (adc > 0) nDigits++;
1087 digits->SetDataUnchecked(iRow,iCol,iTime,adc);
1091 // Create normal digits
1094 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1095 Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1096 // Pad and time coupling
1097 signalAmp *= coupling;
1099 signalAmp = TMath::Max((Double_t) gRandom->Gaus(signalAmp,fPar->GetNoise()),0.0);
1101 signalAmp *= convert;
1102 // Convert to ADC counts. Set the overflow-bit fADCoutRange if the
1103 // signal is larger than fADCinRange
1105 if (signalAmp >= fPar->GetADCinRange()) {
1106 adc = ((Int_t) fPar->GetADCoutRange());
1109 adc = ((Int_t) (signalAmp * (fPar->GetADCoutRange()
1110 / fPar->GetADCinRange())));
1113 outADC[iTime] = adc;
1116 // Apply the tail cancelation via the digital filter
1118 DeConvExp(inADC,outADC,nTimeTotal,fPar->GetTCnexp());
1121 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1122 // Store the amplitude of the digit if above threshold
1123 if (outADC[iTime] > fPar->GetADCthreshold()) {
1125 printf(" iRow = %d, iCol = %d, iTime = %d, adc = %f\n"
1126 ,iRow,iCol,iTime,outADC[iTime]);
1129 digits->SetDataUnchecked(iRow,iCol,iTime,outADC[iTime]);
1140 // Compress the arrays
1142 digits->Compress(1,0);
1143 for (iDict = 0; iDict < kNDict; iDict++) {
1144 dictionary[iDict]->Compress(1,0);
1147 totalSizeDigits += digits->GetSize();
1148 totalSizeDict0 += dictionary[0]->GetSize();
1149 totalSizeDict1 += dictionary[1]->GetSize();
1150 totalSizeDict2 += dictionary[2]->GetSize();
1152 Float_t nPixel = nRowMax * nColMax * nTimeMax;
1154 printf("<AliTRDdigitizer::MakeDigits> ");
1155 printf("Found %d digits in detector %d (%3.0f).\n"
1157 ,100.0 * ((Float_t) nDigits) / nPixel);
1160 if (fCompress) signals->Compress(1,0);
1170 delete signalsArray;
1175 printf("<AliTRDdigitizer::MakeDigits> ");
1176 printf("Total number of analyzed hits = %d\n",countHits);
1178 printf("<AliTRDdigitizer::MakeDigits> ");
1179 printf("Total digits data size = %d, %d, %d, %d\n",totalSizeDigits
1190 //_____________________________________________________________________________
1191 void AliTRDdigitizer::AddSDigitsManager(AliTRDdigitsManager *man)
1194 // Add a digits manager for s-digits to the input list.
1197 fSDigitsManagerList->Add(man);
1201 //_____________________________________________________________________________
1202 void AliTRDdigitizer::DeleteSDigitsManager()
1205 // Removes digits manager from the input list.
1208 fSDigitsManagerList->Delete();
1212 //_____________________________________________________________________________
1213 Bool_t AliTRDdigitizer::ConvertSDigits()
1216 // Converts s-digits to normal digits
1219 // Number of track dictionary arrays
1220 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1222 // Converts number of electrons to fC
1223 const Double_t kEl2fC = 1.602E-19 * 1.0E15;
1231 fPar = new AliTRDparameter("TRDparameter","Standard parameter");
1233 printf("<AliTRDdigitizer::ConvertSDigits> ");
1234 printf("Create the default parameter object\n");
1238 Double_t sDigitsScale = 1.0 / GetSDigitsScale();
1239 Double_t noise = fPar->GetNoise();
1240 Double_t padCoupling = fPar->GetPadCoupling();
1241 Double_t timeCoupling = fPar->GetTimeCoupling();
1242 Double_t chipGain = fPar->GetChipGain();
1243 Double_t coupling = padCoupling * timeCoupling;
1244 Double_t convert = kEl2fC * chipGain;
1245 Double_t adcInRange = fPar->GetADCinRange();
1246 Double_t adcOutRange = fPar->GetADCoutRange();
1247 Int_t adcThreshold = fPar->GetADCthreshold();
1249 AliTRDdataArrayI *digitsIn;
1250 AliTRDdataArrayI *digitsOut;
1251 AliTRDdataArrayI *dictionaryIn[kNDict];
1252 AliTRDdataArrayI *dictionaryOut[kNDict];
1254 // Loop through the detectors
1255 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1258 printf("<AliTRDdigitizer::ConvertSDigits> ");
1259 printf("Convert detector %d to digits.\n",iDet);
1262 Int_t plane = fGeo->GetPlane(iDet);
1263 Int_t sector = fGeo->GetSector(iDet);
1264 Int_t chamber = fGeo->GetChamber(iDet);
1265 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1266 Int_t nColMax = fPar->GetColMax(plane);
1267 Int_t nTimeTotal = fPar->GetTimeTotal();
1269 Double_t *inADC = new Double_t[nTimeTotal];
1270 Double_t *outADC = new Double_t[nTimeTotal];
1272 digitsIn = fSDigitsManager->GetDigits(iDet);
1274 digitsOut = fDigitsManager->GetDigits(iDet);
1275 digitsOut->Allocate(nRowMax,nColMax,nTimeTotal);
1276 for (iDict = 0; iDict < kNDict; iDict++) {
1277 dictionaryIn[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1278 dictionaryIn[iDict]->Expand();
1279 dictionaryOut[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1280 dictionaryOut[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1283 for (iRow = 0; iRow < nRowMax; iRow++ ) {
1284 for (iCol = 0; iCol < nColMax; iCol++ ) {
1286 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1287 Double_t signal = (Double_t) digitsIn->GetDataUnchecked(iRow,iCol,iTime);
1288 signal *= sDigitsScale;
1289 // Pad and time coupling
1292 signal = TMath::Max((Double_t) gRandom->Gaus(signal,noise),0.0);
1295 // Convert to ADC counts. Set the overflow-bit adcOutRange if the
1296 // signal is larger than adcInRange
1298 if (signal >= adcInRange) {
1299 adc = ((Int_t) adcOutRange);
1302 adc = ((Int_t) (signal * (adcOutRange / adcInRange)));
1305 outADC[iTime] = adc;
1308 // Apply the tail cancelation via the digital filter
1310 DeConvExp(inADC,outADC,nTimeTotal,fPar->GetTCnexp());
1313 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1314 // Store the amplitude of the digit if above threshold
1315 if (outADC[iTime] > adcThreshold) {
1316 digitsOut->SetDataUnchecked(iRow,iCol,iTime,outADC[iTime]);
1317 // Copy the dictionary
1318 for (iDict = 0; iDict < kNDict; iDict++) {
1319 Int_t track = dictionaryIn[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1320 dictionaryOut[iDict]->SetDataUnchecked(iRow,iCol,iTime,track);
1329 digitsIn->Compress(1,0);
1330 digitsOut->Compress(1,0);
1331 for (iDict = 0; iDict < kNDict; iDict++) {
1332 dictionaryIn[iDict]->Compress(1,0);
1333 dictionaryOut[iDict]->Compress(1,0);
1346 //_____________________________________________________________________________
1347 Bool_t AliTRDdigitizer::MergeSDigits()
1350 // Merges the input s-digits:
1351 // - The amplitude of the different inputs are summed up.
1352 // - Of the track IDs from the input dictionaries only one is
1353 // kept for each input. This works for maximal 3 different merged inputs.
1356 // Number of track dictionary arrays
1357 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1360 fPar = new AliTRDparameter("TRDparameter","Standard parameter");
1362 printf("<AliTRDdigitizer::MergeSDigits> ");
1363 printf("Create the default parameter object\n");
1370 AliTRDdataArrayI *digitsA;
1371 AliTRDdataArrayI *digitsB;
1372 AliTRDdataArrayI *dictionaryA[kNDict];
1373 AliTRDdataArrayI *dictionaryB[kNDict];
1375 // Get the first s-digits
1376 fSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->First();
1377 if (!fSDigitsManager) return kFALSE;
1379 // Loop through the other sets of s-digits
1380 AliTRDdigitsManager *mergeSDigitsManager;
1381 mergeSDigitsManager = (AliTRDdigitsManager *)
1382 fSDigitsManagerList->After(fSDigitsManager);
1385 if (mergeSDigitsManager) {
1386 printf("<AliTRDdigitizer::MergeSDigits> ");
1387 printf("Merge %d input files.\n",fSDigitsManagerList->GetSize());
1390 printf("<AliTRDdigitizer::MergeSDigits> ");
1391 printf("Only one input file.\n");
1396 while (mergeSDigitsManager) {
1400 // Loop through the detectors
1401 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1403 Int_t plane = fGeo->GetPlane(iDet);
1404 Int_t sector = fGeo->GetSector(iDet);
1405 Int_t chamber = fGeo->GetChamber(iDet);
1406 Int_t nRowMax = fPar->GetRowMax(plane,chamber,sector);
1407 Int_t nColMax = fPar->GetColMax(plane);
1408 Int_t nTimeTotal = fPar->GetTimeTotal();
1410 // Loop through the pixels of one detector and add the signals
1411 digitsA = fSDigitsManager->GetDigits(iDet);
1412 digitsB = mergeSDigitsManager->GetDigits(iDet);
1415 for (iDict = 0; iDict < kNDict; iDict++) {
1416 dictionaryA[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1417 dictionaryB[iDict] = mergeSDigitsManager->GetDictionary(iDet,iDict);
1418 dictionaryA[iDict]->Expand();
1419 dictionaryB[iDict]->Expand();
1422 // Merge only detectors that contain a signal
1423 Bool_t doMerge = kTRUE;
1424 if (fMergeSignalOnly) {
1425 if (digitsA->GetOverThreshold(0) == 0) {
1433 printf("<AliTRDdigitizer::MergeSDigits> ");
1434 printf("Merge detector %d of input no.%d\n",iDet,iMerge+1);
1437 for (Int_t iRow = 0; iRow < nRowMax; iRow++ ) {
1438 for (Int_t iCol = 0; iCol < nColMax; iCol++ ) {
1439 for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {
1441 // Add the amplitudes of the summable digits
1442 Int_t ampA = digitsA->GetDataUnchecked(iRow,iCol,iTime);
1443 Int_t ampB = digitsB->GetDataUnchecked(iRow,iCol,iTime);
1445 digitsA->SetDataUnchecked(iRow,iCol,iTime,ampA);
1447 // Add the mask to the track id if defined.
1448 for (iDict = 0; iDict < kNDict; iDict++) {
1449 Int_t trackB = dictionaryB[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1450 if ((fMasks) && (trackB > 0)) {
1451 for (jDict = 0; jDict < kNDict; jDict++) {
1452 Int_t trackA = dictionaryA[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1454 trackA = trackB + fMasks[iMerge];
1455 dictionaryA[iDict]->SetDataUnchecked(iRow,iCol,iTime,trackA);
1468 digitsA->Compress(1,0);
1469 digitsB->Compress(1,0);
1470 for (iDict = 0; iDict < kNDict; iDict++) {
1471 dictionaryA[iDict]->Compress(1,0);
1472 dictionaryB[iDict]->Compress(1,0);
1478 // The next set of s-digits
1479 mergeSDigitsManager = (AliTRDdigitsManager *)
1480 fSDigitsManagerList->After(mergeSDigitsManager);
1488 //_____________________________________________________________________________
1489 Bool_t AliTRDdigitizer::SDigits2Digits()
1492 // Merges the input s-digits and converts them to normal digits
1495 if (!MergeSDigits()) return kFALSE;
1497 return ConvertSDigits();
1501 //_____________________________________________________________________________
1502 Bool_t AliTRDdigitizer::CheckDetector(Int_t plane, Int_t chamber, Int_t sector)
1505 // Checks whether a detector is enabled
1508 if (fSimpleSim) return kTRUE;
1510 if ((fTRD->GetSensChamber() >= 0) &&
1511 (fTRD->GetSensChamber() != chamber)) return kFALSE;
1512 if ((fTRD->GetSensPlane() >= 0) &&
1513 (fTRD->GetSensPlane() != plane)) return kFALSE;
1514 if ( fTRD->GetSensSector() >= 0) {
1515 Int_t sens1 = fTRD->GetSensSector();
1516 Int_t sens2 = sens1 + fTRD->GetSensSectorRange();
1517 sens2 -= ((Int_t) (sens2 / AliTRDgeometry::Nsect()))
1518 * AliTRDgeometry::Nsect();
1519 if (sens1 < sens2) {
1520 if ((sector < sens1) || (sector >= sens2)) return kFALSE;
1523 if ((sector < sens1) && (sector >= sens2)) return kFALSE;
1531 //_____________________________________________________________________________
1532 Bool_t AliTRDdigitizer::WriteDigits() const
1535 // Writes out the TRD-digits and the dictionaries
1538 // Store the digits and the dictionary in the tree
1539 return fDigitsManager->WriteDigits();
1543 //_____________________________________________________________________________
1544 void AliTRDdigitizer::DeConvExp(Double_t *source, Double_t *target
1545 , Int_t n, Int_t nexp)
1548 // Does the deconvolution by the digital filter.
1550 // Author: Marcus Gutfleisch, KIP Heidelberg
1551 // Optimized for: New TRF from Venelin Angelov, simulated with CADENCE
1552 // Pad-ground capacitance = 25 pF
1553 // Pad-pad cross talk capacitance = 6 pF
1554 // For 10 MHz digitization, corresponding to 20 time bins
1555 // in the drift region
1559 Double_t coefficients[2];
1561 /* initialize (coefficient = alpha, rates = lambda) */
1564 rates[0] = 0.466998;
1566 coefficients[0] = 1.0;
1569 rates[0] = 0.8988162;
1570 coefficients[0] = 0.11392069;
1571 rates[1] = 0.3745688;
1572 coefficients[1] = 0.8860793;
1574 Float_t sumc = coefficients[0]+coefficients[1];
1575 coefficients[0] /= sumc;
1576 coefficients[1] /= sumc;
1580 Double_t reminder[2];
1581 Double_t correction, result;
1583 /* attention: computation order is important */
1585 for ( k = 0; k < nexp; k++ ) reminder[k]=0.0;
1587 for ( i = 0; i < n; i++ ) {
1588 result = ( source[i] - correction ); /* no rescaling */
1591 for ( k = 0; k < nexp; k++ ) reminder[k] = rates[k]
1592 * ( reminder[k] + coefficients[k] * result);
1595 for ( k = 0; k < nexp; k++ ) correction += reminder[k];
1600 //_____________________________________________________________________________
1601 void AliTRDdigitizer::InitOutput(TFile *file, Int_t iEvent)
1604 // Initializes the output branches
1608 fDigitsManager->MakeTreeAndBranches(file,iEvent);