Geant4 Cross Reference |
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: G4BetaPlusDecay.cc // 29 // Author: D.H. Wright (SLAC) // 30 // Date: 14 November 2014 // 31 // Modifications: // 32 // 23.08.2023 V.Ivanchenko // 33 // // 34 //////////////////////////////////////////////////////////////////////////////// 35 36 #include "G4BetaPlusDecay.hh" 37 #include "G4BetaDecayCorrections.hh" 38 #include "G4IonTable.hh" 39 #include "G4ThreeVector.hh" 40 #include "G4LorentzVector.hh" 41 #include "G4DynamicParticle.hh" 42 #include "G4DecayProducts.hh" 43 #include "G4PhysicalConstants.hh" 44 #include "G4SystemOfUnits.hh" 45 #include "G4Positron.hh" 46 #include "G4NeutrinoE.hh" 47 #include "G4RandomDirection.hh" 48 #include "G4BetaSpectrumSampler.hh" 49 #include <iostream> 50 #include <iomanip> 51 52 namespace { 53 const G4double eMass = CLHEP::electron_mass_c2; 54 } 55 56 G4BetaPlusDecay::G4BetaPlusDecay(const G4ParticleDefinition* theParentNucleus, 57 const G4double& branch, const G4double& e0, 58 const G4double& excitationE, 59 const G4Ions::G4FloatLevelBase& flb, 60 const G4BetaDecayType& betaType) 61 : G4NuclearDecay("beta+ decay", BetaPlus, excitationE, flb), 62 maxEnergy(e0/eMass - 2.0), 63 estep(maxEnergy/(G4double)(npti - 1)) 64 { 65 SetParent(theParentNucleus); // Store name of parent nucleus, delete G4MT_parent 66 SetBR(branch); 67 SetNumberOfDaughters(3); 68 69 fPrimaryIon = theParentNucleus; 70 fLepton = G4Positron::Positron(); 71 fNeutrino = G4NeutrinoE::NeutrinoE(); 72 73 G4IonTable* theIonTable = G4ParticleTable::GetParticleTable()->GetIonTable(); 74 G4int daughterZ = theParentNucleus->GetAtomicNumber() - 1; 75 G4int daughterA = theParentNucleus->GetAtomicMass(); 76 fResIon = const_cast<const G4ParticleDefinition*>(theIonTable->GetIon(daughterZ, daughterA, 77 excitationE, flb)); 78 79 parentMass = theParentNucleus->GetPDGMass(); 80 resMass = fResIon->GetPDGMass(); 81 82 SetUpBetaSpectrumSampler(daughterZ, daughterA, betaType); 83 84 SetDaughter(0, fResIon); 85 SetDaughter(1, fLepton); 86 SetDaughter(2, fNeutrino); 87 88 // Fill G4MT_parent with theParentNucleus (stored by SetParent in ctor) 89 CheckAndFillParent(); 90 91 // Fill G4MT_daughters with e+, nu and residual nucleus (stored by SetDaughter) 92 CheckAndFillDaughters(); 93 } 94 95 G4DecayProducts* G4BetaPlusDecay::DecayIt(G4double) 96 { 97 // Set up final state 98 // parentParticle is set at rest here because boost with correct momentum 99 // is done later 100 G4DynamicParticle prim(fPrimaryIon, G4ThreeVector(0,0,1), 0.0); 101 G4DecayProducts* products = new G4DecayProducts(prim); 102 103 // Generate positron isotropic in angle, with energy from stored spectrum 104 const G4double eKE = eMass*G4BetaSpectrumSampler::shoot(npti, cdf, estep); 105 106 G4double eMomentum = std::sqrt(eKE*(eKE + 2.*eMass)); 107 G4ThreeVector dir = G4RandomDirection(); 108 G4DynamicParticle* dp = new G4DynamicParticle(fLepton, dir, eKE); 109 products->PushProducts(dp); 110 /* 111 G4cout << "G4BetaPlusDecay::DecayIt: " << fPrimaryIon->GetParticleName() 112 << " -> " << fResIon->GetParticleName() << " + " << fLepton->GetParticleName() 113 << " + " << fNeutrino->GetParticleName() << " Ee(MeV)=" << eKE 114 << G4endl; 115 */ 116 // 4-momentum of residual ion and neutrino 117 G4LorentzVector lv(-eMomentum*dir.x(), -eMomentum*dir.y(), -eMomentum*dir.z(), 118 parentMass - eKE - eMass); 119 120 const G4double elim = CLHEP::eV; 121 // centrum of mass system 122 G4double M = lv.mag(); 123 G4double edel = M - resMass; 124 // Free energy should be above limit 125 if (edel >= elim) { 126 // neutrino 127 G4double eNu = 0.5*(M - resMass*resMass/M); 128 G4LorentzVector lvnu(eNu*G4RandomDirection(), eNu); 129 lvnu.boost(lv.boostVector()); 130 dir = lvnu.vect().unit(); 131 dp = new G4DynamicParticle(fNeutrino, dir, lvnu.e()); 132 products->PushProducts(dp); 133 134 // residual 135 lv -= lvnu; 136 dir = lv.vect().unit(); 137 G4double ekin = std::max(lv.e() - resMass, 0.0); 138 dp = new G4DynamicParticle(fResIon, dir, ekin); 139 products->PushProducts(dp); 140 141 } else { 142 // neglecting relativistic kinematic and giving all energy to neutrino 143 dp = new G4DynamicParticle(fNeutrino, G4RandomDirection(), elim); 144 products->PushProducts(dp); 145 dp = new G4DynamicParticle(fResIon, G4ThreeVector(0.0,0.0,1.0), 0.0); 146 products->PushProducts(dp); 147 } 148 149 return products; 150 } 151 152 153 void 154 G4BetaPlusDecay::SetUpBetaSpectrumSampler(const G4int& daughterZ, 155 const G4int& daughterA, 156 const G4BetaDecayType& betaType) 157 { 158 cdf[0] = 0.0; 159 160 // Check for cases in which Q < 2Me (e.g. z67.a162) 161 if (maxEnergy > 0.) { 162 G4BetaDecayCorrections corrections(-daughterZ, daughterA); 163 164 // Fill array to store cumulative spectrum 165 G4double ex; // Positron kinetic energy 166 G4double p; // Positron momentum in units of electron mass 167 G4double f; // Spectral shape function 168 G4double f0 = 0.0; 169 G4double sum = 0.0; 170 for (G4int i = 1; i < npti-1; ++i) { 171 ex = estep*i; 172 p = std::sqrt(ex*(ex + 2.)); 173 f = p*(1. + ex)*(maxEnergy - ex)*(maxEnergy - ex); 174 175 // Apply Fermi factor to get allowed shape 176 f *= corrections.FermiFunction(1. + ex); 177 178 // Apply shape factor for forbidden transitions 179 f *= corrections.ShapeFactor(betaType, p, maxEnergy - ex); 180 sum += f + f0; 181 cdf[i] = sum; 182 f0 = f; 183 } 184 cdf[npti-1] = sum + f0; 185 } else { 186 for (G4int i = 0; i < npti; ++i) { cdf[i] = 0.0; } 187 } 188 } 189 190 void G4BetaPlusDecay::DumpNuclearInfo() 191 { 192 G4cout << " G4BetaPlusDecay " << fPrimaryIon->GetParticleName() 193 << " -> " << fResIon->GetParticleName() << " + " << fLepton->GetParticleName() 194 << " + " << fNeutrino->GetParticleName() << " Eemax(MeV)=" 195 << maxEnergy*eMass << " BR=" << GetBR() << "%" << G4endl; 196 } 197