Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4NucleonNuclearCrossSection.hh
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: V.Ivanchenko
33//
34
35#ifndef G4NucleonNuclearCrossSection_h
36#define G4NucleonNuclearCrossSection_h
37
40#include "globals.hh"
41
44
46{
47public:
48
51
52 static const char* Default_Name() {return "BarashenkovNucleonXS";}
53
55 G4int Z, const G4Material* mat) final;
56
58 G4int Z, const G4Material* mat=nullptr) final;
59
60 void BuildPhysicsTable(const G4ParticleDefinition&) final;
61
62 void CrossSectionDescription(std::ostream&) const final;
63
65 G4int Z);
66
67 inline G4double GetTotalXsc() { return fTotalXsc; };
68 inline G4double GetInelasticXsc() { return fInelasticXsc; };
69 inline G4double GetElasticXsc() { return fElasticXsc; };
70
71private:
72
73 void ComputeCrossSections(const G4ParticleDefinition*,
74 G4double kinEnergy, G4int Z);
75
77
78 const G4ParticleDefinition* theProton;
79 const G4ParticleDefinition* theNeutron;
80
81 G4double fTotalXsc;
82 G4double fInelasticXsc;
83 G4double fElasticXsc;
84
85};
86
87inline
89 const G4DynamicParticle* dp, G4int Z)
90{
91 ComputeCrossSections(dp->GetDefinition(), dp->GetKineticEnergy(), Z);
92 return fElasticXsc;
93}
94
95#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4double GetElementCrossSection(const G4DynamicParticle *aParticle, G4int Z, const G4Material *mat=nullptr) final
G4double GetElasticCrossSection(const G4DynamicParticle *aParticle, G4int Z)
void CrossSectionDescription(std::ostream &) const final
G4bool IsElementApplicable(const G4DynamicParticle *aParticle, G4int Z, const G4Material *mat) final
#define inline
Definition: internal.h:104