source: trunk/source/processes/hadronic/models/radioactive_decay/include/G4RadioactiveDecay.hh @ 1340

Last change on this file since 1340 was 1315, checked in by garnier, 14 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

File size: 10.4 KB
Line 
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//
44// 29 February 2000, P R Truscott, DERA UK
45// 0.b.3 release.
46//
47// 13 April 2000, F Lei, DERA UK
48// 0.b.4 release. No change to this file     
49//
50// %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
51////////////////////////////////////////////////////////////////////////////////
52//
53#include "G4ios.hh"
54#include "globals.hh"
55#include "G4VRestDiscreteProcess.hh"
56#include "G4ParticleChangeForRadDecay.hh"
57#include "G4RadioactiveDecaymessenger.hh" 
58
59#include "G4NucleusLimits.hh"
60#include "G4RadioactiveDecayRate.hh"
61#include "G4RadioactiveDecayRateVector.hh"
62#include "G4RIsotopeTable.hh"
63#include <vector>
64
65//class G4UserlimitsForRD;
66
67////////////////////////////////////////////////////////////////////////////////
68//
69class G4RadioactiveDecaymessenger;
70
71typedef std::vector<G4RadioactiveDecayRateVector> G4RadioactiveDecayRateTable;
72typedef std::vector<G4RadioactiveDecayRate> G4RadioactiveDecayRates;
73
74class G4RadioactiveDecay : public G4VRestDiscreteProcess
75{
76  // class description
77 
78  // Implementation of nuclear radioactive decay process.
79  // It simulates the decays of radioactive nuclei, which is submitted to RDM in the form
80  // of G4Ions. Half-liffe and decay schemes are hold in the Radioactivity database
81  // All decay products are submitted back to the particle tracking process through
82  // the G4ParticleChangeForRadDecay object.
83
84  // class description - end
85 
86public: // with description
87
88  G4RadioactiveDecay (const G4String& processName="RadioactiveDecay");
89  ~G4RadioactiveDecay();
90  // constructor and destructor
91  //
92
93  G4bool IsApplicable(const G4ParticleDefinition &);
94  // Returns true if the specified isotope is
95  //  1) defined as "nucleus" and
96  //  2) it is within theNucleusLimit
97  //
98  G4bool IsLoaded (const G4ParticleDefinition &);
99  // Returns true if the decay table of the specified nuclei is ready.
100  //
101  void SelectAVolume(const G4String aVolume);
102  // Select a logical volume in which RDM applies.
103  //
104  void DeselectAVolume(const G4String aVolume);
105  // remove a logical volume from the RDM applied list
106  //
107  void SelectAllVolumes();
108  // Select all logical volumes for the application of RDM.
109  //
110  void DeselectAllVolumes();
111  // Remove all logical volumes from RDM applications.
112  //
113  void SetDecayBias (G4String filename);
114  //   Sets the decay biasing scheme using the data in "filename"
115  //
116  void SetHLThreshold (G4double hl) {halflifethreshold = hl;}
117  // Set the half-life threshold for isomer production
118  //
119  void SetICM (G4bool icm) {applyICM = icm;}
120  // Enable/disable ICM
121  //
122  void SetARM (G4bool arm) {applyARM = arm;}
123  // Enable/disable ARM
124  //
125  void SetSourceTimeProfile (G4String filename) ;
126  //  Sets source exposure function using histograms in "filename"
127  //
128  G4bool IsRateTableReady(const G4ParticleDefinition &);
129  // Returns true if the coefficient and decay time table for all the
130  // descendants of the specified isotope are ready.
131  //
132  // used in VR decay mode only
133  //
134  void AddDecayRateTable(const G4ParticleDefinition &);
135  // Calculates the coefficient and decay time table for all the descendants
136  // of the specified isotope.  Adds the calculated table to the private data
137  // member "theDecayRateTableVector".
138  //
139  //
140  // used in VR decay mode only
141  //
142  void GetDecayRateTable(const G4ParticleDefinition &);
143  // Used to retrieve the coefficient and decay time table for all the
144  // descendants of the specified isotope from "theDecayRateTableVector"
145  // and place it in "theDecayRateTable".
146  //
147  //
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  //
154  //  //
155  // used in VR decay mode only
156  //
157
158  G4DecayTable *LoadDecayTable (G4ParticleDefinition & theParentNucleus);
159  // Load the decay data of isotope theParentNucleus.
160  //
161  inline void  SetVerboseLevel(G4int value) {verboseLevel = value;}
162  // Sets the VerboseLevel which controls duggering display.
163  //
164  inline G4int GetVerboseLevel() const {return verboseLevel;}
165  // Returns the VerboseLevel which controls level of debugging output.
166  //
167  inline void SetNucleusLimits(G4NucleusLimits theNucleusLimits1)
168  {theNucleusLimits = theNucleusLimits1 ;}
169  //  Sets theNucleusLimits which specifies the range of isotopes
170  //  the G4RadioactiveDecay applies.
171  //
172  inline G4NucleusLimits GetNucleusLimits() const
173  {return theNucleusLimits;}
174  //  Returns theNucleusLimits which specifies the range of isotopes
175  //  the G4RadioactiveDecay applies.
176  //
177  inline void SetAnalogueMonteCarlo (G4bool r ) { AnalogueMC  = r; }
178  //   Controls whether G4RadioactiveDecay runs in analogue mode or
179  //   variance reduction mode.
180  inline void SetFBeta (G4bool r ) { FBeta  = r; }
181  //   Controls whether G4RadioactiveDecay uses fast beta simulation mode
182  //
183  inline G4bool IsAnalogueMonteCarlo () {return AnalogueMC;}
184  //   Returns true if the simulation is an analogue Monte Carlo, and false if
185  //   any of the biassing schemes have been selected.
186  //
187  inline void SetBRBias (G4bool r ) { BRBias  = r;
188  SetAnalogueMonteCarlo(0);}
189  //   Sets whether branching ration bias scheme applies.
190  //
191  inline void SetSplitNuclei (G4int r) { NSplit = r; SetAnalogueMonteCarlo(0); }
192  //  Sets the N number for the Nuclei spliting bias scheme
193  //
194  inline G4int GetSplitNuclei () {return NSplit;}
195  //  Returns the N number used for the Nuclei spliting bias scheme
196  //
197public:
198 
199  void BuildPhysicsTable(const G4ParticleDefinition &);
200
201protected:
202
203  G4VParticleChange* DecayIt( const G4Track& theTrack,
204                              const G4Step&  theStep);
205
206  G4DecayProducts* DoDecay(G4ParticleDefinition& theParticleDef);
207
208  G4double GetMeanFreePath(const G4Track& theTrack,
209                           G4double previousStepSize,
210                           G4ForceCondition* condition);
211
212  G4double GetMeanLifeTime(const G4Track& theTrack,
213                           G4ForceCondition* condition);
214
215  G4double GetTaoTime(G4double,G4double);
216
217  G4double GetDecayTime();
218
219  G4int GetDecayTimeBin(const G4double aDecayTime);
220
221private:
222
223  G4RadioactiveDecay(const G4RadioactiveDecay &right);
224  G4RadioactiveDecay & operator=(const G4RadioactiveDecay &right);
225
226private:
227
228  G4RadioactiveDecaymessenger  *theRadioactiveDecaymessenger;
229
230  G4PhysicsTable               *aPhysicsTable;
231
232  G4RIsotopeTable              *theIsotopeTable;
233
234  G4NucleusLimits               theNucleusLimits;
235
236  const G4double                HighestBinValue;
237  const G4double                LowestBinValue;
238
239  const G4int                   TotBin;
240
241  G4bool                        AnalogueMC;
242  G4bool                        BRBias;
243  G4bool                        FBeta;
244  G4int                         NSplit;
245
246  G4double                      halflifethreshold;
247  G4bool                        applyICM;
248  G4bool                        applyARM;
249
250  G4int                         NSourceBin;
251  G4double                      SBin[99];
252  G4double                      SProfile[99];
253
254  G4int                         NDecayBin;
255  G4double                      DBin[99];
256  G4double                      DProfile[99];
257
258  std::vector<G4String>         LoadedNuclei;
259  std::vector<G4String>         ValidVolumes;
260
261  G4RadioactiveDecayRate        theDecayRate;
262  G4RadioactiveDecayRates       theDecayRateVector;
263  G4RadioactiveDecayRateVector  theDecayRateTable;
264  G4RadioactiveDecayRateTable   theDecayRateTableVector;
265
266  static const G4double         levelTolerance;
267
268  // Remainder of life time at rest
269  G4double                      fRemainderLifeTime;
270
271  G4int                         verboseLevel;
272
273
274  // ParticleChange for decay process
275  G4ParticleChangeForRadDecay   fParticleChangeForRadDecay;
276
277  inline G4double AtRestGetPhysicalInteractionLength
278    (const G4Track& track, G4ForceCondition* condition)
279  {fRemainderLifeTime = G4VRestDiscreteProcess::
280    AtRestGetPhysicalInteractionLength(track, condition );
281  return fRemainderLifeTime;}
282
283  inline G4VParticleChange* AtRestDoIt
284    (const G4Track& theTrack, const G4Step& theStep)
285  {return DecayIt(theTrack, theStep);}
286
287  inline G4VParticleChange* PostStepDoIt
288    (const G4Track& theTrack, const G4Step& theStep)
289  {return DecayIt(theTrack, theStep);}
290
291};
292#endif
293
294
295
296
297
298
299
Note: See TracBrowser for help on using the repository browser.