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

#include <G4NeutronHPInelasticData.hh>

+ Inheritance diagram for G4NeutronHPInelasticData:

Public Member Functions

 G4NeutronHPInelasticData ()
 
 ~G4NeutronHPInelasticData ()
 
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 &)
 
void IgnoreOnFlightDopplerBroadening ()
 
void EnableOnFlightDopplerBroadening ()
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
 
virtual G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
const G4StringGetName () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4VCrossSectionDataSet
void SetName (const G4String &)
 
- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 

Detailed Description

Definition at line 50 of file G4NeutronHPInelasticData.hh.

Constructor & Destructor Documentation

◆ G4NeutronHPInelasticData()

G4NeutronHPInelasticData::G4NeutronHPInelasticData ( )

Definition at line 44 of file G4NeutronHPInelasticData.cc.

45:G4VCrossSectionDataSet("NeutronHPInelasticXS")
46{
47
48 SetMinKinEnergy( 0*MeV );
49 SetMaxKinEnergy( 20*MeV );
50
51 ke_cache = 0.0;
52 xs_cache = 0.0;
53 element_cache = NULL;
54 material_cache = NULL;
55
56 onFlightDB = true;
57 theCrossSections = 0;
59}
void BuildPhysicsTable(const G4ParticleDefinition &)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)

◆ ~G4NeutronHPInelasticData()

G4NeutronHPInelasticData::~G4NeutronHPInelasticData ( )

Definition at line 61 of file G4NeutronHPInelasticData.cc.

62{
63 if ( theCrossSections != 0 ) theCrossSections->clearAndDestroy();
64 delete theCrossSections;
65}
void clearAndDestroy()

Member Function Documentation

◆ BuildPhysicsTable()

void G4NeutronHPInelasticData::BuildPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 107 of file G4NeutronHPInelasticData.cc.

108{
109 if(&aP!=G4Neutron::Neutron())
110 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
111
112//080428
113 if ( getenv( "G4NEUTRONHP_NEGLECT_DOPPLER" ) )
114 {
115 G4cout << "Find environment variable of \"G4NEUTRONHP_NEGLECT_DOPPLER\"." << G4endl;
116 G4cout << "On the fly Doppler broadening will be neglect in the cross section calculation of inelastic scattering of neutrons (<20MeV)." << G4endl;
117 onFlightDB = false;
118 }
119
120 size_t numberOfElements = G4Element::GetNumberOfElements();
121// theCrossSections = new G4PhysicsTable( numberOfElements );
122// TKDB
123 //if ( theCrossSections == 0 )
124 //{ theCrossSections = new G4PhysicsTable( numberOfElements ); }
125 if ( theCrossSections == NULL )
126 theCrossSections = new G4PhysicsTable( numberOfElements );
127 else
128 theCrossSections->clearAndDestroy();
129
130 // make a PhysicsVector for each element
131
132 static const G4ElementTable *theElementTable = G4Element::GetElementTable();
133 for( size_t i=0; i<numberOfElements; ++i )
134 {
136 Instance()->MakePhysicsVector((*theElementTable)[i], this);
137 theCrossSections->push_back(physVec);
138 }
139}
std::vector< G4Element * > G4ElementTable
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:406
static const G4ElementTable * GetElementTable()
Definition: G4Element.cc:399
G4PhysicsVector * MakePhysicsVector(G4Element *thE, G4NeutronHPFissionData *theP)
static G4NeutronHPData * Instance()
void push_back(G4PhysicsVector *)

Referenced by G4NeutronHPInelasticData().

◆ DumpPhysicsTable()

void G4NeutronHPInelasticData::DumpPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 141 of file G4NeutronHPInelasticData.cc.

142{
143 if(&aP!=G4Neutron::Neutron())
144 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
145
146//
147// Dump element based cross section
148// range 10e-5 eV to 20 MeV
149// 10 point per decade
150// in barn
151//
152
153 G4cout << G4endl;
154 G4cout << G4endl;
155 G4cout << "Inelastic Cross Section of Neutron HP"<< G4endl;
156 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
157 G4cout << G4endl;
158 G4cout << "Name of Element" << G4endl;
159 G4cout << "Energy[eV] XS[barn]" << G4endl;
160 G4cout << G4endl;
161
162 size_t numberOfElements = G4Element::GetNumberOfElements();
163 static const G4ElementTable *theElementTable = G4Element::GetElementTable();
164
165 for ( size_t i = 0 ; i < numberOfElements ; ++i )
166 {
167
168 G4cout << (*theElementTable)[i]->GetName() << G4endl;
169
170 G4int ie = 0;
171
172 for ( ie = 0 ; ie < 130 ; ie++ )
173 {
174 G4double eKinetic = 1.0e-5 * std::pow ( 10.0 , ie/10.0 ) *eV;
175 G4bool outOfRange = false;
176
177 if ( eKinetic < 20*MeV )
178 {
179 G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
180 }
181
182 }
183
184 G4cout << G4endl;
185 }
186
187 //G4cout << "G4NeutronHPInelasticData::DumpPhysicsTable still to be implemented"<<G4endl;
188}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67

◆ EnableOnFlightDopplerBroadening()

void G4NeutronHPInelasticData::EnableOnFlightDopplerBroadening ( )
inline

Definition at line 82 of file G4NeutronHPInelasticData.hh.

82{ onFlightDB = true; };

◆ GetCrossSection()

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

Definition at line 192 of file G4NeutronHPInelasticData.cc.

194{
195 G4double result = 0;
196 G4bool outOfRange;
197 G4int index = anE->GetIndex();
198
199 // prepare neutron
200 G4double eKinetic = aP->GetKineticEnergy();
201
202 // T. K.
203//if ( getenv( "G4NEUTRONHP_NEGLECT_DOPPLER" ) )
204//080428
205 if ( !onFlightDB )
206 {
207 G4double factor = 1.0;
208 if ( eKinetic < aT * k_Boltzmann )
209 {
210 // below 0.1 eV neutrons
211 // Have to do some, but now just igonre.
212 // Will take care after performance check.
213 // factor = factor * targetV;
214 }
215 return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
216 }
217
218 G4ReactionProduct theNeutron( aP->GetDefinition() );
219 theNeutron.SetMomentum( aP->GetMomentum() );
220 theNeutron.SetKineticEnergy( eKinetic );
221
222 // prepare thermal nucleus
223 G4Nucleus aNuc;
224 G4double eps = 0.0001;
225 G4double theA = anE->GetN();
226 G4double theZ = anE->GetZ();
227 G4double eleMass;
228 eleMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theA+eps), static_cast<G4int>(theZ+eps))
230
231 G4ReactionProduct boosted;
232 G4double aXsection;
233
234 // MC integration loop
235 G4int counter = 0;
236 G4int failCount = 0;
237 G4double buffer = 0;
238 G4int size = G4int(std::max(10., aT/60*kelvin));
239 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
240 G4double neutronVMag = neutronVelocity.mag();
241
242 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.01*buffer)
243 {
244 if(counter) buffer = result/counter;
245 while (counter<size)
246 {
247 counter ++;
248 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
249 boosted.Lorentz(theNeutron, aThermalNuc);
250 G4double theEkin = boosted.GetKineticEnergy();
251 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
252 if(aXsection <0)
253 {
254 if(failCount<1000)
255 {
256 failCount++;
257 counter--;
258 continue;
259 }
260 else
261 {
262 aXsection = 0;
263 }
264 }
265 // velocity correction.
266 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
267 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
268 result += aXsection;
269 }
270 size += size;
271 }
272 result /= counter;
273/*
274 // Checking impact of G4NEUTRONHP_NEGLECT_DOPPLER
275 G4cout << " result " << result << " "
276 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
277 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
278*/
279 return result;
280}
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:134
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition: G4Nucleus.cc:130
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
G4double GetMass() const
#define buffer
Definition: xmlparse.cc:611

Referenced by GetIsoCrossSection().

◆ GetIsoCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 80 of file G4NeutronHPInelasticData.cc.

85{
86 if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
87
88 ke_cache = dp->GetKineticEnergy();
89 element_cache = element;
90 material_cache = material;
91 G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
92 xs_cache = xs;
93 return xs;
94 //return GetCrossSection( dp , element , material->GetTemperature() );
95}
G4double GetTemperature() const
Definition: G4Material.hh:181
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)

◆ IgnoreOnFlightDopplerBroadening()

void G4NeutronHPInelasticData::IgnoreOnFlightDopplerBroadening ( )
inline

Definition at line 81 of file G4NeutronHPInelasticData.hh.

81{ onFlightDB = false; };

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 67 of file G4NeutronHPInelasticData.cc.

71{
72 G4double eKin = dp->GetKineticEnergy();
73 if ( eKin > GetMaxKinEnergy()
74 || eKin < GetMinKinEnergy()
75 || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
76
77 return true;
78}

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