Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4BGGPionInelasticXS Class Reference

#include <G4BGGPionInelasticXS.hh>

Inheritance diagram for G4BGGPionInelasticXS:
Collaboration diagram for G4BGGPionInelasticXS:

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
 
- 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 68 of file G4BGGPionInelasticXS.hh.

Constructor & Destructor Documentation

G4BGGPionInelasticXS::G4BGGPionInelasticXS ( const G4ParticleDefinition p)

Definition at line 56 of file G4BGGPionInelasticXS.cc.

57  : G4VCrossSectionDataSet("Barashenkov-Glauber-Gribov")
58 {
59  verboseLevel = 0;
60  fGlauberEnergy = 91.*GeV;
61  fLowEnergy = 20.*MeV;
62  fSAIDHighEnergyLimit = 2.6*GeV;
63  SetMinKinEnergy(0.0);
64  SetMaxKinEnergy(100*TeV);
65 
66  for (G4int i = 0; i < 93; i++) {
67  theGlauberFac[i] = 0.0;
68  theCoulombFac[i] = 0.0;
69  theA[i] = 1;
70  }
71  fPion = 0;
72  fGlauber = 0;
73  fHadron = 0;
74  fSAID = 0;
75 
76  fG4pow = G4Pow::GetInstance();
77 
78  particle = p;
79  theProton= G4Proton::Proton();
80  isPiplus = false;
81  isInitialized = false;
82 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4VCrossSectionDataSet(const G4String &nam="")
const char * p
Definition: xmltok.h:285
int G4int
Definition: G4Types.hh:78
static constexpr double TeV
Definition: G4SIunits.hh:218
void SetMinKinEnergy(G4double value)
static G4Proton * Proton()
Definition: G4Proton.cc:93
void SetMaxKinEnergy(G4double value)
static constexpr double GeV
Definition: G4SIunits.hh:217
static constexpr double MeV
Definition: G4SIunits.hh:214
G4bool isInitialized()

Here is the call graph for this function:

G4BGGPionInelasticXS::~G4BGGPionInelasticXS ( )
virtual

Definition at line 85 of file G4BGGPionInelasticXS.cc.

86 {
87  delete fSAID;
88  delete fHadron;
89  delete fPion;
90  delete fGlauber;
91 }

Member Function Documentation

void G4BGGPionInelasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 184 of file G4BGGPionInelasticXS.cc.

185 {
186  if(&p == G4PionPlus::PionPlus() || &p == G4PionMinus::PionMinus()) {
187  particle = &p;
188  } else {
189  G4cout << "### G4BGGPionInelasticXS WARNING: is not applicable to "
190  << p.GetParticleName()
191  << G4endl;
192  throw G4HadronicException(__FILE__, __LINE__,
193  "G4BGGPionInelasticXS::BuildPhysicsTable is used for wrong particle");
194  return;
195  }
196 
197  if(isInitialized) { return; }
198  isInitialized = true;
199 
200  fPion = new G4UPiNuclearCrossSection();
201  fGlauber = new G4ComponentGGHadronNucleusXsc();
202  fHadron = new G4HadronNucleonXsc();
203  fSAID = new G4ComponentSAIDTotalXS();
204 
205  fPion->BuildPhysicsTable(*particle);
206  fGlauber->BuildPhysicsTable(*particle);
207 
208  if(particle == G4PionPlus::PionPlus()) { isPiplus = true; }
209 
210  G4ThreeVector mom(0.0,0.0,1.0);
211  G4DynamicParticle dp(particle, mom, fGlauberEnergy);
212 
214 
215  G4double csup, csdn;
216  G4int A;
217 
218  if(verboseLevel > 0) {
219  G4cout << "### G4BGGPionInelasticXS::Initialise for "
220  << particle->GetParticleName()
221  << " isPiplus: " << isPiplus
222  << G4endl;
223  }
224 
225  for(G4int iz=2; iz<93; iz++) {
226 
227  A = G4lrint(nist->GetAtomicMassAmu(iz));
228  theA[iz] = A;
229 
230  csup = fGlauber->GetInelasticGlauberGribov(&dp, iz, A);
231  csdn = fPion->GetInelasticCrossSection(&dp, iz, theA[iz]);
232 
233  theGlauberFac[iz] = csdn/csup;
234  if(verboseLevel > 0) {
235  G4cout << "Z= " << iz << " A= " << A
236  << " factor= " << theGlauberFac[iz] << G4endl;
237  }
238  }
239  dp.SetKineticEnergy(fSAIDHighEnergyLimit);
240  fHadron->GetHadronNucleonXscPDG(&dp, theProton);
241  theCoulombFac[1] = fSAIDHighEnergyLimit*
242  (fSAID->GetInelasticIsotopeCrossSection(particle,fSAIDHighEnergyLimit,1,1)
243  /fHadron->GetInelasticHadronNucleonXsc() - 1);
244 
245  if(isPiplus) {
246  dp.SetKineticEnergy(2*MeV);
247  for(G4int iz=2; iz<93; iz++) {
248  theCoulombFac[iz] = fPion->GetInelasticCrossSection(&dp, iz, theA[iz])
249  /CoulombFactor(2*MeV,iz);
250  }
251 
252  } else {
253  dp.SetKineticEnergy(fLowEnergy);
254  //fHadron->GetHadronNucleonXscNS(&dp, theProton);
255  //theCoulombFac[1] = theGlauberFac[1]*fHadron->GetInelasticHadronNucleonXsc();
256  for(G4int iz=2; iz<93; iz++) {
257  theCoulombFac[iz] = fPion->GetInelasticCrossSection(&dp, iz, theA[iz]);
258  }
259  }
260 }
virtual G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
const char * p
Definition: xmltok.h:285
G4double GetInelasticCrossSection(const G4DynamicParticle *aParticle, G4int Z, G4int A)
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
int G4lrint(double ad)
Definition: templates.hh:163
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
G4double GetAtomicMassAmu(const G4String &symb) const
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
double G4double
Definition: G4Types.hh:76
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4bool isInitialized()
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetInelasticHadronNucleonXsc()

Here is the call graph for this function:

void G4BGGPionInelasticXS::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 289 of file G4BGGPionInelasticXS.cc.

290 {
291  outFile << "The Barashenkov-Glauber-Gribov cross section handles inelastic\n"
292  << "pion scattering from nuclei at all energies. The Barashenkov\n"
293  << "parameterization is used below 91 GeV and the Glauber-Gribov\n"
294  << "parameterization is used above 91 GeV.\n";
295 }
G4double G4BGGPionInelasticXS::GetElementCrossSection ( const G4DynamicParticle dp,
G4int  Z,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 115 of file G4BGGPionInelasticXS.cc.

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 && !isPiplus) {
129  cross = theCoulombFac[Z];
130  } else if(ekin <= 2*MeV && isPiplus) {
131  cross = theCoulombFac[Z]*CoulombFactor(ekin, Z);
132  } else if(ekin > fGlauberEnergy) {
133  cross = theGlauberFac[Z]*fGlauber->GetInelasticGlauberGribov(dp, Z, theA[Z]);
134  } else {
135  cross = fPion->GetInelasticCrossSection(dp, Z, theA[Z]);
136  }
137  }
138  if(verboseLevel > 1) {
139  G4cout << "G4BGGPionInelasticXS::GetCrossSection for "
140  << dp->GetDefinition()->GetParticleName()
141  << " Ekin(GeV)= " << dp->GetKineticEnergy()
142  << " in nucleus Z= " << Z << " A= " << theA[Z]
143  << " XS(b)= " << cross/barn
144  << G4endl;
145  }
146  return cross;
147 }
G4double GetKineticEnergy() const
G4double GetInelasticCrossSection(const G4DynamicParticle *aParticle, G4int Z, G4int A)
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
double G4double
Definition: G4Types.hh:76
static constexpr double barn
Definition: G4SIunits.hh:105

Here is the call graph for this function:

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

157 {
158  // this method should be called only for Z = 1
159 
160  G4double cross = 0.0;
161  G4double ekin = dp->GetKineticEnergy();
162 
163  if(ekin <= fSAIDHighEnergyLimit) {
164  cross = fSAID->GetInelasticIsotopeCrossSection(particle, ekin, 1, 1);
165  } else {
166  fHadron->GetHadronNucleonXscPDG(dp, theProton);
167  cross = (theCoulombFac[1]/ekin + 1)*fHadron->GetInelasticHadronNucleonXsc();
168  }
169  cross *= A;
170 
171  if(verboseLevel > 1) {
172  G4cout << "G4BGGPionInelasticXS::GetCrossSection for "
173  << dp->GetDefinition()->GetParticleName()
174  << " Ekin(GeV)= " << dp->GetKineticEnergy()
175  << " in nucleus Z= " << Z << " A= " << A
176  << " XS(b)= " << cross/barn
177  << G4endl;
178  }
179  return cross;
180 }
virtual G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
G4double GetKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static constexpr double barn
Definition: G4SIunits.hh:105
G4double GetInelasticHadronNucleonXsc()

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 96 of file G4BGGPionInelasticXS.cc.

98 {
99  return true;
100 }
G4bool G4BGGPionInelasticXS::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 G4BGGPionInelasticXS.cc.

108 {
109  return (1 == Z && 2 >= A);
110 }
double A(double temperature)

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