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