Geant4  10.02.p03
G4VAtomDeexcitation Class Referenceabstract

#include <G4VAtomDeexcitation.hh>

Inheritance diagram for G4VAtomDeexcitation:
Collaboration 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=0)=0
 
virtual G4double ComputeShellIonisationCrossSectionPerAtom (const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0)=0
 
void AlongStepDeexcitation (std::vector< G4Track *> &tracks, const G4Step &step, G4double &eLoss, G4int coupleIndex)
 

Private Member Functions

 G4VAtomDeexcitation (G4VAtomDeexcitation &)
 
G4VAtomDeexcitationoperator= (const G4VAtomDeexcitation &right)
 

Private Attributes

G4EmParameterstheParameters
 
const G4ParticleDefinitiongamma
 
G4ProductionCutsTabletheCoupleTable
 
G4int verbose
 
G4String name
 
G4bool isActive
 
G4bool flagAuger
 
G4bool flagAugerCascade
 
G4bool flagPIXE
 
G4bool ignoreCuts
 
std::vector< G4boolactiveZ
 
std::vector< G4boolactiveDeexcitationMedia
 
std::vector< G4boolactiveAugerMedia
 
std::vector< G4boolactivePIXEMedia
 
std::vector< G4StringactiveRegions
 
std::vector< G4booldeRegions
 
std::vector< G4boolAugerRegions
 
std::vector< G4boolPIXERegions
 
std::vector< G4DynamicParticle * > vdyn
 

Static Private Attributes

static G4int pixeIDg = -1
 
static G4int pixeIDe = -1
 

Detailed Description

Definition at line 64 of file G4VAtomDeexcitation.hh.

Constructor & Destructor Documentation

◆ G4VAtomDeexcitation() [1/2]

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

Definition at line 70 of file G4VAtomDeexcitation.cc.

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 }
std::vector< G4DynamicParticle * > vdyn
const G4ParticleDefinition * gamma
G4EmParameters * theParameters
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4ProductionCutsTable * theCoupleTable
static G4int Register(const G4String &)
static G4EmParameters * Instance()
Here is the call graph for this function:

◆ ~G4VAtomDeexcitation()

G4VAtomDeexcitation::~G4VAtomDeexcitation ( )
virtual

Definition at line 87 of file G4VAtomDeexcitation.cc.

88 {}

◆ G4VAtomDeexcitation() [2/2]

G4VAtomDeexcitation::G4VAtomDeexcitation ( G4VAtomDeexcitation )
private

Member Function Documentation

◆ AlongStepDeexcitation()

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

Definition at line 227 of file G4VAtomDeexcitation.cc.

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);
277  G4double bindingEnergy = shell->BindingEnergy();
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) {
312  t->SetCreatorModelIndex(pixeIDg);
313  } else {
314  t->SetCreatorModelIndex(pixeIDe);
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 }
std::vector< G4Element * > G4ElementVector
std::vector< G4DynamicParticle * > vdyn
virtual G4double GetShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0)=0
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
int G4int
Definition: G4Types.hh:78
string material
Definition: eplot.py:19
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
G4double GetKineticEnergy() const
#define G4UniformRand()
Definition: Randomize.hh:97
Float_t Z
const G4ParticleDefinition * gamma
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:206
TString part[npart]
G4ProductionCutsTable * theCoupleTable
G4bool CheckAugerActiveRegion(G4int coupleIndex)
int G4lrint(double ad)
Definition: templates.hh:163
void GenerateParticles(std::vector< G4DynamicParticle *> *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4ParticleDefinition * GetDefinition() const
std::vector< G4bool > activeZ
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
double G4double
Definition: G4Types.hh:76
G4double BindingEnergy() const
#define DBL_MAX
Definition: templates.hh:83
G4AtomicShellEnumerator
std::vector< G4bool > activePIXEMedia
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CheckAugerActiveRegion()

G4bool G4VAtomDeexcitation::CheckAugerActiveRegion ( G4int  coupleIndex)
inline

Definition at line 267 of file G4VAtomDeexcitation.hh.

268 {
269 
270  return (flagAuger || activeAugerMedia[coupleIndex]);
271 }
std::vector< G4bool > activeAugerMedia
Here is the caller graph for this function:

◆ CheckDeexcitationActiveRegion()

G4bool G4VAtomDeexcitation::CheckDeexcitationActiveRegion ( G4int  coupleIndex)
inline

Definition at line 261 of file G4VAtomDeexcitation.hh.

262 {
263  return (isActive || activeDeexcitationMedia[coupleIndex]);
264 }
std::vector< G4bool > activeDeexcitationMedia
Here is the caller graph for this function:

◆ ComputeShellIonisationCrossSectionPerAtom()

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

Implemented in G4UAtomicDeexcitation.

Here is the caller graph for this function:

◆ GenerateParticles() [1/2]

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
Float_t Z
G4ProductionCutsTable * theCoupleTable
G4bool CheckAugerActiveRegion(G4int coupleIndex)
void GenerateParticles(std::vector< G4DynamicParticle *> *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GenerateParticles() [2/2]

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

Implemented in G4UAtomicDeexcitation.

◆ GetAtomicShell()

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

Implemented in G4UAtomicDeexcitation.

Here is the caller graph for this function:

◆ GetListOfActiveAtoms()

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

Definition at line 245 of file G4VAtomDeexcitation.hh.

246 {
247  return activeZ;
248 }
std::vector< G4bool > activeZ

◆ GetName()

const G4String & G4VAtomDeexcitation::GetName ( void  ) const
inline

Definition at line 239 of file G4VAtomDeexcitation.hh.

240 {
241  return name;
242 }
Here is the caller graph for this function:

◆ GetShellIonisationCrossSectionPerAtom()

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

Implemented in G4UAtomicDeexcitation.

Here is the caller graph for this function:

◆ GetVerboseLevel()

G4int G4VAtomDeexcitation::GetVerboseLevel ( ) const
inline

Definition at line 255 of file G4VAtomDeexcitation.hh.

256 {
257  return verbose;
258 }

◆ InitialiseAtomicDeexcitation()

void G4VAtomDeexcitation::InitialiseAtomicDeexcitation ( )

Definition at line 92 of file G4VAtomDeexcitation.cc.

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 }
std::vector< G4String > activeRegions
G4bool DeexcitationIgnoreCut() const
G4bool Auger() const
G4bool Pixe() const
G4ProductionCuts * GetProductionCuts() const
void SetDeexcitationActiveRegion(const G4String &rname, G4bool valDeexcitation, G4bool valAuger, G4bool valPIXE)
int G4int
Definition: G4Types.hh:78
std::vector< G4bool > activeAugerMedia
std::vector< G4bool > AugerRegions
static const G4double reg
static G4RegionStore * GetInstance()
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
const G4String & PIXECrossSectionModel()
Float_t Z
G4EmParameters * theParameters
G4ProductionCutsTable * theCoupleTable
G4bool AugerCascade() const
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::vector< G4bool > activeDeexcitationMedia
int G4lrint(double ad)
Definition: templates.hh:163
G4bool Fluo() const
virtual void InitialiseForNewRun()=0
std::vector< G4bool > deRegions
#define G4endl
Definition: G4ios.hh:61
std::vector< G4bool > activeZ
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
G4ProductionCuts * GetProductionCuts() const
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
const G4String & PIXEElectronCrossSectionModel()
std::vector< G4bool > PIXERegions
std::vector< G4bool > activePIXEMedia
Here is the call graph for this function:
Here is the caller graph for this function:

◆ InitialiseForExtraAtom()

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

Implemented in G4UAtomicDeexcitation.

◆ InitialiseForNewRun()

virtual void G4VAtomDeexcitation::InitialiseForNewRun ( )
pure virtual

Implemented in G4UAtomicDeexcitation.

Here is the caller graph for this function:

◆ IsAugerActive()

G4bool G4VAtomDeexcitation::IsAugerActive ( ) const
inline

Definition at line 210 of file G4VAtomDeexcitation.hh.

211 {
212  return flagAuger;
213 }
Here is the caller graph for this function:

◆ IsAugerCascadeActive()

G4bool G4VAtomDeexcitation::IsAugerCascadeActive ( ) const
inline

Definition at line 222 of file G4VAtomDeexcitation.hh.

223 {
224  return flagAugerCascade;
225 }
Here is the caller graph for this function:

◆ IsFluoActive()

G4bool G4VAtomDeexcitation::IsFluoActive ( ) const
inline

Definition at line 198 of file G4VAtomDeexcitation.hh.

199 {
200  return isActive;
201 }
Here is the caller graph for this function:

◆ IsPIXEActive()

G4bool G4VAtomDeexcitation::IsPIXEActive ( ) const
inline

Definition at line 234 of file G4VAtomDeexcitation.hh.

235 {
236  return flagPIXE;
237 }
Here is the caller graph for this function:

◆ operator=()

G4VAtomDeexcitation& G4VAtomDeexcitation::operator= ( const G4VAtomDeexcitation right)
private

◆ SetAuger()

void G4VAtomDeexcitation::SetAuger ( G4bool  val)
inline

Definition at line 203 of file G4VAtomDeexcitation.hh.

204 {
205  flagAuger = val;
206  if(val) { isActive = true; }
207  theParameters->SetAuger(val);
208 }
void SetAuger(G4bool val)
G4EmParameters * theParameters
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetAugerCascade()

void G4VAtomDeexcitation::SetAugerCascade ( G4bool  val)
inline

Definition at line 215 of file G4VAtomDeexcitation.hh.

216 {
217  flagAugerCascade = val;
218  if(val) { isActive = true; }
220 }
void SetAugerCascade(G4bool val)
G4EmParameters * theParameters
Here is the call graph for this function:

◆ SetDeexcitationActiveRegion()

void G4VAtomDeexcitation::SetDeexcitationActiveRegion ( const G4String rname,
G4bool  valDeexcitation,
G4bool  valAuger,
G4bool  valPIXE 
)

Definition at line 177 of file G4VAtomDeexcitation.cc.

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 }
std::vector< G4String > activeRegions
const G4String & GetName() const
void SetDeexcitationActiveRegion(const G4String &rname, G4bool valDeexcitation, G4bool valAuger, G4bool valPIXE)
int G4int
Definition: G4Types.hh:78
std::vector< G4bool > AugerRegions
static G4RegionStore * GetInstance()
Char_t n[5]
std::vector< G4bool > deRegions
std::vector< G4bool > PIXERegions
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetFluo()

void G4VAtomDeexcitation::SetFluo ( G4bool  val)
inline

Definition at line 192 of file G4VAtomDeexcitation.hh.

193 {
194  isActive = val;
195  theParameters->SetFluo(val);
196 }
G4EmParameters * theParameters
void SetFluo(G4bool val)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetPIXE()

void G4VAtomDeexcitation::SetPIXE ( G4bool  val)
inline

Definition at line 227 of file G4VAtomDeexcitation.hh.

228 {
229  flagPIXE = val;
230  if(val) { isActive = true; }
231  theParameters->SetPixe(val);
232 }
G4EmParameters * theParameters
void SetPixe(G4bool val)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetVerboseLevel()

void G4VAtomDeexcitation::SetVerboseLevel ( G4int  val)
inline

Definition at line 250 of file G4VAtomDeexcitation.hh.

251 {
252  verbose = val;
253 }
Here is the caller graph for this function:

Member Data Documentation

◆ activeAugerMedia

std::vector<G4bool> G4VAtomDeexcitation::activeAugerMedia
private

Definition at line 180 of file G4VAtomDeexcitation.hh.

◆ activeDeexcitationMedia

std::vector<G4bool> G4VAtomDeexcitation::activeDeexcitationMedia
private

Definition at line 179 of file G4VAtomDeexcitation.hh.

◆ activePIXEMedia

std::vector<G4bool> G4VAtomDeexcitation::activePIXEMedia
private

Definition at line 181 of file G4VAtomDeexcitation.hh.

◆ activeRegions

std::vector<G4String> G4VAtomDeexcitation::activeRegions
private

Definition at line 182 of file G4VAtomDeexcitation.hh.

◆ activeZ

std::vector<G4bool> G4VAtomDeexcitation::activeZ
private

Definition at line 178 of file G4VAtomDeexcitation.hh.

◆ AugerRegions

std::vector<G4bool> G4VAtomDeexcitation::AugerRegions
private

Definition at line 184 of file G4VAtomDeexcitation.hh.

◆ deRegions

std::vector<G4bool> G4VAtomDeexcitation::deRegions
private

Definition at line 183 of file G4VAtomDeexcitation.hh.

◆ flagAuger

G4bool G4VAtomDeexcitation::flagAuger
private

Definition at line 174 of file G4VAtomDeexcitation.hh.

◆ flagAugerCascade

G4bool G4VAtomDeexcitation::flagAugerCascade
private

Definition at line 175 of file G4VAtomDeexcitation.hh.

◆ flagPIXE

G4bool G4VAtomDeexcitation::flagPIXE
private

Definition at line 176 of file G4VAtomDeexcitation.hh.

◆ gamma

const G4ParticleDefinition* G4VAtomDeexcitation::gamma
private

Definition at line 168 of file G4VAtomDeexcitation.hh.

◆ ignoreCuts

G4bool G4VAtomDeexcitation::ignoreCuts
private

Definition at line 177 of file G4VAtomDeexcitation.hh.

◆ isActive

G4bool G4VAtomDeexcitation::isActive
private

Definition at line 173 of file G4VAtomDeexcitation.hh.

◆ name

G4String G4VAtomDeexcitation::name
private

Definition at line 172 of file G4VAtomDeexcitation.hh.

◆ pixeIDe

G4int G4VAtomDeexcitation::pixeIDe = -1
staticprivate

Definition at line 189 of file G4VAtomDeexcitation.hh.

◆ pixeIDg

G4int G4VAtomDeexcitation::pixeIDg = -1
staticprivate

Definition at line 188 of file G4VAtomDeexcitation.hh.

◆ PIXERegions

std::vector<G4bool> G4VAtomDeexcitation::PIXERegions
private

Definition at line 185 of file G4VAtomDeexcitation.hh.

◆ theCoupleTable

G4ProductionCutsTable* G4VAtomDeexcitation::theCoupleTable
private

Definition at line 170 of file G4VAtomDeexcitation.hh.

◆ theParameters

G4EmParameters* G4VAtomDeexcitation::theParameters
private

Definition at line 167 of file G4VAtomDeexcitation.hh.

◆ vdyn

std::vector<G4DynamicParticle*> G4VAtomDeexcitation::vdyn
private

Definition at line 186 of file G4VAtomDeexcitation.hh.

◆ verbose

G4int G4VAtomDeexcitation::verbose
private

Definition at line 171 of file G4VAtomDeexcitation.hh.


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