Geant4 9.6.0
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 *)
 
virtual ~G4BGGNucleonInelasticXS ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
 
virtual G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
const G4StringGetName () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4VCrossSectionDataSet
void SetName (const G4String &)
 
- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 

Detailed Description

Definition at line 65 of file G4BGGNucleonInelasticXS.hh.

Constructor & Destructor Documentation

◆ G4BGGNucleonInelasticXS()

G4BGGNucleonInelasticXS::G4BGGNucleonInelasticXS ( const G4ParticleDefinition p)

Definition at line 62 of file G4BGGNucleonInelasticXS.cc.

63 : G4VCrossSectionDataSet("Barashenkov-Glauber")
64{
65 verboseLevel = 0;
66 fGlauberEnergy = 91.*GeV;
67 fLowEnergy = 14.*MeV;
68 fHighEnergy = 5.*GeV;
69 fSAIDHighEnergyLimit = 1.3*GeV;
70 for (G4int i = 0; i < 93; ++i) {
71 theGlauberFac[i] = 0.0;
72 theCoulombFac[i] = 0.0;
73 theA[i] = 1;
74 }
75 fNucleon = 0;
76 fGlauber = 0;
77 fHadron = 0;
78 // fGHEISHA = 0;
79 fSAID = 0;
80 particle = p;
81 theProton= G4Proton::Proton();
82 isProton = false;
83 isInitialized = false;
84}
int G4int
Definition: G4Types.hh:66
static G4Proton * Proton()
Definition: G4Proton.cc:93

◆ ~G4BGGNucleonInelasticXS()

G4BGGNucleonInelasticXS::~G4BGGNucleonInelasticXS ( )
virtual

Definition at line 88 of file G4BGGNucleonInelasticXS.cc.

89{
90 delete fHadron;
91 delete fSAID;
92}

Member Function Documentation

◆ BuildPhysicsTable()

void G4BGGNucleonInelasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 193 of file G4BGGNucleonInelasticXS.cc.

194{
195 if(&p == theProton || &p == G4Neutron::Neutron()) {
196 particle = &p;
197 } else {
198 G4cout << "### G4BGGNucleonInelasticXS WARNING: is not applicable to "
199 << p.GetParticleName()
200 << G4endl;
201 throw G4HadronicException(__FILE__, __LINE__,
202 "G4BGGNucleonElasticXS::BuildPhysicsTable is used for wrong particle");
203 return;
204 }
205
206 if(isInitialized) { return; }
207 isInitialized = true;
208
211
212 fHadron = new G4HadronNucleonXsc();
213 //fGHEISHA = new G4HadronInelasticDataSet();
214 fSAID = new G4ComponentSAIDTotalXS();
215
216 fNucleon->BuildPhysicsTable(*particle);
217 fGlauber->BuildPhysicsTable(*particle);
218
219 if(particle == theProton) {
220 isProton = true;
221 fSAIDHighEnergyLimit = 2*GeV;
222 fHighEnergy = 2*GeV;
223 }
224
225 G4ParticleDefinition* part = const_cast<G4ParticleDefinition*>(particle);
226 G4ThreeVector mom(0.0,0.0,1.0);
227 G4DynamicParticle dp(part, mom, fGlauberEnergy);
228
230 G4int A;
231
232 G4double csup, csdn;
233
234 if(verboseLevel > 0) {
235 G4cout << "### G4BGGNucleonInelasticXS::Initialise for "
236 << particle->GetParticleName() << G4endl;
237 }
238 for(G4int iz=2; iz<93; iz++) {
239
240 A = G4lrint(nist->GetAtomicMassAmu(iz));
241 theA[iz] = A;
242
243 csup = fGlauber->GetInelasticGlauberGribov(&dp, iz, A);
244 csdn = fNucleon->GetElementCrossSection(&dp, iz);
245
246 theGlauberFac[iz] = csdn/csup;
247 if(verboseLevel > 0) {
248 G4cout << "Z= " << iz << " A= " << A
249 << " GlauberFactor= " << theGlauberFac[iz] << G4endl;
250 }
251 }
252 //const G4Material* mat = 0;
253
254 dp.SetKineticEnergy(fSAIDHighEnergyLimit);
255 fHadron->GetHadronNucleonXscNS(&dp, theProton);
256 theCoulombFac[0] = fSAIDHighEnergyLimit*
257 (fSAID->GetInelasticIsotopeCrossSection(particle,fSAIDHighEnergyLimit,1,1)
258 /fHadron->GetInelasticHadronNucleonXsc() - 1);
259
260 //G4cout << "Z=1 E(GeV)= " << fSAIDHighEnergyLimit/GeV
261 // << " xsNS(b)= " << fHadron->GetInelasticHadronNucleonXsc()/barn;
262 fHadron->GetHadronNucleonXscPDG(&dp, theProton);
263 //G4cout << " xsPDG(b)= " << fHadron->GetInelasticHadronNucleonXsc()/barn;
264 //G4cout << " xsSAID(b)= " << fSAID->GetInelasticIsotopeCrossSection(particle,fSAIDHighEnergyLimit,1,1)/barn << G4endl;
265
266 dp.SetKineticEnergy(fHighEnergy);
267 fHadron->GetHadronNucleonXscPDG(&dp, theProton);
269
270 //G4cout << "Z=1 E(GeV)= " << fHighEnergy/GeV
271 // << " xsPDG(b)= " << fHadron->GetInelasticHadronNucleonXsc()/barn;
272
273 fHadron->GetHadronNucleonXscNS(&dp, theProton);
274 theCoulombFac[1] = fHighEnergy*((theCoulombFac[0]/fHighEnergy + 1)
275 *fHadron->GetInelasticHadronNucleonXsc()/x - 1);
276
277 fHadron->GetHadronNucleonXscNS(&dp, theProton);
278 //G4cout << " xsNS(b)= " << fHadron->GetInelasticHadronNucleonXsc()/barn << G4endl;
279
280 if(verboseLevel > 0) {
281 G4cout << "Z=1 A=1" << " CoulombFactor[0]= " << theCoulombFac[0]
282 << " CoulombFactor[1]= " << theCoulombFac[1] << G4endl;
283 }
284 theCoulombFac[2] = 1.0;
285
286 dp.SetKineticEnergy(fLowEnergy);
287 for(G4int iz=3; iz<93; iz++) {
288 theCoulombFac[iz] =
289 fNucleon->GetElementCrossSection(&dp, iz)/CoulombFactor(fLowEnergy, iz);
290
291 if(verboseLevel > 0) {
292 G4cout << "Z= " << iz << " A= " << theA[iz]
293 << " CoulombFactor= " << theCoulombFac[iz] << G4endl;
294 }
295 }
296}
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
virtual G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
G4VCrossSectionDataSet * GetCrossSectionDataSet(const G4String &name, G4bool warning=true)
static G4CrossSectionDataSetRegistry * Instance()
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4ParticleDefinition *)
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
G4double GetInelasticHadronNucleonXsc()
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const
virtual G4double GetElementCrossSection(const G4DynamicParticle *aParticle, G4int Z, const G4Material *mat=0)
const G4String & GetParticleName() const
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
int G4lrint(double ad)
Definition: templates.hh:163

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 341 of file G4BGGNucleonInelasticXS.cc.

342{
343 outFile << "The Barashenkov-Glauber-Gribov cross section calculates inelastic\n"
344 << "scattering of protons and neutrons from nuclei using the\n"
345 << "Barashenkov parameterization below 91 GeV and the Glauber-Gribov\n"
346 << "parameterization above 91 GeV. It uses the G4HadronNucleonXsc\n"
347 << "cross section component for hydrogen targets, and the\n"
348 << "G4GlauberGribovCrossSection component for other targets.\n";
349}

◆ GetElementCrossSection()

G4double G4BGGNucleonInelasticXS::GetElementCrossSection ( const G4DynamicParticle dp,
G4int  Z,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 115 of file G4BGGNucleonInelasticXS.cc.

117{
118 // this method should be called only for Z > 1
119
120 G4double cross = 0.0;
121 G4double ekin = dp->GetKineticEnergy();
122 G4int Z = ZZ;
123 if(1 == Z) {
124 cross = 1.0115*GetIsoCrossSection(dp,1,1);
125 } else if(2 == Z) {
126 if(ekin > fGlauberEnergy) {
127 cross = theGlauberFac[Z]*fGlauber->GetInelasticGlauberGribov(dp, Z, theA[Z]);
128 } else {
129 cross = fNucleon->GetElementCrossSection(dp, Z);
130 }
131
132 } else {
133 if(Z > 92) { Z = 92; }
134
135 if(ekin <= fLowEnergy) {
136 cross = theCoulombFac[Z]*CoulombFactor(ekin, Z);
137 } else if(ekin > fGlauberEnergy) {
138 cross = theGlauberFac[Z]*fGlauber->GetInelasticGlauberGribov(dp, Z, theA[Z]);
139 } else {
140 cross = fNucleon->GetElementCrossSection(dp, Z);
141 }
142 }
143
144 if(verboseLevel > 1) {
145 G4cout << "G4BGGNucleonInelasticXS::GetCrossSection for "
147 << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
148 << " in nucleus Z= " << Z << " A= " << theA[Z]
149 << " XS(b)= " << cross/barn
150 << G4endl;
151 }
152 return cross;
153}
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const

◆ GetIsoCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 158 of file G4BGGNucleonInelasticXS.cc.

163{
164 // this method should be called only for Z = 1
165
166 G4double cross = 0.0;
167 G4double ekin = dp->GetKineticEnergy();
168
169 if(ekin <= fSAIDHighEnergyLimit) {
170 cross = fSAID->GetInelasticIsotopeCrossSection(particle, ekin, 1, 1);
171 } else if(ekin < fHighEnergy) {
172 fHadron->GetHadronNucleonXscNS(dp, theProton);
173 cross = (theCoulombFac[0]/ekin + 1)*fHadron->GetInelasticHadronNucleonXsc();
174 } else {
175 fHadron->GetHadronNucleonXscPDG(dp, theProton);
176 cross = (theCoulombFac[1]/ekin + 1)*fHadron->GetInelasticHadronNucleonXsc();
177 }
178 cross *= A;
179
180 if(verboseLevel > 1) {
181 G4cout << "G4BGGNucleonInelasticXS::GetCrossSection for "
183 << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
184 << " in nucleus Z= " << Z << " A= " << A
185 << " XS(b)= " << cross/barn
186 << G4endl;
187 }
188 return cross;
189}

Referenced by GetElementCrossSection().

◆ IsElementApplicable()

G4bool G4BGGNucleonInelasticXS::IsElementApplicable ( const G4DynamicParticle ,
G4int  Z,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Reimplemented in G4NeutronHPBGGNucleonInelasticXS.

Definition at line 96 of file G4BGGNucleonInelasticXS.cc.

98{
99 return true;
100}

◆ IsIsoApplicable()

G4bool G4BGGNucleonInelasticXS::IsIsoApplicable ( const G4DynamicParticle ,
G4int  Z,
G4int  A,
const G4Element elm = 0,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Reimplemented in G4NeutronHPBGGNucleonInelasticXS.

Definition at line 104 of file G4BGGNucleonInelasticXS.cc.

108{
109 return (1 == Z && 2 >= A);
110}

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