Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BGGNucleonElasticXS Class Referencefinal

#include <G4BGGNucleonElasticXS.hh>

+ Inheritance diagram for G4BGGNucleonElasticXS:

Public Member Functions

 G4BGGNucleonElasticXS (const G4ParticleDefinition *)
 
 ~G4BGGNucleonElasticXS () final
 
G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat) final
 
G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm, const G4Material *mat) final
 
G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat) final
 
G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
 
void BuildPhysicsTable (const G4ParticleDefinition &) final
 
void CrossSectionDescription (std::ostream &) const final
 
G4BGGNucleonElasticXSoperator= (const G4BGGNucleonElasticXS &right)=delete
 
 G4BGGNucleonElasticXS (const G4BGGNucleonElasticXS &)=delete
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
virtual G4double ComputeCrossSectionPerElement (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, const G4Element *, const G4Material *mat=nullptr)
 
virtual G4double ComputeIsoCrossSection (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
 
virtual const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
bool ForAllAtomsAndEnergies () const
 
void SetForAllAtomsAndEnergies (G4bool val)
 
const G4StringGetName () const
 
void SetName (const G4String &nam)
 
G4VCrossSectionDataSetoperator= (const G4VCrossSectionDataSet &right)=delete
 
 G4VCrossSectionDataSet (const G4VCrossSectionDataSet &)=delete
 

Additional Inherited Members

- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel {0}
 
G4String name
 

Detailed Description

Definition at line 62 of file G4BGGNucleonElasticXS.hh.

Constructor & Destructor Documentation

◆ G4BGGNucleonElasticXS() [1/2]

G4BGGNucleonElasticXS::G4BGGNucleonElasticXS ( const G4ParticleDefinition * p)
explicit

Definition at line 59 of file G4BGGNucleonElasticXS.cc.

60 : G4VCrossSectionDataSet("BarashenkovGlauberGribov")
61{
62 verboseLevel = 0;
63 fGlauberEnergy = 91.*GeV;
64 fLowEnergy = 14.0*MeV;
65 fNucleon = new G4NucleonNuclearCrossSection();
66 fGlauber = new G4ComponentGGHadronNucleusXsc();
67 fHadron = new G4HadronNucleonXsc();
68
69 theProton = G4Proton::Proton();
70 isProton = (theProton == p);
72
73 if (0 == theA[0]) { Initialise(); }
74}
static G4Proton * Proton()
Definition G4Proton.cc:90
G4VCrossSectionDataSet(const G4String &nam="")
void SetForAllAtomsAndEnergies(G4bool val)

Referenced by G4BGGNucleonElasticXS(), and operator=().

◆ ~G4BGGNucleonElasticXS()

G4BGGNucleonElasticXS::~G4BGGNucleonElasticXS ( )
final

Definition at line 78 of file G4BGGNucleonElasticXS.cc.

79{
80 delete fHadron;
81}

◆ G4BGGNucleonElasticXS() [2/2]

G4BGGNucleonElasticXS::G4BGGNucleonElasticXS ( const G4BGGNucleonElasticXS & )
delete

Member Function Documentation

◆ BuildPhysicsTable()

void G4BGGNucleonElasticXS::BuildPhysicsTable ( const G4ParticleDefinition & p)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 168 of file G4BGGNucleonElasticXS.cc.

169{
170 if(&p == theProton || &p == G4Neutron::Neutron()) {
171 isProton = (theProton == &p);
172
173 } else {
175 ed << "This BGG cross section is applicable only to nucleons and not to "
176 << p.GetParticleName() << G4endl;
177 G4Exception("G4BGGNucleonElasticXS::BuildPhysicsTable", "had001",
178 FatalException, ed);
179 }
180}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
#define G4endl
Definition G4ios.hh:67
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
const G4String & GetParticleName() const

◆ CrossSectionDescription()

void G4BGGNucleonElasticXS::CrossSectionDescription ( std::ostream & outFile) const
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 246 of file G4BGGNucleonElasticXS.cc.

247{
248 outFile << "The Barashenkov-Glauber-Gribov cross section handles elastic\n"
249 << "scattering of protons and neutrons from nuclei using the\n"
250 << "Barashenkov parameterization below 91 GeV and the Glauber-Gribov\n"
251 << "parameterization above 91 GeV. n";
252}

◆ GetElementCrossSection()

G4double G4BGGNucleonElasticXS::GetElementCrossSection ( const G4DynamicParticle * dp,
G4int Z,
const G4Material * mat )
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 105 of file G4BGGNucleonElasticXS.cc.

107{
108 // this method should be called only for Z > 1
109
110 G4double cross = 0.0;
111 G4int Z = std::min(ZZ, 92);
112 G4double ekin = dp->GetKineticEnergy();
113 if(1 == Z) {
114 cross = 1.0115*GetIsoCrossSection(dp,1,1);
115 } else {
116 if(ekin <= fLowEnergy) {
117 cross = (isProton) ? theCoulombFacP[Z] : theCoulombFacN[Z];
118 cross *= CoulombFactor(ekin, Z);
119 } else if(ekin > fGlauberEnergy) {
120 cross = (isProton) ? theGlauberFacP[Z] : theGlauberFacN[Z];
121 cross *= fGlauber->GetElasticGlauberGribov(dp, Z, theA[Z]);
122 } else {
123 cross = fNucleon->GetElasticCrossSection(dp, Z);
124 }
125 }
126#ifdef G4VERBOSE
127 if (verboseLevel > 1) {
128 G4cout << "G4BGGNucleonElasticXS::GetElementCrossSection for "
130 << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
131 << " in nucleus Z= " << Z << " A= " << theA[Z]
132 << " XS(b)= " << cross/barn
133 << G4endl;
134 }
135#endif
136 return cross;
137}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
G4GLOB_DLL std::ostream G4cout
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const

◆ GetIsoCrossSection()

G4double G4BGGNucleonElasticXS::GetIsoCrossSection ( const G4DynamicParticle * dp,
G4int Z,
G4int A,
const G4Isotope * iso = nullptr,
const G4Element * elm = nullptr,
const G4Material * mat = nullptr )
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 142 of file G4BGGNucleonElasticXS.cc.

147{
148 // this method should be called only for Z = 1
149 fHadron->HadronNucleonXscNS(dp->GetDefinition(), theProton,
150 dp->GetKineticEnergy());
151 G4double cross = A*fHadron->GetElasticHadronNucleonXsc();
152
153#ifdef G4VERBOSE
154 if (verboseLevel > 1) {
155 G4cout << "G4BGGNucleonElasticXS::GetIsoCrossSection for "
157 << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
158 << " in nucleus Z=1 A=" << A
159 << " XS(b)= " << cross/barn
160 << G4endl;
161 }
162#endif
163 return cross;
164}
const G4double A[17]

Referenced by GetElementCrossSection().

◆ IsElementApplicable()

G4bool G4BGGNucleonElasticXS::IsElementApplicable ( const G4DynamicParticle * ,
G4int Z,
const G4Material * mat )
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 86 of file G4BGGNucleonElasticXS.cc.

88{
89 return true;
90}

◆ IsIsoApplicable()

G4bool G4BGGNucleonElasticXS::IsIsoApplicable ( const G4DynamicParticle * ,
G4int Z,
G4int A,
const G4Element * elm,
const G4Material * mat )
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 94 of file G4BGGNucleonElasticXS.cc.

98{
99 return (1 == Z);
100}

◆ operator=()

G4BGGNucleonElasticXS & G4BGGNucleonElasticXS::operator= ( const G4BGGNucleonElasticXS & right)
delete

The documentation for this class was generated from the following files: