Geant4  10.02.p03
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
 

Private Member Functions

G4double CoulombFactor (G4double kinEnergy, G4int Z)
 
G4BGGPionInelasticXSoperator= (const G4BGGPionInelasticXS &right)
 
 G4BGGPionInelasticXS (const G4BGGPionInelasticXS &)
 

Private Attributes

G4double fGlauberEnergy
 
G4double fLowEnergy
 
G4double fSAIDHighEnergyLimit
 
G4double theGlauberFac [93]
 
G4double theCoulombFac [93]
 
G4int theA [93]
 
const G4ParticleDefinitionparticle
 
const G4ParticleDefinitiontheProton
 
G4PowfG4pow
 
G4ComponentGGHadronNucleusXscfGlauber
 
G4UPiNuclearCrossSectionfPion
 
G4HadronNucleonXscfHadron
 
G4ComponentSAIDTotalXSfSAID
 
G4bool isPiplus
 
G4bool isInitialized
 

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() [1/2]

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;
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 
77 
78  particle = p;
80  isPiplus = false;
81  isInitialized = false;
82 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
static const double MeV
Definition: G4SIunits.hh:211
G4VCrossSectionDataSet(const G4String &nam="")
const G4ParticleDefinition * particle
int G4int
Definition: G4Types.hh:78
G4UPiNuclearCrossSection * fPion
G4HadronNucleonXsc * fHadron
void SetMinKinEnergy(G4double value)
G4ComponentSAIDTotalXS * fSAID
static G4Proton * Proton()
Definition: G4Proton.cc:93
static const double GeV
Definition: G4SIunits.hh:214
const G4ParticleDefinition * theProton
G4ComponentGGHadronNucleusXsc * fGlauber
void SetMaxKinEnergy(G4double value)
static const double TeV
Definition: G4SIunits.hh:215
Here is the call graph for this function:

◆ ~G4BGGPionInelasticXS()

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 }
G4UPiNuclearCrossSection * fPion
G4HadronNucleonXsc * fHadron
G4ComponentSAIDTotalXS * fSAID
G4ComponentGGHadronNucleusXsc * fGlauber

◆ G4BGGPionInelasticXS() [2/2]

G4BGGPionInelasticXS::G4BGGPionInelasticXS ( const G4BGGPionInelasticXS )
private

Member Function Documentation

◆ BuildPhysicsTable()

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 
202  fHadron = new G4HadronNucleonXsc();
204 
207 
208  if(particle == G4PionPlus::PionPlus()) { isPiplus = true; }
209 
210  G4ThreeVector mom(0.0,0.0,1.0);
212 
214 
215  G4double csup, csdn;
216  G4int A;
217 
218  if(verboseLevel > 0) {
219  G4cout << "### G4BGGPionInelasticXS::Initialise for "
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);
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 }
static const double MeV
Definition: G4SIunits.hh:211
virtual G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
G4double GetInelasticCrossSection(const G4DynamicParticle *aParticle, G4int Z, G4int A)
const G4ParticleDefinition * particle
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
G4UPiNuclearCrossSection * fPion
G4HadronNucleonXsc * fHadron
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
G4double CoulombFactor(G4double kinEnergy, G4int Z)
G4double iz
Definition: TRTMaterials.hh:39
G4ComponentSAIDTotalXS * fSAID
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
const G4ParticleDefinition * theProton
G4ComponentGGHadronNucleusXsc * fGlauber
int G4lrint(double ad)
Definition: templates.hh:163
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
#define G4endl
Definition: G4ios.hh:61
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
double G4double
Definition: G4Types.hh:76
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetAtomicMassAmu(const G4String &symb) const
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetInelasticHadronNucleonXsc()
Here is the call graph for this function:

◆ CoulombFactor()

G4double G4BGGPionInelasticXS::CoulombFactor ( G4double  kinEnergy,
G4int  Z 
)
private

Definition at line 264 of file G4BGGPionInelasticXS.cc.

265 {
266  G4int A = theA[Z];
267  G4double res= 0.0;
268  if(kinEnergy <= DBL_MIN) { return res; }
269  else if(A < 2) { return kinEnergy*kinEnergy; }
270 
271  G4double elog = fG4pow->log10A(6.7*kinEnergy/GeV);
272  G4double aa = A;
273 
274  // from G4ProtonInelasticCrossSection
275  G4double ff1 = 0.70 - 0.002*aa; // slope of the drop at medium energies.
276  G4double ff2 = 1.00 + 1/aa; // start of the slope.
277  G4double ff3 = 0.8 + 18/aa - 0.002*aa; // stephight
278  res = 1.0 + ff3*(1.0 - (1.0/(1+fG4pow->expA(-8*ff1*(elog + 1.37*ff2)))));
279 
280  ff1 = 1. - 1./aa - 0.001*aa; // slope of the rise
281  ff2 = 1.17 - 2.7/aa-0.0014*aa; // start of the rise
282  res /= (1 + fG4pow->expA(-8.*ff1*(elog + 2*ff2)));
283  return res;
284 }
G4double log10A(G4double A) const
Definition: G4Pow.hh:230
int G4int
Definition: G4Types.hh:78
double A(double temperature)
Float_t Z
static const double GeV
Definition: G4SIunits.hh:214
#define DBL_MIN
Definition: templates.hh:75
double G4double
Definition: G4Types.hh:76
G4double expA(G4double A) const
Definition: G4Pow.hh:235
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CrossSectionDescription()

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 }

◆ GetElementCrossSection()

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) {
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 }
static const double MeV
Definition: G4SIunits.hh:211
G4double GetInelasticCrossSection(const G4DynamicParticle *aParticle, G4int Z, G4int A)
int G4int
Definition: G4Types.hh:78
G4UPiNuclearCrossSection * fPion
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4double GetKineticEnergy() const
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
Float_t Z
G4double CoulombFactor(G4double kinEnergy, G4int Z)
G4ComponentGGHadronNucleusXsc * fGlauber
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
static const double barn
Definition: G4SIunits.hh:104
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ GetIsoCrossSection()

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 {
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)
const G4ParticleDefinition * particle
G4HadronNucleonXsc * fHadron
G4double GetKineticEnergy() const
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
Float_t Z
G4ComponentSAIDTotalXS * fSAID
const G4ParticleDefinition * theProton
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
static const double barn
Definition: G4SIunits.hh:104
double G4double
Definition: G4Types.hh:76
G4double GetInelasticHadronNucleonXsc()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ IsElementApplicable()

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 }

◆ IsIsoApplicable()

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)
Float_t Z

◆ operator=()

G4BGGPionInelasticXS& G4BGGPionInelasticXS::operator= ( const G4BGGPionInelasticXS right)
private

Member Data Documentation

◆ fG4pow

G4Pow* G4BGGPionInelasticXS::fG4pow
private

Definition at line 117 of file G4BGGPionInelasticXS.hh.

◆ fGlauber

G4ComponentGGHadronNucleusXsc* G4BGGPionInelasticXS::fGlauber
private

Definition at line 119 of file G4BGGPionInelasticXS.hh.

◆ fGlauberEnergy

G4double G4BGGPionInelasticXS::fGlauberEnergy
private

Definition at line 107 of file G4BGGPionInelasticXS.hh.

◆ fHadron

G4HadronNucleonXsc* G4BGGPionInelasticXS::fHadron
private

Definition at line 121 of file G4BGGPionInelasticXS.hh.

◆ fLowEnergy

G4double G4BGGPionInelasticXS::fLowEnergy
private

Definition at line 108 of file G4BGGPionInelasticXS.hh.

◆ fPion

G4UPiNuclearCrossSection* G4BGGPionInelasticXS::fPion
private

Definition at line 120 of file G4BGGPionInelasticXS.hh.

◆ fSAID

G4ComponentSAIDTotalXS* G4BGGPionInelasticXS::fSAID
private

Definition at line 122 of file G4BGGPionInelasticXS.hh.

◆ fSAIDHighEnergyLimit

G4double G4BGGPionInelasticXS::fSAIDHighEnergyLimit
private

Definition at line 109 of file G4BGGPionInelasticXS.hh.

◆ isInitialized

G4bool G4BGGPionInelasticXS::isInitialized
private

Definition at line 124 of file G4BGGPionInelasticXS.hh.

◆ isPiplus

G4bool G4BGGPionInelasticXS::isPiplus
private

Definition at line 123 of file G4BGGPionInelasticXS.hh.

◆ particle

const G4ParticleDefinition* G4BGGPionInelasticXS::particle
private

Definition at line 114 of file G4BGGPionInelasticXS.hh.

◆ theA

G4int G4BGGPionInelasticXS::theA[93]
private

Definition at line 112 of file G4BGGPionInelasticXS.hh.

◆ theCoulombFac

G4double G4BGGPionInelasticXS::theCoulombFac[93]
private

Definition at line 111 of file G4BGGPionInelasticXS.hh.

◆ theGlauberFac

G4double G4BGGPionInelasticXS::theGlauberFac[93]
private

Definition at line 110 of file G4BGGPionInelasticXS.hh.

◆ theProton

const G4ParticleDefinition* G4BGGPionInelasticXS::theProton
private

Definition at line 115 of file G4BGGPionInelasticXS.hh.


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