Geant4  10.02.p01
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 92921 2015-09-21 15:06:51Z 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 
69 
71  : verbose(1), name(modname), isActive(false),
72  flagAuger(false),flagAugerCascade(false),
73  flagPIXE(false), ignoreCuts(false)
74 {
76  vdyn.reserve(5);
77  theCoupleTable = nullptr;
78  G4String gg = "gammaPIXE";
79  G4String ee = "e-PIXE";
83 }
84 
85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
86 
88 {}
89 
90 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
91 
93 {
94  // Define list of couples
96  G4int numOfCouples = theCoupleTable->GetTableSize();
97 
98  // needed for unit tests
99  G4int nn = std::max(numOfCouples, 1);
100  activeDeexcitationMedia.resize(nn, false);
101  activeAugerMedia.resize(nn, false);
102  activePIXEMedia.resize(nn, false);
103  activeZ.resize(93, false);
104 
105  // initialisation of flags and options
111 
112  // check if deexcitation is active for the given run
113  if( !isActive ) { return; }
114 
115  // Define list of regions
116  size_t nRegions = deRegions.size();
117 
118  if(0 == nRegions) {
120  nRegions = deRegions.size();
121  }
122 
123  if(0 < verbose) {
124  G4cout << G4endl;
125  G4cout << "### === Deexcitation model " << name
126  << " is activated for " << nRegions;
127  if(1 == nRegions) { G4cout << " region:" << G4endl; }
128  else { G4cout << " regions:" << G4endl;}
129  }
130 
131  // Identify active media
132  G4RegionStore* regionStore = G4RegionStore::GetInstance();
133  for(size_t j=0; j<nRegions; ++j) {
134  const G4Region* reg = regionStore->GetRegion(activeRegions[j], false);
135  if(reg && 0 < numOfCouples) {
136  const G4ProductionCuts* rpcuts = reg->GetProductionCuts();
137  if(0 < verbose) {
138  G4cout << " " << activeRegions[j] << G4endl;
139  }
140  for(G4int i=0; i<numOfCouples; ++i) {
141  const G4MaterialCutsCouple* couple =
143  if (couple->GetProductionCuts() == rpcuts) {
147  }
148  }
149  }
150  }
152  //G4cout << nelm << G4endl;
153  for(G4int k=0; k<nelm; ++k) {
154  G4int Z = G4lrint((*(G4Element::GetElementTable()))[k]->GetZ());
155  if(Z > 5 && Z < 93) {
156  activeZ[Z] = true;
157  //G4cout << "!!! Active de-excitation Z= " << Z << G4endl;
158  }
159  }
160 
161  // Initialise derived class
163 
164  if(0 < verbose && flagPIXE) {
165  G4cout << "### === PIXE model for hadrons: "
167  << G4endl;
168  G4cout << "### === PIXE model for e+-: "
170  << G4endl;
171  }
172 }
173 
174 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
175 
176 void
178  G4bool valDeexcitation,
179  G4bool valAuger,
180  G4bool valPIXE)
181 {
182  // no PIXE in parallel world
183  if(rname == "DefaultRegionForParallelWorld") { return; }
184 
185  G4String ss = rname;
186  //G4cout << "### G4VAtomDeexcitation::SetDeexcitationActiveRegion " << ss
187  // << " " << valDeexcitation << " " << valAuger
188  // << " " << valPIXE << G4endl;
189  if(ss == "world" || ss == "World" || ss == "WORLD") {
190  ss = "DefaultRegionForTheWorld";
191  }
192  size_t n = deRegions.size();
193  if(n > 0) {
194  for(size_t i=0; i<n; ++i) {
195 
196  // Region already exist
197  if(ss == activeRegions[i]) {
198  deRegions[i] = valDeexcitation;
199  AugerRegions[i] = valAuger;
200  PIXERegions[i] = valPIXE;
201  return;
202  }
203  }
204  }
205  // New region
206  activeRegions.push_back(ss);
207  deRegions.push_back(valDeexcitation);
208  AugerRegions.push_back(valAuger);
209  PIXERegions.push_back(valPIXE);
210 
211  // if de-excitation defined fo rthe world volume
212  // it should be active everywhere
213  if(ss == "DefaultRegionForTheWorld") {
215  G4int nn = regions->size();
216  for(G4int i=0; i<nn; ++i) {
217  SetDeexcitationActiveRegion((*regions)[i]->GetName(), valDeexcitation,
218  valAuger, valPIXE);
219 
220  }
221  }
222 }
223 
224 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
225 
226 void
227 G4VAtomDeexcitation::AlongStepDeexcitation(std::vector<G4Track*>& tracks,
228  const G4Step& step,
229  G4double& eLossMax,
230  G4int coupleIndex)
231 {
232  G4double truelength = step.GetStepLength();
233  if(!flagPIXE && !activePIXEMedia[coupleIndex]) { return; }
234  if(eLossMax <= 0.0 || truelength <= 0.0) { return; }
235 
236  // step parameters
237  const G4StepPoint* preStep = step.GetPreStepPoint();
238  G4ThreeVector prePos = preStep->GetPosition();
239  G4ThreeVector delta = step.GetPostStepPoint()->GetPosition() - prePos;
240  G4double preTime = preStep->GetGlobalTime();
241  G4double dt = step.GetPostStepPoint()->GetGlobalTime() - preTime;
242 
243  // particle parameters
244  const G4Track* track = step.GetTrack();
245  const G4ParticleDefinition* part = track->GetDefinition();
246  G4double ekin = preStep->GetKineticEnergy();
247 
248  // media parameters
249  G4double gCut = (*theCoupleTable->GetEnergyCutsVector(0))[coupleIndex];
250  if(ignoreCuts) { gCut = 0.0; }
251  G4double eCut = DBL_MAX;
252  if(CheckAugerActiveRegion(coupleIndex)) {
253  eCut = (*theCoupleTable->GetEnergyCutsVector(1))[coupleIndex];
254  if(ignoreCuts) { eCut = 0.0; }
255  }
256 
257  //G4cout<<"!Sample PIXE gCut(MeV)= "<<gCut<<" eCut(MeV)= "<<eCut
258  // <<" Ekin(MeV)= " << ekin/MeV << G4endl;
259 
260  const G4Material* material = preStep->GetMaterial();
261  const G4ElementVector* theElementVector = material->GetElementVector();
262  const G4double* theAtomNumDensityVector =
263  material->GetVecNbOfAtomsPerVolume();
264  G4int nelm = material->GetNumberOfElements();
265 
266  // loop over deexcitations
267  for(G4int i=0; i<nelm; ++i) {
268  G4int Z = G4lrint((*theElementVector)[i]->GetZ());
269  if(activeZ[Z] && Z < 93) {
270  G4int nshells =
271  std::min(9,(*theElementVector)[i]->GetNbOfAtomicShells());
272  G4double rho = truelength*theAtomNumDensityVector[i];
273  //G4cout<<" Z "<< Z <<" is active x(mm)= " << truelength/mm << G4endl;
274  for(G4int ii=0; ii<nshells; ++ii) {
276  const G4AtomicShell* shell = GetAtomicShell(Z, as);
278 
279  if(gCut > bindingEnergy) { break; }
280 
281  if(eLossMax > bindingEnergy) {
282  G4double sig = rho*
283  GetShellIonisationCrossSectionPerAtom(part, Z, as, ekin, material);
284 
285  // mfp is mean free path in units of step size
286  if(sig > 0.0) {
287  G4double mfp = 1.0/sig;
288  G4double stot = 0.0;
289  //G4cout << " Shell " << ii << " mfp(mm)= " << mfp/mm << G4endl;
290  // sample ionisation points
291  do {
292  stot -= mfp*std::log(G4UniformRand());
293  if( stot > 1.0 || eLossMax < bindingEnergy) { break; }
294  // sample deexcitation
295  vdyn.clear();
296  GenerateParticles(&vdyn, shell, Z, gCut, eCut);
297  G4int nsec = vdyn.size();
298  if(nsec > 0) {
299  G4ThreeVector r = prePos + stot*delta;
300  G4double time = preTime + stot*dt;
301  for(G4int j=0; j<nsec; ++j) {
302  G4DynamicParticle* dp = vdyn[j];
303  G4double e = dp->GetKineticEnergy();
304 
305  // save new secondary if there is enough energy
306  if(eLossMax >= e) {
307  eLossMax -= e;
308  G4Track* t = new G4Track(dp, time, r);
309 
310  // defined secondary type
311  if(dp->GetDefinition() == gamma) {
313  } else {
315  }
316 
317  tracks.push_back(t);
318  } else {
319  delete dp;
320  }
321  }
322  }
323  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
324  } while (stot < 1.0);
325  }
326  }
327  }
328  }
329  }
330  return;
331 }
332 
333 //....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< G4String > activeRegions
std::vector< G4Element * > G4ElementVector
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4double GetStepLength() const
G4bool AugerCascade() const
std::vector< G4DynamicParticle * > vdyn
G4String name
Definition: TRTMaterials.hh:40
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)
G4ParticleDefinition * GetDefinition() const
G4double BindingEnergy() const
G4bool Fluo() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
G4bool DeexcitationIgnoreCut() const
std::vector< G4bool > activeAugerMedia
std::vector< G4bool > AugerRegions
void SetCreatorModelIndex(G4int idx)
static const G4double reg
static G4RegionStore * GetInstance()
G4StepPoint * GetPreStepPoint() const
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:206
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
const G4ThreeVector & GetPosition() const
const G4String & PIXECrossSectionModel()
const G4ParticleDefinition * gamma
bool G4bool
Definition: G4Types.hh:79
Definition: G4Step.hh:76
G4EmParameters * theParameters
const G4int n
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4ProductionCutsTable * theCoupleTable
static G4ProductionCutsTable * GetProductionCutsTable()
G4bool CheckAugerActiveRegion(G4int coupleIndex)
G4bool Auger() const
std::vector< G4bool > activeDeexcitationMedia
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
int G4lrint(double ad)
Definition: templates.hh:163
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static G4int Register(const G4String &)
const G4String & GetName() const
G4StepPoint * GetPostStepPoint() const
static G4EmParameters * Instance()
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
virtual void InitialiseForNewRun()=0
std::vector< G4bool > deRegions
#define G4endl
Definition: G4ios.hh:61
G4double GetGlobalTime() const
G4bool Pixe() const
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
std::vector< G4bool > activeZ
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
G4Track * GetTrack() const
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
G4double bindingEnergy(G4int A, G4int Z)
G4ProductionCuts * GetProductionCuts() const
#define DBL_MAX
Definition: templates.hh:83
G4VAtomDeexcitation(const G4String &modname="Deexcitation")
const G4String & PIXEElectronCrossSectionModel()
std::vector< G4bool > PIXERegions
G4AtomicShellEnumerator
std::vector< G4bool > activePIXEMedia