Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/radioactive_decay/include/G4VRadioactiveDecay.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 //  G4VRadioactiveDecay
 31 //
 32 //  Authors: F. Lei and P.R. Truscott
 33 //  Date:   15 June 2000
 34 //
 35 //  9 August 2017 D.H. Wright (SLAC) introduced multi-threaded mode
 36 //  30 May 2024   V. Ivanchenko rename the class 
 37 //
 38 //  Description: This class is the base for simulation of radioactive
 39 //               decays. It performs alpha, beta, gamma, electron capture,
 40 //               and isomeric transition decays of radioactive nuclei.
 41 //
 42 ///////////////////////////////////////////////////////////////////////////////
 43 
 44 #ifndef G4VRadioactiveDecay_h
 45 #define G4VRadioactiveDecay_h 1
 46 
 47 #include <vector>
 48 #include <map>
 49 #include <CLHEP/Units/SystemOfUnits.h>
 50 
 51 #include "G4ios.hh"
 52 #include "globals.hh"
 53 #include "G4VRestDiscreteProcess.hh"
 54 #include "G4ParticleChangeForRadDecay.hh"
 55 
 56 #include "G4NucleusLimits.hh"
 57 #include "G4ThreeVector.hh"
 58 #include "G4RadioactiveDecayMode.hh"
 59 
 60 class G4Fragment;
 61 class G4RadioactiveDecayMessenger;
 62 class G4PhotonEvaporation;
 63 class G4Ions;
 64 class G4DecayTable;
 65 class G4ITDecay;
 66 
 67 typedef std::map<G4String, G4DecayTable*> DecayTableMap;
 68 
 69 class G4VRadioactiveDecay : public G4VRestDiscreteProcess 
 70 {
 71   // class description
 72 
 73   // Implementation of the radioactive decay process which simulates the
 74   // decays of radioactive nuclei.  These nuclei are submitted to RDM as
 75   // G4Ions.  The required half-lives and decay schemes are retrieved from
 76   // the Radioactivity database which was derived from ENSDF.
 77   // All decay products are submitted back to the particle tracking process
 78   // through the G4ParticleChangeForRadDecay object.
 79   // class description - end 
 80 
 81   public: // with description
 82 
 83     G4VRadioactiveDecay(const G4String& processName="RadioactiveDecay",
 84                         const G4double timeThreshold=-1.0);
 85     ~G4VRadioactiveDecay() override;
 86 
 87     G4bool IsApplicable(const G4ParticleDefinition&) override;
 88     // Return true if the specified isotope is
 89     //  1) defined as "nucleus" and
 90     //  2) it is within theNucleusLimit
 91 
 92     G4VParticleChange* AtRestDoIt(const G4Track& theTrack,
 93                                   const G4Step& theStep) override;
 94 
 95     G4VParticleChange* PostStepDoIt(const G4Track& theTrack,
 96                                     const G4Step& theStep) override;
 97 
 98     void BuildPhysicsTable(const G4ParticleDefinition &) override;
 99 
100     void ProcessDescription(std::ostream& outFile) const override;
101 
102     virtual G4VParticleChange* DecayIt(const G4Track& theTrack,
103                                        const G4Step& theStep);
104 
105     // Return decay table if it exists, if not, load it from file
106     G4DecayTable* GetDecayTable(const G4ParticleDefinition*);
107 
108     // Select a logical volume in which RDM applies
109     void SelectAVolume(const G4String& aVolume);
110 
111     // Remove a logical volume from the RDM applied list
112     void DeselectAVolume(const G4String& aVolume);
113 
114     // Select all logical volumes for the application of RDM
115     void SelectAllVolumes();
116 
117     // Remove all logical volumes from RDM applications
118     void DeselectAllVolumes();
119 
120     // Enable/disable ARM
121     inline void SetARM(G4bool arm) { applyARM = arm; }
122 
123     G4DecayTable* LoadDecayTable(const G4Ions*);
124     // Load the decay data of isotope theParentNucleus
125 
126     void AddUserDecayDataFile(G4int Z, G4int A, const G4String& filename);
127     // Allow the user to replace the radio-active decay data provided in Geant4
128     // by its own data file for a given isotope
129 
130     inline void SetNucleusLimits(G4NucleusLimits theNucleusLimits1)
131       { theNucleusLimits = theNucleusLimits1; }
132     // Sets theNucleusLimits which specifies the range of isotopes
133     // the G4VRadioactiveDecay applies.
134 
135     // Returns theNucleusLimits which specifies the range of isotopes used
136     // by G4VRadioactiveDecay
137     inline G4NucleusLimits GetNucleusLimits() const { return theNucleusLimits; }
138 
139     inline void SetDecayDirection(const G4ThreeVector& theDir) {
140       forceDecayDirection = theDir.unit();
141     }
142 
143     inline const G4ThreeVector& GetDecayDirection() const {
144       return forceDecayDirection; 
145     }
146 
147     inline void SetDecayHalfAngle(G4double halfAngle=0.*CLHEP::deg) {
148       forceDecayHalfAngle = std::min(std::max(0.*CLHEP::deg,halfAngle),180.*CLHEP::deg);
149     }
150 
151     inline G4double GetDecayHalfAngle() const { return forceDecayHalfAngle; }
152 
153     // Force direction (random within half-angle) for "visible" daughters
154     // (applies to electrons, positrons, gammas, neutrons, protons or alphas)
155     inline void SetDecayCollimation(const G4ThreeVector& theDir,
156                                     G4double halfAngle = 0.*CLHEP::deg) {
157       SetDecayDirection(theDir);
158       SetDecayHalfAngle(halfAngle);
159     }
160 
161     // Ignore radioactive decays at rest of nuclides happening
162     // after this (very long) time threshold
163     inline void SetThresholdForVeryLongDecayTime(const G4double inputThreshold) {
164       fThresholdForVeryLongDecayTime = std::max( 0.0, inputThreshold );
165     }
166     inline G4double GetThresholdForVeryLongDecayTime() const {
167       return fThresholdForVeryLongDecayTime;
168     }
169 
170     void StreamInfo(std::ostream& os, const G4String& endline);
171 
172     G4VRadioactiveDecay(const G4VRadioactiveDecay& right) = delete;
173     G4VRadioactiveDecay& operator=(const G4VRadioactiveDecay& right) = delete;
174 
175   protected:
176 
177     G4double GetMeanFreePath(const G4Track& theTrack, G4double previousStepSize,
178                              G4ForceCondition* condition) override;
179 
180     G4double GetMeanLifeTime(const G4Track& theTrack,
181                              G4ForceCondition* condition) override;
182 
183     // sampling of products 
184     void DecayAnalog(const G4Track& theTrack, G4DecayTable*);
185 
186     // sampling products at rest
187     G4DecayProducts* DoDecay(const G4ParticleDefinition&, G4DecayTable*);
188 
189     // Apply directional bias for "visible" daughters (e+-, gamma, n, p, alpha)
190     void CollimateDecay(G4DecayProducts* products);
191     void CollimateDecayProduct(G4DynamicParticle* product);
192     G4ThreeVector ChooseCollimationDirection() const;
193 
194     // ParticleChange for decay process
195     G4ParticleChangeForRadDecay fParticleChangeForRadDecay;
196 
197     G4RadioactiveDecayMessenger* theRadioactiveDecayMessenger;
198     G4PhotonEvaporation* photonEvaporation;
199     G4ITDecay* decayIT;
200 
201     std::vector<G4String> ValidVolumes;
202     G4bool isAllVolumesMode{true};
203 
204     static const G4double levelTolerance;
205 
206     // Library of decay tables
207     static DecayTableMap* master_dkmap;
208 
209   private:
210 
211     G4NucleusLimits theNucleusLimits;
212 
213     G4bool isInitialised{false};
214     G4bool applyARM{true};
215 
216     // Parameters for pre-collimated (biased) decay products
217     G4ThreeVector forceDecayDirection{G4ThreeVector(0., 0., 0.)};
218     G4double forceDecayHalfAngle{0.0};
219     static const G4ThreeVector origin;  // (0,0,0) for convenience
220 
221     // Radioactive decay database directory path 
222     static G4String dirPath;
223 
224     // User define radioactive decay data files replacing some files in the G4RADECAY database
225     static std::map<G4int, G4String>* theUserRDataFiles;
226 
227     // The last RadDecayMode
228     G4RadioactiveDecayMode theRadDecayMode{G4RadioactiveDecayMode::IT};
229  
230     // Ignore radioactive decays at rest of nuclides happening after this (very long) time threshold
231     G4double fThresholdForVeryLongDecayTime;
232 };
233 
234 #endif
235 
236