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
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039 #ifndef G4RadioactiveDecay_h
0040 #define G4RadioactiveDecay_h 1
0041
0042 #include <vector>
0043 #include <map>
0044 #include <CLHEP/Units/SystemOfUnits.h>
0045
0046 #include "G4ios.hh"
0047 #include "globals.hh"
0048 #include "G4VRestDiscreteProcess.hh"
0049 #include "G4ParticleChangeForRadDecay.hh"
0050
0051 #include "G4NucleusLimits.hh"
0052 #include "G4ThreeVector.hh"
0053 #include "G4RadioactiveDecayMode.hh"
0054
0055 class G4Fragment;
0056 class G4RadioactiveDecayMessenger;
0057 class G4PhotonEvaporation;
0058 class G4Ions;
0059 class G4DecayTable;
0060 class G4ITDecay;
0061
0062 typedef std::map<G4String, G4DecayTable*> DecayTableMap;
0063
0064
0065 class G4RadioactiveDecay : public G4VRestDiscreteProcess
0066 {
0067
0068
0069
0070
0071
0072
0073
0074
0075
0076
0077 public:
0078
0079 G4RadioactiveDecay(const G4String& processName="RadioactiveDecay",
0080 const G4double timeThreshold=-1.0);
0081 ~G4RadioactiveDecay() override;
0082
0083 G4bool IsApplicable(const G4ParticleDefinition&) override;
0084
0085
0086
0087
0088 G4VParticleChange* AtRestDoIt(const G4Track& theTrack,
0089 const G4Step& theStep) override;
0090
0091 G4VParticleChange* PostStepDoIt(const G4Track& theTrack,
0092 const G4Step& theStep) override;
0093
0094 void BuildPhysicsTable(const G4ParticleDefinition &) override;
0095
0096 void ProcessDescription(std::ostream& outFile) const override;
0097
0098 virtual G4VParticleChange* DecayIt(const G4Track& theTrack,
0099 const G4Step& theStep);
0100
0101
0102 G4DecayTable* GetDecayTable(const G4ParticleDefinition*);
0103
0104
0105 void SelectAVolume(const G4String& aVolume);
0106
0107
0108 void DeselectAVolume(const G4String& aVolume);
0109
0110
0111 void SelectAllVolumes();
0112
0113
0114 void DeselectAllVolumes();
0115
0116
0117 void SetARM(G4bool arm) {applyARM = arm;}
0118
0119 G4DecayTable* LoadDecayTable(const G4Ions*);
0120
0121
0122 void AddUserDecayDataFile(G4int Z, G4int A, const G4String& filename);
0123
0124
0125
0126 inline void SetNucleusLimits(G4NucleusLimits theNucleusLimits1)
0127 {theNucleusLimits = theNucleusLimits1 ;}
0128
0129
0130
0131
0132
0133 inline G4NucleusLimits GetNucleusLimits() const {return theNucleusLimits;}
0134
0135 inline void SetDecayDirection(const G4ThreeVector& theDir) {
0136 forceDecayDirection = theDir.unit();
0137 }
0138
0139 inline const G4ThreeVector& GetDecayDirection() const {
0140 return forceDecayDirection;
0141 }
0142
0143 inline void SetDecayHalfAngle(G4double halfAngle=0.*CLHEP::deg) {
0144 forceDecayHalfAngle = std::min(std::max(0.*CLHEP::deg,halfAngle),180.*CLHEP::deg);
0145 }
0146
0147 inline G4double GetDecayHalfAngle() const {return forceDecayHalfAngle;}
0148
0149
0150
0151 inline void SetDecayCollimation(const G4ThreeVector& theDir,
0152 G4double halfAngle = 0.*CLHEP::deg) {
0153 SetDecayDirection(theDir);
0154 SetDecayHalfAngle(halfAngle);
0155 }
0156
0157
0158 inline void SetThresholdForVeryLongDecayTime(const G4double inputThreshold) {
0159 fThresholdForVeryLongDecayTime = std::max( 0.0, inputThreshold );
0160 }
0161 inline G4double GetThresholdForVeryLongDecayTime() const {return fThresholdForVeryLongDecayTime;}
0162
0163 void StreamInfo(std::ostream& os, const G4String& endline);
0164
0165 G4RadioactiveDecay(const G4RadioactiveDecay& right) = delete;
0166 G4RadioactiveDecay& operator=(const G4RadioactiveDecay& right) = delete;
0167
0168 protected:
0169
0170 G4double GetMeanFreePath(const G4Track& theTrack, G4double previousStepSize,
0171 G4ForceCondition* condition) override;
0172
0173 G4double GetMeanLifeTime(const G4Track& theTrack,
0174 G4ForceCondition* condition) override;
0175
0176
0177 void DecayAnalog(const G4Track& theTrack, G4DecayTable*);
0178
0179
0180 G4DecayProducts* DoDecay(const G4ParticleDefinition&, G4DecayTable*);
0181
0182
0183 void CollimateDecay(G4DecayProducts* products);
0184 void CollimateDecayProduct(G4DynamicParticle* product);
0185 G4ThreeVector ChooseCollimationDirection() const;
0186
0187
0188 G4ParticleChangeForRadDecay fParticleChangeForRadDecay;
0189
0190 G4RadioactiveDecayMessenger* theRadioactiveDecayMessenger;
0191 G4PhotonEvaporation* photonEvaporation;
0192 G4ITDecay* decayIT;
0193
0194 std::vector<G4String> ValidVolumes;
0195 bool isAllVolumesMode{true};
0196
0197 static const G4double levelTolerance;
0198
0199
0200 static DecayTableMap* master_dkmap;
0201
0202 private:
0203
0204 G4NucleusLimits theNucleusLimits;
0205
0206 G4bool isInitialised{false};
0207 G4bool applyARM{true};
0208
0209
0210 G4ThreeVector forceDecayDirection{G4ThreeVector(0., 0., 0.)};
0211 G4double forceDecayHalfAngle{0.0};
0212 static const G4ThreeVector origin;
0213
0214
0215 static G4String dirPath;
0216
0217
0218 static std::map<G4int, G4String>* theUserRDataFiles;
0219
0220
0221 G4RadioactiveDecayMode theRadDecayMode{G4RadioactiveDecayMode::IT};
0222
0223
0224 G4double fThresholdForVeryLongDecayTime;
0225 };
0226
0227 #endif
0228