Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4VAtomDeexcitation Class Referenceabstract

#include <G4VAtomDeexcitation.hh>

Inheritance diagram for G4VAtomDeexcitation:
G4UAtomicDeexcitation

Public Member Functions

 G4VAtomDeexcitation (const G4String &modname="Deexcitation", const G4String &pixename="")
 
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 SetPIXE (G4bool)
 
G4bool IsPIXEActive () const
 
const G4StringGetName () const
 
void SetPIXECrossSectionModel (const G4String &)
 
const G4StringPIXECrossSectionModel () const
 
void SetPIXEElectronCrossSectionModel (const G4String &)
 
const G4StringPIXEElectronCrossSectionModel () 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)
 

Detailed Description

Definition at line 63 of file G4VAtomDeexcitation.hh.

Constructor & Destructor Documentation

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

Definition at line 70 of file G4VAtomDeexcitation.cc.

References G4Gamma::Gamma(), and G4PhysicsModelCatalog::Register().

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 }
const XML_Char * name
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
string pname
Definition: eplot.py:33
static G4int Register(G4String &)
G4VAtomDeexcitation::~G4VAtomDeexcitation ( )
virtual

Definition at line 86 of file G4VAtomDeexcitation.cc.

87 {}

Member Function Documentation

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

Definition at line 208 of file G4VAtomDeexcitation.cc.

References G4InuclSpecialFunctions::bindingEnergy(), G4AtomicShell::BindingEnergy(), CheckAugerActiveRegion(), DBL_MAX, G4lrint(), G4UniformRand, GenerateParticles(), GetAtomicShell(), G4Track::GetDefinition(), G4DynamicParticle::GetDefinition(), G4Material::GetElementVector(), G4ProductionCutsTable::GetEnergyCutsVector(), G4StepPoint::GetGlobalTime(), G4StepPoint::GetKineticEnergy(), G4DynamicParticle::GetKineticEnergy(), G4StepPoint::GetMaterial(), G4Material::GetNumberOfElements(), G4StepPoint::GetPosition(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), GetShellIonisationCrossSectionPerAtom(), G4Step::GetStepLength(), G4Step::GetTrack(), G4Material::GetVecNbOfAtomsPerVolume(), eplot::material, G4INCL::Math::min(), and G4Track::SetCreatorModelIndex().

Referenced by G4VEnergyLossProcess::AlongStepDoIt().

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 }
G4ParticleDefinition * GetDefinition() const
const std::vector< G4double > * GetEnergyCutsVector(size_t pcIdx) const
std::vector< G4Element * > G4ElementVector
G4double GetKineticEnergy() const
G4double GetStepLength() const
G4Material * GetMaterial() const
virtual G4double GetShellIonisationCrossSectionPerAtom(const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0)=0
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
string material
Definition: eplot.py:19
G4StepPoint * GetPreStepPoint() const
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:204
#define G4UniformRand()
Definition: Randomize.hh:87
const G4ThreeVector & GetPosition() const
G4bool CheckAugerActiveRegion(G4int coupleIndex)
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
G4double GetGlobalTime() const
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
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
G4bool G4VAtomDeexcitation::CheckAugerActiveRegion ( G4int  coupleIndex)
inline

Definition at line 282 of file G4VAtomDeexcitation.hh.

Referenced by AlongStepDeexcitation(), and GenerateParticles().

283 {
284 
285  return (flagAuger || activeAugerMedia[coupleIndex]);
286 }
G4bool G4VAtomDeexcitation::CheckDeexcitationActiveRegion ( G4int  coupleIndex)
inline
virtual G4double G4VAtomDeexcitation::ComputeShellIonisationCrossSectionPerAtom ( const G4ParticleDefinition ,
G4int  Z,
G4AtomicShellEnumerator  shell,
G4double  kinE,
const G4Material mat = 0 
)
pure virtual
void G4VAtomDeexcitation::GenerateParticles ( std::vector< G4DynamicParticle * > *  secVect,
const G4AtomicShell as,
G4int  Z,
G4int  coupleIndex 
)
inline

Definition at line 289 of file G4VAtomDeexcitation.hh.

References CheckAugerActiveRegion(), DBL_MAX, and G4ProductionCutsTable::GetEnergyCutsVector().

Referenced by AlongStepDeexcitation(), G4NuclearDecayChannel::DecayIt(), G4LivermorePolarizedPhotoElectricModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4DNABornIonisationModel::SampleSecondaries(), G4LivermoreComptonModifiedModel::SampleSecondaries(), G4LivermorePhotoElectricModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4PenelopePhotoElectricModel::SampleSecondaries(), G4MicroElecInelasticModel::SampleSecondaries(), G4MuElecInelasticModel::SampleSecondaries(), G4PEEffectFluoModel::SampleSecondaries(), G4PenelopeComptonModel::SampleSecondaries(), G4PenelopeIonisationModel::SampleSecondaries(), and G4LowEPComptonModel::SampleSecondaries().

293 {
294  G4double gCut = DBL_MAX;
295  if (theCoupleTable) {
296  gCut = (*(theCoupleTable->GetEnergyCutsVector(0)))[idx];
297  }
298  if(gCut < as->BindingEnergy()) {
299  G4double eCut = DBL_MAX;
300  if(CheckAugerActiveRegion(idx)) {
301  if (theCoupleTable) {
302  eCut = (*(theCoupleTable->GetEnergyCutsVector(1)))[idx];
303  }
304  }
305  GenerateParticles(v, as, Z, gCut, eCut);
306  }
307 }
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
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
const std::vector< G4bool > & G4VAtomDeexcitation::GetListOfActiveAtoms ( ) const
inline

Definition at line 260 of file G4VAtomDeexcitation.hh.

261 {
262  return activeZ;
263 }
const G4String & G4VAtomDeexcitation::GetName ( void  ) const
inline

Definition at line 230 of file G4VAtomDeexcitation.hh.

231 {
232  return name;
233 }
const XML_Char * name
virtual G4double G4VAtomDeexcitation::GetShellIonisationCrossSectionPerAtom ( const G4ParticleDefinition ,
G4int  Z,
G4AtomicShellEnumerator  shell,
G4double  kinE,
const G4Material mat = 0 
)
pure virtual
G4int G4VAtomDeexcitation::GetVerboseLevel ( ) const
inline

Definition at line 270 of file G4VAtomDeexcitation.hh.

271 {
272  return verbose;
273 }
void G4VAtomDeexcitation::InitialiseAtomicDeexcitation ( )

Definition at line 91 of file G4VAtomDeexcitation.cc.

References G4cout, G4endl, G4lrint(), G4Material::GetElementVector(), G4RegionStore::GetInstance(), G4MaterialCutsCouple::GetMaterial(), G4ProductionCutsTable::GetMaterialCutsCouple(), G4Material::GetNumberOfElements(), G4MaterialCutsCouple::GetProductionCuts(), G4Region::GetProductionCuts(), G4ProductionCutsTable::GetProductionCutsTable(), G4RegionStore::GetRegion(), G4ProductionCutsTable::GetTableSize(), InitialiseForNewRun(), IsPIXEActive(), and SetDeexcitationActiveRegion().

Referenced by G4LossTableManager::BuildPhysicsTable(), G4RadioactiveDecay::BuildPhysicsTable(), and G4LossTableManager::LocalPhysicsTables().

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 }
G4ProductionCuts * GetProductionCuts() const
G4Region * GetRegion(const G4String &name, G4bool verbose=true) const
std::vector< G4Element * > G4ElementVector
G4bool IsPIXEActive() const
void SetDeexcitationActiveRegion(const G4String &rname, G4bool valDeexcitation, G4bool valAuger, G4bool valPIXE)
const XML_Char * name
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
static G4RegionStore * GetInstance()
G4GLOB_DLL std::ostream G4cout
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
int G4lrint(double ad)
Definition: templates.hh:163
virtual void InitialiseForNewRun()=0
#define G4endl
Definition: G4ios.hh:61
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
G4ProductionCuts * GetProductionCuts() const
const G4Material * GetMaterial() const
virtual void G4VAtomDeexcitation::InitialiseForExtraAtom ( G4int  Z)
pure virtual

Implemented in G4UAtomicDeexcitation.

virtual void G4VAtomDeexcitation::InitialiseForNewRun ( )
pure virtual

Implemented in G4UAtomicDeexcitation.

Referenced by InitialiseAtomicDeexcitation().

G4bool G4VAtomDeexcitation::IsAugerActive ( ) const
inline

Definition at line 214 of file G4VAtomDeexcitation.hh.

215 {
216  return flagAuger;
217 }
G4bool G4VAtomDeexcitation::IsFluoActive ( ) const
inline
G4bool G4VAtomDeexcitation::IsPIXEActive ( ) const
inline
const G4String & G4VAtomDeexcitation::PIXECrossSectionModel ( ) const
inline

Definition at line 248 of file G4VAtomDeexcitation.hh.

Referenced by G4UAtomicDeexcitation::InitialiseForNewRun().

249 {
250  return namePIXE;
251 }
const G4String & G4VAtomDeexcitation::PIXEElectronCrossSectionModel ( ) const
inline

Definition at line 254 of file G4VAtomDeexcitation.hh.

Referenced by G4UAtomicDeexcitation::InitialiseForNewRun().

255 {
256  return nameElectronPIXE;
257 }
void G4VAtomDeexcitation::SetAuger ( G4bool  val)
inline

Definition at line 208 of file G4VAtomDeexcitation.hh.

Referenced by G4EmProcessOptions::SetAuger().

209 {
210  flagAuger = val;
211  if(val) { isActive = true; }
212 }
void G4VAtomDeexcitation::SetDeexcitationActiveRegion ( const G4String rname,
G4bool  valDeexcitation,
G4bool  valAuger,
G4bool  valPIXE 
)

Definition at line 173 of file G4VAtomDeexcitation.cc.

References n.

Referenced by InitialiseAtomicDeexcitation(), and G4EmProcessOptions::SetDeexcitationActiveRegion().

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 }
const G4int n
void G4VAtomDeexcitation::SetFluo ( G4bool  val)
inline
void G4VAtomDeexcitation::SetPIXE ( G4bool  val)
inline

Definition at line 219 of file G4VAtomDeexcitation.hh.

Referenced by G4EmProcessOptions::SetPIXE().

220 {
221  flagPIXE = val;
222  if(val) { isActive = true; }
223 }
void G4VAtomDeexcitation::SetPIXECrossSectionModel ( const G4String n)
inline

Definition at line 236 of file G4VAtomDeexcitation.hh.

References n.

Referenced by G4UAtomicDeexcitation::InitialiseForNewRun(), and G4EmProcessOptions::SetPIXECrossSectionModel().

237 {
238  namePIXE = n;
239 }
const G4int n
void G4VAtomDeexcitation::SetPIXEElectronCrossSectionModel ( const G4String n)
inline

Definition at line 242 of file G4VAtomDeexcitation.hh.

References n.

Referenced by G4UAtomicDeexcitation::InitialiseForNewRun(), and G4EmProcessOptions::SetPIXEElectronCrossSectionModel().

243 {
244  nameElectronPIXE = n;
245 }
const G4int n
void G4VAtomDeexcitation::SetVerboseLevel ( G4int  val)
inline

Definition at line 265 of file G4VAtomDeexcitation.hh.

Referenced by G4LossTableManager::SetVerbose().

266 {
267  verbose = val;
268 }

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