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 | // $Id: G4PenelopeComptonModel.cc,v 1.2 2008/12/04 14:11:21 pandola Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-02 $ |
---|
28 | // |
---|
29 | // Author: Luciano Pandola |
---|
30 | // |
---|
31 | // History: |
---|
32 | // -------- |
---|
33 | // 02 Oct 2008 L Pandola Migration from process to model |
---|
34 | // 28 Oct 2008 L Pandola Treat the database data from Penelope according to the |
---|
35 | // original model, namely merging levels below 15 eV in |
---|
36 | // a single one. Still, it is not fully compliant with the |
---|
37 | // original Penelope model, because plasma excitation is not |
---|
38 | // considered. |
---|
39 | // 22 Nov 2008 L Pandola Make unit of measurements explicit for binding energies |
---|
40 | // that are read from the external files. |
---|
41 | // 24 Nov 2008 L Pandola Find a cleaner way to delete vectors. |
---|
42 | // |
---|
43 | |
---|
44 | #include "G4PenelopeComptonModel.hh" |
---|
45 | #include "G4ParticleDefinition.hh" |
---|
46 | #include "G4MaterialCutsCouple.hh" |
---|
47 | #include "G4ProductionCutsTable.hh" |
---|
48 | #include "G4DynamicParticle.hh" |
---|
49 | #include "G4VEMDataSet.hh" |
---|
50 | #include "G4PhysicsTable.hh" |
---|
51 | #include "G4ElementTable.hh" |
---|
52 | #include "G4Element.hh" |
---|
53 | #include "G4PhysicsLogVector.hh" |
---|
54 | #include "G4PenelopeIntegrator.hh" |
---|
55 | #include "G4AtomicTransitionManager.hh" |
---|
56 | #include "G4AtomicDeexcitation.hh" |
---|
57 | #include "G4AtomicShell.hh" |
---|
58 | #include "G4Gamma.hh" |
---|
59 | #include "G4Electron.hh" |
---|
60 | |
---|
61 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
62 | |
---|
63 | |
---|
64 | G4PenelopeComptonModel::G4PenelopeComptonModel(const G4ParticleDefinition*, |
---|
65 | const G4String& nam) |
---|
66 | :G4VEmModel(nam),ionizationEnergy(0),hartreeFunction(0), |
---|
67 | occupationNumber(0),isInitialised(false) |
---|
68 | { |
---|
69 | fIntrinsicLowEnergyLimit = 100.0*eV; |
---|
70 | fIntrinsicHighEnergyLimit = 100.0*GeV; |
---|
71 | SetLowEnergyLimit(fIntrinsicLowEnergyLimit); |
---|
72 | SetHighEnergyLimit(fIntrinsicHighEnergyLimit); |
---|
73 | // |
---|
74 | energyForIntegration = 0.0; |
---|
75 | ZForIntegration = 1; |
---|
76 | |
---|
77 | fUseAtomicDeexcitation = true; |
---|
78 | verboseLevel= 0; |
---|
79 | // Verbosity scale: |
---|
80 | // 0 = nothing |
---|
81 | // 1 = warning for energy non-conservation |
---|
82 | // 2 = details of energy budget |
---|
83 | // 3 = calculation of cross sections, file openings, sampling of atoms |
---|
84 | // 4 = entering in methods |
---|
85 | |
---|
86 | //These vectors do not change when materials or cut change. |
---|
87 | //Therefore I can read it at the constructor |
---|
88 | ionizationEnergy = new std::map<G4int,G4DataVector*>; |
---|
89 | hartreeFunction = new std::map<G4int,G4DataVector*>; |
---|
90 | occupationNumber = new std::map<G4int,G4DataVector*>; |
---|
91 | |
---|
92 | ReadData(); //Read data from file |
---|
93 | |
---|
94 | } |
---|
95 | |
---|
96 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
97 | |
---|
98 | G4PenelopeComptonModel::~G4PenelopeComptonModel() |
---|
99 | { |
---|
100 | std::map <G4int,G4DataVector*>::iterator i; |
---|
101 | for (i=ionizationEnergy->begin();i != ionizationEnergy->end();i++) |
---|
102 | if (i->second) delete i->second; |
---|
103 | for (i=hartreeFunction->begin();i != hartreeFunction->end();i++) |
---|
104 | if (i->second) delete i->second; |
---|
105 | for (i=occupationNumber->begin();i != occupationNumber->end();i++) |
---|
106 | if (i->second) delete i->second; |
---|
107 | |
---|
108 | |
---|
109 | if (ionizationEnergy) |
---|
110 | delete ionizationEnergy; |
---|
111 | if (hartreeFunction) |
---|
112 | delete hartreeFunction; |
---|
113 | if (occupationNumber) |
---|
114 | delete occupationNumber; |
---|
115 | } |
---|
116 | |
---|
117 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
118 | |
---|
119 | void G4PenelopeComptonModel::Initialise(const G4ParticleDefinition* particle, |
---|
120 | const G4DataVector& cuts) |
---|
121 | { |
---|
122 | if (verboseLevel > 3) |
---|
123 | G4cout << "Calling G4PenelopeComptonModel::Initialise()" << G4endl; |
---|
124 | |
---|
125 | InitialiseElementSelectors(particle,cuts); |
---|
126 | if (LowEnergyLimit() < fIntrinsicLowEnergyLimit) |
---|
127 | { |
---|
128 | G4cout << "G4PenelopeComptonModel: low energy limit increased from " << |
---|
129 | LowEnergyLimit()/eV << " eV to " << fIntrinsicLowEnergyLimit/eV << " eV" << G4endl; |
---|
130 | SetLowEnergyLimit(fIntrinsicLowEnergyLimit); |
---|
131 | } |
---|
132 | if (HighEnergyLimit() > fIntrinsicHighEnergyLimit) |
---|
133 | { |
---|
134 | G4cout << "G4PenelopeComptonModel: high energy limit decreased from " << |
---|
135 | HighEnergyLimit()/GeV << " GeV to " << fIntrinsicHighEnergyLimit/GeV << " GeV" << G4endl; |
---|
136 | SetHighEnergyLimit(fIntrinsicHighEnergyLimit); |
---|
137 | } |
---|
138 | |
---|
139 | G4cout << "Penelope Compton model is initialized " << G4endl |
---|
140 | << "Energy range: " |
---|
141 | << LowEnergyLimit() / keV << " keV - " |
---|
142 | << HighEnergyLimit() / GeV << " GeV" |
---|
143 | << G4endl; |
---|
144 | |
---|
145 | if(isInitialised) return; |
---|
146 | |
---|
147 | if(pParticleChange) |
---|
148 | fParticleChange = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange); |
---|
149 | else |
---|
150 | fParticleChange = new G4ParticleChangeForGamma(); |
---|
151 | isInitialised = true; |
---|
152 | } |
---|
153 | |
---|
154 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
155 | |
---|
156 | G4double G4PenelopeComptonModel::ComputeCrossSectionPerAtom( |
---|
157 | const G4ParticleDefinition*, |
---|
158 | G4double energy, |
---|
159 | G4double Z, G4double, |
---|
160 | G4double, G4double) |
---|
161 | { |
---|
162 | // Penelope model to calculate the Compton scattering cross section: |
---|
163 | // D. Brusa et al., Nucl. Instrum. Meth. A 379 (1996) 167 |
---|
164 | // |
---|
165 | // The cross section for Compton scattering is calculated according to the Klein-Nishina |
---|
166 | // formula for energy > 5 MeV. |
---|
167 | // For E < 5 MeV it is used a parametrization for the differential cross-section dSigma/dOmega, |
---|
168 | // which is integrated numerically in cos(theta), G4PenelopeComptonModel::DifferentialCrossSection(). |
---|
169 | // The parametrization includes the J(p) |
---|
170 | // distribution profiles for the atomic shells, that are tabulated from Hartree-Fock calculations |
---|
171 | // from F. Biggs et al., At. Data Nucl. Data Tables 16 (1975) 201 |
---|
172 | // |
---|
173 | if (verboseLevel > 3) |
---|
174 | G4cout << "Calling ComputeCrossSectionPerAtom() of G4PenelopeComptonModel" << G4endl; |
---|
175 | |
---|
176 | G4int iZ = (G4int) Z; |
---|
177 | G4double cs=0.0; |
---|
178 | energyForIntegration=energy; |
---|
179 | ZForIntegration = iZ; |
---|
180 | if (energy< 5*MeV) |
---|
181 | { |
---|
182 | // numerically integrate differential cross section dSigma/dOmega |
---|
183 | G4PenelopeIntegrator<G4PenelopeComptonModel,G4double (G4PenelopeComptonModel::*)(G4double)> |
---|
184 | theIntegrator; |
---|
185 | cs = theIntegrator.Calculate(this,&G4PenelopeComptonModel::DifferentialCrossSection,-1.0,1.0,1e-05); |
---|
186 | } |
---|
187 | else |
---|
188 | { |
---|
189 | // use Klein-Nishina formula |
---|
190 | G4double ki=energy/electron_mass_c2; |
---|
191 | G4double ki3=ki*ki; |
---|
192 | G4double ki2=1.0+2*ki; |
---|
193 | G4double ki1=ki3-ki2-1.0; |
---|
194 | G4double t0=1.0/(ki2); |
---|
195 | G4double csl = 0.5*ki3*t0*t0+ki2*t0+ki1*std::log(t0)-(1.0/t0); |
---|
196 | G4int nosc = occupationNumber->find(iZ)->second->size(); |
---|
197 | for (G4int i=0;i<nosc;i++) |
---|
198 | { |
---|
199 | G4double ionEnergy = (*(ionizationEnergy->find(iZ)->second))[i]; |
---|
200 | G4double tau=(energy-ionEnergy)/energy; |
---|
201 | if (tau > t0) |
---|
202 | { |
---|
203 | G4double csu = 0.5*ki3*tau*tau+ki2*tau+ki1*std::log(tau)-(1.0/tau); |
---|
204 | G4int f = (G4int) (*(occupationNumber->find(iZ)->second))[i]; |
---|
205 | cs = cs + f*(csu-csl); |
---|
206 | } |
---|
207 | } |
---|
208 | cs=pi*classic_electr_radius*classic_electr_radius*cs/(ki*ki3); |
---|
209 | } |
---|
210 | |
---|
211 | if (verboseLevel > 2) |
---|
212 | G4cout << "Compton cross Section at " << energy/keV << " keV for Z=" << Z << |
---|
213 | " = " << cs/barn << " barn" << G4endl; |
---|
214 | return cs; |
---|
215 | } |
---|
216 | |
---|
217 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
218 | |
---|
219 | void G4PenelopeComptonModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect, |
---|
220 | const G4MaterialCutsCouple* couple, |
---|
221 | const G4DynamicParticle* aDynamicGamma, |
---|
222 | G4double, |
---|
223 | G4double) |
---|
224 | { |
---|
225 | |
---|
226 | // Penelope model to sample the Compton scattering final state. |
---|
227 | // D. Brusa et al., Nucl. Instrum. Meth. A 379 (1996) 167 |
---|
228 | // The model determines also the original shell from which the electron is expelled, |
---|
229 | // in order to produce fluorescence de-excitation (from G4DeexcitationManager) |
---|
230 | // |
---|
231 | // The final state for Compton scattering is calculated according to the Klein-Nishina |
---|
232 | // formula for energy > 5 MeV. In this case, the Doppler broadening is negligible and |
---|
233 | // one can assume that the target electron is at rest. |
---|
234 | // For E < 5 MeV it is used the parametrization for the differential cross-section dSigma/dOmega, |
---|
235 | // to sample the scattering angle and the energy of the emerging electron, which is |
---|
236 | // G4PenelopeComptonModel::DifferentialCrossSection(). The rejection method is |
---|
237 | // used to sample cos(theta). The efficiency increases monotonically with photon energy and is |
---|
238 | // nearly independent on the Z; typical values are 35%, 80% and 95% for 1 keV, 1 MeV and 10 MeV, |
---|
239 | // respectively. |
---|
240 | // The parametrization includes the J(p) distribution profiles for the atomic shells, that are |
---|
241 | // tabulated |
---|
242 | // from Hartree-Fock calculations from F. Biggs et al., At. Data Nucl. Data Tables 16 (1975) 201. |
---|
243 | // Doppler broadening is included. |
---|
244 | // |
---|
245 | |
---|
246 | if (verboseLevel > 3) |
---|
247 | G4cout << "Calling SampleSecondaries() of G4PenelopeComptonModel" << G4endl; |
---|
248 | |
---|
249 | G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy(); |
---|
250 | |
---|
251 | if (photonEnergy0 <= LowEnergyLimit()) |
---|
252 | { |
---|
253 | fParticleChange->ProposeTrackStatus(fStopAndKill); |
---|
254 | fParticleChange->SetProposedKineticEnergy(0.); |
---|
255 | fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0); |
---|
256 | return ; |
---|
257 | } |
---|
258 | |
---|
259 | G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection(); |
---|
260 | |
---|
261 | // Select randomly one element in the current material |
---|
262 | if (verboseLevel > 2) |
---|
263 | G4cout << "Going to select element in " << couple->GetMaterial()->GetName() << G4endl; |
---|
264 | // atom can be selected effitiantly if element selectors are initialised |
---|
265 | const G4Element* anElement = SelectRandomAtom(couple,G4Gamma::GammaDefinition(),photonEnergy0); |
---|
266 | G4int Z = (G4int) anElement->GetZ(); |
---|
267 | if (verboseLevel > 2) |
---|
268 | G4cout << "Selected " << anElement->GetName() << G4endl; |
---|
269 | |
---|
270 | const G4int nmax = 64; |
---|
271 | G4double rn[nmax],pac[nmax]; |
---|
272 | |
---|
273 | G4double ki,ki1,ki2,ki3,taumin,a1,a2; |
---|
274 | G4double tau,TST; |
---|
275 | G4double S=0.0; |
---|
276 | G4double epsilon,cosTheta; |
---|
277 | G4double harFunc = 0.0; |
---|
278 | G4int occupNb= 0; |
---|
279 | G4double ionEnergy=0.0; |
---|
280 | G4int nosc = occupationNumber->find(Z)->second->size(); |
---|
281 | G4int iosc = nosc; |
---|
282 | ki = photonEnergy0/electron_mass_c2; |
---|
283 | ki2 = 2*ki+1.0; |
---|
284 | ki3 = ki*ki; |
---|
285 | ki1 = ki3-ki2-1.0; |
---|
286 | taumin = 1.0/ki2; |
---|
287 | a1 = std::log(ki2); |
---|
288 | a2 = a1+2.0*ki*(1.0+ki)/(ki2*ki2); |
---|
289 | //If the incoming photon is above 5 MeV, the quicker approach based on the |
---|
290 | //pure Klein-Nishina formula is used |
---|
291 | if (photonEnergy0 > 5*MeV) |
---|
292 | { |
---|
293 | do{ |
---|
294 | do{ |
---|
295 | if ((a2*G4UniformRand()) < a1) |
---|
296 | { |
---|
297 | tau = std::pow(taumin,G4UniformRand()); |
---|
298 | } |
---|
299 | else |
---|
300 | { |
---|
301 | tau = std::sqrt(1.0+G4UniformRand()*(taumin*taumin-1.0)); |
---|
302 | } |
---|
303 | //rejection function |
---|
304 | TST = (1+tau*(ki1+tau*(ki2+tau*ki3)))/(ki3*tau*(1.0+tau*tau)); |
---|
305 | }while (G4UniformRand()> TST); |
---|
306 | epsilon=tau; |
---|
307 | cosTheta = 1.0 - (1.0-tau)/(ki*tau); |
---|
308 | //Target shell electrons |
---|
309 | TST = Z*G4UniformRand(); |
---|
310 | iosc = nosc; |
---|
311 | S=0.0; |
---|
312 | for (G4int j=0;j<nosc;j++) |
---|
313 | { |
---|
314 | occupNb = (G4int) (*(occupationNumber->find(Z)->second))[j]; |
---|
315 | S = S + occupNb; |
---|
316 | if (S > TST) iosc = j; |
---|
317 | if (S > TST) break; |
---|
318 | } |
---|
319 | ionEnergy = (*(ionizationEnergy->find(Z)->second))[iosc]; |
---|
320 | }while((epsilon*photonEnergy0-photonEnergy0+ionEnergy) >0); |
---|
321 | } |
---|
322 | else //photonEnergy0<5 MeV |
---|
323 | { |
---|
324 | //Incoherent scattering function for theta=PI |
---|
325 | G4double s0=0.0; |
---|
326 | G4double pzomc=0.0,rni=0.0; |
---|
327 | G4double aux=0.0; |
---|
328 | for (G4int i=0;i<nosc;i++){ |
---|
329 | ionEnergy = (*(ionizationEnergy->find(Z)->second))[i]; |
---|
330 | if (photonEnergy0 > ionEnergy) |
---|
331 | { |
---|
332 | G4double aux = photonEnergy0*(photonEnergy0-ionEnergy)*2.0; |
---|
333 | harFunc = (*(hartreeFunction->find(Z)->second))[i]/fine_structure_const; |
---|
334 | occupNb = (G4int) (*(occupationNumber->find(Z)->second))[i]; |
---|
335 | pzomc = harFunc*(aux-electron_mass_c2*ionEnergy)/ |
---|
336 | (electron_mass_c2*std::sqrt(2.0*aux+ionEnergy*ionEnergy)); |
---|
337 | if (pzomc > 0) |
---|
338 | { |
---|
339 | rni = 1.0-0.5*std::exp(0.5-(std::sqrt(0.5)+std::sqrt(2.0)*pzomc)* |
---|
340 | (std::sqrt(0.5)+std::sqrt(2.0)*pzomc)); |
---|
341 | } |
---|
342 | else |
---|
343 | { |
---|
344 | rni = 0.5*std::exp(0.5-(std::sqrt(0.5)-std::sqrt(2.0)*pzomc)* |
---|
345 | (std::sqrt(0.5)-std::sqrt(2.0)*pzomc)); |
---|
346 | } |
---|
347 | s0 = s0 + occupNb*rni; |
---|
348 | } |
---|
349 | } |
---|
350 | |
---|
351 | //Sampling tau |
---|
352 | G4double cdt1; |
---|
353 | do |
---|
354 | { |
---|
355 | if ((G4UniformRand()*a2) < a1) |
---|
356 | { |
---|
357 | tau = std::pow(taumin,G4UniformRand()); |
---|
358 | } |
---|
359 | else |
---|
360 | { |
---|
361 | tau = std::sqrt(1.0+G4UniformRand()*(taumin*taumin-1.0)); |
---|
362 | } |
---|
363 | cdt1 = (1.0-tau)/(ki*tau); |
---|
364 | S=0.0; |
---|
365 | //Incoherent scattering function |
---|
366 | for (G4int i=0;i<nosc;i++){ |
---|
367 | ionEnergy = (*(ionizationEnergy->find(Z)->second))[i]; |
---|
368 | if (photonEnergy0 > ionEnergy) //sum only on excitable levels |
---|
369 | { |
---|
370 | aux = photonEnergy0*(photonEnergy0-ionEnergy)*cdt1; |
---|
371 | harFunc = (*(hartreeFunction->find(Z)->second))[i]/fine_structure_const; |
---|
372 | occupNb = (G4int) (*(occupationNumber->find(Z)->second))[i]; |
---|
373 | pzomc = harFunc*(aux-electron_mass_c2*ionEnergy)/ |
---|
374 | (electron_mass_c2*std::sqrt(2.0*aux+ionEnergy*ionEnergy)); |
---|
375 | if (pzomc > 0) |
---|
376 | { |
---|
377 | rn[i] = 1.0-0.5*std::exp(0.5-(std::sqrt(0.5)+std::sqrt(2.0)*pzomc)* |
---|
378 | (std::sqrt(0.5)+std::sqrt(2.0)*pzomc)); |
---|
379 | } |
---|
380 | else |
---|
381 | { |
---|
382 | rn[i] = 0.5*std::exp(0.5-(std::sqrt(0.5)-std::sqrt(2.0)*pzomc)* |
---|
383 | (std::sqrt(0.5)-std::sqrt(2.0)*pzomc)); |
---|
384 | } |
---|
385 | S = S + occupNb*rn[i]; |
---|
386 | pac[i] = S; |
---|
387 | } |
---|
388 | else |
---|
389 | { |
---|
390 | pac[i] = S-(1e-06); |
---|
391 | } |
---|
392 | } |
---|
393 | //Rejection function |
---|
394 | TST = S*(1.0+tau*(ki1+tau*(ki2+tau*ki3)))/(ki3*tau*(1.0+tau*tau)); |
---|
395 | }while ((G4UniformRand()*s0) > TST); |
---|
396 | //Target electron shell |
---|
397 | cosTheta = 1.0 - cdt1; |
---|
398 | G4double fpzmax=0.0,fpz=0.0; |
---|
399 | G4double A=0.0; |
---|
400 | do |
---|
401 | { |
---|
402 | do |
---|
403 | { |
---|
404 | TST =S*G4UniformRand(); |
---|
405 | iosc=nosc; |
---|
406 | for (G4int i=0;i<nosc;i++){ |
---|
407 | if (pac[i]>TST) iosc = i; |
---|
408 | if (pac[i]>TST) break; |
---|
409 | } |
---|
410 | A = G4UniformRand()*rn[iosc]; |
---|
411 | harFunc = (*(hartreeFunction->find(Z)->second))[iosc]/fine_structure_const; |
---|
412 | occupNb = (G4int) (*(occupationNumber->find(Z)->second))[iosc]; |
---|
413 | if (A < 0.5) { |
---|
414 | pzomc = (std::sqrt(0.5)-std::sqrt(0.5-std::log(2.0*A)))/ |
---|
415 | (std::sqrt(2.0)*harFunc); |
---|
416 | } |
---|
417 | else |
---|
418 | { |
---|
419 | pzomc = (std::sqrt(0.5-std::log(2.0-2.0*A))-std::sqrt(0.5))/ |
---|
420 | (std::sqrt(2.0)*harFunc); |
---|
421 | } |
---|
422 | } while (pzomc < -1); |
---|
423 | // F(EP) rejection |
---|
424 | G4double XQC = 1.0+tau*(tau-2.0*cosTheta); |
---|
425 | G4double AF = std::sqrt(XQC)*(1.0+tau*(tau-cosTheta)/XQC); |
---|
426 | if (AF > 0) { |
---|
427 | fpzmax = 1.0+AF*0.2; |
---|
428 | } |
---|
429 | else |
---|
430 | { |
---|
431 | fpzmax = 1.0-AF*0.2; |
---|
432 | } |
---|
433 | fpz = 1.0+AF*std::max(std::min(pzomc,0.2),-0.2); |
---|
434 | }while ((fpzmax*G4UniformRand())>fpz); |
---|
435 | |
---|
436 | //Energy of the scattered photon |
---|
437 | G4double T = pzomc*pzomc; |
---|
438 | G4double b1 = 1.0-T*tau*tau; |
---|
439 | G4double b2 = 1.0-T*tau*cosTheta; |
---|
440 | if (pzomc > 0.0) |
---|
441 | { |
---|
442 | epsilon = (tau/b1)*(b2+std::sqrt(std::abs(b2*b2-b1*(1.0-T)))); |
---|
443 | } |
---|
444 | else |
---|
445 | { |
---|
446 | epsilon = (tau/b1)*(b2-std::sqrt(std::abs(b2*b2-b1*(1.0-T)))); |
---|
447 | } |
---|
448 | } |
---|
449 | |
---|
450 | |
---|
451 | //Ok, the kinematics has been calculated. |
---|
452 | G4double sinTheta = std::sqrt(1-cosTheta*cosTheta); |
---|
453 | G4double phi = twopi * G4UniformRand() ; |
---|
454 | G4double dirx = sinTheta * std::cos(phi); |
---|
455 | G4double diry = sinTheta * std::sin(phi); |
---|
456 | G4double dirz = cosTheta ; |
---|
457 | |
---|
458 | // Update G4VParticleChange for the scattered photon |
---|
459 | G4ThreeVector photonDirection1(dirx,diry,dirz); |
---|
460 | photonDirection1.rotateUz(photonDirection0); |
---|
461 | fParticleChange->ProposeMomentumDirection(photonDirection1) ; |
---|
462 | G4double photonEnergy1 = epsilon * photonEnergy0; |
---|
463 | |
---|
464 | if (photonEnergy1 > 0.) |
---|
465 | { |
---|
466 | fParticleChange->SetProposedKineticEnergy(photonEnergy1) ; |
---|
467 | } |
---|
468 | else |
---|
469 | { |
---|
470 | fParticleChange->SetProposedKineticEnergy(0.) ; |
---|
471 | fParticleChange->ProposeTrackStatus(fStopAndKill); |
---|
472 | } |
---|
473 | |
---|
474 | |
---|
475 | //Create scattered electron |
---|
476 | G4double diffEnergy = photonEnergy0*(1-epsilon); |
---|
477 | ionEnergy = (*(ionizationEnergy->find(Z)->second))[iosc]; |
---|
478 | G4double Q2 = photonEnergy0*photonEnergy0+photonEnergy1*(photonEnergy1-2.0*photonEnergy0*cosTheta); |
---|
479 | G4double cosThetaE; //scattering angle for the electron |
---|
480 | if (Q2 > 1.0e-12) |
---|
481 | { |
---|
482 | cosThetaE = (photonEnergy0-photonEnergy1*cosTheta)/std::sqrt(Q2); |
---|
483 | } |
---|
484 | else |
---|
485 | { |
---|
486 | cosThetaE = 1.0; |
---|
487 | } |
---|
488 | G4double sinThetaE = std::sqrt(1-cosThetaE*cosThetaE); |
---|
489 | |
---|
490 | //initialize here, then check photons created by Atomic-Deexcitation, and the final state e- |
---|
491 | std::vector<G4DynamicParticle*>* photonVector=0; |
---|
492 | |
---|
493 | const G4AtomicTransitionManager* transitionManager = G4AtomicTransitionManager::Instance(); |
---|
494 | const G4AtomicShell* shell = transitionManager->Shell(Z,iosc); |
---|
495 | G4double bindingEnergy = shell->BindingEnergy(); |
---|
496 | G4int shellId = shell->ShellId(); |
---|
497 | G4double ionEnergyInPenelopeDatabase = ionEnergy; |
---|
498 | ionEnergy = std::max(bindingEnergy,ionEnergyInPenelopeDatabase); //protection against energy non-conservation |
---|
499 | |
---|
500 | G4double eKineticEnergy = diffEnergy - ionEnergy; //subtract the excitation energy. If not emitted by fluorescence, |
---|
501 | //the ionization energy is deposited as local energy deposition |
---|
502 | G4double localEnergyDeposit = ionEnergy; |
---|
503 | G4double energyInFluorescence = 0.; //testing purposes only |
---|
504 | |
---|
505 | if (eKineticEnergy < 0) |
---|
506 | { |
---|
507 | //It means that there was some problem/mismatch between the two databases. Try to make it work |
---|
508 | //In this case available Energy (diffEnergy) < ionEnergy |
---|
509 | //Full residual energy is deposited locally |
---|
510 | localEnergyDeposit = diffEnergy; |
---|
511 | eKineticEnergy = 0.0; |
---|
512 | } |
---|
513 | G4double cutForLowEnergySecondaryPhotons = 250.0*eV; |
---|
514 | |
---|
515 | //Get the cut for electrons |
---|
516 | const G4ProductionCutsTable* theCoupleTable= |
---|
517 | G4ProductionCutsTable::GetProductionCutsTable(); |
---|
518 | size_t indx = couple->GetIndex(); |
---|
519 | G4double cutE = (*(theCoupleTable->GetEnergyCutsVector(1)))[indx]; |
---|
520 | cutE = std::max(cutForLowEnergySecondaryPhotons,cutE); |
---|
521 | |
---|
522 | //the local energy deposit is what remains: part of this may be spent for fluorescence. |
---|
523 | if (fUseAtomicDeexcitation) |
---|
524 | { |
---|
525 | G4int nPhotons=0; |
---|
526 | |
---|
527 | G4double cutg = (*(theCoupleTable->GetEnergyCutsVector(0)))[indx]; |
---|
528 | cutg = std::max(cutForLowEnergySecondaryPhotons,cutg); |
---|
529 | |
---|
530 | G4DynamicParticle* aPhoton; |
---|
531 | G4AtomicDeexcitation deexcitationManager; |
---|
532 | |
---|
533 | if (Z>5 && (localEnergyDeposit > cutg || localEnergyDeposit > cutE)) |
---|
534 | { |
---|
535 | photonVector = deexcitationManager.GenerateParticles(Z,shellId); |
---|
536 | for (size_t k=0;k<photonVector->size();k++){ |
---|
537 | aPhoton = (*photonVector)[k]; |
---|
538 | if (aPhoton) |
---|
539 | { |
---|
540 | G4double itsCut = cutg; |
---|
541 | if (aPhoton->GetDefinition() == G4Electron::Electron()) itsCut = cutE; |
---|
542 | G4double itsEnergy = aPhoton->GetKineticEnergy(); |
---|
543 | if (itsEnergy > itsCut && itsEnergy <= localEnergyDeposit) |
---|
544 | { |
---|
545 | nPhotons++; |
---|
546 | localEnergyDeposit -= itsEnergy; |
---|
547 | energyInFluorescence += itsEnergy; |
---|
548 | } |
---|
549 | else |
---|
550 | { |
---|
551 | delete aPhoton; |
---|
552 | (*photonVector)[k]=0; |
---|
553 | } |
---|
554 | } |
---|
555 | } |
---|
556 | } |
---|
557 | } |
---|
558 | |
---|
559 | //Produce explicitely the electron only if its proposed kinetic energy is |
---|
560 | //above the cut, otherwise add local energy deposition |
---|
561 | G4DynamicParticle* electron = 0; |
---|
562 | if (eKineticEnergy > cutE) |
---|
563 | { |
---|
564 | G4double xEl = sinThetaE * std::cos(phi+pi); |
---|
565 | G4double yEl = sinThetaE * std::sin(phi+pi); |
---|
566 | G4double zEl = cosThetaE; |
---|
567 | G4ThreeVector eDirection(xEl,yEl,zEl); //electron direction |
---|
568 | eDirection.rotateUz(photonDirection0); |
---|
569 | electron = new G4DynamicParticle (G4Electron::Electron(), |
---|
570 | eDirection,eKineticEnergy) ; |
---|
571 | fvect->push_back(electron); |
---|
572 | } |
---|
573 | else |
---|
574 | { |
---|
575 | localEnergyDeposit += eKineticEnergy; |
---|
576 | } |
---|
577 | |
---|
578 | //This block below is executed only if there is at least one secondary photon produced by |
---|
579 | //AtomicDeexcitation |
---|
580 | if (photonVector) |
---|
581 | { |
---|
582 | for (size_t ll=0;ll<photonVector->size();ll++) |
---|
583 | { |
---|
584 | if ((*photonVector)[ll]) |
---|
585 | { |
---|
586 | G4DynamicParticle* aFluorescencePhoton = (*photonVector)[ll]; |
---|
587 | fvect->push_back(aFluorescencePhoton); |
---|
588 | } |
---|
589 | } |
---|
590 | } |
---|
591 | delete photonVector; |
---|
592 | if (localEnergyDeposit < 0) |
---|
593 | { |
---|
594 | G4cout << "WARNING-" |
---|
595 | << "G4PenelopeComptonModel::SampleSecondaries - Negative energy deposit" |
---|
596 | << G4endl; |
---|
597 | localEnergyDeposit=0.; |
---|
598 | } |
---|
599 | fParticleChange->ProposeLocalEnergyDeposit(localEnergyDeposit); |
---|
600 | |
---|
601 | G4double electronEnergy = 0.; |
---|
602 | if (verboseLevel > 1) |
---|
603 | { |
---|
604 | G4cout << "-----------------------------------------------------------" << G4endl; |
---|
605 | G4cout << "Energy balance from G4PenelopeCompton" << G4endl; |
---|
606 | G4cout << "Incoming photon energy: " << photonEnergy0/keV << " keV" << G4endl; |
---|
607 | G4cout << "-----------------------------------------------------------" << G4endl; |
---|
608 | G4cout << "Scattered photon: " << photonEnergy1/keV << " keV" << G4endl; |
---|
609 | if (electron) |
---|
610 | electronEnergy = eKineticEnergy; |
---|
611 | G4cout << "Scattered electron " << electronEnergy/keV << " keV" << G4endl; |
---|
612 | G4cout << "Fluorescence: " << energyInFluorescence/keV << " keV" << G4endl; |
---|
613 | G4cout << "Local energy deposit " << localEnergyDeposit/keV << " keV" << G4endl; |
---|
614 | G4cout << "Total final state: " << (photonEnergy1+electronEnergy+energyInFluorescence+ |
---|
615 | localEnergyDeposit)/keV << |
---|
616 | " keV" << G4endl; |
---|
617 | G4cout << "-----------------------------------------------------------" << G4endl; |
---|
618 | } |
---|
619 | if (verboseLevel > 0) |
---|
620 | { |
---|
621 | G4double energyDiff = std::fabs(photonEnergy1+ |
---|
622 | electronEnergy+energyInFluorescence+ |
---|
623 | localEnergyDeposit-photonEnergy0); |
---|
624 | if (energyDiff > 0.05*keV) |
---|
625 | G4cout << "Warning from G4PenelopeCompton: problem with energy conservation: " << |
---|
626 | (photonEnergy1+electronEnergy+energyInFluorescence+localEnergyDeposit)/keV << |
---|
627 | " keV (final) vs. " << |
---|
628 | photonEnergy0/keV << " keV (initial)" << G4endl; |
---|
629 | } |
---|
630 | } |
---|
631 | |
---|
632 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
633 | |
---|
634 | void G4PenelopeComptonModel::ReadData() |
---|
635 | { |
---|
636 | char* path = getenv("G4LEDATA"); |
---|
637 | if (!path) |
---|
638 | { |
---|
639 | G4String excep = "G4PenelopeComptonModel - G4LEDATA environment variable not set!"; |
---|
640 | G4Exception(excep); |
---|
641 | } |
---|
642 | G4String pathString(path); |
---|
643 | G4String pathFile = pathString + "/penelope/compton-pen.dat"; |
---|
644 | std::ifstream file(pathFile); |
---|
645 | |
---|
646 | if (!file.is_open()) |
---|
647 | { |
---|
648 | G4String excep = "G4PenelopeComptonModel - data file " + pathFile + " not found!"; |
---|
649 | G4Exception(excep); |
---|
650 | } |
---|
651 | |
---|
652 | G4int k1,test,test1; |
---|
653 | G4double a1,a2; |
---|
654 | G4int Z=1,nLevels=0; |
---|
655 | |
---|
656 | if (!ionizationEnergy || !hartreeFunction || !occupationNumber) |
---|
657 | { |
---|
658 | G4String excep = "G4PenelopeComptonModel: problem with reading data from file"; |
---|
659 | G4Exception(excep); |
---|
660 | } |
---|
661 | |
---|
662 | do{ |
---|
663 | G4double harOfElectronsBelowThreshold = 0; |
---|
664 | G4int nbOfElectronsBelowThreshold = 0; |
---|
665 | G4DataVector* occVector = new G4DataVector; |
---|
666 | G4DataVector* harVector = new G4DataVector; |
---|
667 | G4DataVector* bindingEVector = new G4DataVector; |
---|
668 | file >> Z >> nLevels; |
---|
669 | for (G4int h=0;h<nLevels;h++) |
---|
670 | { |
---|
671 | file >> k1 >> a1 >> a2; |
---|
672 | //Make explicit unit of measurements for ionisation energy, which is MeV |
---|
673 | a1 *= MeV; |
---|
674 | if (a1 > 15*eV) |
---|
675 | { |
---|
676 | occVector->push_back((G4double) k1); |
---|
677 | bindingEVector->push_back(a1); |
---|
678 | harVector->push_back(a2); |
---|
679 | } |
---|
680 | else |
---|
681 | { |
---|
682 | nbOfElectronsBelowThreshold += k1; |
---|
683 | harOfElectronsBelowThreshold += k1*a2; |
---|
684 | } |
---|
685 | } |
---|
686 | //Add the "final" level |
---|
687 | if (nbOfElectronsBelowThreshold) |
---|
688 | { |
---|
689 | occVector->push_back(nbOfElectronsBelowThreshold); |
---|
690 | bindingEVector->push_back(0*eV); |
---|
691 | G4double averageHartree = |
---|
692 | harOfElectronsBelowThreshold/((G4double) nbOfElectronsBelowThreshold); |
---|
693 | harVector->push_back(averageHartree); |
---|
694 | } |
---|
695 | //Ok, done for element Z |
---|
696 | occupationNumber->insert(std::make_pair(Z,occVector)); |
---|
697 | ionizationEnergy->insert(std::make_pair(Z,bindingEVector)); |
---|
698 | hartreeFunction->insert(std::make_pair(Z,harVector)); |
---|
699 | file >> test >> test1; //-1 -1 close the data for each Z |
---|
700 | if (test > 0) { |
---|
701 | G4String excep = "G4PenelopeComptonModel - data file corrupted!"; |
---|
702 | G4Exception(excep); |
---|
703 | } |
---|
704 | }while (test != -2); //the very last Z is closed with -2 instead of -1 |
---|
705 | file.close(); |
---|
706 | if (verboseLevel > 2) |
---|
707 | { |
---|
708 | G4cout << "Data from G4PenelopeComptonModel read " << G4endl; |
---|
709 | } |
---|
710 | } |
---|
711 | |
---|
712 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
713 | |
---|
714 | G4double G4PenelopeComptonModel::DifferentialCrossSection(G4double cosTheta) |
---|
715 | { |
---|
716 | // |
---|
717 | // Penelope model for the Compton scattering differential cross section |
---|
718 | // dSigma/dOmega. |
---|
719 | // D. Brusa et al., Nucl. Instrum. Meth. A 379 (1996) 167 |
---|
720 | // The parametrization includes the J(p) distribution profiles for the atomic shells, |
---|
721 | // that are tabulated from Hartree-Fock calculations |
---|
722 | // from F. Biggs et al., At. Data Nucl. Data Tables 16 (1975) 201 |
---|
723 | // |
---|
724 | const G4double k2 = std::sqrt(2.0); |
---|
725 | const G4double k1 = std::sqrt(0.5); |
---|
726 | const G4double k12 = 0.5; |
---|
727 | G4double cdt1 = 1.0-cosTheta; |
---|
728 | G4double energy = energyForIntegration; |
---|
729 | G4int Z = ZForIntegration; |
---|
730 | //energy of Compton line; |
---|
731 | G4double EOEC = 1.0+(energy/electron_mass_c2)*cdt1; |
---|
732 | G4double ECOE = 1.0/EOEC; |
---|
733 | //Incoherent scattering function (analytical profile) |
---|
734 | G4double sia = 0.0; |
---|
735 | G4int nosc = occupationNumber->find(Z)->second->size(); |
---|
736 | for (G4int i=0;i<nosc;i++){ |
---|
737 | G4double ionEnergy = (*(ionizationEnergy->find(Z)->second))[i]; |
---|
738 | //Sum only of those shells for which E>Eion |
---|
739 | if (energy > ionEnergy) |
---|
740 | { |
---|
741 | G4double aux = energy * (energy-ionEnergy)*cdt1; |
---|
742 | G4double Pzimax = |
---|
743 | (aux - electron_mass_c2*ionEnergy)/(electron_mass_c2*std::sqrt(2*aux+ionEnergy*ionEnergy)); |
---|
744 | G4double harFunc = (*(hartreeFunction->find(Z)->second))[i]/fine_structure_const; |
---|
745 | G4int occupNb = (G4int) (*(occupationNumber->find(Z)->second))[i]; |
---|
746 | G4double x = harFunc*Pzimax; |
---|
747 | G4double siap = 0; |
---|
748 | if (x > 0) |
---|
749 | { |
---|
750 | siap = 1.0-0.5*std::exp(k12-(k1+k2*x)*(k1+k2*x)); |
---|
751 | } |
---|
752 | else |
---|
753 | { |
---|
754 | siap = 0.5*std::exp(k12-(k1-k2*x)*(k1-k2*x)); |
---|
755 | } |
---|
756 | sia = sia + occupNb*siap; //sum of all contributions; |
---|
757 | } |
---|
758 | } |
---|
759 | G4double XKN = EOEC+ECOE-1+cosTheta*cosTheta; |
---|
760 | G4double diffCS = pi*classic_electr_radius*classic_electr_radius*ECOE*ECOE*XKN*sia; |
---|
761 | return diffCS; |
---|
762 | } |
---|
763 | |
---|