Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4BGGPionInelasticXS.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id$
27//
28// -------------------------------------------------------------------
29//
30// GEANT4 Class file
31//
32//
33// File name: G4BGGPionInelasticXS
34//
35// Author: Vladimir Ivanchenko
36//
37// Creation date: 01.10.2003
38// Modifications:
39//
40// -------------------------------------------------------------------
41//
42
44#include "G4SystemOfUnits.hh"
47#include "G4HadronNucleonXsc.hh"
48//#include "G4HadronInelasticDataSet.hh"
50
51#include "G4Proton.hh"
52#include "G4PionPlus.hh"
53#include "G4PionMinus.hh"
54#include "G4NistManager.hh"
55
56
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}
82
83
85{
86 delete fGlauber;
87 delete fPion;
88 delete fHadron;
89 delete fSAID;
90}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
93
94G4bool
96 const G4Material*)
97{
98 return true;
99}
100
101//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
102
104 G4int Z, G4int A,
105 const G4Element*,
106 const G4Material*)
107{
108 return (1 == Z && 2 >= A);
109}
110
111//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
112
115 G4int ZZ, const G4Material*)
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}
147
148//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
149
152 G4int Z, G4int A,
153 const G4Isotope*,
154 const G4Element*,
155 const G4Material*)
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}
200
201//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
202
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}
295
296//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
297
298G4double G4BGGPionInelasticXS::CoulombFactor(G4double kinEnergy, G4int Z)
299{
300 G4int A = theA[Z];
301 G4double res= 0.0;
302 if(kinEnergy <= DBL_MIN) { return res; }
303 else if(A < 2) { return kinEnergy*kinEnergy; }
304
305 G4double elog = std::log10(6.7*kinEnergy/GeV);
306 G4double aa = A;
307
308 // from G4ProtonInelasticCrossSection
309 G4double ff1 = 0.70 - 0.002*aa; // slope of the drop at medium energies.
310 G4double ff2 = 1.00 + 1/aa; // start of the slope.
311 G4double ff3 = 0.8 + 18/aa - 0.002*aa; // stephight
312 res = 1.0 + ff3*(1.0 - (1.0/(1+std::exp(-8*ff1*(elog + 1.37*ff2)))));
313
314 ff1 = 1. - 1./aa - 0.001*aa; // slope of the rise
315 ff2 = 1.17 - 2.7/aa-0.0014*aa; // start of the rise
316 res /= (1 + std::exp(-8.*ff1*(elog + 2*ff2)));
317 /*
318 G4double f1 = 8.0 - 8.0/aa - 0.008*aa;
319 G4double f2 = 2.34 - 5.4/aa - 0.0028*aa;
320
321 res = 1.0/(1.0 + std::exp(-f1*(elog + f2)));
322
323 f1 = 5.6 - 0.016*aa;
324 f2 = 1.37 + 1.37/aa;
325 res *= ( 1.0 + (0.8 + 18./aa - 0.002*aa)/(1.0 + std::exp(f1*(elog + f2))));
326 */
327 return res;
328}
329
330//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
331
332void
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}
340
341//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
virtual G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, 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 G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
virtual void CrossSectionDescription(std::ostream &) const
G4BGGPionInelasticXS(const G4ParticleDefinition *)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
void SetKineticEnergy(G4double aEnergy)
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
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double GetInelasticCrossSection(const G4DynamicParticle *aParticle, G4int Z, G4int A)
void BuildPhysicsTable(const G4ParticleDefinition &)
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
int G4lrint(double ad)
Definition: templates.hh:163
#define DBL_MIN
Definition: templates.hh:75