[819] | 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 | // GEANT4 Class file |
---|
| 29 | // |
---|
| 30 | // |
---|
| 31 | // File name: G4hShellCrossSection |
---|
| 32 | // |
---|
| 33 | // Author: S. Dussoni and A. Mantero (Alfonso.Mantero@ge.infn.it) |
---|
| 34 | // |
---|
| 35 | // History: |
---|
| 36 | // ----------- |
---|
| 37 | // 23 Oct 2001 A. Mantero 1st implementation |
---|
| 38 | // 24 Oct 2001 MGP Cleaned up |
---|
| 39 | // 30 Oct 2001 V.Ivanchenko Include formula (53) |
---|
| 40 | // 07 Oct 2002 V.Ivanchenko Fix in formula (53) |
---|
| 41 | // 22 Apr 2004 S.Saliceti Add testFlag variable and GetCrossSection method |
---|
| 42 | // ---------------------------------------------------------------------------- |
---|
| 43 | |
---|
| 44 | #include "globals.hh" |
---|
| 45 | #include "G4hShellCrossSection.hh" |
---|
| 46 | #include "G4AtomicTransitionManager.hh" |
---|
| 47 | #include "G4Electron.hh" |
---|
| 48 | #include "G4Proton.hh" |
---|
| 49 | #include "G4ParticleDefinition.hh" |
---|
| 50 | |
---|
| 51 | G4hShellCrossSection::G4hShellCrossSection() |
---|
| 52 | { } |
---|
| 53 | |
---|
| 54 | |
---|
| 55 | G4hShellCrossSection::~G4hShellCrossSection() |
---|
| 56 | { } |
---|
| 57 | |
---|
| 58 | |
---|
| 59 | |
---|
| 60 | std::vector<G4double> G4hShellCrossSection::GetCrossSection(G4int Z, |
---|
| 61 | G4double incidentEnergy, |
---|
| 62 | G4double mass, |
---|
| 63 | G4double deltaEnergy, |
---|
| 64 | G4bool) const |
---|
| 65 | { |
---|
| 66 | return CalculateCrossSections(Z,incidentEnergy,mass,deltaEnergy,false); |
---|
| 67 | } |
---|
| 68 | |
---|
| 69 | |
---|
| 70 | std::vector<G4double> G4hShellCrossSection::CalculateCrossSections(G4int Z, |
---|
| 71 | G4double incidentEnergy, |
---|
| 72 | G4double hMass, |
---|
| 73 | G4double deltaEnergy, |
---|
| 74 | G4bool testFlag) const |
---|
| 75 | { |
---|
| 76 | // Cross-sections for proton ionization calculated as in |
---|
| 77 | // "M. Gryzinski, Two-Particle Collisions. I. General Relations for |
---|
| 78 | // Collisions in the Laboratory system, Phys.Rev. 138 A305" |
---|
| 79 | // Other reference papers are Gryzinski's "Paper I" and "Paper II" |
---|
| 80 | // V.Ivanchenko add only implementation of the formula (53) |
---|
| 81 | // last factor neglected because it is 1 with a good accuracy |
---|
| 82 | |
---|
| 83 | const G4AtomicTransitionManager* transitionManager = |
---|
| 84 | G4AtomicTransitionManager::Instance(); |
---|
| 85 | |
---|
| 86 | size_t nShells = transitionManager->NumberOfShells(Z); |
---|
| 87 | |
---|
| 88 | // Vector that stores the calculated cross-sections for each shell: |
---|
| 89 | std::vector<G4double> crossSections; |
---|
| 90 | |
---|
| 91 | // In this loop we calculate cross-section for every shell in the atom |
---|
| 92 | for (size_t k=0; k<nShells; k++) |
---|
| 93 | { |
---|
| 94 | G4double bindingEnergy = transitionManager->Shell(Z,k)->BindingEnergy(); |
---|
| 95 | G4double xEnergy = 0.5*bindingEnergy; |
---|
| 96 | G4double y = incidentEnergy*electron_mass_c2/(xEnergy*hMass); |
---|
| 97 | G4double dele = deltaEnergy + xEnergy; |
---|
| 98 | G4double x = electron_mass_c2/dele; |
---|
| 99 | |
---|
| 100 | G4double aCrossSection = (dele/(xEnergy*(1. + 1./y)) |
---|
| 101 | + 4.*std::log(2.7 + std::sqrt(y))/3.) * x*x*x; |
---|
| 102 | |
---|
| 103 | // Fill the vector of cross sections with the value just calculated |
---|
| 104 | crossSections.push_back(aCrossSection); |
---|
| 105 | |
---|
| 106 | if (testFlag) |
---|
| 107 | { |
---|
| 108 | G4cout <<"Element: " <<Z<<" Shell: "<<k<<" Particle Energy(MeV): "<<incidentEnergy/MeV<<G4endl; |
---|
| 109 | G4cout <<"Delta Ray Energy(MeV): "<< deltaEnergy/MeV<<" Binding Energy(MeV): "<<bindingEnergy/MeV<<G4endl; |
---|
| 110 | G4cout <<"Cross Section: "<<aCrossSection/barn<<" barns"<< G4endl; |
---|
| 111 | } |
---|
| 112 | } |
---|
| 113 | |
---|
| 114 | return crossSections; |
---|
| 115 | } |
---|
| 116 | |
---|
| 117 | |
---|
| 118 | |
---|
| 119 | // Normalization of relative cross-sections to 1 |
---|
| 120 | |
---|
| 121 | std::vector<G4double> G4hShellCrossSection::Probabilities( |
---|
| 122 | G4int Z, |
---|
| 123 | G4double incidentEnergy, |
---|
| 124 | G4double hMass, |
---|
| 125 | G4double deltaEnergy |
---|
| 126 | ) const |
---|
| 127 | { |
---|
| 128 | std::vector<G4double> crossSection = CalculateCrossSections(Z, incidentEnergy, hMass, deltaEnergy,true); |
---|
| 129 | G4double nShells = crossSection.size(); |
---|
| 130 | |
---|
| 131 | // Partial and total cross-section used for normalization of crossSections: |
---|
| 132 | G4double totalCrossSection = 0.; |
---|
| 133 | |
---|
| 134 | // Calculation of total cross-section |
---|
| 135 | for (size_t j=0; j<nShells; j++) |
---|
| 136 | { |
---|
| 137 | totalCrossSection += crossSection[j]; // equivalent to totalCrossSection = totalCrossSection + aCrossSection |
---|
| 138 | } |
---|
| 139 | |
---|
| 140 | // Calculation of normal cross section |
---|
| 141 | for(size_t i=0; i<nShells; i++) |
---|
| 142 | { |
---|
| 143 | crossSection[i] = crossSection[i] / totalCrossSection; |
---|
| 144 | } |
---|
| 145 | // Returns the normalized vector |
---|
| 146 | return crossSection; |
---|
| 147 | } |
---|
| 148 | |
---|
| 149 | |
---|
| 150 | |
---|