Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Attributes | Static Protected Attributes
G4NuclearDecayChannel Class Reference

#include <G4NuclearDecayChannel.hh>

Inheritance diagram for G4NuclearDecayChannel:
G4GeneralPhaseSpaceDecay G4VDecayChannel G4AlphaDecayChannel G4BetaMinusDecayChannel G4BetaPlusDecayChannel G4ITDecayChannel G4KshellECDecayChannel G4LshellECDecayChannel G4MshellECDecayChannel

Public Member Functions

 G4NuclearDecayChannel (const G4RadioactiveDecayMode &theMode, G4int Verbose)
 
 G4NuclearDecayChannel (const G4RadioactiveDecayMode &theMode, G4int Verbose, const G4ParticleDefinition *theParentNucleus, G4double theBR, G4double theQtransition, G4int A, G4int Z, G4double theDaughterExcitation)
 
 G4NuclearDecayChannel (const G4RadioactiveDecayMode &theMode, G4int Verbose, const G4ParticleDefinition *theParentNucleus, G4double theBR, G4double theQtransition, G4int A, G4int Z, G4double theDaughterExcitation, const G4String theDaughterName1)
 
 G4NuclearDecayChannel (const G4RadioactiveDecayMode &theMode, G4int Verbose, const G4ParticleDefinition *theParentNucleus, G4double theBR, G4double theFFN, G4bool betaS, G4RandGeneral *randBeta, G4double theQtransition, G4int A, G4int Z, G4double theDaughterExcitation, const G4String theDaughterName1, const G4String theDaughterName2)
 
virtual ~G4NuclearDecayChannel ()
 
G4DecayProductsDecayIt (G4double)
 
void SetHLThreshold (G4double hl)
 
void SetICM (G4bool icm)
 
void SetARM (G4bool arm)
 
G4RadioactiveDecayMode GetDecayMode ()
 
G4double GetDaughterExcitation ()
 
G4ParticleDefinitionGetDaughterNucleus ()
 
- Public Member Functions inherited from G4GeneralPhaseSpaceDecay
 G4GeneralPhaseSpaceDecay (G4int Verbose=1)
 
 G4GeneralPhaseSpaceDecay (const G4String &theParentName, G4double theBR, G4int theNumberOfDaughters, const G4String &theDaughterName1, const G4String &theDaughterName2="", const G4String &theDaughterName3="")
 
 G4GeneralPhaseSpaceDecay (const G4String &theParentName, G4double theParentMass, G4double theBR, G4int theNumberOfDaughters, const G4String &theDaughterName1, const G4String &theDaughterName2="", const G4String &theDaughterName3="")
 
 G4GeneralPhaseSpaceDecay (const G4String &theParentName, G4double theParentMass, G4double theBR, G4int theNumberOfDaughters, const G4String &theDaughterName1, const G4String &theDaughterName2, const G4String &theDaughterName3, const G4double *masses)
 
virtual ~G4GeneralPhaseSpaceDecay ()
 
G4double GetParentMass () const
 
void SetParentMass (const G4double aParentMass)
 
- Public Member Functions inherited from G4VDecayChannel
 G4VDecayChannel (const G4String &aName, G4int Verbose=1)
 
 G4VDecayChannel (const G4String &aName, const G4String &theParentName, G4double theBR, G4int theNumberOfDaughters, const G4String &theDaughterName1, const G4String &theDaughterName2="", const G4String &theDaughterName3="", const G4String &theDaughterName4="")
 
virtual ~G4VDecayChannel ()
 
G4int operator== (const G4VDecayChannel &right) const
 
G4int operator!= (const G4VDecayChannel &right) const
 
G4int operator< (const G4VDecayChannel &right) const
 
const G4StringGetKinematicsName () const
 
G4double GetBR () const
 
G4int GetNumberOfDaughters () const
 
G4ParticleDefinitionGetParent ()
 
G4ParticleDefinitionGetDaughter (G4int anIndex)
 
G4int GetAngularMomentum ()
 
const G4StringGetParentName () const
 
const G4StringGetDaughterName (G4int anIndex) const
 
G4double GetParentMass () const
 
G4double GetDaughterMass (G4int anIndex) const
 
void SetParent (const G4ParticleDefinition *particle_type)
 
void SetParent (const G4String &particle_name)
 
void SetBR (G4double value)
 
void SetNumberOfDaughters (G4int value)
 
void SetDaughter (G4int anIndex, const G4ParticleDefinition *particle_type)
 
void SetDaughter (G4int anIndex, const G4String &particle_name)
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
void DumpInfo ()
 

Protected Attributes

const G4RadioactiveDecayMode decayMode
 
G4double daughterExcitation
 
G4int daughterA
 
G4int daughterZ
 
G4ParticleDefinitiondaughterNucleus
 
const G4double Qtransition
 
G4double halflifethreshold
 
G4bool applyICM
 
G4bool applyARM
 
G4RandGeneralRandomEnergy
 
- Protected Attributes inherited from G4VDecayChannel
G4String kinematics_name
 
G4double rbranch
 
G4int numberOfDaughters
 
G4Stringparent_name
 
G4String ** daughters_name
 
G4ParticleTableparticletable
 
G4int verboseLevel
 
G4ParticleDefinitionG4MT_parent
 
G4ParticleDefinition ** G4MT_daughters
 
G4double G4MT_parent_mass
 
G4doubleG4MT_daughters_mass
 

Static Protected Attributes

static const G4double pTolerance = 0.001
 
static const G4double levelTolerance = 2.0*keV
 
static G4ThreadLocal
G4DynamicParticle
dynamicDaughter = 0
 
- Static Protected Attributes inherited from G4VDecayChannel
static const G4String noName = " "
 

Additional Inherited Members

- Static Public Member Functions inherited from G4GeneralPhaseSpaceDecay
static G4double Pmx (G4double e, G4double p1, G4double p2)
 
- Protected Member Functions inherited from G4GeneralPhaseSpaceDecay
G4DecayProductsOneBodyDecayIt ()
 
G4DecayProductsTwoBodyDecayIt ()
 
G4DecayProductsThreeBodyDecayIt ()
 
G4DecayProductsManyBodyDecayIt ()
 
- Protected Member Functions inherited from G4VDecayChannel
void ClearDaughtersName ()
 
void FillDaughters ()
 
void FillParent ()
 
 G4VDecayChannel ()
 
 G4VDecayChannel (const G4VDecayChannel &)
 
G4VDecayChanneloperator= (const G4VDecayChannel &)
 

Detailed Description

Definition at line 41 of file G4NuclearDecayChannel.hh.

Constructor & Destructor Documentation

G4NuclearDecayChannel::G4NuclearDecayChannel ( const G4RadioactiveDecayMode theMode,
G4int  Verbose 
)
inline

Definition at line 54 of file G4NuclearDecayChannel.hh.

55  : G4GeneralPhaseSpaceDecay(Verbose), decayMode(theMode),Qtransition(0) {;}
const G4RadioactiveDecayMode decayMode
G4NuclearDecayChannel::G4NuclearDecayChannel ( const G4RadioactiveDecayMode theMode,
G4int  Verbose,
const G4ParticleDefinition theParentNucleus,
G4double  theBR,
G4double  theQtransition,
G4int  A,
G4int  Z,
G4double  theDaughterExcitation 
)

Definition at line 98 of file G4NuclearDecayChannel.cc.

References applyARM, applyICM, G4VDecayChannel::FillDaughters(), G4VDecayChannel::FillParent(), G4cout, G4endl, G4VDecayChannel::G4MT_parent_mass, G4ParticleDefinition::GetPDGMass(), G4VDecayChannel::GetVerboseLevel(), halflifethreshold, python.hepunit::nanosecond, G4VDecayChannel::SetBR(), G4VDecayChannel::SetNumberOfDaughters(), and G4VDecayChannel::SetParent().

106  :G4GeneralPhaseSpaceDecay(Verbose), decayMode(theMode),
107  Qtransition(theQtransition), RandomEnergy(0)
108 {
109 #ifdef G4VERBOSE
110  if (GetVerboseLevel() > 1) {
111  G4cout << "G4NuclearDecayChannel constructor for " << G4int(theMode)
112  << G4endl;
113  }
114 #endif
115  SetParent(theParentNucleus);
116  FillParent();
117  G4MT_parent_mass = theParentNucleus->GetPDGMass();
118  SetBR(theBR);
120  FillDaughterNucleus(0, A, Z, theDaughterExcitation);
122  applyICM = true;
123  applyARM = true;
124  FillDaughters();
125 }
void SetBR(G4double value)
G4double G4MT_parent_mass
const G4RadioactiveDecayMode decayMode
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
void SetNumberOfDaughters(G4int value)
int nanosecond
Definition: hepunit.py:82
G4int GetVerboseLevel() const
G4double GetPDGMass() const
void SetParent(const G4ParticleDefinition *particle_type)
#define G4endl
Definition: G4ios.hh:61
G4NuclearDecayChannel::G4NuclearDecayChannel ( const G4RadioactiveDecayMode theMode,
G4int  Verbose,
const G4ParticleDefinition theParentNucleus,
G4double  theBR,
G4double  theQtransition,
G4int  A,
G4int  Z,
G4double  theDaughterExcitation,
const G4String  theDaughterName1 
)

Definition at line 130 of file G4NuclearDecayChannel.cc.

References applyARM, applyICM, G4VDecayChannel::FillDaughters(), G4VDecayChannel::FillParent(), G4cout, G4endl, G4VDecayChannel::G4MT_parent_mass, G4ParticleDefinition::GetPDGMass(), G4VDecayChannel::GetVerboseLevel(), halflifethreshold, python.hepunit::nanosecond, G4VDecayChannel::SetBR(), G4VDecayChannel::SetDaughter(), G4VDecayChannel::SetNumberOfDaughters(), and G4VDecayChannel::SetParent().

139  :G4GeneralPhaseSpaceDecay(Verbose), decayMode(theMode),
140  Qtransition(theQtransition), RandomEnergy(0)
141 {
142 #ifdef G4VERBOSE
143  if (GetVerboseLevel() > 1) {
144  G4cout << "G4NuclearDecayChannel constructor for " << G4int(theMode)
145  << G4endl;
146  }
147 #endif
148  SetParent(theParentNucleus);
149  FillParent();
150  G4MT_parent_mass = theParentNucleus->GetPDGMass();
151  SetBR(theBR);
153  SetDaughter(0, theDaughterName1);
154  FillDaughterNucleus(1, A, Z, theDaughterExcitation);
156  applyICM = true;
157  applyARM = true;
158  FillDaughters();
159 }
void SetBR(G4double value)
G4double G4MT_parent_mass
const G4RadioactiveDecayMode decayMode
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
void SetNumberOfDaughters(G4int value)
int nanosecond
Definition: hepunit.py:82
G4int GetVerboseLevel() const
G4double GetPDGMass() const
void SetParent(const G4ParticleDefinition *particle_type)
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)
#define G4endl
Definition: G4ios.hh:61
G4NuclearDecayChannel::G4NuclearDecayChannel ( const G4RadioactiveDecayMode theMode,
G4int  Verbose,
const G4ParticleDefinition theParentNucleus,
G4double  theBR,
G4double  theFFN,
G4bool  betaS,
G4RandGeneral randBeta,
G4double  theQtransition,
G4int  A,
G4int  Z,
G4double  theDaughterExcitation,
const G4String  theDaughterName1,
const G4String  theDaughterName2 
)

Definition at line 164 of file G4NuclearDecayChannel.cc.

References applyARM, applyICM, G4VDecayChannel::FillDaughters(), G4VDecayChannel::FillParent(), G4cout, G4endl, G4VDecayChannel::G4MT_parent_mass, G4ParticleDefinition::GetPDGMass(), G4VDecayChannel::GetVerboseLevel(), halflifethreshold, python.hepunit::nanosecond, RandomEnergy, G4VDecayChannel::SetBR(), G4VDecayChannel::SetDaughter(), G4VDecayChannel::SetNumberOfDaughters(), and G4VDecayChannel::SetParent().

177  :G4GeneralPhaseSpaceDecay(Verbose), decayMode(theMode),
178  Qtransition(theQtransition)
179 {
180 #ifdef G4VERBOSE
181  if (GetVerboseLevel() > 1) {
182  G4cout << "G4NuclearDecayChannel constructor for " << G4int(theMode)
183  << G4endl;
184  }
185 #endif
186  SetParent(theParentNucleus);
187  FillParent();
188  G4MT_parent_mass = theParentNucleus->GetPDGMass();
189  SetBR (theBR);
191  SetDaughter(0, theDaughterName1);
192  SetDaughter(2, theDaughterName2);
193  FillDaughterNucleus(1, A, Z, theDaughterExcitation);
194  RandomEnergy = randBeta;
195 
197  applyICM = true;
198  applyARM = true;
199  FillDaughters();
200 }
void SetBR(G4double value)
G4double G4MT_parent_mass
const G4RadioactiveDecayMode decayMode
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
void SetNumberOfDaughters(G4int value)
int nanosecond
Definition: hepunit.py:82
G4int GetVerboseLevel() const
G4double GetPDGMass() const
void SetParent(const G4ParticleDefinition *particle_type)
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)
#define G4endl
Definition: G4ios.hh:61
G4NuclearDecayChannel::~G4NuclearDecayChannel ( )
virtual

Definition at line 202 of file G4NuclearDecayChannel.cc.

203 {}

Member Function Documentation

G4DecayProducts * G4NuclearDecayChannel::DecayIt ( G4double  )
virtual

Reimplemented from G4GeneralPhaseSpaceDecay.

Definition at line 238 of file G4NuclearDecayChannel.cc.

References applyARM, applyICM, G4LossTableManager::AtomDeexcitation(), CLHEP::HepLorentzVector::boost(), CLHEP::HepLorentzVector::boostVector(), G4PhotonEvaporation::BreakUp(), daughterA, daughterExcitation, daughterZ, decayMode, G4VDecayChannel::DumpInfo(), dynamicDaughter, CLHEP::HepLorentzVector::e(), FatalException, G4endl, G4Exception(), G4VDecayChannel::G4MT_parent, G4UniformRand, G4VAtomDeexcitation::GenerateParticles(), G4DynamicParticle::Get4Momentum(), G4VAtomDeexcitation::GetAtomicShell(), G4IonTable::GetIon(), G4ParticleTable::GetIonTable(), G4AtomicShells::GetNumberOfShells(), G4ParticleTable::GetParticleTable(), G4ParticleDefinition::GetPDGMass(), G4PhotonEvaporation::GetVacantShellNumber(), G4VDecayChannel::GetVerboseLevel(), G4LossTableManager::Instance(), G4VAtomDeexcitation::IsFluoActive(), JustWarning, python.hepunit::keV, KshellEC, LshellEC, MshellEC, G4VDecayChannel::numberOfDaughters, G4GeneralPhaseSpaceDecay::OneBodyDecayIt(), G4VDecayChannel::parent_name, G4DecayProducts::PopProducts(), G4DecayProducts::PushProducts(), RDM_ERROR, G4PhotonEvaporation::RDMForced(), G4DynamicParticle::Set4Momentum(), CLHEP::HepLorentzVector::setE(), G4PhotonEvaporation::SetICM(), G4GeneralPhaseSpaceDecay::SetParentMass(), G4PhotonEvaporation::SetVerboseLevel(), and G4GeneralPhaseSpaceDecay::TwoBodyDecayIt().

239 {
240  G4double deltaM;
241  if (decayMode == 1) { // beta- decay
242  deltaM = CLHEP::electron_mass_c2;
243  } else if (decayMode == 2) { // beta+ decay
244  deltaM = 2.*CLHEP::electron_mass_c2;
245  } else if (decayMode < 6 && decayMode > 2) { // EC
246  deltaM = -CLHEP::electron_mass_c2;
247  } else { // all others
248  deltaM = 0.0;
249  }
250 
251  // Mass available for decay in rest frame of parent after correcting for
252  // the appropriate number of electron masses and reserving the daughter
253  // excitation energy to be applied later
254  G4double massOfParent = G4MT_parent->GetPDGMass(); // PDG mass includes excitation energy
255  SetParentMass(massOfParent - deltaM - daughterExcitation);
256 
257  // Define a product vector.
258  G4DecayProducts* products = 0;
259 
260  // Depending upon the number of daughters, select the appropriate decay
261  // kinematics scheme.
262  switch (numberOfDaughters) {
263  case 0:
264  {
266  ed << " No daughters defined " << G4endl;
267  G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_005",
268  JustWarning, ed);
269  }
270  break;
271  case 1:
272  products = OneBodyDecayIt();
273  break;
274  case 2:
275  products = TwoBodyDecayIt();
276  break;
277  case 3:
278  products = BetaDecayIt();
279  break;
280  default:
281  {
283  ed << " More than 3 daughters in decay: N = " << numberOfDaughters
284  << G4endl;
285  G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_007",
286  FatalException, ed);
287  }
288  }
289  if (products == 0) {
291  ed << " Parent nucleus " << *parent_name << " was not decayed " << G4endl;
292  G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_008",
293  JustWarning, ed);
294  DumpInfo();
295  } else {
296  // If the decay is to an excited state of the daughter nuclide, the photon
297  // evaporation process must be applied.
298 
299  // Need to hold the shell idex after ICM
300  G4int shellIndex = -1;
301  if (daughterExcitation > 0.0) {
302  // Pop the daughter nucleus off the product vector to get its 4-momentum
303  dynamicDaughter = products->PopProducts();
304  G4LorentzVector daughterMomentum = dynamicDaughter->Get4Momentum();
305  if (dynamicDaughter) delete dynamicDaughter;
306 
307  // Using daughter nucleus, set up a G4Fragment for photon evaporatation
308  if (decayMode == 0) {
309  G4double exe = ((const G4Ions*)(G4MT_parent))->GetExcitationEnergy();
310  daughterMomentum.setE(daughterMomentum.e() + exe);
311  }
312  G4Fragment nucleus(daughterA, daughterZ, daughterMomentum);
313  G4PhotonEvaporation* deexcitation = new G4PhotonEvaporation;
314  deexcitation->SetVerboseLevel(GetVerboseLevel());
315 
316  // switch on/off internal electron conversion
317  deexcitation->SetICM(applyICM);
318 
319  // In IT mode, we need to force the transition
320  if (decayMode == 0) {
321  deexcitation->RDMForced(true);
322  // at this point, de-excitation will occur even if IT state is long-lived (>1ns)
323  // Why does it need to be forced?
324  } else {
325  // Not forced, but decay will still happen if lifetime < 1 ns,
326  // otherwise no gamma decay is performed
327  deexcitation->RDMForced(false);
328  }
329 
330  //////////////////////////////////////////////////////////////////////////
331  // //
332  // Apply photon evaporation if Isomeric Transition is indicated. //
333  // Use G4PhotonEvaporation::BreakUp() which does only one transition. //
334  // This allows IC to be done. //
335  // //
336  //////////////////////////////////////////////////////////////////////////
337 
338  G4IonTable* theIonTable =
340  G4ParticleDefinition* daughterIon = 0;
341 
342  if (decayMode != 0) {
343  daughterIon = theIonTable->GetIon(daughterZ, daughterA, daughterExcitation);
344  } else {
345  // The fragment vector from photon evaporation contains the list of
346  // evaporated gammas, some of which may have been replaced by conversion
347  // electrons. The last element is the residual nucleus.
348  G4FragmentVector* gammas = deexcitation->BreakUp(nucleus); // Evaporate only one photon
349  G4int nGammas = gammas->size() - 1;
350  G4double finalDaughterExcitation =
351  gammas->operator[](nGammas)->GetExcitationEnergy();
352  // Go through each gamma/e- and add it to the decay product. The
353  // angular distribution of the gammas is isotropic, and the residual
354  // nucleus is assumed not to have suffered any recoil as a result of
355  // this de-excitation.
356  G4double Egamma = 0.0;
357  for (G4int ig = 0; ig < nGammas; ig++) {
358  G4DynamicParticle* theGammaRay =
359  new G4DynamicParticle(gammas->operator[](ig)->GetParticleDefinition(),
360  gammas->operator[](ig)->GetMomentum());
361  theGammaRay -> SetProperTime(gammas->operator[](ig)->GetCreationTime());
362  products->PushProducts (theGammaRay);
363  Egamma = (gammas->operator[](ig)->GetMomentum()).e();
364  }
365  if (finalDaughterExcitation <= 1.0*keV) finalDaughterExcitation = 0;
366 
367  // Get new ion with excitation energy reduced by emitted gamma energy
368  daughterIon =
369  theIonTable->GetIon(daughterZ, daughterA, finalDaughterExcitation);
370  daughterMomentum.setE(daughterMomentum.e() - Egamma);
371 
372  // Delete/reset variables associated with the gammas.
373  while (!gammas->empty() ) {
374  delete *(gammas->end()-1);
375  gammas->pop_back();
376  }
377  delete gammas;
378 
379  } // end if decayMode == 0
380 
381  G4ThreeVector const daughterMomentum1(static_cast<const G4LorentzVector> (daughterMomentum));
382  dynamicDaughter = new G4DynamicParticle(daughterIon, daughterMomentum1);
383  products->PushProducts(dynamicDaughter);
384 
385  // retrieve the ICM shell index
386  shellIndex = deexcitation->GetVacantShellNumber();
387 
388  delete deexcitation;
389  } // if daughter excitation > 0
390 
391  // Now take care of the EC products which have to go through the ARM
392  G4int eShell = -1;
393  if (decayMode == 3 || decayMode == 4 || decayMode == 5) {
394  switch (decayMode)
395  {
396  case KshellEC:
397  {
398  eShell = 0; // --> 0 from 1 (f.lei 30/4/2008)
399  }
400  break;
401  case LshellEC:
402  {
403  eShell = G4int(G4UniformRand()*3)+1;
404  }
405  break;
406  case MshellEC:
407  {
408  // limit the shell index to 6 as specified by the ARM (F.Lei 06/05/2010)
409  // eShell = G4int(G4UniformRand()*5)+4;
410  eShell = G4int(G4UniformRand()*3)+4;
411  }
412  break;
413  case RDM_ERROR:
414  default:
415  G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_009",
416  FatalException, "Incorrect decay mode selection");
417  }
418  } else {
419  // For other cases eShell comes from shellIndex resulting from the photo decay
420  // modeled by G4PhotonEvaporation* de-excitation (see above)
421  eShell = shellIndex;
422  }
423 
424  // now apply ARM if it is requested and there is a vaccancy
425  if (applyARM && eShell != -1) {
426  G4int aZ = daughterZ;
427 
428  // V.Ivanchenko migration to new interface to atomic deexcitation
429  // no check on index of G4MaterialCutsCouple, simplified
430  // check on secondary energy Esec < 0.1 keV
431  G4VAtomDeexcitation* atomDeex =
433  if (atomDeex) {
434  if(atomDeex->IsFluoActive() && aZ > 5 && aZ < 100) { // only applies to 5< Z <100
435  if (eShell >= G4AtomicShells::GetNumberOfShells(aZ)){
436  eShell = G4AtomicShells::GetNumberOfShells(aZ)-1;
437  }
439  const G4AtomicShell* shell = atomDeex->GetAtomicShell(aZ, as);
440  std::vector<G4DynamicParticle*> armProducts;
441  const G4double deexLimit = 0.1*keV;
442  atomDeex->GenerateParticles(&armProducts, shell, aZ, deexLimit, deexLimit);
443  size_t narm = armProducts.size();
444  if (narm > 0) {
445  // L.Desorgher: need to initialize dynamicDaughter in some decay
446  // cases (for example Hg194)
447  dynamicDaughter = products->PopProducts();
449  for (size_t i = 0; i<narm; ++i) {
450  G4DynamicParticle* dp = armProducts[i];
451  G4LorentzVector lv = dp->Get4Momentum().boost(bst);
452  dp->Set4Momentum(lv);
453  products->PushProducts(dp);
454  }
455  products->PushProducts(dynamicDaughter);
456  }
457  }
458  }
459  }
460  } // Parent nucleus decayed
461  /*
462 
463  if (atomDeex && aZ > 5 && aZ < 100) { // only applies to 5< Z <100
464  // Retrieve the corresponding identifier and binding energy of the selected shell
465  const G4AtomicTransitionManager* transitionManager = G4AtomicTransitionManager::Instance();
466 
467  //check that eShell is smaller than the max number of shells
468  //this to avoid a warning from transitionManager(otherwise the output is the same)
469  //Correction to Bug 1662. L Desorgher (04/02/2011)
470  if (eShell >= transitionManager->NumberOfShells(aZ)){
471  eShell=transitionManager->NumberOfShells(aZ)-1;
472  }
473 
474  const G4AtomicShell* shell = transitionManager->Shell(aZ, eShell);
475  G4double bindingEnergy = shell->BindingEnergy();
476  G4int shellId = shell->ShellId();
477 
478  G4AtomicDeexcitation* atomDeex = new G4AtomicDeexcitation();
479  //the default is no Auger electron generation.
480  // Switch it on/off here!
481  atomDeex->ActivateAugerElectronProduction(true);
482  std::vector<G4DynamicParticle*>* armProducts = atomDeex->GenerateParticles(aZ,shellId);
483 
484  // pop up the daughter before insertion;
485  // f.lei (30/04/2008) check if the total kinetic energy is less than
486  // the shell binding energy; if true add the difference to the daughter to conserve the energy
487  dynamicDaughter = products->PopProducts();
488  G4double tARMEnergy = 0.0;
489  for (size_t i = 0; i < armProducts->size(); i++) {
490  products->PushProducts ((*armProducts)[i]);
491  tARMEnergy += (*armProducts)[i]->GetKineticEnergy();
492  }
493  if ((bindingEnergy - tARMEnergy) > 0.1*keV){
494  G4double dEnergy = dynamicDaughter->GetKineticEnergy() + (bindingEnergy - tARMEnergy);
495  dynamicDaughter->SetKineticEnergy(dEnergy);
496  }
497  products->PushProducts(dynamicDaughter);
498 
499 #ifdef G4VERBOSE
500  if (GetVerboseLevel()>0)
501  {
502  G4cout <<"G4NuclearDecayChannel::Selected shell number for ARM = " <<shellId <<G4endl;
503  G4cout <<"G4NuclearDecayChannel::ARM products = " <<armProducts->size()<<G4endl;
504  G4cout <<" The binding energy = " << bindingEnergy << G4endl;
505  G4cout <<" Total ARM particle kinetic energy = " << tARMEnergy << G4endl;
506  }
507 #endif
508 
509  delete armProducts;
510  delete atomDeex;
511  }
512  }
513  */
514 
515  return products;
516 }
Hep3Vector boostVector() const
static G4LossTableManager * Instance()
G4bool IsFluoActive() const
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4int PushProducts(G4DynamicParticle *aParticle)
G4ParticleDefinition * G4MT_parent
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:449
const G4RadioactiveDecayMode decayMode
int G4int
Definition: G4Types.hh:78
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
G4IonTable * GetIonTable() const
#define G4UniformRand()
Definition: Randomize.hh:87
Definition: G4Ions.hh:51
HepLorentzVector & boost(double, double, double)
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:65
virtual G4FragmentVector * BreakUp(const G4Fragment &nucleus)
G4String * parent_name
G4LorentzVector Get4Momentum() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void Set4Momentum(const G4LorentzVector &momentum)
G4int GetVerboseLevel() const
G4double GetPDGMass() const
static G4ParticleTable * GetParticleTable()
G4DynamicParticle * PopProducts()
#define G4endl
Definition: G4ios.hh:61
static G4ThreadLocal G4DynamicParticle * dynamicDaughter
G4VAtomDeexcitation * AtomDeexcitation()
double G4double
Definition: G4Types.hh:76
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void SetParentMass(const G4double aParentMass)
G4AtomicShellEnumerator
void SetVerboseLevel(G4int verbose)
static G4int GetNumberOfShells(G4int Z)
G4double G4NuclearDecayChannel::GetDaughterExcitation ( )
inline

Definition at line 98 of file G4NuclearDecayChannel.hh.

References daughterExcitation.

Referenced by G4RadioactiveDecay::AddDecayRateTable().

G4ParticleDefinition* G4NuclearDecayChannel::GetDaughterNucleus ( )
inline

Definition at line 101 of file G4NuclearDecayChannel.hh.

References daughterNucleus.

Referenced by G4RadioactiveDecay::AddDecayRateTable().

101 {return daughterNucleus;}
G4ParticleDefinition * daughterNucleus
G4RadioactiveDecayMode G4NuclearDecayChannel::GetDecayMode ( )
inline

Definition at line 95 of file G4NuclearDecayChannel.hh.

References decayMode.

Referenced by G4RadioactiveDecay::AddDecayRateTable(), and G4RadioactiveDecay::LoadDecayTable().

95 {return decayMode;}
const G4RadioactiveDecayMode decayMode
void G4NuclearDecayChannel::SetARM ( G4bool  arm)
inline

Definition at line 93 of file G4NuclearDecayChannel.hh.

References applyARM.

Referenced by G4RadioactiveDecay::LoadDecayTable().

void G4NuclearDecayChannel::SetHLThreshold ( G4double  hl)
inline

Definition at line 87 of file G4NuclearDecayChannel.hh.

References halflifethreshold.

Referenced by G4RadioactiveDecay::LoadDecayTable().

void G4NuclearDecayChannel::SetICM ( G4bool  icm)
inline

Definition at line 90 of file G4NuclearDecayChannel.hh.

References applyICM.

Referenced by G4RadioactiveDecay::LoadDecayTable().

Field Documentation

G4bool G4NuclearDecayChannel::applyARM
protected

Definition at line 150 of file G4NuclearDecayChannel.hh.

Referenced by DecayIt(), G4NuclearDecayChannel(), and SetARM().

G4bool G4NuclearDecayChannel::applyICM
protected

Definition at line 149 of file G4NuclearDecayChannel.hh.

Referenced by DecayIt(), G4NuclearDecayChannel(), and SetICM().

G4int G4NuclearDecayChannel::daughterA
protected

Definition at line 143 of file G4NuclearDecayChannel.hh.

Referenced by DecayIt().

G4double G4NuclearDecayChannel::daughterExcitation
protected

Definition at line 142 of file G4NuclearDecayChannel.hh.

Referenced by DecayIt(), and GetDaughterExcitation().

G4ParticleDefinition* G4NuclearDecayChannel::daughterNucleus
protected

Definition at line 145 of file G4NuclearDecayChannel.hh.

Referenced by GetDaughterNucleus().

G4int G4NuclearDecayChannel::daughterZ
protected

Definition at line 144 of file G4NuclearDecayChannel.hh.

Referenced by DecayIt().

const G4RadioactiveDecayMode G4NuclearDecayChannel::decayMode
protected

Definition at line 139 of file G4NuclearDecayChannel.hh.

Referenced by DecayIt(), and GetDecayMode().

G4ThreadLocal G4DynamicParticle * G4NuclearDecayChannel::dynamicDaughter = 0
staticprotected

Definition at line 146 of file G4NuclearDecayChannel.hh.

Referenced by DecayIt().

G4double G4NuclearDecayChannel::halflifethreshold
protected

Definition at line 148 of file G4NuclearDecayChannel.hh.

Referenced by G4NuclearDecayChannel(), and SetHLThreshold().

const G4double G4NuclearDecayChannel::levelTolerance = 2.0*keV
staticprotected

Definition at line 141 of file G4NuclearDecayChannel.hh.

const G4double G4NuclearDecayChannel::pTolerance = 0.001
staticprotected

Definition at line 140 of file G4NuclearDecayChannel.hh.

const G4double G4NuclearDecayChannel::Qtransition
protected

Definition at line 147 of file G4NuclearDecayChannel.hh.

G4RandGeneral* G4NuclearDecayChannel::RandomEnergy
protected

Definition at line 151 of file G4NuclearDecayChannel.hh.

Referenced by G4NuclearDecayChannel().


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