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

Last change on this file since 1014 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

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