Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4VAtomDeexcitation.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class header file
31//
32//
33// File name: G4VAtomDeexcitation
34//
35// Author: Alfonso Mantero & Vladimir Ivanchenko
36//
37// Creation date: 30.06.2009
38//
39// Modifications:
40// 15 Mar 2011 ALF stripped G4AtomicShellEnumerator to its own file
41//
42// Class Description:
43//
44// Abstract interface to energy loss models
45
46// -------------------------------------------------------------------
47//
48
49#ifndef G4VAtomDeexcitation_h
50#define G4VAtomDeexcitation_h 1
51
52#include "globals.hh"
53#include "G4AtomicShell.hh"
56#include "G4Track.hh"
57#include <vector>
58
62
64public:
65
66 G4VAtomDeexcitation(const G4String& modname = "Deexcitation",
67 const G4String& pixename = "");
68
69 virtual ~G4VAtomDeexcitation();
70
71 //========== initialization ==========
72
73 // Overall initialisation before new run
75
76 // Initialisation of deexcitation at the beginning of run
77 virtual void InitialiseForNewRun() = 0;
78
79 // Initialisation for a concrete atom
80 // May be called at run time
81 virtual void InitialiseForExtraAtom(G4int Z) = 0;
82
83 // Activation of deexcitation per detector region
84 void SetDeexcitationActiveRegion(const G4String& rname,
85 G4bool valDeexcitation,
86 G4bool valAuger,
87 G4bool valPIXE);
88
89 // Activation of deexcitation
90 inline void SetFluo(G4bool);
91 inline G4bool IsFluoActive() const;
92
93 // Activation of Auger electron production
94 inline void SetAuger(G4bool);
95 inline G4bool IsAugerActive() const;
96
97 // Activation of PIXE simulation
98 inline void SetPIXE(G4bool);
99 // inline void SetPIXEActive(G4bool);
100 inline G4bool IsPIXEActive() const;
101
102 // Deexcitation model name
103 inline const G4String& GetName() const;
104
105 // PIXE model name
106 inline void SetPIXECrossSectionModel(const G4String&);
107 inline const G4String& PIXECrossSectionModel() const;
108
109 // PIXE model name for e+e-
110 inline void SetPIXEElectronCrossSectionModel(const G4String&);
111 inline const G4String& PIXEElectronCrossSectionModel() const;
112
113 // Access to the list of atoms active for deexcitation
114 inline const std::vector<G4bool>& GetListOfActiveAtoms() const;
115
116 // Verbosity level
117 inline void SetVerboseLevel(G4int);
118 inline G4int GetVerboseLevel() const;
119
120 //========== Run time methods ==========
121
122 // Check if deexcitation is active for a given geometry volume
123 inline G4bool CheckDeexcitationActiveRegion(G4int coupleIndex);
124 inline G4bool CheckAugerActiveRegion(G4int coupleIndex);
125
126 // Get atomic shell by shell index, used by discrete processes
127 // (for example, photoelectric), when shell vacancy sampled by the model
128 virtual
130 G4AtomicShellEnumerator shell) = 0;
131
132 // generation of deexcitation for given atom and shell vacancy
133 inline void GenerateParticles(std::vector<G4DynamicParticle*>* secVect,
134 const G4AtomicShell*,
135 G4int Z,
136 G4int coupleIndex);
137
138 // generation of deexcitation for given atom and shell vacancy
139 virtual void GenerateParticles(std::vector<G4DynamicParticle*>* secVect,
140 const G4AtomicShell*,
141 G4int Z,
142 G4double gammaCut,
143 G4double eCut) = 0;
144
145 // access or compute PIXE cross section
146 virtual G4double
148 G4int Z,
150 G4double kinE,
151 const G4Material* mat = 0) = 0;
152
153 // access or compute PIXE cross section
154 virtual G4double
156 G4int Z,
158 G4double kinE,
159 const G4Material* mat = 0) = 0;
160
161 // Sampling of PIXE for ionisation processes
162 void AlongStepDeexcitation(std::vector<G4Track*>& tracks,
163 const G4Step& step,
164 G4double& eLoss,
165 G4int coupleIndex);
166
167private:
168
169 // copy constructor and hide assignment operator
171 G4VAtomDeexcitation & operator=(const G4VAtomDeexcitation &right);
172
173 G4ProductionCutsTable* theCoupleTable;
174 G4double lowestKinEnergy;
175 G4int verbose;
176 G4String name;
177 G4String namePIXE;
178 G4String nameElectronPIXE;
179 G4bool isActive;
180 G4bool flagAuger;
181 G4bool flagPIXE;
182 std::vector<G4bool> activeZ;
183 std::vector<G4bool> activeDeexcitationMedia;
184 std::vector<G4bool> activeAugerMedia;
185 std::vector<G4bool> activePIXEMedia;
186 std::vector<G4String> activeRegions;
187 std::vector<G4bool> deRegions;
188 std::vector<G4bool> AugerRegions;
189 std::vector<G4bool> PIXERegions;
190 std::vector<G4DynamicParticle*> vdyn;
191};
192
194{
195 isActive = val;
196}
197
199{
200 return isActive;
201}
202
204{
205 flagAuger = val;
206 if(val) { isActive = true; }
207}
208
210{
211 return flagAuger;
212}
213
215{
216 flagPIXE = val;
217 if(val) { isActive = true; }
218}
219
221{
222 return flagPIXE;
223}
224
226{
227 return name;
228}
229
230inline
232{
233 namePIXE = n;
234}
235
236inline void
238{
239 nameElectronPIXE = n;
240}
241
242inline
244{
245 return namePIXE;
246}
247
248inline
250{
251 return nameElectronPIXE;
252}
253
254inline const std::vector<G4bool>&
256{
257 return activeZ;
258}
259
261{
262 verbose = val;
263}
264
266{
267 return verbose;
268}
269
270inline G4bool
272{
273 return (isActive || activeDeexcitationMedia[coupleIndex]);
274}
275
276inline G4bool
278{
279
280 return (flagAuger || activeAugerMedia[coupleIndex]);
281}
282
283inline void
284G4VAtomDeexcitation::GenerateParticles(std::vector<G4DynamicParticle*>* v,
285 const G4AtomicShell* as,
286 G4int Z,
287 G4int idx)
288{
289 G4double gCut = DBL_MAX;
290 if (theCoupleTable) {
291 gCut = (*(theCoupleTable->GetEnergyCutsVector(0)))[idx];
292 }
293 if(gCut < as->BindingEnergy()) {
294 G4double eCut = DBL_MAX;
295 if(CheckAugerActiveRegion(idx)) {
296 if (theCoupleTable) {
297 eCut = (*(theCoupleTable->GetEnergyCutsVector(1)))[idx];
298 }
299 }
300 GenerateParticles(v, as, Z, gCut, eCut);
301 }
302}
303
304#endif
305
G4AtomicShellEnumerator
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
Definition: G4Step.hh:78
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
const std::vector< G4bool > & GetListOfActiveAtoms() const
virtual G4double ComputeShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0)=0
virtual G4double GetShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0)=0
G4bool IsAugerActive() const
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
const G4String & GetName() const
G4bool CheckAugerActiveRegion(G4int coupleIndex)
virtual void InitialiseForExtraAtom(G4int Z)=0
const G4String & PIXECrossSectionModel() const
virtual void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4double gammaCut, G4double eCut)=0
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void SetPIXECrossSectionModel(const G4String &)
const G4String & PIXEElectronCrossSectionModel() const
void SetPIXEElectronCrossSectionModel(const G4String &)
void AlongStepDeexcitation(std::vector< G4Track * > &tracks, const G4Step &step, G4double &eLoss, G4int coupleIndex)
void SetDeexcitationActiveRegion(const G4String &rname, G4bool valDeexcitation, G4bool valAuger, G4bool valPIXE)
virtual void InitialiseForNewRun()=0
#define DBL_MAX
Definition: templates.hh:83