]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDv1.cxx
Correct ansi scoping not accepted by HP compilers
[u/mrichter/AliRoot.git] / TRD / AliTRDv1.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /*
17 $Log$
18 Revision 1.13  1999/11/02 17:14:51  fca
19 Correct ansi scoping not accepted by HP compilers
20
21 Revision 1.12  1999/11/02 16:35:56  fca
22 New version of TRD introduced
23
24 Revision 1.11  1999/11/01 20:41:51  fca
25 Added protections against using the wrong version of FRAME
26
27 Revision 1.10  1999/09/29 09:24:35  fca
28 Introduction of the Copyright and cvs Log
29
30 */
31
32 ///////////////////////////////////////////////////////////////////////////////
33 //                                                                           //
34 //  Transition Radiation Detector version 2 -- slow simulator                //
35 //                                                                           //
36 //Begin_Html
37 /*
38 <img src="picts/AliTRDfullClass.gif">
39 */
40 //End_Html
41 //                                                                           //
42 //                                                                           //
43 ///////////////////////////////////////////////////////////////////////////////
44
45 #include <TMath.h>
46 #include <TVector.h>
47 #include <TRandom.h>
48
49 #include "AliTRDv1.h"
50 #include "AliTRDmatrix.h"
51 #include "AliRun.h"
52 #include "AliMC.h"
53 #include "AliConst.h"
54
55 ClassImp(AliTRDv1)
56
57 //_____________________________________________________________________________
58 AliTRDv1::AliTRDv1(const char *name, const char *title) 
59          :AliTRD(name, title) 
60 {
61   //
62   // Standard constructor for Transition Radiation Detector version 2
63   //
64
65   fIdSens        = 0;
66
67   fIdChamber1    = 0;
68   fIdChamber2    = 0;
69   fIdChamber3    = 0;
70
71   fSensSelect    = 0;
72   fSensPlane     = 0;
73   fSensChamber   = 0;
74   fSensSector    = 0;
75
76   fGasGain       = 0;
77   fNoise         = 0;
78   fChipGain      = 0;
79   fADCoutRange   = 0;
80   fADCinRange    = 0;
81   fADCthreshold  = 0;
82
83   fDiffusionT    = 0;
84   fDiffusionL    = 0;
85
86   fClusMaxThresh = 0;
87   fClusSigThresh = 0;
88   fClusMethod    = 0;
89
90   fDeltaE        = NULL;
91
92   SetBufferSize(128000);
93
94 }
95
96 //_____________________________________________________________________________
97 AliTRDv1::~AliTRDv1()
98 {
99
100   if (fDeltaE) delete fDeltaE;
101
102 }
103  
104 //_____________________________________________________________________________
105 void AliTRDv1::CreateGeometry()
106 {
107   //
108   // Create the GEANT geometry for the Transition Radiation Detector - Version 2
109   // This version covers the full azimuth. 
110   //
111
112   // Check that FRAME is there otherwise we have no place where to put the TRD
113   AliModule* FRAME = gAlice->GetModule("FRAME");
114   if (!FRAME) return;
115
116   // Define the chambers
117   AliTRD::CreateGeometry();
118
119 }
120
121 //_____________________________________________________________________________
122 void AliTRDv1::CreateMaterials()
123 {
124   //
125   // Create materials for the Transition Radiation Detector version 2
126   //
127
128   AliTRD::CreateMaterials();
129
130 }
131
132 //_____________________________________________________________________________
133 void AliTRDv1::Diffusion(Float_t driftlength, Float_t *xyz)
134 {
135   //
136   // Applies the diffusion smearing to the position of a single electron
137   //
138
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);
147   }
148   else {
149     xyz[0] = 0.0;
150     xyz[1] = 0.0;
151     xyz[2] = 0.0;
152   }
153
154 }
155
156 //_____________________________________________________________________________
157 void AliTRDv1::Hits2Digits()
158 {
159   //
160   // Creates TRD digits from hits. This procedure includes the following:
161   //      - Diffusion
162   //      - Gas gain including fluctuations
163   //      - Pad-response (simple Gaussian approximation)
164   //      - Electronics noise
165   //      - Electronics gain
166   //      - Digitization
167   //      - ADC threshold
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
170   // Init-function).
171   // To produce digits from a root-file with TRD-hits use the
172   // slowDigitsCreate.C macro.
173   //
174
175   printf("AliTRDv1::Hits2Digits -- Start creating digits\n");
176
177   ///////////////////////////////////////////////////////////////
178   // Parameter 
179   ///////////////////////////////////////////////////////////////
180
181   // Converts number of electrons to fC
182   const Float_t el2fC = 1.602E-19 * 1.0E15; 
183
184   ///////////////////////////////////////////////////////////////
185
186   Int_t nBytes = 0;
187
188   Int_t iRow;
189
190   AliTRDhit *TRDhit;
191
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();
196
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();
200
201   Int_t chamBeg = 0;
202   Int_t chamEnd = kNcham;
203   if (fSensChamber) chamEnd = chamBeg = fSensChamber;
204   Int_t planBeg = 0;
205   Int_t planEnd = kNplan;
206   if (fSensPlane)   planEnd = planBeg = fSensPlane;
207   Int_t sectBeg = 0;
208   Int_t sectEnd = kNsect;
209   if (fSensSector)  sectEnd = sectBeg = fSensSector;
210
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++) {
215
216         Int_t nDigits = 0;
217
218         printf("AliTRDv1::Hits2Digits -- Digitizing chamber %d, plane %d, sector %d\n"
219               ,icham+1,iplan+1,isect+1);
220
221         // Create a detector matrix to keep the signal and track numbers
222         AliTRDmatrix *matrix = new AliTRDmatrix(fRowMax[iplan][icham][isect]
223                                                ,fColMax[iplan]
224                                                ,fTimeMax
225                                                ,isect+1,icham+1,iplan+1);
226
227         // Loop through all entries in the tree
228         for (Int_t iTrack = 0; iTrack < nTrack; iTrack++) {
229
230           gAlice->ResetHits();
231           nBytes += HitTree->GetEvent(iTrack);
232
233           // Get the number of hits in the TRD created by this particle
234           Int_t nHit = fHits->GetEntriesFast();
235
236           // Loop through the TRD hits  
237           for (Int_t iHit = 0; iHit < nHit; iHit++) {
238
239             if (!(TRDhit = (AliTRDhit *) fHits->UncheckedAt(iHit))) 
240               continue;
241
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;        
250
251             if ((sector  != isect+1) ||
252                 (plane   != iplan+1) ||
253                 (chamber != icham+1)) 
254               continue;
255
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);
261             Float_t zRot =  z;
262
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);
268             // The time bucket
269             Int_t timeH = (Int_t) ((xRot - fTime0[iplan]              ) / fTimeBinSize);
270
271             // Array to sum up the signal in a box surrounding the
272             // hit postition
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;
281                 }
282               }
283             }
284
285             // Loop over all electrons of this hit
286             Int_t nEl = (Int_t) q;
287             for (Int_t iEl = 0; iEl < nEl; iEl++) {
288
289               // Apply the diffusion smearing
290               Float_t driftlength = xRot - fTime0[iplan];
291               Float_t xyz[3];
292               xyz[0] = xRot;
293               xyz[1] = yRot;
294               xyz[2] = zRot;
295               Diffusion(driftlength,xyz);
296
297               // At this point absorption effects that depend on the 
298               // driftlength could be taken into account.              
299
300               // The electron position and the distance to the hit position
301               // in pad units
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;
308               // The time bucket
309               Int_t timeE = (Int_t) ((xyz[0] - fTime0[iplan]              ) / fTimeBinSize);
310               Int_t timeD = timeH - timeE;
311
312               // Apply the gas gain including fluctuations
313               Int_t signal = (Int_t) (-fGasGain * TMath::Log(gRandom->Rndm()));
314
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) 
318                            / fColPadSize;
319
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;
328
329             }
330
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++) {
335
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);
339
340                   Float_t signalB = signalSum[iRow][iCol][iTime];
341                   if (signalB > 0.0) {
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");
345                   }
346
347                 }
348               }
349             }
350
351           }
352
353         }
354
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++) {         
359
360               Float_t signalAmp = matrix->GetSignal(iRow,iCol,iTime);
361
362               // Add the noise
363               signalAmp  = TMath::Max(gRandom->Gaus(signalAmp,fNoise),(Float_t) 0.0);
364               // Convert to fC
365               signalAmp *= el2fC;
366               // Convert to mV
367               signalAmp *= fChipGain;
368               // Convert to ADC counts
369               Int_t adc  = (Int_t) (signalAmp * (fADCoutRange / fADCinRange));
370
371               // Apply threshold on ADC value
372               if (adc > fADCthreshold) {
373
374                 Int_t trackSave[3];
375                 for (Int_t ii = 0; ii < 3; ii++) {
376                   trackSave[ii] = matrix->GetTrack(iRow,iCol,iTime,ii);
377                 }
378
379                 Int_t digits[7];
380                 digits[0] = matrix->GetSector();
381                 digits[1] = matrix->GetChamber();
382                 digits[2] = matrix->GetPlane();
383                 digits[3] = iRow;
384                 digits[4] = iCol;
385                 digits[5] = iTime;
386                 digits[6] = adc;
387
388                 // Add this digit to the TClonesArray
389                 AddDigit(trackSave,digits);
390                 nDigits++;
391
392               }
393
394             }
395           }
396         }
397
398         printf("AliTRDv1::Hits2Digits -- Number of digits found: %d\n",nDigits);
399
400         // Clean up
401         delete matrix;
402
403       }
404     }
405   }
406
407   // Fill the digits-tree
408   printf("AliTRDv1::Hits2Digits -- Fill the digits tree\n");
409   DigitsTree->Fill();
410
411 }
412
413 //_____________________________________________________________________________
414 void AliTRDv1::Digits2Clusters()
415 {
416
417   //
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.
422   //
423   
424   Int_t row;
425
426   printf("AliTRDv1::Digits2Clusters -- Start creating clusters\n");
427
428   AliTRDdigit  *TRDdigit;
429   TClonesArray *TRDDigits;
430
431   // Parameters
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)
435
436   const Float_t epsilon    = 0.01;             // iteration limit for unfolding procedure
437
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();
442
443   // Get the pointer to the digits container
444   TRDDigits = Digits();
445
446   Int_t chamBeg = 0;
447   Int_t chamEnd = kNcham;
448   if (fSensChamber) chamEnd = chamBeg = fSensChamber;
449   Int_t planBeg = 0;
450   Int_t planEnd = kNplan;
451   if (fSensPlane)   planEnd = planBeg = fSensPlane;
452   Int_t sectBeg = 0;
453   Int_t sectEnd = kNsect;
454   if (fSensSector)  sectEnd = sectBeg = fSensSector;
455
456   // Import the digit tree
457   gAlice->ResetDigits();
458   Int_t nbytes;
459   nbytes += DigitTree->GetEvent(1);
460
461   // Get the number of digits in the detector
462   Int_t nTRDDigits = TRDDigits->GetEntriesFast();
463
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++) {
468
469         Int_t nClusters = 0;
470         printf("AliTRDv1::Digits2Clusters -- Finding clusters in chamber %d, plane %d, sector %d\n"
471                ,icham+1,iplan+1,isect+1);
472
473         // Create a detector matrix to keep maxima
474         AliTRDmatrix *digitMatrix = new AliTRDmatrix(fRowMax[iplan][icham][isect]
475                                                      ,fColMax[iplan]
476                                                      ,fTimeMax,isect+1
477                                                      ,icham+1,iplan+1);
478         // Create a matrix to contain maximum flags
479         AliTRDmatrix *maximaMatrix = new AliTRDmatrix(fRowMax[iplan][icham][isect]
480                                                       ,fColMax[iplan]
481                                                       ,fTimeMax
482                                                       ,isect+1,icham+1,iplan+1);
483
484         // Loop through all TRD digits
485         for (Int_t iTRDDigits = 0; iTRDDigits < nTRDDigits; iTRDDigits++) {
486
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;
496
497           Int_t   track[3];
498           for (Int_t iTrack = 0; iTrack < 3; iTrack++) {
499             track[iTrack]  = TRDdigit->AliDigit::fTracks[iTrack];
500           }
501
502           if ((sector  != isect+1) ||
503               (plane   != iplan+1) ||
504               (chamber != icham+1))
505             continue;
506
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]);
513               }
514             }
515           }
516
517         }
518
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++) {
523
524               if (digitMatrix->GetSignal(row,col,time) 
525                   < digitMatrix->GetSignal(row,col - 1,time)) {
526                 // really maximum?
527                 if (col > 1) {
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);
532                   }
533                   else maximaMatrix->SetSignal(row,col - 1,time,0);
534                 }
535               }
536
537             }   // time
538           }     // col
539         }       // row
540
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++) {
545
546               if ((maximaMatrix->GetSignal(row,col,time) > 0)
547                   && (digitMatrix->GetSignal(row,col,time) > maxThresh)) {
548
549                 Int_t   clusters[5]        = {0};   // cluster-object data
550
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
556
557                 // setting values
558                 clusters[0] = isect+1;    // = isect ????
559                 clusters[1] = icham+1;    // = ichamber ????
560                 clusters[2] = iplan+1;    // = iplane ????
561                 clusters[3] = time;
562
563                 clusterPads[0] = icham+1;
564                 clusterPads[1] = isect+1;
565                 clusterPads[2] = iplan+1;
566
567                 for (Int_t iPad = 0; iPad < 3; iPad++) {
568                   clusterSignal[iPad] = digitMatrix->GetSignal(row,col-1+iPad,time);
569                 }
570
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);
576                     }
577
578                     // unfold:
579                     ratio = Unfold(epsilon, padSignal);
580
581                     // set signal on overlapping pad to ratio
582                     clusterSignal[2] *= ratio;
583                   }
584                 }
585
586                 switch (clusteringMethod) {
587                 case 1:
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;
593
594                   nClusters++;
595                   break;
596                 case 2:
597                   // method 2: integral gauss fit on 3 pads
598                   TH1F *hPadCharges = new TH1F("hPadCharges", "Charges on center 3 pads"
599                                                             , 5, -1.5, 3.5);
600                   for (Int_t iCol = -1; iCol <= 3; iCol++) {
601                     if (clusterSignal[iCol] < 1) clusterSignal[iCol] = 1;
602                     hPadCharges->Fill(iCol, clusterSignal[iCol]);
603                   }
604                   hPadCharges->Fit("gaus", "IQ", "SAME", -0.5, 2.5);
605                   TF1     *fPadChargeFit = hPadCharges->GetFunction("gaus");
606                   Double_t  colMean = fPadChargeFit->GetParameter(1);
607
608                   clusterPads[3] = row + 0.5;
609                   clusterPads[4] = col - 1.5 + colMean;
610                   clusterPads[5] = time + 0.5;
611
612                   delete hPadCharges;
613
614                   nClusters++;
615                   break;
616                 }
617
618                 Float_t clusterCharge =   clusterSignal[0]
619                                         + clusterSignal[1]
620                                         + clusterSignal[2];
621                 clusters[4] = (Int_t)clusterCharge;
622
623                 Int_t trackSave[3];
624                 for (Int_t iTrack = 0; iTrack < 3; iTrack++) {
625                   trackSave[iTrack] = digitMatrix->GetTrack(row,col,time,iTrack);
626                 }
627
628                 // Calculate cluster position in ALICE refFrame coords
629                 // and set array clusterPos to calculated values
630                 Pads2XYZ(clusterPads, clusterPos);
631
632                 // Add cluster to reconstruction tree
633                 AddCluster(trackSave,clusters,clusterPos);
634
635               }
636
637             }  // time
638           }    // col
639         }      // row
640
641         printf("AliTRDv1::Digits2Clusters -- Number of clusters found: %d\n",nClusters);
642
643         delete digitMatrix;
644         delete maximaMatrix;
645
646       }          // isect
647     }            // iplan
648   }              // icham
649
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");
654   ClusterTree->Fill();
655
656 }
657
658 //_____________________________________________________________________________
659 void AliTRDv1::Init() 
660 {
661   //
662   // Initialise Transition Radiation Detector after geometry has been built.
663   // Includes the default settings of all parameter for the digitization.
664   //
665
666   AliTRD::Init();
667
668   printf("          Slow simulator\n");
669   if (fSensPlane)
670     printf("          Only plane %d is sensitive\n",fSensPlane);
671   if (fSensChamber)   
672     printf("          Only chamber %d is sensitive\n",fSensChamber);
673   if (fSensSector)
674     printf("          Only sector %d is sensitive\n",fSensSector);
675
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);
684
685   // Identifier of the sensitive volume (drift region)
686   fIdSens     = gMC->VolId("UL05");
687
688   // Identifier of the TRD-driftchambers
689   fIdChamber1 = gMC->VolId("UCIO");
690   fIdChamber2 = gMC->VolId("UCIM");
691   fIdChamber3 = gMC->VolId("UCII");
692
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;
700
701   // Transverse and longitudinal diffusion coefficients (Xe/Isobutane)
702   if (!(fDiffusionT))    fDiffusionT    = 0.060;
703   if (!(fDiffusionL))    fDiffusionL    = 0.017;
704
705   // The default parameter for the clustering
706   if (!(fClusMaxThresh)) fClusMaxThresh = 5.0;
707   if (!(fClusSigThresh)) fClusSigThresh = 2.0;
708   if (!(fClusMethod))    fClusMethod    = 1;
709
710   for (Int_t i = 0; i < 80; i++) printf("*");
711   printf("\n");
712
713 }
714
715 //_____________________________________________________________________________
716 Float_t AliTRDv1::PadResponse(Float_t x)
717 {
718   //
719   // The pad response for the chevron pads. 
720   // We use a simple Gaussian approximation which should be good
721   // enough for our purpose.
722   //
723
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;
729
730   Float_t pr = aa * (bb + TMath::Exp(-x*x / (2. * cc2)));
731
732   return (pr);
733
734 }
735
736 //_____________________________________________________________________________
737 void AliTRDv1::SetSensPlane(Int_t iplane)
738 {
739   //
740   // Defines the hit-sensitive plane (1-6)
741   //
742
743   if ((iplane < 0) || (iplane > 6)) {
744     printf("Wrong input value: %d\n",iplane);
745     printf("Use standard setting\n");
746     fSensPlane  = 0;
747     fSensSelect = 0;
748     return;
749   }
750
751   fSensSelect = 1;
752   fSensPlane  = iplane;
753
754 }
755
756 //_____________________________________________________________________________
757 void AliTRDv1::SetSensChamber(Int_t ichamber)
758 {
759   //
760   // Defines the hit-sensitive chamber (1-5)
761   //
762
763   if ((ichamber < 0) || (ichamber > 5)) {
764     printf("Wrong input value: %d\n",ichamber);
765     printf("Use standard setting\n");
766     fSensChamber = 0;
767     fSensSelect  = 0;
768     return;
769   }
770
771   fSensSelect  = 1;
772   fSensChamber = ichamber;
773
774 }
775
776 //_____________________________________________________________________________
777 void AliTRDv1::SetSensSector(Int_t isector)
778 {
779   //
780   // Defines the hit-sensitive sector (1-18)
781   //
782
783   if ((isector < 0) || (isector > 18)) {
784     printf("Wrong input value: %d\n",isector);
785     printf("Use standard setting\n");
786     fSensSector = 0;
787     fSensSelect = 0;
788     return;
789   }
790
791   fSensSelect = 1;
792   fSensSector = isector;
793
794 }
795
796 //_____________________________________________________________________________
797 void AliTRDv1::StepManager()
798 {
799   //
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.
805   //
806
807   Int_t    iIdSens, icSens;
808   Int_t    iIdSpace, icSpace;
809   Int_t    iIdChamber, icChamber;
810   Int_t    vol[3]; 
811   Int_t    iPid;
812
813   Float_t  hits[4];
814   Float_t  random[1];
815   Float_t  charge;
816   Float_t  aMass;
817
818   Double_t pTot;
819   Double_t qTot;
820   Double_t eDelta;
821   Double_t betaGamma, pp;
822
823   TLorentzVector pos, mom;
824   TClonesArray  &lhits = *fHits;
825
826   const Double_t kBig = 1.0E+12;
827
828   // Ionization energy
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;
841
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); 
845
846   // Use only charged tracks 
847   if (( gMC->TrackCharge()       ) &&
848       (!gMC->IsTrackStop()       ) && 
849       (!gMC->IsTrackDisappeared())) {
850
851     // Inside a sensitive volume?
852     iIdSens = gMC->CurrentVolID(icSens);
853     if (iIdSens == fIdSens) { 
854
855       iIdSpace   = gMC->CurrentVolOffID(4,icSpace  );
856       iIdChamber = gMC->CurrentVolOffID(1,icChamber);
857
858       // Calculate the energy of the delta-electrons
859       eDelta = TMath::Exp(fDeltaE->GetRandom()) - kPoti;
860       eDelta = TMath::Max(eDelta,0.0);
861
862       // The number of secondary electrons created
863       qTot = (Double_t) ((Int_t) (eDelta / kWion) + 1);
864
865       // The hit coordinates and charge
866       gMC->TrackPosition(pos);
867       hits[0] = pos[0];
868       hits[1] = pos[1];
869       hits[2] = pos[2];
870       hits[3] = qTot;
871
872       // The sector number
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;
875
876       // The chamber number 
877       //   1: outer left
878       //   2: middle left
879       //   3: inner
880       //   4: middle right
881       //   5: outer right
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)       
887         vol[1] = 3;
888
889       // The plane number
890       vol[2] = icChamber - TMath::Nint((Float_t) (icChamber / 7)) * 6;
891
892       // Check on selected volumes
893       Int_t addthishit = 1;
894       if (fSensSelect) {
895         if ((fSensPlane)   && (vol[2] != fSensPlane  )) addthishit = 0;
896         if ((fSensChamber) && (vol[1] != fSensChamber)) addthishit = 0;
897         if ((fSensSector)  && (vol[0] != fSensSector )) addthishit = 0;
898       }
899
900       // Add this hit
901       if (addthishit) {
902
903         new(lhits[fNhits++]) AliTRDhit(fIshunt,gAlice->CurrentTrack(),vol,hits);
904
905         // The energy loss according to Bethe Bloch
906         gMC->TrackMomentum(mom);
907         pTot = mom.Rho();
908         iPid = gMC->TrackPid();
909         if ( (iPid >  3) ||
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;
917         }
918         // Electrons above 20 Mev/c are at the plateau
919         else {
920           pp = kPrim * kPlateau;
921         }
922       
923         // Calculate the maximum step size for the next tracking step
924         if (pp > 0) {
925           do 
926             gMC->Rndm(random,1);
927           while ((random[0] == 1.) || (random[0] == 0.));
928           gMC->SetMaxStep( - TMath::Log(random[0]) / pp);
929         }
930
931       }
932       else {
933         // set step size to maximal value
934         gMC->SetMaxStep(kBig); 
935       }
936
937     }
938
939   }
940
941 }
942
943 //_____________________________________________________________________________
944 Double_t AliTRDv1::BetheBloch(Double_t bg) 
945 {
946   //
947   // Parametrization of the Bethe-Bloch-curve
948   // The parametrization is the same as for the TPC and is taken from Lehrhaus.
949   //
950
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;
957
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;
965
966   if (bg > 0) {
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);
972   }
973   else
974     return 0;
975
976 }
977
978 //_____________________________________________________________________________
979 Double_t Ermilova(Double_t *x, Double_t *)
980 {
981   //
982   // Calculates the delta-ray energy distribution according to Ermilova.
983   // Logarithmic scale !
984   //
985
986   Double_t energy;
987   Double_t dpos;
988   Double_t dnde;
989
990   Int_t    pos1, pos2;
991
992   const Int_t nV = 31;
993
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
1000                     ,11.5129 };
1001
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
1008                     ,  0.004 };
1009
1010   energy = x[0];
1011
1012   // Find the position 
1013   pos1 = pos2 = 0;
1014   dpos = 0;
1015   do {
1016     dpos = energy - vxe[pos2++];
1017   } 
1018   while (dpos > 0);
1019   pos2--; 
1020   if (pos2 > nV) pos2 = nV;
1021   pos1 = pos2 - 1;
1022
1023   // Differentiate between the sampling points
1024   dnde = (vye[pos1] - vye[pos2]) / (vxe[pos2] - vxe[pos1]);
1025
1026   return dnde;
1027
1028 }
1029
1030 //_____________________________________________________________________________
1031 void AliTRDv1::Pads2XYZ(Float_t *pads, Float_t *pos)
1032 {
1033   // Method to convert pad coordinates (row,col,time)
1034   // into ALICE reference frame coordinates (x,y,z)
1035
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)
1039
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)
1043
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"
1047
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;
1053
1054   /**
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) /
1059   **/
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);
1066
1067   // Setting values
1068   pos[0] = x;
1069   pos[1] = y;
1070   pos[2] = z;
1071
1072 }
1073
1074 //_____________________________________________________________________________
1075 Float_t AliTRDv1::Unfold(Float_t eps, Float_t* padSignal)
1076 {
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.
1081
1082   Int_t   itStep            = 0;      // count iteration steps
1083
1084   Float_t ratio             = 0.5;    // start value for ratio
1085   Float_t prevRatio         = 0;      // store previous ratio
1086
1087   Float_t newLeftSignal[3]  = {0};    // array to store left cluster signal
1088   Float_t newRightSignal[3] = {0};    // array to store right cluster signal
1089
1090   // start iteration:
1091   while ((TMath::Abs(prevRatio - ratio) > eps) && (itStep < 10)) {
1092
1093     itStep++;
1094     prevRatio = ratio;
1095
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]);
1101
1102     // set cluster charge ratio
1103     Float_t ampLeft  = padSignal[1];
1104     Float_t ampRight = padSignal[3];
1105
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);
1110
1111     newRightSignal[0] = ampRight*PadResponse(-1 - maxRight);
1112     newRightSignal[1] = ampRight*PadResponse( 0 - maxRight);
1113     newRightSignal[2] = ampRight*PadResponse( 1 - maxRight);
1114
1115     // calculate new overlapping ratio
1116     ratio = newLeftSignal[2]/(newLeftSignal[2] + newRightSignal[0]);
1117
1118   }
1119
1120   return ratio;
1121
1122 }
1123