Geant4 11.2.2
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")
 
 ~G4KleinNishinaModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4KleinNishinaModeloperator= (const G4KleinNishinaModel &right)=delete
 
 G4KleinNishinaModel (const G4KleinNishinaModel &)=delete
 
- 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 level, const G4ParticleDefinition *, G4double kineticEnergy)
 
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 *, const G4double &length, G4double &eloss)
 
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 FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
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)
 
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)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
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 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 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 *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
void SetLPMFlag (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleDefinitiontheGamma
 
G4ParticleDefinitiontheElectron
 
G4ParticleChangeForGammafParticleChange
 
G4double lowestSecondaryEnergy
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
std::size_t currentCoupleIndex = 0
 
std::size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

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() [1/2]

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:91
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
G4ParticleDefinition * theElectron
G4ParticleChangeForGamma * fParticleChange
G4ParticleDefinition * theGamma
void SetDeexcitationFlag(G4bool val)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

◆ ~G4KleinNishinaModel()

G4KleinNishinaModel::~G4KleinNishinaModel ( )
overridedefault

◆ G4KleinNishinaModel() [2/2]

G4KleinNishinaModel::G4KleinNishinaModel ( const G4KleinNishinaModel & )
delete

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 109 of file G4KleinNishinaModel.cc.

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

90{
91 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
92 if(IsMaster()) { InitialiseElementSelectors(p, cuts); }
93 if(nullptr == fParticleChange) {
95 }
96}
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4bool IsMaster() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)

◆ InitialiseLocal()

void G4KleinNishinaModel::InitialiseLocal ( const G4ParticleDefinition * ,
G4VEmModel * masterModel )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 100 of file G4KleinNishinaModel.cc.

102{
104}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
std::vector< G4EmElementSelector * > * GetElementSelectors()

◆ operator=()

G4KleinNishinaModel & G4KleinNishinaModel::operator= ( const G4KleinNishinaModel & right)
delete

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 158 of file G4KleinNishinaModel.cc.

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

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

◆ lowestSecondaryEnergy

G4double G4KleinNishinaModel::lowestSecondaryEnergy
protected

Definition at line 95 of file G4KleinNishinaModel.hh.

Referenced by G4KleinNishinaModel(), and SampleSecondaries().

◆ theElectron

G4ParticleDefinition* G4KleinNishinaModel::theElectron
protected

Definition at line 93 of file G4KleinNishinaModel.hh.

Referenced by G4KleinNishinaModel(), and SampleSecondaries().

◆ theGamma

G4ParticleDefinition* G4KleinNishinaModel::theGamma
protected

Definition at line 92 of file G4KleinNishinaModel.hh.

Referenced by G4KleinNishinaModel(), and SampleSecondaries().


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