File indexing completed on 2025-09-16 08:56:46
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042 #ifndef G4RadioactiveDecay_h
0043 #define G4RadioactiveDecay_h 1
0044
0045 #include "G4VRadioactiveDecay.hh"
0046
0047 #include <vector>
0048 #include <map>
0049 #include <CLHEP/Units/SystemOfUnits.h>
0050
0051 #include "G4ios.hh"
0052 #include "globals.hh"
0053 #include "G4VRestDiscreteProcess.hh"
0054 #include "G4ParticleChangeForRadDecay.hh"
0055
0056 #include "G4NucleusLimits.hh"
0057 #include "G4RadioactiveDecayRatesToDaughter.hh"
0058 #include "G4RadioactiveDecayChainsFromParent.hh"
0059 #include "G4RadioactivityTable.hh"
0060 #include "G4ThreeVector.hh"
0061 #include "G4Threading.hh"
0062
0063 class G4Fragment;
0064 class G4RadioactivationMessenger;
0065
0066 typedef std::vector<G4RadioactiveDecayChainsFromParent> G4RadioactiveDecayParentChainTable;
0067 typedef std::vector<G4RadioactiveDecayRatesToDaughter> G4RadioactiveDecayRates;
0068 typedef std::map<G4String, G4DecayTable*> DecayTableMap;
0069
0070 class G4RadioactiveDecay : public G4VRadioactiveDecay
0071 {
0072 public:
0073
0074 G4RadioactiveDecay(const G4String& processName="RadioactiveDecay",
0075 const G4double timeThresholdForRadioactiveDecays=-1.0);
0076 ~G4RadioactiveDecay() override;
0077
0078 G4VParticleChange* DecayIt(const G4Track& theTrack,
0079 const G4Step& theStep) override;
0080
0081 void ProcessDescription(std::ostream& outFile) const override;
0082
0083
0084 void SetDecayBias(const G4String& filename);
0085
0086
0087 void SetHLThreshold(G4double hl) {halflifethreshold = hl;}
0088
0089 void SetSourceTimeProfile(const G4String& filename);
0090
0091
0092 G4bool IsRateTableReady(const G4ParticleDefinition &);
0093
0094
0095
0096
0097 void CalculateChainsFromParent(const G4ParticleDefinition&);
0098
0099
0100
0101
0102
0103 void GetChainsFromParent(const G4ParticleDefinition&);
0104
0105
0106
0107
0108
0109 void SetDecayRate(G4int,G4int,G4double, G4int, std::vector<G4double>&,
0110 std::vector<G4double>&);
0111
0112
0113
0114 std::vector<G4RadioactivityTable*>& GetTheRadioactivityTables()
0115 {return theRadioactivityTables;}
0116
0117
0118
0119
0120
0121
0122 inline void SetAnalogueMonteCarlo (G4bool r) {
0123 AnalogueMC = r;
0124 }
0125
0126
0127
0128 inline G4bool IsAnalogueMonteCarlo () {return AnalogueMC;}
0129
0130
0131 inline void SetBRBias(G4bool r) {
0132 BRBias = r;
0133 AnalogueMC = false;
0134 }
0135
0136
0137 inline void SetSplitNuclei(G4int r) {
0138 NSplit = r;
0139 AnalogueMC = false;
0140 }
0141
0142
0143 inline G4int GetSplitNuclei () {return NSplit;}
0144
0145 G4RadioactiveDecay(const G4RadioactiveDecay& right) = delete;
0146 G4RadioactiveDecay& operator=(const G4RadioactiveDecay& right) = delete;
0147
0148 protected:
0149
0150 G4double ConvolveSourceTimeProfile(const G4double, const G4double);
0151 G4double GetDecayTime();
0152 G4int GetDecayTimeBin(const G4double aDecayTime);
0153
0154 G4double GetMeanLifeTime(const G4Track& theTrack,
0155 G4ForceCondition* condition) override;
0156
0157
0158 void AddDeexcitationSpectrumForBiasMode(G4ParticleDefinition* apartDef,
0159 G4double weight,
0160 G4double currenTime,
0161 std::vector<double>& weights_v,
0162 std::vector<double>& times_v,
0163 std::vector<G4DynamicParticle*>& secondaries_v);
0164
0165 private:
0166
0167 G4RadioactivationMessenger* theRadioactivationMessenger;
0168 G4bool AnalogueMC;
0169 G4bool BRBias;
0170 G4int NSplit;
0171
0172 G4double halflifethreshold;
0173
0174 G4int NSourceBin;
0175 G4double SBin[100];
0176 G4double SProfile[100];
0177 G4int NDecayBin;
0178 G4double DBin[100];
0179 G4double DProfile[100];
0180
0181 G4RadioactiveDecayRatesToDaughter ratesToDaughter;
0182 G4RadioactiveDecayRates theDecayRateVector;
0183 G4RadioactiveDecayChainsFromParent chainsFromParent;
0184 G4RadioactiveDecayParentChainTable theParentChainTable;
0185
0186
0187 std::vector<G4RadioactivityTable*> theRadioactivityTables;
0188 G4int decayWindows[100];
0189 };
0190
0191 #endif
0192