Geant4 10.7.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=nullptr, const G4String &nam="Klein-Nishina")
 
virtual ~G4KleinNishinaCompton ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) override
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
- 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 void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
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 GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, 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 *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
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
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleDefinitiontheGamma
 
G4ParticleDefinitiontheElectron
 
G4ParticleChangeForGammafParticleChange
 
G4double lowestSecondaryEnergy
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

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 57 of file G4KleinNishinaCompton.hh.

Constructor & Destructor Documentation

◆ G4KleinNishinaCompton()

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

Definition at line 64 of file G4KleinNishinaCompton.cc.

66 : G4VEmModel(nam)
67{
70 lowestSecondaryEnergy = 100.0*eV;
71 fParticleChange = nullptr;
72}
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
G4ParticleDefinition * theGamma
G4ParticleChangeForGamma * fParticleChange
G4ParticleDefinition * theElectron

◆ ~G4KleinNishinaCompton()

G4KleinNishinaCompton::~G4KleinNishinaCompton ( )
virtual

Definition at line 76 of file G4KleinNishinaCompton.cc.

77{}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Reimplemented in G4PolarizedComptonModel.

Definition at line 100 of file G4KleinNishinaCompton.cc.

105{
106 G4double xSection = 0.0 ;
107 if (GammaEnergy <= LowEnergyLimit()) { return xSection; }
108
109 static const G4double a = 20.0 , b = 230.0 , c = 440.0;
110
111 static const G4double
112 d1= 2.7965e-1*CLHEP::barn, d2=-1.8300e-1*CLHEP::barn,
113 d3= 6.7527 *CLHEP::barn, d4=-1.9798e+1*CLHEP::barn,
114 e1= 1.9756e-5*CLHEP::barn, e2=-1.0205e-2*CLHEP::barn,
115 e3=-7.3913e-2*CLHEP::barn, e4= 2.7079e-2*CLHEP::barn,
116 f1=-3.9178e-7*CLHEP::barn, f2= 6.8241e-5*CLHEP::barn,
117 f3= 6.0480e-5*CLHEP::barn, f4= 3.0274e-4*CLHEP::barn;
118
119 G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
120 p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
121
122 G4double T0 = 15.0*keV;
123 if (Z < 1.5) { T0 = 40.0*keV; }
124
125 G4double X = max(GammaEnergy, T0) / electron_mass_c2;
126 xSection = p1Z*G4Log(1.+2.*X)/X
127 + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
128
129 // modification for low energy. (special case for Hydrogen)
130 if (GammaEnergy < T0) {
131 static const G4double dT0 = keV;
132 X = (T0+dT0) / electron_mass_c2 ;
133 G4double sigma = p1Z*G4Log(1.+2*X)/X
134 + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
135 G4double c1 = -T0*(sigma-xSection)/(xSection*dT0);
136 G4double c2 = 0.150;
137 if (Z > 1.5) { c2 = 0.375-0.0556*G4Log(Z); }
138 G4double y = G4Log(GammaEnergy/T0);
139 xSection *= G4Exp(-y*(c1+c2*y));
140 }
141 // G4cout<<"e= "<< GammaEnergy<<" Z= "<<Z<<" cross= " << xSection << G4endl;
142 return xSection;
143}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
double G4double
Definition: G4Types.hh:83
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
T max(const T t1, const T t2)
brief Return the largest of the two arguments

Referenced by G4PolarizedComptonModel::ComputeCrossSectionPerAtom().

◆ Initialise()

void G4KleinNishinaCompton::Initialise ( const G4ParticleDefinition p,
const G4DataVector cuts 
)
overridevirtual

Implements G4VEmModel.

Definition at line 81 of file G4KleinNishinaCompton.cc.

83{
84 if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
85 if(nullptr == fParticleChange) {
87 }
88}
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
G4bool IsMaster() const
Definition: G4VEmModel.hh:736
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148

◆ InitialiseLocal()

void G4KleinNishinaCompton::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 92 of file G4KleinNishinaCompton.cc.

94{
96}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:842
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:834

◆ SampleSecondaries()

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

Implements G4VEmModel.

Reimplemented in G4HeatedKleinNishinaCompton, and G4PolarizedComptonModel.

Definition at line 147 of file G4KleinNishinaCompton.cc.

153{
154 // The scattered gamma energy is sampled according to Klein - Nishina formula.
155 // The random number techniques of Butcher & Messel are used
156 // (Nuc Phys 20(1960),15).
157 // Note : Effects due to binding of atomic electrons are negliged.
158
159 G4double gamEnergy0 = aDynamicGamma->GetKineticEnergy();
160
161 // do nothing below the threshold
162 if(gamEnergy0 <= LowEnergyLimit()) { return; }
163
164 G4double E0_m = gamEnergy0 / electron_mass_c2 ;
165
166 G4ThreeVector gamDirection0 = aDynamicGamma->GetMomentumDirection();
167
168 //
169 // sample the energy rate of the scattered gamma
170 //
171
172 G4double epsilon, epsilonsq, onecost, sint2, greject ;
173
174 G4double eps0 = 1./(1. + 2.*E0_m);
175 G4double epsilon0sq = eps0*eps0;
176 G4double alpha1 = - G4Log(eps0);
177 G4double alpha2 = alpha1 + 0.5*(1.- epsilon0sq);
178
179 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
180 G4double rndm[3];
181
182 static const G4int nlooplim = 1000;
183 G4int nloop = 0;
184 do {
185 ++nloop;
186 // false interaction if too many iterations
187 if(nloop > nlooplim) { return; }
188
189 // 3 random numbers to sample scattering
190 rndmEngineMod->flatArray(3, rndm);
191
192 if ( alpha1 > alpha2*rndm[0] ) {
193 epsilon = G4Exp(-alpha1*rndm[1]); // eps0**r
194 epsilonsq = epsilon*epsilon;
195
196 } else {
197 epsilonsq = epsilon0sq + (1.- epsilon0sq)*rndm[1];
198 epsilon = sqrt(epsilonsq);
199 };
200
201 onecost = (1.- epsilon)/(epsilon*E0_m);
202 sint2 = onecost*(2.-onecost);
203 greject = 1. - epsilon*sint2/(1.+ epsilonsq);
204
205 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
206 } while (greject < rndm[2]);
207
208 //
209 // scattered gamma angles. ( Z - axis along the parent gamma)
210 //
211
212 if(sint2 < 0.0) { sint2 = 0.0; }
213 G4double cosTeta = 1. - onecost;
214 G4double sinTeta = sqrt (sint2);
215 G4double Phi = twopi * rndmEngineMod->flat();
216
217 //
218 // update G4VParticleChange for the scattered gamma
219 //
220
221 G4ThreeVector gamDirection1(sinTeta*cos(Phi), sinTeta*sin(Phi), cosTeta);
222 gamDirection1.rotateUz(gamDirection0);
223 G4double gamEnergy1 = epsilon*gamEnergy0;
224 G4double edep = 0.0;
225 if(gamEnergy1 > lowestSecondaryEnergy) {
228 } else {
231 edep = gamEnergy1;
232 }
233
234 //
235 // kinematic of the scattered electron
236 //
237
238 G4double eKinEnergy = gamEnergy0 - gamEnergy1;
239
240 if(eKinEnergy > lowestSecondaryEnergy) {
241 G4ThreeVector eDirection = gamEnergy0*gamDirection0 - gamEnergy1*gamDirection1;
242 eDirection = eDirection.unit();
243
244 // create G4DynamicParticle object for the electron.
245 G4DynamicParticle* dp = new G4DynamicParticle(theElectron,eDirection,eKinEnergy);
246 fvect->push_back(dp);
247 } else {
248 edep += eKinEnergy;
249 }
250 // energy balance
251 if(edep > 0.0) {
253 }
254}
double epsilon(double density, double temperature)
@ fStopAndKill
int G4int
Definition: G4Types.hh:85
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)

Member Data Documentation

◆ fParticleChange

◆ lowestSecondaryEnergy

G4double G4KleinNishinaCompton::lowestSecondaryEnergy
protected

◆ theElectron

◆ theGamma

G4ParticleDefinition* G4KleinNishinaCompton::theGamma
protected

Definition at line 89 of file G4KleinNishinaCompton.hh.

Referenced by G4KleinNishinaCompton().


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