101#ifdef G4MULTITHREADED
106G4int& G4RadioactiveDecay::NumberOfInstances()
108 static G4int numberOfInstances = 0;
109 return numberOfInstances;
115 forceDecayDirection(0.,0.,0.), forceDecayHalfAngle(0.*deg), dirPath(
""),
117 fThresholdForVeryLongDecayTime( 1.0e+27*
CLHEP::nanosecond )
121 G4cout <<
"G4RadioactiveDecay constructor: processName = " << processName
143 "Environment variable G4RADIOACTIVEDATA is not set");
146 std::ostringstream os;
147 os << dirPath <<
"/z1.a3";
148 std::ifstream testFile;
149 testFile.open(os.str() );
150 if (!testFile.is_open() )
152 "Environment variable G4RADIOACTIVEDATA is set, but does not point to correct directory");
156 theUserRadioactiveDataFiles.clear();
159#ifdef G4MULTITHREADED
160 G4AutoLock lk(&G4RadioactiveDecay::radioactiveDecayMutex);
161 NumberOfInstances()++;
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";
189 for (DecayTableMap::iterator i =
dkmap->begin(); i !=
dkmap->end(); i++) {
194#ifdef G4MULTITHREADED
195 G4AutoLock lk(&G4RadioactiveDecay::radioactiveDecayMutex);
196 --NumberOfInstances();
197 if(NumberOfInstances()==0)
199 for (DecayTableMap::iterator i = master_dkmap->begin(); i != master_dkmap->end(); i++) {
202 master_dkmap->clear();
212 if (((
const G4Ions*)(&aParticle))->GetExcitationEnergy() > 0.) {
return true;}
221 G4int A = ((
const G4Ions*) (&aParticle))->GetAtomicMass();
222 G4int Z = ((
const G4Ions*) (&aParticle))->GetAtomicNumber();
234 DecayTableMap::iterator table_ptr =
dkmap->find(key);
237 if (table_ptr ==
dkmap->end() ) {
239 if(theDecayTable) (*dkmap)[key] = theDecayTable;
241 theDecayTable = table_ptr->second;
243 return theDecayTable;
251 volume = theLogicalVolumes->
GetVolume(aVolume);
252 if (volume !=
nullptr)
259 G4cout <<
" Radioactive decay applied to " << aVolume <<
G4endl;
264 ed << aVolume <<
" is not a valid logical volume name."
265 <<
" Decay not activated for it."
267 G4Exception(
"G4RadioactiveDecay::SelectAVolume()",
"HAD_RDM_300",
277 volume = theLogicalVolumes->
GetVolume(aVolume);
278 if (volume !=
nullptr)
287 G4cout <<
" G4RadioactiveDecay::DeselectAVolume: " << aVolume
288 <<
" is removed from list " <<
G4endl;
293 ed << aVolume <<
" is not in the list. No action taken." <<
G4endl;
294 G4Exception(
"G4RadioactiveDecay::DeselectAVolume()",
"HAD_RDM_300",
301 ed << aVolume <<
" is not a valid logical volume name. No action taken."
303 G4Exception(
"G4RadioactiveDecay::DeselectAVolume()",
"HAD_RDM_300",
318 for (std::size_t i = 0; i < theLogicalVolumes->size(); ++i){
319 volume = (*theLogicalVolumes)[i];
357 G4cout <<
"G4RadioactiveDecay::GetMeanLifeTime() " <<
G4endl;
359 <<
" GeV, Mass: " << theParticle->
GetMass()/GeV
360 <<
" GeV, Life time: " << theLife/
ns <<
" ns " <<
G4endl;
364 else if (theLife < 0.0) {meanlife =
DBL_MAX;}
365 else {meanlife = theLife;}
368 if (((
const G4Ions*)(theParticleDef))->GetExcitationEnergy() > 0. &&
369 meanlife ==
DBL_MAX) {meanlife = 0.;}
372 G4cout <<
" mean life time: " << meanlife/s <<
" s " <<
G4endl;
394 G4cout <<
"G4RadioactiveDecay::GetMeanFreePath() " <<
G4endl;
396 <<
" GeV, Mass: " << aMass/GeV <<
" GeV, tau: " << tau <<
" ns "
407 }
else if (tau < 0.0) {
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",
419 pathlength = c_light*tau*betaGamma;
425 G4cout <<
"G4Decay::GetMeanFreePath: "
427 <<
" stops, kinetic energy = "
437 G4cout <<
"mean free path: "<< pathlength/m <<
" m" <<
G4endl;
451 if (!isInitialised) {
452 isInitialised =
true;
469G4RadioactiveDecay::StreamInfo(std::ostream& os,
const G4String& endline)
476 G4long prec = os.precision(5);
477 os <<
"======================================================================"
479 os <<
"====== Radioactive Decay Physics Parameters ======="
481 os <<
"======================================================================"
483 os <<
"min MeanLife (from G4NuclideTable) "
484 <<
G4BestUnit(minMeanLife,
"Time") << endline;
485 os <<
"Max life time (from G4DeexPrecoParameters) "
487 os <<
"Internal e- conversion flag "
489 os <<
"Stored internal conversion coefficients "
491 os <<
"Enabled atomic relaxation mode "
492 << applyARM << endline;
493 os <<
"Enable correlated gamma emission "
495 os <<
"Max 2J for sampling of angular correlations "
497 os <<
"Atomic de-excitation enabled "
498 << emparam->
Fluo() << endline;
499 os <<
"Auger electron emission enabled "
500 << emparam->
Auger() << endline;
501 os <<
"Check EM cuts disabled for atomic de-excitation "
503 os <<
"Use Bearden atomic level energies "
505 os <<
"Use ANSTO fluorescence model "
507 os <<
"Threshold for very long decay time at rest "
508 <<
G4BestUnit(fThresholdForVeryLongDecayTime,
"Time") << endline;
509 os <<
"======================================================================"
527 G4int A = ((
const G4Ions*)(&theParentNucleus))->GetAtomicMass();
528 G4int Z = ((
const G4Ions*)(&theParentNucleus))->GetAtomicNumber();
530 G4double levelEnergy = ((
const G4Ions*)(&theParentNucleus))->GetExcitationEnergy();
532 ((
const G4Ions*)(&theParentNucleus))->GetFloatLevelBase();
534#ifdef G4MULTITHREADED
535 G4AutoLock lk(&G4RadioactiveDecay::radioactiveDecayMutex);
538 DecayTableMap::iterator master_table_ptr = master_dkmap->find(key);
540 if (master_table_ptr != master_dkmap->end() ) {
541 return master_table_ptr->second;
546 G4String file = theUserRadioactiveDataFiles[1000*
A+
Z];
549 std::ostringstream os;
550 os << dirPath <<
"/z" <<
Z <<
".a" <<
A <<
'\0';
557 std::ifstream DecaySchemeFile;
558 DecaySchemeFile.open(file);
560 if (DecaySchemeFile.good()) {
564 G4double modeTotalBR[nMode] = {0.0};
566 for (
G4int i = 0; i < nMode; i++) {
570 char inputChars[120]={
' '};
591 while (!complete && !DecaySchemeFile.getline(inputChars, 120).eof()) {
594 G4Exception(
"G4RadioactiveDecay::LoadDecayTable()",
"HAD_RDM_100",
599 inputLine = inputChars;
601 if (inputChars[0] !=
'#' && inputLine.length() != 0) {
602 std::istringstream tmpStream(inputLine);
604 if (inputChars[0] ==
'P') {
607 tmpStream >> recordType >> parentExcitation >> floatingFlag >> dummy;
615 found = (std::abs(parentExcitation*keV - levelEnergy) <
levelTolerance);
616 if (floatingLevel !=
noFloat) {
619 if (!floatMatch) found =
false;
628 if (inputLine.length() < 72) {
629 tmpStream >> theDecayMode >> dummy >> decayModeTotal;
630 switch (theDecayMode) {
636 anITChannel->
SetARM(applyARM);
637 theDecayTable->
Insert(anITChannel);
642 modeTotalBR[
BetaMinus] = decayModeTotal;
break;
644 modeTotalBR[
BetaPlus] = decayModeTotal;
break;
646 modeTotalBR[
KshellEC] = decayModeTotal;
break;
648 modeTotalBR[
LshellEC] = decayModeTotal;
break;
650 modeTotalBR[
MshellEC] = decayModeTotal;
break;
652 modeTotalBR[
NshellEC] = decayModeTotal;
break;
654 modeTotalBR[
Alpha] = decayModeTotal;
break;
656 modeTotalBR[
Proton] = decayModeTotal;
break;
658 modeTotalBR[
Neutron] = decayModeTotal;
break;
660 modeTotalBR[
SpFission] = decayModeTotal;
break;
674 modeTotalBR[
Triton] = decayModeTotal;
break;
678 G4Exception(
"G4RadioactiveDecay::LoadDecayTable()",
"HAD_RDM_000",
683 if (inputLine.length() < 84) {
684 tmpStream >> theDecayMode >> a >> daughterFloatFlag >> b >> c;
687 tmpStream >> theDecayMode >> a >> daughterFloatFlag >> b >> c >> betaType;
697 switch (theDecayMode) {
702 daughterFloatLevel, betaType);
705 theDecayTable->
Insert(aBetaMinusChannel);
714 daughterFloatLevel, betaType);
717 theDecayTable->
Insert(aBetaPlusChannel);
725 new G4ECDecay(&theParentNucleus, b, c*MeV, a*MeV,
729 aKECChannel->
SetARM(applyARM);
730 theDecayTable->
Insert(aKECChannel);
738 new G4ECDecay(&theParentNucleus, b, c*MeV, a*MeV,
742 aLECChannel->
SetARM(applyARM);
743 theDecayTable->
Insert(aLECChannel);
751 new G4ECDecay(&theParentNucleus, b, c*MeV, a*MeV,
755 aMECChannel->
SetARM(applyARM);
756 theDecayTable->
Insert(aMECChannel);
764 new G4ECDecay(&theParentNucleus, b, c*MeV, a*MeV,
768 aNECChannel->
SetARM(applyARM);
769 theDecayTable->
Insert(aNECChannel);
781 theDecayTable->
Insert(anAlphaChannel);
782 modeSumBR[
Alpha] += b;
793 theDecayTable->
Insert(aProtonChannel);
805 theDecayTable->
Insert(aNeutronChannel);
814 new G4SFDecay(&theParentNucleus, b, c*MeV, a*MeV,
816 theDecayTable->
Insert(aSpontFissChannel);
858 theDecayTable->
Insert(aTritonChannel);
866 G4Exception(
"G4RadioactiveDecay::LoadDecayTable()",
"HAD_RDM_000",
884 theNuclearDecayChannel =
static_cast<G4NuclearDecay*
>(theChannel);
887 if (theDecayMode !=
IT) {
888 theBR = theChannel->
GetBR();
889 theChannel->
SetBR(theBR*modeTotalBR[theDecayMode]/modeSumBR[theDecayMode]);
894 DecaySchemeFile.close();
896 if (!found && levelEnergy > 0) {
902 anITChannel->
SetARM(applyARM);
903 theDecayTable->
Insert(anITChannel);
910#ifdef G4MULTITHREADED
913 return theDecayTable;
921 std::ifstream DecaySchemeFile(filename);
922 if (DecaySchemeFile) {
924 theUserRadioactiveDataFiles[ID_ion] = filename;
927 ed << filename <<
" does not exist! " <<
G4endl;
928 G4Exception(
"G4RadioactiveDecay::AddUserDecayDataFile()",
"HAD_RDM_001",
954 G4cout <<
"G4RadioactiveDecay::DecayIt : "
956 <<
" is not selected for the RDM"<<
G4endl;
978 G4cout <<
"G4RadioactiveDecay::DecayIt : "
980 <<
" is not an ion or is outside (Z,A) limits set for the decay. "
981 <<
" Set particle change accordingly. "
996 if (theDecayTable == 0 || theDecayTable->
entries() == 0) {
1001 G4cout <<
"G4RadioactiveDecay::DecayIt : "
1002 <<
"decay table not defined for "
1004 <<
". Set particle change accordingly. "
1115 if (products->
entries() == 1) {
1148 if (temptime < 0.) temptime = 0.;
1149 finalGlobalTime += temptime;
1150 finalLocalTime += temptime;
1163 if ( finalGlobalTime > fThresholdForVeryLongDecayTime ) {
1172 products->
Boost(ParentEnergy, ParentDirection);
1179 G4cout <<
"G4RadioactiveDecay::DecayAnalog: Decay vertex :";
1180 G4cout <<
" Time: " << finalGlobalTime/
ns <<
"[ns]";
1185 G4cout <<
"G4Decay::DecayIt : decay products in Lab. Frame" <<
G4endl;
1191 G4int modelID = modelID_forIT + 10*theRadDecayMode;
1192 const G4int modelID_forAtomicRelaxation =
1194 for (
G4int index = 0; index < numberOfSecondaries; ++index ) {
1200 if ( theRadDecayMode ==
IT && index > 0 ) {
1201 if ( index == numberOfSecondaries-1 ) {
1207 index < numberOfSecondaries-1 ) {
1240 if (theDecayChannel == 0) {
1244 G4Exception(
"G4RadioactiveDecay::DoDecay",
"HAD_RDM_013",
1250 G4cout <<
"G4RadioactiveDecay::DoIt : selected decay channel addr: "
1251 << theDecayChannel <<
G4endl;
1254 theRadDecayMode = (
static_cast<G4NuclearDecay*
>(theDecayChannel))->GetDecayMode();
1269 if (origin == forceDecayDirection)
return;
1270 if (180.*deg == forceDecayHalfAngle)
return;
1271 if (0 == products || 0 == products->
entries())
return;
1291 if (daughterType == electron || daughterType == positron ||
1292 daughterType == neutron || daughterType == gamma ||
1293 daughterType == alpha || daughterType == triton || daughterType == proton)
CollimateDecayProduct(daughter);
1300 G4cout <<
"CollimateDecayProduct for daughter "
1313 if (origin == forceDecayDirection)
return origin;
1314 if (forceDecayHalfAngle == 180.*deg)
return origin;
1319 if (forceDecayHalfAngle > 0.) {
1322 G4double cosMin = std::cos(forceDecayHalfAngle);
1331 G4cout <<
" ChooseCollimationDirection returns " << dir <<
G4endl;
const char * G4FindDataDir(const char *)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
std::map< G4String, G4DecayTable * > DecayTableMap
@ G4RadioactiveDecayModeSize
std::map< G4String, G4DecayTable * > DecayTableMap
#define G4MUTEX_INITIALIZER
G4GLOB_DLL std::ostream G4cout
static G4Alpha * Definition()
G4DynamicParticle * PopProducts()
void Boost(G4double totalEnergy, const G4ThreeVector &momentumDirection)
G4VDecayChannel * GetDecayChannel(G4int index) const
void Insert(G4VDecayChannel *aChannel)
G4VDecayChannel * SelectADecayChannel(G4double parentMass=-1.)
G4bool CorrelatedGamma() const
G4bool StoreICLevelData() const
G4double GetMaxLifeTime() const
G4bool GetInternalConversionFlag() const
void SetMomentumDirection(const G4ThreeVector &aDirection)
const G4ThreeVector & GetMomentumDirection() const
const G4ParticleDefinition * GetParticleDefinition() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetTotalMomentum() const
void SetARM(G4bool onoff)
static G4Electron * Definition()
static G4EmParameters * Instance()
G4bool DeexcitationIgnoreCut() const
static G4Gamma * Definition()
static G4GenericIon * GenericIon()
static G4HadronicParameters * Instance()
void RegisterExtraProcess(G4VProcess *)
void RegisterParticleForExtraProcess(G4VProcess *, const G4ParticleDefinition *)
static G4HadronicProcessStore * Instance()
void SetARM(G4bool onoff)
static G4Ions::G4FloatLevelBase FloatLevelBase(char flbChar)
G4LogicalVolume * GetVolume(const G4String &name, G4bool verbose=true, G4bool reverseSearch=false) const
static G4LogicalVolumeStore * GetInstance()
const G4String & GetName() const
static G4Neutron * Definition()
G4RadioactiveDecayMode GetDecayMode()
G4DeexPrecoParameters * GetParameters()
static G4NuclearLevelData * GetInstance()
G4double GetMeanLifeThreshold()
static G4NuclideTable * GetInstance()
void Initialize(const G4Track &) final
void ProposeLocalTime(G4double t)
G4bool GetPDGStable() const
const G4String & GetParticleType() const
G4double GetPDGMass() const
G4double GetPDGLifeTime() const
const G4String & GetParticleName() const
void RDMForced(G4bool) override
void SetICM(G4bool) override
static G4int GetModelID(const G4int modelIndex)
static G4Positron * Definition()
static G4Proton * Definition()
void SelectAVolume(const G4String &aVolume)
void DecayAnalog(const G4Track &theTrack)
G4RadioactiveDecayMessenger * theRadioactiveDecayMessenger
void BuildPhysicsTable(const G4ParticleDefinition &)
std::vector< G4String > ValidVolumes
G4DecayTable * LoadDecayTable(const G4ParticleDefinition &theParentNucleus)
G4RadioactiveDecay(const G4String &processName="RadioactiveDecay")
G4double GetMeanLifeTime(const G4Track &theTrack, G4ForceCondition *condition)
G4DecayProducts * DoDecay(const G4ParticleDefinition &theParticleDef)
void CollimateDecayProduct(G4DynamicParticle *product)
G4DecayTable * GetDecayTable(const G4ParticleDefinition *)
static const G4double levelTolerance
G4ParticleChangeForRadDecay fParticleChangeForRadDecay
G4int GetVerboseLevel() const
G4double GetMeanFreePath(const G4Track &theTrack, G4double previousStepSize, G4ForceCondition *condition)
G4ThreeVector ChooseCollimationDirection() const
G4VParticleChange * DecayIt(const G4Track &theTrack, const G4Step &theStep)
G4bool IsApplicable(const G4ParticleDefinition &)
virtual void ProcessDescription(std::ostream &outFile) const
void DeselectAllVolumes()
void AddUserDecayDataFile(G4int Z, G4int A, G4String filename)
void DeselectAVolume(const G4String &aVolume)
void CollimateDecay(G4DecayProducts *products)
G4PhotonEvaporation * photonEvaporation
G4TrackStatus GetTrackStatus() const
G4VPhysicalVolume * GetVolume() 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)
static G4Triton * Definition()
void SetBR(G4double value)
virtual G4DecayProducts * DecayIt(G4double parentMass=-1.0)=0
void ProposeTrackStatus(G4TrackStatus status)
void ProposeWeight(G4double finalWeight)
void AddSecondary(G4Track *aSecondary)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4LogicalVolume * GetLogicalVolume() const
void ClearNumberOfInteractionLengthLeft()
void SetProcessSubType(G4int)
G4VParticleChange * pParticleChange
void rstrip(G4String &str, char ch=' ')
Remove trailing characters from string.