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

#include <G4KleinNishinaModel.hh>

+ Inheritance diagram for G4KleinNishinaModel:

Public Member Functions

 G4KleinNishinaModel (const G4String &nam="KleinNishina")
 
virtual ~G4KleinNishinaModel ()
 
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 58 of file G4KleinNishinaModel.hh.

Constructor & Destructor Documentation

◆ G4KleinNishinaModel()

G4KleinNishinaModel::G4KleinNishinaModel ( const G4String nam = "KleinNishina")
explicit

Definition at line 66 of file G4KleinNishinaModel.cc.

67 : G4VEmModel(nam),
68 lv1(0.,0.,0.,0.),
69 lv2(0.,0.,0.,0.),
70 bst(0.,0.,0.)
71{
75 limitFactor = 4;
76 fProbabilities.resize(9,0.0);
78 fParticleChange = nullptr;
79 fAtomDeexcitation = nullptr;
80}
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
G4ParticleDefinition * theElectron
G4ParticleChangeForGamma * fParticleChange
G4ParticleDefinition * theGamma
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:813

◆ ~G4KleinNishinaModel()

G4KleinNishinaModel::~G4KleinNishinaModel ( )
virtual

Definition at line 84 of file G4KleinNishinaModel.cc.

85{}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 110 of file G4KleinNishinaModel.cc.

114{
115 G4double xSection = 0.0 ;
116 if (gammaEnergy <= LowEnergyLimit()) { return xSection; }
117
118 static const G4double a = 20.0 , b = 230.0 , c = 440.0;
119
120static const G4double
121 d1= 2.7965e-1*CLHEP::barn, d2=-1.8300e-1*CLHEP::barn,
122 d3= 6.7527 *CLHEP::barn, d4=-1.9798e+1*CLHEP::barn,
123 e1= 1.9756e-5*CLHEP::barn, e2=-1.0205e-2*CLHEP::barn,
124 e3=-7.3913e-2*CLHEP::barn, e4= 2.7079e-2*CLHEP::barn,
125 f1=-3.9178e-7*CLHEP::barn, f2= 6.8241e-5*CLHEP::barn,
126 f3= 6.0480e-5*CLHEP::barn, f4= 3.0274e-4*CLHEP::barn;
127
128 G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
129 p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
130
131 G4double T0 = 15.0*keV;
132 if (Z < 1.5) { T0 = 40.0*keV; }
133
134 G4double X = max(gammaEnergy, T0) / electron_mass_c2;
135 xSection = p1Z*G4Log(1.+2.*X)/X
136 + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
137
138 // modification for low energy. (special case for Hydrogen)
139 static const G4double dT0 = keV;
140 if (gammaEnergy < T0) {
141 X = (T0+dT0) / electron_mass_c2 ;
142 G4double sigma = p1Z*G4Log(1.+2*X)/X
143 + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
144 G4double c1 = -T0*(sigma-xSection)/(xSection*dT0);
145 G4double c2 = 0.150;
146 if (Z > 1.5) { c2 = 0.375-0.0556*G4Log(Z); }
147 G4double y = G4Log(gammaEnergy/T0);
148 xSection *= G4Exp(-y*(c1+c2*y));
149 }
150
151 if(xSection < 0.0) { xSection = 0.0; }
152 // G4cout << "e= " << GammaEnergy << " Z= " << Z
153 // << " cross= " << xSection << G4endl;
154 return xSection;
155}
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

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 89 of file G4KleinNishinaModel.cc.

91{
92 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
93 if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
94 if(nullptr == fParticleChange) {
96 }
97}
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
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 G4KleinNishinaModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 101 of file G4KleinNishinaModel.cc.

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 159 of file G4KleinNishinaModel.cc.

165{
166 // primary gamma
167 G4double energy = aDynamicGamma->GetKineticEnergy();
168
169 // do nothing below the threshold
170 if(energy <= LowEnergyLimit()) { return; }
171
172 G4ThreeVector direction = aDynamicGamma->GetMomentumDirection();
173
174 // select atom
175 const G4Element* elm = SelectRandomAtom(couple, theGamma, energy);
176
177 // select shell first
178 G4int nShells = elm->GetNbOfAtomicShells();
179 if(nShells > (G4int)fProbabilities.size()) { fProbabilities.resize(nShells); }
180 G4double totprob = 0.0;
181 G4int i;
182 for(i=0; i<nShells; ++i) {
183 //G4double bindingEnergy = elm->GetAtomicShell(i);
184 totprob += elm->GetNbOfShellElectrons(i);
185 //totprob += elm->GetNbOfShellElectrons(i)/(bindingEnergy*bindingEnergy);
186 fProbabilities[i] = totprob;
187 }
188
189 // Loop on sampling
190 static const G4int nlooplim = 1000;
191 G4int nloop = 0;
192
193 G4double bindingEnergy, ePotEnergy, eKinEnergy;
194 G4double gamEnergy0, gamEnergy1;
195
196 CLHEP::HepRandomEngine* rndmEngineMod = G4Random::getTheEngine();
197 G4double rndm[4];
198
199 do {
200 ++nloop;
201
202 // 4 random numbers to select e-
203 rndmEngineMod->flatArray(4, rndm);
204 G4double xprob = totprob*rndm[0];
205
206 // select shell
207 for(i=0; i<nShells; ++i) { if(xprob <= fProbabilities[i]) { break; } }
208
210 lv1.set(0.0,0.0,energy,energy);
211 /*
212 G4cout << "nShells= " << nShells << " i= " << i
213 << " Egamma= " << energy << " Ebind= " << bindingEnergy
214 << G4endl;
215 */
216 // for rest frame of the electron
217 G4double x = -G4Log(rndm[1]);
218 eKinEnergy = bindingEnergy*x;
219 ePotEnergy = bindingEnergy*(1.0 + x);
220
221 // for rest frame of the electron
222 G4double eTotMomentum = sqrt(eKinEnergy*(eKinEnergy + 2*electron_mass_c2));
223 G4double phi = rndm[2]*twopi;
224 G4double costet = 2*rndm[3] - 1;
225 G4double sintet = sqrt((1 - costet)*(1 + costet));
226 lv2.set(eTotMomentum*sintet*cos(phi),eTotMomentum*sintet*sin(phi),
227 eTotMomentum*costet,eKinEnergy + electron_mass_c2);
228 bst = lv2.boostVector();
229 lv1.boost(-bst);
230
231 gamEnergy0 = lv1.e();
232
233 // In the rest frame of the electron
234 // The scattered gamma energy is sampled according to Klein-Nishina formula
235 // The random number techniques of Butcher & Messel are used
236 // (Nuc Phys 20(1960),15).
237 G4double E0_m = gamEnergy0/electron_mass_c2;
238
239 //G4cout << "Nloop= "<< nloop << " Ecm(keV)= " << gamEnergy0/keV << G4endl;
240 //
241 // sample the energy rate of the scattered gamma
242 //
243
244 G4double epsilon, epsilonsq, onecost, sint2, greject ;
245
246 G4double eps0 = 1./(1 + 2*E0_m);
247 G4double epsilon0sq = eps0*eps0;
248 G4double alpha1 = - G4Log(eps0);
249 G4double alpha2 = alpha1 + 0.5*(1 - epsilon0sq);
250
251 do {
252 ++nloop;
253 // false interaction if too many iterations
254 if(nloop > nlooplim) { return; }
255
256 // 3 random numbers to sample scattering
257 rndmEngineMod->flatArray(3, rndm);
258
259 if ( alpha1 > alpha2*rndm[0] ) {
260 epsilon = G4Exp(-alpha1*rndm[1]); // epsilon0**r
261 epsilonsq = epsilon*epsilon;
262
263 } else {
264 epsilonsq = epsilon0sq + (1.- epsilon0sq)*rndm[1];
265 epsilon = sqrt(epsilonsq);
266 }
267
268 onecost = (1.- epsilon)/(epsilon*E0_m);
269 sint2 = onecost*(2.-onecost);
270 greject = 1. - epsilon*sint2/(1.+ epsilonsq);
271
272 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
273 } while (greject < rndm[2]);
274 gamEnergy1 = epsilon*gamEnergy0;
275
276 // before scattering total 4-momentum in e- system
277 lv2.set(0.0,0.0,0.0,electron_mass_c2);
278 lv2 += lv1;
279
280 //
281 // scattered gamma angles. ( Z - axis along the parent gamma)
282 //
283 if(sint2 < 0.0) { sint2 = 0.0; }
284 costet = 1. - onecost;
285 sintet = sqrt(sint2);
286 phi = twopi * rndmEngineMod->flat();
287
288 // e- recoil
289 //
290 // in rest frame of the electron
291 G4ThreeVector gamDir = lv1.vect().unit();
292 G4ThreeVector v = G4ThreeVector(sintet*cos(phi),sintet*sin(phi),costet);
293 v.rotateUz(gamDir);
294 lv1.set(gamEnergy1*v.x(),gamEnergy1*v.y(),gamEnergy1*v.z(),gamEnergy1);
295 lv2 -= lv1;
296 //G4cout<<"Egam(keV)= " << lv1.e()/keV
297 // <<" Ee(keV)= " << (lv2.e()-electron_mass_c2)/keV << G4endl;
298 lv2.boost(bst);
299 eKinEnergy = lv2.e() - electron_mass_c2 - ePotEnergy;
300 //G4cout << "Nloop= " << nloop << " eKinEnergy= " << eKinEnergy << G4endl;
301
302 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
303 } while ( eKinEnergy < 0.0 );
304
305 //
306 // update G4VParticleChange for the scattered gamma
307 //
308
309 lv1.boost(bst);
310 gamEnergy1 = lv1.e();
311 if(gamEnergy1 > lowestSecondaryEnergy) {
312 G4ThreeVector gamDirection1 = lv1.vect().unit();
313 gamDirection1.rotateUz(direction);
315 } else {
317 gamEnergy1 = 0.0;
318 }
320
321 //
322 // kinematic of the scattered electron
323 //
324
325 if(eKinEnergy > lowestSecondaryEnergy) {
326 G4ThreeVector eDirection = lv2.vect().unit();
327 eDirection.rotateUz(direction);
328 G4DynamicParticle* dp =
329 new G4DynamicParticle(theElectron,eDirection,eKinEnergy);
330 fvect->push_back(dp);
331 } else { eKinEnergy = 0.0; }
332
333 G4double edep = energy - gamEnergy1 - eKinEnergy;
334 G4double esec = 0.0;
335
336 // sample deexcitation
337 //
338 if(fAtomDeexcitation) {
339 G4int index = couple->GetIndex();
340 if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
341 G4int Z = elm->GetZasInt();
343 const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
344 G4int nbefore = fvect->size();
345 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
346 G4int nafter = fvect->size();
347 //G4cout << "N1= " << nbefore << " N2= " << nafter << G4endl;
348 for (G4int j=nbefore; j<nafter; ++j) {
349 G4double e = ((*fvect)[j])->GetKineticEnergy();
350 if(esec + e > edep) {
351 // correct energy in order to have energy balance
352 e = edep - esec;
353 ((*fvect)[j])->SetKineticEnergy(e);
354 esec += e;
355 /*
356 G4cout << "### G4KleinNishinaModel Edep(eV)= " << edep/eV
357 << " Esec(eV)= " << esec/eV
358 << " E["<< j << "](eV)= " << e/eV
359 << " N= " << nafter
360 << " Z= " << Z << " shell= " << i
361 << " Ebind(keV)= " << bindingEnergy/keV
362 << " Eshell(keV)= " << shell->BindingEnergy()/keV
363 << G4endl;
364 */
365 // delete the rest of secondaries (should not happens)
366 for (G4int jj=nafter-1; jj>j; --jj) {
367 delete (*fvect)[jj];
368 fvect->pop_back();
369 }
370 break;
371 }
372 esec += e;
373 }
374 edep -= esec;
375 }
376 }
377 if(std::abs(energy - gamEnergy1 - eKinEnergy - esec - edep) > eV) {
378 G4cout << "### G4KleinNishinaModel dE(eV)= "
379 << (energy - gamEnergy1 - eKinEnergy - esec - edep)/eV
380 << " shell= " << i
381 << " E(keV)= " << energy/keV
382 << " Ebind(keV)= " << bindingEnergy/keV
383 << " Eg(keV)= " << gamEnergy1/keV
384 << " Ee(keV)= " << eKinEnergy/keV
385 << " Esec(keV)= " << esec/keV
386 << " Edep(keV)= " << edep/keV
387 << G4endl;
388 }
389 // energy balance
390 if(edep > 0.0) {
392 }
393}
G4AtomicShellEnumerator
double epsilon(double density, double temperature)
CLHEP::Hep3Vector G4ThreeVector
@ fStopAndKill
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
Hep3Vector unit() const
double x() const
double y() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
void set(double x, double y, double z, double t)
virtual double flat()=0
virtual void flatArray(const int size, double *vect)=0
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4int GetNbOfAtomicShells() const
Definition: G4Element.hh:146
G4int GetZasInt() const
Definition: G4Element.hh:131
G4int GetNbOfShellElectrons(G4int index) const
Definition: G4Element.cc:381
G4double GetAtomicShell(G4int index) const
Definition: G4Element.cc:366
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:570
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double energy(const ThreeVector &p, const G4double m)
G4double bindingEnergy(G4int A, G4int Z)

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForGamma* G4KleinNishinaModel::fParticleChange
protected

Definition at line 91 of file G4KleinNishinaModel.hh.

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

◆ lowestSecondaryEnergy

G4double G4KleinNishinaModel::lowestSecondaryEnergy
protected

Definition at line 92 of file G4KleinNishinaModel.hh.

Referenced by G4KleinNishinaModel(), and SampleSecondaries().

◆ theElectron

G4ParticleDefinition* G4KleinNishinaModel::theElectron
protected

Definition at line 90 of file G4KleinNishinaModel.hh.

Referenced by G4KleinNishinaModel(), and SampleSecondaries().

◆ theGamma

G4ParticleDefinition* G4KleinNishinaModel::theGamma
protected

Definition at line 89 of file G4KleinNishinaModel.hh.

Referenced by G4KleinNishinaModel(), and SampleSecondaries().


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