Geant4  10.02.p03
G4KleinNishinaCompton Class Reference

#include <G4KleinNishinaCompton.hh>

Inheritance diagram for G4KleinNishinaCompton:
Collaboration 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 &)
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
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 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, const G4ParticleDefinition *, G4double)
 
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 *> *)
 
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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
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
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
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 SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () 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)
 

Protected Attributes

G4ParticleDefinitiontheGamma
 
G4ParticleDefinitiontheElectron
 
G4ParticleChangeForGamma * fParticleChange
 
G4double lowestSecondaryEnergy
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChange * pParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Private Member Functions

G4KleinNishinaComptonoperator= (const G4KleinNishinaCompton &right)
 
 G4KleinNishinaCompton (const G4KleinNishinaCompton &)
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLoss * GetParticleChangeForLoss ()
 
G4ParticleChangeForGamma * GetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 58 of file G4KleinNishinaCompton.hh.

Constructor & Destructor Documentation

◆ G4KleinNishinaCompton() [1/2]

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

Definition at line 74 of file G4KleinNishinaCompton.cc.

76  : G4VEmModel(nam)
77 {
80  lowestSecondaryEnergy = 100.0*eV;
81  fParticleChange = nullptr;
82 }
G4ParticleChangeForGamma * fParticleChange
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:69
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static const double eV
Definition: G4SIunits.hh:212
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4ParticleDefinition * theElectron
G4ParticleDefinition * theGamma
Here is the call graph for this function:

◆ ~G4KleinNishinaCompton()

G4KleinNishinaCompton::~G4KleinNishinaCompton ( )
virtual

Definition at line 86 of file G4KleinNishinaCompton.cc.

87 {}

◆ G4KleinNishinaCompton() [2/2]

G4KleinNishinaCompton::G4KleinNishinaCompton ( const G4KleinNishinaCompton )
private

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 110 of file G4KleinNishinaCompton.cc.

115 {
116  G4double xSection = 0.0 ;
117  if (GammaEnergy <= LowEnergyLimit()) { return xSection; }
118 
119  static const G4double a = 20.0 , b = 230.0 , c = 440.0;
120 
121  G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
122  p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
123 
124  G4double T0 = 15.0*keV;
125  if (Z < 1.5) { T0 = 40.0*keV; }
126 
127  G4double X = max(GammaEnergy, T0) / electron_mass_c2;
128  xSection = p1Z*G4Log(1.+2.*X)/X
129  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
130 
131  // modification for low energy. (special case for Hydrogen)
132  if (GammaEnergy < T0) {
133  static const G4double dT0 = keV;
134  X = (T0+dT0) / electron_mass_c2 ;
135  G4double sigma = p1Z*G4Log(1.+2*X)/X
136  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
137  G4double c1 = -T0*(sigma-xSection)/(xSection*dT0);
138  G4double c2 = 0.150;
139  if (Z > 1.5) { c2 = 0.375-0.0556*G4Log(Z); }
140  G4double y = G4Log(GammaEnergy/T0);
141  xSection *= G4Exp(-y*(c1+c2*y));
142  }
143  // G4cout<<"e= "<< GammaEnergy<<" Z= "<<Z<<" cross= " << xSection << G4endl;
144  return xSection;
145 }
static const G4double d3
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
static const G4double d1
static const G4double f2
static const G4double e2
static const G4double e4
Float_t X
Double_t y
Float_t Z
static const G4double f4
float electron_mass_c2
Definition: hepunit.py:274
static const G4double f1
static const G4double e1
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static const G4double f3
TCanvas * c2
Definition: plot_hist.C:75
static const double keV
Definition: G4SIunits.hh:213
double G4double
Definition: G4Types.hh:76
static const G4double e3
static const G4double d2
static const G4double dT0
static const G4double d4
Here is the call graph for this function:
Here is the caller graph for this function:

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 91 of file G4KleinNishinaCompton.cc.

93 {
94  if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
95  if(nullptr == fParticleChange) {
97  }
98 }
G4ParticleChangeForGamma * fParticleChange
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148
G4bool IsMaster() const
Definition: G4VEmModel.hh:718
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:134
Here is the call graph for this function:

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 102 of file G4KleinNishinaCompton.cc.

104 {
106 }
void SetElementSelectors(std::vector< G4EmElementSelector *> *)
Definition: G4VEmModel.hh:810
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:802
Here is the call graph for this function:

◆ operator=()

G4KleinNishinaCompton& G4KleinNishinaCompton::operator= ( const G4KleinNishinaCompton right)
private

◆ SampleSecondaries()

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

Implements G4VEmModel.

Reimplemented in G4PolarizedComptonModel, G4HeatedKleinNishinaCompton, and MyKleinNishinaCompton.

Definition at line 149 of file G4KleinNishinaCompton.cc.

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

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForGamma* G4KleinNishinaCompton::fParticleChange
protected

Definition at line 91 of file G4KleinNishinaCompton.hh.

◆ lowestSecondaryEnergy

G4double G4KleinNishinaCompton::lowestSecondaryEnergy
protected

Definition at line 92 of file G4KleinNishinaCompton.hh.

◆ theElectron

G4ParticleDefinition* G4KleinNishinaCompton::theElectron
protected

Definition at line 90 of file G4KleinNishinaCompton.hh.

◆ theGamma

G4ParticleDefinition* G4KleinNishinaCompton::theGamma
protected

Definition at line 89 of file G4KleinNishinaCompton.hh.


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