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 ////////////////////////////////////////////////////////////////////////////
20 // Transition Radiation Detector version 1 -- slow simulator //
22 ////////////////////////////////////////////////////////////////////////////
27 #include <TLorentzVector.h>
31 #include <TVirtualMC.h>
32 #include <TGeoManager.h>
33 #include <TGeoMatrix.h>
34 #include <TGeoPhysicalNode.h>
38 #include "AliTrackReference.h"
41 #include "AliGeomManager.h"
43 #include "AliTRDgeometry.h"
44 #include "AliTRDCommonParam.h"
45 #include "AliTRDhit.h"
46 #include "AliTRDsimTR.h"
51 //_____________________________________________________________________________
56 ,fTypeOfStepManager(0)
65 // Default constructor
70 //_____________________________________________________________________________
71 AliTRDv1::AliTRDv1(const char *name, const char *title)
75 ,fTypeOfStepManager(2)
84 // Standard constructor for Transition Radiation Detector version 1
87 SetBufferSize(128000);
89 if (AliTRDCommonParam::Instance()->IsXenon()) {
90 fWion = 23.53; // Ionization energy XeCO2 (85/15)
92 else if (AliTRDCommonParam::Instance()->IsArgon()) {
93 fWion = 27.21; // Ionization energy ArCO2 (82/18)
96 AliFatal("Wrong gas mixture");
102 //_____________________________________________________________________________
103 AliTRDv1::~AliTRDv1()
106 // AliTRDv1 destructor
126 //_____________________________________________________________________________
127 void AliTRDv1::AddAlignableVolumes() const
130 // Create entries for alignable volumes associating the symbolic volume
131 // name with the corresponding volume path. Needs to be syncronized with
132 // eventual changes in the geometry.
138 TString vpStr = "ALIC_1/B077_1/BSEGMO";
139 TString vpApp1 = "_1/BTRD";
140 TString vpApp2 = "_1";
141 TString vpApp3a = "/UTR1_1/UTS1_1/UTI1_1/UT";
142 TString vpApp3b = "/UTR2_1/UTS2_1/UTI2_1/UT";
143 TString vpApp3c = "/UTR3_1/UTS3_1/UTI3_1/UT";
145 TString snStr = "TRD/sm";
146 TString snApp1 = "/st";
147 TString snApp2 = "/pl";
151 // The symbolic names are: TRD/sm00
155 for (Int_t isector = 0; isector < AliTRDgeometry::Nsector(); isector++) {
164 symName += Form("%02d",isector);
166 gGeoManager->SetAlignableEntry(symName.Data(),volPath.Data());
171 // The readout chambers
172 // The symbolic names are: TRD/sm00/st0/pl0
176 AliGeomManager::ELayerID idTRD1 = AliGeomManager::kTRD1;
179 for (Int_t isector = 0; isector < AliTRDgeometry::Nsector(); isector++) {
181 if (fGeometry->GetSMstatus(isector) == 0) continue;
183 for (Int_t istack = 0; istack < AliTRDgeometry::Nstack(); istack++) {
184 for (Int_t ilayer = 0; ilayer < AliTRDgeometry::Nlayer(); ilayer++) {
186 layer = idTRD1 + ilayer;
187 modUID = AliGeomManager::LayerToVolUIDSafe(layer,isector*5+istack);
189 Int_t idet = AliTRDgeometry::GetDetectorSec(ilayer,istack);
212 volPath += Form("%02d",idet);
216 symName += Form("%02d",isector);
222 TGeoPNEntry *alignableEntry =
223 gGeoManager->SetAlignableEntry(symName.Data(),volPath.Data(),modUID);
225 // Add the tracking to local matrix following the TPC example
226 if (alignableEntry) {
227 TGeoHMatrix *globMatrix = alignableEntry->GetGlobalOrig();
228 Double_t sectorAngle = 20.0 * (isector % 18) + 10.0;
229 TGeoHMatrix *t2lMatrix = new TGeoHMatrix();
230 t2lMatrix->RotateZ(sectorAngle);
231 t2lMatrix->MultiplyLeft(&(globMatrix->Inverse()));
232 alignableEntry->SetMatrix(t2lMatrix);
235 AliError(Form("Alignable entry %s is not valid!",symName.Data()));
244 //_____________________________________________________________________________
245 void AliTRDv1::CreateGeometry()
248 // Create the GEANT geometry for the Transition Radiation Detector - Version 1
249 // This version covers the full azimuth.
252 // Check that FRAME is there otherwise we have no place where to put the TRD
253 AliModule* frame = gAlice->GetModule("FRAME");
255 AliError("TRD needs FRAME to be present\n");
259 // Define the chambers
260 AliTRD::CreateGeometry();
264 //_____________________________________________________________________________
265 void AliTRDv1::CreateMaterials()
268 // Create materials for the Transition Radiation Detector version 1
271 AliTRD::CreateMaterials();
275 //_____________________________________________________________________________
276 void AliTRDv1::CreateTRhit(Int_t det)
279 // Creates an electron cluster from a TR photon.
280 // The photon is assumed to be created a the end of the radiator. The
281 // distance after which it deposits its energy takes into account the
282 // absorbtion of the entrance window and of the gas mixture in drift
286 // Maximum number of TR photons per track
287 const Int_t kNTR = 50;
296 gMC->TrackMomentum(mom);
297 Float_t pTot = mom.Rho();
298 fTR->CreatePhotons(11,pTot,nTR,eTR);
300 AliFatal(Form("Boundary error: nTR = %d, kNTR = %d",nTR,kNTR));
303 // Loop through the TR photons
304 for (Int_t iTR = 0; iTR < nTR; iTR++) {
306 Float_t energyMeV = eTR[iTR] * 0.001;
307 Float_t energyeV = eTR[iTR] * 1000.0;
308 Float_t absLength = 0.0;
311 // Take the absorbtion in the entrance window into account
312 Double_t muMy = fTR->GetMuMy(energyMeV);
313 sigma = muMy * fFoilDensity;
315 absLength = gRandom->Exp(1.0/sigma);
316 if (absLength < AliTRDgeometry::MyThick()) {
324 // The absorbtion cross sections in the drift gas
325 // Gas-mixture (Xe/CO2)
327 if (AliTRDCommonParam::Instance()->IsXenon()) {
328 muNo = fTR->GetMuXe(energyMeV);
330 else if (AliTRDCommonParam::Instance()->IsArgon()) {
331 muNo = fTR->GetMuAr(energyMeV);
333 Double_t muCO = fTR->GetMuCO(energyMeV);
334 sigma = (fGasNobleFraction * muNo + (1.0 - fGasNobleFraction) * muCO)
338 // The distance after which the energy of the TR photon
341 absLength = gRandom->Exp(1.0/sigma);
342 if (absLength > (AliTRDgeometry::DrThick()
343 + AliTRDgeometry::AmThick())) {
351 // The position of the absorbtion
353 gMC->TrackPosition(pos);
354 posHit[0] = pos[0] + mom[0] / pTot * absLength;
355 posHit[1] = pos[1] + mom[1] / pTot * absLength;
356 posHit[2] = pos[2] + mom[2] / pTot * absLength;
359 Int_t q = ((Int_t) (energyeV / fWion));
361 // Add the hit to the array. TR photon hits are marked
362 // by negative charge
363 AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber()
367 ,gMC->TrackTime()*1.0e06
374 //_____________________________________________________________________________
375 void AliTRDv1::Init()
378 // Initialise Transition Radiation Detector after geometry has been built.
383 AliDebug(1,"Slow simulator\n");
385 // Switch on TR simulation as default
387 AliInfo("TR simulation off");
390 fTR = new AliTRDsimTR();
393 // First ionization potential (eV) for the gas mixture
394 const Float_t kPoti = 12.1;
395 // Maximum energy (50 keV);
396 const Float_t kEend = 50000.0;
397 // Ermilova distribution for the delta-ray spectrum
398 Float_t poti = TMath::Log(kPoti);
399 Float_t eEnd = TMath::Log(kEend);
401 // Ermilova distribution for the delta-ray spectrum
402 fDeltaE = new TF1("deltae" ,Ermilova ,poti,eEnd,0);
404 // Geant3 distribution for the delta-ray spectrum
405 fDeltaG = new TF1("deltag",IntSpecGeant,2.421257,28.536469,0);
407 AliDebug(1,"+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++");
411 //_____________________________________________________________________________
412 void AliTRDv1::StepManager()
415 // Slow simulator. Every charged track produces electron cluster as hits
416 // along its path across the drift volume.
419 switch (fTypeOfStepManager) {
421 StepManagerErmilova();
427 StepManagerFixedStep();
430 AliWarning("Not a valid Step Manager.");
435 //_____________________________________________________________________________
436 void AliTRDv1::SelectStepManager(Int_t t)
439 // Selects a step manager type:
442 // 2 - Fixed step size
445 fTypeOfStepManager = t;
446 AliInfo(Form("Step Manager type %d was selected",fTypeOfStepManager));
450 //_____________________________________________________________________________
451 void AliTRDv1::StepManagerGeant()
454 // Slow simulator. Every charged track produces electron cluster as hits
455 // along its path across the drift volume. The step size is set acording
456 // to Bethe-Bloch. The energy distribution of the delta electrons follows
457 // a spectrum taken from Geant3.
459 // Works only for Xe/CO2!!
461 // Version by A. Bercuci
479 Double_t stepSize = 0;
481 Bool_t drRegion = kFALSE;
482 Bool_t amRegion = kFALSE;
489 TString cIdSensDr = "J";
490 TString cIdSensAm = "K";
491 Char_t cIdChamber[3];
499 const Int_t kNlayer = AliTRDgeometry::Nlayer();
500 const Int_t kNstack = AliTRDgeometry::Nstack();
501 const Int_t kNdetsec = kNlayer * kNstack;
503 const Double_t kBig = 1.0e+12; // Infinitely big
504 const Float_t kPTotMaxEl = 0.002; // Maximum momentum for e+ e- g
506 // Minimum energy for the step size adjustment
507 const Float_t kEkinMinStep = 1.0e-5;
508 // energy threshold for production of delta electrons
509 const Float_t kECut = 1.0e4;
510 // Parameters entering the parametrized range for delta electrons
511 const Float_t kRa = 5.37e-4;
512 const Float_t kRb = 0.9815;
513 const Float_t kRc = 3.123e-3;
514 // Gas density -> To be made user adjustable !
515 // [0.85*0.00549+0.15*0.00186 (Xe-CO2 85-15)]
516 const Float_t kRho = 0.004945 ;
518 // Plateau value of the energy-loss for electron in xenon
519 // The averaged value (26/3/99)
520 const Float_t kPlateau = 1.55;
521 // dN1/dx|min for the gas mixture (90% Xe + 10% CO2)
522 const Float_t kPrim = 19.34;
523 // First ionization potential (eV) for the gas mixture (90% Xe + 10% CO2)
524 const Float_t kPoti = 12.1;
526 const Int_t kPdgElectron = 11;
528 // Set the maximum step size to a very large number for all
529 // neutral particles and those outside the driftvolume
530 gMC->SetMaxStep(kBig);
532 // Use only charged tracks
533 if (( gMC->TrackCharge() ) &&
534 (!gMC->IsTrackDisappeared())) {
536 // Inside a sensitive volume?
539 cIdCurrent = gMC->CurrentVolName();
540 if (cIdSensDr == cIdCurrent[1]) {
543 if (cIdSensAm == cIdCurrent[1]) {
546 if (drRegion || amRegion) {
548 // The hit coordinates and charge
549 gMC->TrackPosition(pos);
554 // The sector number (0 - 17), according to standard coordinate system
555 cIdPath = gGeoManager->GetPath();
556 cIdSector[0] = cIdPath[21];
557 cIdSector[1] = cIdPath[22];
558 sector = atoi(cIdSector);
560 // The layer and stack number
561 cIdChamber[0] = cIdCurrent[2];
562 cIdChamber[1] = cIdCurrent[3];
563 Int_t idChamber = (atoi(cIdChamber) % kNdetsec);
564 stack = ((Int_t) idChamber / kNlayer);
565 layer = ((Int_t) idChamber % kNlayer);
567 // The detector number
568 det = fGeometry->GetDetector(layer,stack,sector);
570 // Special hits only in the drift region
572 (gMC->IsTrackEntering())) {
574 // Create a track reference at the entrance of each
575 // chamber that contains the momentum components of the particle
576 gMC->TrackMomentum(mom);
577 AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber(), AliTrackReference::kTRD);
579 // Create the hits from TR photons if electron/positron is
580 // entering the drift volume
582 (TMath::Abs(gMC->TrackPid()) == kPdgElectron)) {
587 else if ((amRegion) &&
588 (gMC->IsTrackExiting())) {
590 // Create a track reference at the exit of each
591 // chamber that contains the momentum components of the particle
592 gMC->TrackMomentum(mom);
593 AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber(), AliTrackReference::kTRD);
597 // Calculate the energy of the delta-electrons
598 // modified by Alex Bercuci (A.Bercuci@gsi.de) on 26.01.06
599 // take into account correlation with the underlying GEANT tracking
601 // http://www-linux.gsi.de/~abercuci/Contributions/TRD/index.html
603 // determine the most significant process (last on the processes list)
604 // which caused this hit
605 gMC->StepProcesses(processes);
606 Int_t nofprocesses = processes.GetSize();
612 pid = processes[nofprocesses-1];
615 // Generate Edep according to GEANT parametrisation
616 eDelta = TMath::Exp(fDeltaG->GetRandom()) - kPoti;
617 eDelta = TMath::Max(eDelta,0.0);
618 Float_t prRange = 0.0;
619 Float_t range = gMC->TrackLength() - fTrackLength0;
620 // merge GEANT tracker information with locally cooked one
621 if (gAlice->GetMCApp()->GetCurrentTrackNumber() == fPrimaryTrackPid) {
623 if (eDelta >= kECut) {
624 prRange = kRa * eDelta * 0.001
625 * (1.0 - kRb / (1.0 + kRc * eDelta * 0.001)) / kRho;
626 if (prRange >= (3.7 - range)) {
632 if (eDelta < kECut) {
636 prRange = kRa * eDelta * 0.001
637 * (1.0 - kRb / (1.0 + kRc * eDelta * 0.001)) / kRho;
638 if (prRange >= ((AliTRDgeometry::DrThick()
639 + AliTRDgeometry::AmThick()) - range)) {
655 // Generate the electron cluster size
658 qTot = ((Int_t) (eDelta / fWion) + 1);
660 // Create a new dEdx hit
661 AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber()
665 ,gMC->TrackTime()*1.0e06
670 // Calculate the maximum step size for the next tracking step
671 // Produce only one hit if Ekin is below cutoff
672 aMass = gMC->TrackMass();
673 if ((gMC->Etot() - aMass) > kEkinMinStep) {
675 // The energy loss according to Bethe Bloch
676 iPdg = TMath::Abs(gMC->TrackPid());
677 if ((iPdg != kPdgElectron) ||
678 ((iPdg == kPdgElectron) &&
679 (pTot < kPTotMaxEl))) {
680 gMC->TrackMomentum(mom);
682 betaGamma = pTot / aMass;
683 pp = BetheBlochGeant(betaGamma);
684 // Take charge > 1 into account
685 charge = gMC->TrackCharge();
686 if (TMath::Abs(charge) > 1) {
687 pp = pp * charge*charge;
691 // Electrons above 20 Mev/c are at the plateau
692 pp = kPrim * kPlateau;
697 nsteps = gRandom->Poisson(pp);
699 stepSize = 1.0 / nsteps;
700 gMC->SetMaxStep(stepSize);
710 //_____________________________________________________________________________
711 void AliTRDv1::StepManagerErmilova()
714 // Slow simulator. Every charged track produces electron cluster as hits
715 // along its path across the drift volume. The step size is set acording
716 // to Bethe-Bloch. The energy distribution of the delta electrons follows
717 // a spectrum taken from Ermilova et al.
719 // Works only for Xe/CO2!!
740 Bool_t drRegion = kFALSE;
741 Bool_t amRegion = kFALSE;
748 TString cIdSensDr = "J";
749 TString cIdSensAm = "K";
750 Char_t cIdChamber[3];
756 const Int_t kNlayer = AliTRDgeometry::Nlayer();
757 const Int_t kNstack = AliTRDgeometry::Nstack();
758 const Int_t kNdetsec = kNlayer * kNstack;
760 const Double_t kBig = 1.0e+12; // Infinitely big
761 const Float_t kPTotMaxEl = 0.002; // Maximum momentum for e+ e- g
763 // Minimum energy for the step size adjustment
764 const Float_t kEkinMinStep = 1.0e-5;
766 // Plateau value of the energy-loss for electron in xenon
767 // The averaged value (26/3/99)
768 const Float_t kPlateau = 1.55;
769 // dN1/dx|min for the gas mixture (90% Xe + 10% CO2)
770 const Float_t kPrim = 48.0;
771 // First ionization potential (eV) for the gas mixture (90% Xe + 10% CO2)
772 const Float_t kPoti = 12.1;
774 const Int_t kPdgElectron = 11;
776 // Set the maximum step size to a very large number for all
777 // neutral particles and those outside the driftvolume
778 gMC->SetMaxStep(kBig);
780 // Use only charged tracks
781 if (( gMC->TrackCharge() ) &&
782 (!gMC->IsTrackDisappeared())) {
784 // Inside a sensitive volume?
787 cIdCurrent = gMC->CurrentVolName();
788 if (cIdSensDr == cIdCurrent[1]) {
791 if (cIdSensAm == cIdCurrent[1]) {
794 if (drRegion || amRegion) {
796 // The hit coordinates and charge
797 gMC->TrackPosition(pos);
802 // The sector number (0 - 17), according to standard coordinate system
803 cIdPath = gGeoManager->GetPath();
804 cIdSector[0] = cIdPath[21];
805 cIdSector[1] = cIdPath[22];
806 sector = atoi(cIdSector);
808 // The plane and chamber number
809 cIdChamber[0] = cIdCurrent[2];
810 cIdChamber[1] = cIdCurrent[3];
811 Int_t idChamber = (atoi(cIdChamber) % kNdetsec);
812 stack = ((Int_t) idChamber / kNlayer);
813 layer = ((Int_t) idChamber % kNlayer);
815 // The detector number
816 det = fGeometry->GetDetector(layer,stack,sector);
818 // Special hits only in the drift region
820 (gMC->IsTrackEntering())) {
822 // Create a track reference at the entrance of each
823 // chamber that contains the momentum components of the particle
824 gMC->TrackMomentum(mom);
825 AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber(), AliTrackReference::kTRD);
827 // Create the hits from TR photons if electron/positron is
828 // entering the drift volume
830 (TMath::Abs(gMC->TrackPid()) == kPdgElectron)) {
835 else if ((amRegion) &&
836 (gMC->IsTrackExiting())) {
838 // Create a track reference at the exit of each
839 // chamber that contains the momentum components of the particle
840 gMC->TrackMomentum(mom);
841 AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber(), AliTrackReference::kTRD);
845 // Calculate the energy of the delta-electrons
846 eDelta = TMath::Exp(fDeltaE->GetRandom()) - kPoti;
847 eDelta = TMath::Max(eDelta,0.0);
849 // Generate the electron cluster size
852 qTot = ((Int_t) (eDelta / fWion) + 1);
854 // Create a new dEdx hit
856 AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber()
860 ,gMC->TrackTime()*1.0e06
864 AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber()
868 ,gMC->TrackTime()*1.0e06
874 // Calculate the maximum step size for the next tracking step
875 // Produce only one hit if Ekin is below cutoff
876 aMass = gMC->TrackMass();
877 if ((gMC->Etot() - aMass) > kEkinMinStep) {
879 // The energy loss according to Bethe Bloch
880 iPdg = TMath::Abs(gMC->TrackPid());
881 if ((iPdg != kPdgElectron) ||
882 ((iPdg == kPdgElectron) &&
883 (pTot < kPTotMaxEl))) {
884 gMC->TrackMomentum(mom);
886 betaGamma = pTot / aMass;
887 pp = kPrim * BetheBloch(betaGamma);
888 // Take charge > 1 into account
889 charge = gMC->TrackCharge();
890 if (TMath::Abs(charge) > 1) {
891 pp = pp * charge*charge;
895 // Electrons above 20 Mev/c are at the plateau
896 pp = kPrim * kPlateau;
901 gMC->GetRandom()->RndmArray(1,random);
903 while ((random[0] == 1.0) ||
905 stepSize = - TMath::Log(random[0]) / pp;
906 gMC->SetMaxStep(stepSize);
917 //_____________________________________________________________________________
918 void AliTRDv1::StepManagerFixedStep()
921 // Slow simulator. Every charged track produces electron cluster as hits
922 // along its path across the drift volume. The step size is fixed in
923 // this version of the step manager.
925 // Works for Xe/CO2 as well as Ar/CO2
929 const Int_t kPdgElectron = 11;
940 Bool_t drRegion = kFALSE;
941 Bool_t amRegion = kFALSE;
948 TString cIdSensDr = "J";
949 TString cIdSensAm = "K";
950 Char_t cIdChamber[3];
956 const Int_t kNlayer = AliTRDgeometry::Nlayer();
957 const Int_t kNstack = AliTRDgeometry::Nstack();
958 const Int_t kNdetsec = kNlayer * kNstack;
960 const Double_t kBig = 1.0e+12;
961 const Float_t kEkinMinStep = 1.0e-5; // Minimum energy for the step size adjustment
963 // Set the maximum step size to a very large number for all
964 // neutral particles and those outside the driftvolume
965 gMC->SetMaxStep(kBig);
967 // If not charged track or already stopped or disappeared, just return.
968 if ((!gMC->TrackCharge()) ||
969 gMC->IsTrackDisappeared()) {
973 // Inside a sensitive volume?
974 cIdCurrent = gMC->CurrentVolName();
976 if (cIdSensDr == cIdCurrent[1]) {
979 if (cIdSensAm == cIdCurrent[1]) {
988 // The hit coordinates and charge
989 gMC->TrackPosition(pos);
994 // The sector number (0 - 17), according to standard coordinate system
995 cIdPath = gGeoManager->GetPath();
996 cIdSector[0] = cIdPath[21];
997 cIdSector[1] = cIdPath[22];
998 sector = atoi(cIdSector);
1000 // The plane and chamber number
1001 cIdChamber[0] = cIdCurrent[2];
1002 cIdChamber[1] = cIdCurrent[3];
1003 Int_t idChamber = (atoi(cIdChamber) % kNdetsec);
1004 stack = ((Int_t) idChamber / kNlayer);
1005 layer = ((Int_t) idChamber % kNlayer);
1007 // The detector number
1008 det = fGeometry->GetDetector(layer,stack,sector);
1010 // 0: InFlight 1:Entering 2:Exiting
1013 // Special hits only in the drift region
1015 (gMC->IsTrackEntering())) {
1017 // Create a track reference at the entrance of each
1018 // chamber that contains the momentum components of the particle
1019 gMC->TrackMomentum(mom);
1020 AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber(), AliTrackReference::kTRD);
1023 // Create the hits from TR photons if electron/positron is
1024 // entering the drift volume
1026 (TMath::Abs(gMC->TrackPid()) == kPdgElectron)) {
1031 else if ((amRegion) &&
1032 (gMC->IsTrackExiting())) {
1034 // Create a track reference at the exit of each
1035 // chamber that contains the momentum components of the particle
1036 gMC->TrackMomentum(mom);
1037 AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber(), AliTrackReference::kTRD);
1042 // Calculate the charge according to GEANT Edep
1043 // Create a new dEdx hit
1044 eDep = TMath::Max(gMC->Edep(),0.0) * 1.0e+09;
1045 qTot = (Int_t) (eDep / fWion);
1048 AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber()
1052 ,gMC->TrackTime()*1.0e06
1056 // Set Maximum Step Size
1057 // Produce only one hit if Ekin is below cutoff
1058 if ((gMC->Etot() - gMC->TrackMass()) < kEkinMinStep) {
1061 gMC->SetMaxStep(fStepSize);
1065 //_____________________________________________________________________________
1066 Double_t AliTRDv1::BetheBloch(Double_t bg)
1069 // Parametrization of the Bethe-Bloch-curve
1070 // The parametrization is the same as for the TPC and is taken from Lehrhaus.
1073 // This parameters have been adjusted to averaged values from GEANT
1074 const Double_t kP1 = 7.17960e-02;
1075 const Double_t kP2 = 8.54196;
1076 const Double_t kP3 = 1.38065e-06;
1077 const Double_t kP4 = 5.30972;
1078 const Double_t kP5 = 2.83798;
1080 // Lower cutoff of the Bethe-Bloch-curve to limit step sizes
1081 const Double_t kBgMin = 0.8;
1082 const Double_t kBBMax = 6.83298;
1085 Double_t yy = bg / TMath::Sqrt(1.0 + bg*bg);
1086 Double_t aa = TMath::Power(yy,kP4);
1087 Double_t bb = TMath::Power((1.0/bg),kP5);
1088 bb = TMath::Log(kP3 + bb);
1089 return ((kP2 - aa - bb) * kP1 / aa);
1097 //_____________________________________________________________________________
1098 Double_t AliTRDv1::BetheBlochGeant(Double_t bg)
1101 // Return dN/dx (number of primary collisions per centimeter)
1102 // for given beta*gamma factor.
1104 // Implemented by K.Oyama according to GEANT 3 parametrization shown in
1105 // A.Andronic's webpage: http://www-alice.gsi.de/trd/papers/dedx/dedx.html
1106 // This must be used as a set with IntSpecGeant.
1111 Double_t arrG[20] = { 1.100000, 1.200000, 1.300000, 1.500000
1112 , 1.800000, 2.000000, 2.500000, 3.000000
1113 , 4.000000, 7.000000, 10.000000, 20.000000
1114 , 40.000000, 70.000000, 100.000000, 300.000000
1115 , 600.000000, 1000.000000, 3000.000000, 10000.000000 };
1117 Double_t arrNC[20] = { 75.009056, 45.508083, 35.299252, 27.116327
1118 , 22.734999, 21.411915, 19.934095, 19.449375
1119 , 19.344431, 20.185553, 21.027925, 22.912676
1120 , 24.933352, 26.504053, 27.387468, 29.566597
1121 , 30.353779, 30.787134, 31.129285, 31.157350 };
1123 // Betagamma to gamma
1124 Double_t g = TMath::Sqrt(1.0 + bg*bg);
1126 // Find the index just before the point we need.
1127 for (i = 0; i < 18; i++) {
1128 if ((arrG[i] < g) &&
1134 // Simple interpolation.
1135 Double_t pp = ((arrNC[i+1] - arrNC[i]) / (arrG[i+1] - arrG[i]))
1136 * (g - arrG[i]) + arrNC[i];
1142 //_____________________________________________________________________________
1143 Double_t Ermilova(Double_t *x, Double_t *)
1146 // Calculates the delta-ray energy distribution according to Ermilova.
1147 // Logarithmic scale !
1157 const Int_t kNv = 31;
1159 Float_t vxe[kNv] = { 2.3026, 2.9957, 3.4012, 3.6889, 3.9120
1160 , 4.0943, 4.2485, 4.3820, 4.4998, 4.6052
1161 , 4.7005, 5.0752, 5.2983, 5.7038, 5.9915
1162 , 6.2146, 6.5221, 6.9078, 7.3132, 7.6009
1163 , 8.0064, 8.5172, 8.6995, 8.9872, 9.2103
1164 , 9.4727, 9.9035, 10.3735, 10.5966, 10.8198
1167 Float_t vye[kNv] = { 80.0, 31.0, 23.3, 21.1, 21.0
1168 , 20.9, 20.8, 20.0, 16.0, 11.0
1169 , 8.0, 6.0, 5.2, 4.6, 4.0
1170 , 3.5, 3.0, 1.4, 0.67, 0.44
1171 , 0.3, 0.18, 0.12, 0.08, 0.056
1172 , 0.04, 0.023, 0.015, 0.011, 0.01
1177 // Find the position
1182 dpos = energy - vxe[pos2++];
1191 // Differentiate between the sampling points
1192 dnde = (vye[pos1] - vye[pos2]) / (vxe[pos2] - vxe[pos1]);
1198 //_____________________________________________________________________________
1199 Double_t IntSpecGeant(Double_t *x, Double_t *)
1202 // Integrated spectrum from Geant3
1205 const Int_t npts = 83;
1206 Double_t arre[npts] = { 2.421257, 2.483278, 2.534301, 2.592230
1207 , 2.672067, 2.813299, 3.015059, 3.216819
1208 , 3.418579, 3.620338, 3.868209, 3.920198
1209 , 3.978284, 4.063923, 4.186264, 4.308605
1210 , 4.430946, 4.553288, 4.724261, 4.837736
1211 , 4.999842, 5.161949, 5.324056, 5.486163
1212 , 5.679688, 5.752998, 5.857728, 5.962457
1213 , 6.067185, 6.171914, 6.315653, 6.393674
1214 , 6.471694, 6.539689, 6.597658, 6.655627
1215 , 6.710957, 6.763648, 6.816338, 6.876198
1216 , 6.943227, 7.010257, 7.106285, 7.252151
1217 , 7.460531, 7.668911, 7.877290, 8.085670
1218 , 8.302979, 8.353585, 8.413120, 8.483500
1219 , 8.541030, 8.592857, 8.668865, 8.820485
1220 , 9.037086, 9.253686, 9.470286, 9.686887
1221 , 9.930838, 9.994655, 10.085822, 10.176990
1222 , 10.268158, 10.359325, 10.503614, 10.627565
1223 , 10.804637, 10.981709, 11.158781, 11.335854
1224 , 11.593397, 11.781165, 12.049404, 12.317644
1225 , 12.585884, 12.854123, 14.278421, 16.975889
1226 , 20.829416, 24.682943, 28.536469 };
1228 Double_t arrdnde[npts] = { 10.960000, 10.960000, 10.359500, 9.811340
1229 , 9.1601500, 8.206670, 6.919630, 5.655430
1230 , 4.6221300, 3.777610, 3.019560, 2.591950
1231 , 2.5414600, 2.712920, 3.327460, 4.928240
1232 , 7.6185300, 10.966700, 12.225800, 8.094750
1233 , 3.3586900, 1.553650, 1.209600, 1.263840
1234 , 1.3241100, 1.312140, 1.255130, 1.165770
1235 , 1.0594500, 0.945450, 0.813231, 0.699837
1236 , 0.6235580, 2.260990, 2.968350, 2.240320
1237 , 1.7988300, 1.553300, 1.432070, 1.535520
1238 , 1.4429900, 1.247990, 1.050750, 0.829549
1239 , 0.5900280, 0.395897, 0.268741, 0.185320
1240 , 0.1292120, 0.103545, 0.0949525, 0.101535
1241 , 0.1276380, 0.134216, 0.123816, 0.104557
1242 , 0.0751843, 0.0521745, 0.0373546, 0.0275391
1243 , 0.0204713, 0.0169234, 0.0154552, 0.0139194
1244 , 0.0125592, 0.0113638, 0.0107354, 0.0102137
1245 , 0.00845984, 0.00683338, 0.00556836, 0.00456874
1246 , 0.0036227, 0.00285991, 0.00226664, 0.00172234
1247 , 0.00131226, 0.00100284, 0.000465492, 7.26607e-05
1248 , 3.63304e-06, 0.0000000, 0.0000000 };
1251 Double_t energy = x[0];
1253 if (energy >= arre[npts-1]) {
1257 for (i = 0; i < npts; i++) {
1258 if (energy < arre[i]) {
1264 AliErrorGeneral("AliTRDv1::IntSpecGeant","Given energy value is too small or zero");