Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4RToEConvForElectron.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4RToEConvForElectron.cc 70745 2013-06-05 10:54:00Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------
31 // GEANT 4 class implementation file/ History:
32 // 5 Oct. 2002, H.Kuirashige : Structure created based on object model
33 // --------------------------------------------------------------
34 
35 #include "G4RToEConvForElectron.hh"
36 #include "G4ParticleDefinition.hh"
37 #include "G4ParticleTable.hh"
38 #include "G4Material.hh"
39 #include "G4PhysicsLogVector.hh"
40 
41 #include "G4ios.hh"
42 #include "G4PhysicalConstants.hh"
43 #include "G4SystemOfUnits.hh"
44 
47  Mass(0.0),
48  Z(-1.),
49  taul(0.0),
50  ionpot(0.0),
51  ionpotlog(-1.0e-10),
52  bremfactor(0.1)
53 {
55  if (theParticle ==0) {
56 #ifdef G4VERBOSE
57  if (GetVerboseLevel()>0) {
58  G4cout << " G4RToEConvForElectron::G4RToEConvForElectron() ";
59  G4cout << " Electron is not defined !!" << G4endl;
60  }
61 #endif
62  } else {
64  }
65 }
66 
68 {
69 }
70 
71 
72 // **********************************************************************
73 // ************************* ComputeLoss ********************************
74 // **********************************************************************
76  G4double KineticEnergy)
77 {
78  const G4double cbr1=0.02, cbr2=-5.7e-5, cbr3=1., cbr4=0.072;
79  const G4double Tlow=10.*keV, Thigh=1.*GeV;
80 
81  // calculate dE/dx for electrons
82  if( std::fabs(AtomicNumber-Z)>0.1 ) {
83  Z = AtomicNumber;
84  taul = Tlow/Mass;
85  ionpot = 1.6e-5*MeV*std::exp(0.9*std::log(Z))/Mass;
86  ionpotlog = std::log(ionpot);
87  }
88 
89 
90  G4double tau = KineticEnergy/Mass;
91  G4double dEdx;
92 
93  if(tau<taul) {
94  G4double t1 = taul+1.;
95  G4double t2 = taul+2.;
96  G4double tsq = taul*taul;
97  G4double beta2 = taul*t2/(t1*t1);
98  G4double f = 1.-beta2+std::log(tsq/2.)
99  +(0.5+0.25*tsq+(1.+2.*taul)*std::log(0.5))/(t1*t1);
100  dEdx = (std::log(2.*taul+4.)-2.*ionpotlog+f)/beta2;
101  dEdx = twopi_mc2_rcl2*Z*dEdx;
102  G4double clow = dEdx*std::sqrt(taul);
103  dEdx = clow/std::sqrt(KineticEnergy/Mass);
104 
105  } else {
106  G4double t1 = tau+1.;
107  G4double t2 = tau+2.;
108  G4double tsq = tau*tau;
109  G4double beta2 = tau*t2/(t1*t1);
110  G4double f = 1.-beta2+std::log(tsq/2.)
111  +(0.5+0.25*tsq+(1.+2.*tau)*std::log(0.5))/(t1*t1);
112  dEdx = (std::log(2.*tau+4.)-2.*ionpotlog+f)/beta2;
113  dEdx = twopi_mc2_rcl2*Z*dEdx;
114 
115  // loss from bremsstrahlung follows
116  G4double cbrem = (cbr1+cbr2*Z)
117  *(cbr3+cbr4*std::log(KineticEnergy/Thigh));
118  cbrem = Z*(Z+1.)*cbrem*tau/beta2;
119 
120  cbrem *= bremfactor ;
121 
122  dEdx += twopi_mc2_rcl2*cbrem;
123  }
124 
125  return dEdx;
126 }
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
virtual G4double ComputeLoss(G4double AtomicNumber, G4double KineticEnergy)
G4GLOB_DLL std::ostream G4cout
G4double GetPDGMass() const
static G4ParticleTable * GetParticleTable()
tuple t1
Definition: plottest35.py:33
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
const G4ParticleDefinition * theParticle