Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4BGGPionInelasticXS.cc
이 파일의 문서화 페이지로 가기
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 110543 2018-05-29 13:38:54Z 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;
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;
83 }
84 
85 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
86 
88 {
89  delete fSAID;
90  delete fHadron;
91  delete fPion;
92  delete fGlauber;
93 }
94 
95 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
96 
97 G4bool
99  const G4Material*)
100 {
101  return true;
102 }
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
105 
107  G4int Z, G4int A,
108  const G4Element*,
109  const G4Material*)
110 {
111  return (1 == Z && 2 >= A);
112 }
113 
114 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
115 
116 G4double
118  G4int ZZ, const G4Material*)
119 {
120  // this method should be called only for Z > 1
121 
122  G4double cross = 0.0;
123  G4double ekin = dp->GetKineticEnergy();
124  G4int Z = ZZ;
125  if(1 == Z) {
126  cross = 1.0115*GetIsoCrossSection(dp,1,1);
127  } else {
128  if(Z > 92) { Z = 92; }
129 
130  if(ekin <= fLowEnergy && !isPiplus) {
131  cross = theCoulombFac[Z];
132  } else if(ekin <= 2*MeV && isPiplus) {
133  cross = theCoulombFac[Z]*CoulombFactor(ekin, Z);
134  } else if(ekin > fGlauberEnergy) {
136  } else {
137  cross = fPion->GetInelasticCrossSection(dp, Z, theA[Z]);
138  }
139  }
140  if(verboseLevel > 1) {
141  G4cout << "G4BGGPionInelasticXS::GetCrossSection for "
142  << dp->GetDefinition()->GetParticleName()
143  << " Ekin(GeV)= " << dp->GetKineticEnergy()
144  << " in nucleus Z= " << Z << " A= " << theA[Z]
145  << " XS(b)= " << cross/barn
146  << G4endl;
147  }
148  return cross;
149 }
150 
151 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
152 
153 G4double
155  G4int Z, G4int A,
156  const G4Isotope*,
157  const G4Element*,
158  const G4Material*)
159 {
160  // this method should be called only for Z = 1
161 
162  G4double cross = 0.0;
163  G4double ekin = dp->GetKineticEnergy();
164 
165  if(ekin <= fSAIDHighEnergyLimit) {
166  cross = fSAID->GetInelasticIsotopeCrossSection(particle, ekin, 1, 1);
167  } else {
169  cross = (theCoulombFac[1]/ekin + 1)*fHadron->GetInelasticHadronNucleonXsc();
170  }
171  cross *= A;
172 
173  if(verboseLevel > 1) {
174  G4cout << "G4BGGPionInelasticXS::GetCrossSection for "
175  << dp->GetDefinition()->GetParticleName()
176  << " Ekin(GeV)= " << dp->GetKineticEnergy()
177  << " in nucleus Z= " << Z << " A= " << A
178  << " XS(b)= " << cross/barn
179  << G4endl;
180  }
181  return cross;
182 }
183 
184 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
185 
187 {
188  if(&p == G4PionPlus::PionPlus() || &p == G4PionMinus::PionMinus()) {
189  particle = &p;
190  } else {
191  G4cout << "### G4BGGPionInelasticXS WARNING: is not applicable to "
192  << p.GetParticleName()
193  << G4endl;
194  throw G4HadronicException(__FILE__, __LINE__,
195  "G4BGGPionInelasticXS::BuildPhysicsTable is used for wrong particle");
196  return;
197  }
198 
199  if(isInitialized) { return; }
200  isInitialized = true;
201 
204  fHadron = new G4HadronNucleonXsc();
206 
209 
210  if(particle == G4PionPlus::PionPlus()) { isPiplus = true; }
211 
212  G4ThreeVector mom(0.0,0.0,1.0);
214 
216 
217  G4double csup, csdn;
218  G4int A;
219 
220  if(verboseLevel > 0) {
221  G4cout << "### G4BGGPionInelasticXS::Initialise for "
223  << " isPiplus: " << isPiplus
224  << G4endl;
225  }
226 
227  for(G4int iz=2; iz<93; iz++) {
228 
229  A = G4lrint(nist->GetAtomicMassAmu(iz));
230  theA[iz] = A;
231 
232  csup = fGlauber->GetInelasticGlauberGribov(&dp, iz, A);
233  csdn = fPion->GetInelasticCrossSection(&dp, iz, theA[iz]);
234 
235  theGlauberFac[iz] = csdn/csup;
236  if(verboseLevel > 0) {
237  G4cout << "Z= " << iz << " A= " << A
238  << " factor= " << theGlauberFac[iz] << G4endl;
239  }
240  }
246 
247  if(isPiplus) {
248  dp.SetKineticEnergy(2*MeV);
249  for(G4int iz=2; iz<93; iz++) {
250  theCoulombFac[iz] = fPion->GetInelasticCrossSection(&dp, iz, theA[iz])
251  /CoulombFactor(2*MeV,iz);
252  }
253 
254  } else {
256  //fHadron->GetHadronNucleonXscNS(&dp, theProton);
257  //theCoulombFac[1] = theGlauberFac[1]*fHadron->GetInelasticHadronNucleonXsc();
258  for(G4int iz=2; iz<93; iz++) {
259  theCoulombFac[iz] = fPion->GetInelasticCrossSection(&dp, iz, theA[iz]);
260  }
261  }
262 }
263 
264 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
265 
267 {
268  G4int A = theA[Z];
269  G4double res= 0.0;
270  if(kinEnergy <= DBL_MIN) { return res; }
271  else if(A < 2) { return kinEnergy*kinEnergy; }
272 
273  G4double elog = fG4pow->log10A(6.7*kinEnergy/GeV);
274  G4double aa = A;
275 
276  // from G4ProtonInelasticCrossSection
277  G4double ff1 = 0.70 - 0.002*aa; // slope of the drop at medium energies.
278  G4double ff2 = 1.00 + 1/aa; // start of the slope.
279  G4double ff3 = 0.8 + 18/aa - 0.002*aa; // stephight
280  res = 1.0 + ff3*(1.0 - (1.0/(1+fG4pow->expA(-8*ff1*(elog + 1.37*ff2)))));
281 
282  ff1 = 1. - 1./aa - 0.001*aa; // slope of the rise
283  ff2 = 1.17 - 2.7/aa-0.0014*aa; // start of the rise
284  res /= (1 + fG4pow->expA(-8.*ff1*(elog + 2*ff2)));
285  return res;
286 }
287 
288 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
289 
290 void
292 {
293  outFile << "The Barashenkov-Glauber-Gribov cross section handles inelastic\n"
294  << "pion scattering from nuclei at all energies. The Barashenkov\n"
295  << "parameterization is used below 91 GeV and the Glauber-Gribov\n"
296  << "parameterization is used above 91 GeV.\n";
297 }
298 
299 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
G4BGGPionInelasticXS(const G4ParticleDefinition *)
virtual void CrossSectionDescription(std::ostream &) const
G4HadronNucleonXsc * fHadron
void SetMinKinEnergy(G4double value)
const G4ParticleDefinition * particle
static constexpr double MeV
Definition: G4SIunits.hh:214
void SetMaxKinEnergy(G4double value)
void SetKineticEnergy(G4double aEnergy)
G4ComponentSAIDTotalXS * fSAID
G4ComponentGGHadronNucleusXsc * fGlauber
#define G4endl
Definition: G4ios.hh:61
const char * p
Definition: xmltok.h:285
const G4ParticleDefinition * theProton
const G4String & GetParticleName() const
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
G4double GetAtomicMassAmu(const G4String &symb) const
Float_t Z
static constexpr double TeV
Definition: G4SIunits.hh:218
G4double expA(G4double A) const
Definition: G4Pow.hh:218
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
static G4Pow * GetInstance()
Definition: G4Pow.cc:57
G4ParticleDefinition * GetDefinition() const
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4double GetInelasticCrossSection(const G4DynamicParticle *aParticle, G4int Z, G4int A)
double A(double temperature)
G4double GetInelasticHadronNucleonXsc()
void SetForAllAtomsAndEnergies(G4bool val)
virtual G4double GetInelasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
G4double CoulombFactor(G4double kinEnergy, G4int Z)
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
G4UPiNuclearCrossSection * fPion
int G4lrint(double ad)
Definition: templates.hh:151
int G4int
Definition: G4Types.hh:78
G4double log10A(G4double A) const
Definition: G4Pow.hh:213
static constexpr double barn
Definition: G4SIunits.hh:105
#define DBL_MIN
Definition: templates.hh:75
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
virtual G4double GetElementCrossSection(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)
static constexpr double GeV
Definition: G4SIunits.hh:217
void BuildPhysicsTable(const G4ParticleDefinition &)
static G4NistManager * Instance()