Geant4-11
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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//
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// //
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 ICM
103 void SetICM(G4bool icm) {applyICM = icm;}
104
105 // Enable/disable ARM
106 void SetARM(G4bool arm) {applyARM = arm;}
107
108 G4DecayTable* LoadDecayTable(const G4ParticleDefinition& theParentNucleus);
109 // Load the decay data of isotope theParentNucleus
110
111 void AddUserDecayDataFile(G4int Z, G4int A,G4String filename);
112 // Allow the user to replace the radio-active decay data provided in Geant4
113 // by its own data file for a given isotope
114
115 inline void SetVerboseLevel(G4int value) {verboseLevel = value;}
116 // Sets the VerboseLevel which controls duggering display
117
118 inline G4int GetVerboseLevel() const {return verboseLevel;}
119 // Returns the VerboseLevel which controls level of debugging output
120
121 inline void SetNucleusLimits(G4NucleusLimits theNucleusLimits1)
122 {theNucleusLimits = theNucleusLimits1 ;}
123 // Sets theNucleusLimits which specifies the range of isotopes
124 // the G4RadioactiveDecay applies.
125
126 // Returns theNucleusLimits which specifies the range of isotopes used
127 // by G4RadioactiveDecay
129
130 inline void SetDecayDirection(const G4ThreeVector& theDir) {
131 forceDecayDirection = theDir.unit();
132 }
133
134 inline const G4ThreeVector& GetDecayDirection() const {
135 return forceDecayDirection;
136 }
137
138 inline void SetDecayHalfAngle(G4double halfAngle=0.*CLHEP::deg) {
140 }
141
143
144 // Force direction (random within half-angle) for "visible" daughters
145 // (applies to electrons, positrons, gammas, neutrons, protons or alphas)
146 inline void SetDecayCollimation(const G4ThreeVector& theDir,
147 G4double halfAngle = 0.*CLHEP::deg) {
148 SetDecayDirection(theDir);
149 SetDecayHalfAngle(halfAngle);
150 }
151
152 // Ignore radioactive decays at rest of nuclides happening after this (very long) time threshold
153 inline void SetThresholdForVeryLongDecayTime(const G4double inputThreshold) {
154 fThresholdForVeryLongDecayTime = std::max( 0.0, inputThreshold );
155 }
157
159
160 G4VParticleChange* DecayIt(const G4Track& theTrack,
161 const G4Step& theStep);
162
163 protected:
164
165 void DecayAnalog(const G4Track& theTrack);
166
167 G4DecayProducts* DoDecay(const G4ParticleDefinition& theParticleDef);
168
169 // Apply directional bias for "visible" daughters (e+-, gamma, n, p, alpha)
170 void CollimateDecay(G4DecayProducts* products);
173
174 G4double GetMeanFreePath(const G4Track& theTrack, G4double previousStepSize,
176
177 G4double GetMeanLifeTime(const G4Track& theTrack,
179
180 // ParticleChange for decay process
182
185
186 std::vector<G4String> ValidVolumes;
188
190
191 // Library of decay tables
193#ifdef G4MULTITHREADED
194 static DecayTableMap* master_dkmap;
195#endif
196
197 private:
198
199 void StreamInfo(std::ostream& os, const G4String& endline);
200
203
205
207
210
211 // Parameters for pre-collimated (biased) decay products
214 static const G4ThreeVector origin; // (0,0,0) for convenience
215
216 // Radioactive decay database directory path
218
219 //User define radioactive decay data files replacing some files in the G4RADECAY database
220 std::map<G4int, G4String> theUserRadioactiveDataFiles;
221
222 //The last RadDecayMode
224
225// // Library of decay tables
226// DecayTableMap* dkmap;
227// #ifdef G4MULTITHREADED
228// static DecayTableMap* master_dkmap;
229// #endif
230
231 // Remainder of life time at rest
234
235 // Ignore radioactive decays at rest of nuclides happening after this (very long) time threshold
237
238 // inline implementations
239 inline
242 {
245 return fRemainderLifeTime;
246 }
247
248 inline
250 const G4Step& theStep)
251 {return DecayIt(theTrack, theStep);}
252
253 inline
255 const G4Step& theStep)
256 {return DecayIt(theTrack, theStep);}
257
258#ifdef G4MULTITHREADED
259 public:
260 static G4Mutex radioactiveDecayMutex;
261 protected:
262 G4int& NumberOfInstances();
263#endif
264};
265
266#endif
267
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
std::map< G4String, G4DecayTable * > DecayTableMap
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
static const G4ThreeVector origin
std::map< G4int, G4String > theUserRadioactiveDataFiles
void SelectAVolume(const G4String &aVolume)
void SetVerboseLevel(G4int value)
void DecayAnalog(const G4Track &theTrack)
G4RadioactiveDecayMessenger * theRadioactiveDecayMessenger
void StreamInfo(std::ostream &os, const G4String &endline)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4VParticleChange * PostStepDoIt(const G4Track &theTrack, const G4Step &theStep)
void SetICM(G4bool icm)
std::vector< G4String > ValidVolumes
G4double GetDecayHalfAngle() const
G4DecayTable * LoadDecayTable(const G4ParticleDefinition &theParentNucleus)
G4RadioactiveDecay(const G4String &processName="RadioactiveDecay")
G4double GetMeanLifeTime(const G4Track &theTrack, G4ForceCondition *condition)
G4DecayProducts * DoDecay(const G4ParticleDefinition &theParticleDef)
G4NucleusLimits GetNucleusLimits() const
void CollimateDecayProduct(G4DynamicParticle *product)
G4NucleusLimits theNucleusLimits
G4DecayTable * GetDecayTable(const G4ParticleDefinition *)
const G4ThreeVector & GetDecayDirection() const
static const G4double levelTolerance
G4ParticleChangeForRadDecay fParticleChangeForRadDecay
G4ThreeVector forceDecayDirection
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)
G4RadioactiveDecay & operator=(const G4RadioactiveDecay &right)
G4bool IsApplicable(const G4ParticleDefinition &)
G4RadioactiveDecayMode theRadDecayMode
void SetARM(G4bool arm)
G4VParticleChange * AtRestDoIt(const G4Track &theTrack, const G4Step &theStep)
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)
G4double fThresholdForVeryLongDecayTime
void DeselectAVolume(const G4String &aVolume)
G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition)
void CollimateDecay(G4DecayProducts *products)
G4RadioactiveDecay(const G4RadioactiveDecay &right)
G4PhotonEvaporation * photonEvaporation
G4double GetThresholdForVeryLongDecayTime() const
Definition: G4Step.hh:62
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
static constexpr double deg
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments