Geant4 11.1.1
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")
 
 ~G4RadioactiveDecay ()
 
virtual void ProcessDescription (std::ostream &outFile) const
 
G4bool IsApplicable (const G4ParticleDefinition &)
 
G4DecayTableGetDecayTable (const G4ParticleDefinition *)
 
void SelectAVolume (const G4String &aVolume)
 
void DeselectAVolume (const G4String &aVolume)
 
void SelectAllVolumes ()
 
void DeselectAllVolumes ()
 
void SetARM (G4bool arm)
 
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 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 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 ()
 
G4VRestDiscreteProcessoperator= (const G4VRestDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
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
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
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 G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
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
 
virtual void ProcessDescription (std::ostream &outfile) 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

void DecayAnalog (const G4Track &theTrack)
 
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)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
virtual G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Protected Attributes

G4ParticleChangeForRadDecay fParticleChangeForRadDecay
 
G4RadioactiveDecayMessengertheRadioactiveDecayMessenger
 
G4PhotonEvaporationphotonEvaporation
 
std::vector< G4StringValidVolumes
 
bool isAllVolumesMode
 
DecayTableMapdkmap
 
- 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*eV
 

Additional Inherited Members

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

Detailed Description

Definition at line 63 of file G4RadioactiveDecay.hh.

Constructor & Destructor Documentation

◆ G4RadioactiveDecay()

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

Definition at line 113 of file G4RadioactiveDecay.cc.

114 : G4VRestDiscreteProcess(processName, fDecay), isInitialised(false),
115 forceDecayDirection(0.,0.,0.), forceDecayHalfAngle(0.*deg), dirPath(""),
116 verboseLevel(1),
117 fThresholdForVeryLongDecayTime( 1.0e+27*CLHEP::nanosecond ) // Longer than twice Universe's age
118{
119#ifdef G4VERBOSE
120 if (GetVerboseLevel() > 1) {
121 G4cout << "G4RadioactiveDecay constructor: processName = " << processName
122 << G4endl;
123 }
124#endif
125
127
130
131 // Set up photon evaporation for use in G4ITDecay
135
136 // DHW G4DeexPrecoParameters* deex = G4NuclearLevelData::GetInstance()->GetParameters();
137 // DHW deex->SetCorrelatedGamma(true);
138
139 // Check data directory
140 const char* path_var = G4FindDataDir("G4RADIOACTIVEDATA");
141 if (!path_var) {
142 G4Exception("G4RadioactiveDecay()","HAD_RDM_200",FatalException,
143 "Environment variable G4RADIOACTIVEDATA is not set");
144 } else {
145 dirPath = path_var; // convert to string
146 std::ostringstream os;
147 os << dirPath << "/z1.a3"; // used as a dummy
148 std::ifstream testFile;
149 testFile.open(os.str() );
150 if (!testFile.is_open() )
151 G4Exception("G4RadioactiveDecay()","HAD_RDM_201",FatalException,
152 "Environment variable G4RADIOACTIVEDATA is set, but does not point to correct directory");
153 }
154
155 // Reset the list of user defined data files
156 theUserRadioactiveDataFiles.clear();
157
158 // Instantiate the map of decay tables
159#ifdef G4MULTITHREADED
160 G4AutoLock lk(&G4RadioactiveDecay::radioactiveDecayMutex);
161 NumberOfInstances()++;
162 if(!master_dkmap) master_dkmap = new DecayTableMap;
163#endif
164 dkmap = new DecayTableMap;
165
166 // Apply default values
167 applyARM = true;
168
169 // RDM applies to all logical volumes by default
170 isAllVolumesMode = true;
173}
const char * G4FindDataDir(const char *)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
@ fRadioactiveDecay
@ fDecay
std::map< G4String, G4DecayTable * > DecayTableMap
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void RegisterExtraProcess(G4VProcess *)
static G4HadronicProcessStore * Instance()
void RDMForced(G4bool) override
void SetICM(G4bool) override
G4RadioactiveDecayMessenger * theRadioactiveDecayMessenger
G4ParticleChangeForRadDecay fParticleChangeForRadDecay
G4int GetVerboseLevel() const
G4PhotonEvaporation * photonEvaporation
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:325

◆ ~G4RadioactiveDecay()

G4RadioactiveDecay::~G4RadioactiveDecay ( )

Definition at line 185 of file G4RadioactiveDecay.cc.

186{
188 delete photonEvaporation;
189 for (DecayTableMap::iterator i = dkmap->begin(); i != dkmap->end(); i++) {
190 delete i->second;
191 }
192 dkmap->clear();
193 delete dkmap;
194#ifdef G4MULTITHREADED
195 G4AutoLock lk(&G4RadioactiveDecay::radioactiveDecayMutex);
196 --NumberOfInstances();
197 if(NumberOfInstances()==0)
198 {
199 for (DecayTableMap::iterator i = master_dkmap->begin(); i != master_dkmap->end(); i++) {
200 delete i->second;
201 }
202 master_dkmap->clear();
203 delete master_dkmap;
204 }
205#endif
206}

Member Function Documentation

◆ AddUserDecayDataFile()

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

Definition at line 917 of file G4RadioactiveDecay.cc.

918{
919 if (Z < 1 || A < 2) G4cout << "Z and A not valid!" << G4endl;
920
921 std::ifstream DecaySchemeFile(filename);
922 if (DecaySchemeFile) {
923 G4int ID_ion = A*1000 + Z;
924 theUserRadioactiveDataFiles[ID_ion] = filename;
925 } else {
927 ed << filename << " does not exist! " << G4endl;
928 G4Exception("G4RadioactiveDecay::AddUserDecayDataFile()", "HAD_RDM_001",
929 FatalException, ed);
930 }
931}
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]

Referenced by G4RadioactiveDecayMessenger::SetNewValue().

◆ BuildPhysicsTable()

void G4RadioactiveDecay::BuildPhysicsTable ( const G4ParticleDefinition )
virtual

Reimplemented from G4VProcess.

Definition at line 449 of file G4RadioactiveDecay.cc.

450{
451 if (!isInitialised) {
452 isInitialised = true;
453#ifdef G4VERBOSE
455 G4Threading::IsMasterThread()) { StreamInfo(G4cout, "\n"); }
456#endif
457 }
460}
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
static G4HadronicParameters * Instance()
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)
G4bool IsMasterThread()
Definition: G4Threading.cc:124

◆ ChooseCollimationDirection()

G4ThreeVector G4RadioactiveDecay::ChooseCollimationDirection ( ) const
protected

Definition at line 1312 of file G4RadioactiveDecay.cc.

1312 {
1313 if (origin == forceDecayDirection) return origin; // Don't do collimation
1314 if (forceDecayHalfAngle == 180.*deg) return origin;
1315
1316 G4ThreeVector dir = forceDecayDirection;
1317
1318 // Return direction offset by random throw
1319 if (forceDecayHalfAngle > 0.) {
1320 // Generate uniform direction around central axis
1321 G4double phi = 2.*pi*G4UniformRand();
1322 G4double cosMin = std::cos(forceDecayHalfAngle);
1323 G4double cosTheta = (1.-cosMin)*G4UniformRand() + cosMin; // [cosMin,1.)
1324
1325 dir.setPhi(dir.phi()+phi);
1326 dir.setTheta(dir.theta()+std::acos(cosTheta));
1327 }
1328
1329#ifdef G4VERBOSE
1330 if (GetVerboseLevel()>1)
1331 G4cout << " ChooseCollimationDirection returns " << dir << G4endl;
1332#endif
1333
1334 return dir;
1335}
double G4double
Definition: G4Types.hh:83
#define G4UniformRand()
Definition: Randomize.hh:52
double phi() const
double theta() const
void setTheta(double)
void setPhi(double)
const G4double pi

Referenced by CollimateDecayProduct().

◆ CollimateDecay()

void G4RadioactiveDecay::CollimateDecay ( G4DecayProducts products)
protected

Definition at line 1267 of file G4RadioactiveDecay.cc.

1267 {
1268
1269 if (origin == forceDecayDirection) return; // No collimation requested
1270 if (180.*deg == forceDecayHalfAngle) return;
1271 if (0 == products || 0 == products->entries()) return;
1272
1273#ifdef G4VERBOSE
1274 if (GetVerboseLevel() > 1) G4cout << "Begin of CollimateDecay..." << G4endl;
1275#endif
1276
1277 // Particles suitable for directional biasing (for if-blocks below)
1281 static const G4ParticleDefinition* gamma = G4Gamma::Definition();
1285
1286 G4ThreeVector newDirection; // Re-use to avoid memory churn
1287 for (G4int i=0; i<products->entries(); i++) {
1288 G4DynamicParticle* daughter = (*products)[i];
1289 const G4ParticleDefinition* daughterType =
1290 daughter->GetParticleDefinition();
1291 if (daughterType == electron || daughterType == positron ||
1292 daughterType == neutron || daughterType == gamma ||
1293 daughterType == alpha || daughterType == triton || daughterType == proton) CollimateDecayProduct(daughter);
1294 }
1295}
static G4Alpha * Definition()
Definition: G4Alpha.cc:48
G4int entries() const
const G4ParticleDefinition * GetParticleDefinition() const
static G4Electron * Definition()
Definition: G4Electron.cc:48
static G4Gamma * Definition()
Definition: G4Gamma.cc:48
static G4Neutron * Definition()
Definition: G4Neutron.cc:53
static G4Positron * Definition()
Definition: G4Positron.cc:48
static G4Proton * Definition()
Definition: G4Proton.cc:48
void CollimateDecayProduct(G4DynamicParticle *product)
static G4Triton * Definition()
Definition: G4Triton.cc:49

Referenced by DoDecay().

◆ CollimateDecayProduct()

void G4RadioactiveDecay::CollimateDecayProduct ( G4DynamicParticle product)
protected

Definition at line 1297 of file G4RadioactiveDecay.cc.

1297 {
1298#ifdef G4VERBOSE
1299 if (GetVerboseLevel() > 1) {
1300 G4cout << "CollimateDecayProduct for daughter "
1301 << daughter->GetParticleDefinition()->GetParticleName() << G4endl;
1302 }
1303#endif
1304
1306 if (origin != collimate) daughter->SetMomentumDirection(collimate);
1307}
G4ThreeVector ChooseCollimationDirection() const

Referenced by CollimateDecay().

◆ DecayAnalog()

void G4RadioactiveDecay::DecayAnalog ( const G4Track theTrack)
protected

Definition at line 1107 of file G4RadioactiveDecay.cc.

1108{
1109 const G4DynamicParticle* theParticle = theTrack.GetDynamicParticle();
1110 const G4ParticleDefinition* theParticleDef = theParticle->GetDefinition();
1111 G4DecayProducts* products = DoDecay(*theParticleDef);
1112
1113 // Check if the product is the same as input and kill the track if
1114 // necessary to prevent infinite loop (11/05/10, F.Lei)
1115 if (products->entries() == 1) {
1120 delete products;
1121 return;
1122 }
1123
1124 G4double energyDeposit = 0.0;
1125 G4double finalGlobalTime = theTrack.GetGlobalTime();
1126 G4double finalLocalTime = theTrack.GetLocalTime();
1127
1128 // Get parent particle information and boost the decay products to the
1129 // laboratory frame
1130
1131 // ParentEnergy used for the boost should be the total energy of the nucleus
1132 // of the parent ion without the energy of the shell electrons
1133 // (correction for bug 1359 by L. Desorgher)
1134 G4double ParentEnergy = theParticle->GetKineticEnergy()
1135 + theParticle->GetParticleDefinition()->GetPDGMass();
1136 G4ThreeVector ParentDirection(theParticle->GetMomentumDirection());
1137
1138 if (theTrack.GetTrackStatus() == fStopButAlive) {
1139 // this condition seems to be always True, further investigation is needed (L.Desorgher)
1140
1141 // The particle is decayed at rest
1142 // Since the time is for the particle at rest, need to add additional time
1143 // lapsed between particle coming to rest and the actual decay. This time
1144 // is sampled with the mean-life of the particle. Need to protect the case
1145 // PDGTime < 0. (F.Lei 11/05/10)
1146 G4double temptime = -std::log(G4UniformRand() ) *
1147 theParticleDef->GetPDGLifeTime();
1148 if (temptime < 0.) temptime = 0.;
1149 finalGlobalTime += temptime;
1150 finalLocalTime += temptime;
1151 energyDeposit += theParticle->GetKineticEnergy();
1152
1153 // Kill the parent particle, and ignore its decay, if it decays later than the
1154 // threshold fThresholdForVeryLongDecayTime (whose default value corresponds
1155 // to more than twice the age of the universe).
1156 // This kind of cut has been introduced (in April 2021) in order to avoid to
1157 // account energy depositions happening after many billions of years in
1158 // ordinary materials used in calorimetry, in particular Tungsten and Lead
1159 // (via their natural unstable, but very long lived, isotopes, such as
1160 // W183, W180 and Pb204).
1161 // Note that the cut is not on the average, mean lifetime, but on the actual
1162 // sampled global decay time.
1163 if ( finalGlobalTime > fThresholdForVeryLongDecayTime ) {
1168 delete products;
1169 return;
1170 }
1171 }
1172 products->Boost(ParentEnergy, ParentDirection);
1173
1174 // Add products in theParticleChangeForRadDecay.
1175 G4int numberOfSecondaries = products->entries();
1177
1178 if (GetVerboseLevel() > 1) {
1179 G4cout << "G4RadioactiveDecay::DecayAnalog: Decay vertex :";
1180 G4cout << " Time: " << finalGlobalTime/ns << "[ns]";
1181 G4cout << " X:" << (theTrack.GetPosition()).x() /cm << "[cm]";
1182 G4cout << " Y:" << (theTrack.GetPosition()).y() /cm << "[cm]";
1183 G4cout << " Z:" << (theTrack.GetPosition()).z() /cm << "[cm]";
1184 G4cout << G4endl;
1185 G4cout << "G4Decay::DecayIt : decay products in Lab. Frame" << G4endl;
1186 products->DumpInfo();
1187 products->IsChecked();
1188 }
1189
1190 const G4int modelID_forIT = G4PhysicsModelCatalog::GetModelID( "model_RDM_IT" );
1191 G4int modelID = modelID_forIT + 10*theRadDecayMode;
1192 const G4int modelID_forAtomicRelaxation =
1193 G4PhysicsModelCatalog::GetModelID( "model_RDM_AtomicRelaxation" );
1194 for ( G4int index = 0; index < numberOfSecondaries; ++index ) {
1195 G4Track* secondary = new G4Track( products->PopProducts(), finalGlobalTime,
1196 theTrack.GetPosition() );
1197 secondary->SetWeight( theTrack.GetWeight() );
1198 secondary->SetCreatorModelID( modelID );
1199 // Change for atomics relaxation
1200 if ( theRadDecayMode == IT && index > 0 ) {
1201 if ( index == numberOfSecondaries-1 ) {
1202 secondary->SetCreatorModelID( modelID_forIT );
1203 } else {
1204 secondary->SetCreatorModelID( modelID_forAtomicRelaxation) ;
1205 }
1206 } else if ( theRadDecayMode >= KshellEC && theRadDecayMode <= NshellEC &&
1207 index < numberOfSecondaries-1 ) {
1208 secondary->SetCreatorModelID( modelID_forAtomicRelaxation );
1209 }
1210 secondary->SetGoodForTrackingFlag();
1211 secondary->SetTouchableHandle( theTrack.GetTouchableHandle() );
1213 }
1214
1215 delete products;
1216
1217 // Kill the parent particle
1221
1222 // Reset NumberOfInteractionLengthLeft.
1224}
@ 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 &theParticleDef)
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()
Definition: G4VProcess.hh:428
#define ns(x)
Definition: xmltok.c:1649

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

◆ DecayIt()

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

Definition at line 940 of file G4RadioactiveDecay.cc.

941{
942 // Initialize G4ParticleChange object, get particle details and decay table
945 const G4DynamicParticle* theParticle = theTrack.GetDynamicParticle();
946 const G4ParticleDefinition* theParticleDef = theParticle->GetDefinition();
947
948 // First check whether RDM applies to the current logical volume
949 if (!isAllVolumesMode) {
950 if (!std::binary_search(ValidVolumes.begin(), ValidVolumes.end(),
951 theTrack.GetVolume()->GetLogicalVolume()->GetName())) {
952#ifdef G4VERBOSE
953 if (GetVerboseLevel()>1) {
954 G4cout <<"G4RadioactiveDecay::DecayIt : "
955 << theTrack.GetVolume()->GetLogicalVolume()->GetName()
956 << " is not selected for the RDM"<< G4endl;
957 G4cout << " There are " << ValidVolumes.size() << " volumes" << G4endl;
958 G4cout << " The Valid volumes are " << G4endl;
959 for (std::size_t i = 0; i< ValidVolumes.size(); ++i)
960 G4cout << ValidVolumes[i] << G4endl;
961 }
962#endif
964
965 // Kill the parent particle.
970 }
971 }
972
973 // Now check if particle is valid for RDM
974 if (!(IsApplicable(*theParticleDef) ) ) {
975 // Particle is not an ion or is outside the nucleuslimits for decay
976#ifdef G4VERBOSE
977 if (GetVerboseLevel() > 1) {
978 G4cout << "G4RadioactiveDecay::DecayIt : "
979 << theParticleDef->GetParticleName()
980 << " is not an ion or is outside (Z,A) limits set for the decay. "
981 << " Set particle change accordingly. "
982 << G4endl;
983 }
984#endif
986
987 // Kill the parent particle
992 }
993
994 G4DecayTable* theDecayTable = GetDecayTable(theParticleDef);
995
996 if (theDecayTable == 0 || theDecayTable->entries() == 0) {
997 // No data in the decay table. Set particle change parameters
998 // to indicate this.
999#ifdef G4VERBOSE
1000 if (GetVerboseLevel() > 1) {
1001 G4cout << "G4RadioactiveDecay::DecayIt : "
1002 << "decay table not defined for "
1003 << theParticleDef->GetParticleName()
1004 << ". Set particle change accordingly. "
1005 << G4endl;
1006 }
1007#endif
1009
1010 // Kill the parent particle.
1015
1016 } else {
1017 // Data found. Try to decay nucleus
1018
1019/*
1020 G4double energyDeposit = 0.0;
1021 G4double finalGlobalTime = theTrack.GetGlobalTime();
1022 G4double finalLocalTime = theTrack.GetLocalTime();
1023 G4int index;
1024 G4ThreeVector currentPosition;
1025 currentPosition = theTrack.GetPosition();
1026
1027 G4DecayProducts* products = DoDecay(*theParticleDef);
1028
1029 // If the product is the same as the input kill the track if
1030 // necessary to prevent infinite loop (11/05/10, F.Lei)
1031 if (products->entries() == 1) {
1032 fParticleChangeForRadDecay.SetNumberOfSecondaries(0);
1033 fParticleChangeForRadDecay.ProposeTrackStatus(fStopAndKill);
1034 fParticleChangeForRadDecay.ProposeLocalEnergyDeposit(0.0);
1035 ClearNumberOfInteractionLengthLeft();
1036 return &fParticleChangeForRadDecay;
1037 }
1038
1039 // Get parent particle information and boost the decay products to the
1040 // laboratory frame based on this information.
1041
1042 // The Parent Energy used for the boost should be the total energy of
1043 // the nucleus of the parent ion without the energy of the shell electrons
1044 // (correction for bug 1359 by L. Desorgher)
1045 G4double ParentEnergy = theParticle->GetKineticEnergy()
1046 + theParticle->GetParticleDefinition()->GetPDGMass();
1047 G4ThreeVector ParentDirection(theParticle->GetMomentumDirection());
1048
1049 if (theTrack.GetTrackStatus() == fStopButAlive) {
1050 // This condition seems to be always True, further investigation is needed
1051 // (L.Desorgher)
1052 // The particle is decayed at rest.
1053 // since the time is still for rest particle in G4 we need to add the
1054 // additional time lapsed between the particle come to rest and the
1055 // actual decay. This time is simply sampled with the mean-life of
1056 // the particle. But we need to protect the case PDGTime < 0.
1057 // (F.Lei 11/05/10)
1058 G4double temptime = -std::log( G4UniformRand())
1059 *theParticleDef->GetPDGLifeTime();
1060 if (temptime < 0.) temptime = 0.;
1061 finalGlobalTime += temptime;
1062 finalLocalTime += temptime;
1063 energyDeposit += theParticle->GetKineticEnergy();
1064 }
1065 products->Boost(ParentEnergy, ParentDirection);
1066
1067 // Add products in theParticleChangeForRadDecay.
1068 G4int numberOfSecondaries = products->entries();
1069 fParticleChangeForRadDecay.SetNumberOfSecondaries(numberOfSecondaries);
1070
1071#ifdef G4VERBOSE
1072 if (GetVerboseLevel()>1) {
1073 G4cout <<"G4RadioactiveDecay::DecayIt : Decay vertex :";
1074 G4cout <<" Time: " <<finalGlobalTime/ns <<"[ns]";
1075 G4cout <<" X:" <<(theTrack.GetPosition()).x() /cm <<"[cm]";
1076 G4cout <<" Y:" <<(theTrack.GetPosition()).y() /cm <<"[cm]";
1077 G4cout <<" Z:" <<(theTrack.GetPosition()).z() /cm <<"[cm]";
1078 G4cout << G4endl;
1079 G4cout <<"G4Decay::DecayIt : decay products in Lab. Frame" <<G4endl;
1080 products->DumpInfo();
1081 products->IsChecked();
1082 }
1083#endif
1084 for (index=0; index < numberOfSecondaries; index++) {
1085 G4Track* secondary = new G4Track(products->PopProducts(),
1086 finalGlobalTime, currentPosition);
1087 secondary->SetGoodForTrackingFlag();
1088 secondary->SetTouchableHandle(theTrack.GetTouchableHandle());
1089 fParticleChangeForRadDecay.AddSecondary(secondary);
1090 }
1091 delete products;
1092
1093 // Kill the parent particle
1094 fParticleChangeForRadDecay.ProposeTrackStatus(fStopAndKill) ;
1095 fParticleChangeForRadDecay.ProposeLocalEnergyDeposit(energyDeposit);
1096 fParticleChangeForRadDecay.ProposeLocalTime(finalLocalTime);
1097 // Reset NumberOfInteractionLengthLeft.
1098 ClearNumberOfInteractionLengthLeft();
1099*/
1100 // Decay without variance reduction
1101 DecayAnalog(theTrack);
1103 }
1104}
G4int entries() const
const G4String & GetName() const
void Initialize(const G4Track &) final
const G4String & GetParticleName() const
void DecayAnalog(const G4Track &theTrack)
std::vector< G4String > ValidVolumes
G4DecayTable * GetDecayTable(const G4ParticleDefinition *)
G4bool IsApplicable(const G4ParticleDefinition &)
G4VPhysicalVolume * GetVolume() const
void ProposeWeight(G4double finalWeight)
G4LogicalVolume * GetLogicalVolume() const

◆ DeselectAllVolumes()

void G4RadioactiveDecay::DeselectAllVolumes ( )

Definition at line 332 of file G4RadioactiveDecay.cc.

333{
334 ValidVolumes.clear();
335 isAllVolumesMode=false;
336#ifdef G4VERBOSE
337 if (GetVerboseLevel() > 1) G4cout << "RDM removed from all volumes" << G4endl;
338#endif
339}

Referenced by G4RadioactiveDecayMessenger::SetNewValue().

◆ DeselectAVolume()

void G4RadioactiveDecay::DeselectAVolume ( const G4String aVolume)

Definition at line 273 of file G4RadioactiveDecay.cc.

274{
276 G4LogicalVolume* volume = nullptr;
277 volume = theLogicalVolumes->GetVolume(aVolume);
278 if (volume != nullptr)
279 {
280 auto location= std::find(ValidVolumes.cbegin(),ValidVolumes.cend(),aVolume);
281 if (location != ValidVolumes.cend() )
282 {
283 ValidVolumes.erase(location);
284 std::sort(ValidVolumes.begin(), ValidVolumes.end());
285 isAllVolumesMode = false;
286 if (GetVerboseLevel() > 0)
287 G4cout << " G4RadioactiveDecay::DeselectAVolume: " << aVolume
288 << " is removed from list " << G4endl;
289 }
290 else
291 {
293 ed << aVolume << " is not in the list. No action taken." << G4endl;
294 G4Exception("G4RadioactiveDecay::DeselectAVolume()", "HAD_RDM_300",
295 JustWarning, ed);
296 }
297 }
298 else
299 {
301 ed << aVolume << " is not a valid logical volume name. No action taken."
302 << G4endl;
303 G4Exception("G4RadioactiveDecay::DeselectAVolume()", "HAD_RDM_300",
304 JustWarning, ed);
305 }
306}
@ 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)
protected

Definition at line 1228 of file G4RadioactiveDecay.cc.

1229{
1230 G4DecayProducts* products = 0;
1231 G4DecayTable* theDecayTable = GetDecayTable(&theParticleDef);
1232
1233 // Choose a decay channel.
1234 // G4DecayTable::SelectADecayChannel checks to see if sum of daughter masses
1235 // exceeds parent mass. Pass it the parent mass + maximum Q value to account
1236 // for difference in mass defect.
1237 G4double parentPlusQ = theParticleDef.GetPDGMass() + 30.*MeV;
1238 G4VDecayChannel* theDecayChannel = theDecayTable->SelectADecayChannel(parentPlusQ);
1239
1240 if (theDecayChannel == 0) {
1241 // Decay channel not found.
1243 ed << " Cannot determine decay channel for " << theParticleDef.GetParticleName() << G4endl;
1244 G4Exception("G4RadioactiveDecay::DoDecay", "HAD_RDM_013",
1245 FatalException, ed);
1246 } else {
1247 // A decay channel has been identified, so execute the DecayIt.
1248#ifdef G4VERBOSE
1249 if (GetVerboseLevel() > 1) {
1250 G4cout << "G4RadioactiveDecay::DoIt : selected decay channel addr: "
1251 << theDecayChannel << G4endl;
1252 }
1253#endif
1254 theRadDecayMode = (static_cast<G4NuclearDecay*>(theDecayChannel))->GetDecayMode();
1255 products = theDecayChannel->DecayIt(theParticleDef.GetPDGMass() );
1256
1257 // Apply directional bias if requested by user
1258 CollimateDecay(products);
1259 }
1260
1261 return products;
1262}
G4VDecayChannel * SelectADecayChannel(G4double parentMass=-1.)
Definition: G4DecayTable.cc:82
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 131 of file G4RadioactiveDecay.hh.

131 {
132 return forceDecayDirection;
133 }

◆ GetDecayHalfAngle()

G4double G4RadioactiveDecay::GetDecayHalfAngle ( ) const
inline

Definition at line 139 of file G4RadioactiveDecay.hh.

139{return forceDecayHalfAngle;}

◆ GetDecayTable()

G4DecayTable * G4RadioactiveDecay::GetDecayTable ( const G4ParticleDefinition aNucleus)

Definition at line 231 of file G4RadioactiveDecay.cc.

232{
233 G4String key = aNucleus->GetParticleName();
234 DecayTableMap::iterator table_ptr = dkmap->find(key);
235
236 G4DecayTable* theDecayTable = 0;
237 if (table_ptr == dkmap->end() ) { // If table not there,
238 theDecayTable = LoadDecayTable(*aNucleus); // load from file and
239 if(theDecayTable) (*dkmap)[key] = theDecayTable; // store in library
240 } else {
241 theDecayTable = table_ptr->second;
242 }
243 return theDecayTable;
244}
G4DecayTable * LoadDecayTable(const G4ParticleDefinition &theParentNucleus)

Referenced by DecayIt(), and DoDecay().

◆ GetMeanFreePath()

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

Implements G4VRestDiscreteProcess.

Definition at line 384 of file G4RadioactiveDecay.cc.

386{
387 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
388 const G4ParticleDefinition* aParticleDef = aParticle->GetDefinition();
389 G4double tau = aParticleDef->GetPDGLifeTime();
390 G4double aMass = aParticle->GetMass();
391
392#ifdef G4VERBOSE
393 if (GetVerboseLevel() > 2) {
394 G4cout << "G4RadioactiveDecay::GetMeanFreePath() " << G4endl;
395 G4cout << " KineticEnergy: " << aParticle->GetKineticEnergy()/GeV
396 << " GeV, Mass: " << aMass/GeV << " GeV, tau: " << tau << " ns "
397 << G4endl;
398 }
399#endif
400 G4double pathlength = DBL_MAX;
401 if (tau != -1) {
402 // Ion can decay
403
404 if (tau < -1000.0) {
405 pathlength = DBL_MIN; // nuclide had very short lifetime or wasn't in table
406
407 } else if (tau < 0.0) {
408 G4cout << aParticleDef->GetParticleName() << " has lifetime " << tau << G4endl;
410 ed << "Ion has negative lifetime " << tau
411 << " but is not stable. Setting mean free path to DBL_MAX" << G4endl;
412 G4Exception("G4RadioactiveDecay::GetMeanFreePath()", "HAD_RDM_011",
413 JustWarning, ed);
414 pathlength = DBL_MAX;
415
416 } else {
417 // Calculate mean free path
418 G4double betaGamma = aParticle->GetTotalMomentum()/aMass;
419 pathlength = c_light*tau*betaGamma;
420
421 if (pathlength < DBL_MIN) {
422 pathlength = DBL_MIN;
423#ifdef G4VERBOSE
424 if (GetVerboseLevel() > 2) {
425 G4cout << "G4Decay::GetMeanFreePath: "
426 << aParticleDef->GetParticleName()
427 << " stops, kinetic energy = "
428 << aParticle->GetKineticEnergy()/keV <<" keV " << G4endl;
429 }
430#endif
431 }
432 }
433 }
434
435#ifdef G4VERBOSE
436 if (GetVerboseLevel() > 2) {
437 G4cout << "mean free path: "<< pathlength/m << " m" << G4endl;
438 }
439#endif
440 return pathlength;
441}
G4double GetMass() const
G4double GetTotalMomentum() const
#define DBL_MIN
Definition: templates.hh:54
#define DBL_MAX
Definition: templates.hh:62

◆ GetMeanLifeTime()

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

Implements G4VRestDiscreteProcess.

Definition at line 348 of file G4RadioactiveDecay.cc.

350{
351 G4double meanlife = 0.;
352 const G4DynamicParticle* theParticle = theTrack.GetDynamicParticle();
353 const G4ParticleDefinition* theParticleDef = theParticle->GetDefinition();
354 G4double theLife = theParticleDef->GetPDGLifeTime();
355#ifdef G4VERBOSE
356 if (GetVerboseLevel() > 2) {
357 G4cout << "G4RadioactiveDecay::GetMeanLifeTime() " << G4endl;
358 G4cout << "KineticEnergy: " << theParticle->GetKineticEnergy()/GeV
359 << " GeV, Mass: " << theParticle->GetMass()/GeV
360 << " GeV, Life time: " << theLife/ns << " ns " << G4endl;
361 }
362#endif
363 if (theParticleDef->GetPDGStable()) {meanlife = DBL_MAX;}
364 else if (theLife < 0.0) {meanlife = DBL_MAX;}
365 else {meanlife = theLife;}
366 // Set meanlife to zero for excited istopes which are not in the
367 // RDM database
368 if (((const G4Ions*)(theParticleDef))->GetExcitationEnergy() > 0. &&
369 meanlife == DBL_MAX) {meanlife = 0.;}
370#ifdef G4VERBOSE
371 if (GetVerboseLevel() > 2)
372 G4cout << " mean life time: " << meanlife/s << " s " << G4endl;
373#endif
374
375 return meanlife;
376}
Definition: G4Ions.hh:52
G4bool GetPDGStable() const

Referenced by G4Radioactivation::GetMeanLifeTime().

◆ GetNucleusLimits()

G4NucleusLimits G4RadioactiveDecay::GetNucleusLimits ( ) const
inline

Definition at line 125 of file G4RadioactiveDecay.hh.

125{return theNucleusLimits;}

◆ GetThresholdForVeryLongDecayTime()

G4double G4RadioactiveDecay::GetThresholdForVeryLongDecayTime ( ) const
inline

Definition at line 153 of file G4RadioactiveDecay.hh.

153{return fThresholdForVeryLongDecayTime;}

◆ GetVerboseLevel()

◆ IsApplicable()

G4bool G4RadioactiveDecay::IsApplicable ( const G4ParticleDefinition aParticle)
virtual

Reimplemented from G4VProcess.

Definition at line 209 of file G4RadioactiveDecay.cc.

210{
211 // All particles other than G4Ions, are rejected by default
212 if (((const G4Ions*)(&aParticle))->GetExcitationEnergy() > 0.) {return true;}
213 if (aParticle.GetParticleName() == "GenericIon") {
214 return true;
215 } else if (!(aParticle.GetParticleType() == "nucleus")
216 || aParticle.GetPDGLifeTime() < 0. ) {
217 return false;
218 }
219
220 // Determine whether the nuclide falls into the correct A and Z range
221 G4int A = ((const G4Ions*) (&aParticle))->GetAtomicMass();
222 G4int Z = ((const G4Ions*) (&aParticle))->GetAtomicNumber();
223
224 if (A > theNucleusLimits.GetAMax() || A < theNucleusLimits.GetAMin())
225 {return false;}
226 else if (Z > theNucleusLimits.GetZMax() || Z < theNucleusLimits.GetZMin())
227 {return false;}
228 return true;
229}
G4int GetZMax() const
G4int GetZMin() const
G4int GetAMin() const
G4int GetAMax() const
const G4String & GetParticleType() const

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

◆ LoadDecayTable()

G4DecayTable * G4RadioactiveDecay::LoadDecayTable ( const G4ParticleDefinition theParentNucleus)

Definition at line 523 of file G4RadioactiveDecay.cc.

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

Referenced by GetDecayTable(), and G4Radioactivation::GetDecayTable1().

◆ ProcessDescription()

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

Reimplemented from G4VProcess.

Reimplemented in G4Radioactivation.

Definition at line 175 of file G4RadioactiveDecay.cc.

176{
177 outFile << "The radioactive decay process (G4RadioactiveDecay) handles the\n"
178 << "alpha, beta+, beta-, electron capture and isomeric transition\n"
179 << "decays of nuclei (G4GenericIon) with masses A > 4.\n"
180 << "The required half-lives and decay schemes are retrieved from\n"
181 << "the RadioactiveDecay database which was derived from ENSDF.\n";
182}

◆ SelectAllVolumes()

void G4RadioactiveDecay::SelectAllVolumes ( )

Definition at line 309 of file G4RadioactiveDecay.cc.

310{
312 G4LogicalVolume* volume = nullptr;
313 ValidVolumes.clear();
314#ifdef G4VERBOSE
315 if (GetVerboseLevel()>1)
316 G4cout << " RDM Applies to all Volumes" << G4endl;
317#endif
318 for (std::size_t i = 0; i < theLogicalVolumes->size(); ++i){
319 volume = (*theLogicalVolumes)[i];
320 ValidVolumes.push_back(volume->GetName());
321#ifdef G4VERBOSE
322 if (GetVerboseLevel()>1)
323 G4cout << " RDM Applies to Volume " << volume->GetName() << G4endl;
324#endif
325 }
326 std::sort(ValidVolumes.begin(), ValidVolumes.end());
327 // sort needed in order to allow binary_search
328 isAllVolumesMode=true;
329}

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

◆ SelectAVolume()

void G4RadioactiveDecay::SelectAVolume ( const G4String aVolume)

Definition at line 247 of file G4RadioactiveDecay.cc.

248{
250 G4LogicalVolume* volume = nullptr;
251 volume = theLogicalVolumes->GetVolume(aVolume);
252 if (volume != nullptr)
253 {
254 ValidVolumes.push_back(aVolume);
255 std::sort(ValidVolumes.begin(), ValidVolumes.end());
256 // sort need for performing binary_search
257
258 if (GetVerboseLevel() > 0)
259 G4cout << " Radioactive decay applied to " << aVolume << G4endl;
260 }
261 else
262 {
264 ed << aVolume << " is not a valid logical volume name."
265 << " Decay not activated for it."
266 << G4endl;
267 G4Exception("G4RadioactiveDecay::SelectAVolume()", "HAD_RDM_300",
268 JustWarning, ed);
269 }
270}

Referenced by G4RadioactiveDecayMessenger::SetNewValue().

◆ SetARM()

void G4RadioactiveDecay::SetARM ( G4bool  arm)
inline

Definition at line 103 of file G4RadioactiveDecay.hh.

103{applyARM = arm;}

Referenced by G4RadioactiveDecayMessenger::SetNewValue().

◆ SetDecayCollimation()

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

Definition at line 143 of file G4RadioactiveDecay.hh.

144 {
145 SetDecayDirection(theDir);
146 SetDecayHalfAngle(halfAngle);
147 }
void SetDecayHalfAngle(G4double halfAngle=0.*CLHEP::deg)
void SetDecayDirection(const G4ThreeVector &theDir)

◆ SetDecayDirection()

void G4RadioactiveDecay::SetDecayDirection ( const G4ThreeVector theDir)
inline

Definition at line 127 of file G4RadioactiveDecay.hh.

127 {
128 forceDecayDirection = theDir.unit();
129 }
Hep3Vector unit() const

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

◆ SetDecayHalfAngle()

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

Definition at line 135 of file G4RadioactiveDecay.hh.

135 {
136 forceDecayHalfAngle = std::min(std::max(0.*CLHEP::deg,halfAngle),180.*CLHEP::deg);
137 }

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

◆ SetNucleusLimits()

void G4RadioactiveDecay::SetNucleusLimits ( G4NucleusLimits  theNucleusLimits1)
inline

Definition at line 118 of file G4RadioactiveDecay.hh.

119 {theNucleusLimits = theNucleusLimits1 ;}

◆ SetThresholdForVeryLongDecayTime()

void G4RadioactiveDecay::SetThresholdForVeryLongDecayTime ( const G4double  inputThreshold)
inline

Definition at line 150 of file G4RadioactiveDecay.hh.

150 {
151 fThresholdForVeryLongDecayTime = std::max( 0.0, inputThreshold );
152 }

Referenced by G4RadioactiveDecayMessenger::SetNewValue().

◆ SetVerboseLevel()

void G4RadioactiveDecay::SetVerboseLevel ( G4int  value)
inline

Definition at line 112 of file G4RadioactiveDecay.hh.

112{verboseLevel = value;}

Referenced by G4RadioactiveDecayMessenger::SetNewValue().

Member Data Documentation

◆ dkmap

DecayTableMap* G4RadioactiveDecay::dkmap
protected

◆ fParticleChangeForRadDecay

G4ParticleChangeForRadDecay G4RadioactiveDecay::fParticleChangeForRadDecay
protected

◆ isAllVolumesMode

bool G4RadioactiveDecay::isAllVolumesMode
protected

◆ levelTolerance

const G4double G4RadioactiveDecay::levelTolerance = 10.0*eV
staticprotected

◆ photonEvaporation

◆ theRadioactiveDecayMessenger

G4RadioactiveDecayMessenger* G4RadioactiveDecay::theRadioactiveDecayMessenger
protected

Definition at line 180 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: