Geant4  10.02.p03
G4ParticleHPElasticData Class Reference

#include <G4ParticleHPElasticData.hh>

Inheritance diagram for G4ParticleHPElasticData:
Collaboration 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=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 G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
const G4StringGetName () const
 

Private Attributes

G4PhysicsTabletheCrossSections
 
G4bool onFlightDB
 
G4bool instanceOfWorker
 
G4double ke_cache
 
G4double xs_cache
 
const G4Elementelement_cache
 
const G4Materialmaterial_cache
 

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 49 of file G4ParticleHPElasticData.cc.

50 :G4VCrossSectionDataSet("NeutronHPElasticXS")
51 {
52  SetMinKinEnergy( 0*MeV );
53  SetMaxKinEnergy( 20*MeV );
54 
55  theCrossSections = 0;
56  onFlightDB = true;
57  instanceOfWorker = false;
59  instanceOfWorker = true;
60  }
61 // BuildPhysicsTable( *G4Neutron::Neutron() );
62 }
static const double MeV
Definition: G4SIunits.hh:211
G4VCrossSectionDataSet(const G4String &nam="")
void SetMinKinEnergy(G4double value)
G4bool IsWorkerThread()
Definition: G4Threading.cc:135
void SetMaxKinEnergy(G4double value)
Here is the call graph for this function:

◆ ~G4ParticleHPElasticData()

G4ParticleHPElasticData::~G4ParticleHPElasticData ( )

Definition at line 64 of file G4ParticleHPElasticData.cc.

65 {
66  if ( theCrossSections != NULL && instanceOfWorker != true ) {
68  delete theCrossSections;
69  theCrossSections = NULL;
70  }
71 }
void clearAndDestroy()
Here is the call graph for this function:

Member Function Documentation

◆ BuildPhysicsTable()

void G4ParticleHPElasticData::BuildPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 113 of file G4ParticleHPElasticData.cc.

114 {
115 
116  if(&aP!=G4Neutron::Neutron())
117  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
118 
119 //080428
120  if ( G4ParticleHPManager::GetInstance()->GetNeglectDoppler() )
121  {
122  G4cout << "Find a flag of \"G4NEUTRONHP_NEGLECT_DOPPLER\"." << G4endl;
123  G4cout << "On the fly Doppler broadening will be neglect in the cross section calculation of elastic scattering of neutrons (<20MeV)." << G4endl;
124  onFlightDB = false;
125  }
126 
127  if ( G4Threading::IsWorkerThread() ) {
129  return;
130  }
131 
132  size_t numberOfElements = G4Element::GetNumberOfElements();
133 // TKDB
134  //if ( theCrossSections == 0 ) theCrossSections = new G4PhysicsTable( numberOfElements );
135  if ( theCrossSections == NULL )
136  theCrossSections = new G4PhysicsTable( numberOfElements );
137  else
139 
140  // make a PhysicsVector for each element
141 
142  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
143  for( size_t i=0; i<numberOfElements; ++i )
144  {
146  Instance(G4Neutron::Neutron())->MakePhysicsVector((*theElementTable)[i], this);
147  theCrossSections->push_back(physVec);
148  }
149 
151 }
static G4ParticleHPManager * GetInstance()
void push_back(G4PhysicsVector *)
#define G4ThreadLocal
Definition: tls.hh:89
void RegisterElasticCrossSections(G4PhysicsTable *val)
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
G4PhysicsTable * GetElasticCrossSections()
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4bool IsWorkerThread()
Definition: G4Threading.cc:135
G4PhysicsVector * MakePhysicsVector(G4Element *thE, G4ParticleHPFissionData *theP)
static G4ParticleHPData * Instance(G4ParticleDefinition *projectile)
#define G4endl
Definition: G4ios.hh:61
std::vector< G4Element * > G4ElementTable
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
void clearAndDestroy()
Here is the call graph for this function:

◆ CrossSectionDescription()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 295 of file G4ParticleHPElasticData.cc.

296 {
297  outFile << "High Precision cross data based on Evaluated Nuclear Data Files (ENDF) for elastic reaction of neutrons below 20MeV\n" ;
298 }
Here is the caller graph for this function:

◆ DumpPhysicsTable()

void G4ParticleHPElasticData::DumpPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 153 of file G4ParticleHPElasticData.cc.

154 {
155  if(&aP!=G4Neutron::Neutron())
156  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
157 
158 //
159 // Dump element based cross section
160 // range 10e-5 eV to 20 MeV
161 // 10 point per decade
162 // in barn
163 //
164 
165  G4cout << G4endl;
166  G4cout << G4endl;
167  G4cout << "Elastic Cross Section of Neutron HP"<< G4endl;
168  G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
169  G4cout << G4endl;
170  G4cout << "Name of Element" << G4endl;
171  G4cout << "Energy[eV] XS[barn]" << G4endl;
172  G4cout << G4endl;
173 
174  size_t numberOfElements = G4Element::GetNumberOfElements();
175  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
176 
177  for ( size_t i = 0 ; i < numberOfElements ; ++i )
178  {
179 
180  G4cout << (*theElementTable)[i]->GetName() << G4endl;
181 
182  G4int ie = 0;
183 
184  for ( ie = 0 ; ie < 130 ; ie++ )
185  {
186  G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA ( 10.0 , ie/10.0 ) *eV;
187  G4bool outOfRange = false;
188 
189  if ( eKinetic < 20*MeV )
190  {
191  G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
192  }
193 
194  }
195 
196  G4cout << G4endl;
197  }
198 
199 
200 // G4cout << "G4ParticleHPElasticData::DumpPhysicsTable still to be implemented"<<G4endl;
201 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
static const double MeV
Definition: G4SIunits.hh:211
#define G4ThreadLocal
Definition: tls.hh:89
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
bool G4bool
Definition: G4Types.hh:79
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static const double eV
Definition: G4SIunits.hh:212
#define G4endl
Definition: G4ios.hh:61
static const double barn
Definition: G4SIunits.hh:104
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
std::vector< G4Element * > G4ElementTable
double G4double
Definition: G4Types.hh:76
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
Here is the call graph for this function:

◆ EnableOnFlightDopplerBroadening()

void G4ParticleHPElasticData::EnableOnFlightDopplerBroadening ( )
inline

Definition at line 84 of file G4ParticleHPElasticData.hh.

Here is the call graph for this function:

◆ GetCrossSection()

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

Definition at line 209 of file G4ParticleHPElasticData.cc.

210 {
211  G4double result = 0;
212  G4bool outOfRange;
213  G4int index = anE->GetIndex();
214 
215  // prepare neutron
216  G4double eKinetic = aP->GetKineticEnergy();
217 
218  if ( !onFlightDB )
219  {
220  //NEGLECT_DOPPLER_B.
221  G4double factor = 1.0;
222  if ( eKinetic < aT * k_Boltzmann )
223  {
224  // below 0.1 eV neutrons
225  // Have to do some, but now just igonre.
226  // Will take care after performance check.
227  // factor = factor * targetV;
228  }
229  return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
230  }
231 
232  G4ReactionProduct theNeutron( aP->GetDefinition() );
233  theNeutron.SetMomentum( aP->GetMomentum() );
234  theNeutron.SetKineticEnergy( eKinetic );
235 
236  // prepare thermal nucleus
237  G4Nucleus aNuc;
238  G4double eps = 0.0001;
239  G4double theA = anE->GetN();
240  G4double theZ = anE->GetZ();
241  G4double eleMass;
242 
243 
244  eleMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) )
246 
247  G4ReactionProduct boosted;
248  G4double aXsection;
249 
250  // MC integration loop
251  G4int counter = 0;
252  G4double buffer = 0;
253  G4int size = G4int(std::max(10., aT/60*kelvin));
254  G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
255  G4double neutronVMag = neutronVelocity.mag();
256 
257  while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.03*buffer) // Loop checking, 11.05.2015, T. Koi
258  {
259  if(counter) buffer = result/counter;
260  while (counter<size) // Loop checking, 11.05.2015, T. Koi
261  {
262  counter ++;
263  G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
264  boosted.Lorentz(theNeutron, aThermalNuc);
265  G4double theEkin = boosted.GetKineticEnergy();
266  aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
267  // velocity correction.
268  G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
269  aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
270  result += aXsection;
271  }
272  size += size;
273  }
274  result /= counter;
275 /*
276  // Checking impact of G4NEUTRONHP_NEGLECT_DOPPLER
277  G4cout << " result " << result << " "
278  << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
279  << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
280 */
281  return result;
282 }
G4double GetMass() const
static G4double GetNuclearMass(const G4double A, const G4double Z)
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
Int_t index
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition: G4Nucleus.cc:143
size_t GetIndex() const
Definition: G4Element.hh:181
void SetMomentum(const G4double x, const G4double y, const G4double z)
#define buffer
Definition: xmlparse.cc:628
static const G4double eps
int G4int
Definition: G4Types.hh:78
G4double GetN() const
Definition: G4Element.hh:134
G4double GetKineticEnergy() const
G4ThreeVector GetMomentum() const
float k_Boltzmann
Definition: hepunit.py:299
bool G4bool
Definition: G4Types.hh:79
double mag() const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static const double kelvin
Definition: G4SIunits.hh:278
static const G4double factor
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
G4double GetZ() const
Definition: G4Element.hh:131
G4ThreeVector GetMomentum() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetIsoCrossSection()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 87 of file G4ParticleHPElasticData.cc.

92 {
93  if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
94 
95  ke_cache = dp->GetKineticEnergy();
96  element_cache = element;
98  G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
99  xs_cache = xs;
100  return xs;
101 }
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)
string material
Definition: eplot.py:19
G4double GetKineticEnergy() const
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ GetVerboseLevel()

G4int G4ParticleHPElasticData::GetVerboseLevel ( ) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 285 of file G4ParticleHPElasticData.cc.

286 {
288 }
static G4ParticleHPManager * GetInstance()
Here is the call graph for this function:
Here is the caller graph for this function:

◆ IgnoreOnFlightDopplerBroadening()

void G4ParticleHPElasticData::IgnoreOnFlightDopplerBroadening ( )
inline

Definition at line 83 of file G4ParticleHPElasticData.hh.

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 73 of file G4ParticleHPElasticData.cc.

77 {
78 
79  G4double eKin = dp->GetKineticEnergy();
80  if ( eKin > GetMaxKinEnergy()
81  || eKin < GetMinKinEnergy()
82  || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
83 
84  return true;
85 }
G4double GetKineticEnergy() const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4ParticleDefinition * GetDefinition() const
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ SetVerboseLevel()

void G4ParticleHPElasticData::SetVerboseLevel ( G4int  newValue)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 291 of file G4ParticleHPElasticData.cc.

292 {
294 }
static G4ParticleHPManager * GetInstance()
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ element_cache

const G4Element* G4ParticleHPElasticData::element_cache
private

Definition at line 98 of file G4ParticleHPElasticData.hh.

◆ instanceOfWorker

G4bool G4ParticleHPElasticData::instanceOfWorker
private

Definition at line 94 of file G4ParticleHPElasticData.hh.

◆ ke_cache

G4double G4ParticleHPElasticData::ke_cache
private

Definition at line 96 of file G4ParticleHPElasticData.hh.

◆ material_cache

const G4Material* G4ParticleHPElasticData::material_cache
private

Definition at line 99 of file G4ParticleHPElasticData.hh.

◆ onFlightDB

G4bool G4ParticleHPElasticData::onFlightDB
private

Definition at line 93 of file G4ParticleHPElasticData.hh.

◆ theCrossSections

G4PhysicsTable* G4ParticleHPElasticData::theCrossSections
private

Definition at line 92 of file G4ParticleHPElasticData.hh.

◆ xs_cache

G4double G4ParticleHPElasticData::xs_cache
private

Definition at line 97 of file G4ParticleHPElasticData.hh.


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