Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4BGGPionElasticXS.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: G4BGGPionElasticXS.cc 110543 2018-05-29 13:38:54Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4BGGPionElasticXS
34 //
35 // Author: Vladimir Ivanchenko
36 //
37 // Creation date: 01.10.2003
38 // Modifications:
39 //
40 // -------------------------------------------------------------------
41 //
42 
43 #include "G4BGGPionElasticXS.hh"
44 #include "G4SystemOfUnits.hh"
47 #include "G4HadronNucleonXsc.hh"
49 #include "G4Proton.hh"
50 #include "G4PionPlus.hh"
51 #include "G4PionMinus.hh"
52 #include "G4NistManager.hh"
53 
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55 
57  : G4VCrossSectionDataSet("Barashenkov-Glauber")
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  particle = 0;
77  isPiplus = false;
78  isInitialized = false;
80 }
81 
82 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
83 
85 {
86  delete fSAID;
87  delete fHadron;
88  delete fPion;
89  delete fGlauber;
90 }
91 
92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
93 
94 G4bool
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 
113 G4double
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) {
128  cross = theCoulombFac[Z];
129  } else if(ekin > fGlauberEnergy) {
130  cross = theGlauberFac[Z]*fGlauber->GetElasticGlauberGribov(dp, Z, theA[Z]);
131  } else {
132  cross = fPion->GetElasticCrossSection(dp, Z, theA[Z]);
133  }
134  }
135  if(verboseLevel > 1) {
136  G4cout << "G4BGGPionElasticXS::GetElementCrossSection for "
137  << dp->GetDefinition()->GetParticleName()
138  << " Ekin(GeV)= " << dp->GetKineticEnergy()
139  << " in nucleus Z= " << Z << " A= " << theA[Z]
140  << " XS(b)= " << cross/barn
141  << G4endl;
142  }
143  return cross;
144 }
145 
146 G4double
148  G4int Z, G4int A,
149  const G4Isotope*,
150  const G4Element*,
151  const G4Material*)
152 {
153  // this method should be called only for Z = 1
154 
155  G4double cross = 0.0;
156  G4double ekin = dp->GetKineticEnergy();
157 
158  if(ekin <= fSAIDHighEnergyLimit) {
159  cross = fSAID->GetElasticIsotopeCrossSection(particle, ekin, 1, 1);
160  } else {
163  }
164  cross *= A;
165  /*
166  if(ekin <= fLowEnergy) {
167  cross = theCoulombFac[1];
168 
169  } else if( A < 2) {
170  fHadron->GetHadronNucleonXscNS(dp, G4Proton::Proton());
171  cross = fHadron->GetElasticHadronNucleonXsc();
172  } else {
173  fHadron->GetHadronNucleonXscNS(dp, G4Proton::Proton());
174  cross = fHadron->GetElasticHadronNucleonXsc();
175  fHadron->GetHadronNucleonXscNS(dp, G4Neutron::Neutron());
176  cross += fHadron->GetElasticHadronNucleonXsc();
177  }
178  */
179  if(verboseLevel > 1) {
180  G4cout << "G4BGGPionElasticXS::GetIsoCrossSection for "
181  << dp->GetDefinition()->GetParticleName()
182  << " Ekin(GeV)= " << dp->GetKineticEnergy()
183  << " in nucleus Z= " << Z << " A= " << A
184  << " XS(b)= " << cross/barn
185  << G4endl;
186  }
187  return cross;
188 }
189 
190 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
191 
193 {
194  if(&p == G4PionPlus::PionPlus() || &p == G4PionMinus::PionMinus()) {
195  particle = &p;
196  } else {
197  G4cout << "### G4BGGPionElasticXS WARNING: is not applicable to "
198  << p.GetParticleName()
199  << G4endl;
200  throw G4HadronicException(__FILE__, __LINE__,
201  "G4BGGPionElasticXS::BuildPhysicsTable is used for wrong particle");
202  return;
203  }
204 
205  if(isInitialized) { return; }
206  isInitialized = true;
207 
210  fHadron = new G4HadronNucleonXsc();
212 
215 
216  if(particle == G4PionPlus::PionPlus()) { isPiplus = true; }
217 
218  G4ThreeVector mom(0.0,0.0,1.0);
220 
222 
223  G4double csup, csdn;
224  G4int A;
225 
226  if(verboseLevel > 0) {
227  G4cout << "### G4BGGPionElasticXS::Initialise for "
229  }
230  for(G4int iz=2; iz<93; iz++) {
231 
232  A = G4lrint(nist->GetAtomicMassAmu(iz));
233  theA[iz] = A;
234 
235  csup = fGlauber->GetElasticGlauberGribov(&dp, iz, A);
236  csdn = fPion->GetElasticCrossSection(&dp, iz, A);
237 
238  theGlauberFac[iz] = csdn/csup;
239  if(verboseLevel > 0) {
240  G4cout << "Z= " << iz << " A= " << A
241  << " factor= " << theGlauberFac[iz] << G4endl;
242  }
243  }
244  /*
245  dp.SetKineticEnergy(fLowEnergy);
246  fHadron->GetHadronNucleonXscNS(&dp, G4Proton::Proton());
247  theCoulombFac[1] = fHadron->GetElasticHadronNucleonXsc();
248  */
251  theCoulombFac[1] =
254 
256  for(G4int iz=2; iz<93; iz++) {
257  theCoulombFac[iz] = fPion->GetElasticCrossSection(&dp, iz, theA[iz]);
258  if(verboseLevel > 0) {
259  G4cout << "Z= " << iz << " A= " << A
260  << " factor= " << theCoulombFac[iz] << G4endl;
261  }
262  }
263 }
264 
265 void
266 G4BGGPionElasticXS::CrossSectionDescription(std::ostream& outFile) const
267 {
268  outFile << "The Barashenkov-Glauber-Gribov cross section handles elastic\n"
269  << "scattering of pions from nuclei at all energies. The\n"
270  << "Barashenkov parameterization is used below 91 GeV and the\n"
271  << "Glauber-Gribov parameterization is used above 91 GeV.\n";
272 }
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
const G4ParticleDefinition * theProton
void SetMinKinEnergy(G4double value)
G4ComponentSAIDTotalXS * fSAID
static constexpr double MeV
Definition: G4SIunits.hh:214
void SetMaxKinEnergy(G4double value)
void SetKineticEnergy(G4double aEnergy)
#define G4endl
Definition: G4ios.hh:61
const char * p
Definition: xmltok.h:285
const G4String & GetParticleName() const
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
G4double GetElasticCrossSection(const G4DynamicParticle *aParticle, G4int Z, G4int A)
const G4ParticleDefinition * particle
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4double GetElasticHadronNucleonXsc()
G4double GetAtomicMassAmu(const G4String &symb) const
G4ComponentGGHadronNucleusXsc * fGlauber
Float_t Z
static constexpr double TeV
Definition: G4SIunits.hh:218
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4ParticleDefinition * GetDefinition() const
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4HadronNucleonXsc * fHadron
double A(double temperature)
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
void SetForAllAtomsAndEnergies(G4bool val)
virtual G4double GetElasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
virtual G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
G4double GetElasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
int G4lrint(double ad)
Definition: templates.hh:151
int G4int
Definition: G4Types.hh:78
static constexpr double barn
Definition: G4SIunits.hh:105
G4BGGPionElasticXS(const G4ParticleDefinition *)
G4double GetKineticEnergy() const
G4GLOB_DLL std::ostream G4cout
G4UPiNuclearCrossSection * fPion
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
static constexpr double GeV
Definition: G4SIunits.hh:217
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
void BuildPhysicsTable(const G4ParticleDefinition &)
virtual void CrossSectionDescription(std::ostream &) const
static G4NistManager * Instance()
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *)