Geant4_10
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: G4VAtomDeexcitation.cc 74376 2013-10-04 08:25:47Z gcosmo $
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 
48 #include "G4VAtomDeexcitation.hh"
49 #include "G4SystemOfUnits.hh"
50 #include "G4ParticleDefinition.hh"
51 #include "G4DynamicParticle.hh"
52 #include "G4Step.hh"
53 #include "G4Region.hh"
54 #include "G4RegionStore.hh"
55 #include "G4MaterialCutsCouple.hh"
56 #include "G4MaterialCutsCouple.hh"
57 #include "G4Material.hh"
58 #include "G4Element.hh"
59 #include "G4ElementVector.hh"
60 #include "Randomize.hh"
61 #include "G4VParticleChange.hh"
62 #include "G4PhysicsModelCatalog.hh"
63 #include "G4Gamma.hh"
64 
65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66 
67 G4int G4VAtomDeexcitation::pixeIDg = -1;
68 G4int G4VAtomDeexcitation::pixeIDe = -1;
69 
71  const G4String& pname)
72  : lowestKinEnergy(keV), verbose(1), name(modname), namePIXE(pname),
73  nameElectronPIXE(""), isActive(false), flagAuger(false), flagPIXE(false)
74 {
75  vdyn.reserve(5);
76  theCoupleTable = 0;
77  G4String gg = "gammaPIXE";
78  G4String ee = "e-PIXE";
79  if(pixeIDg < 0) { pixeIDg = G4PhysicsModelCatalog::Register(gg); }
80  if(pixeIDe < 0) { pixeIDe = G4PhysicsModelCatalog::Register(ee); }
81  gamma = G4Gamma::Gamma();
82 }
83 
84 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
85 
87 {}
88 
89 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
90 
92 {
93  // Define list of couples
95  size_t numOfCouples = theCoupleTable->GetTableSize();
96 
97  // needed for unit tests
98  if(0 == numOfCouples) { numOfCouples = 1; }
99 
100  activeDeexcitationMedia.resize(numOfCouples, false);
101  activeAugerMedia.resize(numOfCouples, false);
102  activePIXEMedia.resize(numOfCouples, false);
103  activeZ.resize(93, false);
104 
105  // check if deexcitation is active for the given run
106  if( !isActive ) { return; }
107 
108  // Define list of regions
109  size_t nRegions = deRegions.size();
110 
111  if(0 == nRegions) {
112  SetDeexcitationActiveRegion("World",isActive,flagAuger,flagPIXE);
113  nRegions = 1;
114  }
115 
116  if(0 < verbose) {
117  G4cout << G4endl;
118  G4cout << "### === Deexcitation model " << name
119  << " is activated for " << nRegions;
120  if(1 == nRegions) { G4cout << " region:" << G4endl; }
121  else { G4cout << " regions:" << G4endl;}
122  }
123 
124  // Identify active media
125  G4RegionStore* regionStore = G4RegionStore::GetInstance();
126  for(size_t j=0; j<nRegions; ++j) {
127  const G4Region* reg = regionStore->GetRegion(activeRegions[j], false);
128  const G4ProductionCuts* rpcuts = reg->GetProductionCuts();
129  if(0 < verbose) {
130  G4cout << " " << activeRegions[j] << G4endl;
131  }
132 
133  for(size_t i=0; i<numOfCouples; ++i) {
134  const G4MaterialCutsCouple* couple =
135  theCoupleTable->GetMaterialCutsCouple(i);
136  if (couple->GetProductionCuts() == rpcuts) {
137  activeDeexcitationMedia[i] = deRegions[j];
138  activeAugerMedia[i] = AugerRegions[j];
139  activePIXEMedia[i] = PIXERegions[j];
140  const G4Material* mat = couple->GetMaterial();
141  const G4ElementVector* theElementVector =
142  mat->GetElementVector();
143  G4int nelm = mat->GetNumberOfElements();
144  if(deRegions[j]) {
145  for(G4int k=0; k<nelm; ++k) {
146  G4int Z = G4lrint(((*theElementVector)[k])->GetZ());
147  if(Z > 5 && Z < 93) {
148  activeZ[Z] = true;
149  //G4cout << "!!! Active de-excitation Z= " << Z << G4endl;
150  }
151  }
152  }
153  }
154  }
155  }
156 
157  // Initialise derived class
159 
160  if(0 < verbose && flagPIXE) {
161  G4cout << "### === PIXE model for hadrons: " << namePIXE
162  << " " << IsPIXEActive()
163  << G4endl;
164  G4cout << "### === PIXE model for e+-: " << nameElectronPIXE
165  << " " << IsPIXEActive()
166  << G4endl;
167  }
168 }
169 
170 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
171 
172 void
174  G4bool valDeexcitation,
175  G4bool valAuger,
176  G4bool valPIXE)
177 {
178  G4String ss = rname;
179  //G4cout << "### G4VAtomDeexcitation::SetDeexcitationActiveRegion " << ss
180  // << " " << valDeexcitation << " " << valAuger
181  // << " " << valPIXE << G4endl;
182  if(ss == "world" || ss == "World" || ss == "WORLD") {
183  ss = "DefaultRegionForTheWorld";
184  }
185  size_t n = deRegions.size();
186  if(n > 0) {
187  for(size_t i=0; i<n; ++i) {
188 
189  // Region already exist
190  if(ss == activeRegions[i]) {
191  deRegions[i] = valDeexcitation;
192  AugerRegions[i] = valAuger;
193  PIXERegions[i] = valPIXE;
194  return;
195  }
196  }
197  }
198  // New region
199  activeRegions.push_back(ss);
200  deRegions.push_back(valDeexcitation);
201  AugerRegions.push_back(valAuger);
202  PIXERegions.push_back(valPIXE);
203 }
204 
205 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
206 
207 void
208 G4VAtomDeexcitation::AlongStepDeexcitation(std::vector<G4Track*>& tracks,
209  const G4Step& step,
210  G4double& eLossMax,
211  G4int coupleIndex)
212 {
213  G4double truelength = step.GetStepLength();
214  if(!flagPIXE && !activePIXEMedia[coupleIndex]) { return; }
215  if(eLossMax <= 0.0 || truelength <= 0.0) { return; }
216 
217  // step parameters
218  const G4StepPoint* preStep = step.GetPreStepPoint();
219  G4ThreeVector prePos = preStep->GetPosition();
220  G4ThreeVector delta = step.GetPostStepPoint()->GetPosition() - prePos;
221  G4double preTime = preStep->GetGlobalTime();
222  G4double dt = step.GetPostStepPoint()->GetGlobalTime() - preTime;
223 
224  // particle parameters
225  const G4Track* track = step.GetTrack();
226  const G4ParticleDefinition* part = track->GetDefinition();
227  G4double ekin = preStep->GetKineticEnergy();
228 
229  // media parameters
230  G4double gCut = (*theCoupleTable->GetEnergyCutsVector(0))[coupleIndex];
231  G4double eCut = DBL_MAX;
232  if(CheckAugerActiveRegion(coupleIndex)) {
233  eCut = (*theCoupleTable->GetEnergyCutsVector(1))[coupleIndex];
234  }
235 
236  //G4cout<<"!Sample PIXE gCut(MeV)= "<<gCut<<" eCut(MeV)= "<<eCut
237  // <<" Ekin(MeV)= " << ekin/MeV << G4endl;
238 
239  const G4Material* material = preStep->GetMaterial();
240  const G4ElementVector* theElementVector = material->GetElementVector();
241  const G4double* theAtomNumDensityVector =
242  material->GetVecNbOfAtomsPerVolume();
243  G4int nelm = material->GetNumberOfElements();
244 
245  // loop over deexcitations
246  for(G4int i=0; i<nelm; ++i) {
247  G4int Z = G4lrint((*theElementVector)[i]->GetZ());
248  if(activeZ[Z] && Z < 93) {
249  G4int nshells =
250  std::min(9,(*theElementVector)[i]->GetNbOfAtomicShells());
251  G4double rho = truelength*theAtomNumDensityVector[i];
252  //G4cout<<" Z "<< Z <<" is active x(mm)= " << truelength/mm << G4endl;
253  for(G4int ii=0; ii<nshells; ++ii) {
255  const G4AtomicShell* shell = GetAtomicShell(Z, as);
257 
258  if(gCut > bindingEnergy) { break; }
259 
260  if(eLossMax > bindingEnergy) {
261  G4double sig = rho*
262  GetShellIonisationCrossSectionPerAtom(part, Z, as, ekin, material);
263 
264  // mfp is mean free path in units of step size
265  if(sig > 0.0) {
266  G4double mfp = 1.0/sig;
267  G4double stot = 0.0;
268  //G4cout << " Shell " << ii << " mfp(mm)= " << mfp/mm << G4endl;
269  // sample ionisation points
270  do {
271  stot -= mfp*std::log(G4UniformRand());
272  if( stot > 1.0 || eLossMax < bindingEnergy) { break; }
273  // sample deexcitation
274  vdyn.clear();
275  GenerateParticles(&vdyn, shell, Z, gCut, eCut);
276  G4int nsec = vdyn.size();
277  if(nsec > 0) {
278  G4ThreeVector r = prePos + stot*delta;
279  G4double time = preTime + stot*dt;
280  for(G4int j=0; j<nsec; ++j) {
281  G4DynamicParticle* dp = vdyn[j];
282  G4double e = dp->GetKineticEnergy();
283 
284  // save new secondary if there is enough energy
285  if(eLossMax >= e) {
286  eLossMax -= e;
287  G4Track* t = new G4Track(dp, time, r);
288 
289  // defined secondary type
290  if(dp->GetDefinition() == gamma) {
291  t->SetCreatorModelIndex(pixeIDg);
292  } else {
293  t->SetCreatorModelIndex(pixeIDe);
294  }
295 
296  tracks.push_back(t);
297  } else {
298  delete dp;
299  }
300  }
301  }
302  } while (stot < 1.0);
303  }
304  }
305  }
306  }
307  }
308  return;
309 }
310 
311 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
G4ParticleDefinition * GetDefinition() const
G4ProductionCuts * GetProductionCuts() const
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
void AlongStepDeexcitation(std::vector< G4Track * > &tracks, const G4Step &step, G4double &eLoss, G4int coupleIndex)
std::vector< G4Element * > G4ElementVector
G4double GetKineticEnergy() const
G4double GetStepLength() const
G4bool IsPIXEActive() const
G4Material * GetMaterial() const
virtual G4double GetShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0)=0
void SetDeexcitationActiveRegion(const G4String &rname, G4bool valDeexcitation, G4bool valAuger, G4bool valPIXE)
const XML_Char * name
Definition: expat.h:151
G4ParticleDefinition * GetDefinition() const
G4double BindingEnergy() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
void SetCreatorModelIndex(G4int idx)
Definition: G4Track.hh:245
Float_t mat
Definition: plot.C:40
static G4RegionStore * GetInstance()
string material
Definition: eplot.py:19
G4StepPoint * GetPreStepPoint() const
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
Char_t n[5]
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetPosition() const
Float_t Z
Definition: plot.C:39
bool G4bool
Definition: G4Types.hh:79
G4VAtomDeexcitation(const G4String &modname="Deexcitation", const G4String &pixename="")
TString part[npart]
Definition: Style.C:32
Definition: G4Step.hh:76
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
string pname
Definition: eplot.py:33
static G4ProductionCutsTable * GetProductionCutsTable()
jump r
Definition: plot.C:36
G4bool CheckAugerActiveRegion(G4int coupleIndex)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
int G4lrint(double ad)
Definition: templates.hh:163
G4StepPoint * GetPostStepPoint() const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
virtual void InitialiseForNewRun()=0
#define G4endl
Definition: G4ios.hh:61
G4double GetGlobalTime() const
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
static G4int Register(G4String &)
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
G4Track * GetTrack() const
G4double bindingEnergy(G4int A, G4int Z)
G4ProductionCuts * GetProductionCuts() const
#define DBL_MAX
Definition: templates.hh:83
G4AtomicShellEnumerator
const G4Material * GetMaterial() const