Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NuclearDecayChannel Class Reference

#include <G4NuclearDecayChannel.hh>

+ Inheritance diagram for G4NuclearDecayChannel:

Public Member Functions

 G4NuclearDecayChannel (const G4RadioactiveDecayMode &, 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, G4bool, CLHEP::RandGeneral *randBeta, G4double theQtransition, G4int A, G4int Z, G4double theDaughterExcitation, const G4String theDaughterName1, const G4String theDaughterName2)
 
 ~G4NuclearDecayChannel ()
 
G4DecayProductsDecayIt (G4double theParentMass)
 
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)
 
virtual G4DecayProductsDecayIt (G4double mass=0.0)
 
- 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
 
virtual G4DecayProductsDecayIt (G4double parentMass=-1.0)=0
 
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

G4RadioactiveDecayMode decayMode
 
G4double daughterExcitation
 
G4int daughterA
 
G4int daughterZ
 
G4ParticleDefinitiondaughterNucleus
 
G4DynamicParticledynamicDaughter
 
G4double Qtransition
 
G4double halflifethreshold
 
G4bool applyICM
 
G4bool applyARM
 
CLHEP::RandGeneralRandomEnergy
 
- Protected Attributes inherited from G4VDecayChannel
G4String kinematics_name
 
G4double rbranch
 
G4int numberOfDaughters
 
G4Stringparent_name
 
G4String ** daughters_name
 
G4ParticleTableparticletable
 
G4ParticleDefinitionparent
 
G4ParticleDefinition ** daughters
 
G4double parent_mass
 
G4doubledaughters_mass
 
G4int verboseLevel
 

Static Protected Attributes

static const G4double pTolerance = 0.001
 
static const G4double levelTolerance = 2.0*keV
 
- 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
 G4VDecayChannel ()
 
 G4VDecayChannel (const G4VDecayChannel &)
 
G4VDecayChanneloperator= (const G4VDecayChannel &)
 
void ClearDaughtersName ()
 
void FillDaughters ()
 
void FillParent ()
 

Detailed Description

Definition at line 42 of file G4NuclearDecayChannel.hh.

Constructor & Destructor Documentation

◆ G4NuclearDecayChannel() [1/4]

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

◆ G4NuclearDecayChannel() [2/4]

G4NuclearDecayChannel::G4NuclearDecayChannel ( const G4RadioactiveDecayMode theMode,
G4int  Verbose,
const G4ParticleDefinition theParentNucleus,
G4double  theBR,
G4double  theQtransition,
G4int  A,
G4int  Z,
G4double  theDaughterExcitation 
)

Definition at line 96 of file G4NuclearDecayChannel.cc.

105 :G4GeneralPhaseSpaceDecay(Verbose), decayMode(theMode), dynamicDaughter(0),
106 RandomEnergy(0)
107{
108#ifdef G4VERBOSE
109 if (GetVerboseLevel()>1)
110 {G4cout <<"G4NuclearDecayChannel constructor for " <<G4int(theMode) <<G4endl;}
111#endif
112 SetParent(theParentNucleus);
113 FillParent();
114 parent_mass = theParentNucleus->GetPDGMass();
115 SetBR (theBR);
117 FillDaughterNucleus (0, A, Z, theDaughterExcitation);
118 Qtransition = theQtransition;
119 halflifethreshold = -1.*second;
120 applyICM = true;
121 applyARM = true;
122}
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4RadioactiveDecayMode decayMode
CLHEP::RandGeneral * RandomEnergy
G4DynamicParticle * dynamicDaughter
void SetBR(G4double value)
G4int GetVerboseLevel() const
void SetNumberOfDaughters(G4int value)
void SetParent(const G4ParticleDefinition *particle_type)

◆ G4NuclearDecayChannel() [3/4]

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 127 of file G4NuclearDecayChannel.cc.

137 :G4GeneralPhaseSpaceDecay(Verbose), decayMode(theMode), dynamicDaughter(0),
138 RandomEnergy(0)
139{
140#ifdef G4VERBOSE
141 if (GetVerboseLevel()>1)
142 {G4cout <<"G4NuclearDecayChannel constructor for " <<G4int(theMode) <<G4endl;}
143#endif
144 SetParent (theParentNucleus);
145 FillParent();
146 parent_mass = theParentNucleus->GetPDGMass();
147 SetBR (theBR);
149 SetDaughter(0, theDaughterName1);
150 FillDaughterNucleus (1, A, Z, theDaughterExcitation);
151 Qtransition = theQtransition;
152 halflifethreshold = -1.*second;
153 applyICM = true;
154 applyARM = true;
155}
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)

◆ G4NuclearDecayChannel() [4/4]

G4NuclearDecayChannel::G4NuclearDecayChannel ( const G4RadioactiveDecayMode theMode,
G4int  Verbose,
const G4ParticleDefinition theParentNucleus,
G4double  theBR,
G4double  ,
G4bool  ,
CLHEP::RandGeneral randBeta,
G4double  theQtransition,
G4int  A,
G4int  Z,
G4double  theDaughterExcitation,
const G4String  theDaughterName1,
const G4String  theDaughterName2 
)

Definition at line 160 of file G4NuclearDecayChannel.cc.

175{
176#ifdef G4VERBOSE
177 if (GetVerboseLevel()>1)
178 {G4cout <<"G4NuclearDecayChannel constructor for " <<G4int(theMode) <<G4endl;}
179#endif
180 SetParent (theParentNucleus);
181 FillParent();
182 parent_mass = theParentNucleus->GetPDGMass();
183 SetBR (theBR);
185 SetDaughter(0, theDaughterName1);
186 SetDaughter(2, theDaughterName2);
187 FillDaughterNucleus(1, A, Z, theDaughterExcitation);
188 RandomEnergy = randBeta;
189 Qtransition = theQtransition;
190 halflifethreshold = -1*second;
191 applyICM = true;
192 applyARM = true;
193}

◆ ~G4NuclearDecayChannel()

G4NuclearDecayChannel::~G4NuclearDecayChannel ( )
inline

Definition at line 83 of file G4NuclearDecayChannel.hh.

83{;}

Member Function Documentation

◆ DecayIt()

G4DecayProducts * G4NuclearDecayChannel::DecayIt ( G4double  theParentMass)
virtual

Reimplemented from G4GeneralPhaseSpaceDecay.

Definition at line 235 of file G4NuclearDecayChannel.cc.

236{
237 // Load the details of the parent and daughter particles if they have not
238 // been defined properly
239
240 if (parent == 0) FillParent();
241 if (daughters == 0) FillDaughters();
242
243 // We want to ensure that the difference between the total
244 // parent and daughter masses equals the energy liberated by the transition.
245
246 theParentMass = 0.0;
247 for( G4int index=0; index < numberOfDaughters; index++)
248 {theParentMass += daughters[index]->GetPDGMass();}
249 theParentMass += Qtransition ;
250 // bug fix for beta+ decay (flei 25/09/01)
251 if (decayMode == 2) theParentMass -= 2*0.511 * MeV;
252 //
253#ifdef G4VERBOSE
254 if (GetVerboseLevel()>1) {
255 G4cout << "G4NuclearDecayChannel::DecayIt "<< G4endl;
256 G4cout << "the decay mass = " << theParentMass << G4endl;
257 }
258#endif
259
260 SetParentMass (theParentMass);
261
262 // Define a product vector.
263 G4DecayProducts* products = 0;
264
265 // Depending upon the number of daughters, select the appropriate decay
266 // kinematics scheme.
267 switch (numberOfDaughters) {
268 case 0:
269 G4cerr << "G4NuclearDecayChannel::DecayIt ";
270 G4cerr << " daughters not defined " <<G4endl;
271 break;
272 case 1:
273 products = OneBodyDecayIt();
274 break;
275 case 2:
276 products = TwoBodyDecayIt();
277 break;
278 case 3:
279 products = BetaDecayIt();
280 break;
281 default:
282 // G4cerr <<"Error in G4NuclearDecayChannel::DecayIt" <<G4endl;
283 // G4cerr <<"Number of daughters in decay = " <<numberOfDaughters <<G4endl;
284 // G4Exception(__FILE__, G4inttostring(__LINE__), FatalException, "G4NuclearDecayChannel::DecayIt");
285 {
287 ed << " More than 3 daughters in decay: N = " << numberOfDaughters
288 << G4endl;
289 G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_007",
290 FatalException, ed);
291 }
292 }
293
294 if (products == 0) {
296 ed << " Parent nucleus " << *parent_name << " was not decayed " << G4endl;
297 G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_008",
298 JustWarning, ed);
299 DumpInfo();
300 } else {
301
302 // If the decay is to an excited state of the daughter nuclide, we need
303 // to apply the photo-evaporation process. This includes the IT decay mode itself.
304
305 // Need to hold the shell idex after ICM
306 G4int shellIndex = -1;
307
308 if (daughterExcitation > 0.0) {
309 // Pop the daughter nucleus off the product vector - we need to retain
310 // the momentum of this particle.
311
312 dynamicDaughter = products->PopProducts();
313 G4LorentzVector daughterMomentum = dynamicDaughter->Get4Momentum();
314 G4ThreeVector const daughterMomentum1(static_cast<const G4LorentzVector> (daughterMomentum));
315
316 // Now define a G4Fragment with the correct A, Z and excitation,
317 // and declare and initialise a G4PhotonEvaporation object
318 G4Fragment nucleus(daughterA, daughterZ, daughterMomentum);
319 G4PhotonEvaporation* deexcitation = new G4PhotonEvaporation;
320 deexcitation->SetVerboseLevel(GetVerboseLevel());
321
322 // switch on/off internal electron conversion
323 deexcitation->SetICM(applyICM);
324
325 // Set the maximum life-time for a level that will be treated
326 // Level with life-time longer than this will be outputed as meta-stable
327 // isotope
328 deexcitation->SetMaxHalfLife(halflifethreshold);
329
330 // but in IT mode, we need to force the transition
331 if (decayMode == 0) {
332 deexcitation->RDMForced(true);
333 } else {
334 deexcitation->RDMForced(false);
335 }
336
337 //////////////////////////////////////////////////////////////////////////
338 // //
339 // Now apply photo-evaporation //
340 // Use BreakUp() so limit to one transition at a time, if ICM is //
341 // requested. //
342 // //
343 //////////////////////////////////////////////////////////////////////////
344 // this change is related to bug #1001 (F.Lei 07/05/2010)
345
346 G4FragmentVector* gammas = 0;
347 if (applyICM) {
348 gammas = deexcitation->BreakUp(nucleus);
349 } else {
350 gammas = deexcitation->BreakItUp(nucleus);
351 }
352
353 // the returned G4FragmentVector contains the residual nuclide
354 // as its last entry.
355 G4int nGammas=gammas->size()-1;
356
357 // Go through each gamma/e- and add it to the decay product. The angular
358 // distribution of the gammas is isotropic, and the residual nucleus is
359 // assumed not to have suffered any recoil as a result of this de-excitation.
360 for (G4int ig = 0; ig < nGammas; ig++) {
361 G4DynamicParticle* theGammaRay =
362 new G4DynamicParticle(gammas->operator[](ig)->GetParticleDefinition(),
363 gammas->operator[](ig)->GetMomentum());
364 theGammaRay -> SetProperTime(gammas->operator[](ig)->GetCreationTime());
365 products->PushProducts (theGammaRay);
366 }
367
368 // now the nucleus
369 G4double finalDaughterExcitation = gammas->operator[](nGammas)->GetExcitationEnergy();
370 // f.lei (03/01/03) this is needed to fix the crach in test18
371 if (finalDaughterExcitation <= 1.0*keV) finalDaughterExcitation = 0;
372
374
375 G4IonTable* theIonTable =
378 theIonTable->GetIon(daughterZ,daughterA,finalDaughterExcitation),
379 daughterMomentum1);
380 products->PushProducts (dynamicDaughter);
381
382 // retrive the ICM shell index
383 shellIndex = deexcitation->GetVacantShellNumber();
384
385 // Delete/reset variables associated with the gammas.
386 while (!gammas->empty()) {
387 delete *(gammas->end()-1);
388 gammas->pop_back();
389 }
390 delete gammas;
391 delete deexcitation;
392 } // if daughter excitation > 0
393
394 // Now take care of the EC products which have to go through the ARM
395 G4int eShell = -1;
396 if (decayMode == 3 || decayMode == 4 || decayMode == 5) {
397 switch (decayMode)
398 {
399 case KshellEC:
400 //
401 {
402 eShell = 0; // --> 0 from 1 (f.lei 30/4/2008)
403 }
404 break;
405 case LshellEC:
406 //
407 {
408 eShell = G4int(G4UniformRand()*3)+1;
409 }
410 break;
411 case MshellEC:
412 //
413 {
414 // limit the shell index to 6 as specified by the ARM (F.Lei 06/05/2010)
415 // eShell = G4int(G4UniformRand()*5)+4;
416 eShell = G4int(G4UniformRand()*3)+4;
417 }
418 break;
419 case ERROR:
420 default:
421 G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_009",
422 FatalException, "Incorrect decay mode selection");
423 }
424 } else {
425 // For other cases eShell comes from shellIndex resulting from the photo decay
426 // modeled by G4PhotonEvaporation* de-excitation (see above)
427 eShell = shellIndex;
428 }
429
430 // now apply ARM if it is requested and there is a vaccancy
431 if (applyARM && eShell != -1) {
432 G4int aZ = daughterZ;
433
434 // V.Ivanchenko migration to new interface to atomic deexcitation
435 // no check on index of G4MaterialCutsCouple, simplified
436 // check on secondary energy Esec < 0.1 keV
437 G4VAtomDeexcitation* atomDeex =
439 if (atomDeex) {
440 if(atomDeex->IsFluoActive() && aZ > 5 && aZ < 100) { // only applies to 5< Z <100
441 if (eShell >= G4AtomicShells::GetNumberOfShells(aZ)){
443 }
445 const G4AtomicShell* shell = atomDeex->GetAtomicShell(aZ, as);
446 std::vector<G4DynamicParticle*> armProducts;
447 const G4double deexLimit = 0.1*keV;
448 atomDeex->GenerateParticles(&armProducts, shell, aZ, deexLimit, deexLimit);
449 size_t narm = armProducts.size();
450 if (narm > 0) {
451 // L.Desorgher: need to initialize dynamicDaughter in some decay
452 // cases (for example Hg194)
453 dynamicDaughter = products->PopProducts();
455 for (size_t i = 0; i<narm; ++i) {
456 G4DynamicParticle* dp = armProducts[i];
457 G4LorentzVector lv = dp->Get4Momentum().boost(bst);
458 dp->Set4Momentum(lv);
459 products->PushProducts(dp);
460 }
461 products->PushProducts(dynamicDaughter);
462 }
463 }
464 }
465 }
466 } // Parent nucleus decayed
467 /*
468
469 if (atomDeex && aZ > 5 && aZ < 100) { // only applies to 5< Z <100
470 // Retrieve the corresponding identifier and binding energy of the selected shell
471 const G4AtomicTransitionManager* transitionManager = G4AtomicTransitionManager::Instance();
472
473 //check that eShell is smaller than the max number of shells
474 //this to avoid a warning from transitionManager(otherwise the output is the same)
475 //Correction to Bug 1662. L Desorgher (04/02/2011)
476 if (eShell >= transitionManager->NumberOfShells(aZ)){
477 eShell=transitionManager->NumberOfShells(aZ)-1;
478 }
479
480 const G4AtomicShell* shell = transitionManager->Shell(aZ, eShell);
481 G4double bindingEnergy = shell->BindingEnergy();
482 G4int shellId = shell->ShellId();
483
484 G4AtomicDeexcitation* atomDeex = new G4AtomicDeexcitation();
485 //the default is no Auger electron generation.
486 // Switch it on/off here!
487 atomDeex->ActivateAugerElectronProduction(true);
488 std::vector<G4DynamicParticle*>* armProducts = atomDeex->GenerateParticles(aZ,shellId);
489
490 // pop up the daughter before insertion;
491 // f.lei (30/04/2008) check if the total kinetic energy is less than
492 // the shell binding energy; if true add the difference to the daughter to conserve the energy
493 dynamicDaughter = products->PopProducts();
494 G4double tARMEnergy = 0.0;
495 for (size_t i = 0; i < armProducts->size(); i++) {
496 products->PushProducts ((*armProducts)[i]);
497 tARMEnergy += (*armProducts)[i]->GetKineticEnergy();
498 }
499 if ((bindingEnergy - tARMEnergy) > 0.1*keV){
500 G4double dEnergy = dynamicDaughter->GetKineticEnergy() + (bindingEnergy - tARMEnergy);
501 dynamicDaughter->SetKineticEnergy(dEnergy);
502 }
503 products->PushProducts(dynamicDaughter);
504
505#ifdef G4VERBOSE
506 if (GetVerboseLevel()>0)
507 {
508 G4cout <<"G4NuclearDecayChannel::Selected shell number for ARM = " <<shellId <<G4endl;
509 G4cout <<"G4NuclearDecayChannel::ARM products = " <<armProducts->size()<<G4endl;
510 G4cout <<" The binding energy = " << bindingEnergy << G4endl;
511 G4cout <<" Total ARM particle kinetic energy = " << tARMEnergy << G4endl;
512 }
513#endif
514
515 delete armProducts;
516 delete atomDeex;
517 }
518 }
519 */
520 return products;
521}
G4AtomicShellEnumerator
@ JustWarning
@ FatalException
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:65
#define ERROR(x)
double G4double
Definition: G4Types.hh:64
G4DLLIMPORT std::ostream G4cerr
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
static G4int GetNumberOfShells(G4int Z)
G4DynamicParticle * PopProducts()
G4int PushProducts(G4DynamicParticle *aParticle)
G4LorentzVector Get4Momentum() const
void Set4Momentum(const G4LorentzVector &momentum)
void SetParentMass(const G4double aParentMass)
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int J=0)
Definition: G4IonTable.cc:267
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
static G4ParticleTable * GetParticleTable()
G4IonTable * GetIonTable()
virtual G4FragmentVector * BreakItUp(const G4Fragment &nucleus)
void SetMaxHalfLife(G4double)
virtual G4FragmentVector * BreakUp(const G4Fragment &nucleus)
void SetVerboseLevel(G4int verbose)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
G4String * parent_name
G4ParticleDefinition * parent
G4ParticleDefinition ** daughters
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76

◆ GetDaughterExcitation()

G4double G4NuclearDecayChannel::GetDaughterExcitation ( )
inline

Definition at line 102 of file G4NuclearDecayChannel.hh.

102{return daughterExcitation;}

Referenced by G4RadioactiveDecay::AddDecayRateTable().

◆ GetDaughterNucleus()

G4ParticleDefinition * G4NuclearDecayChannel::GetDaughterNucleus ( )
inline

Definition at line 105 of file G4NuclearDecayChannel.hh.

105{return daughterNucleus;}
G4ParticleDefinition * daughterNucleus

Referenced by G4RadioactiveDecay::AddDecayRateTable().

◆ GetDecayMode()

G4RadioactiveDecayMode G4NuclearDecayChannel::GetDecayMode ( )
inline

◆ SetARM()

void G4NuclearDecayChannel::SetARM ( G4bool  arm)
inline

Definition at line 96 of file G4NuclearDecayChannel.hh.

96{applyARM = arm;}

Referenced by G4RadioactiveDecay::LoadDecayTable().

◆ SetHLThreshold()

void G4NuclearDecayChannel::SetHLThreshold ( G4double  hl)
inline

Definition at line 90 of file G4NuclearDecayChannel.hh.

Referenced by G4RadioactiveDecay::LoadDecayTable().

◆ SetICM()

void G4NuclearDecayChannel::SetICM ( G4bool  icm)
inline

Definition at line 93 of file G4NuclearDecayChannel.hh.

93{applyICM = icm;}

Referenced by G4RadioactiveDecay::LoadDecayTable().

Member Data Documentation

◆ applyARM

G4bool G4NuclearDecayChannel::applyARM
protected

Definition at line 147 of file G4NuclearDecayChannel.hh.

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

◆ applyICM

G4bool G4NuclearDecayChannel::applyICM
protected

Definition at line 146 of file G4NuclearDecayChannel.hh.

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

◆ daughterA

G4int G4NuclearDecayChannel::daughterA
protected

Definition at line 140 of file G4NuclearDecayChannel.hh.

Referenced by DecayIt().

◆ daughterExcitation

G4double G4NuclearDecayChannel::daughterExcitation
protected

Definition at line 139 of file G4NuclearDecayChannel.hh.

Referenced by DecayIt(), and GetDaughterExcitation().

◆ daughterNucleus

G4ParticleDefinition* G4NuclearDecayChannel::daughterNucleus
protected

Definition at line 142 of file G4NuclearDecayChannel.hh.

Referenced by GetDaughterNucleus().

◆ daughterZ

G4int G4NuclearDecayChannel::daughterZ
protected

Definition at line 141 of file G4NuclearDecayChannel.hh.

Referenced by DecayIt().

◆ decayMode

G4RadioactiveDecayMode G4NuclearDecayChannel::decayMode
protected

Definition at line 136 of file G4NuclearDecayChannel.hh.

Referenced by DecayIt(), and GetDecayMode().

◆ dynamicDaughter

G4DynamicParticle* G4NuclearDecayChannel::dynamicDaughter
protected

Definition at line 143 of file G4NuclearDecayChannel.hh.

Referenced by DecayIt().

◆ halflifethreshold

G4double G4NuclearDecayChannel::halflifethreshold
protected

Definition at line 145 of file G4NuclearDecayChannel.hh.

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

◆ levelTolerance

const G4double G4NuclearDecayChannel::levelTolerance = 2.0*keV
staticprotected

Definition at line 138 of file G4NuclearDecayChannel.hh.

◆ pTolerance

const G4double G4NuclearDecayChannel::pTolerance = 0.001
staticprotected

Definition at line 137 of file G4NuclearDecayChannel.hh.

◆ Qtransition

G4double G4NuclearDecayChannel::Qtransition
protected

Definition at line 144 of file G4NuclearDecayChannel.hh.

Referenced by DecayIt(), and G4NuclearDecayChannel().

◆ RandomEnergy

CLHEP::RandGeneral* G4NuclearDecayChannel::RandomEnergy
protected

Definition at line 148 of file G4NuclearDecayChannel.hh.

Referenced by G4NuclearDecayChannel().


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