Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BGGNucleonElasticXS Class Reference

#include <G4BGGNucleonElasticXS.hh>

+ Inheritance diagram for G4BGGNucleonElasticXS:

Public Member Functions

 G4BGGNucleonElasticXS (const G4ParticleDefinition *)
 
virtual ~G4BGGNucleonElasticXS ()
 
virtual G4bool IsElementApplicable (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)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsElementApplicable (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)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
 
virtual G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
const G4StringGetName () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4VCrossSectionDataSet
void SetName (const G4String &)
 
- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 

Detailed Description

Definition at line 65 of file G4BGGNucleonElasticXS.hh.

Constructor & Destructor Documentation

◆ G4BGGNucleonElasticXS()

G4BGGNucleonElasticXS::G4BGGNucleonElasticXS ( const G4ParticleDefinition p)

Definition at line 56 of file G4BGGNucleonElasticXS.cc.

57 : G4VCrossSectionDataSet("Barashenkov-Glauber")
58{
59 verboseLevel = 0;
60 fGlauberEnergy = 91.*GeV;
61 fLowEnergy = 14.*MeV;
62 fSAIDHighEnergyLimit = 1.3*GeV;
63 for (G4int i = 0; i < 93; ++i) {
64 theGlauberFac[i] = 0.0;
65 theCoulombFac[i] = 0.0;
66 theA[i] = 1;
67 }
68 fNucleon = 0;
69 fGlauber = 0;
70 fHadron = 0;
71 fSAID = 0;
72 particle = p;
73 theProton= G4Proton::Proton();
74 isProton = false;
75 isInitialized = false;
76}
int G4int
Definition: G4Types.hh:66
static G4Proton * Proton()
Definition: G4Proton.cc:93

◆ ~G4BGGNucleonElasticXS()

G4BGGNucleonElasticXS::~G4BGGNucleonElasticXS ( )
virtual

Definition at line 80 of file G4BGGNucleonElasticXS.cc.

81{
82 delete fHadron;
83 delete fSAID;
84}

Member Function Documentation

◆ BuildPhysicsTable()

void G4BGGNucleonElasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 180 of file G4BGGNucleonElasticXS.cc.

181{
182 if(&p == theProton || &p == G4Neutron::Neutron()) {
183 particle = &p;
184
185 } else {
186 G4cout << "### G4BGGNucleonElasticXS WARNING: is not applicable to "
187 << p.GetParticleName()
188 << G4endl;
189 throw G4HadronicException(__FILE__, __LINE__,
190 "G4BGGNucleonElasticXS::BuildPhysicsTable is used for wrong particle");
191 return;
192 }
193
194 if(isInitialized) { return; }
195 isInitialized = true;
196
199
200 fHadron = new G4HadronNucleonXsc();
201 fSAID = new G4ComponentSAIDTotalXS();
202 fNucleon->BuildPhysicsTable(*particle);
203 fGlauber->BuildPhysicsTable(*particle);
204 if(particle == theProton) {
205 isProton = true;
206 fSAIDHighEnergyLimit = 3*GeV;
207 }
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 << "### G4BGGNucleonElasticXS::Initialise for "
220 << particle->GetParticleName() << G4endl;
221 }
222
223 for(G4int iz=2; iz<93; iz++) {
224
225 A = G4lrint(nist->GetAtomicMassAmu(iz));
226 theA[iz] = A;
227
228 csup = fGlauber->GetElasticGlauberGribov(&dp, iz, A);
229 csdn = fNucleon->GetElasticCrossSection(&dp, iz);
230
231 theGlauberFac[iz] = csdn/csup;
232 if(verboseLevel > 0) {
233 G4cout << "Z= " << iz << " A= " << A
234 << " factor= " << theGlauberFac[iz] << G4endl;
235 }
236 }
237 dp.SetKineticEnergy(fSAIDHighEnergyLimit);
238 fHadron->GetHadronNucleonXscPDG(&dp, theProton);
239 theCoulombFac[1] =
240 fSAID->GetElasticIsotopeCrossSection(particle,fSAIDHighEnergyLimit,1,1)
241 /fHadron->GetElasticHadronNucleonXsc();
242 if(verboseLevel > 0) {
243 G4cout << "Z=1 A=1" << " CoulombFactor= " << theCoulombFac[1] << G4endl;
244 }
245
246 dp.SetKineticEnergy(fLowEnergy);
247 for(G4int iz=2; iz<93; iz++) {
248 theCoulombFac[iz] = fNucleon->GetElasticCrossSection(&dp, iz);
249 if(verboseLevel > 0) {
250 G4cout << "Z= " << iz << " A= " << theA[iz]
251 << " factor= " << theCoulombFac[iz] << G4endl;
252 }
253 }
254}
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
virtual G4double GetElasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)
G4VCrossSectionDataSet * GetCrossSectionDataSet(const G4String &name, G4bool warning=true)
static G4CrossSectionDataSetRegistry * Instance()
G4double GetElasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
G4double GetHadronNucleonXscPDG(const G4DynamicParticle *, const G4ParticleDefinition *)
G4double GetElasticHadronNucleonXsc()
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const
G4double GetElasticCrossSection(const G4DynamicParticle *aParticle, G4int Z)
const G4String & GetParticleName() const
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
int G4lrint(double ad)
Definition: templates.hh:163

◆ CrossSectionDescription()

void G4BGGNucleonElasticXS::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 258 of file G4BGGNucleonElasticXS.cc.

259{
260 outFile << "The Barashenkov-Glauber-Gribov cross section handles elastic\n"
261 << "scattering of protons and neutrons from nuclei using the\n"
262 << "Barashenkov parameterization below 91 GeV and the Glauber-Gribov\n"
263 << "parameterization above 91 GeV. n";
264}

◆ GetElementCrossSection()

G4double G4BGGNucleonElasticXS::GetElementCrossSection ( const G4DynamicParticle dp,
G4int  Z,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 108 of file G4BGGNucleonElasticXS.cc.

110{
111 // this method should be called only for Z > 1
112
113 G4double cross = 0.0;
114 G4double ekin = dp->GetKineticEnergy();
115 G4int Z = ZZ;
116 if(1 == Z) {
117 cross = 1.0115*GetIsoCrossSection(dp,1,1);
118 } else {
119 if(Z > 92) { Z = 92; }
120
121 if(ekin <= fLowEnergy) {
122 cross = theCoulombFac[Z];
123
124 } else if(ekin > fGlauberEnergy) {
125 cross = theGlauberFac[Z]*fGlauber->GetElasticGlauberGribov(dp, Z, theA[Z]);
126 } else {
127 cross = fNucleon->GetElasticCrossSection(dp, Z);
128 }
129 }
130
131 if(verboseLevel > 1) {
132 G4cout << "G4BGGNucleonElasticXS::GetElementCrossSection for "
134 << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
135 << " in nucleus Z= " << Z << " A= " << theA[Z]
136 << " XS(b)= " << cross/barn
137 << G4endl;
138 }
139 return cross;
140}
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const

◆ GetIsoCrossSection()

G4double G4BGGNucleonElasticXS::GetIsoCrossSection ( const G4DynamicParticle dp,
G4int  Z,
G4int  A,
const G4Isotope iso = 0,
const G4Element elm = 0,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 145 of file G4BGGNucleonElasticXS.cc.

150{
151 // this method should be called only for Z = 1
152
153 G4double cross = 0.0;
154 G4double ekin = dp->GetKineticEnergy();
155
156 if(ekin <= fSAIDHighEnergyLimit) {
157 cross = fSAID->GetElasticIsotopeCrossSection(particle, ekin, 1, 1);
158 } else {
159 fHadron->GetHadronNucleonXscPDG(dp, theProton);
160 cross = theCoulombFac[1]*fHadron->GetElasticHadronNucleonXsc();
161 }
162 // } else if(ekin <= 20*GeV) {
163 // fHadron->GetHadronNucleonXscNS(dp, G4Proton::Proton());
164 // cross = theGlauberFac[1]*fHadron->GetElasticHadronNucleonXsc();
165 cross *= A;
166
167 if(verboseLevel > 1) {
168 G4cout << "G4BGGNucleonElasticXS::GetIsoCrossSection for "
170 << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
171 << " in nucleus Z= " << Z << " A= " << A
172 << " XS(b)= " << cross/barn
173 << G4endl;
174 }
175 return cross;
176}

Referenced by GetElementCrossSection().

◆ IsElementApplicable()

G4bool G4BGGNucleonElasticXS::IsElementApplicable ( const G4DynamicParticle ,
G4int  Z,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 89 of file G4BGGNucleonElasticXS.cc.

91{
92 return true;
93}

◆ IsIsoApplicable()

G4bool G4BGGNucleonElasticXS::IsIsoApplicable ( const G4DynamicParticle ,
G4int  Z,
G4int  A,
const G4Element elm = 0,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 97 of file G4BGGNucleonElasticXS.cc.

101{
102 return (1 == Z && 2 >= A);
103}

The documentation for this class was generated from the following files: