Geant4-11
G4NucleonNuclearCrossSection.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//
26// author: Vladimir.Grichine@cern.ch
27//
28// Implements data from: Barashenkov V.S., Nucleon-Nucleus Cross Section,
29// Preprint JINR P2-89-770, p. 12, Dubna 1989 (scanned version from KEK)
30// Based on G. Folger version of G4PiNuclearCrossSection class
31//
32// Modified: 16.08.2018 V.Ivanchenko - major update
33//
34
36#include "G4DynamicParticle.hh"
38#include "G4Neutron.hh"
39#include "G4Proton.hh"
41
43
45 : G4VCrossSectionDataSet(Default_Name()),
46 fTotalXsc(0.0), fInelasticXsc(0.0), fElasticXsc(0.0)
47{
51}
52
54//
55
57{}
58
60
62 const G4DynamicParticle*, G4int Z, const G4Material*)
63{
64 return (Z > 1);
65}
66
68
70 const G4DynamicParticle* dp, G4int Z, const G4Material*)
71{
73 return fInelasticXsc;
74}
75
77
79 const G4ParticleDefinition* pd,
80 G4double kinEnergy, G4int Z)
81{
82 fBarash->ComputeCrossSections(pd, kinEnergy, Z);
86}
87
89
91{
93}
94
96
97void
99{
100 outFile << "G4NucleonNuclearCrossSection is a variant of the Barashenkov\n"
101 << "cross section parameterization to be used of protons and\n"
102 << "nucleons on targets heavier than hydrogen. It is intended for\n"
103 << "use as a cross section component and is currently used by\n"
104 << "G4BGGNucleonInelasticXS. It is valid for incident energies up\n"
105 << "to 1 TeV.\n";
106}
107
109
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
void BuildPhysicsTable(const G4ParticleDefinition &) final
void ComputeCrossSections(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z)
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
const G4ParticleDefinition * theNeutron
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4double GetElementCrossSection(const G4DynamicParticle *aParticle, G4int Z, const G4Material *mat=nullptr) final
void ComputeCrossSections(const G4ParticleDefinition *, G4double kinEnergy, G4int Z)
G4ComponentBarNucleonNucleusXsc * fBarash
void CrossSectionDescription(std::ostream &) const final
G4bool IsElementApplicable(const G4DynamicParticle *aParticle, G4int Z, const G4Material *mat) final
const G4ParticleDefinition * theProton
static G4Proton * Proton()
Definition: G4Proton.cc:92