source: trunk/source/processes/electromagnetic/lowenergy/src/G4PenelopeRayleighModel.cc@ 1057

Last change on this file since 1057 was 1055, checked in by garnier, 17 years ago

maj sur la beta de geant 4.9.3

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