Geant4 10.7.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 () 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 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 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 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 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 G4int GetVerboseLevel () 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
 

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 59 of file G4ParticleInelasticXS.hh.

Constructor & Destructor Documentation

◆ G4ParticleInelasticXS() [1/2]

G4ParticleInelasticXS::G4ParticleInelasticXS ( const G4ParticleDefinition part)
explicit

Definition at line 64 of file G4ParticleInelasticXS.cc.

65 : G4VCrossSectionDataSet("G4ParticleInelasticXS"),
66 highEnergyXsection(nullptr),
67 particle(part),
68 index(0),
69 isMaster(false)
70{
71 if(!part) {
72 G4Exception("G4ParticleInelasticXS::G4ParticleInelasticXS(..)","had015",
73 FatalException, "NO particle definition in constructor");
74 } else {
75 verboseLevel = 0;
76 const G4String& particleName = particle->GetParticleName();
77 if(verboseLevel > 1) {
78 G4cout << "G4ParticleInelasticXS::G4ParticleInelasticXS for "
79 << particleName << " on atoms with Z < " << MAXZINELP << G4endl;
80 }
81 if(particleName == "proton") {
82 highEnergyXsection = G4CrossSectionDataSetRegistry::Instance()->GetComponentCrossSection("Glauber-Gribov");
83 if(highEnergyXsection == nullptr) {
84 highEnergyXsection = new G4ComponentGGHadronNucleusXsc();
85 }
86 } else {
87 highEnergyXsection = G4CrossSectionDataSetRegistry::Instance()->GetComponentCrossSection("Glauber-Gribov Nucl-nucl");
88 if(highEnergyXsection == nullptr) {
89 highEnergyXsection = new G4ComponentGGNuclNuclXsc();
90 }
91 if(particleName == "deuteron") index = 1;
92 else if(particleName == "triton") index = 2;
93 else if(particleName == "He3") index = 3;
94 else if(particleName == "alpha") index = 4;
95 else {
97 ed << particleName << " is a wrong particle type";
98 G4Exception("G4ParticleInelasticXS::BuildPhysicsTable(..)","had012",
99 FatalException, ed, "");
100 }
101 }
102 }
104 temp.resize(13,0.0);
105}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
const G4int MAXZINELP
#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 285 of file G4ParticleInelasticXS.cc.

286{
287 if(verboseLevel > 0){
288 G4cout << "G4ParticleInelasticXS::BuildPhysicsTable for "
289 << p.GetParticleName() << G4endl;
290 }
291 if(&p != particle) {
293 ed << p.GetParticleName() << " is a wrong particle type -"
294 << particle->GetParticleName() << " is expected";
295 G4Exception("G4ParticleInelasticXS::BuildPhysicsTable(..)","had012",
296 FatalException, ed, "");
297 return;
298 }
299
300 if(data[index] == nullptr) {
301#ifdef G4MULTITHREADED
302 G4MUTEXLOCK(&particleInelasticXSMutex);
303 if(data[index] == nullptr) {
304#endif
305 isMaster = true;
306 data[index] = new G4ElementData();
307 data[index]->SetName(particle->GetParticleName() + "Inelastic");
308 FindDirectoryPath();
309#ifdef G4MULTITHREADED
310 }
311 G4MUTEXUNLOCK(&particleInelasticXSMutex);
312#endif
313 }
314
315 // it is possible re-initialisation for the new run
316 if(isMaster) {
317
318 // Access to elements
319 auto theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
320 size_t numOfCouples = theCoupleTable->GetTableSize();
321 for(size_t j=0; j<numOfCouples; ++j) {
322 auto mat = theCoupleTable->GetMaterialCutsCouple(j)->GetMaterial();
323 auto elmVec = mat->GetElementVector();
324 size_t numOfElem = mat->GetNumberOfElements();
325 for (size_t ie = 0; ie < numOfElem; ++ie) {
326 G4int Z = std::max(1,std::min(((*elmVec)[ie])->GetZasInt(), MAXZINELP-1));
327 if(nullptr == data[index]->GetElementData(Z)) { Initialise(Z); }
328 }
329 }
330 }
331}
#define G4MUTEXLOCK(mutex)
Definition: G4Threading.hh:251
#define G4MUTEXUNLOCK(mutex)
Definition: G4Threading.hh:254
int G4int
Definition: G4Types.hh:85
void SetName(const G4String &nam)
static G4ProductionCutsTable * GetProductionCutsTable()

Referenced by G4IonProtonCrossSection::BuildPhysicsTable().

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

◆ GetElementCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 139 of file G4ParticleInelasticXS.cc.

142{
143 G4double xs = 0.0;
144 G4double ekin = aParticle->GetKineticEnergy();
145
146 G4int Z = (ZZ >= MAXZINELP) ? MAXZINELP - 1 : ZZ;
147
148 auto pv = GetPhysicsVector(Z);
149 if(nullptr == pv) { return xs; }
150 // G4cout << "G4ParticleInelasticXS::GetCrossSection e= " << ekin
151 // << " Z= " << Z << G4endl;
152
153 if(ekin <= pv->GetMaxEnergy()) {
154 xs = pv->LogVectorValue(ekin, aParticle->GetLogKineticEnergy());
155 } else {
156 xs = coeff[Z][index]*highEnergyXsection->GetInelasticElementCrossSection(particle,
157 ekin, Z, aeff[Z]);
158 }
159
160#ifdef G4VERBOSE
161 if(verboseLevel > 1) {
162 G4cout << "ElmXS: Z= " << Z << " Ekin(MeV)= " << ekin/CLHEP::MeV
163 << " xs(bn)= " << xs/CLHEP::barn << " element data for "
164 << particle->GetParticleName()
165 << " idx= " << index << G4endl;
166 }
167#endif
168 return xs;
169}
double G4double
Definition: G4Types.hh:83
G4double GetLogKineticEnergy() const
G4double GetKineticEnergy() const
G4double GetInelasticElementCrossSection(const G4ParticleDefinition *, G4double kinEnergy, const G4Element *)

Referenced by G4IonProtonCrossSection::GetAlphaCrossSection(), G4IonProtonCrossSection::GetDeuteronCrossSection(), G4IonProtonCrossSection::GetHe3CrossSection(), G4IonProtonCrossSection::GetProtonCrossSection(), and G4IonProtonCrossSection::GetTritonCrossSection().

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

175{
176 return IsoCrossSection(aParticle->GetKineticEnergy(),
177 aParticle->GetLogKineticEnergy(),Z, A);
178}
G4double IsoCrossSection(G4double ekin, G4double logE, G4int Z, G4int A)

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

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

Referenced by GetIsoCrossSection(), G4IonProtonCrossSection::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 234 of file G4ParticleInelasticXS.cc.

236{
237 size_t nIso = anElement->GetNumberOfIsotopes();
238 const G4Isotope* iso = anElement->GetIsotope(0);
239
240 //G4cout << "SelectIsotope NIso= " << nIso << G4endl;
241 if(1 == nIso) { return iso; }
242
243 // more than 1 isotope
244 G4int Z = anElement->GetZasInt();
245 //G4cout << "SelectIsotope Z= " << Z << G4endl;
246
247 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
249 G4double sum = 0.0;
250 size_t j;
251
252 // isotope wise cross section not available
253 if(0 == 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 size_t nn = temp.size();
265 if(nn < nIso) { temp.resize(nIso, 0.); }
266
267 for (j=0; j<nIso; ++j) {
268 //G4cout << j << "-th isotope " << (*isoVector)[j]->GetN()
269 // << " abund= " << abundVector[j] << G4endl;
270 sum += abundVector[j]*IsoCrossSection(kinEnergy, logE, Z,
271 anElement->GetIsotope(j)->GetN());
272 temp[j] = sum;
273 }
274 sum *= q;
275 for (j=0; j<nIso; ++j) {
276 if(temp[j] >= sum) {
277 iso = anElement->GetIsotope(j);
278 break;
279 }
280 }
281 return iso;
282}
#define G4UniformRand()
Definition: Randomize.hh:52
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:166
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:169
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:158
G4int GetZasInt() const
Definition: G4Element.hh:131
G4int GetN() const
Definition: G4Isotope.hh:93

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