Geant4-11
Public Member Functions | Protected Attributes | Private Attributes
G4ParticleHPInelasticData Class Reference

#include <G4ParticleHPInelasticData.hh>

Inheritance diagram for G4ParticleHPInelasticData:
G4VCrossSectionDataSet

Public Member Functions

void BuildPhysicsTable (const G4ParticleDefinition &)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
virtual void CrossSectionDescription (std::ostream &) const
 
void DumpPhysicsTable (const G4ParticleDefinition &)
 
bool ForAllAtomsAndEnergies () const
 
 G4ParticleHPInelasticData (G4ParticleDefinition *projectile=G4Neutron::Neutron())
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, G4double aT)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
 
G4double GetIsoCrossSection (const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *)
 
G4double GetMaxKinEnergy () const
 
G4double GetMinKinEnergy () const
 
const G4StringGetName () const
 
G4ParticleDefinitionGetProjectile ()
 
G4int GetVerboseLevel () const
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
 
G4bool IsIsoApplicable (const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *)
 
virtual const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE)
 
void SetForAllAtomsAndEnergies (G4bool val)
 
void SetMaxKinEnergy (G4double value)
 
void SetMinKinEnergy (G4double value)
 
void SetName (const G4String &nam)
 
void SetVerboseLevel (G4int)
 
 ~G4ParticleHPInelasticData ()
 

Protected Attributes

G4String name
 
G4int verboseLevel
 

Private Attributes

const G4Elementelement_cache
 
G4bool instanceOfWorker
 
G4bool isForAllAtomsAndEnergies
 
G4double ke_cache
 
const G4Materialmaterial_cache
 
G4double maxKinEnergy
 
G4double minKinEnergy
 
G4CrossSectionDataSetRegistryregistry
 
G4PhysicsTabletheCrossSections
 
G4ParticleHPDatatheHPData
 
G4ParticleDefinitiontheProjectile
 
G4double xs_cache
 

Detailed Description

Definition at line 53 of file G4ParticleHPInelasticData.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPInelasticData()

G4ParticleHPInelasticData::G4ParticleHPInelasticData ( G4ParticleDefinition projectile = G4Neutron::Neutron())

Definition at line 47 of file G4ParticleHPInelasticData.cc.

49{
50 const char* dataDirVariable;
51 G4String particleName;
52 if( projectile == G4Neutron::Neutron() ) {
53 dataDirVariable = "G4NEUTRONHPDATA";
54 }else if( projectile == G4Proton::Proton() ) {
55 dataDirVariable = "G4PROTONHPDATA";
56 particleName = "Proton";
57 }else if( projectile == G4Deuteron::Deuteron() ) {
58 dataDirVariable = "G4DEUTERONHPDATA";
59 particleName = "Deuteron";
60 }else if( projectile == G4Triton::Triton() ) {
61 dataDirVariable = "G4TRITONHPDATA";
62 particleName = "Triton";
63 }else if( projectile == G4He3::He3() ) {
64 dataDirVariable = "G4HE3HPDATA";
65 particleName = "He3";
66 }else if( projectile == G4Alpha::Alpha() ) {
67 dataDirVariable = "G4ALPHAHPDATA";
68 particleName = "Alpha";
69 } else {
70 G4String message("G4ParticleHPInelasticData may only be called for neutron, proton, deuteron, triton, He3 or alpha, while it is called for " + projectile->GetParticleName());
71 throw G4HadronicException(__FILE__, __LINE__,message.c_str());
72 }
73 // G4cout << this << " G4ParticleHPInelasticData::G4ParticleHPInelasticData " << projectile->GetParticleName() << " DATADIR " << dataDirVariable << G4endl;//GDEB
74 G4String dataName = projectile->GetParticleName()+"HPInelasticXS";
75 dataName.at(0) = toupper(dataName.at(0)) ;
76 SetName( dataName );
77
78 if ( !std::getenv(dataDirVariable) && !std::getenv( "G4PARTICLEHPDATA" ) ){
79 G4String message("Please setenv G4PARTICLEHPDATA (recommended) or, at least setenv " +
80 G4String(dataDirVariable) + " to point to the " + projectile->GetParticleName() + " cross-section files.");
81 throw G4HadronicException(__FILE__, __LINE__,message.c_str());
82 }
83
84 G4String dirName;
85 if ( std::getenv(dataDirVariable) ) {
86 dirName = std::getenv(dataDirVariable);
87 } else {
88 G4String baseName = std::getenv( "G4PARTICLEHPDATA" );
89 dirName = baseName + "/" + particleName;
90 }
91 #ifdef G4VERBOSE
93 G4cout << "@@@ G4ParticleHPInelasticData instantiated for particle " << projectile->GetParticleName() << " data directory variable is " << dataDirVariable << " pointing to " << dirName << G4endl;
94 }
95 #endif
96
99
101 theProjectile=projectile;
102
103 theHPData = NULL;
104 instanceOfWorker = false;
107 } else {
108 instanceOfWorker = true;
109 }
110 element_cache = NULL;
111 material_cache = NULL;
112 ke_cache = 0.0;
113 xs_cache = 0.0;
114}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4Alpha * Alpha()
Definition: G4Alpha.cc:88
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:93
static G4HadronicParameters * Instance()
static G4He3 * He3()
Definition: G4He3.cc:93
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
const G4String & GetParticleName() const
G4ParticleDefinition * theProjectile
static G4Proton * Proton()
Definition: G4Proton.cc:92
static G4Triton * Triton()
Definition: G4Triton.cc:93
G4VCrossSectionDataSet(const G4String &nam="")
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)
void SetName(const G4String &nam)
static constexpr double MeV
G4bool IsMasterThread()
Definition: G4Threading.cc:124

References G4Alpha::Alpha(), G4Deuteron::Deuteron(), element_cache, G4cout, G4endl, G4ParticleDefinition::GetParticleName(), GetVerboseLevel(), G4He3::He3(), G4HadronicParameters::Instance(), instanceOfWorker, G4Threading::IsMasterThread(), ke_cache, material_cache, CLHEP::MeV, G4Neutron::Neutron(), G4Proton::Proton(), G4VCrossSectionDataSet::SetMaxKinEnergy(), G4VCrossSectionDataSet::SetMinKinEnergy(), G4VCrossSectionDataSet::SetName(), theCrossSections, theHPData, theProjectile, G4Triton::Triton(), and xs_cache.

◆ ~G4ParticleHPInelasticData()

G4ParticleHPInelasticData::~G4ParticleHPInelasticData ( )

Definition at line 116 of file G4ParticleHPInelasticData.cc.

117{
118 if ( theCrossSections != NULL && instanceOfWorker != true ) {
120 delete theCrossSections;
121 theCrossSections = NULL;
122 }
123 if ( theHPData != NULL && instanceOfWorker != true ) {
124 delete theHPData;
125 theHPData = NULL;
126 }
127}
void clearAndDestroy()

References G4PhysicsTable::clearAndDestroy(), instanceOfWorker, theCrossSections, and theHPData.

Member Function Documentation

◆ BuildPhysicsTable()

void G4ParticleHPInelasticData::BuildPhysicsTable ( const G4ParticleDefinition projectile)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 169 of file G4ParticleHPInelasticData.cc.

170{
171 // if(&projectile!=G4Neutron::Neutron())
172 // throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
173
176 return;
177 } else {
178 if ( theHPData == NULL ) theHPData = G4ParticleHPData::Instance( const_cast<G4ParticleDefinition*> ( &projectile ) );
179 }
180
181 size_t numberOfElements = G4Element::GetNumberOfElements();
182// theCrossSections = new G4PhysicsTable( numberOfElements );
183// TKDB
184 //if ( theCrossSections == 0 )
185 //{ theCrossSections = new G4PhysicsTable( numberOfElements ); }
186 if ( theCrossSections == NULL )
187 theCrossSections = new G4PhysicsTable( numberOfElements );
188 else
190
191 // make a PhysicsVector for each element
192
193 //G4ParticleHPData* hpData = new G4ParticleHPData(projectile); //NEW
194 static G4ThreadLocal G4ElementTable *theElementTable = 0 ;
195 if (!theElementTable) theElementTable= G4Element::GetElementTable();
196 for( size_t i=0; i<numberOfElements; ++i )
197 {
198 //NEW G4PhysicsVector* physVec = G4ParticleHPData::
199 //NEW Instance(projectile, dataDirVariable)->MakePhysicsVector((*theElementTable)[i], this);
200 //G4PhysicsVector* physVec = hpData->MakePhysicsVector((*theElementTable)[i], this);
201 G4PhysicsVector* physVec = theHPData->MakePhysicsVector((*theElementTable)[i], this);
202 theCrossSections->push_back(physVec);
203 }
204
206}
std::vector< G4Element * > G4ElementTable
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:397
static size_t GetNumberOfElements()
Definition: G4Element.cc:404
static G4ParticleHPData * Instance(G4ParticleDefinition *projectile)
G4PhysicsVector * MakePhysicsVector(G4Element *thE, G4ParticleHPFissionData *theP)
G4PhysicsTable * GetInelasticCrossSections(const G4ParticleDefinition *)
void RegisterInelasticCrossSections(const G4ParticleDefinition *, G4PhysicsTable *)
static G4ParticleHPManager * GetInstance()
void push_back(G4PhysicsVector *)
G4bool IsWorkerThread()
Definition: G4Threading.cc:123
#define G4ThreadLocal
Definition: tls.hh:77

References G4PhysicsTable::clearAndDestroy(), G4ThreadLocal, G4Element::GetElementTable(), G4ParticleHPManager::GetInelasticCrossSections(), G4ParticleHPManager::GetInstance(), G4Element::GetNumberOfElements(), G4ParticleHPData::Instance(), G4Threading::IsWorkerThread(), G4ParticleHPData::MakePhysicsVector(), G4PhysicsTable::push_back(), G4ParticleHPManager::RegisterInelasticCrossSections(), theCrossSections, and theHPData.

◆ ComputeCrossSection()

G4double G4VCrossSectionDataSet::ComputeCrossSection ( const G4DynamicParticle part,
const G4Element elm,
const G4Material mat = nullptr 
)
inherited

Definition at line 81 of file G4VCrossSectionDataSet.cc.

84{
85 G4int Z = elm->GetZasInt();
86
87 if (IsElementApplicable(part, Z, mat)) {
88 return GetElementCrossSection(part, Z, mat);
89 }
90
91 // isotope-wise cross section making sum over available
92 // isotope cross sections, which may be incomplete, so
93 // the result is corrected
94 size_t nIso = elm->GetNumberOfIsotopes();
95 G4double fact = 0.0;
96 G4double xsec = 0.0;
97
98 // user-defined isotope abundances
99 const G4IsotopeVector* isoVector = elm->GetIsotopeVector();
100 const G4double* abundVector = elm->GetRelativeAbundanceVector();
101
102 for (size_t j=0; j<nIso; ++j) {
103 const G4Isotope* iso = (*isoVector)[j];
104 G4int A = iso->GetN();
105 if(abundVector[j] > 0.0 && IsIsoApplicable(part, Z, A, elm, mat)) {
106 fact += abundVector[j];
107 xsec += abundVector[j]*GetIsoCrossSection(part, Z, A, iso, elm, mat);
108 }
109 }
110 return (fact > 0.0) ? xsec/fact : 0.0;
111}
std::vector< G4Isotope * > G4IsotopeVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
G4double * GetRelativeAbundanceVector() const
Definition: G4Element.hh:167
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
G4int GetZasInt() const
Definition: G4Element.hh:132
G4IsotopeVector * GetIsotopeVector() const
Definition: G4Element.hh:163
G4int GetN() const
Definition: G4Isotope.hh:93
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 G4bool IsIsoApplicable(const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)

References A, G4VCrossSectionDataSet::GetElementCrossSection(), G4VCrossSectionDataSet::GetIsoCrossSection(), G4Element::GetIsotopeVector(), G4Isotope::GetN(), G4Element::GetNumberOfIsotopes(), G4Element::GetRelativeAbundanceVector(), G4Element::GetZasInt(), G4VCrossSectionDataSet::IsElementApplicable(), G4VCrossSectionDataSet::IsIsoApplicable(), and Z.

Referenced by G4VCrossSectionDataSet::GetCrossSection().

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 370 of file G4ParticleHPInelasticData.cc.

371{
372 outFile << "Extension of High Precision cross section for inelastic reaction of proton, deuteron, triton, He3 and alpha below 20MeV\n";
373}

◆ DumpPhysicsTable()

void G4ParticleHPInelasticData::DumpPhysicsTable ( const G4ParticleDefinition projectile)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 208 of file G4ParticleHPInelasticData.cc.

209{
210 if(&projectile!=theProjectile)
211 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use ParticleHP data for a wrong projectile!!!");
212
213 #ifdef G4VERBOSE
214 if ( G4HadronicParameters::Instance()->GetVerboseLevel() == 0 ) return;
215//
216// Dump element based cross section
217// range 10e-5 eV to 20 MeV
218// 10 point per decade
219// in barn
220//
221
222 G4cout << G4endl;
223 G4cout << G4endl;
224 G4cout << "Inelastic Cross Section of Neutron HP"<< G4endl;
225 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
226 G4cout << G4endl;
227 G4cout << "Name of Element" << G4endl;
228 G4cout << "Energy[eV] XS[barn]" << G4endl;
229 G4cout << G4endl;
230
231 size_t numberOfElements = G4Element::GetNumberOfElements();
232 static G4ThreadLocal G4ElementTable *theElementTable = 0 ;
233 if (!theElementTable) theElementTable= G4Element::GetElementTable();
234
235 for ( size_t i = 0 ; i < numberOfElements ; ++i )
236 {
237
238 G4cout << (*theElementTable)[i]->GetName() << G4endl;
239
240 G4int ie = 0;
241
242 for ( ie = 0 ; ie < 130 ; ie++ )
243 {
244 G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA ( 10.0 , ie/10.0 ) *CLHEP::eV;
245 G4bool outOfRange = false;
246
247 if ( eKinetic < 20*CLHEP::MeV )
248 {
249 G4cout << eKinetic/CLHEP::eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/CLHEP::barn << G4endl;
250 }
251
252 }
253
254 G4cout << G4endl;
255 }
256
257 //G4cout << "G4ParticleHPInelasticData::DumpPhysicsTable still to be implemented"<<G4endl;
258 #endif
259}
bool G4bool
Definition: G4Types.hh:86
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:230
static constexpr double barn
Definition: SystemOfUnits.h:86
static constexpr double eV

References CLHEP::barn, CLHEP::eV, G4cout, G4endl, G4ThreadLocal, G4Element::GetElementTable(), G4Pow::GetInstance(), G4Element::GetNumberOfElements(), GetVerboseLevel(), G4HadronicParameters::Instance(), CLHEP::MeV, G4Pow::powA(), and theProjectile.

◆ ForAllAtomsAndEnergies()

bool G4VCrossSectionDataSet::ForAllAtomsAndEnergies ( ) const
inlineinherited

◆ GetCrossSection() [1/2]

G4double G4VCrossSectionDataSet::GetCrossSection ( const G4DynamicParticle dp,
const G4Element elm,
const G4Material mat = nullptr 
)
inlineinherited

Definition at line 187 of file G4VCrossSectionDataSet.hh.

190{
191 return ComputeCrossSection(dp, elm, mat);
192}
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)

References G4VCrossSectionDataSet::ComputeCrossSection().

◆ GetCrossSection() [2/2]

G4double G4ParticleHPInelasticData::GetCrossSection ( const G4DynamicParticle projectile,
const G4Element anE,
G4double  aT 
)

Definition at line 263 of file G4ParticleHPInelasticData.cc.

265{
266 G4double result = 0;
267 G4bool outOfRange;
268 G4int index = anE->GetIndex();
269
270 // prepare neutron
271 G4double eKinetic = projectile->GetKineticEnergy();
272
273 if ( G4ParticleHPManager::GetInstance()->GetNeglectDoppler() )
274 {
275 //NEGLECT_DOPPLER
276 G4double factor = 1.0;
277 if ( eKinetic < aT * CLHEP::k_Boltzmann )
278 {
279 // below 0.1 eV neutrons
280 // Have to do some, but now just igonre.
281 // Will take care after performance check.
282 // factor = factor * targetV;
283 }
284 return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
285 }
286
287 G4ReactionProduct theNeutron( projectile->GetDefinition() );
288 theNeutron.SetMomentum( projectile->GetMomentum() );
289 theNeutron.SetKineticEnergy( eKinetic );
290
291 // prepare thermal nucleus
292 G4Nucleus aNuc;
293 G4double eps = 0.0001;
294 G4double theA = anE->GetN();
295 G4double theZ = anE->GetZ();
296 G4double eleMass;
297 eleMass = G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theA+eps), static_cast<G4int>(theZ+eps) );
298
299 G4ReactionProduct boosted;
300 G4double aXsection;
301
302 // MC integration loop
303 G4int counter = 0;
304 G4int failCount = 0;
305 G4double buffer = 0; G4int size = G4int(std::max(10., aT/60*CLHEP::kelvin));
306 G4ThreeVector neutronVelocity = 1./theProjectile->GetPDGMass()*theNeutron.GetMomentum();
307 G4double neutronVMag = neutronVelocity.mag();
308
309 // G4cout << " G4ParticleHPInelasticData 2 " << size << G4endl;//GDEB
310#ifndef G4PHP_DOPPLER_LOOP_ONCE
311 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer) // Loop checking, 11.05.2015, T. Koi
312 {
313 if(counter) buffer = result/counter;
314 while (counter<size) // Loop checking, 11.05.2015, T. Koi
315 {
316 counter ++;
317#endif
318 //G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus( eleMass/theProjectile->GetPDGMass(), aT );
319 //G4Nucleus::GetThermalNucleus requests normalization of mass in neutron mass
320 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus( eleMass/G4Neutron::Neutron()->GetPDGMass(), aT );
321 boosted.Lorentz(theNeutron, aThermalNuc);
322 G4double theEkin = boosted.GetKineticEnergy();
323 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
324 // G4cout << " G4ParticleHPInelasticData aXsection " << aXsection << " index " << index << " theEkin " << theEkin << " outOfRange " << outOfRange <<G4endl;//GDEB
325 if(aXsection <0)
326 {
327 if(failCount<1000)
328 {
329 failCount++;
330#ifndef G4PHP_DOPPLER_LOOP_ONCE
331 counter--;
332 continue;
333#endif
334 }
335 else
336 {
337 aXsection = 0;
338 }
339 }
340 // velocity correction.
341 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
342 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
343 result += aXsection;
344 }
345#ifndef G4PHP_DOPPLER_LOOP_ONCE
346 size += size;
347 }
348 result /= counter;
349#endif
350
351/*
352 // Checking impact of G4PHP_NEGLECT_DOPPLER
353 G4cout << " result " << result << " "
354 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
355 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
356*/
357// G4cout << this << " G4ParticleHPInelasticData result " << result << G4endl; //GDEB
358
359 return result;
360}
static const G4double eps
double mag() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double GetZ() const
Definition: G4Element.hh:131
size_t GetIndex() const
Definition: G4Element.hh:182
G4double GetN() const
Definition: G4Element.hh:135
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition: G4Nucleus.cc:236
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
G4double GetMass() const
static constexpr double k_Boltzmann
static constexpr double kelvin
T max(const T t1, const T t2)
brief Return the largest of the two arguments
#define buffer
Definition: xmlparse.cc:628

References buffer, eps, G4DynamicParticle::GetDefinition(), G4Element::GetIndex(), G4ParticleHPManager::GetInstance(), G4DynamicParticle::GetKineticEnergy(), G4ReactionProduct::GetKineticEnergy(), G4ReactionProduct::GetMass(), G4DynamicParticle::GetMomentum(), G4ReactionProduct::GetMomentum(), G4Element::GetN(), G4NucleiProperties::GetNuclearMass(), G4ParticleDefinition::GetPDGMass(), G4Nucleus::GetThermalNucleus(), G4Element::GetZ(), CLHEP::k_Boltzmann, CLHEP::kelvin, G4ReactionProduct::Lorentz(), CLHEP::Hep3Vector::mag(), G4INCL::Math::max(), G4Neutron::Neutron(), G4ReactionProduct::SetKineticEnergy(), G4ReactionProduct::SetMomentum(), theCrossSections, and theProjectile.

Referenced by GetIsoCrossSection().

◆ GetElementCrossSection()

G4double G4VCrossSectionDataSet::GetElementCrossSection ( const G4DynamicParticle dynPart,
G4int  Z,
const G4Material mat = nullptr 
)
virtualinherited

Reimplemented in G4EMDissociationCrossSection, G4IonsShenCrossSection, G4NeutrinoElectronCcXsc, G4NeutrinoElectronNcXsc, G4NeutrinoElectronTotXsc, G4NeutronElectronElXsc, G4PhotoNuclearCrossSection, G4NeutronCaptureXS, G4NeutronElasticXS, G4NeutronInelasticXS, G4ElectroNuclearCrossSection, G4BGGNucleonElasticXS, G4BGGPionElasticXS, G4BGGPionInelasticXS, G4BGGNucleonInelasticXS, G4CrossSectionElastic, G4CrossSectionInelastic, G4GammaNuclearXS, G4ParticleInelasticXS, G4ZeroXS, G4NucleonNuclearCrossSection, G4MuNeutrinoNucleusTotXsc, and G4KokoulinMuonNuclearXS.

Definition at line 114 of file G4VCrossSectionDataSet.cc.

117{
119 ed << "GetElementCrossSection is not implemented in <" << name << ">\n"
120 << "Particle: " << dynPart->GetDefinition()->GetParticleName()
121 << " Ekin(MeV)= " << dynPart->GetKineticEnergy()/MeV;
122 if(mat) { ed << " material: " << mat->GetName(); }
123 ed << " target Z= " << Z << G4endl;
124 G4Exception("G4VCrossSectionDataSet::GetElementCrossSection", "had001",
125 FatalException, ed);
126 return 0.0;
127}
@ 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
static constexpr double MeV
Definition: G4SIunits.hh:200
const G4String & GetName() const
Definition: G4Material.hh:173

References FatalException, G4endl, G4Exception(), G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4Material::GetName(), G4ParticleDefinition::GetParticleName(), MeV, G4VCrossSectionDataSet::name, and Z.

Referenced by G4QMDReaction::ApplyYourself(), G4VCrossSectionDataSet::ComputeCrossSection(), G4GammaNuclearXS::GetElementCrossSection(), G4GammaNuclearXS::GetIsoCrossSection(), and G4GammaNuclearXS::Initialise().

◆ GetIsoCrossSection()

G4double G4ParticleHPInelasticData::GetIsoCrossSection ( const G4DynamicParticle dp,
G4int  ,
G4int  ,
const G4Isotope ,
const G4Element element,
const G4Material material 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 142 of file G4ParticleHPInelasticData.cc.

147{
148 if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
149
151 element_cache = element;
153 G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
154 xs_cache = xs;
155 return xs;
156}
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)
string material
Definition: eplot.py:19

References element_cache, GetCrossSection(), G4DynamicParticle::GetKineticEnergy(), ke_cache, eplot::material, material_cache, and xs_cache.

◆ GetMaxKinEnergy()

G4double G4VCrossSectionDataSet::GetMaxKinEnergy ( ) const
inlineinherited

◆ GetMinKinEnergy()

G4double G4VCrossSectionDataSet::GetMinKinEnergy ( ) const
inlineinherited

◆ GetName()

const G4String & G4VCrossSectionDataSet::GetName ( ) const
inlineinherited

◆ GetProjectile()

G4ParticleDefinition * G4ParticleHPInelasticData::GetProjectile ( )
inline

Definition at line 85 of file G4ParticleHPInelasticData.hh.

85{return theProjectile;}

References theProjectile.

◆ GetVerboseLevel()

G4int G4ParticleHPInelasticData::GetVerboseLevel ( ) const
virtual

◆ IsElementApplicable()

G4bool G4VCrossSectionDataSet::IsElementApplicable ( const G4DynamicParticle ,
G4int  Z,
const G4Material mat = nullptr 
)
virtualinherited

◆ IsIsoApplicable()

G4bool G4ParticleHPInelasticData::IsIsoApplicable ( const G4DynamicParticle dp,
G4int  ,
G4int  ,
const G4Element ,
const G4Material  
)
virtual

◆ SelectIsotope()

const G4Isotope * G4VCrossSectionDataSet::SelectIsotope ( const G4Element anElement,
G4double  kinEnergy,
G4double  logE 
)
virtualinherited

Reimplemented in G4GammaNuclearXS, G4NeutronCaptureXS, G4NeutronElasticXS, G4NeutronInelasticXS, and G4ParticleInelasticXS.

Definition at line 149 of file G4VCrossSectionDataSet.cc.

151{
152 size_t nIso = anElement->GetNumberOfIsotopes();
153 const G4Isotope* iso = anElement->GetIsotope(0);
154
155 // more than 1 isotope
156 if(1 < nIso) {
157 const G4double* abundVector = anElement->GetRelativeAbundanceVector();
158 G4double sum = 0.0;
160 for (size_t j=0; j<nIso; ++j) {
161 sum += abundVector[j];
162 if(q <= sum) {
163 iso = anElement->GetIsotope(j);
164 break;
165 }
166 }
167 }
168 return iso;
169}
#define G4UniformRand()
Definition: Randomize.hh:52
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170

References G4UniformRand, G4Element::GetIsotope(), G4Element::GetNumberOfIsotopes(), and G4Element::GetRelativeAbundanceVector().

◆ SetForAllAtomsAndEnergies()

void G4VCrossSectionDataSet::SetForAllAtomsAndEnergies ( G4bool  val)
inlineinherited

◆ SetMaxKinEnergy()

void G4VCrossSectionDataSet::SetMaxKinEnergy ( G4double  value)
inlineinherited

◆ SetMinKinEnergy()

void G4VCrossSectionDataSet::SetMinKinEnergy ( G4double  value)
inlineinherited

◆ SetName()

void G4VCrossSectionDataSet::SetName ( const G4String nam)
inlineinherited

Definition at line 240 of file G4VCrossSectionDataSet.hh.

241{
242 name = nam;
243}

References G4VCrossSectionDataSet::name.

Referenced by G4ParticleHPInelasticData().

◆ SetVerboseLevel()

void G4ParticleHPInelasticData::SetVerboseLevel ( G4int  newValue)
virtual

Field Documentation

◆ element_cache

const G4Element* G4ParticleHPInelasticData::element_cache
private

Definition at line 104 of file G4ParticleHPInelasticData.hh.

Referenced by G4ParticleHPInelasticData(), and GetIsoCrossSection().

◆ instanceOfWorker

G4bool G4ParticleHPInelasticData::instanceOfWorker
private

◆ isForAllAtomsAndEnergies

G4bool G4VCrossSectionDataSet::isForAllAtomsAndEnergies
privateinherited

◆ ke_cache

G4double G4ParticleHPInelasticData::ke_cache
private

Definition at line 102 of file G4ParticleHPInelasticData.hh.

Referenced by G4ParticleHPInelasticData(), and GetIsoCrossSection().

◆ material_cache

const G4Material* G4ParticleHPInelasticData::material_cache
private

Definition at line 105 of file G4ParticleHPInelasticData.hh.

Referenced by G4ParticleHPInelasticData(), and GetIsoCrossSection().

◆ maxKinEnergy

G4double G4VCrossSectionDataSet::maxKinEnergy
privateinherited

◆ minKinEnergy

G4double G4VCrossSectionDataSet::minKinEnergy
privateinherited

◆ name

G4String G4VCrossSectionDataSet::name
protectedinherited

◆ registry

G4CrossSectionDataSetRegistry* G4VCrossSectionDataSet::registry
privateinherited

◆ theCrossSections

G4PhysicsTable* G4ParticleHPInelasticData::theCrossSections
private

◆ theHPData

G4ParticleHPData* G4ParticleHPInelasticData::theHPData
private

◆ theProjectile

G4ParticleDefinition* G4ParticleHPInelasticData::theProjectile
private

◆ verboseLevel

G4int G4VCrossSectionDataSet::verboseLevel
protectedinherited

Definition at line 168 of file G4VCrossSectionDataSet.hh.

Referenced by G4BGGNucleonElasticXS::BuildPhysicsTable(), G4BGGPionElasticXS::BuildPhysicsTable(), G4BGGPionInelasticXS::BuildPhysicsTable(), G4GammaNuclearXS::BuildPhysicsTable(), G4NeutronCaptureXS::BuildPhysicsTable(), G4NeutronElasticXS::BuildPhysicsTable(), G4NeutronInelasticXS::BuildPhysicsTable(), G4ParticleInelasticXS::BuildPhysicsTable(), G4BGGNucleonInelasticXS::BuildPhysicsTable(), G4LENDCrossSection::create_used_target_map(), G4BGGNucleonElasticXS::G4BGGNucleonElasticXS(), G4BGGNucleonInelasticXS::G4BGGNucleonInelasticXS(), G4BGGPionElasticXS::G4BGGPionElasticXS(), G4BGGPionInelasticXS::G4BGGPionInelasticXS(), G4GammaNuclearXS::G4GammaNuclearXS(), G4NeutronCaptureXS::G4NeutronCaptureXS(), G4NeutronElasticXS::G4NeutronElasticXS(), G4NeutronInelasticXS::G4NeutronInelasticXS(), G4ParticleInelasticXS::G4ParticleInelasticXS(), G4NeutronCaptureXS::GetElementCrossSection(), G4NeutronElasticXS::GetElementCrossSection(), G4NeutronInelasticXS::GetElementCrossSection(), G4BGGNucleonElasticXS::GetElementCrossSection(), G4BGGPionElasticXS::GetElementCrossSection(), G4BGGPionInelasticXS::GetElementCrossSection(), G4BGGNucleonInelasticXS::GetElementCrossSection(), G4GammaNuclearXS::GetElementCrossSection(), G4ParticleInelasticXS::GetElementCrossSection(), G4BGGNucleonElasticXS::GetIsoCrossSection(), G4BGGPionElasticXS::GetIsoCrossSection(), G4BGGPionInelasticXS::GetIsoCrossSection(), G4GammaNuclearXS::GetIsoCrossSection(), G4BGGNucleonInelasticXS::GetIsoCrossSection(), G4VCrossSectionDataSet::GetVerboseLevel(), G4NeutronElasticXS::Initialise(), G4ParticleInelasticXS::IsoCrossSection(), G4NeutronCaptureXS::IsoCrossSection(), G4NeutronInelasticXS::IsoCrossSection(), G4GammaNuclearXS::RetrieveVector(), G4NeutronCaptureXS::RetrieveVector(), G4NeutronInelasticXS::RetrieveVector(), G4ParticleInelasticXS::RetrieveVector(), and G4VCrossSectionDataSet::SetVerboseLevel().

◆ xs_cache

G4double G4ParticleHPInelasticData::xs_cache
private

Definition at line 103 of file G4ParticleHPInelasticData.hh.

Referenced by G4ParticleHPInelasticData(), and GetIsoCrossSection().


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