Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4NeutronHPElasticData Class Reference

#include <G4NeutronHPElasticData.hh>

Inheritance diagram for G4NeutronHPElasticData:
G4VCrossSectionDataSet

Public Member Functions

 G4NeutronHPElasticData ()
 
 ~G4NeutronHPElasticData ()
 
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
 
- 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)
 
virtual void CrossSectionDescription (std::ostream &) const
 
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 49 of file G4NeutronHPElasticData.hh.

Constructor & Destructor Documentation

G4NeutronHPElasticData::G4NeutronHPElasticData ( )

Definition at line 46 of file G4NeutronHPElasticData.cc.

References python.hepunit::MeV, G4VCrossSectionDataSet::SetMaxKinEnergy(), and G4VCrossSectionDataSet::SetMinKinEnergy().

47 :G4VCrossSectionDataSet("NeutronHPElasticXS")
48 {
49  SetMinKinEnergy( 0*MeV );
50  SetMaxKinEnergy( 20*MeV );
51 
52  ke_cache = 0.0;
53  xs_cache = 0.0;
54  element_cache = NULL;
55  material_cache = NULL;
56 
57  theCrossSections = 0;
58  onFlightDB = true;
59 // BuildPhysicsTable( *G4Neutron::Neutron() );
60 }
G4VCrossSectionDataSet(const G4String &nam="")
void SetMinKinEnergy(G4double value)
void SetMaxKinEnergy(G4double value)
G4NeutronHPElasticData::~G4NeutronHPElasticData ( )

Definition at line 62 of file G4NeutronHPElasticData.cc.

References G4PhysicsTable::clearAndDestroy().

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

Member Function Documentation

void G4NeutronHPElasticData::BuildPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 108 of file G4NeutronHPElasticData.cc.

References G4PhysicsTable::clearAndDestroy(), G4cout, G4endl, G4ThreadLocal, G4Element::GetElementTable(), G4Element::GetNumberOfElements(), G4NeutronHPData::Instance(), G4NeutronHPData::MakePhysicsVector(), G4Neutron::Neutron(), and G4PhysicsTable::push_back().

109 {
110 
111  if(&aP!=G4Neutron::Neutron())
112  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
113 
114 //080428
115  if ( getenv( "G4NEUTRONHP_NEGLECT_DOPPLER" ) )
116  {
117  G4cout << "Find environment variable of \"G4NEUTRONHP_NEGLECT_DOPPLER\"." << G4endl;
118  G4cout << "On the fly Doppler broadening will be neglect in the cross section calculation of elastic scattering of neutrons (<20MeV)." << G4endl;
119  onFlightDB = false;
120  }
121 
122  size_t numberOfElements = G4Element::GetNumberOfElements();
123 // TKDB
124  //if ( theCrossSections == 0 ) 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 G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) 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 }
G4PhysicsVector * MakePhysicsVector(G4Element *thE, G4NeutronHPFissionData *theP)
void push_back(G4PhysicsVector *)
#define G4ThreadLocal
Definition: tls.hh:52
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:402
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4NeutronHPData * Instance()
#define G4endl
Definition: G4ios.hh:61
std::vector< G4Element * > G4ElementTable
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
void clearAndDestroy()
void G4NeutronHPElasticData::DumpPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 141 of file G4NeutronHPElasticData.cc.

References python.hepunit::barn, python.hepunit::eV, G4cout, G4endl, G4ThreadLocal, G4Element::GetElementTable(), G4Element::GetNumberOfElements(), python.hepunit::MeV, and G4Neutron::Neutron().

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 << "Elastic 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 G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) 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 
188 // G4cout << "G4NeutronHPElasticData::DumpPhysicsTable still to be implemented"<<G4endl;
189 }
#define G4ThreadLocal
Definition: tls.hh:52
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
#define G4endl
Definition: G4ios.hh:61
std::vector< G4Element * > G4ElementTable
double G4double
Definition: G4Types.hh:76
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:395
void G4NeutronHPElasticData::EnableOnFlightDopplerBroadening ( )
inline

Definition at line 82 of file G4NeutronHPElasticData.hh.

82 { onFlightDB = true; };
G4double G4NeutronHPElasticData::GetCrossSection ( const G4DynamicParticle aP,
const G4Element anE,
G4double  aT 
)

Definition at line 197 of file G4NeutronHPElasticData.cc.

References buffer, G4DynamicParticle::GetDefinition(), G4Element::GetIndex(), G4DynamicParticle::GetKineticEnergy(), G4ReactionProduct::GetKineticEnergy(), G4ReactionProduct::GetMass(), G4DynamicParticle::GetMomentum(), G4ReactionProduct::GetMomentum(), G4Element::GetN(), G4NucleiProperties::GetNuclearMass(), G4ParticleDefinition::GetPDGMass(), G4Nucleus::GetThermalNucleus(), G4Element::GetZ(), python.hepunit::k_Boltzmann, python.hepunit::kelvin, G4ReactionProduct::Lorentz(), CLHEP::Hep3Vector::mag(), G4INCL::Math::max(), G4Neutron::Neutron(), and G4ReactionProduct::SetMomentum().

Referenced by GetIsoCrossSection().

198 {
199  G4double result = 0;
200  G4bool outOfRange;
201  G4int index = anE->GetIndex();
202 
203  // prepare neutron
204  G4double eKinetic = aP->GetKineticEnergy();
205 
206  // T. K.
207 // if ( getenv( "G4NEUTRONHP_NEGLECT_DOPPLER" ) )
208 //080428
209  if ( !onFlightDB )
210  {
211  G4double factor = 1.0;
212  if ( eKinetic < aT * k_Boltzmann )
213  {
214  // below 0.1 eV neutrons
215  // Have to do some, but now just igonre.
216  // Will take care after performance check.
217  // factor = factor * targetV;
218  }
219  return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
220  }
221 
222  G4ReactionProduct theNeutron( aP->GetDefinition() );
223  theNeutron.SetMomentum( aP->GetMomentum() );
224  theNeutron.SetKineticEnergy( eKinetic );
225 
226  // prepare thermal nucleus
227  G4Nucleus aNuc;
228  G4double eps = 0.0001;
229  G4double theA = anE->GetN();
230  G4double theZ = anE->GetZ();
231  G4double eleMass;
232 
233 
234  eleMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) )
236 
237  G4ReactionProduct boosted;
238  G4double aXsection;
239 
240  // MC integration loop
241  G4int counter = 0;
242  G4double buffer = 0;
243  G4int size = G4int(std::max(10., aT/60*kelvin));
244  G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
245  G4double neutronVMag = neutronVelocity.mag();
246 
247  while(counter == 0 || std::abs(buffer-result/std::max(1,counter)) > 0.03*buffer)
248  {
249  if(counter) buffer = result/counter;
250  while (counter<size)
251  {
252  counter ++;
253  G4ReactionProduct aThermalNuc = aNuc.GetThermalNucleus(eleMass, aT);
254  boosted.Lorentz(theNeutron, aThermalNuc);
255  G4double theEkin = boosted.GetKineticEnergy();
256  aXsection = (*((*theCrossSections)(index))).GetValue(theEkin, outOfRange);
257  // velocity correction.
258  G4ThreeVector targetVelocity = 1./aThermalNuc.GetMass()*aThermalNuc.GetMomentum();
259  aXsection *= (targetVelocity-neutronVelocity).mag()/neutronVMag;
260  result += aXsection;
261  }
262  size += size;
263  }
264  result /= counter;
265 /*
266  // Checking impact of G4NEUTRONHP_NEGLECT_DOPPLER
267  G4cout << " result " << result << " "
268  << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) << " "
269  << (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) /result << G4endl;
270 */
271  return result;
272 }
static G4double GetNuclearMass(const G4double A, const G4double Z)
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
G4double GetKineticEnergy() const
G4double GetN() const
Definition: G4Element.hh:134
void SetMomentum(const G4double x, const G4double y, const G4double z)
G4double GetZ() const
Definition: G4Element.hh:131
#define buffer
Definition: xmlparse.cc:611
G4ParticleDefinition * GetDefinition() const
G4ReactionProduct GetThermalNucleus(G4double aMass, G4double temp=-1) const
Definition: G4Nucleus.cc:130
int G4int
Definition: G4Types.hh:78
float k_Boltzmann
Definition: hepunit.py:299
bool G4bool
Definition: G4Types.hh:79
size_t GetIndex() const
Definition: G4Element.hh:181
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4double GetKineticEnergy() const
G4double GetPDGMass() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4ThreeVector GetMomentum() const
double G4double
Definition: G4Types.hh:76
double mag() const
G4double GetMass() const
G4ThreeVector GetMomentum() const
G4double G4NeutronHPElasticData::GetIsoCrossSection ( const G4DynamicParticle dp,
G4int  ,
G4int  ,
const G4Isotope ,
const G4Element element,
const G4Material material 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 82 of file G4NeutronHPElasticData.cc.

References GetCrossSection(), G4DynamicParticle::GetKineticEnergy(), and eplot::material.

87 {
88  if ( dp->GetKineticEnergy() == ke_cache && element == element_cache && material == material_cache ) return xs_cache;
89 
90  ke_cache = dp->GetKineticEnergy();
91  element_cache = element;
92  material_cache = material;
93  G4double xs = GetCrossSection( dp , element , material->GetTemperature() );
94  xs_cache = xs;
95  return xs;
96 }
G4double GetKineticEnergy() const
string material
Definition: eplot.py:19
double G4double
Definition: G4Types.hh:76
G4double GetCrossSection(const G4DynamicParticle *, const G4Element *, G4double aT)
G4int G4NeutronHPElasticData::GetVerboseLevel ( ) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 275 of file G4NeutronHPElasticData.cc.

References G4NeutronHPManager::GetInstance(), and G4NeutronHPManager::GetVerboseLevel().

276 {
278 }
static G4NeutronHPManager * GetInstance()
void G4NeutronHPElasticData::IgnoreOnFlightDopplerBroadening ( )
inline

Definition at line 81 of file G4NeutronHPElasticData.hh.

81 { onFlightDB = false; };
G4bool G4NeutronHPElasticData::IsIsoApplicable ( const G4DynamicParticle dp,
G4int  ,
G4int  ,
const G4Element ,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 68 of file G4NeutronHPElasticData.cc.

References G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4VCrossSectionDataSet::GetMaxKinEnergy(), G4VCrossSectionDataSet::GetMinKinEnergy(), and G4Neutron::Neutron().

72 {
73 
74  G4double eKin = dp->GetKineticEnergy();
75  if ( eKin > GetMaxKinEnergy()
76  || eKin < GetMinKinEnergy()
77  || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
78 
79  return true;
80 }
G4double GetKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
double G4double
Definition: G4Types.hh:76
void G4NeutronHPElasticData::SetVerboseLevel ( G4int  newValue)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 281 of file G4NeutronHPElasticData.cc.

References G4NeutronHPManager::GetInstance(), and G4NeutronHPManager::SetVerboseLevel().

282 {
284 }
static G4NeutronHPManager * GetInstance()
void SetVerboseLevel(G4int i)

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