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 &, 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 ()

Protected Attributes

G4RadioactiveDecayMode decayMode
G4double daughterExcitation
G4int daughterA
G4int daughterZ
G4ParticleDefinitiondaughterNucleus
G4DynamicParticledynamicDaughter
G4double Qtransition
G4double halflifethreshold
G4bool applyICM
G4bool applyARM
CLHEP::RandGeneral * RandomEnergy

Static Protected Attributes

static const G4double pTolerance = 0.001
static const G4double levelTolerance = 2.0*keV

Detailed Description

Definition at line 42 of file G4NuclearDecayChannel.hh.


Constructor & Destructor Documentation

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

Definition at line 55 of file G4NuclearDecayChannel.hh.

00055                                                                       :
00056     G4GeneralPhaseSpaceDecay(Verbose) {;}
  // default constructor

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

Definition at line 97 of file G4NuclearDecayChannel.cc.

References applyARM, applyICM, G4VDecayChannel::FillParent(), G4cout, G4endl, G4ParticleDefinition::GetPDGMass(), G4VDecayChannel::GetVerboseLevel(), halflifethreshold, G4VDecayChannel::parent_mass, Qtransition, G4VDecayChannel::SetBR(), G4VDecayChannel::SetNumberOfDaughters(), and G4VDecayChannel::SetParent().

00105  :G4GeneralPhaseSpaceDecay(Verbose), decayMode(theMode), dynamicDaughter(0),
00106   RandomEnergy(0)
00107 {
00108 #ifdef G4VERBOSE
00109   if (GetVerboseLevel()>1)
00110     {G4cout <<"G4NuclearDecayChannel constructor for " <<G4int(theMode) <<G4endl;}
00111 #endif
00112   SetParent(theParentNucleus);
00113   FillParent();
00114   parent_mass = theParentNucleus->GetPDGMass();
00115   SetBR (theBR);
00116   SetNumberOfDaughters (1);
00117   FillDaughterNucleus (0, A, Z, theDaughterExcitation);
00118   Qtransition = theQtransition;
00119   halflifethreshold = -1.*second;
00120   applyICM = true;
00121   applyARM = true;
00122 }

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

References applyARM, applyICM, G4VDecayChannel::FillParent(), G4cout, G4endl, G4ParticleDefinition::GetPDGMass(), G4VDecayChannel::GetVerboseLevel(), halflifethreshold, G4VDecayChannel::parent_mass, Qtransition, G4VDecayChannel::SetBR(), G4VDecayChannel::SetDaughter(), G4VDecayChannel::SetNumberOfDaughters(), and G4VDecayChannel::SetParent().

00137  :G4GeneralPhaseSpaceDecay(Verbose), decayMode(theMode), dynamicDaughter(0),
00138   RandomEnergy(0)
00139 {
00140 #ifdef G4VERBOSE
00141   if (GetVerboseLevel()>1)
00142     {G4cout <<"G4NuclearDecayChannel constructor for " <<G4int(theMode) <<G4endl;}
00143 #endif
00144   SetParent (theParentNucleus);
00145   FillParent();
00146   parent_mass = theParentNucleus->GetPDGMass();
00147   SetBR (theBR);
00148   SetNumberOfDaughters (2);
00149   SetDaughter(0, theDaughterName1);
00150   FillDaughterNucleus (1, A, Z, theDaughterExcitation);
00151   Qtransition = theQtransition;
00152   halflifethreshold = -1.*second;
00153   applyICM = true;
00154   applyARM = true;
00155 }

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

References applyARM, applyICM, G4VDecayChannel::FillParent(), G4cout, G4endl, G4ParticleDefinition::GetPDGMass(), G4VDecayChannel::GetVerboseLevel(), halflifethreshold, G4VDecayChannel::parent_mass, Qtransition, RandomEnergy, G4VDecayChannel::SetBR(), G4VDecayChannel::SetDaughter(), G4VDecayChannel::SetNumberOfDaughters(), and G4VDecayChannel::SetParent().

00174  :G4GeneralPhaseSpaceDecay(Verbose), decayMode(theMode), dynamicDaughter(0)
00175 {
00176 #ifdef G4VERBOSE
00177   if (GetVerboseLevel()>1)
00178     {G4cout <<"G4NuclearDecayChannel constructor for " <<G4int(theMode) <<G4endl;}
00179 #endif
00180   SetParent (theParentNucleus);
00181   FillParent();
00182   parent_mass = theParentNucleus->GetPDGMass();
00183   SetBR (theBR);
00184   SetNumberOfDaughters (3);
00185   SetDaughter(0, theDaughterName1);
00186   SetDaughter(2, theDaughterName2);
00187   FillDaughterNucleus(1, A, Z, theDaughterExcitation);
00188   RandomEnergy = randBeta;
00189   Qtransition = theQtransition;
00190   halflifethreshold = -1*second;
00191   applyICM = true;
00192   applyARM = true;
00193 }

G4NuclearDecayChannel::~G4NuclearDecayChannel (  )  [inline]

Definition at line 83 of file G4NuclearDecayChannel.hh.

00083 {;} 


Member Function Documentation

G4DecayProducts * G4NuclearDecayChannel::DecayIt ( G4double  theParentMass  )  [virtual]

Reimplemented from G4GeneralPhaseSpaceDecay.

Definition at line 235 of file G4NuclearDecayChannel.cc.

References applyARM, applyICM, G4LossTableManager::AtomDeexcitation(), G4PhotonEvaporation::BreakItUp(), G4PhotonEvaporation::BreakUp(), daughterA, daughterExcitation, G4VDecayChannel::daughters, daughterZ, decayMode, G4VDecayChannel::DumpInfo(), dynamicDaughter, ERROR, FatalException, G4VDecayChannel::FillDaughters(), G4VDecayChannel::FillParent(), G4cerr, G4cout, G4endl, G4Exception(), G4UniformRand, G4VAtomDeexcitation::GenerateParticles(), G4DynamicParticle::Get4Momentum(), G4VAtomDeexcitation::GetAtomicShell(), G4IonTable::GetIon(), G4AtomicShells::GetNumberOfShells(), G4ParticleTable::GetParticleTable(), G4PhotonEvaporation::GetVacantShellNumber(), G4VDecayChannel::GetVerboseLevel(), halflifethreshold, G4LossTableManager::Instance(), G4VAtomDeexcitation::IsFluoActive(), JustWarning, KshellEC, LshellEC, MshellEC, G4VDecayChannel::numberOfDaughters, G4GeneralPhaseSpaceDecay::OneBodyDecayIt(), G4VDecayChannel::parent, G4VDecayChannel::parent_name, G4DecayProducts::PopProducts(), G4DecayProducts::PushProducts(), Qtransition, G4PhotonEvaporation::RDMForced(), G4PhotonEvaporation::SetICM(), G4PhotonEvaporation::SetMaxHalfLife(), G4GeneralPhaseSpaceDecay::SetParentMass(), G4PhotonEvaporation::SetVerboseLevel(), and G4GeneralPhaseSpaceDecay::TwoBodyDecayIt().

00236 {
00237   // Load the details of the parent and daughter particles if they have not
00238   // been defined properly
00239 
00240   if (parent == 0) FillParent();
00241   if (daughters == 0) FillDaughters();
00242 
00243   // We want to ensure that the difference between the total
00244   // parent and daughter masses equals the energy liberated by the transition.
00245 
00246   theParentMass = 0.0;
00247   for( G4int index=0; index < numberOfDaughters; index++)
00248     {theParentMass += daughters[index]->GetPDGMass();}
00249   theParentMass += Qtransition  ;
00250   // bug fix for beta+ decay (flei 25/09/01)
00251   if (decayMode == 2) theParentMass -= 2*0.511 * MeV;
00252   //
00253 #ifdef G4VERBOSE
00254   if (GetVerboseLevel()>1) {
00255     G4cout << "G4NuclearDecayChannel::DecayIt "<< G4endl;
00256     G4cout << "the decay mass = " << theParentMass << G4endl;
00257   }
00258 #endif
00259 
00260   SetParentMass (theParentMass);
00261   
00262   // Define a product vector.
00263   G4DecayProducts* products = 0;
00264 
00265   // Depending upon the number of daughters, select the appropriate decay
00266   // kinematics scheme.
00267   switch (numberOfDaughters) {
00268     case 0:
00269       G4cerr << "G4NuclearDecayChannel::DecayIt ";
00270       G4cerr << " daughters not defined " <<G4endl;
00271       break;
00272     case 1:
00273       products = OneBodyDecayIt();
00274       break;
00275     case 2:
00276       products = TwoBodyDecayIt();
00277       break;
00278     case 3:
00279       products = BetaDecayIt();
00280       break;
00281     default:
00282       // G4cerr <<"Error in G4NuclearDecayChannel::DecayIt" <<G4endl;
00283       // G4cerr <<"Number of daughters in decay = " <<numberOfDaughters <<G4endl;
00284       // G4Exception(__FILE__, G4inttostring(__LINE__), FatalException,  "G4NuclearDecayChannel::DecayIt");
00285     {
00286       G4ExceptionDescription ed;
00287       ed << " More than 3 daughters in decay: N = " << numberOfDaughters
00288          << G4endl;
00289       G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_007",
00290                   FatalException, ed);
00291     }
00292   }
00293 
00294   if (products == 0) {
00295     G4ExceptionDescription ed;
00296     ed << " Parent nucleus " << *parent_name << " was not decayed " << G4endl;
00297     G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_008",
00298                 JustWarning, ed);
00299     DumpInfo();
00300   } else {
00301 
00302     // If the decay is to an excited state of the daughter nuclide, we need
00303     // to apply the photo-evaporation process. This includes the IT decay mode itself.
00304 
00305     // Need to hold the shell idex after ICM
00306     G4int shellIndex = -1;
00307 
00308     if (daughterExcitation > 0.0) {
00309       // Pop the daughter nucleus off the product vector - we need to retain
00310       // the momentum of this particle.
00311 
00312       dynamicDaughter = products->PopProducts();
00313       G4LorentzVector daughterMomentum = dynamicDaughter->Get4Momentum();
00314       G4ThreeVector const daughterMomentum1(static_cast<const G4LorentzVector> (daughterMomentum));
00315 
00316       // Now define a G4Fragment with the correct A, Z and excitation,
00317       // and declare and initialise a G4PhotonEvaporation object    
00318       G4Fragment nucleus(daughterA, daughterZ, daughterMomentum);
00319       G4PhotonEvaporation* deexcitation = new G4PhotonEvaporation;
00320       deexcitation->SetVerboseLevel(GetVerboseLevel());
00321 
00322       // switch on/off internal electron conversion
00323       deexcitation->SetICM(applyICM);
00324 
00325       // Set the maximum life-time for a level that will be treated
00326       // Level with life-time longer than this will be outputed as meta-stable
00327       // isotope
00328       deexcitation->SetMaxHalfLife(halflifethreshold);
00329 
00330       // but in IT mode, we need to force the transition 
00331       if (decayMode == 0) {
00332         deexcitation->RDMForced(true);
00333       } else {
00334         deexcitation->RDMForced(false);
00335       }
00336 
00338       //                                                                      //
00339       //  Now apply photo-evaporation                                         //
00340       //  Use BreakUp() so limit to one transition at a time, if ICM is       //
00341       //  requested.                                                          //
00342       //                                                                      //
00344       //  this change is related to bug #1001  (F.Lei 07/05/2010)
00345 
00346       G4FragmentVector* gammas = 0;     
00347       if (applyICM) {
00348         gammas = deexcitation->BreakUp(nucleus);        
00349       } else {
00350         gammas = deexcitation->BreakItUp(nucleus);
00351       }
00352 
00353       // the returned G4FragmentVector contains the residual nuclide
00354       // as its last entry.
00355       G4int nGammas=gammas->size()-1;
00356 
00357       // Go through each gamma/e- and add it to the decay product.  The angular
00358       // distribution of the gammas is isotropic, and the residual nucleus is
00359       // assumed not to have suffered any recoil as a result of this de-excitation.
00360       for (G4int ig = 0; ig < nGammas; ig++) {
00361         G4DynamicParticle* theGammaRay =
00362           new G4DynamicParticle(gammas->operator[](ig)->GetParticleDefinition(),
00363                                 gammas->operator[](ig)->GetMomentum());
00364         theGammaRay -> SetProperTime(gammas->operator[](ig)->GetCreationTime());
00365         products->PushProducts (theGammaRay);
00366       }
00367 
00368       // now the nucleus
00369       G4double finalDaughterExcitation = gammas->operator[](nGammas)->GetExcitationEnergy();
00370       // f.lei (03/01/03) this is needed to fix the crach in test18 
00371       if (finalDaughterExcitation <= 1.0*keV) finalDaughterExcitation = 0;
00372       
00373       if (dynamicDaughter) delete dynamicDaughter;
00374       
00375       G4IonTable* theIonTable =
00376           (G4IonTable*)(G4ParticleTable::GetParticleTable()->GetIonTable());      
00377       dynamicDaughter = new G4DynamicParticle(
00378                   theIonTable->GetIon(daughterZ,daughterA,finalDaughterExcitation),
00379                   daughterMomentum1);
00380       products->PushProducts (dynamicDaughter); 
00381       
00382       // retrive the ICM shell index
00383       shellIndex = deexcitation->GetVacantShellNumber();
00384       
00385       // Delete/reset variables associated with the gammas.
00386       while (!gammas->empty()) {
00387         delete *(gammas->end()-1);
00388         gammas->pop_back();
00389       }
00390       delete gammas;
00391       delete deexcitation;
00392     } // if daughter excitation > 0
00393 
00394     // Now take care of the EC products which have to go through the ARM
00395     G4int eShell = -1;
00396     if (decayMode == 3 || decayMode == 4 || decayMode == 5) {
00397       switch (decayMode)
00398         {
00399         case KshellEC:
00400           //
00401           {
00402             eShell = 0; // --> 0 from 1 (f.lei 30/4/2008)
00403           }
00404           break;
00405         case LshellEC:
00406           //
00407           {
00408             eShell = G4int(G4UniformRand()*3)+1;
00409           }
00410           break;
00411         case MshellEC:
00412           //
00413           {
00414             // limit the shell index to 6 as specified by the ARM (F.Lei 06/05/2010)
00415             // eShell = G4int(G4UniformRand()*5)+4;
00416             eShell = G4int(G4UniformRand()*3)+4;
00417           }
00418           break;
00419         case ERROR:
00420         default:
00421         G4Exception("G4NuclearDecayChannel::DecayIt()", "HAD_RDM_009",
00422                     FatalException, "Incorrect decay mode selection");
00423         }
00424     } else {
00425       // For other cases eShell comes from shellIndex resulting from  the photo decay
00426       // modeled by G4PhotonEvaporation* de-excitation (see above)
00427       eShell = shellIndex;
00428     }
00429 
00430     // now apply ARM if it is requested and there is a vaccancy
00431     if (applyARM && eShell != -1) {
00432       G4int aZ = daughterZ;
00433 
00434       // V.Ivanchenko migration to new interface to atomic deexcitation
00435       // no check on index of G4MaterialCutsCouple, simplified 
00436       // check on secondary energy Esec < 0.1 keV
00437       G4VAtomDeexcitation* atomDeex =
00438                G4LossTableManager::Instance()->AtomDeexcitation();
00439       if (atomDeex) {
00440         if(atomDeex->IsFluoActive() && aZ > 5 && aZ < 100) {  // only applies to 5< Z <100 
00441           if (eShell >= G4AtomicShells::GetNumberOfShells(aZ)){
00442             eShell = G4AtomicShells::GetNumberOfShells(aZ)-1;
00443           }
00444           G4AtomicShellEnumerator as = G4AtomicShellEnumerator(eShell);
00445           const G4AtomicShell* shell = atomDeex->GetAtomicShell(aZ, as);    
00446           std::vector<G4DynamicParticle*> armProducts;
00447           const G4double deexLimit = 0.1*keV;
00448           atomDeex->GenerateParticles(&armProducts, shell, aZ, deexLimit, deexLimit);
00449           size_t narm = armProducts.size();
00450           if (narm > 0) {
00451             // L.Desorgher: need to initialize dynamicDaughter in some decay
00452             // cases (for example Hg194)
00453             dynamicDaughter = products->PopProducts();
00454             G4ThreeVector bst = dynamicDaughter->Get4Momentum().boostVector();
00455             for (size_t i = 0; i<narm; ++i) {
00456               G4DynamicParticle* dp = armProducts[i];
00457               G4LorentzVector lv = dp->Get4Momentum().boost(bst);
00458               dp->Set4Momentum(lv);
00459               products->PushProducts(dp);
00460             }
00461             products->PushProducts(dynamicDaughter);
00462           }
00463         }
00464       }
00465     }
00466   } // Parent nucleus decayed
00467   /*
00468 
00469     if (atomDeex && aZ > 5 && aZ < 100) {  // only applies to 5< Z <100 
00470       // Retrieve the corresponding identifier and binding energy of the selected shell
00471       const G4AtomicTransitionManager* transitionManager = G4AtomicTransitionManager::Instance();
00472 
00473       //check that eShell is smaller than the max number of shells
00474       //this to avoid a warning from transitionManager(otherwise the output is the same)
00475       //Correction to Bug 1662. L Desorgher (04/02/2011)
00476       if (eShell >= transitionManager->NumberOfShells(aZ)){
00477           eShell=transitionManager->NumberOfShells(aZ)-1;
00478       }
00479 
00480       const G4AtomicShell* shell = transitionManager->Shell(aZ, eShell);
00481       G4double bindingEnergy = shell->BindingEnergy();
00482       G4int shellId = shell->ShellId();
00483 
00484       G4AtomicDeexcitation* atomDeex = new G4AtomicDeexcitation();  
00485       //the default is no Auger electron generation. 
00486       // Switch it on/off here! 
00487       atomDeex->ActivateAugerElectronProduction(true);
00488       std::vector<G4DynamicParticle*>* armProducts = atomDeex->GenerateParticles(aZ,shellId);
00489 
00490       // pop up the daughter before insertion; 
00491       // f.lei (30/04/2008) check if the total kinetic energy is less than 
00492       // the shell binding energy; if true add the difference to the daughter to conserve the energy 
00493       dynamicDaughter = products->PopProducts();
00494       G4double tARMEnergy = 0.0; 
00495       for (size_t i = 0;  i < armProducts->size(); i++) {
00496         products->PushProducts ((*armProducts)[i]);
00497         tARMEnergy += (*armProducts)[i]->GetKineticEnergy();
00498       }
00499       if ((bindingEnergy - tARMEnergy) > 0.1*keV){
00500         G4double dEnergy = dynamicDaughter->GetKineticEnergy() + (bindingEnergy - tARMEnergy);
00501         dynamicDaughter->SetKineticEnergy(dEnergy);
00502       }
00503       products->PushProducts(dynamicDaughter); 
00504 
00505 #ifdef G4VERBOSE
00506       if (GetVerboseLevel()>0)
00507         {
00508           G4cout <<"G4NuclearDecayChannel::Selected shell number for ARM =  " <<shellId <<G4endl;
00509           G4cout <<"G4NuclearDecayChannel::ARM products =  " <<armProducts->size()<<G4endl;
00510           G4cout <<"                 The binding energy =  " << bindingEnergy << G4endl;              
00511           G4cout <<"  Total ARM particle kinetic energy =  " << tARMEnergy << G4endl;              
00512         }
00513 #endif   
00514 
00515       delete armProducts;
00516       delete atomDeex;
00517     }
00518   }
00519   */
00520   return products;
00521 }

G4double G4NuclearDecayChannel::GetDaughterExcitation (  )  [inline]

Definition at line 102 of file G4NuclearDecayChannel.hh.

References daughterExcitation.

Referenced by G4RadioactiveDecay::AddDecayRateTable().

00102 {return daughterExcitation;}

G4ParticleDefinition* G4NuclearDecayChannel::GetDaughterNucleus (  )  [inline]

Definition at line 105 of file G4NuclearDecayChannel.hh.

References daughterNucleus.

Referenced by G4RadioactiveDecay::AddDecayRateTable().

00105 {return daughterNucleus;}

G4RadioactiveDecayMode G4NuclearDecayChannel::GetDecayMode (  )  [inline]

Definition at line 99 of file G4NuclearDecayChannel.hh.

References decayMode.

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

00099 {return decayMode;}

void G4NuclearDecayChannel::SetARM ( G4bool  arm  )  [inline]

Definition at line 96 of file G4NuclearDecayChannel.hh.

References applyARM.

Referenced by G4RadioactiveDecay::LoadDecayTable().

00096 {applyARM = arm;}

void G4NuclearDecayChannel::SetHLThreshold ( G4double  hl  )  [inline]

Definition at line 90 of file G4NuclearDecayChannel.hh.

References halflifethreshold.

Referenced by G4RadioactiveDecay::LoadDecayTable().

00090 {halflifethreshold = hl;}

void G4NuclearDecayChannel::SetICM ( G4bool  icm  )  [inline]

Definition at line 93 of file G4NuclearDecayChannel.hh.

References applyICM.

Referenced by G4RadioactiveDecay::LoadDecayTable().

00093 {applyICM = icm;}


Field Documentation

G4bool G4NuclearDecayChannel::applyARM [protected]

Definition at line 147 of file G4NuclearDecayChannel.hh.

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

G4bool G4NuclearDecayChannel::applyICM [protected]

Definition at line 146 of file G4NuclearDecayChannel.hh.

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

G4int G4NuclearDecayChannel::daughterA [protected]

Definition at line 140 of file G4NuclearDecayChannel.hh.

Referenced by DecayIt().

G4double G4NuclearDecayChannel::daughterExcitation [protected]

Definition at line 139 of file G4NuclearDecayChannel.hh.

Referenced by DecayIt(), and GetDaughterExcitation().

G4ParticleDefinition* G4NuclearDecayChannel::daughterNucleus [protected]

Definition at line 142 of file G4NuclearDecayChannel.hh.

Referenced by GetDaughterNucleus().

G4int G4NuclearDecayChannel::daughterZ [protected]

Definition at line 141 of file G4NuclearDecayChannel.hh.

Referenced by DecayIt().

G4RadioactiveDecayMode G4NuclearDecayChannel::decayMode [protected]

Definition at line 136 of file G4NuclearDecayChannel.hh.

Referenced by DecayIt(), and GetDecayMode().

G4DynamicParticle* G4NuclearDecayChannel::dynamicDaughter [protected]

Definition at line 143 of file G4NuclearDecayChannel.hh.

Referenced by DecayIt().

G4double G4NuclearDecayChannel::halflifethreshold [protected]

Definition at line 145 of file G4NuclearDecayChannel.hh.

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

const G4double G4NuclearDecayChannel::levelTolerance = 2.0*keV [static, protected]

Definition at line 138 of file G4NuclearDecayChannel.hh.

const G4double G4NuclearDecayChannel::pTolerance = 0.001 [static, protected]

Definition at line 137 of file G4NuclearDecayChannel.hh.

G4double G4NuclearDecayChannel::Qtransition [protected]

Definition at line 144 of file G4NuclearDecayChannel.hh.

Referenced by DecayIt(), and G4NuclearDecayChannel().

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:
Generated on Mon May 27 17:52:44 2013 for Geant4 by  doxygen 1.4.7