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: G4PenelopeIonisationModel.cc,v 1.10 2009/10/23 09:29:24 pandola Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-03 $ |
---|
28 | // |
---|
29 | // Author: Luciano Pandola |
---|
30 | // |
---|
31 | // History: |
---|
32 | // -------- |
---|
33 | // 26 Nov 2008 L Pandola Migration from process to model |
---|
34 | // 17 Apr 2009 V Ivanchenko Cleanup initialisation and generation of secondaries: |
---|
35 | // - apply internal high-energy limit only in constructor |
---|
36 | // - do not apply low-energy limit (default is 0) |
---|
37 | // - added MinEnergyCut method |
---|
38 | // - do not change track status |
---|
39 | // 19 May 2009 L Pandola Explicitely set to zero pointers deleted in |
---|
40 | // Initialise(), since they might be checked later on |
---|
41 | // 21 Oct 2009 L Pandola Remove un-necessary fUseAtomicDeexcitation flag - now managed by |
---|
42 | // G4VEmModel::DeexcitationFlag() |
---|
43 | // Add ActivateAuger() method |
---|
44 | // |
---|
45 | |
---|
46 | #include "G4PenelopeIonisationModel.hh" |
---|
47 | #include "G4ParticleDefinition.hh" |
---|
48 | #include "G4MaterialCutsCouple.hh" |
---|
49 | #include "G4ProductionCutsTable.hh" |
---|
50 | #include "G4DynamicParticle.hh" |
---|
51 | #include "G4Element.hh" |
---|
52 | #include "G4AtomicTransitionManager.hh" |
---|
53 | #include "G4AtomicDeexcitation.hh" |
---|
54 | #include "G4AtomicShell.hh" |
---|
55 | #include "G4Gamma.hh" |
---|
56 | #include "G4Electron.hh" |
---|
57 | #include "G4Positron.hh" |
---|
58 | #include "G4CrossSectionHandler.hh" |
---|
59 | #include "G4AtomicDeexcitation.hh" |
---|
60 | #include "G4VEMDataSet.hh" |
---|
61 | |
---|
62 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
63 | |
---|
64 | |
---|
65 | G4PenelopeIonisationModel::G4PenelopeIonisationModel(const G4ParticleDefinition*, |
---|
66 | const G4String& nam) |
---|
67 | :G4VEmModel(nam),isInitialised(false), |
---|
68 | kineticEnergy1(0.*eV), |
---|
69 | cosThetaPrimary(1.0), |
---|
70 | energySecondary(0.*eV), |
---|
71 | cosThetaSecondary(0.0), |
---|
72 | iOsc(-1), |
---|
73 | crossSectionHandler(0),ionizationEnergy(0), |
---|
74 | resonanceEnergy(0),occupationNumber(0),shellFlag(0), |
---|
75 | theXSTable(0) |
---|
76 | { |
---|
77 | fIntrinsicLowEnergyLimit = 100.0*eV; |
---|
78 | fIntrinsicHighEnergyLimit = 100.0*GeV; |
---|
79 | // SetLowEnergyLimit(fIntrinsicLowEnergyLimit); |
---|
80 | SetHighEnergyLimit(fIntrinsicHighEnergyLimit); |
---|
81 | // |
---|
82 | // Atomic deexcitation model activated by default |
---|
83 | SetDeexcitationFlag(true); |
---|
84 | verboseLevel= 0; |
---|
85 | |
---|
86 | // Verbosity scale: |
---|
87 | // 0 = nothing |
---|
88 | // 1 = warning for energy non-conservation |
---|
89 | // 2 = details of energy budget |
---|
90 | // 3 = calculation of cross sections, file openings, sampling of atoms |
---|
91 | // 4 = entering in methods |
---|
92 | |
---|
93 | //These vectors do not change when materials or cut change. |
---|
94 | //Therefore I can read it at the constructor |
---|
95 | ionizationEnergy = new std::map<G4int,G4DataVector*>; |
---|
96 | resonanceEnergy = new std::map<G4int,G4DataVector*>; |
---|
97 | occupationNumber = new std::map<G4int,G4DataVector*>; |
---|
98 | shellFlag = new std::map<G4int,G4DataVector*>; |
---|
99 | |
---|
100 | ReadData(); //Read data from file |
---|
101 | |
---|
102 | } |
---|
103 | |
---|
104 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
105 | |
---|
106 | G4PenelopeIonisationModel::~G4PenelopeIonisationModel() |
---|
107 | { |
---|
108 | |
---|
109 | if (crossSectionHandler) delete crossSectionHandler; |
---|
110 | |
---|
111 | if (theXSTable) |
---|
112 | { |
---|
113 | for (size_t i=0; i<theXSTable->size(); i++) |
---|
114 | delete (*theXSTable)[i]; |
---|
115 | delete theXSTable; |
---|
116 | } |
---|
117 | |
---|
118 | |
---|
119 | std::map <G4int,G4DataVector*>::iterator i; |
---|
120 | for (i=ionizationEnergy->begin();i != ionizationEnergy->end();i++) |
---|
121 | if (i->second) delete i->second; |
---|
122 | for (i=resonanceEnergy->begin();i != resonanceEnergy->end();i++) |
---|
123 | if (i->second) delete i->second; |
---|
124 | for (i=occupationNumber->begin();i != occupationNumber->end();i++) |
---|
125 | if (i->second) delete i->second; |
---|
126 | for (i=shellFlag->begin();i != shellFlag->end();i++) |
---|
127 | if (i->second) delete i->second; |
---|
128 | |
---|
129 | if (ionizationEnergy) |
---|
130 | delete ionizationEnergy; |
---|
131 | if (resonanceEnergy) |
---|
132 | delete resonanceEnergy; |
---|
133 | if (occupationNumber) |
---|
134 | delete occupationNumber; |
---|
135 | if (shellFlag) |
---|
136 | delete shellFlag; |
---|
137 | } |
---|
138 | |
---|
139 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
140 | |
---|
141 | void G4PenelopeIonisationModel::Initialise(const G4ParticleDefinition* particle, |
---|
142 | const G4DataVector& cuts) |
---|
143 | { |
---|
144 | if (verboseLevel > 3) |
---|
145 | G4cout << "Calling G4PenelopeIonisationModel::Initialise()" << G4endl; |
---|
146 | |
---|
147 | //Delete and re-initialize the cross section handler |
---|
148 | if (crossSectionHandler) |
---|
149 | { |
---|
150 | crossSectionHandler->Clear(); |
---|
151 | delete crossSectionHandler; |
---|
152 | crossSectionHandler = 0; |
---|
153 | } |
---|
154 | |
---|
155 | if (theXSTable) |
---|
156 | { |
---|
157 | for (size_t i=0; i<theXSTable->size(); i++) |
---|
158 | { |
---|
159 | delete (*theXSTable)[i]; |
---|
160 | (*theXSTable)[i] = 0; |
---|
161 | } |
---|
162 | delete theXSTable; |
---|
163 | theXSTable = 0; |
---|
164 | } |
---|
165 | |
---|
166 | crossSectionHandler = new G4CrossSectionHandler(); |
---|
167 | crossSectionHandler->Clear(); |
---|
168 | G4String crossSectionFile = "NULL"; |
---|
169 | |
---|
170 | if (particle == G4Electron::Electron()) |
---|
171 | crossSectionFile = "penelope/ion-cs-el-"; |
---|
172 | else if (particle == G4Positron::Positron()) |
---|
173 | crossSectionFile = "penelope/ion-cs-po-"; |
---|
174 | crossSectionHandler->LoadData(crossSectionFile); |
---|
175 | //This is used to retrieve cross section values later on |
---|
176 | crossSectionHandler->BuildMeanFreePathForMaterials(); |
---|
177 | |
---|
178 | InitialiseElementSelectors(particle,cuts); |
---|
179 | |
---|
180 | if (verboseLevel > 2) |
---|
181 | G4cout << "Loaded cross section files for PenelopeIonisationModel" << G4endl; |
---|
182 | |
---|
183 | if (verboseLevel > 2) { |
---|
184 | G4cout << "Penelope Ionisation model is initialized " << G4endl |
---|
185 | << "Energy range: " |
---|
186 | << LowEnergyLimit() / keV << " keV - " |
---|
187 | << HighEnergyLimit() / GeV << " GeV" |
---|
188 | << G4endl; |
---|
189 | } |
---|
190 | |
---|
191 | if(isInitialised) return; |
---|
192 | fParticleChange = GetParticleChangeForLoss(); |
---|
193 | isInitialised = true; |
---|
194 | } |
---|
195 | |
---|
196 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
197 | |
---|
198 | G4double G4PenelopeIonisationModel::CrossSectionPerVolume(const G4Material* material, |
---|
199 | const G4ParticleDefinition* theParticle, |
---|
200 | G4double energy, |
---|
201 | G4double cutEnergy, |
---|
202 | G4double) |
---|
203 | { |
---|
204 | // Penelope model to calculate the cross section for inelastic collisions above the |
---|
205 | // threshold. It makes use of the Generalised Oscillator Strength (GOS) model from |
---|
206 | // D. Liljequist, J. Phys. D: Appl. Phys. 16 (1983) 1567 |
---|
207 | // |
---|
208 | // The total cross section (hard+soft) is read from a database file (element per |
---|
209 | // element), while the ratio hard-to-total is calculated analytically by taking |
---|
210 | // into account the atomic oscillators coming into the play for a given threshold. |
---|
211 | // This is done by the method CalculateCrossSectionsRatio(). |
---|
212 | // For incident e- the maximum energy allowed for the delta-rays is energy/2. |
---|
213 | // because particles are undistinghishable. |
---|
214 | // |
---|
215 | // The contribution is splitted in three parts: distant longitudinal collisions, |
---|
216 | // distant transverse collisions and close collisions. Each term is described by |
---|
217 | // its own analytical function. |
---|
218 | // Fermi density correction is calculated analytically according to |
---|
219 | // U. Fano, Ann. Rev. Nucl. Sci. 13 (1963),1 |
---|
220 | // |
---|
221 | if (verboseLevel > 3) |
---|
222 | G4cout << "Calling CrossSectionPerVolume() of G4PenelopeIonisationModel" << G4endl; |
---|
223 | |
---|
224 | SetupForMaterial(theParticle, material, energy); |
---|
225 | |
---|
226 | // VI - should be check at initialisation not in run time |
---|
227 | /* |
---|
228 | if (!crossSectionHandler) |
---|
229 | { |
---|
230 | G4cout << "G4PenelopeIonisationModel::CrossSectionPerVolume" << G4endl; |
---|
231 | G4cout << "The cross section handler is not correctly initialized" << G4endl; |
---|
232 | G4Exception(); |
---|
233 | } |
---|
234 | */ |
---|
235 | if (!theXSTable) |
---|
236 | { |
---|
237 | if (verboseLevel > 2) |
---|
238 | { |
---|
239 | G4cout << "G4PenelopeIonisationModel::CrossSectionPerVolume" << G4endl; |
---|
240 | G4cout << "Going to build Cross Section table " << G4endl; |
---|
241 | } |
---|
242 | theXSTable = new std::vector<G4VEMDataSet*>; |
---|
243 | theXSTable = BuildCrossSectionTable(theParticle); |
---|
244 | } |
---|
245 | |
---|
246 | G4double totalCross = 0.0; |
---|
247 | G4double cross = 0.0; |
---|
248 | const G4ElementVector* theElementVector = material->GetElementVector(); |
---|
249 | const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume(); |
---|
250 | G4double electronVolumeDensity = |
---|
251 | material->GetTotNbOfElectPerVolume(); //electron density |
---|
252 | |
---|
253 | if (verboseLevel > 4) |
---|
254 | G4cout << "Electron volume density of " << material->GetName() << ": " << |
---|
255 | electronVolumeDensity*cm3 << " electrons/cm3" << G4endl; |
---|
256 | |
---|
257 | G4int nelm = material->GetNumberOfElements(); |
---|
258 | for (G4int i=0; i<nelm; i++) |
---|
259 | { |
---|
260 | G4int iZ = (G4int) (*theElementVector)[i]->GetZ(); |
---|
261 | G4double ratio = CalculateCrossSectionsRatio(energy,cutEnergy,iZ,electronVolumeDensity, |
---|
262 | theParticle); |
---|
263 | cross += theAtomNumDensityVector[i]* |
---|
264 | crossSectionHandler->FindValue(iZ,energy)*ratio; |
---|
265 | totalCross += theAtomNumDensityVector[i]* |
---|
266 | crossSectionHandler->FindValue(iZ,energy); |
---|
267 | } |
---|
268 | |
---|
269 | if (verboseLevel > 2) |
---|
270 | { |
---|
271 | G4cout << "G4PenelopeIonisationModel " << G4endl; |
---|
272 | G4cout << "Mean free path for delta emission > " << cutEnergy/keV << " keV at " << |
---|
273 | energy/keV << " keV = " << (1./cross)/mm << " mm" << G4endl; |
---|
274 | G4cout << "Total free path for ionisation (no threshold) at " << |
---|
275 | energy/keV << " keV = " << (1./totalCross)/mm << " mm" << G4endl; |
---|
276 | } |
---|
277 | return cross; |
---|
278 | } |
---|
279 | |
---|
280 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
281 | |
---|
282 | G4double G4PenelopeIonisationModel::ComputeDEDXPerVolume(const G4Material* material, |
---|
283 | const G4ParticleDefinition* theParticle, |
---|
284 | G4double kineticEnergy, |
---|
285 | G4double cutEnergy) |
---|
286 | { |
---|
287 | // Penelope model to calculate the stopping power for soft inelastic collisions |
---|
288 | // below the threshold. It makes use of the Generalised Oscillator Strength (GOS) |
---|
289 | // model from |
---|
290 | // D. Liljequist, J. Phys. D: Appl. Phys. 16 (1983) 1567 |
---|
291 | // |
---|
292 | // The stopping power is calculated analytically using the dsigma/dW cross |
---|
293 | // section from the GOS models, which includes separate contributions from |
---|
294 | // distant longitudinal collisions, distant transverse collisions and |
---|
295 | // close collisions. Only the atomic oscillators that come in the play |
---|
296 | // (according to the threshold) are considered for the calculation. |
---|
297 | // Differential cross sections have a different form for e+ and e-. |
---|
298 | // |
---|
299 | // Fermi density correction is calculated analytically according to |
---|
300 | // U. Fano, Ann. Rev. Nucl. Sci. 13 (1963),1 |
---|
301 | |
---|
302 | if (verboseLevel > 3) |
---|
303 | G4cout << "Calling ComputeDEDX() of G4PenelopeIonisationModel" << G4endl; |
---|
304 | |
---|
305 | G4double sPower = 0.0; |
---|
306 | const G4ElementVector* theElementVector = material->GetElementVector(); |
---|
307 | const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume(); |
---|
308 | G4double electronVolumeDensity = |
---|
309 | material->GetTotNbOfElectPerVolume(); //electron density |
---|
310 | |
---|
311 | //Loop on the elements in the material |
---|
312 | G4int nelm = material->GetNumberOfElements(); |
---|
313 | for (G4int i=0; i<nelm; i++) |
---|
314 | { |
---|
315 | G4int iZ = (G4int) (*theElementVector)[i]->GetZ(); |
---|
316 | |
---|
317 | //Calculate stopping power contribution from each element |
---|
318 | //Constants |
---|
319 | G4double gamma = 1.0+kineticEnergy/electron_mass_c2; |
---|
320 | G4double gamma2 = gamma*gamma; |
---|
321 | G4double beta2 = (gamma2-1.0)/gamma2; |
---|
322 | G4double constant = pi*classic_electr_radius*classic_electr_radius |
---|
323 | *2.0*electron_mass_c2/beta2; |
---|
324 | |
---|
325 | G4double delta = CalculateDeltaFermi(kineticEnergy,iZ,electronVolumeDensity); |
---|
326 | G4int nbOsc = (G4int) resonanceEnergy->find(iZ)->second->size(); |
---|
327 | G4double stoppingPowerForElement = 0.0; |
---|
328 | //Loop on oscillators of element Z |
---|
329 | for (G4int iosc=0;iosc<nbOsc;iosc++) |
---|
330 | { |
---|
331 | G4double S1 = 0.0; |
---|
332 | G4double resEnergy = (*(resonanceEnergy->find(iZ)->second))[iosc]; |
---|
333 | if (theParticle == G4Electron::Electron()) |
---|
334 | S1 = ComputeStoppingPowerForElectrons(kineticEnergy,cutEnergy,delta,resEnergy); |
---|
335 | else if (theParticle == G4Positron::Positron()) |
---|
336 | S1 = ComputeStoppingPowerForPositrons(kineticEnergy,cutEnergy,delta,resEnergy); |
---|
337 | G4double occupNb = (*(occupationNumber->find(iZ)->second))[iosc]; |
---|
338 | stoppingPowerForElement += occupNb*constant*S1; |
---|
339 | } |
---|
340 | |
---|
341 | sPower += stoppingPowerForElement*theAtomNumDensityVector[i]; |
---|
342 | } |
---|
343 | if (verboseLevel > 2) |
---|
344 | { |
---|
345 | G4cout << "G4PenelopeIonisationModel " << G4endl; |
---|
346 | G4cout << "Stopping power < " << cutEnergy/keV << " keV at " << |
---|
347 | kineticEnergy/keV << " keV = " << sPower/(keV/mm) << " keV/mm" << G4endl; |
---|
348 | } |
---|
349 | return sPower; |
---|
350 | } |
---|
351 | |
---|
352 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
353 | |
---|
354 | void G4PenelopeIonisationModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect, |
---|
355 | const G4MaterialCutsCouple* couple, |
---|
356 | const G4DynamicParticle* aDynamicParticle, |
---|
357 | G4double cutE, G4double) |
---|
358 | { |
---|
359 | // Penelope model to sample the final state following an hard inelastic interaction. |
---|
360 | // It makes use of the Generalised Oscillator Strength (GOS) model from |
---|
361 | // D. Liljequist, J. Phys. D: Appl. Phys. 16 (1983) 1567 |
---|
362 | // |
---|
363 | // The GOS model is used to calculate the individual cross sections for all |
---|
364 | // the atomic oscillators coming in the play, taking into account the three |
---|
365 | // contributions (distant longitudinal collisions, distant transverse collisions and |
---|
366 | // close collisions). Then the target shell and the interaction channel are |
---|
367 | // sampled. Final state of the delta-ray (energy, angle) are generated according |
---|
368 | // to the analytical distributions (dSigma/dW) for the selected interaction |
---|
369 | // channels. |
---|
370 | // For e-, the maximum energy for the delta-ray is initialEnergy/2. (because |
---|
371 | // particles are indistinghusbable), while it is the full initialEnergy for |
---|
372 | // e+. |
---|
373 | // The efficiency of the random sampling algorithm (e.g. for close collisions) |
---|
374 | // decreases when initial and cutoff energy increase (e.g. 87% for 10-keV primary |
---|
375 | // and 1 keV threshold, 99% for 10-MeV primary and 10-keV threshold). |
---|
376 | // Differential cross sections have a different form for e+ and e-. |
---|
377 | // |
---|
378 | // WARNING: The model provides an _average_ description of inelastic collisions. |
---|
379 | // Anyway, the energy spectrum associated to distant excitations of a given |
---|
380 | // atomic shell is approximated as a single resonance. The simulated energy spectra |
---|
381 | // show _unphysical_ narrow peaks at energies that are multiple of the shell |
---|
382 | // resonance energies. The spurious speaks are automatically smoothed out after |
---|
383 | // multiple inelastic collisions. |
---|
384 | // |
---|
385 | // The model determines also the original shell from which the delta-ray is expelled, |
---|
386 | // in order to produce fluorescence de-excitation (from G4DeexcitationManager) |
---|
387 | // |
---|
388 | // Fermi density correction is calculated analytically according to |
---|
389 | // U. Fano, Ann. Rev. Nucl. Sci. 13 (1963),1 |
---|
390 | |
---|
391 | if (verboseLevel > 3) |
---|
392 | G4cout << "Calling SamplingSecondaries() of G4PenelopeIonisationModel" << G4endl; |
---|
393 | |
---|
394 | G4double kineticEnergy0 = aDynamicParticle->GetKineticEnergy(); |
---|
395 | const G4ParticleDefinition* theParticle = aDynamicParticle->GetDefinition(); |
---|
396 | |
---|
397 | if (kineticEnergy0 <= fIntrinsicLowEnergyLimit) |
---|
398 | { |
---|
399 | fParticleChange->SetProposedKineticEnergy(0.); |
---|
400 | fParticleChange->ProposeLocalEnergyDeposit(kineticEnergy0); |
---|
401 | return ; |
---|
402 | } |
---|
403 | const G4double electronVolumeDensity = |
---|
404 | couple->GetMaterial()->GetTotNbOfElectPerVolume(); |
---|
405 | G4ParticleMomentum particleDirection0 = aDynamicParticle->GetMomentumDirection(); |
---|
406 | |
---|
407 | //Initialise final-state variables. The proper values will be set by the methods |
---|
408 | // CalculateDiscreteForElectrons() and CalculateDiscreteForPositrons() |
---|
409 | kineticEnergy1=kineticEnergy0; |
---|
410 | cosThetaPrimary=1.0; |
---|
411 | energySecondary=0.0; |
---|
412 | cosThetaSecondary=1.0; |
---|
413 | |
---|
414 | // Select randomly one element in the current material |
---|
415 | if (verboseLevel > 2) |
---|
416 | G4cout << "Going to select element in " << couple->GetMaterial()->GetName() << G4endl; |
---|
417 | |
---|
418 | //Use sampler of G4CrossSectionHandler for now |
---|
419 | // atom can be selected effitiantly if element selectors are initialised |
---|
420 | //const G4Element* anElement = SelectRandomAtom(couple,theParticle,kineticEnergy0); |
---|
421 | |
---|
422 | G4int iZ = SampleRandomAtom(couple,kineticEnergy0); |
---|
423 | |
---|
424 | const G4ProductionCutsTable* theCoupleTable= |
---|
425 | G4ProductionCutsTable::GetProductionCutsTable(); |
---|
426 | size_t indx = couple->GetIndex(); |
---|
427 | G4double cutG = (*(theCoupleTable->GetEnergyCutsVector(0)))[indx]; |
---|
428 | |
---|
429 | if (verboseLevel > 2) |
---|
430 | G4cout << "Selected Z = " << iZ << G4endl; |
---|
431 | |
---|
432 | // The method CalculateDiscreteForXXXX() set the private variables: |
---|
433 | // kineticEnergy1 = energy of the primary electron/positron after the interaction |
---|
434 | // cosThetaPrimary = std::cos(theta) of the primary after the interaction |
---|
435 | // energySecondary = energy of the secondary electron |
---|
436 | // cosThetaSecondary = std::cos(theta) of the secondary |
---|
437 | if (theParticle == G4Electron::Electron()) |
---|
438 | CalculateDiscreteForElectrons(kineticEnergy0,cutE,iZ,electronVolumeDensity); |
---|
439 | else if (theParticle == G4Positron::Positron()) |
---|
440 | CalculateDiscreteForPositrons(kineticEnergy0,cutE,iZ,electronVolumeDensity); |
---|
441 | |
---|
442 | // if (energySecondary == 0) return; |
---|
443 | |
---|
444 | //Update the primary particle |
---|
445 | G4double sint = std::sqrt(1. - cosThetaPrimary*cosThetaPrimary); |
---|
446 | G4double phi = twopi * G4UniformRand(); |
---|
447 | G4double dirx = sint * std::cos(phi); |
---|
448 | G4double diry = sint * std::sin(phi); |
---|
449 | G4double dirz = cosThetaPrimary; |
---|
450 | |
---|
451 | G4ThreeVector electronDirection1(dirx,diry,dirz); |
---|
452 | electronDirection1.rotateUz(particleDirection0); |
---|
453 | |
---|
454 | if (kineticEnergy1 > 0) |
---|
455 | { |
---|
456 | fParticleChange->ProposeMomentumDirection(electronDirection1); |
---|
457 | fParticleChange->SetProposedKineticEnergy(kineticEnergy1); |
---|
458 | } |
---|
459 | else |
---|
460 | { |
---|
461 | fParticleChange->SetProposedKineticEnergy(0.); |
---|
462 | } |
---|
463 | |
---|
464 | //Generate the delta ray |
---|
465 | G4int iosc2 = 0; |
---|
466 | G4double ioniEnergy = 0.0; |
---|
467 | if (iOsc > 0) { |
---|
468 | ioniEnergy=(*(ionizationEnergy->find(iZ)->second))[iOsc]; |
---|
469 | iosc2 = (ionizationEnergy->find(iZ)->second->size()) - iOsc; //they are in reversed order |
---|
470 | } |
---|
471 | |
---|
472 | const G4AtomicTransitionManager* transitionManager = G4AtomicTransitionManager::Instance(); |
---|
473 | G4double bindingEnergy = 0.0; |
---|
474 | G4int shellId = 0; |
---|
475 | if (iOsc > 0) |
---|
476 | { |
---|
477 | const G4AtomicShell* shell = transitionManager->Shell(iZ,iosc2-1); // Modified by Alf |
---|
478 | bindingEnergy = shell->BindingEnergy(); |
---|
479 | shellId = shell->ShellId(); |
---|
480 | } |
---|
481 | |
---|
482 | G4double ionEnergy = bindingEnergy; //energy spent to ionise the atom according to G4dabatase |
---|
483 | G4double eKineticEnergy = energySecondary; |
---|
484 | |
---|
485 | //This is an awful thing: Penelope generates the fluorescence only for L and K shells |
---|
486 | //(i.e. Osc = 1 --> 4). For high-Z, the other shells can be quite relevant. In this case |
---|
487 | //one MUST ensure ''by hand'' the energy conservation. Then there is the other problem that |
---|
488 | //the fluorescence database of Penelope doesn not match that of Geant4. |
---|
489 | |
---|
490 | G4double energyBalance = kineticEnergy0 - kineticEnergy1 - energySecondary; //Penelope Balance |
---|
491 | |
---|
492 | if (std::abs(energyBalance) < 1*eV) |
---|
493 | //in this case Penelope didn't subtract the fluorescence energy: do here by hand |
---|
494 | eKineticEnergy = energySecondary - bindingEnergy; |
---|
495 | else |
---|
496 | //Penelope subtracted the fluorescence, but one has to match the databases |
---|
497 | eKineticEnergy = energySecondary+ioniEnergy-bindingEnergy; |
---|
498 | |
---|
499 | G4double localEnergyDeposit = ionEnergy; |
---|
500 | G4double energyInFluorescence = 0.0*eV; |
---|
501 | |
---|
502 | if(DeexcitationFlag() && iZ > 5) |
---|
503 | { |
---|
504 | if (ionEnergy > cutG || ionEnergy > cutE) |
---|
505 | { |
---|
506 | deexcitationManager.SetCutForSecondaryPhotons(cutG); |
---|
507 | deexcitationManager.SetCutForAugerElectrons(cutE); |
---|
508 | std::vector<G4DynamicParticle*> *photonVector = |
---|
509 | deexcitationManager.GenerateParticles(iZ,shellId); |
---|
510 | //Check for secondaries |
---|
511 | if(photonVector) |
---|
512 | { |
---|
513 | for (size_t k=0;k<photonVector->size();k++) |
---|
514 | { |
---|
515 | G4DynamicParticle* aPhoton = (*photonVector)[k]; |
---|
516 | if (aPhoton) |
---|
517 | { |
---|
518 | G4double itsEnergy = aPhoton->GetKineticEnergy(); |
---|
519 | if (itsEnergy <= localEnergyDeposit) |
---|
520 | { |
---|
521 | if(aPhoton->GetDefinition() == G4Gamma::Gamma()) |
---|
522 | energyInFluorescence += itsEnergy; |
---|
523 | localEnergyDeposit -= itsEnergy; |
---|
524 | fvect->push_back(aPhoton); |
---|
525 | } |
---|
526 | else |
---|
527 | { |
---|
528 | delete aPhoton; |
---|
529 | (*photonVector)[k] = 0; |
---|
530 | } |
---|
531 | } |
---|
532 | } |
---|
533 | delete photonVector; |
---|
534 | } |
---|
535 | } |
---|
536 | } |
---|
537 | |
---|
538 | // Generate the delta ray |
---|
539 | G4double sin2 = std::sqrt(1. - cosThetaSecondary*cosThetaSecondary); |
---|
540 | G4double phi2 = twopi * G4UniformRand(); |
---|
541 | |
---|
542 | G4double xEl = sin2 * std::cos(phi2); |
---|
543 | G4double yEl = sin2 * std::sin(phi2); |
---|
544 | G4double zEl = cosThetaSecondary; |
---|
545 | G4ThreeVector eDirection(xEl,yEl,zEl); //electron direction |
---|
546 | eDirection.rotateUz(particleDirection0); |
---|
547 | |
---|
548 | G4DynamicParticle* deltaElectron = new G4DynamicParticle (G4Electron::Electron(), |
---|
549 | eDirection,eKineticEnergy) ; |
---|
550 | fvect->push_back(deltaElectron); |
---|
551 | |
---|
552 | if (localEnergyDeposit < 0) |
---|
553 | { |
---|
554 | G4cout << "WARNING-" |
---|
555 | << "G4PenelopeIonisationModel::SampleSecondaries - Negative energy deposit" |
---|
556 | << G4endl; |
---|
557 | localEnergyDeposit=0.; |
---|
558 | } |
---|
559 | fParticleChange->ProposeLocalEnergyDeposit(localEnergyDeposit); |
---|
560 | |
---|
561 | if (verboseLevel > 1) |
---|
562 | { |
---|
563 | G4cout << "-----------------------------------------------------------" << G4endl; |
---|
564 | G4cout << "Energy balance from G4PenelopeIonisation" << G4endl; |
---|
565 | G4cout << "Incoming primary energy: " << kineticEnergy0/keV << " keV" << G4endl; |
---|
566 | G4cout << "-----------------------------------------------------------" << G4endl; |
---|
567 | G4cout << "Outgoing primary energy: " << kineticEnergy1/keV << " keV" << G4endl; |
---|
568 | G4cout << "Delta ray " << eKineticEnergy/keV << " keV" << G4endl; |
---|
569 | G4cout << "Fluorescence: " << energyInFluorescence/keV << " keV" << G4endl; |
---|
570 | G4cout << "Local energy deposit " << localEnergyDeposit/keV << " keV" << G4endl; |
---|
571 | G4cout << "Total final state: " << (eKineticEnergy+energyInFluorescence+kineticEnergy1+ |
---|
572 | localEnergyDeposit)/keV << |
---|
573 | " keV" << G4endl; |
---|
574 | G4cout << "-----------------------------------------------------------" << G4endl; |
---|
575 | } |
---|
576 | if (verboseLevel > 0) |
---|
577 | { |
---|
578 | G4double energyDiff = std::fabs(eKineticEnergy+energyInFluorescence+kineticEnergy1+ |
---|
579 | localEnergyDeposit-kineticEnergy0); |
---|
580 | if (energyDiff > 0.05*keV) |
---|
581 | G4cout << "Warning from G4PenelopeIonisation: problem with energy conservation: " << |
---|
582 | (eKineticEnergy+energyInFluorescence+kineticEnergy1+localEnergyDeposit)/keV << |
---|
583 | " keV (final) vs. " << |
---|
584 | kineticEnergy0/keV << " keV (initial)" << G4endl; |
---|
585 | } |
---|
586 | } |
---|
587 | |
---|
588 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
589 | |
---|
590 | void G4PenelopeIonisationModel::ReadData() |
---|
591 | { |
---|
592 | if (verboseLevel > 2) |
---|
593 | { |
---|
594 | G4cout << "Data from G4PenelopeIonisationModel read " << G4endl; |
---|
595 | } |
---|
596 | char* path = getenv("G4LEDATA"); |
---|
597 | if (!path) |
---|
598 | { |
---|
599 | G4String excep = "G4PenelopeIonisationModel - G4LEDATA environment variable not set!"; |
---|
600 | G4Exception(excep); |
---|
601 | } |
---|
602 | G4String pathString(path); |
---|
603 | G4String pathFile = pathString + "/penelope/ion-pen.dat"; |
---|
604 | std::ifstream file(pathFile); |
---|
605 | |
---|
606 | if (!file.is_open()) |
---|
607 | { |
---|
608 | G4String excep = "G4PenelopeIonisationModel - data file " + pathFile + " not found!"; |
---|
609 | G4Exception(excep); |
---|
610 | } |
---|
611 | |
---|
612 | if (!ionizationEnergy || !resonanceEnergy || !occupationNumber || !shellFlag) |
---|
613 | { |
---|
614 | G4String excep = "G4PenelopeIonisationModel: problem with reading data from file"; |
---|
615 | G4Exception(excep); |
---|
616 | } |
---|
617 | |
---|
618 | G4int Z=1,nLevels=0; |
---|
619 | G4int test,test1; |
---|
620 | |
---|
621 | do{ |
---|
622 | G4DataVector* occVector = new G4DataVector; |
---|
623 | G4DataVector* ionEVector = new G4DataVector; |
---|
624 | G4DataVector* resEVector = new G4DataVector; |
---|
625 | G4DataVector* shellIndVector = new G4DataVector; |
---|
626 | file >> Z >> nLevels; |
---|
627 | G4double a1,a2,a3,a4; |
---|
628 | G4int k1,k2,k3; |
---|
629 | for (G4int h=0;h<nLevels;h++) |
---|
630 | { |
---|
631 | //index,occup number,ion energy,res energy,fj0,kz,shell flag |
---|
632 | file >> k1 >> a1 >> a2 >> a3 >> a4 >> k2 >> k3; |
---|
633 | //Make explicit unit of measurements for ionisation and resonance energy, |
---|
634 | // which is MeV |
---|
635 | a2 *= MeV; |
---|
636 | a3 *= MeV; |
---|
637 | // |
---|
638 | occVector->push_back(a1); |
---|
639 | ionEVector->push_back(a2); |
---|
640 | resEVector->push_back(a3); |
---|
641 | shellIndVector->push_back((G4double) k3); |
---|
642 | } |
---|
643 | //Ok, done for element Z |
---|
644 | occupationNumber->insert(std::make_pair(Z,occVector)); |
---|
645 | ionizationEnergy->insert(std::make_pair(Z,ionEVector)); |
---|
646 | resonanceEnergy->insert(std::make_pair(Z,resEVector)); |
---|
647 | shellFlag->insert(std::make_pair(Z,shellIndVector)); |
---|
648 | file >> test >> test1; //-1 -1 close the data for each Z |
---|
649 | if (test > 0) { |
---|
650 | G4String excep = "G4PenelopeIonisationModel - data file corrupted!"; |
---|
651 | G4Exception(excep); |
---|
652 | } |
---|
653 | }while (test != -2); //the very last Z is closed with -2 instead of -1 |
---|
654 | } |
---|
655 | |
---|
656 | |
---|
657 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
658 | |
---|
659 | G4double G4PenelopeIonisationModel::CalculateDeltaFermi(G4double kinEnergy ,G4int Z, |
---|
660 | G4double electronVolumeDensity) |
---|
661 | { |
---|
662 | G4double plasmaEnergyCoefficient = 1.377e-39*(MeV*MeV*m3); //(e*hbar)^2/(epsilon0*electron_mass) |
---|
663 | G4double plasmaEnergySquared = plasmaEnergyCoefficient*electronVolumeDensity; |
---|
664 | // std::sqrt(plasmaEnergySquared) is the plasma energy of the solid (MeV) |
---|
665 | G4double gam = 1.0+kinEnergy/electron_mass_c2; |
---|
666 | G4double gam2=gam*gam; |
---|
667 | G4double delta = 0.0; |
---|
668 | |
---|
669 | //Density effect |
---|
670 | G4double TST = ((G4double) Z)/(gam2*plasmaEnergySquared); |
---|
671 | |
---|
672 | G4double wl2 = 0.0; |
---|
673 | G4double fdel=0.0; |
---|
674 | G4double wr=0; |
---|
675 | size_t nbOsc = resonanceEnergy->find(Z)->second->size(); |
---|
676 | for(size_t i=0;i<nbOsc;i++) |
---|
677 | { |
---|
678 | G4int occupNb = (G4int) (*(occupationNumber->find(Z)->second))[i]; |
---|
679 | wr = (*(resonanceEnergy->find(Z)->second))[i]; |
---|
680 | fdel += occupNb/(wr*wr+wl2); |
---|
681 | } |
---|
682 | if (fdel < TST) return delta; |
---|
683 | G4double help1 = (*(resonanceEnergy->find(Z)->second))[nbOsc-1]; |
---|
684 | wl2 = help1*help1; |
---|
685 | do{ |
---|
686 | wl2=wl2*2.0; |
---|
687 | fdel = 0.0; |
---|
688 | for (size_t ii=0;ii<nbOsc;ii++){ |
---|
689 | G4int occupNb = (G4int) (*(occupationNumber->find(Z)->second))[ii]; |
---|
690 | wr = (*(resonanceEnergy->find(Z)->second))[ii]; |
---|
691 | fdel += occupNb/(wr*wr+wl2); |
---|
692 | } |
---|
693 | }while (fdel > TST); |
---|
694 | G4double wl2l=0.0; |
---|
695 | G4double wl2u = wl2; |
---|
696 | G4double control = 0.0; |
---|
697 | do{ |
---|
698 | wl2=0.5*(wl2l+wl2u); |
---|
699 | fdel = 0.0; |
---|
700 | for (size_t jj=0;jj<nbOsc;jj++){ |
---|
701 | G4int occupNb = (G4int) (*(occupationNumber->find(Z)->second))[jj]; |
---|
702 | wr = (*(resonanceEnergy->find(Z)->second))[jj]; |
---|
703 | fdel += occupNb/(wr*wr+wl2); |
---|
704 | } |
---|
705 | if (fdel > TST) |
---|
706 | wl2l = wl2; |
---|
707 | else |
---|
708 | wl2u = wl2; |
---|
709 | control = wl2u-wl2l-wl2*1e-12; |
---|
710 | }while(control>0); |
---|
711 | |
---|
712 | //Density correction effect |
---|
713 | for (size_t kk=0;kk<nbOsc;kk++){ |
---|
714 | G4int occupNb = (G4int) (*(occupationNumber->find(Z)->second))[kk]; |
---|
715 | wr = (*(resonanceEnergy->find(Z)->second))[kk]; |
---|
716 | delta += occupNb*std::log(1.0+wl2/(wr*wr)); |
---|
717 | } |
---|
718 | delta = (delta/((G4double) Z))-wl2/(gam2*plasmaEnergySquared); |
---|
719 | return delta; |
---|
720 | } |
---|
721 | |
---|
722 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
723 | |
---|
724 | void |
---|
725 | G4PenelopeIonisationModel::CalculateDiscreteForElectrons(G4double kinEnergy, |
---|
726 | G4double cutoffEnergy, |
---|
727 | G4int Z, |
---|
728 | G4double electronVolumeDensity) |
---|
729 | { |
---|
730 | if (verboseLevel > 2) |
---|
731 | G4cout << "Entering in CalculateDiscreteForElectrons() for energy " << |
---|
732 | kinEnergy/keV << " keV " << G4endl; |
---|
733 | |
---|
734 | //Initialization of variables to be calculated in this routine |
---|
735 | kineticEnergy1=kinEnergy; |
---|
736 | cosThetaPrimary=1.0; |
---|
737 | energySecondary=0.0; |
---|
738 | cosThetaSecondary=1.0; |
---|
739 | iOsc=-1; |
---|
740 | |
---|
741 | //constants |
---|
742 | G4double rb=kinEnergy+2.0*electron_mass_c2; |
---|
743 | G4double gamma = 1.0+kinEnergy/electron_mass_c2; |
---|
744 | G4double gamma2 = gamma*gamma; |
---|
745 | G4double beta2 = (gamma2-1.0)/gamma2; |
---|
746 | G4double amol = (gamma-1.0)*(gamma-1.0)/gamma2; |
---|
747 | G4double cps = kinEnergy*rb; |
---|
748 | G4double cp = std::sqrt(cps); |
---|
749 | |
---|
750 | G4double delta = CalculateDeltaFermi(kinEnergy,Z,electronVolumeDensity); |
---|
751 | G4double distantTransvCS0 = std::max(std::log(gamma2)-beta2-delta,0.0); |
---|
752 | |
---|
753 | G4double rl,rl1; |
---|
754 | |
---|
755 | if (cutoffEnergy > kinEnergy) return; //delta rays are not generated |
---|
756 | |
---|
757 | G4DataVector* qm = new G4DataVector(); |
---|
758 | G4DataVector* cumulHardCS = new G4DataVector(); |
---|
759 | std::vector<G4int> *typeOfInteraction = new std::vector<G4int>; |
---|
760 | G4DataVector* nbOfLevel = new G4DataVector(); |
---|
761 | |
---|
762 | //Hard close collisions with outer shells |
---|
763 | G4double wmaxc = 0.5*kinEnergy; |
---|
764 | G4double closeCS0 = 0.0; |
---|
765 | G4double closeCS = 0.0; |
---|
766 | if (cutoffEnergy>0.1*eV) |
---|
767 | { |
---|
768 | rl=cutoffEnergy/kinEnergy; |
---|
769 | rl1=1.0-rl; |
---|
770 | if (rl < 0.5) |
---|
771 | closeCS0 = (amol*(0.5-rl)+(1.0/rl)-(1.0/rl1)+(1.0-amol)*std::log(rl/rl1))/kinEnergy; |
---|
772 | } |
---|
773 | |
---|
774 | // Cross sections for the different oscillators |
---|
775 | |
---|
776 | // totalHardCS contains the cumulative hard interaction cross section for the different |
---|
777 | // excitable levels and the different interaction channels (close, distant, etc.), |
---|
778 | // i.e. |
---|
779 | // cumulHardCS[0] = 0.0 |
---|
780 | // cumulHardCS[1] = 1st excitable level (distant longitudinal only) |
---|
781 | // cumulHardCS[2] = 1st excitable level (distant longitudinal + transverse) |
---|
782 | // cumulHardCS[3] = 1st excitable level (distant longitudinal + transverse + close) |
---|
783 | // cumulHardCS[4] = 1st excitable level (all channels) + 2nd excitable level (distant long only) |
---|
784 | // etc. |
---|
785 | // This is used for sampling the atomic level which is ionised and the channel of the |
---|
786 | // interaction. |
---|
787 | // |
---|
788 | // For each index iFill of the cumulHardCS vector, |
---|
789 | // nbOfLevel[iFill] contains the current excitable atomic level and |
---|
790 | // typeOfInteraction[iFill] contains the current interaction channel, with the legenda: |
---|
791 | // 1 = distant longitudinal interaction |
---|
792 | // 2 = distant transverse interaction |
---|
793 | // 3 = close collision |
---|
794 | // 4 = close collision with outer shells (in this case nbOfLevel < 0 --> no binding energy) |
---|
795 | |
---|
796 | |
---|
797 | G4int nOscil = ionizationEnergy->find(Z)->second->size(); |
---|
798 | G4double totalHardCS = 0.0; |
---|
799 | G4double involvedElectrons = 0.0; |
---|
800 | for (G4int i=0;i<nOscil;i++){ |
---|
801 | G4double wi = (*(resonanceEnergy->find(Z)->second))[i]; |
---|
802 | G4int occupNb = (G4int) (*(occupationNumber->find(Z)->second))[i]; |
---|
803 | |
---|
804 | //Distant excitations |
---|
805 | if (wi>cutoffEnergy && wi<kinEnergy) |
---|
806 | { |
---|
807 | if (wi>(1e-6*kinEnergy)) |
---|
808 | { |
---|
809 | G4double cpp=std::sqrt((kinEnergy-wi)*(kinEnergy-wi+2.0*electron_mass_c2)); |
---|
810 | qm->push_back(std::sqrt((cp-cpp)*(cp-cpp)+electron_mass_c2*electron_mass_c2)-electron_mass_c2); |
---|
811 | } |
---|
812 | else |
---|
813 | { |
---|
814 | qm->push_back((wi*wi)/(beta2*2.0*electron_mass_c2)); |
---|
815 | } |
---|
816 | if ((*qm)[i] < wi) |
---|
817 | { |
---|
818 | |
---|
819 | G4double distantLongitCS = occupNb*std::log(wi*((*qm)[i]+2.0*electron_mass_c2)/ |
---|
820 | ((*qm)[i]*(wi+2.0*electron_mass_c2)))/wi; |
---|
821 | cumulHardCS->push_back(totalHardCS); |
---|
822 | typeOfInteraction->push_back(1); //distant longitudinal |
---|
823 | nbOfLevel->push_back((G4double) i); //only excitable level are counted |
---|
824 | totalHardCS += distantLongitCS; |
---|
825 | |
---|
826 | G4double distantTransvCS = occupNb*distantTransvCS0/wi; |
---|
827 | |
---|
828 | cumulHardCS->push_back(totalHardCS); |
---|
829 | typeOfInteraction->push_back(2); //distant tranverse |
---|
830 | nbOfLevel->push_back((G4double) i); |
---|
831 | totalHardCS += distantTransvCS; |
---|
832 | } |
---|
833 | } |
---|
834 | else |
---|
835 | qm->push_back(wi); |
---|
836 | |
---|
837 | //close collisions |
---|
838 | if(wi < wmaxc) |
---|
839 | { |
---|
840 | if (wi < cutoffEnergy) |
---|
841 | involvedElectrons += occupNb; |
---|
842 | else |
---|
843 | { |
---|
844 | rl=wi/kinEnergy; |
---|
845 | rl1=1.0-rl; |
---|
846 | closeCS = occupNb*(amol*(0.5-rl)+(1.0/rl)-(1.0/rl1)+(1.0-amol)*std::log(rl/rl1))/kinEnergy; |
---|
847 | cumulHardCS->push_back(totalHardCS); |
---|
848 | typeOfInteraction->push_back(3); //close |
---|
849 | nbOfLevel->push_back((G4double) i); |
---|
850 | totalHardCS += closeCS; |
---|
851 | } |
---|
852 | } |
---|
853 | } // loop on the levels |
---|
854 | |
---|
855 | cumulHardCS->push_back(totalHardCS); |
---|
856 | typeOfInteraction->push_back(4); //close interaction with outer shells |
---|
857 | nbOfLevel->push_back(-1.0); |
---|
858 | totalHardCS += involvedElectrons*closeCS0; |
---|
859 | cumulHardCS->push_back(totalHardCS); //this is the final value of the totalHardCS |
---|
860 | |
---|
861 | if (totalHardCS < 1e-30) { |
---|
862 | kineticEnergy1=kinEnergy; |
---|
863 | cosThetaPrimary=1.0; |
---|
864 | energySecondary=0.0; |
---|
865 | cosThetaSecondary=0.0; |
---|
866 | iOsc=-1; |
---|
867 | delete qm; |
---|
868 | delete cumulHardCS; |
---|
869 | delete typeOfInteraction; |
---|
870 | delete nbOfLevel; |
---|
871 | return; |
---|
872 | } |
---|
873 | |
---|
874 | //Testing purposes |
---|
875 | if (verboseLevel > 6) |
---|
876 | { |
---|
877 | for (size_t t=0;t<cumulHardCS->size();t++) |
---|
878 | G4cout << (*cumulHardCS)[t] << " " << (*typeOfInteraction)[t] << |
---|
879 | " " << (*nbOfLevel)[t] << G4endl; |
---|
880 | } |
---|
881 | |
---|
882 | //Selection of the active oscillator on the basis of the cumulative cross sections |
---|
883 | G4double TST = totalHardCS*G4UniformRand(); |
---|
884 | G4int is=0; |
---|
885 | G4int js= nbOfLevel->size(); |
---|
886 | do{ |
---|
887 | G4int it=(is+js)/2; |
---|
888 | if (TST > (*cumulHardCS)[it]) is=it; |
---|
889 | if (TST <= (*cumulHardCS)[it]) js=it; |
---|
890 | }while((js-is) > 1); |
---|
891 | |
---|
892 | G4double UII=0.0; |
---|
893 | G4double rkc=cutoffEnergy/kinEnergy; |
---|
894 | G4double dde; |
---|
895 | G4int kks; |
---|
896 | |
---|
897 | G4int sampledInteraction = (*typeOfInteraction)[is]; |
---|
898 | iOsc = (G4int) (*nbOfLevel)[is]; |
---|
899 | |
---|
900 | if (verboseLevel > 4) |
---|
901 | G4cout << "Chosen interaction #:" << sampledInteraction << " on level " << iOsc << G4endl; |
---|
902 | |
---|
903 | //Generates the final state according to the sampled level and |
---|
904 | //interaction channel |
---|
905 | |
---|
906 | if (sampledInteraction == 1) //Hard distant longitudinal collisions |
---|
907 | { |
---|
908 | dde= (*(resonanceEnergy->find(Z)->second))[iOsc]; |
---|
909 | kineticEnergy1=kinEnergy-dde; |
---|
910 | G4double qs=(*qm)[iOsc]/(1.0+((*qm)[iOsc]/(2.0*electron_mass_c2))); |
---|
911 | G4double q=qs/(std::pow((qs/dde)*(1.0+(0.5*dde/electron_mass_c2)),G4UniformRand())-(0.5*qs/electron_mass_c2)); |
---|
912 | G4double qtrev = q*(q+2.0*electron_mass_c2); |
---|
913 | G4double cpps = kineticEnergy1*(kineticEnergy1+2.0*electron_mass_c2); |
---|
914 | cosThetaPrimary = (cpps+cps-qtrev)/(2.0*cp*std::sqrt(cpps)); |
---|
915 | if (cosThetaPrimary>1.0) cosThetaPrimary=1.0; |
---|
916 | //Energy and emission angle of the delta ray |
---|
917 | kks = (G4int) (*(shellFlag->find(Z)->second))[iOsc]; |
---|
918 | //kks > 4 means that we are in an outer shell |
---|
919 | if (kks>4) |
---|
920 | energySecondary=dde; |
---|
921 | else |
---|
922 | energySecondary=dde-(*(ionizationEnergy->find(Z)->second))[iOsc]; |
---|
923 | |
---|
924 | cosThetaSecondary = 0.5*(dde*(kinEnergy+rb-dde)+qtrev)/std::sqrt(cps*qtrev); |
---|
925 | if (cosThetaSecondary>1.0) cosThetaSecondary=1.0; |
---|
926 | } |
---|
927 | else if (sampledInteraction == 2) //Hard distant transverse collisions |
---|
928 | { |
---|
929 | dde=(*(resonanceEnergy->find(Z)->second))[iOsc]; |
---|
930 | kineticEnergy1=kinEnergy-dde; |
---|
931 | cosThetaPrimary=1.0; |
---|
932 | //Energy and emission angle of the delta ray |
---|
933 | kks = (G4int) (*(shellFlag->find(Z)->second))[iOsc]; |
---|
934 | if (kks>4) |
---|
935 | { |
---|
936 | energySecondary=dde; |
---|
937 | } |
---|
938 | else |
---|
939 | { |
---|
940 | energySecondary=dde-(*(ionizationEnergy->find(Z)->second))[iOsc]; |
---|
941 | } |
---|
942 | cosThetaSecondary = 1.0; |
---|
943 | } |
---|
944 | else if (sampledInteraction == 3 || sampledInteraction == 4) //Close interaction |
---|
945 | { |
---|
946 | if (sampledInteraction == 4) //interaction with inner shells |
---|
947 | { |
---|
948 | UII=0.0; |
---|
949 | rkc = cutoffEnergy/kinEnergy; |
---|
950 | iOsc = -1; |
---|
951 | } |
---|
952 | else |
---|
953 | { |
---|
954 | kks = (G4int) (*(shellFlag->find(Z)->second))[iOsc]; |
---|
955 | if (kks > 4) { |
---|
956 | UII=0.0; |
---|
957 | } |
---|
958 | else |
---|
959 | { |
---|
960 | UII = (*(ionizationEnergy->find(Z)->second))[iOsc]; |
---|
961 | } |
---|
962 | rkc = (*(resonanceEnergy->find(Z)->second))[iOsc]/kinEnergy; |
---|
963 | } |
---|
964 | G4double A = 0.5*amol; |
---|
965 | G4double arkc = A*0.5*rkc; |
---|
966 | G4double phi,rk2,rk,rkf; |
---|
967 | do{ |
---|
968 | G4double fb = (1.0+arkc)*G4UniformRand(); |
---|
969 | if (fb<1.0) |
---|
970 | { |
---|
971 | rk=rkc/(1.0-fb*(1.0-(rkc*2.0))); |
---|
972 | } |
---|
973 | else{ |
---|
974 | rk = rkc+(fb-1.0)*(0.5-rkc)/arkc; |
---|
975 | } |
---|
976 | rk2 = rk*rk; |
---|
977 | rkf = rk/(1.0-rk); |
---|
978 | phi = 1.0+(rkf*rkf)-rkf+amol*(rk2+rkf); |
---|
979 | }while ((G4UniformRand()*(1.0+A*rk2)) > phi); |
---|
980 | //Energy and scattering angle (primary electron); |
---|
981 | kineticEnergy1 = kinEnergy*(1.0-rk); |
---|
982 | cosThetaPrimary = std::sqrt(kineticEnergy1*rb/(kinEnergy*(rb-(rk*kinEnergy)))); |
---|
983 | //Energy and scattering angle of the delta ray |
---|
984 | energySecondary = kinEnergy-kineticEnergy1-UII; |
---|
985 | cosThetaSecondary = std::sqrt(rk*kinEnergy*rb/(kinEnergy*(rk*kinEnergy+2.0*electron_mass_c2))); |
---|
986 | } |
---|
987 | else |
---|
988 | { |
---|
989 | G4String excep = "G4PenelopeIonisationModel - Error in the calculation of the final state"; |
---|
990 | G4Exception(excep); |
---|
991 | } |
---|
992 | |
---|
993 | delete qm; |
---|
994 | delete cumulHardCS; |
---|
995 | |
---|
996 | delete typeOfInteraction; |
---|
997 | delete nbOfLevel; |
---|
998 | |
---|
999 | return; |
---|
1000 | } |
---|
1001 | |
---|
1002 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
1003 | |
---|
1004 | void |
---|
1005 | G4PenelopeIonisationModel::CalculateDiscreteForPositrons(G4double kinEnergy, |
---|
1006 | G4double cutoffEnergy, |
---|
1007 | G4int Z, |
---|
1008 | G4double electronVolumeDensity) |
---|
1009 | { |
---|
1010 | kineticEnergy1=kinEnergy; |
---|
1011 | cosThetaPrimary=1.0; |
---|
1012 | energySecondary=0.0; |
---|
1013 | cosThetaSecondary=1.0; |
---|
1014 | iOsc=-1; |
---|
1015 | //constants |
---|
1016 | G4double rb=kinEnergy+2.0*electron_mass_c2; |
---|
1017 | G4double gamma = 1.0+kinEnergy/electron_mass_c2; |
---|
1018 | G4double gamma2 = gamma*gamma; |
---|
1019 | G4double beta2 = (gamma2-1.0)/gamma2; |
---|
1020 | G4double amol = (gamma-1.0)*(gamma-1.0)/gamma2; |
---|
1021 | G4double cps = kinEnergy*rb; |
---|
1022 | G4double cp = std::sqrt(cps); |
---|
1023 | G4double help = (gamma+1.0)*(gamma+1.0); |
---|
1024 | G4double bha1 = amol*(2.0*help-1.0)/(gamma2-1.0); |
---|
1025 | G4double bha2 = amol*(3.0+1.0/help); |
---|
1026 | G4double bha3 = amol*2.0*gamma*(gamma-1.0)/help; |
---|
1027 | G4double bha4 = amol*(gamma-1.0)*(gamma-1.0)/help; |
---|
1028 | |
---|
1029 | G4double delta = CalculateDeltaFermi(kinEnergy,Z,electronVolumeDensity); |
---|
1030 | G4double distantTransvCS0 = std::max(std::log(gamma2)-beta2-delta,0.0); |
---|
1031 | |
---|
1032 | G4double rl,rl1; |
---|
1033 | |
---|
1034 | if (cutoffEnergy > kinEnergy) return; //delta rays are not generated |
---|
1035 | |
---|
1036 | G4DataVector* qm = new G4DataVector(); |
---|
1037 | G4DataVector* cumulHardCS = new G4DataVector(); |
---|
1038 | std::vector<G4int> *typeOfInteraction = new std::vector<G4int>; |
---|
1039 | G4DataVector* nbOfLevel = new G4DataVector(); |
---|
1040 | |
---|
1041 | |
---|
1042 | //Hard close collisions with outer shells |
---|
1043 | G4double wmaxc = kinEnergy; |
---|
1044 | G4double closeCS0 = 0.0; |
---|
1045 | G4double closeCS = 0.0; |
---|
1046 | if (cutoffEnergy>0.1*eV) |
---|
1047 | { |
---|
1048 | rl=cutoffEnergy/kinEnergy; |
---|
1049 | rl1=1.0-rl; |
---|
1050 | if (rl < 1.0) |
---|
1051 | closeCS0 = (((1.0/rl)-1.0) + bha1*std::log(rl) + bha2*rl1 |
---|
1052 | + (bha3/2.0)*((rl*rl)-1.0) |
---|
1053 | + (bha4/3.0)*(1.0-(rl*rl*rl)))/kinEnergy; |
---|
1054 | } |
---|
1055 | |
---|
1056 | // Cross sections for the different oscillators |
---|
1057 | |
---|
1058 | // totalHardCS contains the cumulative hard interaction cross section for the different |
---|
1059 | // excitable levels and the different interaction channels (close, distant, etc.), |
---|
1060 | // i.e. |
---|
1061 | // cumulHardCS[0] = 0.0 |
---|
1062 | // cumulHardCS[1] = 1st excitable level (distant longitudinal only) |
---|
1063 | // cumulHardCS[2] = 1st excitable level (distant longitudinal + transverse) |
---|
1064 | // cumulHardCS[3] = 1st excitable level (distant longitudinal + transverse + close) |
---|
1065 | // cumulHardCS[4] = 1st excitable level (all channels) + 2nd excitable level (distant long only) |
---|
1066 | // etc. |
---|
1067 | // This is used for sampling the atomic level which is ionised and the channel of the |
---|
1068 | // interaction. |
---|
1069 | // |
---|
1070 | // For each index iFill of the cumulHardCS vector, |
---|
1071 | // nbOfLevel[iFill] contains the current excitable atomic level and |
---|
1072 | // typeOfInteraction[iFill] contains the current interaction channel, with the legenda: |
---|
1073 | // 1 = distant longitudinal interaction |
---|
1074 | // 2 = distant transverse interaction |
---|
1075 | // 3 = close collision |
---|
1076 | // 4 = close collision with outer shells (in this case nbOfLevel < 0 --> no binding energy) |
---|
1077 | |
---|
1078 | |
---|
1079 | G4int nOscil = ionizationEnergy->find(Z)->second->size(); |
---|
1080 | G4double totalHardCS = 0.0; |
---|
1081 | G4double involvedElectrons = 0.0; |
---|
1082 | for (G4int i=0;i<nOscil;i++){ |
---|
1083 | G4double wi = (*(resonanceEnergy->find(Z)->second))[i]; |
---|
1084 | G4int occupNb = (G4int) (*(occupationNumber->find(Z)->second))[i]; |
---|
1085 | //Distant excitations |
---|
1086 | if (wi>cutoffEnergy && wi<kinEnergy) |
---|
1087 | { |
---|
1088 | if (wi>(1e-6*kinEnergy)){ |
---|
1089 | G4double cpp=std::sqrt((kinEnergy-wi)*(kinEnergy-wi+2.0*electron_mass_c2)); |
---|
1090 | qm->push_back(std::sqrt((cp-cpp)*(cp-cpp)+ electron_mass_c2 * electron_mass_c2)-electron_mass_c2); |
---|
1091 | } |
---|
1092 | else |
---|
1093 | { |
---|
1094 | qm->push_back(wi*wi/(beta2+2.0*electron_mass_c2)); |
---|
1095 | } |
---|
1096 | //verificare che quando arriva qui il vettore ha SEMPRE l'i-esimo elemento |
---|
1097 | if ((*qm)[i] < wi) |
---|
1098 | { |
---|
1099 | |
---|
1100 | G4double distantLongitCS = occupNb*std::log(wi*((*qm)[i]+2.0*electron_mass_c2)/ |
---|
1101 | ((*qm)[i]*(wi+2.0*electron_mass_c2)))/wi; |
---|
1102 | cumulHardCS->push_back(totalHardCS); |
---|
1103 | typeOfInteraction->push_back(1); //distant longitudinal |
---|
1104 | nbOfLevel->push_back((G4double) i); //only excitable level are counted |
---|
1105 | totalHardCS += distantLongitCS; |
---|
1106 | |
---|
1107 | G4double distantTransvCS = occupNb*distantTransvCS0/wi; |
---|
1108 | |
---|
1109 | cumulHardCS->push_back(totalHardCS); |
---|
1110 | typeOfInteraction->push_back(2); //distant tranverse |
---|
1111 | nbOfLevel->push_back((G4double) i); |
---|
1112 | totalHardCS += distantTransvCS; |
---|
1113 | } |
---|
1114 | } |
---|
1115 | else |
---|
1116 | qm->push_back(wi); |
---|
1117 | |
---|
1118 | //close collisions |
---|
1119 | if(wi < wmaxc) |
---|
1120 | { |
---|
1121 | if (wi < cutoffEnergy) { |
---|
1122 | involvedElectrons += occupNb; |
---|
1123 | } |
---|
1124 | else |
---|
1125 | { |
---|
1126 | rl=wi/kinEnergy; |
---|
1127 | rl1=1.0-rl; |
---|
1128 | closeCS = occupNb*(((1.0/rl)-1.0)+bha1*std::log(rl)+bha2*rl1 |
---|
1129 | + (bha3/2.0)*((rl*rl)-1.0) |
---|
1130 | + (bha4/3.0)*(1.0-(rl*rl*rl)))/kinEnergy; |
---|
1131 | cumulHardCS->push_back(totalHardCS); |
---|
1132 | typeOfInteraction->push_back(3); //close |
---|
1133 | nbOfLevel->push_back((G4double) i); |
---|
1134 | totalHardCS += closeCS; |
---|
1135 | } |
---|
1136 | } |
---|
1137 | } // loop on the levels |
---|
1138 | |
---|
1139 | cumulHardCS->push_back(totalHardCS); |
---|
1140 | typeOfInteraction->push_back(4); //close interaction with outer shells |
---|
1141 | nbOfLevel->push_back(-1.0); |
---|
1142 | totalHardCS += involvedElectrons*closeCS0; |
---|
1143 | cumulHardCS->push_back(totalHardCS); //this is the final value of the totalHardCS |
---|
1144 | |
---|
1145 | if (totalHardCS < 1e-30) { |
---|
1146 | kineticEnergy1=kinEnergy; |
---|
1147 | cosThetaPrimary=1.0; |
---|
1148 | energySecondary=0.0; |
---|
1149 | cosThetaSecondary=0.0; |
---|
1150 | iOsc=-1; |
---|
1151 | delete qm; |
---|
1152 | delete cumulHardCS; |
---|
1153 | delete typeOfInteraction; |
---|
1154 | delete nbOfLevel; |
---|
1155 | return; |
---|
1156 | } |
---|
1157 | |
---|
1158 | //Selection of the active oscillator on the basis of the cumulative cross sections |
---|
1159 | G4double TST = totalHardCS*G4UniformRand(); |
---|
1160 | G4int is=0; |
---|
1161 | G4int js= nbOfLevel->size(); |
---|
1162 | do{ |
---|
1163 | G4int it=(is+js)/2; |
---|
1164 | if (TST > (*cumulHardCS)[it]) is=it; |
---|
1165 | if (TST <= (*cumulHardCS)[it]) js=it; |
---|
1166 | }while((js-is) > 1); |
---|
1167 | |
---|
1168 | G4double UII=0.0; |
---|
1169 | G4double rkc=cutoffEnergy/kinEnergy; |
---|
1170 | G4double dde; |
---|
1171 | G4int kks; |
---|
1172 | |
---|
1173 | G4int sampledInteraction = (*typeOfInteraction)[is]; |
---|
1174 | iOsc = (G4int) (*nbOfLevel)[is]; |
---|
1175 | |
---|
1176 | //Generates the final state according to the sampled level and |
---|
1177 | //interaction channel |
---|
1178 | |
---|
1179 | if (sampledInteraction == 1) //Hard distant longitudinal collisions |
---|
1180 | { |
---|
1181 | dde= (*(resonanceEnergy->find(Z)->second))[iOsc]; |
---|
1182 | kineticEnergy1=kinEnergy-dde; |
---|
1183 | G4double qs=(*qm)[iOsc]/(1.0+((*qm)[iOsc]/(2.0*electron_mass_c2))); |
---|
1184 | G4double q=qs/(std::pow((qs/dde)*(1.0+(0.5*dde/electron_mass_c2)),G4UniformRand())-(0.5*qs/electron_mass_c2)); |
---|
1185 | G4double qtrev = q*(q+2.0*electron_mass_c2); |
---|
1186 | G4double cpps = kineticEnergy1*(kineticEnergy1+2.0*electron_mass_c2); |
---|
1187 | cosThetaPrimary = (cpps+cps-qtrev)/(2.0*cp*std::sqrt(cpps)); |
---|
1188 | if (cosThetaPrimary>1.0) cosThetaPrimary=1.0; |
---|
1189 | //Energy and emission angle of the delta ray |
---|
1190 | kks = (G4int) (*(shellFlag->find(Z)->second))[iOsc]; |
---|
1191 | if (kks>4) |
---|
1192 | energySecondary=dde; |
---|
1193 | |
---|
1194 | else |
---|
1195 | energySecondary=dde-(*(ionizationEnergy->find(Z)->second))[iOsc]; |
---|
1196 | cosThetaSecondary = 0.5*(dde*(kinEnergy+rb-dde)+qtrev)/std::sqrt(cps*qtrev); |
---|
1197 | if (cosThetaSecondary>1.0) cosThetaSecondary=1.0; |
---|
1198 | } |
---|
1199 | else if (sampledInteraction == 2) //Hard distant transverse collisions |
---|
1200 | { |
---|
1201 | dde=(*(resonanceEnergy->find(Z)->second))[iOsc]; |
---|
1202 | kineticEnergy1=kinEnergy-dde; |
---|
1203 | cosThetaPrimary=1.0; |
---|
1204 | //Energy and emission angle of the delta ray |
---|
1205 | kks = (G4int) (*(shellFlag->find(Z)->second))[iOsc]; |
---|
1206 | if (kks>4) |
---|
1207 | { |
---|
1208 | energySecondary=dde; |
---|
1209 | } |
---|
1210 | else |
---|
1211 | { |
---|
1212 | energySecondary=dde-(*(ionizationEnergy->find(Z)->second))[iOsc]; |
---|
1213 | } |
---|
1214 | cosThetaSecondary = 1.0; |
---|
1215 | } |
---|
1216 | else if (sampledInteraction == 3 || sampledInteraction == 4) //Close interaction |
---|
1217 | { |
---|
1218 | if (sampledInteraction == 4) //interaction with inner shells |
---|
1219 | { |
---|
1220 | UII=0.0; |
---|
1221 | rkc = cutoffEnergy/kinEnergy; |
---|
1222 | iOsc = -1; |
---|
1223 | } |
---|
1224 | else |
---|
1225 | { |
---|
1226 | kks = (G4int) (*(shellFlag->find(Z)->second))[iOsc]; |
---|
1227 | if (kks > 4) { |
---|
1228 | UII=0.0; |
---|
1229 | } |
---|
1230 | else |
---|
1231 | { |
---|
1232 | UII = (*(ionizationEnergy->find(Z)->second))[iOsc]; |
---|
1233 | } |
---|
1234 | rkc = (*(resonanceEnergy->find(Z)->second))[iOsc]/kinEnergy; |
---|
1235 | } |
---|
1236 | G4double phi,rk; |
---|
1237 | do{ |
---|
1238 | rk=rkc/(1.0-G4UniformRand()*(1.0-rkc)); |
---|
1239 | phi = 1.0-rk*(bha1-rk*(bha2-rk*(bha3-bha4*rk))); |
---|
1240 | }while ( G4UniformRand() > phi); |
---|
1241 | //Energy and scattering angle (primary electron); |
---|
1242 | kineticEnergy1 = kinEnergy*(1.0-rk); |
---|
1243 | cosThetaPrimary = std::sqrt(kineticEnergy1*rb/(kinEnergy*(rb-(rk*kinEnergy)))); |
---|
1244 | //Energy and scattering angle of the delta ray |
---|
1245 | energySecondary = kinEnergy-kineticEnergy1-UII; |
---|
1246 | cosThetaSecondary = std::sqrt(rk*kinEnergy*rb/(kinEnergy*(rk*kinEnergy+2.0*electron_mass_c2))); |
---|
1247 | } |
---|
1248 | else |
---|
1249 | { |
---|
1250 | G4String excep = "G4PenelopeIonisationModel - Error in the calculation of the final state"; |
---|
1251 | G4Exception(excep); |
---|
1252 | } |
---|
1253 | |
---|
1254 | delete qm; |
---|
1255 | delete cumulHardCS; |
---|
1256 | delete typeOfInteraction; |
---|
1257 | delete nbOfLevel; |
---|
1258 | |
---|
1259 | return; |
---|
1260 | } |
---|
1261 | |
---|
1262 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
1263 | |
---|
1264 | G4double |
---|
1265 | G4PenelopeIonisationModel::CalculateCrossSectionsRatio(G4double kinEnergy, |
---|
1266 | G4double cutoffEnergy, |
---|
1267 | G4int Z, G4double electronVolumeDensity, |
---|
1268 | const G4ParticleDefinition* theParticle) |
---|
1269 | { |
---|
1270 | //Constants |
---|
1271 | G4double gamma = 1.0+kinEnergy/electron_mass_c2; |
---|
1272 | G4double gamma2 = gamma*gamma; |
---|
1273 | G4double beta2 = (gamma2-1.0)/gamma2; |
---|
1274 | G4double constant = pi*classic_electr_radius*classic_electr_radius*2.0*electron_mass_c2/beta2; |
---|
1275 | G4double delta = CalculateDeltaFermi(kinEnergy,Z,electronVolumeDensity); |
---|
1276 | G4int nbOsc = (G4int) resonanceEnergy->find(Z)->second->size(); |
---|
1277 | G4double softCS = 0.0; |
---|
1278 | G4double hardCS = 0.0; |
---|
1279 | for (G4int i=0;i<nbOsc;i++){ |
---|
1280 | G4double resEnergy = (*(resonanceEnergy->find(Z)->second))[i]; |
---|
1281 | std::pair<G4double,G4double> SoftAndHardXS(0.,0.); |
---|
1282 | if (theParticle == G4Electron::Electron()) |
---|
1283 | SoftAndHardXS = CrossSectionsRatioForElectrons(kinEnergy,resEnergy,delta,cutoffEnergy); |
---|
1284 | else if (theParticle == G4Positron::Positron()) |
---|
1285 | SoftAndHardXS = CrossSectionsRatioForPositrons(kinEnergy,resEnergy,delta,cutoffEnergy); |
---|
1286 | G4double occupNb = (*(occupationNumber->find(Z)->second))[i]; |
---|
1287 | softCS += occupNb*constant*SoftAndHardXS.first; |
---|
1288 | hardCS += occupNb*constant*SoftAndHardXS.second; |
---|
1289 | } |
---|
1290 | G4double ratio = 0.0; |
---|
1291 | |
---|
1292 | if (softCS+hardCS) ratio = (hardCS)/(softCS+hardCS); |
---|
1293 | return ratio; |
---|
1294 | } |
---|
1295 | |
---|
1296 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
1297 | |
---|
1298 | std::pair<G4double,G4double> |
---|
1299 | G4PenelopeIonisationModel::CrossSectionsRatioForElectrons(G4double kineticEnergy, |
---|
1300 | G4double resEnergy, |
---|
1301 | G4double densityCorrection, |
---|
1302 | G4double cutoffEnergy) |
---|
1303 | { |
---|
1304 | std::pair<G4double,G4double> theResult(0.,0.); |
---|
1305 | if (kineticEnergy < resEnergy) return theResult; |
---|
1306 | |
---|
1307 | //Calculate constants |
---|
1308 | G4double gamma = 1.0+kineticEnergy/electron_mass_c2; |
---|
1309 | G4double gamma2 = gamma*gamma; |
---|
1310 | G4double beta2 = (gamma2-1.0)/gamma2; |
---|
1311 | G4double cps = kineticEnergy*(kineticEnergy+2.0*electron_mass_c2); |
---|
1312 | G4double amol = (kineticEnergy/(kineticEnergy+electron_mass_c2)) * (kineticEnergy/(kineticEnergy+electron_mass_c2)) ; |
---|
1313 | G4double hardCont = 0.0; |
---|
1314 | G4double softCont = 0.0; |
---|
1315 | |
---|
1316 | //Distant interactions |
---|
1317 | G4double cp1s = (kineticEnergy-resEnergy)*(kineticEnergy-resEnergy+2.0*electron_mass_c2); |
---|
1318 | G4double cp1 = std::sqrt(cp1s); |
---|
1319 | G4double cp = std::sqrt(cps); |
---|
1320 | G4double sdLong=0.0, sdTrans = 0.0, sdDist=0.0; |
---|
1321 | |
---|
1322 | //Distant longitudinal interactions |
---|
1323 | G4double qm = 0.0; |
---|
1324 | |
---|
1325 | if (resEnergy > kineticEnergy*(1e-6)) |
---|
1326 | { |
---|
1327 | qm = std::sqrt((cp-cp1)*(cp-cp1)+(electron_mass_c2*electron_mass_c2))-electron_mass_c2; |
---|
1328 | } |
---|
1329 | else |
---|
1330 | { |
---|
1331 | qm = resEnergy*resEnergy/(beta2*2.0*electron_mass_c2); |
---|
1332 | qm = qm*(1.0-0.5*qm/electron_mass_c2); |
---|
1333 | } |
---|
1334 | |
---|
1335 | if (qm < resEnergy) |
---|
1336 | { |
---|
1337 | sdLong = std::log(resEnergy*(qm+2.0*electron_mass_c2)/(qm*(resEnergy+2.0*electron_mass_c2))); |
---|
1338 | } |
---|
1339 | else |
---|
1340 | { |
---|
1341 | sdLong = 0.0; |
---|
1342 | } |
---|
1343 | |
---|
1344 | if (sdLong > 0) { |
---|
1345 | sdTrans = std::max(std::log(gamma2)-beta2-densityCorrection,0.0); |
---|
1346 | sdDist = sdTrans + sdLong; |
---|
1347 | if (cutoffEnergy > resEnergy) |
---|
1348 | { |
---|
1349 | softCont = sdDist/resEnergy; |
---|
1350 | } |
---|
1351 | else |
---|
1352 | { |
---|
1353 | hardCont = sdDist/resEnergy; |
---|
1354 | } |
---|
1355 | } |
---|
1356 | |
---|
1357 | // Close collisions (Moeller's cross section) |
---|
1358 | G4double wl = std::max(cutoffEnergy,resEnergy); |
---|
1359 | G4double wu = 0.5*kineticEnergy; |
---|
1360 | |
---|
1361 | if (wl < (wu-1*eV)) |
---|
1362 | { |
---|
1363 | hardCont += (1.0/(kineticEnergy-wu))-(1.0/(kineticEnergy-wl)) |
---|
1364 | - (1.0/wu)+(1.0/wl) |
---|
1365 | + (1.0-amol)*std::log(((kineticEnergy-wu)*wl)/((kineticEnergy-wl)*wu))/kineticEnergy |
---|
1366 | + amol*(wu-wl)/(kineticEnergy*kineticEnergy); |
---|
1367 | wu=wl; |
---|
1368 | } |
---|
1369 | |
---|
1370 | wl = resEnergy; |
---|
1371 | if (wl > (wu-1*eV)) |
---|
1372 | { |
---|
1373 | theResult.first = softCont; |
---|
1374 | theResult.second = hardCont; |
---|
1375 | return theResult; |
---|
1376 | } |
---|
1377 | softCont += (1.0/(kineticEnergy-wu))-(1.0/(kineticEnergy-wl)) |
---|
1378 | - (1.0/wu)+(1.0/wl) |
---|
1379 | + (1.0-amol)*std::log(((kineticEnergy-wu)*wl)/((kineticEnergy-wl)*wu))/kineticEnergy |
---|
1380 | + amol*(wu-wl)/(kineticEnergy*kineticEnergy); |
---|
1381 | theResult.first = softCont; |
---|
1382 | theResult.second = hardCont; |
---|
1383 | return theResult; |
---|
1384 | } |
---|
1385 | |
---|
1386 | |
---|
1387 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
1388 | |
---|
1389 | std::pair<G4double,G4double> |
---|
1390 | G4PenelopeIonisationModel::CrossSectionsRatioForPositrons(G4double kineticEnergy, |
---|
1391 | G4double resEnergy, |
---|
1392 | G4double densityCorrection, |
---|
1393 | G4double cutoffEnergy) |
---|
1394 | { |
---|
1395 | |
---|
1396 | std::pair<G4double,G4double> theResult(0.,0.); |
---|
1397 | if (kineticEnergy < resEnergy) return theResult; |
---|
1398 | |
---|
1399 | //Calculate constants |
---|
1400 | G4double gamma = 1.0+kineticEnergy/electron_mass_c2; |
---|
1401 | G4double gamma2 = gamma*gamma; |
---|
1402 | G4double beta2 = (gamma2-1.0)/gamma2; |
---|
1403 | G4double cps = kineticEnergy*(kineticEnergy+2.0*electron_mass_c2); |
---|
1404 | G4double amol = (kineticEnergy/(kineticEnergy+electron_mass_c2)) * (kineticEnergy/(kineticEnergy+electron_mass_c2)) ; |
---|
1405 | G4double help = (gamma+1.0)*(gamma+1.0); |
---|
1406 | G4double bha1 = amol*(2.0*help-1.0)/(gamma2-1.0); |
---|
1407 | G4double bha2 = amol*(3.0+1.0/help); |
---|
1408 | G4double bha3 = amol*2.0*gamma*(gamma-1.0)/help; |
---|
1409 | G4double bha4 = amol*(gamma-1.0)*(gamma-1.0)/help; |
---|
1410 | G4double hardCont = 0.0; |
---|
1411 | G4double softCont = 0.0; |
---|
1412 | |
---|
1413 | //Distant interactions |
---|
1414 | G4double cp1s = (kineticEnergy-resEnergy)*(kineticEnergy-resEnergy+2.0*electron_mass_c2); |
---|
1415 | G4double cp1 = std::sqrt(cp1s); |
---|
1416 | G4double cp = std::sqrt(cps); |
---|
1417 | G4double sdLong=0.0, sdTrans = 0.0, sdDist=0.0; |
---|
1418 | |
---|
1419 | //Distant longitudinal interactions |
---|
1420 | G4double qm = 0.0; |
---|
1421 | |
---|
1422 | if (resEnergy > kineticEnergy*(1e-6)) |
---|
1423 | { |
---|
1424 | qm = std::sqrt((cp-cp1)*(cp-cp1)+(electron_mass_c2*electron_mass_c2))-electron_mass_c2; |
---|
1425 | } |
---|
1426 | else |
---|
1427 | { |
---|
1428 | qm = resEnergy*resEnergy/(beta2*2.0*electron_mass_c2); |
---|
1429 | qm = qm*(1.0-0.5*qm/electron_mass_c2); |
---|
1430 | } |
---|
1431 | |
---|
1432 | if (qm < resEnergy) |
---|
1433 | { |
---|
1434 | sdLong = std::log(resEnergy*(qm+2.0*electron_mass_c2)/(qm*(resEnergy+2.0*electron_mass_c2))); |
---|
1435 | } |
---|
1436 | else |
---|
1437 | { |
---|
1438 | sdLong = 0.0; |
---|
1439 | } |
---|
1440 | |
---|
1441 | if (sdLong > 0) { |
---|
1442 | sdTrans = std::max(std::log(gamma2)-beta2-densityCorrection,0.0); |
---|
1443 | sdDist = sdTrans + sdLong; |
---|
1444 | if (cutoffEnergy > resEnergy) |
---|
1445 | { |
---|
1446 | softCont = sdDist/resEnergy; |
---|
1447 | } |
---|
1448 | else |
---|
1449 | { |
---|
1450 | hardCont = sdDist/resEnergy; |
---|
1451 | } |
---|
1452 | } |
---|
1453 | |
---|
1454 | |
---|
1455 | // Close collisions (Bhabha's cross section) |
---|
1456 | G4double wl = std::max(cutoffEnergy,resEnergy); |
---|
1457 | G4double wu = kineticEnergy; |
---|
1458 | |
---|
1459 | if (wl < (wu-1*eV)) { |
---|
1460 | hardCont += (1.0/wl)-(1.0/wu)-bha1*std::log(wu/wl)/kineticEnergy |
---|
1461 | + bha2*(wu-wl)/(kineticEnergy*kineticEnergy) |
---|
1462 | -bha3*((wu*wu)-(wl*wl))/(2.0*kineticEnergy*kineticEnergy*kineticEnergy) |
---|
1463 | + bha4*((wu*wu*wu)-(wl*wl*wl))/(3.0*kineticEnergy*kineticEnergy*kineticEnergy*kineticEnergy); |
---|
1464 | wu=wl; |
---|
1465 | } |
---|
1466 | wl = resEnergy; |
---|
1467 | if (wl > (wu-1*eV)) |
---|
1468 | { |
---|
1469 | theResult.first = softCont; |
---|
1470 | theResult.second = hardCont; |
---|
1471 | return theResult; |
---|
1472 | } |
---|
1473 | softCont += (1.0/wl)-(1.0/wu)-bha1*std::log(wu/wl)/kineticEnergy |
---|
1474 | + bha2*(wu-wl)/(kineticEnergy*kineticEnergy) |
---|
1475 | -bha3*((wu*wu)-(wl*wl))/(2.0*kineticEnergy*kineticEnergy*kineticEnergy) |
---|
1476 | + bha4*((wu*wu*wu)-(wl*wl*wl))/(3.0*kineticEnergy*kineticEnergy*kineticEnergy*kineticEnergy); |
---|
1477 | |
---|
1478 | |
---|
1479 | theResult.first = softCont; |
---|
1480 | theResult.second = hardCont; |
---|
1481 | return theResult; |
---|
1482 | |
---|
1483 | } |
---|
1484 | |
---|
1485 | G4double G4PenelopeIonisationModel::ComputeStoppingPowerForElectrons(G4double kinEnergy, |
---|
1486 | G4double cutEnergy, |
---|
1487 | G4double deltaFermi, |
---|
1488 | G4double resEnergy) |
---|
1489 | { |
---|
1490 | //Calculate constants |
---|
1491 | G4double gamma = 1.0+kinEnergy/electron_mass_c2; |
---|
1492 | G4double gamma2 = gamma*gamma; |
---|
1493 | G4double beta2 = (gamma2-1.0)/gamma2; |
---|
1494 | G4double cps = kinEnergy*(kinEnergy+2.0*electron_mass_c2); |
---|
1495 | G4double amol = (gamma-1.0)*(gamma-1.0)/gamma2; |
---|
1496 | G4double sPower = 0.0; |
---|
1497 | if (kinEnergy < resEnergy) return sPower; |
---|
1498 | |
---|
1499 | //Distant interactions |
---|
1500 | G4double cp1s = (kinEnergy-resEnergy)*(kinEnergy-resEnergy+2.0*electron_mass_c2); |
---|
1501 | G4double cp1 = std::sqrt(cp1s); |
---|
1502 | G4double cp = std::sqrt(cps); |
---|
1503 | G4double sdLong=0.0, sdTrans = 0.0, sdDist=0.0; |
---|
1504 | |
---|
1505 | //Distant longitudinal interactions |
---|
1506 | G4double qm = 0.0; |
---|
1507 | |
---|
1508 | if (resEnergy > kinEnergy*(1e-6)) |
---|
1509 | { |
---|
1510 | qm = std::sqrt((cp-cp1)*(cp-cp1)+(electron_mass_c2*electron_mass_c2))-electron_mass_c2; |
---|
1511 | } |
---|
1512 | else |
---|
1513 | { |
---|
1514 | qm = resEnergy*resEnergy/(beta2*2.0*electron_mass_c2); |
---|
1515 | qm = qm*(1.0-0.5*qm/electron_mass_c2); |
---|
1516 | } |
---|
1517 | |
---|
1518 | if (qm < resEnergy) |
---|
1519 | sdLong = std::log(resEnergy*(qm+2.0*electron_mass_c2)/(qm*(resEnergy+2.0*electron_mass_c2))); |
---|
1520 | else |
---|
1521 | sdLong = 0.0; |
---|
1522 | |
---|
1523 | if (sdLong > 0) { |
---|
1524 | sdTrans = std::max(std::log(gamma2)-beta2-deltaFermi,0.0); |
---|
1525 | sdDist = sdTrans + sdLong; |
---|
1526 | if (cutEnergy > resEnergy) sPower = sdDist; |
---|
1527 | } |
---|
1528 | |
---|
1529 | |
---|
1530 | // Close collisions (Moeller's cross section) |
---|
1531 | G4double wl = std::max(cutEnergy,resEnergy); |
---|
1532 | G4double wu = 0.5*kinEnergy; |
---|
1533 | |
---|
1534 | if (wl < (wu-1*eV)) wu=wl; |
---|
1535 | wl = resEnergy; |
---|
1536 | if (wl > (wu-1*eV)) return sPower; |
---|
1537 | sPower += std::log(wu/wl)+(kinEnergy/(kinEnergy-wu))-(kinEnergy/(kinEnergy-wl)) |
---|
1538 | + (2.0 - amol)*std::log((kinEnergy-wu)/(kinEnergy-wl)) |
---|
1539 | + amol*((wu*wu)-(wl*wl))/(2.0*kinEnergy*kinEnergy); |
---|
1540 | return sPower; |
---|
1541 | } |
---|
1542 | |
---|
1543 | |
---|
1544 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
1545 | |
---|
1546 | G4double G4PenelopeIonisationModel::ComputeStoppingPowerForPositrons(G4double kinEnergy, |
---|
1547 | G4double cutEnergy, |
---|
1548 | G4double deltaFermi, |
---|
1549 | G4double resEnergy) |
---|
1550 | { |
---|
1551 | //Calculate constants |
---|
1552 | G4double gamma = 1.0+kinEnergy/electron_mass_c2; |
---|
1553 | G4double gamma2 = gamma*gamma; |
---|
1554 | G4double beta2 = (gamma2-1.0)/gamma2; |
---|
1555 | G4double cps = kinEnergy*(kinEnergy+2.0*electron_mass_c2); |
---|
1556 | G4double amol = (kinEnergy/(kinEnergy+electron_mass_c2)) * (kinEnergy/(kinEnergy+electron_mass_c2)); |
---|
1557 | G4double help = (gamma+1.0)*(gamma+1.0); |
---|
1558 | G4double bha1 = amol*(2.0*help-1.0)/(gamma2-1.0); |
---|
1559 | G4double bha2 = amol*(3.0+1.0/help); |
---|
1560 | G4double bha3 = amol*2.0*gamma*(gamma-1.0)/help; |
---|
1561 | G4double bha4 = amol*(gamma-1.0)*(gamma-1.0)/help; |
---|
1562 | |
---|
1563 | G4double sPower = 0.0; |
---|
1564 | if (kinEnergy < resEnergy) return sPower; |
---|
1565 | |
---|
1566 | //Distant interactions |
---|
1567 | G4double cp1s = (kinEnergy-resEnergy)*(kinEnergy-resEnergy+2.0*electron_mass_c2); |
---|
1568 | G4double cp1 = std::sqrt(cp1s); |
---|
1569 | G4double cp = std::sqrt(cps); |
---|
1570 | G4double sdLong=0.0, sdTrans = 0.0, sdDist=0.0; |
---|
1571 | |
---|
1572 | //Distant longitudinal interactions |
---|
1573 | G4double qm = 0.0; |
---|
1574 | |
---|
1575 | if (resEnergy > kinEnergy*(1e-6)) |
---|
1576 | { |
---|
1577 | qm = std::sqrt((cp-cp1)*(cp-cp1)+(electron_mass_c2*electron_mass_c2))-electron_mass_c2; |
---|
1578 | } |
---|
1579 | else |
---|
1580 | { |
---|
1581 | qm = resEnergy*resEnergy/(beta2*2.0*electron_mass_c2); |
---|
1582 | qm = qm*(1.0-0.5*qm/electron_mass_c2); |
---|
1583 | } |
---|
1584 | |
---|
1585 | if (qm < resEnergy) |
---|
1586 | sdLong = std::log(resEnergy*(qm+2.0*electron_mass_c2)/(qm*(resEnergy+2.0*electron_mass_c2))); |
---|
1587 | else |
---|
1588 | sdLong = 0.0; |
---|
1589 | |
---|
1590 | if (sdLong > 0) { |
---|
1591 | sdTrans = std::max(std::log(gamma2)-beta2-deltaFermi,0.0); |
---|
1592 | sdDist = sdTrans + sdLong; |
---|
1593 | if (cutEnergy > resEnergy) sPower = sdDist; |
---|
1594 | } |
---|
1595 | |
---|
1596 | |
---|
1597 | // Close collisions (Bhabha's cross section) |
---|
1598 | G4double wl = std::max(cutEnergy,resEnergy); |
---|
1599 | G4double wu = kinEnergy; |
---|
1600 | |
---|
1601 | if (wl < (wu-1*eV)) wu=wl; |
---|
1602 | wl = resEnergy; |
---|
1603 | if (wl > (wu-1*eV)) return sPower; |
---|
1604 | sPower += std::log(wu/wl)-bha1*(wu-wl)/kinEnergy |
---|
1605 | + bha2*((wu*wu)-(wl*wl))/(2.0*kinEnergy*kinEnergy) |
---|
1606 | - bha3*((wu*wu*wu)-(wl*wl*wl))/(3.0*kinEnergy*kinEnergy*kinEnergy) |
---|
1607 | + bha4*((wu*wu*wu*wu)-(wl*wl*wl*wl))/(4.0*kinEnergy*kinEnergy*kinEnergy*kinEnergy); |
---|
1608 | return sPower; |
---|
1609 | } |
---|
1610 | |
---|
1611 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo... |
---|
1612 | |
---|
1613 | /* Notice: the methods here above are only temporary. They will become obsolete in a while */ |
---|
1614 | |
---|
1615 | #include "G4VDataSetAlgorithm.hh" |
---|
1616 | #include "G4LinLogLogInterpolation.hh" |
---|
1617 | #include "G4CompositeEMDataSet.hh" |
---|
1618 | #include "G4EMDataSet.hh" |
---|
1619 | |
---|
1620 | std::vector<G4VEMDataSet*>* |
---|
1621 | G4PenelopeIonisationModel::BuildCrossSectionTable(const G4ParticleDefinition* theParticle) |
---|
1622 | { |
---|
1623 | std::vector<G4VEMDataSet*>* set = new std::vector<G4VEMDataSet*>; |
---|
1624 | |
---|
1625 | size_t nOfBins = 200; |
---|
1626 | G4PhysicsLogVector* theLogVector = new G4PhysicsLogVector(LowEnergyLimit(), |
---|
1627 | HighEnergyLimit(), |
---|
1628 | nOfBins); |
---|
1629 | G4DataVector* energies; |
---|
1630 | G4DataVector* cs; |
---|
1631 | |
---|
1632 | const G4ProductionCutsTable* theCoupleTable= |
---|
1633 | G4ProductionCutsTable::GetProductionCutsTable(); |
---|
1634 | size_t numOfCouples = theCoupleTable->GetTableSize(); |
---|
1635 | |
---|
1636 | |
---|
1637 | for (size_t m=0; m<numOfCouples; m++) |
---|
1638 | { |
---|
1639 | const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(m); |
---|
1640 | const G4Material* material= couple->GetMaterial(); |
---|
1641 | const G4ElementVector* elementVector = material->GetElementVector(); |
---|
1642 | const G4double* nAtomsPerVolume = material->GetAtomicNumDensityVector(); |
---|
1643 | G4double electronVolumeDensity = |
---|
1644 | material->GetTotNbOfElectPerVolume(); //electron density |
---|
1645 | G4int nElements = material->GetNumberOfElements(); |
---|
1646 | |
---|
1647 | G4double tcut = (*(theCoupleTable->GetEnergyCutsVector(1)))[couple->GetIndex()]; |
---|
1648 | |
---|
1649 | G4VDataSetAlgorithm* algo = new G4LinLogLogInterpolation(); |
---|
1650 | |
---|
1651 | G4VEMDataSet* setForMat = new G4CompositeEMDataSet(algo,1.,1.); |
---|
1652 | |
---|
1653 | for (G4int i=0; i<nElements; i++) |
---|
1654 | { |
---|
1655 | G4int iZ = (G4int) (*elementVector)[i]->GetZ(); |
---|
1656 | energies = new G4DataVector; |
---|
1657 | cs = new G4DataVector; |
---|
1658 | G4double density = nAtomsPerVolume[i]; |
---|
1659 | for (size_t bin=0; bin<nOfBins; bin++) |
---|
1660 | { |
---|
1661 | G4double e = theLogVector->GetLowEdgeEnergy(bin); |
---|
1662 | energies->push_back(e); |
---|
1663 | G4double value = 0.0; |
---|
1664 | if(e > tcut) |
---|
1665 | { |
---|
1666 | G4double ratio = CalculateCrossSectionsRatio(e,tcut,iZ, |
---|
1667 | electronVolumeDensity, |
---|
1668 | theParticle); |
---|
1669 | value = crossSectionHandler->FindValue(iZ,e)*ratio*density; |
---|
1670 | } |
---|
1671 | cs->push_back(value); |
---|
1672 | } |
---|
1673 | G4VDataSetAlgorithm* algo1 = algo->Clone(); |
---|
1674 | G4VEMDataSet* elSet = new G4EMDataSet(i,energies,cs,algo1,1.,1.); |
---|
1675 | setForMat->AddComponent(elSet); |
---|
1676 | } |
---|
1677 | set->push_back(setForMat); |
---|
1678 | } |
---|
1679 | return set; |
---|
1680 | } |
---|
1681 | |
---|
1682 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo... |
---|
1683 | |
---|
1684 | G4int G4PenelopeIonisationModel::SampleRandomAtom(const G4MaterialCutsCouple* couple, |
---|
1685 | G4double e) const |
---|
1686 | { |
---|
1687 | // Select randomly an element within the material, according to the weight |
---|
1688 | // determined by the cross sections in the data set |
---|
1689 | |
---|
1690 | const G4Material* material = couple->GetMaterial(); |
---|
1691 | G4int nElements = material->GetNumberOfElements(); |
---|
1692 | |
---|
1693 | // Special case: the material consists of one element |
---|
1694 | if (nElements == 1) |
---|
1695 | { |
---|
1696 | G4int Z = (G4int) material->GetZ(); |
---|
1697 | return Z; |
---|
1698 | } |
---|
1699 | |
---|
1700 | // Composite material |
---|
1701 | const G4ElementVector* elementVector = material->GetElementVector(); |
---|
1702 | size_t materialIndex = couple->GetIndex(); |
---|
1703 | |
---|
1704 | G4VEMDataSet* materialSet = (*theXSTable)[materialIndex]; |
---|
1705 | G4double materialCrossSection0 = 0.0; |
---|
1706 | G4DataVector cross; |
---|
1707 | cross.clear(); |
---|
1708 | for ( G4int i=0; i < nElements; i++ ) |
---|
1709 | { |
---|
1710 | G4double cr = materialSet->GetComponent(i)->FindValue(e); |
---|
1711 | materialCrossSection0 += cr; |
---|
1712 | cross.push_back(materialCrossSection0); |
---|
1713 | } |
---|
1714 | |
---|
1715 | G4double random = G4UniformRand() * materialCrossSection0; |
---|
1716 | |
---|
1717 | for (G4int k=0 ; k < nElements ; k++ ) |
---|
1718 | { |
---|
1719 | if (random <= cross[k]) return (G4int) (*elementVector)[k]->GetZ(); |
---|
1720 | } |
---|
1721 | // It should never get here |
---|
1722 | return 0; |
---|
1723 | } |
---|
1724 | |
---|
1725 | |
---|
1726 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
1727 | |
---|
1728 | void G4PenelopeIonisationModel::ActivateAuger(G4bool augerbool) |
---|
1729 | { |
---|
1730 | if (!DeexcitationFlag() && augerbool) |
---|
1731 | { |
---|
1732 | G4cout << "WARNING - G4PenelopeIonisationModel" << G4endl; |
---|
1733 | G4cout << "The use of the Atomic Deexcitation Manager is set to false " << G4endl; |
---|
1734 | G4cout << "Therefore, Auger electrons will be not generated anyway" << G4endl; |
---|
1735 | } |
---|
1736 | deexcitationManager.ActivateAugerElectronProduction(augerbool); |
---|
1737 | if (verboseLevel > 1) |
---|
1738 | G4cout << "Auger production set to " << augerbool << G4endl; |
---|
1739 | } |
---|
1740 | |
---|
1741 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|