Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4BGGPionElasticXS.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//
29// GEANT4 Class file
30//
31//
32// File name: G4BGGPionElasticXS
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 01.10.2003
37// Modifications:
38//
39// -------------------------------------------------------------------
40//
41
42#include "G4BGGPionElasticXS.hh"
43#include "G4SystemOfUnits.hh"
46#include "G4HadronNucleonXsc.hh"
47#include "G4NuclearRadii.hh"
48
49#include "G4Proton.hh"
50#include "G4PionPlus.hh"
51#include "G4PionMinus.hh"
52#include "G4NistManager.hh"
54#include "G4Pow.hh"
55
56G4double G4BGGPionElasticXS::theGlauberFacPiPlus[93] = {0.0};
57G4double G4BGGPionElasticXS::theCoulombFacPiPlus[93] = {0.0};
58G4double G4BGGPionElasticXS::theGlauberFacPiMinus[93] = {0.0};
59G4double G4BGGPionElasticXS::theCoulombFacPiMinus[93] = {0.0};
60G4int G4BGGPionElasticXS::theA[93] = {0};
61
62#ifdef G4MULTITHREADED
63G4Mutex G4BGGPionElasticXS::pionElasticXSMutex = G4MUTEX_INITIALIZER;
64#endif
65
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67
69 : G4VCrossSectionDataSet("BarashenkovGlauberGribov")
70{
71 verboseLevel = 0;
72 fGlauberEnergy = 91.*GeV;
73 fLowEnergy = 20.*MeV;
74 fLowestEnergy = 1.*MeV;
75 SetMinKinEnergy(0.0);
77
78 fPion = nullptr;
79 fGlauber = nullptr;
80 fHadron = nullptr;
81
82 fG4pow = G4Pow::GetInstance();
83
84 theProton= G4Proton::Proton();
85 thePiPlus= G4PionPlus::PionPlus();
86 isPiplus = (p == thePiPlus);
87 isMaster = false;
89}
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
92
94{
95 delete fHadron;
96}
97
98//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
99
100G4bool
102 const G4Material*)
103{
104 return true;
105}
106
107//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
108
110 G4int Z, G4int,
111 const G4Element*, const G4Material*)
112{
113 return (1 == Z);
114}
115
116//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
117
120 G4int ZZ, const G4Material*)
121{
122 // this method should be called only for Z > 1
123 G4double cross = 0.0;
124 G4double ekin = std::max(dp->GetKineticEnergy(), fLowestEnergy);
125 G4int Z = std::min(ZZ, 92);
126 if(1 == Z) {
127 cross = 1.0115*GetIsoCrossSection(dp,1,1);
128 } else {
129 if(ekin <= fLowEnergy) {
130 cross = (isPiplus) ? theCoulombFacPiPlus[Z]*CoulombFactorPiPlus(ekin, Z)
131 : theCoulombFacPiMinus[Z]*FactorPiMinus(ekin);
132 } else if(ekin > fGlauberEnergy) {
133 cross = (isPiplus) ? theGlauberFacPiPlus[Z] : theGlauberFacPiMinus[Z];
134 cross *= fGlauber->GetElasticGlauberGribov(dp, Z, theA[Z]);
135 } else {
136 cross = fPion->GetElasticCrossSection(dp, Z, theA[Z]);
137 }
138 }
139 if(verboseLevel > 1) {
140 G4cout << "G4BGGPionElasticXS::GetElementCrossSection for "
142 << " Ekin(GeV)= " << dp->GetKineticEnergy()
143 << " in nucleus Z= " << Z << " A= " << theA[Z]
144 << " XS(b)= " << cross/barn
145 << G4endl;
146 }
147 return cross;
148}
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 fHadron->HadronNucleonXscNS(dp->GetDefinition(), theProton,
159 dp->GetKineticEnergy());
160 G4double cross = A*fHadron->GetElasticHadronNucleonXsc();
161
162 if(verboseLevel > 1) {
163 G4cout << "G4BGGPionElasticXS::GetIsoCrossSection for "
165 << " Ekin(GeV)= " << dp->GetKineticEnergy()
166 << " in nucleus Z= " << Z << " A= " << A
167 << " XS(b)= " << cross/barn
168 << G4endl;
169 }
170 return cross;
171}
172
173//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
174
176{
177 if(fPion) { return; }
178 if(verboseLevel > 1) {
179 G4cout << "G4BGGPionElasticXS::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("G4BGGPionElasticXS::BuildPhysicsTable", "had001",
189 FatalException, ed);
190 return;
191 }
192
193 fPion = new G4UPiNuclearCrossSection();
194 fGlauber = new G4ComponentGGHadronNucleusXsc();
195 fHadron = new G4HadronNucleonXsc();
196
197 fPion->BuildPhysicsTable(p);
198
199 if(0 == theA[0]) {
200#ifdef G4MULTITHREADED
201 G4MUTEXLOCK(&pionElasticXSMutex);
202 if(0 == theA[0]) {
203#endif
204 isMaster = true;
205#ifdef G4MULTITHREADED
206 }
207 G4MUTEXUNLOCK(&pionElasticXSMutex);
208#endif
209 } else {
210 return;
211 }
212
213 if(isMaster && 0 == theA[0]) {
214
215 theA[0] = theA[1] = 1;
216 G4ThreeVector mom(0.0,0.0,1.0);
217 G4DynamicParticle dp(thePiPlus, mom, fGlauberEnergy);
218
220
221 G4double csup, csdn;
222 for(G4int iz=2; iz<93; ++iz) {
223
224 G4int A = G4lrint(nist->GetAtomicMassAmu(iz));
225 theA[iz] = A;
226
227 csup = fGlauber->GetElasticGlauberGribov(&dp, iz, A);
228 csdn = fPion->GetElasticCrossSection(&dp, iz, A);
229 theGlauberFacPiPlus[iz] = csdn/csup;
230 }
231
233 for(G4int iz=2; iz<93; ++iz) {
234 csup = fGlauber->GetElasticGlauberGribov(&dp, iz, theA[iz]);
235 csdn = fPion->GetElasticCrossSection(&dp, iz, theA[iz]);
236 theGlauberFacPiMinus[iz] = csdn/csup;
237
238 if(verboseLevel > 0) {
239 G4cout << "Z= " << iz << " A= " << theA[iz]
240 << " factorPiPlus= " << theGlauberFacPiPlus[iz]
241 << " factorPiMinus= " << theGlauberFacPiMinus[iz]
242 << G4endl;
243 }
244 }
245 theCoulombFacPiPlus[1] = 1.0;
246 theCoulombFacPiMinus[1]= 1.0;
247 dp.SetKineticEnergy(fLowEnergy);
248 dp.SetDefinition(thePiPlus);
249 for(G4int iz=2; iz<93; ++iz) {
250 theCoulombFacPiPlus[iz] = fPion->GetElasticCrossSection(&dp, iz, theA[iz])
251 /CoulombFactorPiPlus(fLowEnergy, iz);
252 }
254 for(G4int iz=2; iz<93; ++iz) {
255 theCoulombFacPiMinus[iz] = fPion->GetElasticCrossSection(&dp, iz, theA[iz])
256 /FactorPiMinus(fLowEnergy);
257
258 if(verboseLevel > 0) {
259 G4cout << "Z= " << iz << " A= " << theA[iz]
260 << " CoulombFactorPiPlus= " << theCoulombFacPiPlus[iz]
261 << " CoulombFactorPiMinus= " << theCoulombFacPiMinus[iz]
262 << G4endl;
263 }
264 }
265 }
266}
267
268//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
269
270G4double G4BGGPionElasticXS::CoulombFactorPiPlus(G4double kinEnergy, G4int Z)
271{
272 return (kinEnergy > 0.0) ?
273 G4NuclearRadii::CoulombFactor(Z, theA[Z], thePiPlus, kinEnergy) : 0.0;
274}
275
276//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
277
278G4double G4BGGPionElasticXS::FactorPiMinus(G4double kinEnergy)
279{
280 return 1.0/std::sqrt(kinEnergy);
281}
282
283//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
284
285void
287{
288 outFile << "The Barashenkov-Glauber-Gribov cross section handles elastic\n"
289 << "scattering of pions from nuclei at all energies. The\n"
290 << "Barashenkov parameterization is used below 91 GeV and the\n"
291 << "Glauber-Gribov parameterization is used above 91 GeV.\n";
292}
293
294//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
std::mutex G4Mutex
Definition: G4Threading.hh:81
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *) final
G4BGGPionElasticXS(const G4ParticleDefinition *)
G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat) final
void CrossSectionDescription(std::ostream &) const final
G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm, const G4Material *mat) final
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4double GetElasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
void SetKineticEnergy(G4double aEnergy)
G4double GetElasticHadronNucleonXsc() const
G4double HadronNucleonXscNS(const G4ParticleDefinition *theParticle, const G4ParticleDefinition *nucleon, G4double ekin)
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()
Definition: G4PionMinus.cc:97
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:97
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
static G4Proton * Proton()
Definition: G4Proton.cc:92
G4double GetElasticCrossSection(const G4DynamicParticle *aParticle, G4int Z, G4int A) const
void BuildPhysicsTable(const G4ParticleDefinition &) final
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)
void SetForAllAtomsAndEnergies(G4bool val)
int G4lrint(double ad)
Definition: templates.hh:134