Geant4 10.7.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 ()
 
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 ()
 
void SetVerboseLevel (G4int)
 
G4int GetVerboseLevel () const
 
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 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 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 G4int GetVerboseLevel () 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
 

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 51 of file G4ParticleHPElasticData.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPElasticData()

G4ParticleHPElasticData::G4ParticleHPElasticData ( )

Definition at line 50 of file G4ParticleHPElasticData.cc.

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

◆ ~G4ParticleHPElasticData()

G4ParticleHPElasticData::~G4ParticleHPElasticData ( )

Definition at line 69 of file G4ParticleHPElasticData.cc.

70{
71 if ( theCrossSections != NULL && instanceOfWorker != true ) {
72 theCrossSections->clearAndDestroy();
73 delete theCrossSections;
74 theCrossSections = NULL;
75 }
76}
void clearAndDestroy()

Member Function Documentation

◆ BuildPhysicsTable()

void G4ParticleHPElasticData::BuildPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 118 of file G4ParticleHPElasticData.cc.

119{
120
121 if(&aP!=G4Neutron::Neutron())
122 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
123
124//080428
125 if ( G4ParticleHPManager::GetInstance()->GetNeglectDoppler() )
126 {
127 onFlightDB = false;
128 #ifdef G4VERBOSE
130 G4cout << "Find a flag of \"G4NEUTRONHP_NEGLECT_DOPPLER\"." << G4endl;
131 G4cout << "On the fly Doppler broadening will be neglect in the cross section calculation of elastic scattering of neutrons (<20MeV)." << G4endl;
132 }
133 #endif
134 }
135
138 return;
139 }
140
141 size_t numberOfElements = G4Element::GetNumberOfElements();
142// TKDB
143 //if ( theCrossSections == 0 ) theCrossSections = new G4PhysicsTable( numberOfElements );
144 if ( theCrossSections == NULL )
145 theCrossSections = new G4PhysicsTable( numberOfElements );
146 else
147 theCrossSections->clearAndDestroy();
148
149 // make a PhysicsVector for each element
150
151 static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
152 for( size_t i=0; i<numberOfElements; ++i )
153 {
155 Instance(G4Neutron::Neutron())->MakePhysicsVector((*theElementTable)[i], this);
156 theCrossSections->push_back(physVec);
157 }
158
160}
std::vector< G4Element * > G4ElementTable
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:397
static size_t GetNumberOfElements()
Definition: G4Element.cc:404
static G4HadronicParameters * Instance()
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4ParticleHPData * Instance(G4ParticleDefinition *projectile)
G4PhysicsVector * MakePhysicsVector(G4Element *thE, G4ParticleHPFissionData *theP)
G4PhysicsTable * GetElasticCrossSections()
void RegisterElasticCrossSections(G4PhysicsTable *val)
static G4ParticleHPManager * GetInstance()
void push_back(G4PhysicsVector *)
#define G4ThreadLocal
Definition: tls.hh:77

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 307 of file G4ParticleHPElasticData.cc.

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

◆ DumpPhysicsTable()

void G4ParticleHPElasticData::DumpPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 162 of file G4ParticleHPElasticData.cc.

163{
164 if(&aP!=G4Neutron::Neutron())
165 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
166
167 #ifdef G4VERBOSE
168 if ( G4HadronicParameters::Instance()->GetVerboseLevel() == 0 ) return;
169
170//
171// Dump element based cross section
172// range 10e-5 eV to 20 MeV
173// 10 point per decade
174// in barn
175//
176
177 G4cout << G4endl;
178 G4cout << G4endl;
179 G4cout << "Elastic Cross Section of Neutron HP"<< G4endl;
180 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
181 G4cout << G4endl;
182 G4cout << "Name of Element" << G4endl;
183 G4cout << "Energy[eV] XS[barn]" << G4endl;
184 G4cout << G4endl;
185
186 size_t numberOfElements = G4Element::GetNumberOfElements();
187 static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
188
189 for ( size_t i = 0 ; i < numberOfElements ; ++i )
190 {
191
192 G4cout << (*theElementTable)[i]->GetName() << G4endl;
193
194 G4int ie = 0;
195
196 for ( ie = 0 ; ie < 130 ; ie++ )
197 {
198 G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA ( 10.0 , ie/10.0 ) *eV;
199 G4bool outOfRange = false;
200
201 if ( eKinetic < 20*MeV )
202 {
203 G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
204 }
205
206 }
207
208 G4cout << G4endl;
209 }
210
211 //G4cout << "G4ParticleHPElasticData::DumpPhysicsTable still to be implemented"<<G4endl;
212 #endif
213}
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:230

◆ EnableOnFlightDopplerBroadening()

void G4ParticleHPElasticData::EnableOnFlightDopplerBroadening ( )
inline

Definition at line 84 of file G4ParticleHPElasticData.hh.

84{ onFlightDB = true; };

◆ GetCrossSection()

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

Definition at line 220 of file G4ParticleHPElasticData.cc.

222{
223 G4double result = 0;
224 G4bool outOfRange;
225 G4int index = anE->GetIndex();
226
227 // prepare neutron
228 G4double eKinetic = aP->GetKineticEnergy();
229
230 if ( !onFlightDB )
231 {
232 //NEGLECT_DOPPLER_B.
233 G4double factor = 1.0;
234 if ( eKinetic < aT * k_Boltzmann )
235 {
236 // below 0.1 eV neutrons
237 // Have to do some, but now just igonre.
238 // Will take care after performance check.
239 // factor = factor * targetV;
240 }
241 return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
242 }
243
244 G4ReactionProduct theNeutron( aP->GetDefinition() );
245 theNeutron.SetMomentum( aP->GetMomentum() );
246 theNeutron.SetKineticEnergy( eKinetic );
247
248 // prepare thermal nucleus
249 G4Nucleus aNuc;
250 G4double eps = 0.0001;
251 G4double theA = anE->GetN();
252 G4double theZ = anE->GetZ();
253 G4double eleMass;
254
255
256 eleMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) )
258
259 G4ReactionProduct boosted;
260 G4double aXsection;
261
262 // MC integration loop
263 G4int counter = 0;
264 G4double buffer = 0;
265 G4int size = G4int(std::max(10., aT/60*kelvin));
266 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
267 G4double neutronVMag = neutronVelocity.mag();
268
269 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.03*buffer) // Loop checking, 11.05.2015, T. Koi
270 {
271 if(counter) buffer = result/counter;
272 while (counter<size) // Loop checking, 11.05.2015, T. Koi
273 {
274 counter ++;
275 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
276 boosted.Lorentz(theNeutron, aThermalNuc);
277 G4double theEkin = boosted.GetKineticEnergy();
278 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
279 // velocity correction.
280 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
281 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
282 result += aXsection;
283 }
284 size += size;
285 }
286 result /= counter;
287/*
288 // Checking impact of G4NEUTRONHP_NEGLECT_DOPPLER
289 G4cout << " result " << result << " "
290 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
291 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
292*/
293 return result;
294}
double mag() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
G4double GetZ() const
Definition: G4Element.hh:130
size_t GetIndex() const
Definition: G4Element.hh:181
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:143
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
G4double GetMass() const
#define buffer
Definition: xmlparse.cc:628

Referenced by GetIsoCrossSection().

◆ GetIsoCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 92 of file G4ParticleHPElasticData.cc.

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

◆ GetVerboseLevel()

G4int G4ParticleHPElasticData::GetVerboseLevel ( ) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 296 of file G4ParticleHPElasticData.cc.

Referenced by BuildPhysicsTable(), and DumpPhysicsTable().

◆ IgnoreOnFlightDopplerBroadening()

void G4ParticleHPElasticData::IgnoreOnFlightDopplerBroadening ( )
inline

Definition at line 83 of file G4ParticleHPElasticData.hh.

83{ onFlightDB = false; };

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 78 of file G4ParticleHPElasticData.cc.

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

◆ SetVerboseLevel()

void G4ParticleHPElasticData::SetVerboseLevel ( G4int  newValue)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 302 of file G4ParticleHPElasticData.cc.

304{
306}

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