Geant4 11.1.1
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"
53#include "G4Threading.hh"
55
56class G4Fragment;
59
60typedef std::map<G4String, G4DecayTable*> DecayTableMap;
61
62
64{
65 // class description
66
67 // Implementation of the radioactive decay process which simulates the
68 // decays of radioactive nuclei. These nuclei are submitted to RDM as
69 // G4Ions. The required half-lives and decay schemes are retrieved from
70 // the Radioactivity database which was derived from ENSDF.
71 // All decay products are submitted back to the particle tracking process
72 // through the G4ParticleChangeForRadDecay object.
73 // class description - end
74
75 public: // with description
76
77 G4RadioactiveDecay(const G4String& processName="RadioactiveDecay");
79
80 virtual void ProcessDescription(std::ostream& outFile) const;
81
83 // Return true if the specified isotope is
84 // 1) defined as "nucleus" and
85 // 2) it is within theNucleusLimit
86
87 // Return decay table if it exists, if not, load it from file
89
90 // Select a logical volume in which RDM applies
91 void SelectAVolume(const G4String& aVolume);
92
93 // Remove a logical volume from the RDM applied list
94 void DeselectAVolume(const G4String& aVolume);
95
96 // Select all logical volumes for the application of RDM
97 void SelectAllVolumes();
98
99 // Remove all logical volumes from RDM applications
100 void DeselectAllVolumes();
101
102 // Enable/disable ARM
103 void SetARM(G4bool arm) {applyARM = arm;}
104
105 G4DecayTable* LoadDecayTable(const G4ParticleDefinition& theParentNucleus);
106 // Load the decay data of isotope theParentNucleus
107
108 void AddUserDecayDataFile(G4int Z, G4int A,G4String filename);
109 // Allow the user to replace the radio-active decay data provided in Geant4
110 // by its own data file for a given isotope
111
112 inline void SetVerboseLevel(G4int value) {verboseLevel = value;}
113 // Sets the VerboseLevel which controls duggering display
114
115 inline G4int GetVerboseLevel() const {return verboseLevel;}
116 // Returns the VerboseLevel which controls level of debugging output
117
118 inline void SetNucleusLimits(G4NucleusLimits theNucleusLimits1)
119 {theNucleusLimits = theNucleusLimits1 ;}
120 // Sets theNucleusLimits which specifies the range of isotopes
121 // the G4RadioactiveDecay applies.
122
123 // Returns theNucleusLimits which specifies the range of isotopes used
124 // by G4RadioactiveDecay
125 inline G4NucleusLimits GetNucleusLimits() const {return theNucleusLimits;}
126
127 inline void SetDecayDirection(const G4ThreeVector& theDir) {
128 forceDecayDirection = theDir.unit();
129 }
130
131 inline const G4ThreeVector& GetDecayDirection() const {
132 return forceDecayDirection;
133 }
134
135 inline void SetDecayHalfAngle(G4double halfAngle=0.*CLHEP::deg) {
136 forceDecayHalfAngle = std::min(std::max(0.*CLHEP::deg,halfAngle),180.*CLHEP::deg);
137 }
138
139 inline G4double GetDecayHalfAngle() const {return forceDecayHalfAngle;}
140
141 // Force direction (random within half-angle) for "visible" daughters
142 // (applies to electrons, positrons, gammas, neutrons, protons or alphas)
143 inline void SetDecayCollimation(const G4ThreeVector& theDir,
144 G4double halfAngle = 0.*CLHEP::deg) {
145 SetDecayDirection(theDir);
146 SetDecayHalfAngle(halfAngle);
147 }
148
149 // Ignore radioactive decays at rest of nuclides happening after this (very long) time threshold
150 inline void SetThresholdForVeryLongDecayTime(const G4double inputThreshold) {
151 fThresholdForVeryLongDecayTime = std::max( 0.0, inputThreshold );
152 }
153 inline G4double GetThresholdForVeryLongDecayTime() const {return fThresholdForVeryLongDecayTime;}
154
156
157 G4VParticleChange* DecayIt(const G4Track& theTrack,
158 const G4Step& theStep);
159
160 protected:
161
162 void DecayAnalog(const G4Track& theTrack);
163
164 G4DecayProducts* DoDecay(const G4ParticleDefinition& theParticleDef);
165
166 // Apply directional bias for "visible" daughters (e+-, gamma, n, p, alpha)
167 void CollimateDecay(G4DecayProducts* products);
170
171 G4double GetMeanFreePath(const G4Track& theTrack, G4double previousStepSize,
173
174 G4double GetMeanLifeTime(const G4Track& theTrack,
176
177 // ParticleChange for decay process
179
182
183 std::vector<G4String> ValidVolumes;
185
187
188 // Library of decay tables
190#ifdef G4MULTITHREADED
191 static DecayTableMap* master_dkmap;
192#endif
193
194 private:
195
196 void StreamInfo(std::ostream& os, const G4String& endline);
197
199 G4RadioactiveDecay& operator=(const G4RadioactiveDecay& right);
200
201 G4NucleusLimits theNucleusLimits;
202
203 G4bool isInitialised;
204
205 G4bool applyARM;
206
207 // Parameters for pre-collimated (biased) decay products
208 G4ThreeVector forceDecayDirection;
209 G4double forceDecayHalfAngle;
210 static const G4ThreeVector origin; // (0,0,0) for convenience
211
212 // Radioactive decay database directory path
213 G4String dirPath;
214
215 //User define radioactive decay data files replacing some files in the G4RADECAY database
216 std::map<G4int, G4String> theUserRadioactiveDataFiles;
217
218 //The last RadDecayMode
219 G4RadioactiveDecayMode theRadDecayMode;
220
221// // Library of decay tables
222// DecayTableMap* dkmap;
223// #ifdef G4MULTITHREADED
224// static DecayTableMap* master_dkmap;
225// #endif
226
227 // Remainder of life time at rest
228 G4double fRemainderLifeTime;
229 G4int verboseLevel;
230
231 // Ignore radioactive decays at rest of nuclides happening after this (very long) time threshold
232 G4double fThresholdForVeryLongDecayTime;
233
234 // inline implementations
235 inline
236 G4double AtRestGetPhysicalInteractionLength(const G4Track& track,
238 {
239 fRemainderLifeTime =
241 return fRemainderLifeTime;
242 }
243
244 inline
245 G4VParticleChange* AtRestDoIt(const G4Track& theTrack,
246 const G4Step& theStep)
247 {return DecayIt(theTrack, theStep);}
248
249 inline
250 G4VParticleChange* PostStepDoIt(const G4Track& theTrack,
251 const G4Step& theStep)
252 {return DecayIt(theTrack, theStep);}
253
254#ifdef G4MULTITHREADED
255 public:
256 static G4Mutex radioactiveDecayMutex;
257 protected:
258 G4int& NumberOfInstances();
259#endif
260};
261
262#endif
263
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
G4RadioactiveDecayMode
std::map< G4String, G4DecayTable * > DecayTableMap
std::mutex G4Mutex
Definition: G4Threading.hh:81
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
Hep3Vector unit() const
void SelectAVolume(const G4String &aVolume)
void SetVerboseLevel(G4int value)
void DecayAnalog(const G4Track &theTrack)
G4RadioactiveDecayMessenger * theRadioactiveDecayMessenger
void BuildPhysicsTable(const G4ParticleDefinition &)
std::vector< G4String > ValidVolumes
G4double GetDecayHalfAngle() const
G4DecayTable * LoadDecayTable(const G4ParticleDefinition &theParentNucleus)
G4double GetMeanLifeTime(const G4Track &theTrack, G4ForceCondition *condition)
G4DecayProducts * DoDecay(const G4ParticleDefinition &theParticleDef)
G4NucleusLimits GetNucleusLimits() const
void CollimateDecayProduct(G4DynamicParticle *product)
G4DecayTable * GetDecayTable(const G4ParticleDefinition *)
const G4ThreeVector & GetDecayDirection() const
static const G4double levelTolerance
G4ParticleChangeForRadDecay fParticleChangeForRadDecay
G4int GetVerboseLevel() const
G4double GetMeanFreePath(const G4Track &theTrack, G4double previousStepSize, G4ForceCondition *condition)
G4ThreeVector ChooseCollimationDirection() const
G4VParticleChange * DecayIt(const G4Track &theTrack, const G4Step &theStep)
void SetDecayHalfAngle(G4double halfAngle=0.*CLHEP::deg)
void SetDecayDirection(const G4ThreeVector &theDir)
G4bool IsApplicable(const G4ParticleDefinition &)
void SetARM(G4bool arm)
void SetDecayCollimation(const G4ThreeVector &theDir, G4double halfAngle=0.*CLHEP::deg)
void SetThresholdForVeryLongDecayTime(const G4double inputThreshold)
virtual void ProcessDescription(std::ostream &outFile) const
void SetNucleusLimits(G4NucleusLimits theNucleusLimits1)
void AddUserDecayDataFile(G4int Z, G4int A, G4String filename)
void DeselectAVolume(const G4String &aVolume)
void CollimateDecay(G4DecayProducts *products)
G4PhotonEvaporation * photonEvaporation
G4double GetThresholdForVeryLongDecayTime() const
Definition: G4Step.hh:62
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)