]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDv1.cxx
Go back to previous version of simulation
[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 /* $Id$ */
17
18 ////////////////////////////////////////////////////////////////////////////
19 //                                                                        //
20 //  Transition Radiation Detector version 1 -- slow simulator             //
21 //                                                                        //
22 ////////////////////////////////////////////////////////////////////////////
23
24 #include <stdlib.h> 
25
26 #include <TF1.h>
27 #include <TLorentzVector.h>
28 #include <TMath.h>
29 #include <TRandom.h>
30 #include <TVector.h>
31 #include <TVirtualMC.h>
32 #include <TGeoManager.h>
33
34 #include "AliConst.h"
35 #include "AliLog.h"
36 #include "AliMC.h"
37 #include "AliRun.h"
38
39 #include "AliTRDgeometry.h"
40 #include "AliTRDhit.h"
41 #include "AliTRDsim.h"
42 #include "AliTRDv1.h"
43
44 ClassImp(AliTRDv1)
45  
46 //_____________________________________________________________________________
47 AliTRDv1::AliTRDv1()
48   :AliTRD()
49   ,fTRon(kFALSE)
50   ,fTR(NULL)
51   ,fTypeOfStepManager(0)
52   ,fStepSize(0)
53   ,fDeltaE(NULL)
54   ,fDeltaG(NULL)
55   ,fTrackLength0(0)
56   ,fPrimaryTrackPid(0)
57 {
58   //
59   // Default constructor
60   //
61
62 }
63
64 //_____________________________________________________________________________
65 AliTRDv1::AliTRDv1(const char *name, const char *title) 
66   :AliTRD(name,title) 
67   ,fTRon(kTRUE)
68   ,fTR(NULL)
69   ,fTypeOfStepManager(2)
70   ,fStepSize(0.1)
71   ,fDeltaE(NULL)
72   ,fDeltaG(NULL)
73   ,fTrackLength0(0)
74   ,fPrimaryTrackPid(0)
75 {
76   //
77   // Standard constructor for Transition Radiation Detector version 1
78   //
79
80   SetBufferSize(128000);
81
82 }
83
84 //_____________________________________________________________________________
85 AliTRDv1::~AliTRDv1()
86 {
87   //
88   // AliTRDv1 destructor
89   //
90
91   if (fDeltaE) {
92     delete fDeltaE;
93     fDeltaE = 0;
94   }
95
96   if (fDeltaG) {
97     delete fDeltaG;
98     fDeltaG = 0;
99   }
100
101   if (fTR) {
102     delete fTR;
103     fTR     = 0;
104   }
105
106 }
107  
108 //_____________________________________________________________________________
109 void AliTRDv1::AddAlignableVolumes() const
110 {
111   //
112   // Create entries for alignable volumes associating the symbolic volume
113   // name with the corresponding volume path. Needs to be syncronized with
114   // eventual changes in the geometry.
115   //
116
117   TString volPath;
118   TString symName;
119
120   TString vpStr  = "ALIC_1/B077_1/BSEGMO";
121   TString vpApp1 = "_1/BTRD";
122   TString vpApp2 = "_1";
123   TString vpApp3 = "/UTR1_1/UTS1_1/UTI1_1/UT";
124
125   TString snStr  = "TRD/sm";
126   TString snApp1 = "/st";
127   TString snApp2 = "/pl";
128
129   //
130   // The super modules
131   // The symbolic names are: TRD/sm00
132   //                           ...
133   //                         TRD/sm17
134   //
135   for (Int_t isect = 0; isect < AliTRDgeometry::Nsect(); isect++) {
136
137     volPath  = vpStr;
138     volPath += isect;
139     volPath += vpApp1;
140     volPath += isect;
141     volPath += vpApp2;
142
143     symName  = snStr;
144     symName += Form("%02d",isect);
145
146     gGeoManager->SetAlignableEntry(symName.Data(),volPath.Data());
147
148   }
149
150   //
151   // The readout chambers
152   // The symbolic names are: TRD/sm00/st0/pl0
153   //                           ...
154   //                         TRD/sm17/st4/pl5
155   //
156   for (Int_t isect = 0; isect < AliTRDgeometry::Nsect(); isect++) {
157     for (Int_t icham = 0; icham < AliTRDgeometry::Ncham(); icham++) {
158       for (Int_t iplan = 0; iplan < AliTRDgeometry::Nplan(); iplan++) {
159
160         Int_t idet = AliTRDgeometry::GetDetectorSec(iplan,icham);
161
162         volPath  = vpStr;
163         volPath += isect;
164         volPath += vpApp1;
165         volPath += isect;
166         volPath += vpApp2;
167         volPath += vpApp3;
168         volPath += Form("%02d",idet);
169         volPath += vpApp2;
170
171         symName  = snStr;
172         symName += Form("%02d",isect);
173         symName += snApp1;
174         symName += icham;
175         symName += snApp2;
176         symName += iplan;
177
178         gGeoManager->SetAlignableEntry(symName.Data(),volPath.Data());
179
180       }
181     }
182   }
183
184 }
185
186 //_____________________________________________________________________________
187 void AliTRDv1::CreateGeometry()
188 {
189   //
190   // Create the GEANT geometry for the Transition Radiation Detector - Version 1
191   // This version covers the full azimuth. 
192   //
193
194   // Check that FRAME is there otherwise we have no place where to put the TRD
195   AliModule* frame = gAlice->GetModule("FRAME");
196   if (!frame) {
197     AliError("TRD needs FRAME to be present\n");
198     return;
199   }
200
201   // Define the chambers
202   AliTRD::CreateGeometry();
203
204 }
205
206 //_____________________________________________________________________________
207 void AliTRDv1::CreateMaterials()
208 {
209   //
210   // Create materials for the Transition Radiation Detector version 1
211   //
212
213   AliTRD::CreateMaterials();
214
215 }
216
217 //_____________________________________________________________________________
218 void AliTRDv1::CreateTRhit(Int_t det)
219 {
220   //
221   // Creates an electron cluster from a TR photon.
222   // The photon is assumed to be created a the end of the radiator. The 
223   // distance after which it deposits its energy takes into account the 
224   // absorbtion of the entrance window and of the gas mixture in drift
225   // volume.
226   //
227
228   // Ionization energy
229   const Float_t kWion        = 23.53;
230
231   // Maximum number of TR photons per track
232   const Int_t   kNTR         = 50;
233
234   TLorentzVector mom;
235   TLorentzVector pos;
236
237   Float_t eTR[kNTR];
238   Int_t   nTR;
239
240   // Create TR photons
241   gMC->TrackMomentum(mom);
242   Float_t pTot = mom.Rho();
243   fTR->CreatePhotons(11,pTot,nTR,eTR);
244   if (nTR > kNTR) {
245     AliFatal(Form("Boundary error: nTR = %d, kNTR = %d",nTR,kNTR));
246   }
247
248   // Loop through the TR photons
249   for (Int_t iTR = 0; iTR < nTR; iTR++) {
250
251     Float_t energyMeV = eTR[iTR] * 0.001;
252     Float_t energyeV  = eTR[iTR] * 1000.0;
253     Float_t absLength = 0.0;
254     Float_t sigma     = 0.0;
255
256     // Take the absorbtion in the entrance window into account
257     Double_t muMy = fTR->GetMuMy(energyMeV);
258     sigma         = muMy * fFoilDensity;
259     if (sigma > 0.0) {
260       absLength = gRandom->Exp(1.0/sigma);
261       if (absLength < AliTRDgeometry::MyThick()) {
262         continue;
263       }
264     }
265     else {
266       continue;
267     }
268
269     // The absorbtion cross sections in the drift gas
270     // Gas-mixture (Xe/CO2)
271     Double_t muXe = fTR->GetMuXe(energyMeV);
272     Double_t muCO = fTR->GetMuCO(energyMeV);
273     sigma = (0.85 * muXe + 0.15 * muCO) * fGasDensity * fTR->GetTemp();
274
275     // The distance after which the energy of the TR photon
276     // is deposited.
277     if (sigma > 0.0) {
278       absLength = gRandom->Exp(1.0/sigma);
279       if (absLength > (AliTRDgeometry::DrThick()
280                      + AliTRDgeometry::AmThick())) {
281         continue;
282       }
283     }
284     else {
285       continue;
286     }
287
288     // The position of the absorbtion
289     Float_t posHit[3];
290     gMC->TrackPosition(pos);
291     posHit[0] = pos[0] + mom[0] / pTot * absLength;
292     posHit[1] = pos[1] + mom[1] / pTot * absLength;
293     posHit[2] = pos[2] + mom[2] / pTot * absLength;
294
295     // Create the charge 
296     Int_t q = ((Int_t) (energyeV / kWion));
297
298     // Add the hit to the array. TR photon hits are marked 
299     // by negative charge
300     AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber()
301           ,det
302           ,posHit
303           ,-q
304           ,gMC->TrackTime()*1.0e06
305           ,kTRUE);
306
307   }
308
309 }
310
311 //_____________________________________________________________________________
312 void AliTRDv1::Init() 
313 {
314   //
315   // Initialise Transition Radiation Detector after geometry has been built.
316   //
317
318   AliTRD::Init();
319
320   AliDebug(1,"Slow simulator\n");
321
322   // Switch on TR simulation as default
323   if (!fTRon) {
324     AliInfo("TR simulation off");
325   }
326   else {
327     fTR = new AliTRDsim();
328   }
329
330   // First ionization potential (eV) for the gas mixture (90% Xe + 10% CO2)
331   const Float_t kPoti = 12.1;
332   // Maximum energy (50 keV);
333   const Float_t kEend = 50000.0;
334   // Ermilova distribution for the delta-ray spectrum
335   Float_t poti        = TMath::Log(kPoti);
336   Float_t eEnd        = TMath::Log(kEend);
337
338   // Ermilova distribution for the delta-ray spectrum
339   fDeltaE = new TF1("deltae" ,Ermilova ,poti,eEnd,0);
340
341   // Geant3 distribution for the delta-ray spectrum
342   fDeltaG = new TF1("deltag",IntSpecGeant,2.421257,28.536469,0);
343
344   AliDebug(1,"+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++");
345
346 }
347
348 //_____________________________________________________________________________
349 void AliTRDv1::StepManager()
350 {
351   //
352   // Slow simulator. Every charged track produces electron cluster as hits
353   // along its path across the drift volume. 
354   //
355
356   switch (fTypeOfStepManager) {
357    case 0: 
358     StepManagerErmilova();
359     break;  
360    case 1: 
361     StepManagerGeant();  
362     break;  
363    case 2: 
364     StepManagerFixedStep();
365     break;  
366    default: 
367     AliWarning("Not a valid Step Manager.");
368   }
369
370 }
371
372 //_____________________________________________________________________________
373 void AliTRDv1::SelectStepManager(Int_t t)
374 {
375   //
376   // Selects a step manager type:
377   //   0 - Ermilova
378   //   1 - Geant3
379   //   2 - Fixed step size
380   //
381
382   fTypeOfStepManager = t;
383   AliInfo(Form("Step Manager type %d was selected",fTypeOfStepManager));
384
385 }
386
387 //_____________________________________________________________________________
388 void AliTRDv1::StepManagerGeant()
389 {
390   //
391   // Slow simulator. Every charged track produces electron cluster as hits
392   // along its path across the drift volume. The step size is set acording
393   // to Bethe-Bloch. The energy distribution of the delta electrons follows
394   // a spectrum taken from Geant3.
395   //
396   // Version by A. Bercuci
397   //
398
399   Int_t    pla = 0;
400   Int_t    cha = 0;
401   Int_t    sec = 0;
402   Int_t    det = 0;
403   Int_t    iPdg;
404   Int_t    qTot;
405
406   Float_t  hits[3];
407   Float_t  charge;
408   Float_t  aMass;
409
410   Double_t pTot     = 0;
411   Double_t eDelta;
412   Double_t betaGamma;
413   Double_t pp;
414   Double_t stepSize = 0;
415
416   Bool_t   drRegion = kFALSE;
417   Bool_t   amRegion = kFALSE;
418
419   TString  cIdPath;
420   Char_t   cIdSector[3];
421            cIdSector[2]  = 0;
422
423   TString  cIdCurrent;
424   TString  cIdSensDr = "J";
425   TString  cIdSensAm = "K";
426   Char_t   cIdChamber[3];
427            cIdChamber[2] = 0;
428
429   TLorentzVector pos;
430   TLorentzVector mom;
431
432   TArrayI        processes;
433
434   const Int_t    kNplan       = AliTRDgeometry::Nplan();
435   const Int_t    kNcham       = AliTRDgeometry::Ncham();
436   const Int_t    kNdetsec     = kNplan * kNcham;
437
438   const Double_t kBig         = 1.0e+12; // Infinitely big
439   const Float_t  kWion        = 23.53;   // Ionization energy
440   const Float_t  kPTotMaxEl   = 0.002;   // Maximum momentum for e+ e- g
441
442   // Minimum energy for the step size adjustment
443   const Float_t  kEkinMinStep = 1.0e-5;
444   // energy threshold for production of delta electrons
445   const Float_t  kECut        = 1.0e4;
446   // Parameters entering the parametrized range for delta electrons
447   const Float_t  kRa          = 5.37e-4;
448   const Float_t  kRb          = 0.9815;
449   const Float_t  kRc          = 3.123e-3;
450   // Gas density -> To be made user adjustable !
451   // [0.85*0.00549+0.15*0.00186 (Xe-CO2 85-15)]
452   const Float_t  kRho         = 0.004945 ; 
453
454   // Plateau value of the energy-loss for electron in xenon
455   // The averaged value (26/3/99)
456   const Float_t  kPlateau     = 1.55;
457   // dN1/dx|min for the gas mixture (90% Xe + 10% CO2)
458   const Float_t  kPrim        = 19.34;  
459   // First ionization potential (eV) for the gas mixture (90% Xe + 10% CO2)
460   const Float_t  kPoti        = 12.1;
461   // PDG code electron
462   const Int_t    kPdgElectron = 11;  
463
464   // Set the maximum step size to a very large number for all
465   // neutral particles and those outside the driftvolume
466   gMC->SetMaxStep(kBig);
467
468   // Use only charged tracks
469   if (( gMC->TrackCharge()       ) &&
470       (!gMC->IsTrackDisappeared())) {
471
472     // Inside a sensitive volume?
473     drRegion = kFALSE;
474     amRegion = kFALSE;
475     cIdCurrent = gMC->CurrentVolName();
476     if (cIdSensDr == cIdCurrent[1]) {
477       drRegion = kTRUE;
478     }
479     if (cIdSensAm == cIdCurrent[1]) {
480       amRegion = kTRUE;
481     }
482     if (drRegion || amRegion) {
483
484       // The hit coordinates and charge
485       gMC->TrackPosition(pos);
486       hits[0] = pos[0];
487       hits[1] = pos[1];
488       hits[2] = pos[2];
489
490       // The sector number (0 - 17), according to standard coordinate system
491       cIdPath      = gGeoManager->GetPath();
492       cIdSector[0] = cIdPath[21];
493       cIdSector[1] = cIdPath[22];
494       sec = atoi(cIdSector);
495
496       // The plane and chamber number
497       cIdChamber[0]   = cIdCurrent[2];
498       cIdChamber[1]   = cIdCurrent[3];
499       Int_t idChamber = (atoi(cIdChamber) % kNdetsec);
500       cha = kNcham - ((Int_t) idChamber / kNplan) - 1;
501       pla = ((Int_t) idChamber % kNplan);
502
503       // The detector number
504       det = fGeometry->GetDetector(pla,cha,sec);
505
506       // Special hits only in the drift region
507       if      ((drRegion) &&
508                (gMC->IsTrackEntering())) {
509
510         // Create a track reference at the entrance of each
511         // chamber that contains the momentum components of the particle
512         gMC->TrackMomentum(mom);
513         AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber());
514
515         // Create the hits from TR photons if electron/positron is
516         // entering the drift volume
517         if ((fTR) && 
518             (TMath::Abs(gMC->TrackPid()) == kPdgElectron)) {
519           CreateTRhit(det);    
520         }
521
522       }
523       else if ((amRegion) && 
524                (gMC->IsTrackExiting())) {
525
526         // Create a track reference at the exit of each
527         // chamber that contains the momentum components of the particle
528         gMC->TrackMomentum(mom);
529         AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber());
530
531       }
532
533       // Calculate the energy of the delta-electrons
534       // modified by Alex Bercuci (A.Bercuci@gsi.de) on 26.01.06
535       // take into account correlation with the underlying GEANT tracking
536       // mechanism. see
537       // http://www-linux.gsi.de/~abercuci/Contributions/TRD/index.html
538       //
539       // determine the most significant process (last on the processes list)
540       // which caused this hit
541       gMC->StepProcesses(processes);
542       Int_t nofprocesses = processes.GetSize();
543       Int_t pid;
544       if (!nofprocesses) {
545         pid = 0;
546       }
547       else {
548         pid =   processes[nofprocesses-1];              
549       }         
550                 
551       // Generate Edep according to GEANT parametrisation
552       eDelta = TMath::Exp(fDeltaG->GetRandom()) - kPoti;
553       eDelta = TMath::Max(eDelta,0.0);
554       Float_t prRange = 0.0;
555       Float_t range   = gMC->TrackLength() - fTrackLength0;
556       // merge GEANT tracker information with locally cooked one
557       if (gAlice->GetMCApp()->GetCurrentTrackNumber() == fPrimaryTrackPid) {
558         if      (pid == 27) { 
559           if (eDelta >= kECut) {                
560             prRange = kRa * eDelta * 0.001
561                     * (1.0 - kRb / (1.0 + kRc * eDelta * 0.001)) / kRho;
562             if (prRange >= (3.7 - range)) {
563               eDelta *= 0.1;
564             }
565           }
566         } 
567         else if (pid ==  1) {   
568           if (eDelta <  kECut) {
569             eDelta *= 0.5;
570           }
571           else {                
572             prRange = kRa * eDelta * 0.001
573                     * (1.0 - kRb / (1.0 + kRc * eDelta * 0.001)) / kRho;
574             if (prRange >= ((AliTRDgeometry::DrThick()
575                            + AliTRDgeometry::AmThick()) - range)) {
576               eDelta *= 0.05;
577             }
578             else {
579               eDelta *= 0.5;
580             }
581           }
582         } 
583         else {
584           eDelta = 0.0;
585         }       
586       } 
587       else {
588         eDelta = 0.0;
589       }
590
591       // Generate the electron cluster size
592       if (eDelta > 0.0) {
593
594         qTot = ((Int_t) (eDelta / kWion) + 1);
595
596         // Create a new dEdx hit
597         AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber()
598               ,det
599               ,hits
600               ,qTot
601               ,gMC->TrackTime()*1.0e06
602               ,drRegion);
603
604       }
605                         
606       // Calculate the maximum step size for the next tracking step
607       // Produce only one hit if Ekin is below cutoff
608       aMass = gMC->TrackMass();
609       if ((gMC->Etot() - aMass) > kEkinMinStep) {
610
611         // The energy loss according to Bethe Bloch
612         iPdg = TMath::Abs(gMC->TrackPid());
613         if ((iPdg != kPdgElectron) ||
614             ((iPdg == kPdgElectron) && 
615              (pTot  < kPTotMaxEl))) {
616           gMC->TrackMomentum(mom);
617           pTot      = mom.Rho();
618           betaGamma = pTot / aMass;
619           pp        = BetheBlochGeant(betaGamma);
620           // Take charge > 1 into account
621           charge     = gMC->TrackCharge();
622           if (TMath::Abs(charge) > 1) {
623             pp = pp * charge*charge;
624           }
625         } 
626         else { 
627           // Electrons above 20 Mev/c are at the plateau
628           pp = kPrim * kPlateau;
629         }
630
631         Int_t nsteps = 0;
632         do {
633           nsteps = gRandom->Poisson(pp);
634         } while(!nsteps);
635         stepSize = 1.0 / nsteps;
636         gMC->SetMaxStep(stepSize);
637
638       }
639
640     }
641
642   }
643
644 }
645
646 //_____________________________________________________________________________
647 void AliTRDv1::StepManagerErmilova()
648 {
649   //
650   // Slow simulator. Every charged track produces electron cluster as hits 
651   // along its path across the drift volume. The step size is set acording
652   // to Bethe-Bloch. The energy distribution of the delta electrons follows
653   // a spectrum taken from Ermilova et al.
654   //
655
656   Int_t    pla = 0;
657   Int_t    cha = 0;
658   Int_t    sec = 0;
659   Int_t    det = 0;
660   Int_t    iPdg;
661   Int_t    qTot;
662
663   Float_t  hits[3];
664   Double_t random[1];
665   Float_t  charge;
666   Float_t  aMass;
667
668   Double_t pTot     = 0.0;
669   Double_t eDelta;
670   Double_t betaGamma;
671   Double_t pp;
672   Double_t stepSize;
673
674   Bool_t   drRegion = kFALSE;
675   Bool_t   amRegion = kFALSE;
676
677   TString  cIdPath;
678   Char_t   cIdSector[3];
679            cIdSector[2]  = 0;
680
681   TString  cIdCurrent;
682   TString  cIdSensDr = "J";
683   TString  cIdSensAm = "K";
684   Char_t   cIdChamber[3];
685            cIdChamber[2] = 0;
686
687   TLorentzVector pos;
688   TLorentzVector mom;
689
690   const Int_t    kNplan       = AliTRDgeometry::Nplan();
691   const Int_t    kNcham       = AliTRDgeometry::Ncham();
692   const Int_t    kNdetsec     = kNplan * kNcham;
693
694   const Double_t kBig         = 1.0e+12; // Infinitely big
695   const Float_t  kWion        = 23.53;   // Ionization energy
696   const Float_t  kPTotMaxEl   = 0.002;   // Maximum momentum for e+ e- g 
697
698   // Minimum energy for the step size adjustment
699   const Float_t  kEkinMinStep = 1.0e-5;
700
701   // Plateau value of the energy-loss for electron in xenon
702   // The averaged value (26/3/99)
703   const Float_t  kPlateau     = 1.55;
704   // dN1/dx|min for the gas mixture (90% Xe + 10% CO2)
705   const Float_t  kPrim        = 48.0;  
706   // First ionization potential (eV) for the gas mixture (90% Xe + 10% CO2)
707   const Float_t  kPoti        = 12.1;
708   // PDG code electron
709   const Int_t    kPdgElectron = 11;  
710
711   // Set the maximum step size to a very large number for all 
712   // neutral particles and those outside the driftvolume
713   gMC->SetMaxStep(kBig); 
714
715   // Use only charged tracks 
716   if (( gMC->TrackCharge()       ) &&
717       (!gMC->IsTrackDisappeared())) {
718
719     // Inside a sensitive volume?
720     drRegion = kFALSE;
721     amRegion = kFALSE;
722     cIdCurrent = gMC->CurrentVolName();
723     if (cIdSensDr == cIdCurrent[1]) {
724       drRegion = kTRUE;
725     }
726     if (cIdSensAm == cIdCurrent[1]) {
727       amRegion = kTRUE;
728     }
729     if (drRegion || amRegion) {
730
731       // The hit coordinates and charge
732       gMC->TrackPosition(pos);
733       hits[0] = pos[0];
734       hits[1] = pos[1];
735       hits[2] = pos[2];
736
737       // The sector number (0 - 17), according to standard coordinate system
738       cIdPath      = gGeoManager->GetPath();
739       cIdSector[0] = cIdPath[21];
740       cIdSector[1] = cIdPath[22];
741       sec = atoi(cIdSector);
742
743       // The plane and chamber number
744       cIdChamber[0] = cIdCurrent[2];
745       cIdChamber[1] = cIdCurrent[3];
746       Int_t idChamber = (atoi(cIdChamber) % kNdetsec);
747       cha = kNcham - ((Int_t) idChamber / kNplan) - 1;
748       pla = ((Int_t) idChamber % kNplan);
749
750       // The detector number
751       det = fGeometry->GetDetector(pla,cha,sec);
752
753       // Special hits only in the drift region
754       if      ((drRegion) &&
755                (gMC->IsTrackEntering())) {
756
757         // Create a track reference at the entrance of each
758         // chamber that contains the momentum components of the particle
759         gMC->TrackMomentum(mom);
760         AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber());
761
762         // Create the hits from TR photons if electron/positron is
763         // entering the drift volume
764         if ((fTR) && 
765             (TMath::Abs(gMC->TrackPid()) == kPdgElectron)) {
766           CreateTRhit(det);    
767         }
768
769       }
770       else if ((amRegion) && 
771                (gMC->IsTrackExiting())) {
772
773         // Create a track reference at the exit of each
774         // chamber that contains the momentum components of the particle
775         gMC->TrackMomentum(mom);
776         AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber());
777
778       }
779
780       // Calculate the energy of the delta-electrons
781       eDelta = TMath::Exp(fDeltaE->GetRandom()) - kPoti;
782       eDelta = TMath::Max(eDelta,0.0);
783
784       // Generate the electron cluster size
785       if (eDelta > 0.0) {
786
787         qTot = ((Int_t) (eDelta / kWion) + 1);
788
789         // Create a new dEdx hit
790         if (drRegion) {
791           AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber()
792                 ,det
793                 ,hits
794                 ,qTot
795                 ,gMC->TrackTime()*1.0e06
796                 ,kTRUE);
797         }
798         else {
799           AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber()
800                 ,det
801                 ,hits
802                 ,qTot
803                 ,gMC->TrackTime()*1.0e06
804                 ,kFALSE);
805         }
806
807       }
808
809       // Calculate the maximum step size for the next tracking step
810       // Produce only one hit if Ekin is below cutoff 
811       aMass = gMC->TrackMass();
812       if ((gMC->Etot() - aMass) > kEkinMinStep) {
813
814         // The energy loss according to Bethe Bloch
815         iPdg  = TMath::Abs(gMC->TrackPid());
816         if ((iPdg != kPdgElectron) ||
817             ((iPdg == kPdgElectron) && 
818              (pTot  < kPTotMaxEl))) {
819           gMC->TrackMomentum(mom);
820           pTot      = mom.Rho();
821           betaGamma = pTot / aMass;
822           pp        = kPrim * BetheBloch(betaGamma);
823           // Take charge > 1 into account
824           charge = gMC->TrackCharge();
825           if (TMath::Abs(charge) > 1) {
826             pp = pp * charge*charge;
827           }
828         } 
829         else { 
830           // Electrons above 20 Mev/c are at the plateau
831           pp = kPrim * kPlateau;
832         }
833       
834         if (pp > 0.0) {
835           do {
836             gMC->GetRandom()->RndmArray(1,random);
837           }
838           while ((random[0] == 1.0) || 
839                  (random[0] == 0.0));
840           stepSize = - TMath::Log(random[0]) / pp; 
841           gMC->SetMaxStep(stepSize);
842         }
843
844       }
845
846     }
847
848   }
849
850 }
851
852 //_____________________________________________________________________________
853 void AliTRDv1::StepManagerFixedStep()
854 {
855   //
856   // Slow simulator. Every charged track produces electron cluster as hits 
857   // along its path across the drift volume. The step size is fixed in
858   // this version of the step manager.
859   //
860
861   // PDG code electron
862   const Int_t   kPdgElectron = 11;
863
864   Int_t    pla = 0;
865   Int_t    cha = 0;
866   Int_t    sec = 0;
867   Int_t    det = 0;
868   Int_t    qTot;
869
870   Float_t  hits[3];
871   Double_t eDep;
872
873   Bool_t   drRegion = kFALSE;
874   Bool_t   amRegion = kFALSE;
875
876   TString  cIdPath;
877   Char_t   cIdSector[3];
878            cIdSector[2]  = 0;
879
880   TString  cIdCurrent;
881   TString  cIdSensDr = "J";
882   TString  cIdSensAm = "K";
883   Char_t   cIdChamber[3];
884            cIdChamber[2] = 0;
885
886   TLorentzVector pos;
887   TLorentzVector mom;
888
889   const Int_t    kNplan       = AliTRDgeometry::Nplan();
890   const Int_t    kNcham       = AliTRDgeometry::Ncham();
891   const Int_t    kNdetsec     = kNplan * kNcham;
892
893   const Double_t kBig         = 1.0e+12;
894
895   const Float_t  kWion        = 23.53;   // Ionization energy
896   const Float_t  kEkinMinStep = 1.0e-5;  // Minimum energy for the step size adjustment
897
898   // Set the maximum step size to a very large number for all 
899   // neutral particles and those outside the driftvolume
900   gMC->SetMaxStep(kBig); 
901
902   // If not charged track or already stopped or disappeared, just return.
903   if ((!gMC->TrackCharge()) || 
904         gMC->IsTrackDisappeared()) {
905     return;
906   }
907
908   // Inside a sensitive volume?
909   cIdCurrent = gMC->CurrentVolName();
910
911   if (cIdSensDr == cIdCurrent[1]) {
912     drRegion = kTRUE;
913   }
914   if (cIdSensAm == cIdCurrent[1]) {
915     amRegion = kTRUE;
916   }
917
918   if ((!drRegion) && 
919       (!amRegion)) {
920     return;
921   }
922
923   // The hit coordinates and charge
924   gMC->TrackPosition(pos);
925   hits[0] = pos[0];
926   hits[1] = pos[1];
927   hits[2] = pos[2];
928
929   // The sector number (0 - 17), according to standard coordinate system
930   cIdPath      = gGeoManager->GetPath();
931   cIdSector[0] = cIdPath[21];
932   cIdSector[1] = cIdPath[22];
933   sec = atoi(cIdSector);
934
935   // The plane and chamber number
936   cIdChamber[0]   = cIdCurrent[2];
937   cIdChamber[1]   = cIdCurrent[3];
938   Int_t idChamber = (atoi(cIdChamber) % kNdetsec);
939   cha = kNcham - ((Int_t) idChamber / kNplan) - 1;
940   pla = ((Int_t) idChamber % kNplan);
941
942   // The detector number
943   det = fGeometry->GetDetector(pla,cha,sec);
944
945   // 0: InFlight 1:Entering 2:Exiting
946   Int_t trkStat = 0;
947
948   // Special hits only in the drift region
949   if      ((drRegion) &&
950            (gMC->IsTrackEntering())) {
951
952     // Create a track reference at the entrance of each
953     // chamber that contains the momentum components of the particle
954     gMC->TrackMomentum(mom);
955     AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber());
956     trkStat = 1;
957
958     // Create the hits from TR photons if electron/positron is
959     // entering the drift volume
960     if ((fTR) && 
961         (TMath::Abs(gMC->TrackPid()) == kPdgElectron)) {
962       CreateTRhit(det);    
963     }
964
965   }
966   else if ((amRegion) && 
967            (gMC->IsTrackExiting())) {
968
969     // Create a track reference at the exit of each
970     // chamber that contains the momentum components of the particle
971     gMC->TrackMomentum(mom);
972     AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber());
973     trkStat = 2;
974
975   }
976   
977   // Calculate the charge according to GEANT Edep
978   // Create a new dEdx hit
979   eDep = TMath::Max(gMC->Edep(),0.0) * 1.0e+09;
980   qTot = (Int_t) (eDep / kWion);
981   if ((qTot) ||
982       (trkStat)) {
983     AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber()
984           ,det
985           ,hits
986           ,qTot
987           ,gMC->TrackTime()*1.0e06
988           ,drRegion);
989   }
990
991   // Set Maximum Step Size
992   // Produce only one hit if Ekin is below cutoff
993   if ((gMC->Etot() - gMC->TrackMass()) < kEkinMinStep) {
994     return;
995   }
996   gMC->SetMaxStep(fStepSize);
997
998 }
999
1000 //_____________________________________________________________________________
1001 Double_t AliTRDv1::BetheBloch(Double_t bg) 
1002 {
1003   //
1004   // Parametrization of the Bethe-Bloch-curve
1005   // The parametrization is the same as for the TPC and is taken from Lehrhaus.
1006   //
1007
1008   // This parameters have been adjusted to averaged values from GEANT
1009   const Double_t kP1    = 7.17960e-02;
1010   const Double_t kP2    = 8.54196;
1011   const Double_t kP3    = 1.38065e-06;
1012   const Double_t kP4    = 5.30972;
1013   const Double_t kP5    = 2.83798;
1014
1015   // Lower cutoff of the Bethe-Bloch-curve to limit step sizes
1016   const Double_t kBgMin = 0.8;
1017   const Double_t kBBMax = 6.83298;
1018
1019   if (bg > kBgMin) {
1020     Double_t yy = bg / TMath::Sqrt(1.0 + bg*bg);
1021     Double_t aa = TMath::Power(yy,kP4);
1022     Double_t bb = TMath::Power((1.0/bg),kP5);
1023              bb = TMath::Log(kP3 + bb);
1024     return ((kP2 - aa - bb) * kP1 / aa);
1025   }
1026   else {
1027     return kBBMax;
1028   }
1029
1030 }
1031
1032 //_____________________________________________________________________________
1033 Double_t AliTRDv1::BetheBlochGeant(Double_t bg)
1034 {
1035   //
1036   // Return dN/dx (number of primary collisions per centimeter)
1037   // for given beta*gamma factor.
1038   //
1039   // Implemented by K.Oyama according to GEANT 3 parametrization shown in
1040   // A.Andronic's webpage: http://www-alice.gsi.de/trd/papers/dedx/dedx.html
1041   // This must be used as a set with IntSpecGeant.
1042   //
1043
1044   Int_t    i = 0;
1045
1046   Double_t arrG[20]  = {     1.100000,     1.200000,     1.300000,     1.500000
1047                        ,     1.800000,     2.000000,     2.500000,     3.000000
1048                        ,     4.000000,     7.000000,    10.000000,    20.000000
1049                        ,    40.000000,    70.000000,   100.000000,   300.000000
1050                        ,   600.000000,  1000.000000,  3000.000000, 10000.000000 };
1051
1052   Double_t arrNC[20] = {    75.009056,    45.508083,    35.299252,    27.116327
1053                        ,    22.734999,    21.411915,    19.934095,    19.449375
1054                        ,    19.344431,    20.185553,    21.027925,    22.912676
1055                        ,    24.933352,    26.504053,    27.387468,    29.566597
1056                        ,    30.353779,    30.787134,    31.129285,    31.157350 };
1057
1058   // Betagamma to gamma
1059   Double_t g = TMath::Sqrt(1.0 + bg*bg);
1060
1061   // Find the index just before the point we need.
1062   for (i = 0; i < 18; i++) {
1063     if ((arrG[i]   < g) && 
1064         (arrG[i+1] > g)) {
1065       break;
1066     }
1067   }
1068
1069   // Simple interpolation.
1070   Double_t pp = ((arrNC[i+1] - arrNC[i]) / (arrG[i+1]  - arrG[i])) 
1071               * (g - arrG[i]) + arrNC[i];
1072
1073   return pp;
1074
1075 }
1076
1077 //_____________________________________________________________________________
1078 Double_t Ermilova(Double_t *x, Double_t *)
1079 {
1080   //
1081   // Calculates the delta-ray energy distribution according to Ermilova.
1082   // Logarithmic scale !
1083   //
1084
1085   Double_t energy;
1086   Double_t dpos;
1087   Double_t dnde;
1088
1089   Int_t    pos1;
1090   Int_t    pos2;
1091
1092   const Int_t kNv = 31;
1093
1094   Float_t vxe[kNv] = {  2.3026,  2.9957,  3.4012,  3.6889,  3.9120  
1095                      ,  4.0943,  4.2485,  4.3820,  4.4998,  4.6052
1096                      ,  4.7005,  5.0752,  5.2983,  5.7038,  5.9915
1097                      ,  6.2146,  6.5221,  6.9078,  7.3132,  7.6009
1098                      ,  8.0064,  8.5172,  8.6995,  8.9872,  9.2103
1099                      ,  9.4727,  9.9035, 10.3735, 10.5966, 10.8198
1100                      , 11.5129 };
1101
1102   Float_t vye[kNv] = { 80.0,    31.0,    23.3,    21.1,    21.0
1103                      , 20.9,    20.8,    20.0,    16.0,    11.0
1104                      ,  8.0,     6.0,     5.2,     4.6,     4.0
1105                      ,  3.5,     3.0,     1.4,     0.67,    0.44
1106                      ,  0.3,     0.18,    0.12,    0.08,    0.056
1107                      ,  0.04,    0.023,   0.015,   0.011,    0.01
1108                      ,  0.004  };
1109
1110   energy = x[0];
1111
1112   // Find the position 
1113   pos1 = 0;
1114   pos2 = 0;
1115   dpos = 0;
1116   do {
1117     dpos = energy - vxe[pos2++];
1118   } 
1119   while (dpos > 0);
1120   pos2--; 
1121   if (pos2 > kNv) {
1122     pos2 = kNv - 1;
1123   }
1124   pos1 = pos2 - 1;
1125
1126   // Differentiate between the sampling points
1127   dnde = (vye[pos1] - vye[pos2]) / (vxe[pos2] - vxe[pos1]);
1128
1129   return dnde;
1130
1131 }
1132
1133 //_____________________________________________________________________________
1134 Double_t IntSpecGeant(Double_t *x, Double_t *)
1135 {
1136   //
1137   // Integrated spectrum from Geant3
1138   //
1139
1140   const Int_t npts = 83;
1141   Double_t arre[npts]    = {  2.421257,     2.483278,    2.534301,     2.592230
1142                            ,  2.672067,     2.813299,    3.015059,     3.216819
1143                            ,  3.418579,     3.620338,    3.868209,     3.920198
1144                            ,  3.978284,     4.063923,    4.186264,     4.308605
1145                            ,  4.430946,     4.553288,    4.724261,     4.837736
1146                            ,  4.999842,     5.161949,    5.324056,     5.486163
1147                            ,  5.679688,     5.752998,    5.857728,     5.962457
1148                            ,  6.067185,     6.171914,    6.315653,     6.393674
1149                            ,  6.471694,     6.539689,    6.597658,     6.655627
1150                            ,  6.710957,     6.763648,    6.816338,     6.876198
1151                            ,  6.943227,     7.010257,    7.106285,     7.252151
1152                            ,  7.460531,     7.668911,    7.877290,     8.085670
1153                            ,  8.302979,     8.353585,    8.413120,     8.483500
1154                            ,  8.541030,     8.592857,    8.668865,     8.820485
1155                            ,  9.037086,     9.253686,    9.470286,     9.686887
1156                            ,  9.930838,     9.994655,   10.085822,    10.176990
1157                            , 10.268158,    10.359325,   10.503614,    10.627565
1158                            , 10.804637,    10.981709,   11.158781,    11.335854
1159                            , 11.593397,    11.781165,   12.049404,    12.317644
1160                            , 12.585884,    12.854123,   14.278421,    16.975889
1161                            , 20.829416,    24.682943,   28.536469 };
1162
1163   Double_t arrdnde[npts] = { 10.960000,    10.960000,   10.359500,     9.811340
1164                            ,  9.1601500,    8.206670,    6.919630,     5.655430
1165                            ,  4.6221300,    3.777610,    3.019560,     2.591950
1166                            ,  2.5414600,    2.712920,    3.327460,     4.928240
1167                            ,  7.6185300,   10.966700,   12.225800,     8.094750
1168                            ,  3.3586900,    1.553650,    1.209600,     1.263840
1169                            ,  1.3241100,    1.312140,    1.255130,     1.165770
1170                            ,  1.0594500,    0.945450,    0.813231,     0.699837
1171                            ,  0.6235580,    2.260990,    2.968350,     2.240320
1172                            ,  1.7988300,    1.553300,    1.432070,     1.535520
1173                            ,  1.4429900,    1.247990,    1.050750,     0.829549
1174                            ,  0.5900280,    0.395897,    0.268741,     0.185320
1175                            ,  0.1292120,    0.103545,    0.0949525,    0.101535
1176                            ,  0.1276380,    0.134216,    0.123816,     0.104557
1177                            ,  0.0751843,    0.0521745,   0.0373546,    0.0275391
1178                            ,  0.0204713,    0.0169234,   0.0154552,    0.0139194
1179                            ,  0.0125592,    0.0113638,   0.0107354,    0.0102137
1180                            ,  0.00845984,   0.00683338,  0.00556836,   0.00456874
1181                            ,  0.0036227,    0.00285991,  0.00226664,   0.00172234
1182                            ,  0.00131226,   0.00100284,  0.000465492,  7.26607e-05
1183                            ,  3.63304e-06,  0.0000000,   0.0000000   };
1184
1185   Int_t    i;
1186   Double_t energy = x[0];
1187
1188   for (i = 0; i < npts; i++) {
1189     if (energy < arre[i]) {
1190       break;
1191     }
1192   }
1193
1194   if (i == 0) {
1195     AliErrorGeneral("AliTRDv1::IntSpecGeant","Given energy value is too small or zero");
1196   }
1197
1198   return arrdnde[i];
1199
1200 }