Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BGGPionInelasticXS Class Reference

#include <G4BGGPionInelasticXS.hh>

+ Inheritance 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 ()
 
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)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, 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 G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
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::G4BGGPionInelasticXS ( const G4ParticleDefinition p)

Definition at line 57 of file G4BGGPionInelasticXS.cc.

58 : G4VCrossSectionDataSet("Barashenkov-Glauber-Gribov")
59{
60 verboseLevel = 0;
61 fGlauberEnergy = 91.*GeV;
62 fLowEnergy = 20.*MeV;
63 fSAIDHighEnergyLimit = 2.6*GeV;
64 SetMinKinEnergy(0.0);
65 SetMaxKinEnergy(100*TeV);
66
67 for (G4int i = 0; i < 93; i++) {
68 theGlauberFac[i] = 0.0;
69 theCoulombFac[i] = 0.0;
70 theA[i] = 1;
71 }
72 fPion = 0;
73 fGlauber = 0;
74 fHadron = 0;
75 // fGHEISHA = 0;
76 fSAID = 0;
77 particle = p;
78 theProton= G4Proton::Proton();
79 isPiplus = false;
80 isInitialized = false;
81}
int G4int
Definition: G4Types.hh:66
static G4Proton * Proton()
Definition: G4Proton.cc:93
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)

◆ ~G4BGGPionInelasticXS()

G4BGGPionInelasticXS::~G4BGGPionInelasticXS ( )
virtual

Definition at line 84 of file G4BGGPionInelasticXS.cc.

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

Member Function Documentation

◆ BuildPhysicsTable()

void G4BGGPionInelasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 203 of file G4BGGPionInelasticXS.cc.

204{
205 if(&p == G4PionPlus::PionPlus() || &p == G4PionMinus::PionMinus()) {
206 particle = &p;
207 } else {
208 G4cout << "### G4BGGPionInelasticXS WARNING: is not applicable to "
209 << p.GetParticleName()
210 << G4endl;
211 throw G4HadronicException(__FILE__, __LINE__,
212 "G4BGGPionInelasticXS::BuildPhysicsTable is used for wrong particle");
213 return;
214 }
215
216 if(isInitialized) { return; }
217 isInitialized = true;
218
219 fPion = new G4UPiNuclearCrossSection();
220 fGlauber = new G4GlauberGribovCrossSection();
221 fHadron = new G4HadronNucleonXsc();
222 // fGHEISHA = new G4HadronInelasticDataSet();
223 fSAID = new G4ComponentSAIDTotalXS();
224
225 fPion->BuildPhysicsTable(*particle);
226 fGlauber->BuildPhysicsTable(*particle);
227 if(particle == G4PionPlus::PionPlus()) { isPiplus = true; }
228
229 G4ParticleDefinition* part = const_cast<G4ParticleDefinition*>(particle);
230 G4ThreeVector mom(0.0,0.0,1.0);
231 G4DynamicParticle dp(part, mom, fGlauberEnergy);
232
234
235 G4double csup, csdn;
236 G4int A;
237
238 if(verboseLevel > 0) {
239 G4cout << "### G4BGGPionInelasticXS::Initialise for "
240 << particle->GetParticleName()
241 << " isPiplus: " << isPiplus
242 << G4endl;
243 }
244
245 for(G4int iz=2; iz<93; iz++) {
246
247 A = G4lrint(nist->GetAtomicMassAmu(iz));
248 theA[iz] = A;
249
250 csup = fGlauber->GetInelasticGlauberGribov(&dp, iz, A);
251 csdn = fPion->GetInelasticCrossSection(&dp, iz, theA[iz]);
252
253 theGlauberFac[iz] = csdn/csup;
254 if(verboseLevel > 0) {
255 G4cout << "Z= " << iz << " A= " << A
256 << " factor= " << theGlauberFac[iz] << G4endl;
257 }
258 }
259 dp.SetKineticEnergy(fSAIDHighEnergyLimit);
260 fHadron->GetHadronNucleonXscPDG(&dp, theProton);
261 theCoulombFac[1] = fSAIDHighEnergyLimit*
262 (fSAID->GetInelasticIsotopeCrossSection(particle,fSAIDHighEnergyLimit,1,1)
263 /fHadron->GetInelasticHadronNucleonXsc() - 1);
264
265 /*
266 dp.SetKineticEnergy(20*GeV);
267 const G4Material* mat = 0;
268 if(isPiplus) {
269 fHadron->GetHadronNucleonXscPDG(&dp, theProton);
270 theCoulombFac[1] = fHadron->GetInelasticHadronNucleonXsc()
271 /fGHEISHA->GetElementCrossSection(&dp, 1, mat);
272 } else {
273 fHadron->GetHadronNucleonXscPDG(&dp, theProton);
274 theGlauberFac[1] = fHadron->GetInelasticHadronNucleonXsc();
275 fHadron->GetHadronNucleonXscNS(&dp, theProton);
276 theGlauberFac[1] /= fHadron->GetInelasticHadronNucleonXsc();
277 }
278 */
279 if(isPiplus) {
280 dp.SetKineticEnergy(2*MeV);
281 for(G4int iz=2; iz<93; iz++) {
282 theCoulombFac[iz] = fPion->GetInelasticCrossSection(&dp, iz, theA[iz])
283 /CoulombFactor(2*MeV,iz);
284 }
285
286 } else {
287 dp.SetKineticEnergy(fLowEnergy);
288 //fHadron->GetHadronNucleonXscNS(&dp, theProton);
289 //theCoulombFac[1] = theGlauberFac[1]*fHadron->GetInelasticHadronNucleonXsc();
290 for(G4int iz=2; iz<93; iz++) {
291 theCoulombFac[iz] = fPion->GetInelasticCrossSection(&dp, iz, theA[iz]);
292 }
293 }
294}
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
virtual G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
G4double GetInelasticHadronNucleonXsc()
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const
const G4String & GetParticleName() const
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
G4double GetInelasticCrossSection(const G4DynamicParticle *aParticle, G4int Z, G4int A)
void BuildPhysicsTable(const G4ParticleDefinition &)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
int G4lrint(double ad)
Definition: templates.hh:163

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 333 of file G4BGGPionInelasticXS.cc.

334{
335 outFile << "The Barashenkov-Glauber-Gribov cross section handles inelastic\n"
336 << "pion scattering from nuclei at all energies. The Barashenkov\n"
337 << "parameterization is used below 91 GeV and the Glauber-Gribov\n"
338 << "parameterization is used above 91 GeV.\n";
339}

◆ GetElementCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 114 of file G4BGGPionInelasticXS.cc.

116{
117 // this method should be called only for Z > 1
118
119 G4double cross = 0.0;
120 G4double ekin = dp->GetKineticEnergy();
121 G4int Z = ZZ;
122 if(1 == Z) {
123 cross = 1.0115*GetIsoCrossSection(dp,1,1);
124 } else {
125 if(Z > 92) { Z = 92; }
126
127 if(ekin <= fLowEnergy && !isPiplus) {
128 cross = theCoulombFac[Z];
129 } else if(ekin <= 2*MeV && isPiplus) {
130 cross = theCoulombFac[Z]*CoulombFactor(ekin, Z);
131 } else if(ekin > fGlauberEnergy) {
132 cross = theGlauberFac[Z]*fGlauber->GetInelasticGlauberGribov(dp, Z, theA[Z]);
133 } else {
134 cross = fPion->GetInelasticCrossSection(dp, Z, theA[Z]);
135 }
136 }
137 if(verboseLevel > 1) {
138 G4cout << "G4BGGPionInelasticXS::GetCrossSection for "
140 << " Ekin(GeV)= " << dp->GetKineticEnergy()
141 << " in nucleus Z= " << Z << " A= " << theA[Z]
142 << " XS(b)= " << cross/barn
143 << G4endl;
144 }
145 return cross;
146}
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const

◆ 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 151 of file G4BGGPionInelasticXS.cc.

156{
157 // this method should be called only for Z = 1
158
159 G4double cross = 0.0;
160 G4double ekin = dp->GetKineticEnergy();
161
162 if(ekin <= fSAIDHighEnergyLimit) {
163 cross = fSAID->GetInelasticIsotopeCrossSection(particle, ekin, 1, 1);
164 } else {
165 fHadron->GetHadronNucleonXscPDG(dp, theProton);
166 cross = (theCoulombFac[1]/ekin + 1)*fHadron->GetInelasticHadronNucleonXsc();
167 }
168 /*
169 if(isPiplus) {
170 if(ekin <= 20*GeV) {
171 cross = theCoulombFac[1]*fGHEISHA->GetElementCrossSection(dp, 1, mat);
172 } else {
173 fHadron->GetHadronNucleonXscPDG(dp, theProton);
174 cross = fHadron->GetInelasticHadronNucleonXsc();
175 }
176 } else {
177 if(ekin <= fLowEnergy) {
178 cross = theCoulombFac[1];
179 } else if(ekin <= 20*GeV) {
180 fHadron->GetHadronNucleonXscNS(dp, theProton);
181 cross = theGlauberFac[1]*fHadron->GetInelasticHadronNucleonXsc();
182 } else {
183 fHadron->GetHadronNucleonXscPDG(dp, theProton);
184 cross = fHadron->GetInelasticHadronNucleonXsc();
185 }
186 }
187 */
188 cross *= A;
189
190 if(verboseLevel > 1) {
191 G4cout << "G4BGGPionInelasticXS::GetCrossSection for "
193 << " Ekin(GeV)= " << dp->GetKineticEnergy()
194 << " in nucleus Z= " << Z << " A= " << A
195 << " XS(b)= " << cross/barn
196 << G4endl;
197 }
198 return cross;
199}

Referenced by GetElementCrossSection().

◆ IsElementApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 95 of file G4BGGPionInelasticXS.cc.

97{
98 return true;
99}

◆ 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 103 of file G4BGGPionInelasticXS.cc.

107{
108 return (1 == Z && 2 >= A);
109}

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