Geant4  10.02.p03
G4BGGPionElasticXS Class Reference

#include <G4BGGPionElasticXS.hh>

Inheritance diagram for G4BGGPionElasticXS:
Collaboration diagram for G4BGGPionElasticXS:

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
 
- 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

G4BGGPionElasticXSoperator= (const G4BGGPionElasticXS &right)
 
 G4BGGPionElasticXS (const G4BGGPionElasticXS &)
 

Private Attributes

G4double fGlauberEnergy
 
G4double fLowEnergy
 
G4double fSAIDHighEnergyLimit
 
G4double theGlauberFac [93]
 
G4double theCoulombFac [93]
 
G4int theA [93]
 
const G4ParticleDefinitionparticle
 
const G4ParticleDefinitiontheProton
 
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 64 of file G4BGGPionElasticXS.hh.

Constructor & Destructor Documentation

◆ G4BGGPionElasticXS() [1/2]

G4BGGPionElasticXS::G4BGGPionElasticXS ( const G4ParticleDefinition )

Definition at line 56 of file G4BGGPionElasticXS.cc.

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

◆ ~G4BGGPionElasticXS()

G4BGGPionElasticXS::~G4BGGPionElasticXS ( )
virtual

Definition at line 83 of file G4BGGPionElasticXS.cc.

84 {
85  delete fSAID;
86  delete fHadron;
87  delete fPion;
88  delete fGlauber;
89 }
G4ComponentGGHadronNucleusXsc * fGlauber
G4ComponentSAIDTotalXS * fSAID
G4UPiNuclearCrossSection * fPion
G4HadronNucleonXsc * fHadron

◆ G4BGGPionElasticXS() [2/2]

G4BGGPionElasticXS::G4BGGPionElasticXS ( const G4BGGPionElasticXS )
private

Member Function Documentation

◆ BuildPhysicsTable()

void G4BGGPionElasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 191 of file G4BGGPionElasticXS.cc.

192 {
193  if(&p == G4PionPlus::PionPlus() || &p == G4PionMinus::PionMinus()) {
194  particle = &p;
195  } else {
196  G4cout << "### G4BGGPionElasticXS WARNING: is not applicable to "
197  << p.GetParticleName()
198  << G4endl;
199  throw G4HadronicException(__FILE__, __LINE__,
200  "G4BGGPionElasticXS::BuildPhysicsTable is used for wrong particle");
201  return;
202  }
203 
204  if(isInitialized) { return; }
205  isInitialized = true;
206 
209  fHadron = new G4HadronNucleonXsc();
211 
214 
215  if(particle == G4PionPlus::PionPlus()) { isPiplus = true; }
216 
217  G4ThreeVector mom(0.0,0.0,1.0);
219 
221 
222  G4double csup, csdn;
223  G4int A;
224 
225  if(verboseLevel > 0) {
226  G4cout << "### G4BGGPionElasticXS::Initialise for "
228  }
229  for(G4int iz=2; iz<93; iz++) {
230 
231  A = G4lrint(nist->GetAtomicMassAmu(iz));
232  theA[iz] = A;
233 
234  csup = fGlauber->GetElasticGlauberGribov(&dp, iz, A);
235  csdn = fPion->GetElasticCrossSection(&dp, iz, A);
236 
237  theGlauberFac[iz] = csdn/csup;
238  if(verboseLevel > 0) {
239  G4cout << "Z= " << iz << " A= " << A
240  << " factor= " << theGlauberFac[iz] << G4endl;
241  }
242  }
243  /*
244  dp.SetKineticEnergy(fLowEnergy);
245  fHadron->GetHadronNucleonXscNS(&dp, G4Proton::Proton());
246  theCoulombFac[1] = fHadron->GetElasticHadronNucleonXsc();
247  */
248  dp.SetKineticEnergy(fSAIDHighEnergyLimit);
250  theCoulombFac[1] =
253 
254  dp.SetKineticEnergy(fLowEnergy);
255  for(G4int iz=2; iz<93; iz++) {
257  if(verboseLevel > 0) {
258  G4cout << "Z= " << iz << " A= " << A
259  << " factor= " << theCoulombFac[iz] << G4endl;
260  }
261  }
262 }
G4ComponentGGHadronNucleusXsc * fGlauber
G4double GetElasticHadronNucleonXsc()
G4double GetElasticCrossSection(const G4DynamicParticle *aParticle, G4int Z, G4int A)
G4ComponentSAIDTotalXS * fSAID
virtual G4double GetElasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
G4UPiNuclearCrossSection * fPion
G4HadronNucleonXsc * fHadron
int G4int
Definition: G4Types.hh:78
static G4NistManager * Instance()
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
const G4ParticleDefinition * theProton
G4double iz
Definition: TRTMaterials.hh:39
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
G4double GetElasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
const G4ParticleDefinition * particle
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
double G4double
Definition: G4Types.hh:76
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetAtomicMassAmu(const G4String &symb) const
void BuildPhysicsTable(const G4ParticleDefinition &)
Here is the call graph for this function:

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 265 of file G4BGGPionElasticXS.cc.

266 {
267  outFile << "The Barashenkov-Glauber-Gribov cross section handles elastic\n"
268  << "scattering of pions from nuclei at all energies. The\n"
269  << "Barashenkov parameterization is used below 91 GeV and the\n"
270  << "Glauber-Gribov parameterization is used above 91 GeV.\n";
271 }

◆ GetElementCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 113 of file G4BGGPionElasticXS.cc.

115 {
116  // this method should be called only for Z > 1
117 
118  G4double cross = 0.0;
119  G4double ekin = dp->GetKineticEnergy();
120  G4int Z = ZZ;
121  if(1 == Z) {
122  cross = 1.0115*GetIsoCrossSection(dp,1,1);
123  } else {
124  if(Z > 92) { Z = 92; }
125 
126  if(ekin <= fLowEnergy) {
127  cross = theCoulombFac[Z];
128  } else if(ekin > fGlauberEnergy) {
129  cross = theGlauberFac[Z]*fGlauber->GetElasticGlauberGribov(dp, Z, theA[Z]);
130  } else {
131  cross = fPion->GetElasticCrossSection(dp, Z, theA[Z]);
132  }
133  }
134  if(verboseLevel > 1) {
135  G4cout << "G4BGGPionElasticXS::GetElementCrossSection for "
136  << dp->GetDefinition()->GetParticleName()
137  << " Ekin(GeV)= " << dp->GetKineticEnergy()
138  << " in nucleus Z= " << Z << " A= " << theA[Z]
139  << " XS(b)= " << cross/barn
140  << G4endl;
141  }
142  return cross;
143 }
G4ComponentGGHadronNucleusXsc * fGlauber
G4double GetElasticCrossSection(const G4DynamicParticle *aParticle, G4int Z, G4int A)
G4UPiNuclearCrossSection * fPion
int G4int
Definition: G4Types.hh:78
G4double GetKineticEnergy() const
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
Float_t Z
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4double GetElasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
G4ParticleDefinition * GetDefinition() const
#define G4endl
Definition: G4ios.hh:61
static const double barn
Definition: G4SIunits.hh:104
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ GetIsoCrossSection()

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

151 {
152  // this method should be called only for Z = 1
153 
154  G4double cross = 0.0;
155  G4double ekin = dp->GetKineticEnergy();
156 
157  if(ekin <= fSAIDHighEnergyLimit) {
158  cross = fSAID->GetElasticIsotopeCrossSection(particle, ekin, 1, 1);
159  } else {
162  }
163  cross *= A;
164  /*
165  if(ekin <= fLowEnergy) {
166  cross = theCoulombFac[1];
167 
168  } else if( A < 2) {
169  fHadron->GetHadronNucleonXscNS(dp, G4Proton::Proton());
170  cross = fHadron->GetElasticHadronNucleonXsc();
171  } else {
172  fHadron->GetHadronNucleonXscNS(dp, G4Proton::Proton());
173  cross = fHadron->GetElasticHadronNucleonXsc();
174  fHadron->GetHadronNucleonXscNS(dp, G4Neutron::Neutron());
175  cross += fHadron->GetElasticHadronNucleonXsc();
176  }
177  */
178  if(verboseLevel > 1) {
179  G4cout << "G4BGGPionElasticXS::GetIsoCrossSection for "
180  << dp->GetDefinition()->GetParticleName()
181  << " Ekin(GeV)= " << dp->GetKineticEnergy()
182  << " in nucleus Z= " << Z << " A= " << A
183  << " XS(b)= " << cross/barn
184  << G4endl;
185  }
186  return cross;
187 }
G4double GetElasticHadronNucleonXsc()
G4ComponentSAIDTotalXS * fSAID
virtual G4double GetElasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
G4HadronNucleonXsc * fHadron
G4double GetKineticEnergy() const
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
const G4ParticleDefinition * theProton
Float_t Z
const G4ParticleDefinition * particle
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
Here is the call graph for this function:
Here is the caller graph for this function:

◆ IsElementApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 94 of file G4BGGPionElasticXS.cc.

96 {
97  return true;
98 }

◆ IsIsoApplicable()

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.

106 {
107  return (1 == Z && 2 >= A);
108 }
double A(double temperature)
Float_t Z

◆ operator=()

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

Member Data Documentation

◆ fGlauber

G4ComponentGGHadronNucleusXsc* G4BGGPionElasticXS::fGlauber
private

Definition at line 110 of file G4BGGPionElasticXS.hh.

◆ fGlauberEnergy

G4double G4BGGPionElasticXS::fGlauberEnergy
private

Definition at line 101 of file G4BGGPionElasticXS.hh.

◆ fHadron

G4HadronNucleonXsc* G4BGGPionElasticXS::fHadron
private

Definition at line 112 of file G4BGGPionElasticXS.hh.

◆ fLowEnergy

G4double G4BGGPionElasticXS::fLowEnergy
private

Definition at line 102 of file G4BGGPionElasticXS.hh.

◆ fPion

G4UPiNuclearCrossSection* G4BGGPionElasticXS::fPion
private

Definition at line 111 of file G4BGGPionElasticXS.hh.

◆ fSAID

G4ComponentSAIDTotalXS* G4BGGPionElasticXS::fSAID
private

Definition at line 113 of file G4BGGPionElasticXS.hh.

◆ fSAIDHighEnergyLimit

G4double G4BGGPionElasticXS::fSAIDHighEnergyLimit
private

Definition at line 103 of file G4BGGPionElasticXS.hh.

◆ isInitialized

G4bool G4BGGPionElasticXS::isInitialized
private

Definition at line 115 of file G4BGGPionElasticXS.hh.

◆ isPiplus

G4bool G4BGGPionElasticXS::isPiplus
private

Definition at line 114 of file G4BGGPionElasticXS.hh.

◆ particle

const G4ParticleDefinition* G4BGGPionElasticXS::particle
private

Definition at line 108 of file G4BGGPionElasticXS.hh.

◆ theA

G4int G4BGGPionElasticXS::theA[93]
private

Definition at line 106 of file G4BGGPionElasticXS.hh.

◆ theCoulombFac

G4double G4BGGPionElasticXS::theCoulombFac[93]
private

Definition at line 105 of file G4BGGPionElasticXS.hh.

◆ theGlauberFac

G4double G4BGGPionElasticXS::theGlauberFac[93]
private

Definition at line 104 of file G4BGGPionElasticXS.hh.

◆ theProton

const G4ParticleDefinition* G4BGGPionElasticXS::theProton
private

Definition at line 109 of file G4BGGPionElasticXS.hh.


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