Geant4 9.6.0
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#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// %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
55////////////////////////////////////////////////////////////////////////////////
56//
57#include <vector>
58#include <map>
60
61#include "G4ios.hh"
62#include "globals.hh"
65// #include "G4RadioactiveDecaymessenger.hh"
66
67#include "G4NucleusLimits.hh"
70#include "G4RIsotopeTable.hh"
72#include "G4ThreeVector.hh"
73
75
76typedef std::vector<G4RadioactiveDecayRateVector> G4RadioactiveDecayRateTable;
77typedef std::vector<G4RadioactiveDecayRate> G4RadioactiveDecayRates;
78
79
81{
82 // class description
83
84 // Implementation of the radioactive decay process which simulates the
85 // decays of radioactive nuclei. These nuclei are submitted to RDM as
86 // G4Ions. The required half-lives and decay schemes are retrieved from
87 // the Radioactivity database which was derived from ENSDF.
88 // All decay products are submitted back to the particle tracking process
89 // through the G4ParticleChangeForRadDecay object.
90 // class description - end
91
92 public: // with description
93
94 G4RadioactiveDecay(const G4String& processName="RadioactiveDecay");
96
98 // Returns true if the specified isotope is
99 // 1) defined as "nucleus" and
100 // 2) it is within theNucleusLimit
101
103 // Returns true if the decay table of the specified nucleus is ready
104
105 void SelectAVolume(const G4String aVolume);
106 // Select a logical volume in which RDM applies
107
108 void DeselectAVolume(const G4String aVolume);
109 // remove a logical volume from the RDM applied list
110
111 void SelectAllVolumes();
112 // Select all logical volumes for the application of RDM
113
114 void DeselectAllVolumes();
115 // Remove all logical volumes from RDM applications
116
117 void SetDecayBias (G4String filename);
118 // Sets the decay biasing scheme using the data in "filename"
119
120 void SetHLThreshold (G4double hl) {halflifethreshold = hl;}
121 // Set the half-life threshold for isomer production
122
123 void SetICM (G4bool icm) {applyICM = icm;}
124 // Enable/disable ICM
125
126 void SetARM (G4bool arm) {applyARM = arm;}
127 // Enable/disable ARM
128
129 void SetSourceTimeProfile (G4String filename) ;
130 // Sets source exposure function using histograms in "filename"
131
133 // Returns true if the coefficient and decay time table for all the
134 // descendants of the specified isotope are ready.
135 // used in VR decay mode only
136
138 // Calculates the coefficient and decay time table for all the descendents
139 // of the specified isotope. Adds the calculated table to the private data
140 // member "theDecayRateTableVector".
141 // used in VR decay mode only
142
144 // Used to retrieve the coefficient and decay time table for all the
145 // descendants of the specified isotope from "theDecayRateTableVector"
146 // and place it in "theDecayRateTable".
147 // used in VR decay mode only
148
149 void SetDecayRate(G4int,G4int,G4double, G4int, std::vector<G4double>,
150 std::vector<G4double>);
151 // Sets "theDecayRate" with data supplied in the arguements.
152 // used in VR decay mode only
153
154 std::vector<G4RadioactivityTable*> GetTheRadioactivityTables()
155 {return theRadioactivityTables;}
156 // Return vector of G4Radioactivity map - should be used in VR mode only
157
159 // Load the decay data of isotope theParentNucleus
160
161 void AddUserDecayDataFile(G4int Z, G4int A,G4String filename);
162 //Allow the user to replace the radio-active decay data provided in Geant4
163 // by its own data file for a given isotope
164
165
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;
186 if (!AnalogueMC) halflifethreshold = 1e-6*CLHEP::s;
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
194 inline G4bool IsAnalogueMonteCarlo () {return AnalogueMC;}
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) {
222 forceDecayHalfAngle = std::min(std::max(0.*CLHEP::deg,halfAngle),180.*CLHEP::deg);
223 }
224
225 inline G4double GetDecayHalfAngle() const {return forceDecayHalfAngle;}
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, or alphas)
235
237
238 protected:
239
240 G4VParticleChange* DecayIt(const G4Track& theTrack,
241 const G4Step& theStep);
242
244
245 // Apply directional bias for "visible" daughters (e+-, gamma, n, alpha)
246 void CollimateDecay(G4DecayProducts* products);
249
250 G4double GetMeanFreePath(const G4Track& theTrack, G4double previousStepSize,
252
253 G4double GetMeanLifeTime(const G4Track& theTrack,
255
257
259
260 G4int GetDecayTimeBin(const G4double aDecayTime);
261
262 private:
263
265 G4RadioactiveDecay & operator=(const G4RadioactiveDecay &right);
266
267 G4RadioactiveDecaymessenger* theRadioactiveDecaymessenger;
268 G4RIsotopeTable* theIsotopeTable;
269
270 G4NucleusLimits theNucleusLimits;
271
272 G4double HighestValue;
273
274 G4bool isInitialised;
275 G4bool AnalogueMC;
276 G4bool BRBias;
277 G4bool FBeta;
278 G4int NSplit;
279
280 G4double halflifethreshold;
281 G4bool applyICM;
282 G4bool applyARM;
283
284 // Parameters for pre-collimated (biased) decay products
285 G4ThreeVector forceDecayDirection;
286 G4double forceDecayHalfAngle;
287 static const G4ThreeVector origin; // (0,0,0) for convenience
288
289 G4int NSourceBin;
290 G4double SBin[100];
291 G4double SProfile[100];
292 G4int NDecayBin;
293 G4double DBin[100];
294 G4double DProfile[100];
295
296 std::vector<G4String> LoadedNuclei;
297 std::vector<G4String> ValidVolumes;
298 bool isAllVolumesMode;
299
300 G4RadioactiveDecayRate theDecayRate;
301 G4RadioactiveDecayRates theDecayRateVector;
302 G4RadioactiveDecayRateVector theDecayRateTable;
303 G4RadioactiveDecayRateTable theDecayRateTableVector;
304
305 // for the radioactivity tables
306 std::vector<G4RadioactivityTable*> theRadioactivityTables;
307 G4int decayWindows[99];
308 static const G4double levelTolerance;
309
310 //User define radioactive decay data files replacing some files in the G4RADECAY database
311 std::map<G4int, G4String> theUserRadioactiveDataFiles;
312
313
314 // Remainder of life time at rest
315 G4double fRemainderLifeTime;
316 G4int verboseLevel;
317
318
319 // ParticleChange for decay process
320 G4ParticleChangeForRadDecay fParticleChangeForRadDecay;
321
322 // inline implementations
323 inline
324 G4double AtRestGetPhysicalInteractionLength(const G4Track& track,
326 {
327 fRemainderLifeTime =
329 return fRemainderLifeTime;
330 }
331
332 inline
333 G4VParticleChange* AtRestDoIt(const G4Track& theTrack,
334 const G4Step& theStep)
335 {return DecayIt(theTrack, theStep);}
336
337 inline
338 G4VParticleChange* PostStepDoIt(const G4Track& theTrack,
339 const G4Step& theStep)
340 {return DecayIt(theTrack, theStep);}
341};
342
343#endif
344
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
std::vector< G4RadioactiveDecayRateVector > G4RadioactiveDecayRateTable
std::vector< G4RadioactiveDecayRate > G4RadioactiveDecayRates
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
Hep3Vector unit() const
G4DecayProducts * DoDecay(G4ParticleDefinition &theParticleDef)
void DeselectAVolume(const G4String aVolume)
void SetSplitNuclei(G4int r)
void AddDecayRateTable(const G4ParticleDefinition &)
void SetVerboseLevel(G4int value)
void SetSourceTimeProfile(G4String filename)
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetTaoTime(const G4double, const G4double)
std::vector< G4RadioactivityTable * > GetTheRadioactivityTables()
void SetDecayBias(G4String filename)
void SetICM(G4bool icm)
G4double GetDecayHalfAngle() const
void SetBRBias(G4bool r)
G4double GetMeanLifeTime(const G4Track &theTrack, G4ForceCondition *condition)
G4NucleusLimits GetNucleusLimits() const
void CollimateDecayProduct(G4DynamicParticle *product)
const G4ThreeVector & GetDecayDirection() const
void SelectAVolume(const G4String aVolume)
G4int GetVerboseLevel() const
G4double GetMeanFreePath(const G4Track &theTrack, G4double previousStepSize, G4ForceCondition *condition)
G4ThreeVector ChooseCollimationDirection() const
G4int GetDecayTimeBin(const G4double aDecayTime)
G4VParticleChange * DecayIt(const G4Track &theTrack, const G4Step &theStep)
void SetDecayHalfAngle(G4double halfAngle=0.*CLHEP::deg)
void SetAnalogueMonteCarlo(G4bool r)
void SetDecayRate(G4int, G4int, G4double, G4int, std::vector< G4double >, std::vector< G4double >)
void SetDecayDirection(const G4ThreeVector &theDir)
G4bool IsApplicable(const G4ParticleDefinition &)
void SetARM(G4bool arm)
void SetHLThreshold(G4double hl)
void GetDecayRateTable(const G4ParticleDefinition &)
void SetDecayCollimation(const G4ThreeVector &theDir, G4double halfAngle=0.*CLHEP::deg)
G4bool IsLoaded(const G4ParticleDefinition &)
G4bool IsRateTableReady(const G4ParticleDefinition &)
void SetNucleusLimits(G4NucleusLimits theNucleusLimits1)
void AddUserDecayDataFile(G4int Z, G4int A, G4String filename)
G4DecayTable * LoadDecayTable(G4ParticleDefinition &theParentNucleus)
void CollimateDecay(G4DecayProducts *products)
Definition: G4Step.hh:78
virtual G4double AtRestGetPhysicalInteractionLength(const G4Track &, G4ForceCondition *)