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

#include <G4KleinNishinaCompton.hh>

+ Inheritance diagram for G4KleinNishinaCompton:

Public Member Functions

 G4KleinNishinaCompton (const G4ParticleDefinition *p=0, const G4String &nam="Klein-Nishina")
 
virtual ~G4KleinNishinaCompton ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
 
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

G4ParticleDefinitiontheGamma
 
G4ParticleDefinitiontheElectron
 
G4ParticleChangeForGammafParticleChange
 
G4double lowestGammaEnergy
 
- 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 58 of file G4KleinNishinaCompton.hh.

Constructor & Destructor Documentation

◆ G4KleinNishinaCompton()

G4KleinNishinaCompton::G4KleinNishinaCompton ( const G4ParticleDefinition p = 0,
const G4String nam = "Klein-Nishina" 
)

Definition at line 63 of file G4KleinNishinaCompton.cc.

65 : G4VEmModel(nam)
66{
69 lowestGammaEnergy = 1.0*eV;
71}
static G4Electron * Electron()
Definition: G4Electron.cc:94
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4ParticleDefinition * theGamma
G4ParticleChangeForGamma * fParticleChange
G4ParticleDefinition * theElectron

◆ ~G4KleinNishinaCompton()

G4KleinNishinaCompton::~G4KleinNishinaCompton ( )
virtual

Definition at line 75 of file G4KleinNishinaCompton.cc.

76{}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4KleinNishinaCompton::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A,
G4double  cut,
G4double  emax 
)
virtual

Reimplemented from G4VEmModel.

Reimplemented in G4PolarizedComptonModel.

Definition at line 88 of file G4KleinNishinaCompton.cc.

93{
94 G4double xSection = 0.0 ;
95 if ( Z < 0.9999 ) return xSection;
96 if ( GammaEnergy < 0.1*keV ) return xSection;
97 // if ( GammaEnergy > (100.*GeV/Z) ) return xSection;
98
99 static const G4double a = 20.0 , b = 230.0 , c = 440.0;
100
101 static const G4double
102 d1= 2.7965e-1*barn, d2=-1.8300e-1*barn, d3= 6.7527 *barn, d4=-1.9798e+1*barn,
103 e1= 1.9756e-5*barn, e2=-1.0205e-2*barn, e3=-7.3913e-2*barn, e4= 2.7079e-2*barn,
104 f1=-3.9178e-7*barn, f2= 6.8241e-5*barn, f3= 6.0480e-5*barn, f4= 3.0274e-4*barn;
105
106 G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
107 p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
108
109 G4double T0 = 15.0*keV;
110 if (Z < 1.5) T0 = 40.0*keV;
111
112 G4double X = max(GammaEnergy, T0) / electron_mass_c2;
113 xSection = p1Z*std::log(1.+2.*X)/X
114 + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
115
116 // modification for low energy. (special case for Hydrogen)
117 if (GammaEnergy < T0) {
118 G4double dT0 = 1.*keV;
119 X = (T0+dT0) / electron_mass_c2 ;
120 G4double sigma = p1Z*log(1.+2*X)/X
121 + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
122 G4double c1 = -T0*(sigma-xSection)/(xSection*dT0);
123 G4double c2 = 0.150;
124 if (Z > 1.5) c2 = 0.375-0.0556*log(Z);
125 G4double y = log(GammaEnergy/T0);
126 xSection *= exp(-y*(c1+c2*y));
127 }
128 // G4cout << "e= " << GammaEnergy << " Z= " << Z << " cross= " << xSection << G4endl;
129 return xSection;
130}
double G4double
Definition: G4Types.hh:64

Referenced by G4PolarizedComptonModel::ComputeCrossSectionPerAtom().

◆ Initialise()

void G4KleinNishinaCompton::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 80 of file G4KleinNishinaCompton.cc.

82{
84}
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109

◆ SampleSecondaries()

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

Implements G4VEmModel.

Reimplemented in G4PolarizedComptonModel.

Definition at line 134 of file G4KleinNishinaCompton.cc.

139{
140 // The scattered gamma energy is sampled according to Klein - Nishina formula.
141 // The random number techniques of Butcher & Messel are used
142 // (Nuc Phys 20(1960),15).
143 // Note : Effects due to binding of atomic electrons are negliged.
144
145 G4double gamEnergy0 = aDynamicGamma->GetKineticEnergy();
146
147 // extra protection
148 if(gamEnergy0 < lowestGammaEnergy) {
152 return;
153 }
154
155 G4double E0_m = gamEnergy0 / electron_mass_c2 ;
156
157 G4ThreeVector gamDirection0 = aDynamicGamma->GetMomentumDirection();
158
159 //
160 // sample the energy rate of the scattered gamma
161 //
162
163 G4double epsilon, epsilonsq, onecost, sint2, greject ;
164
165 G4double eps0 = 1./(1. + 2.*E0_m);
166 G4double epsilon0sq = eps0*eps0;
167 G4double alpha1 = - log(eps0);
168 G4double alpha2 = 0.5*(1.- epsilon0sq);
169
170 do {
171 if ( alpha1/(alpha1+alpha2) > G4UniformRand() ) {
172 epsilon = exp(-alpha1*G4UniformRand()); // eps0**r
173 epsilonsq = epsilon*epsilon;
174
175 } else {
176 epsilonsq = epsilon0sq + (1.- epsilon0sq)*G4UniformRand();
177 epsilon = sqrt(epsilonsq);
178 };
179
180 onecost = (1.- epsilon)/(epsilon*E0_m);
181 sint2 = onecost*(2.-onecost);
182 greject = 1. - epsilon*sint2/(1.+ epsilonsq);
183
184 } while (greject < G4UniformRand());
185
186 //
187 // scattered gamma angles. ( Z - axis along the parent gamma)
188 //
189
190 if(sint2 < 0.0) { sint2 = 0.0; }
191 G4double cosTeta = 1. - onecost;
192 G4double sinTeta = sqrt (sint2);
193 G4double Phi = twopi * G4UniformRand();
194
195 //
196 // update G4VParticleChange for the scattered gamma
197 //
198
199 G4ThreeVector gamDirection1(sinTeta*cos(Phi), sinTeta*sin(Phi), cosTeta);
200 gamDirection1.rotateUz(gamDirection0);
201 G4double gamEnergy1 = epsilon*gamEnergy0;
202 if(gamEnergy1 > lowestGammaEnergy) {
205 } else {
209 }
210
211 //
212 // kinematic of the scattered electron
213 //
214
215 G4double eKinEnergy = gamEnergy0 - gamEnergy1;
216
217 if(eKinEnergy > DBL_MIN) {
218 G4ThreeVector eDirection = gamEnergy0*gamDirection0 - gamEnergy1*gamDirection1;
219 eDirection = eDirection.unit();
220
221 // create G4DynamicParticle object for the electron.
222 G4DynamicParticle* dp = new G4DynamicParticle(theElectron,eDirection,eKinEnergy);
223 fvect->push_back(dp);
224 }
225}
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
#define DBL_MIN
Definition: templates.hh:75

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForGamma* G4KleinNishinaCompton::fParticleChange
protected

◆ lowestGammaEnergy

G4double G4KleinNishinaCompton::lowestGammaEnergy
protected

◆ theElectron

G4ParticleDefinition* G4KleinNishinaCompton::theElectron
protected

◆ theGamma

G4ParticleDefinition* G4KleinNishinaCompton::theGamma
protected

Definition at line 86 of file G4KleinNishinaCompton.hh.

Referenced by G4KleinNishinaCompton().


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