Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4UAtomicDeexcitation Class Reference

#include <G4UAtomicDeexcitation.hh>

Inheritance diagram for G4UAtomicDeexcitation:
G4VAtomDeexcitation

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

Detailed Description

Definition at line 61 of file G4UAtomicDeexcitation.hh.

Constructor & Destructor Documentation

G4UAtomicDeexcitation::G4UAtomicDeexcitation ( )

Definition at line 75 of file G4UAtomicDeexcitation.cc.

References G4Electron::Electron(), G4LossTableManager::EmCorrections(), G4LossTableManager::Instance(), and G4Positron::Positron().

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

Definition at line 90 of file G4UAtomicDeexcitation.cc.

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

Member Function Documentation

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

Implements G4VAtomDeexcitation.

Definition at line 390 of file G4UAtomicDeexcitation.cc.

References GetShellIonisationCrossSectionPerAtom().

396 {
397  return GetShellIonisationCrossSectionPerAtom(p,Z,shell,kinE,mat);
398 }
virtual G4double GetShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0)
void G4UAtomicDeexcitation::GenerateParticles ( std::vector< G4DynamicParticle * > *  secVect,
const G4AtomicShell atomicShell,
G4int  Z,
G4double  gammaCut,
G4double  eCut 
)
virtual

Implements G4VAtomDeexcitation.

Definition at line 233 of file G4UAtomicDeexcitation.cc.

References G4Exception(), JustWarning, and G4AtomicShell::ShellId().

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

Implements G4VAtomDeexcitation.

Definition at line 228 of file G4UAtomicDeexcitation.cc.

References G4AtomicTransitionManager::Shell().

229 {
230  return transitionManager->Shell(Z, size_t(shell));
231 }
G4AtomicShell * Shell(G4int Z, size_t shellIndex) const
G4double G4UAtomicDeexcitation::GetShellIonisationCrossSectionPerAtom ( const G4ParticleDefinition pdef,
G4int  Z,
G4AtomicShellEnumerator  shell,
G4double  kinE,
const G4Material mat = 0 
)
virtual

Implements G4VAtomDeexcitation.

Definition at line 326 of file G4UAtomicDeexcitation.cc.

References G4VhShellCrossSection::CrossSection(), G4EmCorrections::EffectiveChargeSquareRatio(), python.hepunit::eplus, G4AtomicShells::GetNumberOfShells(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), and python.hepunit::proton_mass_c2.

Referenced by ComputeShellIonisationCrossSectionPerAtom().

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

Implements G4VAtomDeexcitation.

Definition at line 225 of file G4UAtomicDeexcitation.cc.

226 {}
void G4UAtomicDeexcitation::InitialiseForNewRun ( )
virtual

Implements G4VAtomDeexcitation.

Definition at line 97 of file G4UAtomicDeexcitation.cc.

References G4cout, G4endl, G4VhShellCrossSection::GetName(), G4AtomicTransitionManager::Instance(), G4VAtomDeexcitation::IsFluoActive(), G4VAtomDeexcitation::IsPIXEActive(), G4VAtomDeexcitation::PIXECrossSectionModel(), G4VAtomDeexcitation::PIXEElectronCrossSectionModel(), G4VAtomDeexcitation::SetPIXECrossSectionModel(), and G4VAtomDeexcitation::SetPIXEElectronCrossSectionModel().

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

Definition at line 384 of file G4UAtomicDeexcitation.cc.

385 {
386  minElectronEnergy = cut;
387 }
void G4UAtomicDeexcitation::SetCutForSecondaryPhotons ( G4double  cut)

Definition at line 379 of file G4UAtomicDeexcitation.cc.

380 {
381  minGammaEnergy = cut;
382 }

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