Geant4  10.03
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 101248 2016-11-10 08:51:37Z 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), flagAuger(false),
72  flagAugerCascade(false), flagPIXE(false), ignoreCuts(false),
73  isActiveLocked(false), isAugerLocked(false),
74  isAugerCascadeLocked(false), isPIXELocked(false)
75 {
77  vdyn.reserve(5);
78  theCoupleTable = nullptr;
79  G4String gg = "gammaPIXE";
80  G4String ee = "e-PIXE";
84 }
85 
86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
87 
89 {}
90 
91 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
92 
94 {
96 
97  // Define list of couples
99  G4int numOfCouples = theCoupleTable->GetTableSize();
100 
101  // needed for unit tests
102  G4int nn = std::max(numOfCouples, 1);
103  activeDeexcitationMedia.resize(nn, false);
104  activeAugerMedia.resize(nn, false);
105  activePIXEMedia.resize(nn, false);
106  activeZ.resize(93, false);
107 
108  // initialisation of flags and options
109  // normally there is no locksed flags
113  if(!isPIXELocked) { flagPIXE = theParameters->Pixe(); }
115 
116  // Define list of regions
117  size_t nRegions = deRegions.size();
118  // check if deexcitation is active for the given run
119  if(!isActive && 0 == nRegions) { return; }
120 
121  // if no active regions add a world
122  if(0 == nRegions) {
124  nRegions = deRegions.size();
125  }
126 
127  if(0 < verbose) {
128  G4cout << G4endl;
129  G4cout << "### === Deexcitation model " << name
130  << " is activated for " << nRegions;
131  if(1 == nRegions) { G4cout << " region:" << G4endl; }
132  else { G4cout << " regions:" << G4endl;}
133  }
134 
135  // Identify active media
136  G4RegionStore* regionStore = G4RegionStore::GetInstance();
137  for(size_t j=0; j<nRegions; ++j) {
138  const G4Region* reg = regionStore->GetRegion(activeRegions[j], false);
139  if(reg && 0 < numOfCouples) {
140  const G4ProductionCuts* rpcuts = reg->GetProductionCuts();
141  if(0 < verbose) {
142  G4cout << " " << activeRegions[j]
143  << " " << deRegions[j] << " " << AugerRegions[j]
144  << " " << PIXERegions[j] << G4endl;
145  }
146  for(G4int i=0; i<numOfCouples; ++i) {
147  const G4MaterialCutsCouple* couple =
149  if (couple->GetProductionCuts() == rpcuts) {
153  }
154  }
155  }
156  }
158  //G4cout << nelm << G4endl;
159  for(G4int k=0; k<nelm; ++k) {
160  G4int Z = (*(G4Element::GetElementTable()))[k]->GetZasInt();
161  if(Z > 5 && Z < 93) {
162  activeZ[Z] = true;
163  //G4cout << "!!! Active de-excitation Z= " << Z << G4endl;
164  }
165  }
166 
167  // Initialise derived class
169 
170  if(0 < verbose && flagAuger) {
171  G4cout << "### === Auger cascade flag: " << flagAugerCascade
172  << G4endl;
173  }
174  if(0 < verbose) {
175  G4cout << "### === Ignore cuts flag: " << ignoreCuts
176  << G4endl;
177  }
178  if(0 < verbose && flagPIXE) {
179  G4cout << "### === PIXE model for hadrons: "
181  << G4endl;
182  G4cout << "### === PIXE model for e+-: "
184  << G4endl;
185  }
186 }
187 
188 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
189 
190 void
192  G4bool valDeexcitation,
193  G4bool valAuger,
194  G4bool valPIXE)
195 {
196  // no PIXE in parallel world
197  if(rname == "DefaultRegionForParallelWorld") { return; }
198 
199  G4String ss = rname;
200  /*
201  G4cout << "### G4VAtomDeexcitation::SetDeexcitationActiveRegion " << ss
202  << " " << valDeexcitation << " " << valAuger
203  << " " << valPIXE << G4endl;
204  */
205  if(ss == "world" || ss == "World" || ss == "WORLD") {
206  ss = "DefaultRegionForTheWorld";
207  }
208  size_t n = deRegions.size();
209  for(size_t i=0; i<n; ++i) {
210 
211  // Region already exist
212  if(ss == activeRegions[i]) {
213  deRegions[i] = valDeexcitation;
214  AugerRegions[i] = valAuger;
215  PIXERegions[i] = valPIXE;
216  return;
217  }
218  }
219  // New region
220  activeRegions.push_back(ss);
221  deRegions.push_back(valDeexcitation);
222  AugerRegions.push_back(valAuger);
223  PIXERegions.push_back(valPIXE);
224 
225  // if de-excitation defined for the world volume
226  // it should be active for all G4Regions
227  if(ss == "DefaultRegionForTheWorld") {
229  G4int nn = regions->size();
230  for(G4int i=0; i<nn; ++i) {
231  if(ss == (*regions)[i]->GetName()) { continue; }
232  SetDeexcitationActiveRegion((*regions)[i]->GetName(), valDeexcitation,
233  valAuger, valPIXE);
234 
235  }
236  }
237 }
238 
239 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
240 
241 void
242 G4VAtomDeexcitation::AlongStepDeexcitation(std::vector<G4Track*>& tracks,
243  const G4Step& step,
244  G4double& eLossMax,
245  G4int coupleIndex)
246 {
247  G4double truelength = step.GetStepLength();
248  if(!flagPIXE && !activePIXEMedia[coupleIndex]) { return; }
249  if(eLossMax <= 0.0 || truelength <= 0.0) { return; }
250 
251  // step parameters
252  const G4StepPoint* preStep = step.GetPreStepPoint();
253  G4ThreeVector prePos = preStep->GetPosition();
254  G4ThreeVector delta = step.GetPostStepPoint()->GetPosition() - prePos;
255  G4double preTime = preStep->GetGlobalTime();
256  G4double dt = step.GetPostStepPoint()->GetGlobalTime() - preTime;
257 
258  // particle parameters
259  const G4Track* track = step.GetTrack();
260  const G4ParticleDefinition* part = track->GetDefinition();
261  G4double ekin = preStep->GetKineticEnergy();
262 
263  // media parameters
264  G4double gCut = (*theCoupleTable->GetEnergyCutsVector(0))[coupleIndex];
265  if(ignoreCuts) { gCut = 0.0; }
266  G4double eCut = DBL_MAX;
267  if(CheckAugerActiveRegion(coupleIndex)) {
268  eCut = (*theCoupleTable->GetEnergyCutsVector(1))[coupleIndex];
269  if(ignoreCuts) { eCut = 0.0; }
270  }
271 
272  //G4cout<<"!Sample PIXE gCut(MeV)= "<<gCut<<" eCut(MeV)= "<<eCut
273  // <<" Ekin(MeV)= " << ekin/MeV << G4endl;
274 
275  const G4Material* material = preStep->GetMaterial();
276  const G4ElementVector* theElementVector = material->GetElementVector();
277  const G4double* theAtomNumDensityVector =
278  material->GetVecNbOfAtomsPerVolume();
279  G4int nelm = material->GetNumberOfElements();
280 
281  // loop over deexcitations
282  for(G4int i=0; i<nelm; ++i) {
283  G4int Z = (*theElementVector)[i]->GetZasInt();
284  if(activeZ[Z] && Z < 93) {
285  G4int nshells =
286  std::min(9,(*theElementVector)[i]->GetNbOfAtomicShells());
287  G4double rho = truelength*theAtomNumDensityVector[i];
288  //G4cout<<" Z "<< Z <<" is active x(mm)= " << truelength/mm << G4endl;
289  for(G4int ii=0; ii<nshells; ++ii) {
291  const G4AtomicShell* shell = GetAtomicShell(Z, as);
293 
294  if(gCut > bindingEnergy) { break; }
295 
296  if(eLossMax > bindingEnergy) {
297  G4double sig = rho*
298  GetShellIonisationCrossSectionPerAtom(part, Z, as, ekin, material);
299 
300  // mfp is mean free path in units of step size
301  if(sig > 0.0) {
302  G4double mfp = 1.0/sig;
303  G4double stot = 0.0;
304  //G4cout << " Shell " << ii << " mfp(mm)= " << mfp/mm << G4endl;
305  // sample ionisation points
306  do {
307  stot -= mfp*std::log(G4UniformRand());
308  if( stot > 1.0 || eLossMax < bindingEnergy) { break; }
309  // sample deexcitation
310  vdyn.clear();
311  GenerateParticles(&vdyn, shell, Z, gCut, eCut);
312  G4int nsec = vdyn.size();
313  if(nsec > 0) {
314  G4ThreeVector r = prePos + stot*delta;
315  G4double time = preTime + stot*dt;
316  for(G4int j=0; j<nsec; ++j) {
317  G4DynamicParticle* dp = vdyn[j];
318  G4double e = dp->GetKineticEnergy();
319 
320  // save new secondary if there is enough energy
321  if(eLossMax >= e) {
322  eLossMax -= e;
323  G4Track* t = new G4Track(dp, time, r);
324 
325  // defined secondary type
326  if(dp->GetDefinition() == gamma) {
328  } else {
330  }
331 
332  tracks.push_back(t);
333  } else {
334  delete dp;
335  }
336  }
337  }
338  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
339  } while (stot < 1.0);
340  }
341  }
342  }
343  }
344  }
345  return;
346 }
347 
348 //....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
virtual G4double GetShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=nullptr)=0
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
G4double GetStepLength() const
G4bool AugerCascade() const
std::vector< G4DynamicParticle * > vdyn
G4Material * GetMaterial() const
void SetDeexcitationActiveRegion(const G4String &rname, G4bool valDeexcitation, G4bool valAuger, G4bool valPIXE)
G4ParticleDefinition * GetDefinition() const
const char * name(G4int ptype)
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:405
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
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()
void DefineRegParamForDeex(G4VAtomDeexcitation *) const
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:398
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