Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VAtomDeexcitation.cc
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 class file
31//
32//
33// File name: G4VAtomDeexcitation
34//
35// Author: Alfonso Mantero & Vladimir Ivanchenko
36//
37// Creation date: 21.04.2010
38//
39// Modifications:
40//
41// Class Description:
42//
43// Abstract interface to energy loss models
44
45// -------------------------------------------------------------------
46//
47
49#include "G4SystemOfUnits.hh"
51#include "G4DynamicParticle.hh"
52#include "G4Step.hh"
53#include "G4Region.hh"
54#include "G4RegionStore.hh"
57#include "G4Material.hh"
58#include "G4Element.hh"
59#include "G4ElementVector.hh"
60#include "Randomize.hh"
61#include "G4VParticleChange.hh"
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
64
66 const G4String& pname)
67 : lowestKinEnergy(keV), verbose(1), name(modname), namePIXE(pname),
68 nameElectronPIXE(""), isActive(false), flagAuger(false), flagPIXE(false)
69{
70 vdyn.reserve(5);
71 theCoupleTable = 0;
72}
73
74//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
75
77{}
78
79//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80
82{
83 // Define list of couples
85 size_t numOfCouples = theCoupleTable->GetTableSize();
86
87 // needed for unit tests
88 if(0 == numOfCouples) { numOfCouples = 1; }
89
90 activeDeexcitationMedia.resize(numOfCouples, false);
91 activeAugerMedia.resize(numOfCouples, false);
92 activePIXEMedia.resize(numOfCouples, false);
93 activeZ.resize(93, false);
94
95 // check if deexcitation is active for the given run
96 if( !isActive ) { return; }
97
98 // Define list of regions
99 size_t nRegions = deRegions.size();
100
101 if(0 == nRegions) {
102 SetDeexcitationActiveRegion("World",isActive,flagAuger,flagPIXE);
103 nRegions = 1;
104 }
105
106 if(0 < verbose) {
107 G4cout << G4endl;
108 G4cout << "### === Deexcitation model " << name
109 << " is activated for " << nRegions;
110 if(1 == nRegions) { G4cout << " region:" << G4endl; }
111 else { G4cout << " regions:" << G4endl;}
112 }
113
114 // Identify active media
116 for(size_t j=0; j<nRegions; ++j) {
117 const G4Region* reg = regionStore->GetRegion(activeRegions[j], false);
118 const G4ProductionCuts* rpcuts = reg->GetProductionCuts();
119 if(0 < verbose) {
120 G4cout << " " << activeRegions[j] << G4endl;
121 }
122
123 for(size_t i=0; i<numOfCouples; ++i) {
124 const G4MaterialCutsCouple* couple =
125 theCoupleTable->GetMaterialCutsCouple(i);
126 if (couple->GetProductionCuts() == rpcuts) {
127 activeDeexcitationMedia[i] = deRegions[j];
128 activeAugerMedia[i] = AugerRegions[j];
129 activePIXEMedia[i] = PIXERegions[j];
130 const G4Material* mat = couple->GetMaterial();
131 const G4ElementVector* theElementVector =
132 mat->GetElementVector();
133 G4int nelm = mat->GetNumberOfElements();
134 if(deRegions[j]) {
135 for(G4int k=0; k<nelm; ++k) {
136 G4int Z = G4lrint(((*theElementVector)[k])->GetZ());
137 if(Z > 5 && Z < 93) {
138 activeZ[Z] = true;
139 //G4cout << "!!! Active de-excitation Z= " << Z << G4endl;
140 }
141 }
142 }
143 }
144 }
145 }
146
147 // Initialise derived class
149
150 if(0 < verbose && flagPIXE) {
151 G4cout << "### === PIXE model for hadrons: " << namePIXE
152 << " " << IsPIXEActive()
153 << G4endl;
154 G4cout << "### === PIXE model for e+-: " << nameElectronPIXE
155 << " " << IsPIXEActive()
156 << G4endl;
157 }
158}
159
160//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
161
162void
164 G4bool valDeexcitation,
165 G4bool valAuger,
166 G4bool valPIXE)
167{
168 G4String ss = rname;
169 //G4cout << "### G4VAtomDeexcitation::SetDeexcitationActiveRegion " << ss
170 // << " " << valDeexcitation << " " << valAuger
171 // << " " << valPIXE << G4endl;
172 if(ss == "world" || ss == "World" || ss == "WORLD") {
173 ss = "DefaultRegionForTheWorld";
174 }
175 size_t n = deRegions.size();
176 if(n > 0) {
177 for(size_t i=0; i<n; ++i) {
178
179 // Region already exist
180 if(ss == activeRegions[i]) {
181 deRegions[i] = valDeexcitation;
182 AugerRegions[i] = valAuger;
183 PIXERegions[i] = valPIXE;
184 return;
185 }
186 }
187 }
188 // New region
189 activeRegions.push_back(ss);
190 deRegions.push_back(valDeexcitation);
191 AugerRegions.push_back(valAuger);
192 PIXERegions.push_back(valPIXE);
193}
194
195//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
196
197void
198G4VAtomDeexcitation::AlongStepDeexcitation(std::vector<G4Track*>& tracks,
199 const G4Step& step,
200 G4double& eLossMax,
201 G4int coupleIndex)
202{
203 G4double truelength = step.GetStepLength();
204 if(!flagPIXE && !activePIXEMedia[coupleIndex]) { return; }
205 if(eLossMax <= 0.0 || truelength <= 0.0) { return; }
206
207 // step parameters
208 const G4StepPoint* preStep = step.GetPreStepPoint();
209 G4ThreeVector prePos = preStep->GetPosition();
210 G4ThreeVector delta = step.GetPostStepPoint()->GetPosition() - prePos;
211 G4double preTime = preStep->GetGlobalTime();
212 G4double dt = step.GetPostStepPoint()->GetGlobalTime() - preTime;
213
214 // particle parameters
215 const G4Track* track = step.GetTrack();
216 const G4ParticleDefinition* part = track->GetDefinition();
217 G4double ekin = preStep->GetKineticEnergy();
218
219 // media parameters
220 G4double gCut = (*theCoupleTable->GetEnergyCutsVector(0))[coupleIndex];
221 G4double eCut = DBL_MAX;
222 if(CheckAugerActiveRegion(coupleIndex)) {
223 eCut = (*theCoupleTable->GetEnergyCutsVector(1))[coupleIndex];
224 }
225
226 //G4cout<<"!Sample PIXE gCut(MeV)= "<<gCut<<" eCut(MeV)= "<<eCut
227 // <<" Ekin(MeV)= " << ekin/MeV << G4endl;
228
229 const G4Material* material = preStep->GetMaterial();
230 const G4ElementVector* theElementVector = material->GetElementVector();
231 const G4double* theAtomNumDensityVector = material->GetVecNbOfAtomsPerVolume();
232 G4int nelm = material->GetNumberOfElements();
233
234 // loop over deexcitations
235 for(G4int i=0; i<nelm; ++i) {
236 G4int Z = G4lrint((*theElementVector)[i]->GetZ());
237 if(activeZ[Z] && Z < 93) {
238 G4int nshells = std::min(9,(*theElementVector)[i]->GetNbOfAtomicShells());
239 G4double rho = truelength*theAtomNumDensityVector[i];
240 //G4cout << " Z " << Z <<" is active x(mm)= " << truelength/mm << G4endl;
241 for(G4int ii=0; ii<nshells; ++ii) {
243 const G4AtomicShell* shell = GetAtomicShell(Z, as);
244 G4double bindingEnergy = shell->BindingEnergy();
245
246 if(gCut > bindingEnergy) { break; }
247
248 if(eLossMax > bindingEnergy) {
249 G4double sig = rho*
250 GetShellIonisationCrossSectionPerAtom(part, Z, as, ekin, material);
251
252 // mfp is mean free path in units of step size
253 if(sig > 0.0) {
254 G4double mfp = 1.0/sig;
255 G4double stot = 0.0;
256 //G4cout << " Shell " << ii << " mfp(mm)= " << mfp/mm << G4endl;
257 // sample ionisation points
258 do {
259 stot -= mfp*std::log(G4UniformRand());
260 if( stot > 1.0 || eLossMax < bindingEnergy) { break; }
261 // sample deexcitation
262 vdyn.clear();
263 GenerateParticles(&vdyn, shell, Z, gCut, eCut);
264 G4int nsec = vdyn.size();
265 if(nsec > 0) {
266 G4ThreeVector r = prePos + stot*delta;
267 G4double time = preTime + stot*dt;
268 for(G4int j=0; j<nsec; ++j) {
269 G4DynamicParticle* dp = vdyn[j];
270 G4double e = dp->GetKineticEnergy();
271
272 // save new secondary if there is enough energy
273 if(eLossMax >= e) {
274 eLossMax -= e;
275 G4Track* t = new G4Track(dp, time, r);
276 tracks.push_back(t);
277 } else {
278 delete dp;
279 }
280 }
281 }
282 } while (stot < 1.0);
283 }
284 }
285 }
286 }
287 }
288 return;
289}
290
291//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4AtomicShellEnumerator
std::vector< G4Element * > G4ElementVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
G4double BindingEnergy() const
G4double GetKineticEnergy() const
const G4Material * GetMaterial() const
G4ProductionCuts * GetProductionCuts() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:185
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:205
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
static G4ProductionCutsTable * GetProductionCutsTable()
static G4RegionStore * GetInstance()
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
G4ProductionCuts * GetProductionCuts() const
G4double GetGlobalTime() const
G4Material * GetMaterial() const
const G4ThreeVector & GetPosition() const
G4double GetKineticEnergy() const
Definition: G4Step.hh:78
G4Track * GetTrack() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
G4StepPoint * GetPostStepPoint() const
G4ParticleDefinition * GetDefinition() const
virtual G4double GetShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0)=0
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
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(const G4String &modname="Deexcitation", const G4String &pixename="")
void SetDeexcitationActiveRegion(const G4String &rname, G4bool valDeexcitation, G4bool valAuger, G4bool valPIXE)
virtual void InitialiseForNewRun()=0
int G4lrint(double ad)
Definition: templates.hh:163
#define DBL_MAX
Definition: templates.hh:83