G4BGGPionElasticXS Class Reference

#include <G4BGGPionElasticXS.hh>

Inheritance diagram for G4BGGPionElasticXS:

G4VCrossSectionDataSet

Public Member Functions

 G4BGGPionElasticXS (const G4ParticleDefinition *)
virtual ~G4BGGPionElasticXS ()
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *)
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 64 of file G4BGGPionElasticXS.hh.


Constructor & Destructor Documentation

G4BGGPionElasticXS::G4BGGPionElasticXS ( const G4ParticleDefinition  ) 

Definition at line 56 of file G4BGGPionElasticXS.cc.

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

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

G4BGGPionElasticXS::~G4BGGPionElasticXS (  )  [virtual]

Definition at line 83 of file G4BGGPionElasticXS.cc.

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


Member Function Documentation

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 191 of file G4BGGPionElasticXS.cc.

References G4VCrossSectionDataSet::BuildPhysicsTable(), G4UPiNuclearCrossSection::BuildPhysicsTable(), G4cout, G4endl, G4lrint(), G4NistManager::GetAtomicMassAmu(), G4UPiNuclearCrossSection::GetElasticCrossSection(), G4GlauberGribovCrossSection::GetElasticGlauberGribov(), G4HadronNucleonXsc::GetElasticHadronNucleonXsc(), G4ComponentSAIDTotalXS::GetElasticIsotopeCrossSection(), G4HadronNucleonXsc::GetHadronNucleonXscPDG(), G4ParticleDefinition::GetParticleName(), G4NistManager::Instance(), G4PionMinus::PionMinus(), G4PionPlus::PionPlus(), G4DynamicParticle::SetKineticEnergy(), and G4VCrossSectionDataSet::verboseLevel.

00192 {
00193   if(&p == G4PionPlus::PionPlus() || &p == G4PionMinus::PionMinus()) {
00194     particle = &p;
00195   } else {
00196     G4cout << "### G4BGGPionElasticXS WARNING: is not applicable to " 
00197            << p.GetParticleName()
00198            << G4endl;
00199     throw G4HadronicException(__FILE__, __LINE__,
00200           "G4BGGPionElasticXS::BuildPhysicsTable is used for wrong particle");
00201     return;
00202   }
00203 
00204   if(isInitialized) { return; }
00205   isInitialized = true;
00206 
00207   fPion    = new G4UPiNuclearCrossSection();
00208   fGlauber = new G4GlauberGribovCrossSection();
00209   fHadron  = new G4HadronNucleonXsc();
00210   fSAID    = new G4ComponentSAIDTotalXS();
00211   fPion->BuildPhysicsTable(*particle);
00212   fGlauber->BuildPhysicsTable(*particle);
00213   if(particle == G4PionPlus::PionPlus()) { isPiplus = true; }
00214 
00215   G4ParticleDefinition* part = const_cast<G4ParticleDefinition*>(particle);
00216   G4ThreeVector mom(0.0,0.0,1.0);
00217   G4DynamicParticle dp(part, mom, fGlauberEnergy);
00218 
00219   G4NistManager* nist = G4NistManager::Instance();
00220 
00221   G4double csup, csdn;
00222   G4int A;
00223 
00224   if(verboseLevel > 0) {
00225     G4cout << "### G4BGGPionElasticXS::Initialise for "
00226            << particle->GetParticleName() << G4endl;
00227   }
00228   for(G4int iz=2; iz<93; iz++) {
00229 
00230     A = G4lrint(nist->GetAtomicMassAmu(iz));
00231     theA[iz] = A;
00232 
00233     csup = fGlauber->GetElasticGlauberGribov(&dp, iz, A);
00234     csdn = fPion->GetElasticCrossSection(&dp, iz, A);
00235 
00236     theGlauberFac[iz] = csdn/csup;
00237     if(verboseLevel > 0) {
00238       G4cout << "Z= " << iz <<  "  A= " << A 
00239              << " factor= " << theGlauberFac[iz] << G4endl; 
00240     }
00241   }
00242   /*
00243   dp.SetKineticEnergy(fLowEnergy);
00244   fHadron->GetHadronNucleonXscNS(&dp, G4Proton::Proton());
00245   theCoulombFac[1] = fHadron->GetElasticHadronNucleonXsc();
00246   */
00247   dp.SetKineticEnergy(fSAIDHighEnergyLimit);
00248   fHadron->GetHadronNucleonXscPDG(&dp, theProton);
00249   theCoulombFac[1] = 
00250     fSAID->GetElasticIsotopeCrossSection(particle,fSAIDHighEnergyLimit,1,1)
00251     /fHadron->GetElasticHadronNucleonXsc();
00252 
00253   dp.SetKineticEnergy(fLowEnergy);
00254   for(G4int iz=2; iz<93; iz++) {
00255     theCoulombFac[iz] = fPion->GetElasticCrossSection(&dp, iz, theA[iz]);
00256     if(verboseLevel > 0) {
00257       G4cout << "Z= " << iz <<  "  A= " << A 
00258              << " factor= " << theCoulombFac[iz] << G4endl; 
00259     }
00260   }
00261 }

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 264 of file G4BGGPionElasticXS.cc.

00265 {
00266   outFile << "The Barashenkov-Glauber-Gribov cross section handles elastic\n"
00267           << "scattering of pions from nuclei at all energies. The\n"
00268           << "Barashenkov parameterization is used below 91 GeV and the\n"
00269           << "Glauber-Gribov parameterization is used above 91 GeV.\n";
00270 }

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 113 of file G4BGGPionElasticXS.cc.

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

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

G4double G4BGGPionElasticXS::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 146 of file G4BGGPionElasticXS.cc.

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

Referenced by GetElementCrossSection().

00151 {
00152   // this method should be called only for Z = 1
00153 
00154   G4double cross = 0.0;
00155   G4double ekin = dp->GetKineticEnergy();
00156 
00157   if(ekin <= fSAIDHighEnergyLimit) {
00158     cross = fSAID->GetElasticIsotopeCrossSection(particle, ekin, 1, 1);
00159   } else {
00160     fHadron->GetHadronNucleonXscPDG(dp, theProton);
00161     cross = theCoulombFac[1]*fHadron->GetElasticHadronNucleonXsc();
00162   } 
00163   cross *= A;
00164   /*
00165   if(ekin <= fLowEnergy) {
00166     cross = theCoulombFac[1];
00167 
00168   } else if( A < 2) {
00169     fHadron->GetHadronNucleonXscNS(dp, G4Proton::Proton());
00170     cross = fHadron->GetElasticHadronNucleonXsc();
00171   } else {
00172     fHadron->GetHadronNucleonXscNS(dp, G4Proton::Proton());
00173     cross = fHadron->GetElasticHadronNucleonXsc();
00174     fHadron->GetHadronNucleonXscNS(dp, G4Neutron::Neutron());
00175     cross += fHadron->GetElasticHadronNucleonXsc();
00176   }
00177   */
00178   if(verboseLevel > 1) {
00179     G4cout << "G4BGGPionElasticXS::GetIsoCrossSection  for "
00180            << dp->GetDefinition()->GetParticleName()
00181            << "  Ekin(GeV)= " << dp->GetKineticEnergy()
00182            << " in nucleus Z= " << Z << "  A= " << A
00183            << " XS(b)= " << cross/barn 
00184            << G4endl;
00185   }
00186   return cross;
00187 }

G4bool G4BGGPionElasticXS::IsElementApplicable ( const G4DynamicParticle ,
G4int  Z,
const G4Material  
) [virtual]

Reimplemented from G4VCrossSectionDataSet.

Definition at line 94 of file G4BGGPionElasticXS.cc.

00096 {
00097   return true;
00098 }

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 102 of file G4BGGPionElasticXS.cc.

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


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