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

Constructor & Destructor Documentation

◆ G4ParticleHPElasticData()

G4ParticleHPElasticData::G4ParticleHPElasticData ( )

Definition at line 55 of file G4ParticleHPElasticData.cc.

56:G4VCrossSectionDataSet("NeutronHPElasticXS")
57{
58 SetMinKinEnergy( 0*MeV );
59 SetMaxKinEnergy( 20*MeV );
60
61 theCrossSections = 0;
62 instanceOfWorker = false;
64 instanceOfWorker = true;
65 }
66 element_cache = nullptr;
67 material_cache = nullptr;
68 ke_cache = 0.0;
69 xs_cache = 0.0;
70// BuildPhysicsTable( *G4Neutron::Neutron() );
71}
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)
G4bool IsWorkerThread()
Definition: G4Threading.cc:123

◆ ~G4ParticleHPElasticData()

G4ParticleHPElasticData::~G4ParticleHPElasticData ( )

Definition at line 73 of file G4ParticleHPElasticData.cc.

74{
75 if ( theCrossSections != nullptr && instanceOfWorker != true ) {
76 theCrossSections->clearAndDestroy();
77 delete theCrossSections;
78 theCrossSections = nullptr;
79 }
80}
void clearAndDestroy()

Member Function Documentation

◆ BuildPhysicsTable()

void G4ParticleHPElasticData::BuildPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 122 of file G4ParticleHPElasticData.cc.

123{
124
125 if(&aP!=G4Neutron::Neutron())
126 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
127
130 return;
131 }
132
133 std::size_t numberOfElements = G4Element::GetNumberOfElements();
134// TKDB
135 //if ( theCrossSections == 0 ) theCrossSections = new G4PhysicsTable( numberOfElements );
136 if ( theCrossSections == nullptr )
137 theCrossSections = new G4PhysicsTable( numberOfElements );
138 else
139 theCrossSections->clearAndDestroy();
140
141 // make a PhysicsVector for each element
142
143 static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
144 for( std::size_t i=0; i<numberOfElements; ++i )
145 {
147 Instance(G4Neutron::Neutron())->MakePhysicsVector((*theElementTable)[i], this);
148 theCrossSections->push_back(physVec);
149 }
150
152}
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 * 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 290 of file G4ParticleHPElasticData.cc.

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

◆ DumpPhysicsTable()

void G4ParticleHPElasticData::DumpPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 154 of file G4ParticleHPElasticData.cc.

155{
156 if(&aP!=G4Neutron::Neutron())
157 throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
158
159#ifdef G4VERBOSE
160 if ( G4HadronicParameters::Instance()->GetVerboseLevel() == 0 ) return;
161
162//
163// Dump element based cross section
164// range 10e-5 eV to 20 MeV
165// 10 point per decade
166// in barn
167//
168
169 G4cout << G4endl;
170 G4cout << G4endl;
171 G4cout << "Elastic Cross Section of Neutron HP"<< G4endl;
172 G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
173 G4cout << G4endl;
174 G4cout << "Name of Element" << G4endl;
175 G4cout << "Energy[eV] XS[barn]" << G4endl;
176 G4cout << G4endl;
177
178 std::size_t numberOfElements = G4Element::GetNumberOfElements();
179 static G4ThreadLocal G4ElementTable *theElementTable = 0 ;
180 if (!theElementTable) theElementTable= G4Element::GetElementTable();
181
182 for ( std::size_t i = 0 ; i < numberOfElements ; ++i )
183 {
184 G4cout << (*theElementTable)[i]->GetName() << G4endl;
185 G4int ie = 0;
186
187 for ( ie = 0 ; ie < 130 ; ++ie )
188 {
189 G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA ( 10.0 , ie/10.0 ) *eV;
190 G4bool outOfRange = false;
191
192 if ( eKinetic < 20*MeV )
193 {
194 G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
195 }
196 }
197 G4cout << G4endl;
198 }
199#endif
200}
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 G4ParticleHPElasticData::GetCrossSection ( const G4DynamicParticle aP,
const G4Element anE,
G4double  aT 
)

Definition at line 202 of file G4ParticleHPElasticData.cc.

204{
205 G4double result = 0;
206 G4bool outOfRange;
207 G4int index = (G4int)anE->GetIndex();
208
209 // prepare neutron
210 G4double eKinetic = aP->GetKineticEnergy();
211
213 {
214 //NEGLECT_DOPPLER_B.
215 G4double factor = 1.0;
216 if ( eKinetic < aT * k_Boltzmann )
217 {
218 // below 0.1 eV neutrons
219 // Have to do some, but now just igonre.
220 // Will take care after performance check.
221 // factor = factor * targetV;
222 }
223 return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
224 }
225
226 G4ReactionProduct theNeutron( aP->GetDefinition() );
227 theNeutron.SetMomentum( aP->GetMomentum() );
228 theNeutron.SetKineticEnergy( eKinetic );
229
230 // prepare thermal nucleus
231 G4Nucleus aNuc;
232 G4double eps = 0.0001;
233 G4double theA = anE->GetN();
234 G4double theZ = anE->GetZ();
235 G4double eleMass;
236
237
238 eleMass = ( G4NucleiProperties::GetNuclearMass( G4int(theA+eps) , G4int(theZ+eps) )
240
241 G4ReactionProduct boosted;
242 G4double aXsection;
243
244 // MC integration loop
245 G4int counter = 0;
246 G4double buffer = 0;
247 G4int size = G4int(std::max(10., aT/60*kelvin));
248 G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
249 G4double neutronVMag = neutronVelocity.mag();
250
251 while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.03*buffer) // Loop checking, 11.05.2015, T. Koi
252 {
253 if(counter) buffer = result/counter;
254 while (counter<size) // Loop checking, 11.05.2015, T. Koi
255 {
256 counter ++;
257 G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
258 boosted.Lorentz(theNeutron, aThermalNuc);
259 G4double theEkin = boosted.GetKineticEnergy();
260 aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
261 // velocity correction.
262 G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
263 aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
264 result += aXsection;
265 }
266 size += size;
267 }
268 result /= counter;
269/*
270 // Checking impact of G4NEUTRONHP_NEGLECT_DOPPLER
271 G4cout << " result " << result << " "
272 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
273 << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
274*/
275 return result;
276}
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 G4ParticleHPElasticData::GetIsoCrossSection ( const G4DynamicParticle dp,
G4int  ,
G4int  ,
const G4Isotope ,
const G4Element element,
const G4Material material 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 96 of file G4ParticleHPElasticData.cc.

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

◆ GetVerboseLevel()

G4int G4ParticleHPElasticData::GetVerboseLevel ( ) const

Definition at line 278 of file G4ParticleHPElasticData.cc.

Referenced by DumpPhysicsTable().

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 82 of file G4ParticleHPElasticData.cc.

86{
87
88 G4double eKin = dp->GetKineticEnergy();
89 if ( eKin > GetMaxKinEnergy()
90 || eKin < GetMinKinEnergy()
91 || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
92
93 return true;
94}

◆ SetVerboseLevel()

void G4ParticleHPElasticData::SetVerboseLevel ( G4int  newValue)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 284 of file G4ParticleHPElasticData.cc.

286{
288}

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