00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035 #include "G4XrayRayleighModel.hh"
00036 #include "G4PhysicalConstants.hh"
00037 #include "G4SystemOfUnits.hh"
00038
00040
00041 using namespace std;
00042
00043 const G4double G4XrayRayleighModel::fCofA = 2.*pi2*Bohr_radius*Bohr_radius;
00044
00045 const G4double G4XrayRayleighModel::fCofR = 8.*pi*classic_electr_radius*classic_electr_radius/3.;
00046
00048
00049 G4XrayRayleighModel::G4XrayRayleighModel(const G4ParticleDefinition*,
00050 const G4String& nam)
00051 :G4VEmModel(nam),isInitialised(false)
00052 {
00053 fParticleChange = 0;
00054 lowEnergyLimit = 250*eV;
00055 highEnergyLimit = 10.*MeV;
00056 fFormFactor = 0.0;
00057
00058
00059 SetHighEnergyLimit(highEnergyLimit);
00060
00061 verboseLevel= 0;
00062
00063
00064
00065
00066
00067
00068
00069 if(verboseLevel > 0)
00070 {
00071 G4cout << "Xray Rayleigh is constructed " << G4endl
00072 << "Energy range: "
00073 << lowEnergyLimit / eV << " eV - "
00074 << highEnergyLimit / MeV << " MeV"
00075 << G4endl;
00076 }
00077 }
00078
00080
00081 G4XrayRayleighModel::~G4XrayRayleighModel()
00082 {
00083
00084 }
00085
00087
00088 void G4XrayRayleighModel::Initialise(const G4ParticleDefinition* particle,
00089 const G4DataVector& cuts)
00090 {
00091 if (verboseLevel > 3)
00092 {
00093 G4cout << "Calling G4XrayRayleighModel::Initialise()" << G4endl;
00094 }
00095
00096 InitialiseElementSelectors(particle,cuts);
00097
00098
00099 if(isInitialised) return;
00100 fParticleChange = GetParticleChangeForGamma();
00101 isInitialised = true;
00102
00103 }
00104
00106
00107 G4double G4XrayRayleighModel::ComputeCrossSectionPerAtom(
00108 const G4ParticleDefinition*,
00109 G4double gammaEnergy,
00110 G4double Z, G4double,
00111 G4double, G4double)
00112 {
00113 if (verboseLevel > 3)
00114 {
00115 G4cout << "Calling CrossSectionPerAtom() of G4XrayRayleighModel" << G4endl;
00116 }
00117 if (gammaEnergy < lowEnergyLimit || gammaEnergy > highEnergyLimit)
00118 {
00119 return 0.0;
00120 }
00121 G4double k = gammaEnergy/hbarc;
00122 k *= Bohr_radius;
00123 G4double p0 = 0.680654;
00124 G4double p1 = -0.0224188;
00125 G4double lnZ = std::log(Z);
00126
00127 G4double lna = p0 + p1*lnZ;
00128
00129 G4double alpha = std::exp(lna);
00130
00131 G4double fo = std::pow(k, alpha);
00132
00133 p0 = 3.68455;
00134 p1 = -0.464806;
00135 lna = p0 + p1*lnZ;
00136
00137 fo *= 0.01*std::exp(lna);
00138
00139 fFormFactor = fo;
00140
00141 G4double b = 1. + 2.*fo;
00142 G4double b2 = b*b;
00143 G4double b3 = b*b2;
00144
00145 G4double xsc = fCofR*Z*Z/b3;
00146 xsc *= fo*fo + (1. + fo)*(1. + fo);
00147
00148
00149 return xsc;
00150
00151 }
00152
00154
00155 void G4XrayRayleighModel::SampleSecondaries(std::vector<G4DynamicParticle*>* ,
00156 const G4MaterialCutsCouple* couple,
00157 const G4DynamicParticle* aDynamicGamma,
00158 G4double,
00159 G4double)
00160 {
00161 if ( verboseLevel > 3)
00162 {
00163 G4cout << "Calling SampleSecondaries() of G4XrayRayleighModel" << G4endl;
00164 }
00165 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
00166
00167 G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
00168
00169
00170
00171
00172
00173 G4double cosDipole, cosTheta, sinTheta;
00174 G4double c, delta, cofA, signc = 1., a, power = 1./3.;
00175
00176 c = 4. - 8.*G4UniformRand();
00177 a = c;
00178
00179 if( c < 0. )
00180 {
00181 signc = -1.;
00182 a = -c;
00183 }
00184 delta = std::sqrt(a*a+4.);
00185 delta += a;
00186 delta *= 0.5;
00187 cofA = -signc*std::pow(delta, power);
00188 cosDipole = cofA - 1./cofA;
00189
00190
00191 const G4Element* elm = SelectRandomAtom(couple, aDynamicGamma->GetParticleDefinition(), photonEnergy0);
00192 G4double Z = elm->GetZ();
00193
00194 G4double k = photonEnergy0/hbarc;
00195 k *= Bohr_radius;
00196 G4double p0 = 0.680654;
00197 G4double p1 = -0.0224188;
00198 G4double lnZ = std::log(Z);
00199
00200 G4double lna = p0 + p1*lnZ;
00201
00202 G4double alpha = std::exp(lna);
00203
00204 G4double fo = std::pow(k, alpha);
00205
00206 p0 = 3.68455;
00207 p1 = -0.464806;
00208 lna = p0 + p1*lnZ;
00209
00210 fo *= 0.01*pi*std::exp(lna);
00211
00212
00213 G4double beta = fo/(1 + fo);
00214
00215 cosTheta = (cosDipole + beta)/(1. + cosDipole*beta);
00216
00217
00218 if( cosTheta > 1.) cosTheta = 1.;
00219 if( cosTheta < -1.) cosTheta = -1.;
00220
00221 sinTheta = std::sqrt( (1. - cosTheta)*(1. + cosTheta) );
00222
00223
00224
00225 G4double phi = twopi * G4UniformRand() ;
00226 G4double dirX = sinTheta*std::cos(phi);
00227 G4double dirY = sinTheta*std::sin(phi);
00228 G4double dirZ = cosTheta;
00229
00230
00231
00232 G4ThreeVector photonDirection1(dirX, dirY, dirZ);
00233 photonDirection1.rotateUz(photonDirection0);
00234 fParticleChange->ProposeMomentumDirection(photonDirection1);
00235
00236 fParticleChange->SetProposedKineticEnergy(photonEnergy0);
00237 }
00238
00239