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: G4PenelopeRayleighModel.cc,v 1.7 2009/12/21 12:49:01 pandola Exp $ |
---|
27 | // GEANT4 tag $Name: geant4-09-04-beta-01 $ |
---|
28 | // |
---|
29 | // Author: Luciano Pandola |
---|
30 | // |
---|
31 | // History: |
---|
32 | // -------- |
---|
33 | // 14 Oct 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 | // 19 May 2009 L Pandola Explicitely set to zero pointers deleted in |
---|
38 | // PrepareConstants(), since they might be checked later on |
---|
39 | // 18 Dec 2009 L Pandola Added a dummy ComputeCrossSectionPerAtom() method issueing a |
---|
40 | // warning if users try to access atomic cross sections via |
---|
41 | // G4EmCalculator |
---|
42 | // |
---|
43 | |
---|
44 | #include "G4PenelopeRayleighModel.hh" |
---|
45 | #include "G4ParticleDefinition.hh" |
---|
46 | #include "G4MaterialCutsCouple.hh" |
---|
47 | #include "G4ProductionCutsTable.hh" |
---|
48 | #include "G4DynamicParticle.hh" |
---|
49 | #include "G4PhysicsTable.hh" |
---|
50 | #include "G4ElementTable.hh" |
---|
51 | #include "G4Element.hh" |
---|
52 | #include "G4PenelopeIntegrator.hh" |
---|
53 | |
---|
54 | |
---|
55 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
56 | |
---|
57 | |
---|
58 | G4PenelopeRayleighModel::G4PenelopeRayleighModel(const G4ParticleDefinition*, |
---|
59 | const G4String& nam) |
---|
60 | :G4VEmModel(nam),samplingFunction_x(0),samplingFunction_xNoLog(0), |
---|
61 | theMaterial(0), |
---|
62 | isInitialised(false) |
---|
63 | { |
---|
64 | fIntrinsicLowEnergyLimit = 100.0*eV; |
---|
65 | fIntrinsicHighEnergyLimit = 100.0*GeV; |
---|
66 | // SetLowEnergyLimit(fIntrinsicLowEnergyLimit); |
---|
67 | SetHighEnergyLimit(fIntrinsicHighEnergyLimit); |
---|
68 | // |
---|
69 | verboseLevel= 0; |
---|
70 | // Verbosity scale: |
---|
71 | // 0 = nothing |
---|
72 | // 1 = warning for energy non-conservation |
---|
73 | // 2 = details of energy budget |
---|
74 | // 3 = calculation of cross sections, file openings, sampling of atoms |
---|
75 | // 4 = entering in methods |
---|
76 | |
---|
77 | PrepareConstants(); |
---|
78 | } |
---|
79 | |
---|
80 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
81 | |
---|
82 | G4PenelopeRayleighModel::~G4PenelopeRayleighModel() |
---|
83 | { |
---|
84 | std::map <const G4Material*,G4DataVector*>::iterator i; |
---|
85 | for(i=SamplingTable.begin(); i != SamplingTable.end(); i++) { |
---|
86 | delete (*i).second; |
---|
87 | } |
---|
88 | if (samplingFunction_x) delete samplingFunction_x; |
---|
89 | if (samplingFunction_xNoLog) delete samplingFunction_xNoLog; |
---|
90 | } |
---|
91 | |
---|
92 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
93 | |
---|
94 | void G4PenelopeRayleighModel::Initialise(const G4ParticleDefinition* , |
---|
95 | const G4DataVector& ) |
---|
96 | { |
---|
97 | if (verboseLevel > 3) |
---|
98 | G4cout << "Calling G4PenelopeRayleighModel::Initialise()" << G4endl; |
---|
99 | |
---|
100 | |
---|
101 | if (verboseLevel > 0) { |
---|
102 | G4cout << "Penelope Rayleigh model is initialized " << G4endl |
---|
103 | << "Energy range: " |
---|
104 | << LowEnergyLimit() / keV << " keV - " |
---|
105 | << HighEnergyLimit() / GeV << " GeV" |
---|
106 | << G4endl; |
---|
107 | } |
---|
108 | |
---|
109 | if(isInitialised) return; |
---|
110 | fParticleChange = GetParticleChangeForGamma(); |
---|
111 | isInitialised = true; |
---|
112 | } |
---|
113 | |
---|
114 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
115 | |
---|
116 | G4double |
---|
117 | G4PenelopeRayleighModel::CrossSectionPerVolume(const G4Material* material, |
---|
118 | const G4ParticleDefinition* p, |
---|
119 | G4double ekin, |
---|
120 | G4double, |
---|
121 | G4double) |
---|
122 | { |
---|
123 | // Penelope model to calculate the Rayleigh scattering inverse mean |
---|
124 | // free path. |
---|
125 | // |
---|
126 | // The basic method is from |
---|
127 | // M. Born, Atomic physics, Ed. Blackie and Sons (1969) |
---|
128 | // using numerical approximations developed in |
---|
129 | // J. Baro' et al., Radiat. Phys. Chem. 44 (1994) 531 |
---|
130 | // Data used for the form factor used in the calculation (numerical integral of |
---|
131 | // dSigma/dOmega) are been derived by fitting the atomic forn factor tables |
---|
132 | // tabulated in |
---|
133 | // J.H. Hubbel et al., J. Phys. Chem. Ref. Data 4 (1975) 471; erratum ibid. |
---|
134 | // 6 (1977) 615. |
---|
135 | // The numerical integration of the differential cross section dSigma/dOmega, |
---|
136 | // which is implemented in the DifferentialCrossSection() method, is performed |
---|
137 | // by 20-point Gaussian method (managed by G4PenelopeIntegrator). |
---|
138 | // |
---|
139 | |
---|
140 | if (verboseLevel > 3) |
---|
141 | G4cout << "Calling CrossSectionPerVolume() of G4PenelopeRayleighModel" << G4endl; |
---|
142 | SetupForMaterial(p, material, ekin); |
---|
143 | |
---|
144 | //Assign local variable "material" to private member "theMaterial", because |
---|
145 | //this information is necessary to calculate the cross section |
---|
146 | theMaterial = material; |
---|
147 | |
---|
148 | G4int nElements = material->GetNumberOfElements(); |
---|
149 | const G4ElementVector* elementVector = material->GetElementVector(); |
---|
150 | |
---|
151 | G4int maxZ=0; |
---|
152 | for (G4int i=0; i<nElements; i++) |
---|
153 | { |
---|
154 | G4int Z = (G4int) (*elementVector)[i]->GetZ(); |
---|
155 | if (Z>maxZ) |
---|
156 | maxZ = Z; |
---|
157 | } |
---|
158 | |
---|
159 | G4double ec=std::min(ekin,0.5*maxZ); |
---|
160 | factorE = 849.3315*(ec/electron_mass_c2)*(ec/electron_mass_c2); |
---|
161 | G4double cs=0; |
---|
162 | // |
---|
163 | //Integrate the Differential Cross Section dSigma/dCosTheta between -1 and 1. |
---|
164 | // |
---|
165 | G4PenelopeIntegrator<G4PenelopeRayleighModel,G4double(G4PenelopeRayleighModel::*)(G4double)> |
---|
166 | theIntegrator; |
---|
167 | cs = |
---|
168 | theIntegrator.Calculate(this,&G4PenelopeRayleighModel::DifferentialCrossSection, |
---|
169 | -1.0,0.90,1e-06); |
---|
170 | cs += theIntegrator.Calculate(this,&G4PenelopeRayleighModel::DifferentialCrossSection, |
---|
171 | 0.90,0.9999999,1e-06); |
---|
172 | cs = cs*(ec/ekin)*(ec/ekin)*pi*classic_electr_radius*classic_electr_radius; |
---|
173 | // |
---|
174 | // Here cs represents the cross section per molecule for materials defined with chemical |
---|
175 | // formulas and the average cross section per atom in compounds (defined with the mass |
---|
176 | // fraction) |
---|
177 | // |
---|
178 | const G4int* stechiometric = material->GetAtomsVector(); |
---|
179 | //This is the total density of atoms in the material |
---|
180 | G4double atomDensity = material->GetTotNbOfAtomsPerVolume(); |
---|
181 | G4double moleculeDensity = 0; |
---|
182 | |
---|
183 | //Default case: the material is a compound. In this case, cs is the average cross section |
---|
184 | // _per_atom_ and one has to multiply for the atom density. |
---|
185 | G4double cross = atomDensity*cs; |
---|
186 | G4bool isAMolecule = false; |
---|
187 | |
---|
188 | //Alternative case: the material is a molecule. In this case cs is the cross section |
---|
189 | // _per_molecule_ and one has to multiply for the molecule density |
---|
190 | if (stechiometric) |
---|
191 | { |
---|
192 | //Calculate the total number of atoms per molecule |
---|
193 | G4int atomsPerMolecule = 0; |
---|
194 | for (G4int k=0;k<nElements;k++) |
---|
195 | { |
---|
196 | atomsPerMolecule += stechiometric[k]; |
---|
197 | if (verboseLevel > 2) |
---|
198 | { |
---|
199 | G4cout << "Element: " << (G4int) (*elementVector)[k]->GetZ() << " has " << |
---|
200 | stechiometric[k] << " atoms/molecule" << G4endl; |
---|
201 | } |
---|
202 | } |
---|
203 | if (atomsPerMolecule) |
---|
204 | { |
---|
205 | isAMolecule = true; |
---|
206 | if (verboseLevel > 3) |
---|
207 | { |
---|
208 | G4cout << "Material " << material->GetName() << " is a molecule composed by " << |
---|
209 | atomsPerMolecule << " atoms" << G4endl; |
---|
210 | } |
---|
211 | moleculeDensity = atomDensity/((G4double) atomsPerMolecule); |
---|
212 | cross = cs*moleculeDensity; |
---|
213 | } |
---|
214 | } |
---|
215 | |
---|
216 | if (verboseLevel > 2) |
---|
217 | { |
---|
218 | if (isAMolecule) |
---|
219 | { |
---|
220 | G4cout << "Rayleigh cross section at " << ekin/keV << " keV for material " << |
---|
221 | material->GetName() << " (molecule) = " << cs/barn << " barn/molecule." << G4endl; |
---|
222 | G4cout << "Mean free path: " << (1./cross)/mm << " mm" << G4endl; |
---|
223 | } |
---|
224 | else |
---|
225 | { |
---|
226 | G4cout << "Rayleigh cross section at " << ekin/keV << " keV for material " << |
---|
227 | material->GetName() << " (compound) = " << cs/barn << " barn/atom." << G4endl; |
---|
228 | G4cout << "Mean free path: " << (1./cross)/mm << " mm" << G4endl; |
---|
229 | } |
---|
230 | } |
---|
231 | return cross; |
---|
232 | } |
---|
233 | |
---|
234 | |
---|
235 | //This is a dummy method. Never inkoved by the tracking, it just issues |
---|
236 | //a warning if one tries to get Cross Sections per Atom via the |
---|
237 | //G4EmCalculator. |
---|
238 | G4double G4PenelopeRayleighModel::ComputeCrossSectionPerAtom(const G4ParticleDefinition*, |
---|
239 | G4double, |
---|
240 | G4double, |
---|
241 | G4double, |
---|
242 | G4double, |
---|
243 | G4double) |
---|
244 | { |
---|
245 | G4cout << "*** G4PenelopeRayleighModel -- WARNING ***" << G4endl; |
---|
246 | G4cout << "Penelope Rayleigh model does not calculate cross section _per atom_ " << G4endl; |
---|
247 | G4cout << "so the result is always zero. For physics values, please invoke " << G4endl; |
---|
248 | G4cout << "GetCrossSectionPerVolume() or GetMeanFreePath() via the G4EmCalculator" << G4endl; |
---|
249 | return 0; |
---|
250 | } |
---|
251 | |
---|
252 | |
---|
253 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
254 | |
---|
255 | void G4PenelopeRayleighModel::SampleSecondaries(std::vector<G4DynamicParticle*>* , |
---|
256 | const G4MaterialCutsCouple* couple, |
---|
257 | const G4DynamicParticle* aDynamicGamma, |
---|
258 | G4double, |
---|
259 | G4double) |
---|
260 | { |
---|
261 | // Penelope model to sample the Rayleigh scattering final state. |
---|
262 | // |
---|
263 | // The angular deflection of the scattered photon is sampled according to the |
---|
264 | // differential cross section dSigma/dOmega used for the numerical integration, |
---|
265 | // and implemented in the DifferentialCrossSection() method. See comments in |
---|
266 | // method CrossSectionPerVolume() for more details on the original references |
---|
267 | // of the model. |
---|
268 | // |
---|
269 | |
---|
270 | if (verboseLevel > 3) |
---|
271 | G4cout << "Calling SamplingSecondaries() of G4PenelopeRayleighModel" << G4endl; |
---|
272 | |
---|
273 | G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy(); |
---|
274 | |
---|
275 | if (photonEnergy0 <= fIntrinsicLowEnergyLimit) |
---|
276 | { |
---|
277 | fParticleChange->ProposeTrackStatus(fStopAndKill); |
---|
278 | fParticleChange->SetProposedKineticEnergy(0.); |
---|
279 | fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0); |
---|
280 | return ; |
---|
281 | } |
---|
282 | |
---|
283 | G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection(); |
---|
284 | |
---|
285 | // Sampling inizialitation (build internal table) |
---|
286 | theMaterial = couple->GetMaterial(); |
---|
287 | InitialiseSampling(); |
---|
288 | |
---|
289 | G4DataVector* samplingFunction_y = SamplingTable.find(theMaterial)->second; |
---|
290 | |
---|
291 | // Sample the angle of the scattered photon |
---|
292 | G4double x2max = 2.0*std::log(41.2148*photonEnergy0/electron_mass_c2); |
---|
293 | G4int jm=0; |
---|
294 | G4int asize = samplingFunction_x->size(); |
---|
295 | if (x2max < (*samplingFunction_x)[1]) |
---|
296 | jm=0; |
---|
297 | else if(x2max>(*samplingFunction_x)[asize-2]) |
---|
298 | jm=asize-2; |
---|
299 | else |
---|
300 | { |
---|
301 | //spacing in the logTable |
---|
302 | G4double logScalingFactor = (*samplingFunction_x)[1]-(*samplingFunction_x)[0]; |
---|
303 | jm=(G4int) ((x2max-(*samplingFunction_x)[0])/logScalingFactor); |
---|
304 | } |
---|
305 | |
---|
306 | G4double rumax = (*samplingFunction_y)[jm]+((*samplingFunction_y)[jm+1]-(*samplingFunction_y)[jm])* |
---|
307 | (x2max-(*samplingFunction_x)[jm])/((*samplingFunction_x)[jm+1]-(*samplingFunction_x)[jm]); |
---|
308 | G4double cosTheta=0; |
---|
309 | G4double rejectionValue = 0; |
---|
310 | do{ |
---|
311 | G4double ru = rumax + std::log(G4UniformRand()); |
---|
312 | G4int j=0; |
---|
313 | G4int ju=jm+1; |
---|
314 | do{ |
---|
315 | G4int jt=(j+ju)/2; //bipartition |
---|
316 | if (ru > (*samplingFunction_y)[jt]) |
---|
317 | j=jt; |
---|
318 | else |
---|
319 | ju=jt; |
---|
320 | }while ((ju-j)>1); |
---|
321 | G4double x2rat = 0; |
---|
322 | G4double denomin = (*samplingFunction_y)[j+1]-(*samplingFunction_y)[j]; |
---|
323 | if (denomin > 1e-12) |
---|
324 | { |
---|
325 | x2rat = (*samplingFunction_x)[j]+(((*samplingFunction_x)[j+1]-(*samplingFunction_x)[j])* |
---|
326 | (ru-(*samplingFunction_y)[j])/denomin)-x2max; |
---|
327 | } |
---|
328 | else |
---|
329 | { |
---|
330 | x2rat = (*samplingFunction_x)[j]-x2max; |
---|
331 | } |
---|
332 | cosTheta = 1.0-2.0*std::exp(x2rat); |
---|
333 | rejectionValue = 0.5*(1.0+cosTheta*cosTheta); |
---|
334 | }while (G4UniformRand() > rejectionValue); |
---|
335 | |
---|
336 | G4double sinTheta = std::sqrt(1-cosTheta*cosTheta); |
---|
337 | |
---|
338 | // Scattered photon angles. ( Z - axis along the parent photon) |
---|
339 | G4double phi = twopi * G4UniformRand() ; |
---|
340 | G4double dirX = sinTheta*std::cos(phi); |
---|
341 | G4double dirY = sinTheta*std::sin(phi); |
---|
342 | G4double dirZ = cosTheta; |
---|
343 | |
---|
344 | // Update G4VParticleChange for the scattered photon |
---|
345 | G4ThreeVector photonDirection1(dirX, dirY, dirZ); |
---|
346 | |
---|
347 | photonDirection1.rotateUz(photonDirection0); |
---|
348 | fParticleChange->ProposeMomentumDirection(photonDirection1) ; |
---|
349 | fParticleChange->SetProposedKineticEnergy(photonEnergy0) ; |
---|
350 | |
---|
351 | return; |
---|
352 | } |
---|
353 | |
---|
354 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
355 | |
---|
356 | G4double G4PenelopeRayleighModel::DifferentialCrossSection(G4double cosTheta) |
---|
357 | { |
---|
358 | //Differential cross section for Rayleigh scattering |
---|
359 | G4double x2 = factorE*(1-cosTheta); |
---|
360 | G4double gradx = (1+cosTheta*cosTheta)*MolecularFormFactor(x2); |
---|
361 | return gradx; |
---|
362 | } |
---|
363 | |
---|
364 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
365 | |
---|
366 | G4double G4PenelopeRayleighModel::MolecularFormFactor(G4double x2) |
---|
367 | { |
---|
368 | //Squared molecular form factor (additivity rule) |
---|
369 | const G4int ntot=95; |
---|
370 | G4double RA1[ntot] = {0.0e0, 3.9265e+0, 4.3100e+1, 5.2757e+1, 2.5021e+1, |
---|
371 | 1.2211e+1, 9.3229e+0, 3.2455e+0, 2.4197e+0, 1.5985e+0, |
---|
372 | 3.0926e+1, 1.5315e+1, 7.7061e+0, 3.9493e+0, 2.2042e+0, |
---|
373 | 1.9453e+1, 1.9354e+1, 8.0374e+0, 8.3779e+1, 5.7370e+1, |
---|
374 | 5.2310e+1, 4.7514e+1, 4.3785e+1, 4.2048e+1, 3.6972e+1, |
---|
375 | 3.3849e+1, 3.1609e+1, 2.8763e+1, 2.7217e+1, 2.4263e+1, |
---|
376 | 2.2403e+1, 1.8606e+1, 1.5143e+1, 1.4226e+1, 1.1792e+1, |
---|
377 | 9.7574e+0, 1.2796e+1, 1.2854e+1, 1.2368e+1, 1.0208e+1, |
---|
378 | 8.2823e+0, 7.4677e+0, 7.6028e+0, 6.1090e+0, 5.5346e+0, |
---|
379 | 4.2340e+0, 4.0444e+0, 4.2905e+0, 4.7950e+0, 5.1112e+0, |
---|
380 | 5.2407e+0, 5.2153e+0, 5.1639e+0, 4.8814e+0, 5.8054e+0, |
---|
381 | 6.6724e+0, 6.5104e+0, 6.3364e+0, 6.2889e+0, 6.3028e+0, |
---|
382 | 6.3853e+0, 6.3475e+0, 6.5779e+0, 6.8486e+0, 7.0993e+0, |
---|
383 | 7.6122e+0, 7.9681e+0, 8.3481e+0, 6.3875e+0, 8.0042e+0, |
---|
384 | 8.0820e+0, 7.6940e+0, 7.1927e+0, 6.6751e+0, 6.1623e+0, |
---|
385 | 5.8335e+0, 5.5599e+0, 4.6551e+0, 4.4327e+0, 4.7601e+0, |
---|
386 | 5.2872e+0, 5.6084e+0, 5.7680e+0, 5.8041e+0, 5.7566e+0, |
---|
387 | 5.6541e+0, 6.3932e+0, 6.9313e+0, 7.0027e+0, 6.8796e+0, |
---|
388 | 6.4739e+0, 6.2405e+0, 6.0081e+0, 5.5708e+0, 5.3680e+0}; |
---|
389 | |
---|
390 | |
---|
391 | G4double RA2[ntot] = {0.0e0, 1.3426e-1, 9.4875e+1,-1.0896e+2,-4.5494e+1, |
---|
392 | -1.9572e+1,-1.2382e+1,-3.6827e+0,-2.4542e+0,-1.4453e+0, |
---|
393 | 1.3401e+2, 7.9717e+1, 6.2164e+1, 4.0300e+1, 3.1682e+1, |
---|
394 | -1.3639e+1,-1.5950e+1,-5.1523e+0, 1.8351e+2, 1.2205e+2, |
---|
395 | 1.0007e+2, 8.5632e+1, 7.9145e+1, 6.3675e+1, 6.2954e+1, |
---|
396 | 5.6601e+1, 5.4171e+1, 4.8752e+1, 3.8062e+1, 3.9933e+1, |
---|
397 | 4.8343e+1, 4.2137e+1, 3.4617e+1, 2.9430e+1, 2.4010e+1, |
---|
398 | 1.9744e+1, 4.0009e+1, 5.1614e+1, 5.0456e+1, 3.9088e+1, |
---|
399 | 2.6824e+1, 2.2953e+1, 2.4773e+1, 1.6893e+1, 1.4548e+1, |
---|
400 | 9.7226e+0, 1.0192e+1, 1.1153e+1, 1.3188e+1, 1.4733e+1, |
---|
401 | 1.5644e+1, 1.5939e+1, 1.5923e+1, 1.5254e+1, 2.0748e+1, |
---|
402 | 2.6901e+1, 2.7032e+1, 2.4938e+1, 2.1528e+1, 2.0362e+1, |
---|
403 | 1.9474e+1, 1.8238e+1, 1.7898e+1, 1.9174e+1, 1.9023e+1, |
---|
404 | 1.8194e+1, 1.8504e+1, 1.8955e+1, 1.4276e+1, 1.7558e+1, |
---|
405 | 1.8651e+1, 1.7984e+1, 1.6793e+1, 1.5469e+1, 1.4143e+1, |
---|
406 | 1.3149e+1, 1.2255e+1, 9.2352e+0, 8.6067e+0, 9.7460e+0, |
---|
407 | 1.1749e+1, 1.3281e+1, 1.4326e+1, 1.4920e+1, 1.5157e+1, |
---|
408 | 1.5131e+1, 1.9489e+1, 2.3649e+1, 2.4686e+1, 2.4760e+1, |
---|
409 | 2.1519e+1, 2.0099e+1, 1.8746e+1, 1.5943e+1, 1.4880e+1}; |
---|
410 | |
---|
411 | G4double RA3[ntot] = {0.0e0, 2.2648e+0, 1.0579e+3, 8.6177e+2, 2.4422e+2, |
---|
412 | 7.8788e+1, 3.8293e+1, 1.2564e+1, 6.9091e+0, 3.7926e+0, |
---|
413 | 0.0000e+0, 0.0000e+0, 1.6759e-9, 1.3026e+1, 3.0569e+0, |
---|
414 | 1.5521e+2, 1.2815e+2, 4.7378e+1, 9.2802e+2, 4.7508e+2, |
---|
415 | 3.6612e+2, 2.7582e+2, 2.1008e+2, 1.5903e+2, 1.2322e+2, |
---|
416 | 9.2898e+1, 7.1345e+1, 5.1651e+1, 3.8474e+1, 2.7410e+1, |
---|
417 | 1.9126e+1, 1.0889e+1, 5.3479e+0, 8.2223e+0, 5.0837e+0, |
---|
418 | 2.8905e+0, 2.7457e+0, 6.7082e-1, 0.0000e+0, 0.0000e+0, |
---|
419 | 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, |
---|
420 | 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, |
---|
421 | 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, |
---|
422 | 0.0000e+0, 0.0000e+0, 0.0000e+0, 1.7264e-1, 2.7322e-1, |
---|
423 | 3.9444e-1, 4.5648e-1, 6.2286e-1, 7.2468e-1, 8.4296e-1, |
---|
424 | 1.1698e+0, 1.2994e+0, 1.4295e+0, 0.0000e+0, 8.1570e-1, |
---|
425 | 6.9349e-1, 4.9536e-1, 3.1211e-1, 1.5931e-1, 2.9512e-2, |
---|
426 | 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, |
---|
427 | 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, |
---|
428 | 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, |
---|
429 | 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0, 0.0000e+0}; |
---|
430 | |
---|
431 | G4double RA4[ntot] = {1.1055e1,6.3519e0,4.7367e+1, 3.9402e+1, 2.2896e+1, |
---|
432 | 1.3979e+1, 1.0766e+1, 6.5252e+0, 5.1631e+0, 4.0524e+0, |
---|
433 | 2.7145e+1, 1.8724e+1, 1.4782e+1, 1.1608e+1, 9.7750e+0, |
---|
434 | 1.6170e+1, 1.5249e+1, 9.1916e+0, 5.4499e+1, 4.1381e+1, |
---|
435 | 3.7395e+1, 3.3815e+1, 3.1135e+1, 2.8273e+1, 2.6140e+1, |
---|
436 | 2.3948e+1, 2.2406e+1, 2.0484e+1, 1.8453e+1, 1.7386e+1, |
---|
437 | 1.7301e+1, 1.5388e+1, 1.3411e+1, 1.2668e+1, 1.1133e+1, |
---|
438 | 9.8081e+0, 1.3031e+1, 1.4143e+1, 1.3815e+1, 1.2077e+1, |
---|
439 | 1.0033e+1, 9.2549e+0, 9.5338e+0, 7.9076e+0, 7.3263e+0, |
---|
440 | 5.9996e+0, 6.0087e+0, 6.2660e+0, 6.7914e+0, 7.1501e+0, |
---|
441 | 7.3367e+0, 7.3729e+0, 7.3508e+0, 7.1465e+0, 8.2731e+0, |
---|
442 | 9.3745e+0, 9.3508e+0, 8.9897e+0, 8.4566e+0, 8.2690e+0, |
---|
443 | 8.1398e+0, 7.9183e+0, 7.9123e+0, 8.1677e+0, 8.1871e+0, |
---|
444 | 8.1766e+0, 8.2881e+0, 8.4227e+0, 7.0273e+0, 8.0002e+0, |
---|
445 | 8.1440e+0, 7.9104e+0, 7.5685e+0, 7.1970e+0, 6.8184e+0, |
---|
446 | 6.5469e+0, 6.3056e+0, 5.4844e+0, 5.2832e+0, 5.5889e+0, |
---|
447 | 6.0919e+0, 6.4340e+0, 6.6426e+0, 6.7428e+0, 6.7636e+0, |
---|
448 | 6.7281e+0, 7.5729e+0, 8.2808e+0, 8.4400e+0, 8.4220e+0, |
---|
449 | 7.8662e+0, 7.5993e+0, 7.3353e+0, 6.7829e+0, 6.5520e+0}; |
---|
450 | |
---|
451 | G4double RA5[ntot] = {0.0e0, 4.9828e+0, 5.5674e+1, 3.0902e+1, 1.1496e+1, |
---|
452 | 4.8936e+0, 2.5506e+0, 1.2236e+0, 7.4698e-1, 4.7042e-1, |
---|
453 | 4.7809e+0, 4.6315e+0, 4.3677e+0, 4.9269e+0, 2.6033e+0, |
---|
454 | 9.6229e+0, 7.2592e+0, 4.1634e+0, 1.3999e+1, 8.6975e+0, |
---|
455 | 6.9630e+0, 5.4681e+0, 4.2653e+0, 3.2848e+0, 2.7354e+0, |
---|
456 | 2.1617e+0, 1.7030e+0, 1.2826e+0, 9.7080e-1, 7.2227e-1, |
---|
457 | 5.0874e-1, 3.1402e-1, 1.6360e-1, 3.2918e-1, 2.3570e-1, |
---|
458 | 1.5868e-1, 1.5146e-1, 9.7662e-2, 7.3151e-2, 6.4206e-2, |
---|
459 | 4.8945e-2, 4.3189e-2, 4.4368e-2, 3.3976e-2, 3.0466e-2, |
---|
460 | 2.4477e-2, 3.7202e-2, 3.7093e-2, 3.8161e-2, 3.8576e-2, |
---|
461 | 3.8403e-2, 3.7806e-2, 3.4958e-2, 3.6029e-2, 4.3087e-2, |
---|
462 | 4.7069e-2, 4.6452e-2, 4.2486e-2, 4.1517e-2, 4.1691e-2, |
---|
463 | 4.2813e-2, 4.2294e-2, 4.5287e-2, 4.8462e-2, 4.9726e-2, |
---|
464 | 5.5097e-2, 5.6568e-2, 5.8069e-2, 1.2270e-2, 3.8006e-2, |
---|
465 | 3.5048e-2, 3.0050e-2, 2.5069e-2, 2.0485e-2, 1.6151e-2, |
---|
466 | 1.4631e-2, 1.4034e-2, 1.1978e-2, 1.1522e-2, 1.2375e-2, |
---|
467 | 1.3805e-2, 1.4954e-2, 1.5832e-2, 1.6467e-2, 1.6896e-2, |
---|
468 | 1.7166e-2, 1.9954e-2, 2.2497e-2, 2.1942e-2, 2.1965e-2, |
---|
469 | 2.0005e-2, 1.8927e-2, 1.8167e-2, 1.6314e-2, 1.5522e-2}; |
---|
470 | |
---|
471 | G4double x=std::sqrt(x2); |
---|
472 | G4double gradx1=0.0; |
---|
473 | |
---|
474 | G4int nElements = theMaterial->GetNumberOfElements(); |
---|
475 | const G4ElementVector* elementVector = theMaterial->GetElementVector(); |
---|
476 | const G4int* stechiometric = theMaterial->GetAtomsVector(); |
---|
477 | const G4double* vector_of_atoms = theMaterial->GetVecNbOfAtomsPerVolume(); |
---|
478 | const G4double tot_atoms = theMaterial->GetTotNbOfAtomsPerVolume(); |
---|
479 | for (G4int i=0;i<nElements;i++) |
---|
480 | { |
---|
481 | G4int Z = (G4int) (*elementVector)[i]->GetZ(); |
---|
482 | if (Z>ntot) Z=95; |
---|
483 | G4double denomin = 1.+x2*(RA4[Z-1]+x2*RA5[Z-1]); |
---|
484 | G4double fa=Z*(1+x2*(RA1[Z-1]+x*(RA2[Z-1]+x*RA3[Z-1])))/(denomin*denomin); |
---|
485 | if (Z>10 && fa>2.0) |
---|
486 | { |
---|
487 | G4double k1=0.3125; |
---|
488 | G4double k2=2.426311e-02; |
---|
489 | G4double Pa=(Z-k1)*fine_structure_const; |
---|
490 | G4double Pg=std::sqrt(1-(Pa*Pa)); |
---|
491 | G4double Pq=k2*x/Pa; |
---|
492 | G4double fb=std::sin(2.0*Pg*std::atan(Pq))/(Pg*Pq*std::pow((1+Pq*Pq),Pg)); |
---|
493 | fa=std::max(fa,fb); |
---|
494 | } |
---|
495 | if (stechiometric && stechiometric[i]!=0) |
---|
496 | gradx1 += stechiometric[i]*(fa*fa); //sum on the molecule |
---|
497 | else |
---|
498 | gradx1 += (vector_of_atoms[i]/tot_atoms)*(fa*fa); //weighted mean |
---|
499 | } |
---|
500 | return gradx1; |
---|
501 | } |
---|
502 | |
---|
503 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... |
---|
504 | |
---|
505 | void G4PenelopeRayleighModel::InitialiseSampling() |
---|
506 | { |
---|
507 | if (!samplingFunction_x || !samplingFunction_xNoLog) |
---|
508 | { |
---|
509 | G4cout << "G4PenelopeRayleighModel::InitialiseSampling(), something wrong" << G4endl; |
---|
510 | G4cout << "It looks like G4PenelopeRayleighModel::PrepareConstants() has not been called" << G4endl; |
---|
511 | G4Exception(); |
---|
512 | } |
---|
513 | if (!SamplingTable.count(theMaterial)) //material not defined yet |
---|
514 | { |
---|
515 | G4double XlowInte = 0.; |
---|
516 | G4double XhighInte=(*samplingFunction_xNoLog)[0]; |
---|
517 | //material not inizialized yet: initialize it now |
---|
518 | G4DataVector* samplingFunction_y = new G4DataVector(); |
---|
519 | G4PenelopeIntegrator<G4PenelopeRayleighModel,G4double(G4PenelopeRayleighModel::*)(G4double)> |
---|
520 | theIntegrator; |
---|
521 | G4double sum = theIntegrator.Calculate(this,&G4PenelopeRayleighModel::MolecularFormFactor, |
---|
522 | XlowInte,XhighInte,1e-10); |
---|
523 | samplingFunction_y->push_back(sum); |
---|
524 | for (G4int i=1;i<nPoints;i++) |
---|
525 | { |
---|
526 | XlowInte=(*samplingFunction_xNoLog)[i-1]; |
---|
527 | XhighInte=(*samplingFunction_xNoLog)[i]; |
---|
528 | sum += theIntegrator.Calculate(this, |
---|
529 | &G4PenelopeRayleighModel::MolecularFormFactor, |
---|
530 | XlowInte,XhighInte,1e-10); |
---|
531 | samplingFunction_y->push_back(sum); |
---|
532 | } |
---|
533 | for (G4int i=0;i<nPoints;i++) |
---|
534 | (*samplingFunction_y)[i]=std::log((*samplingFunction_y)[i]); |
---|
535 | |
---|
536 | // |
---|
537 | /* |
---|
538 | G4String nnn = theMaterial->GetName()+".ntab"; |
---|
539 | std::ofstream file(nnn); |
---|
540 | for (size_t k=0;k<samplingFunction_x->size();k++) |
---|
541 | file << (*samplingFunction_x)[k] << " " << (*samplingFunction_y)[k] << G4endl; |
---|
542 | file.close(); |
---|
543 | */ |
---|
544 | // |
---|
545 | SamplingTable[theMaterial] = samplingFunction_y; |
---|
546 | } |
---|
547 | } |
---|
548 | |
---|
549 | void G4PenelopeRayleighModel::PrepareConstants() |
---|
550 | { |
---|
551 | if (verboseLevel > 3) |
---|
552 | G4cout << "Calling G4PenelopeRayleighModel::PrepareConstants()" << G4endl; |
---|
553 | nPoints=241; |
---|
554 | Xlow=1e-04; |
---|
555 | Xhigh=1e06; |
---|
556 | if (samplingFunction_x) |
---|
557 | { |
---|
558 | delete samplingFunction_x; |
---|
559 | samplingFunction_x = 0; |
---|
560 | } |
---|
561 | if (samplingFunction_xNoLog) |
---|
562 | { |
---|
563 | delete samplingFunction_xNoLog; |
---|
564 | samplingFunction_xNoLog = 0; |
---|
565 | } |
---|
566 | |
---|
567 | samplingFunction_x = new G4DataVector(); |
---|
568 | samplingFunction_xNoLog = new G4DataVector(); |
---|
569 | |
---|
570 | G4double scalingFactor = std::pow((Xhigh/Xlow),((1/((G4double)(nPoints-1))))); |
---|
571 | G4double logScalingFactor=(1/((G4double)(nPoints-1)))*std::log(Xhigh/Xlow); |
---|
572 | //Logarithmic table between log(Xlow) and log(Xhigh) |
---|
573 | samplingFunction_x->push_back(std::log(Xlow)); |
---|
574 | //Table between Xlow and Xhigh with log spacement: needed for InitialiseSampling() |
---|
575 | samplingFunction_xNoLog->push_back(Xlow); |
---|
576 | |
---|
577 | for (G4int i=1;i<nPoints;i++) |
---|
578 | { |
---|
579 | G4double nextx = (*samplingFunction_x)[i-1]+logScalingFactor; |
---|
580 | samplingFunction_x->push_back(nextx); |
---|
581 | nextx = (*samplingFunction_xNoLog)[i-1]*scalingFactor; |
---|
582 | samplingFunction_xNoLog->push_back(nextx); |
---|
583 | } |
---|
584 | } |
---|