Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4RadioactiveDecay.hh
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26////////////////////////////////////////////////////////////////////////////////
27// //
28// File: G4RadioactiveDecay.hh //
29// Author: D.H. Wright (SLAC) //
30// Date: 9 August 2017 //
31// Description: version the G4RadioactiveDecay process by F. Lei and //
32// P.R. Truscott with biasing and activation calculations //
33// removed to a derived class. It performs alpha, beta, //
34// electron capture and isomeric transition decays of //
35// radioactive nuclei. //
36// //
37////////////////////////////////////////////////////////////////////////////////
38
39#ifndef G4RadioactiveDecay_h
40#define G4RadioactiveDecay_h 1
41
42#include <vector>
43#include <map>
45
46#include "G4ios.hh"
47#include "globals.hh"
50
51#include "G4NucleusLimits.hh"
52#include "G4ThreeVector.hh"
54
55class G4Fragment;
58class G4Ions;
59class G4DecayTable;
60class G4ITDecay;
61
62typedef std::map<G4String, G4DecayTable*> DecayTableMap;
63
64
66{
67 // class description
68
69 // Implementation of the radioactive decay process which simulates the
70 // decays of radioactive nuclei. These nuclei are submitted to RDM as
71 // G4Ions. The required half-lives and decay schemes are retrieved from
72 // the Radioactivity database which was derived from ENSDF.
73 // All decay products are submitted back to the particle tracking process
74 // through the G4ParticleChangeForRadDecay object.
75 // class description - end
76
77 public: // with description
78
79 G4RadioactiveDecay(const G4String& processName="RadioactiveDecay",
80 const G4double timeThreshold=-1.0);
81 ~G4RadioactiveDecay() override;
82
84 // Return true if the specified isotope is
85 // 1) defined as "nucleus" and
86 // 2) it is within theNucleusLimit
87
88 G4VParticleChange* AtRestDoIt(const G4Track& theTrack,
89 const G4Step& theStep) override;
90
91 G4VParticleChange* PostStepDoIt(const G4Track& theTrack,
92 const G4Step& theStep) override;
93
94 void BuildPhysicsTable(const G4ParticleDefinition &) override;
95
96 void ProcessDescription(std::ostream& outFile) const override;
97
98 virtual G4VParticleChange* DecayIt(const G4Track& theTrack,
99 const G4Step& theStep);
100
101 // Return decay table if it exists, if not, load it from file
103
104 // Select a logical volume in which RDM applies
105 void SelectAVolume(const G4String& aVolume);
106
107 // Remove a logical volume from the RDM applied list
108 void DeselectAVolume(const G4String& aVolume);
109
110 // Select all logical volumes for the application of RDM
111 void SelectAllVolumes();
112
113 // Remove all logical volumes from RDM applications
114 void DeselectAllVolumes();
115
116 // Enable/disable ARM
117 void SetARM(G4bool arm) {applyARM = arm;}
118
120 // Load the decay data of isotope theParentNucleus
121
122 void AddUserDecayDataFile(G4int Z, G4int A, const G4String& filename);
123 // Allow the user to replace the radio-active decay data provided in Geant4
124 // by its own data file for a given isotope
125
126 inline void SetNucleusLimits(G4NucleusLimits theNucleusLimits1)
127 {theNucleusLimits = theNucleusLimits1 ;}
128 // Sets theNucleusLimits which specifies the range of isotopes
129 // the G4RadioactiveDecay applies.
130
131 // Returns theNucleusLimits which specifies the range of isotopes used
132 // by G4RadioactiveDecay
133 inline G4NucleusLimits GetNucleusLimits() const {return theNucleusLimits;}
134
135 inline void SetDecayDirection(const G4ThreeVector& theDir) {
136 forceDecayDirection = theDir.unit();
137 }
138
139 inline const G4ThreeVector& GetDecayDirection() const {
140 return forceDecayDirection;
141 }
142
143 inline void SetDecayHalfAngle(G4double halfAngle=0.*CLHEP::deg) {
144 forceDecayHalfAngle = std::min(std::max(0.*CLHEP::deg,halfAngle),180.*CLHEP::deg);
145 }
146
147 inline G4double GetDecayHalfAngle() const {return forceDecayHalfAngle;}
148
149 // Force direction (random within half-angle) for "visible" daughters
150 // (applies to electrons, positrons, gammas, neutrons, protons or alphas)
151 inline void SetDecayCollimation(const G4ThreeVector& theDir,
152 G4double halfAngle = 0.*CLHEP::deg) {
153 SetDecayDirection(theDir);
154 SetDecayHalfAngle(halfAngle);
155 }
156
157 // Ignore radioactive decays at rest of nuclides happening after this (very long) time threshold
158 inline void SetThresholdForVeryLongDecayTime(const G4double inputThreshold) {
159 fThresholdForVeryLongDecayTime = std::max( 0.0, inputThreshold );
160 }
161 inline G4double GetThresholdForVeryLongDecayTime() const {return fThresholdForVeryLongDecayTime;}
162
163 void StreamInfo(std::ostream& os, const G4String& endline);
164
167
168 protected:
169
170 G4double GetMeanFreePath(const G4Track& theTrack, G4double previousStepSize,
171 G4ForceCondition* condition) override;
172
173 G4double GetMeanLifeTime(const G4Track& theTrack,
174 G4ForceCondition* condition) override;
175
176 // sampling of products
177 void DecayAnalog(const G4Track& theTrack, G4DecayTable*);
178
179 // sampling products at rest
181
182 // Apply directional bias for "visible" daughters (e+-, gamma, n, p, alpha)
183 void CollimateDecay(G4DecayProducts* products);
186
187 // ParticleChange for decay process
189
193
194 std::vector<G4String> ValidVolumes;
196
198
199 // Library of decay tables
201
202 private:
203
204 G4NucleusLimits theNucleusLimits;
205
206 G4bool isInitialised{false};
207 G4bool applyARM{true};
208
209 // Parameters for pre-collimated (biased) decay products
210 G4ThreeVector forceDecayDirection{G4ThreeVector(0., 0., 0.)};
211 G4double forceDecayHalfAngle{0.0};
212 static const G4ThreeVector origin; // (0,0,0) for convenience
213
214 // Radioactive decay database directory path
215 static G4String dirPath;
216
217 // User define radioactive decay data files replacing some files in the G4RADECAY database
218 static std::map<G4int, G4String>* theUserRDataFiles;
219
220 // The last RadDecayMode
222
223 // Ignore radioactive decays at rest of nuclides happening after this (very long) time threshold
224 G4double fThresholdForVeryLongDecayTime;
225};
226
227#endif
228
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
std::map< G4String, G4DecayTable * > DecayTableMap
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
Hep3Vector unit() const
G4DecayTable * LoadDecayTable(const G4Ions *)
G4double GetMeanFreePath(const G4Track &theTrack, G4double previousStepSize, G4ForceCondition *condition) override
void AddUserDecayDataFile(G4int Z, G4int A, const G4String &filename)
void SelectAVolume(const G4String &aVolume)
G4RadioactiveDecay & operator=(const G4RadioactiveDecay &right)=delete
G4RadioactiveDecayMessenger * theRadioactiveDecayMessenger
void StreamInfo(std::ostream &os, const G4String &endline)
G4bool IsApplicable(const G4ParticleDefinition &) override
G4VParticleChange * AtRestDoIt(const G4Track &theTrack, const G4Step &theStep) override
static DecayTableMap * master_dkmap
std::vector< G4String > ValidVolumes
G4double GetDecayHalfAngle() const
void ProcessDescription(std::ostream &outFile) const override
G4VParticleChange * PostStepDoIt(const G4Track &theTrack, const G4Step &theStep) override
G4DecayProducts * DoDecay(const G4ParticleDefinition &, G4DecayTable *)
G4NucleusLimits GetNucleusLimits() const
void CollimateDecayProduct(G4DynamicParticle *product)
void DecayAnalog(const G4Track &theTrack, G4DecayTable *)
G4DecayTable * GetDecayTable(const G4ParticleDefinition *)
const G4ThreeVector & GetDecayDirection() const
static const G4double levelTolerance
G4ParticleChangeForRadDecay fParticleChangeForRadDecay
G4ThreeVector ChooseCollimationDirection() const
virtual G4VParticleChange * DecayIt(const G4Track &theTrack, const G4Step &theStep)
void SetDecayHalfAngle(G4double halfAngle=0.*CLHEP::deg)
void SetDecayDirection(const G4ThreeVector &theDir)
void SetARM(G4bool arm)
void SetDecayCollimation(const G4ThreeVector &theDir, G4double halfAngle=0.*CLHEP::deg)
G4double GetMeanLifeTime(const G4Track &theTrack, G4ForceCondition *condition) override
void SetThresholdForVeryLongDecayTime(const G4double inputThreshold)
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4RadioactiveDecay(const G4RadioactiveDecay &right)=delete
void SetNucleusLimits(G4NucleusLimits theNucleusLimits1)
void DeselectAVolume(const G4String &aVolume)
void CollimateDecay(G4DecayProducts *products)
G4PhotonEvaporation * photonEvaporation
G4RadioactiveDecay(const G4String &processName="RadioactiveDecay", const G4double timeThreshold=-1.0)
G4double GetThresholdForVeryLongDecayTime() const