Geant4 11.1.1
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 ()
 
G4bool IsIsoApplicable (const G4DynamicParticle *, G4int, G4int, const G4Element *, const G4Material *)
 
G4double GetIsoCrossSection (const G4DynamicParticle *, G4int, G4int, const G4Isotope *, const G4Element *, const G4Material *)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, G4double aT)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
void DumpPhysicsTable (const G4ParticleDefinition &)
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int)
 
virtual void CrossSectionDescription (std::ostream &) const
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=nullptr, 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 GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, 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)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
virtual void SetVerboseLevel (G4int value)
 
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.

50:G4VCrossSectionDataSet("NeutronHPFissionXS")
51{
52 SetMinKinEnergy( 0*MeV );
53 SetMaxKinEnergy( 20*MeV );
54
55 theCrossSections = 0;
56 instanceOfWorker = false;
58 instanceOfWorker = true;
59 }
60 element_cache = nullptr;
61 material_cache = nullptr;
62 ke_cache = 0.0;
63 xs_cache = 0.0;
64 //BuildPhysicsTable(*G4Neutron::Neutron());
65}
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)
G4bool IsWorkerThread()
Definition: G4Threading.cc:123

◆ ~G4ParticleHPFissionData()

G4ParticleHPFissionData::~G4ParticleHPFissionData ( )

Definition at line 67 of file G4ParticleHPFissionData.cc.

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

Member Function Documentation

◆ BuildPhysicsTable()

void G4ParticleHPFissionData::BuildPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 105 of file G4ParticleHPFissionData.cc.

106{
107 if(&aP!=G4Neutron::Neutron())
108 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
109
112 return;
113 }
114
115 std::size_t numberOfElements = G4Element::GetNumberOfElements();
116 if ( theCrossSections == nullptr )
117 theCrossSections = new G4PhysicsTable( numberOfElements );
118 else
119 theCrossSections->clearAndDestroy();
120
121 // make a PhysicsVector for each element
122
123 static G4ThreadLocal G4ElementTable *theElementTable = nullptr ;
124 if (!theElementTable) theElementTable= G4Element::GetElementTable();
125 for( std::size_t i=0; i<numberOfElements; ++i )
126 {
128 Instance(G4Neutron::Neutron())->MakePhysicsVector((*theElementTable)[i], this);
129 theCrossSections->push_back(physVec);
130 }
131
133}
std::vector< G4Element * > G4ElementTable
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:403
static size_t GetNumberOfElements()
Definition: G4Element.cc:410
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4ParticleHPData * Instance(G4ParticleDefinition *projectile)
G4PhysicsVector * MakePhysicsVector(G4Element *thE, G4ParticleHPFissionData *theP)
G4PhysicsTable * GetFissionCrossSections()
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
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 269 of file G4ParticleHPFissionData.cc.

270{
271 outFile << "High Precision cross data based on Evaluated Nuclear Data Files (ENDF)\n"
272 << "for induced fission reaction of neutrons below 20MeV\n" ;
273}

◆ DumpPhysicsTable()

void G4ParticleHPFissionData::DumpPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 135 of file G4ParticleHPFissionData.cc.

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

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

Referenced by GetIsoCrossSection().

◆ GetIsoCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 89 of file G4ParticleHPFissionData.cc.

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

◆ GetVerboseLevel()

G4int G4ParticleHPFissionData::GetVerboseLevel ( ) const

Definition at line 259 of file G4ParticleHPFissionData.cc.

Referenced by DumpPhysicsTable().

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 76 of file G4ParticleHPFissionData.cc.

80{
81 G4double eKin = dp->GetKineticEnergy();
82 if ( eKin > GetMaxKinEnergy()
83 || eKin < GetMinKinEnergy()
84 || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
85
86 return true;
87}

◆ SetVerboseLevel()

void G4ParticleHPFissionData::SetVerboseLevel ( G4int  newValue)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 264 of file G4ParticleHPFissionData.cc.

265{
267}

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