Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4XrayRayleighModel Class Reference

#include <G4XrayRayleighModel.hh>

+ Inheritance diagram for G4XrayRayleighModel:

Public Member Functions

 G4XrayRayleighModel (const G4ParticleDefinition *p=0, const G4String &nam="XrayRayleigh")
 
virtual ~G4XrayRayleighModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Attributes

G4ParticleChangeForGammafParticleChange
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Detailed Description

Definition at line 42 of file G4XrayRayleighModel.hh.

Constructor & Destructor Documentation

◆ G4XrayRayleighModel()

G4XrayRayleighModel::G4XrayRayleighModel ( const G4ParticleDefinition p = 0,
const G4String nam = "XrayRayleigh" 
)

Definition at line 49 of file G4XrayRayleighModel.cc.

51 :G4VEmModel(nam),isInitialised(false)
52{
54 lowEnergyLimit = 250*eV;
55 highEnergyLimit = 10.*MeV;
56 fFormFactor = 0.0;
57
58 // SetLowEnergyLimit(lowEnergyLimit);
59 SetHighEnergyLimit(highEnergyLimit);
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 if(verboseLevel > 0)
70 {
71 G4cout << "Xray Rayleigh is constructed " << G4endl
72 << "Energy range: "
73 << lowEnergyLimit / eV << " eV - "
74 << highEnergyLimit / MeV << " MeV"
75 << G4endl;
76 }
77}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
G4ParticleChangeForGamma * fParticleChange

◆ ~G4XrayRayleighModel()

G4XrayRayleighModel::~G4XrayRayleighModel ( )
virtual

Definition at line 81 of file G4XrayRayleighModel.cc.

82{
83
84}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4XrayRayleighModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0,
G4double  cut = 0,
G4double  emax = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 107 of file G4XrayRayleighModel.cc.

112{
113 if (verboseLevel > 3)
114 {
115 G4cout << "Calling CrossSectionPerAtom() of G4XrayRayleighModel" << G4endl;
116 }
117 if (gammaEnergy < lowEnergyLimit || gammaEnergy > highEnergyLimit)
118 {
119 return 0.0;
120 }
121 G4double k = gammaEnergy/hbarc;
122 k *= Bohr_radius;
123 G4double p0 = 0.680654;
124 G4double p1 = -0.0224188;
125 G4double lnZ = std::log(Z);
126
127 G4double lna = p0 + p1*lnZ;
128
129 G4double alpha = std::exp(lna);
130
131 G4double fo = std::pow(k, alpha);
132
133 p0 = 3.68455;
134 p1 = -0.464806;
135 lna = p0 + p1*lnZ;
136
137 fo *= 0.01*std::exp(lna);
138
139 fFormFactor = fo;
140
141 G4double b = 1. + 2.*fo;
142 G4double b2 = b*b;
143 G4double b3 = b*b2;
144
145 G4double xsc = fCofR*Z*Z/b3;
146 xsc *= fo*fo + (1. + fo)*(1. + fo);
147
148
149 return xsc;
150
151}
double G4double
Definition: G4Types.hh:64

◆ Initialise()

void G4XrayRayleighModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 88 of file G4XrayRayleighModel.cc.

90{
91 if (verboseLevel > 3)
92 {
93 G4cout << "Calling G4XrayRayleighModel::Initialise()" << G4endl;
94 }
95
96 InitialiseElementSelectors(particle,cuts);
97
98
99 if(isInitialised) return;
101 isInitialised = true;
102
103}
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:123

◆ SampleSecondaries()

void G4XrayRayleighModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 155 of file G4XrayRayleighModel.cc.

160{
161 if ( verboseLevel > 3)
162 {
163 G4cout << "Calling SampleSecondaries() of G4XrayRayleighModel" << G4endl;
164 }
165 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
166
167 G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
168
169
170 // Sample the angle of the scattered photon
171 // according to 1 + cosTheta*cosTheta distribution
172
173 G4double cosDipole, cosTheta, sinTheta;
174 G4double c, delta, cofA, signc = 1., a, power = 1./3.;
175
176 c = 4. - 8.*G4UniformRand();
177 a = c;
178
179 if( c < 0. )
180 {
181 signc = -1.;
182 a = -c;
183 }
184 delta = std::sqrt(a*a+4.);
185 delta += a;
186 delta *= 0.5;
187 cofA = -signc*std::pow(delta, power);
188 cosDipole = cofA - 1./cofA;
189
190 // select atom
191 const G4Element* elm = SelectRandomAtom(couple, aDynamicGamma->GetParticleDefinition(), photonEnergy0);
192 G4double Z = elm->GetZ();
193
194 G4double k = photonEnergy0/hbarc;
195 k *= Bohr_radius;
196 G4double p0 = 0.680654;
197 G4double p1 = -0.0224188;
198 G4double lnZ = std::log(Z);
199
200 G4double lna = p0 + p1*lnZ;
201
202 G4double alpha = std::exp(lna);
203
204 G4double fo = std::pow(k, alpha);
205
206 p0 = 3.68455;
207 p1 = -0.464806;
208 lna = p0 + p1*lnZ;
209
210 fo *= 0.01*pi*std::exp(lna);
211
212
213 G4double beta = fo/(1 + fo);
214
215 cosTheta = (cosDipole + beta)/(1. + cosDipole*beta);
216
217
218 if( cosTheta > 1.) cosTheta = 1.;
219 if( cosTheta < -1.) cosTheta = -1.;
220
221 sinTheta = std::sqrt( (1. - cosTheta)*(1. + cosTheta) );
222
223 // Scattered photon angles. ( Z - axis along the parent photon)
224
225 G4double phi = twopi * G4UniformRand() ;
226 G4double dirX = sinTheta*std::cos(phi);
227 G4double dirY = sinTheta*std::sin(phi);
228 G4double dirZ = cosTheta;
229
230 // Update G4VParticleChange for the scattered photon
231
232 G4ThreeVector photonDirection1(dirX, dirY, dirZ);
233 photonDirection1.rotateUz(photonDirection0);
234 fParticleChange->ProposeMomentumDirection(photonDirection1);
235
237}
#define G4UniformRand()
Definition: Randomize.hh:53
const G4ThreeVector & GetMomentumDirection() const
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
G4double GetZ() const
Definition: G4Element.hh:131
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:459
const G4double pi

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForGamma* G4XrayRayleighModel::fParticleChange
protected

Definition at line 70 of file G4XrayRayleighModel.hh.

Referenced by G4XrayRayleighModel(), Initialise(), and SampleSecondaries().


The documentation for this class was generated from the following files: