Geant4 9.6.0
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 ()
 
virtual void InitialiseForNewRun ()
 
virtual void InitialiseForExtraAtom (G4int Z)
 
void SetCutForSecondaryPhotons (G4double cut)
 
void SetCutForAugerElectrons (G4double cut)
 
virtual const G4AtomicShellGetAtomicShell (G4int Z, G4AtomicShellEnumerator shell)
 
virtual void GenerateParticles (std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4double gammaCut, G4double eCut)
 
virtual G4double GetShellIonisationCrossSectionPerAtom (const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0)
 
virtual G4double ComputeShellIonisationCrossSectionPerAtom (const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0)
 
- Public Member Functions inherited from G4VAtomDeexcitation
 G4VAtomDeexcitation (const G4String &modname="Deexcitation", const G4String &pixename="")
 
virtual ~G4VAtomDeexcitation ()
 
void InitialiseAtomicDeexcitation ()
 
virtual void InitialiseForNewRun ()=0
 
virtual void InitialiseForExtraAtom (G4int Z)=0
 
void SetDeexcitationActiveRegion (const G4String &rname, G4bool valDeexcitation, G4bool valAuger, G4bool valPIXE)
 
void SetFluo (G4bool)
 
G4bool IsFluoActive () const
 
void SetAuger (G4bool)
 
G4bool IsAugerActive () const
 
void SetPIXE (G4bool)
 
G4bool IsPIXEActive () const
 
const G4StringGetName () const
 
void SetPIXECrossSectionModel (const G4String &)
 
const G4StringPIXECrossSectionModel () const
 
void SetPIXEElectronCrossSectionModel (const G4String &)
 
const G4StringPIXEElectronCrossSectionModel () const
 
const std::vector< G4bool > & GetListOfActiveAtoms () const
 
void SetVerboseLevel (G4int)
 
G4int GetVerboseLevel () const
 
G4bool CheckDeexcitationActiveRegion (G4int coupleIndex)
 
G4bool CheckAugerActiveRegion (G4int coupleIndex)
 
virtual const G4AtomicShellGetAtomicShell (G4int Z, G4AtomicShellEnumerator shell)=0
 
void GenerateParticles (std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
 
virtual void GenerateParticles (std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4double gammaCut, G4double eCut)=0
 
virtual G4double GetShellIonisationCrossSectionPerAtom (const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0)=0
 
virtual G4double ComputeShellIonisationCrossSectionPerAtom (const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0)=0
 
void AlongStepDeexcitation (std::vector< G4Track * > &tracks, const G4Step &step, G4double &eLoss, G4int coupleIndex)
 

Detailed Description

Definition at line 62 of file G4UAtomicDeexcitation.hh.

Constructor & Destructor Documentation

◆ G4UAtomicDeexcitation()

G4UAtomicDeexcitation::G4UAtomicDeexcitation ( )

Definition at line 76 of file G4UAtomicDeexcitation.cc.

76 :
77 G4VAtomDeexcitation("UAtomDeexcitation"),
78 minGammaEnergy(DBL_MAX),
79 minElectronEnergy(DBL_MAX),
80 emcorr(0)
81{
82 PIXEshellCS = 0;
83 ePIXEshellCS = 0;
85 theElectron = G4Electron::Electron();
86 thePositron = G4Positron::Positron();
87 transitionManager = 0;
88 anaPIXEshellCS = 0;
89}
static G4Electron * Electron()
Definition: G4Electron.cc:94
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
static G4Positron * Positron()
Definition: G4Positron.cc:94
#define DBL_MAX
Definition: templates.hh:83

◆ ~G4UAtomicDeexcitation()

G4UAtomicDeexcitation::~G4UAtomicDeexcitation ( )
virtual

Definition at line 91 of file G4UAtomicDeexcitation.cc.

92{
93 delete PIXEshellCS;
94 delete anaPIXEshellCS;
95 delete ePIXEshellCS;
96}

Member Function Documentation

◆ ComputeShellIonisationCrossSectionPerAtom()

G4double G4UAtomicDeexcitation::ComputeShellIonisationCrossSectionPerAtom ( const G4ParticleDefinition p,
G4int  Z,
G4AtomicShellEnumerator  shell,
G4double  kinE,
const G4Material mat = 0 
)
virtual

Implements G4VAtomDeexcitation.

Definition at line 389 of file G4UAtomicDeexcitation.cc.

395{
396 return GetShellIonisationCrossSectionPerAtom(p,Z,shell,kinE,mat);
397}
virtual G4double GetShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0)

◆ GenerateParticles()

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

Implements G4VAtomDeexcitation.

Definition at line 234 of file G4UAtomicDeexcitation.cc.

240{
241
242 // Defined initial conditions
243 G4int givenShellId = atomicShell->ShellId();
244 // G4cout << "generating particles for vacancy in shellId: " << givenShellId << G4endl;
245 minGammaEnergy = gammaCut;
246 minElectronEnergy = eCut;
247
248 // V.I. short-cut
249 if(!IsAugerActive()) { minElectronEnergy = DBL_MAX; }
250
251 // generation secondaries
252 G4DynamicParticle* aParticle=0;
253 G4int provShellId = 0;
254 G4int counter = 0;
255
256 // let's check that 5<Z<100
257
258 if (Z>5 && Z<100) {
259
260 // The aim of this loop is to generate more than one fluorecence photon
261 // from the same ionizing event
262 do
263 {
264 if (counter == 0)
265 // First call to GenerateParticles(...):
266 // givenShellId is given by the process
267 {
268 provShellId = SelectTypeOfTransition(Z, givenShellId);
269
270 if ( provShellId >0)
271 {
272 aParticle = GenerateFluorescence(Z,givenShellId,provShellId);
273 //if (aParticle != 0) { G4cout << "****FLUO!_1**** " << aParticle->GetParticleDefinition()->GetParticleType() << " " << aParticle->GetKineticEnergy()/keV << G4endl ;} //debug
274 }
275 else if ( provShellId == -1)
276 {
277 aParticle = GenerateAuger(Z, givenShellId);
278 //if (aParticle != 0) { G4cout << "****AUGER!****" << G4endl;} //debug
279 }
280 else
281 {
282 G4Exception("G4UAtomicDeexcitation::GenerateParticles()","de0002",JustWarning, "Energy deposited locally");
283 }
284 }
285 else
286 // Following calls to GenerateParticles(...):
287 // newShellId is given by GenerateFluorescence(...)
288 {
289 provShellId = SelectTypeOfTransition(Z,newShellId);
290 if (provShellId >0)
291 {
292 aParticle = GenerateFluorescence(Z,newShellId,provShellId);
293 //if (aParticle != 0) { G4cout << "****FLUO!_2****" << aParticle->GetParticleDefinition()->GetParticleType() << " " << aParticle->GetKineticEnergy()/keV << G4endl;} //debug
294 }
295 else if ( provShellId == -1)
296 {
297 aParticle = GenerateAuger(Z, newShellId);
298 //if (aParticle != 0) { G4cout << "****AUGER!****" << G4endl;} //debug
299 }
300 else
301 {
302 G4Exception("G4UAtomicDeexcitation::GenerateParticles()","de0002",JustWarning, "Energy deposited locally");
303 }
304 }
305 counter++;
306 if (aParticle != 0)
307 {
308 vectorOfParticles->push_back(aParticle);
309 //G4cout << "Deexcitation Occurred!" << G4endl; //debug
310 }
311 else {provShellId = -2;}
312 }
313 while (provShellId > -2);
314 }
315 else
316 {
317 G4Exception("G4UAtomicDeexcitation::GenerateParticles()","de0001",JustWarning, "Energy deposited locally");
318 }
319
320 // G4cout << "---------FATTO!---------" << G4endl; //debug
321
322}
@ JustWarning
int G4int
Definition: G4Types.hh:66
G4int ShellId() const
G4bool IsAugerActive() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ GetAtomicShell()

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

Implements G4VAtomDeexcitation.

Definition at line 229 of file G4UAtomicDeexcitation.cc.

230{
231 return transitionManager->Shell(Z, size_t(shell));
232}
G4AtomicShell * Shell(G4int Z, size_t shellIndex) const

◆ GetShellIonisationCrossSectionPerAtom()

G4double G4UAtomicDeexcitation::GetShellIonisationCrossSectionPerAtom ( const G4ParticleDefinition pdef,
G4int  Z,
G4AtomicShellEnumerator  shell,
G4double  kinE,
const G4Material mat = 0 
)
virtual

Implements G4VAtomDeexcitation.

Definition at line 325 of file G4UAtomicDeexcitation.cc.

331{
332
333 // we must put a control on the shell that are passed:
334 // some shells should not pass (line "0" or "2")
335
336 // check atomic number
337 G4double xsec = 0.0;
338 if(Z > 93 || Z < 6 ) { return xsec; } //corrected by alf - Z<6 missing
339 G4int idx = G4int(shellEnum);
340 if(idx >= G4AtomicShells::GetNumberOfShells(Z)) { return xsec; }
341
342 //
343 if(pdef == theElectron || pdef == thePositron) {
344 xsec = ePIXEshellCS->CrossSection(Z,shellEnum,kineticEnergy,0.0,mat);
345 return xsec;
346 }
347
348 G4double mass = pdef->GetPDGMass();
349 G4double escaled = kineticEnergy;
350 G4double q2 = 0.0;
351
352 // scaling to protons
353 if ((pdef->GetParticleName() != "proton" && pdef->GetParticleName() != "alpha" ) )
354 {
355 mass = proton_mass_c2;
356 escaled = kineticEnergy*mass/(pdef->GetPDGMass());
357
358 if(mat) {
359 q2 = emcorr->EffectiveChargeSquareRatio(pdef,mat,kineticEnergy);
360 } else {
361 G4double q = pdef->GetPDGCharge()/eplus;
362 q2 = q*q;
363 }
364 }
365
366 if(PIXEshellCS) { xsec = PIXEshellCS->CrossSection(Z,shellEnum,escaled,mass,mat); }
367 if(xsec < 1e-100) {
368
369 xsec = anaPIXEshellCS->CrossSection(Z,shellEnum,escaled,mass,mat);
370
371 }
372
373 if (q2) {xsec *= q2;}
374
375 return xsec;
376}
double G4double
Definition: G4Types.hh:64
static G4int GetNumberOfShells(G4int Z)
G4double EffectiveChargeSquareRatio(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
G4double GetPDGCharge() const
const G4String & GetParticleName() const
virtual G4double CrossSection(G4int Z, G4AtomicShellEnumerator shell, G4double incidentEnergy, G4double mass, const G4Material *mat)=0

Referenced by ComputeShellIonisationCrossSectionPerAtom().

◆ InitialiseForExtraAtom()

void G4UAtomicDeexcitation::InitialiseForExtraAtom ( G4int  Z)
virtual

Implements G4VAtomDeexcitation.

Definition at line 226 of file G4UAtomicDeexcitation.cc.

227{}

◆ InitialiseForNewRun()

void G4UAtomicDeexcitation::InitialiseForNewRun ( )
virtual

Implements G4VAtomDeexcitation.

Definition at line 98 of file G4UAtomicDeexcitation.cc.

99{
100 if(!IsFluoActive()) { return; }
101 transitionManager = G4AtomicTransitionManager::Instance();
102 if(IsPIXEActive()) {
103 G4cout << G4endl;
104 G4cout << "### === G4UAtomicDeexcitation::InitialiseForNewRun()" << G4endl;
105 anaPIXEshellCS = new G4teoCrossSection("Analytical");
106
107 }
108 else {return;}
109 // initializing PIXE x-section name
110 //
111 if (PIXECrossSectionModel() == "" ||
112 PIXECrossSectionModel() == "Empirical" ||
113 PIXECrossSectionModel() == "empirical")
114 {
115 SetPIXECrossSectionModel("Empirical");
116 }
117 else if (PIXECrossSectionModel() == "ECPSSR_Analytical" ||
118 PIXECrossSectionModel() == "Analytical" ||
119 PIXECrossSectionModel() == "analytical")
120 {
121 SetPIXECrossSectionModel("Analytical");
122 }
123 else if (PIXECrossSectionModel() == "ECPSSR_FormFactor" ||
124 PIXECrossSectionModel() == "ECPSSR_Tabulated" ||
125 PIXECrossSectionModel() == "Analytical_Tabulated")
126 {
127 SetPIXECrossSectionModel("ECPSSR_FormFactor");
128 }
129 else
130 {
131 G4cout << "### G4UAtomicDeexcitation::InitialiseForNewRun WARNING "
132 << G4endl;
133 G4cout << " PIXE cross section name " << PIXECrossSectionModel()
134 << " is unknown, Analytical cross section will be used" << G4endl;
135 SetPIXECrossSectionModel("Analytical");
136 }
137
138 // Check if old model should be deleted
139 if(PIXEshellCS)
140 {
141 if(PIXECrossSectionModel() != PIXEshellCS->GetName())
142 {
143 delete PIXEshellCS;
144 PIXEshellCS = 0;
145 }
146 }
147
148 // Instantiate empirical model
149 if(!PIXEshellCS) {
150 if (PIXECrossSectionModel() == "Empirical")
151 {
152 PIXEshellCS = new G4empCrossSection("Empirical");
153 }
154
155 if (PIXECrossSectionModel() == "ECPSSR_FormFactor")
156 {
157 PIXEshellCS = new G4teoCrossSection("ECPSSR_FormFactor");
158 }
159 }
160
161 // Electron cross section
162 // initializing PIXE x-section name
163 //
164 if (PIXEElectronCrossSectionModel() == "" ||
165 PIXEElectronCrossSectionModel() == "Livermore")
166 {
168 }
169 else if (PIXEElectronCrossSectionModel() == "ProtonAnalytical" ||
170 PIXEElectronCrossSectionModel() == "Analytical" ||
171 PIXEElectronCrossSectionModel() == "analytical")
172 {
173 SetPIXEElectronCrossSectionModel("ProtonAnalytical");
174 }
175 else if (PIXEElectronCrossSectionModel() == "ProtonEmpirical" ||
176 PIXEElectronCrossSectionModel() == "Empirical" ||
177 PIXEElectronCrossSectionModel() == "empirical")
178 {
179 SetPIXEElectronCrossSectionModel("ProtonEmpirical");
180 }
181 else if (PIXEElectronCrossSectionModel() == "Penelope")
183 else
184 {
185 G4cout << "### G4UAtomicDeexcitation::InitialiseForNewRun WARNING "
186 << G4endl;
187 G4cout << " PIXE e- cross section name " << PIXEElectronCrossSectionModel()
188 << " is unknown, PIXE is disabled" << G4endl;
190 }
191
192 // Check if old model should be deleted
193 if(ePIXEshellCS)
194 {
195 if(PIXEElectronCrossSectionModel() != ePIXEshellCS->GetName())
196 {
197 delete ePIXEshellCS;
198 ePIXEshellCS = 0;
199 }
200 }
201
202 // Instantiate empirical model
203 if(!ePIXEshellCS)
204 {
205 if(PIXEElectronCrossSectionModel() == "Empirical")
206 {
207 ePIXEshellCS = new G4empCrossSection("Empirical");
208 }
209
210 else if(PIXEElectronCrossSectionModel() == "Analytical")
211 {
212 ePIXEshellCS = new G4teoCrossSection("Analytical");
213 }
214
215 else if(PIXEElectronCrossSectionModel() == "Livermore")
216 {
217 ePIXEshellCS = new G4LivermoreIonisationCrossSection();
218 }
219 else if (PIXEElectronCrossSectionModel() == "Penelope")
220 {
221 ePIXEshellCS = new G4PenelopeIonisationCrossSection();
222 }
223 }
224}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
static G4AtomicTransitionManager * Instance()
const G4String & PIXECrossSectionModel() const
void SetPIXECrossSectionModel(const G4String &)
const G4String & PIXEElectronCrossSectionModel() const
void SetPIXEElectronCrossSectionModel(const G4String &)
const G4String & GetName() const

◆ SetCutForAugerElectrons()

void G4UAtomicDeexcitation::SetCutForAugerElectrons ( G4double  cut)

Definition at line 383 of file G4UAtomicDeexcitation.cc.

384{
385 minElectronEnergy = cut;
386}

◆ SetCutForSecondaryPhotons()

void G4UAtomicDeexcitation::SetCutForSecondaryPhotons ( G4double  cut)

Definition at line 378 of file G4UAtomicDeexcitation.cc.

379{
380 minGammaEnergy = cut;
381}

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