Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/radioactive_decay/src/G4BetaPlusDecay.cc

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 ]

Diff markup

Differences between /processes/hadronic/models/radioactive_decay/src/G4BetaPlusDecay.cc (Version 11.3.0) and /processes/hadronic/models/radioactive_decay/src/G4BetaPlusDecay.cc (Version 10.7)


  1 //                                                  1 //
  2 // *******************************************      2 // ********************************************************************
  3 // * License and Disclaimer                         3 // * License and Disclaimer                                           *
  4 // *                                                4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of th      5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided      6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License      7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/      8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.           9 // * include a list of copyright holders.                             *
 10 // *                                               10 // *                                                                  *
 11 // * Neither the authors of this software syst     11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing fin     12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warran     13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assum     14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file      15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitatio     16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                               17 // *                                                                  *
 18 // * This  code  implementation is the result      18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboratio     19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distri     20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  ag     21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publicati     22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Sof     23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // *******************************************     24 // ********************************************************************
 25 //                                                 25 //
 26 //////////////////////////////////////////////     26 ////////////////////////////////////////////////////////////////////////////////
 27 //                                                 27 //                                                                            //
 28 //  File:   G4BetaPlusDecay.cc                     28 //  File:   G4BetaPlusDecay.cc                                                //
 29 //  Author: D.H. Wright (SLAC)                     29 //  Author: D.H. Wright (SLAC)                                                //
 30 //  Date:   14 November 2014                       30 //  Date:   14 November 2014                                                  //
 31 //  Modifications:                             << 
 32 //    23.08.2023 V.Ivanchenko                  << 
 33 //                                                 31 //                                                                            //
 34 //////////////////////////////////////////////     32 ////////////////////////////////////////////////////////////////////////////////
 35                                                    33 
 36 #include "G4BetaPlusDecay.hh"                      34 #include "G4BetaPlusDecay.hh"
 37 #include "G4BetaDecayCorrections.hh"               35 #include "G4BetaDecayCorrections.hh"
 38 #include "G4IonTable.hh"                           36 #include "G4IonTable.hh"
 39 #include "G4ThreeVector.hh"                        37 #include "G4ThreeVector.hh"
 40 #include "G4LorentzVector.hh"                  << 
 41 #include "G4DynamicParticle.hh"                    38 #include "G4DynamicParticle.hh"
 42 #include "G4DecayProducts.hh"                      39 #include "G4DecayProducts.hh"
 43 #include "G4PhysicalConstants.hh"                  40 #include "G4PhysicalConstants.hh"
 44 #include "G4SystemOfUnits.hh"                      41 #include "G4SystemOfUnits.hh"
 45 #include "G4Positron.hh"                       << 
 46 #include "G4NeutrinoE.hh"                      << 
 47 #include "G4RandomDirection.hh"                << 
 48 #include "G4BetaSpectrumSampler.hh"            << 
 49 #include <iostream>                                42 #include <iostream>
 50 #include <iomanip>                                 43 #include <iomanip>
 51                                                    44 
 52 namespace {                                    << 
 53   const G4double eMass = CLHEP::electron_mass_ << 
 54 }                                              << 
 55                                                << 
 56 G4BetaPlusDecay::G4BetaPlusDecay(const G4Parti     45 G4BetaPlusDecay::G4BetaPlusDecay(const G4ParticleDefinition* theParentNucleus,
 57                                  const G4doubl     46                                  const G4double& branch, const G4double& e0,
 58                                  const G4doubl     47                                  const G4double& excitationE,
 59                                  const G4Ions:     48                                  const G4Ions::G4FloatLevelBase& flb,
 60                                  const G4BetaD     49                                  const G4BetaDecayType& betaType)
 61  : G4NuclearDecay("beta+ decay", BetaPlus, exc     50  : G4NuclearDecay("beta+ decay", BetaPlus, excitationE, flb),
 62    maxEnergy(e0/eMass - 2.0),                  <<  51    endpointEnergy(e0 - 2.*CLHEP::electron_mass_c2)
 63    estep(maxEnergy/(G4double)(npti - 1))       << 
 64 {                                                  52 {
 65   SetParent(theParentNucleus);  // Store name      53   SetParent(theParentNucleus);  // Store name of parent nucleus, delete G4MT_parent 
 66   SetBR(branch);                                   54   SetBR(branch);
 67   SetNumberOfDaughters(3);                     << 
 68                                                << 
 69   fPrimaryIon = theParentNucleus;              << 
 70   fLepton = G4Positron::Positron();            << 
 71   fNeutrino = G4NeutrinoE::NeutrinoE();        << 
 72                                                    55 
 73   G4IonTable* theIonTable = G4ParticleTable::G <<  56   SetNumberOfDaughters(3);
                                                   >>  57   G4IonTable* theIonTable =
                                                   >>  58     (G4IonTable*)(G4ParticleTable::GetParticleTable()->GetIonTable());
 74   G4int daughterZ = theParentNucleus->GetAtomi     59   G4int daughterZ = theParentNucleus->GetAtomicNumber() - 1;
 75   G4int daughterA = theParentNucleus->GetAtomi <<  60   G4int daughterA = theParentNucleus->GetAtomicMass(); 
 76   fResIon = const_cast<const G4ParticleDefinit <<  61   SetDaughter(0, theIonTable->GetIon(daughterZ, daughterA, excitationE, flb) );
 77                                                <<  62   SetUpBetaSpectrumSampler(daughterZ, daughterA, betaType);
                                                   >>  63   SetDaughter(1, "e+");
                                                   >>  64   SetDaughter(2, "nu_e");
                                                   >>  65 }
 78                                                    66 
 79   parentMass = theParentNucleus->GetPDGMass(); << 
 80   resMass = fResIon->GetPDGMass();             << 
 81                                                    67 
 82   SetUpBetaSpectrumSampler(daughterZ, daughter <<  68 G4BetaPlusDecay::~G4BetaPlusDecay()
                                                   >>  69 {
                                                   >>  70   delete spectrumSampler;
                                                   >>  71 }
 83                                                    72 
 84   SetDaughter(0, fResIon);                     << 
 85   SetDaughter(1, fLepton);                     << 
 86   SetDaughter(2, fNeutrino);                   << 
 87                                                    73 
                                                   >>  74 G4DecayProducts* G4BetaPlusDecay::DecayIt(G4double)
                                                   >>  75 {
 88   // Fill G4MT_parent with theParentNucleus (s     76   // Fill G4MT_parent with theParentNucleus (stored by SetParent in ctor)  
 89   CheckAndFillParent();                            77   CheckAndFillParent();
 90                                                    78 
 91   // Fill G4MT_daughters with e+, nu and resid <<  79   // Fill G4MT_daughters with e-, nu and residual nucleus (stored by SetDaughter)  
 92   CheckAndFillDaughters();                         80   CheckAndFillDaughters();
 93 }                                              << 
 94                                                    81 
 95 G4DecayProducts* G4BetaPlusDecay::DecayIt(G4do <<  82   G4double parentMass = G4MT_parent->GetPDGMass();
 96 {                                              <<  83   G4double eMass = G4MT_daughters[1]->GetPDGMass();
                                                   >>  84   G4double nucleusMass = G4MT_daughters[0]->GetPDGMass();
 97   // Set up final state                            85   // Set up final state
 98   // parentParticle is set at rest here becaus     86   // parentParticle is set at rest here because boost with correct momentum 
 99   // is done later                                 87   // is done later
100   G4DynamicParticle prim(fPrimaryIon, G4ThreeV <<  88   G4DynamicParticle parentParticle(G4MT_parent, G4ThreeVector(0,0,0), 0.0);
101   G4DecayProducts* products = new G4DecayProdu <<  89   G4DecayProducts* products = new G4DecayProducts(parentParticle);
102                                                << 
103   // Generate positron isotropic in angle, wit << 
104   const G4double eKE = eMass*G4BetaSpectrumSam << 
105                                                    90 
106   G4double eMomentum = std::sqrt(eKE*(eKE + 2. <<  91   if (spectrumSampler) {
107   G4ThreeVector dir = G4RandomDirection();     <<  92     // Generate positron isotropic in angle, with energy from stored spectrum
108   G4DynamicParticle* dp = new G4DynamicParticl <<  93     G4double eKE = endpointEnergy*spectrumSampler->shoot(G4Random::getTheEngine() );
109   products->PushProducts(dp);                  <<  94     G4double eMomentum = std::sqrt(eKE*(eKE + 2.*eMass) );
110   /*                                           <<  95 
111   G4cout << "G4BetaPlusDecay::DecayIt: " << fP <<  96     G4double cosTheta = 2.*G4UniformRand() - 1.0;
112    << " -> " << fResIon->GetParticleName() <<  <<  97     G4double sinTheta = std::sqrt(1.0 - cosTheta*cosTheta);
113    << " + " << fNeutrino->GetParticleName() << <<  98     G4double phi = twopi*G4UniformRand()*rad;
114          << G4endl;                            <<  99     G4double sinPhi = std::sin(phi);
115   */                                           << 100     G4double cosPhi = std::cos(phi);
116   // 4-momentum of residual ion and neutrino   << 101 
117   G4LorentzVector lv(-eMomentum*dir.x(), -eMom << 102     G4ParticleMomentum eDirection(sinTheta*cosPhi, sinTheta*sinPhi, cosTheta);
118                      parentMass - eKE - eMass) << 103     G4DynamicParticle* dynamicPositron
119                                                << 104       = new G4DynamicParticle(G4MT_daughters[1], eDirection*eMomentum);
120   const G4double elim = CLHEP::eV;             << 105     products->PushProducts(dynamicPositron);
121   // centrum of mass system                    << 106 
122   G4double M = lv.mag();                       << 107     // Generate neutrino with angle relative to positron, and energy from
123   G4double edel = M - resMass;                 << 108     // energy-momentum conservation using endpoint energy of reaction
124   // Free energy should be above limit         << 109     G4double cosThetaENu = 2.*G4UniformRand() - 1.;
125   if (edel >= elim) {                          << 110     G4double eTE = eMass + eKE;
126     // neutrino                                << 111     G4double nuEnergy = ((endpointEnergy - eKE)*(parentMass + nucleusMass - eTE)
127     G4double eNu = 0.5*(M - resMass*resMass/M) << 112           - eMomentum*eMomentum)/(parentMass - eTE + eMomentum*cosThetaENu)/2.;
128     G4LorentzVector lvnu(eNu*G4RandomDirection << 113 
129     lvnu.boost(lv.boostVector());              << 114     G4double sinThetaENu = std::sqrt(1.0 - cosThetaENu*cosThetaENu);
130     dir = lvnu.vect().unit();                  << 115     phi = twopi*G4UniformRand()*rad;
131     dp = new G4DynamicParticle(fNeutrino, dir, << 116     G4double sinPhiNu = std::sin(phi);
132     products->PushProducts(dp);                << 117     G4double cosPhiNu = std::cos(phi);
133                                                << 118 
134     // residual                                << 119     G4ParticleMomentum nuDirection;
135     lv -= lvnu;                                << 120     nuDirection.setX(sinThetaENu*cosPhiNu*cosTheta*cosPhi -
136     dir = lv.vect().unit();                    << 121                      sinThetaENu*sinPhiNu*sinPhi + cosThetaENu*sinTheta*cosPhi);
137     G4double ekin = std::max(lv.e() - resMass, << 122     nuDirection.setY(sinThetaENu*cosPhiNu*cosTheta*sinPhi +
138     dp = new G4DynamicParticle(fResIon, dir, e << 123                      sinThetaENu*sinPhiNu*cosPhi + cosThetaENu*sinTheta*sinPhi);
139     products->PushProducts(dp);                << 124     nuDirection.setZ(-sinThetaENu*cosPhiNu*sinTheta + cosThetaENu*cosTheta);
                                                   >> 125 
                                                   >> 126     G4DynamicParticle* dynamicNeutrino
                                                   >> 127       = new G4DynamicParticle(G4MT_daughters[2], nuDirection*nuEnergy);
                                                   >> 128     products->PushProducts(dynamicNeutrino);
                                                   >> 129 
                                                   >> 130     // Generate daughter nucleus from sum of positron and neutrino 4-vectors:
                                                   >> 131     // p_D = - p_e - p_nu
                                                   >> 132     G4DynamicParticle* dynamicDaughter =
                                                   >> 133       new G4DynamicParticle(G4MT_daughters[0],
                                                   >> 134                             -eDirection*eMomentum - nuDirection*nuEnergy);
                                                   >> 135     products->PushProducts(dynamicDaughter);
140                                                   136 
141   } else {                                        137   } else {
142     // neglecting relativistic kinematic and g << 138     // positron energy below threshold -> no decay
143     dp = new G4DynamicParticle(fNeutrino, G4Ra << 139     G4DynamicParticle* noDecay =
144     products->PushProducts(dp);                << 140       new G4DynamicParticle(G4MT_parent, G4ThreeVector(0,0,0), 0.0);
145     dp = new G4DynamicParticle(fResIon, G4Thre << 141     products->PushProducts(noDecay);
146     products->PushProducts(dp);                << 
147   }                                               142   }
148                                                   143 
                                                   >> 144   // Check energy conservation against endpoint value, not nuclear masses
                                                   >> 145   /*
                                                   >> 146   G4int nProd = products->entries();
                                                   >> 147   G4DynamicParticle* temp = 0;
                                                   >> 148   G4double Esum = 0.0;
                                                   >> 149   for (G4int i = 0; i < nProd; i++) {
                                                   >> 150     temp = products->operator[](i);
                                                   >> 151     Esum += temp->GetKineticEnergy();
                                                   >> 152   }
                                                   >> 153   G4double eCons = (endpointEnergy - Esum)/keV;
                                                   >> 154   if (eCons > 0.001) G4cout << " Beta+ check: eCons (keV) = " << eCons << G4endl;
                                                   >> 155   */
149   return products;                                156   return products;
150 }                                                 157 }
151                                                   158 
152                                                   159 
153 void                                              160 void
154 G4BetaPlusDecay::SetUpBetaSpectrumSampler(cons    161 G4BetaPlusDecay::SetUpBetaSpectrumSampler(const G4int& daughterZ,
155                                           cons    162                                           const G4int& daughterA,
156                                           cons    163                                           const G4BetaDecayType& betaType)
157 {                                                 164 {
158   cdf[0] = 0.0;                                << 165   G4double e0 = endpointEnergy/CLHEP::electron_mass_c2;
                                                   >> 166   G4BetaDecayCorrections corrections(-daughterZ, daughterA);
                                                   >> 167   spectrumSampler = 0;
159                                                   168 
160   // Check for cases in which Q < 2Me (e.g. z6    169   // Check for cases in which Q < 2Me (e.g. z67.a162) 
161   if (maxEnergy > 0.) {                        << 170   if (e0 > 0.) {
162     G4BetaDecayCorrections corrections(-daught << 171     // Array to store spectrum pdf
163                                                << 172     G4int npti = 100;
164     // Fill array to store cumulative spectrum << 173     G4double* pdf = new G4double[npti];
165     G4double ex;  // Positron kinetic energy   << 174 
166     G4double p;   // Positron momentum in unit << 175     G4double e;  // Total positron energy in units of electron mass
167     G4double f;   // Spectral shape function   << 176     G4double p;  // Positron momentum in units of electron mass
168     G4double f0 = 0.0;                         << 177     G4double f;  // Spectral shap function
169     G4double sum = 0.0;                        << 178     for (G4int ptn = 0; ptn < npti; ptn++) {
170     for (G4int i = 1; i < npti-1; ++i) {       << 179       // Calculate simple phase space
171       ex = estep*i;                            << 180       e = 1. + e0*(ptn + 0.5)/G4double(npti);
172       p = std::sqrt(ex*(ex + 2.));             << 181       p = std::sqrt(e*e - 1.);
173       f = p*(1. + ex)*(maxEnergy - ex)*(maxEne << 182       f = p*e*(e0 - e + 1.)*(e0 - e + 1.);
174                                                   183 
175       // Apply Fermi factor to get allowed sha    184       // Apply Fermi factor to get allowed shape
176       f *= corrections.FermiFunction(1. + ex); << 185       f *= corrections.FermiFunction(e);
177                                                   186 
178       // Apply shape factor for forbidden tran    187       // Apply shape factor for forbidden transitions
179       f *= corrections.ShapeFactor(betaType, p << 188       f *= corrections.ShapeFactor(betaType, p, e0-e+1.);
180       sum += f + f0;                           << 189       pdf[ptn] = f;
181       cdf[i] = sum;                            << 
182       f0 = f;                                  << 
183     }                                             190     }
184     cdf[npti-1] = sum + f0;                    << 191     spectrumSampler = new G4RandGeneral(pdf, npti);
185   } else {                                     << 192     delete[] pdf;
186     for (G4int i = 0; i < npti; ++i) { cdf[i]  << 
187   }                                               193   }
188 }                                                 194 }
189                                                   195 
                                                   >> 196 
190 void G4BetaPlusDecay::DumpNuclearInfo()           197 void G4BetaPlusDecay::DumpNuclearInfo()
191 {                                                 198 {
192   G4cout << " G4BetaPlusDecay  " << fPrimaryIo << 199   G4cout << " G4BetaPlusDecay for parent nucleus " << GetParentName() << G4endl;
193    << " -> " << fResIon->GetParticleName() <<  << 200   G4cout << " decays to " << GetDaughterName(0) << " , " << GetDaughterName(1) 
194    << " + " << fNeutrino->GetParticleName() << << 201          << " and " << GetDaughterName(2) << " with branching ratio " << GetBR()
195    << maxEnergy*eMass << " BR=" << GetBR() <<  << 202          << "% and endpoint energy " << endpointEnergy/keV << " keV " << G4endl;
196 }                                                 203 }
                                                   >> 204 
197                                                   205