G4BGGPionInelasticXS Class Reference

#include <G4BGGPionInelasticXS.hh>

Inheritance diagram for G4BGGPionInelasticXS:

G4VCrossSectionDataSet

Public Member Functions

 G4BGGPionInelasticXS (const G4ParticleDefinition *)
virtual ~G4BGGPionInelasticXS ()
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 68 of file G4BGGPionInelasticXS.hh.


Constructor & Destructor Documentation

G4BGGPionInelasticXS::G4BGGPionInelasticXS ( const G4ParticleDefinition  ) 

Definition at line 57 of file G4BGGPionInelasticXS.cc.

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

00058  : G4VCrossSectionDataSet("Barashenkov-Glauber-Gribov")
00059 {
00060   verboseLevel = 0;
00061   fGlauberEnergy = 91.*GeV;
00062   fLowEnergy = 20.*MeV;
00063   fSAIDHighEnergyLimit = 2.6*GeV;
00064   SetMinKinEnergy(0.0);
00065   SetMaxKinEnergy(100*TeV);
00066 
00067   for (G4int i = 0; i < 93; i++) {
00068     theGlauberFac[i] = 0.0;
00069     theCoulombFac[i] = 0.0;
00070     theA[i] = 1;
00071   }
00072   fPion = 0;
00073   fGlauber = 0;
00074   fHadron  = 0;
00075   //  fGHEISHA = 0;
00076   fSAID    = 0;
00077   particle = p;
00078   theProton= G4Proton::Proton();
00079   isPiplus = false;
00080   isInitialized = false;
00081 }

G4BGGPionInelasticXS::~G4BGGPionInelasticXS (  )  [virtual]

Definition at line 84 of file G4BGGPionInelasticXS.cc.

00085 {
00086   delete fGlauber;
00087   delete fPion;
00088   delete fHadron;
00089   delete fSAID;
00090 }


Member Function Documentation

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 203 of file G4BGGPionInelasticXS.cc.

References G4VCrossSectionDataSet::BuildPhysicsTable(), G4UPiNuclearCrossSection::BuildPhysicsTable(), G4cout, G4endl, G4lrint(), G4NistManager::GetAtomicMassAmu(), G4HadronNucleonXsc::GetHadronNucleonXscPDG(), G4UPiNuclearCrossSection::GetInelasticCrossSection(), G4GlauberGribovCrossSection::GetInelasticGlauberGribov(), G4HadronNucleonXsc::GetInelasticHadronNucleonXsc(), G4ComponentSAIDTotalXS::GetInelasticIsotopeCrossSection(), G4ParticleDefinition::GetParticleName(), G4NistManager::Instance(), G4PionMinus::PionMinus(), G4PionPlus::PionPlus(), G4DynamicParticle::SetKineticEnergy(), and G4VCrossSectionDataSet::verboseLevel.

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

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 333 of file G4BGGPionInelasticXS.cc.

00334 {
00335   outFile << "The Barashenkov-Glauber-Gribov cross section handles inelastic\n"
00336           << "pion scattering from nuclei at all energies.  The Barashenkov\n"
00337           << "parameterization is used below 91 GeV and the Glauber-Gribov\n"
00338           << "parameterization is used above 91 GeV.\n";
00339 }

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 114 of file G4BGGPionInelasticXS.cc.

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

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

G4double G4BGGPionInelasticXS::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 151 of file G4BGGPionInelasticXS.cc.

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

Referenced by GetElementCrossSection().

00156 {
00157   // this method should be called only for Z = 1
00158 
00159   G4double cross = 0.0;
00160   G4double ekin = dp->GetKineticEnergy();
00161 
00162   if(ekin <= fSAIDHighEnergyLimit) {
00163     cross = fSAID->GetInelasticIsotopeCrossSection(particle, ekin, 1, 1);
00164   } else {
00165     fHadron->GetHadronNucleonXscPDG(dp, theProton);
00166     cross = (theCoulombFac[1]/ekin + 1)*fHadron->GetInelasticHadronNucleonXsc();
00167   } 
00168   /*
00169   if(isPiplus) { 
00170     if(ekin <= 20*GeV) {
00171       cross = theCoulombFac[1]*fGHEISHA->GetElementCrossSection(dp, 1, mat);
00172     } else {
00173       fHadron->GetHadronNucleonXscPDG(dp, theProton);
00174       cross = fHadron->GetInelasticHadronNucleonXsc();
00175     }
00176   } else {
00177     if(ekin <= fLowEnergy) {
00178       cross = theCoulombFac[1];
00179     } else if(ekin <= 20*GeV) {
00180       fHadron->GetHadronNucleonXscNS(dp, theProton);
00181       cross = theGlauberFac[1]*fHadron->GetInelasticHadronNucleonXsc();
00182     } else {
00183       fHadron->GetHadronNucleonXscPDG(dp, theProton);
00184       cross = fHadron->GetInelasticHadronNucleonXsc();
00185     }
00186   }
00187   */
00188   cross *= A;
00189 
00190   if(verboseLevel > 1) {
00191     G4cout << "G4BGGPionInelasticXS::GetCrossSection  for "
00192            << dp->GetDefinition()->GetParticleName()
00193            << "  Ekin(GeV)= " << dp->GetKineticEnergy()
00194            << " in nucleus Z= " << Z << "  A= " << A
00195            << " XS(b)= " << cross/barn 
00196            << G4endl;
00197   }
00198   return cross;
00199 }

G4bool G4BGGPionInelasticXS::IsElementApplicable ( const G4DynamicParticle ,
G4int  Z,
const G4Material mat = 0 
) [virtual]

Reimplemented from G4VCrossSectionDataSet.

Definition at line 95 of file G4BGGPionInelasticXS.cc.

00097 {
00098   return true;
00099 }

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 103 of file G4BGGPionInelasticXS.cc.

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


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