Geant4  10.03.p01
 All Classes 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 //
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>
60 
61 #include "G4ios.hh"
62 #include "globals.hh"
65 
66 #include "G4NucleusLimits.hh"
69 #include "G4RadioactivityTable.hh"
70 #include "G4ThreeVector.hh"
71 #include "G4Threading.hh"
72 
73 class G4Fragment;
75 
76 typedef std::vector<G4RadioactiveDecayRateVector> G4RadioactiveDecayRateTable;
77 typedef std::vector<G4RadioactiveDecayRate> G4RadioactiveDecayRates;
78 typedef std::map<G4String, G4DecayTable*> DecayTableMap;
79 
80 
82 {
83  // class description
84 
85  // Implementation of the radioactive decay process which simulates the
86  // decays of radioactive nuclei. These nuclei are submitted to RDM as
87  // G4Ions. The required half-lives and decay schemes are retrieved from
88  // the Radioactivity database which was derived from ENSDF.
89  // All decay products are submitted back to the particle tracking process
90  // through the G4ParticleChangeForRadDecay object.
91  // class description - end
92 
93  public: // with description
94 
95  G4RadioactiveDecay(const G4String& processName="RadioactiveDecay");
97 
98  // Return true if the specified isotope is
99  // 1) defined as "nucleus" and
100  // 2) it is within theNucleusLimit
102 
103  // Return decay table if it exists, if not, load it from file
105 
106  // Select a logical volume in which RDM applies
107  void SelectAVolume(const G4String aVolume);
108 
109  // Remove a logical volume from the RDM applied list
110  void DeselectAVolume(const G4String aVolume);
111 
112  // Select all logical volumes for the application of RDM
113  void SelectAllVolumes();
114 
115  // Remove all logical volumes from RDM applications
116  void DeselectAllVolumes();
117 
118  // Set the decay biasing scheme using the data in "filename"
119  void SetDecayBias(G4String filename);
120 
121  // Set the half-life threshold for isomer production
122  void SetHLThreshold(G4double hl) {halflifethreshold = hl;}
123 
124  // Enable/disable ICM
125  void SetICM(G4bool icm) {applyICM = icm;}
126 
127  // Enable/disable ARM
128  void SetARM(G4bool arm) {applyARM = arm;}
129 
130  // Set source exposure function using histograms in "filename"
131  void SetSourceTimeProfile(G4String filename);
132 
134  // Returns true if the coefficient and decay time table for all the
135  // descendants of the specified isotope are ready.
136  // used in VR decay mode only
137 
139  // Calculates the coefficient and decay time table for all the descendents
140  // of the specified isotope. Adds the calculated table to the private data
141  // member "theDecayRateTableVector".
142  // used in VR decay mode only
143 
145  // Used to retrieve the coefficient and decay time table for all the
146  // descendants of the specified isotope from "theDecayRateTableVector"
147  // and place it in "theDecayRateTable".
148  // used in VR decay mode only
149 
150  void SetDecayRate(G4int,G4int,G4double, G4int, std::vector<G4double>,
151  std::vector<G4double>);
152  // Sets "theDecayRate" with data supplied in the arguements.
153  // used in VR decay mode only
154 
155  std::vector<G4RadioactivityTable*> GetTheRadioactivityTables()
156  {return theRadioactivityTables;}
157  // Return vector of G4Radioactivity map - should be used in VR mode only
158 
159  G4DecayTable* LoadDecayTable(const G4ParticleDefinition& theParentNucleus);
160  // Load the decay data of isotope theParentNucleus
161 
162  void AddUserDecayDataFile(G4int Z, G4int A,G4String filename);
163  // Allow the user to replace the radio-active decay data provided in Geant4
164  // by its own data file for a given isotope
165 
166 
167  inline void SetVerboseLevel(G4int value) {verboseLevel = value;}
168  // Sets the VerboseLevel which controls duggering display
169 
170  inline G4int GetVerboseLevel() const {return verboseLevel;}
171  // Returns the VerboseLevel which controls level of debugging output
172 
173  inline void SetNucleusLimits(G4NucleusLimits theNucleusLimits1)
174  {theNucleusLimits = theNucleusLimits1 ;}
175  // Sets theNucleusLimits which specifies the range of isotopes
176  // the G4RadioactiveDecay applies.
177 
179  {return theNucleusLimits;}
180  // Returns theNucleusLimits which specifies the range of isotopes
181  // the G4RadioactiveDecay applies
182 
183  inline void SetAnalogueMonteCarlo (G4bool r ) {
184  AnalogueMC = r;
185  if (!AnalogueMC) halflifethreshold = 1e-6*CLHEP::s;
186  }
187  // Controls whether G4RadioactiveDecay runs in analogue mode or
188  // variance reduction mode.
189 
190  inline void SetFBeta (G4bool r ) { FBeta = r; }
191  // Controls whether G4RadioactiveDecay uses fast beta simulation mode
192 
193  inline G4bool IsAnalogueMonteCarlo () {return AnalogueMC;}
194  // Returns true if the simulation is an analogue Monte Carlo, and false if
195  // any of the biassing schemes have been selected.
196 
197  inline void SetBRBias (G4bool r) {
198  BRBias = r;
200  }
201  // Sets whether branching ration bias scheme applies.
202 
203  inline void SetSplitNuclei (G4int r) {
204  NSplit = r;
206  }
207  // Sets the N number for the Nuclei spliting bias scheme
208 
209  inline G4int GetSplitNuclei () {return NSplit;}
210  // Returns the N number used for the Nuclei spliting bias scheme
211 
212  inline void SetDecayDirection(const G4ThreeVector& theDir) {
213  forceDecayDirection = theDir.unit();
214  }
215 
216  inline const G4ThreeVector& GetDecayDirection() const {
217  return forceDecayDirection;
218  }
219 
220  inline void SetDecayHalfAngle(G4double halfAngle=0.*CLHEP::deg) {
221  forceDecayHalfAngle = std::min(std::max(0.*CLHEP::deg,halfAngle),180.*CLHEP::deg);
222  }
223 
224  inline G4double GetDecayHalfAngle() const {return forceDecayHalfAngle;}
225 
226  inline void SetDecayCollimation(const G4ThreeVector& theDir,
227  G4double halfAngle = 0.*CLHEP::deg) {
228  SetDecayDirection(theDir);
229  SetDecayHalfAngle(halfAngle);
230  }
231 
232  // Force direction (random within half-angle) for "visible" daughters
233  // (applies to electrons, positrons, gammas, neutrons, protons or alphas)
234 
236 
237  G4VParticleChange* DecayIt(const G4Track& theTrack,
238  const G4Step& theStep);
239 
240  // static G4ThreadLocal G4Fragment* polarizedNucleus;
241 
242  protected:
243 
244  G4DecayProducts* DoDecay(const G4ParticleDefinition& theParticleDef);
245 
246  // Apply directional bias for "visible" daughters (e+-, gamma, n, p, alpha)
247  void CollimateDecay(G4DecayProducts* products);
250 
251  G4double GetMeanFreePath(const G4Track& theTrack, G4double previousStepSize,
253 
254  G4double GetMeanLifeTime(const G4Track& theTrack,
256 
258 
260 
261  G4int GetDecayTimeBin(const G4double aDecayTime);
262 
263  //Add gamma,Xray,conversion,and auger electrons for bias mode
265  G4double weight,
266  G4double currenTime,
267  std::vector<double>& weights_v,
268  std::vector<double>& times_v,
269  std::vector<G4DynamicParticle*>& secondaries_v);
270 
271  private:
272 
274  G4RadioactiveDecay & operator=(const G4RadioactiveDecay &right);
275 
276  G4RadioactiveDecaymessenger* theRadioactiveDecaymessenger;
277 
278  G4NucleusLimits theNucleusLimits;
279 
280  G4bool isInitialised;
281  G4bool AnalogueMC;
282  G4bool BRBias;
283  G4bool FBeta;
284  G4int NSplit;
285 
286  G4double halflifethreshold;
287  G4bool applyICM;
288  G4bool applyARM;
289 
290  // Parameters for pre-collimated (biased) decay products
291  G4ThreeVector forceDecayDirection;
292  G4double forceDecayHalfAngle;
293  static const G4ThreeVector origin; // (0,0,0) for convenience
294 
295  G4int NSourceBin;
296  G4double SBin[100];
297  G4double SProfile[100];
298  G4int NDecayBin;
299  G4double DBin[100];
300  G4double DProfile[100];
301 
302  std::vector<G4String> ValidVolumes;
303  bool isAllVolumesMode;
304 
305  G4RadioactiveDecayRate theDecayRate;
306  G4RadioactiveDecayRates theDecayRateVector;
307  G4RadioactiveDecayRateVector theDecayRateTable;
308  G4RadioactiveDecayRateTable theDecayRateTableVector;
309 
310  // for the radioactivity tables
311  std::vector<G4RadioactivityTable*> theRadioactivityTables;
312  G4int decayWindows[100];
313  static const G4double levelTolerance;
314 
315  //User define radioactive decay data files replacing some files in the G4RADECAY database
316  std::map<G4int, G4String> theUserRadioactiveDataFiles;
317 
318  // Library of decay tables
319  DecayTableMap* dkmap;
320 #ifdef G4MULTITHREADED
321  static DecayTableMap* master_dkmap;
322 #endif
323 
324  // Remainder of life time at rest
325  G4double fRemainderLifeTime;
326  G4int verboseLevel;
327 
328 
329  // ParticleChange for decay process
330  G4ParticleChangeForRadDecay fParticleChangeForRadDecay;
331 
332  // inline implementations
333  inline
334  G4double AtRestGetPhysicalInteractionLength(const G4Track& track,
336  {
337  fRemainderLifeTime =
339  return fRemainderLifeTime;
340  }
341 
342  inline
343  G4VParticleChange* AtRestDoIt(const G4Track& theTrack,
344  const G4Step& theStep)
345  {return DecayIt(theTrack, theStep);}
346 
347  inline
348  G4VParticleChange* PostStepDoIt(const G4Track& theTrack,
349  const G4Step& theStep)
350  {return DecayIt(theTrack, theStep);}
351 
352 #ifdef G4MULTITHREADED
353  public:
354  static G4Mutex radioactiveDecayMutex;
355 #endif
356 };
357 
358 #endif
359 
G4double condition(const G4ErrorSymMatrix &m)
std::vector< G4RadioactiveDecayRateVector > G4RadioactiveDecayRateTable
G4DecayProducts * DoDecay(const G4ParticleDefinition &theParticleDef)
void SelectAVolume(const G4String aVolume)
std::map< G4String, G4DecayTable * > DecayTableMap
G4bool IsRateTableReady(const G4ParticleDefinition &)
G4double ConvolveSourceTimeProfile(const G4double, const G4double)
std::vector< G4RadioactiveDecayRate > G4RadioactiveDecayRates
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)
double A(double temperature)
const XML_Char int const XML_Char * value
Definition: expat.h:331
void SetDecayRate(G4int, G4int, G4double, G4int, std::vector< G4double >, std::vector< G4double >)
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
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)
static constexpr double s
G4DecayTable * GetDecayTable(const G4ParticleDefinition *)
void SetNucleusLimits(G4NucleusLimits theNucleusLimits1)
G4double GetMeanFreePath(const G4Track &theTrack, G4double previousStepSize, G4ForceCondition *condition)
G4int G4Mutex
Definition: G4Threading.hh:173
static constexpr double deg
void GetDecayRateTable(const G4ParticleDefinition &)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
void SetSplitNuclei(G4int r)
G4double GetDecayHalfAngle() const
Hep3Vector unit() 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)
G4RadioactiveDecay(const G4String &processName="RadioactiveDecay")
void SetVerboseLevel(G4int value)
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)
void AddDeexcitationSpectrumForBiasMode(G4ParticleDefinition *apartDef, G4double weight, G4double currenTime, std::vector< double > &weights_v, std::vector< double > &times_v, std::vector< G4DynamicParticle * > &secondaries_v)
void AddUserDecayDataFile(G4int Z, G4int A, G4String filename)
void CollimateDecay(G4DecayProducts *products)
void SetDecayBias(G4String filename)
G4VParticleChange * DecayIt(const G4Track &theTrack, const G4Step &theStep)