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

#include <G4BGGNucleonElasticXS.hh>

+ Inheritance diagram for G4BGGNucleonElasticXS:

Public Member Functions

 G4BGGNucleonElasticXS (const G4ParticleDefinition *)
 
 ~G4BGGNucleonElasticXS () final
 
G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat) final
 
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 GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
 
void BuildPhysicsTable (const G4ParticleDefinition &) final
 
void CrossSectionDescription (std::ostream &) const final
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
virtual G4double ComputeCrossSectionPerElement (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, const G4Element *, const G4Material *mat=nullptr)
 
virtual G4double ComputeIsoCrossSection (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
 
virtual const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
bool ForAllAtomsAndEnergies () const
 
void SetForAllAtomsAndEnergies (G4bool val)
 
const G4StringGetName () const
 
void SetName (const G4String &nam)
 
G4VCrossSectionDataSetoperator= (const G4VCrossSectionDataSet &right)=delete
 
 G4VCrossSectionDataSet (const G4VCrossSectionDataSet &)=delete
 

Additional Inherited Members

- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 
G4String name
 

Detailed Description

Definition at line 63 of file G4BGGNucleonElasticXS.hh.

Constructor & Destructor Documentation

◆ G4BGGNucleonElasticXS()

G4BGGNucleonElasticXS::G4BGGNucleonElasticXS ( const G4ParticleDefinition * p)
explicit

Definition at line 63 of file G4BGGNucleonElasticXS.cc.

64 : G4VCrossSectionDataSet("BarashenkovGlauberGribov")
65{
66 verboseLevel = 0;
67 fGlauberEnergy = 91.*GeV;
68 fLowEnergy = 14.0*MeV;
69 fNucleon = nullptr;
70 fGlauber = nullptr;
71 fHadron = nullptr;
72
73 theProton= G4Proton::Proton();
74 isProton = (theProton == p);
75 isMaster = false;
77}
static G4Proton * Proton()
Definition G4Proton.cc:90
G4VCrossSectionDataSet(const G4String &nam="")
void SetForAllAtomsAndEnergies(G4bool val)

◆ ~G4BGGNucleonElasticXS()

G4BGGNucleonElasticXS::~G4BGGNucleonElasticXS ( )
final

Definition at line 81 of file G4BGGNucleonElasticXS.cc.

82{
83 delete fHadron;
84}

Member Function Documentation

◆ BuildPhysicsTable()

void G4BGGNucleonElasticXS::BuildPhysicsTable ( const G4ParticleDefinition & p)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 167 of file G4BGGNucleonElasticXS.cc.

168{
169 if(nullptr != fNucleon) { return; }
170 if(&p == theProton || &p == G4Neutron::Neutron()) {
171 isProton = (theProton == &p);
172
173 } else {
175 ed << "This BGG cross section is applicable only to nucleons and not to "
176 << p.GetParticleName() << G4endl;
177 G4Exception("G4BGGNucleonElasticXS::BuildPhysicsTable", "had001",
178 FatalException, ed);
179 return;
180 }
181
182 fNucleon = new G4NucleonNuclearCrossSection();
183 fGlauber = new G4ComponentGGHadronNucleusXsc();
184 fHadron = new G4HadronNucleonXsc();
185
186 fNucleon->BuildPhysicsTable(p);
187
188 if(0 == theA[0]) {
189#ifdef G4MULTITHREADED
190 G4MUTEXLOCK(&nucleonElasticXSMutex);
191 if(0 == theA[0]) {
192#endif
193 isMaster = true;
194#ifdef G4MULTITHREADED
195 }
196 G4MUTEXUNLOCK(&nucleonElasticXSMutex);
197#endif
198 } else {
199 return;
200 }
201
202 if(isMaster && 0 == theA[0]) {
203
204 theA[0] = theA[1] = 1;
205 G4ThreeVector mom(0.0,0.0,1.0);
206 G4DynamicParticle dp(theProton, mom, fGlauberEnergy);
207
209 G4double csup, csdn;
210
211 if(verboseLevel > 0) {
212 G4cout << "### G4BGGNucleonElasticXS::Initialise for "
213 << p.GetParticleName() << G4endl;
214 }
215
216 for(G4int iz=2; iz<93; ++iz) {
217 G4int A = G4lrint(nist->GetAtomicMassAmu(iz));
218 theA[iz] = A;
219
220 csup = fGlauber->GetElasticGlauberGribov(&dp, iz, A);
221 csdn = fNucleon->GetElasticCrossSection(&dp, iz);
222 theGlauberFacP[iz] = csdn/csup;
223 }
224
225 dp.SetDefinition(G4Neutron::Neutron());
226 for(G4int iz=2; iz<93; ++iz) {
227 csup = fGlauber->GetElasticGlauberGribov(&dp, iz, theA[iz]);
228 csdn = fNucleon->GetElasticCrossSection(&dp, iz);
229 theGlauberFacN[iz] = csdn/csup;
230
231 if(verboseLevel > 0) {
232 G4cout << "Z= " << iz << " A= " << theA[iz]
233 << " GFactorP= " << theGlauberFacP[iz]
234 << " GFactorN= " << theGlauberFacN[iz] << G4endl;
235 }
236 }
237
238 theCoulombFacP[0] = theCoulombFacP[1] =
239 theCoulombFacN[0] = theCoulombFacN[1] = 1.0;
240 dp.SetDefinition(theProton);
241 dp.SetKineticEnergy(fLowEnergy);
242 for(G4int iz=2; iz<93; ++iz) {
243 theCoulombFacP[iz] = fNucleon->GetElasticCrossSection(&dp, iz)
244 /CoulombFactor(fLowEnergy, iz);
245 }
246 dp.SetDefinition(G4Neutron::Neutron());
247 for(G4int iz=2; iz<93; ++iz) {
248 theCoulombFacN[iz] = fNucleon->GetElasticCrossSection(&dp, iz)
249 /CoulombFactor(fLowEnergy, iz);
250
251 if(verboseLevel > 0) {
252 G4cout << "Z= " << iz << " A= " << theA[iz]
253 << " CFactorP= " << theCoulombFacP[iz]
254 << " CFactorN= " << theCoulombFacN[iz] << G4endl;
255 }
256 }
257 }
258}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
#define G4MUTEXLOCK(mutex)
#define G4MUTEXUNLOCK(mutex)
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
const G4double A[17]
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4double GetElasticGlauberGribov(const G4DynamicParticle *, G4int Z, G4int A)
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const
void BuildPhysicsTable(const G4ParticleDefinition &) final
G4double GetElasticCrossSection(const G4DynamicParticle *aParticle, G4int Z)
const G4String & GetParticleName() const
int G4lrint(double ad)
Definition templates.hh:134

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 273 of file G4BGGNucleonElasticXS.cc.

274{
275 outFile << "The Barashenkov-Glauber-Gribov cross section handles elastic\n"
276 << "scattering of protons and neutrons from nuclei using the\n"
277 << "Barashenkov parameterization below 91 GeV and the Glauber-Gribov\n"
278 << "parameterization above 91 GeV. n";
279}

◆ GetElementCrossSection()

G4double G4BGGNucleonElasticXS::GetElementCrossSection ( const G4DynamicParticle * dp,
G4int Z,
const G4Material * mat )
finalvirtual

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 G4int Z = std::min(ZZ, 92);
115 G4double ekin = dp->GetKineticEnergy();
116 if(1 == Z) {
117 cross = 1.0115*GetIsoCrossSection(dp,1,1);
118 } else {
119 if(ekin <= fLowEnergy) {
120 cross = (isProton) ? theCoulombFacP[Z] : theCoulombFacN[Z];
121 cross *= CoulombFactor(ekin, Z);
122 } else if(ekin > fGlauberEnergy) {
123 cross = (isProton) ? theGlauberFacP[Z] : theGlauberFacN[Z];
124 cross *= fGlauber->GetElasticGlauberGribov(dp, Z, theA[Z]);
125 } else {
126 cross = fNucleon->GetElasticCrossSection(dp, Z);
127 }
128 }
129 if(verboseLevel > 1) {
130 G4cout << "G4BGGNucleonElasticXS::GetElementCrossSection for "
132 << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
133 << " in nucleus Z= " << Z << " A= " << theA[Z]
134 << " XS(b)= " << cross/barn
135 << G4endl;
136 }
137 return cross;
138}
G4double GetIsoCrossSection(const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const

◆ GetIsoCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 143 of file G4BGGNucleonElasticXS.cc.

148{
149 // this method should be called only for Z = 1
150 fHadron->HadronNucleonXscNS(dp->GetDefinition(), theProton,
151 dp->GetKineticEnergy());
152 G4double cross = A*fHadron->GetElasticHadronNucleonXsc();
153
154 if(verboseLevel > 1) {
155 G4cout << "G4BGGNucleonElasticXS::GetIsoCrossSection for "
157 << " Ekin(GeV)= " << dp->GetKineticEnergy()/CLHEP::GeV
158 << " in nucleus Z= " << Z << " A= " << A
159 << " XS(b)= " << cross/barn
160 << G4endl;
161 }
162 return cross;
163}
G4double GetElasticHadronNucleonXsc() const
G4double HadronNucleonXscNS(const G4ParticleDefinition *theParticle, const G4ParticleDefinition *nucleon, G4double ekin)

Referenced by GetElementCrossSection().

◆ IsElementApplicable()

G4bool G4BGGNucleonElasticXS::IsElementApplicable ( const G4DynamicParticle * ,
G4int Z,
const G4Material * mat )
finalvirtual

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,
const G4Material * mat )
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 97 of file G4BGGNucleonElasticXS.cc.

101{
102 return (1 == Z);
103}

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