Geant4 10.7.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")
 
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 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)
 
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=nullptr)=0
 
virtual G4double ComputeShellIonisationCrossSectionPerAtom (const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=nullptr)=0
 
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()

G4UAtomicDeexcitation::G4UAtomicDeexcitation ( )

Definition at line 75 of file G4UAtomicDeexcitation.cc.

75 :
76 G4VAtomDeexcitation("UAtomDeexcitation"),
77 minGammaEnergy(DBL_MAX),
78 minElectronEnergy(DBL_MAX)
79{
80 anaPIXEshellCS = nullptr;
81 PIXEshellCS = nullptr;
82 ePIXEshellCS = nullptr;
84 theElectron = G4Electron::Electron();
85 thePositron = G4Positron::Positron();
86 transitionManager = G4AtomicTransitionManager::Instance();
87}
static G4AtomicTransitionManager * Instance()
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4LossTableManager * Instance()
G4EmCorrections * EmCorrections()
static G4Positron * Positron()
Definition: G4Positron.cc:93
#define DBL_MAX
Definition: templates.hh:62

◆ ~G4UAtomicDeexcitation()

G4UAtomicDeexcitation::~G4UAtomicDeexcitation ( )
virtual

Definition at line 89 of file G4UAtomicDeexcitation.cc.

90{
91 delete anaPIXEshellCS;
92 delete PIXEshellCS;
93 delete ePIXEshellCS;
94}

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 403 of file G4UAtomicDeexcitation.cc.

409{
410 return GetShellIonisationCrossSectionPerAtom(p,Z,shell,kinE,mat);
411}
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 175 of file G4UAtomicDeexcitation.cc.

181{
182
183 // Defined initial conditions
184 G4int givenShellId = atomicShell->ShellId();
185 //G4cout << "generating particles for vacancy in shellId: "
186 // << givenShellId << G4endl; // debug
187 minGammaEnergy = gammaCut;
188 minElectronEnergy = eCut;
189
190 // generation secondaries
191 G4DynamicParticle* aParticle=0;
192 G4int provShellId = 0;
193
194//ORIGINAL METHOD BY ALFONSO MANTERO
196{
197//----------------------------
198 G4int counter = 0;
199
200 // let's check that 5<Z<100
201
202 if (Z>5 && Z<100) {
203
204 // The aim of this loop is to generate more than one fluorecence photon
205 // from the same ionizing event
206 do
207 {
208 if (counter == 0)
209 // First call to GenerateParticles(...):
210 // givenShellId is given by the process
211 {
212 provShellId = SelectTypeOfTransition(Z, givenShellId);
213
214 if ( provShellId >0)
215 {
216 aParticle = GenerateFluorescence(Z,givenShellId,provShellId);
217 //if (aParticle != 0) {
218 // G4cout << "****FLUO!_1**** "
219 // << aParticle->GetParticleDefinition()->GetParticleType()
220 // << " " << aParticle->GetKineticEnergy()/keV << G4endl ;}
221 }
222 else if ( provShellId == -1)
223 {
224 // G4cout << "Try to generate Auger 1" << G4endl;
225 aParticle = GenerateAuger(Z, givenShellId);
226 // if (aParticle != 0) { G4cout << "****AUGER!****" << G4endl;}
227 }
228 else
229 {
230 //G4Exception("G4UAtomicDeexcitation::GenerateParticles()",
231 // "de0002",JustWarning, "Energy deposited locally");
232 }
233 }
234 else
235 // Following calls to GenerateParticles(...):
236 // newShellId is given by GenerateFluorescence(...)
237 {
238 provShellId = SelectTypeOfTransition(Z,newShellId);
239 if (provShellId >0)
240 {
241 aParticle = GenerateFluorescence(Z,newShellId,provShellId);
242 //if (aParticle != 0) { G4cout << "****FLUO!_2****" << aParticle->GetParticleDefinition()->GetParticleType() << " " << aParticle->GetKineticEnergy()/keV << G4endl;} //debug
243 }
244 else if ( provShellId == -1)
245 {
246 // G4cout << "Try to generate Auger 2" << G4endl; //debug
247 aParticle = GenerateAuger(Z, newShellId);
248 // if (aParticle != 0) { G4cout << "****AUGER!****" << G4endl;} //debug
249 }
250 else
251 {
252 //G4Exception("G4UAtomicDeexcitation::GenerateParticles()","de0002",JustWarning, "Energy deposited locally");
253 }
254 }
255 counter++;
256 if (aParticle != 0)
257 {
258 vectorOfParticles->push_back(aParticle);
259 //G4cout << "Deexcitation Occurred!" << G4endl; //debug
260 }
261 else {provShellId = -2;}
262 }
263 while (provShellId > -2);
264 }
265 else
266 {
267 //G4Exception("G4UAtomicDeexcitation::GenerateParticles()","de0001",JustWarning, "Energy deposited locally");
268 }
269
270 //G4cout << "---------FATTO!---------" << G4endl; //debug
271
272} // Auger cascade is not active
273
274//END OF ORIGINAL METHOD BY ALFONSO MANTERO
275//----------------------
276
277// NEW METHOD
278// Auger cascade by Burkhant Suerfu on March 24 2015 (Bugzilla 1727)
279
281{
282//----------------------
283
284 vacancyArray.push_back(givenShellId);
285
286 // let's check that 5<Z<100
287 if (Z<6 || Z>99){
288 //G4Exception("G4UAtomicDeexcitation::GenerateParticles()","de0001",JustWarning, "Energy deposited locally");
289 return;
290 }
291
292 // as long as there is vacancy to be filled by either fluo or auger, stay in the loop.
293 while(!vacancyArray.empty()){
294
295// prepare to process the last element, and then delete it from the vector.
296 givenShellId = vacancyArray[0];
297 provShellId = SelectTypeOfTransition(Z,givenShellId);
298
299 //G4cout<<"\n------ Atom Transition with Z: "<<Z<<"\tbetween current:"
300 // <<givenShellId<<" & target:"<<provShellId<<G4endl;
301 if(provShellId>0){
302 aParticle = GenerateFluorescence(Z,givenShellId,provShellId);
303// if (aParticle != 0) {
304// G4cout << "****FLUO!_1**** "
305// << aParticle->GetParticleDefinition()->GetParticleType()
306// << " " << aParticle->GetKineticEnergy()/keV << G4endl ;}
307 }
308 else if(provShellId == -1){
309 aParticle = GenerateAuger(Z, givenShellId);
310// if (aParticle != 0) { G4cout << "****AUGER!****" <<
311// aParticle->GetParticleDefinition()->GetParticleType()
312// << " " << aParticle->GetKineticEnergy()/keV << G4endl ; }
313// else G4cout<<G4endl;
314 }
315 //else
316 // G4Exception("G4UAtomicDeexcitation::GenerateParticles()","de0002",JustWarning, "Energy deposited locally");
317
318// if a particle is created, put it in the vector of new particles
319 if(aParticle!=0)
320 vectorOfParticles->push_back(aParticle);
321 else{;}
322// one vacancy has been processed. Erase it.
323 vacancyArray.erase(vacancyArray.begin());
324 }
325
326
327//----------------------
328//End of Auger cascade by Burkhant Suerfu on March 24 2015 (Bugzilla 1727)
329
330} // Auger cascade is active
331
332//ENDSI
333}
int G4int
Definition: G4Types.hh:85
G4int ShellId() const
G4bool IsAugerCascadeActive() const

◆ GetAtomicShell()

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

Implements G4VAtomDeexcitation.

Definition at line 170 of file G4UAtomicDeexcitation.cc.

171{
172 return transitionManager->Shell(Z, size_t(shell));
173}
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 336 of file G4UAtomicDeexcitation.cc.

342{
343 // we must put a control on the shell that are passed:
344 // some shells should not pass (line "0" or "2")
345 //G4cout << pdef->GetParticleName() << " Z= " << Z << " Shell= " << shellEnum
346 // << " E= " << kineticEnergy << G4endl;
347
348 // check atomic number
349 G4double xsec = 0.0;
350 if(Z > 93 || Z < 6 ) { return xsec; } //corrected by alf - Z<6 missing
351 G4int idx = G4int(shellEnum);
352 if(idx >= G4AtomicShells::GetNumberOfShells(Z)) { return xsec; }
353
354 //G4cout << pdef->GetParticleName() << " Z= " << Z << " " << PIXEshellCS
355 // << " " << ePIXEshellCS << G4endl;
356 //
357 if(pdef == theElectron || pdef == thePositron) {
358 xsec = ePIXEshellCS->CrossSection(Z,shellEnum,kineticEnergy,0.0,mat);
359 return xsec;
360 }
361
362 G4double mass = pdef->GetPDGMass();
363 G4double escaled = kineticEnergy;
364 G4double q2 = 0.0;
365
366 // scaling to protons
367 if ((pdef->GetParticleName() != "proton" && pdef->GetParticleName() != "alpha" ) )
368 {
369 mass = proton_mass_c2;
370 escaled = kineticEnergy*mass/(pdef->GetPDGMass());
371
372 if(mat) {
373 q2 = emcorr->EffectiveChargeSquareRatio(pdef,mat,kineticEnergy);
374 } else {
375 G4double q = pdef->GetPDGCharge()/eplus;
376 q2 = q*q;
377 }
378 }
379
380 if(PIXEshellCS) { xsec = PIXEshellCS->CrossSection(Z,shellEnum,escaled,mass,mat); }
381 if(xsec < 1e-100) {
382
383 xsec = anaPIXEshellCS->CrossSection(Z,shellEnum,escaled,mass,mat);
384
385 }
386
387 if (q2) {xsec *= q2;}
388
389 return xsec;
390}
double G4double
Definition: G4Types.hh:83
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 166 of file G4UAtomicDeexcitation.cc.

167{}

◆ InitialiseForNewRun()

void G4UAtomicDeexcitation::InitialiseForNewRun ( )
virtual

Implements G4VAtomDeexcitation.

Definition at line 96 of file G4UAtomicDeexcitation.cc.

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

◆ SetCutForAugerElectrons()

void G4UAtomicDeexcitation::SetCutForAugerElectrons ( G4double  cut)

Definition at line 397 of file G4UAtomicDeexcitation.cc.

398{
399 minElectronEnergy = cut;
400}

◆ SetCutForSecondaryPhotons()

void G4UAtomicDeexcitation::SetCutForSecondaryPhotons ( G4double  cut)

Definition at line 392 of file G4UAtomicDeexcitation.cc.

393{
394 minGammaEnergy = cut;
395}

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