Geant4  10.01.p03
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 #ifndef G4RadioactiveDecay_h
27 #define G4RadioactiveDecay_h 1
28 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
29 //
30 // MODULE: G4RadioactiveDecay.hh
31 //
32 // Version: 0.b.4
33 // Date: 14/04/00
34 // Author: F Lei & P R Truscott
35 // Organisation: DERA UK
36 // Customer: ESA/ESTEC, NOORDWIJK
37 // Contract: 12115/96/JG/NL Work Order No. 3
38 //
39 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
40 //
41 // CHANGE HISTORY
42 // --------------
43 // 17 October 2011, L Desorgher - Add the method AddUserDecayDataFile
44 //
45 // 01 June 2011, M. Kelsey -- Add directional biasing interface to allow for
46 // "collimation" of decay daughters.
47 //
48 // 29 February 2000, P R Truscott, DERA UK
49 // 0.b.3 release.
50 //
51 // 13 April 2000, F Lei, DERA UK
52 // 0.b.4 release. No change to this file
53 //
54 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
56 
57 #include <vector>
58 #include <map>
59 #include <CLHEP/Units/SystemOfUnits.h>
60 
61 #include "G4ios.hh"
62 #include "globals.hh"
65 // #include "G4RadioactiveDecaymessenger.hh"
66 
67 #include "G4NucleusLimits.hh"
70 #include "G4RIsotopeTable.hh"
71 #include "G4RadioactivityTable.hh"
72 #include "G4ThreeVector.hh"
73 #include "G4Threading.hh"
74 
76 
77 typedef std::vector<G4RadioactiveDecayRateVector> G4RadioactiveDecayRateTable;
78 typedef std::vector<G4RadioactiveDecayRate> G4RadioactiveDecayRates;
79 typedef std::map<G4String, G4DecayTable*> DecayTableMap;
80 
81 
83 {
84  // class description
85 
86  // Implementation of the radioactive decay process which simulates the
87  // decays of radioactive nuclei. These nuclei are submitted to RDM as
88  // G4Ions. The required half-lives and decay schemes are retrieved from
89  // the Radioactivity database which was derived from ENSDF.
90  // All decay products are submitted back to the particle tracking process
91  // through the G4ParticleChangeForRadDecay object.
92  // class description - end
93 
94  public: // with description
95 
96  G4RadioactiveDecay(const G4String& processName="RadioactiveDecay");
98 
99  // Return true if the specified isotope is
100  // 1) defined as "nucleus" and
101  // 2) it is within theNucleusLimit
103 
104  // Return decay table if it exists, if not, load it from file
106 
107  // Select a logical volume in which RDM applies
108  void SelectAVolume(const G4String aVolume);
109 
110  // Remove a logical volume from the RDM applied list
111  void DeselectAVolume(const G4String aVolume);
112 
113  // Select all logical volumes for the application of RDM
114  void SelectAllVolumes();
115 
116  // Remove all logical volumes from RDM applications
117  void DeselectAllVolumes();
118 
119  // Set the decay biasing scheme using the data in "filename"
120  void SetDecayBias(G4String filename);
121 
122  // Set the half-life threshold for isomer production
124 
125  // Enable/disable ICM
126  void SetICM(G4bool icm) {applyICM = icm;}
127 
128  // Enable/disable ARM
129  void SetARM(G4bool arm) {applyARM = arm;}
130 
131  // Set source exposure function using histograms in "filename"
132  void SetSourceTimeProfile(G4String filename);
133 
135  // Returns true if the coefficient and decay time table for all the
136  // descendants of the specified isotope are ready.
137  // used in VR decay mode only
138 
140  // Calculates the coefficient and decay time table for all the descendents
141  // of the specified isotope. Adds the calculated table to the private data
142  // member "theDecayRateTableVector".
143  // used in VR decay mode only
144 
146  // Used to retrieve the coefficient and decay time table for all the
147  // descendants of the specified isotope from "theDecayRateTableVector"
148  // and place it in "theDecayRateTable".
149  // used in VR decay mode only
150 
151  void SetDecayRate(G4int,G4int,G4double, G4int, std::vector<G4double>,
152  std::vector<G4double>);
153  // Sets "theDecayRate" with data supplied in the arguements.
154  // used in VR decay mode only
155 
156  std::vector<G4RadioactivityTable*> GetTheRadioactivityTables()
157  {return theRadioactivityTables;}
158  // Return vector of G4Radioactivity map - should be used in VR mode only
159 
160  G4DecayTable* LoadDecayTable(const G4ParticleDefinition& theParentNucleus);
161  // Load the decay data of isotope theParentNucleus
162 
163  void AddUserDecayDataFile(G4int Z, G4int A,G4String filename);
164  // Allow the user to replace the radio-active decay data provided in Geant4
165  // by its own data file for a given isotope
166 
167 
168  inline void SetVerboseLevel(G4int value) {verboseLevel = value;}
169  // Sets the VerboseLevel which controls duggering display
170 
171  inline G4int GetVerboseLevel() const {return verboseLevel;}
172  // Returns the VerboseLevel which controls level of debugging output
173 
174  inline void SetNucleusLimits(G4NucleusLimits theNucleusLimits1)
175  {theNucleusLimits = theNucleusLimits1 ;}
176  // Sets theNucleusLimits which specifies the range of isotopes
177  // the G4RadioactiveDecay applies.
178 
180  {return theNucleusLimits;}
181  // Returns theNucleusLimits which specifies the range of isotopes
182  // the G4RadioactiveDecay applies
183 
184  inline void SetAnalogueMonteCarlo (G4bool r ) {
185  AnalogueMC = r;
187  }
188  // Controls whether G4RadioactiveDecay runs in analogue mode or
189  // variance reduction mode.
190 
191  inline void SetFBeta (G4bool r ) { FBeta = r; }
192  // Controls whether G4RadioactiveDecay uses fast beta simulation mode
193 
195  // Returns true if the simulation is an analogue Monte Carlo, and false if
196  // any of the biassing schemes have been selected.
197 
198  inline void SetBRBias (G4bool r) {
199  BRBias = r;
201  }
202  // Sets whether branching ration bias scheme applies.
203 
204  inline void SetSplitNuclei (G4int r) {
205  NSplit = r;
207  }
208  // Sets the N number for the Nuclei spliting bias scheme
209 
210  inline G4int GetSplitNuclei () {return NSplit;}
211  // Returns the N number used for the Nuclei spliting bias scheme
212 
213  inline void SetDecayDirection(const G4ThreeVector& theDir) {
214  forceDecayDirection = theDir.unit();
215  }
216 
217  inline const G4ThreeVector& GetDecayDirection() const {
218  return forceDecayDirection;
219  }
220 
221  inline void SetDecayHalfAngle(G4double halfAngle=0.*CLHEP::deg) {
223  }
224 
226 
227  inline void SetDecayCollimation(const G4ThreeVector& theDir,
228  G4double halfAngle = 0.*CLHEP::deg) {
229  SetDecayDirection(theDir);
230  SetDecayHalfAngle(halfAngle);
231  }
232 
233  // Force direction (random within half-angle) for "visible" daughters
234  // (applies to electrons, positrons, gammas, neutrons, protons or alphas)
235 
237 
238  protected:
239 
240  G4VParticleChange* DecayIt(const G4Track& theTrack,
241  const G4Step& theStep);
242 
243  G4DecayProducts* DoDecay(const G4ParticleDefinition& theParticleDef);
244 
245  // Apply directional bias for "visible" daughters (e+-, gamma, n, p, alpha)
246  void CollimateDecay(G4DecayProducts* products);
249 
250  G4double GetMeanFreePath(const G4Track& theTrack, G4double previousStepSize,
252 
253  G4double GetMeanLifeTime(const G4Track& theTrack,
255 
256  G4double GetTaoTime(const G4double, const G4double);
257 
259 
260  G4int GetDecayTimeBin(const G4double aDecayTime);
261 
262  private:
263 
266 
269 
271 
277 
281 
282  // Parameters for pre-collimated (biased) decay products
285  static const G4ThreeVector origin; // (0,0,0) for convenience
286 
293 
294  std::vector<G4String> ValidVolumes;
296 
301 
302  // for the radioactivity tables
303  std::vector<G4RadioactivityTable*> theRadioactivityTables;
305  static const G4double levelTolerance;
306 
307  //User define radioactive decay data files replacing some files in the G4RADECAY database
308  std::map<G4int, G4String> theUserRadioactiveDataFiles;
309 
310  // Library of decay tables
312 #ifdef G4MULTITHREADED
313  static DecayTableMap* master_dkmap;
314 #endif
315 
316  // Remainder of life time at rest
319 
320 
321  // ParticleChange for decay process
323 
324  // inline implementations
325  inline
328  {
329  fRemainderLifeTime =
331  return fRemainderLifeTime;
332  }
333 
334  inline
336  const G4Step& theStep)
337  {return DecayIt(theTrack, theStep);}
338 
339  inline
341  const G4Step& theStep)
342  {return DecayIt(theTrack, theStep);}
343 
344 #ifdef G4MULTITHREADED
345  public:
346  static G4Mutex radioactiveDecayMutex;
347 #endif
348 };
349 
350 #endif
351 
G4VParticleChange * PostStepDoIt(const G4Track &theTrack, const G4Step &theStep)
G4double condition(const G4ErrorSymMatrix &m)
std::vector< G4RadioactiveDecayRateVector > G4RadioactiveDecayRateTable
G4RadioactiveDecayRateVector theDecayRateTable
G4DecayProducts * DoDecay(const G4ParticleDefinition &theParticleDef)
G4RadioactiveDecay & operator=(const G4RadioactiveDecay &right)
void SelectAVolume(const G4String aVolume)
std::map< G4String, G4DecayTable * > DecayTableMap
CLHEP::Hep3Vector G4ThreeVector
G4bool IsRateTableReady(const G4ParticleDefinition &)
std::vector< G4String > ValidVolumes
std::vector< G4RadioactiveDecayRate > G4RadioactiveDecayRates
G4VParticleChange * AtRestDoIt(const G4Track &theTrack, const G4Step &theStep)
G4RadioactiveDecayRates theDecayRateVector
G4double GetMeanLifeTime(const G4Track &theTrack, G4ForceCondition *condition)
G4int GetVerboseLevel() const
void SetSourceTimeProfile(G4String filename)
G4bool IsApplicable(const G4ParticleDefinition &)
int G4int
Definition: G4Types.hh:78
void SetBRBias(G4bool r)
G4ParticleChangeForRadDecay fParticleChangeForRadDecay
static const double s
Definition: G4SIunits.hh:150
static const G4ThreeVector origin
void SetDecayRate(G4int, G4int, G4double, G4int, std::vector< G4double >, std::vector< G4double >)
static const double deg
Definition: G4SIunits.hh:133
std::vector< G4RadioactivityTable * > GetTheRadioactivityTables()
bool G4bool
Definition: G4Types.hh:79
void CollimateDecayProduct(G4DynamicParticle *product)
const G4ThreeVector & GetDecayDirection() const
void SetARM(G4bool arm)
void SetDecayDirection(const G4ThreeVector &theDir)
G4ThreeVector ChooseCollimationDirection() const
G4int GetDecayTimeBin(const G4double aDecayTime)
Definition: G4Step.hh:76
G4NucleusLimits GetNucleusLimits() const
static const G4double A[nN]
G4double GetTaoTime(const G4double, const G4double)
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
G4NucleusLimits theNucleusLimits
G4DecayTable * GetDecayTable(const G4ParticleDefinition *)
G4double AtRestGetPhysicalInteractionLength(const G4Track &track, G4ForceCondition *condition)
void SetNucleusLimits(G4NucleusLimits theNucleusLimits1)
G4double GetMeanFreePath(const G4Track &theTrack, G4double previousStepSize, G4ForceCondition *condition)
G4RadioactiveDecayRate theDecayRate
static const G4double levelTolerance
G4int G4Mutex
Definition: G4Threading.hh:173
G4RadioactiveDecayRateTable theDecayRateTableVector
std::map< G4int, G4String > theUserRadioactiveDataFiles
void GetDecayRateTable(const G4ParticleDefinition &)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4RadioactiveDecaymessenger * theRadioactiveDecaymessenger
void SetSplitNuclei(G4int r)
G4double GetDecayHalfAngle() const
void AddDecayRateTable(const G4ParticleDefinition &)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4DecayTable * LoadDecayTable(const G4ParticleDefinition &theParentNucleus)
std::vector< G4RadioactivityTable * > theRadioactivityTables
G4RadioactiveDecay(const G4String &processName="RadioactiveDecay")
void SetVerboseLevel(G4int value)
G4ThreeVector forceDecayDirection
void SetDecayCollimation(const G4ThreeVector &theDir, G4double halfAngle=0.*CLHEP::deg)
void SetICM(G4bool icm)
double G4double
Definition: G4Types.hh:76
void SetHLThreshold(G4double hl)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4ForceCondition
void SetAnalogueMonteCarlo(G4bool r)
void SetDecayHalfAngle(G4double halfAngle=0.*CLHEP::deg)
void DeselectAVolume(const G4String aVolume)
G4RIsotopeTable * theIsotopeTable
void AddUserDecayDataFile(G4int Z, G4int A, G4String filename)
void CollimateDecay(G4DecayProducts *products)
void SetDecayBias(G4String filename)
G4VParticleChange * DecayIt(const G4Track &theTrack, const G4Step &theStep)