Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Attributes
G4KleinNishinaCompton Class Reference

#include <G4KleinNishinaCompton.hh>

Inheritance diagram for G4KleinNishinaCompton:
G4VEmModel G4PolarizedComptonModel MyKleinNishinaCompton

Public Member Functions

 G4KleinNishinaCompton (const G4ParticleDefinition *p=0, 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 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 *)
 
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
 

Protected Attributes

G4ParticleDefinitiontheGamma
 
G4ParticleDefinitiontheElectron
 
G4ParticleChangeForGammafParticleChange
 
G4double lowestGammaEnergy
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

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 ( const G4ParticleDefinition p = 0,
const G4String nam = "Klein-Nishina" 
)

Definition at line 70 of file G4KleinNishinaCompton.cc.

References G4Electron::Electron(), python.hepunit::eV, fParticleChange, G4Gamma::Gamma(), lowestGammaEnergy, theElectron, and theGamma.

72  : G4VEmModel(nam)
73 {
76  lowestGammaEnergy = 1.0*eV;
77  fParticleChange = 0;
78 }
G4ParticleChangeForGamma * fParticleChange
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4ParticleDefinition * theElectron
G4ParticleDefinition * theGamma
G4KleinNishinaCompton::~G4KleinNishinaCompton ( )
virtual

Definition at line 82 of file G4KleinNishinaCompton.cc.

83 {}

Member Function Documentation

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

References test::a, test::b, test::c, plottest35::c1, python.hepunit::electron_mass_c2, G4Exp(), G4Log(), python.hepunit::keV, and G4INCL::Math::max().

Referenced by G4PolarizedComptonModel::ComputeCrossSectionPerAtom().

109 {
110  G4double xSection = 0.0 ;
111  if ( Z < 0.9999 ) return xSection;
112  if ( GammaEnergy < 0.1*keV ) return xSection;
113  // if ( GammaEnergy > (100.*GeV/Z) ) return xSection;
114 
115  static const G4double a = 20.0 , b = 230.0 , c = 440.0;
116 
117  G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
118  p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
119 
120  G4double T0 = 15.0*keV;
121  if (Z < 1.5) T0 = 40.0*keV;
122 
123  G4double X = max(GammaEnergy, T0) / electron_mass_c2;
124  xSection = p1Z*G4Log(1.+2.*X)/X
125  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
126 
127  // modification for low energy. (special case for Hydrogen)
128  if (GammaEnergy < T0) {
129  G4double dT0 = 1.*keV;
130  X = (T0+dT0) / electron_mass_c2 ;
131  G4double sigma = p1Z*G4Log(1.+2*X)/X
132  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
133  G4double c1 = -T0*(sigma-xSection)/(xSection*dT0);
134  G4double c2 = 0.150;
135  if (Z > 1.5) c2 = 0.375-0.0556*G4Log(Z);
136  G4double y = G4Log(GammaEnergy/T0);
137  xSection *= G4Exp(-y*(c1+c2*y));
138  }
139  // G4cout << "e= " << GammaEnergy << " Z= " << Z << " cross= " << xSection << G4endl;
140  return xSection;
141 }
float electron_mass_c2
Definition: hepunit.py:274
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
T max(const T t1, const T t2)
brief Return the largest of the two arguments
double G4double
Definition: G4Types.hh:76
tuple c1
Definition: plottest35.py:14
void G4KleinNishinaCompton::Initialise ( const G4ParticleDefinition p,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 87 of file G4KleinNishinaCompton.cc.

References fParticleChange, G4VEmModel::GetParticleChangeForGamma(), G4VEmModel::InitialiseElementSelectors(), and G4VEmModel::IsMaster().

89 {
90  if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
92 }
G4ParticleChangeForGamma * fParticleChange
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:135
G4bool IsMaster() const
Definition: G4VEmModel.hh:676
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
void G4KleinNishinaCompton::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 96 of file G4KleinNishinaCompton.cc.

References G4VEmModel::GetElementSelectors(), and G4VEmModel::SetElementSelectors().

98 {
100 }
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:760
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:768
void G4KleinNishinaCompton::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Reimplemented in G4PolarizedComptonModel, and MyKleinNishinaCompton.

Definition at line 145 of file G4KleinNishinaCompton.cc.

References DBL_MIN, python.hepunit::electron_mass_c2, fParticleChange, fStopAndKill, G4Exp(), G4Log(), G4UniformRand, G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), lowestGammaEnergy, G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), G4VParticleChange::ProposeTrackStatus(), CLHEP::Hep3Vector::rotateUz(), G4ParticleChangeForGamma::SetProposedKineticEnergy(), theElectron, python.hepunit::twopi, and CLHEP::Hep3Vector::unit().

150 {
151  // The scattered gamma energy is sampled according to Klein - Nishina formula.
152  // The random number techniques of Butcher & Messel are used
153  // (Nuc Phys 20(1960),15).
154  // Note : Effects due to binding of atomic electrons are negliged.
155 
156  G4double gamEnergy0 = aDynamicGamma->GetKineticEnergy();
157 
158  // extra protection
159  if(gamEnergy0 < lowestGammaEnergy) {
163  return;
164  }
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 = 0.5*(1.- epsilon0sq);
180 
181  do {
182  if ( alpha1/(alpha1+alpha2) > G4UniformRand() ) {
183  epsilon = G4Exp(-alpha1*G4UniformRand()); // eps0**r
184  epsilonsq = epsilon*epsilon;
185 
186  } else {
187  epsilonsq = epsilon0sq + (1.- epsilon0sq)*G4UniformRand();
188  epsilon = sqrt(epsilonsq);
189  };
190 
191  onecost = (1.- epsilon)/(epsilon*E0_m);
192  sint2 = onecost*(2.-onecost);
193  greject = 1. - epsilon*sint2/(1.+ epsilonsq);
194 
195  } while (greject < G4UniformRand());
196 
197  //
198  // scattered gamma angles. ( Z - axis along the parent gamma)
199  //
200 
201  if(sint2 < 0.0) { sint2 = 0.0; }
202  G4double cosTeta = 1. - onecost;
203  G4double sinTeta = sqrt (sint2);
204  G4double Phi = twopi * G4UniformRand();
205 
206  //
207  // update G4VParticleChange for the scattered gamma
208  //
209 
210  G4ThreeVector gamDirection1(sinTeta*cos(Phi), sinTeta*sin(Phi), cosTeta);
211  gamDirection1.rotateUz(gamDirection0);
212  G4double gamEnergy1 = epsilon*gamEnergy0;
213  if(gamEnergy1 > lowestGammaEnergy) {
216  } else {
220  }
221 
222  //
223  // kinematic of the scattered electron
224  //
225 
226  G4double eKinEnergy = gamEnergy0 - gamEnergy1;
227 
228  if(eKinEnergy > DBL_MIN) {
229  G4ThreeVector eDirection = gamEnergy0*gamDirection0 - gamEnergy1*gamDirection1;
230  eDirection = eDirection.unit();
231 
232  // create G4DynamicParticle object for the electron.
233  G4DynamicParticle* dp = new G4DynamicParticle(theElectron,eDirection,eKinEnergy);
234  fvect->push_back(dp);
235  }
236 }
G4ParticleChangeForGamma * fParticleChange
G4double GetKineticEnergy() const
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
#define G4UniformRand()
Definition: Randomize.hh:87
const G4ThreeVector & GetMomentumDirection() const
float electron_mass_c2
Definition: hepunit.py:274
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:180
Hep3Vector unit() const
#define DBL_MIN
Definition: templates.hh:75
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4ParticleDefinition * theElectron
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)

Field Documentation

G4ParticleChangeForGamma* G4KleinNishinaCompton::fParticleChange
protected
G4double G4KleinNishinaCompton::lowestGammaEnergy
protected
G4ParticleDefinition* G4KleinNishinaCompton::theElectron
protected
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: