Geant4  10.02.p03
G4ParticleHPCaptureData Class Reference

#include <G4ParticleHPCaptureData.hh>

Inheritance diagram for G4ParticleHPCaptureData:
Collaboration diagram for G4ParticleHPCaptureData:

Public Member Functions

 G4ParticleHPCaptureData ()
 
 ~G4ParticleHPCaptureData ()
 
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 ()
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int)
 
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 G4ParticleHPCaptureData.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPCaptureData()

G4ParticleHPCaptureData::G4ParticleHPCaptureData ( )

Definition at line 50 of file G4ParticleHPCaptureData.cc.

51 :G4VCrossSectionDataSet("NeutronHPCaptureXS")
52 {
53  SetMinKinEnergy( 0*MeV );
54  SetMaxKinEnergy( 20*MeV );
55 
56  theCrossSections = 0;
57  onFlightDB = true;
58 
59  instanceOfWorker = false;
61  instanceOfWorker = true;
62  }
63  //BuildPhysicsTable(*G4Neutron::Neutron());
64 }
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:

◆ ~G4ParticleHPCaptureData()

G4ParticleHPCaptureData::~G4ParticleHPCaptureData ( )

Definition at line 66 of file G4ParticleHPCaptureData.cc.

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

Member Function Documentation

◆ BuildPhysicsTable()

void G4ParticleHPCaptureData::BuildPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 114 of file G4ParticleHPCaptureData.cc.

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 capture reaction of neutrons (<20MeV)." << G4endl;
124  onFlightDB = false;
125  }
126 
127  if ( G4Threading::IsWorkerThread() ) {
129  return;
130  }
131 
132  size_t numberOfElements = G4Element::GetNumberOfElements();
133  // G4cout << "CALLED G4ParticleHPCaptureData::BuildPhysicsTable "<<numberOfElements<<G4endl;
134  // TKDB
135  //if ( theCrossSections == 0 ) theCrossSections = new G4PhysicsTable( numberOfElements );
136  if ( theCrossSections == NULL )
137  theCrossSections = new G4PhysicsTable( numberOfElements );
138  else
140 
141  // make a PhysicsVector for each element
142 
143  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
144  for( size_t i=0; i<numberOfElements; ++i )
145  {
146  if(getenv("CaptureDataIndexDebug"))
147  {
148  G4int index_debug = ((*theElementTable)[i])->GetIndex();
149  G4cout << "IndexDebug "<< i <<" "<<index_debug<<G4endl;
150  }
152  Instance(G4Neutron::Neutron())->MakePhysicsVector((*theElementTable)[i], this);
153  theCrossSections->push_back(physVec);
154  }
155 
157 }
static G4ParticleHPManager * GetInstance()
void push_back(G4PhysicsVector *)
void RegisterCaptureCrossSections(G4PhysicsTable *val)
#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
G4PhysicsTable * GetCaptureCrossSections()
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 G4ParticleHPCaptureData::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 292 of file G4ParticleHPCaptureData.cc.

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

◆ DumpPhysicsTable()

void G4ParticleHPCaptureData::DumpPhysicsTable ( const G4ParticleDefinition aP)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 159 of file G4ParticleHPCaptureData.cc.

160 {
161  if(&aP!=G4Neutron::Neutron())
162  throw G4HadronicException(__FILE__, __LINE__, "Attempt to use NeutronHP data for particles other than neutrons!!!");
163 
164 //
165 // Dump element based cross section
166 // range 10e-5 eV to 20 MeV
167 // 10 point per decade
168 // in barn
169 //
170 
171  G4cout << G4endl;
172  G4cout << G4endl;
173  G4cout << "Capture Cross Section of Neutron HP"<< G4endl;
174  G4cout << "(Pointwise cross-section at 0 Kelvin.)" << G4endl;
175  G4cout << G4endl;
176  G4cout << "Name of Element" << G4endl;
177  G4cout << "Energy[eV] XS[barn]" << G4endl;
178  G4cout << G4endl;
179 
180  size_t numberOfElements = G4Element::GetNumberOfElements();
181  static G4ThreadLocal G4ElementTable *theElementTable = 0 ; if (!theElementTable) theElementTable= G4Element::GetElementTable();
182 
183  for ( size_t i = 0 ; i < numberOfElements ; ++i )
184  {
185 
186  G4cout << (*theElementTable)[i]->GetName() << G4endl;
187 
188  G4int ie = 0;
189 
190  for ( ie = 0 ; ie < 130 ; ie++ )
191  {
192  G4double eKinetic = 1.0e-5 * G4Pow::GetInstance()->powA ( 10.0 , ie/10.0 ) *eV;
193  G4bool outOfRange = false;
194 
195  if ( eKinetic < 20*MeV )
196  {
197  G4cout << eKinetic/eV << " " << (*((*theCrossSections)(i))).GetValue(eKinetic, outOfRange)/barn << G4endl;
198  }
199 
200  }
201 
202  G4cout << G4endl;
203  }
204 
205 
206 // G4cout << "G4ParticleHPCaptureData::DumpPhysicsTable still to be implemented"<<G4endl;
207 }
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 G4ParticleHPCaptureData::EnableOnFlightDopplerBroadening ( )
inline

Definition at line 83 of file G4ParticleHPCaptureData.hh.

Here is the call graph for this function:

◆ GetCrossSection()

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

Definition at line 212 of file G4ParticleHPCaptureData.cc.

213 {
214  G4double result = 0;
215  G4bool outOfRange;
216  G4int index = anE->GetIndex();
217 
218  // prepare neutron
219  G4double eKinetic = aP->GetKineticEnergy();
220 
221  if ( !onFlightDB )
222  {
223  //NEGLECT_DOPPLER
224  G4double factor = 1.0;
225  if ( eKinetic < aT * k_Boltzmann )
226  {
227  // below 0.1 eV neutrons
228  // Have to do some, but now just igonre.
229  // Will take care after performance check.
230  // factor = factor * targetV;
231  }
232  return ( (*((*theCrossSections)(index))).GetValue(eKinetic, outOfRange) )* factor;
233  }
234 
235  G4ReactionProduct theNeutron( aP->GetDefinition() );
236  theNeutron.SetMomentum( aP->GetMomentum() );
237  theNeutron.SetKineticEnergy( eKinetic );
238 
239  // prepare thermal nucleus
240  G4Nucleus aNuc;
241  G4double eps = 0.0001;
242  G4double theA = anE->GetN();
243  G4double theZ = anE->GetZ();
244  G4double eleMass;
245  eleMass = G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theA+eps) , static_cast<G4int>(theZ+eps) ) / G4Neutron::Neutron()->GetPDGMass();
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, or luminosity factor...
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 G4ParticleHPCaptureData::GetIsoCrossSection ( const G4DynamicParticle dp,
G4int  ,
G4int  ,
const G4Isotope ,
const G4Element element,
const G4Material material 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 88 of file G4ParticleHPCaptureData.cc.

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

◆ GetVerboseLevel()

G4int G4ParticleHPCaptureData::GetVerboseLevel ( ) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 284 of file G4ParticleHPCaptureData.cc.

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

◆ IgnoreOnFlightDopplerBroadening()

void G4ParticleHPCaptureData::IgnoreOnFlightDopplerBroadening ( )
inline

Definition at line 82 of file G4ParticleHPCaptureData.hh.

◆ IsIsoApplicable()

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 75 of file G4ParticleHPCaptureData.cc.

79 {
80  G4double eKin = dp->GetKineticEnergy();
81  if ( eKin > GetMaxKinEnergy()
82  || eKin < GetMinKinEnergy()
83  || dp->GetDefinition() != G4Neutron::Neutron() ) return false;
84 
85  return true;
86 }
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 G4ParticleHPCaptureData::SetVerboseLevel ( G4int  newValue)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 288 of file G4ParticleHPCaptureData.cc.

289 {
291 }
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* G4ParticleHPCaptureData::element_cache
private

Definition at line 99 of file G4ParticleHPCaptureData.hh.

◆ instanceOfWorker

G4bool G4ParticleHPCaptureData::instanceOfWorker
private

Definition at line 95 of file G4ParticleHPCaptureData.hh.

◆ ke_cache

G4double G4ParticleHPCaptureData::ke_cache
private

Definition at line 97 of file G4ParticleHPCaptureData.hh.

◆ material_cache

const G4Material* G4ParticleHPCaptureData::material_cache
private

Definition at line 100 of file G4ParticleHPCaptureData.hh.

◆ onFlightDB

G4bool G4ParticleHPCaptureData::onFlightDB
private

Definition at line 94 of file G4ParticleHPCaptureData.hh.

◆ theCrossSections

G4PhysicsTable* G4ParticleHPCaptureData::theCrossSections
private

Definition at line 92 of file G4ParticleHPCaptureData.hh.

◆ xs_cache

G4double G4ParticleHPCaptureData::xs_cache
private

Definition at line 98 of file G4ParticleHPCaptureData.hh.


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