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

#include <G4ParticleHPFissionData.hh>

+ Inheritance diagram for G4ParticleHPFissionData:

Public Member Functions

 G4ParticleHPFissionData ()
 
 ~G4ParticleHPFissionData () 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
 
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
 
G4String name
 

Detailed Description

Definition at line 48 of file G4ParticleHPFissionData.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPFissionData()

G4ParticleHPFissionData::G4ParticleHPFissionData ( )

Definition at line 49 of file G4ParticleHPFissionData.cc.

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

◆ ~G4ParticleHPFissionData()

G4ParticleHPFissionData::~G4ParticleHPFissionData ( )
override

Definition at line 66 of file G4ParticleHPFissionData.cc.

67{
68 if (theCrossSections != nullptr && !instanceOfWorker) {
69 theCrossSections->clearAndDestroy();
70 delete theCrossSections;
71 theCrossSections = nullptr;
72 }
73}
void clearAndDestroy()

Member Function Documentation

◆ BuildPhysicsTable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 100 of file G4ParticleHPFissionData.cc.

101{
102 if (&aP != G4Neutron::Neutron())
103 throw G4HadronicException(__FILE__, __LINE__,
104 "Attempt to use NeutronHP data for particles other than neutrons!!!");
105
108 return;
109 }
110
111 std::size_t numberOfElements = G4Element::GetNumberOfElements();
112 if (theCrossSections == nullptr)
113 theCrossSections = new G4PhysicsTable(numberOfElements);
114 else
115 theCrossSections->clearAndDestroy();
116
117 // make a PhysicsVector for each element
118
119 static G4ThreadLocal G4ElementTable* theElementTable = nullptr;
120 if (theElementTable == nullptr) theElementTable = G4Element::GetElementTable();
121 for (std::size_t i = 0; i < numberOfElements; ++i) {
123 ->MakePhysicsVector((*theElementTable)[i], this);
124 theCrossSections->push_back(physVec);
125 }
126
128}
std::vector< G4Element * > G4ElementTable
static G4ElementTable * GetElementTable()
Definition G4Element.cc:389
static size_t GetNumberOfElements()
Definition G4Element.cc:393
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
G4PhysicsVector * MakePhysicsVector(const G4Element *thE, G4ParticleHPFissionData *theP)
static G4ParticleHPData * Instance(G4ParticleDefinition *projectile)
G4PhysicsTable * GetFissionCrossSections() const
static G4ParticleHPManager * GetInstance()
void RegisterFissionCrossSections(G4PhysicsTable *val)
void push_back(G4PhysicsVector *)
#define G4ThreadLocal
Definition tls.hh:77

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 264 of file G4ParticleHPFissionData.cc.

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

◆ DumpPhysicsTable()

void G4ParticleHPFissionData::DumpPhysicsTable ( const G4ParticleDefinition & aP)
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 130 of file G4ParticleHPFissionData.cc.

131{
132 if (&aP != G4Neutron::Neutron())
133 throw G4HadronicException(__FILE__, __LINE__,
134 "Attempt to use NeutronHP data for particles other than neutrons!!!");
135
136#ifdef G4VERBOSE
138
139 //
140 // Dump element based cross section
141 // range 10e-5 eV to 20 MeV
142 // 10 point per decade
143 // in barn
144 //
145 G4cout << G4endl;
146 G4cout << G4endl;
147 G4cout << "Fission Cross Section of Neutron HP" << G4endl;
148 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
149 G4cout << G4endl;
150 G4cout << "Name of Element" << G4endl;
151 G4cout << "Energy[eV] XS[barn]" << G4endl;
152 G4cout << G4endl;
153
154 std::size_t numberOfElements = G4Element::GetNumberOfElements();
155 static G4ThreadLocal G4ElementTable* theElementTable = nullptr;
156 if (theElementTable == nullptr) theElementTable = G4Element::GetElementTable();
157
158 for (std::size_t i = 0; i < numberOfElements; ++i) {
159 G4cout << (*theElementTable)[i]->GetName() << G4endl;
160
161 if ((*((*theCrossSections)(i))).GetVectorLength() == 0) {
162 G4cout << "The cross-section data of the fission of this element is not available." << G4endl;
163 G4cout << G4endl;
164 continue;
165 }
166
167 for (G4int ie = 0; ie < 130; ++ie) {
168 G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA(10.0, ie / 10.0) * eV;
169 G4bool outOfRange = false;
170
171 if (eKinetic < 20 * MeV) {
172 G4cout << eKinetic / eV << " "
173 << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange) / barn << G4endl;
174 }
175 }
176
177 G4cout << G4endl;
178 }
179#endif
180}
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 G4ParticleHPFissionData::GetCrossSection ( const G4DynamicParticle * aP,
const G4Element * anE,
G4double aT )

Definition at line 182 of file G4ParticleHPFissionData.cc.

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 84 of file G4ParticleHPFissionData.cc.

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

◆ GetVerboseLevel()

G4int G4ParticleHPFissionData::GetVerboseLevel ( ) const

Definition at line 254 of file G4ParticleHPFissionData.cc.

Referenced by DumpPhysicsTable().

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 75 of file G4ParticleHPFissionData.cc.

78{
79 G4double eKin = dp->GetKineticEnergy();
80 return eKin <= GetMaxKinEnergy() && eKin >= GetMinKinEnergy()
82}

◆ SetVerboseLevel()

void G4ParticleHPFissionData::SetVerboseLevel ( G4int newValue)
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 259 of file G4ParticleHPFissionData.cc.


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