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

#include <G4RadioactiveDecay.hh>

+ Inheritance diagram for G4RadioactiveDecay:

Public Member Functions

 G4RadioactiveDecay (const G4String &processName="RadioactiveDecay", const G4double timeThreshold=-1.0)
 
 ~G4RadioactiveDecay () override
 
G4bool IsApplicable (const G4ParticleDefinition &) override
 
G4VParticleChangeAtRestDoIt (const G4Track &theTrack, const G4Step &theStep) override
 
G4VParticleChangePostStepDoIt (const G4Track &theTrack, const G4Step &theStep) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
void ProcessDescription (std::ostream &outFile) const override
 
virtual G4VParticleChangeDecayIt (const G4Track &theTrack, const G4Step &theStep)
 
G4DecayTableGetDecayTable (const G4ParticleDefinition *)
 
void SelectAVolume (const G4String &aVolume)
 
void DeselectAVolume (const G4String &aVolume)
 
void SelectAllVolumes ()
 
void DeselectAllVolumes ()
 
void SetARM (G4bool arm)
 
G4DecayTableLoadDecayTable (const G4Ions *)
 
void AddUserDecayDataFile (G4int Z, G4int A, const G4String &filename)
 
void SetNucleusLimits (G4NucleusLimits theNucleusLimits1)
 
G4NucleusLimits GetNucleusLimits () const
 
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 SetThresholdForVeryLongDecayTime (const G4double inputThreshold)
 
G4double GetThresholdForVeryLongDecayTime () const
 
void StreamInfo (std::ostream &os, const G4String &endline)
 
 G4RadioactiveDecay (const G4RadioactiveDecay &right)=delete
 
G4RadioactiveDecayoperator= (const G4RadioactiveDecay &right)=delete
 
- Public Member Functions inherited from G4VRestDiscreteProcess
 G4VRestDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VRestDiscreteProcess (G4VRestDiscreteProcess &)
 
virtual ~G4VRestDiscreteProcess ()
 
G4VRestDiscreteProcessoperator= (const G4VRestDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
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 ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool 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 const G4VProcessGetCreatorProcess () const
 
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

G4double GetMeanFreePath (const G4Track &theTrack, G4double previousStepSize, G4ForceCondition *condition) override
 
G4double GetMeanLifeTime (const G4Track &theTrack, G4ForceCondition *condition) override
 
void DecayAnalog (const G4Track &theTrack, G4DecayTable *)
 
G4DecayProductsDoDecay (const G4ParticleDefinition &, G4DecayTable *)
 
void CollimateDecay (G4DecayProducts *products)
 
void CollimateDecayProduct (G4DynamicParticle *product)
 
G4ThreeVector ChooseCollimationDirection () const
 
- Protected Member Functions inherited from G4VRestDiscreteProcess
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4ParticleChangeForRadDecay fParticleChangeForRadDecay
 
G4RadioactiveDecayMessengertheRadioactiveDecayMessenger
 
G4PhotonEvaporationphotonEvaporation
 
G4ITDecaydecayIT
 
std::vector< G4StringValidVolumes
 
bool isAllVolumesMode {true}
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Static Protected Attributes

static const G4double levelTolerance = 10.0*CLHEP::eV
 
static DecayTableMapmaster_dkmap = nullptr
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 

Detailed Description

Definition at line 65 of file G4RadioactiveDecay.hh.

Constructor & Destructor Documentation

◆ G4RadioactiveDecay() [1/2]

G4RadioactiveDecay::G4RadioactiveDecay ( const G4String & processName = "RadioactiveDecay",
const G4double timeThreshold = -1.0 )

Definition at line 111 of file G4RadioactiveDecay.cc.

113 : G4VRestDiscreteProcess(processName, fDecay),
114 fThresholdForVeryLongDecayTime( 1.0*CLHEP::year )
115{
116 if (GetVerboseLevel() > 1) {
117 G4cout << "G4RadioactiveDecay constructor: processName = " << processName
118 << G4endl;
119 }
120
122
125
126 // Check data directory
127 if (dirPath.empty()) {
128 const char* path_var = G4FindDataDir("G4RADIOACTIVEDATA");
129 if (nullptr == path_var) {
130 G4Exception("G4RadioactiveDecay()", "HAD_RDM_200", FatalException,
131 "Environment variable G4RADIOACTIVEDATA is not set");
132 } else {
133 dirPath = path_var; // convert to string
134 std::ostringstream os;
135 os << dirPath << "/z1.a3"; // used as a dummy
136 std::ifstream testFile;
137 testFile.open(os.str() );
138 if ( !testFile.is_open() )
139 G4Exception("G4RadioactiveDecay()","HAD_RDM_201",FatalException,
140 "Environment variable G4RADIOACTIVEDATA is set, but does not point to correct directory");
141 }
142 }
143 // Set up photon evaporation for use in G4ITDecay
148
149 // Instantiate the map of decay tables
150 if (nullptr == master_dkmap) {
152 }
153 if (nullptr == theUserRDataFiles) {
154 theUserRDataFiles = new std::map<G4int, G4String>;
155 }
156
157 // RDM applies to all logical volumes by default
160
161 // The time threshold for radioactive decays can be set in 3 ways:
162 // 1. Via C++ interface: G4HadronicParameters::Instance()->SetTimeThresholdForRadioactiveDecay(value)
163 // 2. Via the second parameter of the G4RadioactiveDecay constructor
164 // 3. Via UI command: /process/had/rdm/thresholdForVeryLongDecayTime value
165 // If both 1. and 2. are specified (at the moment when the G4RadioactiveDecay constructor is called),
166 // then we take the larger value, to be conservative.
167 // If, later on (after invoking the G4RadioactiveDecay constructor) 3. is specified,
168 // then this value is used (and the eventual values 1. and/or 2. are ignored).
170 if ( timeThreshold > 0.0 || timeThresholdBis > 0.0 ) {
171 if ( timeThreshold > timeThresholdBis ) timeThresholdBis = timeThreshold;
172 fThresholdForVeryLongDecayTime = timeThresholdBis;
173 }
174}
const char * G4FindDataDir(const char *)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
@ fRadioactiveDecay
@ fDecay
std::map< G4String, G4DecayTable * > DecayTableMap
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4HadronicParameters * Instance()
G4double GetTimeThresholdForRadioactiveDecay() const
void RegisterExtraProcess(G4VProcess *)
static G4HadronicProcessStore * Instance()
void RDMForced(G4bool) override
void SetICM(G4bool) override
G4RadioactiveDecayMessenger * theRadioactiveDecayMessenger
static DecayTableMap * master_dkmap
G4ParticleChangeForRadDecay fParticleChangeForRadDecay
G4PhotonEvaporation * photonEvaporation
G4int GetVerboseLevel() const
void SetProcessSubType(G4int)
G4VParticleChange * pParticleChange

◆ ~G4RadioactiveDecay()

G4RadioactiveDecay::~G4RadioactiveDecay ( )
override

Definition at line 201 of file G4RadioactiveDecay.cc.

202{
204 delete photonEvaporation;
205 delete decayIT;
206 if (nullptr != master_dkmap) {
207 G4AutoLock lk(&radioactiveDecayMutex);
208 if (nullptr != master_dkmap) {
209 for (auto const & i : *master_dkmap) {
210 delete i.second;
211 }
212 master_dkmap->clear();
213 delete master_dkmap;
214 master_dkmap = nullptr;
215 }
216 delete theUserRDataFiles;
217 theUserRDataFiles = nullptr;
218 lk.unlock();
219 }
220}

◆ G4RadioactiveDecay() [2/2]

G4RadioactiveDecay::G4RadioactiveDecay ( const G4RadioactiveDecay & right)
delete

Member Function Documentation

◆ AddUserDecayDataFile()

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

Definition at line 898 of file G4RadioactiveDecay.cc.

900{
901 if (Z < 1 || A < 2) G4cout << "Z and A not valid!" << G4endl;
902
903 std::ifstream DecaySchemeFile(filename);
904 if (DecaySchemeFile) {
905 G4int ID_ion = A*1000 + Z;
906 (*theUserRDataFiles)[ID_ion] = filename;
907 } else {
909 ed << filename << " does not exist! " << G4endl;
910 G4Exception("G4RadioactiveDecay::AddUserDecayDataFile()", "HAD_RDM_001",
911 FatalException, ed);
912 }
913}
std::ostringstream G4ExceptionDescription
int G4int
Definition G4Types.hh:85
const G4double A[17]

Referenced by G4RadioactiveDecayMessenger::SetNewValue().

◆ AtRestDoIt()

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

Reimplemented from G4VRestDiscreteProcess.

Definition at line 177 of file G4RadioactiveDecay.cc.

179{
180 return DecayIt(theTrack, theStep);
181}
virtual G4VParticleChange * DecayIt(const G4Track &theTrack, const G4Step &theStep)

◆ BuildPhysicsTable()

void G4RadioactiveDecay::BuildPhysicsTable ( const G4ParticleDefinition & p)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 446 of file G4RadioactiveDecay.cc.

447{
448 if (isInitialised) { return; }
449 isInitialised = true;
451 G4Threading::IsMasterThread() && "GenericIon" == p.GetParticleName()) {
452 StreamInfo(G4cout, "\n");
453 }
457 decayIT->SetARM(applyARM);
458
461}
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)
void PrintInfo(const G4ParticleDefinition *)
void SetARM(G4bool onoff)
Definition G4ITDecay.hh:65
const G4String & GetParticleName() const
void StreamInfo(std::ostream &os, const G4String &endline)
G4bool IsMasterThread()

◆ ChooseCollimationDirection()

G4ThreeVector G4RadioactiveDecay::ChooseCollimationDirection ( ) const
protected

Definition at line 1200 of file G4RadioactiveDecay.cc.

1200 {
1201 if (origin == forceDecayDirection) return origin; // Don't do collimation
1202 if (forceDecayHalfAngle == 180.*deg) return origin;
1203
1204 G4ThreeVector dir = forceDecayDirection;
1205
1206 // Return direction offset by random throw
1207 if (forceDecayHalfAngle > 0.) {
1208 // Generate uniform direction around central axis
1209 G4double phi = 2.*pi*G4UniformRand();
1210 G4double cosMin = std::cos(forceDecayHalfAngle);
1211 G4double cosTheta = (1.-cosMin)*G4UniformRand() + cosMin; // [cosMin,1.)
1212
1213 dir.setPhi(dir.phi()+phi);
1214 dir.setTheta(dir.theta()+std::acos(cosTheta));
1215 }
1216
1217#ifdef G4VERBOSE
1218 if (GetVerboseLevel()>1)
1219 G4cout << " ChooseCollimationDirection returns " << dir << G4endl;
1220#endif
1221
1222 return dir;
1223}
#define G4UniformRand()
Definition Randomize.hh:52
double phi() const
double theta() const
void setTheta(double)
void setPhi(double)

Referenced by CollimateDecayProduct().

◆ CollimateDecay()

void G4RadioactiveDecay::CollimateDecay ( G4DecayProducts * products)
protected

Definition at line 1153 of file G4RadioactiveDecay.cc.

1153 {
1154
1155 if (origin == forceDecayDirection) return; // No collimation requested
1156 if (180.*deg == forceDecayHalfAngle) return;
1157 if (0 == products || 0 == products->entries()) return;
1158
1159#ifdef G4VERBOSE
1160 if (GetVerboseLevel() > 1) G4cout << "Begin of CollimateDecay..." << G4endl;
1161#endif
1162
1163 // Particles suitable for directional biasing (for if-blocks below)
1167 static const G4ParticleDefinition* gamma = G4Gamma::Definition();
1171
1172 G4ThreeVector newDirection; // Re-use to avoid memory churn
1173 for (G4int i=0; i<products->entries(); ++i) {
1174 G4DynamicParticle* daughter = (*products)[i];
1175 const G4ParticleDefinition* daughterType =
1176 daughter->GetParticleDefinition();
1177 if (daughterType == electron || daughterType == positron ||
1178 daughterType == neutron || daughterType == gamma ||
1179 daughterType == alpha || daughterType == triton || daughterType == proton) {
1180 CollimateDecayProduct(daughter);
1181 }
1182 }
1183}
static G4Alpha * Definition()
Definition G4Alpha.cc:43
G4int entries() const
const G4ParticleDefinition * GetParticleDefinition() const
static G4Electron * Definition()
Definition G4Electron.cc:45
static G4Gamma * Definition()
Definition G4Gamma.cc:43
static G4Neutron * Definition()
Definition G4Neutron.cc:50
static G4Positron * Definition()
Definition G4Positron.cc:45
static G4Proton * Definition()
Definition G4Proton.cc:45
void CollimateDecayProduct(G4DynamicParticle *product)
static G4Triton * Definition()
Definition G4Triton.cc:45

Referenced by DoDecay().

◆ CollimateDecayProduct()

void G4RadioactiveDecay::CollimateDecayProduct ( G4DynamicParticle * product)
protected

Definition at line 1185 of file G4RadioactiveDecay.cc.

1185 {
1186#ifdef G4VERBOSE
1187 if (GetVerboseLevel() > 1) {
1188 G4cout << "CollimateDecayProduct for daughter "
1189 << daughter->GetParticleDefinition()->GetParticleName() << G4endl;
1190 }
1191#endif
1192
1194 if (origin != collimate) daughter->SetMomentumDirection(collimate);
1195}
G4ThreeVector ChooseCollimationDirection() const

Referenced by CollimateDecay().

◆ DecayAnalog()

void G4RadioactiveDecay::DecayAnalog ( const G4Track & theTrack,
G4DecayTable * decayTable )
protected

Definition at line 984 of file G4RadioactiveDecay.cc.

986{
987 const G4DynamicParticle* theParticle = theTrack.GetDynamicParticle();
988 const G4ParticleDefinition* theParticleDef = theParticle->GetDefinition();
989 //G4cout << "DecayIt for " << theParticleDef->GetParticleName() << G4endl;
990 G4DecayProducts* products = DoDecay(*theParticleDef, decayTable);
991
992 // Check if the product is the same as input and kill the track if
993 // necessary to prevent infinite loop (11/05/10, F.Lei)
994 if (nullptr == products || products->entries() == 1) {
999 delete products;
1000 return;
1001 }
1002
1003 G4double energyDeposit = 0.0;
1004 G4double finalGlobalTime = theTrack.GetGlobalTime();
1005 G4double finalLocalTime = theTrack.GetLocalTime();
1006
1007 // Get parent particle information and boost the decay products to the
1008 // laboratory frame
1009
1010 // ParentEnergy used for the boost should be the total energy of the nucleus
1011 // of the parent ion without the energy of the shell electrons
1012 // (correction for bug 1359 by L. Desorgher)
1013 G4double ParentEnergy = theParticle->GetKineticEnergy()
1014 + theParticle->GetParticleDefinition()->GetPDGMass();
1015 G4ThreeVector ParentDirection(theParticle->GetMomentumDirection());
1016
1017 if (theTrack.GetTrackStatus() == fStopButAlive) {
1018 // this condition seems to be always True, further investigation is needed (L.Desorgher)
1019
1020 // The particle is decayed at rest
1021 // Since the time is for the particle at rest, need to add additional time
1022 // lapsed between particle coming to rest and the actual decay. This time
1023 // is sampled with the mean-life of the particle. Need to protect the case
1024 // PDGTime < 0. (F.Lei 11/05/10)
1025 G4double temptime = -std::log(G4UniformRand() ) *
1026 theParticleDef->GetPDGLifeTime();
1027 if (temptime < 0.) temptime = 0.;
1028 finalGlobalTime += temptime;
1029 finalLocalTime += temptime;
1030 energyDeposit += theParticle->GetKineticEnergy();
1031
1032 // Kill the parent particle, and ignore its decay, if it decays later than the
1033 // threshold fThresholdForVeryLongDecayTime (whose default value corresponds
1034 // to more than twice the age of the universe).
1035 // This kind of cut has been introduced (in April 2021) in order to avoid to
1036 // account energy depositions happening after many billions of years in
1037 // ordinary materials used in calorimetry, in particular Tungsten and Lead
1038 // (via their natural unstable, but very long lived, isotopes, such as
1039 // W183, W180 and Pb204).
1040 // Note that the cut is not on the average, mean lifetime, but on the actual
1041 // sampled global decay time.
1042 if ( finalGlobalTime > fThresholdForVeryLongDecayTime ) {
1047 delete products;
1048 return;
1049 }
1050 }
1051 products->Boost(ParentEnergy, ParentDirection);
1052
1053 // Add products in theParticleChangeForRadDecay.
1054 G4int numberOfSecondaries = products->entries();
1056
1057 if (GetVerboseLevel() > 1) {
1058 G4cout << "G4RadioactiveDecay::DecayAnalog: Decay vertex :";
1059 G4cout << " Time: " << finalGlobalTime/ns << "[ns]";
1060 G4cout << " X:" << (theTrack.GetPosition()).x() /cm << "[cm]";
1061 G4cout << " Y:" << (theTrack.GetPosition()).y() /cm << "[cm]";
1062 G4cout << " Z:" << (theTrack.GetPosition()).z() /cm << "[cm]";
1063 G4cout << G4endl;
1064 G4cout << "G4Decay::DecayIt : decay products in Lab. Frame" << G4endl;
1065 products->DumpInfo();
1066 products->IsChecked();
1067 }
1068
1069 const G4int modelID_forIT = G4PhysicsModelCatalog::GetModelID( "model_RDM_IT" );
1070 G4int modelID = modelID_forIT + 10*theRadDecayMode;
1071 const G4int modelID_forAtomicRelaxation =
1072 G4PhysicsModelCatalog::GetModelID( "model_RDM_AtomicRelaxation" );
1073 for ( G4int index = 0; index < numberOfSecondaries; ++index ) {
1074 G4Track* secondary = new G4Track( products->PopProducts(), finalGlobalTime,
1075 theTrack.GetPosition() );
1076 secondary->SetWeight( theTrack.GetWeight() );
1077 secondary->SetCreatorModelID( modelID );
1078 // Change for atomics relaxation
1079 if ( theRadDecayMode == IT && index > 0 ) {
1080 if ( index == numberOfSecondaries-1 ) {
1081 secondary->SetCreatorModelID( modelID_forIT );
1082 } else {
1083 secondary->SetCreatorModelID( modelID_forAtomicRelaxation) ;
1084 }
1085 } else if ( theRadDecayMode >= KshellEC && theRadDecayMode <= NshellEC &&
1086 index < numberOfSecondaries-1 ) {
1087 secondary->SetCreatorModelID( modelID_forAtomicRelaxation );
1088 }
1089 secondary->SetGoodForTrackingFlag();
1090 secondary->SetTouchableHandle( theTrack.GetTouchableHandle() );
1092 }
1093
1094 delete products;
1095
1096 // Kill the parent particle
1100
1101 // Reset NumberOfInteractionLengthLeft.
1103}
@ fStopAndKill
@ fStopButAlive
void DumpInfo() const
G4DynamicParticle * PopProducts()
G4bool IsChecked() const
void Boost(G4double totalEnergy, const G4ThreeVector &momentumDirection)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetPDGLifeTime() const
static G4int GetModelID(const G4int modelIndex)
G4DecayProducts * DoDecay(const G4ParticleDefinition &, G4DecayTable *)
G4TrackStatus GetTrackStatus() const
G4double GetWeight() const
void SetWeight(G4double aValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
G4double GetLocalTime() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
void SetCreatorModelID(const G4int id)
void SetGoodForTrackingFlag(G4bool value=true)
void ProposeTrackStatus(G4TrackStatus status)
void AddSecondary(G4Track *aSecondary)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
void ClearNumberOfInteractionLengthLeft()

Referenced by G4Radioactivation::DecayIt(), and DecayIt().

◆ DecayIt()

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

Reimplemented in G4Radioactivation.

Definition at line 922 of file G4RadioactiveDecay.cc.

923{
924 // Initialize G4ParticleChange object, get particle details and decay table
927 const G4DynamicParticle* theParticle = theTrack.GetDynamicParticle();
928 const G4ParticleDefinition* theParticleDef = theParticle->GetDefinition();
929
930 // First check whether RDM applies to the current logical volume
931 if (!isAllVolumesMode) {
932 if (!std::binary_search(ValidVolumes.begin(), ValidVolumes.end(),
933 theTrack.GetVolume()->GetLogicalVolume()->GetName())) {
934#ifdef G4VERBOSE
935 if (GetVerboseLevel()>1) {
936 G4cout <<"G4RadioactiveDecay::DecayIt : "
937 << theTrack.GetVolume()->GetLogicalVolume()->GetName()
938 << " is not selected for the RDM"<< G4endl;
939 G4cout << " There are " << ValidVolumes.size() << " volumes" << G4endl;
940 G4cout << " The Valid volumes are: ";
941 for (auto const & vol : ValidVolumes) { G4cout << vol << " " << G4endl; }
942 G4cout << G4endl;
943 }
944#endif
946
947 // Kill the parent particle.
952 }
953 }
954
955 // Now check if particle is valid for RDM
956 G4DecayTable* theDecayTable = GetDecayTable(theParticleDef);
957 if ( theDecayTable == nullptr || theDecayTable->entries() == 0) {
958 // Particle is not an ion or is outside the nucleuslimits for decay
959#ifdef G4VERBOSE
960 if (GetVerboseLevel() > 1) {
961 G4cout << "G4RadioactiveDecay::DecayIt : "
962 << theParticleDef->GetParticleName()
963 << " is outside (Z,A) limits set for the decay or has no decays."
964 << G4endl;
965 }
966#endif
968
969 // Kill the parent particle
974 }
975 //G4cout << "DecayIt for " << theParticleDef->GetParticleName()
976 // << " isAllVolumesMode:" << isAllVolumesMode
977 // << " decayTable=" << theDecayTable << G4endl;
978
979 // Data found. Decay nucleus without variance reduction.
980 DecayAnalog(theTrack, theDecayTable);
982}
G4int entries() const
const G4String & GetName() const
void Initialize(const G4Track &) final
std::vector< G4String > ValidVolumes
void DecayAnalog(const G4Track &theTrack, G4DecayTable *)
G4DecayTable * GetDecayTable(const G4ParticleDefinition *)
G4VPhysicalVolume * GetVolume() const
void ProposeWeight(G4double finalWeight)
G4LogicalVolume * GetLogicalVolume() const

Referenced by AtRestDoIt(), and PostStepDoIt().

◆ DeselectAllVolumes()

void G4RadioactiveDecay::DeselectAllVolumes ( )

Definition at line 357 of file G4RadioactiveDecay.cc.

358{
359 ValidVolumes.clear();
360 isAllVolumesMode=false;
361#ifdef G4VERBOSE
362 if (GetVerboseLevel() > 1) G4cout << "RDM removed from all volumes" << G4endl;
363#endif
364}

Referenced by G4RadioactiveDecayMessenger::SetNewValue().

◆ DeselectAVolume()

void G4RadioactiveDecay::DeselectAVolume ( const G4String & aVolume)

Definition at line 298 of file G4RadioactiveDecay.cc.

299{
301 G4LogicalVolume* volume = nullptr;
302 volume = theLogicalVolumes->GetVolume(aVolume);
303 if (volume != nullptr)
304 {
305 auto location= std::find(ValidVolumes.cbegin(),ValidVolumes.cend(),aVolume);
306 if (location != ValidVolumes.cend() )
307 {
308 ValidVolumes.erase(location);
309 std::sort(ValidVolumes.begin(), ValidVolumes.end());
310 isAllVolumesMode = false;
311 if (GetVerboseLevel() > 0)
312 G4cout << " G4RadioactiveDecay::DeselectAVolume: " << aVolume
313 << " is removed from list " << G4endl;
314 }
315 else
316 {
318 ed << aVolume << " is not in the list. No action taken." << G4endl;
319 G4Exception("G4RadioactiveDecay::DeselectAVolume()", "HAD_RDM_300",
320 JustWarning, ed);
321 }
322 }
323 else
324 {
326 ed << aVolume << " is not a valid logical volume name. No action taken."
327 << G4endl;
328 G4Exception("G4RadioactiveDecay::DeselectAVolume()", "HAD_RDM_300",
329 JustWarning, ed);
330 }
331}
@ JustWarning
G4LogicalVolume * GetVolume(const G4String &name, G4bool verbose=true, G4bool reverseSearch=false) const
static G4LogicalVolumeStore * GetInstance()

Referenced by G4RadioactiveDecayMessenger::SetNewValue().

◆ DoDecay()

G4DecayProducts * G4RadioactiveDecay::DoDecay ( const G4ParticleDefinition & theParticleDef,
G4DecayTable * theDecayTable )
protected

Definition at line 1107 of file G4RadioactiveDecay.cc.

1109{
1110 G4DecayProducts* products = nullptr;
1111
1112 // Choose a decay channel.
1113 // G4DecayTable::SelectADecayChannel checks to see if sum of daughter masses
1114 // exceeds parent mass. Pass it the parent mass + maximum Q value to account
1115 // for difference in mass defect.
1116 G4double parentPlusQ = theParticleDef.GetPDGMass() + 30.*MeV;
1117 G4VDecayChannel* theDecayChannel = theDecayTable->SelectADecayChannel(parentPlusQ);
1118
1119 if (theDecayChannel == nullptr) {
1120 // Decay channel not found.
1122 ed << " Cannot determine decay channel for " << theParticleDef.GetParticleName() << G4endl;
1123 G4Exception("G4RadioactiveDecay::DoDecay", "HAD_RDM_013",
1124 FatalException, ed);
1125 } else {
1126 // A decay channel has been identified, so execute the DecayIt.
1127#ifdef G4VERBOSE
1128 if (GetVerboseLevel() > 1) {
1129 G4cout << "G4RadioactiveDecay::DoIt : selected decay channel addr: "
1130 << theDecayChannel << G4endl;
1131 }
1132#endif
1133 theRadDecayMode = static_cast<G4NuclearDecay*>(theDecayChannel)->GetDecayMode();
1134
1135 // for IT decay use local G4ITDecay class
1136 if (theRadDecayMode == IT) {
1137 decayIT->SetupDecay(&theParticleDef);
1138 products = decayIT->DecayIt(0.0);
1139 } else {
1140 // for others decayes use shared class
1141 products = theDecayChannel->DecayIt(theParticleDef.GetPDGMass());
1142 }
1143
1144 // Apply directional bias if requested by user
1145 CollimateDecay(products);
1146 }
1147
1148 return products;
1149}
G4VDecayChannel * SelectADecayChannel(G4double parentMass=-1.)
void SetupDecay(const G4ParticleDefinition *)
Definition G4ITDecay.cc:67
G4DecayProducts * DecayIt(G4double) override
Definition G4ITDecay.cc:74
void CollimateDecay(G4DecayProducts *products)
virtual G4DecayProducts * DecayIt(G4double parentMass=-1.0)=0

Referenced by DecayAnalog(), and G4Radioactivation::DecayIt().

◆ GetDecayDirection()

const G4ThreeVector & G4RadioactiveDecay::GetDecayDirection ( ) const
inline

Definition at line 139 of file G4RadioactiveDecay.hh.

139 {
140 return forceDecayDirection;
141 }

◆ GetDecayHalfAngle()

G4double G4RadioactiveDecay::GetDecayHalfAngle ( ) const
inline

Definition at line 147 of file G4RadioactiveDecay.hh.

147{return forceDecayHalfAngle;}

◆ GetDecayTable()

G4DecayTable * G4RadioactiveDecay::GetDecayTable ( const G4ParticleDefinition * aNucleus)

Definition at line 253 of file G4RadioactiveDecay.cc.

254{
255 G4String key = aNucleus->GetParticleName();
256 auto ptr = master_dkmap->find(key);
257
258 G4DecayTable* theDecayTable = nullptr;
259 if ( ptr == master_dkmap->end() ) {
260 // Load new file if table not there
261 const G4Ions* ion = dynamic_cast<const G4Ions*>(aNucleus);
262 if (nullptr != ion) {
263 theDecayTable = LoadDecayTable(ion);
264 }
265 } else {
266 theDecayTable = ptr->second;
267 }
268 return theDecayTable;
269}
G4DecayTable * LoadDecayTable(const G4Ions *)

Referenced by G4Radioactivation::CalculateChainsFromParent(), G4Radioactivation::DecayIt(), and DecayIt().

◆ GetMeanFreePath()

G4double G4RadioactiveDecay::GetMeanFreePath ( const G4Track & theTrack,
G4double previousStepSize,
G4ForceCondition * condition )
overrideprotectedvirtual

Implements G4VRestDiscreteProcess.

Definition at line 416 of file G4RadioactiveDecay.cc.

418{
419 G4double res = DBL_MAX;
420 G4double lifeTime = GetMeanLifeTime(aTrack, fc);
421 if (lifeTime > 0.0 && lifeTime < DBL_MAX) {
422 auto dParticle = aTrack.GetDynamicParticle();
423 res = lifeTime*dParticle->GetTotalEnergy()*aTrack.GetVelocity()/dParticle->GetMass();
424 } else {
425 res = lifeTime;
426 }
427
428#ifdef G4VERBOSE
429 if (GetVerboseLevel() > 2) {
430 G4cout << "G4RadioactiveDecay::GetMeanFreePath() for "
431 << aTrack.GetDefinition()->GetParticleName() << G4endl;
432 G4cout << " kinEnergy(GeV)=" << aTrack.GetKineticEnergy()/CLHEP::GeV
433 << " lifeTime(ns)=" << lifeTime
434 << " mean free path(cm)=" << res/CLHEP::cm << G4endl;
435 }
436#endif
437 return res;
438}
G4double GetMeanLifeTime(const G4Track &theTrack, G4ForceCondition *condition) override
#define DBL_MAX
Definition templates.hh:62

◆ GetMeanLifeTime()

G4double G4RadioactiveDecay::GetMeanLifeTime ( const G4Track & theTrack,
G4ForceCondition * condition )
overrideprotectedvirtual

Implements G4VRestDiscreteProcess.

Definition at line 373 of file G4RadioactiveDecay.cc.

375{
376 G4double meanlife = DBL_MAX;
377 const G4ParticleDefinition* theParticleDef = theTrack.GetParticleDefinition();
378 if (!IsApplicable(*theParticleDef)) { return meanlife; }
379 G4double theLife = theParticleDef->GetPDGLifeTime();
380#ifdef G4VERBOSE
381 if (GetVerboseLevel() > 2) {
382 G4cout << "G4RadioactiveDecay::GetMeanLifeTime() for "
383 << theParticleDef->GetParticleName() << G4endl;
384 G4cout << "KineticEnergy(GeV)=" << theTrack.GetKineticEnergy()/CLHEP::GeV
385 << " Mass(GeV)=" << theParticleDef->GetPDGMass()/CLHEP::GeV
386 << " LifeTime(ns)=" << theLife/CLHEP::ns << G4endl;
387 }
388#endif
389 if (theLife >= 0.0 && theLife <= fThresholdForVeryLongDecayTime) {
390 meanlife = theLife;
391 }
392
393 if (meanlife == DBL_MAX) {
394 const G4Ions* ion = dynamic_cast<const G4Ions*>(theParticleDef);
395 if (nullptr != ion && ion->GetExcitationEnergy() > 0.0) {
396 meanlife = 0.0;
397 }
398 }
399
400#ifdef G4VERBOSE
401 if (GetVerboseLevel() > 2)
402 G4cout << "G4RadioactiveDecay::GetMeanLifeTime: "
403 << meanlife/CLHEP::s << " second " << G4endl;
404#endif
405
406 return meanlife;
407}
G4double GetExcitationEnergy() const
Definition G4Ions.hh:149
G4bool IsApplicable(const G4ParticleDefinition &) override
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const

Referenced by GetMeanFreePath(), and G4Radioactivation::GetMeanLifeTime().

◆ GetNucleusLimits()

G4NucleusLimits G4RadioactiveDecay::GetNucleusLimits ( ) const
inline

Definition at line 133 of file G4RadioactiveDecay.hh.

133{return theNucleusLimits;}

◆ GetThresholdForVeryLongDecayTime()

G4double G4RadioactiveDecay::GetThresholdForVeryLongDecayTime ( ) const
inline

Definition at line 161 of file G4RadioactiveDecay.hh.

161{return fThresholdForVeryLongDecayTime;}

◆ IsApplicable()

G4bool G4RadioactiveDecay::IsApplicable ( const G4ParticleDefinition & aParticle)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 223 of file G4RadioactiveDecay.cc.

224{
225 const G4String& pname = aParticle.GetParticleName();
226 if (pname == "GenericIon" || pname == "triton") { return true; }
227 // All particles other than G4Ions, are rejected by default
228 const G4Ions* p = dynamic_cast<const G4Ions*>(&aParticle);
229 if (nullptr == p) { return false; }
230
231 // excited isomere may decay via gamma evaporation
232 if (p->GetExcitationEnergy() > 0.0) { return true; }
233
234 // Check on life time
235 G4double lifeTime = p->GetPDGLifeTime();
236 if (lifeTime < 0.0 || lifeTime > fThresholdForVeryLongDecayTime) {
237 return false;
238 }
239
240 // Determine whether the nuclide falls into the correct A and Z range
241 G4int A = p->GetAtomicMass();
242 G4int Z = p->GetAtomicNumber();
243
244 if (A > theNucleusLimits.GetAMax() || A < theNucleusLimits.GetAMin() ||
245 Z > theNucleusLimits.GetZMax() || Z < theNucleusLimits.GetZMin()) {
246 return false;
247 }
248
249 return true;
250}
G4int GetZMax() const
G4int GetZMin() const
G4int GetAMin() const
G4int GetAMax() const
G4int GetAtomicNumber() const
G4int GetAtomicMass() const

Referenced by G4Radioactivation::CalculateChainsFromParent(), G4Radioactivation::DecayIt(), and GetMeanLifeTime().

◆ LoadDecayTable()

G4DecayTable * G4RadioactiveDecay::LoadDecayTable ( const G4Ions * theIon)

Definition at line 523 of file G4RadioactiveDecay.cc.

524{
525 G4AutoLock lk(&radioactiveDecayMutex);
526 const G4String key = theIon->GetParticleName();
527 auto dtptr = master_dkmap->find(key);
528 if (dtptr != master_dkmap->end()) {
529 lk.unlock();
530 return dtptr->second;
531 }
532
533 // Generate input data file name using Z and A of the parent nucleus
534 // file containing radioactive decay data.
535 G4int A = theIon->GetAtomicMass();
536 G4int Z = theIon->GetAtomicNumber();
537
538 //G4cout << "LoadDecayTable for " << key << " Z=" << Z << " A=" << A << G4endl;
539
540 G4double levelEnergy = theIon->GetExcitationEnergy();
541 G4Ions::G4FloatLevelBase floatingLevel = theIon->GetFloatLevelBase();
542
543 //Check if data have been provided by the user
544 G4String file;
545 G4int ke = 1000*A + Z;
546 auto ptr = theUserRDataFiles->find(ke);
547 if (ptr != theUserRDataFiles->end()) {
548 file = ptr->second;
549 } else {
550 std::ostringstream os;
551 os << dirPath << "/z" << Z << ".a" << A << '\0';
552 file = os.str();
553 }
554
555 G4DecayTable* theDecayTable = new G4DecayTable();
556 G4bool found(false); // True if energy level matches one in table
557
558 std::ifstream DecaySchemeFile;
559 DecaySchemeFile.open(file);
560
561 if (DecaySchemeFile.good()) {
562 // Initialize variables used for reading in radioactive decay data
563 G4bool floatMatch(false);
564 const G4int nMode = G4RadioactiveDecayModeSize;
565 G4double modeTotalBR[nMode] = {0.0};
566 G4double modeSumBR[nMode];
567 for (G4int i = 0; i < nMode; ++i) {
568 modeSumBR[i] = 0.0;
569 }
570
571 char inputChars[120]={' '};
572 G4String inputLine;
573 G4String recordType("");
574 G4String floatingFlag("");
575 G4String daughterFloatFlag("");
576 G4Ions::G4FloatLevelBase daughterFloatLevel;
577 G4RadioactiveDecayMode theDecayMode;
578 G4double decayModeTotal(0.0);
579 G4double parentExcitation(0.0);
580 G4double a(0.0);
581 G4double b(0.0);
582 G4double c(0.0);
583 G4double dummy(0.0);
584 G4BetaDecayType betaType(allowed);
585
586 // Loop through each data file record until you identify the decay
587 // data relating to the nuclide of concern.
588
589 G4bool complete(false); // bool insures only one set of values read for any
590 // given parent energy level
591 G4int loop = 0;
592 /* Loop checking, 01.09.2015, D.Wright */
593 while (!complete && !DecaySchemeFile.getline(inputChars, 120).eof()) {
594 loop++;
595 if (loop > 100000) {
596 G4Exception("G4RadioactiveDecay::LoadDecayTable()", "HAD_RDM_100",
597 JustWarning, "While loop count exceeded");
598 break;
599 }
600
601 inputLine = inputChars;
602 G4StrUtil::rstrip(inputLine);
603 if (inputChars[0] != '#' && inputLine.length() != 0) {
604 std::istringstream tmpStream(inputLine);
605
606 if (inputChars[0] == 'P') {
607 // Nucleus is a parent type. Check excitation level to see if it
608 // matches that of theParentNucleus
609 tmpStream >> recordType >> parentExcitation >> floatingFlag >> dummy;
610 // "dummy" takes the place of half-life
611 // Now read in from ENSDFSTATE in particle category
612
613 if (found) {
614 complete = true;
615 } else {
616 // Take first level which matches excitation energy regardless of floating level
617 found = (std::abs(parentExcitation*keV - levelEnergy) < levelTolerance);
618 if (floatingLevel != noFloat) {
619 // If floating level specificed, require match of both energy and floating level
620 floatMatch = (floatingLevel == G4Ions::FloatLevelBase(floatingFlag.back()) );
621 if (!floatMatch) found = false;
622 }
623 }
624
625 } else if (found) {
626 // The right part of the radioactive decay data file has been found. Search
627 // through it to determine the mode of decay of the subsequent records.
628
629 // Store for later the total decay probability for each decay mode
630 if (inputLine.length() < 72) {
631 tmpStream >> theDecayMode >> dummy >> decayModeTotal;
632
633 switch (theDecayMode) {
634 case IT:
635 {
636 G4ITDecay* anITChannel = new G4ITDecay(theIon, decayModeTotal, 0.0, 0.0);
637 theDecayTable->Insert(anITChannel);
638 }
639 break;
640 case BetaMinus:
641 modeTotalBR[BetaMinus] = decayModeTotal; break;
642 case BetaPlus:
643 modeTotalBR[BetaPlus] = decayModeTotal; break;
644 case KshellEC:
645 modeTotalBR[KshellEC] = decayModeTotal; break;
646 case LshellEC:
647 modeTotalBR[LshellEC] = decayModeTotal; break;
648 case MshellEC:
649 modeTotalBR[MshellEC] = decayModeTotal; break;
650 case NshellEC:
651 modeTotalBR[NshellEC] = decayModeTotal; break;
652 case Alpha:
653 modeTotalBR[Alpha] = decayModeTotal; break;
654 case Proton:
655 modeTotalBR[Proton] = decayModeTotal; break;
656 case Neutron:
657 modeTotalBR[Neutron] = decayModeTotal; break;
658 case SpFission:
659 modeTotalBR[SpFission] = decayModeTotal; break;
660 case BDProton:
661 /* Not yet implemented */ break;
662 case BDNeutron:
663 /* Not yet implemented */ break;
664 case Beta2Minus:
665 /* Not yet implemented */ break;
666 case Beta2Plus:
667 /* Not yet implemented */ break;
668 case Proton2:
669 /* Not yet implemented */ break;
670 case Neutron2:
671 /* Not yet implemented */ break;
672 case Triton:
673 modeTotalBR[Triton] = decayModeTotal; break;
674 case RDM_ERROR:
675
676 default:
677 G4Exception("G4RadioactiveDecay::LoadDecayTable()", "HAD_RDM_000",
678 FatalException, "Selected decay mode does not exist");
679 } // switch
680
681 } else {
682 if (inputLine.length() < 84) {
683 tmpStream >> theDecayMode >> a >> daughterFloatFlag >> b >> c;
684 betaType = allowed;
685 } else {
686 tmpStream >> theDecayMode >> a >> daughterFloatFlag >> b >> c >> betaType;
687 }
688
689 // Allowed transitions are the default. Forbidden transitions are
690 // indicated in the last column.
691 a /= 1000.;
692 c /= 1000.;
693 b /= 100.;
694 daughterFloatLevel = G4Ions::FloatLevelBase(daughterFloatFlag.back());
695
696 switch (theDecayMode) {
697 case BetaMinus:
698 {
699 G4BetaMinusDecay* aBetaMinusChannel =
700 new G4BetaMinusDecay(theIon, b, c*MeV, a*MeV,
701 daughterFloatLevel, betaType);
702 //aBetaMinusChannel->DumpNuclearInfo();
703 theDecayTable->Insert(aBetaMinusChannel);
704 modeSumBR[BetaMinus] += b;
705 }
706 break;
707
708 case BetaPlus:
709 {
710 G4BetaPlusDecay* aBetaPlusChannel =
711 new G4BetaPlusDecay(theIon, b, c*MeV, a*MeV,
712 daughterFloatLevel, betaType);
713 //aBetaPlusChannel->DumpNuclearInfo();
714 theDecayTable->Insert(aBetaPlusChannel);
715 modeSumBR[BetaPlus] += b;
716 }
717 break;
718
719 case KshellEC: // K-shell electron capture
720 {
721 G4ECDecay* aKECChannel =
722 new G4ECDecay(theIon, b, c*MeV, a*MeV,
723 daughterFloatLevel, KshellEC);
724 //aKECChannel->DumpNuclearInfo();
725 aKECChannel->SetARM(applyARM);
726 theDecayTable->Insert(aKECChannel);
727 modeSumBR[KshellEC] += b;
728 }
729 break;
730
731 case LshellEC: // L-shell electron capture
732 {
733 G4ECDecay* aLECChannel =
734 new G4ECDecay(theIon, b, c*MeV, a*MeV,
735 daughterFloatLevel, LshellEC);
736// aLECChannel->DumpNuclearInfo();
737 aLECChannel->SetARM(applyARM);
738 theDecayTable->Insert(aLECChannel);
739 modeSumBR[LshellEC] += b;
740 }
741 break;
742
743 case MshellEC: // M-shell electron capture
744 {
745 G4ECDecay* aMECChannel =
746 new G4ECDecay(theIon, b, c*MeV, a*MeV,
747 daughterFloatLevel, MshellEC);
748// aMECChannel->DumpNuclearInfo();
749 aMECChannel->SetARM(applyARM);
750 theDecayTable->Insert(aMECChannel);
751 modeSumBR[MshellEC] += b;
752 }
753 break;
754
755 case NshellEC: // N-shell electron capture
756 {
757 G4ECDecay* aNECChannel =
758 new G4ECDecay(theIon, b, c*MeV, a*MeV,
759 daughterFloatLevel, NshellEC);
760// aNECChannel->DumpNuclearInfo();
761 aNECChannel->SetARM(applyARM);
762 theDecayTable->Insert(aNECChannel);
763 modeSumBR[NshellEC] += b;
764 }
765 break;
766
767 case Alpha:
768 {
769 G4AlphaDecay* anAlphaChannel =
770 new G4AlphaDecay(theIon, b, c*MeV, a*MeV,
771 daughterFloatLevel);
772// anAlphaChannel->DumpNuclearInfo();
773 theDecayTable->Insert(anAlphaChannel);
774 modeSumBR[Alpha] += b;
775 }
776 break;
777
778 case Proton:
779 {
780 G4ProtonDecay* aProtonChannel =
781 new G4ProtonDecay(theIon, b, c*MeV, a*MeV,
782 daughterFloatLevel);
783// aProtonChannel->DumpNuclearInfo();
784 theDecayTable->Insert(aProtonChannel);
785 modeSumBR[Proton] += b;
786 }
787 break;
788
789 case Neutron:
790 {
791 G4NeutronDecay* aNeutronChannel =
792 new G4NeutronDecay(theIon, b, c*MeV, a*MeV,
793 daughterFloatLevel);
794// aNeutronChannel->DumpNuclearInfo();
795 theDecayTable->Insert(aNeutronChannel);
796 modeSumBR[Neutron] += b;
797 }
798 break;
799
800 case SpFission:
801 {
802 G4SFDecay* aSpontFissChannel =
803 new G4SFDecay(theIon, b, c*MeV, a*MeV, daughterFloatLevel);
804 theDecayTable->Insert(aSpontFissChannel);
805 modeSumBR[SpFission] += b;
806 }
807 break;
808
809 case BDProton:
810 // Not yet implemented
811 // G4cout << " beta-delayed proton decay, a = " << a << ", b = " << b << ", c = " << c << G4endl;
812 break;
813
814 case BDNeutron:
815 // Not yet implemented
816 // G4cout << " beta-delayed neutron decay, a = " << a << ", b = " << b << ", c = " << c << G4endl;
817 break;
818
819 case Beta2Minus:
820 // Not yet implemented
821 // G4cout << " Double beta- decay, a = " << a << ", b = " << b << ", c = " << c << G4endl;
822 break;
823
824 case Beta2Plus:
825 // Not yet implemented
826 // G4cout << " Double beta+ decay, a = " << a << ", b = " << b << ", c = " << c << G4endl;
827 break;
828
829 case Proton2:
830 // Not yet implemented
831 // G4cout << " Double proton decay, a = " << a << ", b = " << b << ", c = " << c << G4endl;
832 break;
833
834 case Neutron2:
835 // Not yet implemented
836 // G4cout << " Double beta- decay, a = " << a << ", b = " << b << ", c = " << c << G4endl;
837 break;
838
839 case Triton:
840 {
841 G4TritonDecay* aTritonChannel =
842 new G4TritonDecay(theIon, b, c*MeV, a*MeV, daughterFloatLevel);
843 theDecayTable->Insert(aTritonChannel);
844 modeSumBR[Triton] += b;
845 }
846 break;
847
848 case RDM_ERROR:
849
850 default:
851 G4Exception("G4RadioactiveDecay::LoadDecayTable()", "HAD_RDM_000",
852 FatalException, "Selected decay mode does not exist");
853 } // switch
854 } // line < 72
855 } // if char == P
856 } // if char != #
857 } // While
858
859 // Go through the decay table and make sure that the branching ratios are
860 // correctly normalised.
861
862 G4VDecayChannel* theChannel = 0;
863 G4NuclearDecay* theNuclearDecayChannel = 0;
864 G4String mode = "";
865
866 G4double theBR = 0.0;
867 for (G4int i = 0; i < theDecayTable->entries(); ++i) {
868 theChannel = theDecayTable->GetDecayChannel(i);
869 theNuclearDecayChannel = static_cast<G4NuclearDecay*>(theChannel);
870 theDecayMode = theNuclearDecayChannel->GetDecayMode();
871
872 if (theDecayMode != IT) {
873 theBR = theChannel->GetBR();
874 theChannel->SetBR(theBR*modeTotalBR[theDecayMode]/modeSumBR[theDecayMode]);
875 }
876 }
877 } // decay file exists
878
879 DecaySchemeFile.close();
880
881 if (!found && levelEnergy > 0) {
882 // Case where IT cascade for excited isotopes has no entries in RDM database
883 // Decay mode is isomeric transition.
884 G4ITDecay* anITChannel = new G4ITDecay(theIon, 1.0, 0.0, 0.0);
885 theDecayTable->Insert(anITChannel);
886 }
887
888 if (GetVerboseLevel() > 1) {
889 theDecayTable->DumpInfo();
890 }
891
892 // store in master library
893 (*master_dkmap)[theIon->GetParticleName()] = theDecayTable;
894 lk.unlock();
895 return theDecayTable;
896}
G4BetaDecayType
@ allowed
#define noFloat
Definition G4Ions.hh:119
@ G4RadioactiveDecayModeSize
bool G4bool
Definition G4Types.hh:86
G4VDecayChannel * GetDecayChannel(G4int index) const
void Insert(G4VDecayChannel *aChannel)
void DumpInfo() const
void SetARM(G4bool onoff)
Definition G4ECDecay.hh:56
G4Ions::G4FloatLevelBase GetFloatLevelBase() const
Definition G4Ions.hh:159
static G4Ions::G4FloatLevelBase FloatLevelBase(char flbChar)
Definition G4Ions.cc:107
G4FloatLevelBase
Definition G4Ions.hh:80
G4RadioactiveDecayMode GetDecayMode() const
static const G4double levelTolerance
G4double GetBR() const
void SetBR(G4double value)

Referenced by GetDecayTable().

◆ operator=()

G4RadioactiveDecay & G4RadioactiveDecay::operator= ( const G4RadioactiveDecay & right)
delete

◆ PostStepDoIt()

G4VParticleChange * G4RadioactiveDecay::PostStepDoIt ( const G4Track & theTrack,
const G4Step & theStep )
overridevirtual

Reimplemented from G4VRestDiscreteProcess.

Definition at line 184 of file G4RadioactiveDecay.cc.

186{
187 return DecayIt(theTrack, theStep);
188}

◆ ProcessDescription()

void G4RadioactiveDecay::ProcessDescription ( std::ostream & outFile) const
overridevirtual

Reimplemented from G4VProcess.

Definition at line 191 of file G4RadioactiveDecay.cc.

192{
193 outFile << "The radioactive decay process (G4RadioactiveDecay) handles the\n"
194 << "alpha, beta+, beta-, electron capture and isomeric transition\n"
195 << "decays of nuclei (G4GenericIon) with masses A > 4.\n"
196 << "The required half-lives and decay schemes are retrieved from\n"
197 << "the RadioactiveDecay database which was derived from ENSDF.\n";
198}

◆ SelectAllVolumes()

void G4RadioactiveDecay::SelectAllVolumes ( )

Definition at line 334 of file G4RadioactiveDecay.cc.

335{
337 G4LogicalVolume* volume = nullptr;
338 ValidVolumes.clear();
339#ifdef G4VERBOSE
340 if (GetVerboseLevel()>1)
341 G4cout << " RDM Applies to all Volumes" << G4endl;
342#endif
343 for (std::size_t i = 0; i < theLogicalVolumes->size(); ++i){
344 volume = (*theLogicalVolumes)[i];
345 ValidVolumes.push_back(volume->GetName());
346#ifdef G4VERBOSE
347 if (GetVerboseLevel()>1)
348 G4cout << " RDM Applies to Volume " << volume->GetName() << G4endl;
349#endif
350 }
351 std::sort(ValidVolumes.begin(), ValidVolumes.end());
352 // sort needed in order to allow binary_search
353 isAllVolumesMode=true;
354}

Referenced by G4RadioactiveDecay(), and G4RadioactiveDecayMessenger::SetNewValue().

◆ SelectAVolume()

void G4RadioactiveDecay::SelectAVolume ( const G4String & aVolume)

Definition at line 272 of file G4RadioactiveDecay.cc.

273{
275 G4LogicalVolume* volume = nullptr;
276 volume = theLogicalVolumes->GetVolume(aVolume);
277 if (volume != nullptr)
278 {
279 ValidVolumes.push_back(aVolume);
280 std::sort(ValidVolumes.begin(), ValidVolumes.end());
281 // sort need for performing binary_search
282
283 if (GetVerboseLevel() > 0)
284 G4cout << " Radioactive decay applied to " << aVolume << G4endl;
285 }
286 else
287 {
289 ed << aVolume << " is not a valid logical volume name."
290 << " Decay not activated for it."
291 << G4endl;
292 G4Exception("G4RadioactiveDecay::SelectAVolume()", "HAD_RDM_300",
293 JustWarning, ed);
294 }
295}

Referenced by G4RadioactiveDecayMessenger::SetNewValue().

◆ SetARM()

void G4RadioactiveDecay::SetARM ( G4bool arm)
inline

Definition at line 117 of file G4RadioactiveDecay.hh.

117{applyARM = arm;}

Referenced by G4RadioactiveDecayMessenger::SetNewValue().

◆ SetDecayCollimation()

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

Definition at line 151 of file G4RadioactiveDecay.hh.

152 {
153 SetDecayDirection(theDir);
154 SetDecayHalfAngle(halfAngle);
155 }
void SetDecayHalfAngle(G4double halfAngle=0.*CLHEP::deg)
void SetDecayDirection(const G4ThreeVector &theDir)

◆ SetDecayDirection()

void G4RadioactiveDecay::SetDecayDirection ( const G4ThreeVector & theDir)
inline

Definition at line 135 of file G4RadioactiveDecay.hh.

135 {
136 forceDecayDirection = theDir.unit();
137 }
Hep3Vector unit() const

Referenced by SetDecayCollimation(), and G4RadioactiveDecayMessenger::SetNewValue().

◆ SetDecayHalfAngle()

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

Definition at line 143 of file G4RadioactiveDecay.hh.

143 {
144 forceDecayHalfAngle = std::min(std::max(0.*CLHEP::deg,halfAngle),180.*CLHEP::deg);
145 }

Referenced by SetDecayCollimation(), and G4RadioactiveDecayMessenger::SetNewValue().

◆ SetNucleusLimits()

void G4RadioactiveDecay::SetNucleusLimits ( G4NucleusLimits theNucleusLimits1)
inline

Definition at line 126 of file G4RadioactiveDecay.hh.

127 {theNucleusLimits = theNucleusLimits1 ;}

◆ SetThresholdForVeryLongDecayTime()

void G4RadioactiveDecay::SetThresholdForVeryLongDecayTime ( const G4double inputThreshold)
inline

Definition at line 158 of file G4RadioactiveDecay.hh.

158 {
159 fThresholdForVeryLongDecayTime = std::max( 0.0, inputThreshold );
160 }

Referenced by G4RadioactiveDecayMessenger::SetNewValue().

◆ StreamInfo()

void G4RadioactiveDecay::StreamInfo ( std::ostream & os,
const G4String & endline )

Definition at line 470 of file G4RadioactiveDecay.cc.

471{
476
477 G4long prec = os.precision(5);
478 os << "======================================================================"
479 << endline;
480 os << "====== Radioactive Decay Physics Parameters ======="
481 << endline;
482 os << "======================================================================"
483 << endline;
484 os << "min MeanLife (from G4NuclideTable) "
485 << G4BestUnit(minMeanLife, "Time") << endline;
486 os << "Max life time (from G4DeexPrecoParameters) "
487 << G4BestUnit(deex->GetMaxLifeTime(), "Time") << endline;
488 os << "Internal e- conversion flag "
489 << deex->GetInternalConversionFlag() << endline;
490 os << "Stored internal conversion coefficients "
491 << deex->StoreICLevelData() << endline;
492 os << "Enabled atomic relaxation mode "
493 << applyARM << endline;
494 os << "Enable correlated gamma emission "
495 << deex->CorrelatedGamma() << endline;
496 os << "Max 2J for sampling of angular correlations "
497 << deex->GetTwoJMAX() << endline;
498 os << "Atomic de-excitation enabled "
499 << emparam->Fluo() << endline;
500 os << "Auger electron emission enabled "
501 << emparam->Auger() << endline;
502 os << "Check EM cuts disabled for atomic de-excitation "
503 << emparam->DeexcitationIgnoreCut() << endline;
504 os << "Use Bearden atomic level energies "
505 << emparam->BeardenFluoDir() << endline;
506 os << "Use ANSTO fluorescence model "
507 << emparam->ANSTOFluoDir() << endline;
508 os << "Threshold for very long decay time at rest "
509 << G4BestUnit(fThresholdForVeryLongDecayTime, "Time") << endline;
510 os << "======================================================================"
511 << G4endl;
512 os.precision(prec);
513}
#define G4BestUnit(a, b)
long G4long
Definition G4Types.hh:87
G4bool GetInternalConversionFlag() const
static G4EmParameters * Instance()
G4bool Fluo() const
G4bool DeexcitationIgnoreCut() const
G4bool Auger() const
G4bool BeardenFluoDir()
G4DeexPrecoParameters * GetParameters()
static G4NuclearLevelData * GetInstance()
G4double GetMeanLifeThreshold()
static G4NuclideTable * GetInstance()

Referenced by BuildPhysicsTable().

Member Data Documentation

◆ decayIT

◆ fParticleChangeForRadDecay

G4ParticleChangeForRadDecay G4RadioactiveDecay::fParticleChangeForRadDecay
protected

◆ isAllVolumesMode

bool G4RadioactiveDecay::isAllVolumesMode {true}
protected

◆ levelTolerance

const G4double G4RadioactiveDecay::levelTolerance = 10.0*CLHEP::eV
staticprotected

◆ master_dkmap

DecayTableMap * G4RadioactiveDecay::master_dkmap = nullptr
staticprotected

◆ photonEvaporation

G4PhotonEvaporation* G4RadioactiveDecay::photonEvaporation
protected

◆ theRadioactiveDecayMessenger

G4RadioactiveDecayMessenger* G4RadioactiveDecay::theRadioactiveDecayMessenger
protected

Definition at line 190 of file G4RadioactiveDecay.hh.

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

◆ ValidVolumes

std::vector<G4String> G4RadioactiveDecay::ValidVolumes
protected

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