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.26 2001/11/06 17:19:41 cblume
19 Add detailed geometry and simple simulator
21 Revision 1.25 2001/06/27 09:54:44 cblume
22 Moved fField initialization to InitDetector()
24 Revision 1.24 2001/05/21 16:45:47 hristov
25 Last minute changes (C.Blume)
27 Revision 1.23 2001/05/07 08:04:48 cblume
28 New TRF and PRF. Speedup of the code. Digits from amplification region included
30 Revision 1.22 2001/03/30 14:40:14 cblume
31 Update of the digitization parameter
33 Revision 1.21 2001/03/13 09:30:35 cblume
34 Update of digitization. Moved digit branch definition to AliTRD
36 Revision 1.20 2001/02/25 20:19:00 hristov
37 Minor correction: loop variable declared only once for HP, Sun
39 Revision 1.19 2001/02/14 18:22:26 cblume
40 Change in the geometry of the padplane
42 Revision 1.18 2001/01/26 19:56:57 hristov
43 Major upgrade of AliRoot code
45 Revision 1.17 2000/12/08 12:53:27 cblume
46 Change in Copy() function for HP-compiler
48 Revision 1.16 2000/12/07 12:20:46 cblume
49 Go back to array compression. Use sampled PRF to speed up digitization
51 Revision 1.15 2000/11/23 14:34:08 cblume
52 Fixed bug in expansion routine of arrays (initialize buffers properly)
54 Revision 1.14 2000/11/20 08:54:44 cblume
55 Switch off compression as default
57 Revision 1.13 2000/11/10 14:57:52 cblume
58 Changes in the geometry constants for the DEC compiler
60 Revision 1.12 2000/11/01 14:53:20 cblume
61 Merge with TRD-develop
63 Revision 1.1.4.9 2000/10/26 17:00:22 cblume
64 Fixed bug in CheckDetector()
66 Revision 1.1.4.8 2000/10/23 13:41:35 cblume
67 Added protection against Log(0) in the gas gain calulation
69 Revision 1.1.4.7 2000/10/17 02:27:34 cblume
70 Get rid of global constants
72 Revision 1.1.4.6 2000/10/16 01:16:53 cblume
73 Changed timebin 0 to be the one closest to the readout
75 Revision 1.1.4.5 2000/10/15 23:34:29 cblume
76 Faster version of the digitizer
78 Revision 1.1.4.4 2000/10/06 16:49:46 cblume
81 Revision 1.1.4.3 2000/10/04 16:34:58 cblume
82 Replace include files by forward declarations
84 Revision 1.1.4.2 2000/09/22 14:41:10 cblume
85 Bug fix in PRF. Included time response. New structure
87 Revision 1.10 2000/10/05 07:27:53 cblume
88 Changes in the header-files by FCA
90 Revision 1.9 2000/10/02 21:28:19 fca
91 Removal of useless dependecies via forward declarations
93 Revision 1.8 2000/06/09 11:10:07 cblume
94 Compiler warnings and coding conventions, next round
96 Revision 1.7 2000/06/08 18:32:58 cblume
97 Make code compliant to coding conventions
99 Revision 1.6 2000/06/07 16:27:32 cblume
100 Try to remove compiler warnings on Sun and HP
102 Revision 1.5 2000/05/09 16:38:57 cblume
103 Removed PadResponse(). Merge problem
105 Revision 1.4 2000/05/08 15:53:45 cblume
106 Resolved merge conflict
108 Revision 1.3 2000/04/28 14:49:27 cblume
109 Only one declaration of iDict in MakeDigits()
111 Revision 1.1.4.1 2000/05/08 14:42:04 cblume
112 Introduced AliTRDdigitsManager
114 Revision 1.1 2000/02/28 19:00:13 cblume
119 ///////////////////////////////////////////////////////////////////////////////
121 // Creates and handles digits from TRD hits //
122 // Author: C. Blume (C.Blume@gsi.de) //
124 // The following effects are included: //
127 // - Gas gain including fluctuations //
128 // - Pad-response (simple Gaussian approximation) //
129 // - Time-response //
130 // - Electronics noise //
131 // - Electronics gain //
133 // - ADC threshold //
134 // The corresponding parameter can be adjusted via the various //
135 // Set-functions. If these parameters are not explicitly set, default //
136 // values are used (see Init-function). //
137 // As an example on how to use this class to produce digits from hits //
138 // have a look at the macro hits2digits.C //
139 // The production of summable digits is demonstrated in hits2sdigits.C //
140 // and the subsequent conversion of the s-digits into normal digits is //
141 // explained in sdigits2digits.C. //
143 ///////////////////////////////////////////////////////////////////////////////
160 #include "AliTRDhit.h"
161 #include "AliTRDdigitizer.h"
162 #include "AliTRDdataArrayI.h"
163 #include "AliTRDdataArrayF.h"
164 #include "AliTRDsegmentArray.h"
165 #include "AliTRDdigitsManager.h"
166 #include "AliTRDgeometry.h"
168 ClassImp(AliTRDdigitizer)
170 //_____________________________________________________________________________
171 AliTRDdigitizer::AliTRDdigitizer():TNamed()
174 // AliTRDdigitizer default constructor
178 fDigitsManager = NULL;
179 fSDigitsManagerList = NULL;
180 fSDigitsManager = NULL;
202 fDriftVelocity = 0.0;
225 //_____________________________________________________________________________
226 AliTRDdigitizer::AliTRDdigitizer(const Text_t *name, const Text_t *title)
230 // AliTRDdigitizer default constructor
235 fDigitsManager = NULL;
236 fSDigitsManager = NULL;
237 fSDigitsManagerList = NULL;
254 //_____________________________________________________________________________
255 AliTRDdigitizer::AliTRDdigitizer(const AliTRDdigitizer &d)
258 // AliTRDdigitizer copy constructor
261 ((AliTRDdigitizer &) d).Copy(*this);
265 //_____________________________________________________________________________
266 AliTRDdigitizer::~AliTRDdigitizer()
269 // AliTRDdigitizer destructor
278 if (fDigitsManager) {
279 delete fDigitsManager;
280 fDigitsManager = NULL;
283 if (fSDigitsManager) {
284 delete fSDigitsManager;
285 fSDigitsManager = NULL;
288 if (fSDigitsManagerList) {
289 fSDigitsManagerList->Delete();
290 delete fSDigitsManagerList;
291 fSDigitsManagerList = NULL;
295 //_____________________________________________________________________________
296 AliTRDdigitizer &AliTRDdigitizer::operator=(const AliTRDdigitizer &d)
299 // Assignment operator
302 if (this != &d) ((AliTRDdigitizer &) d).Copy(*this);
307 //_____________________________________________________________________________
308 void AliTRDdigitizer::Copy(TObject &d)
316 ((AliTRDdigitizer &) d).fInputFile = NULL;
317 ((AliTRDdigitizer &) d).fSDigitsManagerList = NULL;
318 ((AliTRDdigitizer &) d).fSDigitsManager = NULL;
319 ((AliTRDdigitizer &) d).fDigitsManager = NULL;
320 ((AliTRDdigitizer &) d).fTRD = NULL;
321 ((AliTRDdigitizer &) d).fGeo = NULL;
323 ((AliTRDdigitizer &) d).fEvent = 0;
325 ((AliTRDdigitizer &) d).fGasGain = fGasGain;
326 ((AliTRDdigitizer &) d).fNoise = fNoise;
327 ((AliTRDdigitizer &) d).fChipGain = fChipGain;
328 ((AliTRDdigitizer &) d).fADCoutRange = fADCoutRange;
329 ((AliTRDdigitizer &) d).fADCinRange = fADCinRange;
330 ((AliTRDdigitizer &) d).fADCthreshold = fADCthreshold;
331 ((AliTRDdigitizer &) d).fDiffusionOn = fDiffusionOn;
332 ((AliTRDdigitizer &) d).fDiffusionT = fDiffusionT;
333 ((AliTRDdigitizer &) d).fDiffusionL = fDiffusionL;
334 ((AliTRDdigitizer &) d).fElAttachOn = fElAttachOn;
335 ((AliTRDdigitizer &) d).fElAttachProp = fElAttachProp;
336 ((AliTRDdigitizer &) d).fExBOn = fExBOn;
337 ((AliTRDdigitizer &) d).fOmegaTau = fOmegaTau;
338 ((AliTRDdigitizer &) d).fLorentzFactor = fLorentzFactor;
339 ((AliTRDdigitizer &) d).fDriftVelocity = fDriftVelocity;
340 ((AliTRDdigitizer &) d).fPadCoupling = fPadCoupling;
341 ((AliTRDdigitizer &) d).fTimeCoupling = fTimeCoupling;
342 ((AliTRDdigitizer &) d).fTimeBinWidth = fTimeBinWidth;
343 ((AliTRDdigitizer &) d).fField = fField;
344 ((AliTRDdigitizer &) d).fPRFOn = fPRFOn;
345 ((AliTRDdigitizer &) d).fTRFOn = fTRFOn;
347 ((AliTRDdigitizer &) d).fCompress = fCompress;
348 ((AliTRDdigitizer &) d).fVerbose = fVerbose;
349 ((AliTRDdigitizer &) d).fSDigits = fSDigits;
350 ((AliTRDdigitizer &) d).fSDigitsScale = fSDigitsScale;
352 ((AliTRDdigitizer &) d).fPRFbin = fPRFbin;
353 ((AliTRDdigitizer &) d).fPRFlo = fPRFlo;
354 ((AliTRDdigitizer &) d).fPRFhi = fPRFhi;
355 ((AliTRDdigitizer &) d).fPRFwid = fPRFwid;
356 ((AliTRDdigitizer &) d).fPRFpad = fPRFpad;
357 if (((AliTRDdigitizer &) d).fPRFsmp) delete ((AliTRDdigitizer &) d).fPRFsmp;
358 ((AliTRDdigitizer &) d).fPRFsmp = new Float_t[fPRFbin];
359 for (iBin = 0; iBin < fPRFbin; iBin++) {
360 ((AliTRDdigitizer &) d).fPRFsmp[iBin] = fPRFsmp[iBin];
362 ((AliTRDdigitizer &) d).fTRFbin = fTRFbin;
363 ((AliTRDdigitizer &) d).fTRFlo = fTRFlo;
364 ((AliTRDdigitizer &) d).fTRFhi = fTRFhi;
365 ((AliTRDdigitizer &) d).fTRFwid = fTRFwid;
366 if (((AliTRDdigitizer &) d).fTRFsmp) delete ((AliTRDdigitizer &) d).fTRFsmp;
367 ((AliTRDdigitizer &) d).fTRFsmp = new Float_t[fTRFbin];
368 for (iBin = 0; iBin < fTRFbin; iBin++) {
369 ((AliTRDdigitizer &) d).fTRFsmp[iBin] = fTRFsmp[iBin];
374 //_____________________________________________________________________________
375 Int_t AliTRDdigitizer::Diffusion(Float_t driftlength, Float_t *xyz)
378 // Applies the diffusion smearing to the position of a single electron
381 Float_t driftSqrt = TMath::Sqrt(driftlength);
382 Float_t sigmaT = driftSqrt * fDiffusionT;
383 Float_t sigmaL = driftSqrt * fDiffusionL;
384 xyz[0] = gRandom->Gaus(xyz[0], sigmaL * fLorentzFactor);
385 xyz[1] = gRandom->Gaus(xyz[1], sigmaT * fLorentzFactor);
386 xyz[2] = gRandom->Gaus(xyz[2], sigmaT);
392 //_____________________________________________________________________________
393 Int_t AliTRDdigitizer::ExB(Float_t driftlength, Float_t *xyz)
396 // Applies E x B effects to the position of a single electron
400 xyz[1] = xyz[1] + fOmegaTau * driftlength;
407 //_____________________________________________________________________________
408 Int_t AliTRDdigitizer::PadResponse(Float_t signal, Float_t dist, Float_t *pad)
411 // Applies the pad response
414 Int_t iBin = ((Int_t) (( - dist - fPRFlo) / fPRFwid));
416 Int_t iBin0 = iBin - fPRFpad;
418 Int_t iBin2 = iBin + fPRFpad;
420 if ((iBin0 >= 0) && (iBin2 < fPRFbin)) {
422 pad[0] = signal * fPRFsmp[iBin0];
423 pad[1] = signal * fPRFsmp[iBin1];
424 pad[2] = signal * fPRFsmp[iBin2];
437 //_____________________________________________________________________________
438 Float_t AliTRDdigitizer::TimeResponse(Float_t time)
441 // Applies the preamp shaper time response
444 Int_t iBin = ((Int_t) ((time - fTRFlo) / fTRFwid));
445 if ((iBin >= 0) && (iBin < fTRFbin)) {
446 return fTRFsmp[iBin];
454 //_____________________________________________________________________________
455 void AliTRDdigitizer::Init()
458 // Initializes the digitization procedure with standard values
461 // The default parameter for the digitization
465 fADCoutRange = 1023.; // 10-bit ADC
466 fADCinRange = 1000.; // 1V input range
469 // For the summable digits
470 fSDigitsScale = 100.;
472 // The drift velocity (cm / mus)
473 fDriftVelocity = 1.5;
481 // Propability for electron attachment
485 // The pad response function
488 // The time response function
491 // The pad coupling factor (same number as for the TPC)
494 // The time coupling factor (same number as for the TPC)
499 //_____________________________________________________________________________
500 void AliTRDdigitizer::ReInit()
503 // Reinitializes the digitization procedure after a change in the parameter
507 printf("AliTRDdigitizer::ReInit -- ");
508 printf("No geometry defined. Run InitDetector() first\n");
512 // Calculate the time bin width in ns
513 fTimeBinWidth = fGeo->GetTimeBinSize() / fDriftVelocity * 1000.0;
515 // The range and the binwidth for the sampled TRF
517 // Start 0.2 mus before the signal
518 fTRFlo = -0.2 * fDriftVelocity;
519 // End the maximum driftlength after the signal
520 fTRFhi = AliTRDgeometry::DrThick()
521 + fGeo->GetTimeAfter() * fGeo->GetTimeBinSize();
522 fTRFwid = (fTRFhi - fTRFlo) / ((Float_t) fTRFbin);
524 // Transverse and longitudinal diffusion coefficients (Xe/CO2)
525 fDiffusionT = GetDiffusionT(fDriftVelocity,fField);
526 fDiffusionL = GetDiffusionL(fDriftVelocity,fField);
528 // omega * tau.= tan(Lorentz-angle)
529 fOmegaTau = GetOmegaTau(fDriftVelocity,fField);
531 // The Lorentz factor
533 fLorentzFactor = 1.0 / (1.0 + fOmegaTau*fOmegaTau);
536 fLorentzFactor = 1.0;
541 //_____________________________________________________________________________
542 void AliTRDdigitizer::SampleTRF()
545 // Samples the time response function
546 // It is defined according to Vasiles simulation of the preamp shaper
547 // output and includes the effect of the ion tail (based on Tariqs
548 // Garfield simulation) and a shaping time of 125 ns FWHM
555 const Int_t kNpasa = 200;
556 Float_t time[kNpasa] = { -0.280000, -0.270000, -0.260000, -0.250000
557 , -0.240000, -0.230000, -0.220000, -0.210000
558 , -0.200000, -0.190000, -0.180000, -0.170000
559 , -0.160000, -0.150000, -0.140000, -0.130000
560 , -0.120000, -0.110000, -0.100000, -0.090000
561 , -0.080000, -0.070000, -0.060000, -0.050000
562 , -0.040000, -0.030000, -0.020000, -0.010000
563 , -0.000000, 0.010000, 0.020000, 0.030000
564 , 0.040000, 0.050000, 0.060000, 0.070000
565 , 0.080000, 0.090000, 0.100000, 0.110000
566 , 0.120000, 0.130000, 0.140000, 0.150000
567 , 0.160000, 0.170000, 0.180000, 0.190000
568 , 0.200000, 0.210000, 0.220000, 0.230000
569 , 0.240000, 0.250000, 0.260000, 0.270000
570 , 0.280000, 0.290000, 0.300000, 0.310000
571 , 0.320000, 0.330000, 0.340000, 0.350000
572 , 0.360000, 0.370000, 0.380000, 0.390000
573 , 0.400000, 0.410000, 0.420000, 0.430000
574 , 0.440000, 0.450000, 0.460000, 0.470000
575 , 0.480000, 0.490000, 0.500000, 0.510000
576 , 0.520000, 0.530000, 0.540000, 0.550000
577 , 0.560000, 0.570000, 0.580000, 0.590000
578 , 0.600000, 0.610000, 0.620000, 0.630000
579 , 0.640000, 0.650000, 0.660000, 0.670000
580 , 0.680000, 0.690000, 0.700000, 0.710000
581 , 0.720000, 0.730000, 0.740000, 0.750000
582 , 0.760000, 0.770000, 0.780000, 0.790000
583 , 0.800000, 0.810000, 0.820000, 0.830000
584 , 0.840000, 0.850000, 0.860000, 0.870000
585 , 0.880000, 0.890000, 0.900000, 0.910000
586 , 0.920000, 0.930000, 0.940000, 0.950000
587 , 0.960000, 0.970000, 0.980000, 0.990000
588 , 1.000000, 1.010000, 1.020000, 1.030000
589 , 1.040000, 1.050000, 1.060000, 1.070000
590 , 1.080000, 1.090000, 1.100000, 1.110000
591 , 1.120000, 1.130000, 1.140000, 1.150000
592 , 1.160000, 1.170000, 1.180000, 1.190000
593 , 1.200000, 1.210000, 1.220000, 1.230000
594 , 1.240000, 1.250000, 1.260000, 1.270000
595 , 1.280000, 1.290000, 1.300000, 1.310000
596 , 1.320000, 1.330000, 1.340000, 1.350000
597 , 1.360000, 1.370000, 1.380000, 1.390000
598 , 1.400000, 1.410000, 1.420000, 1.430000
599 , 1.440000, 1.450000, 1.460000, 1.470000
600 , 1.480000, 1.490000, 1.500000, 1.510000
601 , 1.520000, 1.530000, 1.540000, 1.550000
602 , 1.560000, 1.570000, 1.580000, 1.590000
603 , 1.600000, 1.610000, 1.620000, 1.630000
604 , 1.640000, 1.650000, 1.660000, 1.670000
605 , 1.680000, 1.690000, 1.700000, 1.710000 };
607 Float_t signal[kNpasa] = { 0.000000, 0.000000, 0.000000, 0.000000
608 , 0.000000, 0.000000, 0.000000, 0.000000
609 , 0.000000, 0.000000, 0.000000, 0.000000
610 , 0.000000, 0.000000, 0.000000, 0.000098
611 , 0.003071, 0.020056, 0.066053, 0.148346
612 , 0.263120, 0.398496, 0.540226, 0.674436
613 , 0.790977, 0.883083, 0.947744, 0.985714
614 , 0.999248, 0.992105, 0.967669, 0.930827
615 , 0.884586, 0.833083, 0.778571, 0.723684
616 , 0.669173, 0.617293, 0.567669, 0.521805
617 , 0.479699, 0.440977, 0.405639, 0.373985
618 , 0.345526, 0.320038, 0.297256, 0.276917
619 , 0.258797, 0.242632, 0.228195, 0.215301
620 , 0.203759, 0.193383, 0.184023, 0.175564
621 , 0.167895, 0.160940, 0.154549, 0.148722
622 , 0.143308, 0.138346, 0.133722, 0.129398
623 , 0.125376, 0.121617, 0.118045, 0.114699
624 , 0.111541, 0.108571, 0.105714, 0.103008
625 , 0.100414, 0.097970, 0.095602, 0.093346
626 , 0.091165, 0.089060, 0.087068, 0.085150
627 , 0.083308, 0.081541, 0.079812, 0.078158
628 , 0.076541, 0.075000, 0.073496, 0.072068
629 , 0.070677, 0.069286, 0.068008, 0.066729
630 , 0.065489, 0.064286, 0.063120, 0.061992
631 , 0.060902, 0.059850, 0.058797, 0.057820
632 , 0.056842, 0.055902, 0.054962, 0.054060
633 , 0.053158, 0.052293, 0.051466, 0.050639
634 , 0.049850, 0.049060, 0.048308, 0.047556
635 , 0.046842, 0.046128, 0.045451, 0.044774
636 , 0.044098, 0.043459, 0.042820, 0.042218
637 , 0.041617, 0.041015, 0.040451, 0.039887
638 , 0.039323, 0.038797, 0.038271, 0.037744
639 , 0.037237, 0.036744, 0.036259, 0.035786
640 , 0.035323, 0.034872, 0.034429, 0.033996
641 , 0.033575, 0.033162, 0.032756, 0.032361
642 , 0.031974, 0.031594, 0.031222, 0.030857
643 , 0.030496, 0.030143, 0.029793, 0.029451
644 , 0.029109, 0.028774, 0.028444, 0.028113
645 , 0.027793, 0.027477, 0.027165, 0.026861
646 , 0.026564, 0.026271, 0.025981, 0.025699
647 , 0.025421, 0.025147, 0.024880, 0.024613
648 , 0.024353, 0.024094, 0.023842, 0.023590
649 , 0.023346, 0.023102, 0.022865, 0.022628
650 , 0.022398, 0.022173, 0.021951, 0.021733
651 , 0.021519, 0.021308, 0.021098, 0.020891
652 , 0.020688, 0.020485, 0.020286, 0.020090
653 , 0.019895, 0.019707, 0.019519, 0.019335
654 , 0.019150, 0.018974, 0.018797, 0.018624
655 , 0.018451, 0.018282, 0.018113, 0.017947
656 , 0.017782, 0.017617, 0.017455, 0.017297 };
658 if (fTRFsmp) delete fTRFsmp;
659 fTRFsmp = new Float_t[fTRFbin];
661 Float_t loTRF = TMath::Max(fTRFlo / fDriftVelocity,time[0]);
662 Float_t hiTRF = TMath::Min(fTRFhi / fDriftVelocity,time[kNpasa-1]);
663 Float_t binWidth = (hiTRF - loTRF) / ((Float_t) fTRFbin);
665 // Take the linear interpolation
666 for (Int_t iBin = 0; iBin < fTRFbin; iBin++) {
668 Float_t bin = (((Float_t) iBin) + 0.5) * binWidth + loTRF;
672 diff = bin - time[ipos2++];
675 if (ipos2 > kNpasa) ipos2 = kNpasa - 1;
678 fTRFsmp[iBin] = signal[ipos2]
679 + diff * (signal[ipos2] - signal[ipos1])
680 / ( time[ipos2] - time[ipos1]);
686 //_____________________________________________________________________________
687 void AliTRDdigitizer::SamplePRF()
690 // Samples the pad response function
693 const Int_t kPRFbin = 61;
694 Float_t prf[kPRFbin] = { 0.002340, 0.003380, 0.004900, 0.007080, 0.010220
695 , 0.014740, 0.021160, 0.030230, 0.042800, 0.059830
696 , 0.082030, 0.109700, 0.142550, 0.179840, 0.220610
697 , 0.263980, 0.309180, 0.355610, 0.402790, 0.450350
698 , 0.497930, 0.545190, 0.591740, 0.637100, 0.680610
699 , 0.721430, 0.758400, 0.790090, 0.814720, 0.830480
700 , 0.835930, 0.830480, 0.814710, 0.790070, 0.758390
701 , 0.721410, 0.680590, 0.637080, 0.591730, 0.545180
702 , 0.497920, 0.450340, 0.402790, 0.355610, 0.309190
703 , 0.263990, 0.220630, 0.179850, 0.142570, 0.109720
704 , 0.082040, 0.059830, 0.042820, 0.030230, 0.021170
705 , 0.014740, 0.010230, 0.007080, 0.004900, 0.003380
711 fPRFwid = (fPRFhi - fPRFlo) / ((Float_t) fPRFbin);
712 fPRFpad = ((Int_t) (1.0 / fPRFwid));
714 if (fPRFsmp) delete fPRFsmp;
715 fPRFsmp = new Float_t[fPRFbin];
716 for (Int_t iBin = 0; iBin < fPRFbin; iBin++) {
717 fPRFsmp[iBin] = prf[iBin];
722 //_____________________________________________________________________________
723 Bool_t AliTRDdigitizer::Open(const Char_t *name, Int_t nEvent)
726 // Opens a ROOT-file with TRD-hits and reads in the hit-tree
729 // Connect the AliRoot file containing Geometry, Kine, and Hits
730 fInputFile = (TFile*) gROOT->GetListOfFiles()->FindObject(name);
733 printf("AliTRDdigitizer::Open -- ");
734 printf("Open the AliROOT-file %s.\n",name);
736 fInputFile = new TFile(name,"UPDATE");
740 printf("AliTRDdigitizer::Open -- ");
741 printf("%s is already open.\n",name);
745 gAlice = (AliRun*) fInputFile->Get("gAlice");
748 printf("AliTRDdigitizer::Open -- ");
749 printf("AliRun object found on file.\n");
753 printf("AliTRDdigitizer::Open -- ");
754 printf("Could not find AliRun object.\n");
760 // Import the Trees for the event nEvent in the file
761 Int_t nparticles = gAlice->GetEvent(fEvent);
762 if (nparticles <= 0) {
763 printf("AliTRDdigitizer::Open -- ");
764 printf("No entries in the trees for event %d.\n",fEvent);
768 if (InitDetector()) {
777 //_____________________________________________________________________________
778 Bool_t AliTRDdigitizer::InitDetector()
781 // Sets the pointer to the TRD detector and the geometry
784 // Get the pointer to the detector class and check for version 1
785 fTRD = (AliTRD*) gAlice->GetDetector("TRD");
786 if (fTRD->IsVersion() != 1) {
787 printf("AliTRDdigitizer::InitDetector -- ");
788 printf("TRD must be version 1 (slow simulator).\n");
793 fGeo = fTRD->GetGeometry();
795 printf("AliTRDdigitizer::InitDetector -- ");
796 printf("Geometry version %d\n",fGeo->IsVersion());
799 // The magnetic field strength in Tesla
800 fField = 0.2 * gAlice->Field()->Factor();
802 // Create a digits manager
803 fDigitsManager = new AliTRDdigitsManager();
804 fDigitsManager->SetSDigits(fSDigits);
805 fDigitsManager->CreateArrays();
806 fDigitsManager->SetEvent(fEvent);
807 fDigitsManager->SetVerbose(fVerbose);
809 // The list for the input s-digits manager to be merged
810 fSDigitsManagerList = new TList();
818 //_____________________________________________________________________________
819 Bool_t AliTRDdigitizer::MakeBranch(const Char_t *file)
822 // Create the branches for the digits array
825 return fDigitsManager->MakeBranch(file);
829 //_____________________________________________________________________________
830 Bool_t AliTRDdigitizer::MakeDigits()
836 ///////////////////////////////////////////////////////////////
838 ///////////////////////////////////////////////////////////////
840 // Converts number of electrons to fC
841 const Double_t kEl2fC = 1.602E-19 * 1.0E15;
843 ///////////////////////////////////////////////////////////////
845 // Number of pads included in the pad response
846 const Int_t kNpad = 3;
848 // Number of track dictionary arrays
849 const Int_t kNDict = AliTRDdigitsManager::kNDict;
851 // Half the width of the amplification region
852 const Float_t kAmWidth = AliTRDgeometry::AmThick() / 2.;
854 Int_t iRow, iCol, iTime, iPad;
858 Int_t totalSizeDigits = 0;
859 Int_t totalSizeDict0 = 0;
860 Int_t totalSizeDict1 = 0;
861 Int_t totalSizeDict2 = 0;
863 Int_t timeTRDbeg = 0;
864 Int_t timeTRDend = 1;
869 Float_t padSignal[kNpad];
870 Float_t signalOld[kNpad];
872 AliTRDdataArrayF *signals = 0;
873 AliTRDdataArrayI *digits = 0;
874 AliTRDdataArrayI *dictionary[kNDict];
876 // Create a container for the amplitudes
877 AliTRDsegmentArray *signalsArray
878 = new AliTRDsegmentArray("AliTRDdataArrayF",AliTRDgeometry::Ndet());
881 timeTRDbeg = ((Int_t) (-fTRFlo / fGeo->GetTimeBinSize())) - 1;
882 timeTRDend = ((Int_t) ( fTRFhi / fGeo->GetTimeBinSize())) - 1;
884 printf("AliTRDdigitizer::MakeDigits -- ");
885 printf("Sample the TRF between -%d and %d\n",timeTRDbeg,timeTRDend);
889 Float_t elAttachProp = fElAttachProp / 100.;
891 // Create the sampled PRF
894 // Create the sampled TRF
898 printf("AliTRDdigitizer::MakeDigits -- ");
899 printf("No geometry defined\n");
904 printf("AliTRDdigitizer::MakeDigits -- ");
905 printf("Start creating digits.\n");
908 // Get the pointer to the hit tree
909 TTree *HitTree = gAlice->TreeH();
911 // Get the number of entries in the hit tree
912 // (Number of primary particles creating a hit somewhere)
913 Int_t nTrack = (Int_t) HitTree->GetEntries();
915 printf("AliTRDdigitizer::MakeDigits -- ");
916 printf("Found %d primary particles\n",nTrack);
919 Int_t detectorOld = -1;
922 // Loop through all entries in the tree
923 for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
926 nBytes += HitTree->GetEvent(iTrack);
928 // Get the number of hits in the TRD created by this particle
929 Int_t nHit = fTRD->Hits()->GetEntriesFast();
931 printf("AliTRDdigitizer::MakeDigits -- ");
932 printf("Found %d hits for primary particle %d\n",nHit,iTrack);
935 // Loop through the TRD hits
936 for (Int_t iHit = 0; iHit < nHit; iHit++) {
940 AliTRDhit *hit = (AliTRDhit *) fTRD->Hits()->UncheckedAt(iHit);
944 Float_t q = hit->GetCharge();
945 Int_t track = hit->Track();
946 Int_t detector = hit->GetDetector();
947 Int_t plane = fGeo->GetPlane(detector);
948 Int_t sector = fGeo->GetSector(detector);
949 Int_t chamber = fGeo->GetChamber(detector);
951 if (!(CheckDetector(plane,chamber,sector))) continue;
953 Int_t nRowMax = fGeo->GetRowMax(plane,chamber,sector);
954 Int_t nColMax = fGeo->GetColMax(plane);
955 Int_t nTimeMax = fGeo->GetTimeMax();
956 Int_t nTimeBefore = fGeo->GetTimeBefore();
957 Int_t nTimeAfter = fGeo->GetTimeAfter();
958 Int_t nTimeTotal = fGeo->GetTimeTotal();
959 Float_t row0 = fGeo->GetRow0(plane,chamber,sector);
960 Float_t col0 = fGeo->GetCol0(plane);
961 Float_t time0 = fGeo->GetTime0(plane);
962 Float_t rowPadSize = fGeo->GetRowPadSize(plane,chamber,sector);
963 Float_t colPadSize = fGeo->GetColPadSize(plane);
964 Float_t timeBinSize = fGeo->GetTimeBinSize();
965 Float_t divideRow = 1.0 / rowPadSize;
966 Float_t divideCol = 1.0 / colPadSize;
967 Float_t divideTime = 1.0 / timeBinSize;
970 printf("Analyze hit no. %d ",iHit);
971 printf("-----------------------------------------------------------\n");
973 printf("plane = %d, sector = %d, chamber = %d\n"
974 ,plane,sector,chamber);
975 printf("nRowMax = %d, nColMax = %d, nTimeMax = %d\n"
976 ,nRowMax,nColMax,nTimeMax);
977 printf("nTimeBefore = %d, nTimeAfter = %d, nTimeTotal = %d\n"
978 ,nTimeBefore,nTimeAfter,nTimeTotal);
979 printf("row0 = %f, col0 = %f, time0 = %f\n"
981 printf("rowPadSize = %f, colPadSize = %f, timeBinSize = %f\n"
982 ,rowPadSize,colPadSize,timeBinSize);
985 // Don't analyze test hits
986 if (hit->FromTest()) continue;
988 if (detector != detectorOld) {
991 printf("AliTRDdigitizer::MakeDigits -- ");
992 printf("Get new container. New det = %d, Old det = %d\n"
993 ,detector,detectorOld);
995 // Compress the old one if enabled
996 if ((fCompress) && (detectorOld > -1)) {
998 printf("AliTRDdigitizer::MakeDigits -- ");
999 printf("Compress the old container ...");
1001 signals->Compress(1,0);
1002 for (iDict = 0; iDict < kNDict; iDict++) {
1003 dictionary[iDict]->Compress(1,0);
1005 if (fVerbose > 1) printf("done\n");
1007 // Get the new container
1008 signals = (AliTRDdataArrayF *) signalsArray->At(detector);
1009 if (signals->GetNtime() == 0) {
1010 // Allocate a new one if not yet existing
1012 printf("AliTRDdigitizer::MakeDigits -- ");
1013 printf("Allocate a new container ... ");
1015 signals->Allocate(nRowMax,nColMax,nTimeTotal);
1018 // Expand an existing one
1021 printf("AliTRDdigitizer::MakeDigits -- ");
1022 printf("Expand an existing container ... ");
1027 // The same for the dictionary
1028 for (iDict = 0; iDict < kNDict; iDict++) {
1029 dictionary[iDict] = fDigitsManager->GetDictionary(detector,iDict);
1030 if (dictionary[iDict]->GetNtime() == 0) {
1031 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1034 if (fCompress) dictionary[iDict]->Expand();
1037 if (fVerbose > 1) printf("done\n");
1038 detectorOld = detector;
1041 // Rotate the sectors on top of each other
1042 fGeo->Rotate(detector,pos,rot);
1044 // The driftlength. It is negative if the hit is in the
1045 // amplification region.
1046 Float_t driftlength = time0 - rot[0];
1048 // Take also the drift in the amplification region into account
1049 // The drift length is at the moment still the same, regardless of
1050 // the position relativ to the wire. This non-isochronity needs still
1051 // to be implemented.
1052 Float_t driftlengthL = TMath::Abs(driftlength + kAmWidth);
1053 if (fExBOn) driftlengthL /= TMath::Sqrt(fLorentzFactor);
1055 // Loop over all electrons of this hit
1056 // TR photons produce hits with negative charge
1057 Int_t nEl = ((Int_t) TMath::Abs(q));
1058 for (Int_t iEl = 0; iEl < nEl; iEl++) {
1064 // Electron attachment
1066 if (gRandom->Rndm() < (driftlengthL * elAttachProp))
1070 // Apply the diffusion smearing
1072 if (!(Diffusion(driftlengthL,xyz))) continue;
1075 // Apply E x B effects (depends on drift direction)
1077 if (!(ExB(driftlength+kAmWidth,xyz))) continue;
1080 // The electron position after diffusion and ExB in pad coordinates
1081 // The pad row (z-direction)
1082 Int_t rowE = ((Int_t) ((xyz[2] - row0) * divideRow));
1083 if ((rowE < 0) || (rowE >= nRowMax)) continue;
1085 // The pad column (rphi-direction)
1086 Int_t colE = ((Int_t) ((xyz[1] - col0) * divideCol));
1087 if ((colE < 0) || (colE >= nColMax)) continue;
1089 // The time bin (negative for hits in the amplification region)
1090 // In the amplification region the electrons drift from both sides
1091 // to the middle (anode wire plane)
1092 Float_t timeDist = time0 - xyz[0];
1093 Float_t timeOffset = 0;
1097 timeE = ((Int_t) (timeDist * divideTime));
1098 // The distance of the position to the middle of the timebin
1099 timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) - timeDist;
1102 // Difference between half of the amplification gap width and
1103 // the distance to the anode wire
1104 Float_t anodeDist = kAmWidth - TMath::Abs(timeDist + kAmWidth);
1106 timeE = -1 * (((Int_t ) (anodeDist * divideTime)) + 1);
1107 // The distance of the position to the middle of the timebin
1108 timeOffset = ((((Float_t) timeE) + 0.5) * timeBinSize) + anodeDist;
1111 // Apply the gas gain including fluctuations
1112 Float_t ggRndm = 0.0;
1114 ggRndm = gRandom->Rndm();
1115 } while (ggRndm <= 0);
1116 Int_t signal = (Int_t) (-fGasGain * TMath::Log(ggRndm));
1118 // Apply the pad response
1120 // The distance of the electron to the center of the pad
1121 // in units of pad width
1122 Float_t dist = (xyz[1] - col0 - (colE + 0.5) * colPadSize)
1124 if (!(PadResponse(signal,dist,padSignal))) continue;
1128 padSignal[1] = signal;
1132 // Sample the time response inside the drift region
1133 // + additional time bins before and after.
1134 // The sampling is done always in the middle of the time bin
1135 for (Int_t iTimeBin = TMath::Max(timeE-timeTRDbeg, -nTimeBefore)
1136 ;iTimeBin < TMath::Min(timeE+timeTRDend,nTimeMax+nTimeAfter )
1139 // Apply the time response
1140 Float_t timeResponse = 1.0;
1142 Float_t time = (iTimeBin - timeE) * timeBinSize + timeOffset;
1143 timeResponse = TimeResponse(time);
1150 for (iPad = 0; iPad < kNpad; iPad++) {
1152 Int_t colPos = colE + iPad - 1;
1153 if (colPos < 0) continue;
1154 if (colPos >= nColMax) break;
1157 // Note: The time bin number is shifted by nTimeBefore to avoid negative
1158 // time bins. This has to be subtracted later.
1159 Int_t iCurrentTimeBin = iTimeBin + nTimeBefore;
1160 signalOld[iPad] = signals->GetDataUnchecked(rowE,colPos,iCurrentTimeBin);
1161 signalOld[iPad] += padSignal[iPad] * timeResponse;
1162 signals->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,signalOld[iPad]);
1164 // Store the track index in the dictionary
1165 // Note: We store index+1 in order to allow the array to be compressed
1166 if (signalOld[iPad] > 0) {
1167 for (iDict = 0; iDict < kNDict; iDict++) {
1168 Int_t oldTrack = dictionary[iDict]->GetDataUnchecked(rowE
1171 if (oldTrack == track+1) break;
1172 if (oldTrack == 0) {
1173 dictionary[iDict]->SetDataUnchecked(rowE,colPos,iCurrentTimeBin,track+1);
1187 } // All hits finished
1190 printf("AliTRDdigitizer::MakeDigits -- ");
1191 printf("Finished analyzing %d hits\n",countHits);
1194 // The total conversion factor
1195 Float_t convert = kEl2fC * fPadCoupling * fTimeCoupling * fChipGain;
1197 // Loop through all chambers to finalize the digits
1198 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1200 Int_t plane = fGeo->GetPlane(iDet);
1201 Int_t sector = fGeo->GetSector(iDet);
1202 Int_t chamber = fGeo->GetChamber(iDet);
1203 Int_t nRowMax = fGeo->GetRowMax(plane,chamber,sector);
1204 Int_t nColMax = fGeo->GetColMax(plane);
1205 Int_t nTimeMax = fGeo->GetTimeMax();
1206 Int_t nTimeTotal = fGeo->GetTimeTotal();
1209 printf("AliTRDdigitizer::MakeDigits -- ");
1210 printf("Digitization for chamber %d\n",iDet);
1213 // Add a container for the digits of this detector
1214 digits = fDigitsManager->GetDigits(iDet);
1215 // Allocate memory space for the digits buffer
1216 digits->Allocate(nRowMax,nColMax,nTimeTotal);
1218 // Get the signal container
1219 signals = (AliTRDdataArrayF *) signalsArray->At(iDet);
1220 if (signals->GetNtime() == 0) {
1221 // Create missing containers
1222 signals->Allocate(nRowMax,nColMax,nTimeTotal);
1225 // Expand the container if neccessary
1226 if (fCompress) signals->Expand();
1228 // Create the missing dictionary containers
1229 for (iDict = 0; iDict < kNDict; iDict++) {
1230 dictionary[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1231 if (dictionary[iDict]->GetNtime() == 0) {
1232 dictionary[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1238 // Don't create noise in detectors that are switched off
1239 if (CheckDetector(plane,chamber,sector)) {
1241 // Create the digits for this chamber
1242 for (iRow = 0; iRow < nRowMax; iRow++ ) {
1243 for (iCol = 0; iCol < nColMax; iCol++ ) {
1244 for (iTime = 0; iTime < nTimeTotal; iTime++) {
1246 // Create summable digits
1249 Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1250 signalAmp *= fSDigitsScale;
1251 signalAmp = TMath::Min(signalAmp,1.0e9);
1252 Int_t adc = (Int_t) signalAmp;
1254 digits->SetDataUnchecked(iRow,iCol,iTime,adc);
1257 // Create normal digits
1260 Float_t signalAmp = signals->GetDataUnchecked(iRow,iCol,iTime);
1263 signalAmp = TMath::Max((Double_t) gRandom->Gaus(signalAmp,fNoise),0.0);
1265 signalAmp *= convert;
1266 // Convert to ADC counts. Set the overflow-bit fADCoutRange if the
1267 // signal is larger than fADCinRange
1269 if (signalAmp >= fADCinRange) {
1270 adc = ((Int_t) fADCoutRange);
1273 adc = ((Int_t) (signalAmp * (fADCoutRange / fADCinRange)));
1276 // Store the amplitude of the digit if above threshold
1277 if (adc > fADCthreshold) {
1279 printf(" iRow = %d, iCol = %d, iTime = %d\n"
1281 printf(" signal = %f, adc = %d\n",signalAmp,adc);
1284 digits->SetDataUnchecked(iRow,iCol,iTime,adc);
1295 // Compress the arrays
1296 digits->Compress(1,0);
1297 for (iDict = 0; iDict < kNDict; iDict++) {
1298 dictionary[iDict]->Compress(1,0);
1301 totalSizeDigits += digits->GetSize();
1302 totalSizeDict0 += dictionary[0]->GetSize();
1303 totalSizeDict1 += dictionary[1]->GetSize();
1304 totalSizeDict2 += dictionary[2]->GetSize();
1306 Float_t nPixel = nRowMax * nColMax * nTimeMax;
1308 printf("AliTRDdigitizer::MakeDigits -- ");
1309 printf("Found %d digits in detector %d (%3.0f).\n"
1311 ,100.0 * ((Float_t) nDigits) / nPixel);
1314 if (fCompress) signals->Compress(1,0);
1319 printf("AliTRDdigitizer::MakeDigits -- ");
1320 printf("Total number of analyzed hits = %d\n",countHits);
1321 printf("AliTRDdigitizer::MakeDigits -- ");
1322 printf("Total digits data size = %d, %d, %d, %d\n",totalSizeDigits
1332 //_____________________________________________________________________________
1333 void AliTRDdigitizer::AddSDigitsManager(AliTRDdigitsManager *man)
1336 // Add a digits manager for s-digits to the input list.
1339 fSDigitsManagerList->Add(man);
1343 //_____________________________________________________________________________
1344 Bool_t AliTRDdigitizer::ConvertSDigits()
1347 // Converts s-digits to normal digits
1350 // Number of track dictionary arrays
1351 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1353 // Converts number of electrons to fC
1354 const Double_t kEl2fC = 1.602E-19 * 1.0E15;
1362 Double_t sDigitsScale = 1.0 / GetSDigitsScale();
1363 Double_t noise = GetNoise();
1364 Double_t padCoupling = GetPadCoupling();
1365 Double_t timeCoupling = GetTimeCoupling();
1366 Double_t chipGain = GetChipGain();
1367 Double_t convert = kEl2fC * padCoupling * timeCoupling * chipGain;;
1368 Double_t adcInRange = GetADCinRange();
1369 Double_t adcOutRange = GetADCoutRange();
1370 Int_t adcThreshold = GetADCthreshold();
1372 AliTRDdataArrayI *digitsIn;
1373 AliTRDdataArrayI *digitsOut;
1374 AliTRDdataArrayI *dictionaryIn[kNDict];
1375 AliTRDdataArrayI *dictionaryOut[kNDict];
1377 // Loop through the detectors
1378 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1381 printf("AliTRDdigitizer::ConvertSDigits -- ");
1382 printf("Convert detector %d to digits.\n",iDet);
1385 Int_t plane = fGeo->GetPlane(iDet);
1386 Int_t sector = fGeo->GetSector(iDet);
1387 Int_t chamber = fGeo->GetChamber(iDet);
1388 Int_t nRowMax = fGeo->GetRowMax(plane,chamber,sector);
1389 Int_t nColMax = fGeo->GetColMax(plane);
1390 Int_t nTimeTotal = fGeo->GetTimeTotal();
1392 digitsIn = fSDigitsManager->GetDigits(iDet);
1394 digitsOut = fDigitsManager->GetDigits(iDet);
1395 digitsOut->Allocate(nRowMax,nColMax,nTimeTotal);
1396 for (iDict = 0; iDict < kNDict; iDict++) {
1397 dictionaryIn[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1398 dictionaryIn[iDict]->Expand();
1399 dictionaryOut[iDict] = fDigitsManager->GetDictionary(iDet,iDict);
1400 dictionaryOut[iDict]->Allocate(nRowMax,nColMax,nTimeTotal);
1403 for (Int_t iRow = 0; iRow < nRowMax; iRow++ ) {
1404 for (Int_t iCol = 0; iCol < nColMax; iCol++ ) {
1405 for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {
1407 Double_t signal = (Double_t) digitsIn->GetDataUnchecked(iRow,iCol,iTime);
1408 signal *= sDigitsScale;
1410 signal = TMath::Max((Double_t) gRandom->Gaus(signal,noise),0.0);
1413 // Convert to ADC counts. Set the overflow-bit adcOutRange if the
1414 // signal is larger than adcInRange
1416 if (signal >= adcInRange) {
1417 adc = ((Int_t) adcOutRange);
1420 adc = ((Int_t) (signal * (adcOutRange / adcInRange)));
1422 // Store the amplitude of the digit if above threshold
1423 if (adc > adcThreshold) {
1424 digitsOut->SetDataUnchecked(iRow,iCol,iTime,adc);
1426 // Copy the dictionary
1427 for (iDict = 0; iDict < kNDict; iDict++) {
1428 Int_t track = dictionaryIn[iDict]->GetDataUnchecked(iRow,iCol,iTime);
1429 dictionaryOut[iDict]->SetDataUnchecked(iRow,iCol,iTime,track);
1437 digitsIn->Compress(1,0);
1438 digitsOut->Compress(1,0);
1439 for (iDict = 0; iDict < kNDict; iDict++) {
1440 dictionaryIn[iDict]->Compress(1,0);
1441 dictionaryOut[iDict]->Compress(1,0);
1451 //_____________________________________________________________________________
1452 Bool_t AliTRDdigitizer::MergeSDigits()
1455 // Merges the input s-digits:
1456 // - The amplitude of the different inputs are summed up.
1457 // - Of the track IDs from the input dictionaries only one is
1458 // kept for each input. This works for maximal 3 different merged inputs.
1461 // Number of track dictionary arrays
1462 const Int_t kNDict = AliTRDdigitsManager::kNDict;
1466 AliTRDdataArrayI *digitsA;
1467 AliTRDdataArrayI *digitsB;
1468 AliTRDdataArrayI *dictionaryA[kNDict];
1469 AliTRDdataArrayI *dictionaryB[kNDict];
1471 // Get the first s-digits
1472 fSDigitsManager = (AliTRDdigitsManager *) fSDigitsManagerList->First();
1473 if (!fSDigitsManager) return kFALSE;
1475 // Loop through the other sets of s-digits
1476 AliTRDdigitsManager *mergeSDigitsManager;
1477 mergeSDigitsManager = (AliTRDdigitsManager *)
1478 fSDigitsManagerList->After(fSDigitsManager);
1481 if (mergeSDigitsManager) {
1482 printf("AliTRDdigitizer::MergeSDigits -- ");
1483 printf("Merge serveral input files.\n");
1486 printf("AliTRDdigitizer::MergeSDigits -- ");
1487 printf("Only one input file.\n");
1492 while (mergeSDigitsManager) {
1496 // Loop through the detectors
1497 for (Int_t iDet = 0; iDet < AliTRDgeometry::Ndet(); iDet++) {
1499 Int_t plane = fGeo->GetPlane(iDet);
1500 Int_t sector = fGeo->GetSector(iDet);
1501 Int_t chamber = fGeo->GetChamber(iDet);
1502 Int_t nRowMax = fGeo->GetRowMax(plane,chamber,sector);
1503 Int_t nColMax = fGeo->GetColMax(plane);
1504 Int_t nTimeTotal = fGeo->GetTimeTotal();
1506 // Loop through the pixels of one detector and add the signals
1507 digitsA = fSDigitsManager->GetDigits(iDet);
1508 digitsB = mergeSDigitsManager->GetDigits(iDet);
1511 for (iDict = 0; iDict < kNDict; iDict++) {
1512 dictionaryA[iDict] = fSDigitsManager->GetDictionary(iDet,iDict);
1513 dictionaryB[iDict] = mergeSDigitsManager->GetDictionary(iDet,iDict);
1514 dictionaryA[iDict]->Expand();
1515 dictionaryB[iDict]->Expand();
1519 printf("AliTRDdigitizer::MergeSDigits -- ");
1520 printf("Merge detector %d of input no.%d.\n",iDet,iMerge);
1523 for (Int_t iRow = 0; iRow < nRowMax; iRow++ ) {
1524 for (Int_t iCol = 0; iCol < nColMax; iCol++ ) {
1525 for (Int_t iTime = 0; iTime < nTimeTotal; iTime++) {
1527 // Add the amplitudes of the summable digits
1528 Int_t ampA = digitsA->GetDataUnchecked(iRow,iCol,iTime);
1529 Int_t ampB = digitsB->GetDataUnchecked(iRow,iCol,iTime);
1531 digitsA->SetDataUnchecked(iRow,iCol,iTime,ampA);
1533 // Take only one track from each input
1534 Int_t track = dictionaryB[0]->GetDataUnchecked(iRow,iCol,iTime);
1535 if (iMerge < kNDict) {
1536 dictionaryA[iMerge]->SetDataUnchecked(iRow,iCol,iTime,track);
1544 digitsA->Compress(1,0);
1545 digitsB->Compress(1,0);
1546 for (iDict = 0; iDict < kNDict; iDict++) {
1547 dictionaryA[iDict]->Compress(1,0);
1548 dictionaryB[iDict]->Compress(1,0);
1554 // The next set of s-digits
1555 mergeSDigitsManager = (AliTRDdigitsManager *)
1556 fSDigitsManagerList->After(mergeSDigitsManager);
1564 //_____________________________________________________________________________
1565 Bool_t AliTRDdigitizer::SDigits2Digits()
1568 // Merges the input s-digits and converts them to normal digits
1571 if (!MergeSDigits()) return kFALSE;
1573 return ConvertSDigits();
1577 //_____________________________________________________________________________
1578 Bool_t AliTRDdigitizer::CheckDetector(Int_t plane, Int_t chamber, Int_t sector)
1581 // Checks whether a detector is enabled
1584 if ((fTRD->GetSensChamber() >= 0) &&
1585 (fTRD->GetSensChamber() != chamber)) return kFALSE;
1586 if ((fTRD->GetSensPlane() >= 0) &&
1587 (fTRD->GetSensPlane() != plane)) return kFALSE;
1588 if ( fTRD->GetSensSector() >= 0) {
1589 Int_t sens1 = fTRD->GetSensSector();
1590 Int_t sens2 = sens1 + fTRD->GetSensSectorRange();
1591 sens2 -= ((Int_t) (sens2 / AliTRDgeometry::Nsect()))
1592 * AliTRDgeometry::Nsect();
1593 if (sens1 < sens2) {
1594 if ((sector < sens1) || (sector >= sens2)) return kFALSE;
1597 if ((sector < sens1) && (sector >= sens2)) return kFALSE;
1605 //_____________________________________________________________________________
1606 Bool_t AliTRDdigitizer::WriteDigits()
1609 // Writes out the TRD-digits and the dictionaries
1612 // Store the digits and the dictionary in the tree
1613 return fDigitsManager->WriteDigits();
1617 //_____________________________________________________________________________
1618 Float_t AliTRDdigitizer::GetDiffusionL(Float_t vd, Float_t b)
1621 // Returns the longitudinal diffusion coefficient for a given drift
1622 // velocity <vd> and a B-field <b> for Xe/CO2 (15%).
1623 // The values are according to a GARFIELD simulation.
1626 const Int_t kNb = 5;
1627 Float_t p0[kNb] = { 0.007440, 0.007493, 0.007513, 0.007672, 0.007831 };
1628 Float_t p1[kNb] = { 0.019252, 0.018912, 0.018636, 0.018012, 0.017343 };
1629 Float_t p2[kNb] = { -0.005042, -0.004926, -0.004867, -0.004650, -0.004424 };
1630 Float_t p3[kNb] = { 0.000195, 0.000189, 0.000195, 0.000182, 0.000169 };
1632 Int_t ib = ((Int_t) (10 * (b - 0.15)));
1633 ib = TMath::Max( 0,ib);
1634 ib = TMath::Min(kNb,ib);
1636 Float_t diff = p0[ib]
1639 + p3[ib] * vd*vd*vd;
1645 //_____________________________________________________________________________
1646 Float_t AliTRDdigitizer::GetDiffusionT(Float_t vd, Float_t b)
1649 // Returns the transverse diffusion coefficient for a given drift
1650 // velocity <vd> and a B-field <b> for Xe/CO2 (15%).
1651 // The values are according to a GARFIELD simulation.
1654 const Int_t kNb = 5;
1655 Float_t p0[kNb] = { 0.009550, 0.009599, 0.009674, 0.009757, 0.009850 };
1656 Float_t p1[kNb] = { 0.006667, 0.006539, 0.006359, 0.006153, 0.005925 };
1657 Float_t p2[kNb] = { -0.000853, -0.000798, -0.000721, -0.000635, -0.000541 };
1658 Float_t p3[kNb] = { 0.000131, 0.000122, 0.000111, 0.000098, 0.000085 };
1660 Int_t ib = ((Int_t) (10 * (b - 0.15)));
1661 ib = TMath::Max( 0,ib);
1662 ib = TMath::Min(kNb,ib);
1664 Float_t diff = p0[ib]
1667 + p3[ib] * vd*vd*vd;
1673 //_____________________________________________________________________________
1674 Float_t AliTRDdigitizer::GetOmegaTau(Float_t vd, Float_t b)
1677 // Returns omega*tau (tan(Lorentz-angle)) for a given drift velocity <vd>
1678 // and a B-field <b> for Xe/CO2 (15%).
1679 // The values are according to a GARFIELD simulation.
1682 const Int_t kNb = 5;
1683 Float_t p0[kNb] = { 0.004810, 0.007412, 0.010252, 0.013409, 0.016888 };
1684 Float_t p1[kNb] = { 0.054875, 0.081534, 0.107333, 0.131983, 0.155455 };
1685 Float_t p2[kNb] = { -0.008682, -0.012896, -0.016987, -0.020880, -0.024623 };
1686 Float_t p3[kNb] = { 0.000155, 0.000238, 0.000330, 0.000428, 0.000541 };
1688 Int_t ib = ((Int_t) (10 * (b - 0.15)));
1689 ib = TMath::Max( 0,ib);
1690 ib = TMath::Min(kNb,ib);
1692 Float_t alphaL = p0[ib]
1695 + p3[ib] * vd*vd*vd;
1697 return TMath::Tan(alphaL);