Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// -------------------------------------------------------------------
27//
28// GEANT4 Class file
29//
30//
31// File name: G4BGGPionInelasticXS
32//
33// Author: Vladimir Ivanchenko
34//
35// Creation date: 01.10.2003
36// Modifications:
37//
38// -------------------------------------------------------------------
39//
40
42#include "G4SystemOfUnits.hh"
45#include "G4HadronNucleonXsc.hh"
46#include "G4NuclearRadii.hh"
47
48#include "G4Proton.hh"
49#include "G4PionPlus.hh"
50#include "G4PionMinus.hh"
51#include "G4NistManager.hh"
52#include "G4Pow.hh"
54
55G4double G4BGGPionInelasticXS::theGlauberFacPiPlus[93] = {0.0};
56G4double G4BGGPionInelasticXS::theGlauberFacPiMinus[93] = {0.0};
57G4double G4BGGPionInelasticXS::theLowEPiPlus[93] = {0.0};
58G4double G4BGGPionInelasticXS::theLowEPiMinus[93] = {0.0};
59G4int G4BGGPionInelasticXS::theA[93] = {0};
60
62 : G4VCrossSectionDataSet("BarashenkovGlauberGribov")
63{
64 verboseLevel = 0;
65 fGlauberEnergy = 91.*CLHEP::GeV;
66 fLowEnergy = 20.*CLHEP::MeV;
67 fLowestEnergy = 1.*CLHEP::MeV;
68 SetMinKinEnergy(0.0);
70
71 fPion = new G4UPiNuclearCrossSection();
72 fGlauber = new G4ComponentGGHadronNucleusXsc();
73 fHadron = new G4HadronNucleonXsc();
74
75 fG4pow = G4Pow::GetInstance();
76
77 theProton = G4Proton::Proton();
78 thePiPlus = G4PionPlus::PionPlus();
79 isPiplus = (p == thePiPlus);
81
82 if (0 == theA[0]) { Initialise(); }
83}
84
85//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
86
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
93
94G4bool
100
101//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
102
104 G4int Z, G4int,
105 const G4Element*,
106 const G4Material*)
107{
108 return (1 == Z);
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 = std::max(dp->GetKineticEnergy(), fLowestEnergy);
121 G4int Z = std::min(ZZ, 92);
122
123 if(1 == Z) {
124 cross = 1.0115*GetIsoCrossSection(dp,1,1);
125 } else if(ekin < fLowEnergy) {
126 cross = (isPiplus) ? theLowEPiPlus[Z]*CoulombFactorPiPlus(ekin, Z)
127 : theLowEPiMinus[Z]*FactorPiMinus(ekin);
128 } else if(ekin > fGlauberEnergy) {
129 cross = (isPiplus) ? theGlauberFacPiPlus[Z] : theGlauberFacPiMinus[Z];
130 cross *= fGlauber->GetInelasticGlauberGribov(dp, Z, theA[Z]);
131 } else {
132 cross = fPion->GetInelasticCrossSection(dp, Z, theA[Z]);
133 }
134#ifdef G4VERBOSE
135 if(verboseLevel > 1) {
136 G4cout << "G4BGGPionInelasticXS::GetCrossSection for "
138 << " Ekin(GeV)= " << dp->GetKineticEnergy()
139 << " in nucleus Z= " << Z << " A= " << theA[Z]
140 << " XS(b)= " << cross/barn
141 << G4endl;
142 }
143#endif
144 return cross;
145}
146
147//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
148
151 G4int, G4int A,
152 const G4Isotope*,
153 const G4Element*,
154 const G4Material*)
155{
156 // this method should be called only for Z = 1
157 fHadron->HadronNucleonXscNS(dp->GetDefinition(), theProton,
158 dp->GetKineticEnergy());
159 G4double cross = A*fHadron->GetInelasticHadronNucleonXsc();
160
161#ifdef G4VERBOSE
162 if(verboseLevel > 1) {
163 G4cout << "G4BGGPionInelasticXS::GetCrossSection for "
165 << " Ekin(GeV)= " << dp->GetKineticEnergy()
166 << " in nucleus Z=1 A=" << A
167 << " XS(b)= " << cross/barn
168 << G4endl;
169 }
170#endif
171 return cross;
172}
173
174//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
175
177{
178 if(verboseLevel > 1) {
179 G4cout << "G4BGGPionInelasticXS::BuildPhysicsTable for "
180 << p.GetParticleName() << G4endl;
181 }
182 if(&p == G4PionPlus::PionPlus() || &p == G4PionMinus::PionMinus()) {
183 isPiplus = (&p == G4PionPlus::PionPlus());
184 } else {
186 ed << "This BGG cross section is applicable only to pions and not to "
187 << p.GetParticleName() << G4endl;
188 G4Exception("G4BGGPionInelasticXS::BuildPhysicsTable", "had001",
189 FatalException, ed);
190 }
191}
192
193//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
194
195void G4BGGPionInelasticXS::Initialise()
196{
197 theA[0] = theA[1] = 1;
198 G4ThreeVector mom(0.0,0.0,1.0);
199 G4DynamicParticle dp(thePiPlus, mom, fGlauberEnergy);
200
202 G4double csup, csdn;
203
204 for (G4int iz=2; iz<93; ++iz) {
205 G4int A = G4lrint(nist->GetAtomicMassAmu(iz));
206 theA[iz] = A;
207
208 csup = fGlauber->GetInelasticGlauberGribov(&dp, iz, A);
209 csdn = fPion->GetInelasticCrossSection(&dp, iz, A);
210 theGlauberFacPiPlus[iz] = csdn/csup;
211 }
212
213 dp.SetDefinition(G4PionMinus::PionMinus());
214 for (G4int iz=2; iz<93; ++iz) {
215 csup = fGlauber->GetInelasticGlauberGribov(&dp, iz, theA[iz]);
216 csdn = fPion->GetInelasticCrossSection(&dp, iz, theA[iz]);
217 theGlauberFacPiMinus[iz] = csdn/csup;
218
219 if(verboseLevel > 1) {
220 G4cout << "Z= " << iz << " A= " << theA[iz]
221 << " factorPiPlus= " << theGlauberFacPiPlus[iz]
222 << " factorPiMinus= " << theGlauberFacPiMinus[iz]
223 << G4endl;
224 }
225 }
226
227 theLowEPiPlus[1] = theLowEPiMinus[1]= 1.0;
228 dp.SetDefinition(thePiPlus);
229 dp.SetKineticEnergy(fLowEnergy);
230 for (G4int iz=2; iz<93; ++iz) {
231 theLowEPiPlus[iz] = fPion->GetInelasticCrossSection(&dp, iz, theA[iz])
232 /CoulombFactorPiPlus(fLowEnergy, iz);
233 }
234
235 dp.SetDefinition(G4PionMinus::PionMinus());
236 for (G4int iz=2; iz<93; ++iz) {
237 theLowEPiMinus[iz] = fPion->GetInelasticCrossSection(&dp, iz, theA[iz])
238 /FactorPiMinus(fLowEnergy);
239
240 if (verboseLevel > 1) {
241 G4cout << "Z= " << iz << " A= " << theA[iz]
242 << " LowEtorPiPlus= " << theLowEPiPlus[iz]
243 << " LowEtorPiMinus= " << theLowEPiMinus[iz]
244 << G4endl;
245 }
246 }
247}
248
249//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
250
251G4double G4BGGPionInelasticXS::CoulombFactorPiPlus(G4double kinEnergy, G4int Z)
252{
253 return (kinEnergy > 0.0) ?
254 G4NuclearRadii::CoulombFactor(Z, theA[Z], thePiPlus, kinEnergy) : 0.0;
255}
256
257//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
258
259G4double G4BGGPionInelasticXS::FactorPiMinus(G4double kinEnergy)
260{
261 return 1.0/std::sqrt(kinEnergy);
262}
263
264//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
265
266void
268{
269 outFile << "The Barashenkov-Glauber-Gribov cross section handles inelastic\n"
270 << "pion scattering from nuclei at all energies. The Barashenkov\n"
271 << "parameterization is used below 91 GeV and the Glauber-Gribov\n"
272 << "parameterization is used above 91 GeV.\n";
273}
274
275//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void CrossSectionDescription(std::ostream &) const final
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat) final
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
G4BGGPionInelasticXS(const G4ParticleDefinition *)
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm, const G4Material *mat) final
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat) final
G4double GetInelasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4HadronicParameters * Instance()
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const
static G4double CoulombFactor(const G4ParticleDefinition *theParticle, const G4ParticleDefinition *nucleon, G4double ekin)
const G4String & GetParticleName() const
static G4PionMinus * PionMinus()
static G4PionPlus * PionPlus()
Definition G4PionPlus.cc:93
static G4Pow * GetInstance()
Definition G4Pow.cc:41
static G4Proton * Proton()
Definition G4Proton.cc:90
G4double GetInelasticCrossSection(const G4DynamicParticle *aParticle, G4int Z, G4int A) const
G4VCrossSectionDataSet(const G4String &nam="")
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)
void SetForAllAtomsAndEnergies(G4bool val)
int G4lrint(double ad)
Definition templates.hh:134