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

#include <G4LivermoreComptonModel.hh>

+ Inheritance diagram for G4LivermoreComptonModel:

Public Member Functions

 G4LivermoreComptonModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="LivermoreCompton")
 
 ~G4LivermoreComptonModel () override
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
void InitialiseForElement (const G4ParticleDefinition *, G4int Z) override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4LivermoreComptonModeloperator= (const G4LivermoreComptonModel &right)=delete
 
 G4LivermoreComptonModel (const G4LivermoreComptonModel &)=delete
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
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
 

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 *)
 
- 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
 

Detailed Description

Definition at line 46 of file G4LivermoreComptonModel.hh.

Constructor & Destructor Documentation

◆ G4LivermoreComptonModel() [1/2]

G4LivermoreComptonModel::G4LivermoreComptonModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "LivermoreCompton" )
explicit

Definition at line 72 of file G4LivermoreComptonModel.cc.

73 : G4VEmModel(nam)
74{
75 fParticleChange = nullptr;
76 verboseLevel = 1;
77 // Verbosity scale:
78 // 0 = nothing
79 // 1 = warning for energy non-conservation
80 // 2 = details of energy budget
81 // 3 = calculation of cross sections, file openings, sampling of atoms
82 // 4 = entering in methods
83
84 if (verboseLevel > 1) {
85 G4cout << "Livermore Compton model is constructed " << G4endl;
86 }
87
88 // Mark this model as "applicable" for atomic deexcitation
90
91 fParticleChange = nullptr;
92 fAtomDeexcitation = nullptr;
93}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void SetDeexcitationFlag(G4bool val)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

◆ ~G4LivermoreComptonModel()

G4LivermoreComptonModel::~G4LivermoreComptonModel ( )
override

Definition at line 97 of file G4LivermoreComptonModel.cc.

98{
99 if (IsMaster()) {
100 delete shellData;
101 shellData = nullptr;
102 delete profileData;
103 profileData = nullptr;
104 for (G4int i = 0; i <= maxZ; ++i) {
105 if (data[i]) {
106 delete data[i];
107 data[i] = nullptr;
108 }
109 }
110 }
111}
int G4int
Definition G4Types.hh:85
G4bool IsMaster() const

◆ G4LivermoreComptonModel() [2/2]

G4LivermoreComptonModel::G4LivermoreComptonModel ( const G4LivermoreComptonModel & )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4LivermoreComptonModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition * ,
G4double kinEnergy,
G4double Z,
G4double A = 0,
G4double cut = 0,
G4double emax = DBL_MAX )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 235 of file G4LivermoreComptonModel.cc.

238{
239 if (verboseLevel > 3) {
240 G4cout << "G4LivermoreComptonModel::ComputeCrossSectionPerAtom()" << G4endl;
241 }
242 G4double cs = 0.0;
243
244 if (GammaEnergy < LowEnergyLimit()) {
245 return 0.0;
246 }
247
248 G4int intZ = G4lrint(Z);
249 if (intZ < 1 || intZ > maxZ) {
250 return cs;
251 }
252
253 G4PhysicsFreeVector* pv = data[intZ];
254
255 // if element was not initialised
256 // do initialisation safely for MT mode
257 if (pv == nullptr) {
258 InitialiseForElement(nullptr, intZ);
259 pv = data[intZ];
260 if (pv == nullptr) {
261 return cs;
262 }
263 }
264
265 auto n = G4int(pv->GetVectorLength() - 1);
266 G4double e1 = pv->Energy(0);
267 G4double e2 = pv->Energy(n);
268
269 if (GammaEnergy <= e1) {
270 cs = GammaEnergy / (e1 * e1) * pv->Value(e1);
271 }
272 else if (GammaEnergy <= e2) {
273 cs = pv->Value(GammaEnergy) / GammaEnergy;
274 }
275 else if (GammaEnergy > e2) {
276 cs = pv->Value(e2) / GammaEnergy;
277 }
278
279 return cs;
280}
double G4double
Definition G4Types.hh:83
void InitialiseForElement(const G4ParticleDefinition *, G4int Z) override
G4double Energy(const std::size_t index) const
G4double Value(const G4double energy, std::size_t &lastidx) const
std::size_t GetVectorLength() const
G4double LowEnergyLimit() const
int G4lrint(double ad)
Definition templates.hh:134

◆ Initialise()

void G4LivermoreComptonModel::Initialise ( const G4ParticleDefinition * particle,
const G4DataVector & cuts )
overridevirtual

Implements G4VEmModel.

Definition at line 115 of file G4LivermoreComptonModel.cc.

117{
118 if (verboseLevel > 1) {
119 G4cout << "Calling G4LivermoreComptonModel::Initialise()" << G4endl;
120 }
121
122 // Initialise element selector
123 if (IsMaster()) {
124 // Initialise element selector
125 InitialiseElementSelectors(particle, cuts);
126
127 // Access to elements
128 const G4ElementTable* elemTable = G4Element::GetElementTable();
129 size_t numElems = (*elemTable).size();
130 for (size_t ie = 0; ie < numElems; ++ie) {
131 const G4Element* elem = (*elemTable)[ie];
132 const G4int Z = std::min(maxZ, elem->GetZasInt());
133 if (data[Z] == nullptr) {
134 ReadData(Z);
135 }
136 }
137
138 // For Doppler broadening
139 if (shellData == nullptr) {
140 shellData = new G4ShellData();
141 shellData->SetOccupancyData();
142 G4String file = "/doppler/shell-doppler";
143 shellData->LoadData(file);
144 }
145 if (profileData == nullptr) {
146 profileData = new G4DopplerProfile();
147 }
148 }
149
150 if (verboseLevel > 2) {
151 G4cout << "Loaded cross section files" << G4endl;
152 }
153
154 if (verboseLevel > 1) {
155 G4cout << "G4LivermoreComptonModel is initialized " << G4endl
156 << "Energy range: " << LowEnergyLimit() / eV << " eV - " << HighEnergyLimit() / GeV
157 << " GeV" << G4endl;
158 }
159 //
160 if (isInitialised) {
161 return;
162 }
163
164 fParticleChange = GetParticleChangeForGamma();
165 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
166 isInitialised = true;
167}
std::vector< G4Element * > G4ElementTable
#define elem(i, j)
static G4ElementTable * GetElementTable()
Definition G4Element.cc:389
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
void SetOccupancyData()
void LoadData(const G4String &fileName)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double HighEnergyLimit() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)

◆ InitialiseForElement()

void G4LivermoreComptonModel::InitialiseForElement ( const G4ParticleDefinition * ,
G4int Z )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 551 of file G4LivermoreComptonModel.cc.

552{
553 G4AutoLock l(&LivermoreComptonModelMutex);
554 if (data[Z] == nullptr) {
555 ReadData(Z);
556 }
557 l.unlock();
558}

Referenced by ComputeCrossSectionPerAtom().

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 171 of file G4LivermoreComptonModel.cc.

172{
174}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
std::vector< G4EmElementSelector * > * GetElementSelectors()

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 284 of file G4LivermoreComptonModel.cc.

288{
289 // The scattered gamma energy is sampled according to Klein - Nishina
290 // formula then accepted or rejected depending on the Scattering Function
291 // multiplied by factor from Klein - Nishina formula.
292 // Expression of the angular distribution as Klein Nishina
293 // angular and energy distribution and Scattering fuctions is taken from
294 // D. E. Cullen "A simple model of photon transport" Nucl. Instr. Meth.
295 // Phys. Res. B 101 (1995). Method of sampling with form factors is different
296 // data are interpolated while in the article they are fitted.
297 // Reference to the article is from J. Stepanek New Photon, Positron
298 // and Electron Interaction Data for GEANT in Energy Range from 1 eV to 10
299 // TeV (draft).
300 // The random number techniques of Butcher & Messel are used
301 // (Nucl Phys 20(1960),15).
302
303 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
304
305 if (verboseLevel > 3) {
306 G4cout << "G4LivermoreComptonModel::SampleSecondaries() E(MeV)= " << photonEnergy0 / MeV
307 << " in " << couple->GetMaterial()->GetName() << G4endl;
308 }
309
310 // do nothing below the threshold
311 // should never get here because the XS is zero below the limit
312 if (photonEnergy0 < LowEnergyLimit()) return;
313
314 G4double e0m = photonEnergy0 / electron_mass_c2;
315 G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
316
317 // Select randomly one element in the current material
318 const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
319 const G4Element* elm = SelectRandomAtom(couple, particle, photonEnergy0);
320
321 G4int Z = elm->GetZasInt();
322
323 G4double epsilon0Local = 1. / (1. + 2. * e0m);
324 G4double epsilon0Sq = epsilon0Local * epsilon0Local;
325 G4double alpha1 = -G4Log(epsilon0Local);
326 G4double alpha2 = 0.5 * (1. - epsilon0Sq);
327
328 G4double wlPhoton = h_Planck * c_light / photonEnergy0;
329
330 // Sample the energy of the scattered photon
332 G4double epsilonSq;
333 G4double oneCosT;
334 G4double sinT2;
335 G4double gReject;
336
337 if (verboseLevel > 3) {
338 G4cout << "Started loop to sample gamma energy" << G4endl;
339 }
340
341 do {
342 if (alpha1 / (alpha1 + alpha2) > G4UniformRand()) {
343 epsilon = G4Exp(-alpha1 * G4UniformRand());
344 epsilonSq = epsilon * epsilon;
345 }
346 else {
347 epsilonSq = epsilon0Sq + (1. - epsilon0Sq) * G4UniformRand();
348 epsilon = std::sqrt(epsilonSq);
349 }
350
351 oneCosT = (1. - epsilon) / (epsilon * e0m);
352 sinT2 = oneCosT * (2. - oneCosT);
353 G4double x = std::sqrt(oneCosT / 2.) * cm / wlPhoton;
354 G4double scatteringFunction = ComputeScatteringFunction(x, Z);
355 gReject = (1. - epsilon * sinT2 / (1. + epsilonSq)) * scatteringFunction;
356
357 } while (gReject < G4UniformRand() * Z);
358
359 G4double cosTheta = 1. - oneCosT;
360 G4double sinTheta = std::sqrt(sinT2);
361 G4double phi = twopi * G4UniformRand();
362 G4double dirx = sinTheta * std::cos(phi);
363 G4double diry = sinTheta * std::sin(phi);
364 G4double dirz = cosTheta;
365
366 // Doppler broadening - Method based on:
367 // Y. Namito, S. Ban and H. Hirayama,
368 // "Implementation of the Doppler Broadening of a Compton-Scattered Photon
369 // into the EGS4 Code", NIM A 349, pp. 489-494, 1994
370
371 // Maximum number of sampling iterations
372 static G4int maxDopplerIterations = 1000;
373 G4double bindingE = 0.;
374 G4double photonEoriginal = epsilon * photonEnergy0;
375 G4double photonE = -1.;
376 G4int iteration = 0;
377 G4double eMax = photonEnergy0;
378
379 G4int shellIdx = 0;
380
381 if (verboseLevel > 3) {
382 G4cout << "Started loop to sample broading" << G4endl;
383 }
384
385 do {
386 ++iteration;
387 // Select shell based on shell occupancy
388 shellIdx = shellData->SelectRandomShell(Z);
389 bindingE = shellData->BindingEnergy(Z, shellIdx);
390
391 if (verboseLevel > 3) {
392 G4cout << "Shell ID= " << shellIdx << " Ebind(keV)= " << bindingE / keV << G4endl;
393 }
394
395 eMax = photonEnergy0 - bindingE;
396
397 // Randomly sample bound electron momentum
398 // (memento: the data set is in Atomic Units)
399 G4double pSample = profileData->RandomSelectMomentum(Z, shellIdx);
400 if (verboseLevel > 3) {
401 G4cout << "pSample= " << pSample << G4endl;
402 }
403 // Rescale from atomic units
404 G4double pDoppler = pSample * fine_structure_const;
405 G4double pDoppler2 = pDoppler * pDoppler;
406 G4double var2 = 1. + oneCosT * e0m;
407 G4double var3 = var2 * var2 - pDoppler2;
408 G4double var4 = var2 - pDoppler2 * cosTheta;
409 G4double var = var4 * var4 - var3 + pDoppler2 * var3;
410 if (var > 0.) {
411 G4double varSqrt = std::sqrt(var);
412 G4double scale = photonEnergy0 / var3;
413 // Random select either root
414 if (G4UniformRand() < 0.5) {
415 photonE = (var4 - varSqrt) * scale;
416 }
417 else {
418 photonE = (var4 + varSqrt) * scale;
419 }
420 }
421 else {
422 photonE = -1.;
423 }
424 } while (iteration <= maxDopplerIterations && photonE > eMax);
425
426 // End of recalculation of photon energy with Doppler broadening
427 // Revert to original if maximum number of iterations threshold
428 // has been reached
429 if (iteration >= maxDopplerIterations) {
430 photonE = photonEoriginal;
431 bindingE = 0.;
432 }
433
434 // Update G4VParticleChange for the scattered photon
435 G4ThreeVector photonDirection1(dirx, diry, dirz);
436 photonDirection1.rotateUz(photonDirection0);
437 fParticleChange->ProposeMomentumDirection(photonDirection1);
438
439 G4double photonEnergy1 = photonE;
440
441 if (photonEnergy1 > 0.) {
442 fParticleChange->SetProposedKineticEnergy(photonEnergy1);
443 }
444 else {
445 // photon absorbed
446 photonEnergy1 = 0.;
447 fParticleChange->SetProposedKineticEnergy(0.);
448 fParticleChange->ProposeTrackStatus(fStopAndKill);
449 fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0);
450 return;
451 }
452
453 // Kinematics of the scattered electron
454 G4double eKineticEnergy = photonEnergy0 - photonEnergy1 - bindingE;
455
456 // protection against negative final energy: no e- is created
457 if (eKineticEnergy < 0.0) {
458 fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0 - photonEnergy1);
459 return;
460 }
461
462 G4double eTotalEnergy = eKineticEnergy + electron_mass_c2;
463
464 G4double electronE = photonEnergy0 * (1. - epsilon) + electron_mass_c2;
465 G4double electronP2 = electronE * electronE - electron_mass_c2 * electron_mass_c2;
466 G4double sinThetaE = -1.;
467 G4double cosThetaE = 0.;
468 if (electronP2 > 0.) {
469 cosThetaE = (eTotalEnergy + photonEnergy1) * (1. - epsilon) / std::sqrt(electronP2);
470 sinThetaE = -1. * sqrt(1. - cosThetaE * cosThetaE);
471 }
472
473 G4double eDirX = sinThetaE * std::cos(phi);
474 G4double eDirY = sinThetaE * std::sin(phi);
475 G4double eDirZ = cosThetaE;
476
477 G4ThreeVector eDirection(eDirX, eDirY, eDirZ);
478 eDirection.rotateUz(photonDirection0);
479 auto dp = new G4DynamicParticle(G4Electron::Electron(), eDirection, eKineticEnergy);
480 fvect->push_back(dp);
481
482 // sample deexcitation
483 if (verboseLevel > 3) {
484 G4cout << "Started atomic de-excitation " << fAtomDeexcitation << G4endl;
485 }
486
487 if (nullptr != fAtomDeexcitation && iteration < maxDopplerIterations) {
488 G4int index = couple->GetIndex();
489 if (fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
490 size_t nbefore = fvect->size();
491 auto as = G4AtomicShellEnumerator(shellIdx);
492 const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
493 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
494 size_t nafter = fvect->size();
495 if (nafter > nbefore) {
496 for (size_t i = nbefore; i < nafter; ++i) {
497 // Check if there is enough residual energy
498 if (bindingE >= ((*fvect)[i])->GetKineticEnergy()) {
499 // Ok, this is a valid secondary: keep it
500 bindingE -= ((*fvect)[i])->GetKineticEnergy();
501 }
502 else {
503 // Invalid secondary: not enough energy to create it!
504 // Keep its energy in the local deposit
505 delete (*fvect)[i];
506 (*fvect)[i] = nullptr;
507 }
508 }
509 }
510 }
511 }
512 bindingE = std::max(bindingE, 0.0);
513 fParticleChange->ProposeLocalEnergyDeposit(bindingE);
514}
G4double epsilon(G4double density, G4double temperature)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:180
G4double G4Log(G4double x)
Definition G4Log.hh:227
@ fStopAndKill
const G4double alpha2
#define G4UniformRand()
Definition Randomize.hh:52
G4double RandomSelectMomentum(G4int Z, G4int shellIndex) const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition G4Electron.cc:91
G4int GetZasInt() const
Definition G4Element.hh:120
const G4Material * GetMaterial() const
const G4String & GetName() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
G4double BindingEnergy(G4int Z, G4int shellIndex) const
G4int SelectRandomShell(G4int Z) const
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)

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