Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleHPInelasticData Class Reference

#include <G4ParticleHPInelasticData.hh>

+ Inheritance diagram for G4ParticleHPInelasticData:

Public Member Functions

 G4ParticleHPInelasticData (G4ParticleDefinition *projectile=G4Neutron::Neutron())
 
 ~G4ParticleHPInelasticData () override
 
G4bool IsIsoApplicable (const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *) override
 
G4double GetIsoCrossSection (const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *) override
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, G4double aT)
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
void DumpPhysicsTable (const G4ParticleDefinition &) override
 
G4ParticleDefinitionGetProjectile ()
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int) override
 
void CrossSectionDescription (std::ostream &) const override
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, 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 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)
 
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 53 of file G4ParticleHPInelasticData.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPInelasticData()

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

Definition at line 49 of file G4ParticleHPInelasticData.cc.

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

◆ ~G4ParticleHPInelasticData()

G4ParticleHPInelasticData::~G4ParticleHPInelasticData ( )
override

Definition at line 132 of file G4ParticleHPInelasticData.cc.

133{
134 if (theCrossSections != nullptr && !instanceOfWorker) {
135 theCrossSections->clearAndDestroy();
136 delete theCrossSections;
137 theCrossSections = nullptr;
138 }
139 if (theHPData != nullptr && !instanceOfWorker) {
140 delete theHPData;
141 theHPData = nullptr;
142 }
143}

Member Function Documentation

◆ BuildPhysicsTable()

void G4ParticleHPInelasticData::BuildPhysicsTable ( const G4ParticleDefinition & projectile)
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 170 of file G4ParticleHPInelasticData.cc.

171{
173 theCrossSections = G4ParticleHPManager::GetInstance()->GetInelasticCrossSections(&projectile);
174 return;
175 }
176 if (theHPData == nullptr)
177 theHPData = G4ParticleHPData::Instance(const_cast<G4ParticleDefinition*>(&projectile));
178
179 std::size_t numberOfElements = G4Element::GetNumberOfElements();
180 if (theCrossSections == nullptr)
181 theCrossSections = new G4PhysicsTable(numberOfElements);
182 else
183 theCrossSections->clearAndDestroy();
184
185 // make a PhysicsVector for each element
186
187 auto theElementTable = G4Element::GetElementTable();
188 for (std::size_t i = 0; i < numberOfElements; ++i) {
189 G4PhysicsVector* physVec = theHPData->MakePhysicsVector((*theElementTable)[i], this);
190 theCrossSections->push_back(physVec);
191 }
193}
static std::size_t GetNumberOfElements()
Definition G4Element.cc:408
static const G4ElementTable * GetElementTable()
Definition G4Element.cc:401
static G4ParticleHPData * Instance(G4ParticleDefinition *projectile)
G4PhysicsTable * GetInelasticCrossSections(const G4ParticleDefinition *part) const
void RegisterInelasticCrossSections(const G4ParticleDefinition *part, G4PhysicsTable *ptr)
static G4ParticleHPManager * GetInstance()
G4bool IsWorkerThread()

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 334 of file G4ParticleHPInelasticData.cc.

335{
336 outFile << "Extension of High Precision cross section for inelastic reaction of proton, "
337 "deuteron, triton, He3 and alpha below 20MeV\n";
338}

◆ DumpPhysicsTable()

void G4ParticleHPInelasticData::DumpPhysicsTable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 195 of file G4ParticleHPInelasticData.cc.

196{
197#ifdef G4VERBOSE
199
200 // Dump element based cross section
201 // range 10e-5 eV to 20 MeV
202 // 10 point per decade
203 // in barn
204
205 G4cout << G4endl;
206 G4cout << G4endl;
207 G4cout << "Inelastic Cross Section of Neutron HP" << G4endl;
208 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
209 G4cout << G4endl;
210 G4cout << "Name of Element" << G4endl;
211 G4cout << "Energy[eV] XS[barn]" << G4endl;
212 G4cout << G4endl;
213
214 std::size_t numberOfElements = G4Element::GetNumberOfElements();
215 auto theElementTable = G4Element::GetElementTable();
216
217 for (std::size_t i = 0; i < numberOfElements; ++i) {
218 G4cout << (*theElementTable)[i]->GetName() << G4endl;
219
220 G4int ie = 0;
221
222 for (ie = 0; ie < 130; ie++) {
223 G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA(10.0, ie / 10.0) * CLHEP::eV;
224 G4bool outOfRange = false;
225
226 if (eKinetic < 20 * CLHEP::MeV) {
227 G4cout << eKinetic / CLHEP::eV << " "
228 << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange) / CLHEP::barn
229 << G4endl;
230 }
231 }
232 G4cout << G4endl;
233 }
234#endif
235}
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4double powA(G4double A, G4double y) const
Definition G4Pow.hh:230

◆ GetCrossSection()

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

Definition at line 237 of file G4ParticleHPInelasticData.cc.

239{
240 G4double result = 0;
241 G4bool outOfRange;
242 std::size_t index = anE->GetIndex();
243
244 // prepare neutron
245 G4double eKinetic = projectile->GetKineticEnergy();
246
247 if (G4ParticleHPManager::GetInstance()->GetNeglectDoppler()) {
248 // NEGLECT_DOPPLER
249 G4double factor = 1.0;
250 if (eKinetic < aT * CLHEP::k_Boltzmann) {
251 // below 0.1 eV neutrons
252 // Have to do some, but now just igonre.
253 // Will take care after performance check.
254 // factor = factor * targetV;
255 }
256 return ((*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange)) * factor;
257 }
258
259 G4ReactionProduct theNeutron(projectile->GetDefinition());
260 theNeutron.SetMomentum(projectile->GetMomentum());
261 theNeutron.SetKineticEnergy(eKinetic);
262
263 // prepare thermal nucleus
264 G4Nucleus aNuc;
265 G4double eps = 0.0001;
266 G4double theA = anE->GetN();
267 G4double theZ = anE->GetZ();
268 G4double eleMass;
269 eleMass = G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theA + eps),
270 static_cast<G4int>(theZ + eps));
271
272 G4ReactionProduct boosted;
273 G4double aXsection;
274
275 // MC integration loop
276 G4int counter = 0;
277 G4int failCount = 0;
278 G4double buffer = 0;
279 G4int size = G4int(std::max(10., aT / 60 * CLHEP::kelvin));
280 G4ThreeVector neutronVelocity = 1. / theProjectile->GetPDGMass() * theNeutron.GetMomentum();
281 G4double neutronVMag = neutronVelocity.mag();
282
283#ifndef G4PHP_DOPPLER_LOOP_ONCE
284 while (counter == 0
285 || std::abs(buffer - result / std::max(1, counter))
286 > 0.01 * buffer) // Loop checking, 11.05.2015, T. Koi
287 {
288 if (counter != 0) buffer = result / counter;
289 while (counter < size) // Loop checking, 11.05.2015, T. Koi
290 {
291 ++counter;
292#endif
293 G4ReactionProduct aThermalNuc =
294 aNuc.GetThermalNucleus(eleMass / G4Neutron::Neutron()->GetPDGMass(), aT);
295 boosted.Lorentz(theNeutron, aThermalNuc);
296 G4double theEkin = boosted.GetKineticEnergy();
297 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
298 if (aXsection < 0) {
299 if (failCount < 1000) {
300 ++failCount;
301#ifndef G4PHP_DOPPLER_LOOP_ONCE
302 --counter;
303 continue;
304#endif
305 }
306 else {
307 aXsection = 0;
308 }
309 }
310 // velocity correction.
311 G4ThreeVector targetVelocity = 1. / aThermalNuc.GetMass() * aThermalNuc.GetMomentum();
312 aXsection *= (targetVelocity - neutronVelocity).mag() / neutronVMag;
313 result += aXsection;
314 }
315#ifndef G4PHP_DOPPLER_LOOP_ONCE
316 size += size;
317 }
318 result /= counter;
319#endif
320
321 return result;
322}
CLHEP::Hep3Vector G4ThreeVector
double mag() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double GetZ() const
Definition G4Element.hh:119
std::size_t GetIndex() const
Definition G4Element.hh:159
G4double GetN() const
Definition G4Element.hh:123
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition G4Nucleus.cc:248
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
G4double GetMass() const

Referenced by GetIsoCrossSection().

◆ GetIsoCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 154 of file G4ParticleHPInelasticData.cc.

158{
159 if (dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache)
160 return xs_cache;
161
162 ke_cache = dp->GetKineticEnergy();
163 element_cache = element;
164 material_cache = material;
165 G4double xs = GetCrossSection(dp, element, material->GetTemperature());
166 xs_cache = xs;
167 return xs;
168}
G4double GetTemperature() const
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)

◆ GetProjectile()

G4ParticleDefinition * G4ParticleHPInelasticData::GetProjectile ( )
inline

Definition at line 81 of file G4ParticleHPInelasticData.hh.

81{ return theProjectile; }

◆ GetVerboseLevel()

G4int G4ParticleHPInelasticData::GetVerboseLevel ( ) const

Definition at line 324 of file G4ParticleHPInelasticData.cc.

Referenced by DumpPhysicsTable(), and G4ParticleHPInelasticData().

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 145 of file G4ParticleHPInelasticData.cc.

148{
149 G4double eKin = dp->GetKineticEnergy();
150 return eKin <= GetMaxKinEnergy() && eKin >= GetMinKinEnergy()
151 && dp->GetDefinition() == theProjectile;
152}

◆ SetVerboseLevel()

void G4ParticleHPInelasticData::SetVerboseLevel ( G4int newValue)
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 329 of file G4ParticleHPInelasticData.cc.


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