Geant4  10.03
G4RadioactiveDecay Class Reference

#include <G4RadioactiveDecay.hh>

+ Inheritance diagram for G4RadioactiveDecay:
+ Collaboration diagram for G4RadioactiveDecay:

Public Member Functions

 G4RadioactiveDecay (const G4String &processName="RadioactiveDecay")
 
 ~G4RadioactiveDecay ()
 
G4bool IsApplicable (const G4ParticleDefinition &)
 
G4DecayTableGetDecayTable (const G4ParticleDefinition *)
 
void SelectAVolume (const G4String aVolume)
 
void DeselectAVolume (const G4String aVolume)
 
void SelectAllVolumes ()
 
void DeselectAllVolumes ()
 
void SetDecayBias (G4String filename)
 
void SetHLThreshold (G4double hl)
 
void SetICM (G4bool icm)
 
void SetARM (G4bool arm)
 
void SetSourceTimeProfile (G4String filename)
 
G4bool IsRateTableReady (const G4ParticleDefinition &)
 
void AddDecayRateTable (const G4ParticleDefinition &)
 
void GetDecayRateTable (const G4ParticleDefinition &)
 
void SetDecayRate (G4int, G4int, G4double, G4int, std::vector< G4double >, std::vector< G4double >)
 
std::vector
< G4RadioactivityTable * > 
GetTheRadioactivityTables ()
 
G4DecayTableLoadDecayTable (const G4ParticleDefinition &theParentNucleus)
 
void AddUserDecayDataFile (G4int Z, G4int A, G4String filename)
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
void SetNucleusLimits (G4NucleusLimits theNucleusLimits1)
 
G4NucleusLimits GetNucleusLimits () const
 
void SetAnalogueMonteCarlo (G4bool r)
 
void SetFBeta (G4bool r)
 
G4bool IsAnalogueMonteCarlo ()
 
void SetBRBias (G4bool r)
 
void SetSplitNuclei (G4int r)
 
G4int GetSplitNuclei ()
 
void SetDecayDirection (const G4ThreeVector &theDir)
 
const G4ThreeVectorGetDecayDirection () const
 
void SetDecayHalfAngle (G4double halfAngle=0.*CLHEP::deg)
 
G4double GetDecayHalfAngle () const
 
void SetDecayCollimation (const G4ThreeVector &theDir, G4double halfAngle=0.*CLHEP::deg)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4VParticleChangeDecayIt (const G4Track &theTrack, const G4Step &theStep)
 
- Public Member Functions inherited from G4VRestDiscreteProcess
 G4VRestDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VRestDiscreteProcess (G4VRestDiscreteProcess &)
 
virtual ~G4VRestDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

G4DecayProductsDoDecay (const G4ParticleDefinition &theParticleDef)
 
void CollimateDecay (G4DecayProducts *products)
 
void CollimateDecayProduct (G4DynamicParticle *product)
 
G4ThreeVector ChooseCollimationDirection () const
 
G4double GetMeanFreePath (const G4Track &theTrack, G4double previousStepSize, G4ForceCondition *condition)
 
G4double GetMeanLifeTime (const G4Track &theTrack, G4ForceCondition *condition)
 
G4double ConvolveSourceTimeProfile (const G4double, const G4double)
 
G4double GetDecayTime ()
 
G4int GetDecayTimeBin (const G4double aDecayTime)
 
void AddDeexcitationSpectrumForBiasMode (G4ParticleDefinition *apartDef, G4double weight, G4double currenTime, std::vector< double > &weights_v, std::vector< double > &times_v, std::vector< G4DynamicParticle * > &secondaries_v)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Private Member Functions

 G4RadioactiveDecay (const G4RadioactiveDecay &right)
 
G4RadioactiveDecayoperator= (const G4RadioactiveDecay &right)
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)
 
G4VParticleChangeAtRestDoIt (const G4Track &theTrack, const G4Step &theStep)
 
G4VParticleChangePostStepDoIt (const G4Track &theTrack, const G4Step &theStep)
 

Private Attributes

G4RadioactiveDecaymessengertheRadioactiveDecaymessenger
 
G4NucleusLimits theNucleusLimits
 
G4bool isInitialised
 
G4bool AnalogueMC
 
G4bool BRBias
 
G4bool FBeta
 
G4int NSplit
 
G4double halflifethreshold
 
G4bool applyICM
 
G4bool applyARM
 
G4ThreeVector forceDecayDirection
 
G4double forceDecayHalfAngle
 
G4int NSourceBin
 
G4double SBin [100]
 
G4double SProfile [100]
 
G4int NDecayBin
 
G4double DBin [100]
 
G4double DProfile [100]
 
std::vector< G4StringValidVolumes
 
bool isAllVolumesMode
 
G4RadioactiveDecayRate theDecayRate
 
G4RadioactiveDecayRates theDecayRateVector
 
G4RadioactiveDecayRateVector theDecayRateTable
 
G4RadioactiveDecayRateTable theDecayRateTableVector
 
std::vector
< G4RadioactivityTable * > 
theRadioactivityTables
 
G4int decayWindows [100]
 
std::map< G4int, G4StringtheUserRadioactiveDataFiles
 
DecayTableMapdkmap
 
G4double fRemainderLifeTime
 
G4int verboseLevel
 
G4ParticleChangeForRadDecay fParticleChangeForRadDecay
 

Static Private Attributes

static const G4ThreeVector origin
 
static const G4double levelTolerance = 0.1*keV
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 80 of file G4RadioactiveDecay.hh.

Constructor & Destructor Documentation

G4RadioactiveDecay::G4RadioactiveDecay ( const G4String processName = "RadioactiveDecay")
G4RadioactiveDecay::~G4RadioactiveDecay ( )

Definition at line 212 of file G4RadioactiveDecay.cc.

References dkmap, and theRadioactiveDecaymessenger.

G4RadioactiveDecay::G4RadioactiveDecay ( const G4RadioactiveDecay right)
private

Member Function Documentation

void G4RadioactiveDecay::AddDeexcitationSpectrumForBiasMode ( G4ParticleDefinition apartDef,
G4double  weight,
G4double  currenTime,
std::vector< double > &  weights_v,
std::vector< double > &  times_v,
std::vector< G4DynamicParticle * > &  secondaries_v 
)
protected

Definition at line 2014 of file G4RadioactiveDecay.cc.

References G4ITDecay::DecayIt(), G4DecayProducts::entries(), G4ParticleDefinition::GetBaryonNumber(), G4DynamicParticle::GetDefinition(), G4ParticleDefinition::GetPDGLifeTime(), and G4DecayProducts::PopProducts().

Referenced by DecayIt().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void G4RadioactiveDecay::AddUserDecayDataFile ( G4int  Z,
G4int  A,
G4String  filename 
)

Definition at line 1107 of file G4RadioactiveDecay.cc.

References G4cout, G4endl, and theUserRadioactiveDataFiles.

Referenced by G4RadioactiveDecaymessenger::SetNewValue().

+ Here is the caller graph for this function:

G4VParticleChange* G4RadioactiveDecay::AtRestDoIt ( const G4Track theTrack,
const G4Step theStep 
)
inlineprivatevirtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 340 of file G4RadioactiveDecay.hh.

References DecayIt().

+ Here is the call graph for this function:

G4double G4RadioactiveDecay::AtRestGetPhysicalInteractionLength ( const G4Track track,
G4ForceCondition condition 
)
inlineprivatevirtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 331 of file G4RadioactiveDecay.hh.

References G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength(), and fRemainderLifeTime.

+ Here is the call graph for this function:

void G4RadioactiveDecay::BuildPhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 718 of file G4RadioactiveDecay.cc.

References G4LossTableManager::AtomDeexcitation(), FatalException, G4Exception(), G4LossTableManager::Instance(), and isInitialised.

+ Here is the call graph for this function:

G4ThreeVector G4RadioactiveDecay::ChooseCollimationDirection ( ) const
protected

Definition at line 1988 of file G4RadioactiveDecay.cc.

References deg, forceDecayDirection, forceDecayHalfAngle, G4cout, G4endl, G4UniformRand, GetVerboseLevel(), origin, and pi.

Referenced by CollimateDecayProduct().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void G4RadioactiveDecay::CollimateDecay ( G4DecayProducts products)
protected

Definition at line 1945 of file G4RadioactiveDecay.cc.

References alpha, CollimateDecayProduct(), G4Electron::Definition(), G4Alpha::Definition(), G4Proton::Definition(), G4Positron::Definition(), G4Neutron::Definition(), G4Gamma::Definition(), deg, G4InuclParticleNames::electron, G4DecayProducts::entries(), forceDecayDirection, forceDecayHalfAngle, G4cout, G4endl, G4DynamicParticle::GetParticleDefinition(), GetVerboseLevel(), neutron, origin, G4InuclParticleNames::positron, and G4InuclParticleNames::proton.

Referenced by DoDecay().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void G4RadioactiveDecay::CollimateDecayProduct ( G4DynamicParticle product)
protected

Definition at line 1973 of file G4RadioactiveDecay.cc.

References ChooseCollimationDirection(), G4cout, G4endl, G4DynamicParticle::GetParticleDefinition(), G4ParticleDefinition::GetParticleName(), GetVerboseLevel(), origin, and G4DynamicParticle::SetMomentumDirection().

Referenced by CollimateDecay().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

G4double G4RadioactiveDecay::ConvolveSourceTimeProfile ( const G4double  t,
const G4double  tau 
)
protected

Definition at line 390 of file G4RadioactiveDecay.cc.

References G4cout, G4endl, G4Exception(), GetVerboseLevel(), JustWarning, L, NSourceBin, SBin, and SProfile.

Referenced by DecayIt().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

G4VParticleChange * G4RadioactiveDecay::DecayIt ( const G4Track theTrack,
const G4Step theStep 
)

Definition at line 1561 of file G4RadioactiveDecay.cc.

References AddDecayRateTable(), AddDeexcitationSpectrumForBiasMode(), G4ParticleChangeForRadDecay::AddSecondary(), AnalogueMC, G4DecayProducts::Boost(), BRBias, G4VProcess::ClearNumberOfInteractionLengthLeft(), cm, ConvolveSourceTimeProfile(), DBin, decayWindows, DoDecay(), DProfile, G4DecayTable::DumpInfo(), G4DecayProducts::DumpInfo(), G4DecayTable::entries(), G4DecayProducts::entries(), fParticleChangeForRadDecay, fStopAndKill, fStopButAlive, G4cerr, G4cout, G4endl, G4UniformRand, G4ParticleDefinition::GetBaryonNumber(), G4DecayTable::GetDecayChannel(), GetDecayRateTable(), GetDecayTable(), GetDecayTime(), GetDecayTimeBin(), G4DynamicParticle::GetDefinition(), G4Track::GetDynamicParticle(), G4Track::GetGlobalTime(), G4IonTable::GetIon(), G4ParticleTable::GetIonTable(), G4DynamicParticle::GetKineticEnergy(), G4Track::GetLocalTime(), G4VPhysicalVolume::GetLogicalVolume(), G4DynamicParticle::GetMomentumDirection(), G4LogicalVolume::GetName(), G4DynamicParticle::GetParticleDefinition(), G4ParticleDefinition::GetParticleName(), G4ParticleTable::GetParticleTable(), G4ParticleDefinition::GetPDGLifeTime(), G4ParticleDefinition::GetPDGMass(), G4Track::GetPosition(), G4Track::GetTouchableHandle(), G4Track::GetTrackStatus(), GetVerboseLevel(), G4Track::GetVolume(), G4Track::GetWeight(), G4ParticleChangeForDecay::Initialize(), isAllVolumesMode, IsApplicable(), G4DecayProducts::IsChecked(), IsRateTableReady(), n, ns, NSplit, G4DecayProducts::PopProducts(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForDecay::ProposeLocalTime(), G4VParticleChange::ProposeTrackStatus(), s, G4Track::SetGoodForTrackingFlag(), G4VParticleChange::SetNumberOfSecondaries(), G4Track::SetTouchableHandle(), G4Track::SetWeight(), theDecayRateVector, theRadioactivityTables, and ValidVolumes.

Referenced by AtRestDoIt(), and PostStepDoIt().

+ Here is the caller graph for this function:

void G4RadioactiveDecay::DeselectAllVolumes ( )

Definition at line 340 of file G4RadioactiveDecay.cc.

References G4cout, G4endl, GetVerboseLevel(), isAllVolumesMode, and ValidVolumes.

+ Here is the call graph for this function:

void G4RadioactiveDecay::DeselectAVolume ( const G4String  aVolume)

Definition at line 285 of file G4RadioactiveDecay.cc.

References G4cerr, G4cout, G4endl, G4LogicalVolumeStore::GetInstance(), G4LogicalVolume::GetName(), GetVerboseLevel(), isAllVolumesMode, and ValidVolumes.

+ Here is the call graph for this function:

G4DecayProducts * G4RadioactiveDecay::DoDecay ( const G4ParticleDefinition theParticleDef)
protected

Definition at line 1903 of file G4RadioactiveDecay.cc.

References CollimateDecay(), G4VDecayChannel::DecayIt(), FatalException, G4cerr, G4cout, G4endl, G4Exception(), GetDecayTable(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGMass(), GetVerboseLevel(), MeV, and G4DecayTable::SelectADecayChannel().

Referenced by DecayIt().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

const G4ThreeVector& G4RadioactiveDecay::GetDecayDirection ( ) const
inline

Definition at line 215 of file G4RadioactiveDecay.hh.

References forceDecayDirection.

G4double G4RadioactiveDecay::GetDecayHalfAngle ( ) const
inline

Definition at line 223 of file G4RadioactiveDecay.hh.

References forceDecayHalfAngle.

void G4RadioactiveDecay::GetDecayRateTable ( const G4ParticleDefinition aParticle)

Definition at line 366 of file G4RadioactiveDecay.cc.

References G4cout, G4endl, G4ParticleDefinition::GetParticleName(), GetVerboseLevel(), theDecayRateTableVector, and theDecayRateVector.

Referenced by DecayIt().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

G4DecayTable * G4RadioactiveDecay::GetDecayTable ( const G4ParticleDefinition aNucleus)

Definition at line 245 of file G4RadioactiveDecay.cc.

References dkmap, G4ParticleDefinition::GetParticleName(), and LoadDecayTable().

Referenced by AddDecayRateTable(), DecayIt(), and DoDecay().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

G4double G4RadioactiveDecay::GetDecayTime ( )
protected

Definition at line 558 of file G4RadioactiveDecay.cc.

References DBin, DProfile, G4cout, G4endl, G4Exception(), G4UniformRand, GetVerboseLevel(), JustWarning, and s.

Referenced by DecayIt().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

G4int G4RadioactiveDecay::GetDecayTimeBin ( const G4double  aDecayTime)
protected

Definition at line 586 of file G4RadioactiveDecay.cc.

References DBin, G4endl, G4Exception(), and JustWarning.

Referenced by DecayIt().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

G4double G4RadioactiveDecay::GetMeanFreePath ( const G4Track theTrack,
G4double  previousStepSize,
G4ForceCondition condition 
)
protectedvirtual
G4double G4RadioactiveDecay::GetMeanLifeTime ( const G4Track theTrack,
G4ForceCondition condition 
)
protectedvirtual

Implements G4VRestDiscreteProcess.

Definition at line 611 of file G4RadioactiveDecay.cc.

References AnalogueMC, DBL_MAX, G4cout, G4endl, G4DynamicParticle::GetDefinition(), G4Track::GetDynamicParticle(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMass(), G4ParticleDefinition::GetPDGLifeTime(), G4ParticleDefinition::GetPDGStable(), GetVerboseLevel(), GeV, ns, and s.

+ Here is the call graph for this function:

G4NucleusLimits G4RadioactiveDecay::GetNucleusLimits ( ) const
inline

Definition at line 177 of file G4RadioactiveDecay.hh.

References theNucleusLimits.

G4int G4RadioactiveDecay::GetSplitNuclei ( )
inline

Definition at line 208 of file G4RadioactiveDecay.hh.

References NSplit.

std::vector<G4RadioactivityTable*> G4RadioactiveDecay::GetTheRadioactivityTables ( )
inline

Definition at line 154 of file G4RadioactiveDecay.hh.

References theRadioactivityTables.

Referenced by Run::WriteActivity().

+ Here is the caller graph for this function:

G4int G4RadioactiveDecay::GetVerboseLevel ( ) const
inline
G4bool G4RadioactiveDecay::IsAnalogueMonteCarlo ( )
inline

Definition at line 192 of file G4RadioactiveDecay.hh.

References AnalogueMC.

Referenced by Run::WriteActivity().

+ Here is the caller graph for this function:

G4bool G4RadioactiveDecay::IsApplicable ( const G4ParticleDefinition aParticle)
virtual

Reimplemented from G4VProcess.

Definition at line 223 of file G4RadioactiveDecay.cc.

References A(), G4NucleusLimits::GetAMax(), G4NucleusLimits::GetAMin(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetParticleType(), G4ParticleDefinition::GetPDGLifeTime(), G4NucleusLimits::GetZMax(), G4NucleusLimits::GetZMin(), and theNucleusLimits.

Referenced by AddDecayRateTable(), and DecayIt().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

G4bool G4RadioactiveDecay::IsRateTableReady ( const G4ParticleDefinition aParticle)

Definition at line 351 of file G4RadioactiveDecay.cc.

References G4ParticleDefinition::GetParticleName(), and theDecayRateTableVector.

Referenced by DecayIt().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

G4DecayTable * G4RadioactiveDecay::LoadDecayTable ( const G4ParticleDefinition theParentNucleus)

Definition at line 748 of file G4RadioactiveDecay.cc.

References a, A(), allowed, Alpha, applyARM, BDNeutron, BDProton, Beta2Minus, Beta2Plus, BetaMinus, BetaPlus, G4DecayTable::DumpInfo(), G4DecayTable::entries(), FatalException, G4Ions::FloatLevelBase(), G4cout, G4endl, G4Exception(), G4VDecayChannel::GetBR(), G4DecayTable::GetDecayChannel(), G4NuclearDecay::GetDecayMode(), G4ParticleDefinition::GetParticleName(), GetVerboseLevel(), halflifethreshold, G4DecayTable::Insert(), IT, JustWarning, keV, KshellEC, levelTolerance, LshellEC, MeV, MshellEC, Neutron, Neutron2, noFloat, Proton, Proton2, RDM_ERROR, G4ECDecay::SetARM(), G4ITDecay::SetARM(), G4VDecayChannel::SetBR(), G4NuclearDecay::SetHLThreshold(), SpFission, G4String::strip(), and theUserRadioactiveDataFiles.

Referenced by GetDecayTable().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

G4RadioactiveDecay& G4RadioactiveDecay::operator= ( const G4RadioactiveDecay right)
private
G4VParticleChange* G4RadioactiveDecay::PostStepDoIt ( const G4Track theTrack,
const G4Step theStep 
)
inlineprivatevirtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 345 of file G4RadioactiveDecay.hh.

References DecayIt().

+ Here is the call graph for this function:

void G4RadioactiveDecay::SelectAllVolumes ( )

Definition at line 316 of file G4RadioactiveDecay.cc.

References G4cout, G4endl, G4LogicalVolumeStore::GetInstance(), G4LogicalVolume::GetName(), GetVerboseLevel(), isAllVolumesMode, and ValidVolumes.

Referenced by G4RadioactiveDecay().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void G4RadioactiveDecay::SelectAVolume ( const G4String  aVolume)

Definition at line 262 of file G4RadioactiveDecay.cc.

References G4cerr, G4cout, G4endl, G4LogicalVolumeStore::GetInstance(), G4LogicalVolume::GetName(), GetVerboseLevel(), and ValidVolumes.

+ Here is the call graph for this function:

void G4RadioactiveDecay::SetAnalogueMonteCarlo ( G4bool  r)
inline

Definition at line 182 of file G4RadioactiveDecay.hh.

References AnalogueMC, halflifethreshold, and s.

Referenced by SetBRBias(), SetDecayBias(), SetSourceTimeProfile(), and SetSplitNuclei().

+ Here is the caller graph for this function:

void G4RadioactiveDecay::SetARM ( G4bool  arm)
inline

Definition at line 127 of file G4RadioactiveDecay.hh.

References applyARM.

Referenced by PhysicsList::AddRadioactiveDecay().

+ Here is the caller graph for this function:

void G4RadioactiveDecay::SetBRBias ( G4bool  r)
inline

Definition at line 196 of file G4RadioactiveDecay.hh.

References BRBias, and SetAnalogueMonteCarlo().

+ Here is the call graph for this function:

void G4RadioactiveDecay::SetDecayBias ( G4String  filename)

Definition at line 1500 of file G4RadioactiveDecay.cc.

References DBin, decayWindows, DProfile, FatalException, G4cout, G4endl, G4Exception(), GetVerboseLevel(), JustWarning, NDecayBin, s, SetAnalogueMonteCarlo(), and theRadioactivityTables.

+ Here is the call graph for this function:

void G4RadioactiveDecay::SetDecayCollimation ( const G4ThreeVector theDir,
G4double  halfAngle = 0.*CLHEP::deg 
)
inline

Definition at line 225 of file G4RadioactiveDecay.hh.

References SetDecayDirection(), and SetDecayHalfAngle().

+ Here is the call graph for this function:

void G4RadioactiveDecay::SetDecayDirection ( const G4ThreeVector theDir)
inline

Definition at line 211 of file G4RadioactiveDecay.hh.

References forceDecayDirection.

Referenced by SetDecayCollimation().

+ Here is the caller graph for this function:

void G4RadioactiveDecay::SetDecayHalfAngle ( G4double  halfAngle = 0.*CLHEP::deg)
inline

Definition at line 219 of file G4RadioactiveDecay.hh.

References deg, forceDecayHalfAngle, G4INCL::Math::max(), and G4INCL::Math::min().

Referenced by SetDecayCollimation().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void G4RadioactiveDecay::SetDecayRate ( G4int  theZ,
G4int  theA,
G4double  theE,
G4int  theG,
std::vector< G4double theCoefficients,
std::vector< G4double theTaos 
)

Definition at line 1122 of file G4RadioactiveDecay.cc.

References G4RadioactiveDecayRate::SetA(), G4RadioactiveDecayRate::SetDecayRateC(), G4RadioactiveDecayRate::SetE(), G4RadioactiveDecayRate::SetGeneration(), G4RadioactiveDecayRate::SetTaos(), G4RadioactiveDecayRate::SetZ(), and theDecayRate.

Referenced by AddDecayRateTable().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void G4RadioactiveDecay::SetFBeta ( G4bool  r)
inline

Definition at line 189 of file G4RadioactiveDecay.hh.

References FBeta.

void G4RadioactiveDecay::SetHLThreshold ( G4double  hl)
inline

Definition at line 121 of file G4RadioactiveDecay.hh.

References halflifethreshold.

void G4RadioactiveDecay::SetICM ( G4bool  icm)
inline

Definition at line 124 of file G4RadioactiveDecay.hh.

References applyICM.

void G4RadioactiveDecay::SetNucleusLimits ( G4NucleusLimits  theNucleusLimits1)
inline

Definition at line 172 of file G4RadioactiveDecay.hh.

References theNucleusLimits.

void G4RadioactiveDecay::SetSourceTimeProfile ( G4String  filename)

Definition at line 1450 of file G4RadioactiveDecay.cc.

References FatalException, G4cout, G4endl, G4Exception(), GetVerboseLevel(), JustWarning, NSourceBin, s, SBin, SetAnalogueMonteCarlo(), and SProfile.

+ Here is the call graph for this function:

void G4RadioactiveDecay::SetSplitNuclei ( G4int  r)
inline

Definition at line 202 of file G4RadioactiveDecay.hh.

References NSplit, and SetAnalogueMonteCarlo().

+ Here is the call graph for this function:

void G4RadioactiveDecay::SetVerboseLevel ( G4int  value)
inline

Definition at line 166 of file G4RadioactiveDecay.hh.

References verboseLevel.

Member Data Documentation

G4bool G4RadioactiveDecay::AnalogueMC
private
G4bool G4RadioactiveDecay::applyARM
private

Definition at line 285 of file G4RadioactiveDecay.hh.

Referenced by G4RadioactiveDecay(), LoadDecayTable(), and SetARM().

G4bool G4RadioactiveDecay::applyICM
private

Definition at line 284 of file G4RadioactiveDecay.hh.

Referenced by G4RadioactiveDecay(), and SetICM().

G4bool G4RadioactiveDecay::BRBias
private

Definition at line 279 of file G4RadioactiveDecay.hh.

Referenced by DecayIt(), G4RadioactiveDecay(), and SetBRBias().

G4double G4RadioactiveDecay::DBin[100]
private
G4int G4RadioactiveDecay::decayWindows[100]
private

Definition at line 309 of file G4RadioactiveDecay.hh.

Referenced by DecayIt(), G4RadioactiveDecay(), and SetDecayBias().

DecayTableMap* G4RadioactiveDecay::dkmap
private

Definition at line 316 of file G4RadioactiveDecay.hh.

Referenced by G4RadioactiveDecay(), GetDecayTable(), and ~G4RadioactiveDecay().

G4double G4RadioactiveDecay::DProfile[100]
private

Definition at line 297 of file G4RadioactiveDecay.hh.

Referenced by DecayIt(), G4RadioactiveDecay(), GetDecayTime(), and SetDecayBias().

G4bool G4RadioactiveDecay::FBeta
private

Definition at line 280 of file G4RadioactiveDecay.hh.

Referenced by G4RadioactiveDecay(), and SetFBeta().

G4ThreeVector G4RadioactiveDecay::forceDecayDirection
private
G4double G4RadioactiveDecay::forceDecayHalfAngle
private
G4ParticleChangeForRadDecay G4RadioactiveDecay::fParticleChangeForRadDecay
private

Definition at line 327 of file G4RadioactiveDecay.hh.

Referenced by DecayIt(), and G4RadioactiveDecay().

G4double G4RadioactiveDecay::fRemainderLifeTime
private

Definition at line 322 of file G4RadioactiveDecay.hh.

Referenced by AtRestGetPhysicalInteractionLength().

G4double G4RadioactiveDecay::halflifethreshold
private
bool G4RadioactiveDecay::isAllVolumesMode
private
G4bool G4RadioactiveDecay::isInitialised
private

Definition at line 277 of file G4RadioactiveDecay.hh.

Referenced by BuildPhysicsTable().

const G4double G4RadioactiveDecay::levelTolerance = 0.1*keV
staticprivate

Definition at line 310 of file G4RadioactiveDecay.hh.

Referenced by AddDecayRateTable(), and LoadDecayTable().

G4int G4RadioactiveDecay::NDecayBin
private

Definition at line 295 of file G4RadioactiveDecay.hh.

Referenced by G4RadioactiveDecay(), and SetDecayBias().

G4int G4RadioactiveDecay::NSourceBin
private
G4int G4RadioactiveDecay::NSplit
private

Definition at line 281 of file G4RadioactiveDecay.hh.

Referenced by DecayIt(), G4RadioactiveDecay(), GetSplitNuclei(), and SetSplitNuclei().

const G4ThreeVector G4RadioactiveDecay::origin
staticprivate
G4double G4RadioactiveDecay::SBin[100]
private
G4double G4RadioactiveDecay::SProfile[100]
private
G4RadioactiveDecayRate G4RadioactiveDecay::theDecayRate
private

Definition at line 302 of file G4RadioactiveDecay.hh.

Referenced by AddDecayRateTable(), and SetDecayRate().

G4RadioactiveDecayRateVector G4RadioactiveDecay::theDecayRateTable
private

Definition at line 304 of file G4RadioactiveDecay.hh.

Referenced by AddDecayRateTable().

G4RadioactiveDecayRateTable G4RadioactiveDecay::theDecayRateTableVector
private

Definition at line 305 of file G4RadioactiveDecay.hh.

Referenced by AddDecayRateTable(), GetDecayRateTable(), and IsRateTableReady().

G4RadioactiveDecayRates G4RadioactiveDecay::theDecayRateVector
private

Definition at line 303 of file G4RadioactiveDecay.hh.

Referenced by AddDecayRateTable(), DecayIt(), and GetDecayRateTable().

G4NucleusLimits G4RadioactiveDecay::theNucleusLimits
private

Definition at line 275 of file G4RadioactiveDecay.hh.

Referenced by GetNucleusLimits(), IsApplicable(), and SetNucleusLimits().

G4RadioactiveDecaymessenger* G4RadioactiveDecay::theRadioactiveDecaymessenger
private

Definition at line 273 of file G4RadioactiveDecay.hh.

Referenced by G4RadioactiveDecay(), and ~G4RadioactiveDecay().

std::vector<G4RadioactivityTable*> G4RadioactiveDecay::theRadioactivityTables
private
std::map<G4int, G4String> G4RadioactiveDecay::theUserRadioactiveDataFiles
private

Definition at line 313 of file G4RadioactiveDecay.hh.

Referenced by AddUserDecayDataFile(), G4RadioactiveDecay(), and LoadDecayTable().

std::vector<G4String> G4RadioactiveDecay::ValidVolumes
private
G4int G4RadioactiveDecay::verboseLevel
private

Definition at line 323 of file G4RadioactiveDecay.hh.

Referenced by GetVerboseLevel(), and SetVerboseLevel().


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