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,v 1.6 2009/04/02 02:43:42 kurasige Exp $ |
---|
28 | // GEANT4 tag $Name: geant4-09-03-beta-cand-01 $ |
---|
29 | // |
---|
30 | // |
---|
31 | // -------------------------------------------------------------- |
---|
32 | // GEANT 4 class implementation file/ History: |
---|
33 | // 5 Oct. 2002, H.Kuirashige : Structure created based on object model |
---|
34 | // -------------------------------------------------------------- |
---|
35 | |
---|
36 | #include "G4RToEConvForElectron.hh" |
---|
37 | #include "G4ParticleDefinition.hh" |
---|
38 | #include "G4ParticleTable.hh" |
---|
39 | #include "G4Material.hh" |
---|
40 | #include "G4PhysicsLogVector.hh" |
---|
41 | |
---|
42 | #include "G4ios.hh" |
---|
43 | |
---|
44 | G4RToEConvForElectron::G4RToEConvForElectron() : G4VRangeToEnergyConverter() |
---|
45 | { |
---|
46 | theParticle = G4ParticleTable::GetParticleTable()->FindParticle("e-"); |
---|
47 | if (theParticle ==0) { |
---|
48 | #ifdef G4VERBOSE |
---|
49 | if (GetVerboseLevel()>0) { |
---|
50 | G4cout << " G4RToEConvForElectron::G4RToEConvForElectron() "; |
---|
51 | G4cout << " Electron is not defined !!" << G4endl; |
---|
52 | } |
---|
53 | #endif |
---|
54 | } |
---|
55 | } |
---|
56 | |
---|
57 | G4RToEConvForElectron::~G4RToEConvForElectron() |
---|
58 | { |
---|
59 | } |
---|
60 | |
---|
61 | |
---|
62 | // ********************************************************************** |
---|
63 | // ************************* ComputeLoss ******************************** |
---|
64 | // ********************************************************************** |
---|
65 | G4double G4RToEConvForElectron::ComputeLoss(G4double AtomicNumber, |
---|
66 | G4double KineticEnergy) const |
---|
67 | { |
---|
68 | static G4double Z; |
---|
69 | static G4double taul, ionpot, ionpotlog; |
---|
70 | const G4double cbr1=0.02, cbr2=-5.7e-5, cbr3=1., cbr4=0.072; |
---|
71 | const G4double Tlow=10.*keV, Thigh=1.*GeV; |
---|
72 | static G4double bremfactor= 0.1 ; |
---|
73 | |
---|
74 | G4double Mass = theParticle->GetPDGMass(); |
---|
75 | // calculate dE/dx for electrons |
---|
76 | if( std::fabs(AtomicNumber-Z)>0.1 ) { |
---|
77 | Z = AtomicNumber; |
---|
78 | taul = Tlow/Mass; |
---|
79 | ionpot = 1.6e-5*MeV*std::exp(0.9*std::log(Z))/Mass; |
---|
80 | ionpotlog = std::log(ionpot); |
---|
81 | } |
---|
82 | |
---|
83 | |
---|
84 | G4double tau = KineticEnergy/Mass; |
---|
85 | G4double dEdx; |
---|
86 | |
---|
87 | if(tau<taul) { |
---|
88 | G4double t1 = taul+1.; |
---|
89 | G4double t2 = taul+2.; |
---|
90 | G4double tsq = taul*taul; |
---|
91 | G4double beta2 = taul*t2/(t1*t1); |
---|
92 | G4double f = 1.-beta2+std::log(tsq/2.) |
---|
93 | +(0.5+0.25*tsq+(1.+2.*taul)*std::log(0.5))/(t1*t1); |
---|
94 | dEdx = (std::log(2.*taul+4.)-2.*ionpotlog+f)/beta2; |
---|
95 | dEdx = twopi_mc2_rcl2*Z*dEdx; |
---|
96 | G4double clow = dEdx*std::sqrt(taul); |
---|
97 | dEdx = clow/std::sqrt(KineticEnergy/Mass); |
---|
98 | |
---|
99 | } else { |
---|
100 | G4double t1 = tau+1.; |
---|
101 | G4double t2 = tau+2.; |
---|
102 | G4double tsq = tau*tau; |
---|
103 | G4double beta2 = tau*t2/(t1*t1); |
---|
104 | G4double f = 1.-beta2+std::log(tsq/2.) |
---|
105 | +(0.5+0.25*tsq+(1.+2.*tau)*std::log(0.5))/(t1*t1); |
---|
106 | dEdx = (std::log(2.*tau+4.)-2.*ionpotlog+f)/beta2; |
---|
107 | dEdx = twopi_mc2_rcl2*Z*dEdx; |
---|
108 | |
---|
109 | // loss from bremsstrahlung follows |
---|
110 | G4double cbrem = (cbr1+cbr2*Z) |
---|
111 | *(cbr3+cbr4*std::log(KineticEnergy/Thigh)); |
---|
112 | cbrem = Z*(Z+1.)*cbrem*tau/beta2; |
---|
113 | |
---|
114 | cbrem *= bremfactor ; |
---|
115 | |
---|
116 | dEdx += twopi_mc2_rcl2*cbrem; |
---|
117 | } |
---|
118 | |
---|
119 | return dEdx; |
---|
120 | } |
---|
121 | |
---|
122 | |
---|
123 | void G4RToEConvForElectron::BuildRangeVector(const G4Material* aMaterial, |
---|
124 | G4double maxEnergy, |
---|
125 | G4double aMass, |
---|
126 | G4PhysicsLogVector* rangeVector) |
---|
127 | { |
---|
128 | // create range vector for a material |
---|
129 | const G4double tlim = 10.*keV; |
---|
130 | const G4int maxnbint = 100; |
---|
131 | |
---|
132 | const G4ElementVector* elementVector = aMaterial->GetElementVector(); |
---|
133 | const G4double* atomicNumDensityVector = aMaterial->GetAtomicNumDensityVector(); |
---|
134 | G4int NumEl = aMaterial->GetNumberOfElements(); |
---|
135 | |
---|
136 | // calculate parameters of the low energy part first |
---|
137 | size_t i; |
---|
138 | G4double loss=0.; |
---|
139 | for (i=0; i<size_t(NumEl); i++) { |
---|
140 | G4bool isOut; |
---|
141 | G4int IndEl = (*elementVector)[i]->GetIndex(); |
---|
142 | loss += atomicNumDensityVector[i]* |
---|
143 | (*theLossTable)[IndEl]->GetValue(tlim,isOut); |
---|
144 | } |
---|
145 | G4double taulim = tlim/aMass; |
---|
146 | G4double clim = std::sqrt(taulim)*loss; |
---|
147 | G4double taumax = maxEnergy/aMass; |
---|
148 | |
---|
149 | // now the range vector can be filled |
---|
150 | for ( i=0; i<size_t(TotBin); i++) { |
---|
151 | G4double LowEdgeEnergy = rangeVector->GetLowEdgeEnergy(i); |
---|
152 | G4double tau = LowEdgeEnergy/aMass; |
---|
153 | |
---|
154 | if ( tau <= taulim ) { |
---|
155 | G4double Value = 2.*aMass*tau*std::sqrt(tau)/(3.*clim); |
---|
156 | rangeVector->PutValue(i,Value); |
---|
157 | } else { |
---|
158 | G4double rangelim = 2.*aMass*taulim*std::sqrt(taulim)/(3.*clim); |
---|
159 | G4double ltaulow = std::log(taulim); |
---|
160 | G4double ltauhigh = std::log(tau); |
---|
161 | G4double ltaumax = std::log(taumax); |
---|
162 | G4int nbin = G4int(maxnbint*(ltauhigh-ltaulow)/(ltaumax-ltaulow)); |
---|
163 | if( nbin < 1 ) nbin = 1; |
---|
164 | G4double Value = RangeLogSimpson( NumEl, elementVector, |
---|
165 | atomicNumDensityVector, aMass, |
---|
166 | ltaulow, ltauhigh, nbin) |
---|
167 | + rangelim; |
---|
168 | rangeVector->PutValue(i,Value); |
---|
169 | } |
---|
170 | } |
---|
171 | } |
---|