source: trunk/source/processes/electromagnetic/lowenergy/src/G4PenelopeRayleigh.cc @ 991

Last change on this file since 991 was 991, checked in by garnier, 15 years ago

update

File size: 19.5 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// --------------------------------------------------------------------
27//
28// $Id: G4PenelopeRayleigh.cc,v 1.15 2007/09/03 09:43:14 pandola Exp $
29// GEANT4 tag $Name: geant4-09-02 $
30//
31// Author: L. Pandola (luciano.pandola@cern.ch)
32//
33// History:
34// --------
35// 14 Feb 2003   MG Pia       Corrected compilation errors and warnings
36//                            from SUN
37// 10 Mar 2003 V.Ivanchenko   Remove CutPerMaterial warning
38// 12 Mar 2003 L.Pandola      Code "cleaned" - Cuts per region
39// 17 Mar 2004 L.Pandola      Removed unnecessary calls to std::pow(a,b)
40// 18 Mar 2004 M.Mendenhall   Introduced SamplingTable (performance improvement)
41// 03 Sep 2007 L.Pandola      Bug fix for the filling of physics table for
42//                            compounds defined by the mass fraction (bug #965)
43// --------------------------------------------------------------------
44
45#include "G4PenelopeRayleigh.hh"
46#include "Randomize.hh"
47#include "G4ParticleDefinition.hh"
48#include "G4Track.hh"
49#include "G4Step.hh"
50#include "G4ForceCondition.hh"
51#include "G4Gamma.hh"
52#include "G4Electron.hh"
53#include "G4DynamicParticle.hh"
54#include "G4VParticleChange.hh"
55#include "G4ThreeVector.hh"
56#include "G4VCrossSectionHandler.hh"
57#include "G4CrossSectionHandler.hh"
58#include "G4VEMDataSet.hh"
59#include "G4EMDataSet.hh"
60#include "G4CompositeEMDataSet.hh"
61#include "G4VDataSetAlgorithm.hh"
62#include "G4LogLogInterpolation.hh"
63#include "G4PenelopeIntegrator.hh"
64#include "G4MaterialCutsCouple.hh"
65
66G4PenelopeRayleigh::G4PenelopeRayleigh(const G4String& processName)
67  : G4VDiscreteProcess(processName),
68    lowEnergyLimit(250*eV),             
69    highEnergyLimit(100*GeV),
70    samplingConstant(0.0),
71    nBins(200),
72    intrinsicLowEnergyLimit(10*eV),
73    intrinsicHighEnergyLimit(100*GeV)
74
75{
76  if (lowEnergyLimit < intrinsicLowEnergyLimit ||
77      highEnergyLimit > intrinsicHighEnergyLimit)
78    {
79      G4Exception("G4PenelopeRayleigh::G4PenelopeRayleigh - energy limit outside intrinsic process validity range");
80    }
81
82  samplingFunction_x = 0;
83  samplingFunction_y = 0;
84  meanFreePathTable = 0;
85  material = 0;
86
87   if (verboseLevel > 0)
88     {
89       G4cout << GetProcessName() << " is created " << G4endl
90              << "Energy range: "
91              << lowEnergyLimit / keV << " keV - "
92              << highEnergyLimit / GeV << " GeV"
93              << G4endl;
94     }
95}
96
97G4PenelopeRayleigh::~G4PenelopeRayleigh()
98{
99        delete meanFreePathTable;
100       
101        SamplingTablePair::iterator i;
102        for(i=SamplingTables.begin(); i != SamplingTables.end(); i++) {
103                delete (*i).second.first;
104                delete (*i).second.second;
105        }
106}
107
108
109void G4PenelopeRayleigh::BuildPhysicsTable(const G4ParticleDefinition& )
110{
111
112  G4DataVector energyVector;
113  G4double dBin = std::log10(highEnergyLimit/lowEnergyLimit)/nBins;
114  for (G4int i=0;i<nBins;i++)
115    {
116      energyVector.push_back(std::pow(10.,std::log10(lowEnergyLimit)+i*dBin));
117    }
118
119  const G4MaterialTable* materialTable = G4Material::GetMaterialTable();
120  G4int nMaterials = G4Material::GetNumberOfMaterials();
121
122  size_t nOfBins = energyVector.size();
123  size_t bin=0;
124
125  G4VDataSetAlgorithm* algo = new G4LogLogInterpolation();
126  G4VEMDataSet* materialSet = new G4CompositeEMDataSet(algo,1.,1.);
127  std::vector<G4VEMDataSet*> matCrossSections;
128
129  G4int m;
130  for (m=0; m<nMaterials; m++)
131    {
132      G4DataVector* energies = new G4DataVector;
133      G4DataVector* data = new G4DataVector;
134      material = (*materialTable)[m];
135
136      G4int nElements = material->GetNumberOfElements();
137      const G4ElementVector* elementVector = material->GetElementVector();
138      const G4double k1=849.3315;
139
140
141      G4int IZZ=0;
142      G4int iright=0;
143      for (G4int i=0; i<nElements; i++) {
144        G4int Z = (G4int) (*elementVector)[i]->GetZ();
145        if (Z>IZZ){
146          IZZ = Z;
147          iright=i;
148        }
149      }
150      for (bin=0; bin<nOfBins; bin++)
151        {
152          energies->push_back(energyVector[bin]);
153          G4double ec=std::min(energyVector[bin],0.5*IZZ);
154          G4double energyRatio = ec/electron_mass_c2;
155          facte = k1*energyRatio*energyRatio;
156          G4double cs=0;
157          G4PenelopeIntegrator<G4PenelopeRayleigh,G4double(G4PenelopeRayleigh::*)(G4double)> theIntegrator;
158          cs =
159            theIntegrator.Calculate(this,&G4PenelopeRayleigh::DifferentialCrossSection,-1.0,0.90,1e-06);
160          cs += theIntegrator.Calculate(this,&G4PenelopeRayleigh::DifferentialCrossSection,0.90,0.9999999,1e-06);
161          cs = cs*(ec/energyVector[bin])*(ec/energyVector[bin])*pi*classic_electr_radius*classic_electr_radius;
162          const G4double* vector_of_atoms = material->GetVecNbOfAtomsPerVolume();
163          const G4int* stechiometric = material->GetAtomsVector();
164          G4double density = vector_of_atoms[iright]; //non-bound molecules (default)
165          if (stechiometric)
166            {
167              if (stechiometric[iright])
168                density = vector_of_atoms[iright]/stechiometric[iright]; //number of molecules per volume
169            }
170          G4double cross = density*cs;
171          data->push_back(cross);
172        }
173      G4VEMDataSet* elSet = new G4EMDataSet(0,energies,data,algo); 
174      G4VEMDataSet* setForMat = new G4CompositeEMDataSet(algo);
175      setForMat->AddComponent(elSet); 
176      matCrossSections.push_back(setForMat);
177    }
178  G4double matCS = 0.0;
179
180  for (m=0; m<nMaterials; m++)
181    {
182      G4DataVector* energies = new G4DataVector;
183      G4DataVector* data = new G4DataVector;
184      for (bin=0;bin<nOfBins;bin++){
185        energies->push_back(energyVector[bin]);
186        matCS = (matCrossSections[m]->GetComponent(0))->FindValue(energyVector[bin]);
187        if (matCS > 0.){
188          data->push_back(1./matCS);
189        }
190        else
191          {
192            data->push_back(DBL_MAX);
193          }
194      }
195      G4VEMDataSet* dataSet = new G4EMDataSet(m,energies,data,algo,1.,1.); 
196      materialSet->AddComponent(dataSet);
197    }
198  meanFreePathTable = materialSet;
199}
200
201G4VParticleChange* G4PenelopeRayleigh::PostStepDoIt(const G4Track& aTrack,
202                                                     const G4Step& aStep)
203{
204  aParticleChange.Initialize(aTrack);
205
206  const G4DynamicParticle* incidentPhoton = aTrack.GetDynamicParticle();
207  G4double photonEnergy0 = incidentPhoton->GetKineticEnergy();
208
209  if (photonEnergy0 <= lowEnergyLimit)
210    {
211      aParticleChange.ProposeTrackStatus(fStopAndKill);
212      aParticleChange.ProposeEnergy(0.);
213      aParticleChange.ProposeLocalEnergyDeposit(photonEnergy0);
214      return G4VDiscreteProcess::PostStepDoIt(aTrack,aStep);
215    }
216
217
218  G4ParticleMomentum photonDirection0 = incidentPhoton->GetMomentumDirection();
219  const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
220  material = couple->GetMaterial();
221 
222  // Sampling inizialitation (build internal table)
223 
224  InizialiseSampling();
225
226  // Sample the angle of the scattered photon
227  const G4double xpar=41.2148;
228  G4double x2max = 2.0*std::log(xpar*photonEnergy0/electron_mass_c2);
229  G4int jm;
230  G4int asize = samplingFunction_x->size();
231  if (x2max<(*samplingFunction_x)[1]) 
232    {
233      jm=0;
234    }
235  else if(x2max>(*samplingFunction_x)[asize-2])
236    {
237    jm=asize-2;
238    }
239  else 
240    {
241      jm=(G4int) ((x2max-(*samplingFunction_x)[0])/samplingConstant);
242    }
243  G4double rumax = (*samplingFunction_y)[jm]+((*samplingFunction_y)[jm+1]-(*samplingFunction_y)[jm])*
244    (x2max-(*samplingFunction_x)[jm])/((*samplingFunction_x)[jm+1]-(*samplingFunction_x)[jm]); 
245  G4int j,ju,jt;
246  G4double ru,denomin,x2rat;
247  G4double CDT,G,rand;
248  do{
249    ru = rumax + std::log(G4UniformRand());
250    j=0;
251    ju=jm+1;
252    do{
253     jt=(j+ju)/2; //bipartition
254     if (ru > (*samplingFunction_y)[jt])
255      {
256        j=jt;
257      }
258     else
259      {
260        ju=jt;
261      }
262    }while ((ju-j)>1);
263    denomin = (*samplingFunction_y)[j+1]-(*samplingFunction_y)[j];
264    if (denomin > 1e-12) 
265     {
266       x2rat = (*samplingFunction_x)[j]+(((*samplingFunction_x)[j+1]-(*samplingFunction_x)[j])*
267                                         (ru-(*samplingFunction_y)[j])/denomin)-x2max;
268     }
269    else
270     {
271       x2rat = (*samplingFunction_x)[j]-x2max;
272     }
273    CDT = 1.0-2.0*std::exp(x2rat);
274    G = 0.5*(1.0+CDT*CDT);
275    rand = G4UniformRand();
276   }while (rand>G);
277 
278  G4double cosTheta = CDT;
279  G4double sinTheta = std::sqrt(1-cosTheta*cosTheta);
280 
281
282
283
284  // Scattered photon angles. ( Z - axis along the parent photon)
285  G4double phi = twopi * G4UniformRand() ;
286  G4double dirX = sinTheta*std::cos(phi);
287  G4double dirY = sinTheta*std::sin(phi);
288  G4double dirZ = cosTheta;
289 
290  // Update G4VParticleChange for the scattered photon
291  G4ThreeVector photonDirection1(dirX, dirY, dirZ);
292
293  photonDirection1.rotateUz(photonDirection0);
294  aParticleChange.ProposeEnergy(photonEnergy0);
295  aParticleChange.ProposeMomentumDirection(photonDirection1);
296 
297  aParticleChange.SetNumberOfSecondaries(0);
298
299  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
300}
301
302G4bool G4PenelopeRayleigh::IsApplicable(const G4ParticleDefinition& particle)
303{
304  return ( &particle == G4Gamma::Gamma() ); 
305}
306
307G4double G4PenelopeRayleigh::GetMeanFreePath(const G4Track& track,
308                                             G4double, // previousStepSize
309                                              G4ForceCondition*)
310{
311  const G4DynamicParticle* photon = track.GetDynamicParticle();
312  G4double energy = photon->GetKineticEnergy(); 
313  material = track.GetMaterial();
314  size_t materialIndex = material->GetIndex();
315
316  G4double meanFreePath;
317  if (energy > highEnergyLimit) meanFreePath = meanFreePathTable->FindValue(highEnergyLimit,materialIndex);
318  else if (energy < lowEnergyLimit) meanFreePath = DBL_MAX;
319  else
320    {
321      meanFreePath = meanFreePathTable->FindValue(energy,materialIndex);
322    }
323  return meanFreePath;
324}
325
326void G4PenelopeRayleigh::InizialiseSampling()
327{
328  SamplingTablePair::iterator theTable=SamplingTables.find(material);
329  const G4int points=241;
330  G4double Xlow = 0.;
331  G4double Xhigh=1e-04;
332  G4double fact = std::pow((1e06/Xhigh),(1/240.0));
333  samplingConstant=std::log(fact);
334 
335  if (theTable==SamplingTables.end()) { //material not inizialized yet
336    samplingFunction_x = new G4DataVector();
337    samplingFunction_y = new G4DataVector();
338    G4double sum = 0.0;
339    G4PenelopeIntegrator<G4PenelopeRayleigh,G4double(G4PenelopeRayleigh::*)(G4double)> theIntegrator;
340    sum = theIntegrator.Calculate(this,&G4PenelopeRayleigh::MolecularFormFactor,
341                                  Xlow,Xhigh,1e-10); 
342    samplingFunction_x->push_back(Xhigh);
343    samplingFunction_y->push_back(sum);
344    G4int i;
345    for (i=1;i<points;i++){
346      Xlow=Xhigh;
347      Xhigh=Xhigh*fact;
348      sum = theIntegrator.Calculate(this,
349                                    &G4PenelopeRayleigh::MolecularFormFactor,
350                                    Xlow,Xhigh,1e-10);
351      samplingFunction_x->push_back(Xhigh);
352      samplingFunction_y->push_back(sum+(*samplingFunction_y)[i-1]);
353    }
354    for (i=0;i<points;i++){
355      (*samplingFunction_x)[i]=std::log((*samplingFunction_x)[i]);
356      (*samplingFunction_y)[i]=std::log((*samplingFunction_y)[i]);
357    }
358    SamplingTables[material] = std::pair<G4DataVector*,G4DataVector*> (samplingFunction_x,samplingFunction_y);
359  }
360  else { //material already inizialized
361    samplingFunction_x=(*theTable).second.first;
362    samplingFunction_y=(*theTable).second.second;
363  }
364}
365
366
367G4double G4PenelopeRayleigh::MolecularFormFactor(G4double y)
368{
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(y);
472  G4double gradx1=0.0;
473  G4double fa=0.0;
474
475  G4int nElements = material->GetNumberOfElements();
476  const G4ElementVector* elementVector = material->GetElementVector();
477  const G4int* stechiometric = material->GetAtomsVector();
478  const G4double* vector_of_atoms = material->GetVecNbOfAtomsPerVolume();
479  const G4double tot_atoms = material->GetTotNbOfAtomsPerVolume();
480  for (G4int i=0;i<nElements;i++){
481    G4int Z = (G4int) (*elementVector)[i]->GetZ();
482    if (Z>ntot) Z=95;
483    G4double denomin = 1+y*(RA4[Z-1]+y*RA5[Z-1]);
484    fa=Z*(1+y*(RA1[Z-1]+x*(RA2[Z-1]+x*RA3[Z-1])))/(denomin*denomin);
485    G4bool a = ((Z>10) && (fa>2.0));
486    if (a) 
487      {
488        G4double Pa,Pg,Pq,fb;
489        G4double k1=0.3125;
490        G4double k2=2.426311e-02;
491        Pa=(Z-k1)*fine_structure_const;
492        Pg=std::sqrt(1-(Pa*Pa));
493        Pq=k2*x/Pa;
494        fb=std::sin(2*Pg*std::atan(Pq))/(Pg*Pq*std::pow((1+Pq*Pq),Pg));
495        fa=std::max(fa,fb);
496      }
497    if (stechiometric && stechiometric[i]!=0)
498      {
499        gradx1=gradx1+stechiometric[i]*(fa*fa); //sum on the molecule
500      }
501    else
502      {
503        gradx1=gradx1+(vector_of_atoms[i]/tot_atoms)*(fa*fa); //weighted mean
504      }
505  }
506  return gradx1;
507}
508
509
510G4double G4PenelopeRayleigh::DifferentialCrossSection(G4double x)
511{
512  G4double x2=facte*(1-x);
513  G4double gradx = (1+x*x)*MolecularFormFactor(x2);
514  return gradx;
515}
516
517
Note: See TracBrowser for help on using the repository browser.