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>
38 #include "AliTRDgeometry.h"
39 #include "AliTRDhit.h"
40 #include "AliTRDsim.h"
45 //_____________________________________________________________________________
50 ,fTypeOfStepManager(0)
58 // Default constructor
63 //_____________________________________________________________________________
64 AliTRDv1::AliTRDv1(const char *name, const char *title)
68 ,fTypeOfStepManager(1)
76 // Standard constructor for Transition Radiation Detector version 1
79 SetBufferSize(128000);
83 //_____________________________________________________________________________
84 AliTRDv1::AliTRDv1(const AliTRDv1 &trd)
88 ,fTypeOfStepManager(trd.fTypeOfStepManager)
89 ,fStepSize(trd.fStepSize)
92 ,fTrackLength0(trd.fTrackLength0)
93 ,fPrimaryTrackPid(trd.fPrimaryTrackPid)
99 fDeltaE->Copy(*((AliTRDv1 &) trd).fDeltaE);
100 fDeltaG->Copy(*((AliTRDv1 &) trd).fDeltaG);
101 fTR->Copy(*((AliTRDv1 &) trd).fTR);
105 //_____________________________________________________________________________
106 AliTRDv1::~AliTRDv1()
109 // AliTRDv1 destructor
129 //_____________________________________________________________________________
130 AliTRDv1 &AliTRDv1::operator=(const AliTRDv1 &trd)
133 // Assignment operator
136 if (this != &trd) ((AliTRDv1 &) trd).Copy(*this);
142 //_____________________________________________________________________________
143 void AliTRDv1::Copy(TObject &trd) const
149 ((AliTRDv1 &) trd).fTypeOfStepManager = fTypeOfStepManager;
150 ((AliTRDv1 &) trd).fStepSize = fStepSize;
151 ((AliTRDv1 &) trd).fTRon = fTRon;
152 ((AliTRDv1 &) trd).fTrackLength0 = fTrackLength0;
153 ((AliTRDv1 &) trd).fPrimaryTrackPid = fPrimaryTrackPid;
155 fDeltaE->Copy(*((AliTRDv1 &) trd).fDeltaE);
156 fDeltaG->Copy(*((AliTRDv1 &) trd).fDeltaG);
157 fTR->Copy(*((AliTRDv1 &) trd).fTR);
161 //_____________________________________________________________________________
162 void AliTRDv1::CreateGeometry()
165 // Create the GEANT geometry for the Transition Radiation Detector - Version 1
166 // This version covers the full azimuth.
169 // Check that FRAME is there otherwise we have no place where to put the TRD
170 AliModule* frame = gAlice->GetModule("FRAME");
172 AliError("TRD needs FRAME to be present\n");
176 // Define the chambers
177 AliTRD::CreateGeometry();
181 //_____________________________________________________________________________
182 void AliTRDv1::CreateMaterials()
185 // Create materials for the Transition Radiation Detector version 1
188 AliTRD::CreateMaterials();
192 //_____________________________________________________________________________
193 void AliTRDv1::CreateTRhit(Int_t det)
196 // Creates an electron cluster from a TR photon.
197 // The photon is assumed to be created a the end of the radiator. The
198 // distance after which it deposits its energy takes into account the
199 // absorbtion of the entrance window and of the gas mixture in drift
204 const Int_t kPdgElectron = 11;
207 const Float_t kWion = 23.53;
209 // Maximum number of TR photons per track
210 const Int_t kNTR = 50;
215 // Create TR at the entrance of the chamber
216 if (gMC->IsTrackEntering()) {
218 // Create TR only for electrons
219 Int_t iPdg = gMC->TrackPid();
220 if (TMath::Abs(iPdg) != kPdgElectron) {
228 gMC->TrackMomentum(mom);
229 Float_t pTot = mom.Rho();
230 fTR->CreatePhotons(iPdg,pTot,nTR,eTR);
232 AliFatal(Form("Boundary error: nTR = %d, kNTR = %d",nTR,kNTR));
235 // Loop through the TR photons
236 for (Int_t iTR = 0; iTR < nTR; iTR++) {
238 Float_t energyMeV = eTR[iTR] * 0.001;
239 Float_t energyeV = eTR[iTR] * 1000.0;
240 Float_t absLength = 0.0;
243 // Take the absorbtion in the entrance window into account
244 Double_t muMy = fTR->GetMuMy(energyMeV);
245 sigma = muMy * fFoilDensity;
247 absLength = gRandom->Exp(1.0/sigma);
248 if (absLength < AliTRDgeometry::MyThick()) {
256 // The absorbtion cross sections in the drift gas
257 // Gas-mixture (Xe/CO2)
258 Double_t muXe = fTR->GetMuXe(energyMeV);
259 Double_t muCO = fTR->GetMuCO(energyMeV);
260 sigma = (0.85 * muXe + 0.15 * muCO) * fGasDensity * fTR->GetTemp();
262 // The distance after which the energy of the TR photon
265 absLength = gRandom->Exp(1.0/sigma);
266 if (absLength > (AliTRDgeometry::DrThick()
267 + AliTRDgeometry::AmThick())) {
275 // The position of the absorbtion
277 gMC->TrackPosition(pos);
278 posHit[0] = pos[0] + mom[0] / pTot * absLength;
279 posHit[1] = pos[1] + mom[1] / pTot * absLength;
280 posHit[2] = pos[2] + mom[2] / pTot * absLength;
283 Int_t q = ((Int_t) (energyeV / kWion));
285 // Add the hit to the array. TR photon hits are marked
286 // by negative charge
287 AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber()
299 //_____________________________________________________________________________
300 void AliTRDv1::Init()
303 // Initialise Transition Radiation Detector after geometry has been built.
308 AliDebug(1,"Slow simulator\n");
310 // Switch on TR simulation as default
312 AliInfo("TR simulation off");
315 fTR = new AliTRDsim();
318 // First ionization potential (eV) for the gas mixture (90% Xe + 10% CO2)
319 const Float_t kPoti = 12.1;
320 // Maximum energy (50 keV);
321 const Float_t kEend = 50000.0;
322 // Ermilova distribution for the delta-ray spectrum
323 Float_t poti = TMath::Log(kPoti);
324 Float_t eEnd = TMath::Log(kEend);
326 // Ermilova distribution for the delta-ray spectrum
327 fDeltaE = new TF1("deltae" ,Ermilova ,poti,eEnd,0);
329 // Geant3 distribution for the delta-ray spectrum
330 fDeltaG = new TF1("deltag",IntSpecGeant,2.421257,28.536469,0);
332 AliDebug(1,"+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++");
336 //_____________________________________________________________________________
337 void AliTRDv1::StepManager()
340 // Slow simulator. Every charged track produces electron cluster as hits
341 // along its path across the drift volume.
344 switch (fTypeOfStepManager) {
346 StepManagerErmilova();
352 StepManagerFixedStep();
355 AliWarning("Not a valid Step Manager.");
360 //_____________________________________________________________________________
361 void AliTRDv1::SelectStepManager(Int_t t)
364 // Selects a step manager type:
367 // 2 - Fixed step size
370 fTypeOfStepManager = t;
371 AliInfo(Form("Step Manager type %d was selected",fTypeOfStepManager));
375 //_____________________________________________________________________________
376 void AliTRDv1::StepManagerGeant()
379 // Slow simulator. Every charged track produces electron cluster as hits
380 // along its path across the drift volume. The step size is set acording
381 // to Bethe-Bloch. The energy distribution of the delta electrons follows
382 // a spectrum taken from Geant3.
384 // Version by A. Bercuci
402 Double_t stepSize = 0;
404 Bool_t drRegion = kFALSE;
405 Bool_t amRegion = kFALSE;
408 TString cIdSensDr = "J";
409 TString cIdSensAm = "K";
410 Char_t cIdChamber[3];
418 const Int_t kNplan = AliTRDgeometry::Nplan();
419 const Int_t kNcham = AliTRDgeometry::Ncham();
420 const Int_t kNdetsec = kNplan * kNcham;
422 const Double_t kBig = 1.0e+12; // Infinitely big
423 const Float_t kWion = 23.53; // Ionization energy
424 const Float_t kPTotMaxEl = 0.002; // Maximum momentum for e+ e- g
426 // Minimum energy for the step size adjustment
427 const Float_t kEkinMinStep = 1.0e-5;
428 // energy threshold for production of delta electrons
429 const Float_t kECut = 1.0e4;
430 // Parameters entering the parametrized range for delta electrons
431 const Float_t kRa = 5.37e-4;
432 const Float_t kRb = 0.9815;
433 const Float_t kRc = 3.123e-3;
434 // Gas density -> To be made user adjustable !
435 // [0.85*0.00549+0.15*0.00186 (Xe-CO2 85-15)]
436 const Float_t kRho = 0.004945 ;
438 // Plateau value of the energy-loss for electron in xenon
439 // The averaged value (26/3/99)
440 const Float_t kPlateau = 1.55;
441 // dN1/dx|min for the gas mixture (90% Xe + 10% CO2)
442 const Float_t kPrim = 19.34;
443 // First ionization potential (eV) for the gas mixture (90% Xe + 10% CO2)
444 const Float_t kPoti = 12.1;
446 const Int_t kPdgElectron = 11;
448 // Set the maximum step size to a very large number for all
449 // neutral particles and those outside the driftvolume
450 gMC->SetMaxStep(kBig);
452 // Use only charged tracks
453 if (( gMC->TrackCharge() ) &&
454 (!gMC->IsTrackDisappeared())) {
456 // Inside a sensitive volume?
459 cIdCurrent = gMC->CurrentVolName();
460 if (cIdSensDr == cIdCurrent[1]) {
463 if (cIdSensAm == cIdCurrent[1]) {
466 if (drRegion || amRegion) {
468 // The hit coordinates and charge
469 gMC->TrackPosition(pos);
474 // The sector number (0 - 17)
475 // The numbering goes clockwise and starts at y = 0
476 Float_t phi = kRaddeg*TMath::ATan2(pos[0],pos[1]);
483 sec = ((Int_t) (phi / 20.0));
485 // The plane and chamber number
486 cIdChamber[0] = cIdCurrent[2];
487 cIdChamber[1] = cIdCurrent[3];
488 Int_t idChamber = (atoi(cIdChamber) % kNdetsec);
489 cha = kNcham - ((Int_t) idChamber / kNplan) - 1;
490 pla = ((Int_t) idChamber % kNplan);
492 // Check on selected volumes
493 Int_t addthishit = 1;
498 // The detector number
499 det = fGeometry->GetDetector(pla,cha,sec);
501 // Special hits only in the drift region
504 // Create a track reference at the entrance and
505 // exit of each chamber that contain the
506 // momentum components of the particle
507 if (gMC->IsTrackEntering() ||
508 gMC->IsTrackExiting()) {
509 gMC->TrackMomentum(mom);
510 AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber());
513 if (gMC->IsTrackEntering() &&
514 !gMC->IsNewTrack()) {
515 // determine if hit belong to primary track
516 fPrimaryTrackPid = gAlice->GetMCApp()->GetCurrentTrackNumber();
517 // determine track length when entering the detector
518 fTrackLength0 = gMC->TrackLength();
521 // Create the hits from TR photons
522 if (fTR) CreateTRhit(det);
526 // Calculate the energy of the delta-electrons
527 // modified by Alex Bercuci (A.Bercuci@gsi.de) on 26.01.06
528 // take into account correlation with the underlying GEANT tracking
530 // http://www-linux.gsi.de/~abercuci/Contributions/TRD/index.html
532 // determine the most significant process (last on the processes list)
533 // which caused this hit
534 gMC->StepProcesses(processes);
535 Int_t nofprocesses = processes.GetSize();
541 pid = processes[nofprocesses-1];
544 // generate Edep according to GEANT parametrisation
545 eDelta = TMath::Exp(fDeltaG->GetRandom()) - kPoti;
546 eDelta = TMath::Max(eDelta,0.0);
547 Float_t prRange = 0.0;
548 Float_t range = gMC->TrackLength() - fTrackLength0;
549 // merge GEANT tracker information with locally cooked one
550 if (gAlice->GetMCApp()->GetCurrentTrackNumber() == fPrimaryTrackPid) {
552 if (eDelta >= kECut) {
553 prRange = kRa * eDelta * 0.001
554 * (1.0 - kRb / (1.0 + kRc * eDelta * 0.001)) / kRho;
555 if (prRange >= (3.7 - range)) {
561 if (eDelta < kECut) {
565 prRange = kRa * eDelta * 0.001
566 * (1.0 - kRb / (1.0 + kRc * eDelta * 0.001)) / kRho;
567 if (prRange >= ((AliTRDgeometry::DrThick()
568 + AliTRDgeometry::AmThick()) - range)) {
584 // Generate the electron cluster size
589 qTot = ((Int_t) (eDelta / kWion) + 1);
592 // Create a new dEdx hit
593 AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber()
599 // Calculate the maximum step size for the next tracking step
600 // Produce only one hit if Ekin is below cutoff
601 aMass = gMC->TrackMass();
602 if ((gMC->Etot() - aMass) > kEkinMinStep) {
604 // The energy loss according to Bethe Bloch
605 iPdg = TMath::Abs(gMC->TrackPid());
606 if ((iPdg != kPdgElectron) ||
607 ((iPdg == kPdgElectron) &&
608 (pTot < kPTotMaxEl))) {
609 gMC->TrackMomentum(mom);
611 betaGamma = pTot / aMass;
612 pp = BetheBlochGeant(betaGamma);
613 // Take charge > 1 into account
614 charge = gMC->TrackCharge();
615 if (TMath::Abs(charge) > 1) {
616 pp = pp * charge*charge;
620 // Electrons above 20 Mev/c are at the plateau
621 pp = kPrim * kPlateau;
626 nsteps = gRandom->Poisson(pp);
628 stepSize = 1.0 / nsteps;
629 gMC->SetMaxStep(stepSize);
641 //_____________________________________________________________________________
642 void AliTRDv1::StepManagerErmilova()
645 // Slow simulator. Every charged track produces electron cluster as hits
646 // along its path across the drift volume. The step size is set acording
647 // to Bethe-Bloch. The energy distribution of the delta electrons follows
648 // a spectrum taken from Ermilova et al.
669 Bool_t drRegion = kFALSE;
670 Bool_t amRegion = kFALSE;
673 TString cIdSensDr = "J";
674 TString cIdSensAm = "K";
675 Char_t cIdChamber[3];
681 const Int_t kNplan = AliTRDgeometry::Nplan();
682 const Int_t kNcham = AliTRDgeometry::Ncham();
683 const Int_t kNdetsec = kNplan * kNcham;
685 const Double_t kBig = 1.0e+12; // Infinitely big
686 const Float_t kWion = 23.53; // Ionization energy
687 const Float_t kPTotMaxEl = 0.002; // Maximum momentum for e+ e- g
689 // Minimum energy for the step size adjustment
690 const Float_t kEkinMinStep = 1.0e-5;
692 // Plateau value of the energy-loss for electron in xenon
693 // The averaged value (26/3/99)
694 const Float_t kPlateau = 1.55;
695 // dN1/dx|min for the gas mixture (90% Xe + 10% CO2)
696 const Float_t kPrim = 48.0;
697 // First ionization potential (eV) for the gas mixture (90% Xe + 10% CO2)
698 const Float_t kPoti = 12.1;
700 const Int_t kPdgElectron = 11;
702 // Set the maximum step size to a very large number for all
703 // neutral particles and those outside the driftvolume
704 gMC->SetMaxStep(kBig);
706 // Use only charged tracks
707 if (( gMC->TrackCharge() ) &&
708 (!gMC->IsTrackDisappeared())) {
710 // Inside a sensitive volume?
713 cIdCurrent = gMC->CurrentVolName();
714 if (cIdSensDr == cIdCurrent[1]) {
717 if (cIdSensAm == cIdCurrent[1]) {
720 if (drRegion || amRegion) {
722 // The hit coordinates and charge
723 gMC->TrackPosition(pos);
728 // The sector number (0 - 17)
729 // The numbering goes clockwise and starts at y = 0
730 Float_t phi = kRaddeg*TMath::ATan2(pos[0],pos[1]);
737 sec = ((Int_t) (phi / 20.0));
739 // The plane and chamber number
740 cIdChamber[0] = cIdCurrent[2];
741 cIdChamber[1] = cIdCurrent[3];
742 Int_t idChamber = (atoi(cIdChamber) % kNdetsec);
743 cha = kNcham - ((Int_t) idChamber / kNplan) - 1;
744 pla = ((Int_t) idChamber % kNplan);
746 // Check on selected volumes
747 Int_t addthishit = 1;
752 // The detector number
753 det = fGeometry->GetDetector(pla,cha,sec);
755 // Special hits only in the drift region
758 // Create a track reference at the entrance and
759 // exit of each chamber that contain the
760 // momentum components of the particle
761 if (gMC->IsTrackEntering() ||
762 gMC->IsTrackExiting()) {
763 gMC->TrackMomentum(mom);
764 AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber());
766 // Create the hits from TR photons
773 // Calculate the energy of the delta-electrons
774 eDelta = TMath::Exp(fDeltaE->GetRandom()) - kPoti;
775 eDelta = TMath::Max(eDelta,0.0);
776 // Generate the electron cluster size
781 qTot = ((Int_t) (eDelta / kWion) + 1);
784 // Create a new dEdx hit
786 AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber()
793 AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber()
800 // Calculate the maximum step size for the next tracking step
801 // Produce only one hit if Ekin is below cutoff
802 aMass = gMC->TrackMass();
803 if ((gMC->Etot() - aMass) > kEkinMinStep) {
805 // The energy loss according to Bethe Bloch
806 iPdg = TMath::Abs(gMC->TrackPid());
807 if ((iPdg != kPdgElectron) ||
808 ((iPdg == kPdgElectron) &&
809 (pTot < kPTotMaxEl))) {
810 gMC->TrackMomentum(mom);
812 betaGamma = pTot / aMass;
813 pp = kPrim * BetheBloch(betaGamma);
814 // Take charge > 1 into account
815 charge = gMC->TrackCharge();
816 if (TMath::Abs(charge) > 1) {
817 pp = pp * charge*charge;
821 // Electrons above 20 Mev/c are at the plateau
822 pp = kPrim * kPlateau;
827 gMC->GetRandom()->RndmArray(1,random);
829 while ((random[0] == 1.0) ||
831 stepSize = - TMath::Log(random[0]) / pp;
832 gMC->SetMaxStep(stepSize);
845 //_____________________________________________________________________________
846 void AliTRDv1::StepManagerFixedStep()
849 // Slow simulator. Every charged track produces electron cluster as hits
850 // along its path across the drift volume. The step size is fixed in
851 // this version of the step manager.
863 Bool_t drRegion = kFALSE;
864 Bool_t amRegion = kFALSE;
867 TString cIdSensDr = "J";
868 TString cIdSensAm = "K";
869 Char_t cIdChamber[3];
875 const Int_t kNplan = AliTRDgeometry::Nplan();
876 const Int_t kNcham = AliTRDgeometry::Ncham();
877 const Int_t kNdetsec = kNplan * kNcham;
879 const Double_t kBig = 1.0e+12;
881 const Float_t kWion = 23.53; // Ionization energy
882 const Float_t kEkinMinStep = 1.0e-5; // Minimum energy for the step size adjustment
884 // Set the maximum step size to a very large number for all
885 // neutral particles and those outside the driftvolume
886 gMC->SetMaxStep(kBig);
888 // If not charged track or already stopped or disappeared, just return.
889 if ((!gMC->TrackCharge()) ||
890 gMC->IsTrackDisappeared()) return;
892 // Inside a sensitive volume?
893 cIdCurrent = gMC->CurrentVolName();
895 if (cIdSensDr == cIdCurrent[1]) drRegion = kTRUE;
896 if (cIdSensAm == cIdCurrent[1]) amRegion = kTRUE;
903 // The hit coordinates and charge
904 gMC->TrackPosition(pos);
909 // The sector number (0 - 17)
910 // The numbering goes clockwise and starts at y = 0
911 Float_t phi = kRaddeg*TMath::ATan2(pos[0],pos[1]);
918 sec = ((Int_t) (phi / 20.0));
920 // The plane and chamber number
921 cIdChamber[0] = cIdCurrent[2];
922 cIdChamber[1] = cIdCurrent[3];
923 Int_t idChamber = (atoi(cIdChamber) % kNdetsec);
924 cha = kNcham - ((Int_t) idChamber / kNplan) - 1;
925 pla = ((Int_t) idChamber % kNplan);
927 // Check on selected volumes
928 Int_t addthishit = 1;
934 // The detector number
935 det = fGeometry->GetDetector(pla,cha,sec);
937 // 0: InFlight 1:Entering 2:Exiting
940 // Special hits only in the drift region
943 // Create a track reference at the entrance and exit of each
944 // chamber that contain the momentum components of the particle
946 if (gMC->IsTrackEntering()) {
947 gMC->TrackMomentum(mom);
948 AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber());
951 if (gMC->IsTrackExiting()) {
952 gMC->TrackMomentum(mom);
953 AddTrackReference(gAlice->GetMCApp()->GetCurrentTrackNumber());
957 // Create the hits from TR photons
964 // Calculate the charge according to GEANT Edep
965 // Create a new dEdx hit
966 eDep = TMath::Max(gMC->Edep(),0.0) * 1.0e+09;
967 qTot = (Int_t) (eDep / kWion);
968 AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber()
974 // Set Maximum Step Size
975 // Produce only one hit if Ekin is below cutoff
976 if ((gMC->Etot() - gMC->TrackMass()) < kEkinMinStep) {
979 gMC->SetMaxStep(fStepSize);
983 //_____________________________________________________________________________
984 Double_t AliTRDv1::BetheBloch(Double_t bg)
987 // Parametrization of the Bethe-Bloch-curve
988 // The parametrization is the same as for the TPC and is taken from Lehrhaus.
991 // This parameters have been adjusted to averaged values from GEANT
992 const Double_t kP1 = 7.17960e-02;
993 const Double_t kP2 = 8.54196;
994 const Double_t kP3 = 1.38065e-06;
995 const Double_t kP4 = 5.30972;
996 const Double_t kP5 = 2.83798;
998 // Lower cutoff of the Bethe-Bloch-curve to limit step sizes
999 const Double_t kBgMin = 0.8;
1000 const Double_t kBBMax = 6.83298;
1003 Double_t yy = bg / TMath::Sqrt(1.0 + bg*bg);
1004 Double_t aa = TMath::Power(yy,kP4);
1005 Double_t bb = TMath::Power((1.0/bg),kP5);
1006 bb = TMath::Log(kP3 + bb);
1007 return ((kP2 - aa - bb) * kP1 / aa);
1015 //_____________________________________________________________________________
1016 Double_t AliTRDv1::BetheBlochGeant(Double_t bg)
1019 // Return dN/dx (number of primary collisions per centimeter)
1020 // for given beta*gamma factor.
1022 // Implemented by K.Oyama according to GEANT 3 parametrization shown in
1023 // A.Andronic's webpage: http://www-alice.gsi.de/trd/papers/dedx/dedx.html
1024 // This must be used as a set with IntSpecGeant.
1029 Double_t arrG[20] = { 1.100000, 1.200000, 1.300000, 1.500000
1030 , 1.800000, 2.000000, 2.500000, 3.000000
1031 , 4.000000, 7.000000, 10.000000, 20.000000
1032 , 40.000000, 70.000000, 100.000000, 300.000000
1033 , 600.000000, 1000.000000, 3000.000000, 10000.000000 };
1035 Double_t arrNC[20] = { 75.009056, 45.508083, 35.299252, 27.116327
1036 , 22.734999, 21.411915, 19.934095, 19.449375
1037 , 19.344431, 20.185553, 21.027925, 22.912676
1038 , 24.933352, 26.504053, 27.387468, 29.566597
1039 , 30.353779, 30.787134, 31.129285, 31.157350 };
1041 // Betagamma to gamma
1042 Double_t g = TMath::Sqrt(1.0 + bg*bg);
1044 // Find the index just before the point we need.
1045 for (i = 0; i < 18; i++) {
1046 if ((arrG[i] < g) &&
1052 // Simple interpolation.
1053 Double_t pp = ((arrNC[i+1] - arrNC[i]) / (arrG[i+1] - arrG[i]))
1054 * (g - arrG[i]) + arrNC[i];
1060 //_____________________________________________________________________________
1061 Double_t Ermilova(Double_t *x, Double_t *)
1064 // Calculates the delta-ray energy distribution according to Ermilova.
1065 // Logarithmic scale !
1075 const Int_t kNv = 31;
1077 Float_t vxe[kNv] = { 2.3026, 2.9957, 3.4012, 3.6889, 3.9120
1078 , 4.0943, 4.2485, 4.3820, 4.4998, 4.6052
1079 , 4.7005, 5.0752, 5.2983, 5.7038, 5.9915
1080 , 6.2146, 6.5221, 6.9078, 7.3132, 7.6009
1081 , 8.0064, 8.5172, 8.6995, 8.9872, 9.2103
1082 , 9.4727, 9.9035, 10.3735, 10.5966, 10.8198
1085 Float_t vye[kNv] = { 80.0, 31.0, 23.3, 21.1, 21.0
1086 , 20.9, 20.8, 20.0, 16.0, 11.0
1087 , 8.0, 6.0, 5.2, 4.6, 4.0
1088 , 3.5, 3.0, 1.4, 0.67, 0.44
1089 , 0.3, 0.18, 0.12, 0.08, 0.056
1090 , 0.04, 0.023, 0.015, 0.011, 0.01
1095 // Find the position
1100 dpos = energy - vxe[pos2++];
1109 // Differentiate between the sampling points
1110 dnde = (vye[pos1] - vye[pos2]) / (vxe[pos2] - vxe[pos1]);
1116 //_____________________________________________________________________________
1117 Double_t IntSpecGeant(Double_t *x, Double_t *)
1120 // Integrated spectrum from Geant3
1123 const Int_t npts = 83;
1124 Double_t arre[npts] = { 2.421257, 2.483278, 2.534301, 2.592230
1125 , 2.672067, 2.813299, 3.015059, 3.216819
1126 , 3.418579, 3.620338, 3.868209, 3.920198
1127 , 3.978284, 4.063923, 4.186264, 4.308605
1128 , 4.430946, 4.553288, 4.724261, 4.837736
1129 , 4.999842, 5.161949, 5.324056, 5.486163
1130 , 5.679688, 5.752998, 5.857728, 5.962457
1131 , 6.067185, 6.171914, 6.315653, 6.393674
1132 , 6.471694, 6.539689, 6.597658, 6.655627
1133 , 6.710957, 6.763648, 6.816338, 6.876198
1134 , 6.943227, 7.010257, 7.106285, 7.252151
1135 , 7.460531, 7.668911, 7.877290, 8.085670
1136 , 8.302979, 8.353585, 8.413120, 8.483500
1137 , 8.541030, 8.592857, 8.668865, 8.820485
1138 , 9.037086, 9.253686, 9.470286, 9.686887
1139 , 9.930838, 9.994655, 10.085822, 10.176990
1140 , 10.268158, 10.359325, 10.503614, 10.627565
1141 , 10.804637, 10.981709, 11.158781, 11.335854
1142 , 11.593397, 11.781165, 12.049404, 12.317644
1143 , 12.585884, 12.854123, 14.278421, 16.975889
1144 , 20.829416, 24.682943, 28.536469 };
1146 Double_t arrdnde[npts] = { 10.960000, 10.960000, 10.359500, 9.811340
1147 , 9.1601500, 8.206670, 6.919630, 5.655430
1148 , 4.6221300, 3.777610, 3.019560, 2.591950
1149 , 2.5414600, 2.712920, 3.327460, 4.928240
1150 , 7.6185300, 10.966700, 12.225800, 8.094750
1151 , 3.3586900, 1.553650, 1.209600, 1.263840
1152 , 1.3241100, 1.312140, 1.255130, 1.165770
1153 , 1.0594500, 0.945450, 0.813231, 0.699837
1154 , 0.6235580, 2.260990, 2.968350, 2.240320
1155 , 1.7988300, 1.553300, 1.432070, 1.535520
1156 , 1.4429900, 1.247990, 1.050750, 0.829549
1157 , 0.5900280, 0.395897, 0.268741, 0.185320
1158 , 0.1292120, 0.103545, 0.0949525, 0.101535
1159 , 0.1276380, 0.134216, 0.123816, 0.104557
1160 , 0.0751843, 0.0521745, 0.0373546, 0.0275391
1161 , 0.0204713, 0.0169234, 0.0154552, 0.0139194
1162 , 0.0125592, 0.0113638, 0.0107354, 0.0102137
1163 , 0.00845984, 0.00683338, 0.00556836, 0.00456874
1164 , 0.0036227, 0.00285991, 0.00226664, 0.00172234
1165 , 0.00131226, 0.00100284, 0.000465492, 7.26607e-05
1166 , 3.63304e-06, 0.0000000, 0.0000000 };
1169 Double_t energy = x[0];
1171 for (i = 0; i < npts; i++) {
1172 if (energy < arre[i]) {
1178 AliErrorGeneral("AliTRDv1::IntSpecGeant","Given energy value is too small or zero");