Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 //
27 // //
28 // File: G4ITDecay.cc //
29 // Author: D.H. Wright (SLAC) //
30 // Date: 14 November 2014 //
31 // //
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"
40 #include "G4PhotonEvaporation.hh"
41 #include "G4RadioactiveDecay.hh"
42 #include "G4VAtomDeexcitation.hh"
43 #include "G4AtomicShells.hh"
44 #include "G4Electron.hh"
45 #include "G4LossTableManager.hh"
46 #include "G4Fragment.hh"
47 #include "G4SystemOfUnits.hh"
48 #include "G4PhysicalConstants.hh"
49 
50 
52  const G4double& branch, const G4double& Qvalue,
53  const G4double& excitationE)
54  : G4NuclearDecay("IT decay", IT, excitationE, noFloat), transitionQ(Qvalue),
55  applyARM(true), daughterNucleus(nullptr), photoEvap(0)
56 {
57  SetParent(theParentNucleus); // Store name of parent nucleus, delete G4MT_parent
58  SetBR(branch);
59 
60  parentZ = theParentNucleus->GetAtomicNumber();
61  parentA = theParentNucleus->GetAtomicMass();
62 
64  G4IonTable* theIonTable =
66  SetDaughter(0, theIonTable->GetIon(parentZ, parentA, excitationE, noFloat) );
67 
68  // Let G4PhotonEvaporation do the decay
69  photoEvap = new G4PhotonEvaporation;
70  photoEvap->RDMForced(true);
71  photoEvap->SetICM(true);
72 }
73 
74 
76 {
77  delete photoEvap;
78 }
79 
80 
82 {
83  // Fill G4MT_parent with theParentNucleus (stored by SetParent in ctor)
85 
86  // Set up final state
87  // parentParticle is set at rest here because boost with correct momentum
88  // is done later
90  G4ThreeVector(0.,0.,0.) );
91  G4DynamicParticle parentParticle(G4MT_parent, atRest);
92  G4DecayProducts* products = new G4DecayProducts(parentParticle);
93 
94  // Let G4PhotonEvaporation do the decay
95  G4Fragment parentNucleus(parentA, parentZ, atRest);
96  // G4cout << " START of G4ITDecay::DecayIt: " << G4endl;
97 
98  // Check if the old nuclear polarization can be used
99  G4NuclearPolarization* nucPol(nullptr);
100  G4int oldZ(0), oldA(0);
101  G4double oldMass(0.0);
102  if(daughterNucleus) {
103  nucPol = daughterNucleus->GetNuclearPolarization();
104  if(nucPol) {
105  oldZ = daughterNucleus->GetZ_asInt();
106  oldA = daughterNucleus->GetA_asInt();
107  oldMass = daughterNucleus->GetGroundStateMass() +
108  daughterNucleus->GetExcitationEnergy();
109  }
110  }
111  static const G4double mlimit = 10*CLHEP::eV;
112  if (nucPol && oldZ == parentZ && oldA == parentA
113  && std::abs(atRest.e() - oldMass) < mlimit) {
114  // Continue with existing chain
115  parentNucleus = *daughterNucleus;
116  parentNucleus.SetMomentum(atRest);
117  }
118 
119  // G4cout << " BEFORE TRANSITION fragment = " << G4endl;
120  // G4cout << parentNucleus << G4endl;
121 
122  G4Fragment* eOrGamma = photoEvap->EmittedFragment(&parentNucleus);
123 
124  // G4cout << " AFTER TRANSITION " << G4endl;
125  // G4cout << parentNucleus << G4endl;
126 
127  // Check if IT chain has ended with excited isomere
128  // Take care for deletion of cached G4Fragment
129  if (parentNucleus.GetExcitationEnergy() > mlimit) {
130  // G4cout << " IT chain is continue " << G4endl;
131  delete daughterNucleus;
132  daughterNucleus = new G4Fragment(parentNucleus);
133  } else if(daughterNucleus) {
134  // G4cout << " End of IT chain " << G4endl;
135  delete daughterNucleus;
136  daughterNucleus = nullptr;
137  }
138 
139  // Modified nuclide is returned as dynDaughter
140  G4IonTable* theIonTable =
142  G4ParticleDefinition* daughterIon =
143  theIonTable->GetIon(parentZ, parentA, parentNucleus.GetExcitationEnergy(),
145  G4DynamicParticle* dynDaughter = new G4DynamicParticle(daughterIon,
146  parentNucleus.GetMomentum());
147 
148  if (eOrGamma) {
149  G4DynamicParticle* eOrGammaDyn =
151  eOrGamma->GetMomentum() );
152  eOrGammaDyn->SetProperTime(eOrGamma->GetCreationTime() );
153  products->PushProducts(eOrGammaDyn);
154  delete eOrGamma;
155 
156  // Now do atomic relaxation if e- is emitted
157  if (applyARM) {
158  G4int shellIndex = photoEvap->GetVacantShellNumber();
159  if (shellIndex > -1) {
160  G4VAtomDeexcitation* atomDeex =
162  if (atomDeex->IsFluoActive() && parentZ > 5 && parentZ < 100) {
163  G4int nShells = G4AtomicShells::GetNumberOfShells(parentZ);
164  if (shellIndex >= nShells) shellIndex = nShells;
166  const G4AtomicShell* shell = atomDeex->GetAtomicShell(parentZ, as);
167  std::vector<G4DynamicParticle*> armProducts;
168 
169  // VI, SI
170  // Allows fixing of Bugzilla 1727
171  G4double deexLimit = 0.1*keV;
172  if (G4EmParameters::Instance()->DeexcitationIgnoreCut()) deexLimit =0.;
173  //
174 
175  atomDeex->GenerateParticles(&armProducts, shell, parentZ, deexLimit,
176  deexLimit);
177  G4double productEnergy = 0.;
178  for (G4int i = 0; i < G4int(armProducts.size()); i++)
179  productEnergy += armProducts[i]->GetKineticEnergy();
180 
181  G4double deficit = shell->BindingEnergy() - productEnergy;
182  if (deficit > 0.0) {
183  // Add a dummy electron to make up extra energy
184  G4double cosTh = 1.-2.*G4UniformRand();
185  G4double sinTh = std::sqrt(1.- cosTh*cosTh);
186  G4double phi = twopi*G4UniformRand();
187 
188  G4ThreeVector electronDirection(sinTh*std::sin(phi),
189  sinTh*std::cos(phi), cosTh);
190  G4DynamicParticle* extra =
191  new G4DynamicParticle(G4Electron::Electron(), electronDirection,
192  deficit);
193  armProducts.push_back(extra);
194  }
195 
196  G4int nArm = armProducts.size();
197  if (nArm > 0) {
198  G4ThreeVector bst = dynDaughter->Get4Momentum().boostVector();
199  for (G4int i = 0; i < nArm; ++i) {
200  G4DynamicParticle* dp = armProducts[i];
201  G4LorentzVector lv = dp->Get4Momentum().boost(bst);
202  dp->Set4Momentum(lv);
203  products->PushProducts(dp);
204  }
205  }
206  }
207  }
208  } // if ARM on
209  } // eOrGamma
210 
211  products->PushProducts(dynDaughter);
212 
213  // Energy conservation check
214  /*
215  G4int newSize = products->entries();
216  G4DynamicParticle* temp = 0;
217  G4double KEsum = 0.0;
218  for (G4int i = 0; i < newSize; i++) {
219  temp = products->operator[](i);
220  KEsum += temp->GetKineticEnergy();
221  }
222  G4double eCons = G4MT_parent->GetPDGMass() - dynDaughter->GetMass() - KEsum;
223  G4cout << " IT check: Ediff (keV) = " << eCons/keV << G4endl;
224  */
225 // delete photoEvap;
226  // G4cout << " END G4ITDecay::DecayIt " << G4endl;
227  return products;
228 }
229 
230 
232 {
233  G4cout << " G4ITDecay for parent nucleus " << GetParentName() << G4endl;
234  G4cout << " decays to " << GetDaughterName(0)
235  << " + gammas (or electrons), with branching ratio " << GetBR()
236  << "% and Q value " << transitionQ << G4endl;
237 }
238 
G4double GetBR() const
Hep3Vector boostVector() const
G4int GetFloatingLevelNumber() const
Definition: G4Fragment.hh:427
void SetBR(G4double value)
static G4LossTableManager * Instance()
G4bool IsFluoActive() const
CLHEP::Hep3Vector G4ThreeVector
G4ITDecay(const G4ParticleDefinition *theParentNucleus, const G4double &theBR, const G4double &Qvalue, const G4double &excitation)
Definition: G4ITDecay.cc:51
virtual void RDMForced(G4bool)
const G4String & GetParentName() const
G4int PushProducts(G4DynamicParticle *aParticle)
G4ParticleDefinition * G4MT_parent
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:503
virtual G4Fragment * EmittedFragment(G4Fragment *theNucleus) final
const G4ParticleDefinition * GetParticleDefinition() const
Definition: G4Fragment.hh:438
virtual G4DecayProducts * DecayIt(G4double)
Definition: G4ITDecay.cc:81
virtual void DumpNuclearInfo()
Definition: G4ITDecay.cc:231
G4double BindingEnergy() const
int G4int
Definition: G4Types.hh:78
G4int GetAtomicNumber() const
static constexpr double twopi
Definition: G4SIunits.hh:76
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
G4IonTable * GetIonTable() const
static G4Ions::G4FloatLevelBase FloatLevelBase(char flbChar)
Definition: G4Ions.hh:189
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4int GetA_asInt() const
Definition: G4Fragment.hh:266
G4double GetCreationTime() const
Definition: G4Fragment.hh:448
void SetNumberOfDaughters(G4int value)
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:307
HepLorentzVector & boost(double, double, double)
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:312
virtual ~G4ITDecay()
Definition: G4ITDecay.cc:75
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:288
G4int GetAtomicMass() const
G4LorentzVector Get4Momentum() const
const G4String & GetDaughterName(G4int anIndex) const
G4int GetVacantShellNumber() const
static constexpr double eV
void Set4Momentum(const G4LorentzVector &momentum)
G4double GetPDGMass() const
static G4ParticleTable * GetParticleTable()
void SetParent(const G4ParticleDefinition *particle_type)
static G4EmParameters * Instance()
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)
G4int GetZ_asInt() const
Definition: G4Fragment.hh:271
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
G4VAtomDeexcitation * AtomDeexcitation()
double G4double
Definition: G4Types.hh:76
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void SetProperTime(G4double)
static constexpr double keV
Definition: G4SIunits.hh:216
#define noFloat
Definition: G4Ions.hh:118
G4AtomicShellEnumerator
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:283
G4NuclearPolarization * GetNuclearPolarization() const
Definition: G4Fragment.hh:458
static G4int GetNumberOfShells(G4int Z)