Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Radioactivation.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 G4Radioactivation_h
27#define G4Radioactivation_h 1
28
29#include "G4RadioactiveDecay.hh"
30
31#include <vector>
32#include <map>
34
35#include "G4ios.hh"
36#include "globals.hh"
39
40#include "G4NucleusLimits.hh"
44#include "G4ThreeVector.hh"
45#include "G4Threading.hh"
46
47class G4Fragment;
49
50typedef std::vector<G4RadioactiveDecayChainsFromParent> G4RadioactiveDecayParentChainTable;
51typedef std::vector<G4RadioactiveDecayRatesToDaughter> G4RadioactiveDecayRates;
52typedef std::map<G4String, G4DecayTable*> DecayTableMap;
53
55{
56 public: // with description
57
58 G4Radioactivation(const G4String& processName="Radioactivation",
59 const G4double timeThresholdForRadioactiveDecays=-1.0);
60 ~G4Radioactivation() override;
61
62 G4VParticleChange* DecayIt(const G4Track& theTrack,
63 const G4Step& theStep) override;
64
65 void ProcessDescription(std::ostream& outFile) const override;
66
67 // Set the decay biasing scheme using the data in "filename"
68 void SetDecayBias(const G4String& filename);
69
70 // Set the half-life threshold for isomer production
71 void SetHLThreshold(G4double hl) {halflifethreshold = hl;}
72
73 void SetSourceTimeProfile(const G4String& filename);
74 // Set source exposure function using histograms in "filename"
75
77 // Returns true if the coefficient and decay time table for all the
78 // descendants of the specified isotope are ready.
79 // used in VR decay mode only
80
82 // Calculates the coefficient and decay time table for all the descendents
83 // of the specified isotope. Adds the calculated table to the private data
84 // member "theParentChainTable".
85 // used in VR decay mode only
86
88 // Used to retrieve the coefficient and decay time table for all the
89 // descendants of the specified isotope from "theParentChainTable"
90 // and place it in "chainsFromParent".
91 // used in VR decay mode only
92
93 void SetDecayRate(G4int,G4int,G4double, G4int, std::vector<G4double>&,
94 std::vector<G4double>&);
95 // Sets "theDecayRate" with data supplied in the arguements.
96 // used in VR decay mode only
97
98 std::vector<G4RadioactivityTable*>& GetTheRadioactivityTables()
99 {return theRadioactivityTables;}
100 // Return vector of G4Radioactivity map - should be used in VR mode only
101
102
103 // Controls whether G4Radioactivation runs in analogue mode or
104 // variance reduction mode. SetBRBias, SetSplitNuclei and
105 // SetSourceTimeProfile all turn off analogue mode and use VR mode
107 AnalogueMC = r;
108 }
109
110 // Returns true if the simulation is an analogue Monte Carlo, and false if
111 // any of the biassing schemes have been selected.
112 inline G4bool IsAnalogueMonteCarlo () {return AnalogueMC;}
113
114 // Sets whether branching ration bias scheme applies.
115 inline void SetBRBias(G4bool r) {
116 BRBias = r;
117 AnalogueMC = false;
118 }
119
120 // Sets the number of times a nucleus will decay when biased
121 inline void SetSplitNuclei(G4int r) {
122 NSplit = r;
123 AnalogueMC = false;
124 }
125
126 // Returns the nuclear splitting number
127 inline G4int GetSplitNuclei () {return NSplit;}
128
129 G4Radioactivation(const G4Radioactivation& right) = delete;
131
132 protected:
133
136 G4int GetDecayTimeBin(const G4double aDecayTime);
137
138 G4double GetMeanLifeTime(const G4Track& theTrack,
139 G4ForceCondition* condition) override;
140
141 //Add gamma,Xray,conversion,and auger electrons for bias mode
143 G4double weight,
144 G4double currenTime,
145 std::vector<double>& weights_v,
146 std::vector<double>& times_v,
147 std::vector<G4DynamicParticle*>& secondaries_v);
148
149 private:
150
151 G4RadioactivationMessenger* theRadioactivationMessenger;
152 G4bool AnalogueMC;
153 G4bool BRBias;
154 G4int NSplit;
155
156 G4double halflifethreshold;
157
158 G4int NSourceBin;
159 G4double SBin[100];
160 G4double SProfile[100];
161 G4int NDecayBin;
162 G4double DBin[100];
163 G4double DProfile[100];
164
165 G4RadioactiveDecayRatesToDaughter ratesToDaughter;
166 G4RadioactiveDecayRates theDecayRateVector;
167 G4RadioactiveDecayChainsFromParent chainsFromParent;
168 G4RadioactiveDecayParentChainTable theParentChainTable;
169
170 // for the radioactivity tables
171 std::vector<G4RadioactivityTable*> theRadioactivityTables;
172 G4int decayWindows[100];
173};
174
175#endif
176
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
std::vector< G4RadioactivityTable * > & GetTheRadioactivityTables()
void SetAnalogueMonteCarlo(G4bool r)
G4int GetDecayTimeBin(const G4double aDecayTime)
void AddDeexcitationSpectrumForBiasMode(G4ParticleDefinition *apartDef, G4double weight, G4double currenTime, std::vector< double > &weights_v, std::vector< double > &times_v, std::vector< G4DynamicParticle * > &secondaries_v)
G4bool IsRateTableReady(const G4ParticleDefinition &)
void SetDecayBias(const G4String &filename)
G4double GetMeanLifeTime(const G4Track &theTrack, G4ForceCondition *condition) override
G4double ConvolveSourceTimeProfile(const G4double, const G4double)
void CalculateChainsFromParent(const G4ParticleDefinition &)
G4VParticleChange * DecayIt(const G4Track &theTrack, const G4Step &theStep) override
G4Radioactivation(const G4String &processName="Radioactivation", const G4double timeThresholdForRadioactiveDecays=-1.0)
void SetSplitNuclei(G4int r)
void SetDecayRate(G4int, G4int, G4double, G4int, std::vector< G4double > &, std::vector< G4double > &)
void ProcessDescription(std::ostream &outFile) const override
void SetSourceTimeProfile(const G4String &filename)
G4Radioactivation(const G4Radioactivation &right)=delete
void GetChainsFromParent(const G4ParticleDefinition &)
void SetHLThreshold(G4double hl)
void SetBRBias(G4bool r)
G4Radioactivation & operator=(const G4Radioactivation &right)=delete