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.13 1999/11/02 17:14:51 fca
19 Correct ansi scoping not accepted by HP compilers
21 Revision 1.12 1999/11/02 16:35:56 fca
22 New version of TRD introduced
24 Revision 1.11 1999/11/01 20:41:51 fca
25 Added protections against using the wrong version of FRAME
27 Revision 1.10 1999/09/29 09:24:35 fca
28 Introduction of the Copyright and cvs Log
32 ///////////////////////////////////////////////////////////////////////////////
34 // Transition Radiation Detector version 2 -- slow simulator //
38 <img src="picts/AliTRDfullClass.gif">
43 ///////////////////////////////////////////////////////////////////////////////
50 #include "AliTRDmatrix.h"
57 //_____________________________________________________________________________
58 AliTRDv1::AliTRDv1(const char *name, const char *title)
62 // Standard constructor for Transition Radiation Detector version 2
92 SetBufferSize(128000);
96 //_____________________________________________________________________________
100 if (fDeltaE) delete fDeltaE;
104 //_____________________________________________________________________________
105 void AliTRDv1::CreateGeometry()
108 // Create the GEANT geometry for the Transition Radiation Detector - Version 2
109 // This version covers the full azimuth.
112 // Check that FRAME is there otherwise we have no place where to put the TRD
113 AliModule* FRAME = gAlice->GetModule("FRAME");
116 // Define the chambers
117 AliTRD::CreateGeometry();
121 //_____________________________________________________________________________
122 void AliTRDv1::CreateMaterials()
125 // Create materials for the Transition Radiation Detector version 2
128 AliTRD::CreateMaterials();
132 //_____________________________________________________________________________
133 void AliTRDv1::Diffusion(Float_t driftlength, Float_t *xyz)
136 // Applies the diffusion smearing to the position of a single electron
139 if ((driftlength > 0) &&
140 (driftlength < kDrThick)) {
141 Float_t driftSqrt = TMath::Sqrt(driftlength);
142 Float_t sigmaT = driftSqrt * fDiffusionT;
143 Float_t sigmaL = driftSqrt * fDiffusionL;
144 xyz[0] = gRandom->Gaus(xyz[0], sigmaL);
145 xyz[1] = gRandom->Gaus(xyz[1], sigmaT);
146 xyz[2] = gRandom->Gaus(xyz[2], sigmaT);
156 //_____________________________________________________________________________
157 void AliTRDv1::Hits2Digits()
160 // Creates TRD digits from hits. This procedure includes the following:
162 // - Gas gain including fluctuations
163 // - Pad-response (simple Gaussian approximation)
164 // - Electronics noise
165 // - Electronics gain
168 // The corresponding parameter can be adjusted via the various Set-functions.
169 // If these parameters are not explicitly set, default values are used (see
171 // To produce digits from a root-file with TRD-hits use the
172 // slowDigitsCreate.C macro.
175 printf("AliTRDv1::Hits2Digits -- Start creating digits\n");
177 ///////////////////////////////////////////////////////////////
179 ///////////////////////////////////////////////////////////////
181 // Converts number of electrons to fC
182 const Float_t el2fC = 1.602E-19 * 1.0E15;
184 ///////////////////////////////////////////////////////////////
192 // Get the pointer to the hit tree
193 TTree *HitTree = gAlice->TreeH();
194 // Get the pointer to the digits tree
195 TTree *DigitsTree = gAlice->TreeD();
197 // Get the number of entries in the hit tree
198 // (Number of primary particles creating a hit somewhere)
199 Int_t nTrack = (Int_t) HitTree->GetEntries();
202 Int_t chamEnd = kNcham;
203 if (fSensChamber) chamEnd = chamBeg = fSensChamber;
205 Int_t planEnd = kNplan;
206 if (fSensPlane) planEnd = planBeg = fSensPlane;
208 Int_t sectEnd = kNsect;
209 if (fSensSector) sectEnd = sectBeg = fSensSector;
211 // Loop through all the chambers
212 for (Int_t icham = chamBeg; icham < chamEnd; icham++) {
213 for (Int_t iplan = planBeg; iplan < planEnd; iplan++) {
214 for (Int_t isect = sectBeg; isect < sectEnd; isect++) {
218 printf("AliTRDv1::Hits2Digits -- Digitizing chamber %d, plane %d, sector %d\n"
219 ,icham+1,iplan+1,isect+1);
221 // Create a detector matrix to keep the signal and track numbers
222 AliTRDmatrix *matrix = new AliTRDmatrix(fRowMax[iplan][icham][isect]
225 ,isect+1,icham+1,iplan+1);
227 // Loop through all entries in the tree
228 for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
231 nBytes += HitTree->GetEvent(iTrack);
233 // Get the number of hits in the TRD created by this particle
234 Int_t nHit = fHits->GetEntriesFast();
236 // Loop through the TRD hits
237 for (Int_t iHit = 0; iHit < nHit; iHit++) {
239 if (!(TRDhit = (AliTRDhit *) fHits->UncheckedAt(iHit)))
242 Float_t x = TRDhit->fX;
243 Float_t y = TRDhit->fY;
244 Float_t z = TRDhit->fZ;
245 Float_t q = TRDhit->fQ;
246 Int_t track = TRDhit->fTrack;
247 Int_t plane = TRDhit->fPlane;
248 Int_t sector = TRDhit->fSector;
249 Int_t chamber = TRDhit->fChamber;
251 if ((sector != isect+1) ||
252 (plane != iplan+1) ||
253 (chamber != icham+1))
256 // Rotate the sectors on top of each other
257 Float_t phi = 2.0 * kPI / (Float_t) kNsect
258 * ((Float_t) sector - 0.5);
259 Float_t xRot = -x * TMath::Cos(phi) + y * TMath::Sin(phi);
260 Float_t yRot = x * TMath::Sin(phi) + y * TMath::Cos(phi);
263 // The hit position in pad coordinates (center pad)
264 // The pad row (z-direction)
265 Int_t rowH = (Int_t) ((zRot - fRow0[iplan][icham][isect]) / fRowPadSize);
266 // The pad column (rphi-direction)
267 Int_t colH = (Int_t) ((yRot - fCol0[iplan] ) / fColPadSize);
269 Int_t timeH = (Int_t) ((xRot - fTime0[iplan] ) / fTimeBinSize);
271 // Array to sum up the signal in a box surrounding the
273 const Int_t timeBox = 5;
274 const Int_t colBox = 7;
275 const Int_t rowBox = 5;
276 Float_t signalSum[rowBox][colBox][timeBox];
277 for (iRow = 0; iRow < rowBox; iRow++ ) {
278 for (Int_t iCol = 0; iCol < colBox; iCol++ ) {
279 for (Int_t iTime = 0; iTime < timeBox; iTime++) {
280 signalSum[iRow][iCol][iTime] = 0;
285 // Loop over all electrons of this hit
286 Int_t nEl = (Int_t) q;
287 for (Int_t iEl = 0; iEl < nEl; iEl++) {
289 // Apply the diffusion smearing
290 Float_t driftlength = xRot - fTime0[iplan];
295 Diffusion(driftlength,xyz);
297 // At this point absorption effects that depend on the
298 // driftlength could be taken into account.
300 // The electron position and the distance to the hit position
302 // The pad row (z-direction)
303 Int_t rowE = (Int_t) ((xyz[2] - fRow0[iplan][icham][isect]) / fRowPadSize);
304 Int_t rowD = rowH - rowE;
305 // The pad column (rphi-direction)
306 Int_t colE = (Int_t) ((xyz[1] - fCol0[iplan] ) / fColPadSize);
307 Int_t colD = colH - colE;
309 Int_t timeE = (Int_t) ((xyz[0] - fTime0[iplan] ) / fTimeBinSize);
310 Int_t timeD = timeH - timeE;
312 // Apply the gas gain including fluctuations
313 Int_t signal = (Int_t) (-fGasGain * TMath::Log(gRandom->Rndm()));
315 // The distance of the electron to the center of the pad
316 // in units of pad width
317 Float_t dist = (xyz[1] - fCol0[iplan] - (colE + 0.5) * fColPadSize)
320 // Sum up the signal in the different pixels
321 // and apply the pad response
322 Int_t rowIdx = rowD + (Int_t) ( rowBox / 2);
323 Int_t colIdx = colD + (Int_t) ( colBox / 2);
324 Int_t timeIdx = timeD + (Int_t) (timeBox / 2);
325 signalSum[rowIdx][colIdx-1][timeIdx] += PadResponse(dist-1.) * signal;
326 signalSum[rowIdx][colIdx ][timeIdx] += PadResponse(dist ) * signal;
327 signalSum[rowIdx][colIdx+1][timeIdx] += PadResponse(dist+1.) * signal;
331 // Add the padcluster to the detector matrix
332 for (iRow = 0; iRow < rowBox; iRow++ ) {
333 for (Int_t iCol = 0; iCol < colBox; iCol++ ) {
334 for (Int_t iTime = 0; iTime < timeBox; iTime++) {
336 Int_t rowB = rowH + iRow - (Int_t) ( rowBox / 2);
337 Int_t colB = colH + iCol - (Int_t) ( colBox / 2);
338 Int_t timeB = timeH + iTime - (Int_t) (timeBox / 2);
340 Float_t signalB = signalSum[iRow][iCol][iTime];
342 matrix->AddSignal(rowB,colB,timeB,signalB);
343 if (!(matrix->AddTrack(rowB,colB,timeB,track)))
344 printf(" More than three tracks in a pixel!\n");
355 // Create the hits for this chamber
356 for (Int_t iRow = 0; iRow < fRowMax[iplan][icham][isect]; iRow++ ) {
357 for (Int_t iCol = 0; iCol < fColMax[iplan] ; iCol++ ) {
358 for (Int_t iTime = 0; iTime < fTimeMax ; iTime++) {
360 Float_t signalAmp = matrix->GetSignal(iRow,iCol,iTime);
363 signalAmp = TMath::Max(gRandom->Gaus(signalAmp,fNoise),(Float_t) 0.0);
367 signalAmp *= fChipGain;
368 // Convert to ADC counts
369 Int_t adc = (Int_t) (signalAmp * (fADCoutRange / fADCinRange));
371 // Apply threshold on ADC value
372 if (adc > fADCthreshold) {
375 for (Int_t ii = 0; ii < 3; ii++) {
376 trackSave[ii] = matrix->GetTrack(iRow,iCol,iTime,ii);
380 digits[0] = matrix->GetSector();
381 digits[1] = matrix->GetChamber();
382 digits[2] = matrix->GetPlane();
388 // Add this digit to the TClonesArray
389 AddDigit(trackSave,digits);
398 printf("AliTRDv1::Hits2Digits -- Number of digits found: %d\n",nDigits);
407 // Fill the digits-tree
408 printf("AliTRDv1::Hits2Digits -- Fill the digits tree\n");
413 //_____________________________________________________________________________
414 void AliTRDv1::Digits2Clusters()
418 // Method to convert AliTRDdigits created by AliTRDv1::Hits2Digits()
419 // into AliTRDclusters
420 // To produce cluster from a root-file with TRD-digits use the
421 // slowClusterCreate.C macro.
426 printf("AliTRDv1::Digits2Clusters -- Start creating clusters\n");
428 AliTRDdigit *TRDdigit;
429 TClonesArray *TRDDigits;
432 Float_t maxThresh = fClusMaxThresh; // threshold value for maximum
433 Float_t signalThresh = fClusSigThresh; // threshold value for digit signal
434 Int_t clusteringMethod = fClusMethod; // clustering method option (for testing)
436 const Float_t epsilon = 0.01; // iteration limit for unfolding procedure
438 // Get the pointer to the digits tree
439 TTree *DigitTree = gAlice->TreeD();
440 // Get the pointer to the cluster tree
441 TTree *ClusterTree = gAlice->TreeD();
443 // Get the pointer to the digits container
444 TRDDigits = Digits();
447 Int_t chamEnd = kNcham;
448 if (fSensChamber) chamEnd = chamBeg = fSensChamber;
450 Int_t planEnd = kNplan;
451 if (fSensPlane) planEnd = planBeg = fSensPlane;
453 Int_t sectEnd = kNsect;
454 if (fSensSector) sectEnd = sectBeg = fSensSector;
456 // Import the digit tree
457 gAlice->ResetDigits();
459 nbytes += DigitTree->GetEvent(1);
461 // Get the number of digits in the detector
462 Int_t nTRDDigits = TRDDigits->GetEntriesFast();
464 // *** Start clustering *** in every chamber
465 for (Int_t icham = chamBeg; icham < chamEnd; icham++) {
466 for (Int_t iplan = planBeg; iplan < planEnd; iplan++) {
467 for (Int_t isect = sectBeg; isect < sectEnd; isect++) {
470 printf("AliTRDv1::Digits2Clusters -- Finding clusters in chamber %d, plane %d, sector %d\n"
471 ,icham+1,iplan+1,isect+1);
473 // Create a detector matrix to keep maxima
474 AliTRDmatrix *digitMatrix = new AliTRDmatrix(fRowMax[iplan][icham][isect]
478 // Create a matrix to contain maximum flags
479 AliTRDmatrix *maximaMatrix = new AliTRDmatrix(fRowMax[iplan][icham][isect]
482 ,isect+1,icham+1,iplan+1);
484 // Loop through all TRD digits
485 for (Int_t iTRDDigits = 0; iTRDDigits < nTRDDigits; iTRDDigits++) {
487 // Get the information for this digit
488 TRDdigit = (AliTRDdigit*) TRDDigits->UncheckedAt(iTRDDigits);
489 Int_t signal = TRDdigit->fSignal;
490 Int_t sector = TRDdigit->fSector;
491 Int_t chamber = TRDdigit->fChamber;
492 Int_t plane = TRDdigit->fPlane;
493 Int_t row = TRDdigit->fRow;
494 Int_t col = TRDdigit->fCol;
495 Int_t time = TRDdigit->fTime;
498 for (Int_t iTrack = 0; iTrack < 3; iTrack++) {
499 track[iTrack] = TRDdigit->AliDigit::fTracks[iTrack];
502 if ((sector != isect+1) ||
503 (plane != iplan+1) ||
504 (chamber != icham+1))
507 // Fill the detector matrix
508 if (signal > signalThresh) {
509 digitMatrix->SetSignal(row,col,time,signal);
510 for (Int_t iTrack = 0; iTrack < 3; iTrack++) {
511 if (track[iTrack] > 0) {
512 digitMatrix->AddTrack(row,col,time,track[iTrack]);
519 // Loop chamber and find maxima in digitMatrix
520 for (row = 0; row < fRowMax[iplan][icham][isect]; row++) {
521 for (Int_t col = 1; col < fColMax[iplan] ; col++) {
522 for (Int_t time = 0; time < fTimeMax ; time++) {
524 if (digitMatrix->GetSignal(row,col,time)
525 < digitMatrix->GetSignal(row,col - 1,time)) {
528 if (digitMatrix->GetSignal(row,col - 2,time)
529 < digitMatrix->GetSignal(row,col - 1,time)) {
530 // yes, so set maximum flag
531 maximaMatrix->SetSignal(row,col - 1,time,1);
533 else maximaMatrix->SetSignal(row,col - 1,time,0);
541 // now check maxima and calculate cluster position
542 for (row = 0; row < fRowMax[iplan][icham][isect]; row++) {
543 for (Int_t col = 1; col < fColMax[iplan] ; col++) {
544 for (Int_t time = 0; time < fTimeMax ; time++) {
546 if ((maximaMatrix->GetSignal(row,col,time) > 0)
547 && (digitMatrix->GetSignal(row,col,time) > maxThresh)) {
549 Int_t clusters[5] = {0}; // cluster-object data
551 Float_t ratio = 0; // ratio resulting from unfolding
552 Float_t padSignal[5] = {0}; // signals on max and neighbouring pads
553 Float_t clusterSignal[3] = {0}; // signals from cluster
554 Float_t clusterPos[3] = {0}; // cluster in ALICE refFrame coords
555 Float_t clusterPads[6] = {0}; // cluster pad info
558 clusters[0] = isect+1; // = isect ????
559 clusters[1] = icham+1; // = ichamber ????
560 clusters[2] = iplan+1; // = iplane ????
563 clusterPads[0] = icham+1;
564 clusterPads[1] = isect+1;
565 clusterPads[2] = iplan+1;
567 for (Int_t iPad = 0; iPad < 3; iPad++) {
568 clusterSignal[iPad] = digitMatrix->GetSignal(row,col-1+iPad,time);
571 // neighbouring maximum on right side?
572 if (col < fColMax[iplan] - 2) {
573 if (maximaMatrix->GetSignal(row,col + 2,time) > 0) {
574 for (Int_t iPad = 0; iPad < 5; iPad++) {
575 padSignal[iPad] = digitMatrix->GetSignal(row,col-1+iPad,time);
579 ratio = Unfold(epsilon, padSignal);
581 // set signal on overlapping pad to ratio
582 clusterSignal[2] *= ratio;
586 switch (clusteringMethod) {
588 // method 1: simply center of mass
589 clusterPads[3] = row + 0.5;
590 clusterPads[4] = col - 0.5 + (clusterSignal[2] - clusterSignal[0]) /
591 (clusterSignal[1] + clusterSignal[2] + clusterSignal[3]);
592 clusterPads[5] = time + 0.5;
597 // method 2: integral gauss fit on 3 pads
598 TH1F *hPadCharges = new TH1F("hPadCharges", "Charges on center 3 pads"
600 for (Int_t iCol = -1; iCol <= 3; iCol++) {
601 if (clusterSignal[iCol] < 1) clusterSignal[iCol] = 1;
602 hPadCharges->Fill(iCol, clusterSignal[iCol]);
604 hPadCharges->Fit("gaus", "IQ", "SAME", -0.5, 2.5);
605 TF1 *fPadChargeFit = hPadCharges->GetFunction("gaus");
606 Double_t colMean = fPadChargeFit->GetParameter(1);
608 clusterPads[3] = row + 0.5;
609 clusterPads[4] = col - 1.5 + colMean;
610 clusterPads[5] = time + 0.5;
618 Float_t clusterCharge = clusterSignal[0]
621 clusters[4] = (Int_t)clusterCharge;
624 for (Int_t iTrack = 0; iTrack < 3; iTrack++) {
625 trackSave[iTrack] = digitMatrix->GetTrack(row,col,time,iTrack);
628 // Calculate cluster position in ALICE refFrame coords
629 // and set array clusterPos to calculated values
630 Pads2XYZ(clusterPads, clusterPos);
632 // Add cluster to reconstruction tree
633 AddCluster(trackSave,clusters,clusterPos);
641 printf("AliTRDv1::Digits2Clusters -- Number of clusters found: %d\n",nClusters);
650 // Fill the cluster-tree
651 printf("AliTRDv1::Digits2Clusters -- Total number of clusters found: %d\n"
652 ,fClusters->GetEntries());
653 printf("AliTRDv1::Digits2Clusters -- Fill the cluster tree\n");
658 //_____________________________________________________________________________
659 void AliTRDv1::Init()
662 // Initialise Transition Radiation Detector after geometry has been built.
663 // Includes the default settings of all parameter for the digitization.
668 printf(" Slow simulator\n");
670 printf(" Only plane %d is sensitive\n",fSensPlane);
672 printf(" Only chamber %d is sensitive\n",fSensChamber);
674 printf(" Only sector %d is sensitive\n",fSensSector);
676 // First ionization potential (eV) for the gas mixture (90% Xe + 10% CO2)
677 const Float_t kPoti = 12.1;
678 // Maximum energy (50 keV);
679 const Float_t kEend = 50000.0;
680 // Ermilova distribution for the delta-ray spectrum
681 Float_t Poti = TMath::Log(kPoti);
682 Float_t Eend = TMath::Log(kEend);
683 fDeltaE = new TF1("deltae",Ermilova,Poti,Eend,0);
685 // Identifier of the sensitive volume (drift region)
686 fIdSens = gMC->VolId("UL05");
688 // Identifier of the TRD-driftchambers
689 fIdChamber1 = gMC->VolId("UCIO");
690 fIdChamber2 = gMC->VolId("UCIM");
691 fIdChamber3 = gMC->VolId("UCII");
693 // The default parameter for the digitization
694 if (!(fGasGain)) fGasGain = 2.0E3;
695 if (!(fNoise)) fNoise = 3000.;
696 if (!(fChipGain)) fChipGain = 10.;
697 if (!(fADCoutRange)) fADCoutRange = 255.;
698 if (!(fADCinRange)) fADCinRange = 2000.;
699 if (!(fADCthreshold)) fADCthreshold = 1;
701 // Transverse and longitudinal diffusion coefficients (Xe/Isobutane)
702 if (!(fDiffusionT)) fDiffusionT = 0.060;
703 if (!(fDiffusionL)) fDiffusionL = 0.017;
705 // The default parameter for the clustering
706 if (!(fClusMaxThresh)) fClusMaxThresh = 5.0;
707 if (!(fClusSigThresh)) fClusSigThresh = 2.0;
708 if (!(fClusMethod)) fClusMethod = 1;
710 for (Int_t i = 0; i < 80; i++) printf("*");
715 //_____________________________________________________________________________
716 Float_t AliTRDv1::PadResponse(Float_t x)
719 // The pad response for the chevron pads.
720 // We use a simple Gaussian approximation which should be good
721 // enough for our purpose.
724 // The parameters for the response function
725 const Float_t aa = 0.8872;
726 const Float_t bb = -0.00573;
727 const Float_t cc = 0.454;
728 const Float_t cc2 = cc*cc;
730 Float_t pr = aa * (bb + TMath::Exp(-x*x / (2. * cc2)));
736 //_____________________________________________________________________________
737 void AliTRDv1::SetSensPlane(Int_t iplane)
740 // Defines the hit-sensitive plane (1-6)
743 if ((iplane < 0) || (iplane > 6)) {
744 printf("Wrong input value: %d\n",iplane);
745 printf("Use standard setting\n");
756 //_____________________________________________________________________________
757 void AliTRDv1::SetSensChamber(Int_t ichamber)
760 // Defines the hit-sensitive chamber (1-5)
763 if ((ichamber < 0) || (ichamber > 5)) {
764 printf("Wrong input value: %d\n",ichamber);
765 printf("Use standard setting\n");
772 fSensChamber = ichamber;
776 //_____________________________________________________________________________
777 void AliTRDv1::SetSensSector(Int_t isector)
780 // Defines the hit-sensitive sector (1-18)
783 if ((isector < 0) || (isector > 18)) {
784 printf("Wrong input value: %d\n",isector);
785 printf("Use standard setting\n");
792 fSensSector = isector;
796 //_____________________________________________________________________________
797 void AliTRDv1::StepManager()
800 // Called at every step in the Transition Radiation Detector version 2.
801 // Slow simulator. Every charged track produces electron cluster as hits
802 // along its path across the drift volume. The step size is set acording
803 // to Bethe-Bloch. The energy distribution of the delta electrons follows
804 // a spectrum taken from Ermilova et al.
807 Int_t iIdSens, icSens;
808 Int_t iIdSpace, icSpace;
809 Int_t iIdChamber, icChamber;
821 Double_t betaGamma, pp;
823 TLorentzVector pos, mom;
824 TClonesArray &lhits = *fHits;
826 const Double_t kBig = 1.0E+12;
829 const Float_t kWion = 22.04;
830 // Maximum energy for e+ e- g for the step-size calculation
831 const Float_t kPTotMax = 0.002;
832 // Plateau value of the energy-loss for electron in xenon
833 // taken from: Allison + Comb, Ann. Rev. Nucl. Sci. (1980), 30, 253
834 //const Double_t kPlateau = 1.70;
835 // the averaged value (26/3/99)
836 const Float_t kPlateau = 1.55;
837 // dN1/dx|min for the gas mixture (90% Xe + 10% CO2)
838 const Float_t kPrim = 48.0;
839 // First ionization potential (eV) for the gas mixture (90% Xe + 10% CO2)
840 const Float_t kPoti = 12.1;
842 // Set the maximum step size to a very large number for all
843 // neutral particles and those outside the driftvolume
844 gMC->SetMaxStep(kBig);
846 // Use only charged tracks
847 if (( gMC->TrackCharge() ) &&
848 (!gMC->IsTrackStop() ) &&
849 (!gMC->IsTrackDisappeared())) {
851 // Inside a sensitive volume?
852 iIdSens = gMC->CurrentVolID(icSens);
853 if (iIdSens == fIdSens) {
855 iIdSpace = gMC->CurrentVolOffID(4,icSpace );
856 iIdChamber = gMC->CurrentVolOffID(1,icChamber);
858 // Calculate the energy of the delta-electrons
859 eDelta = TMath::Exp(fDeltaE->GetRandom()) - kPoti;
860 eDelta = TMath::Max(eDelta,0.0);
862 // The number of secondary electrons created
863 qTot = (Double_t) ((Int_t) (eDelta / kWion) + 1);
865 // The hit coordinates and charge
866 gMC->TrackPosition(pos);
873 Float_t phi = pos[1] != 0 ? TMath::Atan2(pos[0],pos[1]) : (pos[0] > 0 ? 180. : 0.);
874 vol[0] = ((Int_t) (phi / 20)) + 1;
876 // The chamber number
882 if (iIdChamber == fIdChamber1)
883 vol[1] = (hits[2] < 0 ? 1 : 5);
884 else if (iIdChamber == fIdChamber2)
885 vol[1] = (hits[2] < 0 ? 2 : 4);
886 else if (iIdChamber == fIdChamber3)
890 vol[2] = icChamber - TMath::Nint((Float_t) (icChamber / 7)) * 6;
892 // Check on selected volumes
893 Int_t addthishit = 1;
895 if ((fSensPlane) && (vol[2] != fSensPlane )) addthishit = 0;
896 if ((fSensChamber) && (vol[1] != fSensChamber)) addthishit = 0;
897 if ((fSensSector) && (vol[0] != fSensSector )) addthishit = 0;
903 new(lhits[fNhits++]) AliTRDhit(fIshunt,gAlice->CurrentTrack(),vol,hits);
905 // The energy loss according to Bethe Bloch
906 gMC->TrackMomentum(mom);
908 iPid = gMC->TrackPid();
910 ((iPid <= 3) && (pTot < kPTotMax))) {
911 aMass = gMC->TrackMass();
912 betaGamma = pTot / aMass;
913 pp = kPrim * BetheBloch(betaGamma);
914 // Take charge > 1 into account
915 charge = gMC->TrackCharge();
916 if (TMath::Abs(charge) > 1) pp = pp * charge*charge;
918 // Electrons above 20 Mev/c are at the plateau
920 pp = kPrim * kPlateau;
923 // Calculate the maximum step size for the next tracking step
927 while ((random[0] == 1.) || (random[0] == 0.));
928 gMC->SetMaxStep( - TMath::Log(random[0]) / pp);
933 // set step size to maximal value
934 gMC->SetMaxStep(kBig);
943 //_____________________________________________________________________________
944 Double_t AliTRDv1::BetheBloch(Double_t bg)
947 // Parametrization of the Bethe-Bloch-curve
948 // The parametrization is the same as for the TPC and is taken from Lehrhaus.
951 // This parameters have been adjusted to averaged values from GEANT
952 const Double_t kP1 = 7.17960e-02;
953 const Double_t kP2 = 8.54196;
954 const Double_t kP3 = 1.38065e-06;
955 const Double_t kP4 = 5.30972;
956 const Double_t kP5 = 2.83798;
958 // This parameters have been adjusted to Xe-data found in:
959 // Allison & Cobb, Ann. Rev. Nucl. Sci. (1980), 30, 253
960 //const Double_t kP1 = 0.76176E-1;
961 //const Double_t kP2 = 10.632;
962 //const Double_t kP3 = 3.17983E-6;
963 //const Double_t kP4 = 1.8631;
964 //const Double_t kP5 = 1.9479;
967 Double_t yy = bg / TMath::Sqrt(1. + bg*bg);
968 Double_t aa = TMath::Power(yy,kP4);
969 Double_t bb = TMath::Power((1./bg),kP5);
970 bb = TMath::Log(kP3 + bb);
971 return ((kP2 - aa - bb)*kP1 / aa);
978 //_____________________________________________________________________________
979 Double_t Ermilova(Double_t *x, Double_t *)
982 // Calculates the delta-ray energy distribution according to Ermilova.
983 // Logarithmic scale !
994 Float_t vxe[nV] = { 2.3026, 2.9957, 3.4012, 3.6889, 3.9120
995 , 4.0943, 4.2485, 4.3820, 4.4998, 4.6052
996 , 4.7005, 5.0752, 5.2983, 5.7038, 5.9915
997 , 6.2146, 6.5221, 6.9078, 7.3132, 7.6009
998 , 8.0064, 8.5172, 8.6995, 8.9872, 9.2103
999 , 9.4727, 9.9035,10.3735,10.5966,10.8198
1002 Float_t vye[nV] = { 80.0 , 31.0 , 23.3 , 21.1 , 21.0
1003 , 20.9 , 20.8 , 20.0 , 16.0 , 11.0
1004 , 8.0 , 6.0 , 5.2 , 4.6 , 4.0
1005 , 3.5 , 3.0 , 1.4 , 0.67 , 0.44
1006 , 0.3 , 0.18 , 0.12 , 0.08 , 0.056
1007 , 0.04 , 0.023, 0.015, 0.011, 0.01
1012 // Find the position
1016 dpos = energy - vxe[pos2++];
1020 if (pos2 > nV) pos2 = nV;
1023 // Differentiate between the sampling points
1024 dnde = (vye[pos1] - vye[pos2]) / (vxe[pos2] - vxe[pos1]);
1030 //_____________________________________________________________________________
1031 void AliTRDv1::Pads2XYZ(Float_t *pads, Float_t *pos)
1033 // Method to convert pad coordinates (row,col,time)
1034 // into ALICE reference frame coordinates (x,y,z)
1036 Int_t chamber = (Int_t) pads[0]; // chamber info (1-5)
1037 Int_t sector = (Int_t) pads[1]; // sector info (1-18)
1038 Int_t plane = (Int_t) pads[2]; // plane info (1-6)
1040 Int_t icham = chamber - 1; // chamber info (0-4)
1041 Int_t isect = sector - 1; // sector info (0-17)
1042 Int_t iplan = plane - 1; // plane info (0-5)
1044 Float_t padRow = pads[3]; // Pad Row position
1045 Float_t padCol = pads[4]; // Pad Column position
1046 Float_t timeSlice = pads[5]; // Time "position"
1048 // calculate (x,y) position in rotated chamber
1049 Float_t yRot = fCol0[iplan] + padCol * fColPadSize;
1050 Float_t xRot = fTime0[iplan] + timeSlice * fTimeBinSize;
1051 // calculate z-position:
1052 Float_t z = fRow0[iplan][icham][isect] + padRow * fRowPadSize;
1055 rotate chamber back to original position
1056 1. mirror at y-axis, 2. rotate back to position (-phi)
1057 / cos(phi) -sin(phi) \ / -1 0 \ / -cos(phi) -sin(phi) \
1058 \ sin(phi) cos(phi) / * \ 0 1 / = \ -sin(phi) cos(phi) /
1060 //Float_t phi = 2*kPI / kNsect * ((Float_t) sector - 0.5);
1061 //Float_t x = -xRot * TMath::Cos(phi) - yRot * TMath::Sin(phi);
1062 //Float_t y = -xRot * TMath::Sin(phi) + yRot * TMath::Cos(phi);
1063 Float_t phi = 2*kPI / kNsect * ((Float_t) sector - 0.5);
1064 Float_t x = -xRot * TMath::Cos(phi) + yRot * TMath::Sin(phi);
1065 Float_t y = xRot * TMath::Sin(phi) + yRot * TMath::Cos(phi);
1074 //_____________________________________________________________________________
1075 Float_t AliTRDv1::Unfold(Float_t eps, Float_t* padSignal)
1077 // Method to unfold neighbouring maxima.
1078 // The charge ratio on the overlapping pad is calculated
1079 // until there is no more change within the range given by eps.
1080 // The resulting ratio is then returned to the calling method.
1082 Int_t itStep = 0; // count iteration steps
1084 Float_t ratio = 0.5; // start value for ratio
1085 Float_t prevRatio = 0; // store previous ratio
1087 Float_t newLeftSignal[3] = {0}; // array to store left cluster signal
1088 Float_t newRightSignal[3] = {0}; // array to store right cluster signal
1091 while ((TMath::Abs(prevRatio - ratio) > eps) && (itStep < 10)) {
1096 // cluster position according to charge ratio
1097 Float_t maxLeft = (ratio*padSignal[2] - padSignal[0]) /
1098 (padSignal[0] + padSignal[1] + ratio*padSignal[2]);
1099 Float_t maxRight = (padSignal[4] - (1-ratio)*padSignal[2]) /
1100 ((1-ratio)*padSignal[2] + padSignal[3] + padSignal[4]);
1102 // set cluster charge ratio
1103 Float_t ampLeft = padSignal[1];
1104 Float_t ampRight = padSignal[3];
1106 // apply pad response to parameters
1107 newLeftSignal[0] = ampLeft*PadResponse(-1 - maxLeft);
1108 newLeftSignal[1] = ampLeft*PadResponse( 0 - maxLeft);
1109 newLeftSignal[2] = ampLeft*PadResponse( 1 - maxLeft);
1111 newRightSignal[0] = ampRight*PadResponse(-1 - maxRight);
1112 newRightSignal[1] = ampRight*PadResponse( 0 - maxRight);
1113 newRightSignal[2] = ampRight*PadResponse( 1 - maxRight);
1115 // calculate new overlapping ratio
1116 ratio = newLeftSignal[2]/(newLeftSignal[2] + newRightSignal[0]);