Geant4  10.02.p02
G4AlphaDecay.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: G4AlphaDecay.cc //
29 // Author: D.H. Wright (SLAC) //
30 // Date: 20 October 2014 //
31 // //
33 
34 #include "G4AlphaDecay.hh"
35 #include "G4IonTable.hh"
36 #include "Randomize.hh"
37 #include "G4ThreeVector.hh"
38 #include "G4DynamicParticle.hh"
39 #include "G4DecayProducts.hh"
40 #include "G4PhysicalConstants.hh"
41 #include "G4SystemOfUnits.hh"
42 #include <iostream>
43 #include <iomanip>
44 
46  const G4double& branch, const G4double& Qvalue,
47  const G4double& excitationE)
48  : G4NuclearDecay("alpha decay", Alpha, excitationE), transitionQ(Qvalue)
49 {
50  SetParent(theParentNucleus); // Store name of parent nucleus, delete G4MT_parent
51  SetBR(branch);
52 
54  G4IonTable* theIonTable =
56  G4int daughterZ = theParentNucleus->GetAtomicNumber() - 2;
57  G4int daughterA = theParentNucleus->GetAtomicMass() - 4;
58  SetDaughter(0, theIonTable->GetIon(daughterZ, daughterA, excitationE) );
59  SetDaughter(1, "alpha");
60 }
61 
62 
64 {}
65 
66 
68 {
69  // Fill G4MT_parent with theParentNucleus (stored by SetParent in ctor)
71 
72  // Fill G4MT_daughters with alpha and residual nucleus (stored by SetDaughter)
74 
75  G4double alphaMass = G4MT_daughters[1]->GetPDGMass();
76  // Excitation energy included in PDG mass
77  G4double nucleusMass = G4MT_daughters[0]->GetPDGMass();
78 
79  // Q value was calculated from atomic masses.
80  // Use it to get correct alpha energy.
81  G4double cmMomentum = std::sqrt(transitionQ*(transitionQ + 2.*alphaMass)*
82  (transitionQ + 2.*nucleusMass)*
83  (transitionQ + 2.*alphaMass + 2.*nucleusMass) )/
84  (transitionQ + alphaMass + nucleusMass)/2.;
85 
86  // Set up final state
87  // parentParticle is set at rest here because boost with correct momentum
88  // is done later
89  G4DynamicParticle parentParticle(G4MT_parent, G4ThreeVector(0,0,0), 0.0);
90  G4DecayProducts* products = new G4DecayProducts(parentParticle);
91 
92  G4double costheta = 2.*G4UniformRand()-1.0;
93  G4double sintheta = std::sqrt(1.0 - costheta*costheta);
95  G4ThreeVector direction(sintheta*std::cos(phi),sintheta*std::sin(phi),
96  costheta);
97 
98  G4double KE = std::sqrt(cmMomentum*cmMomentum + alphaMass*alphaMass)
99  - alphaMass;
100  G4DynamicParticle* daughterparticle =
101  new G4DynamicParticle(G4MT_daughters[1], direction, KE, alphaMass);
102  products->PushProducts(daughterparticle);
103 
104  KE = std::sqrt(cmMomentum*cmMomentum + nucleusMass*nucleusMass) - nucleusMass;
105  daughterparticle =
106  new G4DynamicParticle(G4MT_daughters[0], -1.0*direction, KE, nucleusMass);
107  products->PushProducts(daughterparticle);
108 
109  // Energy conservation check
110  // For alpha decays, do final energy check against reaction Q value
111  // which is well-measured using atomic mass differences. Nuclear masses
112  // should not be used since they are not usually directly measured and we
113  // always decay atoms and not fully stripped nuclei.
114  /*
115  G4int nProd = products->entries();
116  G4DynamicParticle* temp = 0;
117  G4double Esum = 0.0;
118  for (G4int i = 0; i < nProd; i++) {
119  temp = products->operator[](i);
120  Esum += temp->GetKineticEnergy();
121  }
122  G4double eCons = (transitionQ - Esum)/keV;
123  if (eCons > 1.e-07) G4cout << " Alpha decay check: Ediff (keV) = " << eCons << G4endl;
124  */
125  return products;
126 }
127 
128 
130 {
131  G4cout << " G4AlphaDecay for parent nucleus " << GetParentName() << G4endl;
132  G4cout << " decays to " << GetDaughterName(0) << " + " << GetDaughterName(1)
133  << " with branching ratio " << GetBR() << "% and Q value "
134  << transitionQ << G4endl;
135 }
136 
void CheckAndFillDaughters()
const G4double transitionQ
Definition: G4AlphaDecay.hh:56
G4double GetBR() const
void SetBR(G4double value)
CLHEP::Hep3Vector G4ThreeVector
const G4String & GetParentName() const
G4int PushProducts(G4DynamicParticle *aParticle)
G4ParticleDefinition * G4MT_parent
virtual G4DecayProducts * DecayIt(G4double)
Definition: G4AlphaDecay.cc:67
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:491
G4ParticleDefinition ** G4MT_daughters
virtual ~G4AlphaDecay()
Definition: G4AlphaDecay.cc:63
G4AlphaDecay(const G4ParticleDefinition *theParentNucleus, const G4double &theBR, const G4double &Qvalue, const G4double &excitation)
Definition: G4AlphaDecay.cc:45
int G4int
Definition: G4Types.hh:78
G4int GetAtomicNumber() const
G4IonTable * GetIonTable() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
virtual void DumpNuclearInfo()
void SetNumberOfDaughters(G4int value)
static const double twopi
Definition: G4SIunits.hh:75
G4int GetAtomicMass() const
const G4String & GetDaughterName(G4int anIndex) const
static const double rad
Definition: G4SIunits.hh:148
G4double GetPDGMass() const
static G4ParticleTable * GetParticleTable()
void SetParent(const G4ParticleDefinition *particle_type)
void SetDaughter(G4int anIndex, const G4ParticleDefinition *particle_type)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76