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

Last change on this file since 1347 was 1347, checked in by garnier, 13 years ago

geant4 tag 9.4

File size: 10.5 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 "G4RadioactivityTable.hh"
64
65#include <vector>
66
67
68//class G4UserlimitsForRD;
69
70////////////////////////////////////////////////////////////////////////////////
71//
72class G4RadioactiveDecaymessenger;
73
74typedef std::vector<G4RadioactiveDecayRateVector> G4RadioactiveDecayRateTable;
75typedef std::vector<G4RadioactiveDecayRate> G4RadioactiveDecayRates;
76
77class G4RadioactiveDecay : public G4VRestDiscreteProcess
78{
79        // class description
80
81        // Implementation of nuclear radioactive decay process.
82        // It simulates the decays of radioactive nuclei, which is submitted to RDM in the form
83        // of G4Ions. Half-liffe and decay schemes are hold in the Radioactivity database
84        // All decay products are submitted back to the particle tracking process through
85        // the G4ParticleChangeForRadDecay object.
86
87        // class description - end
88
89public: // with description
90
91        G4RadioactiveDecay (const G4String& processName="RadioactiveDecay");
92        ~G4RadioactiveDecay();
93        // constructor and destructor
94        //
95
96        G4bool IsApplicable(const G4ParticleDefinition &);
97        // Returns true if the specified isotope is
98        //  1) defined as "nucleus" and
99        //  2) it is within theNucleusLimit
100        //
101        G4bool IsLoaded (const G4ParticleDefinition &);
102        // Returns true if the decay table of the specified nuclei is ready.
103        //
104        void SelectAVolume(const G4String aVolume);
105        // Select a logical volume in which RDM applies.
106        //
107        void DeselectAVolume(const G4String aVolume);
108        // remove a logical volume from the RDM applied list
109        //
110        void SelectAllVolumes();
111        // Select all logical volumes for the application of RDM.
112        //
113        void DeselectAllVolumes();
114        // Remove all logical volumes from RDM applications.
115        //
116        void SetDecayBias (G4String filename);
117        //   Sets the decay biasing scheme using the data in "filename"
118        //
119        void SetHLThreshold (G4double hl) {halflifethreshold = hl;}
120        // Set the half-life threshold for isomer production
121        //
122        void SetICM (G4bool icm) {applyICM = icm;}
123        // Enable/disable ICM
124        //
125        void SetARM (G4bool arm) {applyARM = arm;}
126        // Enable/disable ARM
127        //
128        void SetSourceTimeProfile (G4String filename) ;
129        //  Sets source exposure function using histograms in "filename"
130        //
131        G4bool IsRateTableReady(const G4ParticleDefinition &);
132        // Returns true if the coefficient and decay time table for all the
133        // descendants of the specified isotope are ready.
134        //
135        // used in VR decay mode only
136        //
137        void AddDecayRateTable(const G4ParticleDefinition &);
138        // Calculates the coefficient and decay time table for all the descendants
139        // of the specified isotope.  Adds the calculated table to the private data
140        // member "theDecayRateTableVector".
141        //
142        //
143        // used in VR decay mode only
144        //
145        void GetDecayRateTable(const G4ParticleDefinition &);
146        // Used to retrieve the coefficient and decay time table for all the
147        // descendants of the specified isotope from "theDecayRateTableVector"
148        // and place it in "theDecayRateTable".
149        //
150        //
151        // used in VR decay mode only
152        //
153        void SetDecayRate(G4int,G4int,G4double, G4int, std::vector<G4double>,
154                std::vector<G4double>);
155        // Sets "theDecayRate" with data supplied in the arguements.
156        //
157        //  //
158        // used in VR decay mode only
159        //
160
161        std::vector<G4RadioactivityTable*> GetTheRadioactivityTables() {return theRadioactivityTables;}
162        // return the vector of G4Radioactivity map
163        // should be used in VR mode only
164
165        G4DecayTable *LoadDecayTable (G4ParticleDefinition & theParentNucleus);
166        // Load the decay data of isotope theParentNucleus.
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        //
179        inline G4NucleusLimits GetNucleusLimits() const
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*s;
187        }
188        //   Controls whether G4RadioactiveDecay runs in analogue mode or
189        //   variance reduction mode.
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 ) { BRBias  = r;
198        SetAnalogueMonteCarlo(0);}
199        //   Sets whether branching ration bias scheme applies.
200        //
201        inline void SetSplitNuclei (G4int r) { NSplit = r; SetAnalogueMonteCarlo(0); }
202        //  Sets the N number for the Nuclei spliting bias scheme
203        //
204        inline G4int GetSplitNuclei () {return NSplit;}
205        //  Returns the N number used for the Nuclei spliting bias scheme
206        //
207public:
208
209        void BuildPhysicsTable(const G4ParticleDefinition &);
210
211protected:
212
213        G4VParticleChange* DecayIt( const G4Track& theTrack,
214                const G4Step&  theStep);
215
216        G4DecayProducts* DoDecay(G4ParticleDefinition& theParticleDef);
217
218        G4double GetMeanFreePath(const G4Track& theTrack,
219                G4double previousStepSize,
220                G4ForceCondition* condition);
221
222        G4double GetMeanLifeTime(const G4Track& theTrack,
223                G4ForceCondition* condition);
224
225        G4double GetTaoTime(G4double,G4double);
226
227        G4double GetDecayTime();
228
229        G4int GetDecayTimeBin(const G4double aDecayTime);
230
231private:
232
233        G4RadioactiveDecay(const G4RadioactiveDecay &right);
234        G4RadioactiveDecay & operator=(const G4RadioactiveDecay &right);
235
236private:
237
238        G4RadioactiveDecaymessenger  *theRadioactiveDecaymessenger;
239
240        G4PhysicsTable               *aPhysicsTable;
241
242        G4RIsotopeTable              *theIsotopeTable;
243
244        G4NucleusLimits               theNucleusLimits;
245
246        const G4double                HighestBinValue;
247        const G4double                LowestBinValue;
248
249        const G4int                   TotBin;
250
251        G4bool                        AnalogueMC;
252        G4bool                        BRBias;
253        G4bool                        FBeta;
254        G4int                         NSplit;
255
256        G4double                      halflifethreshold;
257        G4bool                        applyICM;
258        G4bool                        applyARM;
259
260        G4int                         NSourceBin;
261        G4double                      SBin[99];
262        G4double                      SProfile[99];
263
264        G4int                         NDecayBin;
265        G4double                      DBin[99];
266        G4double                      DProfile[99];
267
268        std::vector<G4String>         LoadedNuclei;
269        std::vector<G4String>         ValidVolumes;
270
271        G4RadioactiveDecayRate        theDecayRate;
272        G4RadioactiveDecayRates       theDecayRateVector;
273        G4RadioactiveDecayRateVector  theDecayRateTable;
274        G4RadioactiveDecayRateTable   theDecayRateTableVector;
275
276        // for the radioactivity tables
277        std::vector<G4RadioactivityTable*>      theRadioactivityTables;
278        G4int                         decayWindows[99];
279
280        //
281        static const G4double             levelTolerance;
282
283        // Remainder of life time at rest
284        G4double                      fRemainderLifeTime;
285
286        G4int                         verboseLevel;
287
288
289        // ParticleChange for decay process
290        G4ParticleChangeForRadDecay   fParticleChangeForRadDecay;
291
292        inline G4double AtRestGetPhysicalInteractionLength
293                (const G4Track& track, G4ForceCondition* condition)
294        {fRemainderLifeTime = G4VRestDiscreteProcess::
295        AtRestGetPhysicalInteractionLength(track, condition );
296        return fRemainderLifeTime;}
297
298        inline G4VParticleChange* AtRestDoIt
299                (const G4Track& theTrack, const G4Step& theStep)
300        {return DecayIt(theTrack, theStep);}
301
302        inline G4VParticleChange* PostStepDoIt
303                (const G4Track& theTrack, const G4Step& theStep)
304        {return DecayIt(theTrack, theStep);}
305
306};
307#endif
308
309
310
311
312
313
314
Note: See TracBrowser for help on using the repository browser.