Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ITDecay.cc
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////////////////////////////////////////////////////////////////////////////////
27// //
28// File: G4ITDecay.cc //
29// Author: D.H. Wright (SLAC) //
30// Date: 14 November 2014 //
31// //
32////////////////////////////////////////////////////////////////////////////////
33
34#include "G4ITDecay.hh"
35#include "G4IonTable.hh"
36#include "G4ThreeVector.hh"
37#include "G4LorentzVector.hh"
38#include "G4DynamicParticle.hh"
39#include "G4DecayProducts.hh"
42#include "G4AtomicShells.hh"
43#include "G4Electron.hh"
44#include "G4LossTableManager.hh"
45#include "G4Fragment.hh"
46#include "G4SystemOfUnits.hh"
48
49
51 : G4NuclearDecay("IT decay", IT, 0.0, noFloat), photonEvaporation(ptr)
52{}
53
55 const G4double& branch, const G4double&,
56 const G4double& excitationE)
57 : G4NuclearDecay("IT decay", IT, excitationE, noFloat)
58{
59 SetParent(theParentNucleus); // Store name of parent nucleus, delete G4MT_parent
60 SetBR(branch);
62 SetDaughter(0, theParentNucleus);
63
64 SetupDecay(theParentNucleus);
65}
66
67void G4ITDecay::SetupDecay(const G4ParticleDefinition* theParentNucleus)
68{
69 theParent = theParentNucleus;
70 parentZ = theParentNucleus->GetAtomicNumber();
71 parentA = theParentNucleus->GetAtomicMass();
72}
73
75{
76 // Set up final state
77 // parentParticle is set at rest here because boost with correct momentum
78 // is done later
79 G4LorentzVector atRest(theParent->GetPDGMass(), G4ThreeVector(0.,0.,0.) );
80 G4DynamicParticle parentParticle(theParent, atRest);
81 G4DecayProducts* products = new G4DecayProducts(parentParticle);
82
83 // Let G4PhotonEvaporation do the decay
84 G4Fragment parentNucleus(parentA, parentZ, atRest);
85
86 // one emission, parent nucleaus become less excited
87 G4Fragment* eOrGamma = photonEvaporation->EmittedFragment(&parentNucleus);
88
89 // Modified nuclide is returned as dynDaughter
90 auto theIonTable = G4ParticleTable::GetParticleTable()->GetIonTable();
91 G4ParticleDefinition* daughterIon =
92 theIonTable->GetIon(parentZ, parentA, parentNucleus.GetExcitationEnergy(),
94 G4DynamicParticle* dynDaughter = new G4DynamicParticle(daughterIon,
95 parentNucleus.GetMomentum());
96
97 if (nullptr != eOrGamma) {
98 G4DynamicParticle* eOrGammaDyn =
100 eOrGamma->GetMomentum() );
101 eOrGammaDyn->SetProperTime(eOrGamma->GetCreationTime() );
102 products->PushProducts(eOrGammaDyn);
103 delete eOrGamma;
104
105 // Now do atomic relaxation if e- is emitted
106 if (applyARM) {
107 G4int shellIndex = photonEvaporation->GetVacantShellNumber();
108 if (shellIndex > -1) {
109 G4VAtomDeexcitation* atomDeex =
111 if (atomDeex->IsFluoActive() && parentZ > 5 && parentZ < 105) {
112 G4int nShells = G4AtomicShells::GetNumberOfShells(parentZ);
113 if (shellIndex >= nShells) shellIndex = nShells;
115 const G4AtomicShell* shell = atomDeex->GetAtomicShell(parentZ, as);
116 std::vector<G4DynamicParticle*> armProducts;
117
118 // VI, SI
119 // Allows fixing of Bugzilla 1727
120 G4double deexLimit = 0.1*keV;
121 if (G4EmParameters::Instance()->DeexcitationIgnoreCut()) deexLimit =0.;
122 //
123
124 atomDeex->GenerateParticles(&armProducts, shell, parentZ, deexLimit,
125 deexLimit);
126 G4double productEnergy = 0.;
127 for (G4int i = 0; i < G4int(armProducts.size()); i++)
128 productEnergy += armProducts[i]->GetKineticEnergy();
129
130 G4double deficit = shell->BindingEnergy() - productEnergy;
131 if (deficit > 0.0) {
132 // Add a dummy electron to make up extra energy
133 G4double cosTh = 1.-2.*G4UniformRand();
134 G4double sinTh = std::sqrt(1.- cosTh*cosTh);
135 G4double phi = twopi*G4UniformRand();
136
137 G4ThreeVector electronDirection(sinTh*std::sin(phi),
138 sinTh*std::cos(phi), cosTh);
139 G4DynamicParticle* extra =
140 new G4DynamicParticle(G4Electron::Electron(), electronDirection,
141 deficit);
142 armProducts.push_back(extra);
143 }
144
145 std::size_t nArm = armProducts.size();
146 if (nArm > 0) {
147 G4ThreeVector bst = dynDaughter->Get4Momentum().boostVector();
148 for (std::size_t i = 0; i < nArm; ++i) {
149 G4DynamicParticle* dp = armProducts[i];
150 G4LorentzVector lv = dp->Get4Momentum().boost(bst);
151 dp->Set4Momentum(lv);
152 products->PushProducts(dp);
153 }
154 }
155 }
156 }
157 } // if ARM on
158 } // eOrGamma
159
160 products->PushProducts(dynDaughter);
161
162 // Energy conservation check
163 /*
164 G4int newSize = products->entries();
165 G4DynamicParticle* temp = 0;
166 G4double KEsum = 0.0;
167 for (G4int i = 0; i < newSize; i++) {
168 temp = products->operator[](i);
169 KEsum += temp->GetKineticEnergy();
170 }
171 G4double eCons = G4MT_parent->GetPDGMass() - dynDaughter->GetMass() - KEsum;
172 G4cout << " IT check: Ediff (keV) = " << eCons/keV << G4endl;
173 */
174 return products;
175}
176
177
179{
180 if (theParent != nullptr) {
181 G4cout << " G4ITDecay for parent nucleus " << theParent->GetParticleName() << G4endl;
182 }
183}
184
#define noFloat
Definition G4Ions.hh:119
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
G4double BindingEnergy() const
static G4int GetNumberOfShells(G4int Z)
G4int PushProducts(G4DynamicParticle *aParticle)
void SetProperTime(G4double)
G4LorentzVector Get4Momentum() const
void Set4Momentum(const G4LorentzVector &momentum)
static G4Electron * Electron()
Definition G4Electron.cc:91
static G4EmParameters * Instance()
G4double GetExcitationEnergy() const
const G4LorentzVector & GetMomentum() const
G4double GetCreationTime() const
G4int GetFloatingLevelNumber() const
const G4ParticleDefinition * GetParticleDefinition() const
void SetupDecay(const G4ParticleDefinition *)
Definition G4ITDecay.cc:67
void DumpNuclearInfo() override
Definition G4ITDecay.cc:178
G4ITDecay(G4PhotonEvaporation *aPhotonEvap)
Definition G4ITDecay.cc:50
G4DecayProducts * DecayIt(G4double) override
Definition G4ITDecay.cc:74
static G4Ions::G4FloatLevelBase FloatLevelBase(char flbChar)
Definition G4Ions.cc:107
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
G4int GetAtomicNumber() const
G4int GetAtomicMass() const
const G4String & GetParticleName() const
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
G4int GetVacantShellNumber() const
G4Fragment * EmittedFragment(G4Fragment *theNucleus) override
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void SetBR(G4double value)
void SetNumberOfDaughters(G4int value)
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)
void SetParent(const G4ParticleDefinition *particle_type)