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