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

Detailed Description

Definition at line 65 of file G4BGGNucleonInelasticXS.hh.


Constructor & Destructor Documentation

G4BGGNucleonInelasticXS::G4BGGNucleonInelasticXS ( const G4ParticleDefinition  ) 

Definition at line 62 of file G4BGGNucleonInelasticXS.cc.

References G4Proton::Proton(), and G4VCrossSectionDataSet::verboseLevel.

00063  : G4VCrossSectionDataSet("Barashenkov-Glauber")
00064 {
00065   verboseLevel = 0;
00066   fGlauberEnergy = 91.*GeV;
00067   fLowEnergy = 14.*MeV;
00068   fHighEnergy = 5.*GeV;
00069   fSAIDHighEnergyLimit = 1.3*GeV;
00070   for (G4int i = 0; i < 93; ++i) {
00071     theGlauberFac[i] = 0.0;
00072     theCoulombFac[i] = 0.0;
00073     theA[i] = 1;
00074   }
00075   fNucleon = 0;
00076   fGlauber = 0;
00077   fHadron  = 0;
00078   //  fGHEISHA = 0;
00079   fSAID    = 0;
00080   particle = p;
00081   theProton= G4Proton::Proton();
00082   isProton = false;
00083   isInitialized = false;
00084 }

G4BGGNucleonInelasticXS::~G4BGGNucleonInelasticXS (  )  [virtual]

Definition at line 88 of file G4BGGNucleonInelasticXS.cc.

00089 {
00090   delete fHadron;
00091   delete fSAID;
00092 }


Member Function Documentation

void G4BGGNucleonInelasticXS::BuildPhysicsTable ( const G4ParticleDefinition  )  [virtual]

Reimplemented from G4VCrossSectionDataSet.

Definition at line 193 of file G4BGGNucleonInelasticXS.cc.

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

00194 {
00195   if(&p == theProton || &p == G4Neutron::Neutron()) {
00196     particle = &p;
00197   } else {
00198     G4cout << "### G4BGGNucleonInelasticXS WARNING: is not applicable to " 
00199            << p.GetParticleName()
00200            << G4endl;
00201     throw G4HadronicException(__FILE__, __LINE__,
00202           "G4BGGNucleonElasticXS::BuildPhysicsTable is used for wrong particle");
00203     return;
00204   }
00205 
00206   if(isInitialized) { return; }
00207   isInitialized = true;
00208 
00209     fNucleon = (G4NucleonNuclearCrossSection*)G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4NucleonNuclearCrossSection::Default_Name());
00210     fGlauber = (G4GlauberGribovCrossSection*)G4CrossSectionDataSetRegistry::Instance()->GetCrossSectionDataSet(G4GlauberGribovCrossSection::Default_Name());
00211 
00212     fHadron  = new G4HadronNucleonXsc();
00213   //fGHEISHA = new G4HadronInelasticDataSet();
00214   fSAID    = new G4ComponentSAIDTotalXS();
00215 
00216   fNucleon->BuildPhysicsTable(*particle);
00217   fGlauber->BuildPhysicsTable(*particle);
00218 
00219   if(particle == theProton) { 
00220     isProton = true; 
00221     fSAIDHighEnergyLimit = 2*GeV;
00222     fHighEnergy = 2*GeV;
00223   }
00224 
00225   G4ParticleDefinition* part = const_cast<G4ParticleDefinition*>(particle);
00226   G4ThreeVector mom(0.0,0.0,1.0);
00227   G4DynamicParticle dp(part, mom, fGlauberEnergy);
00228 
00229   G4NistManager* nist = G4NistManager::Instance();
00230   G4int A;
00231 
00232   G4double csup, csdn;
00233 
00234   if(verboseLevel > 0) {
00235     G4cout << "### G4BGGNucleonInelasticXS::Initialise for "
00236            << particle->GetParticleName() << G4endl;
00237   }
00238   for(G4int iz=2; iz<93; iz++) {
00239 
00240     A = G4lrint(nist->GetAtomicMassAmu(iz));
00241     theA[iz] = A;
00242 
00243     csup = fGlauber->GetInelasticGlauberGribov(&dp, iz, A);
00244     csdn = fNucleon->GetElementCrossSection(&dp, iz);
00245 
00246     theGlauberFac[iz] = csdn/csup;
00247     if(verboseLevel > 0) {
00248       G4cout << "Z= " << iz <<  "  A= " << A 
00249              << " GlauberFactor= " << theGlauberFac[iz] << G4endl; 
00250     }
00251   }
00252   //const G4Material* mat = 0;
00253 
00254   dp.SetKineticEnergy(fSAIDHighEnergyLimit);
00255   fHadron->GetHadronNucleonXscNS(&dp, theProton);
00256   theCoulombFac[0] = fSAIDHighEnergyLimit*
00257     (fSAID->GetInelasticIsotopeCrossSection(particle,fSAIDHighEnergyLimit,1,1)
00258      /fHadron->GetInelasticHadronNucleonXsc() - 1);
00259   
00260   //G4cout << "Z=1 E(GeV)= " << fSAIDHighEnergyLimit/GeV
00261   //     << "  xsNS(b)= " << fHadron->GetInelasticHadronNucleonXsc()/barn;  
00262   fHadron->GetHadronNucleonXscPDG(&dp, theProton);
00263   //G4cout << "  xsPDG(b)= " << fHadron->GetInelasticHadronNucleonXsc()/barn;
00264   //G4cout << "  xsSAID(b)= " << fSAID->GetInelasticIsotopeCrossSection(particle,fSAIDHighEnergyLimit,1,1)/barn << G4endl;
00265 
00266   dp.SetKineticEnergy(fHighEnergy);
00267   fHadron->GetHadronNucleonXscPDG(&dp, theProton);
00268   G4double x = fHadron->GetInelasticHadronNucleonXsc();
00269 
00270   //G4cout << "Z=1 E(GeV)= " << fHighEnergy/GeV
00271   //     << "  xsPDG(b)= " << fHadron->GetInelasticHadronNucleonXsc()/barn;  
00272 
00273   fHadron->GetHadronNucleonXscNS(&dp, theProton);
00274   theCoulombFac[1] = fHighEnergy*((theCoulombFac[0]/fHighEnergy + 1)
00275                                   *fHadron->GetInelasticHadronNucleonXsc()/x - 1);
00276 
00277   fHadron->GetHadronNucleonXscNS(&dp, theProton);
00278   //G4cout << "  xsNS(b)= " << fHadron->GetInelasticHadronNucleonXsc()/barn << G4endl;
00279 
00280   if(verboseLevel > 0) {
00281     G4cout << "Z=1   A=1" << " CoulombFactor[0]= " << theCoulombFac[0]
00282            << " CoulombFactor[1]= " << theCoulombFac[1] << G4endl; 
00283   }
00284   theCoulombFac[2] = 1.0;
00285      
00286   dp.SetKineticEnergy(fLowEnergy);
00287   for(G4int iz=3; iz<93; iz++) {
00288     theCoulombFac[iz] = 
00289       fNucleon->GetElementCrossSection(&dp, iz)/CoulombFactor(fLowEnergy, iz);
00290 
00291     if(verboseLevel > 0) {
00292       G4cout << "Z= " << iz <<  "  A= " << theA[iz] 
00293              << " CoulombFactor= " << theCoulombFac[iz] << G4endl; 
00294     }
00295   }
00296 }

void G4BGGNucleonInelasticXS::CrossSectionDescription ( std::ostream &   )  const [virtual]

Reimplemented from G4VCrossSectionDataSet.

Definition at line 341 of file G4BGGNucleonInelasticXS.cc.

00342 {
00343   outFile << "The Barashenkov-Glauber-Gribov cross section calculates inelastic\n"
00344           << "scattering of protons and neutrons from nuclei using the\n"
00345           << "Barashenkov parameterization below 91 GeV and the Glauber-Gribov\n"
00346           << "parameterization above 91 GeV.  It uses the G4HadronNucleonXsc\n"
00347           << "cross section component for hydrogen targets, and the\n"
00348           << "G4GlauberGribovCrossSection component for other targets.\n";
00349 }

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 115 of file G4BGGNucleonInelasticXS.cc.

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

00117 {
00118   // this method should be called only for Z > 1
00119 
00120   G4double cross = 0.0;
00121   G4double ekin = dp->GetKineticEnergy();
00122   G4int Z = ZZ;
00123   if(1 == Z) {
00124     cross = 1.0115*GetIsoCrossSection(dp,1,1);
00125   } else if(2 == Z) {
00126     if(ekin > fGlauberEnergy) {
00127       cross = theGlauberFac[Z]*fGlauber->GetInelasticGlauberGribov(dp, Z, theA[Z]);
00128     } else {
00129       cross = fNucleon->GetElementCrossSection(dp, Z);
00130     }
00131 
00132   } else {
00133     if(Z > 92) { Z = 92; }
00134 
00135     if(ekin <= fLowEnergy) {
00136       cross = theCoulombFac[Z]*CoulombFactor(ekin, Z);
00137     } else if(ekin > fGlauberEnergy) {
00138       cross = theGlauberFac[Z]*fGlauber->GetInelasticGlauberGribov(dp, Z, theA[Z]);
00139     } else {
00140       cross = fNucleon->GetElementCrossSection(dp, Z);
00141     }
00142   }
00143 
00144   if(verboseLevel > 1) {
00145     G4cout << "G4BGGNucleonInelasticXS::GetCrossSection  for "
00146            << dp->GetDefinition()->GetParticleName()
00147            << "  Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
00148            << " in nucleus Z= " << Z << "  A= " << theA[Z]
00149            << " XS(b)= " << cross/barn 
00150            << G4endl;
00151   }
00152   return cross;
00153 }

G4double G4BGGNucleonInelasticXS::GetIsoCrossSection ( const G4DynamicParticle ,
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.

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

Referenced by GetElementCrossSection().

00163 {
00164   // this method should be called only for Z = 1
00165 
00166   G4double cross = 0.0;
00167   G4double ekin = dp->GetKineticEnergy();
00168 
00169   if(ekin <= fSAIDHighEnergyLimit) {
00170     cross = fSAID->GetInelasticIsotopeCrossSection(particle, ekin, 1, 1);
00171   } else if(ekin < fHighEnergy) {
00172     fHadron->GetHadronNucleonXscNS(dp, theProton);
00173     cross = (theCoulombFac[0]/ekin + 1)*fHadron->GetInelasticHadronNucleonXsc();
00174   } else {
00175     fHadron->GetHadronNucleonXscPDG(dp, theProton);
00176     cross = (theCoulombFac[1]/ekin + 1)*fHadron->GetInelasticHadronNucleonXsc();
00177   } 
00178   cross *= A; 
00179 
00180   if(verboseLevel > 1) {
00181     G4cout << "G4BGGNucleonInelasticXS::GetCrossSection  for "
00182            << dp->GetDefinition()->GetParticleName()
00183            << "  Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
00184            << " in nucleus Z= " << Z << "  A= " << A
00185            << " XS(b)= " << cross/barn 
00186            << G4endl;
00187   }
00188   return cross;
00189 }

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.

00098 {
00099   return true;
00100 }

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.

00108 {
00109   return (1 == Z && 2 >= A);
00110 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:51:32 2013 for Geant4 by  doxygen 1.4.7