Geant4 11.3.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////////////////////////////////////////////////////////////////////////////////
27//
28// GEANT4 Class header file
29//
30// G4RadioactiveDecay
31//
32// Authors: D.H. Wright (SLAC) on base of F. Lei and P.R. Truscott code
33// Date: 9 August 2017
34//
35// 30 May 2024 V. Ivanchenko rename the class
36//
37// Description: This class is the extension of simulation of radioactive
38// decays allowing several biasing modes.
39//
40///////////////////////////////////////////////////////////////////////////////
41
42#ifndef G4RadioactiveDecay_h
43#define G4RadioactiveDecay_h 1
44
46
47#include <vector>
48#include <map>
50
51#include "G4ios.hh"
52#include "globals.hh"
55
56#include "G4NucleusLimits.hh"
60#include "G4ThreeVector.hh"
61#include "G4Threading.hh"
62
63class G4Fragment;
65
66typedef std::vector<G4RadioactiveDecayChainsFromParent> G4RadioactiveDecayParentChainTable;
67typedef std::vector<G4RadioactiveDecayRatesToDaughter> G4RadioactiveDecayRates;
68typedef std::map<G4String, G4DecayTable*> DecayTableMap;
69
71{
72 public: // with description
73
74 G4RadioactiveDecay(const G4String& processName="RadioactiveDecay",
75 const G4double timeThresholdForRadioactiveDecays=-1.0);
76 ~G4RadioactiveDecay() override;
77
78 G4VParticleChange* DecayIt(const G4Track& theTrack,
79 const G4Step& theStep) override;
80
81 void ProcessDescription(std::ostream& outFile) const override;
82
83 // Set the decay biasing scheme using the data in "filename"
84 void SetDecayBias(const G4String& filename);
85
86 // Set the half-life threshold for isomer production
87 void SetHLThreshold(G4double hl) {halflifethreshold = hl;}
88
89 void SetSourceTimeProfile(const G4String& filename);
90 // Set source exposure function using histograms in "filename"
91
93 // Returns true if the coefficient and decay time table for all the
94 // descendants of the specified isotope are ready.
95 // used in VR decay mode only
96
98 // Calculates the coefficient and decay time table for all the descendents
99 // of the specified isotope. Adds the calculated table to the private data
100 // member "theParentChainTable".
101 // used in VR decay mode only
102
104 // Used to retrieve the coefficient and decay time table for all the
105 // descendants of the specified isotope from "theParentChainTable"
106 // and place it in "chainsFromParent".
107 // used in VR decay mode only
108
109 void SetDecayRate(G4int,G4int,G4double, G4int, std::vector<G4double>&,
110 std::vector<G4double>&);
111 // Sets "theDecayRate" with data supplied in the arguements.
112 // used in VR decay mode only
113
114 std::vector<G4RadioactivityTable*>& GetTheRadioactivityTables()
115 {return theRadioactivityTables;}
116 // Return vector of G4Radioactivity map - should be used in VR mode only
117
118
119 // Controls whether G4RadioactiveDecay runs in analogue mode or
120 // variance reduction mode. SetBRBias, SetSplitNuclei and
121 // SetSourceTimeProfile all turn off analogue mode and use VR mode
123 AnalogueMC = r;
124 }
125
126 // Returns true if the simulation is an analogue Monte Carlo, and false if
127 // any of the biassing schemes have been selected.
128 inline G4bool IsAnalogueMonteCarlo () {return AnalogueMC;}
129
130 // Sets whether branching ration bias scheme applies.
131 inline void SetBRBias(G4bool r) {
132 BRBias = r;
133 AnalogueMC = false;
134 }
135
136 // Sets the number of times a nucleus will decay when biased
137 inline void SetSplitNuclei(G4int r) {
138 NSplit = r;
139 AnalogueMC = false;
140 }
141
142 // Returns the nuclear splitting number
143 inline G4int GetSplitNuclei () {return NSplit;}
144
147
148 protected:
149
152 G4int GetDecayTimeBin(const G4double aDecayTime);
153
154 G4double GetMeanLifeTime(const G4Track& theTrack,
155 G4ForceCondition* condition) override;
156
157 //Add gamma,Xray,conversion,and auger electrons for bias mode
159 G4double weight,
160 G4double currenTime,
161 std::vector<double>& weights_v,
162 std::vector<double>& times_v,
163 std::vector<G4DynamicParticle*>& secondaries_v);
164
165 private:
166
167 G4RadioactivationMessenger* theRadioactivationMessenger;
168 G4bool AnalogueMC;
169 G4bool BRBias;
170 G4int NSplit;
171
172 G4double halflifethreshold;
173
174 G4int NSourceBin;
175 G4double SBin[100];
176 G4double SProfile[100];
177 G4int NDecayBin;
178 G4double DBin[100];
179 G4double DProfile[100];
180
181 G4RadioactiveDecayRatesToDaughter ratesToDaughter;
182 G4RadioactiveDecayRates theDecayRateVector;
183 G4RadioactiveDecayChainsFromParent chainsFromParent;
184 G4RadioactiveDecayParentChainTable theParentChainTable;
185
186 // for the radioactivity tables
187 std::vector<G4RadioactivityTable*> theRadioactivityTables;
188 G4int decayWindows[100];
189};
190
191#endif
192
G4double condition(const G4ErrorSymMatrix &m)
G4ForceCondition
std::vector< G4RadioactiveDecayRatesToDaughter > G4RadioactiveDecayRates
std::map< G4String, G4DecayTable * > DecayTableMap
std::vector< G4RadioactiveDecayChainsFromParent > G4RadioactiveDecayParentChainTable
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4VParticleChange * DecayIt(const G4Track &theTrack, const G4Step &theStep) override
G4RadioactiveDecay & operator=(const G4RadioactiveDecay &right)=delete
void SetDecayBias(const G4String &filename)
void ProcessDescription(std::ostream &outFile) const override
std::vector< G4RadioactivityTable * > & GetTheRadioactivityTables()
G4int GetDecayTimeBin(const G4double aDecayTime)
void SetAnalogueMonteCarlo(G4bool r)
void CalculateChainsFromParent(const G4ParticleDefinition &)
G4RadioactiveDecay(const G4String &processName="RadioactiveDecay", const G4double timeThresholdForRadioactiveDecays=-1.0)
void SetHLThreshold(G4double hl)
G4double GetMeanLifeTime(const G4Track &theTrack, G4ForceCondition *condition) override
void SetSourceTimeProfile(const G4String &filename)
G4RadioactiveDecay(const G4RadioactiveDecay &right)=delete
G4bool IsRateTableReady(const G4ParticleDefinition &)
G4double ConvolveSourceTimeProfile(const G4double, const G4double)
void AddDeexcitationSpectrumForBiasMode(G4ParticleDefinition *apartDef, G4double weight, G4double currenTime, std::vector< double > &weights_v, std::vector< double > &times_v, std::vector< G4DynamicParticle * > &secondaries_v)
void SetDecayRate(G4int, G4int, G4double, G4int, std::vector< G4double > &, std::vector< G4double > &)
void GetChainsFromParent(const G4ParticleDefinition &)
G4VRadioactiveDecay(const G4String &processName="RadioactiveDecay", const G4double timeThreshold=-1.0)