Geant4 11.3.0
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 () override=default
 
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 ()
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
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 {0}
 
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 69 of file G4ParticleInelasticXS.cc.

70 : G4VCrossSectionDataSet("G4ParticleInelasticXS"),
71 particle(part),
72 elimit(20*CLHEP::MeV)
73{
74 if (nullptr == part) {
75 G4Exception("G4ParticleInelasticXS::G4ParticleInelasticXS(..)","had015",
76 FatalException, "NO particle definition in constructor");
77 } else {
78 verboseLevel = 0;
79 const G4String& particleName = particle->GetParticleName();
80 if(verboseLevel > 1) {
81 G4cout << "G4ParticleInelasticXS::G4ParticleInelasticXS for "
82 << particleName << " on atoms with Z < " << MAXZINELP << G4endl;
83 }
85 if (particleName == "proton") {
86 highEnergyXsection = xsr->GetComponentCrossSection("Glauber-Gribov");
87 if(highEnergyXsection == nullptr) {
88 highEnergyXsection = new G4ComponentGGHadronNucleusXsc();
89 }
90 } else {
91 highEnergyXsection =
92 xsr->GetComponentCrossSection("Glauber-Gribov Nucl-nucl");
93 if(highEnergyXsection == nullptr) {
94 highEnergyXsection = new G4ComponentGGNuclNuclXsc();
95 }
96 for (index=1; index<5; ++index) {
97 if (particleName == pname[index]) { break; }
98 }
99 index = std::min(index, 4);
100 if (1 < index) { SetMaxKinEnergy(25.6*CLHEP::PeV); }
101 }
102 }
104 if (gDataDirectory.empty()) {
106 }
107 G4String ss = pname[index] + "ParticleXS";
108 SetName(ss);
109 if (data[index] == nullptr) {
110 data[index] = new G4ElementData(MAXZINELP);
111 data[index]->SetName(pname[index] + "PartInel");
112 }
113}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4CrossSectionDataSetRegistry * Instance()
static G4HadronicParameters * Instance()
const G4String & GetDirPARTICLEXS() const
G4VCrossSectionDataSet(const G4String &nam="")
void SetMaxKinEnergy(G4double value)
void SetForAllAtomsAndEnergies(G4bool val)
void SetName(const G4String &nam)

Referenced by G4ParticleInelasticXS(), and operator=().

◆ ~G4ParticleInelasticXS()

G4ParticleInelasticXS::~G4ParticleInelasticXS ( )
overridedefault

◆ 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 296 of file G4ParticleInelasticXS.cc.

297{
298 if (verboseLevel > 0){
299 G4cout << "G4ParticleInelasticXS::BuildPhysicsTable for "
300 << p.GetParticleName() << G4endl;
301 }
302 if (&p != particle) {
304 ed << p.GetParticleName() << " is a wrong particle type -"
305 << particle->GetParticleName() << " is expected";
306 G4Exception("G4ParticleInelasticXS::BuildPhysicsTable(..)","had012",
307 FatalException, ed, "");
308 return;
309 }
310
311 // it is possible re-initialisation for the new run
313
314 // prepare isotope selection
315 std::size_t nIso = temp.size();
316
317 // Access to elements
318 for ( auto const & elm : *table ) {
319 std::size_t n = elm->GetNumberOfIsotopes();
320 if (n > nIso) { nIso = n; }
321 G4int Z = std::min( elm->GetZasInt(), MAXZINELP-1);
322 if ( nullptr == (data[index])->GetElementData(Z) ) {
323 Initialise(Z);
324 }
325 }
326 temp.resize(nIso, 0.0);
327}
std::vector< G4Element * > G4ElementTable
std::ostringstream G4ExceptionDescription
int G4int
Definition G4Types.hh:85
static const G4ElementTable * GetElementTable()
Definition G4Element.cc:401
const G4String & GetParticleName() const

Referenced by G4InterfaceToXS::G4InterfaceToXS(), and ~G4ParticleInelasticXS().

◆ ComputeCrossSectionPerElement()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 148 of file G4ParticleInelasticXS.cc.

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

Referenced by ~G4ParticleInelasticXS().

◆ 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 183 of file G4ParticleInelasticXS.cc.

187{
188 return IsoCrossSection(ekin, loge, Z, A);
189}
const G4double A[17]
G4double IsoCrossSection(G4double ekin, G4double logE, G4int Z, G4int A)

Referenced by ~G4ParticleInelasticXS().

◆ 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}

Referenced by ~G4ParticleInelasticXS().

◆ ElementCrossSection()

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

Definition at line 156 of file G4ParticleInelasticXS.cc.

157{
158 G4int Z = (ZZ >= MAXZINELP) ? MAXZINELP - 1 : ZZ;
159
160 // element data is always valid pointer by construction of XS
161 auto pv = GetPhysicsVector(Z);
162
163 // set to null x-section below lowest energy in the table
164 G4double xs = 0.0;
165 if (ekin > pv->Energy(0)) {
166 xs = (ekin <= pv->GetMaxEnergy()) ? pv->LogVectorValue(ekin, loge)
167 : coeff[Z][index]*highEnergyXsection->GetInelasticElementCrossSection(particle,
168 ekin, Z, aeff[Z]);
169 }
170
171#ifdef G4VERBOSE
172 if(verboseLevel > 1) {
173 G4cout << "ElmXS: Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
174 << " xs(bn)= " << xs/CLHEP::barn << " element data for "
175 << particle->GetParticleName()
176 << " idx= " << index << G4endl;
177 }
178#endif
179 return xs;
180}
double G4double
Definition G4Types.hh:83

Referenced by ComputeCrossSectionPerElement(), GetElementCrossSection(), and ~G4ParticleInelasticXS().

◆ 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(),
144 aParticle->GetLogKineticEnergy(), Z);
145}
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const

Referenced by ~G4ParticleInelasticXS().

◆ 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 192 of file G4ParticleInelasticXS.cc.

195{
196 return IsoCrossSection(aParticle->GetKineticEnergy(),
197 aParticle->GetLogKineticEnergy(), Z, A);
198}

Referenced by ~G4ParticleInelasticXS().

◆ 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}

Referenced by ~G4ParticleInelasticXS().

◆ IsoCrossSection()

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

Definition at line 201 of file G4ParticleInelasticXS.cc.

203{
204 G4double xs = 0.0;
205 G4int Z = (ZZ >= MAXZINELP) ? MAXZINELP - 1 : ZZ;
206
207 // needed here to gurantee upload data for Z
208 auto pv = GetPhysicsVector(Z);
209
210 // compute isotope cross section if applicable
211 if (ekin <= elimit && data[index]->GetNumberOfComponents(Z) > 0) {
212 auto pviso = data[index]->GetComponentDataByID(Z, A);
213 if (pviso != nullptr && ekin > pviso->Energy(0)) {
214 xs = pviso->LogVectorValue(ekin, logE);
215#ifdef G4VERBOSE
216 if(verboseLevel > 1) {
217 G4cout << "G4ParticleInelasticXS::IsoXS: for "
218 << particle->GetParticleName() << " Ekin(MeV)= "
219 << ekin/CLHEP::MeV << " xs(b)= " << xs/CLHEP::barn
220 << " Z= " << Z << " A= " << A
221 << " idx= " << index << G4endl;
222 }
223#endif
224 return xs;
225 }
226 }
227 // use element x-section
228 if (ekin > pv->Energy(0)) {
229 xs = (ekin <= pv->GetMaxEnergy()) ? pv->LogVectorValue(ekin, logE) :
230 coeff[Z][index] *
231 highEnergyXsection->GetInelasticElementCrossSection(particle, ekin, Z, aeff[Z])
232 * A/aeff[Z];
233 }
234#ifdef G4VERBOSE
235 if(verboseLevel > 1) {
236 G4cout << "IsoXS for " << particle->GetParticleName()
237 << " Target Z= " << Z << " A= " << A
238 << " Ekin(MeV)= " << ekin/CLHEP::MeV
239 << " xs(bn)= " << xs/CLHEP::barn
240 << " idx= " << index << G4endl;
241 }
242#endif
243 return xs;
244}

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

◆ 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 246 of file G4ParticleInelasticXS.cc.

248{
249 G4int nIso = (G4int)anElement->GetNumberOfIsotopes();
250 const G4Isotope* iso = anElement->GetIsotope(0);
251
252 if (1 == nIso) { return iso; }
253
254 // more than 1 isotope
255 G4int Z = anElement->GetZasInt();
256
257 // initialisation for given Z
258 GetPhysicsVector(Z);
259
260 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
262 G4double sum = 0.0;
263 G4int j;
264
265 // isotope wise cross section not available
266 if (Z >= MAXZINELP || 0 == data[index]->GetNumberOfComponents(Z)) {
267 for (j=0; j<nIso; ++j) {
268 sum += abundVector[j];
269 if(q <= sum) {
270 iso = anElement->GetIsotope(j);
271 break;
272 }
273 }
274 return iso;
275 }
276
277 G4int nn = (G4int)temp.size();
278 if (nn < nIso) { temp.resize(nIso, 0.); }
279
280 for (j=0; j<nIso; ++j) {
281 sum += abundVector[j]*IsoCrossSection(kinEnergy, logE, Z,
282 anElement->GetIsotope(j)->GetN());
283 temp[j] = sum;
284 }
285 sum *= q;
286 for (j=0; j<nIso; ++j) {
287 if (temp[j] >= sum) {
288 iso = anElement->GetIsotope(j);
289 break;
290 }
291 }
292 return iso;
293}
#define G4UniformRand()
Definition Randomize.hh:52
G4double * GetRelativeAbundanceVector() const
Definition G4Element.hh:149
std::size_t GetNumberOfIsotopes() const
Definition G4Element.hh:143
const G4Isotope * GetIsotope(G4int iso) const
Definition G4Element.hh:151
G4int GetN() const
Definition G4Isotope.hh:83

Referenced by ~G4ParticleInelasticXS().


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