G4RadioactiveDecay.hh

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 #ifndef G4RadioactiveDecay_h
00027 #define G4RadioactiveDecay_h 1
00028 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00029 //
00030 // MODULE:              G4RadioactiveDecay.hh
00031 //
00032 // Version:             0.b.4
00033 // Date:                14/04/00
00034 // Author:              F Lei & P R Truscott
00035 // Organisation:        DERA UK
00036 // Customer:            ESA/ESTEC, NOORDWIJK
00037 // Contract:            12115/96/JG/NL Work Order No. 3
00038 //
00039 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00040 //
00041 // CHANGE HISTORY
00042 // --------------
00043 // 17 October 2011, L Desorgher - Add the method AddUserDecayDataFile
00044 //
00045 // 01 June 2011, M. Kelsey -- Add directional biasing interface to allow for
00046 //              "collimation" of decay daughters.
00047 //
00048 // 29 February 2000, P R Truscott, DERA UK
00049 // 0.b.3 release.
00050 //
00051 // 13 April 2000, F Lei, DERA UK
00052 // 0.b.4 release. No change to this file     
00053 //
00054 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
00056 //
00057 #include <vector>
00058 #include <map>
00059 #include <CLHEP/Units/SystemOfUnits.h>
00060 
00061 #include "G4ios.hh"
00062 #include "globals.hh"
00063 #include "G4VRestDiscreteProcess.hh"
00064 #include "G4ParticleChangeForRadDecay.hh"
00065 // #include "G4RadioactiveDecaymessenger.hh"  
00066 
00067 #include "G4NucleusLimits.hh"
00068 #include "G4RadioactiveDecayRate.hh"
00069 #include "G4RadioactiveDecayRateVector.hh"
00070 #include "G4RIsotopeTable.hh"
00071 #include "G4RadioactivityTable.hh"
00072 #include "G4ThreeVector.hh"
00073 
00074 class G4RadioactiveDecaymessenger;
00075 
00076 typedef std::vector<G4RadioactiveDecayRateVector> G4RadioactiveDecayRateTable;
00077 typedef std::vector<G4RadioactiveDecayRate> G4RadioactiveDecayRates;
00078 
00079 
00080 class G4RadioactiveDecay : public G4VRestDiscreteProcess 
00081 {
00082   // class description
00083 
00084   // Implementation of the radioactive decay process which simulates the
00085   // decays of radioactive nuclei.  These nuclei are submitted to RDM as
00086   // G4Ions.  The required half-lives and decay schemes are retrieved from
00087   // the Radioactivity database which was derived from ENSDF.
00088   // All decay products are submitted back to the particle tracking process
00089   // through the G4ParticleChangeForRadDecay object.
00090   // class description - end 
00091 
00092   public: // with description
00093 
00094     G4RadioactiveDecay(const G4String& processName="RadioactiveDecay");
00095     ~G4RadioactiveDecay();
00096 
00097     G4bool IsApplicable(const G4ParticleDefinition&);
00098     // Returns true if the specified isotope is
00099     //  1) defined as "nucleus" and
00100     //  2) it is within theNucleusLimit
00101 
00102     G4bool IsLoaded(const G4ParticleDefinition &);
00103     // Returns true if the decay table of the specified nucleus is ready
00104 
00105     void SelectAVolume(const G4String aVolume);
00106     // Select a logical volume in which RDM applies
00107 
00108     void DeselectAVolume(const G4String aVolume);
00109     // remove a logical volume from the RDM applied list
00110 
00111     void SelectAllVolumes();
00112     // Select all logical volumes for the application of RDM
00113 
00114     void DeselectAllVolumes();
00115     // Remove all logical volumes from RDM applications
00116 
00117     void SetDecayBias (G4String filename);
00118     // Sets the decay biasing scheme using the data in "filename"
00119 
00120     void SetHLThreshold (G4double hl) {halflifethreshold = hl;}
00121     // Set the half-life threshold for isomer production
00122 
00123     void SetICM (G4bool icm) {applyICM = icm;}
00124     // Enable/disable ICM 
00125 
00126     void SetARM (G4bool arm) {applyARM = arm;}
00127     // Enable/disable ARM
00128 
00129     void SetSourceTimeProfile (G4String filename) ;
00130     // Sets source exposure function using histograms in "filename"
00131 
00132     G4bool IsRateTableReady(const G4ParticleDefinition &);
00133     // Returns true if the coefficient and decay time table for all the
00134     // descendants of the specified isotope are ready.
00135     // used in VR decay mode only
00136 
00137     void AddDecayRateTable(const G4ParticleDefinition&);
00138     // Calculates the coefficient and decay time table for all the descendents
00139     // of the specified isotope.  Adds the calculated table to the private data
00140     // member "theDecayRateTableVector".
00141     // used in VR decay mode only 
00142 
00143     void GetDecayRateTable(const G4ParticleDefinition&);
00144     // Used to retrieve the coefficient and decay time table for all the
00145     // descendants of the specified isotope from "theDecayRateTableVector"
00146     // and place it in "theDecayRateTable".
00147     // used in VR decay mode only 
00148 
00149     void SetDecayRate(G4int,G4int,G4double, G4int, std::vector<G4double>,
00150                       std::vector<G4double>);
00151     // Sets "theDecayRate" with data supplied in the arguements.
00152     // used in VR decay mode only 
00153 
00154     std::vector<G4RadioactivityTable*> GetTheRadioactivityTables()
00155        {return theRadioactivityTables;}
00156     // Return vector of G4Radioactivity map - should be used in VR mode only
00157 
00158     G4DecayTable *LoadDecayTable (G4ParticleDefinition & theParentNucleus);
00159     // Load the decay data of isotope theParentNucleus
00160 
00161     void AddUserDecayDataFile(G4int Z, G4int A,G4String filename);
00162     //Allow the user to replace the radio-active decay data provided in Geant4
00163     // by its own data file for a given isotope
00164 
00165 
00166 
00167 
00168     inline void  SetVerboseLevel(G4int value) {verboseLevel = value;}
00169     // Sets the VerboseLevel which controls duggering display
00170 
00171     inline G4int GetVerboseLevel() const {return verboseLevel;}
00172     // Returns the VerboseLevel which controls level of debugging output
00173 
00174     inline void SetNucleusLimits(G4NucleusLimits theNucleusLimits1)
00175       {theNucleusLimits = theNucleusLimits1 ;}
00176     // Sets theNucleusLimits which specifies the range of isotopes
00177     // the G4RadioactiveDecay applies.
00178 
00179     inline G4NucleusLimits GetNucleusLimits() const
00180       {return theNucleusLimits;}
00181     // Returns theNucleusLimits which specifies the range of isotopes
00182     // the G4RadioactiveDecay applies
00183 
00184     inline void SetAnalogueMonteCarlo (G4bool r ) { 
00185       AnalogueMC  = r; 
00186       if (!AnalogueMC) halflifethreshold = 1e-6*CLHEP::s;
00187     }
00188     // Controls whether G4RadioactiveDecay runs in analogue mode or
00189     // variance reduction mode.
00190 
00191     inline void SetFBeta (G4bool r ) { FBeta  = r; }
00192     // Controls whether G4RadioactiveDecay uses fast beta simulation mode
00193 
00194     inline G4bool IsAnalogueMonteCarlo () {return AnalogueMC;}
00195     // Returns true if the simulation is an analogue Monte Carlo, and false if
00196     // any of the biassing schemes have been selected.
00197 
00198     inline void SetBRBias (G4bool r) {
00199       BRBias = r;
00200       SetAnalogueMonteCarlo(0);
00201      }
00202      // Sets whether branching ration bias scheme applies.
00203 
00204     inline void SetSplitNuclei (G4int r) {
00205       NSplit = r;
00206       SetAnalogueMonteCarlo(0);
00207     }
00208     // Sets the N number for the Nuclei spliting bias scheme
00209 
00210     inline G4int GetSplitNuclei () {return NSplit;}
00211     //  Returns the N number used for the Nuclei spliting bias scheme
00212 
00213     inline void SetDecayDirection(const G4ThreeVector& theDir) {
00214       forceDecayDirection = theDir.unit();
00215     }
00216 
00217     inline const G4ThreeVector& GetDecayDirection() const {
00218       return forceDecayDirection; 
00219     }
00220 
00221     inline void SetDecayHalfAngle(G4double halfAngle=0.*CLHEP::deg) {
00222       forceDecayHalfAngle = std::min(std::max(0.*CLHEP::deg,halfAngle),180.*CLHEP::deg);
00223     }
00224 
00225     inline G4double GetDecayHalfAngle() const {return forceDecayHalfAngle;}
00226 
00227     inline void SetDecayCollimation(const G4ThreeVector& theDir,
00228                                     G4double halfAngle = 0.*CLHEP::deg) {
00229       SetDecayDirection(theDir);
00230       SetDecayHalfAngle(halfAngle);
00231     }
00232 
00233     // Force direction (random within half-angle) for "visible" daughters
00234     // (applies to electrons, positrons, gammas, neutrons, or alphas)
00235 
00236     void BuildPhysicsTable(const G4ParticleDefinition &);
00237 
00238   protected:
00239 
00240     G4VParticleChange* DecayIt(const G4Track& theTrack,
00241                                const G4Step&  theStep);
00242 
00243     G4DecayProducts* DoDecay(G4ParticleDefinition& theParticleDef);
00244 
00245     // Apply directional bias for "visible" daughters (e+-, gamma, n, alpha)
00246     void CollimateDecay(G4DecayProducts* products);
00247     void CollimateDecayProduct(G4DynamicParticle* product);
00248     G4ThreeVector ChooseCollimationDirection() const;
00249 
00250     G4double GetMeanFreePath(const G4Track& theTrack, G4double previousStepSize,
00251                              G4ForceCondition* condition);
00252 
00253     G4double GetMeanLifeTime(const G4Track& theTrack,
00254                              G4ForceCondition* condition);
00255 
00256     G4double GetTaoTime(const G4double,const G4double);
00257 
00258     G4double GetDecayTime();
00259 
00260     G4int GetDecayTimeBin(const G4double aDecayTime);
00261 
00262   private:
00263 
00264     G4RadioactiveDecay(const G4RadioactiveDecay &right);
00265     G4RadioactiveDecay & operator=(const G4RadioactiveDecay &right);
00266 
00267     G4RadioactiveDecaymessenger* theRadioactiveDecaymessenger;
00268     G4RIsotopeTable* theIsotopeTable;
00269 
00270     G4NucleusLimits theNucleusLimits;
00271 
00272     G4double HighestValue;
00273 
00274     G4bool isInitialised;
00275     G4bool AnalogueMC;
00276     G4bool BRBias;
00277     G4bool FBeta;
00278     G4int NSplit;
00279 
00280     G4double halflifethreshold;
00281     G4bool applyICM;
00282     G4bool applyARM;
00283 
00284     // Parameters for pre-collimated (biased) decay products
00285     G4ThreeVector forceDecayDirection;
00286     G4double      forceDecayHalfAngle;
00287     static const G4ThreeVector origin;  // (0,0,0) for convenience
00288 
00289     G4int NSourceBin;
00290     G4double SBin[100];
00291     G4double SProfile[100];
00292     G4int NDecayBin;
00293     G4double DBin[100];
00294     G4double DProfile[100];
00295 
00296     std::vector<G4String> LoadedNuclei;
00297     std::vector<G4String> ValidVolumes;
00298     bool isAllVolumesMode;
00299 
00300     G4RadioactiveDecayRate theDecayRate;
00301     G4RadioactiveDecayRates theDecayRateVector;
00302     G4RadioactiveDecayRateVector theDecayRateTable;
00303     G4RadioactiveDecayRateTable theDecayRateTableVector;
00304 
00305     // for the radioactivity tables
00306     std::vector<G4RadioactivityTable*> theRadioactivityTables;
00307     G4int decayWindows[99];
00308     static const G4double levelTolerance;
00309 
00310     //User define radioactive decay data files replacing some files in the G4RADECAY database
00311     std::map<G4int, G4String> theUserRadioactiveDataFiles;
00312 
00313 
00314     // Remainder of life time at rest
00315     G4double fRemainderLifeTime;
00316     G4int verboseLevel;
00317 
00318 
00319     // ParticleChange for decay process
00320     G4ParticleChangeForRadDecay fParticleChangeForRadDecay;
00321 
00322     // inline implementations 
00323     inline
00324     G4double AtRestGetPhysicalInteractionLength(const G4Track& track,
00325                                                 G4ForceCondition* condition)
00326     {
00327       fRemainderLifeTime =
00328         G4VRestDiscreteProcess::AtRestGetPhysicalInteractionLength(track, condition);
00329       return fRemainderLifeTime;
00330     }
00331 
00332     inline
00333     G4VParticleChange* AtRestDoIt(const G4Track& theTrack,
00334                                   const G4Step& theStep)
00335       {return DecayIt(theTrack, theStep);}
00336 
00337     inline
00338     G4VParticleChange* PostStepDoIt(const G4Track& theTrack,
00339                                     const G4Step& theStep)
00340       {return DecayIt(theTrack, theStep);}
00341 };
00342 
00343 #endif
00344 

Generated on Mon May 27 17:49:41 2013 for Geant4 by  doxygen 1.4.7