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

#include <G4UAtomicDeexcitation.hh>

+ Inheritance diagram for G4UAtomicDeexcitation:

Public Member Functions

 G4UAtomicDeexcitation ()
 
virtual ~G4UAtomicDeexcitation ()
 
void InitialiseForNewRun () override
 initialisation methods
 
void InitialiseForExtraAtom (G4int Z) override
 
void SetCutForSecondaryPhotons (G4double cut)
 Set threshold energy for fluorescence.
 
void SetCutForAugerElectrons (G4double cut)
 Set threshold energy for Auger electron production.
 
const G4AtomicShellGetAtomicShell (G4int Z, G4AtomicShellEnumerator shell) override
 
void GenerateParticles (std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4double gammaCut, G4double eCut) override
 generation of deexcitation for given atom, shell vacancy and cuts
 
G4double GetShellIonisationCrossSectionPerAtom (const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=nullptr) override
 access or compute PIXE cross section
 
G4double ComputeShellIonisationCrossSectionPerAtom (const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=nullptr) override
 access or compute PIXE cross section
 
 G4UAtomicDeexcitation (G4UAtomicDeexcitation &)=delete
 
G4UAtomicDeexcitationoperator= (const G4UAtomicDeexcitation &right)=delete
 
- Public Member Functions inherited from G4VAtomDeexcitation
 G4VAtomDeexcitation (const G4String &modname="Deexcitation")
 
virtual ~G4VAtomDeexcitation ()
 
void InitialiseAtomicDeexcitation ()
 
void SetDeexcitationActiveRegion (const G4String &rname, G4bool valDeexcitation, G4bool valAuger, G4bool valPIXE)
 
void SetFluo (G4bool)
 
G4bool IsFluoActive () const
 
void SetAuger (G4bool)
 
G4bool IsAugerActive () const
 
void SetAugerCascade (G4bool)
 
G4bool IsAugerCascadeActive () const
 
void SetPIXE (G4bool)
 
G4bool IsPIXEActive () const
 
const G4StringGetName () const
 
const std::vector< G4bool > & GetListOfActiveAtoms () const
 
void SetVerboseLevel (G4int)
 
G4int GetVerboseLevel () const
 
G4bool CheckDeexcitationActiveRegion (G4int coupleIndex)
 
G4bool CheckAugerActiveRegion (G4int coupleIndex)
 
void GenerateParticles (std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
 
void AlongStepDeexcitation (std::vector< G4Track * > &tracks, const G4Step &step, G4double &eLoss, G4int coupleIndex)
 
 G4VAtomDeexcitation (G4VAtomDeexcitation &)=delete
 
G4VAtomDeexcitationoperator= (const G4VAtomDeexcitation &right)=delete
 

Detailed Description

Definition at line 60 of file G4UAtomicDeexcitation.hh.

Constructor & Destructor Documentation

◆ G4UAtomicDeexcitation() [1/2]

G4UAtomicDeexcitation::G4UAtomicDeexcitation ( )
explicit

Definition at line 77 of file G4UAtomicDeexcitation.cc.

77 :
78 G4VAtomDeexcitation("UAtomDeexcitation"),
79 minGammaEnergy(DBL_MAX),
80 minElectronEnergy(DBL_MAX),
81 newShellId(-1)
82{
83 anaPIXEshellCS = nullptr;
84 PIXEshellCS = nullptr;
85 ePIXEshellCS = nullptr;
87 theElectron = G4Electron::Electron();
88 thePositron = G4Positron::Positron();
89 transitionManager = G4AtomicTransitionManager::Instance();
90}
static G4AtomicTransitionManager * Instance()
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
static G4Positron * Positron()
Definition G4Positron.cc:90
G4VAtomDeexcitation(const G4String &modname="Deexcitation")
#define DBL_MAX
Definition templates.hh:62

◆ ~G4UAtomicDeexcitation()

G4UAtomicDeexcitation::~G4UAtomicDeexcitation ( )
virtual

Definition at line 94 of file G4UAtomicDeexcitation.cc.

95{
96 delete anaPIXEshellCS;
97 delete PIXEshellCS;
98 delete ePIXEshellCS;
99}

◆ G4UAtomicDeexcitation() [2/2]

G4UAtomicDeexcitation::G4UAtomicDeexcitation ( G4UAtomicDeexcitation & )
delete

Member Function Documentation

◆ ComputeShellIonisationCrossSectionPerAtom()

G4double G4UAtomicDeexcitation::ComputeShellIonisationCrossSectionPerAtom ( const G4ParticleDefinition * p,
G4int Z,
G4AtomicShellEnumerator shell,
G4double kinE,
const G4Material * mat = nullptr )
overridevirtual

access or compute PIXE cross section

Implements G4VAtomDeexcitation.

Definition at line 369 of file G4UAtomicDeexcitation.cc.

375{
376 return GetShellIonisationCrossSectionPerAtom(p,Z,shell,kinE,mat);
377}
G4double GetShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=nullptr) override
access or compute PIXE cross section

◆ GenerateParticles()

void G4UAtomicDeexcitation::GenerateParticles ( std::vector< G4DynamicParticle * > * secVect,
const G4AtomicShell * atomicShell,
G4int Z,
G4double gammaCut,
G4double eCut )
overridevirtual

generation of deexcitation for given atom, shell vacancy and cuts

Implements G4VAtomDeexcitation.

Definition at line 186 of file G4UAtomicDeexcitation.cc.

192{
193 // Defined initial conditions
194 G4int givenShellId = atomicShell->ShellId();
195 minGammaEnergy = gammaCut;
196 minElectronEnergy = eCut;
197
198 // generation secondaries
199 G4DynamicParticle* aParticle=0;
200 G4int provShellId = 0;
201
202 //ORIGINAL METHOD BY ALFONSO MANTERO
204 {
205 //----------------------------
206 G4int counter = 0;
207
208 // limits of the EPDL data
209 if (Z>5 && Z<105) {
210
211 // The aim of this loop is to generate more than one fluorecence photon
212 // from the same ionizing event
213 do
214 {
215 if (counter == 0)
216 // First call to GenerateParticles(...):
217 // givenShellId is given by the process
218 {
219 provShellId = SelectTypeOfTransition(Z, givenShellId);
220
221 if (provShellId >0)
222 {
223 aParticle =
224 GenerateFluorescence(Z, givenShellId, provShellId);
225 }
226 else if (provShellId == -1)
227 {
228 aParticle = GenerateAuger(Z, givenShellId);
229 }
230 }
231 else
232 // Following calls to GenerateParticles(...):
233 // newShellId is given by GenerateFluorescence(...)
234 {
235 provShellId = SelectTypeOfTransition(Z,newShellId);
236 if (provShellId >0)
237 {
238 aParticle = GenerateFluorescence(Z,newShellId,provShellId);
239 }
240 else if ( provShellId == -1)
241 {
242 aParticle = GenerateAuger(Z, newShellId);
243 }
244 }
245 ++counter;
246 if (aParticle != 0)
247 {
248 vectorOfParticles->push_back(aParticle);
249 }
250 else {provShellId = -2;}
251 }
252 while (provShellId > -2);
253 }
254 } // Auger cascade is not active
255
256 //END OF ORIGINAL METHOD BY ALFONSO MANTERO
257 //----------------------
258
259 // NEW METHOD
260 // Auger cascade by Burkhant Suerfu on March 24 2015 (Bugzilla 1727)
262 {
263 //----------------------
264 vacancyArray.push_back(givenShellId);
265
266 // let's check that 5<Z<100
267 if (Z<6 || Z>104){
268 return;
269 }
270
271 // as long as there is vacancy to be filled by either fluo or auger, stay in the loop.
272 while(!vacancyArray.empty()){
273 // prepare to process the last element, and then delete it from the vector.
274 givenShellId = vacancyArray[0];
275 provShellId = SelectTypeOfTransition(Z,givenShellId);
276
277 //G4cout<<"\n------ Atom Transition with Z: "<<Z<<"\tbetween current:"
278 // <<givenShellId<<" & target:"<<provShellId<<G4endl;
279 if(provShellId>0){
280 aParticle = GenerateFluorescence(Z,givenShellId,provShellId);
281 }
282 else if(provShellId == -1){
283 aParticle = GenerateAuger(Z, givenShellId);
284 }
285 // if a particle is created, put it in the vector of new particles
286 if(aParticle!=0)
287 vectorOfParticles->push_back(aParticle);
288
289 // one vacancy has been processed. Erase it.
290 vacancyArray.erase(vacancyArray.begin());
291 }
292 //----------------------
293 //End of Auger cascade by Burkhant Suerfu on March 24 2015 (Bugzilla 1727)
294
295 } // Auger cascade is active
296}
int G4int
Definition G4Types.hh:85
G4int ShellId() const
G4bool IsAugerCascadeActive() const

◆ GetAtomicShell()

const G4AtomicShell * G4UAtomicDeexcitation::GetAtomicShell ( G4int Z,
G4AtomicShellEnumerator shell )
overridevirtual

Get atomic shell by shell index, used by discrete processes (for example, photoelectric), when shell vacancy sampled by the model

Implements G4VAtomDeexcitation.

Definition at line 179 of file G4UAtomicDeexcitation.cc.

180{
181 return transitionManager->Shell(Z, (std::size_t)shell);
182}
G4AtomicShell * Shell(G4int Z, size_t shellIndex) const

◆ GetShellIonisationCrossSectionPerAtom()

G4double G4UAtomicDeexcitation::GetShellIonisationCrossSectionPerAtom ( const G4ParticleDefinition * pdef,
G4int Z,
G4AtomicShellEnumerator shell,
G4double kinE,
const G4Material * mat = nullptr )
overridevirtual

access or compute PIXE cross section

Implements G4VAtomDeexcitation.

Definition at line 301 of file G4UAtomicDeexcitation.cc.

307{
308 // we must put a control on the shell that are passed:
309 // some shells should not pass (line "0" or "2")
310
311 // check atomic number
312 G4double xsec = 0.0;
313 if(Z > 93 || Z < 6 ) { return xsec; } //corrected by alf - Z<6 missing
314 G4int idx = G4int(shellEnum);
315 if(idx >= G4AtomicShells::GetNumberOfShells(Z)) { return xsec; }
316
317 if(pdef == theElectron || pdef == thePositron) {
318 xsec = ePIXEshellCS->CrossSection(Z,shellEnum,kineticEnergy,0.0,mat);
319 return xsec;
320 }
321
322 G4double mass = pdef->GetPDGMass();
323 G4double escaled = kineticEnergy;
324 G4double q2 = 0.0;
325
326 // scaling to protons for all particles excluding protons and alpha
327 G4int pdg = pdef->GetPDGEncoding();
328 if (pdg != 2212 && pdg != 1000020040)
329 {
330 mass = proton_mass_c2;
331 escaled = kineticEnergy*mass/(pdef->GetPDGMass());
332
333 if(mat) {
334 q2 = emcorr->EffectiveChargeSquareRatio(pdef,mat,kineticEnergy);
335 } else {
336 G4double q = pdef->GetPDGCharge()/eplus;
337 q2 = q*q;
338 }
339 }
340
341 if(PIXEshellCS) {
342 xsec = PIXEshellCS->CrossSection(Z,shellEnum,escaled,mass,mat);
343 }
344 if(xsec < 1e-100) {
345 xsec = anaPIXEshellCS->CrossSection(Z,shellEnum,escaled,mass,mat);
346 }
347
348 if (q2) {xsec *= q2;}
349
350 return xsec;
351}
double G4double
Definition G4Types.hh:83
static G4int GetNumberOfShells(G4int Z)
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, const G4double kineticEnergy)
virtual G4double CrossSection(G4int Z, G4AtomicShellEnumerator shell, G4double incidentEnergy, G4double mass, const G4Material *mat)=0

Referenced by ComputeShellIonisationCrossSectionPerAtom().

◆ InitialiseForExtraAtom()

void G4UAtomicDeexcitation::InitialiseForExtraAtom ( G4int Z)
overridevirtual

Implements G4VAtomDeexcitation.

Definition at line 173 of file G4UAtomicDeexcitation.cc.

174{}

◆ InitialiseForNewRun()

void G4UAtomicDeexcitation::InitialiseForNewRun ( )
overridevirtual

initialisation methods

Implements G4VAtomDeexcitation.

Definition at line 103 of file G4UAtomicDeexcitation.cc.

104{
105 if(!IsFluoActive()) { return; }
106 transitionManager->Initialise();
107 if(!IsPIXEActive()) { return; }
108
109 if(!anaPIXEshellCS) {
110 anaPIXEshellCS = new G4teoCrossSection("ECPSSR_Analytical");
111 }
112 G4cout << G4endl;
113 G4cout << "### === G4UAtomicDeexcitation::InitialiseForNewRun()" << G4endl;
114
116 G4String namePIXExsModel = param->PIXECrossSectionModel();
117 G4String namePIXExsElectronModel = param->PIXEElectronCrossSectionModel();
118
119 // Check if old cross section for p/ion should be deleted
120 if(PIXEshellCS && namePIXExsModel != PIXEshellCS->GetName())
121 {
122 delete PIXEshellCS;
123 PIXEshellCS = nullptr;
124 }
125
126 // Instantiate new proton/ion cross section
127 if(!PIXEshellCS) {
128 if (namePIXExsModel == "ECPSSR_FormFactor")
129 {
130 PIXEshellCS = new G4teoCrossSection(namePIXExsModel);
131 }
132 else if(namePIXExsModel == "ECPSSR_ANSTO")
133 {
134 PIXEshellCS = new G4teoCrossSection(namePIXExsModel);
135 }
136 else if(namePIXExsModel == "Empirical")
137 {
138 PIXEshellCS = new G4empCrossSection(namePIXExsModel);
139 }
140 }
141
142 // Check if old cross section for e+- should be deleted
143 if(ePIXEshellCS && namePIXExsElectronModel != ePIXEshellCS->GetName())
144 {
145 delete ePIXEshellCS;
146 ePIXEshellCS = nullptr;
147 }
148
149 // Instantiate new e+- cross section
150 if(nullptr == ePIXEshellCS)
151 {
152 if(namePIXExsElectronModel == "Empirical")
153 {
154 ePIXEshellCS = new G4empCrossSection("Empirical");
155 }
156 else if(namePIXExsElectronModel == "ECPSSR_Analytical")
157 {
158 ePIXEshellCS = new G4teoCrossSection("ECPSSR_Analytical");
159 }
160 else if (namePIXExsElectronModel == "Penelope")
161 {
162 ePIXEshellCS = new G4PenelopeIonisationCrossSection();
163 }
164 else
165 {
166 ePIXEshellCS = new G4LivermoreIonisationCrossSection();
167 }
168 }
169}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void Initialise()
needs to be called once from other code before start of run
static G4EmParameters * Instance()
const G4String & PIXECrossSectionModel()
const G4String & PIXEElectronCrossSectionModel()
const G4String & GetName() const

◆ operator=()

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

◆ SetCutForAugerElectrons()

void G4UAtomicDeexcitation::SetCutForAugerElectrons ( G4double cut)

Set threshold energy for Auger electron production.

Definition at line 362 of file G4UAtomicDeexcitation.cc.

363{
364 minElectronEnergy = cut;
365}

◆ SetCutForSecondaryPhotons()

void G4UAtomicDeexcitation::SetCutForSecondaryPhotons ( G4double cut)

Set threshold energy for fluorescence.

Definition at line 355 of file G4UAtomicDeexcitation.cc.

356{
357 minGammaEnergy = cut;
358}

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