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

#include <G4PenelopeIonisationModel.hh>

+ Inheritance diagram for G4PenelopeIonisationModel:

Public Member Functions

 G4PenelopeIonisationModel (const G4ParticleDefinition *p=0, const G4String &processName="PenIoni")
 
virtual ~G4PenelopeIonisationModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double, G4double, G4double, G4double, G4double)
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *theParticle, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
void SetVerbosityLevel (G4int lev)
 
G4int GetVerbosityLevel ()
 
- 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 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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
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
 
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 ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Attributes

G4ParticleChangeForLossfParticleChange
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

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 64 of file G4PenelopeIonisationModel.hh.

Constructor & Destructor Documentation

◆ G4PenelopeIonisationModel()

G4PenelopeIonisationModel::G4PenelopeIonisationModel ( const G4ParticleDefinition p = 0,
const G4String processName = "PenIoni" 
)

Definition at line 70 of file G4PenelopeIonisationModel.cc.

72 :G4VEmModel(nam),fParticleChange(0),isInitialised(false),
73 fAtomDeexcitation(0),kineticEnergy1(0.*eV),
74 cosThetaPrimary(1.0),energySecondary(0.*eV),
75 cosThetaSecondary(0.0),targetOscillator(-1),
76 theCrossSectionHandler(0)
77{
78 fIntrinsicLowEnergyLimit = 100.0*eV;
79 fIntrinsicHighEnergyLimit = 100.0*GeV;
80 // SetLowEnergyLimit(fIntrinsicLowEnergyLimit);
81 SetHighEnergyLimit(fIntrinsicHighEnergyLimit);
82 nBins = 200;
83 //
85 //
86 verboseLevel= 0;
87
88 // Verbosity scale:
89 // 0 = nothing
90 // 1 = warning for energy non-conservation
91 // 2 = details of energy budget
92 // 3 = calculation of cross sections, file openings, sampling of atoms
93 // 4 = entering in methods
94
95 // Atomic deexcitation model activated by default
97}
G4ParticleChangeForLoss * fParticleChange
static G4PenelopeOscillatorManager * GetOscillatorManager()
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:641

◆ ~G4PenelopeIonisationModel()

G4PenelopeIonisationModel::~G4PenelopeIonisationModel ( )
virtual

Definition at line 101 of file G4PenelopeIonisationModel.cc.

102{
103 if (theCrossSectionHandler)
104 delete theCrossSectionHandler;
105}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4PenelopeIonisationModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  ,
G4double  ,
G4double  ,
G4double  ,
G4double   
)
virtual

Reimplemented from G4VEmModel.

Definition at line 225 of file G4PenelopeIonisationModel.cc.

231{
232 G4cout << "*** G4PenelopeIonisationModel -- WARNING ***" << G4endl;
233 G4cout << "Penelope Ionisation model v2008 does not calculate cross section _per atom_ " << G4endl;
234 G4cout << "so the result is always zero. For physics values, please invoke " << G4endl;
235 G4cout << "GetCrossSectionPerVolume() or GetMeanFreePath() via the G4EmCalculator" << G4endl;
236 return 0;
237}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout

◆ ComputeDEDXPerVolume()

G4double G4PenelopeIonisationModel::ComputeDEDXPerVolume ( const G4Material material,
const G4ParticleDefinition theParticle,
G4double  kineticEnergy,
G4double  cutEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 241 of file G4PenelopeIonisationModel.cc.

245{
246 // Penelope model v2008 to calculate the stopping power for soft inelastic collisions
247 // below the threshold. It makes use of the Generalised Oscillator Strength (GOS)
248 // model from
249 // D. Liljequist, J. Phys. D: Appl. Phys. 16 (1983) 1567
250 //
251 // The stopping power is calculated analytically using the dsigma/dW cross
252 // section from the GOS models, which includes separate contributions from
253 // distant longitudinal collisions, distant transverse collisions and
254 // close collisions. Only the atomic oscillators that come in the play
255 // (according to the threshold) are considered for the calculation.
256 // Differential cross sections have a different form for e+ and e-.
257 //
258 // Fermi density correction is calculated analytically according to
259 // U. Fano, Ann. Rev. Nucl. Sci. 13 (1963),1
260
261 if (verboseLevel > 3)
262 G4cout << "Calling ComputeDEDX() of G4PenelopeIonisationModel" << G4endl;
263
264
265 G4PenelopeCrossSection* theXS =
266 theCrossSectionHandler->GetCrossSectionTableForCouple(theParticle,material,
267 cutEnergy);
268 G4double sPowerPerMolecule = 0.0;
269 if (theXS)
270 sPowerPerMolecule = theXS->GetSoftStoppingPower(kineticEnergy);
271
272
273 G4double atomDensity = material->GetTotNbOfAtomsPerVolume();
274 G4double atPerMol = oscManager->GetAtomsPerMolecule(material);
275
276 G4double moleculeDensity = 0.;
277 if (atPerMol)
278 moleculeDensity = atomDensity/atPerMol;
279
280 G4double sPowerPerVolume = sPowerPerMolecule*moleculeDensity;
281
282 if (verboseLevel > 2)
283 {
284 G4cout << "G4PenelopeIonisationModel " << G4endl;
285 G4cout << "Stopping power < " << cutEnergy/keV << " keV at " <<
286 kineticEnergy/keV << " keV = " <<
287 sPowerPerVolume/(keV/mm) << " keV/mm" << G4endl;
288 }
289 return sPowerPerVolume;
290}
double G4double
Definition: G4Types.hh:64
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:208
G4double GetSoftStoppingPower(G4double energy)
Returns the total stopping power due to soft collisions.
G4PenelopeCrossSection * GetCrossSectionTableForCouple(const G4ParticleDefinition *, const G4Material *, G4double cut)
G4double GetAtomsPerMolecule(const G4Material *)
Returns the total number of atoms per molecule.

◆ CrossSectionPerVolume()

G4double G4PenelopeIonisationModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition theParticle,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 155 of file G4PenelopeIonisationModel.cc.

161{
162 // Penelope model v2008 to calculate the cross section for inelastic collisions above the
163 // threshold. It makes use of the Generalised Oscillator Strength (GOS) model from
164 // D. Liljequist, J. Phys. D: Appl. Phys. 16 (1983) 1567
165 //
166 // The total cross section is calculated analytically by taking
167 // into account the atomic oscillators coming into the play for a given threshold.
168 //
169 // For incident e- the maximum energy allowed for the delta-rays is energy/2.
170 // because particles are undistinghishable.
171 //
172 // The contribution is splitted in three parts: distant longitudinal collisions,
173 // distant transverse collisions and close collisions. Each term is described by
174 // its own analytical function.
175 // Fermi density correction is calculated analytically according to
176 // U. Fano, Ann. Rev. Nucl. Sci. 13 (1963),1
177 //
178 if (verboseLevel > 3)
179 G4cout << "Calling CrossSectionPerVolume() of G4PenelopeIonisationModel" << G4endl;
180
181 SetupForMaterial(theParticle, material, energy);
182
183 G4double totalCross = 0.0;
184 G4double crossPerMolecule = 0.;
185
186 G4PenelopeCrossSection* theXS =
187 theCrossSectionHandler->GetCrossSectionTableForCouple(theParticle,
188 material,
189 cutEnergy);
190
191 if (theXS)
192 crossPerMolecule = theXS->GetHardCrossSection(energy);
193
194 G4double atomDensity = material->GetTotNbOfAtomsPerVolume();
195 G4double atPerMol = oscManager->GetAtomsPerMolecule(material);
196
197 if (verboseLevel > 3)
198 G4cout << "Material " << material->GetName() << " has " << atPerMol <<
199 "atoms per molecule" << G4endl;
200
201 G4double moleculeDensity = 0.;
202 if (atPerMol)
203 moleculeDensity = atomDensity/atPerMol;
204
205 G4double crossPerVolume = crossPerMolecule*moleculeDensity;
206
207 if (verboseLevel > 2)
208 {
209 G4cout << "G4PenelopeIonisationModel " << G4endl;
210 G4cout << "Mean free path for delta emission > " << cutEnergy/keV << " keV at " <<
211 energy/keV << " keV = " << (1./crossPerVolume)/mm << " mm" << G4endl;
212 if (theXS)
213 totalCross = (theXS->GetTotalCrossSection(energy))*moleculeDensity;
214 G4cout << "Total free path for ionisation (no threshold) at " <<
215 energy/keV << " keV = " << (1./totalCross)/mm << " mm" << G4endl;
216 }
217 return crossPerVolume;
218}
const G4String & GetName() const
Definition: G4Material.hh:177
G4double GetHardCrossSection(G4double energy)
Returns hard cross section at the given energy.
G4double GetTotalCrossSection(G4double energy)
Returns total cross section at the given energy.
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:311

◆ GetVerbosityLevel()

G4int G4PenelopeIonisationModel::GetVerbosityLevel ( )
inline

Definition at line 108 of file G4PenelopeIonisationModel.hh.

108{return verboseLevel;};

◆ Initialise()

void G4PenelopeIonisationModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 109 of file G4PenelopeIonisationModel.cc.

111{
112 if (verboseLevel > 3)
113 G4cout << "Calling G4PenelopeIonisationModel::Initialise()" << G4endl;
114
115 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
116 //Issue warning if the AtomicDeexcitation has not been declared
117 if (!fAtomDeexcitation)
118 {
119 G4cout << G4endl;
120 G4cout << "WARNING from G4PenelopeIonisationModel " << G4endl;
121 G4cout << "Atomic de-excitation module is not instantiated, so there will not be ";
122 G4cout << "any fluorescence/Auger emission." << G4endl;
123 G4cout << "Please make sure this is intended" << G4endl;
124 }
125
126 //Set the number of bins for the tables. 20 points per decade
127 nBins = (size_t) (20*std::log10(HighEnergyLimit()/LowEnergyLimit()));
128 nBins = std::max(nBins,(size_t)100);
129
130 //Clear and re-build the tables
131 if (theCrossSectionHandler)
132 {
133 delete theCrossSectionHandler;
134 theCrossSectionHandler = 0;
135 }
136 theCrossSectionHandler = new G4PenelopeIonisationXSHandler(nBins);
137 theCrossSectionHandler->SetVerboseLevel(verboseLevel);
138
139 if (verboseLevel > 2) {
140 G4cout << "Penelope Ionisation model v2008 is initialized " << G4endl
141 << "Energy range: "
142 << LowEnergyLimit() / keV << " keV - "
143 << HighEnergyLimit() / GeV << " GeV. Using "
144 << nBins << " bins."
145 << G4endl;
146 }
147
148 if(isInitialised) return;
150 isInitialised = true;
151}
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
void SetVerboseLevel(G4int vl)
Setter for the verbosity level.
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:95

◆ MinEnergyCut()

G4double G4PenelopeIonisationModel::MinEnergyCut ( const G4ParticleDefinition ,
const G4MaterialCutsCouple  
)
virtual

Definition at line 294 of file G4PenelopeIonisationModel.cc.

296{
297 return fIntrinsicLowEnergyLimit;
298}

◆ SampleSecondaries()

void G4PenelopeIonisationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicParticle,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 302 of file G4PenelopeIonisationModel.cc.

306{
307 // Penelope model v2008 to sample the final state following an hard inelastic interaction.
308 // It makes use of the Generalised Oscillator Strength (GOS) model from
309 // D. Liljequist, J. Phys. D: Appl. Phys. 16 (1983) 1567
310 //
311 // The GOS model is used to calculate the individual cross sections for all
312 // the atomic oscillators coming in the play, taking into account the three
313 // contributions (distant longitudinal collisions, distant transverse collisions and
314 // close collisions). Then the target shell and the interaction channel are
315 // sampled. Final state of the delta-ray (energy, angle) are generated according
316 // to the analytical distributions (dSigma/dW) for the selected interaction
317 // channels.
318 // For e-, the maximum energy for the delta-ray is initialEnergy/2. (because
319 // particles are indistinghusbable), while it is the full initialEnergy for
320 // e+.
321 // The efficiency of the random sampling algorithm (e.g. for close collisions)
322 // decreases when initial and cutoff energy increase (e.g. 87% for 10-keV primary
323 // and 1 keV threshold, 99% for 10-MeV primary and 10-keV threshold).
324 // Differential cross sections have a different form for e+ and e-.
325 //
326 // WARNING: The model provides an _average_ description of inelastic collisions.
327 // Anyway, the energy spectrum associated to distant excitations of a given
328 // atomic shell is approximated as a single resonance. The simulated energy spectra
329 // show _unphysical_ narrow peaks at energies that are multiple of the shell
330 // resonance energies. The spurious speaks are automatically smoothed out after
331 // multiple inelastic collisions.
332 //
333 // The model determines also the original shell from which the delta-ray is expelled,
334 // in order to produce fluorescence de-excitation (from G4DeexcitationManager)
335 //
336 // Fermi density correction is calculated analytically according to
337 // U. Fano, Ann. Rev. Nucl. Sci. 13 (1963),1
338
339 if (verboseLevel > 3)
340 G4cout << "Calling SamplingSecondaries() of G4PenelopeIonisationModel" << G4endl;
341
342 G4double kineticEnergy0 = aDynamicParticle->GetKineticEnergy();
343 const G4ParticleDefinition* theParticle = aDynamicParticle->GetDefinition();
344
345 if (kineticEnergy0 <= fIntrinsicLowEnergyLimit)
346 {
349 return ;
350 }
351
352 const G4Material* material = couple->GetMaterial();
353 G4PenelopeOscillatorTable* theTable = oscManager->GetOscillatorTableIonisation(material);
354
355 G4ParticleMomentum particleDirection0 = aDynamicParticle->GetMomentumDirection();
356
357 //Initialise final-state variables. The proper values will be set by the methods
358 // SampleFinalStateElectron() and SampleFinalStatePositron()
359 kineticEnergy1=kineticEnergy0;
360 cosThetaPrimary=1.0;
361 energySecondary=0.0;
362 cosThetaSecondary=1.0;
363 targetOscillator = -1;
364
365 if (theParticle == G4Electron::Electron())
366 SampleFinalStateElectron(material,cutE,kineticEnergy0);
367 else if (theParticle == G4Positron::Positron())
368 SampleFinalStatePositron(material,cutE,kineticEnergy0);
369 else
370 {
372 ed << "Invalid particle " << theParticle->GetParticleName() << G4endl;
373 G4Exception("G4PenelopeIonisationModel::SamplingSecondaries()",
374 "em0001",FatalException,ed);
375
376 }
377 if (energySecondary == 0) return;
378
379 if (verboseLevel > 3)
380 {
381 G4cout << "G4PenelopeIonisationModel::SamplingSecondaries() for " <<
382 theParticle->GetParticleName() << G4endl;
383 G4cout << "Final eKin = " << kineticEnergy1 << " keV" << G4endl;
384 G4cout << "Final cosTheta = " << cosThetaPrimary << G4endl;
385 G4cout << "Delta-ray eKin = " << energySecondary << " keV" << G4endl;
386 G4cout << "Delta-ray cosTheta = " << cosThetaSecondary << G4endl;
387 G4cout << "Oscillator: " << targetOscillator << G4endl;
388 }
389
390 //Update the primary particle
391 G4double sint = std::sqrt(1. - cosThetaPrimary*cosThetaPrimary);
392 G4double phiPrimary = twopi * G4UniformRand();
393 G4double dirx = sint * std::cos(phiPrimary);
394 G4double diry = sint * std::sin(phiPrimary);
395 G4double dirz = cosThetaPrimary;
396
397 G4ThreeVector electronDirection1(dirx,diry,dirz);
398 electronDirection1.rotateUz(particleDirection0);
399
400 if (kineticEnergy1 > 0)
401 {
402 fParticleChange->ProposeMomentumDirection(electronDirection1);
404 }
405 else
407
408
409 //Generate the delta ray
410 G4double ionEnergyInPenelopeDatabase =
411 (*theTable)[targetOscillator]->GetIonisationEnergy();
412
413 //Now, try to handle fluorescence
414 //Notice: merged levels are indicated with Z=0 and flag=30
415 G4int shFlag = (*theTable)[targetOscillator]->GetShellFlag();
416 G4int Z = (G4int) (*theTable)[targetOscillator]->GetParentZ();
417
418 //initialize here, then check photons created by Atomic-Deexcitation, and the final state e-
420 G4double bindingEnergy = 0.*eV;
421 //G4int shellId = 0;
422
423 const G4AtomicShell* shell = 0;
424 //Real level
425 if (Z > 0 && shFlag<30)
426 {
427 shell = transitionManager->Shell(Z,shFlag-1);
428 bindingEnergy = shell->BindingEnergy();
429 //shellId = shell->ShellId();
430 }
431
432 //correct the energySecondary to account for the fact that the Penelope
433 //database of ionisation energies is in general (slightly) different
434 //from the fluorescence database used in Geant4.
435 energySecondary += ionEnergyInPenelopeDatabase-bindingEnergy;
436
437 G4double localEnergyDeposit = bindingEnergy;
438 //testing purposes only
439 G4double energyInFluorescence = 0;
440 G4double energyInAuger = 0;
441
442 if (energySecondary < 0)
443 {
444 //It means that there was some problem/mismatch between the two databases.
445 //In this case, the available energy is ok to excite the level according
446 //to the Penelope database, but not according to the Geant4 database
447 //Full residual energy is deposited locally
448 localEnergyDeposit += energySecondary;
449 energySecondary = 0.0;
450 }
451
452 //Notice: shell might be NULL (invalid!) if shFlag=30. Must be protected
453 //Now, take care of fluorescence, if required
454 if (fAtomDeexcitation && shell)
455 {
456 G4int index = couple->GetIndex();
457 if (fAtomDeexcitation->CheckDeexcitationActiveRegion(index))
458 {
459 size_t nBefore = fvect->size();
460 fAtomDeexcitation->GenerateParticles(fvect,shell,Z,index);
461 size_t nAfter = fvect->size();
462
463 if (nAfter > nBefore) //actual production of fluorescence
464 {
465 for (size_t j=nBefore;j<nAfter;j++) //loop on products
466 {
467 G4double itsEnergy = ((*fvect)[j])->GetKineticEnergy();
468 localEnergyDeposit -= itsEnergy;
469 if (((*fvect)[j])->GetParticleDefinition() == G4Gamma::Definition())
470 energyInFluorescence += itsEnergy;
471 else if (((*fvect)[j])->GetParticleDefinition() == G4Electron::Definition())
472 energyInAuger += itsEnergy;
473 }
474 }
475 }
476 }
477
478 // Generate the delta ray --> to be done only if above cut
479 if (energySecondary > cutE)
480 {
481 G4DynamicParticle* electron = 0;
482 G4double sinThetaE = std::sqrt(1.-cosThetaSecondary*cosThetaSecondary);
483 G4double phiEl = phiPrimary+pi; //pi with respect to the primary electron/positron
484 G4double xEl = sinThetaE * std::cos(phiEl);
485 G4double yEl = sinThetaE * std::sin(phiEl);
486 G4double zEl = cosThetaSecondary;
487 G4ThreeVector eDirection(xEl,yEl,zEl); //electron direction
488 eDirection.rotateUz(particleDirection0);
489 electron = new G4DynamicParticle (G4Electron::Electron(),
490 eDirection,energySecondary) ;
491 fvect->push_back(electron);
492 }
493 else
494 {
495 localEnergyDeposit += energySecondary;
496 energySecondary = 0;
497 }
498
499 if (localEnergyDeposit < 0)
500 {
501 G4cout << "WARNING-"
502 << "G4PenelopeIonisationModel::SampleSecondaries - Negative energy deposit"
503 << G4endl;
504 localEnergyDeposit=0.;
505 }
506 fParticleChange->ProposeLocalEnergyDeposit(localEnergyDeposit);
507
508 if (verboseLevel > 1)
509 {
510 G4cout << "-----------------------------------------------------------" << G4endl;
511 G4cout << "Energy balance from G4PenelopeIonisation" << G4endl;
512 G4cout << "Incoming primary energy: " << kineticEnergy0/keV << " keV" << G4endl;
513 G4cout << "-----------------------------------------------------------" << G4endl;
514 G4cout << "Outgoing primary energy: " << kineticEnergy1/keV << " keV" << G4endl;
515 G4cout << "Delta ray " << energySecondary/keV << " keV" << G4endl;
516 if (energyInFluorescence)
517 G4cout << "Fluorescence x-rays: " << energyInFluorescence/keV << " keV" << G4endl;
518 if (energyInAuger)
519 G4cout << "Auger electrons: " << energyInAuger/keV << " keV" << G4endl;
520 G4cout << "Local energy deposit " << localEnergyDeposit/keV << " keV" << G4endl;
521 G4cout << "Total final state: " << (energySecondary+energyInFluorescence+kineticEnergy1+
522 localEnergyDeposit+energyInAuger)/keV <<
523 " keV" << G4endl;
524 G4cout << "-----------------------------------------------------------" << G4endl;
525 }
526
527 if (verboseLevel > 0)
528 {
529 G4double energyDiff = std::fabs(energySecondary+energyInFluorescence+kineticEnergy1+
530 localEnergyDeposit+energyInAuger-kineticEnergy0);
531 if (energyDiff > 0.05*keV)
532 G4cout << "Warning from G4PenelopeIonisation: problem with energy conservation: " <<
533 (energySecondary+energyInFluorescence+kineticEnergy1+localEnergyDeposit+energyInAuger)/keV <<
534 " keV (final) vs. " <<
535 kineticEnergy0/keV << " keV (initial)" << G4endl;
536 }
537
538}
@ FatalException
std::vector< G4PenelopeOscillator * > G4PenelopeOscillatorTable
int G4int
Definition: G4Types.hh:66
#define G4UniformRand()
Definition: Randomize.hh:53
G4double BindingEnergy() const
G4AtomicShell * Shell(G4int Z, size_t shellIndex) const
static G4AtomicTransitionManager * Instance()
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Definition()
Definition: G4Electron.cc:49
static G4Electron * Electron()
Definition: G4Electron.cc:94
static G4Gamma * Definition()
Definition: G4Gamma.cc:49
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4String & GetParticleName() const
G4PenelopeOscillatorTable * GetOscillatorTableIonisation(const G4Material *)
static G4Positron * Positron()
Definition: G4Positron.cc:94
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
const G4double pi
G4double bindingEnergy(G4int A, G4int Z)

◆ SetVerbosityLevel()

void G4PenelopeIonisationModel::SetVerbosityLevel ( G4int  lev)
inline

Definition at line 107 of file G4PenelopeIonisationModel.hh.

107{verboseLevel = lev;};

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForLoss* G4PenelopeIonisationModel::fParticleChange
protected

Definition at line 111 of file G4PenelopeIonisationModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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