Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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: G4BGGPionInelasticXS.cc 70848 2013-06-06 12:00:02Z gcosmo $
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 
43 #include "G4BGGPionInelasticXS.hh"
44 #include "G4SystemOfUnits.hh"
47 #include "G4HadronNucleonXsc.hh"
49 
50 #include "G4Proton.hh"
51 #include "G4PionPlus.hh"
52 #include "G4PionMinus.hh"
53 #include "G4NistManager.hh"
54 #include "G4Pow.hh"
55 
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 }
83 
84 
86 {
87  delete fGlauber;
88  delete fPion;
89  delete fHadron;
90  delete fSAID;
91 }
92 
93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
94 
95 G4bool
97  const G4Material*)
98 {
99  return true;
100 }
101 
102 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
103 
105  G4int Z, G4int A,
106  const G4Element*,
107  const G4Material*)
108 {
109  return (1 == Z && 2 >= A);
110 }
111 
112 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
113 
114 G4double
116  G4int ZZ, const G4Material*)
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 }
148 
149 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
150 
151 G4double
153  G4int Z, G4int A,
154  const G4Isotope*,
155  const G4Element*,
156  const G4Material*)
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 }
181 
182 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
183 
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 G4GlauberGribovCrossSection();
202  fHadron = new G4HadronNucleonXsc();
203  fSAID = new G4ComponentSAIDTotalXS();
204 
205  fPion->BuildPhysicsTable(*particle);
206  fGlauber->BuildPhysicsTable(*particle);
207  if(particle == G4PionPlus::PionPlus()) { isPiplus = true; }
208 
209  G4ParticleDefinition* part = const_cast<G4ParticleDefinition*>(particle);
210  G4ThreeVector mom(0.0,0.0,1.0);
211  G4DynamicParticle dp(part, 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 }
261 
262 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
263 
264 G4double G4BGGPionInelasticXS::CoulombFactor(G4double kinEnergy, G4int Z)
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 }
285 
286 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
287 
288 void
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 }
296 
297 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
static G4Pow * GetInstance()
Definition: G4Pow.cc:53
virtual G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
G4double GetKineticEnergy() const
std::ofstream outFile
Definition: GammaRayTel.cc:68
G4double expA(G4double A) const
Definition: G4Pow.hh:238
G4double log10A(G4double A) const
Definition: G4Pow.hh:233
const char * p
Definition: xmltok.h:285
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetInelasticCrossSection(const G4DynamicParticle *aParticle, G4int Z, G4int A)
virtual void CrossSectionDescription(std::ostream &) const
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
static G4NistManager * Instance()
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
virtual G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
void SetMinKinEnergy(G4double value)
bool G4bool
Definition: G4Types.hh:79
G4double iz
Definition: TRTMaterials.hh:39
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
G4BGGPionInelasticXS(const G4ParticleDefinition *)
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
void SetKineticEnergy(G4double aEnergy)
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
void SetMaxKinEnergy(G4double value)
int G4lrint(double ad)
Definition: templates.hh:163
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
#define DBL_MIN
Definition: templates.hh:75
G4double GetAtomicMassAmu(const G4String &symb) const
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetInelasticHadronNucleonXsc()