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

#include <G4BGGNucleonElasticXS.hh>

Inheritance diagram for G4BGGNucleonElasticXS:
G4VCrossSectionDataSet

Public Member Functions

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

Constructor & Destructor Documentation

G4BGGNucleonElasticXS::G4BGGNucleonElasticXS ( const G4ParticleDefinition p)

Definition at line 60 of file G4BGGNucleonElasticXS.cc.

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

61  : G4VCrossSectionDataSet("Barashenkov-Glauber")
62 {
63  verboseLevel = 0;
64  fGlauberEnergy = 91.*GeV;
65  fPDGEnergy = 5*GeV;
66  fLowEnergy = 14.*MeV;
67  fSAIDLowEnergyLimit = 1*MeV;
68  fSAIDHighEnergyLimit = 1.3*GeV;
69  fLowestXSection = millibarn;
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  fSAID = 0;
79  particle = p;
80  theProton= G4Proton::Proton();
81  isProton = false;
82  isInitialized = false;
83 }
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()
G4BGGNucleonElasticXS::~G4BGGNucleonElasticXS ( )
virtual

Definition at line 87 of file G4BGGNucleonElasticXS.cc.

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

Member Function Documentation

void G4BGGNucleonElasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 194 of file G4BGGNucleonElasticXS.cc.

References G4VCrossSectionDataSet::BuildPhysicsTable(), G4NucleonNuclearCrossSection::Default_Name(), G4GlauberGribovCrossSection::Default_Name(), G4cout, G4endl, G4lrint(), G4NistManager::GetAtomicMassAmu(), G4CrossSectionDataSetRegistry::GetCrossSectionDataSet(), G4NucleonNuclearCrossSection::GetElasticCrossSection(), G4GlauberGribovCrossSection::GetElasticGlauberGribov(), G4HadronNucleonXsc::GetElasticHadronNucleonXsc(), G4ComponentSAIDTotalXS::GetElasticIsotopeCrossSection(), G4HadronNucleonXsc::GetHadronNucleonXscPDG(), G4ParticleDefinition::GetParticleName(), python.hepunit::GeV, G4CrossSectionDataSetRegistry::Instance(), G4NistManager::Instance(), iz, G4Neutron::Neutron(), G4DynamicParticle::SetKineticEnergy(), and G4VCrossSectionDataSet::verboseLevel.

195 {
196  if(&p == theProton || &p == G4Neutron::Neutron()) {
197  particle = &p;
198 
199  } else {
200  G4cout << "### G4BGGNucleonElasticXS WARNING: is not applicable to "
201  << p.GetParticleName()
202  << G4endl;
203  throw G4HadronicException(__FILE__, __LINE__,
204  "G4BGGNucleonElasticXS::BuildPhysicsTable is used for wrong particle");
205  return;
206  }
207 
208  if(isInitialized) { return; }
209  isInitialized = true;
210 
213 
214  fHadron = new G4HadronNucleonXsc();
215  fSAID = new G4ComponentSAIDTotalXS();
216  fNucleon->BuildPhysicsTable(*particle);
217  fGlauber->BuildPhysicsTable(*particle);
218  if(particle == theProton) {
219  isProton = true;
220  fSAIDHighEnergyLimit = 3*GeV;
221  }
222 
223  G4ParticleDefinition* part = const_cast<G4ParticleDefinition*>(particle);
224  G4ThreeVector mom(0.0,0.0,1.0);
225  G4DynamicParticle dp(part, mom, fGlauberEnergy);
226 
228 
229  G4double csup, csdn;
230  G4int A;
231 
232  if(verboseLevel > 0) {
233  G4cout << "### G4BGGNucleonElasticXS::Initialise for "
234  << particle->GetParticleName() << G4endl;
235  }
236 
237  for(G4int iz=2; iz<93; iz++) {
238 
239  A = G4lrint(nist->GetAtomicMassAmu(iz));
240  theA[iz] = A;
241 
242  csup = fGlauber->GetElasticGlauberGribov(&dp, iz, A);
243  csdn = fNucleon->GetElasticCrossSection(&dp, iz);
244 
245  theGlauberFac[iz] = csdn/csup;
246  if(verboseLevel > 0) {
247  G4cout << "Z= " << iz << " A= " << A
248  << " factor= " << theGlauberFac[iz] << G4endl;
249  }
250  }
251 
252  theCoulombFac[0] =
253  fSAID->GetElasticIsotopeCrossSection(particle,fSAIDLowEnergyLimit,1,1)
254  /CoulombFactor(fSAIDLowEnergyLimit, 1);
255 
256  dp.SetKineticEnergy(fPDGEnergy);
257  fHadron->GetHadronNucleonXscPDG(&dp, theProton);
258  theCoulombFac[1] = fHadron->GetElasticHadronNucleonXsc();
259 
260  if(verboseLevel > 0) {
261  G4cout << "Z=1 A=1 " << " factor0= " << theCoulombFac[0]
262  << " factor1= " << theCoulombFac[1]
263  << G4endl;
264  }
265 
266  dp.SetKineticEnergy(fLowEnergy);
267  for(G4int iz=2; iz<93; iz++) {
268  theCoulombFac[iz] =
269  fNucleon->GetElasticCrossSection(&dp, iz)/CoulombFactor(fLowEnergy, iz);
270  if(verboseLevel > 0) {
271  G4cout << "Z= " << iz << " A= " << theA[iz]
272  << " factor= " << theCoulombFac[iz] << G4endl;
273  }
274  }
275 }
G4VCrossSectionDataSet * GetCrossSectionDataSet(const G4String &name, G4bool warning=true)
G4double GetElasticHadronNucleonXsc()
const char * p
Definition: xmltok.h:285
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual G4double GetElasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
const G4String & GetParticleName() const
G4double GetElasticCrossSection(const G4DynamicParticle *aParticle, G4int Z)
G4GLOB_DLL std::ostream G4cout
G4double iz
Definition: TRTMaterials.hh:39
static G4CrossSectionDataSetRegistry * Instance()
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
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 GetElasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
void G4BGGNucleonElasticXS::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 306 of file G4BGGNucleonElasticXS.cc.

307 {
308  outFile << "The Barashenkov-Glauber-Gribov cross section handles elastic\n"
309  << "scattering of protons and neutrons from nuclei using the\n"
310  << "Barashenkov parameterization below 91 GeV and the Glauber-Gribov\n"
311  << "parameterization above 91 GeV. n";
312 }
std::ofstream outFile
Definition: GammaRayTel.cc:68
G4double G4BGGNucleonElasticXS::GetElementCrossSection ( const G4DynamicParticle dp,
G4int  Z,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 115 of file G4BGGNucleonElasticXS.cc.

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

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 {
126  if(Z > 92) { Z = 92; }
127 
128  if(ekin <= fLowEnergy) {
129  cross = theCoulombFac[Z]*CoulombFactor(ekin, Z);
130 
131  } else if(ekin > fGlauberEnergy) {
132  cross = theGlauberFac[Z]*fGlauber->GetElasticGlauberGribov(dp, Z, theA[Z]);
133  } else {
134  cross = fNucleon->GetElasticCrossSection(dp, Z);
135  }
136  }
137 
138  if(verboseLevel > 1) {
139  G4cout << "G4BGGNucleonElasticXS::GetElementCrossSection for "
140  << dp->GetDefinition()->GetParticleName()
141  << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
142  << " in nucleus Z= " << Z << " A= " << theA[Z]
143  << " XS(b)= " << cross/barn
144  << G4endl;
145  }
146  //AR-04Dec2013 if(cross <= fLowestXSection) { cross = 0.0; }
147  return cross;
148 }
G4double GetKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
const G4String & GetParticleName() const
G4double GetElasticCrossSection(const G4DynamicParticle *aParticle, G4int Z)
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double GetElasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
G4double G4BGGNucleonElasticXS::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 153 of file G4BGGNucleonElasticXS.cc.

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

Referenced by GetElementCrossSection().

158 {
159  // this method should be called only for Z = 1
160 
161  G4double cross = 0.0;
162  G4double ekin = dp->GetKineticEnergy();
163 
164  if(ekin <= fSAIDLowEnergyLimit) {
165  cross = theCoulombFac[0]*CoulombFactor(ekin, 1);
166  } else if(ekin <= fSAIDHighEnergyLimit) {
167  cross = fSAID->GetElasticIsotopeCrossSection(particle, ekin, 1, 1);
168  } else if(ekin <= fPDGEnergy) {
169  G4double cross1 =
170  fSAID->GetElasticIsotopeCrossSection(particle, fPDGEnergy, 1, 1);
171  G4double cross2 = theCoulombFac[1];
172  cross = cross1 + (cross2 - cross1)*(ekin - fSAIDHighEnergyLimit)
173  /(fPDGEnergy - fSAIDHighEnergyLimit);
174  } else {
175  fHadron->GetHadronNucleonXscPDG(dp, theProton);
176  cross = fHadron->GetElasticHadronNucleonXsc();
177  }
178  cross *= A;
179 
180  if(verboseLevel > 1) {
181  G4cout << "G4BGGNucleonElasticXS::GetIsoCrossSection for "
182  << dp->GetDefinition()->GetParticleName()
183  << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
184  << " in nucleus Z= " << Z << " A= " << A
185  << " XS(b)= " << cross/barn
186  << G4endl;
187  }
188  //AR-04Dec2013 if(cross <= fLowestXSection) { cross = 0.0; }
189  return cross;
190 }
G4double GetElasticHadronNucleonXsc()
G4double GetKineticEnergy() const
virtual G4double GetElasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
G4ParticleDefinition * GetDefinition() const
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4bool G4BGGNucleonElasticXS::IsElementApplicable ( const G4DynamicParticle ,
G4int  Z,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 96 of file G4BGGNucleonElasticXS.cc.

98 {
99  return true;
100 }
G4bool G4BGGNucleonElasticXS::IsIsoApplicable ( const G4DynamicParticle ,
G4int  Z,
G4int  A,
const G4Element elm = 0,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 104 of file G4BGGNucleonElasticXS.cc.

108 {
109  return (1 == Z);
110 }
void G4BGGNucleonElasticXS::SetLowestCrossSection ( G4double  val)
inline

Definition at line 126 of file G4BGGNucleonElasticXS.hh.

127 {
128  fLowestXSection = val;
129 }

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