Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VAtomDeexcitation Class Referenceabstract

#include <G4VAtomDeexcitation.hh>

Inheritance diagram for G4VAtomDeexcitation:

Public Member Functions

 G4VAtomDeexcitation (const G4String &modname="Deexcitation")
 
virtual ~G4VAtomDeexcitation ()
 
void InitialiseAtomicDeexcitation ()
 
virtual void InitialiseForNewRun ()=0
 
virtual void InitialiseForExtraAtom (G4int Z)=0
 
void SetDeexcitationActiveRegion (const G4String &rname, G4bool valDeexcitation, G4bool valAuger, G4bool valPIXE)
 
void SetFluo (G4bool)
 
G4bool IsFluoActive () const
 
void SetAuger (G4bool)
 
G4bool IsAugerActive () const
 
void SetAugerCascade (G4bool)
 
G4bool IsAugerCascadeActive () const
 
void SetPIXE (G4bool)
 
G4bool IsPIXEActive () const
 
const G4StringGetName () const
 
const std::vector< G4bool > & GetListOfActiveAtoms () const
 
void SetVerboseLevel (G4int)
 
G4int GetVerboseLevel () const
 
G4bool CheckDeexcitationActiveRegion (G4int coupleIndex)
 
G4bool CheckAugerActiveRegion (G4int coupleIndex)
 
virtual const G4AtomicShellGetAtomicShell (G4int Z, G4AtomicShellEnumerator shell)=0
 
void GenerateParticles (std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
 
virtual void GenerateParticles (std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4double gammaCut, G4double eCut)=0
 
virtual G4double GetShellIonisationCrossSectionPerAtom (const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=nullptr)=0
 
virtual G4double ComputeShellIonisationCrossSectionPerAtom (const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=nullptr)=0
 
void AlongStepDeexcitation (std::vector< G4Track * > &tracks, const G4Step &step, G4double &eLoss, G4int coupleIndex)
 

Detailed Description

Definition at line 64 of file G4VAtomDeexcitation.hh.

Constructor & Destructor Documentation

G4VAtomDeexcitation::G4VAtomDeexcitation ( const G4String modname = "Deexcitation")
explicit

Definition at line 70 of file G4VAtomDeexcitation.cc.

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 {
76  theParameters = G4EmParameters::Instance();
77  vdyn.reserve(5);
78  theCoupleTable = nullptr;
79  G4String gg = "gammaPIXE";
80  G4String ee = "e-PIXE";
81  if(pixeIDg < 0) { pixeIDg = G4PhysicsModelCatalog::Register(gg); }
82  if(pixeIDe < 0) { pixeIDe = G4PhysicsModelCatalog::Register(ee); }
83  gamma = G4Gamma::Gamma();
84 }
const XML_Char * name
Definition: expat.h:151
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4int Register(const G4String &)
static G4EmParameters * Instance()

Here is the call graph for this function:

G4VAtomDeexcitation::~G4VAtomDeexcitation ( )
virtual

Definition at line 88 of file G4VAtomDeexcitation.cc.

89 {}

Member Function Documentation

void G4VAtomDeexcitation::AlongStepDeexcitation ( std::vector< G4Track * > &  tracks,
const G4Step step,
G4double eLoss,
G4int  coupleIndex 
)

Definition at line 242 of file G4VAtomDeexcitation.cc.

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) {
327  t->SetCreatorModelIndex(pixeIDg);
328  } else {
329  t->SetCreatorModelIndex(pixeIDe);
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 }
G4ParticleDefinition * GetDefinition() const
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
std::vector< G4Element * > G4ElementVector
virtual G4double GetShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=nullptr)=0
G4double GetKineticEnergy() const
G4double GetStepLength() const
G4Material * GetMaterial() const
G4ParticleDefinition * GetDefinition() const
G4double BindingEnergy() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
void SetCreatorModelIndex(G4int idx)
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
const G4ThreeVector & GetPosition() const
G4bool CheckAugerActiveRegion(G4int coupleIndex)
G4StepPoint * GetPostStepPoint() const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double GetGlobalTime() const
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
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)
#define DBL_MAX
Definition: templates.hh:83
G4AtomicShellEnumerator

Here is the call graph for this function:

Here is the caller graph for this function:

G4bool G4VAtomDeexcitation::CheckAugerActiveRegion ( G4int  coupleIndex)
inline

Definition at line 268 of file G4VAtomDeexcitation.hh.

269 {
270  return (activeAugerMedia[coupleIndex]);
271 }

Here is the caller graph for this function:

G4bool G4VAtomDeexcitation::CheckDeexcitationActiveRegion ( G4int  coupleIndex)
inline

Definition at line 262 of file G4VAtomDeexcitation.hh.

263 {
264  return (activeDeexcitationMedia[coupleIndex]);
265 }

Here is the caller graph for this function:

virtual G4double G4VAtomDeexcitation::ComputeShellIonisationCrossSectionPerAtom ( const G4ParticleDefinition ,
G4int  Z,
G4AtomicShellEnumerator  shell,
G4double  kinE,
const G4Material mat = nullptr 
)
pure virtual

Implemented in G4UAtomicDeexcitation.

Here is the caller graph for this function:

void G4VAtomDeexcitation::GenerateParticles ( std::vector< G4DynamicParticle * > *  secVect,
const G4AtomicShell as,
G4int  Z,
G4int  coupleIndex 
)
inline

Definition at line 274 of file G4VAtomDeexcitation.hh.

278 {
279  G4double gCut = DBL_MAX;
280  if(ignoreCuts) {
281  gCut = 0.0;
282  } else if (theCoupleTable) {
283  gCut = (*(theCoupleTable->GetEnergyCutsVector(0)))[idx];
284  }
285  if(gCut < as->BindingEnergy()) {
286  G4double eCut = DBL_MAX;
287  if(CheckAugerActiveRegion(idx)) {
288  if(ignoreCuts) {
289  eCut = 0.0;
290  } else if (theCoupleTable) {
291  eCut = (*(theCoupleTable->GetEnergyCutsVector(1)))[idx];
292  }
293  }
294  GenerateParticles(v, as, Z, gCut, eCut);
295  }
296 }
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
G4bool CheckAugerActiveRegion(G4int coupleIndex)
double G4double
Definition: G4Types.hh:76
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
#define DBL_MAX
Definition: templates.hh:83

Here is the call graph for this function:

Here is the caller graph for this function:

virtual void G4VAtomDeexcitation::GenerateParticles ( std::vector< G4DynamicParticle * > *  secVect,
const G4AtomicShell ,
G4int  Z,
G4double  gammaCut,
G4double  eCut 
)
pure virtual

Implemented in G4UAtomicDeexcitation.

virtual const G4AtomicShell* G4VAtomDeexcitation::GetAtomicShell ( G4int  Z,
G4AtomicShellEnumerator  shell 
)
pure virtual

Implemented in G4UAtomicDeexcitation.

Here is the caller graph for this function:

const std::vector< G4bool > & G4VAtomDeexcitation::GetListOfActiveAtoms ( ) const
inline

Definition at line 246 of file G4VAtomDeexcitation.hh.

247 {
248  return activeZ;
249 }
const G4String & G4VAtomDeexcitation::GetName ( ) const
inline

Definition at line 240 of file G4VAtomDeexcitation.hh.

241 {
242  return name;
243 }
const XML_Char * name
Definition: expat.h:151

Here is the caller graph for this function:

virtual G4double G4VAtomDeexcitation::GetShellIonisationCrossSectionPerAtom ( const G4ParticleDefinition ,
G4int  Z,
G4AtomicShellEnumerator  shell,
G4double  kinE,
const G4Material mat = nullptr 
)
pure virtual

Implemented in G4UAtomicDeexcitation.

Here is the caller graph for this function:

G4int G4VAtomDeexcitation::GetVerboseLevel ( ) const
inline

Definition at line 256 of file G4VAtomDeexcitation.hh.

257 {
258  return verbose;
259 }
void G4VAtomDeexcitation::InitialiseAtomicDeexcitation ( )

Definition at line 93 of file G4VAtomDeexcitation.cc.

94 {
95  theParameters->DefineRegParamForDeex(this);
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
110  if(!isActiveLocked) { isActive = theParameters->Fluo(); }
111  if(!isAugerLocked) { flagAuger = theParameters->Auger(); }
112  if(!isAugerCascadeLocked) { flagAugerCascade = theParameters->AugerCascade(); }
113  if(!isPIXELocked) { flagPIXE = theParameters->Pixe(); }
114  ignoreCuts = theParameters->DeexcitationIgnoreCut();
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) {
123  SetDeexcitationActiveRegion("World",isActive,flagAuger,flagPIXE);
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 =
148  theCoupleTable->GetMaterialCutsCouple(i);
149  if (couple->GetProductionCuts() == rpcuts) {
150  activeDeexcitationMedia[i] = deRegions[j];
151  activeAugerMedia[i] = AugerRegions[j];
152  activePIXEMedia[i] = PIXERegions[j];
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: "
180  << theParameters->PIXECrossSectionModel()
181  << G4endl;
182  G4cout << "### === PIXE model for e+-: "
183  << theParameters->PIXEElectronCrossSectionModel()
184  << G4endl;
185  }
186 }
const XML_Char * name
Definition: expat.h:151
G4ProductionCuts * GetProductionCuts() const
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
G4bool AugerCascade() const
void SetDeexcitationActiveRegion(const G4String &rname, G4bool valDeexcitation, G4bool valAuger, G4bool valPIXE)
G4bool Fluo() const
int G4int
Definition: G4Types.hh:78
G4bool DeexcitationIgnoreCut() const
static const G4double reg
static G4RegionStore * GetInstance()
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:405
const G4String & PIXECrossSectionModel()
static G4ProductionCutsTable * GetProductionCutsTable()
G4bool Auger() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
void DefineRegParamForDeex(G4VAtomDeexcitation *) const
virtual void InitialiseForNewRun()=0
#define G4endl
Definition: G4ios.hh:61
G4bool Pixe() const
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:398
G4ProductionCuts * GetProductionCuts() const
const G4String & PIXEElectronCrossSectionModel()

Here is the call graph for this function:

Here is the caller graph for this function:

virtual void G4VAtomDeexcitation::InitialiseForExtraAtom ( G4int  Z)
pure virtual

Implemented in G4UAtomicDeexcitation.

virtual void G4VAtomDeexcitation::InitialiseForNewRun ( )
pure virtual

Implemented in G4UAtomicDeexcitation.

Here is the caller graph for this function:

G4bool G4VAtomDeexcitation::IsAugerActive ( ) const
inline

Definition at line 215 of file G4VAtomDeexcitation.hh.

216 {
217  return flagAuger;
218 }
G4bool G4VAtomDeexcitation::IsAugerCascadeActive ( ) const
inline

Definition at line 225 of file G4VAtomDeexcitation.hh.

226 {
227  return flagAugerCascade;
228 }

Here is the caller graph for this function:

G4bool G4VAtomDeexcitation::IsFluoActive ( ) const
inline

Definition at line 205 of file G4VAtomDeexcitation.hh.

206 {
207  return isActive;
208 }

Here is the caller graph for this function:

G4bool G4VAtomDeexcitation::IsPIXEActive ( ) const
inline

Definition at line 235 of file G4VAtomDeexcitation.hh.

236 {
237  return flagPIXE;
238 }

Here is the caller graph for this function:

void G4VAtomDeexcitation::SetAuger ( G4bool  val)
inline

Definition at line 210 of file G4VAtomDeexcitation.hh.

211 {
212  if(!isAugerLocked) { flagAuger = val; isAugerLocked = true; }
213 }
void G4VAtomDeexcitation::SetAugerCascade ( G4bool  val)
inline

Definition at line 220 of file G4VAtomDeexcitation.hh.

221 {
222  if(!isAugerCascadeLocked) { flagAugerCascade = val; isAugerCascadeLocked = true; }
223 }
void G4VAtomDeexcitation::SetDeexcitationActiveRegion ( const G4String rname,
G4bool  valDeexcitation,
G4bool  valAuger,
G4bool  valPIXE 
)

Definition at line 191 of file G4VAtomDeexcitation.cc.

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 }
void SetDeexcitationActiveRegion(const G4String &rname, G4bool valDeexcitation, G4bool valAuger, G4bool valPIXE)
int G4int
Definition: G4Types.hh:78
static G4RegionStore * GetInstance()
const G4String & GetName() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4VAtomDeexcitation::SetFluo ( G4bool  val)
inline

Definition at line 200 of file G4VAtomDeexcitation.hh.

201 {
202  if(!isActiveLocked) { isActive = val; isActiveLocked = true; }
203 }
void G4VAtomDeexcitation::SetPIXE ( G4bool  val)
inline

Definition at line 230 of file G4VAtomDeexcitation.hh.

231 {
232  if(!isPIXELocked) { flagPIXE = val; isPIXELocked = true; }
233 }
void G4VAtomDeexcitation::SetVerboseLevel ( G4int  val)
inline

Definition at line 251 of file G4VAtomDeexcitation.hh.

252 {
253  verbose = val;
254 }

The documentation for this class was generated from the following files: