Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BGGNucleonInelasticXS Class Reference

#include <G4BGGNucleonInelasticXS.hh>

+ Inheritance diagram for G4BGGNucleonInelasticXS:

Public Member Functions

 G4BGGNucleonInelasticXS (const G4ParticleDefinition *)
 
 ~G4BGGNucleonInelasticXS () override
 
G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat) override
 
G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm, const G4Material *mat) override
 
G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat) override
 
G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
void CrossSectionDescription (std::ostream &) const override
 
- 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
 
G4String name
 

Detailed Description

Definition at line 63 of file G4BGGNucleonInelasticXS.hh.

Constructor & Destructor Documentation

◆ G4BGGNucleonInelasticXS()

G4BGGNucleonInelasticXS::G4BGGNucleonInelasticXS ( const G4ParticleDefinition * p)
explicit

Definition at line 74 of file G4BGGNucleonInelasticXS.cc.

75 : G4VCrossSectionDataSet("BarashenkovGlauberGribov")
76{
77 verboseLevel = 0;
78 fGlauberEnergy = 91.*GeV;
79 fLowEnergy = 14.*MeV;
80
81 fNucleon = nullptr;
82 fGlauber = nullptr;
83 fHadron = nullptr;
84
85 theProton= G4Proton::Proton();
86 isProton = (theProton == p);
87 isMaster = false;
89}
static G4Proton * Proton()
Definition G4Proton.cc:90
G4VCrossSectionDataSet(const G4String &nam="")
void SetForAllAtomsAndEnergies(G4bool val)

◆ ~G4BGGNucleonInelasticXS()

G4BGGNucleonInelasticXS::~G4BGGNucleonInelasticXS ( )
override

Definition at line 93 of file G4BGGNucleonInelasticXS.cc.

94{
95 delete fHadron;
96}

Member Function Documentation

◆ BuildPhysicsTable()

void G4BGGNucleonInelasticXS::BuildPhysicsTable ( const G4ParticleDefinition & p)
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 177 of file G4BGGNucleonInelasticXS.cc.

178{
179 if(nullptr != fNucleon) { return; }
180 if(&p == theProton || &p == G4Neutron::Neutron()) {
181 isProton = (theProton == &p);
182 } else {
184 ed << "This BGG cross section is applicable only to nucleons and not to "
185 << p.GetParticleName() << G4endl;
186 G4Exception("G4BGGNucleonInelasticXS::BuildPhysicsTable", "had001",
187 FatalException, ed);
188 return;
189 }
190
191 fNucleon = new G4NucleonNuclearCrossSection();
192 fGlauber = new G4ComponentGGHadronNucleusXsc();
193 fHadron = new G4HadronNucleonXsc();
194
195 fNucleon->BuildPhysicsTable(p);
196
197 if(0 == theA[0]) {
198#ifdef G4MULTITHREADED
199 G4MUTEXLOCK(&nucleonInelasticXSMutex);
200 if(0 == theA[0]) {
201#endif
202 isMaster = true;
203#ifdef G4MULTITHREADED
204 }
205 G4MUTEXUNLOCK(&nucleonInelasticXSMutex);
206#endif
207 } else {
208 return;
209 }
210
211 if(isMaster && 0 == theA[0]) {
212
213 theA[0] = theA[1] = 1;
214 G4ThreeVector mom(0.0,0.0,1.0);
215 G4DynamicParticle dp(theProton, mom, fGlauberEnergy);
216
218 G4double csup, csdn;
219
220 if(verboseLevel > 0) {
221 G4cout << "### G4BGGNucleonInelasticXS::Initialise for "
222 << p.GetParticleName() << G4endl;
223 }
224 for(G4int iz=2; iz<93; ++iz) {
225
226 G4int A = G4lrint(nist->GetAtomicMassAmu(iz));
227 theA[iz] = A;
228
229 csup = fGlauber->GetInelasticGlauberGribov(&dp, iz, A);
230 csdn = fNucleon->GetElementCrossSection(&dp, iz);
231 theGlauberFacP[iz] = csdn/csup;
232 }
233
234 dp.SetDefinition(G4Neutron::Neutron());
235 for(G4int iz=2; iz<93; ++iz) {
236 csup = fGlauber->GetInelasticGlauberGribov(&dp, iz, theA[iz]);
237 csdn = fNucleon->GetElementCrossSection(&dp, iz);
238 theGlauberFacN[iz] = csdn/csup;
239
240 if(verboseLevel > 0) {
241 G4cout << "Z= " << iz << " A= " << theA[iz]
242 << " GFactorP= " << theGlauberFacP[iz]
243 << " GFactorN= " << theGlauberFacN[iz] << G4endl;
244 }
245 }
246 theCoulombFacP[1] = theCoulombFacN[1] = 1.0;
247 dp.SetDefinition(theProton);
248 dp.SetKineticEnergy(fLowEnergy);
249 for(G4int iz=2; iz<93; ++iz) {
250 theCoulombFacP[iz] = fNucleon->GetElementCrossSection(&dp, iz)
251 /CoulombFactor(fLowEnergy, iz);
252 }
253 dp.SetDefinition(G4Neutron::Neutron());
254 for(G4int iz=2; iz<93; ++iz) {
255 theCoulombFacN[iz] = fNucleon->GetElementCrossSection(&dp, iz)
256 /CoulombFactor(fLowEnergy, iz);
257
258 if(verboseLevel > 0) {
259 G4cout << "Z= " << iz << " A= " << theA[iz]
260 << " CFactorP= " << theCoulombFacP[iz]
261 << " CFactorN= " << theCoulombFacN[iz] << G4endl;
262 }
263 }
264 }
265}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
#define G4MUTEXLOCK(mutex)
#define G4MUTEXUNLOCK(mutex)
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4double GetElementCrossSection(const G4DynamicParticle *aParticle, G4int Z, const G4Material *mat=nullptr) final
const G4String & GetParticleName() const
int G4lrint(double ad)
Definition templates.hh:134

◆ CrossSectionDescription()

void G4BGGNucleonInelasticXS::CrossSectionDescription ( std::ostream & outFile) const
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 310 of file G4BGGNucleonInelasticXS.cc.

311{
312 outFile << "The Barashenkov-Glauber-Gribov cross section calculates inelastic\n"
313 << "scattering of protons and neutrons from nuclei using the\n"
314 << "Barashenkov parameterization below 91 GeV and the Glauber-Gribov\n"
315 << "parameterization above 91 GeV. It uses the G4HadronNucleonXsc\n"
316 << "cross section component for hydrogen targets, and the\n"
317 << "G4ComponentGGHadronNucleusXsc component for other targets.\n";
318}

◆ GetElementCrossSection()

G4double G4BGGNucleonInelasticXS::GetElementCrossSection ( const G4DynamicParticle * dp,
G4int Z,
const G4Material * mat )
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 119 of file G4BGGNucleonInelasticXS.cc.

121{
122 G4double cross = 0.0;
123 G4double ekin = dp->GetKineticEnergy();
124 G4int Z = std::min(ZZ, 92);
125 if(1 == Z) {
126 cross = 1.0115*GetIsoCrossSection(dp,1,1);
127 } else {
128 if(ekin <= fLowEnergy) {
129 cross = (isProton) ? theCoulombFacP[Z] : theCoulombFacN[Z];
130 cross *= CoulombFactor(ekin, Z);
131 } else if(ekin > fGlauberEnergy) {
132 cross = (isProton) ? theGlauberFacP[Z] : theGlauberFacN[Z];
133 cross *= fGlauber->GetInelasticGlauberGribov(dp, Z, theA[Z]);
134 } else {
135 cross = fNucleon->GetElementCrossSection(dp, Z);
136 }
137 }
138
139 if(verboseLevel > 1) {
140 G4cout << "G4BGGNucleonInelasticXS::GetCrossSection for "
142 << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
143 << " in nucleus Z= " << Z << " A= " << theA[Z]
144 << " XS(b)= " << cross/barn
145 << G4endl;
146 }
147 return cross;
148}
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) override
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const

◆ GetIsoCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 153 of file G4BGGNucleonInelasticXS.cc.

158{
159 // this method should be called only for Z = 1
160 fHadron->HadronNucleonXscNS(dp->GetDefinition(), theProton,
161 dp->GetKineticEnergy());
162 G4double cross = A*fHadron->GetInelasticHadronNucleonXsc();
163
164 if(verboseLevel > 1) {
165 G4cout << "G4BGGNucleonInelasticXS::GetIsoCrossSection for "
167 << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
168 << " in nucleus Z= " << Z << " A= " << theA[Z]
169 << " XS(b)= " << cross/barn
170 << G4endl;
171 }
172 return cross;
173}
G4double GetInelasticHadronNucleonXsc() const
G4double HadronNucleonXscNS(const G4ParticleDefinition *theParticle, const G4ParticleDefinition *nucleon, G4double ekin)

Referenced by GetElementCrossSection().

◆ IsElementApplicable()

G4bool G4BGGNucleonInelasticXS::IsElementApplicable ( const G4DynamicParticle * ,
G4int Z,
const G4Material * mat )
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Reimplemented in G4ParticleHPBGGNucleonInelasticXS.

Definition at line 100 of file G4BGGNucleonInelasticXS.cc.

102{
103 return true;
104}

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Reimplemented in G4ParticleHPBGGNucleonInelasticXS.

Definition at line 108 of file G4BGGNucleonInelasticXS.cc.

112{
113 return (1 == Z);
114}

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