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

#include <G4ParticleInelasticXS.hh>

+ Inheritance diagram for G4ParticleInelasticXS:

Public Member Functions

 G4ParticleInelasticXS (const G4ParticleDefinition *)
 
 ~G4ParticleInelasticXS () final
 
G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr) final
 
G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=nullptr, const G4Material *mat=nullptr) final
 
G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr) final
 
G4double ComputeCrossSectionPerElement (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, const G4Element *, const G4Material *) final
 
G4double ComputeIsoCrossSection (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, G4int Z, G4int A, const G4Isotope *iso, const G4Element *elm, 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
 
const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE) final
 
void BuildPhysicsTable (const G4ParticleDefinition &) final
 
void CrossSectionDescription (std::ostream &) const final
 
G4double ElementCrossSection (G4double kinEnergy, G4double loge, G4int Z)
 
G4double IsoCrossSection (G4double ekin, G4double logE, G4int Z, G4int A)
 
G4ParticleInelasticXSoperator= (const G4ParticleInelasticXS &right)=delete
 
 G4ParticleInelasticXS (const G4ParticleInelasticXS &)=delete
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
 
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 GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, 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 BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
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 56 of file G4ParticleInelasticXS.hh.

Constructor & Destructor Documentation

◆ G4ParticleInelasticXS() [1/2]

G4ParticleInelasticXS::G4ParticleInelasticXS ( const G4ParticleDefinition part)
explicit

Definition at line 67 of file G4ParticleInelasticXS.cc.

68 : G4VCrossSectionDataSet("G4ParticleInelasticXS"),
69 particle(part),
70 elimit(20*CLHEP::MeV)
71{
72 if(nullptr == part) {
73 G4Exception("G4ParticleInelasticXS::G4ParticleInelasticXS(..)","had015",
74 FatalException, "NO particle definition in constructor");
75 } else {
76 verboseLevel = 0;
77 const G4String& particleName = particle->GetParticleName();
78 if(verboseLevel > 1) {
79 G4cout << "G4ParticleInelasticXS::G4ParticleInelasticXS for "
80 << particleName << " on atoms with Z < " << MAXZINELP << G4endl;
81 }
82 if(particleName == "proton") {
83 highEnergyXsection = G4CrossSectionDataSetRegistry::Instance()->GetComponentCrossSection("Glauber-Gribov");
84 if(highEnergyXsection == nullptr) {
85 highEnergyXsection = new G4ComponentGGHadronNucleusXsc();
86 }
87 } else {
88 highEnergyXsection = G4CrossSectionDataSetRegistry::Instance()->GetComponentCrossSection("Glauber-Gribov Nucl-nucl");
89 if(highEnergyXsection == nullptr) {
90 highEnergyXsection = new G4ComponentGGNuclNuclXsc();
91 }
92 if(particleName == "deuteron") index = 1;
93 else if(particleName == "triton") index = 2;
94 else if(particleName == "He3") index = 3;
95 else if(particleName == "alpha") index = 4;
96 else {
98 ed << particleName << " is a wrong particle type";
99 G4Exception("G4ParticleInelasticXS::BuildPhysicsTable(..)","had012",
100 FatalException, ed, "");
101 }
102 }
103 }
105}
@ 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 G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4VComponentCrossSection * GetComponentCrossSection(const G4String &name)
static G4CrossSectionDataSetRegistry * Instance()
const G4String & GetParticleName() const
void SetForAllAtomsAndEnergies(G4bool val)

◆ ~G4ParticleInelasticXS()

G4ParticleInelasticXS::~G4ParticleInelasticXS ( )
final

Definition at line 107 of file G4ParticleInelasticXS.cc.

108{
109 if(isMaster) {
110 delete data[index];
111 data[index] = nullptr;
112 }
113}

◆ G4ParticleInelasticXS() [2/2]

G4ParticleInelasticXS::G4ParticleInelasticXS ( const G4ParticleInelasticXS )
delete

Member Function Documentation

◆ BuildPhysicsTable()

void G4ParticleInelasticXS::BuildPhysicsTable ( const G4ParticleDefinition p)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 283 of file G4ParticleInelasticXS.cc.

284{
285 if(verboseLevel > 0){
286 G4cout << "G4ParticleInelasticXS::BuildPhysicsTable for "
287 << p.GetParticleName() << G4endl;
288 }
289 if(&p != particle) {
291 ed << p.GetParticleName() << " is a wrong particle type -"
292 << particle->GetParticleName() << " is expected";
293 G4Exception("G4ParticleInelasticXS::BuildPhysicsTable(..)","had012",
294 FatalException, ed, "");
295 return;
296 }
297
298 G4int fact = (p.GetParticleName() == "proton") ? 1 : 256;
299 SetMaxKinEnergy(G4HadronicParameters::Instance()->GetMaxEnergy() * fact);
300
301 if(data[index] == nullptr) {
302 G4AutoLock l(&pInelasticXSMutex);
303 if(data[index] == nullptr) {
304 isMaster = true;
305 data[index] = new G4ElementData();
306 data[index]->SetName(particle->GetParticleName() + "Inelastic");
307 FindDirectoryPath();
308 }
309 l.unlock();
310 }
311
312 // it is possible re-initialisation for the new run
314 if(isMaster) {
315
316 // Access to elements
317 for ( auto & elm : *table ) {
318 G4int Z = std::max( 1, std::min( elm->GetZasInt(), MAXZINELP-1) );
319 if ( nullptr == (data[index])->GetElementData(Z) ) { Initialise(Z); }
320 }
321 }
322 // prepare isotope selection
323 std::size_t nIso = temp.size();
324 for ( auto & elm : *table ) {
325 std::size_t n = elm->GetNumberOfIsotopes();
326 if(n > nIso) { nIso = n; }
327 }
328 temp.resize(nIso, 0.0);
329}
std::vector< G4Element * > G4ElementTable
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
void SetName(const G4String &nam)
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
static G4HadronicParameters * Instance()
void SetMaxKinEnergy(G4double value)

◆ ComputeCrossSectionPerElement()

G4double G4ParticleInelasticXS::ComputeCrossSectionPerElement ( G4double  kinEnergy,
G4double  loge,
const G4ParticleDefinition ,
const G4Element elm,
const G4Material  
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 147 of file G4ParticleInelasticXS.cc.

151{
152 return ElementCrossSection(ekin, loge, elm->GetZasInt());
153}
G4int GetZasInt() const
Definition: G4Element.hh:132
G4double ElementCrossSection(G4double kinEnergy, G4double loge, G4int Z)

◆ ComputeIsoCrossSection()

G4double G4ParticleInelasticXS::ComputeIsoCrossSection ( G4double  kinEnergy,
G4double  loge,
const G4ParticleDefinition ,
G4int  Z,
G4int  A,
const G4Isotope iso,
const G4Element elm,
const G4Material mat 
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 176 of file G4ParticleInelasticXS.cc.

180{
181 return IsoCrossSection(ekin, loge, Z, A);
182}
const G4double A[17]
G4double IsoCrossSection(G4double ekin, G4double logE, G4int Z, G4int A)

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 115 of file G4ParticleInelasticXS.cc.

116{
117 outFile << "G4ParticleInelasticXS calculates n, p, d, t, he3, he4 inelastic\n"
118 << "cross section on nuclei using data from the high precision\n"
119 << "neutron database. These data are simplified and smoothed over\n"
120 << "the resonance region in order to reduce CPU time.\n"
121 << "For high energy Glauber-Gribov cross section model is used\n";
122}

◆ ElementCrossSection()

G4double G4ParticleInelasticXS::ElementCrossSection ( G4double  kinEnergy,
G4double  loge,
G4int  Z 
)

Definition at line 155 of file G4ParticleInelasticXS.cc.

156{
157 G4int Z = (ZZ >= MAXZINELP) ? MAXZINELP - 1 : ZZ;
158 auto pv = GetPhysicsVector(Z);
159
160 G4double xs = (ekin <= pv->GetMaxEnergy()) ? pv->LogVectorValue(ekin, loge)
161 : coeff[Z][index]*highEnergyXsection->GetInelasticElementCrossSection(particle,
162 ekin, Z, aeff[Z]);
163
164#ifdef G4VERBOSE
165 if(verboseLevel > 1) {
166 G4cout << "ElmXS: Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
167 << " xs(bn)= " << xs/CLHEP::barn << " element data for "
168 << particle->GetParticleName()
169 << " idx= " << index << G4endl;
170 }
171#endif
172 return xs;
173}
double G4double
Definition: G4Types.hh:83
G4double GetInelasticElementCrossSection(const G4ParticleDefinition *, G4double kinEnergy, const G4Element *)

Referenced by ComputeCrossSectionPerElement(), and GetElementCrossSection().

◆ GetElementCrossSection()

G4double G4ParticleInelasticXS::GetElementCrossSection ( const G4DynamicParticle aParticle,
G4int  Z,
const G4Material mat = nullptr 
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 140 of file G4ParticleInelasticXS.cc.

142{
143 return ElementCrossSection(aParticle->GetKineticEnergy(), aParticle->GetLogKineticEnergy(), Z);
144}
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const

◆ GetIsoCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 185 of file G4ParticleInelasticXS.cc.

188{
189 return IsoCrossSection(aParticle->GetKineticEnergy(),
190 aParticle->GetLogKineticEnergy(), Z, A);
191}

◆ IsElementApplicable()

G4bool G4ParticleInelasticXS::IsElementApplicable ( const G4DynamicParticle ,
G4int  Z,
const G4Material mat = nullptr 
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 125 of file G4ParticleInelasticXS.cc.

127{
128 return true;
129}

◆ IsIsoApplicable()

G4bool G4ParticleInelasticXS::IsIsoApplicable ( const G4DynamicParticle ,
G4int  Z,
G4int  A,
const G4Element elm = nullptr,
const G4Material mat = nullptr 
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 132 of file G4ParticleInelasticXS.cc.

135{
136 return true;
137}

◆ IsoCrossSection()

G4double G4ParticleInelasticXS::IsoCrossSection ( G4double  ekin,
G4double  logE,
G4int  Z,
G4int  A 
)

Definition at line 194 of file G4ParticleInelasticXS.cc.

196{
197 G4double xs = 0.0;
198 G4int Z = (ZZ >= MAXZINELP) ? MAXZINELP - 1 : ZZ;
199 auto pv = GetPhysicsVector(Z);
200
201 // compute isotope cross section if applicable
202 if(ekin <= elimit && amin[Z] < amax[Z] && A >= amin[Z] && A <= amax[Z]) {
203 auto pviso = data[index]->GetComponentDataByIndex(Z, A - amin[Z]);
204 if(pviso != nullptr) {
205 xs = pviso->LogVectorValue(ekin, logE);
206#ifdef G4VERBOSE
207 if(verboseLevel > 1) {
208 G4cout << "G4ParticleInelasticXS::IsoXS: for "
209 << particle->GetParticleName() << " Ekin(MeV)= "
210 << ekin/CLHEP::MeV << " xs(b)= " << xs/CLHEP::barn
211 << " Z= " << Z << " A= " << A
212 << " idx= " << index << G4endl;
213 }
214#endif
215 return xs;
216 }
217 }
218 // use element x-section
219 xs = (ekin <= pv->GetMaxEnergy()) ? pv->LogVectorValue(ekin, logE)
220 : coeff[Z][index] *
221 highEnergyXsection->GetInelasticElementCrossSection(particle,
222 ekin, Z, aeff[Z]);
223 xs *= A/aeff[Z];
224#ifdef G4VERBOSE
225 if(verboseLevel > 1) {
226 G4cout << "IsoXS for " << particle->GetParticleName()
227 << " Target Z= " << Z << " A= " << A
228 << " Ekin(MeV)= " << ekin/CLHEP::MeV
229 << " xs(bn)= " << xs/CLHEP::barn
230 << " idx= " << index << G4endl;
231 }
232#endif
233 return xs;
234}
G4PhysicsVector * GetComponentDataByIndex(G4int Z, G4int idx)
G4double LogVectorValue(const G4double energy, const G4double theLogEnergy) const

Referenced by ComputeIsoCrossSection(), GetIsoCrossSection(), and SelectIsotope().

◆ operator=()

G4ParticleInelasticXS & G4ParticleInelasticXS::operator= ( const G4ParticleInelasticXS right)
delete

◆ SelectIsotope()

const G4Isotope * G4ParticleInelasticXS::SelectIsotope ( const G4Element anElement,
G4double  kinEnergy,
G4double  logE 
)
finalvirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 236 of file G4ParticleInelasticXS.cc.

238{
239 G4int nIso = (G4int)anElement->GetNumberOfIsotopes();
240 const G4Isotope* iso = anElement->GetIsotope(0);
241
242 if(1 == nIso) { return iso; }
243
244 // more than 1 isotope
245 G4int Z = anElement->GetZasInt();
246
247 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
249 G4double sum = 0.0;
250 G4int j;
251
252 // isotope wise cross section not available
253 if(amax[Z] == amin[Z] || Z >= MAXZINELP) {
254 for (j=0; j<nIso; ++j) {
255 sum += abundVector[j];
256 if(q <= sum) {
257 iso = anElement->GetIsotope(j);
258 break;
259 }
260 }
261 return iso;
262 }
263
264 G4int nn = (G4int)temp.size();
265 if(nn < nIso) { temp.resize(nIso, 0.); }
266
267 for (j=0; j<nIso; ++j) {
268 sum += abundVector[j]*IsoCrossSection(kinEnergy, logE, Z,
269 anElement->GetIsotope(j)->GetN());
270 temp[j] = sum;
271 }
272 sum *= q;
273 for (j=0; j<nIso; ++j) {
274 if(temp[j] >= sum) {
275 iso = anElement->GetIsotope(j);
276 break;
277 }
278 }
279 return iso;
280}
#define G4UniformRand()
Definition: Randomize.hh:52
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
G4int GetN() const
Definition: G4Isotope.hh:93

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