Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4BGGNucleonInelasticXS Class Reference

#include <G4BGGNucleonInelasticXS.hh>

Inheritance diagram for G4BGGNucleonInelasticXS:
G4VCrossSectionDataSet G4NeutronHPBGGNucleonInelasticXS

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
 
void SetLowestCrossSection (G4double val)
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
virtual G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual G4int GetVerboseLevel () const
 
virtual 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 64 of file G4BGGNucleonInelasticXS.hh.

Constructor & Destructor Documentation

G4BGGNucleonInelasticXS::G4BGGNucleonInelasticXS ( const G4ParticleDefinition p)

Definition at line 65 of file G4BGGNucleonInelasticXS.cc.

References python.hepunit::GeV, python.hepunit::MeV, python.hepunit::millibarn, G4Proton::Proton(), and G4VCrossSectionDataSet::verboseLevel.

66  : G4VCrossSectionDataSet("Barashenkov-Glauber")
67 {
68  verboseLevel = 0;
69  fGlauberEnergy = 91.*GeV;
70  fLowEnergy = 14.*MeV;
71  fHighEnergy = 5.*GeV;
72  fSAIDHighEnergyLimit = 1.3*GeV;
73  fLowestXSection = millibarn;
74  for (G4int i = 0; i < 93; ++i) {
75  theGlauberFac[i] = 0.0;
76  theCoulombFac[i] = 0.0;
77  theA[i] = 1;
78  }
79  fNucleon = 0;
80  fGlauber = 0;
81  fHadron = 0;
82  fSAID = 0;
83 
84  particle = p;
85  theProton= G4Proton::Proton();
86  isProton = false;
87  isInitialized = false;
88 }
G4VCrossSectionDataSet(const G4String &nam="")
const char * p
Definition: xmltok.h:285
int G4int
Definition: G4Types.hh:78
int millibarn
Definition: hepunit.py:40
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4bool isInitialized()
G4BGGNucleonInelasticXS::~G4BGGNucleonInelasticXS ( )
virtual

Definition at line 92 of file G4BGGNucleonInelasticXS.cc.

93 {
94  delete fHadron;
95  delete fSAID;
96 }

Member Function Documentation

void G4BGGNucleonInelasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 199 of file G4BGGNucleonInelasticXS.cc.

References G4VCrossSectionDataSet::BuildPhysicsTable(), G4NucleonNuclearCrossSection::Default_Name(), G4GlauberGribovCrossSection::Default_Name(), G4cout, G4endl, G4lrint(), G4NistManager::GetAtomicMassAmu(), G4CrossSectionDataSetRegistry::GetCrossSectionDataSet(), G4NucleonNuclearCrossSection::GetElementCrossSection(), G4HadronNucleonXsc::GetHadronNucleonXscNS(), G4HadronNucleonXsc::GetHadronNucleonXscPDG(), G4GlauberGribovCrossSection::GetInelasticGlauberGribov(), G4HadronNucleonXsc::GetInelasticHadronNucleonXsc(), G4ComponentSAIDTotalXS::GetInelasticIsotopeCrossSection(), G4ParticleDefinition::GetParticleName(), python.hepunit::GeV, G4CrossSectionDataSetRegistry::Instance(), G4NistManager::Instance(), iz, G4Neutron::Neutron(), G4DynamicParticle::SetKineticEnergy(), G4VCrossSectionDataSet::verboseLevel, and test::x.

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 346 of file G4BGGNucleonInelasticXS.cc.

347 {
348  outFile << "The Barashenkov-Glauber-Gribov cross section calculates inelastic\n"
349  << "scattering of protons and neutrons from nuclei using the\n"
350  << "Barashenkov parameterization below 91 GeV and the Glauber-Gribov\n"
351  << "parameterization above 91 GeV. It uses the G4HadronNucleonXsc\n"
352  << "cross section component for hydrogen targets, and the\n"
353  << "G4GlauberGribovCrossSection component for other targets.\n";
354 }
std::ofstream outFile
Definition: GammaRayTel.cc:68
G4double G4BGGNucleonInelasticXS::GetElementCrossSection ( const G4DynamicParticle dp,
G4int  Z,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 119 of file G4BGGNucleonInelasticXS.cc.

References python.hepunit::barn, G4cout, G4endl, G4DynamicParticle::GetDefinition(), G4NucleonNuclearCrossSection::GetElementCrossSection(), G4GlauberGribovCrossSection::GetInelasticGlauberGribov(), GetIsoCrossSection(), G4DynamicParticle::GetKineticEnergy(), G4ParticleDefinition::GetParticleName(), and G4VCrossSectionDataSet::verboseLevel.

121 {
122  // this method should be called only for Z > 1
123 
124  G4double cross = 0.0;
125  G4double ekin = dp->GetKineticEnergy();
126  G4int Z = ZZ;
127  if(1 == Z) {
128  cross = 1.0115*GetIsoCrossSection(dp,1,1);
129  } else if(2 == Z) {
130  if(ekin > fGlauberEnergy) {
131  cross = theGlauberFac[Z]*fGlauber->GetInelasticGlauberGribov(dp, Z, theA[Z]);
132  } else {
133  cross = fNucleon->GetElementCrossSection(dp, Z);
134  }
135 
136  } else {
137  if(Z > 92) { Z = 92; }
138 
139  if(ekin <= fLowEnergy) {
140  cross = theCoulombFac[Z]*CoulombFactor(ekin, Z);
141  } else if(ekin > fGlauberEnergy) {
142  cross = theGlauberFac[Z]*fGlauber->GetInelasticGlauberGribov(dp, Z, theA[Z]);
143  } else {
144  cross = fNucleon->GetElementCrossSection(dp, Z);
145  }
146  }
147 
148  if(verboseLevel > 1) {
149  G4cout << "G4BGGNucleonInelasticXS::GetCrossSection for "
150  << dp->GetDefinition()->GetParticleName()
151  << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
152  << " in nucleus Z= " << Z << " A= " << theA[Z]
153  << " XS(b)= " << cross/barn
154  << G4endl;
155  }
156  //AR-18Dec2013 if(cross <= fLowestXSection) { cross = 0.0; }
157  return cross;
158 }
G4double GetKineticEnergy() const
virtual G4double GetElementCrossSection(const G4DynamicParticle *aParticle, G4int Z, const G4Material *mat=0)
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
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 163 of file G4BGGNucleonInelasticXS.cc.

References python.hepunit::barn, G4cout, G4endl, G4DynamicParticle::GetDefinition(), G4HadronNucleonXsc::GetHadronNucleonXscNS(), G4HadronNucleonXsc::GetHadronNucleonXscPDG(), G4HadronNucleonXsc::GetInelasticHadronNucleonXsc(), G4ComponentSAIDTotalXS::GetInelasticIsotopeCrossSection(), G4DynamicParticle::GetKineticEnergy(), G4ParticleDefinition::GetParticleName(), and G4VCrossSectionDataSet::verboseLevel.

Referenced by GetElementCrossSection().

168 {
169  // this method should be called only for Z = 1
170 
171  G4double cross = 0.0;
172  G4double ekin = dp->GetKineticEnergy();
173 
174  if(ekin <= fSAIDHighEnergyLimit) {
175  cross = fSAID->GetInelasticIsotopeCrossSection(particle, ekin, 1, 1);
176  } else if(ekin < fHighEnergy) {
177  fHadron->GetHadronNucleonXscNS(dp, theProton);
178  cross = (theCoulombFac[0]/ekin + 1)*fHadron->GetInelasticHadronNucleonXsc();
179  } else {
180  fHadron->GetHadronNucleonXscPDG(dp, theProton);
181  cross = (theCoulombFac[1]/ekin + 1)*fHadron->GetInelasticHadronNucleonXsc();
182  }
183  cross *= A;
184 
185  if(verboseLevel > 1) {
186  G4cout << "G4BGGNucleonInelasticXS::GetCrossSection for "
187  << dp->GetDefinition()->GetParticleName()
188  << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
189  << " in nucleus Z= " << Z << " A= " << A
190  << " XS(b)= " << cross/barn
191  << G4endl;
192  }
193  //AR-18Dec2013 if(cross <= fLowestXSection) { cross = 0.0; }
194  return cross;
195 }
virtual G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
G4double GetKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
const G4String & GetParticleName() const
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4ParticleDefinition *)
G4GLOB_DLL std::ostream G4cout
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double GetInelasticHadronNucleonXsc()
G4bool G4BGGNucleonInelasticXS::IsElementApplicable ( const G4DynamicParticle ,
G4int  Z,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Reimplemented in G4NeutronHPBGGNucleonInelasticXS.

Definition at line 100 of file G4BGGNucleonInelasticXS.cc.

102 {
103  return true;
104 }
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 108 of file G4BGGNucleonInelasticXS.cc.

112 {
113  return (1 == Z);
114 }
void G4BGGNucleonInelasticXS::SetLowestCrossSection ( G4double  val)
inline

Definition at line 125 of file G4BGGNucleonInelasticXS.hh.

126 {
127  fLowestXSection = val;
128 }

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