Geant4 11.3.0
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 {0}
 
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}
G4VCrossSectionDataSet(const G4String &nam="")
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)
G4bool IsWorkerThread()

◆ ~G4ParticleHPFissionData()

G4ParticleHPFissionData::~G4ParticleHPFissionData ( )
override

Definition at line 65 of file G4ParticleHPFissionData.cc.

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

Member Function Documentation

◆ BuildPhysicsTable()

void G4ParticleHPFissionData::BuildPhysicsTable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 99 of file G4ParticleHPFissionData.cc.

100{
103 return;
104 }
105
106 std::size_t numberOfElements = G4Element::GetNumberOfElements();
107 if (theCrossSections == nullptr)
108 theCrossSections = new G4PhysicsTable(numberOfElements);
109 else
110 theCrossSections->clearAndDestroy();
111
112 // make a PhysicsVector for each element
113
114 auto theElementTable = G4Element::GetElementTable();
115 for (std::size_t i = 0; i < numberOfElements; ++i) {
116 G4PhysicsVector* physVec = G4ParticleHPData::Instance(G4Neutron::Neutron())
117 ->MakePhysicsVector((*theElementTable)[i], this);
118 theCrossSections->push_back(physVec);
119 }
120
122}
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 * GetFissionCrossSections() const
static G4ParticleHPManager * GetInstance()
void RegisterFissionCrossSections(G4PhysicsTable *val)

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 253 of file G4ParticleHPFissionData.cc.

254{
255 outFile << "High Precision cross data based on Evaluated Nuclear Data Files (ENDF)\n"
256 << "for induced fission reaction of neutrons below 20MeV\n";
257}

◆ DumpPhysicsTable()

void G4ParticleHPFissionData::DumpPhysicsTable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 124 of file G4ParticleHPFissionData.cc.

125{
126#ifdef G4VERBOSE
128
129 //
130 // Dump element based cross section
131 // range 10e-5 eV to 20 MeV
132 // 10 point per decade
133 // in barn
134 //
135 G4cout << G4endl;
136 G4cout << G4endl;
137 G4cout << "Fission Cross Section of Neutron HP" << G4endl;
138 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
139 G4cout << G4endl;
140 G4cout << "Name of Element" << G4endl;
141 G4cout << "Energy[eV] XS[barn]" << G4endl;
142 G4cout << G4endl;
143
144 std::size_t numberOfElements = G4Element::GetNumberOfElements();
145 auto theElementTable = G4Element::GetElementTable();
146
147 for (std::size_t i = 0; i < numberOfElements; ++i) {
148 G4cout << (*theElementTable)[i]->GetName() << G4endl;
149
150 if ((*((*theCrossSections)(i))).GetVectorLength() == 0) {
151 G4cout << "The cross-section data of the fission of this element is not available." << G4endl;
152 G4cout << G4endl;
153 continue;
154 }
155
156 for (G4int ie = 0; ie < 130; ++ie) {
157 G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA(10.0, ie / 10.0) * eV;
158 G4bool outOfRange = false;
159
160 if (eKinetic < 20 * MeV) {
161 G4cout << eKinetic / eV << " "
162 << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange) / barn << G4endl;
163 }
164 }
165
166 G4cout << G4endl;
167 }
168#endif
169}
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 171 of file G4ParticleHPFissionData.cc.

173{
174 G4double result = 0;
175 if (anE->GetZ() < 88) return result;
176 G4bool outOfRange;
177 auto index = (G4int)anE->GetIndex();
178
179 if (((*theCrossSections)(index))->GetVectorLength() == 0) return result;
180
181 // prepare neutron
182 G4double eKinetic = aP->GetKineticEnergy();
183 G4ReactionProduct theNeutronRP(aP->GetDefinition());
184 theNeutronRP.SetMomentum(aP->GetMomentum());
185 theNeutronRP.SetKineticEnergy(eKinetic);
186
187 if (G4ParticleHPManager::GetInstance()->GetNeglectDoppler()) {
188 // NEGLECT_DOPPLER
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 // prepare thermal nucleus
200 G4Nucleus aNuc;
201 G4double eps = 0.0001;
202 G4double theA = anE->GetN();
203 G4double theZ = anE->GetZ();
204 G4double eleMass;
205 eleMass = (G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theA + eps),
206 static_cast<G4int>(theZ + eps)))
208
209 G4ReactionProduct boosted;
210 G4double aXsection;
211
212 // MC integration loop
213 G4int counter = 0;
214 G4double buffer = 0;
215 G4int size = G4int(std::max(10., aT / 60 * kelvin));
216 G4ThreeVector neutronVelocity =
217 1. / G4Neutron::Neutron()->GetPDGMass() * theNeutronRP.GetMomentum();
218 G4double neutronVMag = neutronVelocity.mag();
219
220 while (counter == 0
221 || std::abs(buffer - result / std::max(1, counter))
222 > 0.01 * buffer) // Loop checking, 11.05.2015, T. Koi
223 {
224 if (counter != 0) buffer = result / counter;
225 while (counter < size) // Loop checking, 11.05.2015, T. Koi
226 {
227 counter++;
228 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
229 boosted.Lorentz(theNeutronRP, aThermalNuc);
230 G4double theEkin = boosted.GetKineticEnergy();
231 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
232 // velocity correction.
233 G4ThreeVector targetVelocity = 1. / aThermalNuc.GetMass() * aThermalNuc.GetMomentum();
234 aXsection *= (targetVelocity - neutronVelocity).mag() / neutronVMag;
235 result += aXsection;
236 }
237 size += size;
238 }
239 result /= counter;
240 return result;
241}
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 G4ParticleHPFissionData::GetIsoCrossSection ( const G4DynamicParticle * dp,
G4int ,
G4int ,
const G4Isotope * ,
const G4Element * element,
const G4Material * material )
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 83 of file G4ParticleHPFissionData.cc.

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

◆ GetVerboseLevel()

G4int G4ParticleHPFissionData::GetVerboseLevel ( ) const

Definition at line 243 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 74 of file G4ParticleHPFissionData.cc.

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

◆ SetVerboseLevel()

void G4ParticleHPFissionData::SetVerboseLevel ( G4int newValue)
overridevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 248 of file G4ParticleHPFissionData.cc.


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