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

#include <G4ParticleHPElasticData.hh>

+ Inheritance diagram for G4ParticleHPElasticData:

Public Member Functions

 G4ParticleHPElasticData ()
 
 ~G4ParticleHPElasticData () 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
 
void SetVerboseLevel (G4int) override
 
G4int GetVerboseLevel () const
 
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 51 of file G4ParticleHPElasticData.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPElasticData()

G4ParticleHPElasticData::G4ParticleHPElasticData ( )

Definition at line 53 of file G4ParticleHPElasticData.cc.

53 : G4VCrossSectionDataSet("NeutronHPElasticXS")
54{
55 SetMinKinEnergy(0 * MeV);
56 SetMaxKinEnergy(20 * MeV);
57
58 theCrossSections = nullptr;
59 instanceOfWorker = false;
61 instanceOfWorker = true;
62 }
63 element_cache = nullptr;
64 material_cache = nullptr;
65 ke_cache = 0.0;
66 xs_cache = 0.0;
67 // BuildPhysicsTable( *G4Neutron::Neutron() );
68}
G4VCrossSectionDataSet(const G4String &nam="")
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)
G4bool IsWorkerThread()

◆ ~G4ParticleHPElasticData()

G4ParticleHPElasticData::~G4ParticleHPElasticData ( )
override

Definition at line 70 of file G4ParticleHPElasticData.cc.

71{
72 if (theCrossSections != nullptr && !instanceOfWorker) {
73 theCrossSections->clearAndDestroy();
74 delete theCrossSections;
75 theCrossSections = nullptr;
76 }
77}

Member Function Documentation

◆ BuildPhysicsTable()

void G4ParticleHPElasticData::BuildPhysicsTable ( const G4ParticleDefinition & aP)
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 104 of file G4ParticleHPElasticData.cc.

105{
106 if (&aP != G4Neutron::Neutron())
107 throw G4HadronicException(__FILE__, __LINE__,
108 "Attempt to use NeutronHP data for particles other than neutrons!!!");
109
112 return;
113 }
114
115 std::size_t numberOfElements = G4Element::GetNumberOfElements();
116 // TKDB
117 // if ( theCrossSections == 0 ) theCrossSections = new G4PhysicsTable( numberOfElements );
118 if (theCrossSections == nullptr)
119 theCrossSections = new G4PhysicsTable(numberOfElements);
120 else
121 theCrossSections->clearAndDestroy();
122
123 // make a PhysicsVector for each element
124
125 auto theElementTable = G4Element::GetElementTable();
126 for (std::size_t i = 0; i < numberOfElements; ++i) {
127 G4PhysicsVector* physVec = G4ParticleHPData::Instance(G4Neutron::Neutron())
128 ->MakePhysicsVector((*theElementTable)[i], this);
129 theCrossSections->push_back(physVec);
130 }
131
133}
static std::size_t GetNumberOfElements()
Definition G4Element.cc:408
static const G4ElementTable * GetElementTable()
Definition G4Element.cc:401
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
G4PhysicsVector * MakePhysicsVector(const G4Element *thE, G4ParticleHPFissionData *theP)
static G4ParticleHPData * Instance(G4ParticleDefinition *projectile)
G4PhysicsTable * GetElasticCrossSections() const
void RegisterElasticCrossSections(G4PhysicsTable *val)
static G4ParticleHPManager * GetInstance()

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 263 of file G4ParticleHPElasticData.cc.

264{
265 outFile << "High Precision cross data based on Evaluated Nuclear Data Files (ENDF) for elastic "
266 "reaction of neutrons below 20MeV\n";
267}

◆ DumpPhysicsTable()

void G4ParticleHPElasticData::DumpPhysicsTable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 135 of file G4ParticleHPElasticData.cc.

136{
137#ifdef G4VERBOSE
139
140 //
141 // Dump element based cross section
142 // range 10e-5 eV to 20 MeV
143 // 10 point per decade
144 // in barn
145 //
146
147 G4cout << G4endl;
148 G4cout << G4endl;
149 G4cout << "Elastic Cross Section of Neutron HP" << G4endl;
150 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
151 G4cout << G4endl;
152 G4cout << "Name of Element" << G4endl;
153 G4cout << "Energy[eV] XS[barn]" << G4endl;
154 G4cout << G4endl;
155
156 std::size_t numberOfElements = G4Element::GetNumberOfElements();
157 auto theElementTable = G4Element::GetElementTable();
158
159 for (std::size_t i = 0; i < numberOfElements; ++i) {
160 G4cout << (*theElementTable)[i]->GetName() << G4endl;
161 G4int ie = 0;
162
163 for (ie = 0; ie < 130; ++ie) {
164 G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA(10.0, ie / 10.0) * eV;
165 G4bool outOfRange = false;
166
167 if (eKinetic < 20 * MeV) {
168 G4cout << eKinetic / eV << " "
169 << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange) / barn << G4endl;
170 }
171 }
172 G4cout << G4endl;
173 }
174#endif
175}
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4HadronicParameters * Instance()
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4double powA(G4double A, G4double y) const
Definition G4Pow.hh:230

◆ GetCrossSection()

G4double G4ParticleHPElasticData::GetCrossSection ( const G4DynamicParticle * aP,
const G4Element * anE,
G4double aT )

Definition at line 177 of file G4ParticleHPElasticData.cc.

179{
180 G4double result = 0;
181 G4bool outOfRange;
182 auto index = (G4int)anE->GetIndex();
183
184 // prepare neutron
185 G4double eKinetic = aP->GetKineticEnergy();
186
187 if (G4ParticleHPManager::GetInstance()->GetNeglectDoppler()) {
188 // NEGLECT_DOPPLER_B.
189 G4double factor = 1.0;
190 if (eKinetic < aT * k_Boltzmann) {
191 // below 0.1 eV neutrons
192 // Have to do some, but now just igonre.
193 // Will take care after performance check.
194 // factor = factor * targetV;
195 }
196 return ((*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange)) * factor;
197 }
198
199 G4ReactionProduct theNeutron(aP->GetDefinition());
200 theNeutron.SetMomentum(aP->GetMomentum());
201 theNeutron.SetKineticEnergy(eKinetic);
202
203 // prepare thermal nucleus
204 G4Nucleus aNuc;
205 G4double eps = 0.0001;
206 G4double theA = anE->GetN();
207 G4double theZ = anE->GetZ();
208 G4double eleMass;
209
210 eleMass = (G4NucleiProperties::GetNuclearMass(G4int(theA + eps), G4int(theZ + eps)))
212
213 G4ReactionProduct boosted;
214 G4double aXsection;
215
216 // MC integration loop
217 G4int counter = 0;
218 G4double buffer = 0;
219 G4int size = G4int(std::max(10., aT / 60 * kelvin));
220 G4ThreeVector neutronVelocity =
221 1. / G4Neutron::Neutron()->GetPDGMass() * theNeutron.GetMomentum();
222 G4double neutronVMag = neutronVelocity.mag();
223
224 while (counter == 0
225 || std::abs(buffer - result / std::max(1, counter))
226 > 0.03 * buffer) // Loop checking, 11.05.2015, T. Koi
227 {
228 if (counter != 0) buffer = result / counter;
229 while (counter < size) // Loop checking, 11.05.2015, T. Koi
230 {
231 counter++;
232 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
233 boosted.Lorentz(theNeutron, aThermalNuc);
234 G4double theEkin = boosted.GetKineticEnergy();
235 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
236 // velocity correction.
237 G4ThreeVector targetVelocity = 1. / aThermalNuc.GetMass() * aThermalNuc.GetMomentum();
238 aXsection *= (targetVelocity - neutronVelocity).mag() / neutronVMag;
239 result += aXsection;
240 }
241 size += size;
242 }
243 result /= counter;
244 /*
245 // Checking impact of G4NEUTRONHP_NEGLECT_DOPPLER
246 G4cout << " result " << result << " "
247 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
248 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
249 */
250 return result;
251}
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 G4ParticleHPElasticData::GetIsoCrossSection ( const G4DynamicParticle * dp,
G4int ,
G4int ,
const G4Isotope * ,
const G4Element * element,
const G4Material * material )
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 88 of file G4ParticleHPElasticData.cc.

92{
93 if (dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache)
94 return xs_cache;
95
96 ke_cache = dp->GetKineticEnergy();
97 element_cache = element;
98 material_cache = material;
99 G4double xs = GetCrossSection(dp, element, material->GetTemperature());
100 xs_cache = xs;
101 return xs;
102}
G4double GetTemperature() const
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)

◆ GetVerboseLevel()

G4int G4ParticleHPElasticData::GetVerboseLevel ( ) const

Definition at line 253 of file G4ParticleHPElasticData.cc.

Referenced by DumpPhysicsTable().

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 79 of file G4ParticleHPElasticData.cc.

82{
83 G4double eKin = dp->GetKineticEnergy();
84 return eKin <= GetMaxKinEnergy() && eKin >= GetMinKinEnergy()
86}

◆ SetVerboseLevel()

void G4ParticleHPElasticData::SetVerboseLevel ( G4int newValue)
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 258 of file G4ParticleHPElasticData.cc.


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